Logo Coherent WaveBurst  
Reference Guide
Logo
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
wmdc_config_wnb.C
Go to the documentation of this file.
1 {
2  #include <vector>
3 
4  // --------------------------------------------------------
5  // define job parameters
6  // --------------------------------------------------------
7  TString frDir = "frames";
8  TString frLabel = "TestWNB";
9 
10  // define frame segment list
11  TString segmentList = "Segments/S6a_LS-segs.txt";
12 
13  // select segments
14  int jobmin = 210; // start segment
15  int jobmax = 312; // end segment
16  int jobstep = 10; // frames per job
17 
18  // --------------------------------------------------------
19  // define network
20  // --------------------------------------------------------
21  int nIFO=3;
22  TString ifo[nIFO]={"L1","H1","V1"};
23  CWB::mdc MDC(nIFO,ifo);
24 
25  // --------------------------------------------------------
26  // define waveforms
27  // --------------------------------------------------------
28  int nWF=8;
29  double F[nWF] = {100. ,100. ,100 ,1000. ,1000. ,1000. ,1000. ,1000.};
30  double B[nWF] = {10. ,100. ,100 ,100. ,100. ,1000. ,1000. ,1000.};
31  double D[nWF] = {0.1 ,0.1 ,0.01 ,0.1 ,0.01 ,0.1 ,0.01 ,0.001};
32 
33  char wf_name[256];
35  vector<mdcpar> par; par.resize(3);
36  for (int n=0;n<nWF;n++) {
37  par[0].name="frequency"; par[0].value=F[n];
38  par[1].name="bandwidth"; par[1].value=B[n];
39  par[2].name="duration"; par[2].value=D[n];
40  for(int m=0;m<15;m++) MDC.AddWaveform(MDC_WNB, par);
41  }
42  MDC.Print();
43 
44  // --------------------------------------------------------
45  // define injection parameters
46  // --------------------------------------------------------
47  MDC.SetInjHrss(2.5e-21);
48  MDC.SetInjRate(0.01);
49  MDC.SetInjJitter(10.0);
50 
51  // --------------------------------------------------------
52  // define sky distribution
53  // --------------------------------------------------------
54  vector<mdcpar> par; par.resize(3);
55  par[0].name="entries";par[0].value=100000; // pool of events
56  par[1].name="rho_min";par[1].value=1; // min rho // Mpc
57  par[2].name="rho_max";par[2].value=100; // max rho // Mpc
59 
60  //vector<mdcpar> par; par.resize(1);
61  //par[0].name="entries";par[0].value=100000; // number of entries
62  //MDC.SetSkyDistribution(MDC_MNGD,par);
63 
64  //vector<mdcpar> par; par.resize(1);
65  //par[0].name="distance_thr";par[0].value=100; // distance max = 100 Mpc
66  //MDC.SetSkyDistribution(MDC_GWGC,"../data/GWGCCatalog_Rev1d2.txt",par);
67 
68  //vector<mdcpar> par; par.resize(3);
69  //par[0].name="theta";par[0].value=30; // theta
70  //par[1].name="phi"; par[1].value=60; // phi
71  //par[2].name="psi"; par[2].value=90; // psi
72  //MDC.SetSkyDistribution(MDC_EARTH_FIX,par);
73 }
int nWF
Definition: mdc.hh:189
void SetInjRate(double inj_rate=MDC_INJ_RATE)
Definition: mdc.hh:279
TString frLabel
int n
Definition: cwb_net.C:10
double B[nWF]
TString segmentList
TString("c")
CWB::mdc MDC(nIFO, ifo)
waveform wf
vector< mdcpar > par
double F[nWF]
int m
Definition: cwb_net.C:10
int jobmin
mdcid AddWaveform(MDC_TYPE mdc_type, vector< mdcpar > par, TString uname="")
Definition: mdc.cc:445
void Print(int level=0)
Definition: mdc.cc:2707
void SetSkyDistribution(MDC_DISTRIBUTION sky_distribution, vector< mdcpar > par, int seed=0, bool add=false)
Definition: mdc.cc:3415
Definition: mdc.hh:216
TString frDir[NIFO_MAX]
void SetInjJitter(double inj_jitter=MDC_INJ_JITTER)
Definition: mdc.hh:285
double e
int jobstep
int jobmax
char wf_name[256]
double D[nWF]
int nIFO
void SetInjHrss(double inj_hrss=MDC_INJ_HRSS)
Definition: mdc.hh:275
TString ifo[nIFO]
Definition: mdc.hh:127