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