KaliVeda
Toolkit for HIC analysis
Loading...
Searching...
No Matches
Global Variables

A global variable is an analysis tool for condensing the information in a multibody event into one or a few characteristic values. A simple example is the event multiplicity (the number of particles in each event), which can be used to characterize heavy-ion collision events in terms of violence or centrality. The global variable classes below can be used with any event described by a class derived from KVEvent, containing particles described by a class which inherits from KVNucleus.

There are several base classes used for implementing global variables:

  • KVVarGlob implements the basic interface;
  • KVVarGlob1 is a specialised base class for global variables which calculate a single representative value;
  • KVVarGlobMean is a further specialisation for variables which calculate the mean value and variance of some quantity;
  • KVVGSum is a very flexible class which can be easily configured to calculate user-defined global variables;
  • KVVGObjectSum templated class calculating sums of objects;
  • KVVGVectorSum templated class calculating sums of TVector3 physics vectors;
  • KVGVList handles lists of global variables.

Many specific implementations are provided - see the list here: Global Variables.

Using global variables

Each global variable class implements the calculation of a given observable for any set of nuclei; this is done independently of the selection of the nuclei in the set, or in what reference frame kinematics should be calculated for either selection purposes or for the calculation of the observable. The way the observable is calculated is defined by the implementation of the KVVarGlob::fill(const KVNucleus*) (for 1-body observables: see below) and KVVarGlob::Calculate() methods.

Global variables can be of different types:

Typical usage (pseudo-code):

SomeVarGlob VG; // KVVarGlob daughter class implementing 1-body global variable
VG.SetSelection( [particle selection criteria] );
VG.SetFrame( [reference frame for kinematics] );
VG.Init();
// Analysis loop
while( [loop over events] )
{
while( [loop over particles in event] )
{
// calculate contribution of particle to variable, only if it satisfies
// the selection criteria set with method SetSelection() (if any)
VG.Fill( [particle] );
}
VG.Calculate(); // perform any necessary calculations
// accessing the result
auto vg = VG.GetValue(); // retrieve unique or principal value of variable
auto vh = VG.GetValue(2); // retrieve value with index 2 from multi-valued variable (implementation-dependent)
auto vi = VG.GetValue("Toto"); // retrieve value named "Toto" from multi-valued variable (implementation-dependent)
auto vv = VG.GetValueVector(); // retrieve vector containing all values of multi-valued variable
for(auto V : vv) // print all values
{
std::cout << V << std::endl;
}
VG.Reset(); // reinitialise prior to analysis of next event
}
winID h TVirtualViewer3D vv

Note that although the KVVarGlob::Fill() method is called for all particles, only those which satsify the conditions given to KVVarGlob::SetSelection() will be used to calculate the global variable.

Example of use

Here is a more realistic example which you can follow by entering each line of code on the kaliveda command line. Consider an event containing different nuclei, which we can build up like so (see Nuclei and Events):

e.AddNucleus()->SetZ(12);
auto n = e.AddNucleus();
n->SetZAandE(1,2,136.4);
n->SetTheta(65.4);
e.ls();
KVEvent with 2 particles :
------------------------------------
KVNucleus Z=12 A=24 E*=0
KVParticle mass=22335.8 Theta=0 Phi=0 KE=0 Vpar=0
KVNucleus Z=1 A=2 E*=0
KVParticle mass=1875.61 Theta=65.4 Phi=0 KE=136.4 Vpar=4.51675
#define e(i)
Abstract base class container for multi-particle events.
Definition KVEvent.h:67
An event container for KVNucleus objects.
Description of properties and kinematics of atomic nuclei.
Definition KVNucleus.h:126
Base class for relativistic kinematics of massive particles.
Definition KVParticle.h:396
const Int_t n
constexpr Double_t E()

Now we'll characterise our event using two global variables: the largest atomic number \(Z\) of all nuclei in the event, and the total transverse energy of light charged particles (LCP: \(Z\leq 2\)):

// declare and initialise global variables
KVZmax zmax;
zmax.Init();
et12.Init();
Sum of transverse kinetic energies of LCP ( ) in event.
Definition KVEtransLCP.h:17
void Init()
Definition KVVGSum.cpp:52
Global variable used to sort particles in order of decreasing atomic number
Definition KVZmax.h:35
void Init()
Definition KVZmax.h:63

At this stage we can print some informations on the configuration of each variable, for example:

et12.Print();
GLOBVAR: KVVarGlobMean class=KVEtransLCP [One-body variable]
- Single-valued variable:
- Kinematics calculated in default frame
- List of options:
KVNameValueList::Options : (0x7f740ece9180)
<mode=sum>
<method=GetEtran>
- Particle selection criteria:
Z<3
- Number of branches in TTree: 1
Abstract base class for global variables which calculate a mean value.
void Print(Option_t *="") const

Here we can see that KVEtransLCP works by summing the result of KVParticle::GetEtran() (non-relativistic transverse energy) for each particle in the event selected according to \(Z<3\), i.e. LCP. It is in fact a specialised version of KVEtrans (which calculates the total transverse energy of all particles, without selection), which itself is a specialised version of KVVGSum, a generic global variable for calculating multiplicities, summed or averaged quantities.

Now in order to fill our variables, we perform a loop over the nuclei in the event:

for(auto& n : e){
zmax.Fill(&n);
et12.Fill(&n);
}
void Fill(const KVNucleus *c)
Definition KVVarGlob.h:406

and after calculating the final values we find:

et12.Calculate();
et12.GetValue() // no ';' in order to see value in kaliveda terminal!
(double) 112.76329
// check total transverse energy of LCP
e.GetParticle(2)->GetEtran()
(double) 112.76329
zmax.Calculate();
zmax.GetValue()
(double) 12.000000
// KVZmax also provides access to event nuclei sorted by Z
zmax.GetZmax(0)->Print();
KVNucleus Z=12 A=24 E*=0
KVParticle mass=22335.8 Theta=0 Phi=0 KE=0 Vpar=0
virtual void Print(Option_t *t="") const
Display nucleus parameters.
Double_t GetValue(void) const
Definition KVVarGlob.h:443
KVNucleus * GetZmax(Int_t i) const
Definition KVZmax.h:77
void Calculate()
Definition KVZmax.h:67

As we will see below, analyses involving several or many global variables are made much easier using a list of global variables, KVGVList.

Options, parameters, reference frames, particle selection, etc.

Particle selection

The selection of particles which are taken into account is handled by the variable itself by calling method KVVarGlob::SetSelection(). This method takes as argument an object of class KVParticleCondition (typedef for KVTemplateParticleCondition<KVNucleus>):

KVZtot ztot;
ztot.SetSelection({"60<tcm<90", [](const KVNucleus* n){ return KVValueRange<double>(60,90).Contains(n->GetTheta()); }});
Range of values specified by minimum, maximum.
Bool_t Contains(ValueType V) const
void SetSelection(const KVParticleCondition &sel)
Definition KVVarGlob.h:600
Sum of atomic numbers in event.
Definition KVZtot.h:17

This variable will now calculate the total charge \(Z_{tot}\) of all nuclei with polar angles \(60^o\leq\theta\leq90^o\).

Kinematical reference frames

To change the reference frame used by the variable to calculate kinematical properties of particles (including those used for particle selection), call method KVVarGlob::SetFrame() (see KVTemplateEvent::SetFrame() and KVParticle::GetFrame() for how to define and access different frames):

KVZtot ztot;
ztot.SetFrame("cm");
ztot.SetSelection({"60<tcm<90", [](const KVNucleus* n){ return KVValueRange<double>(60,90).Contains(n->GetTheta()); }});
void SetFrame(const Char_t *ref)
Definition KVVarGlob.h:505

Assuming that a kinematical reference frame called CM has been defined for the particles of the event when the global variable is filled, this variable will now calculate the total charge \(Z_{tot}\) of all nuclei with polar angles in the CM frame \(60^o\leq\theta_{CM}\leq90^o\)

Options and parameters

In order to give greater flexibility to global variable classes without the need to add member variables and the associated Get/Set methods, we provide methods to handle generic 'options' and 'parameters' for all variables.

An 'option' is a name-value pair, the value is a character string. Methods to use are:

void SetOption(const Char_t* option, const Char_t* value)
Bool_t IsOptionGiven(const Char_t* option)
KVString& GetOptionString(const Char_t* option) const
void UnsetOption(const Char_t* opt)
bool Bool_t
char Char_t
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition KVString.h:73

A 'parameter' is a name-value pair, the value is a double-precision float value. Methods to use are:

void SetParameter(const Char_t* par, Double_t value)
Bool_t IsParameterSet(const Char_t* par)
Double_t GetParameter(const Char_t* par)
void UnsetParameter(const Char_t* par)
double Double_t

Normalization

A normalization factor can be applied to the value(s) calculated by any global variable, set with one of the following methods:

void SetParameter(const Char_t *par, Double_t value)
Definition KVVarGlob.h:549
virtual void SetNormalization(Double_t norm)
Definition KVVarGlob.h:564

The value calculated for each event will be divided by the normalization factor.

Tagging groups of particles with global variables

Each global variable can be used to 'tag' the particles which are used for its calculation, i.e. those which satisfy the selection criteria set with KVVarGlob::SetSelection(). Each particle will be added to a group which can be subsequently used for iteration (see Nuclei and Events). By default, the group name will be the same as the global variable, but can be changed by giving an argument to the KVVarGlob::SetDefineGroup() method:

// the following variable calculates the multiplicity of 4He nuclei and defines a group
// with name "alpha" containing only those nuclei
KVMult alpha_mult("alpha");
alpha_mult.SetSelection({"4He", [](const KVNucleus* n){ return n->IsIsotope(2,4); });
alpha_mult.SetDefineGroup();
Multiplicity of all nuclei in event (including )
Definition KVMult.h:15

Using lists 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
// etc. etc. ...
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
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:

// Analysis loop
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
// ... any further processing of event, filling histograms, TTree, etc.
}
}

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 ls(Option_t *option="") const override
void ls(Option_t *option="") const override
virtual TObjArray * GetListOfBranches()
#define I(x, y, z)
TLine l

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");
// etc. etc. ...
TTree my_tree;
my_list.MakeBranches(&my_tree); // Note: implicit call to KVGVList::Init()
// Analysis loop
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
}
}
virtual Int_t Fill()

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
void SetEventSelection(const EventSelector &f)
Definition KVVarGlob.h:699
v

(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() )
{
// ... do further analysis, mtot is >=4
}
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_", ...); }

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

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().

User-defined global variable classes

If you want to implement a new global variable, there are two possibilities.

First, check that what you want to do is not possible with e.g. the KVVGSum class, or by changing the selection criteria and/or kinematical reference frame of an existing global variable: there is no need to write a new class if you want to calculate the total charge of fragments with \(5\leq Z\leq 15\) emitted in polar angle range \(10^o\leq\theta\leq 35^o\) in the centre-of-mass frame: just use a suitably modified KVZtot:

KVZtot my_ztot("zsum_z5to15_t10to35");
my_ztot.SetFrame("CM");
my_ztot.SetSelection("5<=Z<=15", [](const KVNucleus* n){ return (n->GetZ()>4 && n->GetZ()<16); });
my_ztot.AddSelection("10<=theta<=35", [](const KVNucleus* n){ return (n->GetTheta()>=10 && n->GetTheta()<=35); });

‍Note that in the second selection, no reference is made to the "CM" reference frame: the lambda function will receive a pointer to the kinematics of the nucleus in the "CM" frame, as defined by the call to KVVarGlob::SetFrame().

If not, you can write your own global variable class. You need to decide first of all which base class to use (looking at the inheritance of existing global variables may help you make your choice). If your global variable will calculate the mean (and possibly also the variance) of some quantity, the base class to use is KVVarGlobMean. If your global variable will calculate a single value, use KVVarGlob1. If neither of the previous two cases apply, use the most general base class, KVVarGlob.

Next, decide if your global variable is of type 1-body (calculation based on properties of individual particles of the event, independently of the others), 2-body (calculation based on correlations between pairs of particles), or N-body (multibody correlations) (looking at existing classes may help to decide).

Then, use one of the following methods on the command line in order to generate skeleton '.h' and '.cpp' files for your class:

KVVarGlob::MakeClass("MyClassName", "The description of my class", type)
KVVarGlob1::MakeClass("MyClassName", "The description of my class", type)
KVVarGlobMean::MakeClass("MyClassName", "The description of my class", type)
static void MakeClass(const Char_t *classname, const Char_t *classdesc, int type=KVVarGlob::kOneBody)
static void MakeClass(const Char_t *classname, const Char_t *classdesc, int type=KVVarGlob::kOneBody)
static void MakeClass(const Char_t *classname, const Char_t *classdesc, int type=kOneBody)
Definition KVVarGlob.cpp:45

with 'type' equal to one of KVVarGlob::kOneBody, KVVarGlob::kTwoBody, KVVarGlob::kNBody.

When you have modified the skeleton according to your needs (see comments in generated code for help), you can test the compilation of your class in an interactive session:

kaliveda[0] .L MyClassName.cpp+

Using your own global variables in a data analysis class

You need to set up your environment correctly in order to tell ROOT where to find the sources for your global variables (this is essential for analysis tasks to run in batch mode), and in order to add your variables to the list available for KaliVeda analysis.

First modify (or create if it doesn't exist) your .rootrc file, adding/modifying the following line:

+Unix.*.Root.MacroPath: /full/path/to/directory/with/source/files
+ACLiC.IncludePaths: -I/full/path/to/directory/with/source/files

You can use environment variables in this definition, as long as you enclose them in ‘’$()', e.g.'/root'`. If you have several source directories, you can put them together:

+Unix.*.Root.MacroPath: /first/path:/second/path:/another/path
+ACLiC.IncludePaths: -I/first/path -I/second/path -I/another/path

After relaunching ROOT/KaliVeda, you will now be able to compile your class(es) even if you launch in a different directory.

Next, modify (or create if it doesn't exist) your .kvrootrc file, adding for each of your global variables a line like this:

+Plugin.KVVarGlob: MyClassName MyClassName MyClassName.cpp+ "MyClassName()"

Then you can use your global variable in an analysis in exactly the same way as other variables:

[in my_analyis.cpp]
auto myvg = AddGV("MyClassName", "myVarGlob");

You should try to avoid adding any new methods to your class (i.e. not defined by the interface of KVVarGlob). If you need additional parameters (named double values) and/or options (named TString values) for your global variables, use the methods provided by KVVarGlob (see Options, parameters, reference frames, particle selection, etc.).

If there is no alternative but to add a new method, you can get away with it if you do:

[in my_analyis.cpp]
#include "MyClassName.h"
...
auto myvg = dynamical_cast<MyClassName*>(AddGV("MyClassName", "myVarGlob"));