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