KaliVeda
Toolkit for HIC analysis
Loading...
Searching...
No Matches
KVSimDirFilterAnalyser.cpp
1//Created by KVClassFactory on Thu Mar 23 14:02:40 2017
2//Author: John Frankland,,,
3
4#include "KVSimDirFilterAnalyser.h"
5#include <KVClassFactory.h>
6#include <KVDataSetManager.h>
7#include <KVExpDB.h>
8#include <KVTriggerConditions.h>
9#include "KVMultiDetArray.h"
10#include "KVNucleusEvent.h"
11
13
14
15
17
20{
21 // Default constructor
22}
23
24
25
26
29
36
37
38
41
43{
44 // Use options passed to KVEventSelector to initialise dataset used for filter
45
46 if (!gDataSetManager) {
47 gDataSetManager = new KVDataSetManager;
48 gDataSetManager->Init();
49 }
50 if (!gDataSet) gDataSetManager->GetDataSet(fAnalysisClass->GetOpt("DataSet"))->cd();
51}
52
53
54
58
60{
61 // Use options passed to KVEventSelector to initialise kinematics for reaction
62 // and build experimental set-up
63
64 KVString system = fAnalysisClass->GetOpt("System");
65 fSystem = gExpDB->GetSystem(system);
66 if (fSystem) {
68 }
69 else {
70 fKinematics = new KV2Body(system);
71 fKinematics->SetBit(kCanDelete); // clean up in destructor
73 }
74 KVString run = fAnalysisClass->GetOpt("Run");
75 fRun = gExpDB->GetDBRun(run.Atoi());
78 // put detector array in simulation mode
79 gMultiDetArray->SetSimMode();
80}
81
82
83
88
90{
91 // Called by KVEventSelector just after reading new event and just before calling user analysis
92 //
93 // Copy calculated energy losses of all particles into the detectors
94
95 for (auto& n : EventIterator(fAnalysisClass->GetEvent())) {
96 for (auto& p : *n.GetParameters()) {
97 auto d = gMultiDetArray->GetDetector(p.GetName());
98 if (d) {
99 // sum up energy losses in each detector
100 d->SetEnergyLoss(p.GetDouble() + d->GetEnergyLoss());
101 }
102 }
103 }
104}
105
106
107
110
112{
113 // Generate a new filtered analysis selector class
114
115 KVClassFactory cf(kvsname, "Analysis of filtered simulated events", "",
116 kTRUE, "ROOT6FilteredEventAnalysisTemplate");
117 cf.AddImplIncludeFile("KVReconstructedNucleus.h");
118 cf.AddImplIncludeFile("KVBatchSystem.h");
119
120 cf.GenerateCode();
121}
122
123
124
141
143{
144 // When called from the InitRun() method of a user's analysis class, this method will ensure that only data
145 // compatible with the experimental trigger will be provided for analysis in the user's Analysis() method.
146 //
147 // This will be done by searching for a KVTriggerConditions plugin class defined for the currently-analysed
148 // dataset, defined like so:
149 //
150 //~~~~
151 //+Plugin.KVTriggerConditions: [dataset] [classname] [libname] "[default constructor]()"
152 //~~~~
153 //
154 // An object of the plugin class will be instantiated, and then its overridden
155 // KVTriggerConditions::SetTriggerConditionsForRun() method will be called with 3 arguments:
156 // - a pointer to the user's analysis class derived from KVEventSelector;
157 // - the number of the run currently being analysed
158 // - the 3rd optional argument will be set to kTRUE to indicate that we are analysing filtered simulated data
159
160 TPluginHandler* ph = KVBase::LoadPlugin("KVTriggerConditions", gDataSet->GetName());
161 if (!ph) {
162 Info("SetTriggerConditionsForRun",
163 "No definition of trigger conditions available for dataset %s",
164 gDataSet->GetName());
165 return;
166 }
167 std::unique_ptr<KVTriggerConditions> trig((KVTriggerConditions*) ph->ExecPlugin(0));
168 trig->SetTriggerConditionsForRun(fAnalysisClass, run, kTRUE);
169}
170
171//____________________________________________________________________________//
172
173
#define d(i)
char Char_t
constexpr Bool_t kTRUE
winID h TVirtualViewer3D TVirtualGLPainter p
Class for iterating over nuclei in events accessed through base pointer/reference.
Relativistic binary kinematics calculator.
Definition KV2Body.h:166
void CalculateKinematics()
Definition KV2Body.cpp:677
static TPluginHandler * LoadPlugin(const Char_t *base, const Char_t *uri="0")
Definition KVBase.cpp:793
Factory class for generating skeleton files for new classes.
void GenerateCode()
Generate header and implementation file for currently-defined class.
void AddImplIncludeFile(const Char_t *filename)
KV2Body * GetKinematics()
Manage all datasets contained in a given data repository.
virtual Bool_t Init(KVDataRepository *=0)
KVDataSet * GetDataSet(Int_t) const
Return pointer to DataSet using index in list of all datasets, index>=0.
void cd() const
KVEvent * GetEvent() const
virtual void SetCurrentRun(KVDBRun *)
TString GetOpt(const Char_t *option) const
virtual KVDBSystem * GetSystem(const Char_t *system) const
Definition KVExpDB.h:85
KVDBRun * GetDBRun(Int_t number) const
Definition KVExpDB.h:76
KVDetector * GetDetector(const Char_t *name) const
Return detector in this structure with given name.
static KVMultiDetArray * MakeMultiDetector(const Char_t *dataset_name, Int_t run=-1, TString classname="KVMultiDetArray")
virtual void SetSimMode(Bool_t on=kTRUE)
Class piloting analyses of simulated data.
Manage user analysis of filtered simulation data.
virtual ~KVSimDirFilterAnalyser()
Destructor.
static void Make(const Char_t *kvsname="MyFilteredAnalysis")
Generate a new filtered analysis selector class.
KVDBRun * fRun
currently analysed run
KVEventSelector * fAnalysisClass
user analysis class
void preInitAnalysis()
Use options passed to KVEventSelector to initialise dataset used for filter.
KVDBSystem * fSystem
currently analysed system
KV2Body * fKinematics
kinematics of reaction
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition KVString.h:73
Set trigger conditions for analysis of reconstructed data.
const char * GetName() const override
virtual const char * GetName() const
void SetBit(UInt_t f)
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
virtual void Info(const char *method, const char *msgfmt,...) const
Longptr_t ExecPlugin(int nargs)
Int_t Atoi() const
const Int_t n
ClassImp(TPyArg)