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().name << "\"." << 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  // open and add to TChain all required files
94  // we force the opening of the files to avoid problems with xrootd which sometimes
95  // seems to have a little difficulty
96  for (auto& run : GetRunList()) {
97  TString fullPathToRunfile = gDataSet->GetFullPathToRunfile(GetDataType(), run);
98  cout << "Opening file " << fullPathToRunfile << endl;
99  TFile* f = gDataSet->OpenRunfile<TFile>(GetDataType(), run);
100  cout << "Adding file " << fullPathToRunfile;
101  cout << " to the TChain." << endl;
102  dynamic_cast<TChain*>(theChain)->Add(fullPathToRunfile);
103  if (f && !f->IsZombie()) {
104  // update run infos in available runs file if necessary
106  if (ARF->InfosNeedUpdate(run, gSystem->BaseName(fullPathToRunfile))) {
107  if (!((TTree*)f->Get("ReconEvents"))) {
108  Error("SubmitTask", "No tree named ReconEvents is present in the current file");
109  delete theChain;
110  return;
111  }
112  TEnv* treeInfos = (TEnv*)((TTree*)f->Get("ReconEvents"))->GetUserInfo()->FindObject("TEnv");
113  if (treeInfos) {
114  TString kvversion = treeInfos->GetValue("KVBase::GetKVVersion()", "");
115  TString username = treeInfos->GetValue("gSystem->GetUserInfo()->fUser", "");
116  if (kvversion != "") ARF->UpdateInfos(run, gSystem->BaseName(fullPathToRunfile), kvversion, username);
117  }
118  else {
119  Info("SubmitTask", "No TEnv object associated to the tree");
120  }
121  }
122  }
123  }
124  TotalEntriesToRead = theChain->GetEntries();
125  TString option = Form("EventsReadInterval=%lld,", GetAnalysisTask()->GetStatusUpdateInterval());
126  option += Form("FullRunList=%s", GetFullRunList().GetList().Data());
127 
128  // Add any user-defined options
129  if (GetUserClassOptions() != "") {
130  option += ",";
131  option += GetUserClassOptions();
132  }
133 
134  TObject* new_selector = GetInstanceOfUserClass();
135 
136  if (!new_selector || !new_selector->InheritsFrom("TSelector")) {
137  cout << "The selector \"" << GetUserClassName() << "\" is not valid." << endl;
138  cout << "Process aborted." << endl;
139  SafeDelete(new_selector);
140  }
141  else {
142  SafeDelete(new_selector);
143  Info("SubmitTask", "Beginning TChain::Process...");
144 #ifdef WITH_CPP11
145  if (GetProofMode() != KVDataAnalyser::EProofMode::None) dynamic_cast<TChain*>(theChain)->SetProof(kTRUE);
146 #else
147  if (GetProofMode() != KVDataAnalyser::None) dynamic_cast<TChain*>(theChain)->SetProof(kTRUE);
148 #endif
149  TString analysis_class;
150  if (GetAnalysisTask()->WithUserClass())
151  analysis_class.Form("%s%s", GetUserClass().full_path_imp().Data(), GetACliCMode());
152  else analysis_class = GetUserClassName();
153 
154  if (GetNbEventToRead()) {
155  theChain->Process(analysis_class, option.Data(), GetNbEventToRead());
156  }
157  else {
158  theChain->Process(analysis_class, option.Data());
159  }
160  }
161  delete theChain;
162  fSelector = nullptr; //deleted by TChain/TTreePlayer
163 }
164 
165 
166 
182 
184 {
185  // Called by currently-processed TSelector when a new file in the TChain is opened.
186  //
187  // We call gMultiDetArray->SetParameters for the current run.
188  // Whether or not only physics parameters are set, or the full set of calibrations and identifications
189  // for each detector/identification telescope is determined by the environment variable
190  //~~~
191  // [dataset].ReconAnalysis.WithCalibInfos: [yes/no]
192  //~~~
193  // This can be overridden for any individual analysis by setting the analysis class option
194  //~~~~
195  // WithCalibInfos=yes
196  //~~~~
197  //
198  // Infos on currently read file/tree are printed.
199 
200  Bool_t physics_parameters_only = !gDataSet->GetDataSetEnv("ReconAnalysis.WithCalibInfos", kTRUE);
201  if (fSelector->IsOptGiven("WithCalibInfos"))
202  physics_parameters_only = (fSelector->GetOpt("WithCalibInfos") != "yes");
203 
204  Int_t run = GetRunNumberFromFileName(theChain->GetCurrentFile()->GetName());
206  KVMultiDetArray::MakeMultiDetector(GetDataSet()->GetName(), run);
207 
208  KVDBRun* CurrentRun = gExpDB->GetDBRun(run);
209  SetCurrentRun(CurrentRun);
210  fSelector->SetCurrentRun(CurrentRun);
211 
212  cout << endl << " =================== New Run =================== " <<
213  endl << endl;
214 
215  CurrentRun->Print();
216  if (CurrentRun->GetSystem()) {
217  if (CurrentRun->GetSystem()->GetKinematics())
218  CurrentRun->GetSystem()->GetKinematics()->Print();
219  }
220 
221  cout << endl << " ================================================= " <<
222  endl << endl;
223 
224  PrintTreeInfos();
225  Info("preInitRun", "Data written with series %s, release %d", GetDataSeries().Data(),
226  GetDataReleaseNumber());
227  fRustines.InitializePatchList(GetDataSet()->GetName(), GetDataType(), run, GetDataSeries(),
228  GetDataReleaseNumber(), theChain->GetCurrentFile()->GetStreamerInfoCache());
229  fRustines.Print();
230 }
231 
232 
233 
236 
238 {
239  // apply any required patches to data
240  if (fRustines.HasActivePatches()) fRustines.Apply(fSelector->GetEvent());
241 }
242 
243 
244 
246 
248 {
249  return (TEnv*)theChain->GetTree()->GetUserInfo()->FindObject("TEnv");
250 }
251 
252 
253 
269 
271 {
272  // When called from the InitRun() method of a user's analysis class, this method will ensure that only data
273  // compatible with the experimental trigger will be provided for analysis in the user's Analysis() method.
274  //
275  // This will be done by searching for a KVTriggerConditions plugin class defined for the currently-analysed
276  // dataset, defined like so:
277  //
278  //~~~~
279  //+Plugin.KVTriggerConditions: [dataset] [classname] [libname] "[default constructor]()"
280  //~~~~
281  //
282  // An object of the plugin class will be instantiated, and then its overridden
283  // KVTriggerConditions::SetTriggerConditionsForRun() method will be called with 2 arguments:
284  // - a pointer to the user's analysis class derived from KVReconEventSelector;
285  // - the number of the run currently being analysed
286 
287  TPluginHandler* ph = KVBase::LoadPlugin("KVTriggerConditions", GetDataSet()->GetName());
288  if (!ph) {
289  Info("SetTriggerConditionsForRun",
290  "No definition of trigger conditions available for dataset %s",
291  GetDataSet()->GetName());
292  return;
293  }
294  std::unique_ptr<KVTriggerConditions> trig((KVTriggerConditions*) ph->ExecPlugin(0));
295  trig->SetTriggerConditionsForRun(fSelector, run);
296 }
297 
298 
299 
302 
304 {
305  // Print informations on currently analysed TTree
306  TEnv* treeInfos = GetReconDataTreeInfos();
307  if (!treeInfos) return;
308  cout << endl << "----------------------------------------------------------------------------------------------------" << endl;
309  cout << "INFORMATIONS ON VERSION OF KALIVEDA USED TO GENERATE FILE:" << endl << endl;
310  fDataVersion = treeInfos->GetValue("KVBase::GetKVVersion()", "(unknown)");
311  cout << "version = " << fDataVersion << endl ;
312  cout << "build date = " << treeInfos->GetValue("KVBase::GetKVBuildDate()", "(unknown)") << endl ;
313  cout << "source directory = " << treeInfos->GetValue("KVBase::GetKVSourceDir()", "(unknown)") << endl ;
314  cout << "KVROOT = " << treeInfos->GetValue("KVBase::GetKVRoot()", "(unknown)") << endl ;
315  if (strcmp(treeInfos->GetValue("KVBase::bzrBranchNick()", "(unknown)"), "(unknown)")) {
316  cout << "BZR branch name = " << treeInfos->GetValue("KVBase::bzrBranchNick()", "(unknown)") << endl ;
317  cout << "BZR revision #" << treeInfos->GetValue("KVBase::bzrRevisionNumber()", "(unknown)") << endl ;
318  cout << "BZR revision ID = " << treeInfos->GetValue("KVBase::bzrRevisionId()", "(unknown)") << endl ;
319  cout << "BZR revision date = " << treeInfos->GetValue("KVBase::bzrRevisionDate()", "(unknown)") << endl ;
320  }
321  else if (strcmp(treeInfos->GetValue("KVBase::gitBranch()", "(unknown)"), "(unknown)")) {
322  cout << "git branch = " << treeInfos->GetValue("KVBase::gitBranch()", "(unknown)") << endl;
323  cout << "git commit = " << treeInfos->GetValue("KVBase::gitCommit()", "(unknown)") << endl;
324  }
325  cout << endl << "INFORMATIONS ON GENERATION OF FILE:" << endl << endl;
326  cout << "Generated by : " << treeInfos->GetValue("gSystem->GetUserInfo()->fUser", "(unknown)") << endl ;
327  cout << "Analysis task : " << treeInfos->GetValue("AnalysisTask", "(unknown)") << endl ;
328  cout << "Job name : " << treeInfos->GetValue("BatchSystem.JobName", "(unknown)") << endl ;
329  cout << "Job submitted from : " << treeInfos->GetValue("LaunchDirectory", "(unknown)") << endl ;
330  cout << "Runs : " << treeInfos->GetValue("Runs", "(unknown)") << endl ;
331  cout << "Number of events requested : " << treeInfos->GetValue("NbToRead", "(unknown)") << endl ;
332  cout << endl << "----------------------------------------------------------------------------------------------------" << endl;
333 
334  // if possible, parse fDataVersion into series and release number
335  // e.g. if fDataVersion = "1.8.10":
336  // => fDataSeries = "1.8" fDataReleaseNum = 10
337  // e.g. if fDataVersion = "1.11/00":
338  // => fDataSeries = "1.11" fDataReleaseNum = 0
339  Int_t a, b, c;
340  if (fDataVersion != "(unknown)") {
341  if (fDataVersion.Contains("/")) {
342  if (sscanf(fDataVersion.Data(), "%d.%d/%d", &a, &b, &c) == 3) {
343  fDataSeries.Form("%d.%d", a, b);
344  fDataReleaseNum = c;
345  }
346  }
347  else {
348  if (sscanf(fDataVersion.Data(), "%d.%d.%d", &a, &b, &c) == 3) {
349  fDataSeries.Form("%d.%d", a, b);
350  fDataReleaseNum = c;
351  }
352  }
353  }
354  else {
355  fDataSeries = "";
356  fDataReleaseNum = -1;
357  }
358 }
359 
360 
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 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)
static TPluginHandler * LoadPlugin(const Char_t *base, const Char_t *uri="0")
Definition: KVBase.cpp:772
Description of an experimental run in database ,,.
Definition: KVDBRun.h:40
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 Bool_t CheckTaskVariables()
void Reset() 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
KVDBRun * GetDBRun(Int_t number) const
Definition: KVExpDB.h:89
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,...)
Add
TArc a
ClassImp(TPyArg)