KaliVeda
Toolkit for HIC analysis
ExampleRawAnalysis.cpp

Example of an analysis class for raw (i.e. unreconstructed) data

This is the analysis class generated by default by KaliVedaGUI for raw data analysis.

#ifndef __EXAMPLERAWANALYSIS_H
#define __EXAMPLERAWANALYSIS_H
#include "KVRawDataAnalyser.h"
#include <vector>
#include <string>
class ExampleRawAnalysis : public KVRawDataAnalyser {
struct raw_data_t {
int n_sig;
std::vector<int> det_index;
std::vector<std::string> det_type;
std::vector<std::string> sig_type;
std::vector<float> sig_val;
void clear()
{
n_sig = 0;
det_index.clear();
det_type.clear();
sig_type.clear();
sig_val.clear();
}
};
raw_data_t raw_data;
public:
ExampleRawAnalysis() {}
virtual ~ExampleRawAnalysis() {}
void InitAnalysis();
void InitRun();
void EndRun() {}
void EndAnalysis() {}
ClassDef(ExampleRawAnalysis, 1) //Analysis of raw data
};
#endif
int Int_t
bool Bool_t
#define ClassDef(name, id)
Abstract base class for user analysis of raw data.
virtual Bool_t Analysis()=0
virtual void InitRun()=0
virtual void EndAnalysis()=0
virtual void InitAnalysis()=0
virtual void EndRun()=0
DisplacementVector3D< CoordSystem, U > Mult(const Matrix &m, const DisplacementVector3D< CoordSystem, U > &v)
#include "ExampleRawAnalysis.h"
#include "KVMultiDetArray.h"
ClassImp(ExampleRawAnalysis)
// This class is derived from the KaliVeda class KVRawDataAnalyser.
// It is to be used for analysis of raw (i.e. unreconstructed) data.
//
// The following member functions are called in turn:
//
// InitAnalysis(): called at the very beginning of the analysis
// InitRun(): called everytime a run starts
// Analysis(): called for each event
// EndRun(): called everytime a run ends
// EndAnalysis(): called after all requested data has been read
//
// Modify these methods as you wish in order to create your analysis class.
void ExampleRawAnalysis::InitAnalysis()
{
// Declaration of histograms, trees, etc.
// Called at the beginning of the analysis
// The examples given are compatible with interactive and batch analyses.
/*** DECLARING SOME HISTOGRAMS ***/
AddHisto<TH1F>("Mult", "Number of fired detectors in each event", 1000, -.5, 999.5);
/*** USING A TREE ***/
// Here we fill a TTree every event with
// - the number of fired detector signals
// - for each signal:
// + detector index
// + detector type
// + signal type
// + signal value
auto tt = AddTree("raw_data", "");
tt->Branch("n_sig", &raw_data.n_sig);
tt->Branch("det_index", &raw_data.det_index);
tt->Branch("det_type", &raw_data.det_type);
tt->Branch("sig_type", &raw_data.sig_type);
tt->Branch("sig_val", &raw_data.sig_val);
/*** DEFINE WHERE TO SAVE THE RESULTS ***/
// This filename will be used for interactive jobs.
// When running in batch mode, this will automatically use the job name.
SetJobOutputFileName("ExampleRawAnalysis_results.root");
}
//____________________________________________________________________________//
void ExampleRawAnalysis::InitRun()
{
// Initialisation performed at beginning of each runfile. The experimental setup (KVMultiDetArray
// or KVExpSetUp object) is initialised with the parameters of the run just before calling this method.
//
// Some useful global pointers/methods:
//
// - gMultiDetArray : global pointer to experiment setup description;
// - gExpDB : global pointer to database (KVExpDB) for current dataset;
// - gDataSet : global pointer to current dataset (KVDataSet)
// - GetRunNumber() : returns current run number;
// - GetRunIndexNumber() : returns full run/index of current runfile;
// - GetCurrentRun() : returns pointer to KVDBRun object describing current run in database;
// - GetCurrentRunFile() : return pointer to KVDBRunFile object describing current runfile in database;
// - GetSystem() : return pointer to KVDBSystem object describing current reaction in database.
//
// See parent classes KVRawDataAnalyser, KVDataSetAnalyser and KVDataAnalyser for more.
Info("InitRun", "Beginning analysis of runfile %s containing %llu events", GetRunIndexNumber().as_string().Data(), GetCurrentRun()->GetEvents());
GetCurrentRunFile()->ls();
// Set tree title to name of system being analysed
GetTree("raw_data")->SetTitle(GetSystem()->GetName());
}
//____________________________________________________________________________//
Bool_t ExampleRawAnalysis::Analysis()
{
//Analysis method called for each event
//
// - GetEventNumber() : returns current event number;
// - GetRunFileReader() : returns object used to read data (KVRawDataReader child class);
// - gMultiDetArray->HandledRawData() : returns kTRUE if interesting data was read from event;
//
// Processing will stop if this method returns kFALSE
if (gMultiDetArray->HandledRawData()) {
raw_data.clear();
// get number of fired detectors
Mult = gMultiDetArray->GetFiredDetectors()->GetEntries();
// Note - if experimental set up is a combination of two or more detector arrays,
// and you are interested in only one of them, you could use
// gMultiDetArray->GetArray("[name]")->GetFiredDetectors()
FillHisto("Mult", Mult);
// get list of fired signals
auto sigs = gMultiDetArray->GetFiredSignals();
KVDetectorSignal* det_sig;
TIter it(sigs);
while ((det_sig = (KVDetectorSignal*)it())) {
++raw_data.n_sig;
raw_data.det_index.push_back(det_sig->GetDetector()->GetIndex());
raw_data.det_type.push_back(det_sig->GetDetector()->GetLabel());
raw_data.sig_type.push_back(det_sig->GetType());
raw_data.sig_val.push_back(det_sig->GetValue());
}
}
return kTRUE;
}
constexpr Bool_t kTRUE
Base class for output signal data produced by a detector.
const Char_t * GetType() const override
virtual Double_t GetValue(const KVNameValueList &params="") const
const KVDetector * GetDetector() const
const KVSeqCollection * GetFiredSignals() const
Bool_t HandledRawData() const
const KVSeqCollection * GetFiredDetectors() const
virtual Int_t GetEntries() const
void Info(const char *location, const char *fmt,...)
void FillTree(TTree &myTree, const RooDataSet &data)
auto * tt
ClassImp(TPyArg)