KaliVeda
Toolkit for HIC analysis
Loading...
Searching...
No Matches
KVEventFiltering.cpp
1//Created by KVClassFactory on Thu Jun 21 17:01:26 2012
2//Author: John Frankland,,,
3
4#include "KVEventFiltering.h"
5#include "KVMultiDetArray.h"
6#include "KV2Body.h"
7#include "KVDBSystem.h"
8#include "KVDBRun.h"
9#include "KVDataSet.h"
10#include "KVDataSetManager.h"
11#include "KVGeoNavigator.h"
12
14
15
16
17
18
20
22{
23 // Default constructor
24 fTransformKinematics = kTRUE;
25 fNewFrame = "";
26 fRotate = kTRUE;
27#ifdef WITH_GEMINI
28 fGemini = kFALSE;
29 fGemDecayPerEvent = 1;
30 fGemAddRotEner = kFALSE;
31#endif
32#ifdef DEBUG_FILTER
33 memory_check = KVClassMonitor::GetInstance();
34 SetEventsReadInterval(100);
35#endif
36}
37
38
39
40
47
49{
50 // Copy constructor
51 // This ctor is used to make a copy of an existing object (for example
52 // when a method returns an object), and it is always a good idea to
53 // implement it.
54 // If your class allocates memory in its constructor(s) then it is ESSENTIAL :-)
55
57 fNewFrame = "";
58 obj.Copy(*this);
59}
60
61
62
65
67{
68 // Destructor
69#ifdef DEBUG_FILTER
70 delete memory_check;
71#endif
72}
73
74
75
76
84
86{
87 // This method copies the current state of 'this' object into 'obj'
88 // You should add here any member variables, for example:
89 // (supposing a member variable KVEventFiltering::fToto)
90 // CastedObj.fToto = fToto;
91 // or
92 // CastedObj.SetToto( GetToto() );
93
96 CastedObj.fRotate = fRotate;
97#ifdef WITH_GEMINI
98 CastedObj.fGemini = fGemini;
99 CastedObj.fGemDecayPerEvent = fGemDecayPerEvent;
100#endif
101}
102
103
104
110
111void KVEventFiltering::RandomRotation(KVEvent* to_rotate, const TString& frame_name) const
112{
113 // do random phi rotation around z-axis
114 // if frame_name is given, apply rotation to that frame
115 //
116 // store phi rotation angle [radians] in event parameter "RANDOM_PHI"
117 TRotation r;
119 r.RotateZ(phi);
120 if (frame_name != "") to_rotate->SetFrame("rotated_frame", frame_name, r);
121 else to_rotate->SetFrame("rotated_frame", r);
122 to_rotate->SetParameter("RANDOM_PHI", phi);
123}
124
125
126
136
138{
139 // Event-by-event filtering of simulated data.
140 // If needed (fTransformKinematics = kTRUE), kinematics of event are transformed
141 // to laboratory frame using C.M. velocity calculated in InitAnalysis().
142 // Detection of particles in event is simulated with KVMultiDetArray::DetectEvent,
143 // then the reconstructed detected event is treated by the same identification and calibration
144 // procedures as for experimental data.
145
146 // few event counter print at the beginning to be sure the process started properly
147 // because in case GEMINI decay is used, it sometimes stops (randomly) after few events
148 if ((fEVN <= 10)) Info("Analysis", "%d events processed", (int)fEVN);
149 else if ((fEVN <= 1000) && !(fEVN % 100)) Info("Analysis", "%d events processed", (int)fEVN);
150
151#ifdef DEBUG_FILTER
152 if (fEVN == 2500) memory_check->SetInitStatistics();
153 else if (fEVN > 4995) memory_check->CompareToInit();
154#endif
156#ifdef WITH_GEMINI
157 Int_t iterations = 1;
158 if (fGemini) iterations = fGemDecayPerEvent;
159#endif
160 if (to_be_detected->GetNumber()) to_be_detected->SetParameter("SIMEVENT_NUMBER", (int)to_be_detected->GetNumber());
161 to_be_detected->SetParameter("SIMEVENT_TREE_ENTRY", (int)fTreeEntry);
162#ifdef WITH_GEMINI
163 do {
164 if (fGemini) {
165 try {
167 }
168 catch (...) {
169 continue;
170 }
171 //Copy any parameters associated with simulated event into the Gemini-decayed event
174 }
175#endif
177 if (fNewFrame == "proj") to_be_detected->SetFrame("lab", fProjVelocity);
178 else to_be_detected->SetFrame("lab", fCMVelocity);
179 if (fRotate) {
181 gMultiDetArray->DetectEvent(to_be_detected, fReconEvent, "rotated_frame");
182 }
183 else {
184 gMultiDetArray->DetectEvent(to_be_detected, fReconEvent, "lab");
185 }
186 }
187 else {
188 if (fRotate) {
190 gMultiDetArray->DetectEvent(to_be_detected, fReconEvent, "rotated_frame");
191 }
192 else {
193 gMultiDetArray->DetectEvent(to_be_detected, fReconEvent);
194 }
195 }
198 FillTree();
199#ifdef WITH_GEMINI
200 }
201 while (fGemini && --iterations);
202#endif
203
204
205 return kTRUE;
206}
207
208
209
211
213{
214 gEnv->SetValue(Form("%s.HasCalibIdentInfos", GetOpt("Dataset").Data()), fIdCalMode.Data());
215}
216
217
218
220
224
225
226
242
244{
245 // Select required dataset for filtering (option "Dataset")
246 // Build the associated multidetector geometry.
247 // Calculate C.M. velocity associated with required experimental collision
248 // kinematics (option "System").
249 //
250 // Set the parameters of the detectors according to the required run
251 // if given (option "Run"), or the first run of the given system otherwise.
252 // If ROOT/TGeo geometry is required (option "Geometry"="ROOT"),
253 // build the TGeometry representation of the detector array.
254 //
255 // Open file for filtered data (see KVEventFiltering::OpenOutputFile), which will
256 // be stored in a TTree with name 'ReconstructedEvents', in a branch with name
257 // 'ReconEvent'. The class used for reconstructed events depends on the dataset,
258 // it is given by KVDataSet::GetReconstructedEventClassName().
259
260 if (fDisableCreateTreeFile) return; //when running with PROOF, only execute on workers
261
262 TString dataset = GetOpt("Dataset").Data();
263 if (!gDataSetManager) {
264 gDataSetManager = new KVDataSetManager;
265 gDataSetManager->Init();
266 }
267 gDataSetManager->GetDataSet(dataset)->cd();
268
270 if (IsOptGiven("System")) system = GetOpt("System").Data();
271 KVDBSystem* sys = (gExpDB && system != "" ? gExpDB->GetSystem(system) : nullptr);
272
273 KV2Body* tb = nullptr;
274
276 if (sys) tb = sys->GetKinematics();
277 else if (system != "") {
278 tb = new KV2Body(system);
279 tb->CalculateKinematics();
281 }
282
283 fCMVelocity = (tb ? tb->GetCMVelocity() : TVector3(0, 0, 0));
284 fCMVelocity *= -1.0;
285
286 fProjVelocity = (tb ? tb->GetNucleus(1)->GetVelocity() : TVector3(0, 0, 0));
287 fProjVelocity *= -1.0;
288
289 if (justcreated)
290 delete tb;
291
292 Int_t run = 0;
293 if (IsOptGiven("Run")) {
294 run = GetOpt("Run").Atoi();
295 Info("InitAnalysis", "Run given in options = %d", run);
296 }
297 if (!run) {
298 if (sys) {
299 run = ((KVDBRun*)sys->GetRuns()->First())->GetNumber();
300 Info("InitAnalysis", "Using first run for system = %d", run);
301 }
302 else {
303 Info("InitAnalysis", "No run information");
304 run = -1;
305 }
306 }
307
308 fIdCalMode = gEnv->GetValue(Form("%s.HasCalibIdentInfos", dataset.Data()), "no");
309
310 TString filt = GetOpt("Filter").Data();
311 if (filt == "GeoThresh") {
312 // in geo+thresholds mode, we ignore any calibration/identification parameters for the dataset,
313 // i.e. all detectors and all identification ntelescopes are supposed to work
314 gEnv->SetValue(Form("%s.HasCalibIdentInfos", dataset.Data()), "no");
315 }
316
318 if (run == -1) gMultiDetArray->InitializeIDTelescopes();
319 gMultiDetArray->SetSimMode();
320
321 TString geo = GetOpt("Geometry").Data();
322 if (geo == "ROOT") {
323 gMultiDetArray->CheckROOTGeometry();
324 Info("InitAnalysis", "Filtering with ROOT geometry");
325 Info("InitAnalysis", "Navigator detector name format = %s", gMultiDetArray->GetNavigator()->GetDetectorNameFormat());
326 }
327 else {
328 gMultiDetArray->SetROOTGeometry(kFALSE);
329 Info("InitAnalysis", "Filtering with KaliVeda geometry");
330 }
331
332 if (filt == "Geo") {
334 Info("InitAnalysis", "Geometric filter");
335 }
336 else if (filt == "GeoThresh") {
338 Info("InitAnalysis", "Geometry + thresholds filter");
339 }
340 else if (filt == "Full") {
342 Info("InitAnalysis", "Geometry + thresholds filter + implemented identifications and calibrations.");
343 }
344 TString kine = GetOpt("Kinematics").Data();
345 if (kine == "lab") {
347 Info("InitAnalysis", "Simulation is in laboratory/detector frame");
348 }
349 else {
350 fNewFrame = kine;
351 Info("InitAnalysis", "Simulation will be transformed to laboratory/detector frame");
352 }
353
354 if (IsOptGiven("PhiRot")) {
355 if (GetOpt("PhiRot") == "no") fRotate = kFALSE;
356 }
357 if (fRotate) Info("InitAnalysis", "Random phi rotation around beam axis performed for each event");
358#ifdef WITH_GEMINI
359 if (IsOptGiven("Gemini")) {
360 if (GetOpt("Gemini") == "yes") fGemini = kTRUE;
362 if (IsOptGiven("GemDecayPerEvent")) fGemDecayPerEvent = GetOpt("GemDecayPerEvent").Atoi();
363 else fGemDecayPerEvent = 1;
364 if (IsOptGiven("GemAddRotEner")) {
365 if (GetOpt("GemAddRotEner") == "yes") fGemAddRotEner = kTRUE;
366 }
367 }
368 if (fGemini) {
369 Info("InitAnalysis", "Statistical decay with Gemini++ for each event before detection");
370 if (fGemDecayPerEvent > 1) Info("InitAnalysis", " -- %d decays per primary event", fGemDecayPerEvent);
371 if (fGemAddRotEner) Info("InitAnalysis", " -- Rotational energy will be added to excitation energy");
372 }
373#endif
374 if (filt == "Full" && !gDataSet->HasCalibIdentInfos()) {
375 // for old data without implementation of calibration/identification
376 // we set status of identifications according to experimental informations in IdentificationBilan.dat
377 // this "turns off" any telescopes which were not working during the experimental run
378
379 // beforehand we have to "turn on" all identification telescopes by initializing them
380 // as this will not have been done yet
381 gMultiDetArray->InitializeIDTelescopes();
382 // this will also set informations on (simulated) mass identification capabilities for
383 // certain telescopes
384
386 if (sys && KVBase::SearchKVFile("IdentificationBilan.dat", fullpath, gDataSet->GetName())) {
387 Info("InitAnalysis", "Setting identification bilan...");
388 TString sysname = sys->GetBatchName();
390 TIter next(gMultiDetArray->GetListOfIDTelescopes());
392 while ((idt = (KVIDTelescope*)next())) {
394 }
395 }
396 }
397 gMultiDetArray->PrintStatusOfIDTelescopes();
398
399 OpenOutputFile(sys, run);
400 TTree* t{nullptr};
401 if (sys) t = AddTree("ReconstructedEvents", Form("%s filtered with %s (%s)", GetOpt("SimTitle").Data(), gMultiDetArray->GetTitle(), sys->GetName()));
402 else t = AddTree("ReconstructedEvents", Form("%s filtered with %s", GetOpt("SimTitle").Data(), gMultiDetArray->GetTitle()));
403
405 fReconEvent = (KVReconstructedEvent*)TClass::GetClass(reconevclass)->New();
406 KVEvent::MakeEventBranch(t, "ReconEvent", fReconEvent);
407}
408
409
410
412
414{
415 fEVN = 0;
416#ifdef DEBUG_FILTER
418#endif
419}
420
421
422
445
447{
448 // Open ROOT file for new filtered events TTree.
449 // The file will be written in the directory given as option "OutputDir".
450 // The filename is built up from the original simulation filename and the values
451 // of various options:
452 //
453 // [simfile]_[Gemini]_geo=[geometry]_filt=[filter-type]_[dataset]_[system]_run=[run-number].root
454 //
455 // In addition, informations on the filtered data are stored in the file as
456 // TNamed objects. These can be read by KVSimDir::AnalyseFile:
457 //
458 // KEY: TNamed System;1 title=[full system name]
459 // KEY: TNamed Dataset;1 title=[dataset name]
460 // KEY: TNamed Run;1 title=[run-number]
461 // KEY: TNamed Geometry;1 title=[geometry-type]
462 // KEY: TNamed Filter;1 title=[filter-type]
463 // KEY: TNamed Origin;1 title=[name of simulation file]
464 // KEY: TNamed RandomPhi;1 title=[yes/no, random rotation about beam axis]
465 // KEY: TNamed Gemini++;1 title=[yes/no, Gemini++ decay before detection]
466 // KEY: TNamed GemDecayPerEvent;1 title=[number of Gemini++ decays per primary event]
467 // KEY: TNamed GemAddRotEner;1 title=[Enable or not the addition of the rotational energy to the excitation energy]
468 //
469 TString basefile = GetOpt("SimFileName");
470 basefile.Remove(basefile.Index(".root"), 5);
472#ifdef WITH_GEMINI
473 if (fGemini) {
474 outfile += "_Gemini";
476 if (fGemAddRotEner) outfile += "AddedRotEner";
477 }
478#endif
479 outfile += "_geo=";
480 outfile += GetOpt("Geometry");
481 outfile += "_filt=";
482 outfile += GetOpt("Filter");
483 outfile += "_";
484 outfile += gDataSet->GetName();
485
486 if (S) {
487 outfile += "_";
488 outfile += S->GetBatchName();
489 }
490 else if (GetOpt("System")) {
491 TString tmp = GetOpt("System");
492 tmp.ReplaceAll(" ", "");
493 tmp.ReplaceAll("@", "");
494 tmp.ReplaceAll("MeV/A", "");
495 tmp.ReplaceAll("+", "");
496 outfile += "_";
497 outfile += tmp.Data();
498 }
499 if (run > 0) {
500 outfile += "_run=";
501 outfile += Form("%d", run);
502 }
503 outfile += ".root";
504
507
509
511 writeFile->cd();
512 if (S) {
513 TNamed* system = new TNamed("System", S->GetName());
514 system->Write();
515 }
516 else if (GetOpt("System"))(new TNamed("System", GetOpt("System").Data()))->Write();
517
518 (new TNamed("Dataset", gDataSet->GetName()))->Write();
519 if (run > 0)(new TNamed("Run", Form("%d", run)))->Write();
520 if (gMultiDetArray->IsROOTGeometry()) {
521 (new TNamed("Geometry", "ROOT"))->Write();
522 }
523 else {
524 (new TNamed("Geometry", "KV"))->Write();
525 }
526 (new TNamed("Filter", GetOpt("Filter").Data()))->Write();
527 (new TNamed("Origin", (basefile + ".root").Data()))->Write();
528 (new TNamed("RandomPhi", (fRotate ? "yes" : "no")))->Write();
529#ifdef WITH_GEMINI
530 (new TNamed("Gemini++", (fGemini ? "yes" : "no")))->Write();
531 (new TNamed("GemDecayPerEvent", Form("%d", fGemDecayPerEvent)))->Write();
532 (new TNamed("GemAddRotEner", (fGemAddRotEner ? "yes" : "no")))->Write();
533#endif
534 curdir->cd();
535}
536
537
int Int_t
ROOT::R::TRInterface & r
bool Bool_t
constexpr Bool_t kFALSE
double Double_t
constexpr Bool_t kTRUE
#define gDirectory
R__EXTERN TEnv * gEnv
R__EXTERN TRandom * gRandom
char * Form(const char *fmt,...)
void AssignAndDelete(TString &target, char *tobedeleted)
R__EXTERN TSystem * gSystem
Relativistic binary kinematics calculator.
Definition KV2Body.h:166
virtual void SetNumber(UInt_t num)
Definition KVBase.h:216
static Bool_t SearchKVFile(const Char_t *name, TString &fullpath, const Char_t *kvsubdir="")
Definition KVBase.cpp:538
static KVClassMonitor * GetInstance()
Return pointer to unique instance of class monitor class.
Description of an experimental run in database ,,.
Definition KVDBRun.h:36
Database class used to store information on different colliding systems studied during an experiment....
Definition KVDBSystem.h:52
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.
virtual const Char_t * GetReconstructedEventClassName() const
Definition KVDataSet.h:217
Bool_t HasCalibIdentInfos() const
Definition KVDataSet.h:232
void cd() const
Filter simulated events with multidetector response.
TString fNewFrame
allow the definition of a specific frame
KVSimEvent fGemEvent
event after decay with Gemini
void OpenOutputFile(KVDBSystem *, Int_t)
Bool_t fTransformKinematics
=kTRUE if simulation not in lab frame
void RandomRotation(KVEvent *to_rotate, const TString &frame_name="") const
Bool_t fGemini
true if Gemini++ decay should be performed before detection [default: no]
void Copy(TObject &) const
TString fIdCalMode
original exp setup hasIDandCalib to be reset in case of modifications
KVReconstructedEvent * fReconEvent
KVEventFiltering()
Default constructor.
Bool_t fGemAddRotEner
true if rotational energy has to be added to excitation energy [default: no]
Int_t fGemDecayPerEvent
number of Gemini++ decays to be performed for each event [default:1]
virtual ~KVEventFiltering()
Destructor.
Long64_t fEVN
event number counter
Bool_t fRotate
true if random phi rotation should be applied [default: yes]
General purpose analysis base class for TTree containing KVEvent objects.
Bool_t fDisableCreateTreeFile
used with PROOF
KVEvent * GetEvent() const
Long64_t fTreeEntry
current tree entry number
void AddTree(TTree *tree)
void SetEventsReadInterval(Long64_t N)
void SetJobOutputFileName(const TString &filename)
void AddHisto(TH1 *histo)
Bool_t IsOptGiven(const Char_t *option)
TString GetOpt(const Char_t *option) const
void FillTree(const Char_t *tree_name="")
Abstract base class container for multi-particle events.
Definition KVEvent.h:67
static void MakeEventBranch(TTree *tree, const TString &branchname, T &event, Int_t bufsize=10000000)
Definition KVEvent.h:210
KVNameValueList * GetParameters() const
Definition KVEvent.h:179
virtual KVDBSystem * GetSystem(const Char_t *system) const
Definition KVExpDB.h:85
void DecayEvent(const KVSimEvent *, KVSimEvent *, bool=true)
Definition KVGemini.cpp:173
const Char_t * GetDetectorNameFormat() const
Base class for all detectors or associations of detectors in array which can identify charged particl...
void CheckIdentificationBilan(const TString &system)
Set status of ID Telescope for given system.
static void OpenIdentificationBilan(const TString &path)
Open IdentificationBilan.dat file with given path.
Bool_t IsROOTGeometry() const
virtual void DetectEvent(KVEvent *event, KVReconstructedEvent *rec_event, const Char_t *detection_frame="")
void SetFilterType(Int_t t)
KVGeoNavigator * GetNavigator() const
KVSeqCollection * GetListOfIDTelescopes() const
static KVMultiDetArray * MakeMultiDetector(const Char_t *dataset_name, Int_t run=-1, TString classname="KVMultiDetArray")
virtual void InitializeIDTelescopes()
virtual void SetROOTGeometry(Bool_t on=kTRUE)
virtual void SetSimMode(Bool_t on=kTRUE)
void Copy(TObject &nvl) const
Event containing KVReconstructedNucleus nuclei reconstructed from hits in detectors.
Container class for simulated nuclei, KVSimNucleus.
Definition KVSimEvent.h:22
void SetFrame(const Char_t *frame, const KVFrameTransform &ft)
void SetFrameName(const KVString &name)
const char * GetName() const override
void * New(ENewType defConstructor=kClassNew, Bool_t quiet=kFALSE) const
Bool_t cd() override
virtual const char * GetValue(const char *name, const char *dflt) const
virtual void SetValue(const char *name, const char *value, EEnvLevel level=kEnvChange, const char *type=nullptr)
const char * GetName() const override
const char * GetTitle() const override
virtual Int_t Write(const char *name=nullptr, Int_t option=0, Int_t bufsize=0)
virtual void Copy(TObject &object) const
virtual void Info(const char *method, const char *msgfmt,...) const
virtual Double_t Uniform(Double_t x1, Double_t x2)
Int_t Atoi() const
const char * Data() const
virtual char * ConcatFileName(const char *dir, const char *name)
RooArgSet S(Args_t &&... args)
constexpr Double_t TwoPi()
ClassImp(TPyArg)