KaliVeda  1.13/01
Heavy-Ion Analysis Toolkit
KVGVList.h
Go to the documentation of this file.
1 #ifndef KVGVList_h
2 #define KVGVList_h
3 #include "KVVarGlob.h"
4 #include "KVUniqueNameList.h"
5 #include "KVParticleCondition.h"
6 #include "TTree.h"
7 #include "KVEventClassifier.h"
8 
122 class KVGVList: public KVUniqueNameList {
123 
125  std::vector<Double_t> fBranchVar;
126  std::vector<Int_t> fIBranchVar;
130 
133  {
134  TString _s(s);
135  _s.ReplaceAll("+", "_");
136  _s.ReplaceAll("*", "_");
137  _s.ReplaceAll("-", "_");
138  _s.ReplaceAll("/", "_");
139  _s.ReplaceAll("=", "_");
140  return _s;
141  }
142  KVVarGlob* prepareGVforAdding(const Char_t* class_name, const Char_t* name);
143 
144 protected:
145  void init_KVGVList(void);
149  void Fill(const KVNucleus* c);
150  void Fill2(const KVNucleus* c1, const KVNucleus* c2);
151  void FillN(const KVEvent* e);
152  void Calculate();
153  void Calculate2();
154  void CalculateN();
155 
156 public:
157  KVGVList(const KVString& name = "default", const KVParticleCondition& selection = KVParticleCondition());
158  KVGVList(const KVGVList& a);
159 
160  KVVarGlob* AddGV(const Char_t* class_name, const Char_t* name);
161  KVVarGlob* AddGVFirst(const Char_t* class_name, const Char_t* name);
162 
163  void Init(void);
164  void Reset(void);
165 
167 
168  KVVarGlob* GetGV(const Char_t* nom) const; //find global variable with name 'nom'
169 
171  KVVarGlob* GetGVType(const Char_t* class_name)
172  {
173  return (KVVarGlob*)FindObjectByClass(class_name);
174  }
175  virtual void Add(TObject* obj) ;
176  virtual void AddFirst(TObject* obj) ;
177 
180  {
181  return (fVG1.GetEntries() > 0);
182  }
183 
186  {
187  return (fVG2.GetEntries() > 0);
188  }
189 
192  {
193  return (fVGN.GetEntries() > 0);
194  }
195 
196  void MakeBranches(TTree*);
197  void FillBranches();
198 
199  bool AbortEventAnalysis() const
200  {
201  return fAbortEventAnalysis;
202  }
203  KVEventClassifier* AddEventClassifier(const TString& varname, const TString& value = "");
204 
205  ClassDef(KVGVList, 3) // List of global variables
206 };
207 #endif
int Int_t
KVTemplateParticleCondition< KVNucleus > KVParticleCondition
char Char_t
bool Bool_t
#define ClassDef(name, id)
Simple class for sorting events according to global variables.
Abstract base class container for multi-particle events.
Definition: KVEvent.h:66
Manage a list of global variables.
Definition: KVGVList.h:122
Bool_t Has1BodyVariables()
returns kTRUE if list contains 1-body variables
Definition: KVGVList.h:179
void FillBranches()
Definition: KVGVList.cpp:459
KVParticleCondition fSelection
used to select particles to iterate over in event
Definition: KVGVList.h:124
KVVarGlob * AddGVFirst(const Char_t *class_name, const Char_t *name)
Definition: KVGVList.cpp:734
Int_t fNbIBranch
Definition: KVGVList.h:128
void Reset(void)
Reset all variables before treating an event.
Definition: KVGVList.cpp:88
void Init(void)
Definition: KVGVList.cpp:67
std::vector< Double_t > fBranchVar
used for automatic creation & filling of TTree branches
Definition: KVGVList.h:125
Bool_t HasNBodyVariables()
returns kTRUE if list contains N-body variables
Definition: KVGVList.h:191
TList fVG2
two-body variables
Definition: KVGVList.h:147
TString NameSanitizer(const Char_t *s) const
Definition: KVGVList.h:132
Int_t fNbBranch
Definition: KVGVList.h:127
virtual void AddFirst(TObject *obj)
Definition: KVGVList.cpp:316
void Calculate2()
Calculate all 2-body observables after filling.
Definition: KVGVList.cpp:163
void CalculateN()
Calculate all N-body observables after filling.
Definition: KVGVList.cpp:182
void CalculateGlobalVariables(KVEvent *e)
Definition: KVGVList.cpp:211
bool fAbortEventAnalysis
set to false if a global variable fails its own event selection criterion
Definition: KVGVList.h:129
void FillN(const KVEvent *e)
Calls KVVarGlob::FillN(KVEvent*) method of all N-body variables in the list.
Definition: KVGVList.cpp:133
std::vector< Int_t > fIBranchVar
used for automatic creation & filling of TTree branches
Definition: KVGVList.h:126
KVEventClassifier * AddEventClassifier(const TString &varname, const TString &value="")
Definition: KVGVList.cpp:555
TList fVGN
N-body variables.
Definition: KVGVList.h:148
void Calculate()
Calculate all 1-body observables after filling.
Definition: KVGVList.cpp:144
KVVarGlob * prepareGVforAdding(const Char_t *class_name, const Char_t *name)
Definition: KVGVList.cpp:690
TList fVG1
one-body variables
Definition: KVGVList.h:146
KVVarGlob * AddGV(const Char_t *class_name, const Char_t *name)
Definition: KVGVList.cpp:648
virtual void Add(TObject *obj)
Definition: KVGVList.cpp:279
void Fill(const KVNucleus *c)
Definition: KVGVList.cpp:102
void init_KVGVList(void)
Definition: KVGVList.cpp:15
KVVarGlob * GetGVType(const Char_t *class_name)
Returns first global variable in list with given class.
Definition: KVGVList.h:171
void MakeBranches(TTree *)
Definition: KVGVList.cpp:361
KVVarGlob * GetGV(const Char_t *nom) const
Return pointer to global variable in list with name 'nom'.
Definition: KVGVList.cpp:262
Bool_t Has2BodyVariables()
returns kTRUE if list contains 2-body variables
Definition: KVGVList.h:185
bool AbortEventAnalysis() const
Definition: KVGVList.h:199
KVGVList(const KVString &name="default", const KVParticleCondition &selection=KVParticleCondition())
Definition: KVGVList.cpp:32
void Fill2(const KVNucleus *c1, const KVNucleus *c2)
Definition: KVGVList.cpp:119
Description of properties and kinematics of atomic nuclei.
Definition: KVNucleus.h:125
TObject * FindObjectByClass(const Char_t *) const
Return (first) object in embedded list with given class.
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition: KVString.h:72
Optimised list in which named objects can only be placed once.
Base class for all global variable implementations.
Definition: KVVarGlob.h:231
virtual Int_t GetEntries() const
TString & ReplaceAll(const char *s1, const char *s2)
const long double s
Definition: KVUnits.h:94