KaliVeda
Toolkit for HIC analysis
KVGVList.cpp
1 // Author: Daniel Cussol
2 //
3 // 17/02/2004
4 #include "KVGVList.h"
5 #include "KVNucleusEvent.h"
6 
7 #include <KVUnownedList.h>
8 
10 
11 
12 //_________________________________________________________________
13 
14 
16 void KVGVList::init_KVGVList(void)
17 {
18  fNbBranch = 0;
19  fNbIBranch = 0;
20 }
21 
22 
23 
32 
34  fSelection(selection)
35 {
36  // Create a list of global variables.
37  //
38  // If given, the KVParticleCondition will be used to select the particles to iterate
39  // over in each event in order to calculate the global variables.
40  //
41  // By default, if no selection is defined, only particles which are "OK" (i.e. for which
42  // KVParticle::IsOK() returns true) will be iterated over - this will change in the future!
43 
44  SetName(name);
45  init_KVGVList();
46 }
47 
48 
49 
50 
52 
54 {
55  init_KVGVList();
56  a.Copy(*this);
57 }
58 
59 
60 
65 
66 void KVGVList::Init(void)
67 {
68  // Initialisation of all global variables in list
69  //
70  // As this method may be called several times we ensure that variables are only initialised once
71 
72  this->R__FOR_EACH(KVVarGlob, ListInit)();
73 #ifdef KVGVLIST_OPTIMIZE_GVLIST
74  optimise_variable_lists();
75 #endif
76 }
77 
78 
79 
82 
83 void KVGVList::Reset(void)
84 {
85  // Reset all variables before treating an event
86  this->R__FOR_EACH(KVVarGlob, Reset)();
87  fAbortEventAnalysis = false;
88 }
89 
90 
91 
96 
97 void KVGVList::Fill(const KVNucleus* c)
98 {
99  // Calls KVVarGlob::Fill(KVNucleus*) method of all one-body variables in the list
100  // for all KVNucleus satisfying the KVParticleCondition given to
101  // KVVarGlob::SetSelection() (if no selection given, all nuclei are used).
102 
103  fVG1.R__FOR_EACH(KVVarGlob, Fill)(c);
104 }
105 
106 
107 
108 
113 
114 void KVGVList::Fill2(const KVNucleus* c1, const KVNucleus* c2)
115 {
116  // Calls KVVarGlob::Fill(KVNucleus*,KVNucleus*) method of all two-body variables in the list
117  // for all pairs of KVNucleus (c1,c2) satisfying the KVParticleCondition given to
118  // KVVarGlob::SetSelection() (if no selection given, all nuclei are used).
119 
120  fVG2.R__FOR_EACH(KVVarGlob, Fill2)(c1, c2);
121 }
122 
123 
124 
127 
128 void KVGVList::FillN(const KVEvent* r)
129 {
130  // Calls KVVarGlob::FillN(KVEvent*) method of all N-body variables in the list
131  fVGN.R__FOR_EACH(KVVarGlob, FillN)(r);
132 }
133 
134 
135 
138 
140 {
141  // Calculate all 1-body observables after filling
142  TIter it(&fVG1);
143  KVVarGlob* vg;
144  while ((vg = (KVVarGlob*)it())) {
145  vg->Calculate();
146  if (!vg->TestEventSelection()) {
147  fAbortEventAnalysis = true;
148  break;
149  }
150  }
151 }
152 
153 
154 
157 
159 {
160  // Calculate all 2-body observables after filling
161  TIter it(&fVG2);
162  KVVarGlob* vg;
163  while ((vg = (KVVarGlob*)it())) {
164  vg->Calculate();
165  if (!vg->TestEventSelection()) {
166  fAbortEventAnalysis = true;
167  break;
168  }
169  }
170 }
171 
172 
173 
176 
178 {
179  // Calculate all N-body observables after filling
180  TIter it(&fVGN);
181  KVVarGlob* vg;
182  while ((vg = (KVVarGlob*)it())) {
183  vg->Calculate();
184  if (!vg->TestEventSelection()) {
185  fAbortEventAnalysis = true;
186  break;
187  }
188  }
189 }
190 
191 
192 
205 
207 {
208  // This method will calculate all global variables defined in the list for the event 'e'.
209  //
210  // This will iterate over all particles in the event which correspond to the KVParticleCondition fSelection;
211  // if none has been defined, only particles for which KVParticle::IsOK() returns true will be used.
212  //
213  // Note that for 2-body variables, the Fill2 method will be called for each distinct pair of particles
214  // in the event, plus each pair made up of a particle and itself.
215  //
216  // For each variable for which an event selection condition was set (see KVVarGlob::SetEventSelection())
217  // the condition is tested as soon as the variable is calculated. If the condition is not satisfied,
218  // calculation of the other variables is abandonded and method AbortEventAnalysis() returns kTRUE.
219 
220  Reset();
221 
222 #ifdef KVGVLIST_OPTIMIZE_GVLIST
223  TIter nxt_list(&fOptimList);
224  optigvlist* current_list;
225  while (current_list = (optigvlist*)nxt_list()) {
226  for (auto it = EventIterator(e, fSelection).begin(); it != KVNucleusEvent::Iterator::End(); ++it) {
227  current_list->iterate_one_body([&it](KVVarGlob * vg) {
228  if (vg->IsGlobalVariable())
229  vg->Fill(it.get_const_pointer());
230  }
231  );
232  if (current_list->has_two_body) {
233  for (auto it2 = it; it2 != KVNucleusEvent::Iterator::End(); ++it2) {
234  // we use every distinct pair of particles (including identical pairs) in the event
235  current_list->iterate_two_body([&it, &it2](KVVarGlob * vg) {
236  vg->Fill2(it.get_const_pointer(), it2.get_const_pointer());
237  }
238  );
239  }
240  }
241  }
242  if (current_list->has_N_body) {
243  current_list->iterate_N_body([e](KVVarGlob * vg) {
244  vg->FillN(e);
245  }
246  );
247  }
248  // call calculate for global variables and for KVEventClassifier objects, test event selection
249  current_list->iterate_one_body_with_test([&](KVVarGlob * vg) {
250  vg->Calculate();
251  if ((fAbortEventAnalysis = !vg->TestEventSelection())) {
252  return false;
253  }
254  vg->DefineNewFrame(e);
255  return true;
256  });
257  if (fAbortEventAnalysis) return;
258  if (current_list->has_two_body) {
259  current_list->iterate_two_body_with_test([&](KVVarGlob * vg) {
260  vg->Calculate();
261  if ((fAbortEventAnalysis = !vg->TestEventSelection())) {
262  return false;
263  }
264  vg->DefineNewFrame(e);
265  return true;
266  });
267  if (fAbortEventAnalysis) return;
268  }
269  if (current_list->has_N_body) {
270  current_list->iterate_N_body_with_test([&](KVVarGlob * vg) {
271  vg->Calculate();
272  if ((fAbortEventAnalysis = !vg->TestEventSelection())) {
273  return false;
274  }
275  vg->DefineNewFrame(e);
276  return true;
277  });
278  if (fAbortEventAnalysis) return;
279  }
280  }
281 #else
282  TIter it(this);
283  KVVarGlob* vg;
284 
285  while ((vg = (KVVarGlob*)it())) {
286 
287  if (vg->IsGlobalVariable()) {
288  // call Fill methods for global variables
289  if (vg->IsNBody())
290  vg->FillN(e);
291  else {
292  for (auto it = EventIterator(e, fSelection).begin(); it != KVNucleusEvent::Iterator::End(); ++it) {
293  if (vg->IsTwoBody()) {
294  for (auto it2 = it; it2 != KVNucleusEvent::Iterator::End(); ++it2) {
295  // we use every distinct pair of particles (including identical pairs) in the event
296  vg->Fill2(it.get_const_pointer(), it2.get_const_pointer());
297  }
298  }
299  else
300  vg->Fill(it.get_const_pointer());
301  }
302  }
303  }
304  // call calculate for global variables and for KVEventClassifier objects
305  vg->Calculate();
306  if ((fAbortEventAnalysis = !vg->TestEventSelection())) {
307  return;
308  }
309  vg->DefineNewFrame(e);
310  }
311 #endif
312 }
313 
314 
315 
318 
319 KVVarGlob* KVGVList::GetGV(const Char_t* nom) const
320 {
321  //Return pointer to global variable in list with name 'nom'
322 
323  return (KVVarGlob*) FindObject(nom);
324 }
325 
326 
327 
335 
337 {
338  // Overrides KVUniqueNameList::Add(TObject*) so that global variable pointers are sorted
339  // between the 3 lists used for 1-body, 2-body & N-body variables.
340  //
341  // We also print a warning in case the user tries to add a global variable with the
342  // same name as one already in the list. In the case we retain only the first global variable,
343  // any others with the same name are ignored
344 
345  KVUniqueNameList::Add(obj); // add object to main list, check duplicates
346  if (!ObjectAdded()) {
347  Warning("Add", "You tried to add a global variable with the same name as one already in the list");
348  Warning("Add", "Only the first variable added to the list will be used: name=%s class=%s",
349  GetGV(obj->GetName())->GetName(), GetGV(obj->GetName())->ClassName());
350  Warning("Add", "The following global variable (the one you tried to add) will be ignored:");
351  printf("\n");
352  obj->Print();
353  return;
354  }
355  if (obj->InheritsFrom("KVVarGlob")) {
356  // put global variable pointer in appropriate list
357  KVVarGlob* vg = (KVVarGlob*)obj;
358  if (vg->IsOneBody()) fVG1.Add(vg);
359  else if (vg->IsTwoBody()) fVG2.Add(vg);
360  else if (vg->IsNBody()) fVGN.Add(vg);
361  }
362 }
363 
364 
372 
374 {
375  // Overrides KVUniqueNameList::AddFirst(TObject*) so that global variable pointers are sorted
376  // between the 3 lists used for 1-body, 2-body & N-body variables.
377  //
378  // We also print a warning in case the user tries to add a global variable with the
379  // same name as one already in the list. In the case we retain only the first global variable,
380  // any others with the same name are ignored
381 
382  KVUniqueNameList::AddFirst(obj); // add object to main list, check duplicates
383  if (!ObjectAdded()) {
384  Warning("AddFirst", "You tried to add a global variable with the same name as one already in the list");
385  Warning("AddFirst", "Only the first variable added to the list will be used: name=%s class=%s",
386  GetGV(obj->GetName())->GetName(), GetGV(obj->GetName())->ClassName());
387  Warning("AddFirst", "The following global variable (the one you tried to add) will be ignored:");
388  printf("\n");
389  obj->Print();
390  return;
391  }
392  if (obj->InheritsFrom("KVVarGlob")) {
393  // put global variable pointer in appropriate list
394  KVVarGlob* vg = (KVVarGlob*)obj;
395  if (vg->IsOneBody()) fVG1.Add(vg);
396  else if (vg->IsTwoBody()) fVG2.Add(vg);
397  else if (vg->IsNBody()) fVGN.Add(vg);
398  }
399 }
400 
401 
402 
403 
417 
419 {
420  // Create a branch in the TTree for each global variable in the list.
421  // A leaf with the name of each global variable will be created to hold the
422  // value of the variable (result of KVVarGlob::GetValue() method).
423  // For multi-valued global variables we add a branch for each value with name
424  //
425  // `GVname.ValueName`
426  //
427  // To avoid problems with TTree::Draw, 'GVName' will be sanitized i.e. any
428  // mathematical symbols will be replaced by '_'.
429  //
430  // Any variable for which KVVarGlob::SetMaxNumBranches() was called with argument `0`
431  // will not be added to the TTree.
432 
433  if (!tree) return;
434 
435  // Make sure all variables are initialised before proceeding
436  Init();
437 
438  fNbBranch = 0;
439  fNbIBranch = 0;
440  fBranchVar.clear();
441  fIBranchVar.clear();
442 
443  // calculate the number of variables which will be stored in each vector
444  // the space has to be pre-reserved not added by successive "push_back"s, otherwise references to vector elements
445  // are invalidated every time the vector reallocates space
446 
447  TIter next(this);
448  KVVarGlob* ob;
449  while ((ob = (KVVarGlob*)next())) {
450  if (ob->GetNumberOfBranches()) { //skip variables for which KVVarGlob::SetMaxNumBranches(0) was called
451  if (ob->GetNumberOfValues() > 1) {
452  // multi-valued variable
453  for (int i = 0; i < ob->GetNumberOfBranches(); i++) {
454  if (ob->GetValueType(i) == 'I') ++fNbIBranch;
455  else ++fNbBranch;
456  }
457  }
458  else {
459  if (ob->GetValueType(0) == 'I') ++fNbIBranch;
460  else ++fNbBranch;
461  }
462  }
463  }
464 
465  fBranchVar.reserve(fNbBranch);
466  fIBranchVar.reserve(fNbIBranch);
467  fNbBranch = 0;
468  fNbIBranch = 0;
469 
470  next.Reset();
471 
472  while ((ob = (KVVarGlob*)next())) {
473  if (ob->GetNumberOfBranches()) { //skip variables for which KVVarGlob::SetMaxNumBranches(0) was called
474  TString sane_varname = NameSanitizer(ob->GetName());
475  if (ob->GetNumberOfValues() > 1) {
476  // multi-valued variable
477  for (int i = 0; i < ob->GetNumberOfBranches(); i++) {
478  // replace any nasty mathematical symbols which could pose problems
479  // in names of TTree leaves/branches
480  TString sane_name(ob->GetValueName(i));
481  sane_name.ReplaceAll("*", "star");
482  if (ob->GetValueType(i) == 'I') {
483  tree->Branch(Form("%s.%s", sane_varname.Data(), sane_name.Data()), &fIBranchVar[ fNbIBranch ], Form("%s.%s/I", sane_varname.Data(), sane_name.Data()));
484  ++fNbIBranch;
485  }
486  else {
487  tree->Branch(Form("%s.%s", sane_varname.Data(), sane_name.Data()), &fBranchVar[ fNbBranch ], Form("%s.%s/D", sane_varname.Data(), sane_name.Data()));
488  ++fNbBranch;
489  }
490  }
491  }
492  else {
493  if (ob->GetValueType(0) == 'I') {
494  tree->Branch(sane_varname, &fIBranchVar[ fNbIBranch ], Form("%s/I", sane_varname.Data()));
495  ++fNbIBranch;
496  }
497  else {
498  tree->Branch(sane_varname, &fBranchVar[ fNbBranch ], Form("%s/D", sane_varname.Data()));
499  ++fNbBranch;
500  }
501  }
502  }
503  }
504 }
505 
506 
507 
508 
515 
517 {
518  // Use this method ONLY if you first use MakeBranches(TTree*) in order to
519  // automatically create branches for your global variables.
520  // Call this method for each event in order to put the values of the variables
521  // in the branches ready for TTree::Fill to be called (note that TTree::Fill is not
522  // called in this method: you should do it after this).
523 
524  if (!fNbBranch && !fNbIBranch) return; // MakeBranches has not been called
525 
526  int INT_index = 0;
527  int FLT_index = 0;
528  TIter next(this);
529  KVVarGlob* ob;
530  while ((ob = (KVVarGlob*)next())) {
531  if (ob->GetNumberOfBranches()) { //skip variables for which KVVarGlob::SetMaxNumBranches(0) was called
532 
533  if (ob->GetNumberOfValues() > 1) {
534  // multi-valued variable
535  for (int j = 0; j < ob->GetNumberOfBranches(); j++) {
536  if (ob->GetValueType(j) == 'I') {
537  fIBranchVar[ INT_index ] = (Int_t)ob->GetValue(j);
538  ++INT_index;
539  }
540  else {
541  fBranchVar[ FLT_index ] = ob->GetValue(j);
542  ++FLT_index;
543  }
544  }
545  }
546  else {
547  if (ob->GetValueType(0) == 'I') {
548  fIBranchVar[ INT_index ] = (Int_t)ob->GetValue();
549  ++INT_index;
550  }
551  else {
552  fBranchVar[ FLT_index ] = ob->GetValue();
553  ++FLT_index;
554  }
555  }
556 
557  }
558  }
559 }
560 
561 
562 
611 
613 {
614  // Add an event classification object to the list, based on the named global variable
615  // (which must already be in the list).
616  //
617  // \param[in] varname name of global variable previously added to list
618  // \param[in] value [optional] for multi-valued variables, you can specify which value to use by name, or a mathematical expression
619  // involving one or more of the available values
620  //
621  // Returns a pointer to the object, in order to add either cuts or bins like so:
622  //
623  // #### Using cuts
624  //~~~~{.cpp}
625  // mtot_cuts->AddCut(5);
626  // mtot_cuts->AddCut(10);
627  // mtot_cuts->AddCut(15);
628  // mtot_cuts->AddCut(20);
629  // mtot_cuts->AddCut(25);
630  // mtot_cuts->AddCut(30);
631  //~~~~
632  // This will class events according to:
633  // | mtot | mtot_EC |
634  // |------------|---------|
635  // | <5 | 0 |
636  // | [5, 9] | 1 |
637  // | [10, 14] | 2 |
638  // | [15, 19] | 3 |
639  // | [20, 24] | 4 |
640  // | [25, 29] | 5 |
641  // | >=30 | 6 |
642  //
643  // `mtot_EC` is the default name for an event classification based on `mtot` and will be
644  // used for the branch added to the user's analysis TTree by method MakeBranches().
645  //
646  // #### Using bins
647  // ~~~~{.cpp}
648  // mtot_cuts->AddBin(5,10);
649  // mtot_cuts->AddBin(15,20);
650  // mtot_cuts->AddBin(25,30);
651  //~~~~
652  //This will class events according to the 2 bins with the following numbers:
653  // | mtot | mtot_EC |
654  // |----------|-------------|
655  // | [5, 9] | 0 |
656  // | [15, 19] | 1 |
657  // | [25, 29] | 2 |
658  // | other | -1 |
659  //
660  // Note that in this case, any value outside of a defined bin is unclassified.
661 
662  KVVarGlob* gv = GetGV(varname);
663  if (!gv) {
664  Warning("AddEventClassifier", "Variable %s not found in list. No classification possible.", varname.Data());
665  }
667  Add(ec);
668  return ec;
669 }
670 
671 
672 
704 
705 KVVarGlob* KVGVList::AddGV(const Char_t* class_name, const Char_t* name)
706 {
707  //Add a global variable to the list.
708  //
709  //`"class_name"` must be the name of a valid class inheriting from KVVarGlob, e.g. any of the default global
710  //variable classes defined as part of the standard KaliVeda package (see the [Global Variables module](group__GlobalVariables.html)),
711  //or the name of a user-defined class (see below).
712  //
713  //`"name"` is a unique name for the new global variable object which will be created and added to the
714  //list of global variables. This name can later be used to retrieve the object (see GetGV()).
715  //
716  //Returns pointer to new global variable object in case more than the usual default initialisation is necessary.
717  //
718  //### User-defined global variables
719  //The user may use her own global variables, without having to add them to the main libraries.
720  //If the given class name is not known, it is assumed to be a user-defined class and we attempt to compile and load
721  //the class from the user's source code. For this to work, the user must:
722  //
723  // 1. add to the ROOT macro path the directory where her class's source code is kept, e.g. in `$HOME/.rootrc` add the following line:
724  //
725  //~~~~~~~~~~~~~~~
726  // +Unix.*.Root.MacroPath: $(HOME)/myVarGlobs
727  //~~~~~~~~~~~~~~~
728  //
729  // 2. for each user-defined class, add a line to $HOME/.kvrootrc to define a "plugin". E.g. for a class called MyNewVarGlob,
730  //
731  //~~~~~~~~~~~~~~~
732  // +Plugin.KVVarGlob: MyNewVarGlob MyNewVarGlob MyNewVarGlob.cpp+ "MyNewVarGlob(const char*)"
733  //~~~~~~~~~~~~~~~
734  //
735  // It is assumed that `MyNewVarGlob.h` and `MyNewVarGlob.cpp` will be found in `$HOME/myVarGlobs` (in this example).
736  // The constructor taking a single character string argument (name of the variable) must be defined in the class.
737 
738  KVVarGlob* vg = prepareGVforAdding(class_name, name);
739  if (vg) Add(vg);
740  return vg;
741 }
742 
743 
744 
746 
747 KVVarGlob* KVGVList::prepareGVforAdding(const Char_t* class_name, const Char_t* name)
748 {
749  KVVarGlob* vg = nullptr;
750  TClass* clas = TClass::GetClass(class_name);
751  if (!clas) {
752  //class not in dictionary - user-defined class ? Look for plugin.
753  TPluginHandler* ph = KVBase::LoadPlugin("KVVarGlob", class_name);
754  if (!ph) {
755  //not found
756  Error("AddGV(const Char_t*,const Char_t*)",
757  "Called with class_name=%s.\nClass is unknown: not in standard libraries, and plugin (user-defined class) not found",
758  class_name);
759  return 0;
760  }
761  else {
762  vg = (KVVarGlob*) ph->ExecPlugin(1, name);
763  }
764  }
765  else if (!clas->InheritsFrom("KVVarGlob")) {
766  Error("AddGV(const Char_t*,const Char_t*)",
767  "%s is not a valid class deriving from KVVarGlob.",
768  class_name);
769  return nullptr;
770  }
771  else {
772  // need to use plugins to create even known classes, in order to call ctor with name
773  TPluginHandler* ph = KVBase::LoadPlugin("KVVarGlob", class_name);
774  if (!ph) {
775  // define plugin handler for known class
776  gPluginMgr->AddHandler("KVVarGlob", class_name, class_name, "KVMultiDetglobvars", Form("%s(const char*)", class_name));
777  ph = KVBase::LoadPlugin("KVVarGlob", class_name);
778  }
779  vg = (KVVarGlob*) ph->ExecPlugin(1, name);
780  }
781  return vg;
782 }
783 
784 
785 #ifdef KVGVLIST_OPTIMIZE_GVLIST
786 
793 
794 void KVGVList::optimise_variable_lists()
795 {
796  // Called when list is initialised (after all variables have been added)
797  //
798  // Variables are sorted into sublists so that:
799  // - event-selecting variables are treated together
800  // - any variable defining a new frame is calculated before all following variables
801 
802  if (fOptimList.GetEntries()) return; //already done - don't repeate
803 
804  optigvlist* current_list = new optigvlist;
805  fOptimList.Add(current_list);
806 
807  TIter next_var(this);
808  KVVarGlob* vg, *last_vg(nullptr);
809  while ((vg = (KVVarGlob*)next_var())) {
810  if (last_vg) {
811  if ((last_vg->IsSelectingEvents() && !vg->IsSelectingEvents()) || (last_vg->IsDefiningNewFrame())) {
812  // 1. we have come to the end of a series of event selectors: start a new list
813  // 2. vg's which define new frames have to be calculated before any others
814  current_list = new optigvlist;
815  fOptimList.Add(current_list);
816  }
817  }
818  if (vg->IsNBody()) {
819  current_list->N_body.Add(vg);
820  current_list->has_N_body = true;
821  }
822  else if (vg->IsTwoBody()) {
823  current_list->two_body.Add(vg);
824  current_list->has_two_body = true;
825  }
826  else
827  current_list->one_body.Add(vg);
828  last_vg = vg;
829  }
830 
831  Info("optimise_variable_lists", "After optimization, global variables have been separated into %d lists as follows:", fOptimList.GetEntries());
832  fOptimList.ls();
833 }
834 
835 #endif
836 
837 
842 
843 KVVarGlob* KVGVList::AddGVFirst(const Char_t* class_name, const Char_t* name)
844 {
845  // Add a global variable at the beginning of the list.
846  //
847  // See AddGV() for details.
848 
849  KVVarGlob* vg = prepareGVforAdding(class_name, name);
850  if (vg) AddFirst(vg);
851  return vg;
852 }
853 
854 
int Int_t
ROOT::R::TRInterface & r
#define c(i)
#define e(i)
char Char_t
#define R__FOR_EACH(type, proc)
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
char name[80]
R__EXTERN TPluginManager * gPluginMgr
char * Form(const char *fmt,...)
Class for iterating over nuclei in events accessed through base pointer/reference.
static TPluginHandler * LoadPlugin(const Char_t *base, const Char_t *uri="0")
Definition: KVBase.cpp:793
Simple class for sorting events according to global variables.
Abstract base class container for multi-particle events.
Definition: KVEvent.h:67
#define KVGVLIST_OPTIMIZE_GVLIST
Definition: KVGVList.h:227
void FillBranches()
Definition: KVGVList.cpp:516
KVParticleCondition fSelection
used to select particles to iterate over in event
Definition: KVGVList.h:229
KVVarGlob * AddGVFirst(const Char_t *class_name, const Char_t *name)
Definition: KVGVList.cpp:843
Int_t fNbIBranch
Definition: KVGVList.h:233
void Reset(void)
Reset all variables before treating an event.
Definition: KVGVList.cpp:83
void Init(void)
Definition: KVGVList.cpp:66
std::vector< Double_t > fBranchVar
used for automatic creation & filling of TTree branches
Definition: KVGVList.h:230
TList fVG2
two-body variables
Definition: KVGVList.h:316
TString NameSanitizer(const Char_t *s) const
Definition: KVGVList.h:237
Int_t fNbBranch
Definition: KVGVList.h:232
virtual void AddFirst(TObject *obj)
Definition: KVGVList.cpp:373
void Calculate2()
Calculate all 2-body observables after filling.
Definition: KVGVList.cpp:158
void CalculateN()
Calculate all N-body observables after filling.
Definition: KVGVList.cpp:177
void CalculateGlobalVariables(KVEvent *e)
Definition: KVGVList.cpp:206
KVGVList(const KVString &name="default", const KVParticleCondition &selection=KVParticleCondition{"OK", [](const KVNucleus *n) { return n->IsOK();}})
Definition: KVGVList.cpp:33
bool fAbortEventAnalysis
set to false if a global variable fails its own event selection criterion
Definition: KVGVList.h:234
void FillN(const KVEvent *e)
Calls KVVarGlob::FillN(KVEvent*) method of all N-body variables in the list.
Definition: KVGVList.cpp:128
std::vector< Int_t > fIBranchVar
used for automatic creation & filling of TTree branches
Definition: KVGVList.h:231
KVEventClassifier * AddEventClassifier(const TString &varname, const TString &value="")
Definition: KVGVList.cpp:612
TList fVGN
N-body variables.
Definition: KVGVList.h:317
void Calculate()
Calculate all 1-body observables after filling.
Definition: KVGVList.cpp:139
KVVarGlob * prepareGVforAdding(const Char_t *class_name, const Char_t *name)
Definition: KVGVList.cpp:747
TList fVG1
one-body variables
Definition: KVGVList.h:315
KVVarGlob * AddGV(const Char_t *class_name, const Char_t *name)
Definition: KVGVList.cpp:705
virtual void Add(TObject *obj)
Definition: KVGVList.cpp:336
void Fill(const KVNucleus *c)
Definition: KVGVList.cpp:97
void init_KVGVList(void)
Definition: KVGVList.cpp:16
void MakeBranches(TTree *)
Definition: KVGVList.cpp:418
KVVarGlob * GetGV(const Char_t *nom) const
Return pointer to global variable in list with name 'nom'.
Definition: KVGVList.cpp:319
void Fill2(const KVNucleus *c1, const KVNucleus *c2)
Definition: KVGVList.cpp:114
Description of properties and kinematics of atomic nuclei.
Definition: KVNucleus.h:126
virtual TObject * FindObject(const char *name) const
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition: KVString.h:73
Optimised list in which named objects can only be placed once.
Bool_t ObjectAdded() const
virtual void Add(TObject *obj)
virtual void AddFirst(TObject *obj)
Base class for all global variable implementations.
Definition: KVVarGlob.h:233
virtual void Calculate()=0
Bool_t IsOneBody() const
Definition: KVVarGlob.h:369
virtual Int_t GetNumberOfValues() const
Definition: KVVarGlob.h:638
Double_t GetValue(void) const
Definition: KVVarGlob.h:443
void Fill(const KVNucleus *c)
Definition: KVVarGlob.h:406
virtual void FillN(const KVEvent *)
Definition: KVVarGlob.h:430
Int_t GetNumberOfBranches() const
Definition: KVVarGlob.h:644
bool IsSelectingEvents() const
Definition: KVVarGlob.h:729
void Fill2(const KVNucleus *n1, const KVNucleus *n2)
Definition: KVVarGlob.h:420
virtual Bool_t IsGlobalVariable() const
Definition: KVVarGlob.h:387
virtual TString GetValueName(Int_t i) const
Definition: KVVarGlob.h:654
virtual Char_t GetValueType(Int_t) const
Definition: KVVarGlob.h:668
void DefineNewFrame(KVEvent *e) const
Definition: KVVarGlob.h:745
bool TestEventSelection() const
Definition: KVVarGlob.h:719
Bool_t IsNBody() const
Definition: KVVarGlob.h:381
Bool_t IsTwoBody() const
Definition: KVVarGlob.h:375
void Copy(TObject &arc) const override
static TClass * GetClass(Bool_t load=kTRUE, Bool_t silent=kFALSE)
Bool_t InheritsFrom(const char *cl) const override
void SetName(const char *name)
TIter begin() const
void Reset()
void Add(TObject *obj) override
const char * GetName() const override
virtual const char * GetName() const
virtual const char * ClassName() const
virtual void Warning(const char *method, const char *msgfmt,...) const
virtual Bool_t InheritsFrom(const char *classname) const
virtual void Error(const char *method, const char *msgfmt,...) const
virtual void Print(Option_t *option="") const
virtual void Info(const char *method, const char *msgfmt,...) const
Longptr_t ExecPlugin(int nargs)
void AddHandler(const char *base, const char *regexp, const char *className, const char *pluginName, const char *ctor=nullptr, const char *origin=nullptr)
const char * Data() const
TString & ReplaceAll(const char *s1, const char *s2)
return c1
return c2
TArc a
ClassImp(TPyArg)