KaliVeda  1.12/06
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 "TTree.h"
6 #include "KVEventClassifier.h"
7 
121 class KVGVList: public KVUniqueNameList {
122 
123  std::vector<Double_t> fBranchVar;
124  std::vector<Int_t> fIBranchVar;
128 
131  {
132  TString _s(s);
133  _s.ReplaceAll("+", "_");
134  _s.ReplaceAll("*", "_");
135  _s.ReplaceAll("-", "_");
136  _s.ReplaceAll("/", "_");
137  _s.ReplaceAll("=", "_");
138  return _s;
139  }
140  KVVarGlob* prepareGVforAdding(const Char_t* class_name, const Char_t* name);
141 
142 protected:
143  void init_KVGVList(void);
147  void Fill(const KVNucleus* c);
148  void Fill2(const KVNucleus* c1, const KVNucleus* c2);
149  void FillN(const KVEvent* e);
150  void Calculate();
151  void Calculate2();
152  void CalculateN();
153 
154 public:
155  KVGVList(void); // constructeur par defaut
156  KVGVList(const KVGVList& a); // constructeur par Copy
157 
158  virtual ~ KVGVList(void) {}
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); // methode d'initialisation des variables globales
164  void Reset(void); // Remise a zero avant le
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  }
204 
205  ClassDef(KVGVList, 3) // List of global variables
206 };
207 #endif
int Int_t
char Char_t
bool Bool_t
#define ClassDef(name, id)
Simple class for sorting events according to global variables.
Base class container for multi-particle events.
Definition: KVEvent.h:176
Manage a list of global variables.
Definition: KVGVList.h:121
Bool_t Has1BodyVariables()
returns kTRUE if list contains 1-body variables
Definition: KVGVList.h:179
void FillBranches()
Definition: KVGVList.cpp:434
KVVarGlob * AddGVFirst(const Char_t *class_name, const Char_t *name)
Definition: KVGVList.cpp:701
Int_t fNbIBranch
Definition: KVGVList.h:126
void Reset(void)
Reset all variables before treating an event.
Definition: KVGVList.cpp:63
void Init(void)
Definition: KVGVList.cpp:49
std::vector< Double_t > fBranchVar
used for automatic creation & filling of TTree branches
Definition: KVGVList.h:123
Bool_t HasNBodyVariables()
returns kTRUE if list contains N-body variables
Definition: KVGVList.h:191
TList fVG2
two-body variables
Definition: KVGVList.h:145
TString NameSanitizer(const Char_t *s) const
Definition: KVGVList.h:130
Int_t fNbBranch
Definition: KVGVList.h:125
virtual void AddFirst(TObject *obj)
Definition: KVGVList.cpp:291
void Calculate2()
Calculate all 2-body observables after filling.
Definition: KVGVList.cpp:140
void CalculateN()
Calculate all N-body observables after filling.
Definition: KVGVList.cpp:161
void CalculateGlobalVariables(KVEvent *e)
Definition: KVGVList.cpp:189
bool fAbortEventAnalysis
set to false if a global variable fails its own event selection criterion
Definition: KVGVList.h:127
void FillN(const KVEvent *e)
Calls KVVarGlob::FillN(KVEvent*) method of all N-body variables in the list.
Definition: KVGVList.cpp:108
std::vector< Int_t > fIBranchVar
used for automatic creation & filling of TTree branches
Definition: KVGVList.h:124
TList fVGN
N-body variables.
Definition: KVGVList.h:146
void Calculate()
Calculate all 1-body observables after filling.
Definition: KVGVList.cpp:119
KVVarGlob * prepareGVforAdding(const Char_t *class_name, const Char_t *name)
Definition: KVGVList.cpp:657
TList fVG1
one-body variables
Definition: KVGVList.h:144
KVVarGlob * AddGV(const Char_t *class_name, const Char_t *name)
Definition: KVGVList.cpp:615
virtual void Add(TObject *obj)
Definition: KVGVList.cpp:254
void Fill(const KVNucleus *c)
Definition: KVGVList.cpp:77
void init_KVGVList(void)
Definition: KVGVList.cpp:16
KVVarGlob * GetGVType(const Char_t *class_name)
Returns first global variable in list with given class.
Definition: KVGVList.h:171
KVEventClassifier * AddEventClassifier(const TString &varname)
Definition: KVGVList.cpp:526
void MakeBranches(TTree *)
Definition: KVGVList.cpp:336
virtual ~ KVGVList(void)
Definition: KVGVList.h:158
KVVarGlob * GetGV(const Char_t *nom) const
Return pointer to global variable in list with name 'nom'.
Definition: KVGVList.cpp:237
Bool_t Has2BodyVariables()
returns kTRUE if list contains 2-body variables
Definition: KVGVList.h:185
KVGVList(void)
Definition: KVGVList.cpp:26
bool AbortEventAnalysis() const
Definition: KVGVList.h:199
void Fill2(const KVNucleus *c1, const KVNucleus *c2)
Definition: KVGVList.cpp:94
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.
Optimised list in which named objects can only be placed once.
Base class for all global variable implementations.
Definition: KVVarGlob.h:217
virtual Int_t GetEntries() const
TString & ReplaceAll(const char *s1, const char *s2)
const long double s
Definition: KVUnits.h:94