Logo Coherent WaveBurst  
Reference Guide
Logo
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
CWB_Plugin_Sim4.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 <vector>
17 
18 size_t SetFactors(wavearray<double> &wi, std::vector<double>* pT, std::vector<double>* pF, double dT);
19 
20 void
22 
23  // Plugin to generate simulated MDC injected 'on the fly' using simulation=4
24 
25  cout << endl;
26  cout << "-----> CWB_Plugin_Sim4.C" << endl;
27  cout << "ifo " << ifo.Data() << endl;
28  cout << "type " << type << endl;
29  cout << endl;
30 
32 
33  // import istage,jstage
34  int gISTAGE=-1; IMPORT(int,gISTAGE)
35  int gJSTAGE=-1; IMPORT(int,gJSTAGE)
36  cout << "-----> CWB_Plugin_Sim4.C -> " << " gISTAGE : " << gISTAGE << " gJSTAGE : " << gJSTAGE << endl;
37 
38  if(type==CWB_PLUGIN_CONFIG) {
39  cfg->dataPlugin=false; // disable read data from frames
40  cfg->mdcPlugin=true; // disable read mdc from frames
41  }
42 
43  if(type==CWB_PLUGIN_MDC) {
44 
45  // import ifactor
46  int gIFACTOR=-1; IMPORT(int,gIFACTOR)
47  cout << "-----> CWB_Plugin_Sim4.C -> " << " gIFACTOR : " << gIFACTOR << endl;
48 
49  // export to CINT net,cfg
50  char cmd[128];
51  sprintf(cmd,"network* net = (network*)%p;",net);
52  gROOT->ProcessLine(cmd);
53  sprintf(cmd,"CWB::config* cfg = (CWB::config*)%p;",cfg);
54  gROOT->ProcessLine(cmd);
55 
56  CWB::mdc MDC(net);
57 
58  // ---------------------------------
59  // read plugin config
60  // ---------------------------------
61 
62  cfg->configPlugin.Exec();
63 
64  // ---------------------------------
65  // set list of mdc waveforms
66  // ---------------------------------
67 
68  // import MDC form plugin config
69  IMPORT(CWB::mdc,MDC)
70  MDC.Print();
71 
72  // import mdcHRSS form plugin config
73  std::vector<double> mdcHRSS;
74  IMPORT(std::vector<double>,mdcHRSS)
75  int M = mdcHRSS.size();
76  for(int i=0;i<M;i++) cout << i << " mdcHRSS = " << mdcHRSS[i] << endl;
77 
78  // ---------------------------------
79  // get mdc data
80  // ---------------------------------
81 
82  MDC.Get(*x,ifo);
83 
84  // compute normalization factors
85  int K = MDC.mdcList.size();
86  int L = MDC.mdcType.size();
87  std::vector<double> mdcFactor(K);
88  double inj_hrss = MDC.GetInjHrss();
89  for(int i=0;i<K;i++) {
90  for(int j=0;j<L;j++) {
91  if(TString(MDC.mdcList[i]).Contains(MDC.mdcType[j])) mdcFactor[i]=mdcHRSS[j]/inj_hrss;
92  }
93  }
94  for(int i=0;i<K;i++) cout << i << " mdcFactor = " << mdcFactor[i] << endl;
95 
96  // normalize x data
97  SetFactors(*x, &MDC.mdcTime, &mdcFactor, cfg->iwindow/2.);
98 
99  // rescale amplitudes stored in the mdcList
100  for(int k=0; k<(int)K; k++) {
101  int ilog[5] = {1,3,12,13,14};
102  for(int l=0;l<5;l++) {
103  double mfactor = l<2 ? mdcFactor[k] : mdcFactor[k]*mdcFactor[k];
104  TString slog = TB.GetMDCLog(MDC.mdcList[k], ilog[l]);
105  MDC.mdcList[k]=TB.SetMDCLog(MDC.mdcList[k], ilog[l], mfactor*slog.Atof());
106  }
107  }
108 
109  // ---------------------------------
110  // set mdc list in the network class
111  // ---------------------------------
112 
113  if(ifo.CompareTo(net->ifoName[0])==0) {
114 
115  net->mdcList.clear();
116  net->mdcType.clear();
117  net->mdcTime.clear();
118  net->mdcList=MDC.mdcList;
119  net->mdcType=MDC.mdcType;
120  net->mdcTime=MDC.mdcTime;
121  }
122 
123  // for multi stage analysis store mdc infos to job root file
124  if((gJSTAGE!=CWB_STAGE_FULL) && (ifo.CompareTo(net->ifoName[0])==0)) {
125 
126  int ifactor;
127  std::vector<std::string>* mdcList = new vector<std::string>; // list of injections
128  std::vector<std::string>* mdcType = new vector<std::string>; // list of injection types
129  std::vector<double>* mdcTime = new vector<double>; // gps time of selected injections
130 
131  int inj_size = (int)MDC.mdcList.size();
132 
133  jfile->cd();
134 
135  TTree* mdc_tree = (TTree*)jfile->Get("injections");
136  if(mdc_tree==NULL) {
137  mdc_tree = new TTree("injections", "injections");
138  mdc_tree->Branch("ifactor", &ifactor, "ifactor/I");
139  mdc_tree->Branch("mdcList", &mdcList, inj_size, 0);
140  mdc_tree->Branch("mdcType", &mdcType, inj_size, 0);
141  mdc_tree->Branch("mdcTime", &mdcTime, inj_size, 0);
142  } else {
143  mdc_tree->SetBranchAddress("ifactor", &ifactor);
144  mdc_tree->SetBranchAddress("mdcList", &mdcList);
145  mdc_tree->SetBranchAddress("mdcType", &mdcType);
146  mdc_tree->SetBranchAddress("mdcTime", &mdcTime);
147  }
148 
149  ifactor=gIFACTOR;
150 
151  *mdcList=MDC.mdcList;
152  *mdcType=MDC.mdcType;
153  *mdcTime=MDC.mdcTime;
154 
155  mdc_tree->Fill();
156 
157  jfile->Write();
158  }
159 
160  cout.precision(14);
161  for(int k=0;k<(int)MDC.mdcList.size();k++) cout << k << " mdcList " << MDC.mdcList[k] << endl;
162  for(int k=0;k<(int)MDC.mdcTime.size();k++) cout << k << " mdcTime " << MDC.mdcTime[k] << endl;
163  for(int k=0;k<(int)MDC.mdcType.size();k++) cout << k << " mdcType " << MDC.mdcType[k] << endl;
164 
165  }
166 
167  // for multi stage analysis read mdc infos to job root file
168  if((gJSTAGE!=CWB_STAGE_FULL) && (type==CWB_PLUGIN_ILIKELIHOOD)) {
169 
170  // import ifactor
171  int gIFACTOR=-1; IMPORT(int,gIFACTOR)
172  cout << "-----> CWB_Plugin_Sim4.C -> " << " gIFACTOR : " << gIFACTOR << endl;
173 
174  int ifactor;
175  std::vector<std::string>* mdcList = new vector<std::string>; // list of injections
176  std::vector<std::string>* mdcType = new vector<std::string>; // list of injection types
177  std::vector<double>* mdcTime = new vector<double>; // gps time of selected injections
178 
179  jfile->cd();
180  TTree* mdc_tree = (TTree*)jfile->Get("injections");
181  if(mdc_tree==NULL) {
182  cout << "-----> CWB_Plugin_Sim4.C : Error-> " << " injection tree not present in job file" << endl;
183  gSystem->Exit(1);
184  } else {
185  mdc_tree->SetBranchAddress("ifactor", &ifactor);
186  mdc_tree->SetBranchAddress("mdcList", &mdcList);
187  mdc_tree->SetBranchAddress("mdcType", &mdcType);
188  mdc_tree->SetBranchAddress("mdcTime", &mdcTime);
189  }
190 
191  int tree_size = (int)mdc_tree->GetEntries();
192 
193  bool ifactor_check=false;
194  for(int i=0;i<tree_size;i++) {
195  mdc_tree->GetEntry(i);
196  if(ifactor==gIFACTOR) {
197  ifactor_check=true;
198  net->mdcList.clear();
199  net->mdcType.clear();
200  net->mdcTime.clear();
201  net->mdcList=*mdcList;
202  net->mdcType=*mdcType;
203  net->mdcTime=*mdcTime;
204  }
205  }
206  if(!ifactor_check) {
207  cout << "-----> CWB_Plugin_Sim4.C : Error-> " << " mdc infos not present in the injection tree " << endl;
208  gSystem->Exit(1);
209  }
210 
211  cout.precision(14);
212  for(int k=0;k<(int)net->mdcList.size();k++) cout << k << " mdcList " << net->mdcList[k] << endl;
213  for(int k=0;k<(int)net->mdcTime.size();k++) cout << k << " mdcTime " << net->mdcTime[k] << endl;
214  for(int k=0;k<(int)net->mdcType.size();k++) cout << k << " mdcType " << net->mdcType[k] << endl;
215  }
216 
217  return;
218 }
219 
220 // modify input signals (wi) at times pT according the factor pF
221 size_t
222 SetFactors(wavearray<double> &wi, std::vector<double>* pT, std::vector<double>* pF, double dT) {
223 
224  int j,nstop,nstrt;
225  size_t k;
226  double F,T;
227  size_t K = pT->size();
228  size_t N = wi.size();
229  size_t M = 0;
230  double rate = wi.rate(); // simulation rate
231 
232  dT = fabs(dT);
233 
234  wavearray<double> w; w=wi;
235 
236  // isolate injections
237  for(k=0; k<K; k++) {
238 
239  F = (*pF)[k];
240  T = (*pT)[k] - w.start();
241 
242  nstrt = int((T - dT)*rate);
243  nstop = int((T + dT)*rate);
244  if(nstrt<=0) nstrt = 0;
245  if(nstop>=int(N)) nstop = N;
246  if(nstop<=0) continue; // outside of the segment
247  if(nstrt>=int(N)) continue; // outside of the segment
248 
249  for(j=nstrt; j<nstop; j++) {
250  w.data[j]*=F;
251  }
252 
253  M++;
254  }
255  wi = w;
256  return M;
257 }
258 
std::vector< char * > ifoName
Definition: network.hh:591
CWB::config * cfg
Definition: TestCWB_Plugin.C:5
double iwindow
Definition: config.hh:182
virtual size_t size() const
Definition: wavearray.hh:127
double GetInjHrss()
Definition: mdc.hh:278
std::vector< double > mdcHRSS(K)
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
virtual void rate(double r)
Definition: wavearray.hh:123
bool dataPlugin
Definition: config.hh:346
TString("c")
std::vector< std::string > mdcType
Definition: mdc.hh:358
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
std::vector< std::string > mdcType
Definition: network.hh:595
CWB::mdc * MDC
#define M
Definition: UniqSLagsList.C:3
virtual void start(double s)
Definition: wavearray.hh:119
int j
Definition: cwb_net.C:10
i drho i
std::vector< double > mdcTime
Definition: network.hh:596
#define N
CWB::Toolbox TB
Definition: ComputeSNR.C:5
void Print(int level=0)
Definition: mdc.cc:2707
char ifo[NIFO_MAX][8]
network ** net
NOISE_MDC_SIMULATION.
wavearray< double > w
Definition: Test1.C:27
size_t SetFactors(wavearray< double > &wi, std::vector< double > *pT, std::vector< double > *pF, double dT)
jfile
Definition: cwb_job_obj.C:25
static TString GetMDCLog(TString log, int pos)
Definition: Toolbox.cc:2258
Definition: mdc.hh:216
#define IMPORT(TYPE, VAR)
Definition: cwb.hh:51
i() int(T_cor *100))
std::vector< std::string > mdcList
Definition: network.hh:594
int gIFACTOR
int k
double F
double T
Definition: testWDM_4.C:11
double fabs(const Complex &x)
Definition: numpy.cc:37
int ifactor
static TString SetMDCLog(TString log, int pos, TString val)
Definition: Toolbox.cc:2218
int l
Definition: cbc_plots.C:434
sprintf(tfres,"(1/%g)x(%g) (sec)x(Hz)", 2 *df, df)
void CWB_Plugin(TFile *jfile, CWB::config *cfg, network *net, WSeries< double > *x, TString ifo, int type)
COHERENCE.
DataType_t * data
Definition: wavearray.hh:301
double dT
Definition: testWDM_5.C:12
std::vector< double > mdcTime
Definition: mdc.hh:360