Logo coherent WaveBurst  
Library Reference Guide
Logo
cwb.cc
Go to the documentation of this file.
1 /*
2 # Copyright (C) 2019 Gabriele Vedovato, Sergey Klimenko
3 #
4 # This program is free software: you can redistribute it and/or modify
5 # it under the terms of the GNU General Public License as published by
6 # the Free Software Foundation, either version 3 of the License, or
7 # (at your option) any later version.
8 #
9 # This program is distributed in the hope that it will be useful,
10 # but WITHOUT ANY WARRANTY; without even the implied warranty of
11 # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12 # GNU General Public License for more details.
13 #
14 # You should have received a copy of the GNU General Public License
15 # along with this program. If not, see <https://www.gnu.org/licenses/>.
16 */
17 
18 
19 #include "cwb.hh"
20 #include "mdc.hh"
21 #include "time.hh"
22 #include "Math/Polar3D.h"
23 #include "TMD5.h"
24 
25 #define EXIT(ERR) gSystem->Exit(ERR) // better exit handling for ROOT stuff
26 
27 using namespace ROOT::Math;
28 
29 ClassImp(cwb)
30 
31 cwb::cwb(CWB_STAGE jstage) {
32 //
33 // Default Constructor
34 //
35 // used only for cwb class streaming
36 //
37 
38  this->istage=CWB_STAGE_FULL;
39  this->jstage=jstage;
40  this->runID=0;
41  this->singleDetector=false;
42  this->lagBuffer.Set(0);
43  SetupStage(jstage);
44  Init();
45 }
46 
48 //
49 // Constructor
50 //
51 // In this method the cwb configuration is loaded
52 // The default configuration is defined in the $CWB_PARAMETERS_FILE file
53 // The user can provide a customized configuration using the fName or xName parameters
54 //
55 // fName : if fName ends with *.C the parameter is used as cwb user configuration file name
56 // The custom config must contains only the parameters which differ from the defaults
57 // and must not be declared with types (Ex : fLow=64; instead of double fLow=64;)
58 // if fName="" only the default cwb configuration file is used ($CWB_PARAMETERS_FILE)
59 // if fName ends with *.root the parameter is used as job file name
60 // the configuration is read from the job file (the previous processed stage)
61 // istage is initialized with the previous processed stage
62 // xName : this parameter is used as auxiliary cwb configuration file when fName ends with *.root
63 // the xName file is used to change the configuration stored in the job file
64 //
65 // jstage : is the final stage to be processed
66 //
67 
68  this->runID=0;
69  this->singleDetector=false;
70  this->lagBuffer.Set(0);
71  this->istage=CWB_STAGE_FULL;
72  if(fName.EndsWith(".root")) {
73  TFile* ifile = new TFile(fName);
74  if(!ifile->IsOpen()) {
75  cout << "cwb::cwb - Error opening root file : " << fName.Data() << endl;
76  EXIT(1);
77  }
78  // read config object
79  if(ifile->Get("config")!=NULL) {
80  // read config object
81  CWB::config* pcfg = (CWB::config*)ifile->Get("config");
82  cfg = *pcfg;
83  // read cwb object
84  cwb* CWB = (cwb*)ifile->Get("cwb");
85  if(CWB==NULL) {
86  cout << "cwb::cwb - Error : cwb is not contained in root file " << fName.Data() << endl;
87  EXIT(1);
88  }
89  // read runID
90  this->runID = CWB->runID;
91  // read singleDetector
92  this->singleDetector = CWB->singleDetector;
93  // read initial stage
94  this->istage = CWB->jstage;
95  if(jstage<=istage) {
96  cout << "cwb::cwb - Error : the stage(job file) " << GetStageString(istage).Data() << " >= "
97  << " input stage " << GetStageString(jstage).Data() << endl;EXIT(1);
98  EXIT(1);
99  }
100  // read frame structures
101  for(int i=0;i<2*NIFO_MAX;i++) {
102  this->nfrFiles[i] = CWB->nfrFiles[i];
103  this->fr[i] = CWB->fr[i];
104  this->FRF[i] = CWB->FRF[i];
105  }
106  // restore lag setting
107  lagBuffer = CWB->lagBuffer;
108  lagBuffer.Set(lagBuffer.GetSize()+1);
109  lagBuffer[lagBuffer.GetSize()-1]=0; // set end of string
110  lagMode[0]=CWB->lagMode[0];lagMode[1]=0; // restore lagMode
111  // restore detSegs
112  detSegs = CWB->detSegs;
113 
114  delete CWB;
115  // set local environment
116  if(gSystem->Getenv("HOME_WAT_FILTERS")==NULL) {
117  cout << "cwb::cwb - Error : environment HOME_WAT_FILTERS is not defined!!!" << endl;EXIT(1);
118  } else {
119  strcpy(cfg.filter_dir,TString(gSystem->Getenv("HOME_WAT_FILTERS")).Data());
120  }
121  // export config to cint
122  cfg.Export();
123  // set auxiliary configuration
124  if(xName!="") {
125  CWB::config icfg = cfg; // cfg from input job file
126  cfg.Import(xName); // cfg from input job file + auxiliary configuration
127  if(cfg.simulation==4) { // initialize array factors
128  if(int(cfg.factors[0])<=0) cfg.factors[0]=1; // set 1 if <=0
129  int ioffset = int(cfg.factors[0]); // extract offset
130  for(int i=ioffset;i<ioffset+cfg.nfactor;i++) cfg.factors[i]=i; // assign integer values
131  }
132  // check if auxiliary parameters are compatible with the previous stage config
133  if(cfg.simulation!=icfg.simulation) { // simulation
134  cout << "cwb::cwb - Error : aux simulation " << cfg.simulation
135  << " != in cfg simulation " << icfg.simulation << endl; EXIT(1);}
136  if(cfg.nfactor!=icfg.nfactor) { // nfactor
137  cout << "cwb::cwb - Error : aux nfactor " << cfg.nfactor
138  << " != in cfg nfactor " << icfg.nfactor << endl; EXIT(1);
139  } else { // factors
140  for(int i=0;i<=cfg.nfactor;i++) if(cfg.factors[i]!=icfg.factors[i]) {
141  cout << "cwb::cwb - Error : aux factors["<<i<<"]="<<cfg.factors[i]
142  << " != in cfg factors["<<i<<"]=="<<icfg.factors[i]<<endl;EXIT(1);}
143  }
144  if(cfg.l_low!=icfg.l_low) { // l_low
145  cout << "cwb::cwb - Error : aux l_low " << cfg.l_low
146  << " != in cfg l_low " << icfg.l_low << endl; EXIT(1);}
147  if(cfg.l_high!=icfg.l_high) { // l_high
148  cout << "cwb::cwb - Error : aux l_high " << cfg.l_high
149  << " != in cfg l_high " << icfg.l_high << endl; EXIT(1);}
150  if(cfg.fLow!=icfg.fLow) { // fLow
151  cout << "cwb::cwb - Error : aux fLow " << cfg.fLow
152  << " != in cfg fLow " << icfg.fLow << endl; EXIT(1);}
153  if(cfg.fHigh!=icfg.fHigh) { // fHigh
154  cout << "cwb::cwb - Error : aux fHigh " << cfg.fHigh
155  << " != in cfg fHigh " << icfg.fHigh << endl; EXIT(1);}
156  if(cfg.healpix!=icfg.healpix) { // healpix
157  cout << "cwb::cwb - Error : aux healpix " << cfg.healpix
158  << " != in cfg healpix " << icfg.healpix << endl; EXIT(1);}
159  }
160  } else {
161  cout << "cwb::cwb - Error : config is not contained in root file " << fName.Data() << endl;
162  EXIT(1);
163  }
164  ifile->Close();
165  iname=fName;
166  } else if(fName.EndsWith(".C")) {
167  if(gSystem->Getenv("CWB_PARAMETERS_FILE")==NULL) {
168  cout << "cwb::cwb - Error : environment CWB_PARAMETERS_FILE is not defined!!!" << endl;
169  EXIT(1);
170  } else {
171  // load default configuration $CWB_PARAMETERS_FILE
172  cfg.Import("$CWB_PARAMETERS_FILE");
173  }
174  cfg.Import(fName);
175  iname="";
176  } else if(fName=="") {
177  if(gSystem->Getenv("CWB_PARAMETERS_FILE")==NULL) {
178  cout << "cwb::cwb - Error : environment CWB_PARAMETERS_FILE is not defined!!!" << endl;
179  EXIT(1);
180  } else {
181  // load default configuration $CWB_PARAMETERS_FILE
182  cfg.Import("$CWB_PARAMETERS_FILE");
183  }
184  iname="";
185  } else {
186  cout << "cwb::cwb - Error : bad input file extension [.C, .root] " << fName.Data() << endl;
187  EXIT(1);
188  }
189  this->jstage=jstage;
190  SetupStage(jstage);
191  Init();
192  if(fName!="") cfg.Check(); // check parameter's consistency
193 }
194 
196 //
197 // Constructor
198 //
199 // use the config object as configuration
200 //
201 // jstage : is the final stage to be processed
202 //
203 
204  this->cfg = cfg;
205  this->istage=CWB_STAGE_FULL;
206  this->jstage=jstage;
207  this->runID=0;
208  this->singleDetector=false;
209  this->lagBuffer.Set(0);
210  SetupStage(jstage);
211  Init();
212  cfg.Check(); // check parameter's consistency
213 }
214 
216 //
217 // Destructor
218 //
219 
220  if(mdc!=NULL) delete mdc;
221  //if(netburst!=NULL) delete netburst; // commented because is already deleted by "delete froot"
222  if(history!=NULL) delete history;
223  if(jfile!=NULL) delete jfile;
224  if(froot!=NULL) delete froot;
225  if(singleDetector && pD[1]!=NULL) {pD[1]->IWFP.clear();pD[1]->RWFP.clear();}
226  for(int n=0;n<NIFO_MAX;n++) if(pD[n]!=NULL) delete pD[n];
227 }
228 
229 void
231 //
232 // Reset/Initialize the pipeline parameters
233 //
234 
235  // There is a bug in ROOT > 5.32.04 & < 5.34.00
236  // GarbageCollector don't remove the deleted objects from root file
237  // Tree saved to root in the Coherence stage could be corrupted
238 #ifdef _USE_ROOT6
239  //cout << "cwb::Init - Warning !!! : This is a testing version for ROOT6" << endl;
240 #else
241  if(gROOT->GetVersionInt()>53204 && gROOT->GetVersionInt()<53400) {
242  cout << "cwb::Init - Error : cWB analysis don't works with ROOT version > 5.32.04 && version < 5.34.00" << endl;
243  cout << "You are running version : ROOT " << gROOT->GetVersion() << endl << endl;
244  EXIT(1);
245  }
246 #endif
247  if((cfg.simulation==4)&&(jstage==CWB_STAGE_STRAIN)&&(cfg.nfactor>1)) {
248  // NOTE : in the ReadData it is necessary to save MDC data for each factor
249  cout << "cwb::Init - Error : stage STRAIN not implemented with simulation=4" << endl;
250  EXIT(1);
251  }
252 
254  cout<<"cwb::Init - Error : jobfOptions=CWB_JOBF_SAVE_TRGFILE is allowed only with CWB_STAGE_FULL!!!"<<endl;
255  EXIT(1);
256  }
257 
258  if(cfg.fResample>0) { // RESAMPLING
260  } else {
262  }
263 
264  // check if rate is an integer at resolution level = l_high
265  if(rateANA-(1<<cfg.l_high)*(rateANA/(1<<cfg.l_high))!=0) {
266  cout << "cwb::Init - Error : data rate : " << rateANA
267  << " is not a multiple of 2^l_high : " << (1<<cfg.l_high) << endl;
268  EXIT(1);
269  }
270 
271  froot=NULL;
272  jfile=NULL;
273  mdc=NULL;
274  netburst=NULL;
275  history=NULL;
276  for(int n=0;n<NIFO_MAX;n++) pD[n]=NULL;
277  return;
278 }
279 
280 void
281 cwb::run(int runID) {
282 //
283 // The method used to start the analysis
284 //
285 // runID : is the job ID number, this is used in InitJob method to identify the
286 // the time range to be analyzed
287 //
288 // These are the main actions performed by this method
289 //
290 // - LoadPlugin
291 // - InitNetwork
292 // - InitHistory
293 // - InitJob
294 // - Loop over Factors
295 // - ReadData
296 // - DataConditioning
297 // - Coherence
298 // - SuperCluster
299 // - Likelihood
300 // - Save Recontructed Parameters
301 // - Save Job File (only for multi stage analysis)
302 //
303 
304  lags = 0;
305  double factor = 1.0; // strain factor
306  int ioffset = 0; // ifactor offset
307  watchJob.Start(); // start job benchmark
308  watchStage.Start(); // start stage benchmark
309 
310  bplugin = (TString(cfg.plugin.GetName())!="") ? true : false; // user plugin
311 
312  unsigned int Pid = gSystem->GetPid(); // used to tag in a unique way the temporary files
313 
314  double DATA_RATE = 0.;
315 
316  char tmpFile[1024];
317  char outFile[1024];
318  char endFile[1024];
319  char outDump[1024];
320  char endDump[1024];
321  char out_CED[1024];
322  char end_CED[1024];
323  char command[1024];
324  int ecommand=0;
325  FileStat_t fstemp;
326  char cmd[1024];
327 
328  if(cfg.nIFO==0) {cout << "cwb::cwb - Error : no detector is presents in the configuration" << endl;EXIT(1);}
329 
330  this->nIFO=cfg.nIFO;
331  for(int n=0;n<nIFO;n++) {
332  if(strlen(cfg.ifo[n])>0) strcpy(this->ifo[n],cfg.ifo[n]); // built in detector
333  else strcpy(this->ifo[n],cfg.detParms[n].name); // user define detector
334  }
335 
336  if(this->runID==0) { // this->runID>0 for multistage analysis
337  if(runID>0) {
338  this->runID=runID;
339  } else {
340  TString srunID = TString(gSystem->Getenv("CWB_JOBID"));
341  this->runID = srunID.Atoi();
342  }
343  }
344 
345  sprintf(cfg.work_dir,gSystem->WorkingDirectory());
346  sprintf(cfg.data_label,gSystem->BaseName(cfg.work_dir));
347 
348  cout<<"job ID : "<<this->runID<<endl;
349  cout<<"output : "<<cfg.output_dir<<endl;
350  cout<<"label : "<<cfg.data_label<<endl;
351  cout<<"nodedir : "<<cfg.nodedir<<endl;
352  cout<<"Pid : "<<Pid<<endl;
353 
354  // create log, nodedir & output directories
355  // Note : this step is necessary in multi stage analysis when pipeline
356  // start from an intermediate stage from a non structured working dir
357  sprintf(cmd,"mkdir -p %s",cfg.log_dir); gSystem->Exec(cmd);
358  sprintf(cmd,"mkdir -p %s",cfg.nodedir); gSystem->Exec(cmd);
359  sprintf(cmd,"mkdir -p %s",cfg.output_dir); gSystem->Exec(cmd);
360 
361  // export to CINT istage,jstage (used by plugins)
362  sprintf(cmd,"gISTAGE = (CWB_STAGE)%d;",istage); EXPORT(CWB_STAGE,gISTAGE,cmd)
363  sprintf(cmd,"gJSTAGE = (CWB_STAGE)%d;",jstage); EXPORT(CWB_STAGE,gJSTAGE,cmd)
364  // export to CINT ifactor (used by plugins)
365  sprintf(cmd,"gIFACTOR = -1;"); EXPORT(int,gIFACTOR,cmd) // init to gIFACTOR=-1
366 
367  // compile & load user plugin
368  if(bplugin) LoadPlugin(cfg.plugin,cfg.configPlugin);
369 
370  if(bplugin) {CWB_Plugin(NULL,&cfg,&NET,NULL,"",CWB_PLUGIN_CONFIG);SetupStage(jstage);}
371 
372  // ---------------------------------------
373  // init network
374  // ---------------------------------------
375 
376  iname=="" ? InitNetwork() : InitNetwork(iname);
377 
378  PrintAnalysis();
379 
380  // ---------------------------------------
381  // init history
382  // ---------------------------------------
383 
384  InitHistory();
385 
386  // ---------------------------------------
387  // init job
388  // ---------------------------------------
389 
390  double mdcShift = iname=="" ? InitJob() : InitJob(iname);
391 
392  // temporary job file
393  sprintf(jname,"%s/job_%d_%s_%d_%d.root",cfg.nodedir,int(Tb),cfg.data_label,this->runID,Pid);
394  cout<<"temporary job file : " << jname<<endl;
395 
396  if(jstage==CWB_STAGE_INIT) goto JOB_END;
397 
398  // ---------------------------------------
399  // read data (simulation!=4)
400  // ---------------------------------------
401 
402  if(cfg.simulation!=4) {
403  DATA_RATE = ReadData(mdcShift,0);
404  if(jstage==CWB_STAGE_STRAIN) goto JOB_END;
405  }
406 
407  // ---------------------------------------
408  // if simulation!=0, loop on the injection strain factors
409  // ---------------------------------------
410 
411  if(!cfg.simulation) cfg.nfactor = 1;
412 
413  ioffset = (cfg.simulation==4) ? int(cfg.factors[0]) : 0; // ifactor offset
414 
415  for(int ifactor=ioffset; ifactor<ioffset+cfg.nfactor; ifactor++) {
416 
417  int ceddir = 0; // flag if ced directory exists
418 
419  // export to CINT ifactor (used by plugins)
420  sprintf(cmd,"gIFACTOR = %d;",ifactor); EXPORT(int,gIFACTOR,cmd)
421 
422  // ---------------------------------------
423  // declaration of output files
424  // ---------------------------------------
425 
426  if(cfg.simulation) {
427  factor=cfg.factors[ifactor];
428  cout<<endl<< "---> Start processing factor["<<ifactor<< "]="<<factor<<endl<< endl;
429  char sfactor[32];
430  if(cfg.simulation==3) {
431  if(factor<0) sprintf(sfactor,"n%g",fabs(factor));
432  if(factor==0) sprintf(sfactor,"z%g",factor);
433  if(factor>0) sprintf(sfactor,"p%g",factor);
434  } else sprintf(sfactor,"%g",factor);
435  char sim_label[512];
436  sprintf(sim_label,"%d_%d_%s_%s_job%d",int(Tb),int(dT),cfg.data_label,sfactor,this->runID);
437 
438  sprintf(outFile,"%s/wave_%s_%d.root",cfg.nodedir,sim_label,Pid);
439  sprintf(endFile,"%s/wave_%s.root",cfg.output_dir,sim_label);
440  sprintf(tmpFile,"%s/wave_%s_%d.root.tmp",cfg.nodedir,sim_label,Pid);
441  sprintf(outDump,"%s/wave_%s_%d.txt",cfg.nodedir,sim_label,Pid);
442  sprintf(endDump,"%s/wave_%s.txt",cfg.output_dir,sim_label);
443  sprintf(out_CED,"%s/ced_%s_%d",cfg.nodedir,sim_label,Pid);
444  sprintf(end_CED,"%s/ced_%s",cfg.output_dir,sim_label);
445 
446  if(!gSystem->GetPathInfo(endFile,fstemp)) {
447  printf("The file %s already exists - skip\n",endFile);
448  fflush(stdout);
449  TFile rf(endFile);
450  if(!rf.IsZombie()) continue;
451  }
452  }
453  else {
454  char prod_label[512];
455  sprintf(prod_label,"%d_%d_%s_slag%d_lag%lu_%lu_job%d",
456  int(Tb),int(dT),cfg.data_label,slagID,cfg.lagOff,cfg.lagSize,this->runID);
457 
458  sprintf(outFile,"%s/wave_%s_%d.root",cfg.nodedir,prod_label,Pid);
459  sprintf(endFile,"%s/wave_%s.root",cfg.output_dir,prod_label);
460  sprintf(tmpFile,"%s/wave_%s_%d.root.tmp",cfg.nodedir,prod_label,Pid);
461  sprintf(outDump,"%s/wave_%s_%d.txt",cfg.nodedir,prod_label,Pid);
462  sprintf(endDump,"%s/wave_%s.txt",cfg.output_dir,prod_label);
463  sprintf(out_CED,"%s/ced_%s_%d",cfg.nodedir,prod_label,Pid);
464  sprintf(end_CED,"%s/ced_%s",cfg.output_dir,prod_label);
465  }
466 
467  // check if out_CED & end_CED already exist
468  if(cfg.cedDump) {
469  Long_t id,size=0,flags,mt;
470  if (!gSystem->GetPathInfo(out_CED,&id,&size,&flags,&mt)) {
471  cout << "cwb::run - Warning !!! - Dir \"" << out_CED << "\" already exist" << endl;
472  EXIT(0);
473  }
474  if (!gSystem->GetPathInfo(end_CED,&id,&size,&flags,&mt)) {
475  cout << "cwb::run - Warning !!! - Dir \"" << end_CED << "\" already exist" << endl;
476  EXIT(0);
477  }
478  }
479 
480  // remove outDump file
481  Long_t xid,xsize,xflags,xmt;
482  int xestat = gSystem->GetPathInfo(outDump,&xid,&xsize,&xflags,&xmt);
483  if (xestat==0) {
484  sprintf(command,"/bin/rm %s",outDump);
485  ecommand=gSystem->Exec(command);
486  if(ecommand) {cout << "cwb::cwb - Warning -> " << command << endl;}
487  }
488 
489  cout<<"output file on the node : "<<outFile<<endl;
490  cout<<"final output file name : "<<endFile<<endl;
491  cout<<"temporary output file : "<<tmpFile<<endl;
492 
493  // ---------------------------------------
494  // read mdc (simulation==4)
495  // ---------------------------------------
496 
497  if(cfg.simulation==4) {
498  DATA_RATE = ReadData(mdcShift,ifactor);
499  if(TString(cfg.analysis)=="2G" && jstage==CWB_STAGE_STRAIN) continue;
500  }
501 
502  // ---------------------------------------
503  // data conditioning
504  // ---------------------------------------
505 
506  DataConditioning(ifactor);
507  if(TString(cfg.analysis)=="2G" && jstage==CWB_STAGE_CSTRAIN) continue;
508 
509  // ---------------------------------------
510  // create output root file
511  // initialization of tree structures
512  // create prod/sim/lag directories
513  // ---------------------------------------
514 
515  froot = new TFile(tmpFile, "RECREATE");
516  if(froot==NULL) {cout << "cwb::cwb - Error opening root file : " << tmpFile << endl;EXIT(1);}
517  TTree* net_tree = !cfg.outPlugin ? netburst->setTree() : NULL; // outPlugin=true -> out tree is disabled
518  TTree* live_tree= live.setTree();
519  TTree* mdc_tree=NULL;
520  TTree* var_tree=NULL;
521  TTree* noise_tree=NULL;
522 
523  if(cfg.simulation) {
524  mdc_tree = mdc->setTree();
525  } else {
526  var_tree = wavevar.setTree();
527  noise_tree = noiserms.setTree();
528  }
529 
530  // ---------------------------------------
531  // start of the coherent search
532  // ---------------------------------------
533 
534  size_t mlagSize=cfg.lagOff+cfg.lagSize;
535  size_t mlagOff=cfg.lagOff;
536  size_t mlagStep=cfg.mlagStep;
537 
538  if(mlagStep==0) mlagStep=cfg.lagSize; // if mlagStep=0 -> standard lag analysis
539 
540  for(size_t mlag=mlagOff;mlag<mlagSize;mlag+=mlagStep) { // multilag loop
541 
542  cfg.lagOff = mlag;
543  cfg.lagSize = cfg.lagOff+mlagStep<=mlagSize ? mlagStep : mlagSize-cfg.lagOff;
544  if(cfg.lagSize==0) continue;
545 
546  cout << "lagSize : " << cfg.lagSize << " lagOff : " << cfg.lagOff << endl;
547 
548  std::vector<double> livTime;
549  if(iname!="") livTime=NET.livTime; // save livTime
550  if(!cfg.simulation) { // setup lags
552  lagBuffer.GetArray(),lagMode,cfg.lagSite);
553  cout<<"lag step: "<<cfg.lagStep<<endl;
554  cout<<"number of time lags: "<<lags<<endl;
555  }
556  else if(!lags) lags = NET.setTimeShifts();
557  if(iname!="") NET.livTime=livTime; // restore livTime
558 
559  // ---------------------------------------
560  // coherence analysis
561  // ---------------------------------------
562 
563  Coherence(ifactor);
564  if(TString(cfg.analysis)=="2G" && jstage==CWB_STAGE_COHERENCE) continue;
565 
566  // ---------------------------------------
567  // supercluster analysis
568  // ---------------------------------------
569 
570  SuperCluster(ifactor);
571  cout<<endl;
572  if(TString(cfg.analysis)=="1G") cout<<"events in the buffer: "<<NET.events()<<"\n";
573  if(TString(cfg.analysis)=="2G" && jstage==CWB_STAGE_SUPERCLUSTER) continue;
574 
575  // ---------------------------------------
576  // likelihood
577  // ---------------------------------------
578 
579  ceddir=Likelihood(ifactor, out_CED, netburst, net_tree, outDump);
580  cout<<"\nSearch done\n";
581  if(!cfg.cedDump) cout<<"reconstructed events: "<<NET.events()<<"\n";
582  if(cfg.simulation) NET.printwc(0);
583 
584  // ---------------------------------------
585  // save data into root file
586  // ---------------------------------------
587 
588  PrintStageInfo(CWB_STAGE_SAVE,"Data Save");
589  froot->cd();
590 
591  live.output(live_tree,&NET,slagShift,detSegs);
592 
593  if(cfg.simulation) {
594  double ofactor=0;
595  if(cfg.simulation==4) ofactor=-factor;
596  else if(cfg.simulation==3) ofactor=-ifactor;
597  else ofactor=factor;
598  if(TString(cfg.analysis)=="1G") {
599  if(cfg.dump) netburst->dopen(outDump,const_cast<char*>("w"));
600  netburst->output(net_tree,&NET,ofactor);
601  if(cfg.dump) netburst->dclose();
602  }
603  mdc->output(mdc_tree,&NET,ofactor);
604  } else {
605  if(TString(cfg.analysis)=="1G") {
606  if(cfg.dump) netburst->dopen(outDump,const_cast<char*>("w"));
607  netburst->output(net_tree,&NET);
608  if(cfg.dump) netburst->dclose();
609  }
610  for(int i=0; i<nIFO; i++) {
611  if(cfg.outfOptions&CWB_OUTF_SAVE_VAR) // write var tree
612  wavevar.output(var_tree,&v[i],i+1,cfg.segEdge);
613  if(cfg.outfOptions&CWB_OUTF_SAVE_NOISE) // write noise tree
614  noiserms.output(noise_tree,&pD[i]->nRMS,i+1,DATA_RATE/2);
615  }
616  }
617 
618  // save log info to final output root file
619  PrintStageInfo(CWB_STAGE_FINISH,"Job Finished",false,true,endFile);
620  history->AddLog( const_cast<char*>("FULL"), const_cast<char*>("STOP JOB"));
621  history->Write("history");
622 
623  froot->Write();
624 
625  } // end mlag loop
626 
627  froot->Close();
628 
629  // ---------------------------------
630  // move output files to end dirs
631  // ---------------------------------
632 
633  if((TString(cfg.analysis)=="2G") && (jstage!=CWB_STAGE_FULL && jstage!=CWB_STAGE_LIKELIHOOD)) {
634  // delete temporary output file
635  sprintf(command,"/bin/rm %s", tmpFile);
636  ecommand=gSystem->Exec(command);
637  if(ecommand) {cout << "cwb::cwb - Warning -> " << command << endl;}
638  continue;
639  }
640  sprintf(command,"/bin/mv %s %s", tmpFile, outFile);
641  if(!cfg.cedDump || (cfg.cedDump && cfg.online)) Exec(command);
642  if(cfg.cedDump && !cfg.online) sprintf(command,"/bin/rm %s",tmpFile);
643  else sprintf(command,"/bin/mv %s %s",outFile,endFile);
644  Exec(command);
645  if(cfg.cedDump && !cfg.online) sprintf(command,"/bin/rm %s",outDump);
646  else sprintf(command,"/bin/mv %s %s",outDump,endDump);
647  if(cfg.dump) {
648  // if file outDump exists then it is moved to endDump
649  if(!gSystem->GetPathInfo(outDump,&xid,&xsize,&xflags,&xmt)) Exec(command);
650  }
651  xestat = gSystem->GetPathInfo(end_CED,&xid,&xsize,&xflags,&xmt);
652  if (xestat==0) {
653  sprintf(command,"/bin/mv %s/* %s/.",out_CED,end_CED);
654  } else {
655  sprintf(command,"/bin/mv %s %s",out_CED,end_CED);
656  }
657  if(cfg.cedDump && ceddir && !(jobfOptions&CWB_JOBF_SAVE_CED)) Exec(command);
658  }
659 
660  JOB_END:
661 
662  sprintf(cmd,"gIFACTOR = -1;"); EXPORT(int,gIFACTOR,cmd) // set gIFACTOR=-1
663 
664  // -------------------------------------------------------------------
665  // write cwb, history, config, network objects into temporary job file
666  // -------------------------------------------------------------------
667 
669  // set options to save trigger file when full stage is used
671  //this->cwb::jobfOptions = CWB_JOBF_SAVE_DISABLE;
673  // clear network skymap filled in Likelihood stage (save space)
674  NET. nSensitivity=0; NET. nAlignment=0; NET. nCorrelation=0;
675  NET. nLikelihood=0; NET. nNullEnergy=0; NET. nPenalty=0;
676  NET. nCorrEnergy=0; NET. nNetIndex=0; NET. nDisbalance=0;
677  NET. nSkyStat=0; NET. nEllipticity=0; NET. nPolarisation=0;
678  NET. nProbability=0;
679  }
680 
681  jfile = new TFile(jname,"UPDATE");
682  if(jfile==NULL||!jfile->IsOpen())
683  {cout << "cwb::cwb - Error opening root file : " << jname << endl;EXIT(1);}
684  jfile->cd();
685  if(jobfOptions&CWB_JOBF_SAVE_CONFIG) { // write config object
686  TString tempName = cfg.configPlugin.GetName(); // save temporary macro name
687  TString tempTitle = cfg.configPlugin.GetTitle(); // save temporary macro title
688  // get original macro name
689  TObjString* objn = cfg.configPlugin.GetLineWith("//#?CONFIG_PLUGIN_NAME ");
690  TString origName = objn ? objn->GetString() : "";
691  cfg.configPlugin.SetName(origName.ReplaceAll("//#?CONFIG_PLUGIN_NAME ",""));
692  // get original macro title
693  TObjString* objt = cfg.configPlugin.GetLineWith("//#?CONFIG_PLUGIN_TITLE ");
694  TString origTitle = objt ? objt->GetString() : "";
695  cfg.configPlugin.SetTitle(origTitle.ReplaceAll("//#?CONFIG_PLUGIN_TITLE ",""));
696  cfg.Write("config");
697  cfg.configPlugin.SetName(tempName); // restore temporary macro name
698  cfg.configPlugin.SetTitle(tempTitle); // restore temporary macro title
699  }
701  // in 2G analysis the network::getNetworkPixels do "pixeLHood = *pTF;"
702  // the Forwad operation after WDM::setTDFilter include in pTF the TD structures
703  // since WDM and SymmObjArray streamers are incomplete because of its pointers
704  // the read/save to/from root jfile is not a safe operation
705  // a workaround is to set to 0 the TD structures
706  ((WDM<double>*)NET.pixeLHood.pWavelet)->T0.Resize(0);
707  ((WDM<double>*)NET.pixeLHood.pWavelet)->Tx.Resize(0);
708  }
709  if(jobfOptions&CWB_JOBF_SAVE_NETWORK) NET.Write("network"); // write network object
710  if(jobfOptions&CWB_JOBF_SAVE_JNET_MACRO) { // write cwb_jnet macro
711  // check if cwb_jnet.C macro exists
712  TString cwb_jnet_name = gSystem->ExpandPathName("$CWB_MACROS/cwb_jnet.C");
713  TB.checkFile(cwb_jnet_name);
714  TMacro cwb_jnet(cwb_jnet_name);
715  cwb_jnet.Write("cwb_jnet");
716  }
717  PrintStageInfo(CWB_STAGE_FINISH,"Job Finished",false,true);
718  history->AddLog( const_cast<char*>("FULL"), const_cast<char*>("STOP JOB"));
719  if(jobfOptions&CWB_JOBF_SAVE_HISTORY) history->Write("history");
720  if(jobfOptions&CWB_JOBF_SAVE_CWB) this->cwb::Write("cwb");
721  if(bplugin) {
722  wavearray<double> x; // temporary time series
723  for(int i=0; i<nIFO; i++) {
724  x.rate(cfg.inRate); x.start(FRF[i].start); x.resize(int(x.rate()*(FRF[i].stop-FRF[i].start)));
726  }
727  x.resize(0);
728  }
729  jfile->Close();
730 
731  // ---------------------------------
732  // clean-up temporary job file
733  // ---------------------------------
734 
735  if(jobfOptions) {
736  // get rid of duplicated trees created by TFileMergers
737  TFile *f = TFile::Open(jname,"UPDATE");
738  TDirectory* d;
739  for(int ifactor=0; ifactor<cfg.nfactor; ifactor++) {
740  for(int j=0; j<(int)NET.nLag; j++) {
741  netcluster* pwc = NET.getwc(j);
742  if(pwc==NULL) continue;
743  int cycle = cfg.simulation ? ifactor : Long_t(pwc->shift);
744  char trName[64];sprintf(trName,"clusters-cycle:%d;2",cycle);
745  d = (TDirectory*)f->Get("coherence;1");
746  if(d!=NULL) if(d->Get(trName)!=NULL) d->Delete(trName);
747  d = (TDirectory*)f->Get("supercluster;1");
748  if(d!=NULL) if(d->Get(trName)!=NULL) d->Delete(trName);
749  d = (TDirectory*)f->Get("likelihood;1");
750  if(d!=NULL) if(d->Get(trName)!=NULL) d->Delete(trName);
751  }
752  }
753  f->Close();
754 
755  // assign final file name according to the stage name
756  char jlabel[512]; sprintf(jlabel,"%d_%d_%s",int(Tb),int(dT),cfg.data_label);
757  char wlabel[512]; sprintf(wlabel,"wave_%s",jlabel);
758  TString jLABEL=jlabel;
759  if((jstage!=CWB_STAGE_FULL)&&(jstage!=CWB_STAGE_LIKELIHOOD))
760  sprintf(endFile,"%s/wave_%s_job%d.root",cfg.output_dir,jlabel,this->runID);
761  TString ojname = endFile;
762  if(TString(cfg.analysis)=="1G") ojname.ReplaceAll(wlabel,TString("job_")+jLABEL);
763  else { // 2G
764  if(jstage==CWB_STAGE_FULL) ojname.ReplaceAll(wlabel,TString("job_")+jLABEL);
765  if(jstage==CWB_STAGE_INIT) ojname.ReplaceAll(wlabel,TString("init_")+jLABEL);
766  if(jstage==CWB_STAGE_STRAIN) ojname.ReplaceAll(wlabel,TString("strain_")+jLABEL);
767  if(jstage==CWB_STAGE_CSTRAIN) ojname.ReplaceAll(wlabel,TString("cstrain_")+jLABEL);
768  if(jstage==CWB_STAGE_COHERENCE) ojname.ReplaceAll(wlabel,TString("coherence_")+jLABEL);
769  if(jstage==CWB_STAGE_SUPERCLUSTER) ojname.ReplaceAll(wlabel,TString("supercluster_")+jLABEL);
770  if(jstage==CWB_STAGE_LIKELIHOOD) ojname.ReplaceAll(wlabel,TString("job_")+jLABEL);
771  }
772 
773  // if CWB_JOBF_SAVE_NODE then 2G to store the intermediate stage job files to nodedir
774  // and creates a symbolic link to the final working dir output file name
776  TString olname = ojname; // symbolic link path destination
777  ojname.ReplaceAll(TString(cfg.output_dir)+TString("/"),TString(cfg.nodedir)+TString("/"));
778  // symbolic link path origin
779  if(ojname.BeginsWith("/")) {
780  sprintf(command,"/bin/ln -sf %s %s",ojname.Data(),olname.Data());
781  } else {
782  sprintf(command,"/bin/ln -sf ../%s %s",ojname.Data(),olname.Data());
783  }
784  cout << command << endl;
785  ecommand=gSystem->Exec(command);
786  if(ecommand) {cout << "cwb::cwb - Warning -> " << command << endl;}
787  }
788 
789  sprintf(endFile,"%s",ojname.Data()); // used by PrintStageInfo
790 
791  // move temporary job file to final position
792  sprintf(command,"/bin/mv %s %s",jname,ojname.Data());
793  ecommand=gSystem->Exec(command);
794  if(ecommand) {cout << "cwb::cwb - Warning -> " << command << endl;}
795 
796  } else {
797 
798  // remove temporary job file
799  sprintf(command,"/bin/rm %s",jname);
800  ecommand=gSystem->Exec(command);
801  if(ecommand) {cout << "cwb::cwb - Warning -> " << command << endl;}
802  }
803 
804  // remove temporary configPlugin file
805  if(TString(cfg.configPlugin.GetTitle())!="")
806  gSystem->Exec(TString("rm ")+cfg.configPlugin.GetTitle());
807 
808  // ---------------------------------
809  // print final infos
810  // ---------------------------------
811 
812  cout << endl << endl;
813  PrintStageInfo(CWB_STAGE_FINISH,"Job Finished",true,false,endFile);
814  int job_data_size_sec = dT;
815  double job_speed_factor = double(job_data_size_sec)/double(watchJob.RealTime());
816  cout << endl;
817  printf("Job Speed Factor - %2.2fX\n",job_speed_factor);
818  cout << endl;
819 
820  watchJob.Stop();watchJob.Reset();
821 
822  // ONLINE : create empty file "finished" in the output directory
823  if(cfg.online) {
824  sprintf(command,"touch %s/finished",cfg.output_dir);
825  gSystem->Exec(command);
826  }
827 
828  return;
829 }
830 
831 void
833 //
834 // Init the Network object : set detectors and network parameters
835 //
836 // read the network object from the job file (the previous processed stage)
837 //
838 
839  PrintStageInfo(CWB_STAGE_INIT,"cwb::InitNetwork from file");
840 
841  // read network object from job file
842  TFile* ifile = new TFile(fName);
843  if(ifile==NULL) {cout << "cwb::InitNetwork - Error opening root file : " << fName.Data() << endl;EXIT(1);}
844  network* pnet = (network*)ifile->Get("network");
845  if(pnet!=NULL) {
846  // copy network object into local NET object
847  NET = *pnet;
848  delete pnet;
849  } else {
850  cout << "cwb::InitNetwork - Error : net is not contained in root file " << fName.Data() << endl;
851  EXIT(1);
852  }
853  ifile->Close();
854  // save detector pointers into local structures pD
855  for(int n=0; n<nIFO; n++) pD[n] = NET.getifo(n);
856  // set the original sampling rate
857  for(int i=0; i<nIFO; i++) pD[i]->rate = cfg.fResample>0 ? cfg.fResample : cfg.inRate;
858 
859  // restore skymaps (to save space the network skymaps was not saved into job file)
860  if(cfg.healpix) NET.setSkyMaps(int(cfg.healpix));
862  NET.setAntenna();
863 
864  // restore network parameters
867  NET.Edge = cfg.segEdge;
868  NET.netCC = cfg.netCC;
869  NET.netRHO = cfg.netRHO;
870  NET.EFEC = cfg.EFEC;
872  NET.nSky = cfg.nSky;
874  NET.setRunID(runID);
876  NET.optim=cfg.optim;
878 
879  // restore sky mask
880  if(TString(cfg.analysis)=="1G") {
882  else if(cfg.mask<0 && strlen(cfg.skyMaskFile)>0)
883  SetSkyMask(&NET,&cfg,cfg.skyMaskFile,'e');
884  } else {
885  if(strlen(cfg.skyMaskFile)>0) SetSkyMask(&NET,&cfg,cfg.skyMaskFile,'e');
886  }
887  if(strlen(cfg.skyMaskCCFile)>0) SetSkyMask(&NET,&cfg,cfg.skyMaskCCFile,'c');
888 
889  // restore detector names defined in the network class
890  // there is a issue in network class (no TClass for char*)
891  // detectors names are not saved properly
892  NET.ifoName.clear();
893  for(int n=0; n<nIFO; n++) NET.ifoName.push_back(pD[n]->Name);
894 
895  // execute user plugin
896  if(bplugin) CWB_Plugin(NULL,&cfg,&NET,NULL,"",CWB_PLUGIN_NETWORK);
897 
898  // create injection (injected events) and netevent (reconstructed events) objects
899  mdc = new injection(nIFO);
900  netburst = new netevent(nIFO,cfg.Psave);
901 
902  // print system infos
903  gSystem->Exec("/bin/date");
904  gSystem->Exec("/bin/hostname");
905 
906  return;
907 }
908 
909 void
911 //
912 // Init the Network object :
913 //
914 // - create object & add detector objects to network
915 // - set network parameters : search type, regulators, thresholds, sky map grid resolution
916 // - set the user sky mask : enable/disable the sky map tiles (earth coordinates)
917 // - set the user celestial sky mask : enable/disable the sky map tiles (celestial coordinates)
918 //
919 
920  PrintStageInfo(CWB_STAGE_INIT,"cwb::InitNetwork");
921 
922  // ---------------------------------------
923  // Check single detector mode
924  // ---------------------------------------
925 
926  // if nIFO=1 the analysis is done as a network of 2 equal detectors
927  if(nIFO==1) {
928  // ------> cwb in Sigle Detector Mode !!!
930  nIFO = cfg.nIFO;
931  strcpy(ifo[1],ifo[0]);
932  singleDetector=true;
933  } else singleDetector=false;
934 
935  // ---------------------------------------
936  // init network
937  // ---------------------------------------
938 
939  // when ifo!="" check if user detectors defined in detParams has been initialized
940  for(int i=0; i<nIFO; i++) {
941  if(strlen(cfg.ifo[i])==0 && strlen(cfg.detParms[i].name)==0) {
942  cout << "cwb::InitNetwork - Error : user detector name at position "
943  << i << " is not defined (detParms)" << endl;
944  EXIT(1);
945  }
946  }
947 
948  // create detector objects
949  for(int i=0; i<nIFO; i++) {
950  if(strlen(cfg.ifo[i])>0) pD[i] = new detector(cfg.ifo[i]); // built in detector
951  else pD[i] = new detector(cfg.detParms[i]); // user define detector
952  }
953  // set the original sampling rate
954  for(int i=0; i<nIFO; i++) pD[i]->rate = cfg.fResample>0 ? cfg.fResample : cfg.inRate;
955  // add detector object to network
956  for(int i=0; i<nIFO; i++) NET.add(pD[i]);
957 
958  // set network skymaps
959  if(cfg.healpix) NET.setSkyMaps(int(cfg.healpix));
961  NET.setAntenna();
962 
963  // restore network parameters
966  NET.Edge = cfg.segEdge;
967  NET.netCC = cfg.netCC;
968  NET.netRHO = cfg.netRHO;
969  NET.EFEC = cfg.EFEC;
971  NET.nSky = cfg.nSky;
973  NET.setRunID(runID);
975  NET.optim=cfg.optim;
977 
978  // set sky mask
979  if(TString(cfg.analysis)=="1G") {
981  else if(cfg.mask<0 && strlen(cfg.skyMaskFile)>0)
982  SetSkyMask(&NET,&cfg,cfg.skyMaskFile,'e');
983  } else {
984  if(strlen(cfg.skyMaskFile)>0) SetSkyMask(&NET,&cfg,cfg.skyMaskFile,'e');
985  }
986  if(strlen(cfg.skyMaskCCFile)>0) SetSkyMask(&NET,&cfg,cfg.skyMaskCCFile,'c');
987 
988  // execute user plugin
989  if(bplugin) CWB_Plugin(NULL,&cfg,&NET,NULL,"",CWB_PLUGIN_NETWORK);
990 
991  // create injection (injected events) and netevent (reconstructed events) objects
992  mdc = new injection(nIFO);
993  netburst = new netevent(nIFO,cfg.Psave);
994 
995  // print system infos
996  gSystem->Exec("/bin/date");
997  gSystem->Exec("/bin/hostname");
998 
999  return;
1000 }
1001 
1003 //
1004 // Init History : it is used for bookkeping
1005 //
1006 // for each of these stages :
1007 // (FULL,INIT,STRAIN,CSTRAIN,COHERENCE,SUPERCLUSTER,LIKELIHOOD)
1008 // the following informations are saved
1009 //
1010 // CWB_ENV : cwb enviromental variable
1011 // CWB_ENV_MD5 : the MD5 of the CWB_ENV string
1012 // WATVERSION : WAT version
1013 // XIFO : XIFO (max detectors's number)
1014 // GITVERSION : git version
1015 // GITBRANCH : git branch
1016 // WORKDIR : working directory path
1017 // FRLIBVERSION : framelib version
1018 // ROOTVERSION : ROOT version
1019 // LALVERSION : LAL version
1020 // DATALABEL : label used to tag the analysis
1021 // CMDLINE : command line used to start the job
1022 // ROOTLOGON : the cwb rootlogon.C
1023 // ROOTLOGON_MD5 : the MD5 of the ROOTLOGON_MD5 string
1024 // PARAMETERS : the user_parameters.C file
1025 // PARAMETERS_MD5 : the MD5 of the PARAMETERS_MD5 string
1026 // CWB_CONFIG_URL : cWB/config url
1027 // CWB_CONFIG_PATH : cWB/config path
1028 // CWB_CONFIG_BRANCH : cWB/config branch
1029 // CWB_CONFIG_TAG : cWB/config tag
1030 // CWB_CONFIG_HASH : cWB/config hash code
1031 // CWB_CONFIG_DIFF : cWB/config difference (modified)
1032 //
1033 
1034  PrintStageInfo(CWB_STAGE_INIT,"cwb::InitHistory");
1035 
1036  // check if configuration files exist
1037  TB.checkFile(gSystem->Getenv("CWB_ROOTLOGON_FILE"));
1038 
1039  // declare stages
1040  const char* STAGE_NAMES[7] = {"FULL","INIT","STRAIN","CSTRAIN","COHERENCE","SUPERCLUSTER","LIKELIHOOD"};
1041  // declare types
1042  const char* TYPE_NAMES[22] = {"CWB_ENV","WATVERSION","XIFO","GITVERSION","GITBRANCH","WORKDIR",
1043  "FRLIBVERSION","ROOTVERSION","LALVERSION",
1044  "DATALABEL","CMDLINE","ROOTLOGON","PARAMETERS",
1045  "CWB_ENV_MD5","ROOTLOGON_MD5","PARAMETERS_MD5",
1046  "CWB_CONFIG_URL","CWB_CONFIG_PATH","CWB_CONFIG_BRANCH","CWB_CONFIG_TAG","CWB_CONFIG_HASH","CWB_CONFIG_DIFF"};
1047  // create history object
1048  history = new CWB::History(const_cast<char**>(STAGE_NAMES), 7, const_cast<char**>(TYPE_NAMES), 22);
1050 
1051  // If any, add history from previous processed stage
1052  if(iname!="") {
1053  TFile* ifile = new TFile(iname);
1054  if(ifile==NULL) {cout << "cwb::InitHistory - Error opening root file : " << iname.Data() << endl;EXIT(1);}
1055  // read network object
1056  if(ifile->Get("history")!=NULL) {
1057  // read history object
1058  CWB::History* hst = (CWB::History*)ifile->Get("history");
1059  TList* stageList = hst->GetStageNames();
1060  TList* typeList = hst->GetTypeNames();
1061  for(int i=0;i<stageList->GetSize();i++) {
1062  TObjString* stageString = (TObjString*)stageList->At(i);
1063  for(int j=0;j<typeList->GetSize();j++) {
1064  TObjString* typeString = (TObjString*)typeList->At(j);
1065  char* histData = hst->GetHistory(const_cast<char*>(stageString->GetString().Data()),
1066  const_cast<char*>(typeString->GetString().Data()));
1067  if(histData!=NULL) { // add previous history stage to current history
1068  history->AddHistory(const_cast<char*>(stageString->GetString().Data()),
1069  const_cast<char*>(typeString->GetString().Data()),histData);
1070  delete histData;
1071  }
1072  }
1073  }
1074  delete stageList;
1075  delete typeList;
1076  // read logs
1077  int log_size = hst->GetLogSize(const_cast<char*>("FULL"));
1078  for(int i=0;i<log_size;i++) {
1079  TString log = hst->GetLog(const_cast<char*>("FULL"),i);
1080  // Replace the final stage (STG:8) of the input job file with (STG:jstage)
1081  char stg_label[16];sprintf(stg_label,"STG:%d",jstage);
1082  if(log.Contains("STG:8")) log.ReplaceAll("STG:8",stg_label);
1083  // Add previous history logs to current history
1084  history->AddLog(const_cast<char*>("FULL"), const_cast<char*>(log.Data()));
1085  }
1086  delete hst;
1087  } else {
1088  cout << "cwb::InitHistory - Error : history is not contained in root file " << iname.Data() << endl;
1089  EXIT(1);
1090  }
1091  ifile->Close();
1092  }
1093 
1094  // add current history stage
1095  TString jStageString = GetStageString(jstage).ReplaceAll("CWB_STAGE_","");
1096  char jStage[256];sprintf(jStage,jStageString.Data());
1097 
1098  // get cwb enviromental variable and save to history
1099  char* cwbBuffer = TB.getEnvCWB();
1100  if(cwbBuffer!=NULL) {
1101  history->AddHistory(jStage, const_cast<char*>("CWB_ENV"), cwbBuffer);
1102  TMD5 md5;md5.Update((UChar_t*)cwbBuffer,strlen(cwbBuffer));md5.Final();
1103  history->AddHistory(jStage, const_cast<char*>("CWB_ENV_MD5"), const_cast<char*>(md5.AsString()));
1104  delete [] cwbBuffer;
1105  }
1106 
1107  // save library versions
1108  CWB::mdc MDC;
1109  char framelib_version[32]; sprintf(framelib_version,"%f",FRAMELIB_VERSION);
1110  history->AddHistory(jStage, const_cast<char*>("WATVERSION"), watversion('s'));
1111  history->AddHistory(jStage, const_cast<char*>("XIFO"), watversion('i'));
1112  history->AddHistory(jStage, const_cast<char*>("GITVERSION"), watversion('r'));
1113  history->AddHistory(jStage, const_cast<char*>("GITBRANCH"), watversion('b'));
1114  history->AddHistory(jStage, const_cast<char*>("FRLIBVERSION"), framelib_version);
1115  history->AddHistory(jStage, const_cast<char*>("ROOTVERSION"), const_cast<char*>(gROOT->GetVersion()));
1116  history->AddHistory(jStage, const_cast<char*>("LALVERSION"), const_cast<char*>(GetLALVersion().Data()));
1117 
1118  // save work directory and data label
1119  history->AddHistory(jStage, const_cast<char*>("WORKDIR"), cfg.work_dir);
1120  history->AddHistory(jStage, const_cast<char*>("DATALABEL"), cfg.data_label);
1121 
1122  // save command line used to start the application
1123  char cmd_line[2048]="";
1124  int cmd_line_len=0;
1125  for(int i=0;i<gApplication->Argc();i++) cmd_line_len+=strlen(gApplication->Argv(i));
1126  if(cmd_line_len>2047)
1127  {cout << "cwb::InitHistory - command line too long : " << cmd_line_len << endl;EXIT(1);}
1128  for(int i=0;i<gApplication->Argc();i++) sprintf(cmd_line,"%s %s",cmd_line,gApplication->Argv(i));
1129  history->AddHistory(jStage, const_cast<char*>("CMDLINE"), cmd_line);
1130 
1131  // save rootlogon script (ROOT initialization script)
1132  char* rootlogonBuffer = TB.readFile(gSystem->Getenv("CWB_ROOTLOGON_FILE"));
1133  if(rootlogonBuffer!=NULL) {
1134  history->AddHistory(jStage, const_cast<char*>("ROOTLOGON"), rootlogonBuffer);
1135  TMD5 md5;md5.Update((UChar_t*)rootlogonBuffer,strlen(rootlogonBuffer));md5.Final();
1136  history->AddHistory(jStage, const_cast<char*>("ROOTLOGON_MD5"), const_cast<char*>(md5.AsString()));
1137  delete [] rootlogonBuffer;
1138  }
1139 
1140  // save git cWB/config infos
1141  history->AddHistory(jStage, const_cast<char*>("CWB_CONFIG_URL"), const_cast<char*>(GetGitInfos("url","$CWB_CONFIG").Data()));
1142  history->AddHistory(jStage, const_cast<char*>("CWB_CONFIG_PATH"), const_cast<char*>(GetGitInfos("path","$CWB_CONFIG").Data()));
1143  history->AddHistory(jStage, const_cast<char*>("CWB_CONFIG_BRANCH"), const_cast<char*>(GetGitInfos("branch","$CWB_CONFIG").Data()));
1144  history->AddHistory(jStage, const_cast<char*>("CWB_CONFIG_TAG"), const_cast<char*>(GetGitInfos("tag","$CWB_CONFIG").Data()));
1145  history->AddHistory(jStage, const_cast<char*>("CWB_CONFIG_HASH"), const_cast<char*>(GetGitInfos("hash","$CWB_CONFIG").Data()));
1146  history->AddHistory(jStage, const_cast<char*>("CWB_CONFIG_DIFF"), const_cast<char*>(GetGitInfos("diff","$CWB_CONFIG")!="" ? "M" : ""));
1147 
1148  // save configuration file
1149  char tmpFile[1024];
1150  unsigned int Pid = gSystem->GetPid(); // used to tag in a unique way the temporary files
1151  sprintf(tmpFile,"%s/CWB_Config_%s_%d_job%d.XXXXXX",cfg.nodedir,cfg.data_label,Pid,this->runID);
1152  if(mkstemp(tmpFile)==-1) { // create temporary file
1153  cout << "cwb::InitHistory - mkstemp error in creating tmp file : " << tmpFile << endl;
1154  EXIT(1);
1155  }
1156  char nodedir[1024];strcpy(nodedir,cfg.nodedir); // clean nodedir before dump to history
1157  strcpy(cfg.nodedir,""); // this makes config the same for all jobs
1158  TString tempTitle = cfg.configPlugin.GetTitle(); // save temporary macro title
1159  // get original macro title
1160  TObjString* objt = cfg.configPlugin.GetLineWith("//#?CONFIG_PLUGIN_TITLE ");
1161  TString origTitle = objt ? objt->GetString() : "";
1162  cfg.configPlugin.SetTitle(origTitle.ReplaceAll("//#?CONFIG_PLUGIN_TITLE ",""));
1163  cfg.Print(tmpFile); // write config to temporary file
1164  cfg.configPlugin.SetTitle(tempTitle); // restore temporary macro title
1165  strcpy(cfg.nodedir,nodedir); // restore nodedir
1166  char* parametersBuffer = TB.readFile(tmpFile); // read config from tmp file
1167  gSystem->Exec((TString("/bin/rm ")+TString(tmpFile)).Data()); // delete temporary file
1168  if(parametersBuffer!=NULL) {
1169  history->AddHistory(jStage, const_cast<char*>("PARAMETERS"), parametersBuffer);
1170  TMD5 md5;md5.Update((UChar_t*)parametersBuffer,strlen(parametersBuffer));md5.Final();
1171  history->AddHistory(jStage, const_cast<char*>("PARAMETERS_MD5"), const_cast<char*>(md5.AsString()));
1172  delete [] parametersBuffer;
1173  }
1174 
1175  // save START JOB to Log history
1176  history->AddLog(const_cast<char*>("FULL"), const_cast<char*>("START JOB"));
1177 
1178  return;
1179 }
1180 
1182 //
1183 // Init the Job
1184 //
1185 // read the Job configuration parameters from the job file (the previous processed stage)
1186 //
1187 
1188  PrintStageInfo(CWB_STAGE_INIT,"cwb::InitJob from file");
1189 
1190  double mdcShift=0.;
1191 
1192  // read cwb object from job file
1193  TFile* ifile = new TFile(fName);
1194  if(ifile==NULL||!ifile->IsOpen())
1195  {cout << "cwb::InitJob - Error opening root file : " << fName.Data() << endl;EXIT(1);}
1196  ifile->cd();
1197  cwb* CWB = (cwb*)ifile->Get("cwb");
1198  if(CWB==NULL) {
1199  cout << "cwb::InitJob - Error : cwb is not contained in root file " << fName.Data() << endl;
1200  EXIT(1);
1201  }
1202 
1203  // restore cwb segment parameters from saved cwb object
1204  jobID = CWB->jobID;
1205  Tb = CWB->Tb;
1206  dT = CWB->dT;
1207  Te = CWB->Te;
1208  slagID = CWB->slagID;
1209  for(int n=0;n<=(int)NET.ifoListSize();n++) segID[n]=CWB->segID[n]; // restore segID
1210  for(int n=0;n<=(int)NET.ifoListSize();n++) slagShift[n]=CWB->slagShift[n]; // restore slags
1211  netburst->setSLags(slagShift); // set slags into netevent
1212  ifile->Close();
1213 
1214  // restore TFmap start,rate (size=0) , needed by likelihood to get the initial segment time
1215  for(int n=0;n<(int)NET.ifoListSize();n++) {
1217  NET.getifo(n)->getTFmap()->rate(rateANA);
1218  }
1219 
1220  // set a fake TFmap for the first detector, needed to restore setTimeShifts
1221  if(TString(cfg.analysis)=="2G" && jstage==CWB_STAGE_LIKELIHOOD) {
1222  NET.getifo(0)->getTFmap()->resize(rateANA*(dT+2*cfg.segEdge));
1223  for(int n=0;n<NET.ifoListSize();n++) { // used for ced
1224  NET.getifo(n)->getTFmap()->setlow(cfg.fLow);
1226  }
1227  }
1228 
1229  // init simulation structures for simulation=4
1230  if(cfg.simulation==4) { // initialize array factors
1231  // for simulation=4 factors are used as tags for the output root files
1232  if(fmod(cfg.factors[0],1)) { // must be integer
1233  cout<<"cwb::InitJob : when simulation=4 factors[0] is the offset and must be integer>=0"<<endl;
1234  EXIT(1);
1235  }
1236  if(cfg.factors[0]<0) {
1237  cout<<"cwb::InitJob : when simulation=4 factors[0] is the offset and must be integer>=0"<<endl;
1238  EXIT(1);
1239  }
1240  if(int(cfg.factors[0])<=0) cfg.factors[0]=1; // set 1 if <=0
1241  int ioffset = int(cfg.factors[0]); // extract offset
1242  for(int i=ioffset;i<ioffset+cfg.nfactor;i++) cfg.factors[i]=i; // assign integer values
1243  }
1244 
1245  return mdcShift;
1246 }
1247 
1248 double cwb::InitJob() {
1249 //
1250 // Init the Job structures
1251 //
1252 // Selection of the period to be analyzed by this job (job segment)
1253 // - Read the data quality segment (DQ) files (CAT 0/1/2/4) : CWB::config::DFQ
1254 // CAT 0/1/4 are used to define where the data can be analyzed
1255 // CAT 2 are used to veto bad periods
1256 //
1257 // The job segment is computed according to the segment parameters : CWB::config::segXXX, CWB::config::slagXXX
1258 // Each job segment includes at beginning and at the end a scratch period which length is : CWB::config::segEdge
1259 // There are 2 types of procedures used to select the job segment :
1260 // - 1) the circular time-shifts (LAG) are performed within non shifted detector segments (same as 1G)
1261 // the length of the job segment is optimized to get the maximum livetime.
1262 // The minimum segment len is : CWB::config::segMLS , the maximum segment length is CWB::config::segLen
1263 // - 2) the circular time-shifts (LAG) are performed within shifted detector segments (SLAG)
1264 // all job segments have the same length : CWB::config::segLen
1265 //
1266 // Initialization of frame file structures
1267 // - Read noise/MDC frame lists
1268 // noise/MDC can be read from files or from the plugin (generated 'On The Fly')
1269 // - Read MDC log file (only for simulation)
1270 // Is the list of the injection parameters (inj time, hrss, source location, ...)
1271 //
1272 
1273  PrintStageInfo(CWB_STAGE_INIT,"cwb::InitJob");
1274 
1275  double mdcShift=0.;
1276 
1277  vector<TString> ifos(nIFO);
1278  for(int n=0;n<nIFO;n++) ifos[n]=ifo[n];
1279 
1280  if((cfg.simulation)&&((cfg.slagSize>1)||(cfg.lagSize!=1))) {
1281  cout << "Error : slagSize<=1 & lagSize==1 in simulation mode !!!" << endl;EXIT(1);
1282  }
1283 
1284  // -----------------------------------------------------------------------
1285  // init cat1 segments
1286  //
1287  // job segments are computed using as segEdge the value segEdge+segOverlap
1288  // NB! it is used only to compute the job segments
1289  //
1290  // |x| = segEdge : |+| = segOverlap : |-| = segLen
1291  //
1292  // |xxx|-------------------|++++++|xxx|
1293  // |xxx|-------------------|++++++|xxx|
1294  // -----------------------------------------------------------------------
1295 
1296  if(cfg.slagSize>0) { // SLAG Segments
1297 
1298  cout << endl << "START SLAG Init ..." << endl << endl;
1299 
1300  // Get zero lag merged dq cat0 & cat1 lists
1301  // The cat1List list is used only for slagSegs && slagJobList
1303 
1304  // Get number/list of available super lag jobs
1305  // Compute the available segments with length=segLen contained between the interval [min,max]
1306  // Where min,max are the minimum and macimum times in the cat1List list
1307  // The start time of each segment is forced to be a multiple of segLen
1308  vector<waveSegment> slagJobList=TB.getSlagJobList(cat1List, int(cfg.segLen));
1309  int slagSegs=slagJobList.size();
1310 
1311  // Get super lag list : slagList
1312  // Is the list of available segment shifts according to the slag configuration parameters
1313  vector<slag> slagList=TB.getSlagList(nIFO, cfg.slagSize, slagSegs, cfg.slagOff,
1315 
1316  // Extract SLAG from the the slagList according the input runID
1317  slag SLAG=TB.getSlag(slagList,runID);
1318  if(SLAG.jobId!=runID) {cout << "jobID " << runID << " not found in the slag list !!!" << endl;EXIT(1);}
1319  slagID = SLAG.slagId[0]; // slag identification number
1320  jobID = SLAG.jobId; // for slag -> jobID=runID
1321  for(int n=0; n<nIFO; n++) segID[n]=SLAG.segId[n]; // segment shifts
1322  cout << "SuperLag=" << slagID << " jobID=" << jobID;
1323  for(int n=0; n<nIFO; n++) cout << " segID[" << ifo[n] << "]=" << segID[n];cout << endl;
1324 
1325  // Apply slag shifts to the DQF structures (data quality)
1326  for(int n=0;n<nIFO;n++) ifos[n]=ifo[n];
1327  TB.setSlagShifts(SLAG, ifos, cfg.segLen, cfg.nDQF, cfg.DQF);
1328 
1329  // Get merged dq cat0 & cat1 lists after slag shifts
1331 
1332  // Extract ifo's slag segments ranges contained in the selected slag
1333  // Extract max length segment after the application of cat0+cat1
1334  // Final segment length must be gt [cfg.segMLS+2*(cfg.segEdge+cfg.segOverlap)] & gt [cfg.segMLS]
1336  cout << endl;
1337  cout << endl << "Segment type = SLAG" << endl;
1338 
1339  cout << endl << "END SLAG Init ..." << endl << endl;
1340 
1341  } else { // Segments are not shifted (standard shifts)
1342 
1343  jobID = runID;
1344  slagID = 0;
1345  for(int n=0;n<nIFO;n++) segID[n]=0;
1346 
1347  // Get merged dq cat0 & cat1 lists
1349 
1350  // Extract ifo's slag segments ranges
1351  // Extract max length segment after the application of cat0+cat1
1352  // Final segment length must be gt [cfg.segMLS+2*(cfg.segEdge+cfg.segOverlap)] & gt [cfg.segMLS]
1354  cout << endl;
1355  cout << "Segment type = LAG" << endl;
1356  }
1357  cout << "segLen = " << cfg.segLen << " sec" << endl;
1358  cout << "segMLS = " << cfg.segMLS << " sec" << endl;
1359  cout << "segOverlap = " << cfg.segOverlap << " sec" << endl;
1360  cout << endl;
1361 
1362  // set slags in to the netevent class
1363  for(int n=0; n<nIFO; n++) slagShift[n]=(segID[n]-segID[0])*cfg.segLen;
1364  slagShift[nIFO]=slagID;
1365  netburst->setSLags(slagShift);
1366 
1367  // Check if seg length is compatible with WDM parity (only for 2G)
1368  // This condition is necessary to avoid mixing between odd
1369  // and even pixels when circular buffer is used for lag shift
1370  // The MRAcatalog distinguish odd and even pixels
1371  // If not compatible then the length is modified according the requirements
1372  if(TString(cfg.analysis)=="2G") {
1373  int rate_min = rateANA>>cfg.l_high;
1374  for(int i=0;i<(int)detSegs.size();i++) {
1375  double length = detSegs[i].stop - detSegs[i].start;
1376  if(int(length*rate_min+0.001)&1) detSegs[i].stop-=1;
1377  }
1378  }
1379 
1380  // add segOverlap to the detector's segments stop for this job
1381  for(int i=0;i<nIFO;i++) detSegs[i].stop+=cfg.segOverlap;
1382  // print detector's segments for this job
1383  cout.precision(14);
1384  for(int i=0;i<nIFO;i++) {
1385  cout << "detSegs_dq1[" << ifo[i] << "] GPS range : "
1386  << detSegs[i].start << "-" << detSegs[i].stop << endl;
1387  }
1388  cout << endl;
1389  if(detSegs.size()==0) {cout << "no segments found for this job, job terminated !!!" << endl;EXIT(1);}
1390 
1391  // --------------------------------------------------------------------------
1392  // Fill Frame Files structures (FRF)
1393  // if (dataPlugin=false) noise data is read from frame files else from plugin
1394  // if (mdcPlugin=false) mdc data is read from frame files else from plugin
1395  // --------------------------------------------------------------------------
1396 
1397  // init FRF structures for noise data
1398  for(int n=0; n<nIFO; n++) {
1399  if(!cfg.dataPlugin) { // noise data from frames
1400  if(TString(cfg.frFiles[n])=="") {
1401  cout << "cwb::InitJob - Error : noise frame files list name is not defined" << endl;
1402  EXIT(1);
1403  }
1405  detSegs[n].stop -cfg.dataShift[n]+2*cfg.segEdge);
1406  fr[n].open(cfg.frFiles[n]); // read frame files list
1407  fr[n].setVerbose();
1409  fr[n].setSRIndex(TMath::Nint(TMath::Log2(cfg.inRate))); // force data to be resampled to inRate
1410  nfrFiles[n]=fr[n].getNfiles(); // number of frame files used for this job
1411  cout << ifo[n] << " -> nfrFiles : " << nfrFiles[n] << endl;
1412  FRF[n] = fr[n].getFrList(int(detSegs[n].start-cfg.dataShift[n]),
1413  int(detSegs[n].stop-cfg.dataShift[n]), int(cfg.segEdge));
1414  } else { // noise data from the user plugin
1417  FRF[n].length = FRF[n].stop-FRF[n].start;
1418  }
1419  }
1420 
1421  // init FRF structures for mdc data
1422  if(cfg.simulation) {
1423  for(int n=nIFO; n<2*nIFO; n++) {
1424  // mdc are declared in the array frFiles[nIFO,2*nIFO-1]
1425  // or within a single frame file frFiles[nIFO]
1426  // in such case the frFiles[nIFO+1,2*nIFO-1] are a copy of frFiles[nIFO]
1427  if(TString(cfg.frFiles[n])=="" && n>nIFO) {FRF[n]=FRF[nIFO];continue;}
1428  if(!cfg.mdcPlugin) { // mdc from frames
1429  // set frame file list
1430  if(TString(cfg.frFiles[n])=="") {
1431  cout << "cwb::InitJob - Error : MDC frame files list name is not defined" << endl;
1432  EXIT(1);
1433  }
1434  if(cfg.mdc_shift.startMDC>=0) { // cfg.mdc_shift.startMDC<0 -> automatically mdc time shift
1435  fr[n].setTimeRange(detSegs[n-nIFO].start-2*cfg.segEdge,detSegs[n-nIFO].stop+2*cfg.segEdge);
1436  }
1437  fr[n].open(cfg.frFiles[n]);
1438  fr[n].setVerbose();
1440  fr[n].setSRIndex(TMath::Nint(TMath::Log2(cfg.inRate))); // force data to be resampled to inRate
1441  nfrFiles[n]=fr[n].getNfiles();
1442  cout << "MDC " << " -> nfrFiles : " << nfrFiles[n] << endl;
1443  if(cfg.mdc_shift.startMDC<0) { // the mdc shift is automatically selected
1444  // read mdc range from mdc frl files
1445  waveSegment mdc_range = fr[n].getFrRange();
1446  cout << "mdc_range : " << mdc_range.start << " " << mdc_range.stop << endl;
1447  cfg.mdc_shift.startMDC=mdc_range.start;
1448  cfg.mdc_shift.stopMDC=mdc_range.stop;
1449  }
1450  mdcShift = TB.getMDCShift(cfg.mdc_shift, detSegs[n-nIFO].start);
1451  cout << "mdcShift : " << mdcShift << endl;
1452  FRF[n] = fr[n].getFrList(int(detSegs[n-nIFO].start-mdcShift),
1453  int(detSegs[n-nIFO].stop-mdcShift), int(cfg.segEdge));
1454  } else { // mdc from user plugin
1455  FRF[n]=FRF[n-nIFO];
1456  }
1457  }
1458 
1459  // init simulation structures for simulation=4
1460  if(cfg.simulation==4) { // initialize array factors
1461  // for simulation=4 factors are used as tags for the output root files
1462  if(fmod(cfg.factors[0],1)) { // must be integer
1463  cout<<"cwb::InitJob : when simulation=4 factors[0] is the offset and must be integer>=0"<<endl;
1464  EXIT(1);
1465  }
1466  if(cfg.factors[0]<0) {
1467  cout<<"cwb::InitJob : when simulation=4 factors[0] is the offset and must be integer>=0"<<endl;
1468  EXIT(1);
1469  }
1470  if(int(cfg.factors[0])<=0) cfg.factors[0]=1; // set 1 if <=0
1471  int ioffset = int(cfg.factors[0]); // extract offset
1472  for(int i=ioffset;i<ioffset+cfg.nfactor;i++) cfg.factors[i]=i; // assign integer values
1473  }
1474  }
1475 
1476  // ------------------------------
1477  // init cat2 segments
1478  // ------------------------------
1479 
1480  // get shifted merged dq cat0+cat1+cat2 list
1482  // store cat2List into network class
1484 
1485  // check if seg+cat2 data length is >= segTHR
1486  vector<waveSegment> detSegs_dq2;
1487  detSegs_dq2.push_back(detSegs[0]);
1488  detSegs_dq2 = TB.mergeSegLists(detSegs_dq2,cat2List);
1489  for(int i=0;i<(int)detSegs_dq2.size();i++) {
1490  cout << "detSegs_dq2[" << i << "] GPS range : "
1491  << detSegs_dq2[i].start << "-" << detSegs_dq2[i].stop << endl;
1492  }
1493  cout << endl;
1494  double detSegs_ctime = TB.getTimeSegList(detSegs_dq2);
1495  cout << "live time after cat 2 : " << detSegs_ctime << " sec" << endl;
1496  if(detSegs_ctime<cfg.segTHR) {
1497  cout << "cwb::InitJob : job segment live time after cat2 < segTHR="
1498  << cfg.segTHR << " sec, job terminated !!!" << endl;
1499  cout << endl << "To remove this check set segTHR=0 in the parameter file" << endl << endl;
1500  EXIT(1);
1501  }
1502 
1503  Tb=detSegs[0].start; // start seg time (excluded segEdge)
1504  Te=detSegs[0].stop; // stop seg time (excluded segEdge)
1505  dT = Te-Tb; // length seg time (excluded segEdge)
1506 
1507  // --------------------------------------------------------------------------------
1508  // read input mdc log file when data loaded from frame files (cfg.mdcPlugin==false)
1509  // --------------------------------------------------------------------------------
1510  if(cfg.simulation && !cfg.mdcPlugin) { // read MDC log file
1511  int i=NET.readMDClog(cfg.injectionList,double(long(Tb))-mdcShift);
1512  printf("GPS: %16.6f saved, injections: %d\n",double(long(Tb)),i);
1513  if(!cfg.mdcPlugin) TB.shiftBurstMDCLog(NET.mdcList, ifos, mdcShift);
1514  for(int i=0;i<(int)NET.mdcTime.size();i++) NET.mdcTime[i]+=mdcShift;
1515 
1516  // check if seg+cat2+inj data length is not zero
1517  vector<waveSegment> mdcSegs(NET.mdcTime.size());
1518  for(int k=0;k<(int)NET.mdcTime.size();k++) {
1519  mdcSegs[k].start=NET.mdcTime[k]-cfg.iwindow/2.;
1520  mdcSegs[k].stop=NET.mdcTime[k]+cfg.iwindow/2.;
1521  }
1522  vector<waveSegment> mdcSegs_dq2 = TB.mergeSegLists(detSegs_dq2,mdcSegs);
1523  double mdcSegs_ctime = TB.getTimeSegList(mdcSegs_dq2);
1524  cout << "live time in zero lag after cat2+inj : " << mdcSegs_ctime << endl;
1525  if(mdcSegs_ctime==0) {
1526  cout << "cwb::InitJob : job segment with zero cat2+inj "
1527  << "live time in zero lag, job terminated !!!" << endl;
1528  cout << "Warning : MDC data frames could be zero !!!" << endl;
1529  EXIT(1);
1530  }
1531  }
1532 
1533  // ------------------------------------------------------------
1534  // store the contents of the user lagFile to string (lagBuffer)
1535  // (used in multistages analysis to restore lags)
1536  // ------------------------------------------------------------
1537  lagBuffer.Set(0);
1538  lagMode[0]=cfg.lagMode[0];lagMode[1]=0;
1539  if(!cfg.simulation && cfg.lagFile!=NULL) {
1540  if(TString(cfg.lagMode)=="r") { // read lags from file
1541  ifstream in;
1542  in.open(cfg.lagFile, ios::in);
1543  if(!in.good()) {
1544  cout << "cwb::InitJob : Error Opening File : " << cfg.lagFile << endl;
1545  EXIT(1);
1546  }
1547  // get length of file:
1548  in.seekg (0, in.end);
1549  int len = in.tellg();
1550  in.seekg (0, in.beg);
1551  lagBuffer.Set(len);
1552  in.read(lagBuffer.GetArray(),len);// store the contents of the file into the lagBuffer
1553  lagMode[0]='s';lagMode[1]=0; // change lagMode : setTimeShifts read lags from string
1554  if(!in) {
1555  cout << "cwb::InitJob : Error Reading File : " << cfg.lagFile << endl;
1556  EXIT(1);
1557  }
1558  in.close();
1559  // set end of string
1560  lagBuffer.Set(lagBuffer.GetSize()+1);
1561  lagBuffer[lagBuffer.GetSize()-1]=0;
1562  } else { // write lags to file
1563  lagBuffer.Set(strlen(cfg.lagFile),cfg.lagFile);
1564  lagMode[0]='w';lagMode[1]=0;
1565  }
1566  }
1567 
1568  // execute user plugin
1569  if(bplugin) {
1570  wavearray<double> x; // temporary time series
1571  for(int i=0; i<nIFO; i++) {
1572  x.rate(cfg.inRate); x.start(FRF[i].start); x.resize(int(x.rate()*(FRF[i].stop-FRF[i].start)));
1574  }
1575  x.resize(0);
1576  }
1577 
1578  return mdcShift;
1579 }
1580 
1581 double
1583 //
1584 // Read unconditioned Noise/MDC from input job file (the previous processed stage)
1585 // Data are saved into current job file
1586 //
1587 
1588  PrintStageInfo(CWB_STAGE_STRAIN,"cwb::ReadData from file");
1589 
1590  // open in update mode the current job file (jfile)
1591  jfile = new TFile(jname,"UPDATE");
1592  if(jfile==NULL||!jfile->IsOpen())
1593  {cout << "cwb::ReadData - Error opening root file : " << jname << endl;EXIT(1);}
1594  TDirectory* cdstrain = (TDirectory*)jfile->Get("strain");
1595  if(cdstrain==NULL) cdstrain = jfile->mkdir("strain");
1596 
1597  // read unconditioned Noise/MDC from input job file (ifile)
1598  // write unconditioned Noise/MDC to current job file (jfile)
1599  TFile* ifile = new TFile(fName);
1600  if(ifile==NULL) {cout << "cwb::ReadData - Error opening root file : " << fName << endl;EXIT(1);}
1601 
1602  WSeries<double>* pws=NULL;
1603  for(int i=0; i<nIFO; i++) {
1604  ifile->cd();
1605  // read strain from ifile
1606  pws = (WSeries<double>*)ifile->Get(TString("strain/")+pD[i]->Name);
1607  if(pws==NULL)
1608  {cout << "cwb::ReadData - Error : data not present in file : " << fName << endl;EXIT(1);}
1609  // write strain to jfile
1610  pTF[i] = pD[i]->getTFmap();
1611  cdstrain->cd();pws->Write(ifo[i]);
1612  if(cfg.simulation) {
1613  ifile->cd();
1614  // read mdc from ifile
1615  pws = (WSeries<double>*)ifile->Get(TString("mdc/")+ifo[i]);
1616  if(pws==NULL)
1617  {cout << "cwb::ReadData - Error : mdc not present in file : " << fName << endl;EXIT(1);}
1618  // write mdc to jfile
1619  TDirectory* cdmdc = (TDirectory*)jfile->Get("mdc");
1620  if(cdmdc==NULL) cdmdc = jfile->mkdir("mdc");
1621  cdmdc->cd();pws->Write(ifo[i]);
1622  }
1623  if(singleDetector) break;
1624  }
1625 
1626  // close ifile/jfile
1627  ifile->Close();
1628  jfile->Close();
1629  return pws->rate();
1630 }
1631 
1632 void
1634 //
1635 // Read conditioned/whitened Noise/MDC from input job file (the previous processed stage)
1636 // Data are saved into current job file
1637 //
1638 
1639  PrintStageInfo(CWB_STAGE_CSTRAIN,"cwb::DataConditioning from file");
1640 
1641  char cdcstrain_name[32];sprintf(cdcstrain_name,"cstrain-f%d",ifactor);
1642 
1643  // open in update mode the current job file (jfile)
1644  jfile=NULL;
1645  TDirectory* jcdcstrain = NULL;
1647  jfile = new TFile(jname,"UPDATE");
1648  if(jfile==NULL||!jfile->IsOpen())
1649  {cout << "cwb::DataConditioning - Error opening root file : " << jname << endl;EXIT(1);}
1650  jcdcstrain=jfile->mkdir(cdcstrain_name);
1651  }
1652 
1653  // read conditioned data from input job file (ifile)
1654  // write conditioned data to current job file (jfile)
1655  TFile* ifile = new TFile(fName);
1656  if(ifile==NULL)
1657  {cout << "cwb::DataConditioning - Error opening root file : " << fName << endl;EXIT(1);}
1658 
1659  WSeries<double>* pws=NULL;
1660  for(int i=0; i<nIFO; i++) {
1661  ifile->cd();
1662  // read whitened data from ifile
1663  pws = (WSeries<double>*)ifile->Get(TString(cdcstrain_name)+TString("/")+pD[i]->Name);
1664  if(pws==NULL)
1665  {cout << "cwb::DataConditioning - Error : data not present, job terminated!!!" << endl;EXIT(1);}
1666  // write whitened data to jfile
1667  if(jfile!=NULL) {jfile->cd();jcdcstrain->cd();pws->Write(ifo[i]);}
1668  pTF[i] = pD[i]->getTFmap();
1669  *pTF[i]=*pws;
1670  delete pws;
1671  if(singleDetector) break;
1672  }
1673 
1674  // close ifile/jfile
1675  if(jfile!=NULL) jfile->Close();
1676  ifile->Close();
1677  return;
1678 }
1679 
1680 void
1682 //
1683 // Base class virtual Coherence method
1684 // It is implemented in cwb1G/cwb2G
1685 //
1686 
1687  PrintStageInfo(CWB_STAGE_COHERENCE,"cwb::Coherence from file");
1688 
1689  return;
1690 }
1691 
1692 void
1694 //
1695 // Base class virtual SuperCluster method
1696 // It is implemented in cwb1G/cwb2G
1697 //
1698 
1699  PrintStageInfo(CWB_STAGE_SUPERCLUSTER,"cwb::SuperCluster from file");
1700 
1701  return;
1702 }
1703 
1704 void
1705 cwb::PrintAnalysis(bool stageInfos) {
1706 //
1707 // Print Analysis Setup
1708 // input : stageInfos : if true stage infos are printed
1709 //
1710 // output format (Example) :
1711 //
1712 // jobID : 1
1713 //
1714 // Analysis : 2G
1715 // stage : CSTRAIN
1716 //
1717 // detectors : L1 H1 V1
1718 // search : un-modeled(r)-MRA-Standard
1719 // simulation : 2
1720 // maximum time delay between detectors : 0.0264483
1721 // maximum time delay difference : 0.0264483
1722 // HEALPix order : 7
1723 // levelR, rateANA, fLow, fHigh : 4, 1024, 32, 512
1724 // bpp, Acore : 0.001, 1.5
1725 // netRHO and netCC : 5, 0.5
1726 // clustering TFgap, Tgap, Fgap : 6.0, 0.0, 0.0
1727 // pattern : 0
1728 // subnetcut subnet, subcut : 0.6, 0.33
1729 // regulator delta, gamma : 0.05, 0.5
1730 // precision csize, order : 50, 5
1731 //
1732 
1733  if(stageInfos) PrintStageInfo(CWB_STAGE_INIT,"cwb::PrintAnalysis");
1734 
1735  cout << " jobID : " << runID << endl;
1736  cout << endl;
1737  cout << " Analysis : " << cfg.analysis << endl;
1738  cout << " stage : ";
1739  if(TString(cfg.analysis)=="1G") cout << "FULL" << endl;
1740  else { // 2G
1741  if(jstage==CWB_STAGE_FULL) cout << "FULL" << endl;
1742  if(jstage==CWB_STAGE_INIT) cout << "INIT" << endl;
1743  if(jstage==CWB_STAGE_STRAIN) cout << "STRAIN" << endl;
1744  if(jstage==CWB_STAGE_CSTRAIN) cout << "CSTRAIN" << endl;
1745  if(jstage==CWB_STAGE_COHERENCE) cout << "COHERENCE" << endl;
1746  if(jstage==CWB_STAGE_SUPERCLUSTER) cout << "SUPERCLUSTER" << endl;
1747  if(jstage==CWB_STAGE_LIKELIHOOD) cout << "LIKELIHOOD" << endl;
1748  }
1749 
1750  cout << "\n detectors : ";
1751  for(int i=0; i<cfg.nIFO; i++) cout<<ifo[i]<<" ";
1752  cout<<endl;
1753 
1754  cout << " search : ";
1755  if(TString(cfg.analysis)=="1G") {
1756  if(cfg.search=='E' || cfg.search=='E') cout<<"un-modeled (Energy)";
1757  if(cfg.search=='b' || cfg.search=='B') cout<<"un-modeled (single stream)";
1758  if(cfg.search=='r' || cfg.search=='R') cout<<"un-modeled (dual stream)";
1759  if(cfg.search=='i' || cfg.search=='I') cout<<"elliptical polarisation";
1760  if(cfg.search=='g' || cfg.search=='G') cout<<"circular polarisation";
1761  if(cfg.search=='s' || cfg.search=='S') cout<<"linear polarisation";
1762  cout << "(" << cfg.search << ")" << endl;
1763  } else if(TString(cfg.analysis)=="2G") {
1764  char _search = std::tolower(cfg.search);
1765  if(_search=='r') cout<<"un-modeled";
1766  if(_search=='i') cout<<"iota-wave";
1767  if(_search=='p') cout<<"psi-wave";
1768  if(_search=='l' || _search=='s') cout<<"linear polarisation";
1769  if(_search=='c' || _search=='g') cout<<"circular polarisation";
1770  if(_search=='e' || _search=='b') cout<<"elliptical polarisation";
1771  cout << "(" << cfg.search << ")";
1772  if(cfg.optim) cout<<"-SRA";
1773  else cout<<"-MRA";
1774  if(cfg.pattern==0) cout<<"-Standard";
1775  if((cfg.pattern!=0 && cfg.pattern<0)) cout<<"-Mixed";
1776  if((cfg.pattern!=0 && cfg.pattern>0)) cout<<"-Packet";
1777  cout << endl;
1778  } else {
1779  cout<<"\n unavailable analysis type !!!"<<endl;EXIT(1);
1780  }
1781  cout << " simulation : " << cfg.simulation << endl;;
1782 
1783  mTau=NET.getDelay(const_cast<char*>("MAX")); // maximum time delay
1784  dTau=NET.getDelay(const_cast<char*>("")); // time delay difference
1785  cout<<"maximum time delay between detectors : "<<mTau<<endl;
1786  cout<<" maximum time delay difference : "<<dTau<<endl;
1787  if(cfg.healpix) {
1788  cout<<" HEALPix order : "<<cfg.healpix<<endl;
1789  } else {
1790  cout<<" skymap angular resolution : "<<cfg.angle<<endl;
1791  cout<<" skymap size in polar angle : "<<cfg.Theta1<<", "<<cfg.Theta2<<endl;
1792  cout<<" skymap size in azimuthal angle : "<<cfg.Phi1<<", "<<cfg.Phi2<<endl;
1793  }
1794  cout<<" levelR, rateANA, fLow, fHigh : "<<cfg.levelR<<", "<<rateANA<<", "
1795  <<cfg.fLow<<", "<<cfg.fHigh<<endl;
1796  cout<<" bpp, Acore : "<<cfg.bpp<<", "<<cfg.Acore<<endl;
1797  cout<<" netRHO and netCC : "<<cfg.netRHO<<", "<<cfg.netCC<<endl;
1798  if(TString(cfg.analysis)=="1G") {
1799  cout<<" regulator delta, gamma : "<<cfg.delta<<", "<<NET.gamma<<endl;
1800  } else { // 2G
1801  cout<<" clustering TFgap, Tgap, Fgap : "<<cfg.TFgap<<", "<<cfg.Tgap<<", "<<cfg.Fgap<<endl;
1802 
1803  cout<<" pattern : "<<cfg.pattern<<" -> ";
1804  int pattern = abs(cfg.pattern);
1805  if(pattern== 0) cout<<"('*' single pixel standard)"<<endl;
1806  if(pattern== 1) cout<<"('3|' packet)"<<endl;
1807  if(pattern== 2) cout<<"('3-' packet)"<<endl;
1808  if(pattern== 3) cout<<"('3/' packet - chirp)"<<endl;
1809  if(pattern== 4) cout<<"('3\\' packet - ringdown)"<<endl;
1810  if(pattern== 5) cout<<"('5/' packet - chirp)"<<endl;
1811  if(pattern== 6) cout<<"('5\\' packet - ringdown)"<<endl;
1812  if(pattern== 7) cout<<"('3+' packet)"<<endl;
1813  if(pattern== 8) cout<<"('3x' cross packet)"<<endl;
1814  if(pattern== 9) cout<<"('9p' 9-pixel square packet)"<<endl;
1815  if(pattern > 9) cout<<"('*' single pixel packet)"<<endl;
1816 
1817  cout<<" subnet, subcut : "<<cfg.subnet<<", "<<cfg.subcut<<endl;
1818  cout<<" regulator delta, gamma : "<<cfg.delta<<", "<<NET.gamma<<endl;
1819  }
1820  if(cfg.precision!=0) {
1821  int precision = int(fabs(cfg.precision));
1822  int csize = precision%65536; // get number of pixels threshold per level
1823  int order = (precision-csize)/65536; // get resampled order
1824  cout<<" precision csize, order : "<<csize<<", "<<order<<endl;
1825  }
1826  if((TString(cfg.analysis)=="1G")&&(cfg.mask>0.)) {
1827  cout<<" earth sky mask applied ; "<<cfg.skyMaskFile<<endl;
1828  cout<<" mask ; "<<cfg.mask<<endl;
1829  }
1830  if((TString(cfg.analysis)=="1G")&&(cfg.mask<0.)&&(strlen(cfg.skyMaskFile)>0))
1831  cout<<" earth sky mask applied ; "<<cfg.skyMaskFile<<endl;
1832  if((TString(cfg.analysis)=="2G")&&(strlen(cfg.skyMaskFile)>0))
1833  cout<<" earth sky mask applied ; "<<cfg.skyMaskFile<<endl;
1834  if(strlen(cfg.skyMaskCCFile)>0)
1835  cout<<" celestial sky mask applied ; "<<cfg.skyMaskCCFile<<endl;
1836  cout<<endl;
1837 
1838  // print ifo infos
1839  //cout << " nIFO : " << nIFO << endl;
1840  //NET.print();
1841  //for(int n=0; n<nIFO; n++) pD[n]->print();
1842 
1843  return;
1844 }
1845 
1846 size_t
1847 cwb::GetProcInfo(bool mvirtual) {
1848 //
1849 // Print to screen System Infos
1850 //
1851 // mvirtual : true/false -> {virtual memory)/(resident memory)
1852 //
1853 // output format (Example) :
1854 //
1855 // Mon Feb 3 17:34:12 UTC 2014
1856 // Memory - virtual : 596 (mb) rss : 315 (mb)
1857 //
1858 //
1859 
1860  ProcInfo_t info;
1861  gSystem->GetProcInfo(&info);
1862 
1863  cout << "Memory - virtual : " << info.fMemVirtual / 1024
1864  << " (mb) rss : " << info.fMemResident / 1024 << " (mb)" << endl;
1865  return mvirtual ? info.fMemVirtual/1024 : info.fMemResident/1024;
1866 }
1867 
1868 void
1870 //
1871 // Compile & load the plugin
1872 // Load the configPlugin
1873 //
1874 // plugin.GetName() : is the plugin macro name
1875 // plugin.GetTitle() : is the macro path
1876 // if ends with .C the plugin code is compiled
1877 // if ends with .so then it is used as path for the pre-compiled lib
1878 // configPlugin.GetName() : is the configPlugin macro name
1879 //
1880 
1881  if(TString(plugin.GetName())=="") return;
1882 
1883  unsigned int Pid = gSystem->GetPid(); // used to tag in a unique way the temporary files
1884 
1885  TString pluginPath = plugin.GetTitle();
1886 
1887  int error,check;
1888 
1889  // check if configPlugin exists and is readable
1890  if(TString(configPlugin.GetName())!="") {
1891  // save original name/title into the macro file, is used in the InitHistory
1892  configPlugin.AddLine(TString::Format("//#?CONFIG_PLUGIN_NAME %s",configPlugin.GetName()));
1893  configPlugin.AddLine(TString::Format("//#?CONFIG_PLUGIN_TITLE %s",configPlugin.GetTitle()));
1894 
1895  // set a unique name for macro, it is used by configPlugin.Exec() to create a temporary file
1896  // extract file name
1897 
1898  char configPluginName[1024];
1899  if((gROOT->GetVersionInt()>=53400)&&(gROOT->GetVersionInt()<=53407)) {
1900  // starting from ROOT 5.34.00 to 5.34.07 TMacro::Exec() creates the temporary macro file under /tmp
1901  // the temporary configPlugin name must be defined without the directory name
1902  sprintf(configPluginName,"CWB_Plugin_Config_%s_%d_job%d.XXXXXX",
1903  cfg.data_label,Pid,this->runID);
1904 
1905  } else {
1906  sprintf(configPluginName,"%s/CWB_Plugin_Config_%s_%d_job%d.XXXXXX",
1907  cfg.nodedir,cfg.data_label,Pid,this->runID);
1908  }
1909  if(mkstemp(configPluginName)==-1) { // create temporary file, must be deleted only at the end of run
1910  cout << "cwb::LoadPlugin - Error : mkstemp error in creating tmp file : " << configPluginName << endl;
1911  EXIT(1);
1912  }
1913  configPlugin.SetName("CWB_PluginConfig");
1914  configPlugin.SetTitle(configPluginName);
1915  }
1916 
1917  // if the plugin path extension ends with 'so' the pre-compiled plugin is used
1918  check=0;
1919  if(pluginPath.EndsWith(".so")) {
1920  // check if the plugin compilation date is up to date
1921  // the plugin compilation date must be > wat compilation date
1922  FileStat_t fStatus;
1923  int estat = gSystem->GetPathInfo(pluginPath.Data(),fStatus);
1924  if(!estat) {
1925  wat::Time plugin_comp_date(double(fStatus.fMtime)); // plugin compilation date (unix)
1926  plugin_comp_date.UnixToGps(); // plugin compilation date (gps)
1927  wat::Time wat_comp_date(watversion('T')); // wat compilation date
1928  if(plugin_comp_date.GetGPS()<wat_comp_date.GetGPS()) {
1929  TString pluginSrc = pluginPath;
1930  pluginSrc.ReplaceAll("_C.so",".C");
1931  cout << endl;
1932  cout << "cwb::LoadPlugin : The plugin compilation date is not up to date !!! " << endl;
1933  cout << " plugin compilation date : " << plugin_comp_date.GetDateString() << endl;
1934  cout << " wat compilation date : " << wat_comp_date.GetDateString() << endl;
1935  cout << endl;
1936  cout << "Recompile the plugin : " << endl;
1937  cout << "cwb_compile " << pluginSrc << endl << endl;
1938  EXIT(1);
1939  }
1940  }
1941  cout << endl << "cwb::LoadPlugin - Load pre-compiled plugin ..." << endl << endl;
1942  check = gROOT->LoadMacro(pluginPath.Data(),&error);
1943  if(check==-1) {
1944  cout << "cwb::LoadPlugin : Load pre-compiled Plugin Failed !!! " << endl;
1945  cout << "cwb::LoadPlugin : The plugin is compiled 'On-The-Fly'!!! " << endl;
1946  }
1947  }
1948 
1949  // if the plugin path extension don't ends with 'so' or check!=0
1950  // then a copy is saved on tmp dir and the plugin is compiled 'On-The-Fly'
1951  if(!(pluginPath.EndsWith(".so"))||(check!=0)) {
1952  cout << endl << "cwb::LoadPlugin - compile and load plugin ..." << endl << endl;
1953  // create unique temporary file
1954  char tmpFile[1024]="";
1955  sprintf(tmpFile,"%s/CWB_Plugin_%s_%d_job%d.XXXXXX",cfg.nodedir,cfg.data_label,Pid,this->runID);
1956  if(mkstemp(tmpFile)==-1) { // create temporary file
1957  cout << "cwb::LoadPlugin - mkstemp error in creating tmp file : " << tmpFile << endl;
1958  EXIT(1);
1959  }
1960  char pluginSrc[1024]="";
1961  sprintf(pluginSrc,"%s.C",tmpFile); // must added .C to be compilable
1962 
1963  CWB::Toolbox::addCWBFlags(plugin,pluginSrc); // get macro + cWB compiler flags
1964 
1965  // compile & load plugin
1966  int success = gSystem->CompileMacro(TString(pluginSrc).Data(),"f");
1967  gSystem->Exec(TString("/bin/rm ")+tmpFile); // remove temporary file
1968  if(!success) {
1969  cout << "cwb::LoadPlugin : Plugin Compilation Failed !!! " << endl;
1970  EXIT(1);
1971  }
1972 
1973  // check if the compiled filename exists and is readable
1974  TString pluginShr = pluginSrc;
1975  pluginShr.ReplaceAll(".C","_C.so");
1976  check = gROOT->LoadMacro(pluginShr.Data(),&error,true);
1977  if(check==-1) {
1978  cout << "cwb::LoadPlugin : Plugin Compilation Failed !!! " << endl;
1979  EXIT(1);
1980  }
1981 
1982  // removed temporary plugin files & shared object created by ACLiC
1983  TString pluginTmp = pluginSrc;
1984  char tmpStr[1024];
1985 
1986  sprintf(tmpStr, "/bin/rm -f %s", pluginTmp.Data());
1987  system(tmpStr);
1988  pluginTmp.ReplaceAll(".C","_C.so");
1989  sprintf(tmpStr, "/bin/rm -f %s", pluginTmp.Data());
1990  system(tmpStr);
1991  pluginTmp.ReplaceAll("_C.so","_C.d");
1992  sprintf(tmpStr, "/bin/rm -f %s", pluginTmp.Data());
1993  system(tmpStr);
1994  }
1995 
1996  return;
1997 }
1998 
1999 void
2001 //
2002 // convert job_elapsed_time to (hh:mm:ss) format and print it
2003 //
2004 // job_elapsed_time : time (seconds)
2005 //
2006 // info : info string added to (hh:mm:ss)
2007 //
2008 
2009  int job_elapsed_hour = int(job_elapsed_time/3600);
2010  int job_elapsed_min = int((job_elapsed_time-3600*job_elapsed_hour)/60);
2011  int job_elapsed_sec = int(job_elapsed_time-3600*job_elapsed_hour-60*job_elapsed_min);
2012  char buf[1024];
2013  sprintf(buf,"%s %02d:%02d:%02d (hh:mm:ss)\n",info.Data(),job_elapsed_hour,job_elapsed_min,job_elapsed_sec);
2014  cout << buf;
2015 
2016  return;
2017 }
2018 
2019 void
2021 //
2022 // Print stage infos to screen/history
2023 //
2024 // stage : stage
2025 // comment : comment
2026 // out : true -> print to standard output
2027 // log : true -> save to history
2028 // fname : fname : job/output file name
2029 //
2030 
2031  TString info = GetStageInfo(stage,comment,fname);
2032  // print to standard output
2033  if(out) cout << info.Data() << endl;
2034  // save to history
2035  if(log&&history) history->AddLog(const_cast<char*>("FULL"), const_cast<char*>(info.Data()));
2036 }
2037 
2038 TString
2040 //
2041 // return string filled using comment and internal analysis status parameters
2042 //
2043 // stage : stage
2044 // comment : comment
2045 // fname : fname : job/output file name
2046 //
2047 // The return string has this format (Example):
2048 //
2049 // --------------------------------------------------------------------
2050 // cwb::PrintAnalysis
2051 // --------------------------------------------------------------------
2052 // UTC - 2014-02-03 21:17:17 UTC Mon
2053 // Job Elapsed Time - 00:00:09 (hh:mm:ss)
2054 // Stage Elapsed Time - 00:00:01 (hh:mm:ss)
2055 // Memory - virtual : 472 (mb) rss : 187 (mb)
2056 // Job File Size - 0 (bytes) : 0 (kb) : 0 (mb)
2057 // --------------------------------------------------------------------
2058 // GPS:1075497453-JOB:1-STG:1-FCT:-1-JET:9-SET:1-MEM:472-JFS:0
2059 // --------------------------------------------------------------------
2060 //
2061 // where :
2062 //
2063 // GPS : gps time (sec)
2064 // JOB : job number
2065 // STG : stage number (defined in cwb.hh)
2066 // FCT : factors index array
2067 // JET : Job Elapsed Time (sec)
2068 // SET : Stage Elapsed Time (sec)
2069 // MEM : Memory usage (MB)
2070 // JFS : Job File Size - temporary file (bytes)
2071 //
2072 
2073  size_t GPS; // GPS Time
2074  size_t JOB=runID; // Job ID
2075  size_t STG=stage; // stage ID
2076  int FCT; // factor ID
2077  size_t JET; // Job Elapsed Time
2078  size_t SET; // Stage Elapsed Time
2079  size_t MEM; // Virtual Memory
2080  size_t JFS; // Job File Size
2081 
2082  // get ifactor
2083  int gIFACTOR; IMPORT(int,gIFACTOR)
2084  FCT = gIFACTOR;
2085 
2086  // redirect cout to buffer and return contents
2087  std::stringstream buffer;
2088  std::streambuf * old = std::cout.rdbuf(buffer.rdbuf());
2089 
2090  // get stage info
2091  cout << endl;
2092  cout << "--------------------------------------------------------------------" << endl;
2093  cout << comment.Data();
2094  if(cfg.simulation&&(FCT>=0)) {
2095  cout << " - factor[" << FCT << "]=" << cfg.factors[FCT] << endl;
2096  } else cout<<endl;
2097  cout << "--------------------------------------------------------------------" << endl;
2098  cout << "UTC - "; cout.flush();
2099  wat::Time date("now"); cout << date.GetDateString() << endl; GPS=date.GetGPS();
2100  watchJob.Stop();
2101  JET = watchJob.RealTime();
2102  PrintElapsedTime(watchJob.RealTime(),"Job Elapsed Time - ");
2103  watchJob.Continue();
2104  watchStage.Stop();
2105  SET = watchStage.RealTime();
2106  PrintElapsedTime(watchStage.RealTime(),"Stage Elapsed Time - ");
2107  watchStage.Reset();
2108  watchStage.Start();
2109  MEM=GetProcInfo();
2110  Long_t xid,xsize,xflags,xmt;
2111  TString xname = fname=="" ? jname : fname;
2112  int xestat = gSystem->GetPathInfo(xname.Data(),&xid,&xsize,&xflags,&xmt);
2113  if(xestat==0) JFS=xsize; else JFS=0;
2114  cout << "Job File Size - " << JFS << " (bytes) : "
2115  << int(JFS/1024.) << " (kb) : " << int(JFS/1024./1024) << " (mb)" << endl;
2116  cout << "--------------------------------------------------------------------" << endl;
2117  cout << "GPS:" << GPS << "-JOB:" << JOB << "-STG:" << STG << "-FCT:" << FCT
2118  << "-JET:" << JET << "-SET:" << SET << "-MEM:" << MEM << "-JFS:" << JFS << endl;
2119  cout << "--------------------------------------------------------------------" << endl;
2120  cout << endl;
2121 
2122  // restore cout
2123  std::cout.rdbuf(old);
2124 
2125  return buffer.str();
2126 }
2127 
2128 void
2130 //
2131 // Print to screen/history the comment string
2132 //
2133 // stage : stage
2134 // comment : comment to be showed/saved
2135 // info : extra user info
2136 // out : true -> print to standard output
2137 // log : true -> save to history
2138 //
2139 
2140  TString ainfo = GetAnalysisInfo(stage,comment,info);
2141  // print to standard output
2142  if(out) cout << ainfo.Data() << endl;
2143  // save to history
2144  if(log&&history) history->AddLog(const_cast<char*>("FULL"), const_cast<char*>(ainfo.Data()));
2145 }
2146 
2147 TString
2149 //
2150 // return string filled using comment/info and internal analysis status parameters
2151 //
2152 // stage : stage
2153 // comment : comment
2154 // info : extra user info
2155 //
2156 // The return string has this format (Example):
2157 //
2158 // --------------------------------------------------------------------
2159 // GPS:1075497453-JOB:1-STG:1-FCT:-1-JET:9-SET:1-MEM:472-JFS:0
2160 // --------------------------------------------------------------------
2161 //
2162 // where :
2163 //
2164 // GPS : gps time (sec)
2165 // JOB : job number
2166 // STG : stage number (defined in cwb.hh)
2167 // FCT : factors index array
2168 // JET : Job Elapsed Time (sec)
2169 // SET : Stage Elapsed Time (sec)
2170 // MEM : Memory usage (MB)
2171 // JFS : Job File Size - temporary file (bytes)
2172 //
2173 
2174  size_t GPS; // GPS Time
2175  size_t JOB=runID; // Job ID
2176  size_t STG=stage; // stage ID
2177  int FCT; // factor ID
2178 
2179  // get ifactor
2180  TGlobal* global = gROOT->GetGlobal("gIFACTOR",true);
2181  if(global!=NULL) FCT = *(int*)global->GetAddress(); else FCT=-1;
2182 
2183  // redirect cout to buffer and return contents
2184  std::stringstream buffer;
2185  std::streambuf * old = std::cout.rdbuf(buffer.rdbuf());
2186 
2187  // get stage info
2188  cout << endl;
2189  cout << "--------------------------------------------------------------------" << endl;
2190  cout << comment.Data();
2191  if(cfg.simulation&&(FCT>=0)) {
2192  cout << " - factor[" << FCT << "]=" << cfg.factors[FCT] << endl;
2193  } else cout<<endl;
2194  cout << "--------------------------------------------------------------------" << endl;
2195  cout << "GPS:" << GPS << "-JOB:" << JOB << "-STG:" << STG << "-FCT:" << FCT << info << endl;
2196  cout << "--------------------------------------------------------------------" << endl;
2197  cout << endl;
2198 
2199  // restore cout
2200  std::cout.rdbuf(old);
2201 
2202  return buffer.str();
2203 }
2204 
2205 TString
2207 //
2208 // return stage string format
2209 //
2210 // jstage : is the stage number defined the enum CWB_STAGE
2211 //
2212 
2213  switch(jstage) {
2214  case CWB_STAGE_FULL :
2215  return("CWB_STAGE_FULL");
2216  break;
2217  case CWB_STAGE_INIT :
2218  return("CWB_STAGE_INIT");
2219  break;
2220  case CWB_STAGE_STRAIN :
2221  return("CWB_STAGE_STRAIN");
2222  break;
2223  case CWB_STAGE_CSTRAIN :
2224  return("CWB_STAGE_CSTRAIN");
2225  break;
2226  case CWB_STAGE_COHERENCE :
2227  return("CWB_STAGE_COHERENCE");
2228  break;
2229  case CWB_STAGE_SUPERCLUSTER :
2230  return("CWB_STAGE_SUPERCLUSTER");
2231  break;
2232  case CWB_STAGE_LIKELIHOOD :
2233  return("CWB_STAGE_LIKELIHOOD");
2234  break;
2235  default :
2236  return("");
2237  break;
2238  }
2239 }
2240 
2241 void
2243 //
2244 // set save options for job file which is used to pass informations between stages
2245 //
2246 // jstage : is the stage number defined the enum CWB_STAGE
2247 //
2248 
2250 
2251  switch(jstage) {
2252  case CWB_STAGE_INIT :
2257  break;
2258  case CWB_STAGE_STRAIN :
2265  break;
2266  case CWB_STAGE_CSTRAIN :
2272  break;
2273  case CWB_STAGE_COHERENCE :
2280  break;
2281  case CWB_STAGE_SUPERCLUSTER :
2287  break;
2288  case CWB_STAGE_LIKELIHOOD :
2290  break;
2291  case CWB_STAGE_FULL :
2292  // if option CWB_JOBF_SAVE_TRGFILE and jstage=CWB_STAGE_FULL then
2293  // all necessary objects are saved in the trigger file
2294  if(!(jobfOptions&CWB_JOBF_SAVE_TRGFILE)) break;
2301  break;
2302  default :
2303  break;
2304  }
2305 
2306  return;
2307 }
2308 
2309 void
2310 cwb::Exec(char* command, int maxtry, bool verbose) {
2311 //
2312 // Execute System Command
2313 //
2314 // command : system command
2315 // maxtry : number of retry before to exit with error (default=3)
2316 // verbose : show the retry infos (default=true)
2317 //
2318 
2319  int ntry=0;
2320  int ecommand=1;
2321  while(ecommand&&(ntry<maxtry)) {
2322  ecommand=gSystem->Exec(command);
2323  if(ecommand) gSystem->Sleep(int(gRandom->Uniform(10000,30000))); // wait [10,30] sec
2324  ntry++;
2325  if(verbose) {
2326  if(ecommand) {
2327  cout << command << endl;
2328  cout << "cwb::Exec - NTRY " << ntry << endl;
2329  cout.flush();
2330  }
2331  }
2332  }
2333  if(ecommand) {cout << "cwb::Exec - Error -> " << command << endl;EXIT(1);}
2334  return;
2335 }
2336 
2337 void
2338 cwb::FileGarbageCollector(TString ifName, TString ofName, vector<TString> delObjList) {
2339 //
2340 // for ROOT < 5.34.05
2341 // the objects marked as deleted do not free the disk space
2342 // this method is used to free the deleted object space contained in the root file
2343 //
2344 // for ROOT >= 5.34.05
2345 // All objected contained in ifName and are not in the delObjList are copied to ofName
2346 // This method is used to cleanup the ROOT file ifName
2347 //
2348 // ifName : input file name
2349 // ofName : output file name (cleaned) - if ofName="" -> ofName=ifName
2350 // delObjList : list of object which must removed from the ifName (only for ROOT >= 5.34.05)
2351 //
2352 
2353  bool replace = false;
2354  if(ofName=="") {
2355  ofName = ifName;
2356  ofName.ReplaceAll(".root","_tmp.root");
2357  replace = true;
2358  }
2359 
2360  if(delObjList.size()==0) {
2361  if(!replace) {
2362  char command[1024];
2363  sprintf(command,"/bin/mv %s %s",ofName.Data(),ifName.Data());
2364  gSystem->Exec(command);
2365  } else return;
2366  }
2367 
2368  if(gROOT->GetVersionInt()<53405) {
2369  // delete the objects contained in the delObjList
2370  TFile* kfile = new TFile(ifName, "UPDATE");
2371  for(int i=0;i<delObjList.size();i++) kfile->Delete(delObjList[i]+";1");
2372  kfile->Close();
2373  }
2374 
2375  gErrorIgnoreLevel=kBreak; // disable root kInfo & kWarning & kError messages
2376  TFileMerger M(false);
2377  M.AddFile(ifName,false);
2378  M.OutputFile(ofName);
2379  if(gROOT->GetVersionInt()<53400) {
2380  // TFileMerger free the deleted object space
2381  if(!M.Merge()) {
2382  cout << "cwb::FileGarbageCollector : Error - Merge failed !!!" << endl;EXIT(1);
2383  }
2384  }
2385  if(gROOT->GetVersionInt()>=53400 && gROOT->GetVersionInt()<53405) {
2386  // starting from ROOT 5.34.00 TFileMerger::Merge() merge also the deleted object
2387  // The following is a workaround to fix this issue
2388  if(!M.PartialMerge(TFileMerger::kAllIncremental | TFileMerger::kRegular)) {
2389  cout << "cwb::FileGarbageCollector : Error - Merge failed !!!" << endl;EXIT(1);
2390  }
2391  }
2392 #if ROOT_VERSION_CODE >= ROOT_VERSION(5,34,5)
2393  if(gROOT->GetVersionInt()>=53405) {
2394  // add files which must not be copied from ifName to ofName
2395  for(int i=0;i<delObjList.size();i++) M.AddObjectNames(delObjList[i]);
2396  // Must add new merging flag on top of the the default ones
2397  Int_t default_mode = TFileMerger::kAll | TFileMerger::kIncremental;
2398  Int_t mode = default_mode | TFileMerger::kSkipListed;
2399  if(!M.PartialMerge(mode)) {
2400  cout << "cwb::FileGarbageCollector : Error - Merge failed !!!" << endl;EXIT(1);
2401  }
2402  M.Reset();
2403  }
2404 #endif
2405  gErrorIgnoreLevel=kUnset; // re-enable root kInfo & kWarning messages
2406 
2407  if(replace) { // replace ifName with ofName
2408  char command[1024];
2409  sprintf(command,"/bin/mv %s %s",ofName.Data(),ifName.Data());
2410  gSystem->Exec(command);
2411  }
2412 
2413  return;
2414 }
2415 
2416 void
2417 cwb::MakeSkyMask(skymap& SkyMask, double theta, double phi, double radius) {
2418 //
2419 // make a sky mask
2420 // is used to define what are the sky locations that are analyzed
2421 //
2422 // theta,phi,radius : used to define the Celestial SkyMap
2423 // define a circle centered in (theta,phi) and radius=radius
2424 // theta : [-90,90], phi : [0,360], radius : degrees
2425 //
2426 // SkyMask : output sky celestial mask
2427 // inside the circle is filled with 1 otherwise with 0
2428 //
2429 
2430  int L = SkyMask.size();
2431  int healpix = SkyMask.getOrder();
2432  // check input parameters
2433  if(fabs(theta)>90 || (phi<0 || phi>360) || radius<=0 || L<=0) {
2434  cout << "cwb::MakeSkyMask : wrong input parameters !!! " << endl;
2435  if(fabs(theta)>90) cout << theta << " theta must be in the range [-90,90]" << endl;
2436  if(phi<0 || phi>360) cout << phi << " phi must be in the range [0,360]" << endl;
2437  if(radius<=0) cout << radius << " radius must be > 0" << endl;
2438  if(L<=0) cout << L << " SkyMask size must be > 0" << endl;
2439  EXIT(1);
2440  }
2441 
2442  if (!gROOT->GetClass("Polar3DVector")) gSystem->Load("libMathCore");
2443 
2444  // compute the minimun available radius
2445  // must be greater than the side of a pixel
2446  if(healpix) {
2447  int npix = 12*(int)pow(4.,(double)healpix);
2448  double sphere_solid_angle = 4*TMath::Pi()*pow(180./TMath::Pi(),2.);
2449  double skyres = sphere_solid_angle/npix;
2450  if(radius < sqrt(skyres)) radius = sqrt(skyres);
2451  } else {
2452  if(radius < SkyMask.sms) radius = SkyMask.sms;
2453  }
2454 
2455  double ph,th;
2456  GeographicToCwb(phi,theta,ph,th); // Geographic to CWB coordinates
2457 
2458  Polar3DVector ov1(1, th*TMath::Pi()/180, ph*TMath::Pi()/180);
2459 
2460  int nset=0;
2461  for (int l=0;l<L;l++) {
2462  double phi = SkyMask.getPhi(l);
2463  double theta = SkyMask.getTheta(l);
2464 
2465  Polar3DVector ov2(1, theta*TMath::Pi()/180, phi*TMath::Pi()/180 );
2466  double Dot = ov1.Dot(ov2);
2467  double dOmega = 180.*TMath::ACos(Dot)/TMath::Pi();
2468  //cout << "dOmega : " << dOmega << endl;
2469 
2470  if(dOmega<=radius) {SkyMask.set(l,1);nset++;} else SkyMask.set(l,0);
2471  }
2472 
2473  if(!nset) {
2474  cout << "cwb::MakeSkyMask : no sky positions setted !!! " << endl;
2475  cout << "check input mask parameters : theta = "
2476  << theta << " phi = " << phi << " radius : " << radius << endl << endl;
2477  EXIT(1);
2478  }
2479 
2480  return;
2481 }
2482 
2483 int
2484 cwb::SetSkyMask(network* net, CWB::config* cfg, char* options, char skycoord, double skyres) {
2485 //
2486 // set earth/celestial sky mask
2487 // is used to define what are the sky locations that are analyzed
2488 // by default all sky is used
2489 //
2490 // net : pointer to the network object
2491 // cfg : pointer to the cwb config object
2492 // options : used to define the earth/celestial SkyMap
2493 // option A :
2494 // skycoord='e'
2495 // --theta THETA --phi PHI --radius RADIUS
2496 // define a circle centered in (THETA,PHI) and radius=RADIUS
2497 // THETA : [-90,90], PHI : [0,360], RADIUS : degrees
2498 // skycoord='c'
2499 // --theta DEC --phi RA --radius RADIUS
2500 // define a circle centered in (DEC,RA) and radius=RADIUS
2501 // DEC : [-90,90], RA : [0,360], RADIUS : degrees
2502 // option B:
2503 // file name
2504 // format : two columns ascii file -> [sky_index value]
2505 // sky_index : is the sky grid index
2506 // value : if !=0 the index sky location is used for the analysis
2507 // skycoord : sky coordinates : 'e'=earth, 'c'=celestial
2508 // skyres : sky resolution : def=-1 -> use the value defined in parameters (angle,healpix)
2509 //
2510 // return 0 if successful
2511 
2512  if(skycoord!='e' && skycoord!='c') {
2513  cout << "cwb::SetSkyMask - Error : wrong input sky coordinates "
2514  << " must be 'e'/'c' earth/celestial" << endl;;
2515  EXIT(1);
2516  }
2517 
2518  if(strlen(options)>0) { // set SkyMask
2519  if(!TString(options).Contains("--")) { // input parameter is the skyMask file
2520  if(skyres>=0) return 1;
2521  int ret = net->setSkyMask(options,skycoord);
2522  if(ret==0) {
2523  cout << endl;
2524  cout << "cwb::SetSkyMask - Error : skyMask file"
2525  << " not exist or it has a wrong format" << endl;
2526  cout << endl;
2527  cout << " format : two columns ascii file -> [sky_index value]" << endl;
2528  cout << " sky_index : is the sky grid index" << endl;
2529  cout << " value : if !=0 the index sky location is used for the analysis" << endl;
2530  cout << endl;
2531  EXIT(1);
2532  }
2533  if(skycoord=='e') net->setIndexMode(0);
2534  return 0;
2535  }
2536 
2537  double theta=-1000;
2538  TString THETA = TB.getParameter(options,"--theta"); // get theta
2539  if(THETA.IsFloat()) theta=THETA.Atof();
2540  double phi=-1000;
2541  TString PHI = TB.getParameter(options,"--phi"); // get phi
2542  if(PHI.IsFloat()) phi=PHI.Atof();
2543  double radius=-1000;
2544  TString RADIUS = TB.getParameter(options,"--radius"); // get Radius
2545  if(RADIUS.IsFloat()) radius=RADIUS.Atof();
2546 
2547  if(theta==-1000 || phi==-1000 || radius==-1000) { // input parameter are the skyMask params
2548  cout << endl << "cwb::SetSkyMask - Error : wrong input skyMask params" << endl << endl;
2549  cout << "wrong input options : " << options << endl;
2550  if(skycoord=='e')
2551  cout<<"options must be : --theta THETA --phi PHI --radius RADIUS"<<endl<<endl;
2552  if(skycoord=='c')
2553  cout<<"options must be : --theta DEC --phi AR --radius RADIUS"<<endl<<endl;
2554  if(fabs(theta)>90) cout << theta << " theta must be in the range [-90,90]" << endl;
2555  if(phi<0 || phi>360) cout << phi << " phi must be in the range [0,360]" << endl;
2556  if(radius<=0) cout << radius << " radius must be > 0" << endl;
2557  cout << endl;
2558  EXIT(1);
2559  } else { // create & set SkyMask
2560  skymap* SkyMask=NULL;
2561  if(skyres<0) skyres = cfg->healpix ? cfg->healpix : cfg->angle;
2562  if(cfg->healpix) SkyMask = new skymap(int(skyres));
2563  else SkyMask = new skymap(skyres,cfg->Theta1,cfg->Theta2,cfg->Phi1,cfg->Phi2);
2564  MakeSkyMask(*SkyMask, theta, phi, radius);
2565  net->setSkyMask(*SkyMask,skycoord);
2566  if(skycoord=='e') net->setIndexMode(0);
2567  if(SkyMask!=NULL) delete SkyMask;
2568  }
2569  }
2570  return 0;
2571 }
2572 
2573 vector<frfile>
2575 //
2576 // return the FRF frame files list
2577 // input : ifoID : ifo iD , if ifoID=-1 then all ifo FRF are returned
2578 //
2579 
2580  vector<frfile> frlist;
2581 
2582  // fill frlist
2583  for(int i=0;i<2*nIFO;i++) {
2584  if((ifoID==-1)||(i==ifoID)||(i==(ifoID+nIFO))) frlist.push_back(FRF[i]);
2585  }
2586 
2587  return frlist;
2588 }
2589 
2590 vector<frfile>
2592 //
2593 // return the FRF frame files list
2594 // input : ifo : ifo name
2595 //
2596 
2597  // get ifo id
2598  int ifoID=-1;
2599  for(int n=0;n<nIFO;n++) if(ifo.CompareTo(this->ifo[n])==0) {ifoID=n;break;}
2600 
2601  if(ifoID==-1) {
2602  cout << "cwb::GetFrList - Error : requested ifo " << ifo
2603  << " not present !!!" << endl;;
2604  EXIT(1);
2605  }
2606 
2607  return GetFrList(ifoID);
2608 }
2609 
std::vector< char * > ifoName
Definition: network.hh:609
CWB_JOBF_OPTIONS jobfOptions
Definition: config.hh:291
void sethigh(double f)
Definition: wseries.hh:132
detector * getifo(size_t n)
param: detector index
Definition: network.hh:436
CWB::config * cfg
char analysis[8]
Definition: config.hh:117
TString outDump
int jobID
category 2 data quality list
Definition: cwb.hh:291
double Te
Definition: cwb.hh:282
TString GetLALVersion(TString options="")
Definition: Toolfun.hh:1004
double iwindow
Definition: config.hh:200
double netCC
Definition: network.hh:596
TString GetGitInfos(TString option="path", TString igit_path="$CWB_CONFIG")
Definition: Toolfun.hh:345
virtual void resize(unsigned int)
Definition: wseries.cc:901
CWB::frame fr[2 *NIFO_MAX]
Definition: cwb.hh:251
void PrintElapsedTime(int job_elapsed_time, TString info)
Definition: cwb.cc:2000
TString GetAnalysisInfo(CWB_STAGE stage, TString comment, TString info)
Definition: cwb.cc:2148
double startMDC
Definition: Toolbox.hh:106
cwb(CWB_STAGE jstage=CWB_STAGE_FULL)
Definition: cwb.cc:31
size_t nLag
Definition: network.hh:573
static size_t GetProcInfo(bool mvirtual=true)
Definition: cwb.cc:1847
double precision
Definition: config.hh:286
void printwc(size_t)
Definition: network.cc:2658
double start
Definition: network.hh:55
TString ofName
int job_elapsed_min
Definition: cwb_net.C:738
INT_4S GetGPS()
Definition: time.hh:125
TString GetDateString()
Definition: time.cc:461
WSeries< double > pixeLHood
Definition: network.hh:634
bool optim
Definition: config.hh:122
CWB::History * hst
void Export(TString fname="")
Definition: config.cc:406
double M
Definition: DrawEBHH.C:13
void Init()
Definition: ChirpMass.C:284
void cwb_jnet(TString jName="", TString uName="")
Definition: cwb_jnet.C:22
void FileGarbageCollector(TString ifName, TString ofName="", vector< TString > delObjList=vector< TString >())
Definition: cwb.cc:2338
int stop
Definition: Toolbox.hh:94
void Print(Option_t *option="")
Definition: config.cc:737
size_t add(detector *)
param: detector structure return number of detectors in the network
Definition: network.cc:2559
bool singleDetector
used for the stage stuff
Definition: cwb.hh:277
TString live
TMacro configPlugin
Definition: config.hh:362
size_t readMDClog(char *, double=0., int=11, int=12)
param: MDC log file param: approximate gps time
Definition: network.cc:3370
void setAntenna(detector *)
param: detector (use theta, phi index array)
Definition: network.cc:2846
bool mdcPlugin
Definition: config.hh:365
Definition: ced.hh:42
int nfrFiles[2 *NIFO_MAX]
Definition: cwb.hh:252
virtual void rate(double r)
Definition: wavearray.hh:141
float factor
#define EXIT(ERR)
Definition: cwb.cc:25
char name[32]
Definition: detector.hh:50
char skyMaskFile[1024]
Definition: config.hh:308
int error
Definition: cwb_compile.C:43
WSeries< float > v[nIFO]
Definition: cwb_net.C:80
CWB SetSkyMask(net, cfg, cfg->skyMaskCCFile, 'c')
double job_speed_factor
Definition: cwb_net.C:741
bool dataPlugin
Definition: config.hh:364
int n
Definition: cwb_net.C:28
void Init()
Definition: cwb.cc:230
double dT
Definition: cwb.hh:283
void CWB_Plugin(TFile *jfile, CWB::config *, network *, WSeries< double > *, TString, int)
COHERENCE.
vector< frfile > GetFrList(int ifoID=-1)
Definition: cwb.cc:2574
char * watversion(char c='s')
Definition: watversion.hh:21
TString("c")
char comment[1024]
size_t setIndexMode(size_t=0)
Definition: network.cc:8105
int nset
Definition: cwb_mkeff.C:75
double Tb
Definition: cwb.hh:281
virtual double InitJob()
Definition: cwb.cc:1248
int start
Definition: Toolbox.hh:93
double stopMDC
Definition: Toolbox.hh:107
ofstream out
Definition: cwb_merge.C:214
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<< "\ layers : "<< nLAYERS<< "\ dF(hz) : "<< dF<< "\ 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
void Check()
Definition: config.cc:1411
std::vector< wavearray< double > * > RWFP
Definition: detector.hh:382
double fLow
Definition: config.hh:140
float theta
int pattern
Definition: config.hh:241
bool eDisbalance
Definition: network.hh:598
netcluster * pwc
Definition: cwb_job_obj.C:38
TString mdc[4]
CWB::frame fr(FRLIST_NAME)
TH2F * ph
CWB::Toolbox TB
vector< int > segId
Definition: Toolbox.hh:102
CWB_STAGE jstage
Definition: cwb.hh:248
bool cedDump
Definition: config.hh:297
void setlow(double f)
Definition: wseries.hh:125
int slagID
Definition: cwb.hh:291
size_t setSkyMask(double f, char *fname)
Definition: network.cc:3213
Long_t flags
CWB_JOBF_OPTIONS jobfOptions
TString GetStageInfo(CWB_STAGE stage, TString comment, TString fname="")
Definition: cwb.cc:2039
virtual double ReadData(double mdcShift, int ifactor)
Definition: cwb.hh:201
dqfile DQF[DQF_MAX]
Definition: config.hh:349
CWB::mdc * MDC
size_t * slagSite
Definition: config.hh:185
double getTheta(size_t i)
Definition: skymap.hh:224
double Te
double netRHO
Definition: config.hh:147
CWB::History * history
char * lagFile
Definition: config.hh:172
Long_t size
int slagMin
Definition: config.hh:182
static double getTimeSegList(vector< waveSegment > list)
Definition: Toolbox.cc:611
char command[1024]
Definition: cwb_compile.C:44
virtual void start(double s)
Definition: wavearray.hh:137
vector< waveSegment > detSegs_dq2
Definition: cwb_net.C:297
int j
Definition: cwb_net.C:28
i drho i
std::vector< double > mdcTime
Definition: network.hh:614
void Import(TString umacro="")
Definition: config.cc:352
char nodedir[1024]
Definition: config.hh:352
double segTHR
Definition: config.hh:163
static bool checkFile(TString fName, bool question=false, TString message="")
Definition: Toolbox.cc:4670
void setRunID(size_t n)
param: run
Definition: network.hh:455
double Phi1
Definition: config.hh:279
int jobID
Definition: cwb_net.C:195
cout<< "SNR "<< xsnr<< endl;wavearray< double > f
Definition: ComputeSNR.C:75
bool outPlugin
Definition: config.hh:369
int jobId
Definition: Toolbox.hh:100
double segOverlap
Definition: config.hh:165
char ifo[NIFO_MAX][8]
double Edge
Definition: network.hh:578
int slagOff
Definition: config.hh:184
network ** net
NOISE_MDC_SIMULATION.
size_t ifoListSize()
Definition: network.hh:431
waveSegment getFrRange()
Definition: frame.hh:107
bool online
Definition: config.hh:118
static char * getEnvCWB()
Definition: Toolbox.cc:5383
double mdcShift
Definition: cwb_net.C:196
void constraint(double d=1., double g=0.0001)
param: constraint parameter, p=0 - no constraint
Definition: network.hh:463
char * rootlogonBuffer
Definition: cwb_net.C:158
int Psave
Definition: config.hh:193
vector< waveSegment > detSegs
time delay difference
Definition: cwb.hh:288
bool singleDetector
bool optim
Definition: network.hh:590
size_t mode
double segEdge
Definition: config.hh:164
#define nIFO
void UnixToGps()
Definition: time.hh:231
virtual ~cwb()
Definition: cwb.cc:215
void setTimeRange(int xstart=0, int xstop=0)
Definition: frame.hh:148
float slagShift[20]
Definition: cwb.hh:292
detectorParams detParms[NIFO_MAX]
Definition: config.hh:128
double Theta2
Definition: config.hh:278
tlive_fix Write()
float phi
double precision
Definition: network.hh:593
size_t fResample
Definition: config.hh:142
TGlobal * global
Definition: config.cc:23
static char * readFile(TString ifName)
Definition: Toolbox.cc:2835
jfile
Definition: cwb_job_obj.C:43
char end_CED[1024]
Definition: cwb_job_obj.C:35
TTree * net_tree
char filter_dir[1024]
Definition: config.hh:305
double shift
Definition: netcluster.hh:382
size_t lagOff
Definition: config.hh:170
Definition: mdc.hh:248
char lagMode[2]
Definition: test_config1.C:57
std::vector< double > livTime
Definition: network.hh:611
static TString GetStageString(CWB_STAGE jstage)
Definition: cwb.cc:2206
int simulation
Definition: config.hh:199
int l_high
Definition: config.hh:156
int length
Definition: Toolbox.hh:95
int nDQF
Definition: config.hh:348
double Fgap
Definition: config.hh:136
int pattern
void PrintAnalysisInfo(CWB_STAGE stage, TString comment, TString info, bool out=true, bool log=true)
Definition: cwb.cc:2129
void SetupStage(CWB_STAGE jstage)
Definition: cwb.cc:2242
CWB_OUTF_OPTIONS outfOptions
Definition: config.hh:292
#define IMPORT(TYPE, VAR)
Definition: cwb.hh:69
int segID[20]
Definition: cwb.hh:291
char search
Definition: config.hh:121
size_t * lagSite
Definition: config.hh:174
vector< slag > slagList
char * GetLog(char *Stage, int index)
Definition: History.cc:304
i() int(T_cor *100))
int GetLogSize(char *Stage)
Definition: History.cc:298
char output_dir[1024]
Definition: config.hh:318
network NET
Definition: cwb_dump_inj.C:30
double Pi
int getOrder()
Definition: skymap.hh:314
virtual void run(int runID=0)
Definition: cwb.cc:281
int ifactor
std::vector< std::string > mdcList
Definition: network.hh:612
const int NIFO_MAX
Definition: wat.hh:22
bool log
Definition: WaveMDC.C:41
double getDelay(const char *c="")
Definition: network.cc:2818
char lagMode[1]
Definition: cwb.hh:294
int slagMax
Definition: config.hh:183
WSeries< double > pTF[nRES]
Definition: revMonster.cc:8
static vector< waveSegment > getSlagJobList(vector< waveSegment > ilist, int seglen=600)
Definition: Toolbox.cc:1791
printf("total live time: non-zero lags = %10.1f \, liveTot)
int l
double delta
Definition: config.hh:252
int mlagStep
double sms
Definition: skymap.hh:320
char fname[1024]
bool eDisbalance
Definition: config.hh:270
double dataShift[NIFO_MAX]
Definition: config.hh:205
void setVerbose(bool verbose=true)
Definition: frame.hh:137
#define RADIUS
Definition: Wavelet.hh:49
double precision
int slagSize
Definition: config.hh:181
frfile FRF[nIFO+1]
Definition: cwb_net.C:269
double Acore
Definition: config.hh:143
void PrintStageInfo(CWB_STAGE stage, TString comment, bool out=true, bool log=true, TString fname="")
Definition: cwb.cc:2020
TList * GetStageNames()
Definition: History.cc:427
void Exec(char *command, int maxtry=3, bool verbose=true)
Definition: cwb.cc:2310
int k
int frRetryTime
Definition: config.hh:344
int pattern
Definition: network.hh:601
int getNfiles()
Definition: frame.hh:110
std::vector< waveSegment > segList
Definition: network.hh:616
char log_dir[1024]
Definition: config.hh:323
int nfactor
Definition: config.hh:201
virtual void DataConditioning(int ifactor)
Definition: cwb.hh:205
std::vector< wavearray< double > * > IWFP
Definition: detector.hh:379
void PrintElapsedTime(int job_elapsed_time, double cpu_time, TString info)
size_t healpix
Definition: config.hh:282
char refIFO[4]
Definition: config.hh:125
int nfrFiles[nIFO+1]
Definition: cwb_net.C:192
vector< waveSegment > cat1List
size_t lagSize
Definition: config.hh:168
void setSRIndex(int srIndex)
Definition: frame.hh:114
void setDelay(const char *="L1")
Definition: network.cc:2767
Definition: skymap.hh:63
long nSky
Definition: config.hh:299
cout<< endl;cout<< "Unfinished Jobs : "<< cnt<< "/"<< jobList.size()<< endl;cout<< endl;sprintf(dagfile,"%s/%s.dag.recovery.%d", condor_dir, data_label, iversion);cout<< "To submit condor recovered jobs, type :"<< endl;cout<< "cwb_condor submit "<< dagfile<< endl;} cout<< endl;if(gSystem->Getenv("_USE_LSF")!=NULL) { TString cwb_stage_label="supercluster_";if(cwb_stage_input=="FULL") cwb_stage_label="wave_";if(cwb_stage_input=="INIT") cwb_stage_label="init_";if(cwb_stage_input=="STRAIN") cwb_stage_label="strain_";if(cwb_stage_input=="CSTRAIN") cwb_stage_label="cstrain_";if(cwb_stage_input=="COHERENCE") cwb_stage_label="coherence_";if(cwb_stage_input=="SUPERCLUSTER") cwb_stage_label="supercluster_";if(cwb_stage_input=="LIKELIHOOD") cwb_stage_label="wave_";int jobID=1;TString exec_cmd=TString::Format("export file_n_st=""$(ls %s*_job%i.root)""", cwb_stage_label.Data(), jobID);gSystem-> Exec(exec_cmd)
TFile * ifile
frfile FRF[2 *NIFO_MAX]
Definition: cwb.hh:253
TList * GetTypeNames()
Definition: History.cc:440
char lagMode[2]
Definition: config.hh:173
size_t rateANA
Definition: test_config1.C:21
char injectionList[1024]
Definition: config.hh:307
double netCC
Definition: config.hh:148
bool EFEC
Definition: config.hh:274
TFile * froot
void setAcore(double a)
Definition: network.hh:401
int npix
static double getMDCShift(mdcshift mshift, double time)
Definition: Toolbox.cc:2296
double subnet
Definition: config.hh:247
mdcshift mdc_shift
Definition: config.hh:211
static slag getSlag(vector< slag > slagList, int jobid)
Definition: Toolbox.cc:1844
size_t events()
Definition: network.hh:329
static int shiftBurstMDCLog(std::vector< std::string > &mdcList, vector< TString > ifos, double mdc_shift)
Definition: Toolbox.cc:2174
vector< waveSegment > detSegs
Definition: cwb_net.C:193
static void setSlagShifts(slag SLAG, vector< TString > ifos, double segLen, int nDQF, dqfile *DQF)
Definition: Toolbox.cc:2017
int l_low
Definition: config.hh:155
char options[256]
#define EXPORT(TYPE, VAR, CMD)
Definition: config.cc:92
static vector< slag > getSlagList(size_t nIFO, size_t slagSize, int slagSegs, int slagOff, size_t slagMin, size_t slagMax, size_t *slagSite, char *slagFile=NULL)
Definition: Toolbox.cc:809
double getPhi(size_t i)
Definition: skymap.hh:164
static vector< waveSegment > readSegList(dqfile DQF)
Definition: Toolbox.cc:409
char cmd_line[512]
Definition: cwb_net.C:154
int job_elapsed_time
Definition: cwb_net.C:736
double Tb
static TString getParameter(TString options, TString param="")
Definition: Toolbox.cc:6727
ifstream in
char * slagFile
Definition: config.hh:186
double fHigh
Definition: config.hh:141
void PrintAnalysis(bool stageInfos=true)
Definition: cwb.cc:1705
int mlagStep
Definition: config.hh:178
WSeries< double > * getTFmap()
param: no parameters
Definition: detector.hh:179
static void MakeSkyMask(skymap &SkyMask, double theta, double phi, double radius)
Definition: cwb.cc:2417
char Name[16]
Definition: detector.hh:327
virtual void InitHistory()
Definition: cwb.cc:1002
TMacro plugin
Definition: config.hh:361
char ifo[NIFO_MAX][8]
Definition: config.hh:124
vector< int > slagId
Definition: Toolbox.hh:101
double fabs(const Complex &x)
Definition: numpy.cc:55
char work_dir[1024]
Definition: config.hh:314
int job_data_size_sec
Definition: cwb_net.C:740
int SetSkyMask(network *net, CWB::config *cfg, char *options, char skycoord, double skyres=-1)
Definition: cwb.cc:2484
int gIFACTOR
double gamma
Definition: network.hh:592
double segLen
Definition: config.hh:161
#define GPS
double mTau
int slagID
Definition: cwb_net.C:195
static vector< waveSegment > getSegList(int jobId, int nIFO, double segLen, double segMLS, double segEdge, vector< waveSegment > dqList)
Definition: Toolbox.cc:2799
char cmd[1024]
int estat
void AddLog(char *Stage, char *Log, TDatime *Time=NULL)
Definition: History.cc:187
void LoadPlugin(TMacro &plugin, TMacro &configPlugin)
Definition: cwb.cc:1869
strcpy(RunLabel, RUN_LABEL)
Definition: Toolbox.hh:99
cout<< "total cat1 livetime : "<< int(cat1_time)<< " sec "<< cat1_time/3600.<< " h "<< cat1_time/86400.<< " day"<< endl;cout<< endl;vector< waveSegment > cat2List
Definition: cwb_dump_job.C:40
netcluster * getwc(size_t n)
param: delay index
Definition: network.hh:439
sprintf(tfres,"(1/%g)x(%g) (sec)x(Hz)", 2 *df, df)
double subcut
Definition: config.hh:248
double Phi2
Definition: config.hh:280
Long_t mt
frfile getFrList(int istart, int istop, int segEdge)
Definition: frame.cc:527
TMacro plugin
void set(size_t i, double a)
param: sky index param: value to set
Definition: skymap.hh:122
static vector< waveSegment > mergeSegLists(vector< waveSegment > &ilist1, vector< waveSegment > &ilist2)
Definition: Toolbox.cc:350
char skyMaskCCFile[1024]
Definition: config.hh:309
int nIFO
Definition: config.hh:120
Definition: cwb.hh:136
double mask
Definition: config.hh:281
double Tgap
Definition: config.hh:134
char nodedir[1024]
Definition: test_config1.C:187
SortOrderType SetSortOrder(SortOrderType SortOrder)
Definition: History.cc:589
long nSky
Definition: network.hh:574
TString jname
Long_t id
double detSegs_ctime
Definition: cwb_net.C:303
enum WAVETYPE m_WaveType
Definition: Wavelet.hh:106
double factors[FACTORS_MAX]
Definition: config.hh:202
int job_elapsed_hour
Definition: cwb_net.C:737
TMacro configPlugin
int segID[20]
Definition: cwb_net.C:195
WaveDWT< DataType_t > * pWavelet
Definition: wseries.hh:456
double dT
Definition: testWDM_5.C:12
void GeographicToCwb(double ilongitude, double ilatitude, double &olongitude, double &olatitude)
Definition: skycoord.hh:421
size_t lagMax
Definition: config.hh:171
char * GetHistory(char *StageName, char *Type)
Definition: History.cc:273
TArrayC lagBuffer
Definition: cwb.hh:295
char data_label[1024]
Definition: config.hh:332
void setRetryTime(int frRetryTime=60)
Definition: frame.hh:141
cout<< fr.getNfiles()<< endl;std::vector< frfile > frlist
virtual void SuperCluster(int ifactor)
Definition: cwb.hh:213
void open(TString ioFile, TString chName="", Option_t *option="", bool onDisk=false, TString label=".gwf", unsigned int mode=0)
Definition: frame.cc:230
int runID
Definition: cwb.hh:246
double lagStep
Definition: config.hh:169
virtual void Coherence(int ifactor)
Definition: cwb.hh:209
size_t size()
Definition: skymap.hh:136
int job_elapsed_sec
Definition: cwb_net.C:739
char fName[256]
virtual void resize(unsigned int)
Definition: wavearray.cc:463
double TFgap
Definition: config.hh:138
double bpp
Definition: config.hh:133
int check
int setTimeShifts(size_t=1, double=1., size_t=0, size_t=0, const char *=NULL, const char *="w", size_t *=NULL)
param number of time lags param time shift step in seconds param first lag ID param maximum lag ID pa...
Definition: network.cc:7321
double length
Definition: TestBandPass.C:18
double Theta1
Definition: config.hh:277
double stop
Definition: network.hh:56
double netRHO
Definition: network.hh:597
size_t inRate
Definition: config.hh:132
double gamma
Definition: config.hh:258
void GetProcInfo(TString str="")
Definition: Toolfun.hh:241
double segMLS
Definition: config.hh:162
void AddHistory(char *Stage, char *Type, char *History, TDatime *Time=NULL)
Definition: History.cc:225
int levelR
Definition: config.hh:152
bool EFEC
Definition: network.hh:587
char frFiles[2 *NIFO_MAX][1024]
Definition: config.hh:342
detector ** pD
virtual void InitNetwork()
Definition: cwb.cc:910
CWB_STAGE
Definition: cwb.hh:122
void SetSingleDetectorMode()
Definition: config.cc:1352
TString ifos[60]
double angle
Definition: config.hh:276
void setSkyMaps(double, double=0., double=180., double=0., double=360.)
param: sky map granularity step, degrees param: theta begin, degrees param: theta end...
Definition: network.cc:2687
static TString addCWBFlags(TMacro macro, TString ofname="")
Definition: Toolbox.cc:4830
size_t healpix
bool dump
Definition: config.hh:295