KaliVeda
Toolkit for HIC analysis
KVEventFiltering Class Reference

Detailed Description

Filter simulated events with multidetector response.

Use this KVEventSelector on simulated data TTrees containing a branch with KVSimEvent-derived objects:

my_tree->Process("KVEventFiltering", "[options]");

where "[options]" is the list of options in the form "BranchName=toto,Dataset=titi,System=tata, ...".

The following options must be given:

  • SimFileName= name of file containing simulated events
  • SimTitle= description of simulation
  • BranchName= name of branch containing simulated events
  • Dataset= name of experimental dataset (defines multidetector array to use etc.)
  • System= name of experimental dataset system (defines collision kinematics etc.)
  • Geometry= type of geometry, either 'KV' (KaliVeda geometry) or 'ROOT' (ROOT TGeometry package)
  • Filter= type of filter, either Geo (geometric filter), GeoThresh (geometry + detector thresholds), or Full (full simulation of detector response, including experimental identification and calibration routines)
  • OutputDir= directory path in which to write filtered data file
  • Kinematics= kinematical frame for simulation, either cm, projectile or lab. if cm/projectile, we use the c.m./projectile velocity of the selected System to transform events into the detector (laboratory) frame. if lab no transformation is performed: simulated events are already in laboratory frame.

The following are optional options:

  • Run: run number to use for detector status, setup, parameters, etc. if not given, first run of the given experimental system is used.
  • PhiRot: by default, a random rotation around the beam axis will be performed before simulating detection of the event. If you don't want this to happen, give option PhiRot=no If used, filtered events will have a parameter RANDOM_PHI with the applied rotation (in radians)
  • Gemini: if option Gemini=yes, then each event will be "decayed" with Gemini++, if KaliVeda has been compiled with Gemini++ support. See below for extra information on Gemini decay stored in particle parameter lists.
  • GemDecayPerEvent: if option Gemini=yes then by default 1 Gemini++ decay will be performed for each event. you can change this by giving a value for this option

The filtered data will be written in the directory given as option "OutputDir". The filename is built up from the original simulation filename and the values of various options:

[simfile]_geo=[geometry]_filt=[filter-type]_[dataset]_[system]_run=[run-number].root

The data will be stored in a TTree with name ReconstructedEvents, in a branch with name ReconEvent. The class used for reconstructed events depends on the dataset, it is given by KVDataSet::GetReconstructedEventClassName() (usually KVReconstructedEvent).

Each filtered event will have some or all of the following parameters in its list:

  • SIMEVENT_TREE_ENTRY = index of simulated event in TTree we are reading
  • SIMEVENT_NUMBER = event number of simulated event in TTree we are reading (if defined i.e. if !=0)
  • RANDOM_PHI = rotation around beam axis [in radians]

When using Gemini++ to decay simulated events before filtering, each particle will have the following parameter defined:

  • GEMINI_PARENT_INDEX = index of parent nucleus in simulated event

The combination of SIMEVENT_TREE_ENTRY and GEMINI_PARENT_INDEX allows to connect primary fragments in the simulated events with their detected decay products in the filtered event.

Example: suppose recev is a filtered event and sim_tree is the TTree containing the original simulation, while simev is a pointer to a KVSimEvent object connected to the appropriate branch in sim_tree. Then we can retrieve the parent nucleus responsible for producing a particular nucleus in recev like so:

sim_tree.GetEntry( recev.GetParameters()->GetIntValue("SIMEVENT_TREE_ENTRY" ) );
KVSimNucleus* parent_nuc = (KVSimNucleus*)simev->GetParticle( recev.GetNucleus(15)->GetParameters()->GetIntValue("GEMINI_PARENT_INDEX") ); // parent nucleus of 15th nucleus in filtered event 'recev'
Nucleus in a simulated event.
Definition: KVSimNucleus.h:32

Definition at line 97 of file KVEventFiltering.h.

#include <KVEventFiltering.h>

Inheritance diagram for KVEventFiltering:

Public Member Functions

 KVEventFiltering ()
 Default constructor. More...
 
 KVEventFiltering (const KVEventFiltering &)
 
virtual ~KVEventFiltering ()
 Destructor. More...
 
Bool_t Analysis ()
 
void Copy (TObject &) const
 
void EndAnalysis ()
 
void EndRun ()
 
void InitAnalysis ()
 
void InitRun ()
 
void OpenOutputFile (KVDBSystem *, Int_t)
 
- Public Member Functions inherited from KVEventSelector
 KVEventSelector (TTree *=0)
 
KVVarGlobAddGV (const Char_t *class_name, const Char_t *name)
 
void AddGV (KVVarGlob *vg)
 
KVGVListAddGVList (const KVString &list_name, const KVParticleCondition &selection=KVParticleCondition())
 
template<typename HistoType , typename... Args>
HistoType * AddHisto (Args &&... args)
 
void AddHisto (TH1 *histo)
 
TTree * AddTree (const TString &name, const TString &title="")
 
void AddTree (TTree *tree)
 
virtual void Begin (TTree *tree)
 
Bool_t CreateTreeFile (const Char_t *filename="")
 
void FillHisto (const Char_t *sname, const Char_t *label, Double_t weight=1)
 
void FillHisto (const Char_t *sname, Double_t one, Double_t two=1, Double_t three=1, Double_t four=1)
 
void FillTree (const Char_t *tree_name="")
 
const Char_t * GetBranchName () const
 
virtual Int_t GetEntry (Long64_t entry, Int_t getall=0)
 
KVEventGetEvent () const
 
Int_t GetEventNumber () const
 
KVEventGetFriendEvent () const
 
Int_t GetFriendTreeEntry (Long64_t entry, Int_t getall=0)
 
KVVarGlobGetGV (const Char_t *name) const
 
KVGVListGetGVList (const KVString &list_name="default")
 
const KVGVListGetGVList (const KVString &list_name="default") const
 
TH1 * GetHisto (const Char_t *name) const
 
const KVHashListGetHistoList () const
 
TString GetOpt (const Char_t *option) const
 
virtual TList * GetOutputList () const
 
TTree * GetTree (const Char_t *name) const
 
const KVHashListGetTreeList () const
 GetTreeList. More...
 
virtual void Init (TTree *tree)
 
Bool_t IsOptGiven (const Char_t *option)
 
virtual Bool_t Notify ()
 
virtual Bool_t Process (Long64_t entry)
 
virtual void SaveHistos (const Char_t *="", Option_t *="recreate", Bool_t=kFALSE)
 
virtual void SetAdditionalBranchAddress ()
 
void SetBranchName (const Char_t *n)
 
virtual void SetCurrentRun (KVDBRun *)
 
void SetEventsReadInterval (Long64_t N)
 
virtual void SetInputList (TList *input)
 
void SetJobOutputFileName (const TString &filename)
 
virtual void SetObject (TObject *obj)
 
void SetOpt (const Char_t *option, const Char_t *value)
 
void SetParticleConditions (const KVParticleCondition &cond)
 
void SetTriggerConditionsForRun (int run)
 
virtual void SlaveBegin (TTree *tree)
 
virtual void SlaveTerminate ()
 
virtual void Terminate ()
 
void UnsetOpt (const Char_t *opt)
 
virtual Int_t Version () const
 

Public Attributes

TVector3 fCMVelocity
 
TString fNewFrame
 allow the definition of a specific frame More...
 
TVector3 fProjVelocity
 
KVReconstructedEventfReconEvent
 
Bool_t fTransformKinematics
 =kTRUE if simulation not in lab frame More...
 

Private Member Functions

void RandomRotation (KVEvent *to_rotate, const TString &frame_name="") const
 

Private Attributes

Long64_t fEVN
 event number counter More...
 
Bool_t fGemAddRotEner
 true if rotational energy has to be added to excitation energy [default: no] More...
 
Int_t fGemDecayPerEvent
 number of Gemini++ decays to be performed for each event [default:1] More...
 
KVSimEvent fGemEvent
 event after decay with Gemini More...
 
Bool_t fGemini
 true if Gemini++ decay should be performed before detection [default: no] More...
 
TString fIdCalMode
 original exp setup hasIDandCalib to be reset in case of modifications More...
 
Bool_t fRotate
 true if random phi rotation should be applied [default: yes] More...
 
KVGemini GEM
 

Constructor & Destructor Documentation

◆ KVEventFiltering() [1/2]

KVEventFiltering::KVEventFiltering ( )

Default constructor.

Definition at line 21 of file KVEventFiltering.cpp.

◆ KVEventFiltering() [2/2]

KVEventFiltering::KVEventFiltering ( const KVEventFiltering obj)

Copy constructor This ctor is used to make a copy of an existing object (for example when a method returns an object), and it is always a good idea to implement it. If your class allocates memory in its constructor(s) then it is ESSENTIAL :-)

Definition at line 48 of file KVEventFiltering.cpp.

◆ ~KVEventFiltering()

KVEventFiltering::~KVEventFiltering ( )
virtual

Destructor.

Definition at line 66 of file KVEventFiltering.cpp.

Member Function Documentation

◆ Analysis()

Bool_t KVEventFiltering::Analysis ( void  )
virtual

Event-by-event filtering of simulated data. If needed (fTransformKinematics = kTRUE), kinematics of event are transformed to laboratory frame using C.M. velocity calculated in InitAnalysis(). Detection of particles in event is simulated with KVMultiDetArray::DetectEvent, then the reconstructed detected event is treated by the same identification and calibration procedures as for experimental data. few event counter print at the beginning to be sure the process started properly because in case GEMINI decay is used, it sometimes stops (randomly) after few events

Reimplemented from KVEventSelector.

Definition at line 137 of file KVEventFiltering.cpp.

◆ Copy()

void KVEventFiltering::Copy ( TObject &  obj) const

This method copies the current state of 'this' object into 'obj' You should add here any member variables, for example: (supposing a member variable KVEventFiltering::fToto) CastedObj.fToto = fToto; or CastedObj.SetToto( GetToto() );

Definition at line 85 of file KVEventFiltering.cpp.

◆ EndAnalysis()

void KVEventFiltering::EndAnalysis ( )
virtual

Reimplemented from KVEventSelector.

Definition at line 212 of file KVEventFiltering.cpp.

◆ EndRun()

void KVEventFiltering::EndRun ( )
virtual

Reimplemented from KVEventSelector.

Definition at line 221 of file KVEventFiltering.cpp.

◆ InitAnalysis()

void KVEventFiltering::InitAnalysis ( )
virtual

Select required dataset for filtering (option "Dataset") Build the associated multidetector geometry. Calculate C.M. velocity associated with required experimental collision kinematics (option "System").

Set the parameters of the detectors according to the required run if given (option "Run"), or the first run of the given system otherwise. If ROOT/TGeo geometry is required (option "Geometry"="ROOT"), build the TGeometry representation of the detector array.

Open file for filtered data (see KVEventFiltering::OpenOutputFile), which will be stored in a TTree with name 'ReconstructedEvents', in a branch with name 'ReconEvent'. The class used for reconstructed events depends on the dataset, it is given by KVDataSet::GetReconstructedEventClassName().

Reimplemented from KVEventSelector.

Definition at line 243 of file KVEventFiltering.cpp.

◆ InitRun()

void KVEventFiltering::InitRun ( void  )
virtual

Reimplemented from KVEventSelector.

Definition at line 413 of file KVEventFiltering.cpp.

◆ OpenOutputFile()

void KVEventFiltering::OpenOutputFile ( KVDBSystem S,
Int_t  run 
)

Open ROOT file for new filtered events TTree. The file will be written in the directory given as option "OutputDir". The filename is built up from the original simulation filename and the values of various options:

[simfile]_[Gemini]_geo=[geometry]_filt=[filter-type]_[dataset]_[system]_run=[run-number].root

In addition, informations on the filtered data are stored in the file as TNamed objects. These can be read by KVSimDir::AnalyseFile:

KEY: TNamed System;1 title=[full system name] KEY: TNamed Dataset;1 title=[dataset name] KEY: TNamed Run;1 title=[run-number] KEY: TNamed Geometry;1 title=[geometry-type] KEY: TNamed Filter;1 title=[filter-type] KEY: TNamed Origin;1 title=[name of simulation file] KEY: TNamed RandomPhi;1 title=[yes/no, random rotation about beam axis] KEY: TNamed Gemini++;1 title=[yes/no, Gemini++ decay before detection] KEY: TNamed GemDecayPerEvent;1 title=[number of Gemini++ decays per primary event] KEY: TNamed GemAddRotEner;1 title=[Enable or not the addition of the rotational energy to the excitation energy]

Definition at line 446 of file KVEventFiltering.cpp.

◆ RandomRotation()

void KVEventFiltering::RandomRotation ( KVEvent to_rotate,
const TString &  frame_name = "" 
) const
private

do random phi rotation around z-axis if frame_name is given, apply rotation to that frame

store phi rotation angle [radians] in event parameter "RANDOM_PHI"

Definition at line 111 of file KVEventFiltering.cpp.

Member Data Documentation

◆ fCMVelocity

TVector3 KVEventFiltering::fCMVelocity

Definition at line 128 of file KVEventFiltering.h.

◆ fEVN

Long64_t KVEventFiltering::fEVN
private

event number counter

Definition at line 101 of file KVEventFiltering.h.

◆ fGemAddRotEner

Bool_t KVEventFiltering::fGemAddRotEner
private

true if rotational energy has to be added to excitation energy [default: no]

Definition at line 105 of file KVEventFiltering.h.

◆ fGemDecayPerEvent

Int_t KVEventFiltering::fGemDecayPerEvent
private

number of Gemini++ decays to be performed for each event [default:1]

Definition at line 106 of file KVEventFiltering.h.

◆ fGemEvent

KVSimEvent KVEventFiltering::fGemEvent
private

event after decay with Gemini

Definition at line 107 of file KVEventFiltering.h.

◆ fGemini

Bool_t KVEventFiltering::fGemini
private

true if Gemini++ decay should be performed before detection [default: no]

Definition at line 104 of file KVEventFiltering.h.

◆ fIdCalMode

TString KVEventFiltering::fIdCalMode
private

original exp setup hasIDandCalib to be reset in case of modifications

Definition at line 111 of file KVEventFiltering.h.

◆ fNewFrame

TString KVEventFiltering::fNewFrame

allow the definition of a specific frame

Definition at line 131 of file KVEventFiltering.h.

◆ fProjVelocity

TVector3 KVEventFiltering::fProjVelocity

Definition at line 129 of file KVEventFiltering.h.

◆ fReconEvent

KVReconstructedEvent* KVEventFiltering::fReconEvent

Definition at line 127 of file KVEventFiltering.h.

◆ fRotate

Bool_t KVEventFiltering::fRotate
private

true if random phi rotation should be applied [default: yes]

Definition at line 102 of file KVEventFiltering.h.

◆ fTransformKinematics

Bool_t KVEventFiltering::fTransformKinematics

=kTRUE if simulation not in lab frame

Definition at line 130 of file KVEventFiltering.h.

◆ GEM

KVGemini KVEventFiltering::GEM
private

Definition at line 108 of file KVEventFiltering.h.