KaliVeda
Toolkit for HIC analysis
KVINDRAReconDataAnalyser.cpp
1 /*
2 $Id: KVINDRAReconDataAnalyser.cpp,v 1.4 2007/11/15 14:59:45 franklan Exp $
3 $Revision: 1.4 $
4 $Date: 2007/11/15 14:59:45 $
5 */
6 
7 //Created by KVClassFactory on Wed Apr 5 23:50:04 2006
8 //Author: John Frankland
9 
10 #include "KVBase.h"
11 #include "KVINDRAReconDataAnalyser.h"
12 #include "KVINDRADBRun.h"
13 #include "KVINDRADB.h"
14 #include "KVDataAnalysisTask.h"
15 #include "KVDataSet.h"
16 #include "TChain.h"
17 #include "TObjString.h"
18 #include "TChain.h"
19 #include "KVAvailableRunsFile.h"
20 #include "KVINDRA.h"
21 #include "TProof.h"
22 #include "KVINDRATriggerConditions.h"
23 #include "KVINDRAEventSelector.h"
24 
25 #include <KVReconEventSelector.h>
26 
27 using namespace std;
28 
30 
31 
32 
35  : fSelector(nullptr), fOldSelector(nullptr), theChain(nullptr),
36  theRawData(nullptr), theGeneData(nullptr),
37  ParVal(nullptr), ParNum(nullptr), parList(nullptr)
38 {
39 }
40 
41 
42 
45 
47 {
48  //Reset task variables
50  theChain = nullptr;
51  theRawData = nullptr;
52  theGeneData = nullptr;
53  ParVal = nullptr;
54  ParNum = nullptr;
55  fSelector = nullptr;
56  fOldSelector = nullptr;
57  TotalEntriesToRead = 0;
58 }
59 
60 
61 
62 
63 
64 
67 
69 {
70  // Checks the task variables
71 
73 
74  cout << "============> Analysis summary <=============" << endl;
75  cout << "Analysis of runs " << GetRunList().
76  GetList() << " with the class ";
77  cout << "\"" << GetUserClassName() << "\"." << endl;
78  if (GetNbEventToRead()) {
79  cout << GetNbEventToRead() << " events will be processed." << endl;
80  }
81  else {
82  cout << "All events will be processed." << endl;
83  }
84  cout << "=============================================" << endl;
85 
86  return kTRUE;
87 }
88 
89 
90 
91 
96 
98 {
99  //Run the interactive analysis
100 
101  //make the chosen dataset the active dataset ( = gDataSet; note this also opens database
102  //and positions gDataBase & gIndraDB).
103  GetDataSet()->cd();
104  fSelector = nullptr;
105  fOldSelector = nullptr;
106 
107  theChain = new TChain("ReconstructedEvents");
108  theChain->SetDirectory(0); // we handle delete
109 
110  // open and add to TChain all required files
111  // we force the opening of the files to avoid problems with xrootd which sometimes
112  // seems to have a little difficulty
113  for (auto& run : GetRunList()) {
114  TString fullPathToRunfile = gDataSet->GetFullPathToRunfile(GetDataType(), run);
115  cout << "Opening file " << fullPathToRunfile << endl;
116  TFile* f = gDataSet->OpenRunfile<TFile>(GetDataType(), run);
117  cout << "Adding file " << fullPathToRunfile;
118  cout << " to the TChain." << endl;
119  dynamic_cast<TChain*>(theChain)->Add(fullPathToRunfile);
120  if (f && !f->IsZombie()) {
121  // update run infos in available runs file if necessary
123  if (ARF->InfosNeedUpdate(run, gSystem->BaseName(fullPathToRunfile))) {
124  if (!((TTree*)f->Get("ReconstructedEvents"))) {
125  Error("SubmitTask", "No tree named ReconstructedEvents is present in the current file");
126  delete theChain;
127  return;
128  }
129  TEnv* treeInfos = (TEnv*)((TTree*)f->Get("ReconstructedEvents"))->GetUserInfo()->FindObject("TEnv");
130  if (treeInfos) {
131  TString kvversion = treeInfos->GetValue("KVBase::GetKVVersion()", "");
132  TString username = treeInfos->GetValue("gSystem->GetUserInfo()->fUser", "");
133  if (kvversion != "") ARF->UpdateInfos(run, gSystem->BaseName(fullPathToRunfile), kvversion, username);
134  }
135  else {
136  Info("SubmitTask", "No TEnv object associated to the tree");
137  }
138  }
139  }
140  }
141  TotalEntriesToRead = theChain->GetEntries();
142  TString option("");
143 
144  // Add the total run list in option
145  if (!(option.IsNull())) option += ",";
146  option += Form("FullRunList=%s", GetFullRunList().GetList().Data());
147 
148  // Add any user-defined options
149  if (GetUserClassOptions() != "") {
150  option += ",";
151  option += GetUserClassOptions();
152  }
153 
154  set_up_analyser_for_task(this);
155 
156  // for backwards compatibility, we allow user class to inherit from
157  // KVOldINDRASelector instead of KVINDRAEventSelector
158  // for forwards compatibility, we allow it to inherit from KVReconEventSelector too!
159  TObject* new_selector = GetInstanceOfUserClass("KVOldINDRASelector,KVReconEventSelector");
160 
161  if (!new_selector || !new_selector->InheritsFrom("TSelector")) {
162  cout << "The selector \"" << GetUserClassName() << "\" is not valid." << endl;
163  cout << "Process aborted." << endl;
164  SafeDelete(new_selector);
165  }
166  else {
167  SafeDelete(new_selector);
168  Info("SubmitTask", "Beginning TChain::Process...");
169 #ifdef WITH_CPP11
170  if (GetProofMode() != KVDataAnalyser::EProofMode::None) dynamic_cast<TChain*>(theChain)->SetProof(kTRUE);
171 #else
172  if (GetProofMode() != KVDataAnalyser::None) dynamic_cast<TChain*>(theChain)->SetProof(kTRUE);
173 #endif
174  TString analysis_class;
175  if (GetAnalysisTask()->WithUserClass()) analysis_class.Form("%s%s", GetUserClass().full_path_imp().Data(), GetACliCMode());
176  else analysis_class = GetUserClassName();
177 
178  if (GetNbEventToRead()) {
179  theChain->Process(analysis_class, option.Data(), GetNbEventToRead());
180  }
181  else {
182  theChain->Process(analysis_class, option.Data());
183  }
184  }
185  delete theChain;
186  fSelector = nullptr; //deleted by TChain/TTreePlayer
187  fOldSelector = nullptr; //deleted by TChain/TTreePlayer
188 }
189 
190 
191 
192 
198 
200 {
201  //Save (in the TEnv fBatchEnv) all necessary information on analysis task which can be used to execute it later
202  //(i.e. when batch processing system executes the job).
203  //If save=kTRUE (default), write the information in a file whose name is given by ".jobname"
204  //where 'jobname' is the name of the job as given to the batch system.
205 
207  // backwards-compatible fix for old KVSelector analysis classes
208  // for forwards compatibility, we allow it to inherit from KVReconEventSelector too!
209  GetBatchInfoFile()->SetValue("UserClassAlternativeBaseClass", "KVOldINDRASelector,KVReconEventSelector");
210  if (save) GetBatchInfoFile()->SaveLevel(kEnvUser);
211 }
212 
213 
214 
215 
225 
227 {
228  // Called by currently-processed KVEventSelector before user's InitAnalysis() method.
229  // We build the multidetector for the current dataset in case informations on
230  // detector are needed e.g. to define histograms in InitAnalysis().
231  // Note that at this stage we are not analysing a given run, so the parameters
232  // of the array are not set (they will be set in preInitRun()).
233  //
234  // Note for PROOF: as this will be called both on master and on slave workers,
235  // in order to reduce memory footprint we only build INDRA on the slaves
236 
237  if (!gProof || !gProof->IsMaster()) {
238  if (!gIndra) KVMultiDetArray::MakeMultiDetector(GetDataSet()->GetName());
239  }
240  // in case data is being analysed with a "vanilla" KVReconEventSelector, we need
241  // to set the correct branch name for reconstructed INDRA data which is "INDRAReconEvent"
242  // instead of the default "ReconEvent"
243  if (fSelector) fSelector->SetBranchName("INDRAReconEvent");
244  else Warning("preInitAnalysis", "could not set branch name correctly");
245 }
246 
247 
248 
250 
252 {
253  if (fSelector) {
254  if (fSelector->InheritsFrom("KVINDRAEventSelector"))
255  dynamic_cast<KVINDRAEventSelector*>(fSelector)->SetCurrentRun(CurrentRun);
256  // WARNING - horrible kludge...
257  else if (fSelector->InheritsFrom("KVReconEventSelector"))
258  dynamic_cast<KVReconEventSelector*>(fSelector)->SetCurrentRun(CurrentRun);
259  }
260  else {
261  fOldSelector->SetCurrentRun(CurrentRun);
262  }
263 }
264 
265 
266 
282 
284 {
285  // Called by currently-processed TSelector when a new file in the TChain is opened.
286  // We call gIndra->SetParameters for the current run: whether only physics parameters will be set,
287  // or the full set of identification/calibration parameters depends on the value of the
288  // environment variable
289  //~~~
290  // [dataset].ReconAnalysis.WithCalibInfos: [yes/no]
291  //~~~
292  //This can be overridden for any individual analysis by setting the analysis class option
293  //~~~~
294  // WithCalibInfos=yes
295  //~~~~
296  // We connect the acquisition parameter objects to the branches of the raw data tree.
297  // Infos on currently read file/tree are printed.
298  // Any required data patches ("rustines") are initialized.
299 
300  Bool_t physics_parameters_only = !gDataSet->GetDataSetEnv("ReconAnalysis.WithCalibInfos", kTRUE);
301  if (fSelector->IsOptGiven("WithCalibInfos"))
302  physics_parameters_only = (fSelector->GetOpt("WithCalibInfos") != "yes");
303 
304  Int_t run = GetRunNumberFromFileName(theChain->GetCurrentFile()->GetName());
305  gIndra->SetParameters(run, physics_parameters_only);
306  KVINDRADBRun* CurrentRun = gIndraDB->GetRun(run);
307  SetCurrentRun(CurrentRun);
308  SetSelectorCurrentRun(CurrentRun);
309  cout << endl << " =================== New Run =================== " <<
310  endl << endl;
311 
312  CurrentRun->Print();
313  if (CurrentRun->GetSystem()) {
314  SetSystem(CurrentRun->GetSystem());
315  if (CurrentRun->GetSystem()->GetKinematics())
316  CurrentRun->GetSystem()->GetKinematics()->Print();
317  }
318 
319  cout << endl << " ================================================= " <<
320  endl << endl;
321 
322  ConnectRawDataTree();
323  ConnectGeneDataTree();
324  PrintTreeInfos();
325  Info("preInitRun", "Data written with series %s, release %d", GetDataSeries().Data(),
326  GetDataReleaseNumber());
327  fRustines.InitializePatchList(GetDataSet()->GetName(), GetDataType(), run, GetDataSeries(),
328  GetDataReleaseNumber(), theChain->GetCurrentFile()->GetStreamerInfoCache());
329  fRustines.Print();
330 }
331 
332 
333 
335 
337 {
338  Long64_t rawEntry = (fSelector ? fSelector->GetEventNumber() - 1
339  : fOldSelector->GetEventNumber() - 1);
340 
341  return rawEntry;
342 }
343 
344 
345 
347 
349 {
350  return (fSelector ? dynamic_cast<KVReconstructedEvent*>(fSelector->GetEvent()) : fOldSelector->GetEvent());
351 }
352 
353 
354 #ifdef USING_ROOT6
355 
357 
359 {
361  trig.SetTriggerConditionsForRun(fSelector, run);
362 }
363 
364 #endif
365 
366 
370 
372 {
373  // Read and set raw data for the current reconstructed event
374  // Any required data patches ("rustines") are applied.
375 
376  if (!theRawData) return;
377  // all recon events are numbered 1, 2, ... : therefore entry number is N-1
378  Long64_t rawEntry = GetRawEntryNumber();
379  theRawData->GetEntry(rawEntry);
381  for (int i = 0; i < NbParFired; i++)
382  gIndra->handle_ebyedat_raw_data_parameter((*parList)[ParNum[i]]->GetName(), ParVal[i]);
383 
384  // as rustines often depend on a knowledge of the original raw data,
385  // we apply them after it has been read in
386  KVINDRAReconEvent* event = (KVINDRAReconEvent*)GetReconstructedEvent();
387  if (fRustines.HasActivePatches()) fRustines.Apply(event);
388 }
389 
390 
391 
398 
400 {
401  // Called by preInitRun().
402  // When starting to read a new run (=new file), we look for the TTree "RawData" in the
403  // current file (it should have been created by KVINDRARawDataReconstructor).
404  // If found, it will be used by ReadRawData() to set the values of all acquisition parameters
405  // for each event.
406 
407  theRawData = (TTree*)theChain->GetCurrentFile()->Get("RawData");
408  if (!theRawData) {
409  Warning("ConnectRawDataTree", "RawData tree not found in file; raw data parameters of detectors will not be available in analysis");
410  return;
411  }
412  else
413  Info("ConnectRawDataTree", "Found RawData tree in file");
414  Int_t maxNopar = theRawData->GetMaximum("NbParFired");
415  if (ParVal) delete [] ParVal;
416  if (ParNum) delete [] ParNum;
417  ParVal = new UShort_t[maxNopar];
418  ParNum = new UInt_t[maxNopar];
419  parList = (TObjArray*)theRawData->GetUserInfo()->FindObject("ParameterList");
420  theRawData->SetBranchAddress("NbParFired", &NbParFired);
421  theRawData->SetBranchAddress("ParNum", ParNum);
422  theRawData->SetBranchAddress("ParVal", ParVal);
423  Info("ConnectRawDataTree", "Connected raw data parameters");
424  Entry = 0;
425 }
426 
427 
428 
433 
435 {
436  // Called by preInitRun().
437  // When starting to read a new run (=new file), we look for the TTree "GeneData" in the
438  // current file (it should have been created by KVINDRARawDataReconstructor).
439 
440  theGeneData = (TTree*)theChain->GetCurrentFile()->Get("GeneData");
441  if (!theGeneData) {
442  cout << " --> No pulser & laser data for this run !!!" << endl << endl;
443  }
444  else {
445  cout << " --> Pulser & laser data tree contains " << theGeneData->GetEntries()
446  << " events" << endl << endl;
447  }
448 }
449 
450 
451 
453 
455 {
456  return (TEnv*)theChain->GetTree()->GetUserInfo()->FindObject("TEnv");
457 }
458 
459 
460 
463 
465 {
466  // Print informations on currently analysed TTree
467  TEnv* treeInfos = GetReconDataTreeInfos();
468  if (!treeInfos) return;
469  cout << endl << "----------------------------------------------------------------------------------------------------" << endl;
470  cout << "INFORMATIONS ON VERSION OF KALIVEDA USED TO GENERATE FILE:" << endl << endl;
471  fDataVersion = treeInfos->GetValue("KVBase::GetKVVersion()", "(unknown)");
472  cout << "version = " << fDataVersion << endl ;
473  cout << "build date = " << treeInfos->GetValue("KVBase::GetKVBuildDate()", "(unknown)") << endl ;
474  cout << "source directory = " << treeInfos->GetValue("KVBase::GetKVSourceDir()", "(unknown)") << endl ;
475  cout << "KVROOT = " << treeInfos->GetValue("KVBase::GetKVRoot()", "(unknown)") << endl ;
476  cout << "BZR branch name = " << treeInfos->GetValue("KVBase::bzrBranchNick()", "(unknown)") << endl ;
477  cout << "BZR revision #" << treeInfos->GetValue("KVBase::bzrRevisionNumber()", "(unknown)") << endl ;
478  cout << "BZR revision ID = " << treeInfos->GetValue("KVBase::bzrRevisionId()", "(unknown)") << endl ;
479  cout << "BZR revision date = " << treeInfos->GetValue("KVBase::bzrRevisionDate()", "(unknown)") << endl ;
480  cout << endl << "INFORMATIONS ON GENERATION OF FILE:" << endl << endl;
481  cout << "Generated by : " << treeInfos->GetValue("gSystem->GetUserInfo()->fUser", "(unknown)") << endl ;
482  cout << "Analysis task : " << treeInfos->GetValue("AnalysisTask", "(unknown)") << endl ;
483  cout << "Job name : " << treeInfos->GetValue("BatchSystem.JobName", "(unknown)") << endl ;
484  cout << "Job submitted from : " << treeInfos->GetValue("LaunchDirectory", "(unknown)") << endl ;
485  cout << "Runs : " << treeInfos->GetValue("Runs", "(unknown)") << endl ;
486  cout << "Number of events requested : " << treeInfos->GetValue("NbToRead", "(unknown)") << endl ;
487  cout << endl << "----------------------------------------------------------------------------------------------------" << endl;
488 
489  // if possible, parse fDataVersion into series and release number
490  // e.g. if fDataVersion = "1.8.10":
491  // => fDataSeries = "1.8" fDataReleaseNum = 10
492  Int_t a, b, c;
493  if (fDataVersion != "(unknown)") {
494  if (sscanf(fDataVersion.Data(), "%d.%d.%d", &a, &b, &c) == 3) {
495  fDataSeries.Form("%d.%d", a, b);
496  fDataReleaseNum = c;
497  }
498  }
499  else {
500  fDataSeries = "";
501  fDataReleaseNum = -1;
502  }
503 }
504 
505 
506 
508 
510 {
511  if (theRawData) theRawData->CloneTree(-1, "fast"); //copy raw data tree to file
512  if (theGeneData) theGeneData->CloneTree(-1, "fast"); //copy pulser & laser (gene) tree to file
513 }
514 
515 
516 
523 
525 {
526  // Overrides KVDataAnalyser method, for backwards & forwards compatibility.
527  //
528  // Old user analysis classes may derive from KVOldINDRASelector (aka KVSelector) instead of KVINDRAEventSelector.
529  //
530  // It is also possible to analyse data with a "vanilla" KVReconEventSelector analysis class.
531  return KVDataAnalyser::CheckIfUserClassIsValid("KVOldINDRASelector,KVReconEventSelector");
532 }
533 
534 
int Int_t
unsigned int UInt_t
#define SafeDelete(p)
#define f(i)
#define c(i)
bool Bool_t
unsigned short UShort_t
constexpr Bool_t kFALSE
constexpr Bool_t kTRUE
kEnvUser
Option_t Option_t option
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t b
R__EXTERN TProof * gProof
char * Form(const char *fmt,...)
R__EXTERN TSystem * gSystem
void Print(Option_t *opt="") const override
Definition: KV2Body.cpp:813
Handles lists of available runs for different datasets and types of data.
Bool_t InfosNeedUpdate(const run_index_t &run, const KVString &filename)
void UpdateInfos(const run_index_t &run, const KVString &filename, const KVString &kvversion, const KVString &username)
void Print(Option_t *option="") const override
Definition: KVDBRun.cpp:61
KVDBSystem * GetSystem() const
Definition: KVDBRun.cpp:237
KV2Body * GetKinematics()
Definition: KVDBSystem.cpp:80
virtual void Reset()
virtual Bool_t CheckIfUserClassIsValid(const KVString &alternative_base_class="")
virtual Bool_t CheckTaskVariables()
void WriteBatchEnvFile(const TString &jobname, Bool_t save=kTRUE) override
TString GetFullPathToRunfile(const KVString &type, const run_index_t &run) const
Definition: KVDataSet.cpp:931
KVAvailableRunsFile * GetAvailableRunsFile(const Char_t *type) const
Definition: KVDataSet.cpp:50
KVString GetDataSetEnv(const Char_t *type, const Char_t *defval="") const
Definition: KVDataSet.cpp:784
FileType * OpenRunfile(const KVString &type, const run_index_t &run)
Definition: KVDataSet.h:164
Database entry for each run of an INDRA experiment.
Definition: KVINDRADBRun.h:30
KVINDRADBRun * GetRun(Int_t run) const
Definition: KVINDRADB.h:133
Base class for analysis of reconstructed INDRA events.
Manage analysis of reconstructed INDRA data.
void SetTriggerConditionsForRun(int) override
void SetSelectorCurrentRun(KVINDRADBRun *CurrentRun)
void WriteBatchEnvFile(const TString &, Bool_t sav=kTRUE) override
void PrintTreeInfos()
Print informations on currently analysed TTree.
KVReconstructedEvent * GetReconstructedEvent()
Bool_t CheckIfUserClassIsValid(const KVString &="") override
Bool_t CheckTaskVariables(void) override
Checks the task variables.
void Reset() override
Reset task variables.
Event reconstructed from energy losses in INDRA multidetector.
Set trigger conditions for analysis of reconstructed INDRA data.
void SetTriggerConditionsForRun(KVEventSelector *, Int_t, Bool_t=kFALSE) override
void handle_ebyedat_raw_data_parameter(const char *param_name, uint16_t val)
Definition: KVINDRA.cpp:417
void prepare_to_handle_new_raw_data()
reset acquisition parameters etc. before reading new raw data event
virtual void SetParameters(UInt_t n, Bool_t physics_parameters_only=kFALSE)
static KVMultiDetArray * MakeMultiDetector(const Char_t *dataset_name, Int_t run=-1, TString classname="KVMultiDetArray")
Base class for user analysis of reconstructed data.
void SetCurrentRun(KVDBRun *r) override
Event containing KVReconstructedNucleus nuclei reconstructed from hits in detectors.
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition: KVString.h:73
virtual void SetProof(Bool_t on=kTRUE, Bool_t refresh=kFALSE, Bool_t gettreeheader=kFALSE)
virtual const char * GetValue(const char *name, const char *dflt) const
virtual TObject * FindObject(const char *name) const
virtual Bool_t InheritsFrom(const char *classname) const
Bool_t IsMaster() const
void Form(const char *fmt,...)
virtual const char * BaseName(const char *pathname)
long long Long64_t
BinData::ErrorType GetDataType(const TGraph *gr, DataOptions &fitOpt)
void Error(const char *location, const char *fmt,...)
void Info(const char *location, const char *fmt,...)
void Warning(const char *location, const char *fmt,...)
Add
TArc a
ClassImp(TPyArg)