BmnRoot
Loading...
Searching...
No Matches
run_reco_bmn.cxx
Go to the documentation of this file.
1#if !defined(__CLING__) || defined(__MAKECLING__)
2// ROOT includes
3#include "TKey.h"
4#include "TList.h"
5#include "TPRegexp.h"
6#include "TROOT.h"
7#include "TRandom3.h"
8#include "TStopwatch.h"
9#include "TString.h"
10#include "TSystem.h"
11
12// Fair includes
13#include "FairCave.h"
14#include "FairField.h"
15#include "FairFileSource.h"
16#include "FairMagnet.h"
17#include "FairParAsciiFileIo.h"
18#include "FairParRootFileIo.h"
19#include "FairPipe.h"
20#include "FairRootFileSink.h"
21#include "FairRunAna.h"
22#include "FairRuntimeDb.h"
23#include "FairTask.h"
24#include "FairTrackParP.h"
25
26// BM@N includes
27#include "BmnBeamTracking.h"
28#include "BmnCSCHitMaker.h"
29#include "BmnCombineVertexFinder.h"
30#include "BmnDchHitProducer.h"
31#include "BmnDchTrackFinder.h"
32#include "BmnEventSelector.h"
33#include "BmnFHCalReconstructor.h"
34#include "BmnFieldMap.h"
35#include "BmnFileSource.h"
36#include "BmnFillDstTask.h"
37#include "BmnFillTrigInfoDst.h"
38#include "BmnFunctionSet.h"
39#include "BmnGemResiduals.h"
40#include "BmnGemStripHitMaker.h"
41#include "BmnGemTrack.h"
42#include "BmnGlobalTracking.h"
43#include "BmnHodoReconstructor.h"
44#include "BmnInnerTrackingRun7.h"
45#include "BmnKFPrimaryVertexFinder.h"
46#include "BmnNewFieldMap.h"
47// #include "BmnLANDHitProducer.h"
48#include "BmnMatchRecoToMC.h"
49#include "BmnMwpcHitFinder.h"
50#include "BmnMwpcHitProducer.h"
51#include "BmnMwpcTrackFinder.h"
52#include "BmnNdetReconstructor.h"
53#include "BmnPid.h"
54#include "BmnScWallReconstructor.h"
55#include "BmnSiBTHitMaker.h"
56#include "BmnSiResiduals.h"
57#include "BmnSiliconHitMaker.h"
58#include "BmnStsMatchTracks.h"
59#include "BmnStsVectorFinder.h"
60#include "BmnToCbmHitConverter.h"
61#include "BmnTof1HitProducer.h"
62#include "BmnTof701HitProducer.h"
63#include "BmnTofHitProducer.h"
64#include "BmnVacWallReconstructor.h"
65#include "BmnVertexFinder.h"
66#include "BmnZdcAnalyzer.h"
67#include "CbmFindPrimaryVertex.h"
68#include "CbmKF.h"
69#include "CbmL1.h"
70#include "CbmL1StsTrackFinder.h"
71#include "CbmPVFinderKF.h"
72#include "CbmPrimaryVertexFinder.h"
73#include "CbmStsFindTracks.h"
74#include "CbmStsTrackFinder.h"
75#include "CentralityClusterizer.h"
76#include "DstEventHeader.h"
77#include "MpdGetNumEvents.h"
78#include "UniRun.h"
79
80#include <iostream>
81using namespace std;
82#endif
83
84// Macro for reconstruction of simulated and experimental BM@N events
85//
86// inFileName - input file with data (MC or exp. data)
87// dstFileName - output file with reconstructed data (full DST format)
88// nStartEvent - number of first event to process (starts with zero), default: 0
89// nEvents - number of events to process after nStartEvent, 0 - all events of the given file after nStartEvent to
90R__ADD_INCLUDE_PATH($VMCWORKDIR)
97
98void run_reco_bmn(TString inFileName = "$VMCWORKDIR/macro/run/bmnsim.root",
99 TString dstFileName = "$VMCWORKDIR/macro/run/bmndst.root",
100 Int_t nStartEvent = 0,
101 Int_t nEvents = 10,
102 enumTracking tracking = VectorFinder)
103{
104 gDebug = 0; // ROOT debug option
105 // Verbosity level (0 = summary only, 1 = event-level, 2 = track-level, 3 = full debug)
106 Int_t iVerbose = 0;
107 // Log level for LOG(severity): "fatal", "error", "alarm", "important", "warn", "state", "info", "detail", "debug",
108 // "debug[1-4]", "trace", "nolog"
109 FairLogger::GetLogger()->SetLogScreenLevel("info");
110
111 // ----- Timer ---------------------------------------------------------
112 TStopwatch timer;
113 timer.Start();
114
115 // check input file exists
116 if (!BmnFunctionSet::CheckFileExist(inFileName, 1))
117 exit(-1);
118
119 // ----- Reconstruction run --------------------------------------------
120 FairRunAna* fRunAna = new FairRunAna();
121 fRunAna->SetEventHeader(new DstEventHeader());
122
123 Bool_t isTarget = kTRUE; // flag for tracking (run with target or not)
124 Bool_t isExp = (BmnFunctionSet::isSimulationFile(inFileName)
125 == 0); // flag for FairRunAna tasks whether simulation or experimental data are used
126
127 // Declare input source as simulation file or experimental data
128 FairSource* fFileSource;
129
130 // -1 means use of the BM@N-setup when processing MC-input
131 // DO NOT change it manually!
132 Int_t run_period = 8, run_number = -1;
133 Double_t fieldScale = 0.;
134 if (!isExp) // for simulation files
135 fFileSource = new FairFileSource(inFileName);
136 else // for experimental files
137 {
138 // set source as raw root data file (without additional directories)
139 fFileSource = new BmnFileSource(inFileName, run_period, run_number);
140
141 // get geometry for run
142 gRandom->SetSeed(0);
143 TString geoFileName = Form("current_geo_file_%d.root", UInt_t(gRandom->Integer(UINT32_MAX)));
144 Int_t res_code = UniRun::ReadGeometryFile(run_period, run_number, (char*)geoFileName.Data());
145 if (res_code != 0) {
146 cout << "ERROR: could not read geometry file from the database" << endl;
147 exit(-3);
148 }
149
150 // get gGeoManager from ROOT file (if required)
151 TFile* geoFile = new TFile(geoFileName, "READ");
152 if (!geoFile->IsOpen()) {
153 cout << "ERROR: could not open ROOT file with geometry: " + geoFileName << endl;
154 exit(-4);
155 }
156 TList* keyList = geoFile->GetListOfKeys();
157 TIter next(keyList);
158 TKey* key = (TKey*)next();
159 TString className(key->GetClassName());
160 if (className.BeginsWith("TGeoManager"))
161 key->ReadObj();
162 else {
163 cout << "ERROR: TGeoManager is not top element in geometry file " + geoFileName << endl;
164 exit(-5);
165 }
166
167 // set magnet field with factor corresponding to the given run
168 UniRun* pCurrentRun = UniRun::GetRun(run_period, run_number);
169 if (pCurrentRun == 0)
170 exit(-6);
171 Double_t* field_voltage = pCurrentRun->GetFieldVoltage();
172 if (field_voltage == NULL) {
173 cout << "ERROR: no field voltage was found for run " << run_period << ":" << run_number << endl;
174 exit(-7);
175 }
176 Double_t map_current = (run_period == 8) ? 112.0 : (run_period == 7) ? 55.87 : 0.0;
177 if (*field_voltage < 10) {
178 fieldScale = 0;
179 } else
180 fieldScale = (*field_voltage) / map_current;
181
182 BmnFieldMap* magField = new BmnNewFieldMap("FieldMap_1900_extrap_noPed.root");
183 magField->SetScale(fieldScale);
184 magField->Init();
185 fRunAna->SetField(magField);
186 isExp = kTRUE;
187 TString targ = "-";
188 if (pCurrentRun->GetTargetParticle() == "") {
189 isTarget = kFALSE;
190 } else {
191 targ = pCurrentRun->GetTargetParticle();
192 isTarget = kTRUE;
193 }
194 TString beam = pCurrentRun->GetBeamParticle();
195 if (beam == "")
196 beam = "-";
197 cout << "\n\n|||||||||||||||| EXPERIMENTAL RUN SUMMARY ||||||||||||||||" << endl;
198 cout << "||\t\t\t\t\t\t\t||" << endl;
199 cout << "||\t\tPeriod:\t\t" << run_period << "\t\t\t||" << endl;
200 cout << "||\t\tNumber:\t\t" << run_number << "\t\t\t||" << endl;
201 cout << "||\t\tBeam:\t\t" << beam << "\t\t\t||" << endl;
202 cout << "||\t\tTarget:\t\t" << targ << "\t\t\t||" << endl;
203 cout << "||\t\tField scale:\t" << setprecision(4) << fieldScale << "\t\t\t||" << endl;
204 cout << "||\t\t\t\t\t\t\t||" << endl;
205 cout << "||||||||||||||||||||||||||||||||||||||||||||||||||||||||||\n\n" << endl;
206 remove(geoFileName.Data());
207 }
208
209 fRunAna->SetSource(fFileSource);
210
211 if (dstFileName == "")
212 dstFileName = "bmndst.root";
213 // if directory for the output file does not exist, then create
214 if (BmnFunctionSet::CreateDirectoryTree(dstFileName, 1) < 0)
215 exit(-2);
216 fRunAna->SetSink(new FairRootFileSink(dstFileName));
217 fRunAna->SetGenerateRunInfo(false);
218
219 // if nEvents is equal 0 then all events of the given file starting with "nStartEvent" should be processed
220 if (nEvents == 0)
221 nEvents = MpdGetNumEvents::GetNumROOTEvents((char*)inFileName.Data()) - nStartEvent;
222
223 // Digitisation files.
224 // Add TObjectString file names to a TList which is passed as input to the FairParAsciiFileIo.
225 // The FairParAsciiFileIo will create on fly a concatenated input parameter file, which is then used during the
226 // reconstruction.
227 TList* parFileNameList = new TList();
228
229 // ====================================================================== //
230 // === Pileup Event Classification === //
231 // ====================================================================== //
232 if (run_period >= 8 && isExp) {
233 BmnEventSelector* eventSelector = new BmnEventSelector();
234 eventSelector->SetRunId(run_number);
235 eventSelector->SetInputFileName("$VMCWORKDIR/input/eventSelector_calib_run8.root");
236 fRunAna->AddTask(eventSelector);
237 }
238
239 // ====================================================================== //
240 // === Fill basic trigger information === //
241 // ====================================================================== //
242 BmnFillTrigInfoDst* trig = new BmnFillTrigInfoDst(isExp);
243 fRunAna->AddTask(trig);
244
245 // ====================================================================== //
246 // === Hits in front of the target === //
247 // ====================================================================== //
248 if (run_period >= 8) {
249 if (isExp) {
250 BmnSiBTHitMaker* sibtHM = new BmnSiBTHitMaker(run_period, run_number, isExp);
251 fRunAna->AddTask(sibtHM);
252 }
253 } else {
254 if (!isExp) {
256 fRunAna->AddTask(mwpcHP);
257 }
258 BmnMwpcHitFinder* mwpcHM = new BmnMwpcHitFinder(isExp, run_period, run_number);
259 fRunAna->AddTask(mwpcHM);
260 }
261
262 // ====================================================================== //
263 // === Silicon hit finder === //
264 // ====================================================================== //
265 BmnSiliconHitMaker* siliconHM = new BmnSiliconHitMaker(run_period, run_number, isExp);
266 if (isExp)
267 siliconHM->SetHitMatching(kFALSE);
268 fRunAna->AddTask(siliconHM);
269
270 // ====================================================================== //
271 // === GEM hit finder === //
272 // ====================================================================== //
273 BmnGemStripHitMaker* gemHM = new BmnGemStripHitMaker(run_period, run_number, isExp);
274 if (isExp)
275 gemHM->SetHitMatching(kFALSE);
276 fRunAna->AddTask(gemHM);
277
278 // ====================================================================== //
279 // === CSC hit finder === //
280 // ====================================================================== //
281 BmnCSCHitMaker* cscHM = new BmnCSCHitMaker(run_period, run_number, isExp);
282 if (!isExp)
283 cscHM->SetCurrentConfig(BmnCSCConfiguration::FullCSC_Run8); // set explicitly
284 if (isExp)
285 cscHM->SetHitMatching(kFALSE);
286 fRunAna->AddTask(cscHM);
287
288 // ====================================================================== //
289 // === TOF400 hit finder === //
290 // ====================================================================== //
291 BmnTof1HitProducer* tof1HP = new BmnTof1HitProducer("TOF1", !isExp, iVerbose, kFALSE);
292 tof1HP->SetPeriodRun(run_period, run_number);
293 fRunAna->AddTask(tof1HP);
294
295 // ====================================================================== //
296 // === TOF700 hit finder === //
297 // ====================================================================== //
298 BmnTof701HitProducer* tof701HP = new BmnTof701HitProducer("TOF701", !isExp, iVerbose, kFALSE);
299 tof701HP->SetPeriodRun(run_period, run_number);
300 fRunAna->AddTask(tof701HP);
301
302 // ====================================================================== //
303 // === DCH hit finder === //
304 // ====================================================================== //
305 // if (!isExp) {
306 // BmnDchHitProducer* dchHP = new BmnDchHitProducer();
307 // fRunAna->AddTask(dchHP);
308 // }
309
310 // ====================================================================== //
311 // === FHCAL/ZDC === //
312 // ====================================================================== //
313 bool IsGlobalCoordinates = true;
314 if (run_period >= 8) {
315 // TODO check isExp
316 BmnFHCalReconstructor* fhcalReco = new BmnFHCalReconstructor(isExp, IsGlobalCoordinates);
317 fhcalReco->SetRecoCutsFile("FHCAL_reco_cuts_period8.txt");
318 fRunAna->AddTask(fhcalReco);
319 } else {
320 BmnZdcAnalyzer* zdcAna = new BmnZdcAnalyzer();
321 fRunAna->AddTask(zdcAna);
322 }
323
324 // ====================================================================== //
325 // === ScWall === //
326 // ====================================================================== //
327 IsGlobalCoordinates = true;
328 BmnScWallReconstructor* scwallReco = new BmnScWallReconstructor(isExp, IsGlobalCoordinates);
329 scwallReco->SetRecoCutsFile("SCWALL_reco_cuts_period8.txt");
330 fRunAna->AddTask(scwallReco);
331
332 // ====================================================================== //
333 // === Hodo === //
334 // ====================================================================== //
335
336 BmnHodoReconstructor* hodoReco = new BmnHodoReconstructor("", isExp);
337 hodoReco->SetRecoCutsFile("HODO_reco_cuts_period8.txt");
338 fRunAna->AddTask(hodoReco);
339
340 // ====================================================================== //
341 // === Ndet === //
342 // ====================================================================== //
343 IsGlobalCoordinates = true;
344 BmnNdetReconstructor* ndetReco = new BmnNdetReconstructor(isExp, IsGlobalCoordinates);
345 ndetReco->SetRecoCutsFile("NDET_reco_cuts_period8.txt");
346 fRunAna->AddTask(ndetReco);
347
348 // ====================================================================== //
349 // === Vacuum Wall === //
350 // ====================================================================== //
351
353 fRunAna->AddTask(vacReco);
354
355 // ====================================================================== //
356 // === Inner Tracking === //
357 // ====================================================================== //
358 TString innerTrackBranchName; // use different track container
359 if (tracking == L1Tracking)
360 innerTrackBranchName = "StsTrack";
361 else
362 innerTrackBranchName = "StsVector";
363 if (run_period >= 8) {
364 BmnToCbmHitConverter* hitConverter = new BmnToCbmHitConverter(iVerbose, isExp);
365 // AZ-150224 hitConverter->SetFixedErrors(0.05, 0.25, 0.02, 0.65);
366 hitConverter->SetFixedErrors(0.08 / sqrt(12) * 1.2, 0.1234, 0.01 / sqrt(12), 0.1234); // AZ-150224
367 // hitConverter->EnableAlignment(kFALSE);
368 fRunAna->AddTask(hitConverter);
369
370 CbmKF* kalman = new CbmKF();
371 fRunAna->AddTask(kalman);
372
373 CbmL1* l1 = new CbmL1();
374 TString stsMatBudgetFile = "";
375 l1->SetMaterialBudgetFileName(stsMatBudgetFile);
376 fRunAna->AddTask(l1);
377
378 if (tracking == L1Tracking) {
379 FairTask* stsFindTracks = new CbmStsFindTracks(iVerbose, "CbmL1StsTrackFinder");
380 fRunAna->AddTask(stsFindTracks);
381 }
382 if (tracking == VectorFinder) {
384 vf->SetMatBudgetFileName("$VMCWORKDIR/parameters/sts_matbudget_run8.root");
385 fRunAna->AddTask(vf);
386 }
387
388 BmnStsMatchTracks* stsMatchTracks = new BmnStsMatchTracks(iVerbose);
389 stsMatchTracks->SetTrackBranch(innerTrackBranchName);
390 fRunAna->AddTask(stsMatchTracks);
391
392 } else {
393 BmnInnerTrackingRun7* innerTF = new BmnInnerTrackingRun7(run_number, isTarget);
394 innerTF->SetFiltration(isExp); // we use filtration for experimental data only now
395 fRunAna->AddTask(innerTF);
396 }
397
398 // ====================================================================== //
399 // === Beam Tracking === //
400 // ====================================================================== //
401 if (run_period >= 8) {
402 if (isExp) {
403 BmnBeamTracking* beamTF = new BmnBeamTracking(run_period);
404 fRunAna->AddTask(beamTF);
405 }
406 } else {
407 BmnMwpcTrackFinder* mwpcTF = new BmnMwpcTrackFinder(isExp, run_period, run_number);
408 fRunAna->AddTask(mwpcTF);
409 }
410
411 // ====================================================================== //
412 // === Tracking (DCH) === //
413 // ====================================================================== //
414 // BmnDchTrackFinder* dchTF = new BmnDchTrackFinder(run_period, run_number, isExp);
415 // dchTF->SetTransferFunction("transfer_func.txt");
416 // fRunAna->AddTask(dchTF);
417
418 // ====================================================================== //
419 // === Primary vertex finding === //
420 // ====================================================================== //
421 // if (run_period >= 8) {
422 CbmPrimaryVertexFinder* pvFinder = new CbmPVFinderKF();
423 CbmFindPrimaryVertex* findVertex = new CbmFindPrimaryVertex(pvFinder);
424 findVertex->SetTrackBranch(innerTrackBranchName);
425 fRunAna->AddTask(findVertex);
426
427 // AZ-150224
429 fRunAna->AddTask(pvMPD);
430
432 fRunAna->AddTask(combPV);
433 //
434 // // } else {
435 // BmnVertexFinder* gemVF = new BmnVertexFinder(run_period);
436 // fRunAna->AddTask(gemVF);
437 //}
438
439 // ====================================================================== //
440 // === Global Tracking === //
441 // ====================================================================== //
442 Bool_t doAlign = kTRUE;
443 if (!isExp)
444 doAlign = kFALSE;
445 BmnGlobalTracking* glTF = new BmnGlobalTracking(isExp, kFALSE /*doAlign*/);
446 glTF->SetRunNumber(run_number);
447 glTF->SetInnerTracksBranchName(innerTrackBranchName);
448 fRunAna->AddTask(glTF);
449
450 // BmnVertexFinder* gemVF = new BmnVertexFinder(run_period);
451 // fRunAna->AddTask(gemVF);
452
453 // ====================================================================== //
454 // === Matching global track to MC track procedure === //
455 // ====================================================================== //
456 if (!isExp) {
457 BmnMatchRecoToMC* mcMatching = new BmnMatchRecoToMC();
458 mcMatching->SetInnerTracksBranchName(innerTrackBranchName);
459 fRunAna->AddTask(mcMatching);
460 }
461
462 // ====================================================================== //
463 // === PID procedure === //
464 // ====================================================================== //
465 BmnPid* pid_analyser = new BmnPid();
466 fRunAna->AddTask(pid_analyser);
467
468 // ====================================================================== //
469 // === Residual analysis === //
470 // ====================================================================== //
471 if (run_period < 8) {
472 BmnResiduals* res = new BmnResiduals(run_period, run_number);
473 fRunAna->AddTask(res);
474 }
475
476 // ====================================================================== //
477 // === Centrality === //
478 // ====================================================================== //
479 CentralityClusterizer* cclusterizer = new CentralityClusterizer("XeCsI@3.8AGeV_MBT_7clusters_period8.root");
480 fRunAna->AddTask(cclusterizer);
481
482 // Fill DST Event Header
483 BmnFillDstTask* dst_task = new BmnFillDstTask(nStartEvent, nEvents, isExp, run_period, run_number);
484 fRunAna->AddTask(dst_task);
485
486 // ----- Parameter database (input) -------------------------------------
487 FairRuntimeDb* rtdb = fRunAna->GetRuntimeDb();
488 FairParRootFileIo* parIo1 = new FairParRootFileIo();
489 FairParAsciiFileIo* parIo2 = new FairParAsciiFileIo();
490 parIo1->open(inFileName.Data());
491 parIo2->open(parFileNameList, "in");
492 rtdb->setFirstInput(parIo1);
493 rtdb->setSecondInput(parIo2);
494 // -------------------------------------------------------------------------
495
496 // ----- Initialize FaiRunAna --------------------------------------------
497 fRunAna->GetMainTask()->SetVerbose(iVerbose);
498 fRunAna->Init();
499 // -------------------------------------------------------------------------
500
501 // ----- Parameter database (output) ------------------------------------
502 FairParRootFileIo* parOut = new FairParRootFileIo(kTRUE);
503 // parOut->open(dstFileName);
504 parOut->open(gFile);
505 rtdb->setOutput(parOut);
506 rtdb->writeContainer(rtdb->getContainer("FairGeoParSet"), rtdb->getCurrentRun());
507 rtdb->saveOutput();
508 // rtdb->print();
509 // -------------------------------------------------------------------------
510
511 // ----- Run FairRunAna -------------------------------------------------
512 fRunAna->Run(nStartEvent, nStartEvent + nEvents);
513 // -------------------------------------------------------------------------
514
515 // ----- Cleanup --------------------------------------------------------
516 delete fRunAna;
517 delete pvFinder;
518 delete parFileNameList;
519 delete parIo1;
520 delete parIo2;
521 delete parOut;
522 // Global singletons destruction
524 delete TDatabasePDG::Instance();
525
526 // ----- Finish --------------------------------------------------------
527 timer.Stop();
528 Double_t rtime = timer.RealTime();
529 Double_t ctime = timer.CpuTime();
530 cout << endl << endl;
531 cout << "Macro finished successfully." << endl; // marker of successful execution
532 cout << "Input file is " + inFileName << endl;
533 cout << "Output file is " + dstFileName << endl;
534 cout << "Real time " << rtime << " s, CPU time " << ctime << " s" << endl;
535 cout << endl;
536 // ------------------------------------------------------------------------
537}
538
539int main(int argc, char** arg)
540{
541 run_reco_bmn();
542}
friend F32vec4 sqrt(const F32vec4 &a)
Definition P4_F32vec4.h:34
void SetCurrentConfig(BmnCSCConfiguration::CSC_CONFIG config)
void SetHitMatching(Bool_t opt=kTRUE)
void SetInputFileName(const char *name)
void SetRunId(int runId)
Class for BmnFHCalEvent reconstruction (creation) from BmnFHCalDigi {Data} or BmnFHCalDigit {Sim}.
void SetRecoCutsFile(TString reco_cuts_file)
void SetScale(Double_t factor)
Definition BmnFieldMap.h:55
virtual void Init()
static int CreateDirectoryTree(TString &fileName, int iVerbose=0, EAccessMode mode=kWritePermission)
static int isSimulationFile(TString fileName="")
static int CheckFileExist(TString &fileName, int iVerbose=0, EAccessMode mode=kFileExists)
void SetHitMatching(Bool_t opt=kTRUE)
void SetRunNumber(Int_t r)
void SetInnerTracksBranchName(TString name)
Class for BmnHodoEvent reconstruction (creation) from BmnHodoDigi {Data} or BmnHodoDigit {Sim}.
void SetRecoCutsFile(const std::string &cuts)
void SetInnerTracksBranchName(TString name)
Class for BmnNdetEvent reconstruction (creation) from BmnNdetDigi {Data} or BmnNdetDigit {Sim}.
void SetRecoCutsFile(TString reco_cuts_file)
Class for BmnScWallEvent reconstruction (creation) from BmnScWallDigi {Data} or BmnScWallDigit {Sim}.
void SetRecoCutsFile(TString reco_cuts_file)
void SetHitMatching(Bool_t opt=kTRUE)
void SetTrackBranch(TString trBranch="StsTrack")
void SetMatBudgetFileName(TString s)
void SetFixedErrors(Float_t dXgem=0.015, Float_t dYgem=0.058, Float_t dXsil=0.003, Float_t dYsil=0.021)
void SetPeriodRun(Int_t p, Int_t r)
Class for BmnVacWall points info transition on reconstruction stage.
void SetTrackBranch(TString trBranch="StsTrack")
Definition CbmKF.h:28
Definition CbmL1.h:49
void SetMaterialBudgetFileName(TString s)
Definition CbmL1.h:78
static CbmStsDigiScheme * Instance(int version=1)
Class for event centrality determination through clusterization.
static Int_t GetNumROOTEvents(const char *filename, const char *treename="", int iVerbose=2)
static UniRun * GetRun(int period_number, int run_number)
get run from the database
Definition UniRun.cxx:177
static int ReadGeometryFile(int period_number, int run_number, char *geo_file_path, bool usePrevGeometryIfMissing=false)
Definition UniRun.cxx:1105
TString GetBeamParticle()
get beam particle of the current run
Definition UniRun.h:113
double * GetFieldVoltage()
get field voltage of the current run
Definition UniRun.h:125
TString GetTargetParticle()
get target particle of the current run
Definition UniRun.h:115
STL namespace.
void run_reco_bmn(TString inFileName="$VMCWORKDIR/macro/run/bmnsim.root", TString dstFileName="$VMCWORKDIR/macro/run/bmndst.root", Int_t nStartEvent=0, Int_t nEvents=10, enumTracking tracking=VectorFinder)
enumTracking
@ L1Tracking
@ VectorFinder
@ CellAuto
int main()