Logo Coherent WaveBurst  
Reference Guide
Logo
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
DrawSGC.C
Go to the documentation of this file.
2 
3 void DrawSGC() {
4  //
5  // Show how to use mdc class to get & draw SGC waveforms
6  // Author : Gabriele Vedovato
7 
8  #include <vector>
9  #define N_IFO 3
10 
11  TString ifo[N_IFO]={"L1","H1","V1"};
12  MDC = new CWB::mdc(N_IFO,ifo);
13 
14  char wf_name[256];
15  waveform wf;
16  vector<mdcpar> par;
17 
18  par.resize(2);
19 
20  par[0].name="frequency"; par[0].value=100.;
21  par[1].name="Q"; par[1].value=8.9;
22  MDC->AddWaveform(MDC_SGC, par);
23 
24  par[0].name="frequency"; par[0].value=153.;
25  par[1].name="Q"; par[1].value=8.9;
26  MDC->AddWaveform(MDC_SGC, par);
27 
28  par[0].name="frequency"; par[0].value=1053.;
29  par[1].name="Q"; par[1].value=9;
30  MDC->AddWaveform(MDC_SGC, par);
31 
32  par[0].name="frequency"; par[0].value=1304.;
33  par[1].name="Q"; par[1].value=9;
34  MDC->AddWaveform(MDC_SGC, par);
35 
36  MDC->Print();
37 
38  // --------------------------------------------------------
39  // define injection parameters
40  // --------------------------------------------------------
41  MDC->SetInjHrss(2.5e-21);
42  MDC->SetInjRate(0.02);
43  MDC->SetInjJitter(10.0);
44 
45  // --------------------------------------------------------
46  // define sky distribution
47  // --------------------------------------------------------
48  // 934361014 ra : 3.4909999370575 dec : -27.687000274658 radius : 1.3999999761581
49  float DEC = -27.687000274658;
50  float RA = 3.4909999370575;
51  float seed = 0;
52  par.resize(3);
53  par[0].name="theta"; par[0].value=DEC;
54  par[1].name="phi"; par[1].value=RA;
55  par[2].name="entries"; par[2].value=1000;
57 
58  // Get the first waveform hp,hx components
59  wf = MDC->GetWaveform("SGC1304Q9",0);
60 
61  cout << "size : " << wf.hp.size() << " rate : " << wf.hp.rate() << " start : " << (int)wf.hp.start() << endl;
62  wf.hp.start(0); // set start to 0 (needed by draw Method)
63  wf.hx.start(0);
64 
65  //MDC->Draw(wf.hp,MDC_TIME,"ALP ZOOM");
66  //MDC->Draw(wf.hx,MDC_TIME,"same",kRed);
67 
68  //MDC->Draw(wf.hp,MDC_FFT,"ALP ZOOM"); // draw hp in frequency domain
69  //MDC->Draw(wf.hp,MDC_TF,"ZOOM"); // draw hp in time-frequency domain;
70 
71 }
virtual size_t size() const
Definition: wavearray.hh:127
#define N_IFO
Definition: mdc.hh:189
virtual void rate(double r)
Definition: wavearray.hh:123
void SetInjRate(double inj_rate=MDC_INJ_RATE)
Definition: mdc.hh:279
TString("c")
TString mdc[4]
wavearray< double > hp
Definition: mdc.hh:196
waveform wf
virtual void start(double s)
Definition: wavearray.hh:119
mdcid AddWaveform(MDC_TYPE mdc_type, vector< mdcpar > par, TString uname="")
Definition: mdc.cc:445
void Print(int level=0)
Definition: mdc.cc:2707
char ifo[NIFO_MAX][8]
void SetSkyDistribution(MDC_DISTRIBUTION sky_distribution, vector< mdcpar > par, int seed=0, bool add=false)
Definition: mdc.cc:3415
Definition: mdc.hh:216
i() int(T_cor *100))
CWB::mdc * MDC
Definition: DrawSGC.C:1
Definition: mdc.hh:121
void SetInjJitter(double inj_jitter=MDC_INJ_JITTER)
Definition: mdc.hh:285
vector< mdcpar > par
double e
wavearray< double > hx
Definition: mdc.hh:197
void DrawSGC()
Definition: DrawSGC.C:3
waveform GetWaveform(int ID, int id=0)
Definition: mdc.cc:2512
char wf_name[256]
void SetInjHrss(double inj_hrss=MDC_INJ_HRSS)
Definition: mdc.hh:275
TString name
Definition: mdc.hh:178