KaliVeda  1.12/06
Heavy-Ion Analysis Toolkit
KVEventSelector.h
Go to the documentation of this file.
1 #ifndef KVEventSelector_h
2 #define KVEventSelector_h
3 
4 #include <TROOT.h>
5 #include <TChain.h>
6 #include <TFile.h>
7 #include <TSelector.h>
8 #include "KVEvent.h"
9 #include "KVGVList.h"
10 #include "KVString.h"
11 #include "KVParticleCondition.h"
12 #include "TDatime.h"
13 #include "KVHashList.h"
14 #include <TH3.h>
15 #include <TH2.h>
16 #include <TProfile2D.h>
17 #include "KVNameValueList.h"
18 #include "TProofOutputFile.h"
19 #include "KVDataAnalyser.h"
20 
153 class KVEventSelector : public TSelector {
154 
155 protected :
158 
162 
165 
168 
170 
173 
175 
179 
182 
185 
188 
190 
191  void add_histo(TH1* histo);
193  {
196 
197  if (fDisableCreateTreeFile) return;
198  ltree->Add(tree);
199  }
200  void FillTH1(TH1* h1, Double_t one, Double_t two);
201  void FillTProfile(TProfile* h1, Double_t one, Double_t two, Double_t three);
202  void FillTH2(TH2* h2, Double_t one, Double_t two, Double_t three);
203  void FillTProfile2D(TProfile2D* h2, Double_t one, Double_t two, Double_t three, Double_t four);
204  void FillTH3(TH3* h3, Double_t one, Double_t two, Double_t three, Double_t four);
205 
206  void SetUpAuxEventChain();
207  void SetCombinedOutputFile(const TString& filename)
208  {
223  fCombinedOutputFile = filename;
224  }
225 public:
228  Bool_t CreateTreeFile(const Char_t* filename = "");
229 
230  virtual void ParseOptions();
231 
232  KVEventSelector(TTree* /*tree*/ = 0) : fChain(0), fAuxChain(0), fBranchName("data"), fFirstEvent(kTRUE),
234  {
235  lhisto = new KVHashList();
236  ltree = new KVHashList();
237  }
239  {
240  lhisto->Clear();
241  delete lhisto;
242  lhisto = 0;
243  ltree->Clear();
244  delete ltree;
245  ltree = 0;
246  }
248  {
250  }
251  virtual Int_t Version() const
252  {
253  return 3;
254  }
255  virtual void Begin(TTree* tree);
256  virtual void SlaveBegin(TTree* tree);
257  virtual void Init(TTree* tree);
258  void InitFriendTree(TTree* tree, const TString& branchname);
259  virtual Bool_t Notify();
260  virtual Bool_t Process(Long64_t entry);
261  virtual void CheckEndOfRun();
262  virtual Int_t GetEntry(Long64_t entry, Int_t getall = 0)
263  {
264  return fChain ? fChain->GetTree()->GetEntry(entry, getall) : 0;
265  }
267  {
268  return fAuxChain ? fAuxChain->GetTree()->GetEntry(entry, getall) : 0;
269  }
271  {
272  return AuxEvent;
273  }
274  virtual void SetObject(TObject* obj)
275  {
276  fObject = obj;
277  }
278  virtual void SetInputList(TList* input)
279  {
280  fInput = input;
281  }
282  virtual TList* GetOutputList() const
283  {
284  return fOutput;
285  }
286  virtual void SlaveTerminate();
287  virtual void Terminate();
288 
289  void SetBranchName(const Char_t* n)
290  {
291  fBranchName = n;
292  }
293  const Char_t* GetBranchName() const
294  {
295  return fBranchName.Data();
296  }
297  virtual void SetAnalysisFrame()
298  {
304  }
305  virtual void SetCurrentRun(KVDBRun*) {}
306 
307  KVEvent* GetEvent() const
308  {
309  return Event;
310  }
311  void SetEvent(KVEvent* e)
312  {
313  Event = e;
314  }
316  {
317  if (!GetEvent()) {
318  Error("GetEventNumber", "No event defined!!!");
319  return -1;
320  }
321  return GetEvent()->GetNumber();
322  }
323 
324  /* user entry points */
325  virtual void InitAnalysis()
326  {
327  AbstractMethod("InitAnalysis");
328  }
329  virtual void InitRun()
330  {
331  AbstractMethod("InitRun");
332  }
333  virtual Bool_t Analysis()
334  {
336  return kTRUE;
337  }
338  virtual void EndRun()
339  {
340  AbstractMethod("EndRun");
341  }
342  virtual void EndAnalysis()
343  {
344  AbstractMethod("EndAnalysis");
345  }
347  {
349  return &gvlist;
350  }
351  const KVGVList* GetGVList(void) const
352  {
354  return &gvlist;
355  }
356  void AddGV(KVVarGlob* vg)
357  {
360  if (!vg)
361  Error("AddGV(KVVarGlob*)", "KVVarGlob pointer is null");
362  else
363  GetGVList()->Add(vg);
364  }
365  KVVarGlob* AddGV(const Char_t* class_name, const Char_t* name);
366  KVVarGlob* GetGV(const Char_t* name) const
367  {
371 
372  return GetGVList()->GetGV(name);
373  }
374  virtual void RecalculateGlobalVariables();
376  {
378  return (fTreeEntry + 1 == fChain->GetTree()->GetEntries());
379  }
380 
381  void SetParticleConditions(const KVParticleCondition&, const KVString& = "");
383  {
384  fPartName = t;
385  }
386 
387  void AddHisto(TH1* histo)
388  {
389  Deprecate(Form("You should use e.g. 'auto h = AddHisto<%s>(\"%s\",\"%s\",%d,...);' to add histograms to your analysis.",
390  histo->ClassName(), histo->GetName(), histo->GetTitle(), histo->GetNbinsX()));
391  add_histo(histo);
392  }
393 
394  template<typename HistoType, typename... Args>
395  HistoType* AddHisto(Args&& ... args)
396  {
409 
410  auto h = new HistoType(std::forward<Args>(args)...);
411  add_histo(h);
412  return h;
413  }
414 
416  {
417  Deprecate(Form("You should use e.g. 'auto t = AddTree(\"%s\", \"%s\");' to add a TTree to your analysis",
418  tree->GetName(), tree->GetTitle()));
419  add_tree(tree);
420  }
421  TTree* AddTree(const TString& name, const TString& title = "", Int_t splitLevel = 99, TDirectory* = gDirectory);
422 
423  void FillHisto(const Char_t* sname, Double_t one, Double_t two = 1, Double_t three = 1, Double_t four = 1);
424  void FillHisto(const Char_t* sname, const Char_t* label, Double_t weight = 1);
425  void FillTree(const Char_t* sname = "");
426 
427  const KVHashList* GetHistoList() const;
428  const KVHashList* GetTreeList() const;
429 
430  TH1* GetHisto(const Char_t* name) const;
431  TTree* GetTree(const Char_t* name) const;
432 
433  virtual void SaveHistos(const Char_t* filename = "", Option_t* option = "recreate", Bool_t onlyfilled = kFALSE);
434 
435  virtual void SetOpt(const Char_t* option, const Char_t* value);
436  virtual Bool_t IsOptGiven(const Char_t* option);
437  virtual TString GetOpt(const Char_t* option) const;
438  virtual void UnsetOpt(const Char_t* opt);
439 
441  {
444  }
445 
446  void SetJobOutputFileName(const TString& filename)
447  {
454 
455 #ifdef WITH_CPP11
456  if (KVDataAnalyser::IsRunningBatchAnalysis() && (gDataAnalyser->GetProofMode() == KVDataAnalyser::EProofMode::None))
457 #else
459 #endif
461 
462  else
463  SetCombinedOutputFile(filename);
464  }
465 #ifdef USING_ROOT6
466  void SetTriggerConditionsForRun(int);
467 #endif
468  ClassDef(KVEventSelector, 0)//General purpose analysis class for TTrees containing KVEvent objects
469 };
470 
475 #define AddVar(var,type) Branch(dadastr(var), &var, didixstr(duduvartype(var,type)))
476 #define AddVarBranch(var,branch,type) Branch(dadastr(branch), &var, didixstr(duduvartype(branch,type)))
477 #define duduvartype(var,type) var/type
478 #define didixstr(s) dadastr(s)
479 #define dadastr(s) #s
480 
481 #endif
int Int_t
KVDataAnalyser * gDataAnalyser
#define e(i)
char Char_t
const Bool_t kFALSE
bool Bool_t
double Double_t
const Bool_t kTRUE
const char Option_t
#define ClassDef(name, id)
#define N
char * Form(const char *fmt,...)
TVector3 GetCMVelocity() const
Return vector velocity of centre of mass of reaction (units: cm/ns)
Definition: KV2Body.cpp:573
UInt_t GetNumber() const
Definition: KVBase.h:213
virtual const Char_t * GetJobName() const
Description of an experimental run in database ,.
Definition: KVDBRun.h:35
virtual const KVBatchSystem * GetBatchSystem()
virtual const KV2Body * GetKinematics() const
EProofMode GetProofMode() const
static Bool_t IsRunningBatchAnalysis()
General purpose analysis class for TTree containing KVEvent objects.
virtual void SaveHistos(const Char_t *filename="", Option_t *option="recreate", Bool_t onlyfilled=kFALSE)
virtual void SetInputList(TList *input)
void FillTProfile2D(TProfile2D *h2, Double_t one, Double_t two, Double_t three, Double_t four)
virtual Int_t Version() const
KVHashList * ltree
!
void AddGV(KVVarGlob *vg)
const KVHashList * GetTreeList() const
return the list of created trees
KVVarGlob * GetGV(const Char_t *name) const
Bool_t fDisableCreateTreeFile
used with PROOF
KVEventSelector(TTree *=0)
Bool_t fFirstEvent
set to kFALSE after first event is read
Int_t GetFriendTreeEntry(Long64_t entry, Int_t getall=0)
TH1 * GetHisto(const Char_t *name) const
TProofOutputFile * mergeFile
for merging with PROOF
void SetParticleConditions(const KVParticleCondition &, const KVString &="")
const Char_t * GetBranchName() const
void FillTree(const Char_t *sname="")
virtual void InitAnalysis()
void InitFriendTree(TTree *tree, const TString &branchname)
KVGVList gvlist
List of global variables.
virtual void SetOpt(const Char_t *option, const Char_t *value)
Set a value for an option.
virtual void EndRun()
TTree * fAuxChain
[optional] pointer to another TTree or TChain which may be used during analysis
const KVHashList * GetHistoList() const
return the list of created trees
void SetCombinedOutputFile(const TString &filename)
virtual Bool_t Analysis()
virtual Int_t GetEntry(Long64_t entry, Int_t getall=0)
KVEvent * Event
Declaration of leaf types.
virtual void SetAnalysisFrame()
TTree * GetTree(const Char_t *name) const
return the tree named tree_name
Long64_t fEventsRead
number of events read
virtual void CheckEndOfRun()
Testing whether EndRun() should be called.
void FillTH2(TH2 *h2, Double_t one, Double_t two, Double_t three)
virtual void Begin(TTree *tree)
Long64_t fTreeEntry
current tree entry number
TTree * fChain
pointer to the analyzed TTree or TChain
const KVGVList * GetGVList(void) const
Bool_t CreateTreeFile(const Char_t *filename="")
KVString fCombinedOutputFile
optional name for single results file with trees and histos
void AddTree(TTree *tree)
KVString fBranchName
name of branch which contains events to analyse
virtual void SlaveBegin(TTree *tree)
virtual void SlaveTerminate()
virtual void Terminate()
virtual Bool_t Notify()
void SetEventsReadInterval(Long64_t N)
void SetJobOutputFileName(const TString &filename)
void add_histo(TH1 *histo)
void AddHisto(TH1 *histo)
void FillHisto(const Char_t *sname, Double_t one, Double_t two=1, Double_t three=1, Double_t four=1)
Bool_t AtEndOfRun(void)
TBranch * b_Event
List of branches.
void SetBranchName(const Char_t *n)
virtual void SetAdditionalBranchAddress()
virtual Bool_t IsOptGiven(const Char_t *option)
Returns kTRUE if the option 'opt' has been set.
void FillTProfile(TProfile *h1, Double_t one, Double_t two, Double_t three)
virtual void SetCurrentRun(KVDBRun *)
KVEvent * GetFriendEvent() const
KVParticleCondition fPartCond
(optional) conditions for selecting particles
virtual ~KVEventSelector()
KVEvent * GetEvent() const
void SetTriggerConditionsForRun(int)
void FillTH1(TH1 *h1, Double_t one, Double_t two)
virtual void SetObject(TObject *obj)
virtual void Init(TTree *tree)
void SetParticleConditionsParticleClassName(const KVString &t)
void add_tree(TTree *tree)
virtual void RecalculateGlobalVariables()
virtual Bool_t Process(Long64_t entry)
Bool_t fNotifyCalled
avoid multiple calls to Notify/InitRun
virtual TString GetOpt(const Char_t *option) const
KVHashList * lhisto
!
KVGVList * GetGVList(void)
virtual void EndAnalysis()
virtual void ParseOptions()
void SetEvent(KVEvent *e)
KVEvent * AuxEvent
[optional] events in fAuxChain
Long64_t fEventsReadInterval
interval at which to print number of events read
KVNameValueList fOptionList
parsed list of options given to TTree::Process
HistoType * AddHisto(Args &&... args)
virtual void InitRun()
Int_t GetEventNumber() const
void FillTH3(TH3 *h3, Double_t one, Double_t two, Double_t three, Double_t four)
virtual TList * GetOutputList() const
KVString fPartName
(optional) classname for upcasting in KVParticleCondition::Optimize
virtual void UnsetOpt(const Char_t *opt)
Removes the option 'opt' from the internal lists, as if it had never been set.
Base class container for multi-particle events.
Definition: KVEvent.h:176
void SetFrame(const Char_t *frame, const KVFrameTransform &ft)
Definition: KVEvent.cpp:762
Manage a list of global variables.
Definition: KVGVList.h:121
virtual void Add(TObject *obj)
Definition: KVGVList.cpp:254
KVVarGlob * GetGV(const Char_t *nom) const
Return pointer to global variable in list with name 'nom'.
Definition: KVGVList.cpp:237
Extended version of ROOT THashList.
Definition: KVHashList.h:28
Handles lists of named parameters with different types, a list of KVNamedParameter objects.
Handles particle selection criteria for data analysis classes ,.
virtual void Clear(Option_t *option="")
virtual void Add(TObject *obj)
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition: KVString.h:72
Base class for all global variable implementations.
Definition: KVVarGlob.h:217
virtual Int_t GetNbinsX() const
virtual const char * GetName() const
virtual const char * GetTitle() const
void AbstractMethod(const char *method) const
virtual const char * ClassName() const
virtual void Error(const char *method, const char *msgfmt,...) const
TList * fInput
TSelectorList * fOutput
TObject * fObject
const char * Data() const
virtual Int_t GetEntry(Long64_t entry, Int_t getall=0)
virtual TTree * GetTree() const
virtual Long64_t GetEntries() const
long long Long64_t
const Int_t n
TH1 * h