KaliVeda  1.13/01
Heavy-Ion Analysis Toolkit
KVEvent.h
Go to the documentation of this file.
1 #ifndef KVBASEEVENT_H
2 #define KVBASEEVENT_H
3 
4 #include <KVBase.h>
5 #include <TTree.h>
6 #include "KVNameValueList.h"
7 #include "TClonesArray.h"
8 class KVFrameTransform;
9 class KVParticle;
10 class KVNucleus;
11 
12 #include <TH1.h>
13 #include <iterator>
14 
15 class KVIntegerList;
16 
66 class KVEvent: public KVBase {
67 
68 protected:
69 
72 
73 #ifdef __WITHOUT_TCA_CONSTRUCTED_AT
74  TObject* ConstructedAt(Int_t idx)
75  {
88 
89  TObject* obj = (*fParticles)[idx];
90  if (obj && obj->TestBit(TObject::kNotDeleted)) {
91  return obj;
92  }
93  return (fParticles->GetClass()) ? static_cast<TObject*>(fParticles->GetClass()->New(obj)) : 0;
94  }
96  TObject* ConstructedAt(Int_t idx, Option_t* clear_options)
97  {
109 
110  TObject* obj = (*fParticles)[idx];
111  if (obj && obj->TestBit(TObject::kNotDeleted)) {
112  obj->Clear(clear_options);
113  return obj;
114  }
115  return (fParticles->GetClass()) ? static_cast<TObject*>(fParticles->GetClass()->New(obj)) : 0;
116  }
117 #endif
118 
119 public:
120  KVEvent(const TClass* particle_class, Int_t mult = 50)
121  : fParticles(new TClonesArray(particle_class, mult)),
122  fParameters("EventParameters", "Parameters associated with an event")
123  {
124  CustomStreamer();
125  }
126  virtual ~KVEvent()
127  {
130 
131  fParticles->Delete();
133  }
134 
135  void Copy(TObject& obj) const;
136 
138  {
139  return fParticles;
140  }
141  virtual Int_t GetMult(Option_t* opt = "") const
142  {
151 
152  (void)opt;
153  return fParticles->GetEntriesFast();
154  }
155  virtual KVParticle* GetNextParticle(Option_t* = "") const = 0;
156  virtual void ResetGetNextParticle() const = 0;
157  virtual KVParticle* GetParticle(Int_t npart) const = 0;
158  virtual KVParticle* AddParticle() = 0;
159 
160  KVNucleus* GetNextNucleus(Option_t* opt = "") const;
161  void ResetGetNextNucleus() const
162  {
165  }
166  KVNucleus* GetNucleus(Int_t npart) const;
168 
169  virtual void SetFrame(const Char_t*, const KVFrameTransform&) = 0;
170  virtual void SetFrame(const Char_t*, const Char_t*, const KVFrameTransform&) = 0;
171  virtual Bool_t IsOK() const
172  {
178 
179  return (GetMult("ok") >= GetMinimumOKMultiplicity());
180  }
182  {
185  SetParameter("MIN_OK_MULT", x);
186  }
188  {
192  Int_t x = GetParameters()->GetIntValue("MIN_OK_MULT");
193  if (x == -1) return 1;
194  return x;
195  }
196  virtual void MergeEventFragments(TCollection* events, Option_t* opt = "");
197 
199  {
201  }
202 
204  {
205  return (KVNameValueList*)&fParameters;
206  }
207 
208  virtual Double_t GetSum(const Char_t*, Option_t* = "") = 0;
209  virtual Double_t GetChannelQValue() const = 0;
210  virtual KVString GetPartitionName() = 0;
211  virtual void SetFrameName(const KVString& name) = 0;
212  virtual void ChangeDefaultFrame(const Char_t*, const Char_t* = "") = 0;
213  const Char_t* GetFrameName() const
214  {
217 
218  return (GetParameters()->HasStringParameter("defaultFrame") ?
219  GetParameters()->GetStringValue("defaultFrame") : "");
220  }
221  template<typename ValType> void SetParameter(const Char_t* name, ValType value) const
222  {
226 
227  GetParameters()->SetValue(name, value);
228  }
229  virtual void GetMasses(std::vector<Double_t>&) = 0;
230 
232 #define KVEVENT_MAKE_EVENT_BRANCH_NO_VOID_PTR 1
233  template<typename T>
234  static void MakeEventBranch(TTree* tree, const TString& branchname, T& event, Int_t bufsize = 10000000)
235  {
249 
250  tree->Branch(branchname, event->ClassName(), &event, bufsize, 0)->SetAutoDelete(kFALSE);
251  }
252  static KVEvent* Factory(const char* plugin)
253  {
255 
256  TPluginHandler* ph = LoadPlugin("KVEvent", plugin);
257  if (ph) {
258  return (KVEvent*)ph->ExecPlugin(0);
259  }
260  return nullptr;
261  }
262  void Clear(Option_t* opt = "")
263  {
268 
269  if (strcmp(opt, "")) { // pass options to particle class Clear() method
270  TString Opt = Form("C+%s", opt);
271  fParticles->Clear(Opt);
272  }
273  else
274  fParticles->Clear("C");
275  fParameters.Clear();
277  }
278 
279  ClassDef(KVEvent, 4)
280 };
281 
282 #endif // KVBASEEVENT_H
int Int_t
typedef void(GLAPIENTRYP _GLUfuncptr)(void)
#define SafeDelete(p)
char Char_t
const Bool_t kFALSE
bool Bool_t
double Double_t
const char Option_t
#define ClassDef(name, id)
char * Form(const char *fmt,...)
Base class for KaliVeda framework.
Definition: KVBase.h:141
static TPluginHandler * LoadPlugin(const Char_t *base, const Char_t *uri="0")
Definition: KVBase.cpp:793
Abstract base class container for multi-particle events.
Definition: KVEvent.h:66
KVNameValueList * GetParameters() const
Definition: KVEvent.h:203
virtual Double_t GetSum(const Char_t *, Option_t *="")=0
void SetMinimumOKMultiplicity(Int_t x)
Definition: KVEvent.h:181
virtual KVParticle * GetParticle(Int_t npart) const =0
void Copy(TObject &obj) const
Definition: KVEvent.cpp:123
virtual Bool_t IsOK() const
Definition: KVEvent.h:171
KVNameValueList fParameters
general-purpose list of parameters
Definition: KVEvent.h:71
KVNucleus * GetNucleus(Int_t npart) const
Definition: KVEvent.cpp:92
void Clear(Option_t *opt="")
Clear object properties : name, type/title, number, label.
Definition: KVEvent.h:262
virtual Double_t GetChannelQValue() const =0
static KVEvent * Factory(const char *plugin)
Definition: KVEvent.h:252
virtual void SetFrameName(const KVString &name)=0
KVNucleus * AddNucleus()
Definition: KVEvent.cpp:109
virtual void ResetGetNextParticle() const =0
static void MakeEventBranch(TTree *tree, const TString &branchname, T &event, Int_t bufsize=10000000)
Definition: KVEvent.h:234
virtual void SetFrame(const Char_t *, const Char_t *, const KVFrameTransform &)=0
const Char_t * GetFrameName() const
Definition: KVEvent.h:213
Int_t GetMinimumOKMultiplicity() const
Definition: KVEvent.h:187
virtual KVString GetPartitionName()=0
virtual ~KVEvent()
Definition: KVEvent.h:126
TClonesArray * fParticles
array of particles in event
Definition: KVEvent.h:70
virtual KVParticle * AddParticle()=0
virtual KVParticle * GetNextParticle(Option_t *="") const =0
virtual void SetFrame(const Char_t *, const KVFrameTransform &)=0
const TClonesArray * GetParticleArray() const
Definition: KVEvent.h:137
KVEvent(const TClass *particle_class, Int_t mult=50)
Definition: KVEvent.h:120
virtual void GetMasses(std::vector< Double_t > &)=0
KVNucleus * GetNextNucleus(Option_t *opt="") const
Definition: KVEvent.cpp:54
virtual void ChangeDefaultFrame(const Char_t *, const Char_t *="")=0
virtual Int_t GetMult(Option_t *opt="") const
Definition: KVEvent.h:141
void SetParameter(const Char_t *name, ValType value) const
Definition: KVEvent.h:221
void ResetGetNextNucleus() const
Definition: KVEvent.h:161
void CustomStreamer()
Definition: KVEvent.h:198
virtual void MergeEventFragments(TCollection *events, Option_t *opt="")
Definition: KVEvent.cpp:154
Utility class for kinematical transformations of KVParticle class.
Handle a list of positive integers (partition)
Definition: KVIntegerList.h:68
Handles lists of named parameters with different types, a list of KVNamedParameter objects.
Int_t GetIntValue(const Char_t *name) const
void SetValue(const Char_t *name, value_type value)
virtual void Clear(Option_t *opt="")
Description of properties and kinematics of atomic nuclei.
Definition: KVNucleus.h:125
Base class for relativistic kinematics of massive particles.
Definition: KVParticle.h:398
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition: KVString.h:72
void * New(ENewType defConstructor=kClassNew, Bool_t quiet=kFALSE) const
virtual void Delete(Option_t *option="")
void BypassStreamer(Bool_t bypass=kTRUE)
virtual void Clear(Option_t *option="")
TClass * GetClass() const
Int_t GetEntriesFast() const
virtual void Clear(Option_t *="")
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Longptr_t ExecPlugin(int nargs, const T &... params)
Double_t x[n]