KaliVeda
Toolkit for HIC analysis
KVReconDataAnalyser.cpp
1 //Created by KVClassFactory on Thu Jul 19 2018
2 //Author: John Frankland
3 
4 #include "KVBase.h"
5 #include "KVReconDataAnalyser.h"
6 #include "KVDBRun.h"
7 #include "KVMultiDetArray.h"
8 #include "KVDataAnalysisTask.h"
9 #include "KVDataSet.h"
10 #include "TChain.h"
11 #include "TObjString.h"
12 #include "TChain.h"
13 #include "KVAvailableRunsFile.h"
14 #include "TProof.h"
15 
16 #include <KVClassFactory.h>
17 #include <KVTriggerConditions.h>
18 
19 using namespace std;
20 
22 
23 
24 
27  : fSelector(nullptr), theChain(nullptr)
28 {
29 }
30 
31 
32 
35 
37 {
38  //Reset task variables
40  theChain = nullptr;
41  fSelector = nullptr;
42  TotalEntriesToRead = 0;
43 }
44 
45 
46 
47 
48 
49 
52 
54 {
55  // Checks the task variables
56 
58 
59  cout << "============> Analysis summary <=============" << endl;
60  cout << "Analysis of runs " << GetRunList().GetList() << " with the class ";
61  cout << "\"" << GetUserClass() << "\"." << endl;
62  if (GetNbEventToRead()) {
63  cout << GetNbEventToRead() << " events will be processed." << endl;
64  }
65  else {
66  cout << "All events will be processed." << endl;
67  }
68  cout << "=============================================" << endl;
69 
70  return kTRUE;
71 }
72 
73 
74 
75 
80 
82 {
83  //Run the interactive analysis
84 
85  //make the chosen dataset the active dataset ( = gDataSet; note this also opens database
86  //and positions gDataBase & gExpDB).
87  GetDataSet()->cd();
88  fSelector = nullptr;
89 
90  theChain = new TChain("ReconEvents");
91  theChain->SetDirectory(0); // we handle delete
92 
93  GetRunList().Begin();
94  Int_t run;
95 
96  // open and add to TChain all required files
97  // we force the opening of the files to avoid problems with xrootd which sometimes
98  // seems to have a little difficulty
99  while (!GetRunList().End()) {
100  run = GetRunList().Next();
101  TString fullPathToRunfile = gDataSet->GetFullPathToRunfile(GetDataType(), run);
102  cout << "Opening file " << fullPathToRunfile << endl;
103  TFile* f = gDataSet->OpenRunfile<TFile>(GetDataType(), run);
104  cout << "Adding file " << fullPathToRunfile;
105  cout << " to the TChain." << endl;
106  dynamic_cast<TChain*>(theChain)->Add(fullPathToRunfile);
107  if (f && !f->IsZombie()) {
108  // update run infos in available runs file if necessary
110  if (ARF->InfosNeedUpdate(run, gSystem->BaseName(fullPathToRunfile))) {
111  if (!((TTree*)f->Get("ReconEvents"))) {
112  Error("SubmitTask", "No tree named ReconEvents is present in the current file");
113  delete theChain;
114  return;
115  }
116  TEnv* treeInfos = (TEnv*)((TTree*)f->Get("ReconEvents"))->GetUserInfo()->FindObject("TEnv");
117  if (treeInfos) {
118  TString kvversion = treeInfos->GetValue("KVBase::GetKVVersion()", "");
119  TString username = treeInfos->GetValue("gSystem->GetUserInfo()->fUser", "");
120  if (kvversion != "") ARF->UpdateInfos(run, gSystem->BaseName(fullPathToRunfile), kvversion, username);
121  }
122  else {
123  Info("SubmitTask", "No TEnv object associated to the tree");
124  }
125  }
126  }
127  }
128  TotalEntriesToRead = theChain->GetEntries();
129  TString option = Form("EventsReadInterval=%lld,", GetAnalysisTask()->GetStatusUpdateInterval());
130  option += Form("FullRunList=%s", GetFullRunList().GetList());
131 
132  // Add any user-defined options
133  if (GetUserClassOptions() != "") {
134  option += ",";
135  option += GetUserClassOptions();
136  }
137 
138  TObject* new_selector = GetInstanceOfUserClass();
139 
140  if (!new_selector || !new_selector->InheritsFrom("TSelector")) {
141  cout << "The selector \"" << GetUserClass() << "\" is not valid." << endl;
142  cout << "Process aborted." << endl;
143  SafeDelete(new_selector);
144  }
145  else {
146  SafeDelete(new_selector);
147  Info("SubmitTask", "Beginning TChain::Process...");
148 #ifdef WITH_CPP11
149  if (GetProofMode() != KVDataAnalyser::EProofMode::None) dynamic_cast<TChain*>(theChain)->SetProof(kTRUE);
150 #else
151  if (GetProofMode() != KVDataAnalyser::None) dynamic_cast<TChain*>(theChain)->SetProof(kTRUE);
152 #endif
153  TString analysis_class;
154  if (GetAnalysisTask()->WithUserClass()) analysis_class.Form("%s%s", GetUserClassImp().Data(), GetACliCMode());
155  else analysis_class = GetUserClass();
156 
157  if (GetNbEventToRead()) {
158  theChain->Process(analysis_class, option.Data(), GetNbEventToRead());
159  }
160  else {
161  theChain->Process(analysis_class, option.Data());
162  }
163  }
164  delete theChain;
165  fSelector = nullptr; //deleted by TChain/TTreePlayer
166 }
167 
168 
169 
185 
187 {
188  // Called by currently-processed TSelector when a new file in the TChain is opened.
189  //
190  // We call gMultiDetArray->SetParameters for the current run.
191  // Whether or not only physics parameters are set, or the full set of calibrations and identifications
192  // for each detector/identification telescope is determined by the environment variable
193  //~~~
194  // [dataset].ReconAnalysis.WithCalibInfos: [yes/no]
195  //~~~
196  // This can be overridden for any individual analysis by setting the analysis class option
197  //~~~~
198  // WithCalibInfos=yes
199  //~~~~
200  //
201  // Infos on currently read file/tree are printed.
202 
203  Bool_t physics_parameters_only = !gDataSet->GetDataSetEnv("ReconAnalysis.WithCalibInfos", kTRUE);
204  if (fSelector->IsOptGiven("WithCalibInfos"))
205  physics_parameters_only = (fSelector->GetOpt("WithCalibInfos") != "yes");
206 
207  Int_t run = GetRunNumberFromFileName(theChain->GetCurrentFile()->GetName());
209  KVMultiDetArray::MakeMultiDetector(GetDataSet()->GetName(), run);
210 
211  KVDBRun* CurrentRun = gExpDB->GetDBRun(run);
212  SetCurrentRun(CurrentRun);
213  fSelector->SetCurrentRun(CurrentRun);
214 
215  cout << endl << " =================== New Run =================== " <<
216  endl << endl;
217 
218  CurrentRun->Print();
219  if (CurrentRun->GetSystem()) {
220  if (CurrentRun->GetSystem()->GetKinematics())
221  CurrentRun->GetSystem()->GetKinematics()->Print();
222  }
223 
224  cout << endl << " ================================================= " <<
225  endl << endl;
226 
227  PrintTreeInfos();
228  Info("preInitRun", "Data written with series %s, release %d", GetDataSeries().Data(),
229  GetDataReleaseNumber());
230  fRustines.InitializePatchList(GetDataSet()->GetName(), GetDataType(), run, GetDataSeries(),
231  GetDataReleaseNumber(), theChain->GetCurrentFile()->GetStreamerInfoCache());
232  fRustines.Print();
233 }
234 
235 
236 
239 
241 {
242  // apply any required patches to data
243  if (fRustines.HasActivePatches()) fRustines.Apply(fSelector->GetEvent());
244 }
245 
246 
247 
249 
251 {
252  return (TEnv*)theChain->GetTree()->GetUserInfo()->FindObject("TEnv");
253 }
254 
255 
256 
272 
274 {
275  // When called from the InitRun() method of a user's analysis class, this method will ensure that only data
276  // compatible with the experimental trigger will be provided for analysis in the user's Analysis() method.
277  //
278  // This will be done by searching for a KVTriggerConditions plugin class defined for the currently-analysed
279  // dataset, defined like so:
280  //
281  //~~~~
282  //+Plugin.KVTriggerConditions: [dataset] [classname] [libname] "[default constructor]()"
283  //~~~~
284  //
285  // An object of the plugin class will be instantiated, and then its overridden
286  // KVTriggerConditions::SetTriggerConditionsForRun() method will be called with 2 arguments:
287  // - a pointer to the user's analysis class derived from KVReconEventSelector;
288  // - the number of the run currently being analysed
289 
290  TPluginHandler* ph = KVBase::LoadPlugin("KVTriggerConditions", GetDataSet()->GetName());
291  if (!ph) {
292  Info("SetTriggerConditionsForRun",
293  "No definition of trigger conditions available for dataset %s",
294  GetDataSet()->GetName());
295  return;
296  }
297  std::unique_ptr<KVTriggerConditions> trig((KVTriggerConditions*) ph->ExecPlugin(0));
298  trig->SetTriggerConditionsForRun(fSelector, run);
299 }
300 
301 
302 
305 
307 {
308  // Print informations on currently analysed TTree
309  TEnv* treeInfos = GetReconDataTreeInfos();
310  if (!treeInfos) return;
311  cout << endl << "----------------------------------------------------------------------------------------------------" << endl;
312  cout << "INFORMATIONS ON VERSION OF KALIVEDA USED TO GENERATE FILE:" << endl << endl;
313  fDataVersion = treeInfos->GetValue("KVBase::GetKVVersion()", "(unknown)");
314  cout << "version = " << fDataVersion << endl ;
315  cout << "build date = " << treeInfos->GetValue("KVBase::GetKVBuildDate()", "(unknown)") << endl ;
316  cout << "source directory = " << treeInfos->GetValue("KVBase::GetKVSourceDir()", "(unknown)") << endl ;
317  cout << "KVROOT = " << treeInfos->GetValue("KVBase::GetKVRoot()", "(unknown)") << endl ;
318  if (strcmp(treeInfos->GetValue("KVBase::bzrBranchNick()", "(unknown)"), "(unknown)")) {
319  cout << "BZR branch name = " << treeInfos->GetValue("KVBase::bzrBranchNick()", "(unknown)") << endl ;
320  cout << "BZR revision #" << treeInfos->GetValue("KVBase::bzrRevisionNumber()", "(unknown)") << endl ;
321  cout << "BZR revision ID = " << treeInfos->GetValue("KVBase::bzrRevisionId()", "(unknown)") << endl ;
322  cout << "BZR revision date = " << treeInfos->GetValue("KVBase::bzrRevisionDate()", "(unknown)") << endl ;
323  }
324  else if (strcmp(treeInfos->GetValue("KVBase::gitBranch()", "(unknown)"), "(unknown)")) {
325  cout << "git branch = " << treeInfos->GetValue("KVBase::gitBranch()", "(unknown)") << endl;
326  cout << "git commit = " << treeInfos->GetValue("KVBase::gitCommit()", "(unknown)") << endl;
327  }
328  cout << endl << "INFORMATIONS ON GENERATION OF FILE:" << endl << endl;
329  cout << "Generated by : " << treeInfos->GetValue("gSystem->GetUserInfo()->fUser", "(unknown)") << endl ;
330  cout << "Analysis task : " << treeInfos->GetValue("AnalysisTask", "(unknown)") << endl ;
331  cout << "Job name : " << treeInfos->GetValue("BatchSystem.JobName", "(unknown)") << endl ;
332  cout << "Job submitted from : " << treeInfos->GetValue("LaunchDirectory", "(unknown)") << endl ;
333  cout << "Runs : " << treeInfos->GetValue("Runs", "(unknown)") << endl ;
334  cout << "Number of events requested : " << treeInfos->GetValue("NbToRead", "(unknown)") << endl ;
335  cout << endl << "----------------------------------------------------------------------------------------------------" << endl;
336 
337  // if possible, parse fDataVersion into series and release number
338  // e.g. if fDataVersion = "1.8.10":
339  // => fDataSeries = "1.8" fDataReleaseNum = 10
340  // e.g. if fDataVersion = "1.11/00":
341  // => fDataSeries = "1.11" fDataReleaseNum = 0
342  Int_t a, b, c;
343  if (fDataVersion != "(unknown)") {
344  if (fDataVersion.Contains("/")) {
345  if (sscanf(fDataVersion.Data(), "%d.%d/%d", &a, &b, &c) == 3) {
346  fDataSeries.Form("%d.%d", a, b);
347  fDataReleaseNum = c;
348  }
349  }
350  else {
351  if (sscanf(fDataVersion.Data(), "%d.%d.%d", &a, &b, &c) == 3) {
352  fDataSeries.Form("%d.%d", a, b);
353  fDataReleaseNum = c;
354  }
355  }
356  }
357  else {
358  fDataSeries = "";
359  fDataReleaseNum = -1;
360  }
361 }
362 
363 
int Int_t
#define SafeDelete(p)
#define f(i)
#define c(i)
bool Bool_t
constexpr Bool_t kFALSE
constexpr Bool_t kTRUE
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
char * Form(const char *fmt,...)
R__EXTERN TSystem * gSystem
void Print(Option_t *opt="") const
Definition: KV2Body.cpp:813
Handles lists of available runs for different datasets and types of data.
virtual void UpdateInfos(Int_t run, const Char_t *filename, const Char_t *kvversion, const Char_t *username)
virtual Bool_t InfosNeedUpdate(Int_t run, const Char_t *filename)
static TPluginHandler * LoadPlugin(const Char_t *base, const Char_t *uri="0")
Definition: KVBase.cpp:793
Description of an experimental run in database ,,.
Definition: KVDBRun.h:36
KVDBSystem * GetSystem() const
Definition: KVDBRun.cpp:242
virtual void Print(Option_t *option="") const
Definition: KVDBRun.cpp:69
KV2Body * GetKinematics()
Definition: KVDBSystem.cpp:80
virtual Bool_t CheckTaskVariables()
KVAvailableRunsFile * GetAvailableRunsFile(const Char_t *type) const
Definition: KVDataSet.cpp:50
FileType * OpenRunfile(const Char_t *type, Int_t run)
Definition: KVDataSet.h:167
const Char_t * GetDataSetEnv(const Char_t *type, const Char_t *defval="") const
Definition: KVDataSet.cpp:767
TString GetFullPathToRunfile(const Char_t *type, Int_t run) const
Definition: KVDataSet.cpp:897
KVDBRun * GetDBRun(Int_t number) const
Definition: KVExpDB.h:76
static void SetMakeMultiDetectorPhysicsParametersOnly(Bool_t on=kTRUE)
static KVMultiDetArray * MakeMultiDetector(const Char_t *dataset_name, Int_t run=-1, TString classname="KVMultiDetArray")
Manages user analysis of reconstructed experimental data.
void preAnalysis() override
apply any required patches to data
TEnv * GetReconDataTreeInfos() const
void PrintTreeInfos()
Print informations on currently analysed TTree.
Bool_t CheckTaskVariables(void) override
Checks the task variables.
virtual void SetTriggerConditionsForRun(int) override
void SubmitTask() override
void Reset() override
Reset task variables.
Set trigger conditions for analysis of reconstructed data.
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
Longptr_t ExecPlugin(int nargs)
void Form(const char *fmt,...)
virtual const char * BaseName(const char *pathname)
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 End()
Add
TArc a
ClassImp(TPyArg)