KaliVeda
Toolkit for HIC analysis
KVEvent.cpp
1 #include "KVEvent.h"
2 #include "KVNucleus.h"
3 
5 
6 
7 
12 void KVEvent::Streamer(TBuffer& R__b)
13 {
14  // Customised Streamer for KVEvent.
15  // This is just the automatic Streamer with the addition of a call to the Clear()
16  // method before reading a new object (avoid memory leaks with lists of parameters).
17 
18  if (R__b.IsReading()) {
19  Clear();
20  R__b.ReadClassBuffer(KVEvent::Class(), this);
21  }
22  else {
23  R__b.WriteClassBuffer(KVEvent::Class(), this);
24  }
25 }
26 
27 
28 
54 
56 {
57  // \warning Only use with events containing objects derived from KVNucleus
58  //
59  // Use this method to iterate over the list of nuclei in the event
60  // After the last particle GetNextNucleus() returns a null pointer and
61  // resets itself ready for a new iteration over the particle list.
62  //
63  // If opt="" all particles are included in the iteration.
64  // If opt="ok" or "OK" only nuclei whose KVNucleus::IsOK() method returns kTRUE are included.
65  //
66  // Any other value of opt is interpreted as a (case-insensitive) particle group name: only
67  // particles with BelongsToGroup(opt) returning kTRUE are included.
68  //
69  // If you want to start from the beginning again before getting to the end
70  // of the list, especially if you want to change the selection criteria,
71  // call method ResetGetNextNucleus() before continuing.
72  //
73  // If you interrupt an iteration before the end, then start another iteration
74  // without calling ResetGetNextNucleus(), even if you change the argument of
75  // the call to GetNextNucleus(), you will repeat exactly the same iteration
76  // as the previous one.
77  //
78  // \warning Only one iteration at a time over the event can be performed
79  // using this method. If you want/need to perform several i.e. nested
80  // iterations, see KVTemplateEvent::EventIterator
81 
82  return dynamic_cast<KVNucleus*>(GetNextParticle(opt));
83 }
84 
85 
86 
92 
94 {
95  // \warning Only use with events containing objects derived from KVNucleus
96  //
97  // \param[in] npart index of particle in event, which is a non-zero value from 1 to the value returned by GetMult()
98  // \returns pointer to the particle with index npart
99 
100  return dynamic_cast<KVNucleus*>(GetParticle(npart));
101 }
102 
103 
104 
109 
111 {
112  // Add a particle to the event
113  //
114  // \returns pointer to new particle if it inherits from KVNucleus, nullptr if not
115  return dynamic_cast<KVNucleus*>(AddParticle());
116 }
117 
118 
119 
123 
124 void KVEvent::Copy(TObject& obj) const
125 {
126  // Copy this event into the object referenced by obj,
127  // assumed to be at least derived from KVEvent.
128  KVBase::Copy(obj);
130  Int_t MTOT = fParticles->GetEntriesFast();
131  for (Int_t nn = 0; nn < MTOT; nn += 1) {
132  GetParticle(nn + 1)->Copy(*((KVEvent&) obj).AddParticle());
133  }
134 }
135 
136 
137 
154 
156 {
157  // Merge all events in the list into one event (this one)
158  //
159  // We also merge/sum the parameter lists of the events
160  //
161  // First we clear this event, then we fill it with copies of each particle in each event
162  // in the list.
163  //
164  // If option "opt" is given, it is given as argument to each call to
165  // KVEvent::Clear() - this option is then passed on to the KVParticle::Clear()
166  // method of each particle in each event.
167  //
168  // \param[in] events A list of events to merge
169  // \param[in] opt Optional argument transmitted to KVEvent::Clear()
170  //
171  // \note the events in the list will be empty and useless after this!
172 
173  Clear(opt);
174  TIter it(events);
175  KVEvent* e;
176  while ((e = (KVEvent*)it())) {
177  e->ResetGetNextParticle();
178  while (auto n = e->GetNextParticle()) {
179  n->Copy(*AddParticle());
180  }
181  GetParameters()->Merge(*(e->GetParameters()));
182  e->Clear(opt);
183  }
184 }
185 
186 
187 
int Int_t
#define e(i)
const char Option_t
void Copy(TObject &) const override
Make a copy of this object.
Definition: KVBase.cpp:373
Abstract base class container for multi-particle events.
Definition: KVEvent.h:67
KVNameValueList * GetParameters() const
Definition: KVEvent.h:179
virtual KVParticle * GetParticle(Int_t npart) const =0
KVNameValueList fParameters
general-purpose list of parameters
Definition: KVEvent.h:72
KVNucleus * GetNucleus(Int_t npart) const
Definition: KVEvent.cpp:93
KVNucleus * AddNucleus()
Definition: KVEvent.cpp:110
TClonesArray * fParticles
array of particles in event
Definition: KVEvent.h:71
virtual KVParticle * AddParticle()=0
virtual KVParticle * GetNextParticle(Option_t *="") const =0
KVNucleus * GetNextNucleus(Option_t *opt="") const
Definition: KVEvent.cpp:55
void Clear(Option_t *opt="") override
Definition: KVEvent.h:238
void Copy(TObject &obj) const override
Definition: KVEvent.cpp:124
virtual void MergeEventFragments(TCollection *events, Option_t *opt="")
Definition: KVEvent.cpp:155
void Copy(TObject &nvl) const override
void Merge(const KVNameValueList &)
Description of properties and kinematics of atomic nuclei.
Definition: KVNucleus.h:123
void Copy(TObject &) const override
Definition: KVParticle.cpp:286
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=nullptr)=0
Bool_t IsReading() const
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
void Streamer(TBuffer &) override
static TClass * Class()
Int_t GetEntriesFast() const
const Int_t n
ClassImp(TPyArg)