KaliVeda  1.12/06
Heavy-Ion Analysis Toolkit
KVEvent.h
Go to the documentation of this file.
1 /***************************************************************************
2  kvevent.h - description
3  -------------------
4  begin : Sun May 19 2002
5  copyright : (C) 2002 by J.D. Frankland
6  email : frankland@ganil.fr
7 
8 $Id: KVEvent.h,v 1.29 2008/12/17 11:23:12 ebonnet Exp $
9  ***************************************************************************/
10 
11 /***************************************************************************
12  * *
13  * This program is free software; you can redistribute it and/or modify *
14  * it under the terms of the GNU General Public License as published by *
15  * the Free Software Foundation; either version 2 of the License, or *
16  * (at your option) any later version. *
17  * *
18  ***************************************************************************/
19 
20 #ifndef KVEVENT_H
21 #define KVEVENT_H
22 
23 #define KVEVENT_PART_INDEX_OOB "Particle index %d out of bounds [1,%d]"
24 
25 #include "TTree.h"
26 #include "TVector3.h"
27 #include "TClonesArray.h"
28 #include "KVNucleus.h"
29 #include "KVBase.h"
30 #include "KVConfig.h"
31 #include "TRotation.h"
32 #include "TLorentzRotation.h"
33 #include "KVParticleCondition.h"
34 #include "KVNameValueList.h"
35 #include "TMethodCall.h"
36 
37 #include <TH1.h>
38 #include <iterator>
39 
40 class KVIntegerList;
41 
176 class KVEvent: public KVBase {
177 
178 protected:
179 
182 #ifdef __WITHOUT_TCA_CONSTRUCTED_AT
183  TObject* ConstructedAt(Int_t idx);
184  TObject* ConstructedAt(Int_t idx, Option_t* clear_options);
185 #endif
186 public:
187 
188  class Iterator {
189  public:
190  typedef std::forward_iterator_tag iterator_category;
192  typedef std::ptrdiff_t difference_type;
193  typedef KVNucleus* pointer;
195 
196  enum Type { // type of iterator
197  Null, // null value
198  All, // include all particles
199  OK, // include particles which are "OK"
200  Group // include particles belonging to group
201  };
202 
203  private:
205  Type fType;//iterator type
206  mutable Bool_t fIterating;//=kTRUE when iteration in progress
207  TString fGroup;//groupname for group iterations
209  {
212 
213  switch (fType) {
214  case OK:
215  return current()->IsOK();
216  break;
217  case Group:
218  return current()->BelongsToGroup(fGroup);
219  break;
220  case All:
221  default:
222  return kTRUE;
223  }
224  return kTRUE;
225  }
227  {
230  return reinterpret_cast<KVNucleus*>(*fIter);
231  }
232  public:
234  : fIter(static_cast<TIterator*>(nullptr)),
235 #ifdef WITH_CPP11
236  fType(Type::Null),
237 #else
238  fType(Null),
239 #endif
241  fGroup()
242  {}
243  Iterator(const Iterator& i)
244  : fIter(i.fIter),
245  fType(i.fType),
247  fGroup(i.fGroup)
248  {}
249 
250 #ifdef WITH_CPP11
251  Iterator(const KVEvent* e, Type t = Type::All, const TString& grp = "")
252 #else
253  Iterator(const KVEvent* e, Type t = All, const TString& grp = "")
254 #endif
255  : fIter(e->fParticles), fType(t), fIterating(kTRUE), fGroup(grp)
256  {
265 
267  fIter.Begin();
268  while ((current() != nullptr) && !AcceptableIteration()) ++fIter;
269  }
270 
271 #ifdef WITH_CPP11
272  Iterator(const KVEvent& e, Type t = Type::All, const TString& grp = "")
273 #else
274  Iterator(const KVEvent& e, Type t = All, const TString& grp = "")
275 #endif
276  : fIter(e.fParticles), fType(t), fIterating(kTRUE), fGroup(grp)
277  {
286 
288  fIter.Begin();
289  while ((current() != nullptr) && !AcceptableIteration()) ++fIter;
290  }
291 
293  {
295 
296  return *(current());
297  }
298  template<typename PointerType>
299  PointerType* get_pointer() const
300  {
301  return dynamic_cast<PointerType*>(current());
302  }
303  template<typename ReferenceType>
304  ReferenceType& get_reference() const
305  {
306  return dynamic_cast<ReferenceType&>(*current());
307  }
308  Bool_t operator!= (const Iterator& it) const
309  {
311  return current() != it.current();
312  }
313  Bool_t operator== (const Iterator& it) const
314  {
316  return current() == it.current();
317  }
319  {
322  if (current() == nullptr) fIterating = kFALSE;
323  ++fIter;
324  while (current() != nullptr && !AcceptableIteration()) ++fIter;
325  return *this;
326  }
328  {
331  Iterator tmp(*this);
332  operator++();
333  return tmp;
334  }
336  {
338  if (this != &rhs) { // check self-assignment based on address of object
339  fIter = rhs.fIter;
340  fType = rhs.fType;
341  fGroup = rhs.fGroup;
342  fIterating = rhs.fIterating;
343  }
344  return *this;
345  }
346 
347  static Iterator End()
348  {
349  return Iterator();
350  }
351  virtual ~Iterator() {}
352 #ifdef WITH_CPP11
353  void Reset(Type t = Type::Null, TString grp = "")
354 #else
355  void Reset(Type t = Null, TString grp = "")
356 #endif
357  {
361 #ifdef WITH_CPP11
362  if (t != Type::Null) {
363 #else
364  if (t != Null) {
365 #endif
366  fType = t;
367  fGroup = grp;
368  }
369  fIter.Begin();
370  fIterating = kTRUE;
371  while ((current() != nullptr) && !AcceptableIteration()) ++fIter;
372  }
374  {
376  return fIterating;
377  }
378  void SetIsIterating(Bool_t on = kTRUE)
379  {
381  fIterating = on;
382  }
383 
384  ClassDef(Iterator, 0) //Iterator class for KVEvent
385  };
386 protected:
387  mutable Iterator fIter;
388 
389 public:
392  {
393  return (KVNameValueList*)&fParameters;
394  }
395 
396  KVEvent(Int_t mult = 50, const char* classname = "KVNucleus");
397  virtual ~ KVEvent();
398 
399 #if ROOT_VERSION_CODE >= ROOT_VERSION(3,4,0)
400  virtual void Copy(TObject& obj) const;
401 #else
402  virtual void Copy(TObject& obj);
403 #endif
404 
406  KVNucleus* GetParticle(Int_t npart) const;
407  virtual Int_t GetMult(Option_t* opt = "") const;
408  Int_t GetMultiplicity(Int_t Z, Int_t A = 0, Option_t* opt = "");
409  void GetMultiplicities(Int_t mult[], const TString& species, Option_t* opt = "");
410  Double_t GetSum(const Char_t* KVNucleus_method, Option_t* opt = "");
411  Double_t GetSum(const Char_t* KVNucleus_method, const Char_t* method_prototype, const Char_t* args, Option_t* opt = "");
412  void FillHisto(TH1* h, const Char_t* KVNucleus_method, Option_t* opt = "");
413  void FillHisto(TH1* h, const Char_t* KVNucleus_method, const Char_t* method_prototype, const Char_t* args, Option_t* opt = "");
414  virtual void Clear(Option_t* opt = "");
415  virtual void Print(Option_t* t = "") const;
416  virtual void ls(Option_t* t = "") const
417  {
418  Print(t);
419  }
420  KVNucleus* GetParticleWithName(const Char_t* name) const;
421  KVNucleus* GetParticle(const Char_t* group_name) const;
422 
423  Iterator begin() const
424  {
426  return Iterator(this);
427  }
428  Iterator end() const
429  {
431  return Iterator::End();
432  }
433 
434  KVNucleus* GetNextParticle(Option_t* opt = "") const;
435  void ResetGetNextParticle() const;
436 
437  void ResetEnergies();
438 
439  virtual Bool_t IsOK();
442 
444  {
446  }
447 
448  void DefineGroup(const Char_t* groupname, const Char_t* from = "");
449  void DefineGroup(const Char_t* groupname, KVParticleCondition* cond, const Char_t* from = "");
450 
451  void SetFrame(const Char_t* frame, const KVFrameTransform& ft);
452  void SetFrame(const Char_t* newframe, const Char_t* oldframe, const KVFrameTransform& ft);
453  void ChangeFrame(const KVFrameTransform&, const KVString& = "");
454  void ChangeDefaultFrame(const Char_t*, const Char_t* defname = "");
455  void UpdateAllFrames();
456 
457  virtual void FillIntegerList(KVIntegerList*, Option_t* opt);
458 
459  virtual void GetMasses(std::vector<Double_t>&);
460  virtual void GetGSMasses(std::vector<Double_t>&);
461  Double_t GetChannelQValue() const;
463  const Char_t* GetPartitionName();
464 
466 #define KVEVENT_MAKE_EVENT_BRANCH_NO_VOID_PTR 1
467  template<typename T>
468  static void MakeEventBranch(TTree* tree, const TString& branchname, const TString& classname, T& event, Int_t bufsize = 10000000)
469  {
477 
478  tree->Branch(branchname, classname, &event, bufsize, 0)->SetAutoDelete(kFALSE);
479  }
480 
481  virtual void MergeEventFragments(TCollection*, Option_t* opt = "");
482  static KVEvent* Factory(const char*);
483  void SetFrameName(const KVString&);
484  template<typename ValType> void SetParameter(const Char_t* name, ValType value) const
485  {
486  GetParameters()->SetValue(name, value);
487  }
488  const Char_t* GetFrameName() const
489  {
492 
493  return (GetParameters()->HasStringParameter("defaultFrame") ?
494  GetParameters()->GetStringValue("defaultFrame") : "");
495  }
496 
497  ClassDef(KVEvent, 4) //Base class for all types of multiparticle event
498 };
499 
502 #ifdef WITH_CPP11
504 #else
506 #endif
507  : it(event, t, grp)
508  {}
510  {
511  return it;
512  }
514  {
515  return KVEvent::Iterator::End();
516  }
517 };
518 
521 #ifdef WITH_CPP11
522  EventIterator(event, KVEvent::Iterator::Type::OK)
523 #else
524  EventIterator(event, KVEvent::Iterator::OK)
525 #endif
526  {}
527 };
528 
531 #ifdef WITH_CPP11
532  EventIterator(event, KVEvent::Iterator::Type::Group, grp)
533 #else
534  EventIterator(event, KVEvent::Iterator::Group, grp)
535 #endif
536  {}
537 };
538 #endif
int Int_t
#define e(i)
char Char_t
const Bool_t kFALSE
bool Bool_t
double Double_t
const Bool_t kTRUE
const char Option_t
#define ClassDef(name, id)
Base class for KaliVeda framework.
Definition: KVBase.h:135
void SetIsIterating(Bool_t on=kTRUE)
Definition: KVEvent.h:378
KVNucleus value_type
Definition: KVEvent.h:191
void Reset(Type t=Type::Null, TString grp="")
Definition: KVEvent.h:353
TIter fIter
internal iterator used by GetNextParticle()
Definition: KVEvent.h:204
Bool_t operator!=(const Iterator &it) const
Definition: KVEvent.h:308
Iterator & operator=(const Iterator &rhs)
Definition: KVEvent.h:335
Bool_t IsIterating() const
Definition: KVEvent.h:373
Iterator(const KVEvent &e, Type t=Type::All, const TString &grp="")
Definition: KVEvent.h:272
KVNucleus & operator*() const
Definition: KVEvent.h:292
std::ptrdiff_t difference_type
Definition: KVEvent.h:192
Bool_t fIterating
Definition: KVEvent.h:206
Iterator(const Iterator &i)
Definition: KVEvent.h:243
static Iterator End()
Definition: KVEvent.h:347
std::forward_iterator_tag iterator_category
Definition: KVEvent.h:190
KVNucleus & reference
Definition: KVEvent.h:194
virtual ~Iterator()
Definition: KVEvent.h:351
Bool_t AcceptableIteration()
Definition: KVEvent.h:208
TString fGroup
Definition: KVEvent.h:207
Iterator(const KVEvent *e, Type t=Type::All, const TString &grp="")
Definition: KVEvent.h:251
const Iterator & operator++()
Definition: KVEvent.h:318
Bool_t operator==(const Iterator &it) const
Definition: KVEvent.h:313
KVNucleus * current() const
Definition: KVEvent.h:226
KVNucleus * pointer
Definition: KVEvent.h:193
ReferenceType & get_reference() const
Definition: KVEvent.h:304
PointerType * get_pointer() const
Definition: KVEvent.h:299
Base class container for multi-particle events.
Definition: KVEvent.h:176
KVNameValueList * GetParameters() const
Definition: KVEvent.h:391
Double_t GetGSChannelQValue() const
Definition: KVEvent.cpp:1087
Iterator begin() const
Definition: KVEvent.h:423
virtual void FillIntegerList(KVIntegerList *, Option_t *opt)
Definition: KVEvent.cpp:982
void ChangeDefaultFrame(const Char_t *, const Char_t *defname="")
Definition: KVEvent.cpp:839
void SetFrameName(const KVString &)
Definition: KVEvent.cpp:1215
virtual void Copy(TObject &obj) const
Copy this to obj.
Definition: KVEvent.cpp:116
KVNameValueList fParameters
general-purpose list of parameters
Definition: KVEvent.h:181
virtual void Clear(Option_t *opt="")
Definition: KVEvent.cpp:200
Iterator GetNextParticleIterator(Option_t *opt) const
Definition: KVEvent.cpp:47
Int_t GetMultiplicity(Int_t Z, Int_t A=0, Option_t *opt="")
Definition: KVEvent.cpp:486
KVNucleus * GetParticle(Int_t npart) const
Definition: KVEvent.cpp:137
Iterator end() const
Definition: KVEvent.h:428
Iterator fIter
internal iterator used by GetNextParticle()
Definition: KVEvent.h:387
KVEvent(Int_t mult=50, const char *classname="KVNucleus")
Definition: KVEvent.cpp:77
virtual Bool_t IsOK()
Definition: KVEvent.cpp:636
const Char_t * GetFrameName() const
Definition: KVEvent.h:488
Int_t GetMinimumOKMultiplicity() const
Definition: KVEvent.cpp:666
static void MakeEventBranch(TTree *tree, const TString &branchname, const TString &classname, T &event, Int_t bufsize=10000000)
Definition: KVEvent.h:468
Double_t GetChannelQValue() const
Definition: KVEvent.cpp:1044
void SetMinimumOKMultiplicity(Int_t)
Definition: KVEvent.cpp:652
virtual ~ KVEvent()
void DefineGroup(const Char_t *groupname, const Char_t *from="")
Definition: KVEvent.cpp:712
TClonesArray * fParticles
array of particles in event
Definition: KVEvent.h:180
KVNucleus * AddParticle()
Definition: KVEvent.cpp:166
void ResetEnergies()
Definition: KVEvent.cpp:688
KVNucleus * GetParticleWithName(const Char_t *name) const
Definition: KVEvent.cpp:242
virtual void MergeEventFragments(TCollection *, Option_t *opt="")
Definition: KVEvent.cpp:1164
void GetMultiplicities(Int_t mult[], const TString &species, Option_t *opt="")
Definition: KVEvent.cpp:513
virtual Int_t GetMult(Option_t *opt="") const
Definition: KVEvent.cpp:278
virtual void Print(Option_t *t="") const
Definition: KVEvent.cpp:224
virtual void ls(Option_t *t="") const
Definition: KVEvent.h:416
KVNucleus * GetNextParticle(Option_t *opt="") const
Definition: KVEvent.cpp:564
void ChangeFrame(const KVFrameTransform &, const KVString &="")
Definition: KVEvent.cpp:815
void SetParameter(const Char_t *name, ValType value) const
Definition: KVEvent.h:484
void UpdateAllFrames()
Definition: KVEvent.cpp:861
void CustomStreamer()
Definition: KVEvent.h:443
void ResetGetNextParticle() const
Definition: KVEvent.cpp:618
Double_t GetSum(const Char_t *KVNucleus_method, Option_t *opt="")
Definition: KVEvent.cpp:310
const Char_t * GetPartitionName()
Definition: KVEvent.cpp:1123
virtual void GetMasses(std::vector< Double_t > &)
Fill vector with mass of each nucleus of event (in MeV) [note: this is the mass including any excitat...
Definition: KVEvent.cpp:1001
void FillHisto(TH1 *h, const Char_t *KVNucleus_method, Option_t *opt="")
Definition: KVEvent.cpp:403
virtual void GetGSMasses(std::vector< Double_t > &)
Fill vector with ground state mass of each nucleus of event (in MeV).
Definition: KVEvent.cpp:1016
void SetFrame(const Char_t *frame, const KVFrameTransform &ft)
Definition: KVEvent.cpp:762
static KVEvent * Factory(const char *)
Create and return pointer to new event of class given by plugin.
Definition: KVEvent.cpp:1194
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.
void SetValue(const Char_t *name, value_type value)
Description of properties and kinematics of atomic nuclei.
Definition: KVNucleus.h:125
Handles particle selection criteria for data analysis classes ,.
virtual Bool_t IsOK()
Definition: KVParticle.cpp:312
Bool_t BelongsToGroup(const Char_t *groupname) const
Definition: KVParticle.cpp:597
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition: KVString.h:72
void BypassStreamer(Bool_t bypass=kTRUE)
TIter & Begin()
Type
auto All(const RVec< T > &v) -> decltype(v[0]==false)
EventIterator(KVEvent &event, KVEvent::Iterator::Type t=KVEvent::Iterator::Type::All, const TString &grp="")
Definition: KVEvent.h:503
KVEvent::Iterator it
Definition: KVEvent.h:501
KVEvent::Iterator begin() const
Definition: KVEvent.h:509
KVEvent::Iterator end() const
Definition: KVEvent.h:513
GroupEventIterator(KVEvent &event, const TString &grp)
Definition: KVEvent.h:530
OKEventIterator(KVEvent &event)
Definition: KVEvent.h:520