KaliVeda  1.12/06
Heavy-Ion Analysis Toolkit
KVEvent.cpp
Go to the documentation of this file.
1 /***************************************************************************
2 $Id: KVEvent.cpp,v 1.41 2008/12/17 11:23:12 ebonnet Exp $
3  kvevent.cpp - description
4  -------------------
5  begin : Sun May 19 2002
6  copyright : (C) 2002 by J.D. Frankland
7  email : frankland@ganil.fr
8  ***************************************************************************/
9 
10 /***************************************************************************
11  * *
12  * This program is free software; you can redistribute it and/or modify *
13  * it under the terms of the GNU General Public License as published by *
14  * the Free Software Foundation; either version 2 of the License, or *
15  * (at your option) any later version. *
16  * *
17  ***************************************************************************/
18 
19 #include "Riostream.h"
20 #include "TClass.h"
21 #include "KVEvent.h"
22 #include "TStreamerInfo.h"
23 #include "KVParticleCondition.h"
24 #include "TClass.h"
25 #include "KVIntegerList.h"
26 
27 using namespace std;
28 
30 
32 /*
33 <h2>KVEvent</h2>
34 
35 */
37 
38 
39 
46 
48 {
49  // Provide correct iterator using same options as for GetNextParticle() method:
50  //
51  // - if opt="" (default) => iterator over all particles
52  // - if opt="ok"/"OK" => iterator for all "OK" particles
53  // - if opt!="" && opt!="ok"/"OK" => iterator for all particles in group with name given by opt (case-insensitive)
54 
55  TString Opt(opt);
56  Opt.ToUpper();
57  Bool_t ok_iter = (Opt == "OK");
58  Bool_t grp_iter = (!ok_iter && Opt.Length());
59  KVEvent& event = const_cast<KVEvent&>(*this);
60  if (ok_iter) return OKEventIterator(event).begin();
61  else if (grp_iter) return GroupEventIterator(event, Opt).begin();
62  return EventIterator(event).begin();
63 }
64 
65 
66 
76 
77 KVEvent::KVEvent(Int_t mult, const char* classname)
78  : fParameters("EventParameters", "Parameters associated with an event")
79 {
80  //Initialise KVEvent to hold mult events of "classname" objects
81  //(the class must inherit from KVNucleus).
82  //"mult" is the approximate maximum multiplicity of the events (if too
83  //small, extra space will be allocated automatically by TClonesArray).
84  //
85  // Default argument :
86  // classname = "KVNucleus"
87  //
88 
89  fParticles = new TClonesArray(classname, mult);
90  CustomStreamer();//force use of KVEvent::Streamer function for reading/writing
91 }
92 
93 
94 
95 
99 
100 KVEvent::~KVEvent()
101 {
102  //Destructor. Destroys all objects stored in TClonesArray and releases
103  //allocated memory.
104 
105  fParticles->Delete();
107 }
108 
109 
110 
111 #if ROOT_VERSION_CODE >= ROOT_VERSION(3,4,0)
112 
115 
116 void KVEvent::Copy(TObject& obj) const
117 #else
118 void KVEvent::Copy(TObject& obj)
119 #endif
120 {
121  //Copy this to obj
122  KVBase::Copy(obj);
123  fParameters.Copy(((KVEvent&)obj).fParameters);
124  for (Int_t nn = 0; nn < fParticles->GetEntriesFast(); nn += 1) {
125  //printf("avant=%s - ",GetParticle(nn+1)->GetName());
126  GetParticle(nn + 1)->Copy(*((KVEvent&) obj).AddParticle());
127  //printf("apres=%s\n",((KVEvent &) obj).GetParticle(nn+1)->GetName());
128  }
129 }
130 
131 
132 
136 
138 {
139  //Access to event member with index npart (1<=npart<=GetMult() : error if out of bounds)
140  //This method may be overridden in event classes derived from KVEvent.
141 
142  if (npart < 1 || npart > fParticles->GetEntriesFast()) {
143  Error("GetParticle", KVEVENT_PART_INDEX_OOB, npart,
145  return 0;
146  }
147 
148  return (KVNucleus*)((*fParticles)[npart - 1]);
149 }
150 
151 
152 
153 
165 
167 {
168  // Method used for building an event particle by particle.
169  // DO NOT USE FOR READING EVENTS - use GetParticle(Int_t npart)!!
170  //
171  // This method increases the multiplicity fMult by one and "creates" a new particle with index (fMult-1).
172  // In actual fact a new object is only created if needed i.e. if the new multiplicity is greater
173  // than previously. Particle objects in the array are reused from one event to another.
174  // The default constructor for the class corresponding to the "new" particle will only be called
175  // once during its lifetime (i.e. if N events are generated, the particle ctor will be called only once,
176  // not N times). Once created, in subsequent events we just call the particle's Clear() method
177  // in order to reset its internal variables ready for a new event.
178 
179  Int_t mult = GetMult();
180 #ifdef __WITHOUT_TCA_CONSTRUCTED_AT
181  KVNucleus* tmp = (KVNucleus*) ConstructedAt(mult, "C");
182 #else
183  KVNucleus* tmp = (KVNucleus*) fParticles->ConstructedAt(mult, "C");
184 #endif
185  if (!tmp) {
186  Error("AddParticle", "Allocation failure, Mult=%d", mult);
187  return 0;
188  }
189  return tmp;
190 }
191 
192 
193 
194 
199 
201 {
202  // Reset the event to zero ready for new event.
203  // Any option string is passed on to the Clear() method of the particle
204  // objects in the TClonesArray fParticles.
205 
206  if (strcmp(opt, "")) {
207  // pass options to particle class Clear() method
208  TString Opt = Form("C+%s", opt);
209  fParticles->Clear(Opt);
210  }
211  else
212  fParticles->Clear("C");
213  fParameters.Clear();
215 }
216 
217 
218 
219 
223 
224 void KVEvent::Print(Option_t* t) const
225 {
226  //Print a list of all particles in the event with some characteristics.
227  //Optional argument t can be used to select particles (="ok", "groupname", ...)
228 
229  cout << "\nKVEvent with " << GetMult(t) << " particles :" << endl;
230  cout << "------------------------------------" << endl;
231  fParameters.Print();
232  for (Iterator it = GetNextParticleIterator(t); it != end(); ++it)(*it).Print();
233 }
234 
235 
236 
237 
241 
243 {
244  //Find particle using its name (SetName()/GetName() methods)
245  //In case more than one particle has the same name, the first one found is returned.
246 
247  KVNucleus* tmp = (KVNucleus*)fParticles->FindObject(name);
248  return tmp;
249 }
250 
251 
252 
253 
256 
257 KVNucleus* KVEvent::GetParticle(const Char_t* group_name) const
258 {
259  // Find first particle in event belonging to group with name "group_name"
260 
261  Iterator it = GetNextParticleIterator(group_name);
262  KVNucleus* tmp = it.get_pointer<KVNucleus>();
263  if (!tmp) Warning("GetParticle", "Particle not found: %s", group_name);
264  return tmp;
265 }
266 
267 
268 
269 
277 
279 {
280  //Returns multiplicity (number of particles) of event.
281  //If opt = "" (default), returns number of particles in TClonesArray* fParticles
282  // i.e. the value of fParticles->GetEntriesFast() (we assume there are no gaps
283  // in the list)
284  //If opt = "ok" only particles with IsOK()==kTRUE are included.
285  //If opt = "name" only particles belonging to group "name" are included.
286 
287  Int_t fMultOK = 0;
288  if (TString(opt) == "") {
289  //get total multiplicity of event
290  return fParticles->GetEntriesFast();
291  }
292  else {
293  for (Iterator it = GetNextParticleIterator(opt); it != end(); ++it) ++fMultOK;
294  }
295  return fMultOK;
296 }
297 
298 
299 
300 
309 
310 Double_t KVEvent::GetSum(const Char_t* KVNucleus_method, Option_t* opt)
311 {
312  //Returns sum over particles of the observable given by the indicated KVNucleus_method
313  //for example
314  //if the method is called this way GetSum("GetZ"), it returns the sum of the charge
315  //of particles in the current event
316  //
317  //If opt = "ok" only particles with IsOK()==kTRUE are considered.
318  //If opt = "name" only particles belonging to group "name" are considered.
319 
320  Double_t fSum = 0;
321  TMethodCall mt;
322  mt.InitWithPrototype(KVNucleus::Class(), KVNucleus_method, "");
323 
324  if (mt.IsValid()) {
326  if (mt.ReturnType() == TMethodCall::kLong) {
327  Long_t ret;
328  for (; it != end(); ++it) {
329  KVNucleus* tmp = it.get_pointer<KVNucleus>();
330  mt.Execute(tmp, "", ret);
331  fSum += ret;
332  }
333  }
334  else if (mt.ReturnType() == TMethodCall::kDouble) {
335  Double_t ret;
336  for (; it != end(); ++it) {
337  KVNucleus* tmp = it.get_pointer<KVNucleus>();
338  mt.Execute(tmp, "", ret);
339  fSum += ret;
340  }
341  }
342  }
343 
344  return fSum;
345 }
346 
347 
348 
356 
357 Double_t KVEvent::GetSum(const Char_t* KVNucleus_method, const Char_t* method_prototype, const Char_t* args, Option_t* opt)
358 {
359  //Returns sum over particles of the observable given by the indicated KVNucleus_method
360  // with given prototype (e.g. method_prototype="int,int") and argument values
361  // e.g. args="2,4")
362  //
363  //If opt = "ok" only particles with IsOK()==kTRUE are considered.
364  //If opt = "name" only particles belonging to group "name" are considered.
365 
366  Double_t fSum = 0;
367  TMethodCall mt;
368  mt.InitWithPrototype(KVNucleus::Class(), KVNucleus_method, method_prototype);
369 
370  if (mt.IsValid()) {
372  if (mt.ReturnType() == TMethodCall::kLong) {
373  Long_t ret;
374  for (; it != end(); ++it) {
375  KVNucleus* tmp = it.get_pointer<KVNucleus>();
376  mt.Execute(tmp, args, ret);
377  fSum += ret;
378  }
379  }
380  else if (mt.ReturnType() == TMethodCall::kDouble) {
381  Double_t ret;
382  for (; it != end(); ++it) {
383  KVNucleus* tmp = it.get_pointer<KVNucleus>();
384  mt.Execute(tmp, args, ret);
385  fSum += ret;
386  }
387  }
388  }
389 
390  return fSum;
391 }
392 
393 
394 
402 
403 void KVEvent::FillHisto(TH1* h, const Char_t* KVNucleus_method, Option_t* opt)
404 {
405  // Fill histogram with values of the observable given by the indicated KVNucleus_method.
406  // For example: if the method is called this way - FillHisto(h,"GetZ") - it fills histogram
407  // with the charge of all particles in the current event.
408  //
409  // If opt = "ok" only particles with IsOK()==kTRUE are considered.
410  // If opt = "name" only particles belonging to group "name" are considered.
411 
412  TMethodCall mt;
413  mt.InitWithPrototype(KVNucleus::Class(), KVNucleus_method, "");
414 
415  if (mt.IsValid()) {
417  if (mt.ReturnType() == TMethodCall::kLong) {
418  Long_t ret;
419  for (; it != end(); ++it) {
420  KVNucleus* tmp = it.get_pointer<KVNucleus>();
421  mt.Execute(tmp, "", ret);
422  h->Fill((Double_t)ret);
423  }
424  }
425  else if (mt.ReturnType() == TMethodCall::kDouble) {
426  Double_t ret;
427  for (; it != end(); ++it) {
428  KVNucleus* tmp = it.get_pointer<KVNucleus>();
429  mt.Execute(tmp, "", ret);
430  h->Fill(ret);
431  }
432  }
433  }
434 }
435 
436 
437 
444 
445 void KVEvent::FillHisto(TH1* h, const Char_t* KVNucleus_method, const Char_t* method_prototype, const Char_t* args, Option_t* opt)
446 {
447  // Fill histogram with values of given method with given prototype (e.g. method_prototype="int,int") and argument values
448  // e.g. args="2,4") for each particle in event.
449  //
450  // If opt = "ok" only particles with IsOK()==kTRUE are considered.
451  // If opt = "name" only particles belonging to group "name" are considered.
452 
453  TMethodCall mt;
454  mt.InitWithPrototype(KVNucleus::Class(), KVNucleus_method, method_prototype);
455 
456  if (mt.IsValid()) {
458  if (mt.ReturnType() == TMethodCall::kLong) {
459  Long_t ret;
460  for (; it != end(); ++it) {
461  KVNucleus* tmp = it.get_pointer<KVNucleus>();
462  mt.Execute(tmp, args, ret);
463  h->Fill((Double_t)ret);
464  }
465  }
466  else if (mt.ReturnType() == TMethodCall::kDouble) {
467  Double_t ret;
468  for (; it != end(); ++it) {
469  KVNucleus* tmp = it.get_pointer<KVNucleus>();
470  mt.Execute(tmp, args, ret);
471  h->Fill(ret);
472  }
473  }
474  }
475 }
476 
477 
478 
485 
487 {
488  // Calculate the multiplicity of nuclei given Z (if A not given)
489  // or of nuclei with given Z & A (if given)
490  //
491  //If opt = "ok" only particles with IsOK()==kTRUE are considered.
492  //If opt = "name" only particles belonging to group "name" are considered.
493 
494  if (A > 0) return (Int_t)GetSum("IsIsotope", "int,int", Form("%d,%d", Z, A), opt);
495  return (Int_t)GetSum("IsElement", "int", Form("%d", Z), opt);
496 }
497 
498 
499 
512 
513 void KVEvent::GetMultiplicities(Int_t mult[], const TString& species, Option_t* opt)
514 {
515  // Fill array mult[] with the number of each nuclear species in the
516  // comma-separated list in this event. Make sure that mult[] is
517  // large enough for the list.
518  //
519  // Example:
520  // Int_t mult[4];
521  // event.GetMultiplicities(mult, "1n,1H,2H,3H");
522  //
523  // N.B. the species name must correspond to that given by KVNucleus::GetSymbol
524  //
525  // If given, "opt" will be used to select particles ("OK" or groupname)
526 
527  unique_ptr<TObjArray> spec(species.Tokenize(", "));// remove any spaces
528  Int_t nspec = spec->GetEntries();
529  memset(mult, 0, nspec * sizeof(Int_t)); // set multiplicities to zero
530  for (Iterator it = GetNextParticleIterator(opt); it != end(); ++it) {
531  for (int i = 0; i < nspec; i++) {
532  if (((TObjString*)(*spec)[i])->String() == (*it).GetSymbol()) mult[i] += 1;
533  }
534  }
535 }
536 
537 
539 
540 
563 
565 {
566  // Use this method to iterate over the list of particles in the event
567  // After the last particle GetNextParticle() returns a null pointer and
568  // resets itself ready for a new iteration over the particle list.
569  //
570  // If opt="" all particles are included in the iteration.
571  // If opt="ok" or "OK" only particles whose IsOK() method returns kTRUE are included.
572  //
573  // Any other value of opt is interpreted as a (case-insensitive) particle group name: only
574  // particles with BelongsToGroup(opt) returning kTRUE are included.
575  //
576  // If you want to start from the beginning again before getting to the end
577  // of the list, especially if you want to change the selection criteria,
578  // call method ResetGetNextParticle() before continuing.
579  // If you interrupt an iteration before the end, then start another iteration
580  // without calling ResetGetNextParticle(), even if you change the argument of
581  // the call to GetNextParticle(), you will repeat exactly the same iteration
582  // as the previous one.
583  //
584  // WARNING: Only one iteration at a time over the event can be performed
585  // using this method. If you want/need to perform several i.e. nested
586  // iterations, use the KVEvent::Iterator class
587 
588  TString Opt(opt);
589  Opt.ToUpper();
590 
591  // continue iteration
592  if (fIter.IsIterating()) return &(*(fIter++));
593 
594  // start new iteration
595  Bool_t ok_iter = (Opt == "OK");
596  Bool_t grp_iter = (!ok_iter && Opt.Length());
597 
598 #ifdef WITH_CPP11
599  if (ok_iter) fIter = Iterator(this, Iterator::Type::OK);
600  else if (grp_iter) fIter = Iterator(this, Iterator::Type::Group, Opt);
601  else fIter = Iterator(this, Iterator::Type::All);
602 #else
603  if (ok_iter) fIter = Iterator(this, Iterator::OK);
604  else if (grp_iter) fIter = Iterator(this, Iterator::Group, Opt);
605  else fIter = Iterator(this, Iterator::All);
606 #endif
607  return &(*(fIter++));
608 }
609 
610 
611 
612 
617 
619 {
620  // Reset iteration over event particles so that next call
621  // to GetNextParticle will begin a new iteration (possibly with
622  // different criteria).
623 
625 }
626 
627 
628 
629 
635 
637 {
638  // Returns kTRUE if the event is OK for analysis.
639  // This means there must be at least MOKmin particles with IsOK()=kTRUE,
640  // where MOKmin is set by calling SetMinimumOKMultiplicity(Int_t)
641  // (value stored in parameter MIN_OK_MULT)
642 
643  return (GetMult("ok") >= GetMinimumOKMultiplicity());
644 }
645 
646 
647 
651 
653 {
654  // Set minimum number of particles with IsOK()=kTRUE in event for
655  // it to be considered 'good' for analysis
656  SetParameter("MIN_OK_MULT", x);
657 }
658 
659 
660 
665 
667 {
668  // Get minimum number of particles with IsOK()=kTRUE in event for
669  // it to be considered 'good' for analysis
670  // NB: if no minimum has been set, we return 1
671  Int_t x = GetParameters()->GetIntValue("MIN_OK_MULT");
672  if (x == -1) return 1;
673  return x;
674 }
675 
676 
677 
678 
687 
689 {
690  //Used for simulated events after "detection" by some multidetector array.
691  //
692  //The passage of the event's particles through the different absorbers modifies their
693  //kinetic energies, indeed all those which are correctly identified by the detector
694  //actually stop. Calling this method will reset all the particles' energies to their
695  //initial value i.e. before they entered the first absorber.
696  //Particles which have not encountered any absorbers/detectors are left as they are.
697 
698  for (Iterator it = begin(); it != end(); ++it) {
699  (*it).ResetEnergy();
700  }
701 }
702 
703 
704 
705 
711 
712 void KVEvent::DefineGroup(const Char_t* groupname, const Char_t* from)
713 {
714  //In the same philosophy as KVINDRAReconEvent::AcceptIDCodes() method
715  // allow to affiliate a group name to particles of the event
716  // if "from" is not null, a test of previously stored group name
717  // such as "OK" is checked
718 
719  for (Iterator it = GetNextParticleIterator(from); it != end(); ++it) {
720  (*it).AddGroup(groupname);
721  }
722 }
723 
724 
725 
726 
735 
736 void KVEvent::DefineGroup(const Char_t* groupname, KVParticleCondition* cond, const Char_t* from)
737 {
738  //In the same philosophy as KVINDRAReconEvent::AcceptIDCodes() method
739  // allow to affiliate using a KVParticleCondition a group name to particles of the event
740  // if "from" is not null, a test of previously stored group name
741  // such as "OK" is checked
742  // the method used in KVParticleCondition has to be compatible with the KVNucleus
743  // concerned class.
744  // This can be done using first KVParticleCondition::SetParticleClassName(const Char_t* cl)
745 
746  for (Iterator it = GetNextParticleIterator(from); it != end(); ++it) {
747  (*it).AddGroup(groupname, cond);
748  }
749 }
750 
751 
752 
753 
754 
761 
762 void KVEvent::SetFrame(const Char_t* frame, const KVFrameTransform& ft)
763 {
764  //Define a Lorentz-boosted and/or rotated frame for all "ok" particles in the event.
765  //See KVParticle::SetFrame() for details.
766  //
767  //In order to access the kinematics in the boosted frame, use the GetFrame() method of the
768  //individual particles (see KVParticle::GetFrame()).
769 
770  for (Iterator it = GetNextParticleIterator("ok"); it != end(); ++it) {
771  (*it).SetFrame(frame, ft);
772  }
773 }
774 
775 
776 
777 
788 
789 void KVEvent::SetFrame(const Char_t* newframe, const Char_t* oldframe, const KVFrameTransform& ft)
790 {
791  //Define a Lorentz-boosted frame "newframe" for all "ok" particles in the event.
792  //The transformation is applied to the particle coordinates in the existing frame "oldframe"
793  //
794  //See KVParticle::SetFrame() for details.
795  //
796  //In order to access the kinematics in the boosted frame, use the GetFrame() method of the
797  //individual particles in either of these ways :
798  // KVParticle* newframe = particle->GetFrame("newframe");
799  // KVParticle* newframe = particle->GetFrame("oldframe")->GetFrame("newframe");
800 
801  for (Iterator it = GetNextParticleIterator("ok"); it != end(); ++it) {
802  (*it).SetFrame(newframe, oldframe, ft);
803  }
804 }
805 
806 
807 
814 
815 void KVEvent::ChangeFrame(const KVFrameTransform& ft, const KVString& name)
816 {
817  //Permanently change the reference frame used for particle kinematics in the event.
818  //The transformation is applied to all "ok" particles in the event.
819  //You can optionally set the name of this new default kinematical reference frame.
820  //
821  //See KVParticle::ChangeFrame() and KVParticle::SetFrame() for details.
822 
823 
824  for (Iterator it = GetNextParticleIterator("ok"); it != end(); ++it) {
825  (*it).ChangeFrame(ft, name);
826  }
827  if (name != "") SetParameter("defaultFrame", name);
828 }
829 
830 
831 
838 
839 void KVEvent::ChangeDefaultFrame(const Char_t* newdef, const Char_t* defname)
840 {
841  // Make existing reference frame 'newdef' the new default frame for particle kinematics.
842  // The current default frame will then be accessible from the list of frames
843  // using its name (previously set with SetFrameName). You can change this name with 'defname'.
844  //
845  //See KVParticle::ChangeDefaultFrame() and KVParticle::SetFrame() for details.
846 
847  for (Iterator it = GetNextParticleIterator("ok"); it != end(); ++it) {
848  (*it).ChangeDefaultFrame(newdef, defname);
849  }
850  SetParameter("defaultFrame", newdef);
851 }
852 
853 
854 
860 
862 {
863  // If the kinematics of particles in their default reference frame have been modified,
864  // call this method to update the kinematics in all defined reference frames.
865  //
866  //See KVParticle::UpdateAllFrames() for details.
867 
868  for (Iterator it = GetNextParticleIterator("ok"); it != end(); ++it) {
869  (*it).UpdateAllFrames();
870  }
871 }
872 
873 
874 
875 
880 
881 void KVEvent::Streamer(TBuffer& R__b)
882 {
883  // Customised Streamer for KVEvent.
884  // This is just the automatic Streamer with the addition of a call to the Clear()
885  // method before reading a new object (avoid memory leaks with lists of parameters).
886 
887  if (R__b.IsReading()) {
888  Clear();
889  R__b.ReadClassBuffer(KVEvent::Class(), this);
890  }
891  else {
892  R__b.WriteClassBuffer(KVEvent::Class(), this);
893  }
894 }
895 
896 
897 
898 #ifdef __WITHOUT_TCA_CONSTRUCTED_AT
899 
913 
914 TObject* KVEvent::ConstructedAt(Int_t idx)
915 {
916  // Get an object at index 'idx' that is guaranteed to have been constructed.
917  // It might be either a freshly allocated object or one that had already been
918  // allocated (and assumingly used). In the later case, it is the callers
919  // responsability to insure that the object is returned to a known state,
920  // usually by calling the Clear method on the TClonesArray.
921  //
922  // Tests to see if the destructor has been called on the object.
923  // If so, or if the object has never been constructed the class constructor is called using
924  // New(). If not, return a pointer to the correct memory location.
925  // This explicitly to deal with TObject classes that allocate memory
926  // which will be reset (but not deallocated) in their Clear()
927  // functions.
928 
929  TObject* obj = (*fParticles)[idx];
930  if (obj && obj->TestBit(TObject::kNotDeleted)) {
931  return obj;
932  }
933  return (fParticles->GetClass()) ? static_cast<TObject*>(fParticles->GetClass()->New(obj)) : 0;
934 }
935 
936 
949 
950 TObject* KVEvent::ConstructedAt(Int_t idx, Option_t* clear_options)
951 {
952  // Get an object at index 'idx' that is guaranteed to have been constructed.
953  // It might be either a freshly allocated object or one that had already been
954  // allocated (and assumingly used). In the later case, the function Clear
955  // will be called and passed the value of 'clear_options'
956  //
957  // Tests to see if the destructor has been called on the object.
958  // If so, or if the object has never been constructed the class constructor is called using
959  // New(). If not, return a pointer to the correct memory location.
960  // This explicitly to deal with TObject classes that allocate memory
961  // which will be reset (but not deallocated) in their Clear()
962  // functions.
963 
964  TObject* obj = (*fParticles)[idx];
965  if (obj && obj->TestBit(TObject::kNotDeleted)) {
966  obj->Clear(clear_options);
967  return obj;
968  }
969  return (fParticles->GetClass()) ? static_cast<TObject*>(fParticles->GetClass()->New(obj)) : 0;
970 }
971 
972 #endif
973 
974 
981 
983 {
984  // Clear & fill the KVIntegerList with the contents of this event,
985  // the option will be passed to GetNextParticle(opt).
986  // IntegerList is then 'Update()'d.
987  // (This method was originally KVIntegerList::Fill(KVEvent*,Option_t*),
988  // it was moved here in order to make KVIntegerList a base class)
989 
990  IL->Clear();
991  for (Iterator it = GetNextParticleIterator(opt); it != end(); ++it) IL->Add((*it).GetZ());
992  IL->SetPopulation(1);
993  IL->CheckForUpdate();
994 }
995 
996 
997 
1000 
1001 void KVEvent::GetMasses(std::vector<Double_t>& mass)
1002 {
1003  // Fill vector with mass of each nucleus of event (in MeV) [note: this is the mass including any excitation energy, not ground state]
1004 
1005  mass.clear();
1006  mass.reserve(GetMult());
1007  int i = 0;
1008  for (Iterator it = begin(); it != end(); ++it) mass[i++] = (*it).GetMass();
1009 }
1010 
1011 
1012 
1015 
1016 void KVEvent::GetGSMasses(std::vector<Double_t>& mass)
1017 {
1018  // Fill vector with ground state mass of each nucleus of event (in MeV).
1019 
1020  mass.clear();
1021  mass.reserve(GetMult());
1022  int i = 0;
1023  for (Iterator it = begin(); it != end(); ++it) mass[i++] = (*it).GetMassGS();
1024 }
1025 
1026 
1027 
1043 
1045 {
1046  // Calculate the Q-value [MeV] for this event as if all nuclei were produced by
1047  // the decay of an initial compound nucleus containing the sum of all nuclei
1048  // in the event, i.e.
1049  // A -> a1 + a2 + a3 + ...
1050  // We take into account any excitation energy of the nuclei of the event
1051  // (see GetGSChannelQValue() for an alternative), i.e. we calculate
1052  // Q = M(A) - ( m(a1) + m(a2) + m(a3) + ... )
1053  // where
1054  // M(X) = ground state mass of X
1055  // m(X) = M(X) + E*(X)
1056  // If Q<0, the excitation energy of the initial compound nucleus, A,
1057  // would have to be at least equal to (-Q) in order for the decay to occur.
1058  // i.e. decay is possible if
1059  // E*(A) > -Q
1060 
1061  Double_t sumM = 0;
1062  KVNucleus CN;
1063  Int_t M = GetMult();
1064  for (int i = 1; i <= M; i++) {
1065  sumM += GetParticle(i)->GetMass();
1066  CN += *(GetParticle(i));
1067  }
1068  return CN.GetMassGS() - sumM;
1069 }
1070 
1071 
1072 
1086 
1088 {
1089  // Calculate the Q-value [MeV] for this event as if all nuclei were produced by
1090  // the decay of an initial compound nucleus containing the sum of all nuclei
1091  // in the event, i.e.
1092  // A -> a1 + a2 + a3 + ...
1093  // i.e. we calculate
1094  // Q = M(A) - ( M(a1) + M(a2) + M(a3) + ... )
1095  // where
1096  // M(X) = ground state mass of X
1097  // If Q<0, the excitation energy of the initial compound nucleus, A,
1098  // would have to be at least equal to (-Q) in order for the decay to occur.
1099  // i.e. decay is possible if
1100  // E*(A) > -Q
1101 
1102  Double_t sumM = 0;
1103  KVNucleus CN;
1104  Int_t M = GetMult();
1105  for (int i = 1; i <= M; i++) {
1106  sumM += GetParticle(i)->GetMassGS();
1107  CN += *(GetParticle(i));
1108  }
1109  return CN.GetMassGS() - sumM;
1110 }
1111 
1112 
1113 
1122 
1124 {
1125  //
1126  //return list of isotopes of the event with the format :
1127  // symbol1(population1) symbol2(population2) ....
1128  // if population==1, it is not indicated :
1129  // Example :
1130  // 15C 12C(2) 4He 3He 1H(4) 1n(3)
1131  //
1132  fParticles->Sort();
1133  static KVString partition;
1134 
1135  KVNameValueList nvl;
1136  partition = "";
1137  for (Iterator it = begin(); it != end(); ++it) {
1138  TString st = (*it).GetSymbol();
1139  Int_t pop = TMath::Max(nvl.GetIntValue(st.Data()), 0);
1140  pop += 1;
1141  nvl.SetValue(st.Data(), pop);
1142  }
1143  for (Int_t ii = 0; ii < nvl.GetEntries(); ii += 1) {
1144  Int_t pop = nvl.GetIntValue(ii);
1145  if (pop == 1) partition += nvl.GetNameAt(ii);
1146  else partition += Form("%s(%d)", nvl.GetNameAt(ii), pop);
1147  if (ii < nvl.GetEntries() - 1) partition += " ";
1148  }
1149  return partition.Data();
1150 }
1151 
1152 
1153 
1163 
1165 {
1166  // Merge all events in the list into one event (this one)
1167  // We also merge/sum the parameter lists of the events
1168  // First we clear this event, then we fill it with copies of each particle in each event
1169  // in the list.
1170  // If option "opt" is given, it is given as argument to each call to
1171  // KVEvent::Clear() - this option is then passed on to the Clear()
1172  // method of each particle in each event.
1173  // NOTE: the events in the list will be empty and useless after this!
1174 
1175  Clear(opt);
1176  TIter it(events);
1177  KVEvent* e;
1178  while ((e = (KVEvent*)it())) {
1179  KVNucleus* n;
1180  e->ResetGetNextParticle();
1181  while ((n = e->GetNextParticle())) {
1182  n->Copy(*AddParticle());
1183  }
1184  GetParameters()->Merge(*(e->GetParameters()));
1185  e->Clear(opt);
1186  }
1187 }
1188 
1189 
1190 
1193 
1194 KVEvent* KVEvent::Factory(const char* plugin)
1195 {
1196  // Create and return pointer to new event of class given by plugin
1197 
1198  TPluginHandler* ph = LoadPlugin("KVEvent", plugin);
1199  if (ph) {
1200  return (KVEvent*)ph->ExecPlugin(0);
1201  }
1202  return nullptr;
1203 }
1204 
1205 
1206 
1214 
1216 {
1217  // Set name of default frame for all particles in event
1218  // After using this method, calls to
1219  // KVParticle::GetFrame(name)
1220  // will return the address of the particle in question, i.e.
1221  // its default kinematics
1222  // The default frame name is stored as a parameter "defaultFrame"
1223 
1224 #ifdef WITH_CPP11
1225  for (KVEvent::Iterator it = std::begin(*this); it != std::end(*this); ++it)
1226 #else
1227  for (KVEvent::Iterator it = begin(); it != end(); ++it)
1228 #endif
1229  {
1230  (*it).SetFrameName(name);
1231  }
1232  SetParameter("defaultFrame", name);
1233 }
1234 
1235 
1237 
1238 
int Int_t
long Long_t
#define KVEVENT_PART_INDEX_OOB
Definition: KVEvent.h:23
KVIonRangeTableMaterial * M
ClassImp(KVPartitionList) void KVPartitionList
Initialisation.
#define SafeDelete(p)
#define e(i)
char Char_t
const Bool_t kFALSE
bool Bool_t
double Double_t
const char Option_t
char * Form(const char *fmt,...)
virtual void Copy(TObject &) const
Make a copy of this object.
Definition: KVBase.cpp:397
static TPluginHandler * LoadPlugin(const Char_t *base, const Char_t *uri="0")
Definition: KVBase.cpp:756
void SetIsIterating(Bool_t on=kTRUE)
Definition: KVEvent.h:378
Bool_t IsIterating() const
Definition: KVEvent.h:373
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
Int_t GetMinimumOKMultiplicity() const
Definition: KVEvent.cpp:666
Double_t GetChannelQValue() const
Definition: KVEvent.cpp:1044
void SetMinimumOKMultiplicity(Int_t)
Definition: KVEvent.cpp:652
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
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
void Add(TArrayI *tab)
void Fill(Double_t* tab,Int_t mult);
void Clear(Option_t *option="")
Classe dérivée de TNamed, Reinitialisation de l'object.
void SetPopulation(Int_t pop)
Initialise la population à "pop".
Handles lists of named parameters with different types, a list of KVNamedParameter objects.
virtual void Print(Option_t *opt="") const
Int_t GetIntValue(const Char_t *name) const
void SetValue(const Char_t *name, value_type value)
virtual void Clear(Option_t *opt="")
const Char_t * GetNameAt(Int_t idx) const
Int_t GetEntries() const
void Merge(const KVNameValueList &)
Description of properties and kinematics of atomic nuclei.
Definition: KVNucleus.h:125
Double_t GetMassGS() const
Definition: KVNucleus.h:292
Handles particle selection criteria for data analysis classes ,.
Double_t GetMass() const
Definition: KVParticle.h:531
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition: KVString.h:72
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=0)=0
Bool_t IsReading() const
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
void * New(ENewType defConstructor=kClassNew, Bool_t quiet=kFALSE) const
virtual void Delete(Option_t *option="")
virtual void Clear(Option_t *option="")
TObject * ConstructedAt(Int_t idx)
virtual void Sort(Int_t upto=kMaxInt)
TClass * GetClass() const
virtual Int_t Fill(const char *name, Double_t w)
EReturnType ReturnType()
static const EReturnType kLong
void InitWithPrototype(const char *function, const char *proto, ROOT::EFunctionMatchMode mode=ROOT::kConversionMatch)
Bool_t IsValid() const
static const EReturnType kDouble
void Execute()
Int_t GetEntriesFast() const
virtual TObject * FindObject(const char *name) const
friend friend class TClonesArray
virtual void Clear(Option_t *="")
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
virtual void Warning(const char *method, const char *msgfmt,...) const
virtual void Error(const char *method, const char *msgfmt,...) const
Longptr_t ExecPlugin(int nargs, const T &... params)
Ssiz_t Length() const
void ToUpper()
TObjArray * Tokenize(const TString &delim) const
const char * Data() const
auto All(const RVec< T > &v) -> decltype(v[0]==false)
Double_t x[n]
const Int_t n
TH1 * h
Double_t Max(Double_t a, Double_t b)
KVEvent::Iterator begin() const
Definition: KVEvent.h:509