Logo Coherent WaveBurst  
Reference Guide
Logo
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
CWB_Plugin_SkyProb.C
Go to the documentation of this file.
1 #define XIFO 4
2 
3 #pragma GCC system_header
4 
5 #include "cwb.hh"
6 #include "config.hh"
7 #include "network.hh"
8 #include "wavearray.hh"
9 #include "TString.h"
10 #include "TObjArray.h"
11 #include "TObjString.h"
12 #include "TRandom.h"
13 #include "TComplex.h"
14 #include "TMath.h"
15 #include "mdc.hh"
16 #include "watplot.hh"
17 #include <vector>
18 
19 //#define WATPLOT // enable event monster plots
20 //#define MDC_OTF // enable MDC On The Fly
21 //#define NOISE_OTF // enable NOISE On The Fly
22 
23 void
25 //!NOISE_MDC_SIMULATION
26 // Plugin to save SkyProb into the output wave root file
27 
28  cout << endl;
29  cout << "-----> CWB_Plugin_SkyProb.C" << endl;
30  cout << "ifo " << ifo.Data() << endl;
31  cout << "type " << type << endl;
32  cout << endl;
33 
34  if(type==CWB_PLUGIN_CONFIG) {
35 #ifdef NOISE_OTF
36  cfg->dataPlugin=true; // disable read data from frames
37 #endif
38 #ifdef MDC_OTF
39  cfg->mdcPlugin=true; // disable read mdc from frames
40 #endif
41  cfg->outPlugin=true; // disable built-in output root file
42  cfg->Psave = 0; // disable save Skymap probability to root file
43  }
44 
45 
46 #ifdef NOISE_OTF
47  if(type==CWB_PLUGIN_DATA) {
48 
50 
51  int seed;
52  if(ifo.CompareTo("L1")==0) seed=1000;
53  if(ifo.CompareTo("H1")==0) seed=2000;
54  if(ifo.CompareTo("V1")==0) seed=3000;
55  if(ifo.CompareTo("J1")==0) seed=4000;
56  if(ifo.CompareTo("A2")==0) seed=5000;
57  if(ifo.CompareTo("Y2")==0) seed=6000;
58  if(ifo.CompareTo("Y3")==0) seed=7000;
59 
60  TString fName;
61  if(ifo.CompareTo("L1")==0) fName="plugins/strains/advLIGO_NSNS_Opt_8khz_one_side.txt";
62  if(ifo.CompareTo("H1")==0) fName="plugins/strains/advLIGO_NSNS_Opt_8khz_one_side.txt";
63  if(ifo.CompareTo("V1")==0) fName="plugins/strains/advVIRGO_sensitivity_12May09_8khz_one_side.txt";
64  if(ifo.CompareTo("J1")==0) fName="plugins/strains/LCGT_sensitivity_8khz_one_side.txt";
65  if(ifo.CompareTo("A2")==0) fName="plugins/strains/advLIGO_NSNS_Opt_8khz_one_side.txt";
66  if(ifo.CompareTo("Y2")==0) fName="plugins/strains/advLIGO_NSNS_Opt_8khz_one_side.txt";
67  if(ifo.CompareTo("Y3")==0) fName="plugins/strains/advLIGO_NSNS_Opt_8khz_one_side.txt";
68 
69  int size=x->size();
70  double start=x->start();
71  TB.getSimNoise(*x, fName, seed, NET->nRun);
72  x->resize(size);
73  x->start(start);
74  }
75 #endif
76 
77 #ifdef MDC_OTF
78  if(type==CWB_PLUGIN_MDC) {
79 
80  char cmd[128];
81  sprintf(cmd,"network* net = (network*)%p;",NET);
82  gROOT->ProcessLine(cmd);
83 
84  CWB::mdc MDC(NET);
85 
86  // ---------------------------------
87  // read plugin config
88  // ---------------------------------
89 
90  cfg->configPlugin.Exec();
91 
92  // ---------------------------------
93  // set list of mdc waveforms
94  // ---------------------------------
95 
96  IMPORT(CWB::mdc,MDC)
97  MDC.Print();
98 
99  // ---------------------------------
100  // get mdc data
101  // ---------------------------------
102 
103  MDC.Get(*x,ifo);
104 
105  // ---------------------------------
106  // set mdc list in the network class
107  // ---------------------------------
108 
109  if(ifo.CompareTo(NET->ifoName[0])==0) {
110  NET->mdcList.clear();
111  NET->mdcType.clear();
112  NET->mdcTime.clear();
113  NET->mdcList=MDC.mdcList;
114  NET->mdcType=MDC.mdcType;
115  NET->mdcTime=MDC.mdcTime;
116  }
117 
118  cout.precision(14);
119  for(int k=0;k<(int)NET->mdcList.size();k++) cout << k << " mdcList " << MDC.mdcList[k] << endl;
120  for(int k=0;k<(int)NET->mdcTime.size();k++) cout << k << " mdcTime " << MDC.mdcTime[k] << endl;
121  for(int k=0;k<(int)NET->mdcType.size();k++) cout << k << " mdcType " << MDC.mdcType[k] << endl;
122  }
123 #endif
124 
125  if(type==CWB_PLUGIN_OLIKELIHOOD) {
126 
127  if(TString(cfg->analysis)!="2G") {
128  cout << "CWB_Plugin_SkyProb.C -> CWB_PLUGIN_OLIKELIHOOD implemented only for 2G" << endl;
129  gSystem->Exit(1);
130  }
131 
132 #ifdef WATPLOT
133  watplot SMS(const_cast<char*>("sms"),200,20,800,400);
134 #endif
135 
136  // import ifactor
137  int gIFACTOR=-1; IMPORT(int,gIFACTOR)
138  cout << "-----> CWB_Plugin_SkyProb.C -> " << " gIFACTOR : " << gIFACTOR << endl;
139 
140  // import slagShift
141  float* gSLAGSHIFT=NULL; IMPORT(float*,gSLAGSHIFT)
142 
143  int nIFO = NET->ifoListSize(); // number of detectors
144  int K = NET->nLag; // number of time lag
145  netevent* EVT;
147  double factor = cfg->simulation==3||cfg->simulation==4 ? -gIFACTOR : cfg->factors[gIFACTOR];
148  int rate = 0; // select all resolutions
149 
150  // search output root file in the system list
151  TFile* froot = NULL;
152  TList *files = (TList*)gROOT->GetListOfFiles();
153  TString outDump="";
154  if (files) {
155  TIter next(files);
156  TSystemFile *file;
157  TString fname;
158  bool check=false;
159  while ((file=(TSystemFile*)next())) {
160  fname = file->GetName();
161  // set output root file as the current file
162  if(fname.Contains("wave_")) {
163  froot=(TFile*)file;froot->cd();
164  outDump=fname;
165  outDump.ReplaceAll(".root.tmp",".txt");
166  cout << "output file name : " << fname << endl;
167  }
168  }
169  if(!froot) {
170  cout << "CWB_Plugin_SkyProb.C : Error - output root file not found" << endl;
171  gSystem->Exit(1);
172  }
173  } else {
174  cout << "CWB_Plugin_SkyProb.C : Error - output root file not found" << endl;
175  gSystem->Exit(1);
176  }
177 
178  TTree* net_tree = (TTree *) froot->Get("waveburst");
179  if(net_tree!=NULL) {
180  EVT = new netevent(net_tree,nIFO);
181  } else {
182  EVT = new netevent(nIFO);
183  net_tree = EVT->setTree();
184  }
185  EVT->setSLags(gSLAGSHIFT); // set slags into netevent
186 
187  skymap* Psm = new skymap(int(0));
188  net_tree->SetBranchAddress("Psm",&Psm);
189 
190  for(int k=0; k<K; k++) { // loop over the lags
191 
192  id = NET->getwc(k)->get(const_cast<char*>("ID"), 0, 'L', rate);
193 
194  for(int j=0; j<(int)id.size(); j++) { // loop over cluster index
195 
196  int ID = size_t(id.data[j]+0.5);
197 
198  if(NET->getwc(k)->sCuts[ID-1]!=-1) continue; // skip rejected/processed clusters
199 
200  *Psm = NET->nProbability;
201 
202  double ofactor=0;
203  if(cfg->simulation==4) ofactor=-factor;
204  else if(cfg->simulation==3) ofactor=-gIFACTOR;
205  else ofactor=factor;
206  if(cfg->dump) EVT->dopen(outDump.Data(),const_cast<char*>("a"),false);
207  EVT->output2G(net_tree,NET,ID,k,ofactor); // get reconstructed parameters
208  if(cfg->dump) EVT->dclose();
209  if(!cfg->cedDump) NET->getwc(k)->sCuts[ID-1]=1; // mark as processed
210 
211 #ifdef WATPLOT // monster event display
212  SMS.canvas->cd();
213  char fname2[1024];
214  sprintf(fname2, "skyprob_%d.png",ID);
215  cout << "write " << fname2 << endl;
216  SMS.plot(NET->nProbability, 0);
217  SMS.canvas->Print(fname2);
218  SMS.clear();
219 #endif
220  }
221  }
222 
223  jfile->cd();
224 
225  delete Psm;
226  if(EVT) delete EVT;
227  }
228 
229  return;
230 }
231 
std::vector< char * > ifoName
Definition: network.hh:591
CWB::config * cfg
Definition: TestCWB_Plugin.C:5
char analysis[8]
Definition: config.hh:99
TString outDump
virtual void resize(unsigned int)
Definition: wseries.cc:883
virtual size_t size() const
Definition: wavearray.hh:127
size_t nLag
Definition: network.hh:555
void setSLags(float *slag)
Definition: netevent.cc:404
TMacro configPlugin
Definition: config.hh:344
char cmd[1024]
TString Get(wavearray< double > &x, TString ifo)
Definition: mdc.cc:1502
bool mdcPlugin
Definition: config.hh:347
std::vector< std::string > mdcList
Definition: mdc.hh:357
bool dataPlugin
Definition: config.hh:346
cout<< "skymap size : "<< L<< endl;for(int l=0;l< L;l++) sm.set(l, l);sm > const_cast< char * >("skymap.dat")
TTree * setTree()
Definition: netevent.cc:412
TString("c")
std::vector< std::string > mdcType
Definition: mdc.hh:358
size_t nRun
Definition: network.hh:554
int ID
Definition: TestMDC.C:70
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
wavearray< double > get(char *name, size_t index=0, char atype='R', int type=1, bool=true)
param: string with parameter name param: index in the amplitude array, which define detector param: c...
Definition: netcluster.cc:2188
bool cedDump
Definition: config.hh:279
std::vector< std::string > mdcType
Definition: network.hh:595
#define SMS
Definition: FrDisplay.cc:108
CWB::mdc * MDC
Long_t size
virtual void start(double s)
Definition: wavearray.hh:119
int j
Definition: cwb_net.C:10
std::vector< double > mdcTime
Definition: network.hh:596
void plot(wavearray< double > &, char *=NULL, int=1, double=0., double=0., bool=false, float=0., float=0., bool=false, float=0., bool=false)
Definition: watplot.cc:132
bool outPlugin
Definition: config.hh:351
CWB::Toolbox TB
Definition: ComputeSNR.C:5
void Print(int level=0)
Definition: mdc.cc:2707
char ifo[NIFO_MAX][8]
void dopen(const char *fname, char *mode, bool header=true)
Definition: netevent.hh:291
size_t ifoListSize()
Definition: network.hh:413
int Psave
Definition: config.hh:175
void clear()
Definition: watplot.hh:40
#define nIFO
TCanvas * canvas
Definition: watplot.hh:174
double factor
jfile
Definition: cwb_job_obj.C:25
void output2G(TTree *, network *, size_t, int, double)
Definition: netevent.cc:687
static void getSimNoise(wavearray< double > &u, TString fName, int seed, int run)
Definition: Toolbox.cc:4809
TTree * net_tree
Definition: mdc.hh:216
int simulation
Definition: config.hh:181
void CWB_Plugin(TFile *jfile, CWB::config *cfg, network *NET, WSeries< double > *x, TString ifo, int type)
COHERENCE.
#define IMPORT(TYPE, VAR)
Definition: cwb.hh:51
i() int(T_cor *100))
network NET
Definition: cwb_dump_inj.C:12
std::vector< std::string > mdcList
Definition: network.hh:594
int gIFACTOR
TIter next(twave->GetListOfBranches())
char fname[1024]
int k
Definition: skymap.hh:45
TFile * froot
netevent EVT(itree, nifo)
std::vector< int > sCuts
Definition: netcluster.hh:374
void dclose()
Definition: netevent.hh:300
string file
Definition: cwb_online.py:385
netcluster * getwc(size_t n)
param: delay index
Definition: network.hh:421
sprintf(tfres,"(1/%g)x(%g) (sec)x(Hz)", 2 *df, df)
Long_t id
double factors[FACTORS_MAX]
Definition: config.hh:184
std::vector< double > mdcTime
Definition: mdc.hh:360
list files
Definition: cwb_online.py:529
char fName[256]
int check
skymap nProbability
Definition: network.hh:613
bool dump
Definition: config.hh:277