KaliVeda
Toolkit for HIC analysis
KVRawDataAnalyser.cpp
1 //Created by KVClassFactory on Thu Sep 24 11:07:45 2009
2 //Author: John Frankland,,,
3 
4 #include "KVRawDataAnalyser.h"
5 #include "KVMultiDetArray.h"
6 #include "KVClassFactory.h"
7 #include "KVDataSet.h"
8 #include "TSystem.h"
9 #include "TROOT.h"
10 
11 using namespace std;
12 
13 ClassImp(KVRawDataAnalyser)
14 
15 
16 
17 
18 
22 {
23  // Default constructor
24  fRunFile = 0;
25  TotalEntriesToRead = 0;
26  fTreeFile = nullptr;
27 }
28 
29 
30 
33 
35 {
36  // Destructor
37 }
38 
39 
40 
53 
55 {
56  // Perform treatment of a given run
57  // Before processing each run, after opening the associated file, user's InitRun() method is called.
58  // After each run, user's EndRun() is called.
59  // For each event of each run, user's Analysis() method is called just after calling
60  // gMultiDetArray->HandleRawDataEvent(fRunfile). In the Analysis() method the user can
61  // test gMultiDetArray->HandledRawData() to know if some pertinent data was found in
62  // the event or not.
63  //
64  // For further customisation, the pre/post-methods are called just before and just after
65  // each of these methods (preInitRun(), postAnalysis(), etc. etc.)
66 
67  //Open data file
68  KVString raw_file = gDataSet->GetFullPathToRunfile(GetDataType().Data(), fRunNumber);
69  fRunFile = gDataSet->OpenRunfile<KVRawDataReader>(GetDataType().Data(), fRunNumber);
70  if ((!fRunFile) || fRunFile->IsZombie()) {
71  //skip run if file cannot be opened
72  if (fRunFile) delete fRunFile;
73  return;
74  }
75 
76  //warning! real number of run may be different from that deduced from file name
77  //we get the real run number from the data and use it to name any new files
78  // Not possible for INDRAFAZIA MFM data (we use 'fake' run numbers)
79  // Is this still necessary? (which dataset was concerned? camp5?)
80 // Int_t newrun = fRunFile->GetRunNumberReadFromFile();
81 // if (newrun && newrun != fRunNumber) {
82 // cout << " *** WARNING *** run number read from file = " << newrun << endl;
83 // fRunNumber = newrun;
84 // }
85 
86  KVMultiDetArray::MakeMultiDetector(gDataSet->GetName(), fRunNumber);
87  fCurrentRun = gExpDB->GetDBRun(fRunNumber);
88 
89  // perform any initialisations which may be required, dependent on the
90  // format of the data, the multidetector, etc.
91  gMultiDetArray->InitialiseRawDataReading(fRunFile);
92 
93  fEventNumber = 1; //event number
94 
95  Long64_t nevents = GetNbEventToRead();
96  if (nevents <= 0) {
97  nevents = 1000000000;
98  cout << endl << "Reading all events from file " << raw_file.Data() << endl;
99  }
100  else {
101  cout << endl << "Reading " << nevents << " events from file " << raw_file.Data() << endl;
102  }
103 
104  cout << "Starting analysis of run " << fRunNumber << " on : ";
105  TDatime now;
106  cout << now.AsString() << endl << endl;
107 
108  fCurrentRun->Print();
109 
110  preInitRun();
111  //call user's beginning of run
112  InitRun();
113  postInitRun();
114 
115  //loop over events in file
116  while ((nevents-- ? fRunFile->GetNextEvent() : kFALSE) && !AbortProcessingLoop()) {
117 
118  gMultiDetArray->HandleRawDataEvent(fRunFile);
119 
120  // skip events that any active patches reject
121  if(fRustines.HasActivePatches())
122  if(fRustines.SkipEvent(gMultiDetArray))
123  continue;
124 
125  preAnalysis();
126  //call user's analysis. stop if returns kFALSE.
127  if (!Analysis()) break;
128  postAnalysis();
129 
130  if (CheckStatusUpdateInterval(fEventNumber)) DoStatusUpdate(fEventNumber);
131 
132  fEventNumber += 1;
133  }
134 
135  cout << "Ending analysis of run " << fRunNumber << " on : ";
136  TDatime now2;
137  cout << now2.AsString() << endl << endl;
138  cout << endl << "Finished reading " << fEventNumber - 1 << " events from file " << raw_file.Data() << endl << endl;
139 
140  preEndRun();
141  //call user's end of run function
142  EndRun();
143  postEndRun();
144 
145  delete fRunFile;
146 }
147 
148 
149 
150 
162 
164 {
165  // Perform analysis of chosen runs
166  // Before beginning the loop over the runs, the user's InitAnalysis() method is called.
167  // After completing the analysis of all runs, the user's EndAnalysis() method is called.
168  //
169  // Further customisation of the event loop is possible by overriding the methods
170  // preInitAnalysis()
171  // postInitAnalysis()
172  // preEndAnalysis()
173  // postEndAnalysis()
174  // which are executed respectively just before and just after the methods.
175 
176  if (gDataSet != GetDataSet()) GetDataSet()->cd();
177 
178  // prepare any user-supplied options for the analysis
179  fOptionList.ParseOptions(GetUserClassOptions());
180 
181  preInitAnalysis();
182  //call user's initialisation
183  InitAnalysis();
184  postInitAnalysis();
185 
186  CalculateTotalEventsToRead();
187 
188  //loop over runs
189  GetRunList().Begin();
190  while (!GetRunList().End() && !AbortProcessingLoop()) {
191  fRunNumber = GetRunList().Next();
192  fRustines.InitializePatchList(gDataSet->GetName(),fRunNumber);
193  fRustines.Print();
194  ProcessRun();
195  }
196 
197  if (fCombinedOutputFile != "") {
198  Info("Terminate", "combine = %s", fCombinedOutputFile.Data());
199  // combine histograms and trees from analysis into one file
200  TString file1, file2;
201  file1.Form("HistoFileFrom%s.root", ClassName());
202  file2.Form("TreeFileFrom%s.root", ClassName());
203  if (fHistoList.GetEntries()) {
204  if (fTreeFile) {
205  Info("Terminate", "both");
206  SaveHistos();
207  fTreeFile->Write();
208  delete fTreeFile;
209  KVBase::CombineFiles(file1, file2, fCombinedOutputFile, kFALSE);
210  }
211  else {
212  // no trees - just rename histo file
213  Info("Terminate", "histo");
214  SaveHistos(fCombinedOutputFile);
215  }
216  }
217  else if (fTreeFile) {
218  // no histos - just rename tree file
219  Info("Terminate", "tree");
220  fTreeFile->Write();
221  delete fTreeFile;
222  gSystem->Rename(file2, fCombinedOutputFile);
223  }
224  else Info("Terminate", "none");
225  }
226  else {
227  SaveHistos();
228  if (fTreeFile) {
229  fTreeFile->Write();
230  delete fTreeFile;
231  }
232  }
233 
234  preEndAnalysis();
235  //call user's end of analysis
236  EndAnalysis();
237  postEndAnalysis();
238 }
239 
240 
241 
242 
245 
247 {
248  //loop over runs and calculate total events
249  TotalEntriesToRead = 0;
250  GetRunList().Begin();
251  while (!GetRunList().End()) {
252  Int_t r = GetRunList().Next();
253  TotalEntriesToRead += gExpDB->GetDBRun(r)->GetEvents();
254  }
255 }
256 
257 
258 
263 
264 void KVRawDataAnalyser::AddTree(TTree* tree)
265 {
266  // Add a user analysis results TTree to list which will be automatically
267  // written to disk at end of analysis.
268  // User must call CreateTreeFile(const Char_t*) before calling this method
269 
270  if (!fTreeFile) {
271  Error("AddTree", "You must call CreateTreeFile(const Char_t*) before using this method!");
272  return;
273  }
274  tree->SetDirectory(fTreeFile);
275  tree->AutoSave();
276  fTreeList.Add(tree);
277 }
278 
279 
280 
286 
287 void KVRawDataAnalyser::FillHisto(const Char_t* histo_name, Double_t one, Double_t two, Double_t three, Double_t four)
288 {
289  //Find in the list, if there is an histogram named "histo_name"
290  //If not print an error message
291  //If yes redirect to the right method according to its closest mother class
292  //to fill it
293 
294  TH1* h1 = 0;
295  if ((h1 = GetHisto(histo_name))) {
296  if (h1->InheritsFrom("TH3"))
297  FillTH3((TH3*)h1, one, two, three, four);
298  else if (h1->InheritsFrom("TProfile2D"))
299  FillTProfile2D((TProfile2D*)h1, one, two, three, four);
300  else if (h1->InheritsFrom("TH2"))
301  FillTH2((TH2*)h1, one, two, three);
302  else if (h1->InheritsFrom("TProfile"))
303  FillTProfile((TProfile*)h1, one, two, three);
304  else if (h1->InheritsFrom("TH1"))
305  FillTH1(h1, one, two);
306  else
307  Warning("FillHisto", "%s -> Classe non prevue ...", fHistoList.FindObject(histo_name)->ClassName());
308  }
309  else {
310  Warning("FillHisto", "%s introuvable", histo_name);
311  }
312 }
313 
314 
315 
318 
319 void KVRawDataAnalyser::FillHisto(const Char_t* histo_name, const Char_t* label, Double_t weight)
320 {
321  // Fill 1D histogram with named bins
322  TH1* h1 = 0;
323  if ((h1 = GetHisto(histo_name))) {
324  h1->Fill(label, weight);
325  }
326  else {
327  Warning("FillHisto", "%s introuvable", histo_name);
328  }
329 }
330 
331 
332 
339 
340 void KVRawDataAnalyser::FillTree(const Char_t* tree_name)
341 {
342  //Filltree method, the tree named tree_name
343  //has to be declared with AddTTree(TTree*) method
344  //
345  //if no sname="", all trees in the list is filled
346  //
347  if (!strcmp(tree_name, "")) {
348  fTreeList.Execute("Fill", "");
349  }
350  else {
351  TTree* tt = 0;
352  if ((tt = GetTree(tree_name))) {
353  tt->Fill();
354  }
355  else {
356  Warning("FillTree", "%s introuvable", tree_name);
357  }
358  }
359 }
360 
361 
362 
376 
377 void KVRawDataAnalyser::SaveHistos(const Char_t* filename, Option_t* option, Bool_t onlyfilled)
378 {
379  // Write in file all histograms declared with AddHisto(TH1*)
380  //
381  // If no filename is specified, set default name : HistoFileFrom[name_of_class].root
382  //
383  // If a filename is specified, search in gROOT->GetListOfFiles() if
384  // this file has been already opened
385  // - if yes write in it
386  // - if not, create it with the corresponding option, write in it
387  // and close it just after
388  //
389  // onlyfilled flag allow to write all (onlyfilled=kFALSE, default)
390  // or only histograms (onlyfilled=kTRUE) those have been filled
391 
392  if (!fHistoList.GetEntries()) return;
393 
394  TString histo_file_name = "";
395  if (!strcmp(filename, ""))
396  histo_file_name.Form("HistoFileFrom%s.root", ClassName());
397  else
398  histo_file_name = filename;
399 
400  Bool_t justopened = kFALSE;
401 
402  TFile* file = 0;
403  TDirectory* pwd = gDirectory;
404  //if filename correspond to an already opened file, write in it
405  //if not open/create it, depending on the option ("recreate" by default)
406  //and write in it
407  if (!(file = (TFile*)gROOT->GetListOfFiles()->FindObject(histo_file_name.Data()))) {
408  file = new TFile(histo_file_name.Data(), option);
409  justopened = kTRUE;
410  }
411  file->cd();
412  TIter next(GetHistoList());
413  TObject* obj = 0;
414  while ((obj = next())) {
415  if (obj->InheritsFrom("TH1")) {
416  if (onlyfilled) {
417  if (((TH1*)obj)->GetEntries() > 0) {
418  obj->Write();
419  }
420  }
421  else {
422  obj->Write();
423  }
424  }
425  }
426  if (justopened)
427  file->Close();
428  pwd->cd();
429 }
430 
431 
432 
436 
437 Bool_t KVRawDataAnalyser::CreateTreeFile(const Char_t* filename)
438 {
439  // This method must be called before creating any user TTree in InitAnalysis().
440  // If no filename is given, default name="TreeFileFrom[name of analysis class].root"
441 
442  TString tree_file_name;
443  if (!strcmp(filename, ""))
444  tree_file_name.Form("TreeFileFrom%s.root", ClassName());
445  else
446  tree_file_name = filename;
447 
448  TDirectory* savedir = gDirectory;
449  fTreeFile = new TFile(tree_file_name, "RECREATE");
450  savedir->cd();
451 
452  return kTRUE;
453 }
454 
455 
456 
457 
460 
461 void KVRawDataAnalyser::Make(const Char_t* kvsname)
462 {
463  //Automatic generation of derived class for raw data analysis
464 
465  KVClassFactory cf(kvsname, "Analysis of raw data", "",
466  kTRUE, "RawAnalysisTemplate");
467 
468  cf.AddImplIncludeFile("KVMultiDetArray.h");
469  cf.GenerateCode();
470 }
471 
472 
473 
476 
478 {
479  // Method called to abort analysis during processing of a run
480 }
481 
482 
static void CombineFiles(const Char_t *file1, const Char_t *file2, const Char_t *newfilename, Bool_t keep=kTRUE)
Definition: KVBase.cpp:1527
Factory class for generating skeleton files for new classes.
void GenerateCode()
Generate header and implementation file for currently-defined class.
void AddImplIncludeFile(const Char_t *filename)
ULong64_t GetEvents() const
Definition: KVDBRun.h:134
FileType * OpenRunfile(const Char_t *type, Int_t run)
Definition: KVDataSet.h:167
TString GetFullPathToRunfile(const Char_t *type, Int_t run) const
Definition: KVDataSet.cpp:897
KVDBRun * GetDBRun(Int_t number) const
Definition: KVExpDB.h:76
virtual void InitialiseRawDataReading(KVRawDataReader *)
static KVMultiDetArray * MakeMultiDetector(const Char_t *dataset_name, Int_t run=-1, TString classname="KVMultiDetArray")
virtual Bool_t HandleRawDataEvent(KVRawDataReader *)
Abstract base class for user analysis of raw data.
virtual void SaveHistos(const Char_t *filename="", Option_t *option="recreate", Bool_t onlyfilled=kFALSE)
Bool_t CreateTreeFile(const Char_t *filename="")
void FillTree(const Char_t *sname="")
void FillHisto(const Char_t *sname, Double_t one, Double_t two=1, Double_t three=1, Double_t four=1)
virtual ~KVRawDataAnalyser()
Destructor.
static void Make(const Char_t *kvsname="MyOwnRawDataAnalyser")
Automatic generation of derived class for raw data analysis.
void CalculateTotalEventsToRead()
loop over runs and calculate total events
void AddTree(TTree *tree)
void AbortDuringRunProcessing()
Method called to abort analysis during processing of a run.
Abstract base class for reading raw (DAQ) data.
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition: KVString.h:73