KaliVeda
Toolkit for HIC analysis
KVEvent.h
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 
67 class KVEvent: public KVBase {
68 
69 protected:
70 
73 
74 #ifdef __WITHOUT_TCA_CONSTRUCTED_AT
75  TObject* ConstructedAt(Int_t idx)
76  {
89 
90  TObject* obj = (*fParticles)[idx];
91  if (obj && obj->TestBit(TObject::kNotDeleted)) {
92  return obj;
93  }
94  return (fParticles->GetClass()) ? static_cast<TObject*>(fParticles->GetClass()->New(obj)) : 0;
95  }
97  TObject* ConstructedAt(Int_t idx, Option_t* clear_options)
98  {
110 
111  TObject* obj = (*fParticles)[idx];
112  if (obj && obj->TestBit(TObject::kNotDeleted)) {
113  obj->Clear(clear_options);
114  return obj;
115  }
116  return (fParticles->GetClass()) ? static_cast<TObject*>(fParticles->GetClass()->New(obj)) : 0;
117  }
118 #endif
119 
120 public:
121  KVEvent(const TClass* particle_class, Int_t mult = 50)
122  : fParticles(new TClonesArray(particle_class, mult)),
123  fParameters("EventParameters", "Parameters associated with an event")
124  {
125  CustomStreamer();
126  }
127  virtual ~KVEvent()
128  {
131 
132  fParticles->Delete();
134  }
135 
136  void Copy(TObject& obj) const;
137 
139  {
140  return fParticles;
141  }
142  virtual Int_t GetMult(Option_t* opt = "") const
143  {
152 
153  (void)opt;
154  return fParticles->GetEntriesFast();
155  }
156  virtual KVParticle* GetNextParticle(Option_t* = "") const = 0;
157  virtual void ResetGetNextParticle() const = 0;
158  virtual KVParticle* GetParticle(Int_t npart) const = 0;
159  virtual KVParticle* AddParticle() = 0;
160 
161  KVNucleus* GetNextNucleus(Option_t* opt = "") const;
162  void ResetGetNextNucleus() const
163  {
166  }
167  KVNucleus* GetNucleus(Int_t npart) const;
169 
170  virtual void SetFrame(const Char_t*, const KVFrameTransform&) = 0;
171  virtual void SetFrame(const Char_t*, const Char_t*, const KVFrameTransform&) = 0;
172  virtual void MergeEventFragments(TCollection* events, Option_t* opt = "");
173 
175  {
177  }
178 
180  {
181  return (KVNameValueList*)&fParameters;
182  }
183 
184  virtual Double_t GetSum(const Char_t*, Option_t* = "") = 0;
185  virtual Double_t GetChannelQValue() const = 0;
186  virtual KVString GetPartitionName() = 0;
187  virtual void SetFrameName(const KVString& name) = 0;
188  virtual void ChangeDefaultFrame(const Char_t*, const Char_t* = "") = 0;
189  const Char_t* GetFrameName() const
190  {
193 
194  return (GetParameters()->HasStringParameter("defaultFrame") ?
195  GetParameters()->GetStringValue("defaultFrame") : "");
196  }
197  template<typename ValType> void SetParameter(const Char_t* name, ValType value) const
198  {
202 
204  }
205  virtual void GetMasses(std::vector<Double_t>&) = 0;
206 
208 #define KVEVENT_MAKE_EVENT_BRANCH_NO_VOID_PTR 1
209  template<typename T>
210  static void MakeEventBranch(TTree* tree, const TString& branchname, T& event, Int_t bufsize = 10000000)
211  {
225 
226  tree->Branch(branchname, event->ClassName(), &event, bufsize, 0)->SetAutoDelete(kFALSE);
227  }
228  static KVEvent* Factory(const char* plugin)
229  {
231 
232  TPluginHandler* ph = LoadPlugin("KVEvent", plugin);
233  if (ph) {
234  return (KVEvent*)ph->ExecPlugin(0);
235  }
236  return nullptr;
237  }
238  void Clear(Option_t* opt = "")
239  {
244 
245  if (strcmp(opt, "")) { // pass options to particle class Clear() method
246  TString Opt = Form("C+%s", opt);
247  fParticles->Clear(Opt);
248  }
249  else
250  fParticles->Clear("C");
251  fParameters.Clear();
253  }
254 
255  ClassDef(KVEvent, 4)
256 };
257 
258 #endif // KVBASEEVENT_H
int Int_t
#define SafeDelete(p)
char Char_t
constexpr Bool_t kFALSE
double Double_t
const char Option_t
#define ClassDef(name, id)
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
char name[80]
char * Form(const char *fmt,...)
Base class for KaliVeda framework.
Definition: KVBase.h:142
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:67
KVNameValueList * GetParameters() const
Definition: KVEvent.h:179
virtual Double_t GetSum(const Char_t *, Option_t *="")=0
virtual KVParticle * GetParticle(Int_t npart) const =0
void Copy(TObject &obj) const
Definition: KVEvent.cpp:123
KVNameValueList fParameters
general-purpose list of parameters
Definition: KVEvent.h:72
KVNucleus * GetNucleus(Int_t npart) const
Definition: KVEvent.cpp:92
void Clear(Option_t *opt="")
Definition: KVEvent.h:238
virtual Double_t GetChannelQValue() const =0
static KVEvent * Factory(const char *plugin)
Definition: KVEvent.h:228
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:210
virtual void SetFrame(const Char_t *, const Char_t *, const KVFrameTransform &)=0
const Char_t * GetFrameName() const
Definition: KVEvent.h:189
virtual KVString GetPartitionName()=0
virtual ~KVEvent()
Definition: KVEvent.h:127
TClonesArray * fParticles
array of particles in event
Definition: KVEvent.h:71
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:138
KVEvent(const TClass *particle_class, Int_t mult=50)
Definition: KVEvent.h:121
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:142
void SetParameter(const Char_t *name, ValType value) const
Definition: KVEvent.h:197
void ResetGetNextNucleus() const
Definition: KVEvent.h:162
void CustomStreamer()
Definition: KVEvent.h:174
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:69
Handles lists of named parameters with different types, a list of KVNamedParameter objects.
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:126
Base class for relativistic kinematics of massive particles.
Definition: KVParticle.h:396
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition: KVString.h:73
void * New(ENewType defConstructor=kClassNew, Bool_t quiet=kFALSE) const
void BypassStreamer(Bool_t bypass=kTRUE)
TClass * GetClass() const
void Clear(Option_t *option="") override
void Delete(Option_t *option="") override
Int_t GetEntriesFast() const
virtual void Clear(Option_t *="")
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Longptr_t ExecPlugin(int nargs)
void(off) SmallVectorTemplateBase< T