Logo Coherent WaveBurst  
Reference Guide
Logo
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
frame.cc
Go to the documentation of this file.
1 #include "frame.hh"
2 #include "Meyer.hh"
3 
4 ////////////////////////////////////////////////////////////////////////////////
5 /* BEGIN_HTML
6 The frame class is a partial c++ wrapper to the c frame library<br>
7 framelib : (http://lappweb.in2p3.fr/virgo/FrameL/)<br>
8 it is designed to implements only the functions used by cWB<br>
9 <br>
10 the following are some examples which illustrate how to use the class
11 
12 <h3><a name="example1">Example 1</a></h3>
13 This example shows how to read a frame from frame file list<br>
14 using the frfile structure
15 <br><br>
16 <pre>
17  root[] #define FRLIST_NAME "input/H1_LDAS_C02_L2.frl" // frame files list name
18  root[] #define CHANNEL_NAME "H1:LDAS-STRAIN" // channel to be read
19  root[] double start = 942449664; // start time
20  root[] double stop = start+600; // stop time
21  root[] wavearray<double> x; // array
22  root[] CWB::frame fr(FRLIST_NAME); // open frame file list
23  root[] int nfiles=fr.getNfiles(); // get number of files in frame list
24  root[] cout << "nfiles : " << nfiles << endl;
25  root[] frfile FRF = fr.getFrList(start, stop, 0); // fill the frfile structure
26  root[] fr.readFrames(FRF,CHANNEL_NAME,x); // read data in x array
27 </pre>
28 <h3><a name="example2">Example 2</a></h3>
29 This example shows how to read a frame from frame file list using the read operator ">>"
30 <br><br>
31 <pre>
32  root[] #define FRLIST_NAME "input/H1_LDAS_C02_L2.frl" // frame files list name
33  root[] #define CHANNEL_NAME "H1:LDAS-STRAIN" // channel to be read
34  root[] wavearray<double> x; // array
35  root[] x.start(942449664); // start time
36  root[] x.stop(x.start()+600); // stop time
37  root[] frame ifr(FRLIST_NAME,CHANNEL_NAME,"READ"); // open frame file list
38  root[] ifr >> x; // read data in x array
39  root[] cout << x.start() << " " << x.rate() << " " << x.size() << endl; // print data info
40 </pre>
41 <h3><a name="example3">Example 3</a></h3>
42 This example shows how to write data to frame file using the write operator "<<"
43 <br><br>
44 <pre>
45  root[] #define CHANNEL_NAME "H1:LDAS-STRAIN" // channel to be read
46  root[] #define FRNAME "TEST" // output frame name
47  root[] wavearray<double> x; // array
48  root[] x.start(942449664); // set start time
49  root[] x.stop(x.start()+600); // set stop time
50  root[] x.rate(16384); // set rate (Hz)
51  root[] for(int i=0;i&lt;x.size();i++) x[i]=i; // fill array
52  root[] frame ofr("TEST.gwf","CHNAME","WRITE"); // open output frame file
53  root[] ofr.setFrName(FRNAME); // set frame name
54  root[] ofr << x; // read data in x array
55  root[] ofr.close(); // close output frame file
56 </pre>
57 <h3><a name="example4">Example 4</a></h3>
58 This example shows how to resample read data
59 <br><br>
60 <pre>
61  root[] #define FRLIST_NAME "input/H1_LDAS_C02_L2.frl" // frame files list name
62  root[] #define CHANNEL_NAME "H1:LDAS-STRAIN" // channel to be read
63  root[] wavearray<double> x; // array
64  root[] x.start(942449664); // start time
65  root[] x.stop(x.start()+600); // stop time
66  root[] frame ifr(FRLIST_NAME,CHANNEL_NAME,"READ"); // open frame file list
67  root[] ifr.setSRIndex(10); // read data are resampled @ 1024 Hz
68  root[] ifr >> x; // read data in x array
69  root[] cout << x.start() << " " << x.rate() << " " << x.size() << endl; // print data info
70 </pre>
71 
72 
73 END_HTML */
74 ////////////////////////////////////////////////////////////////////////////////
75 
76 ClassImp(CWB::frame)
77 
78 //______________________________________________________________________________
79 CWB::frame::frame() : frtree_List(NULL), chName(""), frName("Frame"),
80  rfName(""), frFile(NULL), fOption(""), srIndex(0),
81  verbose(false), frRetryTime(60), xstart(0.), xstop(0.) {
82 // default constructor
83 // initialize only the frame class parameters
84 //
85 }
86 
87 //______________________________________________________________________________
88 CWB::frame::frame(TString ioFile, TString chName, Option_t* option, bool onDisk, TString label, unsigned int mode) :
89  frtree_List(NULL), chName(""), frName("Frame"), rfName(""),
90  frFile(NULL), fOption(""), srIndex(0), verbose(false), frRetryTime(60), xstart(0.), xstop(0.) {
91 //
92 // initialize frame class parameters & open frame file list
93 //
94 // see the 'CWB::frame::open' method description
95 //
96  open(ioFile, chName, option, onDisk, label, mode);
97 }
98 
99 //______________________________________________________________________________
101 //
102 // destructor
103 //
104 
105  // remove temporary root file
106  if(rfName!="") gSystem->Exec(TString("rm "+rfName).Data());
107  close();
108 }
109 
110 //______________________________________________________________________________
112 //
113 // operator : read frame into wavearray x
114 //
115 
116  fr.readFrames(x);
117  return x;
118 }
119 
120 //______________________________________________________________________________
121 CWB::frame& operator << (CWB::frame& fr, wavearray<double>& x) {
122 //
123 // operator : write wavearray x to frame
124 //
125 
126  if(fr.getOption()=="READ") {
127  cout << "CWB::frame::operator(<<) : Not allowed in READ mode" << endl;
128  exit(1);
129  }
130 
131  if(fr.getChName()=="") {
132  cout << "CWB::frame::operator(<<) : channel not defined";
133  exit(1);
134  }
135 
136  if(fr.getFrName()=="") {
137  cout << "CWB::frame::operator(<<) : frame name not defined" << endl;
138  exit(1);
139  }
140 
142  return fr;
143 }
144 
145 //______________________________________________________________________________
147 //
148 // operator : write wavearray x to frame
149 //
150 
151  fr << x;
152  return fr;
153 }
154 
155 //______________________________________________________________________________
156 void
158 //
159 // write in frame the data contained in the wavearray x
160 // frame is saved in the frFile (object must be opened in write mode)
161 //
162 // x : in - data array which contains the data
163 // chName : in - name of channel to be extracted
164 // frName : in - frame name
165 //
166 
167  if (frFile==NULL) {
168  cout << "CWB::frame::writeFrame : output frame file close" << endl;
169  exit(1);
170  }
171  if(x.rate()<=0) {
172  cout << "CWB::frame::writeFrame : input rate must be > 0" << endl;
173  exit(1);
174  }
175  if(x.size()==0) {
176  cout << "CWB::frame::writeFrame : input size must be > 0" << endl;
177  exit(1);
178  }
179 
180  /*----------------------- Create a new frame ---------*/
181  FrameH* frFrame = FrameNew(const_cast<char*>(frName.Data()));
182  frFrame->frame = 0;
183  frFrame->run = -1;
184  frFrame->dt = x.size()/x.rate();
185  frFrame->GTimeS = x.start();
186  frFrame->GTimeN = 0;
187 
188  cout << "Size (sec) " << x.size()/x.rate() << endl;
189  FrProcData* proc = FrProcDataNew(frFrame,const_cast<char*>(chName.Data()),x.rate(),x.size(),-64);
190  if(proc == NULL) {
191  cout << "CWB::frame::writeFrame : Cannot create FrProcData" << endl;
192  exit(1);
193  }
194  proc->timeOffset = 0;
195  proc->tRange = frFrame->dt;
196  proc->type = 1; // Time Serie
197 
198  for (int i=0;i<(int)proc->data->nData;i++) proc->data->dataD[i] = x[i];
199 
200  int err=FrameWrite(frFrame,frFile);
201  if (err) {
202  cout << "CWB::frame::writeFrame : Error writing frame" << endl;
203  exit(1);
204  }
205  FrameFree(frFrame);
206 
207  return;
208 }
209 
210 //______________________________________________________________________________
211 void
212 CWB::frame::open(TString ioFile, TString chName, Option_t* option, bool onDisk, TString label, unsigned int mode) {
213 //
214 // open frame file in read or write mode
215 //
216 // ioFile : in - if option="READ" -> ioFile = input frame file list name
217 // out - if option="WRITE" -> ioFile = output frame file name
218 // chName : in - name of the data channel (def="")
219 // option : in - READ/WRITE mode (def="")
220 // onDisk : in - true -> tmp tree on root file
221 // false (def) -> tmp tree is stored in ram
222 // true must be used only if the input file frame list is huge
223 // label : in - this is a string used to select files listed in the input frame file list
224 // def = ".gwf"
225 // mode : in - if mode=0 (default) the frame file list is a list of frame file paths
226 //
227 // has the following format :
228 // DIR_NAME/XXXXX-GPS-LEN.gwf where :
229 // DIR_NAME : is the directory which contains the frame file.
230 // The strings 'file://localhost' or 'framefile=' at the
231 // beginning of the path will be automatically removed
232 // GPS : is the start frame time in sec
233 // LEN : is the frame length time in sec
234 //
235 // if mode=1 the frame file list is a list of line with the following format :
236 //
237 // FILE_PATH GPS LEN
238 //
239 
240  if(rfName!="") gSystem->Exec(TString("rm "+rfName).Data());
241  close();
242 
243  this->chName=chName;
244  this->frName="Frame";
245 
246  // if onDisk store tree on root file otherwise it is stored in ram
247  if(onDisk) {
248  // create temporary file
249  gRandom->SetSeed(0);
250  int rnID = int(gRandom->Rndm(13)*1.e9);
251  UserGroup_t* uinfo = gSystem->GetUserInfo();
252  TString uname = uinfo->fUser;
253  // create temporary dir
254  gSystem->Exec(TString("mkdir -p /dev/shm/")+uname);
255  char fName[1024]="";
256  sprintf(fName,"/dev/shm/%s/cwbframe_%d.root",uname.Data(),rnID);
257  rfName=fName;
258  }
259 
260  fOption = option;
261  fOption.ToUpper();
262  if((fOption!="READ")&&(fOption!="WRITE")) fOption="READ";
263 
264  if(fOption=="READ") nfiles = frl2FrTree(ioFile, rfName, label, mode);
265  if(fOption=="WRITE") {
266  if(!fNameCheck(ioFile)) exit(1);
267  frFile = FrFileONew(const_cast<char*>(ioFile.Data()),1); // gzip compression
268  if(frFile==NULL) {
269  cout << "CWB::frame::frame : Error opening file " << ioFile.Data() << endl;
270  exit(1);
271  }
272  }
273  return;
274 }
275 
276 //______________________________________________________________________________
277 void
279 //
280 // close frame file and reset parameters
281 //
282 
283  if(frtree_List!=NULL) frtree_List->Reset();
284  frtree_List=NULL;
285  if (frFile!=NULL) FrFileOEnd(frFile);
286  frFile=NULL;
287  fOption="";
288  chName="";
289  frName="";
290 }
291 
292 //______________________________________________________________________________
293 int
294 CWB::frame::frl2FrTree(TString iFile, TString rfName, TString label, unsigned int mode) {
295 //
296 // convert input frame file list to a tree
297 // note : this method can be used only in READ mode
298 //
299 // iFile : in - input frame file list name
300 // rfName : in - name of auxiliary file name used to store the temporary tree
301 // if rfName=="" the tree is saved & sorted to the private frtree_List
302 // label : in - this is a string used to select files listed in the input frame file list
303 // def = ".gwf"
304 // mode : in - 0/1 (def=0)
305 //
306 
307  if(fOption!="READ") {
308  cout << "CWB::frame::frl2FrTree : allowed only in READ mode" << endl;
309  exit(1);
310  }
311 
312  TFile* ofile = NULL;
313  if(rfName.Sizeof()>1) ofile=new TFile(rfName.Data(),"RECREATE");
314 
315  TTree* frtree = new TTree("frl","frl");
316  double gps;
317  double length;
318  double start;
319  double stop;
320  char atlas_path[1024];
321  frtree->Branch("gps",&gps,"gps/D");
322  frtree->Branch("length",&length,"length/D");
323  frtree->Branch("start",&start,"start/D");
324  frtree->Branch("stop",&stop,"stop/D");
325  frtree->Branch("path",atlas_path,"path/C");
326 
327  if(rfName.Sizeof()>1) frtree->SetDirectory(ofile);
328  else frtree->SetDirectory(0);
329 
330  char ifile_name[1024];
331  sprintf(ifile_name,"%s",iFile.Data());
332  cout << "CWB::frame::frl2FrTree - " << ifile_name << endl;
333  FILE *in;
334  char s[512];
335  char f[512];
336  TString idir_name;
337  cout.precision(14);
338  int cnt=0;
339  // read frame file list
340  if( (in=fopen(ifile_name,"r"))!=NULL ) {
341  while(fgets(s,512,in) != NULL) {
342  sprintf(f,"%s",s);
343  //cout << f << endl;
344  TString line(f);
345  // remove header created with gw_data_find
346  line.ReplaceAll("framefile=","");
347  line.ReplaceAll("file://localhost","");
348  line.ReplaceAll("gsiftp://ldr.aei.uni-hannover.de:15000","");
349  //cout << line.Data() << endl;
350  if (line.Contains(label)) {
351  // select the first token of the string
352  TObjArray* token0 = TString(line).Tokenize(TString(" "));
353  TObjString* sline = (TObjString*)token0->At(0);
354  TString path = sline->GetString().Data();
355 
356  //cout << f << endl;
357  sprintf(atlas_path,"%s",path.Data());
358 
359  if(mode==0) {
360  TObjArray* token1 = TString(path).Tokenize(TString("/"));
361  TObjString* path_tok1 = (TObjString*)token1->At(token1->GetEntries()-1);
362  //cout << path_tok1->GetString().Data() << endl;
363 
364  TObjArray* token2 = TString(path_tok1->GetString()).Tokenize(TString("-"));
365  TObjString* gps_tok = (TObjString*)token2->At(token2->GetEntries()-2);
366  TString sgps = gps_tok->GetString().Data();
367  //cout << sgps.Data() << endl;
368  gps = sgps.Atof();
369 
370  TObjString* length_tok = (TObjString*)token2->At(token2->GetEntries()-1);
371  TString slength = length_tok->GetString().Data();
372  slength.ReplaceAll(label,"");
373  //cout << slength.Data() << endl;
374  length = slength.Atof();
375 
376  delete token1;
377  delete token2;
378  } else {
379  TObjArray* token1 = TString(line).Tokenize(TString(" "));
380 
381  if(token1->GetEntries()<3) {
382  cout << "CWB::frame::frl2FrTree : bad format " << ifile_name << endl;
383  exit(1);
384  }
385 
386  TObjString* gps_tok = (TObjString*)token1->At(1);
387  TString sgps = gps_tok->GetString().Data();
388  //cout << sgps.Data() << endl;
389  gps = sgps.Atof();
390 
391  TObjString* length_tok = (TObjString*)token1->At(2);
392  TString slength = length_tok->GetString().Data();
393  //cout << slength.Data() << endl;
394  length = slength.Atof();
395 
396  delete token1;
397  }
398  delete token0;
399 
400  start = gps;
401 
402  stop = start+length;
403 
404  // skip file if it is out of time range
405  if(xstart && stop<xstart) continue;
406  if(xstop && start>xstop) continue;
407 
408  frtree->Fill();
409 
410  //cout << "Frame info : " << gps << " " << atlas_path << " " << length << " " << start << " " << stop << endl;
411  //if (++cnt%100==0) cout << cnt << endl;
412  ++cnt;
413  }
414  }
415  } else {
416  cout << "CWB::frame::frl2FrTree : Error opening file " << ifile_name << endl;
417  exit(1);
418  }
419 
420  if(ofile!=NULL) {
421  frtree->Write();
422  ofile->Close();
423  } else { // copy tree to internal frtree
424  if(frtree_List!=NULL) delete frtree_List;
425  frtree_List = (TTree*)frtree->CloneTree();
426  sortFrTree();
427  delete frtree;
428  }
429  return cnt;
430 }
431 
432 //______________________________________________________________________________
433 int
435 //
436 // sort the tree contained in iFile according the GPS time and save it to the auxiliary file rfName
437 // note : this method can be used only in READ mode
438 //
439 // iFile : in - name of input root file containing tree to be sorted
440 // rfName : in - name of auxiliary file name used to store the temporary tree
441 // if rfName=="" the tree is saved & sorted to the private frtree_List
442 //
443 
444  if(fOption!="READ") {
445  cout << "CWB::frame::sortFrTree : allowed only in READ mode" << endl;
446  exit(1);
447  }
448  TFile f(iFile.Data());
449  TTree *tree = (TTree*)f.Get(LST_TREE_NAME);
450  Int_t nentries = (Int_t)tree->GetEntries();
451  tree->Draw("gps","","goff");
452  Int_t *index = new Int_t[nentries];
453  TMath::Sort(nentries,tree->GetV1(),index,false);
454 
455  TFile f2(rfName.Data(),"recreate");
456  TTree *tsorted = (TTree*)tree->CloneTree(0);
457  for (Int_t i=0;i<nentries;i++) {
458  //if (i%1000==0) cout << i << endl;
459  tree->GetEntry(index[i]);
460  tsorted->Fill();
461  }
462  tsorted->Write();
463  f2.Close();
464  delete [] index;
465 
466  return nentries;
467 }
468 
469 //______________________________________________________________________________
470 int
472 //
473 // sort the the auxiliary tree according the GPS time
474 // note : this method can be used only in READ mode
475 //
476 
477  if(fOption!="READ") {
478  cout << "CWB::frame::sortFrTree : allowed only in READ mode" << endl;
479  exit(1);
480  }
481  Int_t nentries = (Int_t)frtree_List->GetEntries();
482  frtree_List->Draw("gps","","goff");
483  double* gps = frtree_List->GetV1();
484  Int_t *index = new Int_t[nentries];
485  TMath::Sort(nentries,gps,index,false);
486 
487  double pgps = -1;
488  TTree *sorted_frtree_List = (TTree*)frtree_List->CloneTree(0);
489  for (Int_t i=0;i<nentries;i++) {
490  //if (i%1000==0) cout << i << endl;
491  int j = index[i]; // j = sorted index
492  if(gps[j]==pgps) continue; // skip duplicated entries
493  pgps=gps[j];
494  frtree_List->GetEntry(j);
495  sorted_frtree_List->Fill();
496  }
497  delete [] index;
498 
499  // copy sorted tree to the original not sorted tree
500  delete frtree_List;
501  frtree_List = (TTree*)sorted_frtree_List->CloneTree();
502  delete sorted_frtree_List;
503 
504  return nentries;
505 }
506 
507 //______________________________________________________________________________
508 frfile
509 CWB::frame::getFrList(int istart, int istop, int segEdge) {
510 //
511 // return the list of frfile structures of the frame files
512 // cointained in the range [istart,istop]
513 //
514 // istart : in - start time sec
515 // istop : in - stop time sec
516 // segEdge : in - wavelet boundary offset [sec]
517 //
518 
519  frfile frf;
520  if(rfName!="") {
521  frf = getFrList(rfName, istart, istop, segEdge);
522  } else {
523  frf = getFrList(istart, istop, segEdge, NULL);
524  }
525  return frf;
526 }
527 
528 //______________________________________________________________________________
529 frfile
530 CWB::frame::getFrList(TString rfName, int istart, int istop, int segEdge) {
531 //
532 // return the frame list of frames contained in the range [start-segEdge,stop-segEdge]
533 // note : this method can be used only in READ mode
534 //
535 // rfName : in - name of auxiliary file name used to store the temporary tree
536 // use the tree contained in rfName root file
537 // istart : in - start time sec
538 // istop : in - stop time sec
539 // segEdge : in - wavelet boundary offset [sec]
540 //
541 
542  if(fOption!="READ") {
543  cout << "CWB::frame::getFrList : allowed only in READ mode" << endl;
544  exit(1);
545  }
546 
547  frfile frf;
548 
549  TFile *ifile = TFile::Open(rfName.Data());
550  if (ifile==NULL) {cout << "No " << rfName.Data() << " file !!!" << endl;exit(1);}
551  TTree* itree = (TTree *) gROOT->FindObject(LST_TREE_NAME);
552  if (itree==NULL) {cout << "No ffl tree name match !!!" << endl;exit(1);}
553  frf=getFrList(istart, istop, segEdge, itree);
554  delete itree;
555  delete ifile;
556 
557  return frf;
558 }
559 
560 //______________________________________________________________________________
561 vector<frfile>
562 CWB::frame::getFrList(int istart, int istop) {
563 //
564 // return the list of frfile structures of the frame files
565 // cointained in the range [istart,istop]
566 // note : this method can be used only in READ mode
567 //
568 // istart : in - start time sec
569 // istop : in - stop time sec
570 //
571 
572  if(fOption!="READ") {
573  cout << "CWB::frame::getFrList : allowed only in READ mode" << endl;
574  exit(1);
575  }
576 
577  vector<frfile> frlist;
578  frfile frf;
579 
580  TTree* itree=frtree_List;
581  if(itree==NULL) {cout << "CWB::frame::getFrList - Error : No tree defined !!!" << endl;exit(1);}
582 
583  int size = itree->GetEntries();
584  //cout << size << endl;
585 
586  double gps;
587  double length;
588  char path[1024];
589  itree->SetBranchAddress("gps",&gps);
590  itree->SetBranchAddress("length",&length);
591  itree->SetBranchAddress("path",path);
592 
593  double start = istart!=0 ? istart : 0;
594  double stop = istop!=0 ? istop : 2000000000;
595  char cut[1024];
596  sprintf(cut,"(gps>=%f && gps<=%f)", start,stop);
597  //cout << cut << endl;
598  itree->SetEstimate(size);
599  itree->Draw("Entry$",cut,"goff");
600  int isel = itree->GetSelectedRows();
601  double* entry = itree->GetV1();
602 
603  for (int i=0;i<isel;i++) {
604  itree->GetEntry(i);
605  frf.start=gps;
606  frf.stop=gps+length;
607  frf.length=length;
608  frf.file.push_back(path);
609  frlist.push_back(frf);
610  frf.file.clear();
611  }
612 
613  return frlist;
614 }
615 
616 //______________________________________________________________________________
617 frfile
618 CWB::frame::getFrList(int istart, int istop, int segEdge, TTree* itree) {
619 //
620 // return the frame list of frames contained in the range [start-segEdge,stop-segEdge]
621 // note : this method can be used only in READ mode
622 //
623 // istart : in - start time sec
624 // istop : in - stop time sec
625 // segEdge : in - wavelet boundary offset [sec]
626 // itree : in - if itree=NULL use the private frtree_List
627 //
628 
629  if(fOption!="READ") {
630  cout << "CWB::frame::getFrList : allowed only in READ mode" << endl;
631  exit(1);
632  }
633 
634  frfile frf;
635 
636  if(itree==NULL) itree=frtree_List;
637  if(itree==NULL) {cout << "CWB::frame::getFrList - Error : No tree defined !!!" << endl;exit(1);}
638 
639  int size = itree->GetEntries();
640  //cout << size << endl;
641 
642  double gps;
643  double length;
644  char path[1024];
645  itree->SetBranchAddress("gps",&gps);
646  itree->SetBranchAddress("length",&length);
647  itree->SetBranchAddress("path",path);
648 
649  double start=istart-segEdge;
650  double stop= istop+segEdge;
651  char cut[1024];
652  sprintf(cut,"((%f>=gps && %f<=gps+length) || (%f>=gps && %f<=gps+length) || (%f<=gps && %f>=gps+length))",
653  start,start,stop,stop,start,stop);
654  //cout << cut << endl;
655  itree->SetEstimate(size);
656  itree->Draw("Entry$",cut,"goff");
657  int isel = itree->GetSelectedRows();
658  double* entry = itree->GetV1();
659 
660 
661  // CHECK
662  if (isel==0) {
663  cout << endl;
664  cout << "CWB::frame::getFrList - Error : No files matched in the range "
665  << start << ":" << stop << " !!!" << endl << endl;
666  sprintf(cut,"gps>=%f && gps<=%f",start-500,stop+500);
667  cout << "Dump frame segments available in the input list in the range : " << cut << endl << endl;
668  itree->SetScanField(0);
669  itree->Scan("Entry$:gps:length",cut,"goff");
670  exit(1);
671  }
672  if (isel==1) {
673  itree->GetEntry(entry[0]);
674  if (start < gps) {
675  cout << endl;
676  cout << "CWB::frame::getFrList - Error : " << " No files matched !!!" << endl;
677  cout << "start buffer (" << start << ") < begin gps frame (" << gps << ")" << endl;
678  cout << "check input frame list" << endl << endl;
679  exit(1);
680  }
681  if (stop > gps+length) {
682  cout << endl;
683  cout << "CWB::frame::getFrList - Error : " << " No files matched !!!" << endl;
684  cout << "stop buffer (" << stop << ") > end gps frame (" << gps+length << ")" << endl;
685  cout << "check input frame list" << endl << endl;
686  exit(1);
687  }
688  }
689  if (isel>1) {
690  itree->GetEntry(entry[0]);
691  if (start < gps) {
692  cout << endl;
693  cout << "CWB::frame::getFrList - Error : " << " No files matched !!!" << endl;
694  cout << "start buffer (" << start << ") < begin gps frame (" << gps << ")" << endl;
695  cout << "check input frame list" << endl << endl;
696  exit(1);
697  }
698  itree->GetEntry(entry[isel-1]);
699  if (stop > gps+length) {
700  cout << endl;
701  cout << "CWB::frame::getFrList - Error : " << " No files matched !!!" << endl;
702  cout << "stop buffer (" << stop << ") > end gps frame (" << gps+length << ")" << endl;
703  cout << "check input frame list" << endl << endl;
704  exit(1);
705  }
706  // check if frames are contiguous
707  itree->GetEntry(entry[0]);
708  double stop0=gps+length;
709  for(int i=1;i<isel;i++) {
710  itree->GetEntry(entry[i]);
711  if(gps!=stop0) {
712  cout << "CWB::frame::getFrList - Error : " << " Requested frames are not contiguous !!!" << endl;
713  cout << "start : " << gps << " is != from the end of the previous segment " << stop0 << endl;
714  exit(1);
715  }
716  stop0=gps+length;
717  }
718  }
719 
720  itree->GetEntry(entry[isel-1]);
721 /*
722  cout << isel << endl;
723  cout << int(length) << endl;
724  cout << int(start) << endl;
725  cout << int(stop) << endl;
726 */
727  frf.start=start;
728  frf.stop=stop;
729  frf.length=length;
730  frf.file.clear();
731  for (int i=0;i<isel;i++) {
732  itree->GetEntry(entry[i]);
733  frf.file.push_back(path);
734  }
735 
736  return frf;
737 }
738 
739 //______________________________________________________________________________
742 //
743 // return the begin and end range of the frame list
744 // note : this method can be used only in READ mode
745 //
746 // itree : in - if itree=NULL use the private frtree_List
747 //
748 
749  if(fOption!="READ") {
750  cout << "CWB::frame::getFrRange : allowed only in READ mode" << endl;
751  exit(1);
752  }
753 
754  waveSegment SEG={0,0.,0.};
755 
756  if(itree==NULL) itree=frtree_List;
757  if(itree==NULL) {cout << "CWB::frame::getFrRange - No tree defined !!!" << endl;exit(1);}
758 
759  int size = itree->GetEntries();
760 
761  double gps;
762  double length;
763  itree->SetBranchAddress("gps",&gps);
764  itree->SetBranchAddress("length",&length);
765 
766  double min=4000000000.;
767  double max=0.;
768  for (int i=0;i<size;i++) {
769  itree->GetEntry(i);
770  if(min>gps) min=gps;
771  if(max<gps+length) max=gps+length;
772  }
773 
774  SEG.start=min;
775  SEG.stop=max;
776 
777  return SEG;
778 }
779 
780 //______________________________________________________________________________
781 void
783 //
784 // read in w the channel data contained in the frames
785 // the range of data to be read is defined in the wavearray object
786 // user must initialize the object with :
787 // - w.start(GPS)
788 // - w.stop(GPS+LEN)
789 //
790 // note : this method can be used only in READ mode
791 //
792 // w : out - data array which contains the data read
793 //
794 
795  if(getOption()!="READ") {
796  cout << "CWB::frame::readFrame : allowed only in READ mode" << endl;
797  exit(1);
798  }
799  if(getChName()=="") {
800  cout << "CWB::frame::readFrame : channel not defined" << endl;
801  exit(1);
802  }
803  frfile frf = getFrList(w.start(), w.stop(), 0);
804  readFrames(frf,const_cast<char*>(getChName().Data()),w);
805  return;
806 }
807 
808 //______________________________________________________________________________
809 void
810 CWB::frame::readFrames(char* filename, char *channel, wavearray<double> &w) {
811 //
812 // read in w the channel data contained in the frames defined in the filename
813 // this method is obsolete, was used by the wat-5.3.9 pipeline infrastructure
814 // the input file must contains the following informations :
815 //
816 // - nfiles // number of files to be read
817 // - frame length // dummy, not used anymore
818 // - start time // (sec)
819 // - stop time // (sec)
820 // - rate // dummy, not used anymore
821 // - frame file path // this field is repeated nfiles times
822 //
823 // note : this method can be used only in READ mode
824 //
825 // filename : in - name of frame file
826 // channel : in - name of channel to be extracted
827 // w : out - data array which contains the data read
828 //
829 
830  if(fOption!="READ") {
831  cout << "CWB::frame::readFrames : allowed only in READ mode" << endl;
832  exit(1);
833  }
834 
835  int buffersize=1024;
836  char *buffer=new char[buffersize];
837 
838  double rf_start, rf_stop;
839  double frlen;
840  int nframes;
841  double rf_rate;
842 
843  vector<frfile> frlist;
844  frfile frf;
845 
846  FILE *framelist=fopen(filename,"r");
847 
848  fgets(buffer,buffersize,framelist);
849  nframes=(int)atoi(buffer);
850 
851  if(nframes<=0) {
852  fprintf(stderr,"nframes=%d\n",nframes);
853  exit(1);
854  }
855 
856  fgets(buffer,buffersize,framelist);
857  frlen=(double)atof(buffer);//dummy, not used anymore
858  fgets(buffer,buffersize,framelist);
859  rf_start=(double)atof(buffer);
860  fgets(buffer,buffersize,framelist);
861  rf_stop=(double)atof(buffer);
862  fgets(buffer,buffersize,framelist);
863  rf_rate=(double)atof(buffer);//dummy, not used anymore
864 
865  if(verbose) {
866  fprintf(stdout,"-------------\n");
867  fprintf(stdout,"filename=%s channel=%s nframes=%d frlen=%f start=%f stop=%f rate=%f \n",
868  filename, channel, nframes, frlen, rf_start, rf_stop, rf_rate);
869  }
870 
871  frf.start=rf_start;
872  frf.stop=rf_stop;
873  frf.length=frlen;
874  for(int i=0;i<nframes;i++) {
875  fgets(buffer,buffersize,framelist);
876  if(verbose) fprintf(stdout,"framefile=%s\n",buffer);
877  frf.file.push_back(TString(buffer));
878  }
879  delete [] buffer;
880  fclose(framelist);
881 
882  readFrames(frf, channel, w);
883 
884  return;
885 }
886 
887 //______________________________________________________________________________
888 void
890 //
891 // read in w the channel data contained in the frames defined in the frf structure
892 // the range of data to be read is defined in the frf structure
893 //
894 // note : this method can be used only in READ mode
895 //
896 // frf : in - frame file structure
897 // channel : in - name of channel to be extracted
898 // w : out - data array which contains the data read
899 //
900 
901  if(fOption!="READ") {
902  cout << "CWB::frame::readFrames : allowed only in READ mode" << endl;
903  exit(1);
904  }
905 
906  FrVect *v;
907  FrFile *ifile;
908 
909  double rf_start, rf_stop, fstart, fend;
910  double frlen;
911  int nframes;
912  double rf_rate;
913  unsigned int wcounter=0;
914  double rf_time=0;
915  double rf_step=-1.0;
916 
917  nframes=frf.file.size();
918  rf_start=frf.start;
919  rf_stop=frf.stop;
920  frlen=frf.length;
921 
922  // BIG FRAMES : when frlen>8000 (Virgo VSR1-V2) FrFileIGetV crash (GV)
923  bool bigFrame = false;
924  if (frlen>1024) {
925  if(verbose) cout << "CWB::frame::readFrames - bigFrame !!! : " << frlen << endl;
926  bigFrame = true;
927  }
928 
929  double fend0=0;
930  for(int i=0;i<nframes;i++) {
931 
932  if(verbose) fprintf(stdout,"CWB::frame::readFrames - framefile=%s channel=%s\n",frf.file[i].Data(),channel);
933 
934  // Read Loop
935  int ntry=0;
936  ifile=NULL;
937  while(ifile==NULL&&ntry<3) {
938  ifile=FrFileINew(const_cast<char*>(frf.file[i].Data()));
939  if(ifile!=NULL || frRetryTime==0) break;
940  ntry++;
941  cout << "CWB::frame::readFrames - Failed to read file - "
942  << "retry after " << ntry*frRetryTime << " sec ..."
943  << "\t[" << ntry << "/3]" << endl;
944  if(ifile==NULL) gSystem->Sleep(1000*ntry*frRetryTime); // wait
945  }
946 
947  if(ifile==NULL) {
948  fprintf(stderr,"CWB::frame::readFrames - Cannot read file %s\n",frf.file[i].Data());
949  fflush(stderr);
950  fflush(stdout);
951  exit(1);
952  // break;
953  }
954 
955  if(ifile->error) {
956  fprintf(stderr,"CWB::frame::readFrames - Cannot read file %s, error=%d\n",frf.file[i].Data(),ifile->error);
957  fflush(stderr);
958  fflush(stdout);
959  exit(1);
960  // break;
961  }
962  fstart=(double)FrFileITStart(ifile);
963  fend=(double)FrFileITEnd(ifile);
964 
965  // check if frames are contiguous
966  if(fend0!=0 && fstart!=fend0) {
967  cout << "CWB::frame::readFrames - Error : " << " Requested frames are not contiguous !!!" << endl;
968  cout << "start : " << fstart << " is != from the end of the previous segment " << fend0 << endl;
969  exit(1);
970  }
971  fend0=fend;
972 
973  // BIG FRAMES CHECK (GV)
974  if (bigFrame) {
975  if (fstart<rf_start) fstart=rf_start;
976  if (fend>rf_stop) fend=rf_stop;
977  }
978  frlen=fend-fstart;
979  if(verbose) cout << "frlen : " << frlen << endl;
980  if (frlen<=0) continue;
981 
982  rf_time=fstart;
983 
984  v=FrFileIGetV(ifile,channel, fstart, frlen);
985 
986  // BIG FRAMES CHECK (GV)
987  if (bigFrame&&(v!=NULL)) {
988  double time_offset = fstart-v->GTime;
989  if (time_offset!=0) {
990  fprintf(stdout,"CWB::frame::readFrames - BigFrame Check Mismatch : v->GTime=%.8f != fstart=%.8f \n",v->GTime,fstart);
991  FrVectFree(v);
992  v=FrFileIGetV(ifile,channel, fstart+time_offset, frlen);
993  }
994  }
995 
996  if(v==NULL) {
997  fprintf(stderr,"CWB::frame::readFrames - FrFileIGetV failed: channel %s not present\n",channel);
998  exit(1);
999  }
1000  if((v->type!=FR_VECT_4R && v->type!=FR_VECT_8R)&&
1001  (v->type!=FR_VECT_2S && v->type!=FR_VECT_4S && v->type!=FR_VECT_8S)&&
1002  (v->type!=FR_VECT_2U && v->type!=FR_VECT_4U && v->type!=FR_VECT_8U)) {
1003  fprintf(stderr,"CWB::frame::readFrames - Wrong vector type %d\n",v->type);fflush(stderr);
1004  exit(1);
1005  }
1006  rf_step=v->dx[0];
1007  rf_rate=1./rf_step;
1008  w.resize(int((rf_stop-rf_start)*rf_rate));
1009  w.rate(rf_rate);
1010  w.start(rf_start);
1011  w.stop(rf_stop);
1012  long samples=long(frlen*rf_rate);
1013 
1014  if(verbose) {
1015  fprintf(stdout,"CWB::frame::readFrames - samples=%ld v->nDim=%d v->nData=%ld v->type=%d 8R=%d v->GTime=%.8f v->localtime=%d rate=%.1f step=%.8f\n",
1016  samples,(int)v->nDim,v->nData,v->type,FR_VECT_8R,v->GTime, v->localTime, rf_rate, rf_step);
1017  for(long j=0;j<v->nDim;j++) {
1018  fprintf(stdout,"v->nx[%ld]=%ld v->dx[%ld]=%.8f v->startX[%ld]=%.8f\n",
1019  j,v->nx[j],j,v->dx[j],j,v->startX[j]);
1020  }
1021  }
1022 
1023  // BIG FRAMES CHECK (GV)
1024  if (bigFrame&&(v->GTime!=fstart)) {fprintf(stdout,"CWB::frame::readFrames - Error : v->GTime=%.8f != fstart=%.8f \n",v->GTime,fstart);exit(1);}
1025 
1026  for(long j=0;j<samples;j++) {
1027  rf_time=fstart+j*rf_step;
1028  if(rf_time>=rf_start && rf_time<rf_stop) {
1029  if(wcounter<w.size()) {
1030  if(v->type==FR_VECT_2S) w.data[wcounter]=double(v->dataS[j]);
1031  if(v->type==FR_VECT_4S) w.data[wcounter]=double(v->dataI[j]);
1032  if(v->type==FR_VECT_8S) w.data[wcounter]=double(v->dataL[j]);
1033  if(v->type==FR_VECT_4R) w.data[wcounter]=double(v->dataF[j]);
1034  if(v->type==FR_VECT_8R) w.data[wcounter]=v->dataD[j];
1035  if(v->type==FR_VECT_2U) w.data[wcounter]=double(v->dataUS[j]);
1036  if(v->type==FR_VECT_4U) w.data[wcounter]=double(v->dataUI[j]);
1037  if(v->type==FR_VECT_8U) w.data[wcounter]=double(v->dataUL[j]);
1038  wcounter++;
1039  } else {
1040  fprintf(stderr,"w overflow\n");
1041  exit(1);
1042  }
1043  }
1044  }
1045 
1046  FrVectFree(v);
1047  FrFileIEnd(ifile);
1048  }
1049  if(verbose) fprintf(stdout,"CWB::frame::readFrames - wcounter=%d w.size()=%d\n",(int)wcounter,(int)w.size());
1050 
1051  if(w.rate()!=(1<<srIndex)) {
1052  w.Resample(1<<srIndex); // resampling
1053  if(verbose) {
1054  fprintf(stdout,"--------------------------------\n");
1055  fprintf(stdout,"After resampling: rate=%f, size=%d, start=%f\n", w.rate(),(int)w.size(),w.start());
1056  fprintf(stdout,"--------------------------------\n");
1057  }
1058  }
1059 
1060  return;
1061 }
1062 
1063 //______________________________________________________________________________
1064 int
1066 //
1067 // dump to file ofName the list of frames defined in the frf structure
1068 // note : this method can be used only in READ mode
1069 // - list file used for the old cWB analysis
1070 // - output file format
1071 //
1072 // size : number of frame files
1073 // length : stop-start
1074 // start : start time
1075 // stop : stop time
1076 // sRate : rate
1077 // file1 :
1078 // ..... :
1079 // fileN :
1080 //
1081 // frf : in - frfile structure
1082 // ofName : in - name of the output file
1083 // sRate : in - sample rate of frame data
1084 //
1085 
1086  if(fOption!="READ") {
1087  cout << "CWB::frame::dumpFrList : allowed only in READ mode" << endl;
1088  exit(1);
1089  }
1090 
1091  char ofile_name[1024];
1092  sprintf(ofile_name,"%s",ofName.Data());
1093  cout << ofile_name << endl;
1094  ofstream out;
1095  out.open(ofile_name,ios::out);
1096  if (!out.good()) {cout << "Error Opening File : " << ofName.Data() << endl;exit(1);}
1097  out.precision(14);
1098 
1099  out << frf.file.size() << endl;
1100  out << frf.length << endl;
1101  out << frf.start << endl;
1102  out << frf.stop << endl;
1103  out << int(sRate) << endl;
1104  for (int i=0;i<(int)frf.file.size();i++) {
1105  out << frf.file[i];
1106  }
1107 
1108  out.close();
1109 
1110  return 0;
1111 }
1112 
1113 //______________________________________________________________________________
1114 bool
1116 //
1117 // check if file has the correct format
1118 // xxxx-GPS-LENGTH.yyy
1119 //
1120 
1121  // remove file extension
1122  TObjArray* token0 = fName.Tokenize(TString("."));
1123  TObjString* body_tok = (TObjString*)token0->At(0);
1124  TString body = body_tok->GetString().Data();
1125 
1126  TObjArray* token = body.Tokenize(TString("-"));
1127 
1128  if(token->GetEntries()<2) {
1129  cout << "CWB::frame::fNameCheck : File Name Format Error - " << fName.Data() << endl;
1130  cout << "GPS & LENGTH must be integers : xxxx-GPS-LENGTH.yyy" << endl;
1131  return false;
1132  }
1133  TObjString* gps_tok = (TObjString*)token->At(token->GetEntries()-2);
1134  TString sgps = gps_tok->GetString().Data();
1135 
1136  TObjString* length_tok = (TObjString*)token->At(token->GetEntries()-1);
1137  TString slength = length_tok->GetString().Data();
1138 
1139  if((slength.IsDigit()) && (sgps.IsDigit())) return true;
1140 
1141  cout << "CWB::frame::fNameCheck : File Name Format Error - " << fName.Data() << endl;
1142  cout << "GPS & LENGTH must be integers : xxxx-GPS-LENGTH.yyy" << endl;
1143  return false;
1144 }
TString getOption()
Definition: frame.hh:114
TTree * tree
Definition: TimeSortTree.C:20
char cut[512]
double sRate
Definition: TestFrame5.C:14
virtual size_t size() const
Definition: wavearray.hh:127
double start
Definition: network.hh:37
TString ofName
~frame()
Definition: frame.cc:100
int nfiles
int stop
Definition: Toolbox.hh:76
tuple f
Definition: cwb_online.py:91
fprintf(stdout,"start=%f duration=%f rate=%f\n", x.start(), x.size()/x.rate(), x.rate())
Definition: ced.hh:24
double min(double x, double y)
Definition: eBBH.cc:13
virtual void rate(double r)
Definition: wavearray.hh:123
wavearray< double > & operator>>(CWB::frame &fr, wavearray< double > &x)
Definition: frame.cc:111
TString getFrName()
Definition: frame.hh:111
Int_t nentries
Definition: TimeSortTree.C:24
WSeries< float > v[nIFO]
Definition: cwb_net.C:62
cout<< "skymap size : "<< L<< endl;for(int l=0;l< L;l++) sm.set(l, l);sm > const_cast< char * >("skymap.dat")
TString("c")
int start
Definition: Toolbox.hh:75
cout<< endl;cout<< "ts size = "<< ts.size()<< " ts rate = "<< ts.rate()<< endl;tf.Forward(ts, wdm);int levels=tf.getLevel();cout<< "tf size = "<< tf.size()<< endl;double dF=tf.resolution();double dT=1./(2 *dF);cout<< "rate(hz) : "<< RATE<< "\t layers : "<< nLAYERS<< "\t dF(hz) : "<< dF<< "\t dT(ms) : "<< dT *1000.<< endl;int itime=TIME_PIXEL_INDEX;int ifreq=FREQ_PIXEL_INDEX;int index=(levels+1)*itime+ifreq;double time=itime *dT;double freq=(ifreq >0)?ifreq *dF:dF/4;cout<< endl;cout<< "PIXEL TIME = "<< time<< " sec "<< endl;cout<< "PIXEL FREQ = "<< freq<< " Hz "<< endl;cout<< endl;wavearray< double > x
void close()
Definition: frame.cc:278
CWB::frame fr(FRLIST_NAME)
frame()
Definition: frame.cc:79
void writeFrame(wavearray< double > x, TString frName, TString chName)
Definition: frame.cc:157
Long_t size
int dumpFrList(frfile frf, TString ofName, double sRate=16384.)
Definition: frame.cc:1065
#define LST_TREE_NAME
Definition: CBCTool.hh:61
char ifile_name[512]
virtual void start(double s)
Definition: wavearray.hh:119
char ofile[512]
double segEdge
Definition: test_config1.C:49
int j
Definition: cwb_net.C:10
i drho i
waveSegment getFrRange()
Definition: frame.hh:89
size_t mode
wavearray< double > w
Definition: Test1.C:27
ofstream out
Definition: cwb_merge.C:196
TString chName[NIFO_MAX]
TString uname
int frl2FrTree(TString iFile, TString rfName="", TString label=".gwf", unsigned int mode=0)
Definition: frame.cc:294
waveSegment SEG
int length
Definition: Toolbox.hh:77
i() int(T_cor *100))
TString label
Definition: MergeTrees.C:21
TF1 * f2
Definition: cbc_plots.C:1710
void Resample(const wavearray< DataType_t > &, double, int=6)
Definition: wavearray.cc:622
bool fNameCheck(TString fName)
Definition: frame.cc:1115
frfile frf
UserGroup_t * uinfo
Definition: cwb_frdisplay.C:73
TObjArray * token
double * entry
Definition: cwb_setcuts.C:206
TFile * ifile
TString getChName()
Definition: frame.hh:105
s s
Definition: cwb_net.C:137
TString frName[NIFO_MAX]
double gps
TLine * line
Definition: compare_bkg.C:482
ifstream in
wavearray< int > index
virtual void stop(double s)
Definition: wavearray.hh:121
void readFrames(char *filename, char *channel, wavearray< double > &w)
Definition: frame.cc:810
int cnt
sprintf(tfres,"(1/%g)x(%g) (sec)x(Hz)", 2 *df, df)
frfile getFrList(int istart, int istop, int segEdge)
Definition: frame.cc:509
TTree * itree
DataType_t * data
Definition: wavearray.hh:301
vector< TString > file
Definition: Toolbox.hh:78
in close()
int rnID
fclose(ftrig)
cout<< fr.getNfiles()<< endl;std::vector< frfile > frlist
int sortFrTree()
Definition: frame.cc:471
void open(TString ioFile, TString chName="", Option_t *option="", bool onDisk=false, TString label=".gwf", unsigned int mode=0)
Definition: frame.cc:212
char fName[256]
virtual void resize(unsigned int)
Definition: wavearray.cc:445
double length
Definition: TestBandPass.C:18
double stop
Definition: network.hh:38
exit(0)
TTree * tsorted
Definition: TimeSortTree.C:34