KaliVeda
Toolkit for HIC analysis
Loading...
Searching...
No Matches
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
11using namespace std;
12
14
15
16
17
18
20
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
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
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
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
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 preAnalysis();
120 //call user's analysis. stop if returns kFALSE.
121 if (!Analysis()) break;
122 postAnalysis();
123
125
126 fEventNumber += 1;
127 }
128
129 cout << "Ending analysis of run " << fRunNumber << " on : ";
131 cout << now2.AsString() << endl << endl;
132 cout << endl << "Finished reading " << fEventNumber - 1 << " events from file " << raw_file.Data() << endl << endl;
133
134 preEndRun();
135 //call user's end of run function
136 EndRun();
137 postEndRun();
138
139 delete fRunFile;
140}
141
142
143
144
156
158{
159 // Perform analysis of chosen runs
160 // Before beginning the loop over the runs, the user's InitAnalysis() method is called.
161 // After completing the analysis of all runs, the user's EndAnalysis() method is called.
162 //
163 // Further customisation of the event loop is possible by overriding the methods
164 // preInitAnalysis()
165 // postInitAnalysis()
166 // preEndAnalysis()
167 // postEndAnalysis()
168 // which are executed respectively just before and just after the methods.
169
170 if (gDataSet != GetDataSet()) GetDataSet()->cd();
171
172 // prepare any user-supplied options for the analysis
174
176 //call user's initialisation
177 InitAnalysis();
179
181
182 //loop over runs
183 GetRunList().Begin();
184 while (!GetRunList().End() && !AbortProcessingLoop()) {
186 ProcessRun();
187 }
188
189 if (fCombinedOutputFile != "") {
190 Info("Terminate", "combine = %s", fCombinedOutputFile.Data());
191 // combine histograms and trees from analysis into one file
193 file1.Form("HistoFileFrom%s.root", ClassName());
194 file2.Form("TreeFileFrom%s.root", ClassName());
195 if (fHistoList.GetEntries()) {
196 if (fTreeFile) {
197 Info("Terminate", "both");
198 SaveHistos();
199 fTreeFile->Write();
200 delete fTreeFile;
202 }
203 else {
204 // no trees - just rename histo file
205 Info("Terminate", "histo");
207 }
208 }
209 else if (fTreeFile) {
210 // no histos - just rename tree file
211 Info("Terminate", "tree");
212 fTreeFile->Write();
213 delete fTreeFile;
215 }
216 else Info("Terminate", "none");
217 }
218 else {
219 SaveHistos();
220 if (fTreeFile) {
221 fTreeFile->Write();
222 delete fTreeFile;
223 }
224 }
225
227 //call user's end of analysis
228 EndAnalysis();
230}
231
232
233
234
237
239{
240 //loop over runs and calculate total events
242 GetRunList().Begin();
243 while (!GetRunList().End()) {
244 Int_t r = GetRunList().Next();
245 TotalEntriesToRead += gExpDB->GetDBRun(r)->GetEvents();
246 }
247}
248
249
250
255
257{
258 // Add a user analysis results TTree to list which will be automatically
259 // written to disk at end of analysis.
260 // User must call CreateTreeFile(const Char_t*) before calling this method
261
262 if (!fTreeFile) {
263 Error("AddTree", "You must call CreateTreeFile(const Char_t*) before using this method!");
264 return;
265 }
266 tree->SetDirectory(fTreeFile);
267 tree->AutoSave();
269}
270
271
272
278
279void KVRawDataAnalyser::FillHisto(const Char_t* histo_name, Double_t one, Double_t two, Double_t three, Double_t four)
280{
281 //Find in the list, if there is an histogram named "histo_name"
282 //If not print an error message
283 //If yes redirect to the right method according to its closest mother class
284 //to fill it
285
286 TH1* h1 = 0;
287 if ((h1 = GetHisto(histo_name))) {
288 if (h1->InheritsFrom("TH3"))
289 FillTH3((TH3*)h1, one, two, three, four);
290 else if (h1->InheritsFrom("TProfile2D"))
292 else if (h1->InheritsFrom("TH2"))
293 FillTH2((TH2*)h1, one, two, three);
294 else if (h1->InheritsFrom("TProfile"))
296 else if (h1->InheritsFrom("TH1"))
297 FillTH1(h1, one, two);
298 else
299 Warning("FillHisto", "%s -> Classe non prevue ...", fHistoList.FindObject(histo_name)->ClassName());
300 }
301 else {
302 Warning("FillHisto", "%s introuvable", histo_name);
303 }
304}
305
306
307
310
311void KVRawDataAnalyser::FillHisto(const Char_t* histo_name, const Char_t* label, Double_t weight)
312{
313 // Fill 1D histogram with named bins
314 TH1* h1 = 0;
315 if ((h1 = GetHisto(histo_name))) {
316 h1->Fill(label, weight);
317 }
318 else {
319 Warning("FillHisto", "%s introuvable", histo_name);
320 }
321}
322
323
324
331
333{
334 //Filltree method, the tree named tree_name
335 //has to be declared with AddTTree(TTree*) method
336 //
337 //if no sname="", all trees in the list is filled
338 //
339 if (!strcmp(tree_name, "")) {
340 fTreeList.Execute("Fill", "");
341 }
342 else {
343 TTree* tt = 0;
344 if ((tt = GetTree(tree_name))) {
345 tt->Fill();
346 }
347 else {
348 Warning("FillTree", "%s introuvable", tree_name);
349 }
350 }
351}
352
353
354
368
369void KVRawDataAnalyser::SaveHistos(const Char_t* filename, Option_t* option, Bool_t onlyfilled)
370{
371 // Write in file all histograms declared with AddHisto(TH1*)
372 //
373 // If no filename is specified, set default name : HistoFileFrom[name_of_class].root
374 //
375 // If a filename is specified, search in gROOT->GetListOfFiles() if
376 // this file has been already opened
377 // - if yes write in it
378 // - if not, create it with the corresponding option, write in it
379 // and close it just after
380 //
381 // onlyfilled flag allow to write all (onlyfilled=kFALSE, default)
382 // or only histograms (onlyfilled=kTRUE) those have been filled
383
384 if (!fHistoList.GetEntries()) return;
385
387 if (!strcmp(filename, ""))
388 histo_file_name.Form("HistoFileFrom%s.root", ClassName());
389 else
391
393
394 TFile* file = 0;
395 TDirectory* pwd = gDirectory;
396 //if filename correspond to an already opened file, write in it
397 //if not open/create it, depending on the option ("recreate" by default)
398 //and write in it
399 if (!(file = (TFile*)gROOT->GetListOfFiles()->FindObject(histo_file_name.Data()))) {
400 file = new TFile(histo_file_name.Data(), option);
402 }
403 file->cd();
404 TIter next(GetHistoList());
405 TObject* obj = 0;
406 while ((obj = next())) {
407 if (obj->InheritsFrom("TH1")) {
408 if (onlyfilled) {
409 if (((TH1*)obj)->GetEntries() > 0) {
410 obj->Write();
411 }
412 }
413 else {
414 obj->Write();
415 }
416 }
417 }
418 if (justopened)
419 file->Close();
420 pwd->cd();
421}
422
423
424
428
430{
431 // This method must be called before creating any user TTree in InitAnalysis().
432 // If no filename is given, default name="TreeFileFrom[name of analysis class].root"
433
434 TString tree_file_name;
435 if (!strcmp(filename, ""))
436 tree_file_name.Form("TreeFileFrom%s.root", ClassName());
437 else
438 tree_file_name = filename;
439
441 fTreeFile = new TFile(tree_file_name, "RECREATE");
442 savedir->cd();
443
444 return kTRUE;
445}
446
447
448
449
452
453void KVRawDataAnalyser::Make(const Char_t* kvsname)
454{
455 //Automatic generation of derived class for raw data analysis
456
457 KVClassFactory cf(kvsname, "Analysis of raw data", "",
458 kTRUE, "RawAnalysisTemplate");
459
460 cf.AddImplIncludeFile("KVMultiDetArray.h");
461 cf.GenerateCode();
462}
463
464
465
468
470{
471 // Method called to abort analysis during processing of a run
472}
473
474
int Int_t
ROOT::R::TRInterface & r
bool Bool_t
char Char_t
constexpr Bool_t kFALSE
double Double_t
constexpr Bool_t kTRUE
const char Option_t
#define gDirectory
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 Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char filename
#define gROOT
R__EXTERN TSystem * gSystem
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
virtual void Print(Option_t *option="") const
Definition KVDBRun.cpp:69
const KVString & GetUserClassOptions() const
virtual void postEndRun()
void DoStatusUpdate(Long64_t nevents) const
Print infos on events treated, disk usage, memory usage.
virtual void preInitAnalysis()
virtual void postAnalysis()
virtual Bool_t CheckStatusUpdateInterval(Long64_t nevents) const
virtual void preAnalysis()
const KVString & GetDataType() const
virtual void preEndRun()
virtual void postInitRun()
static Bool_t AbortProcessingLoop()
Long64_t GetNbEventToRead(void) const
virtual void postEndAnalysis()
virtual void preEndAnalysis()
virtual void preInitRun()
virtual void postInitAnalysis()
const KVNumberList & GetRunList() const
const KVDataSet * GetDataSet() const
FileType * OpenRunfile(const Char_t *type, Int_t run)
Definition KVDataSet.h:167
TString GetFullPathToRunfile(const Char_t *type, Int_t run) const
void cd() const
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 *)
void Begin(void) const
Int_t Next(void) const
Abstract base class for user analysis of raw data.
virtual Bool_t Analysis()=0
const KVHashList * GetHistoList() const
KVHashList fHistoList
list of histograms of user analysis
KVRawDataReader * fRunFile
currently analysed run file
void FillTProfile(TProfile *h1, Double_t one, Double_t two, Double_t three)
virtual void InitRun()=0
Int_t fRunNumber
run number of current file
virtual void SaveHistos(const Char_t *filename="", Option_t *option="recreate", Bool_t onlyfilled=kFALSE)
void FillTProfile2D(TProfile2D *h2, Double_t one, Double_t two, Double_t three, Double_t four)
KVHashList fTreeList
list of trees of user analysis
KVString fCombinedOutputFile
optional name for single results file with trees and histos
Long64_t fEventNumber
event number in current run
virtual void EndAnalysis()=0
Bool_t CreateTreeFile(const Char_t *filename="")
void FillTree(const Char_t *sname="")
TTree * GetTree(const Char_t *name) const
void FillTH3(TH3 *h3, Double_t one, Double_t two, Double_t three, Double_t four)
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 FillTH1(TH1 *h1, Double_t one, Double_t two)
void AddHisto(TH1 *histo)
TH1 * GetHisto(const Char_t *name) const
virtual void InitAnalysis()=0
void CalculateTotalEventsToRead()
loop over runs and calculate total events
KVDBRun * fCurrentRun
poiner to current run
KVUserAnalysisOptionList fOptionList
list of options set by user for analysis
void FillTH2(TH2 *h2, Double_t one, Double_t two, Double_t three)
virtual void EndRun()=0
void AddTree(TTree *tree)
void AbortDuringRunProcessing()
Method called to abort analysis during processing of a run.
Abstract base class for reading raw (DAQ) data.
virtual Bool_t GetNextEvent()=0
virtual TObject * FindObject(const char *name) const
virtual void Execute(const char *method, const char *params, Int_t *error=0)
virtual void Add(TObject *obj)
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition KVString.h:73
void ParseOptions(const KVString &optlist)
virtual Int_t GetEntries() const
const char * AsString() const
virtual Bool_t cd()
Int_t Write(const char *name=nullptr, Int_t opt=0, Int_t bufsiz=0) const override
virtual Int_t Fill(const char *name, Double_t w)
const char * GetName() const override
virtual const char * ClassName() const
virtual void Warning(const char *method, const char *msgfmt,...) const
R__ALWAYS_INLINE Bool_t IsZombie() const
virtual Int_t Write(const char *name=nullptr, Int_t option=0, Int_t bufsize=0)
virtual Bool_t InheritsFrom(const char *classname) const
virtual void Error(const char *method, const char *msgfmt,...) const
virtual void Info(const char *method, const char *msgfmt,...) const
const char * Data() const
void Form(const char *fmt,...)
virtual int Rename(const char *from, const char *to)
long long Long64_t
TH1F * h1
str tree_name
void End()
auto * tt
ClassImp(TPyArg)