Logo Coherent WaveBurst  
Reference Guide
Logo
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
DrawSGE.C
Go to the documentation of this file.
2 
3 void DrawSGE() {
4  //
5  // Show how to use mdc class to get & draw SG elliptical waveforms
6  // Author : Gabriele Vedovato
7 
8  #define N_IFO 3
9 
10  TString ifo[N_IFO]={"L1","H1","V1"};
11  MDC = new CWB::mdc(N_IFO,ifo);
12 
13  char wf_name[256];
14  waveform wf;
15  vector<mdcpar> par;
16 
17  par.resize(2);
18 
19  par[0].name="frequency"; par[0].value=100.;
20  par[1].name="Q"; par[1].value=8.9;
21  MDC->AddWaveform(MDC_SGE, par);
22 
23  MDC->Print();
24 
25  // --------------------------------------------------------
26  // define injection parameters
27  // --------------------------------------------------------
28  MDC->SetInjHrss(2.5e-21);
29  MDC->SetInjRate(0.02);
30  MDC->SetInjJitter(10.0);
31 
32  // -----------------------------------------------------------------------------------
33  // define sky distribution : fixed earth/celestial direction & ellipticity : iota = 45
34  // -----------------------------------------------------------------------------------
35  float DEC = -30.0; // [-90:90]
36  float RA = 80.0; // [0:360]
37  float PSI = 0.0; // [0:180]
38  float IOTA = 45; // 0/90 -> circular/linear
39  float GPS = 10000000;
40  float seed = 0;
41  par.resize(6);
42  par[0].name="theta"; par[0].value=DEC;
43  par[1].name="phi"; par[1].value=RA;
44  par[2].name="psi"; par[2].value=PSI;
45  par[3].name="entries"; par[3].value=1000;
46  par[4].name="iota"; par[4].value=IOTA; // ellipticity [0:180] deg
47  // if<0 || >180 -> random
48  par[5].name="gps"; par[5].value=GPS; // fixed injection time
49 
50  //MDC->SetSkyDistribution(MDC_CELESTIAL_FIX,par,seed); // theta,phi are in celestial coordinates
51  MDC->SetSkyDistribution(MDC_EARTH_FIX,par,seed); // theta,phi are in earth coordinates
52 
53  // Get the first waveform hp,hx components
54  wf = MDC->GetWaveform("SGE100Q8d9",0);
55  if((wf.hp.size()==0)||(wf.hx.size()==0)) {
56  cout << "Error : Waveform not present in the MDC set !!!" << endl;
57  gSystem->Exit(1);
58  }
59 
60  cout << "hp : size : " << wf.hp.size() << " rate : " << wf.hp.rate() << " start : " << (int)wf.hp.start() << endl;
61  cout << "hx : size : " << wf.hx.size() << " rate : " << wf.hx.rate() << " start : " << (int)wf.hx.start() << endl;
62  wf.hp.start(0); // set start to 0 (needed by draw Method)
63  wf.hx.start(0);
64 
65  // uncomment the following lines to draw hp,hx components
66 
67  //MDC->Draw(wf.hp,MDC_TIME,"ALP ZOOM");
68  //MDC->Draw(wf.hx,MDC_TIME,"same",kRed);
69 
70  //MDC->Draw(wf.hp,MDC_FFT,"ALP ZOOM"); // draw hp in frequency domain
71  //MDC->Draw(wf.hp,MDC_TF,"ZOOM"); // draw hp in time-frequency domain;
72 
73  // build elliptical waveform
74 /*
75  double deg2rad = TMath::Pi()/180.;
76  double ePlus = (1+cos(IOTA*deg2rad)*cos(IOTA*deg2rad))/2;
77  double eCross = cos(IOTA*deg2rad);
78  cout << "ePlus : " << ePlus << " eCross : " << eCross << endl;
79  wavearray<double> ef=wf.hp;
80  for(int i=0;i<ef.size();i++) ef[i]=ePlus*wf.hp[i]+eCross*wf.hx[i];
81  MDC->DrawTime(ef,"ALP ZOOM");
82 */
83 
84  // --------------------------------------------------------
85  // get and plot the waveforms on L1,H1,V1
86  // print the mdc parameters
87  // --------------------------------------------------------
88 
89  // get mdc buffer draw waveform
90  wavearray<double> x(4*16384);
91  x.rate(16384);
92 
93  // draw V1
94  x.start(GPS-2);
95  TString log = MDC->Get(x,ifo[0]); // get log
96  x.start(0);
97  MDC->DrawTime(x,"ALP ZOOM");
98 
99  // get central time
100  double To = MDC->GetCentralTime(x);
101  cout << "Central Time :" << To << " (sec)" << endl;
102 
103  // get central frequency
104  double Fo = MDC->GetCentralFrequency(x);
105  cout << "Central Frequency : " << Fo << " (Hz)" << endl;
106 
107  // draw H1
108  x.start(GPS-2);
109  MDC->Get(x,ifo[1]);
110  x.start(0);
111  MDC->DrawTime(x,"SAME",kRed);
112 
113  // draw L1
114  x.start(GPS-2);
115  MDC->Get(x,ifo[2]);
116  x.start(0);
117  MDC->DrawTime(x,"SAME",kBlue);
118 
119  // print mdc parameters
120  cout << endl;
121  cout << log << endl;
122 }
virtual size_t size() const
Definition: wavearray.hh:127
Definition: mdc.hh:189
TString Get(wavearray< double > &x, TString ifo)
Definition: mdc.cc:1502
virtual void rate(double r)
Definition: wavearray.hh:123
void SetInjRate(double inj_rate=MDC_INJ_RATE)
Definition: mdc.hh:279
static double GetCentralTime(wavearray< double > x)
Definition: mdc.cc:2761
TString("c")
cout<< endl;cout<< "ts size = "<< ts.size()<< " ts rate = "<< ts.rate()<< endl;tf.Forward(ts, wdm);int levels=tf.getLevel();cout<< "tf size = "<< tf.size()<< endl;double dF=tf.resolution();double dT=1./(2 *dF);cout<< "rate(hz) : "<< RATE<< "\t layers : "<< nLAYERS<< "\t dF(hz) : "<< dF<< "\t dT(ms) : "<< dT *1000.<< endl;int itime=TIME_PIXEL_INDEX;int ifreq=FREQ_PIXEL_INDEX;int index=(levels+1)*itime+ifreq;double time=itime *dT;double freq=(ifreq >0)?ifreq *dF:dF/4;cout<< endl;cout<< "PIXEL TIME = "<< time<< " sec "<< endl;cout<< "PIXEL FREQ = "<< freq<< " Hz "<< endl;cout<< endl;wavearray< double > x
TString mdc[4]
wavearray< double > hp
Definition: mdc.hh:196
waveform wf
virtual void start(double s)
Definition: wavearray.hh:119
watplot * DrawTime(wavearray< double > &x, TString options="ALP", Color_t color=kBlack)
Definition: mdc.cc:2381
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))
void DrawSGE()
Definition: DrawSGE.C:3
bool log
Definition: WaveMDC.C:41
Definition: mdc.hh:122
void SetInjJitter(double inj_jitter=MDC_INJ_JITTER)
Definition: mdc.hh:285
vector< mdcpar > par
double e
wavearray< double > hx
Definition: mdc.hh:197
#define GPS
waveform GetWaveform(int ID, int id=0)
Definition: mdc.cc:2512
CWB::mdc * MDC
Definition: DrawSGE.C:1
char wf_name[256]
void SetInjHrss(double inj_hrss=MDC_INJ_HRSS)
Definition: mdc.hh:275
static double GetCentralFrequency(wavearray< double > x)
Definition: mdc.cc:2803
#define N_IFO