KaliVeda
Toolkit for HIC analysis
KVGVList Class Reference

Detailed Description

#define KVGVLIST_OPTIMIZE_GVLIST

Manage a list of global variables.

The KVGVList class handles lists of global variables. Schematic usage is as follows:

KVVGList my_list;
my_list.Add( new SomeVarGlob("var1") ); // add variable using explicit call to constructor
my_list.AddGV("SomeOtherVarGlob", "var2"); // add variable using class name
my_list.Init(); // initialise all variables

Particle selection criteria can be set for a list which will be applied to all variables in the list, using KVParticleCondition (KVTemplateParticleCondition<KVNucleus>):

KVGVList sel_list("only_LCP", {"LCP",[](const KVNucleus* n){ return n->GetZ()<3; }});
sel_list.AddGV("KVZtot", "ztot_lcp"); // will calculate total charge of LCP
sel_list.AddGV("KVEtrans", "et12"); // will calculate total transverse energy of LCP
#define KVGVLIST_OPTIMIZE_GVLIST
Definition: KVGVList.h:227
KVVarGlob * AddGV(const Char_t *class_name, const Char_t *name)
Definition: KVGVList.cpp:705
Description of properties and kinematics of atomic nuclei.
Definition: KVNucleus.h:126
Warning
default selection for KVGVList using default constructor (as above) is to only consider particles for which KVParticle::IsOK() returns true. This is a historical anomaly which can be very troublesome when analysing experimental data. Changing this behaviour would break many existing analysis classes!

Calculation of all global variables in the list in an event loop then proceeds schematically as follows:

while( [loop over events] )
{
my_list.CalculateGlobalVariables( [event] ); // calculate contribution of each particle to each variable
if( !my_list.AbortEventAnalysis() ) // in case cuts for event selection were set - see below
{
auto valueOfvar1 = my_list.GetGV("var1")->GetValue(); // retrieve value of "var1" for event
}
}

1-body, 2-body and \(N\)-body variables are all calculated by KVGVList::CalculateGlobalVariables(). You should be aware that variables are calculated in the same order as they are added to the list with KVGVList::AddGV(). This is important to bear in mind as:

  • any variable for whom event selection criteria are defined (see Event selection) will stop further processing of the list (and set KVGVList::AbortEventAnalysis() to true) as soon as its value is calculated and fails the criteria;
  • any variable VAR2 which uses a kinematical frame defined by another global variable VAR1 in the list (see Definition of new kinematical frames) must be added to the list after VAR1: the new kinematical frame defined by VAR1 will be applied to all particles of the event as soon as VAR1 has been calculated.

Automatic TTree branch creation and filling

One or more branches can be added to a TTree and filled with the values of all or some of the global variables in the list using method KVGVList::MakeBranches(TTree*). For each single-valued variable a branch with name

[varName]

will be added. For multi-valued variables (those for which KVVarGlob::GetNumberOfValues() returns >1) we add

[varName].[value0_name]
[varName].[value1_name]
...
[varName].[valueN_name]

branches for each named value defined for the variable (see KVFlowTensor for an example of a multi-valued variable). The total number of branches added can be modified by KVVarGlob::SetMaxNumBranches(): calling this method with argument 0 will prevent any branch being added for the variable in question:

l.AddGV("KVFlowTensor","tensor");
l.GetGV("tensor")->GetNumberOfValues()
(int) 9 // KVFlowTensor calculates 9 different values!
l.GetGV("tensor")->GetNumberOfBranches()
(int) 3 // by default, the number of branches is limited to 3,
/ corresponding to the first 3 values:
/ FlowAngle, Sphericity, Coplanarity
TTree a_tree;
l.MakeBranches(&a_tree);
a_tree.GetListOfBranches()->ls();
OBJ: TObjArray TObjArray An array of objects : 0
OBJ: TBranch tensor.FlowAngle tensor.FlowAngle/D : 0 at: 0x562ae3f73270
OBJ: TBranch tensor.Sphericity tensor.Sphericity/D : 0 at: 0x562ae4076000
OBJ: TBranch tensor.Coplanarity tensor.Coplanarity/D : 0 at: 0x562ae3f5ca90
l.GetGV("tensor")->SetMaxNumBranches(9); // add all values to TTree
TTree another_tree;
l.MakeBranches(&another_tree);
another_tree.GetListOfBranches()->ls();
OBJ: TObjArray TObjArray An array of objects : 0
OBJ: TBranch tensor.FlowAngle tensor.FlowAngle/D : 0 at: 0x562ae41fda50
OBJ: TBranch tensor.Sphericity tensor.Sphericity/D : 0 at: 0x562ae4277a60
OBJ: TBranch tensor.Coplanarity tensor.Coplanarity/D : 0 at: 0x562ae412f270
OBJ: TBranch tensor.NumberParts tensor.NumberParts/I : 0 at: 0x562ae41fb700
OBJ: TBranch tensor.KinFlowRatio13 tensor.KinFlowRatio13/D : 0 at: 0x562ae41cc0d0
OBJ: TBranch tensor.KinFlowRatio23 tensor.KinFlowRatio23/D : 0 at: 0x562ae4285a10
OBJ: TBranch tensor.PhiReacPlane tensor.PhiReacPlane/D : 0 at: 0x562ae41e2d20
OBJ: TBranch tensor.SqueezeAngle tensor.SqueezeAngle/D : 0 at: 0x562ae41d7c20
OBJ: TBranch tensor.SqueezeRatio tensor.SqueezeRatio/D : 0 at: 0x562ae426c4e0
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
virtual Int_t GetNumberOfValues() const
Definition: KVVarGlob.h:638
void SetMaxNumBranches(Int_t n)
Definition: KVVarGlob.h:683
Int_t GetNumberOfBranches() const
Definition: KVVarGlob.h:644

The event-by-event treatment including filling a TTree with the values of all global variables for each event then becomes particularly simple:

KVVGList my_list;
my_list.AddGV("SomeVarGlob", "var1");
my_list.AddGV("SomeOtherVarGlob", "var2");
TTree my_tree;
my_list.MakeBranches(&my_tree); // Note: implicit call to KVGVList::Init()
while( [loop over events] )
{
my_list.CalculateGlobalVariables( [event] ); // calculate contribution of each particle to each variable
if( !my_list.AbortEventAnalysis() ) // in case cuts for event selection were set - see below
{
my_list.FillBranches(); // read values of all variables
my_tree.Fill(); // write data in TTree
}
}

Event selection

Conditions ('cuts') can be set on each variable which decide whether or not to retain an event for analysis. If any variable in the list fails the test, processing of the list is immediately abandoned.

Selection criteria are set using lambda expressions. In this example, the variable "mtot" must have a value of at least 4 for the event to be retained:

KVGVList vglist;
auto mtot = vglist.AddGV("KVMult","mtot");
mtot->SetEventSelection([](const KVVarGlob* v){ return v->GetValue()>=4; });
Base class for all global variable implementations.
Definition: KVVarGlob.h:233
Double_t GetValue(void) const
Definition: KVVarGlob.h:443
void SetEventSelection(const EventSelector &f)
Definition: KVVarGlob.h:699

(the base class KVVarGlob pointer passed as argument to the lambda is the global variable itself, i.e. mtot in the present example).

Any event selection criterion is tested as soon as each variable has been calculated. If the test fails, no further variables are calculated and the KVGVList goes into 'abort event' mode:

KVEvent* event_to_analyse;
vglist.CalculateGlobalVariables( event_to_analyse );
if( !vglist.AbortEventAnalysis() )
{
}
Abstract base class container for multi-particle events.
Definition: KVEvent.h:67
void CalculateGlobalVariables(KVEvent *e)
Definition: KVGVList.cpp:206
bool AbortEventAnalysis() const
Definition: KVGVList.h:372

Definition of new kinematical frames

Variables in the list can be used to define new kinematical frames which can in turn be used by any variables which occur after them in the list. In order to do so, call method KVVarGlob::SetNewFrameDefinition() for the variable(s) in question with a lambda function having the following signature:

[](KVEvent* e, const KVVarGlob* vg){ e->SetFrame("_frame_name_", ...); }
virtual void SetFrame(const Char_t *, const KVFrameTransform &)=0

When called (e.g. by KVGVList::CalculateGlobalVariables()), the KVVarGlob pointer gives access to the global variable used to define the new frame.

As an example of use, imagine that KVZmax is used to find the heaviest (largest Z) fragment in the forward CM hemisphere, then the velocity of this fragment is used to define a "QP_FRAME" in order to calculate the KVFlowTensor in this frame:

KVGVList vglist;
auto vg = vglist.AddGV("KVZmax", "zmax");
vg->SetFrame("CM");
vg->SetSelection( {"V>0", [](const KVNucleus* n){ return n->GetVpar()>0; }} );
vg->SetNewFrameDefinition(
[](KVEvent* e, const KVVarGlob* v){
e->SetFrame("QP_FRAME", static_cast<const KVZmax*>(v)->GetZmax(0)->GetVelocity());
});
vg = AddGV("KVFlowTensor", "qp_tensor");
vg->SetFrame("QP_FRAME"); // frame will have been defined before tensor is filled
void SetFrame(const Char_t *ref)
Definition: KVVarGlob.h:505
Global variable used to sort particles in order of decreasing atomic number
Definition: KVZmax.h:35

Event sorting

Event classifier objects can be defined for any global variable in the list using method KVGVList::AddEventClassifier():

KVGVList vglist;
vglist.AddGV("KVMult", "mtot");
auto mt_cuts = vglist.AddEventClassifier("mtot");
KVEventClassifier * AddEventClassifier(const TString &varname, const TString &value="")
Definition: KVGVList.cpp:612

which returns a pointer to the event classifier object (KVEventClassifier), in order to set up either cuts or bins.

mtot_EC is the default name for an event classification based on mtot and will be used for the branch added to the user's analysis TTree by method KVGVList::MakeBranches().

Examples
globvars_kvzmax.C.

Definition at line 227 of file KVGVList.h.

#include <KVGVList.h>

Inheritance diagram for KVGVList:

Public Member Functions

 KVGVList (const KVGVList &a)
 
 KVGVList (const KVString &name="default", const KVParticleCondition &selection=KVParticleCondition{"OK", [](const KVNucleus *n) { return n->IsOK();}})
 
bool AbortEventAnalysis () const
 
virtual void Add (TObject *obj)
 
KVEventClassifierAddEventClassifier (const TString &varname, const TString &value="")
 
virtual void AddFirst (TObject *obj)
 
KVVarGlobAddGV (const Char_t *class_name, const Char_t *name)
 
KVVarGlobAddGVFirst (const Char_t *class_name, const Char_t *name)
 
void CalculateGlobalVariables (KVEvent *e)
 
void FillBranches ()
 
KVVarGlobGetGV (const Char_t *nom) const
 Return pointer to global variable in list with name 'nom'. More...
 
KVVarGlobGetGVType (const Char_t *class_name)
 Returns first global variable in list with given class. More...
 
Bool_t Has1BodyVariables ()
 returns kTRUE if list contains 1-body variables More...
 
Bool_t Has2BodyVariables ()
 returns kTRUE if list contains 2-body variables More...
 
Bool_t HasNBodyVariables ()
 returns kTRUE if list contains N-body variables More...
 
void Init (void)
 
void MakeBranches (TTree *)
 
void Reset (void)
 Reset all variables before treating an event. More...
 
- Public Member Functions inherited from KVUniqueNameList
 KVUniqueNameList (Bool_t R=kFALSE)
 Default constructor. More...
 
virtual ~KVUniqueNameList ()
 Destructor. More...
 
virtual void AddAfter (const TObject *after, TObject *obj)
 
virtual void AddAt (TObject *obj, Int_t idx)
 
virtual void AddBefore (const TObject *before, TObject *obj)
 
virtual void AddLast (TObject *obj)
 
Bool_t ObjectAdded () const
 
void ReplaceObjects (Bool_t yes=kTRUE)
 
- Public Member Functions inherited from KVHashList
 KVHashList (Int_t capacity=TCollection::kInitHashTableCapacity, Int_t rehash=2)
 
virtual ~KVHashList ()
 Destructor. More...
 
Float_t AverageCollisions () const
 
template<typename T >
Bool_t ContainsObjectWithName (const T &o)
 
const TList * GetListForObject (const char *name) const
 
const TList * GetListForObject (const TObject *obj) const
 
void Rehash (Int_t newCapacity=0)
 
void Sort (Bool_t order=kSortAscending)
 
- Public Member Functions inherited from KVSeqCollection
 KVSeqCollection ()
 Default constructor. More...
 
 KVSeqCollection (const Char_t *collection_classname)
 
 KVSeqCollection (const KVSeqCollection &)
 
virtual ~KVSeqCollection ()
 
virtual TObject * After (const TObject *obj) const
 
virtual TObject * At (Int_t idx) const
 
virtual TObject * Before (const TObject *obj) const
 
virtual void Clear (Option_t *option="")
 
const Char_t * CollectionClassName () const
 
virtual void Copy (TObject &obj) const
 
virtual void Delete (Option_t *option="")
 
virtual void Execute (const char *method, const char *params, Int_t *error=0)
 
virtual void Execute (TMethod *method, TObjArray *params, Int_t *error=0)
 
virtual TObject * FindObject (const char *name) const
 
virtual TObject * FindObject (const TObject *obj) const
 
virtual TObject * FindObjectAny (const Char_t *att, const Char_t *keys, Bool_t contains_all=kFALSE, Bool_t case_sensitive=kTRUE) const
 
TObject * FindObjectByClass (const Char_t *) const
 Return (first) object in embedded list with given class. More...
 
TObject * FindObjectByClass (const TClass *) const
 Return (first) object in embedded list with given class. More...
 
virtual TObject * FindObjectByLabel (const Char_t *) const
 
virtual TObject * FindObjectByName (const Char_t *name) const
 
virtual TObject * FindObjectByNumber (UInt_t num) const
 
virtual TObject * FindObjectByTitle (const Char_t *) const
 Will return object with given title (value of TObject::GetTitle() method). More...
 
virtual TObject * FindObjectByType (const Char_t *) const
 
virtual TObject * FindObjectWithMethod (const Char_t *retvalue, const Char_t *method) const
 
virtual TObject * FindObjectWithNameAndType (const Char_t *name, const Char_t *type) const
 
virtual TObject * First () const
 
template<typename T >
T * get_object (const TString &name) const
 
TSeqCollection * GetCollection () const
 
virtual TObject ** GetObjectRef (const TObject *obj) const
 Return reference to object. More...
 
virtual Int_t GetSize () const
 
KVSeqCollectionGetSubListWithClass (const Char_t *class_name) const
 
KVSeqCollectionGetSubListWithClass (const TClass *_class) const
 
KVSeqCollectionGetSubListWithLabel (const Char_t *retvalue) const
 
KVSeqCollectionGetSubListWithMethod (const Char_t *retvalue, const Char_t *method) const
 
KVSeqCollectionGetSubListWithName (const Char_t *retvalue) const
 
KVSeqCollectionGetSubListWithType (const Char_t *retvalue) const
 
virtual Bool_t IsCleanup () const
 
virtual Bool_t IsSendingModifiedSignals () const
 
virtual Bool_t IsSortable () const
 
virtual Bool_t IsSorted () const
 
virtual TObject * Last () const
 
virtual TIterator * MakeIterator (Bool_t dir=kIterForward) const
 Make and return iterator for the list. More...
 
virtual void Modified ()
 
KVSeqCollectionoperator= (const KVSeqCollection &)
 
virtual void RecursiveRemove (TObject *obj)
 
virtual TObject * Remove (TObject *obj)
 Remove object from list. More...
 
virtual void SendModifiedSignals (Bool_t yes=kTRUE)
 
virtual void SetCleanup (Bool_t enable=kTRUE)
 
virtual void SetOwner (Bool_t enable=kTRUE)
 

Private Member Functions

void Calculate ()
 Calculate all 1-body observables after filling. More...
 
void Calculate2 ()
 Calculate all 2-body observables after filling. More...
 
void CalculateN ()
 Calculate all N-body observables after filling. More...
 
void Fill (const KVNucleus *c)
 
void Fill2 (const KVNucleus *c1, const KVNucleus *c2)
 
void FillN (const KVEvent *e)
 Calls KVVarGlob::FillN(KVEvent*) method of all N-body variables in the list. More...
 
void init_KVGVList (void)
 
TString NameSanitizer (const Char_t *s) const
 
KVVarGlobprepareGVforAdding (const Char_t *class_name, const Char_t *name)
 

Private Attributes

bool fAbortEventAnalysis
 set to false if a global variable fails its own event selection criterion More...
 
std::vector< Double_t > fBranchVar
 used for automatic creation & filling of TTree branches More...
 
std::vector< Int_t > fIBranchVar
 used for automatic creation & filling of TTree branches More...
 
Int_t fNbBranch
 
Int_t fNbIBranch
 
KVParticleCondition fSelection
 used to select particles to iterate over in event More...
 
TList fVG1
 one-body variables More...
 
TList fVG2
 two-body variables More...
 
TList fVGN
 N-body variables. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from KVSeqCollection
static KVSeqCollectionMakeListFromFile (TFile *file)
 
static KVSeqCollectionMakeListFromFileWithClass (TFile *file, const Char_t *class_name)
 
static KVSeqCollectionMakeListFromFileWithClass (TFile *file, const TClass *_class)
 
static KVSeqCollectionMakeListFromFileWithMethod (TFile *file, const Char_t *retvalue, const Char_t *method)
 
static void RehashCleanupList ()
 

Constructor & Destructor Documentation

◆ KVGVList() [1/2]

KVGVList::KVGVList ( const KVString name = "default",
const KVParticleCondition selection = KVParticleCondition{"OK", [](const KVNucleus * n) { return n->IsOK(); }} 
)

Create a list of global variables.

If given, the KVParticleCondition will be used to select the particles to iterate over in each event in order to calculate the global variables.

By default, if no selection is defined, only particles which are "OK" (i.e. for which KVParticle::IsOK() returns true) will be iterated over - this will change in the future!

Definition at line 33 of file KVGVList.cpp.

◆ KVGVList() [2/2]

KVGVList::KVGVList ( const KVGVList a)

Definition at line 53 of file KVGVList.cpp.

Member Function Documentation

◆ AbortEventAnalysis()

bool KVGVList::AbortEventAnalysis ( ) const
inline

Definition at line 372 of file KVGVList.h.

◆ Add()

void KVGVList::Add ( TObject *  obj)
virtual

Overrides KVUniqueNameList::Add(TObject*) so that global variable pointers are sorted between the 3 lists used for 1-body, 2-body & N-body variables.

We also print a warning in case the user tries to add a global variable with the same name as one already in the list. In the case we retain only the first global variable, any others with the same name are ignored

Reimplemented from KVUniqueNameList.

Examples
globvars_kvzmax.C.

Definition at line 336 of file KVGVList.cpp.

◆ AddEventClassifier()

KVEventClassifier * KVGVList::AddEventClassifier ( const TString &  varname,
const TString &  value = "" 
)

Add an event classification object to the list, based on the named global variable (which must already be in the list).

Parameters
[in]varnamename of global variable previously added to list
[in]value[optional] for multi-valued variables, you can specify which value to use by name, or a mathematical expression involving one or more of the available values

Returns a pointer to the object, in order to add either cuts or bins like so:

Using cuts

mtot_cuts->AddCut(5);
mtot_cuts->AddCut(10);
mtot_cuts->AddCut(15);
mtot_cuts->AddCut(20);
mtot_cuts->AddCut(25);
mtot_cuts->AddCut(30);

This will class events according to:

mtot mtot_EC
<5 0
[5, 9] 1
[10, 14] 2
[15, 19] 3
[20, 24] 4
[25, 29] 5
>=30 6

mtot_EC is the default name for an event classification based on mtot and will be used for the branch added to the user's analysis TTree by method MakeBranches().

Using bins

mtot_cuts->AddBin(5,10);
mtot_cuts->AddBin(15,20);
mtot_cuts->AddBin(25,30);

This will class events according to the 2 bins with the following numbers:

mtot mtot_EC
[5, 9] 0
[15, 19] 1
[25, 29] 2
other -1

Note that in this case, any value outside of a defined bin is unclassified.

Definition at line 612 of file KVGVList.cpp.

◆ AddFirst()

void KVGVList::AddFirst ( TObject *  obj)
virtual

Overrides KVUniqueNameList::AddFirst(TObject*) so that global variable pointers are sorted between the 3 lists used for 1-body, 2-body & N-body variables.

We also print a warning in case the user tries to add a global variable with the same name as one already in the list. In the case we retain only the first global variable, any others with the same name are ignored

Reimplemented from KVUniqueNameList.

Definition at line 373 of file KVGVList.cpp.

◆ AddGV()

KVVarGlob * KVGVList::AddGV ( const Char_t *  class_name,
const Char_t *  name 
)

Add a global variable to the list.

"class_name" must be the name of a valid class inheriting from KVVarGlob, e.g. any of the default global variable classes defined as part of the standard KaliVeda package (see the Global Variables module), or the name of a user-defined class (see below).

"name" is a unique name for the new global variable object which will be created and added to the list of global variables. This name can later be used to retrieve the object (see GetGV()).

Returns pointer to new global variable object in case more than the usual default initialisation is necessary.

User-defined global variables

The user may use her own global variables, without having to add them to the main libraries. If the given class name is not known, it is assumed to be a user-defined class and we attempt to compile and load the class from the user's source code. For this to work, the user must:

  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:
+Unix.*.Root.MacroPath: $(HOME)/myVarGlobs
  1. for each user-defined class, add a line to $HOME/.kvrootrc to define a "plugin". E.g. for a class called MyNewVarGlob,
+Plugin.KVVarGlob: MyNewVarGlob MyNewVarGlob MyNewVarGlob.cpp+ "MyNewVarGlob(const char*)"

It is assumed that MyNewVarGlob.h and MyNewVarGlob.cpp will be found in $HOME/myVarGlobs (in this example). The constructor taking a single character string argument (name of the variable) must be defined in the class.

Definition at line 705 of file KVGVList.cpp.

◆ AddGVFirst()

KVVarGlob * KVGVList::AddGVFirst ( const Char_t *  class_name,
const Char_t *  name 
)

Add a global variable at the beginning of the list.

See AddGV() for details.

Definition at line 843 of file KVGVList.cpp.

◆ Calculate()

void KVGVList::Calculate ( void  )
private

Calculate all 1-body observables after filling.

Definition at line 139 of file KVGVList.cpp.

◆ Calculate2()

void KVGVList::Calculate2 ( )
private

Calculate all 2-body observables after filling.

Definition at line 158 of file KVGVList.cpp.

◆ CalculateGlobalVariables()

void KVGVList::CalculateGlobalVariables ( KVEvent e)

This method will calculate all global variables defined in the list for the event 'e'.

This will iterate over all particles in the event which correspond to the KVParticleCondition fSelection; if none has been defined, only particles for which KVParticle::IsOK() returns true will be used.

Note that for 2-body variables, the Fill2 method will be called for each distinct pair of particles in the event, plus each pair made up of a particle and itself.

For each variable for which an event selection condition was set (see KVVarGlob::SetEventSelection()) the condition is tested as soon as the variable is calculated. If the condition is not satisfied, calculation of the other variables is abandonded and method AbortEventAnalysis() returns kTRUE.

Examples
globvars_kvzmax.C.

Definition at line 206 of file KVGVList.cpp.

◆ CalculateN()

void KVGVList::CalculateN ( )
private

Calculate all N-body observables after filling.

Definition at line 177 of file KVGVList.cpp.

◆ Fill()

void KVGVList::Fill ( const KVNucleus c)
private

Calls KVVarGlob::Fill(KVNucleus*) method of all one-body variables in the list for all KVNucleus satisfying the KVParticleCondition given to KVVarGlob::SetSelection() (if no selection given, all nuclei are used).

Definition at line 97 of file KVGVList.cpp.

◆ Fill2()

void KVGVList::Fill2 ( const KVNucleus c1,
const KVNucleus c2 
)
private

Calls KVVarGlob::Fill(KVNucleus*,KVNucleus*) method of all two-body variables in the list for all pairs of KVNucleus (c1,c2) satisfying the KVParticleCondition given to KVVarGlob::SetSelection() (if no selection given, all nuclei are used).

Definition at line 114 of file KVGVList.cpp.

◆ FillBranches()

void KVGVList::FillBranches ( )

Use this method ONLY if you first use MakeBranches(TTree*) in order to automatically create branches for your global variables. Call this method for each event in order to put the values of the variables in the branches ready for TTree::Fill to be called (note that TTree::Fill is not called in this method: you should do it after this).

Examples
globvars_kvzmax.C.

Definition at line 516 of file KVGVList.cpp.

◆ FillN()

void KVGVList::FillN ( const KVEvent e)
private

Calls KVVarGlob::FillN(KVEvent*) method of all N-body variables in the list.

Definition at line 128 of file KVGVList.cpp.

◆ GetGV()

KVVarGlob * KVGVList::GetGV ( const Char_t *  nom) const

Return pointer to global variable in list with name 'nom'.

Definition at line 319 of file KVGVList.cpp.

◆ GetGVType()

KVVarGlob* KVGVList::GetGVType ( const Char_t *  class_name)
inline

Returns first global variable in list with given class.

Definition at line 344 of file KVGVList.h.

◆ Has1BodyVariables()

Bool_t KVGVList::Has1BodyVariables ( )
inline

returns kTRUE if list contains 1-body variables

Definition at line 352 of file KVGVList.h.

◆ Has2BodyVariables()

Bool_t KVGVList::Has2BodyVariables ( )
inline

returns kTRUE if list contains 2-body variables

Definition at line 358 of file KVGVList.h.

◆ HasNBodyVariables()

Bool_t KVGVList::HasNBodyVariables ( )
inline

returns kTRUE if list contains N-body variables

Definition at line 364 of file KVGVList.h.

◆ Init()

void KVGVList::Init ( void  )

Initialisation of all global variables in list

As this method may be called several times we ensure that variables are only initialised once

Examples
globvars_kvzmax.C.

Definition at line 66 of file KVGVList.cpp.

◆ init_KVGVList()

void KVGVList::init_KVGVList ( void  )
private

Definition at line 16 of file KVGVList.cpp.

◆ MakeBranches()

void KVGVList::MakeBranches ( TTree *  tree)

Create a branch in the TTree for each global variable in the list. A leaf with the name of each global variable will be created to hold the value of the variable (result of KVVarGlob::GetValue() method). For multi-valued global variables we add a branch for each value with name

GVname.ValueName

To avoid problems with TTree::Draw, 'GVName' will be sanitized i.e. any mathematical symbols will be replaced by '_'.

Any variable for which KVVarGlob::SetMaxNumBranches() was called with argument 0 will not be added to the TTree.

Examples
globvars_kvzmax.C.

Definition at line 418 of file KVGVList.cpp.

◆ NameSanitizer()

TString KVGVList::NameSanitizer ( const Char_t *  s) const
inlineprivate

Definition at line 237 of file KVGVList.h.

◆ prepareGVforAdding()

KVVarGlob * KVGVList::prepareGVforAdding ( const Char_t *  class_name,
const Char_t *  name 
)
private

Definition at line 747 of file KVGVList.cpp.

◆ Reset()

void KVGVList::Reset ( void  )

Reset all variables before treating an event.

Examples
globvars_kvzmax.C.

Definition at line 83 of file KVGVList.cpp.

Member Data Documentation

◆ fAbortEventAnalysis

bool KVGVList::fAbortEventAnalysis
private

set to false if a global variable fails its own event selection criterion

Definition at line 234 of file KVGVList.h.

◆ fBranchVar

std::vector<Double_t> KVGVList::fBranchVar
private

used for automatic creation & filling of TTree branches

Definition at line 230 of file KVGVList.h.

◆ fIBranchVar

std::vector<Int_t> KVGVList::fIBranchVar
private

used for automatic creation & filling of TTree branches

Definition at line 231 of file KVGVList.h.

◆ fNbBranch

Int_t KVGVList::fNbBranch
private

Definition at line 232 of file KVGVList.h.

◆ fNbIBranch

Int_t KVGVList::fNbIBranch
private

Definition at line 233 of file KVGVList.h.

◆ fSelection

KVParticleCondition KVGVList::fSelection
private

used to select particles to iterate over in event

Definition at line 229 of file KVGVList.h.

◆ fVG1

TList KVGVList::fVG1
private

one-body variables

Definition at line 315 of file KVGVList.h.

◆ fVG2

TList KVGVList::fVG2
private

two-body variables

Definition at line 316 of file KVGVList.h.

◆ fVGN

TList KVGVList::fVGN
private

N-body variables.

Definition at line 317 of file KVGVList.h.