Logo Coherent WaveBurst  
Reference Guide
Logo
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
mdc.hh
Go to the documentation of this file.
1 /**********************************************************
2  * Package: mdc Class Library
3  * File name: mdc.hh
4  * Author: Gabriele Vedovato (vedovato@lnl.infn.it)
5  **********************************************************/
6 
7 
8 #ifndef MDC_HH
9 #define MDC_HH
10 
11 #include "TObjString.h"
12 #include "TObjArray.h"
13 #include "TString.h"
14 #include "TROOT.h"
15 #include "TMath.h"
16 #include "TSystem.h"
17 #include "TLatex.h"
18 #include "TTree.h"
19 #include "TFile.h"
20 #include "TNamed.h"
21 #include "TGlobal.h"
22 #include "TRandom.h"
23 #include "TComplex.h"
24 #include "TCanvas.h"
25 #include "TGraph.h"
26 #include "TMultiGraph.h"
27 #include "TColor.h"
28 #include "TF1.h"
29 #include "TF3.h"
30 #include "TMacro.h"
31 #include "Math/VectorUtil.h"
32 
33 #include <string>
34 #include <iostream>
35 #include <fstream>
36 #include <stdlib.h>
37 #include <math.h>
38 #include <ctype.h>
39 #include <vector>
40 
41 #include "wavecomplex.hh"
42 #include "wavearray.hh"
43 #include "watplot.hh"
44 #include "skymap.hh"
45 #include "network.hh"
46 #include "detector.hh"
47 #include "wat.hh"
48 
49 #include "History.hh"
50 #include "Toolbox.hh"
51 #include "config.hh"
52 #include "time.hh"
53 #include "gskymap.hh"
54 #include "STFT.hh"
55 
56 #include "FrameL.h"
57 
58 #ifndef __CINT__
59 #ifdef _USE_LAL
60 
61 #define LAL_USE_OLD_COMPLEX_STRUCTS
62 #define restrict // this is to getrid of restrict which is defined only in c99
63 
64 #include <lal/LIGOLwXMLInspiralRead.h>
65 #include <lal/LALConfig.h>
66 #include <lal/LALStdio.h>
67 #include <lal/LALStdlib.h>
68 #include <lal/LALError.h>
69 #include <lal/LALDatatypes.h>
70 #include <lal/LIGOMetadataUtils.h>
71 #include <lal/LIGOMetadataTables.h>
72 #include <lal/AVFactories.h>
73 #include <lal/NRWaveIO.h>
74 #include <lal/NRWaveInject.h>
75 #include <lal/Inject.h>
76 #include <lal/FileIO.h>
77 #include <lal/Units.h>
78 #include <lal/FrequencySeries.h>
79 #include <lal/TimeSeries.h>
80 #include <lal/TimeFreqFFT.h>
81 #include <lal/VectorOps.h>
82 #include <lal/LALDetectors.h>
83 #include <lal/FindChirp.h>
84 #include <lal/Random.h>
85 #include <lal/LALNoiseModels.h>
86 #include <lal/Date.h>
87 #include <lal/LALStatusMacros.h>
88 #include <lal/LALSimulation.h>
89 #include <lal/LALInspiral.h>
90 #include <lal/LALSimInspiralWaveformFlags.h>
91 #if LAL_VERSION_MAJOR > 6 || (LAL_VERSION_MAJOR == 6 && \
92  (LAL_VERSION_MINOR > 14 || (LAL_VERSION_MINOR == 14 && \
93  LAL_VERSION_MICRO >= 0 ))) // LAL_VERSION >= 6.14.0
94 #include <lal/SnglBurstUtils.h>
95 #else
96 #include <lal/LIGOMetadataBurstUtils.h>
97 #endif
98 #include <lal/LIGOLwXMLBurstRead.h>
99 #include <lal/GenerateBurst.h>
100 #include <lal/LIGOLwXML.h>
101 
102 #endif
103 #endif
104 
105 #ifdef _USE_EBBH
106 #include "eBBH.hh"
107 #endif
108 
109 #define MDC_INJ_RATE 0.01 // sec^-1
110 #define MDC_INJ_JITTER 10. // sec
111 #define MDC_INJ_LENGTH 1. // sec
112 #define MDC_INJ_HRSS 2.5e-21
113 
114 #define MDC_SAMPLE_RATE 16384. // sample/sec
115 
116 #define EPZOOM 0.9999999
117 
118 enum MDC_TYPE {
119  MDC_SG = 0,
120  MDC_SGL = 0,
121  MDC_SGC = 1,
122  MDC_SGE = 2,
123  MDC_RD = 3,
124  MDC_RDL = 3,
125  MDC_RDC = 4,
126  MDC_RDE = 5,
127  MDC_WNB = 6,
128  MDC_GA = 7,
129  MDC_EBBH = 8,
133  MDC_CG = 12,
134  MDC_CGL = 12,
135  MDC_CGC = 13,
136  MDC_CGE = 14,
138  MDC_USER = 16 // must be declared as the last one
139 };
140 
144 };
145 
147  MDC_RANDOM = -1000,
150  MDC_MNGD = 2,
151  MDC_GWGC = 3,
155 };
156 
158  switch (n) {
159  case MDC_RANDOM: return "MDC_RANDOM";
160  case MDC_EARTH_FIX: return "MDC_EARTH_FIX";
161  case MDC_CELESTIAL_FIX: return "MDC_CELESTIAL_FIX";
162  case MDC_MNGD: return "MDC_MNGD";
163  case MDC_GWGC: return "MDC_GWGC";
164  case MDC_LOGFILE: return "MDC_LOGFILE";
165  case MDC_CUSTOM: return "MDC_CUSTOM";
166  case MDC_XMLFILE: return "MDC_XMLFILE";
167  default: return "MDC_RANDOM";
168  }
169 }
170 
171 enum MDC_DRAW {
172  MDC_TIME = 0,
173  MDC_FFT = 1,
174  MDC_TF = 2
175 };
176 
177 struct mdcid {
178  TString name; // name of waveform
179  int ID; // major id of waveform list
180  int id; // minor id of waveform list (Ex : the list WNB with different random waveforms)
181 };
182 
183 struct mdcpar {
185  double value;
187 };
188 
189 struct waveform { // waveform structure
190  MDC_TYPE type; // mdc type
191  bool status; // status flag used for internal mdc processing
192  TString name; // name of waveform
193  TString hpPath; // path of user define hp file
194  TString hxPath; // path of user define hx file
195  vector<mdcpar> par; // waveform parameters
196  wavearray<double> hp; // waveform hp component
197  wavearray<double> hx; // waveform hx component;
198  vector<waveform> list; // list of waveforms belonging to the same waveform class
199 };
200 
201 struct source {
202  double gps; // waveform source structure
203  double theta; // latitude (degrees)
204  double phi; // longitude (degrees)
205  double psi; // polarization (degrees)
206  double rho; // distance
207  double iota; // elliptical inclination angle (degrees)
208  double hrss; // sqrt(hp^2+hx^2)
209  int ID; // major id of waveform list
210  int id; // minor id of waveform list (Ex : the list WNB with different random waveforms)
211  waveform wf; // waveform
212 };
213 
214 namespace CWB {
215 
216 class mdc : public TNamed {
217 
218 public:
219 
220  mdc();
221  mdc(int nIFO, TString* ifo);
222  mdc(int nIFO, detector** pD);
223  mdc(network* net);
224  mdc(const CWB::mdc& value);
225  ~mdc();
226 
227  // operators
228  mdc& operator = (const mdc&);
229 
231 
232  mdcid AddWaveform(MDC_TYPE mdc_type, vector<mdcpar> par, TString uname=""); // built-in waveform
233  void AddWaveform(TString mdc_name, TString hp_fName, TString hx_fName,
234  double srate, vector<mdcpar> par=vector<mdcpar>());
235  void AddWaveform(TString mdc_name, TString hp_fName,
236  double srate, vector<mdcpar> par=vector<mdcpar>());
237  void AddWaveform(TString mdc_name, TString hp_fName, TString hx_fName);
238  void AddWaveform(TString mdc_name, TString hp_fName);
240 #ifdef _USE_LAL
241 #ifndef __CINT__
242  // LAL SimBurst waveforms
243  mdcid AddWaveform(MDC_TYPE mdc_type, SimBurst* sim_burst, vector<mdcpar> par, TString uname="");
244 #endif
245 #endif
246 
247  // read ascii waveform with fixed sample rate (1 column [h])
248  void ReadWaveform(wavearray<double>& x, TString fName, double srate);
249  // read ascii waveform with variable sample rate (2 columns [t,h])
250  void ReadWaveform(wavearray<double>& x, TString fName);
251 
252  // get size of wfList
253  inline size_t wfListSize() { return wfList.size(); }
254 
257 
258  void SetSkyDistribution(MDC_DISTRIBUTION sky_distribution, vector<mdcpar> par, int seed=0, bool add=false);
259  void SetSkyDistribution(MDC_DISTRIBUTION sky_distribution, TString fName, vector<mdcpar> par, int seed=0, bool add=false);
260  void SetSkyDistribution(MDC_DISTRIBUTION sky_distribution, TString fName, int seed, bool add=false);
263  vector<mdcpar> GetSkyParms() {return sky_parms;}
264  void DrawSkyDistribution(TString name = "skymap", TString projection = "",
265  TString coordinate = "Geographic", double resolution = 2, bool background=true);
266 
267  void GetSourceCoordinates(double& theta, double& phi, double& psi, double& rho,
268  double& iota, double& hrss, int& ID, int& id);
269  void GetSourceCoordinates(double gps, double& theta, double& phi, double& psi, double& rho,
270  double& iota, double& hrss, int& ID, int& id);
271 
274  double GetInjLength() {return inj_length;}
276  // if inj_hrss=0 the hrss @10Kpc is the one which is defined by hp,hc waveforms
277  this->inj_hrss = inj_hrss>=0 ? inj_hrss : MDC_INJ_HRSS;}
278  double GetInjHrss() {return inj_hrss;}
280  this->inj_rate = inj_rate>0 ? inj_rate : MDC_INJ_RATE;}
281  double GetInjRate() {return inj_rate;}
282  void SetInjOffset(double inj_offset=0) {
283  this->inj_offset = inj_offset>0 ? inj_offset : 0;}
284  double GetInjOffset() {return inj_offset;}
287  double GetInjJitter() {return inj_jitter;}
288  void SetSourceListSeed(unsigned int srcList_seed) {
289  this->srcList_seed = srcList_seed;}
290  double GetSourceListSeed() {return srcList_seed;}
291  double GetSampleRate() {return MDC_SAMPLE_RATE;}
292 
293 #ifdef _USE_LAL
295  TString GetInspiral() {return inspOptions;}
296  wavearray<double> GetInspiral(TString pol, int gps_start_time, int gps_end_time);
297  TString GetInspiralOption(TString inspOption);
298 #endif
299 
300  void Dump(TString fname, int ID, int id, TString polarization);
303  void Dump(TString dname); // dump all waveforms in dname dir
304 
305  watplot* Draw(TString name, int id=0, TString polarization="hp", MDC_DRAW type=MDC_TIME,
306  TString options = "ALP", Color_t color=kBlack);
307  watplot* Draw(int ID, int id=0, TString polarization="hp", MDC_DRAW type=MDC_TIME,
308  TString options = "ALP", Color_t color=kBlack);
310  TString options = "ALP", Color_t color=kBlack);
311  watplot* Draw(TString ifo, double gpsStart, double gpsEnd, int id,
312  MDC_DRAW type=MDC_TIME, TString options = "ALP", Color_t color=kBlack);
313 
314  void DumpLog(TString fName, TString label="", bool append=false);
315  void DumpLogHeader(TString fName, TString label="", int size=0);
316 
317  waveform GetWaveform(int ID, int id=0);
318  waveform GetWaveform(TString name, int id=0);
320  void GetWaveform(waveform& wf);
321 
322  static double GetCentralTime(wavearray<double> x);
323  static double GetCentralTime(waveform wf);
324  static double GetCentralFrequency(wavearray<double> x);
325  static double GetCentralFrequency(waveform wf);
326  static double GetTimeRange(wavearray<double> x, double& tMin, double& tMax, double efraction = EPZOOM);
327  static void TimeShift(wavearray<double>& x, double tShift=0.);
328  static void PhaseShift(wavearray<double>& x, double pShift=0.);
329 
330  double GetAntennaPattern(TString ifo, double phi, double theta, double psi=0., TString polarization="hp");
331  double GetDelay(TString ifo1, TString ifo2, double phi, double theta);
332 
333  void Print(int level=0); // *MENU*
334  virtual void Browse(TBrowser *b) {Print();}
335 
336  wavearray<double> GetWNB(double frequency, double bandwidth, double duration, int seed=0, bool mode=0);
337  wavearray<double> GetRD(double frequency, double tau, double iota, bool polarization=0);
338  wavearray<double> GetSGQ(double frequency, double Q);
339  wavearray<double> GetCGQ(double frequency, double Q);
341 
342  static void AddGauss(wavearray<double> &td, double v, double u=0.);
343  static void AddExp(wavearray<double> &td, double v, int M);
344  static void AddSGBurst(wavearray<double> &td, double a, double f, double s, double d=0.);
345  static void AddCGBurst(wavearray<double> &td, double a, double f, double s, double d=0.);
346  static void AddWGNoise(wavearray<double> &td, double a, double s);
347 
349  size_t length=1000, bool log=false, vector<TString> chName=vector<TString>());
350 
351  network* GetNetwork() {return net;}
352 
353  vector<waveform> GetWaveformList() {return wfList;}
354 
355  vector<waveform> wfList;
356 
357  std::vector<std::string> mdcList; // list of injections
358  std::vector<std::string> mdcType; // list of injection types
359  std::vector<std::string> xmlType; // list of xml injection types (used with MDC_XMLFILE)
360  std::vector<double> mdcTime; // gps time of selected injections
361  std::vector<std::string> mdcName; // list of injection names used with MDC_LOG
362  std::vector<source> srcList; // list of source MDC parameters
363 
364  double GetPar(TString name, vector<mdcpar> par, bool& error);
365  TString GetParString(TString name, vector<mdcpar> par, bool& error);
366 
367  waveform GetSourceWaveform(int& ID, int& id);
368  vector<source> GetSourceList(double start, double stop);
369  TString GetBurstLog(source src, double FrameGPS, double SimHpHp, double SimHcHc, double SimHpHc);
370 
371  watplot* DrawTime(wavearray<double>& x, TString options = "ALP", Color_t color=kBlack);
372  watplot* DrawFFT(wavearray<double>& x, TString options = "ALP", Color_t color=kBlack);
373  void DrawTF(wavearray<double>& x, TString options = "");
374 
376 
377 #ifdef _USE_LAL
378  TString GetInspiral(wavearray<double>& x, TString ifo);
379  TString GenInspiralXML(int gps_start_time, int gps_end_time, bool rmFile=false);
380  TString GetInspName() {return inspName;}
381  TString GetWaveName() {return waveName;}
382 #ifndef __CINT__
383  TString GetInspiralLog(TString inspName, double FrameGPS, SimInspiralTable *thisInj);
384  void FindChirpInjectSignals(wavearray<double>& w, SimInspiralTable *injections, TString ifo);
385 #endif
386  void CreateBurstXML(TString fName, vector<mdcpar> xml_parms);
387  void FillBurstXML(bool precision=false);
388  void CloseBurstXML();
389  static TString GetBurstNameLAL(TString options);
390 #endif
391  static double e2cosi(double e);
392  static double cosi2e(double cosi);
393 
394  TString GetTemporaryFileName(TString tag="mdc", TString ext="txt", TString dir="/tmp", bool mkdir=false);
395 
396  void SetZoom(double epzoom=EPZOOM) {if(epzoom<0||epzoom>1) epzoom=EPZOOM; this->epzoom=epzoom;}
397  double GetZoom() {return epzoom;}
398 
399  gskymap* GetGSkyMap() {return psp;}
400  CWB::STFT* GetSTFT() {return stft;}
401  watplot* GetWatPlot() {return pts;}
402 
403 private:
404 
405  void Init(int seed=0);
406  void exit(int err) {gSystem->Exit(err);} // overwrite exit system function
407 
408  network* net; //!
410  TTree* inj_tree; //!
411 
416  vector<mdcpar> sky_parms;
417 
418  std::vector<std::string> nameList;
419  std::vector<float> thList;
420  std::vector<float> phList;
421  std::vector<float> psiList;
422  std::vector<float> rhoList;
423  std::vector<float> iotaList;
424  std::vector<double> hrssList;
425  std::vector<double> gpsList;
426  std::vector<int> IDList;
427  std::vector<int> idList;
428 
429  double inj_rate;
430  double inj_offset;
431  double inj_jitter;
432  double inj_hrss;
433  double inj_length; //length sec
434  unsigned int srcList_seed; // value added to the seed used in GetSourceList (default is 0)
435 
436  // inspiral parameters
441  TString waveName; //- PN approximant to be used in computing the waveform
442  // (see enum Approximant in LALSimInspiral.h)
443 
445 
446  double epzoom; // used to zoom plots : energy percentage displayed
447 
448  gskymap* psp; //!
450  watplot* pts; //!
451 
452 #ifdef _USE_LAL
453 #ifndef __CINT__
454  // used by the CreateBurstXML,FillBurstXML,CloseBurstXML functions
455  ProcessTable *xml_process_table_head; //!
456  ProcessTable *xml_process; //!
457  ProcessParamsTable *xml_process_params_table_head; //!
458  SearchSummaryTable *xml_search_summary_table_head; //!
459  SearchSummaryTable *xml_search_summary; //!
460  TimeSlide *xml_time_slide_table_head; //!
461  SimBurst *xml_sim_burst_table_head; //!
462  SimBurst **xml_sim_burst; //!
463 #endif
464 #endif
465 
466  ClassDef(mdc,7)
467 };
468 
469 } // end namespace
470 
471 #endif
MDC_DISTRIBUTION
Definition: mdc.hh:146
MDC_DRAW
Definition: mdc.hh:171
void Dump(TString fname, int ID, int id, TString polarization)
Definition: mdc.cc:4365
int id
Definition: mdc.hh:180
double rho
#define EPZOOM
Definition: mdc.hh:116
int ID
Definition: mdc.hh:209
waveform GetSourceWaveform(int &ID, int &id)
Definition: mdc.cc:2011
TString sky_file
Definition: mdc.hh:414
Definition: mdc.hh:172
double GetInjHrss()
Definition: mdc.hh:278
#define MDC_INJ_RATE
Definition: mdc.hh:109
TString inspDIR
Definition: mdc.hh:438
vector< mdcpar > GetSkyParms()
Definition: mdc.hh:263
watplot * Draw(TString name, int id=0, TString polarization="hp", MDC_DRAW type=MDC_TIME, TString options="ALP", Color_t color=kBlack)
Definition: mdc.cc:2288
double duration
par[0] value
Definition: mdc.hh:189
bool status
Definition: mdc.hh:191
vector< mdcpar > sky_parms
Definition: mdc.hh:416
tuple f
Definition: cwb_online.py:91
static double GetTimeRange(wavearray< double > x, double &tMin, double &tMax, double efraction=EPZOOM)
Definition: mdc.cc:2831
std::vector< float > psiList
Definition: mdc.hh:421
TString Get(wavearray< double > &x, TString ifo)
Definition: mdc.cc:1502
TString inspOptions
Definition: mdc.hh:440
Definition: ced.hh:24
std::vector< std::string > mdcList
Definition: mdc.hh:357
#define MDC_INJ_LENGTH
Definition: mdc.hh:111
void SetInjRate(double inj_rate=MDC_INJ_RATE)
Definition: mdc.hh:279
wavearray< double > a(hp.size())
static double GetCentralTime(wavearray< double > x)
Definition: mdc.cc:2761
WSeries< float > v[nIFO]
Definition: cwb_net.C:62
Definition: mdc.hh:128
Definition: mdc.hh:177
par[0] name
int n
Definition: cwb_net.C:10
CWB::STFT * GetSTFT()
Definition: mdc.hh:400
Definition: mdc.hh:173
void SetSourceListSeed(unsigned int srcList_seed)
Definition: mdc.hh:288
Definition: mdc.hh:151
TString("c")
std::vector< std::string > mdcType
Definition: mdc.hh:358
wavearray< double > GetSGQ(double frequency, double Q)
Definition: mdc.cc:2976
int ID
Definition: TestMDC.C:70
static void PhaseShift(wavearray< double > &x, double pShift=0.)
Definition: mdc.cc:2926
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
double frequency
unsigned int srcList_seed
Definition: mdc.hh:434
void GetSourceCoordinates(double &theta, double &phi, double &psi, double &rho, double &iota, double &hrss, int &ID, int &id)
Definition: mdc.cc:1960
double GetInjRate()
Definition: mdc.hh:281
float theta
network * net
Definition: mdc.hh:408
std::vector< double > gpsList
Definition: mdc.hh:425
void SetInjLength(double inj_length=MDC_INJ_LENGTH)
Definition: mdc.hh:272
TString waveName
Definition: mdc.hh:441
static void AddExp(wavearray< double > &td, double v, int M)
Definition: mdc.cc:3221
MDC_TYPE
Definition: mdc.hh:118
double epzoom
Definition: mdc.hh:446
wavearray< double > GetGA(double duration)
Definition: mdc.cc:3101
std::vector< float > phList
Definition: mdc.hh:420
wavearray< double > hp
Definition: mdc.hh:196
double inj_rate
Definition: mdc.hh:429
waveform wf
int polarization
std::vector< int > IDList
Definition: mdc.hh:426
vector< waveform > GetWaveformList()
Definition: mdc.hh:353
Long_t size
double GetZoom()
Definition: mdc.hh:397
#define M
Definition: UniqSLagsList.C:3
TString GetSkyFile()
Definition: mdc.hh:262
TString name
Definition: mdc.hh:184
double GetDelay(TString ifo1, TString ifo2, double phi, double theta)
Definition: mdc.cc:4494
watplot * pts
Definition: mdc.hh:450
watplot * DrawTime(wavearray< double > &x, TString options="ALP", Color_t color=kBlack)
Definition: mdc.cc:2381
mdc & operator=(const mdc &)
Definition: mdc.cc:280
mdcid AddWaveform(MDC_TYPE mdc_type, vector< mdcpar > par, TString uname="")
Definition: mdc.cc:445
void Print(int level=0)
Definition: mdc.cc:2707
Definition: mdc.hh:133
virtual void Browse(TBrowser *b)
Definition: mdc.hh:334
double value
Definition: mdc.hh:185
char ifo[NIFO_MAX][8]
std::vector< std::string > nameList
Definition: mdc.hh:418
MDC_COORDINATES
Definition: mdc.hh:141
double gps
Definition: mdc.hh:202
double GetSampleRate()
Definition: mdc.hh:291
void SetSkyDistribution(MDC_DISTRIBUTION sky_distribution, vector< mdcpar > par, int seed=0, bool add=false)
Definition: mdc.cc:3415
TString GetTemporaryFileName(TString tag="mdc", TString ext="txt", TString dir="/tmp", bool mkdir=false)
Definition: mdc.cc:5661
Definition: mdc.hh:126
double GetInjJitter()
Definition: mdc.hh:287
#define MDC_INJ_JITTER
Definition: mdc.hh:110
size_t mode
double GetPar(TString name, vector< mdcpar > par, bool &error)
Definition: mdc.cc:405
Definition: mdc.hh:174
int GetWaveformID(TString name)
Definition: mdc.cc:2689
TString inspXML
Definition: mdc.hh:437
wavearray< double > w
Definition: Test1.C:27
nc append(pix)
#define nIFO
std::vector< std::string > xmlType
Definition: mdc.hh:359
void DrawSkyDistribution(TString name="skymap", TString projection="", TString coordinate="Geographic", double resolution=2, bool background=true)
Definition: mdc.cc:4135
wavearray< double > GetCGQ(double frequency, double Q)
Definition: mdc.cc:2999
float phi
TString chName[NIFO_MAX]
size_t wfListSize()
Definition: mdc.hh:253
Definition: mdc.hh:125
static void AddGauss(wavearray< double > &td, double v, double u=0.)
Definition: mdc.cc:3201
Definition: mdc.hh:138
double rho
Definition: mdc.hh:206
TString name
Definition: mdc.hh:192
double hrss
Definition: TestMDC.C:70
float psi
void SetInjOffset(double inj_offset=0)
Definition: mdc.hh:282
TString hxPath
Definition: mdc.hh:194
CWB::STFT * stft
Definition: mdc.hh:449
TString uname
int ID
Definition: mdc.hh:179
Definition: mdc.hh:216
TString WriteFrameFile(TString frDir, TString frLabel, size_t gps, size_t length=1000, bool log=false, vector< TString > chName=vector< TString >())
Definition: mdc.cc:4254
Definition: mdc.hh:123
gskymap * GetGSkyMap()
Definition: mdc.hh:399
void DumpLogHeader(TString fName, TString label="", int size=0)
Definition: mdc.cc:4538
#define MDC_SAMPLE_RATE
Definition: mdc.hh:114
double hrss
Definition: mdc.hh:208
int id
Definition: mdc.hh:210
TString label
Definition: MergeTrees.C:21
bool log
Definition: WaveMDC.C:41
vector< mdcpar > par
Definition: mdc.hh:195
cout<< "Injected signals: "<< mdc.GetEntries()<< endl;cout<< "Injected signals in histogram factor_events_inj: "<< NEVTS<< endl;float myifar, ecor, m1, m2, netcc[3], neted, penalty;float rho[2];float chirp[6];float range[2];float frequency[2];float iSNR[3], sSNR[3];sim.SetBranchAddress("mass", mass);sim.SetBranchAddress("factor",&factor);sim.SetBranchAddress("range", range);sim.SetBranchAddress("chirp", chirp);sim.SetBranchAddress("rho", rho);sim.SetBranchAddress("netcc", netcc);sim.SetBranchAddress("neted",&neted);sim.SetBranchAddress("ifar",&myifar);sim.SetBranchAddress("ecor",&ecor);sim.SetBranchAddress("penalty",&penalty);sim.SetBranchAddress("time", mytime);sim.SetBranchAddress("iSNR", iSNR);sim.SetBranchAddress("sSNR", sSNR);sim.SetBranchAddress("spin", spin);sim.SetBranchAddress("frequency", frequency);float **volume=new float *[NBINS_mass1];float **volume_first_shell=new float *[NBINS_mass1];float **radius=new float *[NBINS_mass1];float **error_volume=new float *[NBINS_mass1];float **error_volume_first_shell=new float *[NBINS_mass1];float **error_radius=new float *[NBINS_mass1];for(int i=0;i< NBINS_mass1;i++){volume[i]=new float[NBINS_mass2];volume_first_shell[i]=new float[NBINS_mass2];radius[i]=new float[NBINS_mass2];error_volume[i]=new float[NBINS_mass2];error_volume_first_shell[i]=new float[NBINS_mass2];error_radius[i]=new float[NBINS_mass2];for(int j=0;j< NBINS_mass2;j++){volume[i][j]=0.;volume_first_shell[i][j]=0.;radius[i][j]=0.;error_volume[i][j]=0.;error_volume_first_shell[i][j]=0.;error_radius[i][j]=0.;}}float **spin_mtot_volume=new float *[NBINS_MTOT+1];float **spin_mtot_radius=new float *[NBINS_MTOT+1];float **error_spin_mtot_volume=new float *[NBINS_MTOT+1];float **error_spin_mtot_radius=new float *[NBINS_MTOT+1];for(int i=0;i< NBINS_MTOT+1;i++){spin_mtot_volume[i]=new float[NBINS_SPIN+1];spin_mtot_radius[i]=new float[NBINS_SPIN+1];error_spin_mtot_volume[i]=new float[NBINS_SPIN+1];error_spin_mtot_radius[i]=new float[NBINS_SPIN+1];for(int j=0;j< NBINS_SPIN+1;j++){spin_mtot_volume[i][j]=0.;error_spin_mtot_volume[i][j]=0.;spin_mtot_radius[i][j]=0.;error_spin_mtot_radius[i][j]=0.;}}char fname[1024];sprintf(fname,"%s/recovered_signals.txt", netdir);ofstream fev;fev.open(fname, std::ofstream::out);sprintf(line,"#GPS@L1 FAR[Hz] eFAR[Hz] Pval ""ePval factor rho frequency iSNR sSNR \n");fev<< line<< endl;ofstream *fev_single=new ofstream[nfactor];for(int l=1;l< nfactor+1;l++){sprintf(fname,"%s/recovered_signals_%d.txt", netdir, l);fev_single[l-1].open(fname, std::ofstream::out);fev_single[l-1]<< line<< endl;}double Vrho[RHO_NBINS], eVrho[RHO_NBINS], Rrho[RHO_NBINS], eRrho[RHO_NBINS], Trho[RHO_NBINS];for(int i=0;i< RHO_NBINS;i++){Vrho[i]=0.;eVrho[i]=0.;Rrho[i]=0.;eRrho[i]=0.;Trho[i]=RHO_MIN+i *RHO_BIN;}double dV, dV1, dV_spin_mtot, nevts, internal_volume;int nT;int countv=0;int cnt=0;int cnt2=0;int cntfreq=0;bool bcut=false;double liveTot=sim.GetMaximum("ifar");double BKG_LIVETIME_yr=liveTot/CYS;double BKG_LIVETIME_Myr=BKG_LIVETIME_yr/(1.e6);cout.precision(14);cout<< "Total live time ---> background
const char * DistributionToString(MDC_DISTRIBUTION n)
Definition: mdc.hh:157
TString frDir[NIFO_MAX]
TString GetBurst(wavearray< double > &x, TString ifo)
Definition: mdc.cc:1529
waveform wf
Definition: mdc.hh:211
static double tau
Definition: geodesics.cc:8
char fname[1024]
Definition: mdc.hh:183
wavearray< double > GetWNB(double frequency, double bandwidth, double duration, int seed=0, bool mode=0)
Definition: mdc.cc:3024
double precision
Definition: mdc.hh:120
Definition: mdc.hh:121
double GetInjOffset()
Definition: mdc.hh:284
MDC_COORDINATES GetCoordinatesSystem()
Definition: mdc.hh:256
Definition: mdc.hh:122
r add(tfmap,"hchannel")
void SetInjJitter(double inj_jitter=MDC_INJ_JITTER)
Definition: mdc.hh:285
Definition: skymap.hh:45
std::vector< int > idList
Definition: mdc.hh:427
vector< mdcpar > par
TString inspName
Definition: mdc.hh:439
double e
void DrawTF(wavearray< double > &x, TString options="")
Definition: mdc.cc:2473
watplot * DrawFFT(wavearray< double > &x, TString options="ALP", Color_t color=kBlack)
Definition: mdc.cc:2427
Definition: mdc.hh:134
TString svalue
Definition: mdc.hh:186
vector< waveform > wfList
Definition: mdc.hh:355
char tag[256]
Definition: cwb_merge.C:74
double inj_jitter
Definition: mdc.hh:431
static double cosi2e(double cosi)
Definition: mdc.cc:5011
static void AddCGBurst(wavearray< double > &td, double a, double f, double s, double d=0.)
Definition: mdc.cc:3292
double inj_offset
Definition: mdc.hh:430
double iota
Definition: mdc.hh:207
TString GetBurstLog(source src, double FrameGPS, double SimHpHp, double SimHcHc, double SimHpHc)
Definition: mdc.cc:2151
network * GetNetwork()
Definition: mdc.hh:351
static void TimeShift(wavearray< double > &x, double tShift=0.)
Definition: mdc.cc:2874
double GetSourceListSeed()
Definition: mdc.hh:290
s s
Definition: cwb_net.C:137
injection * inj
Definition: mdc.hh:409
MDC_COORDINATES mdc_coordinates
Definition: mdc.hh:412
static void AddWGNoise(wavearray< double > &td, double a, double s)
Definition: mdc.cc:3334
char options[256]
double inj_length
Definition: mdc.hh:433
Definition: mdc.hh:119
Definition: mdc.hh:150
mdc()
Definition: mdc.cc:178
wavearray< double > hx
Definition: mdc.hh:197
skymap sm
Definition: mdc.hh:444
double psi
Definition: mdc.hh:205
double gps
std::vector< std::string > mdcName
Definition: mdc.hh:361
watplot * GetWatPlot()
Definition: mdc.hh:401
TString xml_filename
Definition: mdc.hh:415
MDC_DISTRIBUTION sky_distribution
Definition: mdc.hh:413
std::vector< float > rhoList
Definition: mdc.hh:422
vector< source > GetSourceList(double start, double stop)
Definition: mdc.cc:2026
MDC_DISTRIBUTION GetSkyDistribution()
Definition: mdc.hh:261
double inj_hrss
Definition: mdc.hh:432
double Q
double theta
Definition: mdc.hh:203
waveform GetWaveform(int ID, int id=0)
Definition: mdc.cc:2512
void ReadWaveform(wavearray< double > &x, TString fName, double srate)
Definition: mdc.cc:1848
static double e2cosi(double e)
Definition: mdc.cc:4971
void SetZoom(double epzoom=EPZOOM)
Definition: mdc.hh:396
~mdc()
Definition: mdc.cc:260
Definition: mdc.hh:124
TString GetParString(TString name, vector< mdcpar > par, bool &error)
Definition: mdc.cc:425
std::vector< double > hrssList
Definition: mdc.hh:424
std::vector< float > iotaList
Definition: mdc.hh:423
double GetInjLength()
Definition: mdc.hh:274
void DumpLog(TString fName, TString label="", bool append=false)
Definition: mdc.cc:5039
void Init(int seed=0)
Definition: mdc.cc:354
void SetInjHrss(double inj_hrss=MDC_INJ_HRSS)
Definition: mdc.hh:275
std::vector< source > srcList
Definition: mdc.hh:362
static double GetCentralFrequency(wavearray< double > x)
Definition: mdc.cc:2803
double GetAntennaPattern(TString ifo, double phi, double theta, double psi=0., TString polarization="hp")
Definition: mdc.cc:4466
std::vector< double > mdcTime
Definition: mdc.hh:360
double phi
Definition: mdc.hh:204
Definition: mdc.hh:201
gskymap * psp
Definition: mdc.hh:448
wavearray< double > GetRD(double frequency, double tau, double iota, bool polarization=0)
Definition: mdc.cc:3139
#define MDC_INJ_HRSS
Definition: mdc.hh:112
char fName[256]
void exit(int err)
Definition: mdc.hh:406
Definition: mdc.hh:127
double length
Definition: TestBandPass.C:18
void SetCoordinatesSystem(MDC_COORDINATES mdc_coordinates=MDC_EARTH)
Definition: mdc.hh:255
vector< waveform > list
Definition: mdc.hh:198
MDC_TYPE type
Definition: mdc.hh:190
TString name
Definition: mdc.hh:178
TString frLabel[NIFO_MAX]
Definition: mdc.hh:136
TTree * inj_tree
Definition: mdc.hh:410
TString hpPath
Definition: mdc.hh:193
Definition: mdc.hh:129
Definition: mdc.hh:135
detector ** pD
static void AddSGBurst(wavearray< double > &td, double a, double f, double s, double d=0.)
Definition: mdc.cc:3250
std::vector< float > thList
Definition: mdc.hh:419
MDC SetInspiral("inspNameTEST", inspOptions)