KaliVeda
Toolkit for HIC analysis
Data Analysis

Here are some guidelines for analysing data. Data can either be raw experimental data, reconstructed experimental data, or simulated data, either filtered with a given experimental setup (in which case it is again reconstructed data) or not.

Each data analysis requires a user analysis class to define the tasks to be performed, histograms/trees to fill etc. The user's analysis class is piloted by a class derived from KVDataAnalyser which handles opening the data files requested by the user and setting up the necessary environment for the analysis to proceed.

Two graphical interfaces are provided for analysing data: KaliVedaGUI and kaliveda-sim. Both interfaces are capable of generating example analysis classes for the data you want to analyse. This is the recommended method for creating a new analysis class, which you can then modify according to your needs.

N.B. In the case of experimental data, the exact inheritance required for your analysis class to work may depend on the dataset being analysed. This is handled automatically when generating a new class with KaliVedaGUI

Different analysis tasks for experimental data

When analysing experimental data with KaliVedaGUI, different analysis tasks appear in the drop-down menu depending on the types of available data. Here is a brief explanation of each task.

Analysis tasks for <tt>raw</tt> data

raw data is data written during the experiment, i.e. the files produced by the DAQ system. Several analysis tasks are possible:

  • Analysis of raw data (base class: KVRawDataAnalyser)
    • This is the most basic analysis task. The user's Analysis() method (see below) is called after each event is read from the DAQ file. It is the user's responsibility to test whether any useful information was actually present in the event by testing gMultiDetArray->HandledRawData() before proceeding with any analysis (in case the data being read comes from several independent detector arrays, this method returns kTRUE if at least one of the arrays in the set-up found data it could handle; you can test individual arrays with gMultiDetArray->GetArray([name])->HandledRawData() - see class KVExpSetUp). Values of acquisition parameters etc. will be set in the detectors of the array(s) which fired in the event. Other useful information may be stored in the KVNameValueList returned by gMultiDetArray->GetReconParameters() (or the equivalent for one of the individual arrays in case of a multi-array setup).
  • Event reconstruction from raw data (base class: KVRawDataReconstructor)
    • This task will generate new recon data files from raw data (no user analysis is possible). Event reconstruction is handled by the KVEventReconstructor class which manages a set of KVGroupReconstructor-derived objects which may be specialised for particle reconstruction in the groups of each multidetector array (a group corresponds to the largest subset of detectors which can be treated independently of all others in the array - see Using Detector Array Geometries). The reconstructed events are written in a ROOT TTree in a new file. The generated data can be analysed using the reconstructed data analysis tasks presented below.
  • Analysis of reconstructed raw data (base class: KVReconRawDataAnalyser)
    • Similar to analysis of raw data, except that event reconstruction is performed just before the user's Analysis() method is called. It is the user's responsibility to test whether any useful information was actually present in the event by testing gMultiDetArray->HandledRawData() : if it returns kTRUE then the reconstructed event can be accessed using method GetReconstructedEvent().

Analysis tasks for <tt>recon</tt> (reconstructed) data

Reconstructed data is contained in ROOT TTrees generated by the Event reconstruction from raw data task presented above. Specifically, for each DAQ event a KVReconstructedEvent is filled with KVReconstructedNucleus objects which are deduced from the detector hits and the array geometry.

Generating a new analysis class

  • With [KaliVedaGUI]:
    • select the data you want to analyse, then in the drop-down User Class list select [NEW]. You will be asked to provide a name for your new analysis class, then the newly-generated template files will be opened in your favourite text editor.
  • With [kaliveda-sim]:
    • click on either New simulated analysis class or New filtered analysis class, depending on which kind of data you want to analyse. You will be asked to provide a name for your new analysis class, which you can then open in your favourite text editor (it isn't automatic).

Structure of the analysis class

Whatever kind of data you wish to analyse, you will now have a more-or-less complete example of an analysis class contained in the files MyClass.h and MyClass.cpp (assuming you replied MyClass when asked for a name for the class). The exact heritage of this class will depend on the kind of data being analysed, but in all cases it will have the following 5 methods:

void MyClass::InitAnalysis();
void MyClass::InitRun();
Bool_t MyClass::Analysis();
void MyClass::EndRun();
void MyClass::EndAnalysis();
bool Bool_t

InitAnalysis() and EndAnalysis() are executed once only at the beginning and end, respectively, of the analysis. InitRun() and EndRun() are called at the beginning and end, respectively, of each data file which is being analysed. Analysis() is called for each and every event read from the files. Each of these methods in your MyClass.cpp file should contain useful comments and/or examples of use.

It should be noted that for all except raw data analysis (i.e. for all data stored in ROOT trees), the analysis classes inherit from KVEventSelector which allows them to be used with PROOF on multi-core machines. It is important to read the guidelines given in the KVEventSelector class reference for this to work.

The user's analysis class can acces the KVDataAnalyser object which is running the analysis through the gDataAnalyser global pointer (see below).

Defining histograms, trees, global variables: <tt>InitAnalysis()</tt>

This is the place to create/open any ROOT files, histograms, trees etc. you require in order to store the output of your analysis. Also the place to define any Global Variables needed for the analysis: see the dedicated chapter for more details. Also see the examples given and comments in the example analysis classes.

If you are analysing reconstructed experimental/simulated data, do not try to access the experimental set-up (gMultiDetArray, etc.) associated with the data you want to analyse here - it will not be defined until InitRun().

Defining particle selection criteria: <tt>InitRun()</tt>

Called once at the beginning of each new data file. In case of analysing experimental reconstructed data you can access the experimental set-up (gMultiDetArray, system being analysed, run number, etc.) associated with the data here if you need to:

GetCurrentRun()->GetNumber(); // the number of the run analysed
GetCurrentRun()->GetSystem(); // returns KVDBSystem* pointer to collision system being analysed

For reconstructed data of either experimental or filtered simulated type, you can access the kinematics of the collision system being analysed like so:

gDataAnalyser->GetKinematics(); // returns KV2Body* pointer to collision kinematics
virtual const KV2Body * GetKinematics() const

Reconstructed particle selection

When analysing reconstructed data, InitRun() is the place to modify or fine-tune if you require the selections which will be applied to the particles of each event in order to include them in or exclude them from the analysis. These are the particles which will be used for the calculation of any global variables declared in InitAnalysis(), and those which will be labelled as "OK" (for example when iterating over the event).

Identification & calibration codes

The first selection applied is to consider "OK" only those particles which have acceptable identification codes (KVReconstructedNucleus::GetIDCode()) and calibration codes (KVReconstructedNucleus::GetECode()). For each multidetector array, default values of the acceptable codes are defined using configuration variables such as

[dataset].[name].ReconstructedNuclei.AcceptIDCodes: 2-4,6
[dataset].[name].ReconstructedNuclei.AcceptECodes: 1-2

where [name] is the name of the array and optionally the values can be dataset dependent.

To change the default acceptable codes for particles reconstructed by a given array to be used in your analysis, use

gMultiDetArray->AcceptIDCodes(const KVNumberList&)
gMultiDetArray->AcceptECodes(const KVNumberList&)
void AcceptECodes(const KVNumberList &codelist)
void AcceptIDCodes(const KVNumberList &codelist)
Strings used to represent a set of ranges of values.
Definition: KVNumberList.h:85

using the same format for the lists of codes as in the configuration variable example given above. If analysing data for a set-up consisting of several different detector arrays, you should use

gMultiDetArray->GetArray("[name]")->AcceptIDCodes(const KVNumberList&)
gMultiDetArray->GetArray("[name]")->AcceptECodes(const KVNumberList&)
virtual KVMultiDetArray * GetArray(const Char_t *) const

replacing [name] by the appropriate array name.

Customised particle selections

For any other selections which need to be applied to exclude 'bad' particles from the analysis, you can use the KVParticleCondition class to define any number of selections based on particle properties which will be applied to all particles which passed the identification/calibration code selections.

Note that you cannot 'retrieve' a particle which has been excluded due to bad identification or calibration codes using a KVParticleCondition selection

These selections are then added to the analysis selection criteria using the SetParticleConditions() method:

#ifdef USING_ROOT6
SetParticleConditions({"Z>0", [](const KVNucleus *n){ return n->GetZ()>0; }});
#else
SetParticleConditions("_NUC_->GetZ()>0");
#endif
Description of properties and kinematics of atomic nuclei.
Definition: KVNucleus.h:126
const Int_t n

See the example analysis classes for more examples. See the class reference for KVParticleCondition for usage.

Event-by-event analysis: <tt>Analysis()</tt>

This method is called once for each event to analyse. Implement the event-by-event analysis here.

During reconstructed data analysis, when this method is called

  • all particles have been labelled OK or not according to the criteria set in InitRun();
  • for all OK particles, the CM centre-of-mass reference frame is defined;
  • all global variables added in InitAnalysis() have been calculated using all OK particles

To iterate over all or a subset of the particles of the event, see here.

<tt>EndRun()</tt>, <tt>EndAnalysis()</tt>

There is usually no need to implement these methods.