SND@LHC Software
Loading...
Searching...
No Matches
MCEventBuilder Class Reference

#include <MCEventBuilder.h>

Inheritance diagram for MCEventBuilder:
Collaboration diagram for MCEventBuilder:

Public Member Functions

 MCEventBuilder (const std::string &outputFileName, bool saveOnlyFirst25)
 
 ~MCEventBuilder ()
 
virtual InitStatus Init ()
 
virtual void Exec (Option_t *opt)
 
virtual void FinishTask ()
 

Private Member Functions

std::vector< int > OrderedIds (const std::vector< double > &times, double firstTime) const
 
void ProcessEvent ()
 
bool FastNoiseFilterMu_Hits (TClonesArray *muArray)
 
bool FastNoiseFilterMu_Boards (TClonesArray *muArray)
 
bool FastNoiseFilterScifi_Hits (TClonesArray *scifiArray, const std::map< Int_t, std::map< Int_t, std::array< float, 2 > > > &siPMFibres)
 
bool FastNoiseFilterScifi_Boards (TClonesArray *scifiArray, const std::map< Int_t, std::map< Int_t, std::array< float, 2 > > > &siPMFibres)
 
bool AdvancedNoiseFilterScifi (TClonesArray *scifiArray, const std::map< Int_t, std::map< Int_t, std::array< float, 2 > > > &siPMFibres)
 
bool AdvancedNoiseFilterMu (TClonesArray *muArray)
 

Private Attributes

bool fSaveFirst25nsOnly
 
FairMCEventHeader * fInHeader
 
TClonesArray * fInMufiArray
 
TClonesArray * fInSciFiArray
 
TClonesArray * fInMCTrackArray
 
std::string fOutputFileName
 
TFile * fOutFile
 
TTree * fOutTree
 
FairMCEventHeader * fOutHeader
 
TClonesArray * fOutMufiArray
 
TClonesArray * fOutSciFiArray
 
TClonesArray * fOutMCTrackArray
 

Detailed Description

Definition at line 19 of file MCEventBuilder.h.

Constructor & Destructor Documentation

◆ MCEventBuilder()

MCEventBuilder::MCEventBuilder ( const std::string &  outputFileName,
bool  saveOnlyFirst25 
)

Definition at line 61 of file MCEventBuilder.cxx.

63 : FairTask("MCEventBuilder"),
64 fOutputFileName(outputFileName),
65 fSaveFirst25nsOnly(saveFirst25nsOnly),
66 fOutFile(nullptr),
67 fOutTree(nullptr),
68 fInMufiArray(nullptr),
69 fInSciFiArray(nullptr),
70 fInMCTrackArray(nullptr),
71 fInHeader(nullptr),
72 fOutMufiArray(nullptr),
73 fOutSciFiArray(nullptr),
74 fOutMCTrackArray(nullptr),
75 fOutHeader(nullptr)
76{}
FairMCEventHeader * fOutHeader
TClonesArray * fInMCTrackArray
TClonesArray * fInMufiArray
TClonesArray * fInSciFiArray
FairMCEventHeader * fInHeader
TClonesArray * fOutSciFiArray
TClonesArray * fOutMCTrackArray
TClonesArray * fOutMufiArray
std::string fOutputFileName

◆ ~MCEventBuilder()

MCEventBuilder::~MCEventBuilder ( )

Definition at line 78 of file MCEventBuilder.cxx.

78{}

Member Function Documentation

◆ AdvancedNoiseFilterMu()

bool MCEventBuilder::AdvancedNoiseFilterMu ( TClonesArray *  muArray)
private

Definition at line 483 of file MCEventBuilder.cxx.

483 {
484 std::set<int> veto_planes, us_planes, ds_planes;
485 for (int i = 0; i < muArray->GetEntriesFast(); ++i) {
486 auto* p = static_cast<MuFilterPoint*>(muArray->At(i));
487 int detID = p->GetDetectorID();
488 int system = detID / 10000;
489 if (system == 1) {
490 veto_planes.insert(detID / 1000);
491 if (veto_planes.size() >= min_veto_planes) {
492 return true;
493 }
494 } else if (system == 2) {
495 us_planes.insert(detID / 1000);
496 if (us_planes.size() >= min_us_planes) {
497 return true;
498 }
499 } else if (system == 3) {
500 ds_planes.insert(detID / 1000);
501 if (ds_planes.size() >= min_ds_planes) {
502 return true;
503 }
504 }
505 }
506 return false;
507}
int i
Definition ShipAna.py:86

◆ AdvancedNoiseFilterScifi()

bool MCEventBuilder::AdvancedNoiseFilterScifi ( TClonesArray *  scifiArray,
const std::map< Int_t, std::map< Int_t, std::array< float, 2 > > > &  siPMFibres 
)
private

Definition at line 560 of file MCEventBuilder.cxx.

563 {
564 std::set<int> UniquePlanes;
565 for (int i = 0; i < scifiArray->GetEntriesFast(); ++i) {
566 auto* p = static_cast<ScifiPoint*>(scifiArray->At(i));
567 int point_detID = p->GetDetectorID();
568 int locFibreID = point_detID % 100000;
569
570 for (const auto& sipmChan : siPMFibresMap.at(locFibreID)) {
571 int channel = sipmChan.first;
572 int detID_geo = (point_detID / 100000) * 100000 + channel;
573 UniquePlanes.insert(detID_geo / 100000);
574
575 if (UniquePlanes.size() >= min_scifi_planes) {
576 return true;
577 }
578 }
579 }
580 return false;
581 }

◆ Exec()

void MCEventBuilder::Exec ( Option_t *  opt)
virtual

Definition at line 142 of file MCEventBuilder.cxx.

142 {
143 ProcessEvent();
144}

◆ FastNoiseFilterMu_Boards()

bool MCEventBuilder::FastNoiseFilterMu_Boards ( TClonesArray *  muArray)
private

Definition at line 447 of file MCEventBuilder.cxx.

447 {
448 //The digits used here to assign IDs to systems are set to match the DAQ board IDs from the 2022 board_mapping.json file
449 std::set<int> us, ds;
450 for (int i = 0; i < muArray->GetEntriesFast(); ++i) {
451 auto* p = static_cast<MuFilterPoint*>(muArray->At(i));
452 int detID = p->GetDetectorID();
453 int system = detID / 10000;
454 int TypeofPlane = (detID / 1000) % 10;
455 if (system == 2) {
456 if (TypeofPlane == 1 || TypeofPlane == 2) {
457 us.insert(7);
458 } else if (TypeofPlane == 3 || TypeofPlane == 4) {
459 us.insert(60);
460 } else if (TypeofPlane == 5) {
461 us.insert(52);
462 }
463 if (us.size() >= min_us_boards) {
464 return true;
465 }
466 } else if (system == 3) {
467 if (TypeofPlane == 1) {
468 ds.insert(41);
469 } else if (TypeofPlane == 2) {
470 ds.insert(35);
471 } else if (TypeofPlane == 3 || TypeofPlane == 4) {
472 ds.insert(55);
473 }
474 if (ds.size() >= min_ds_boards) {
475 return true;
476 }
477 }
478 }
479 return false;
480}

◆ FastNoiseFilterMu_Hits()

bool MCEventBuilder::FastNoiseFilterMu_Hits ( TClonesArray *  muArray)
private

Definition at line 425 of file MCEventBuilder.cxx.

425 {
426 std::set<int> veto, ds;
427 for (int i = 0; i < muArray->GetEntriesFast(); ++i) {
428 auto* p = static_cast<MuFilterPoint*>(muArray->At(i));
429 int detID = p->GetDetectorID();
430 int system = detID / 10000;
431 if (system == 1) {
432 veto.insert(detID);
433 if (veto.size() >= min_veto_hits) {
434 return true;
435 }
436 }
437 else if (system == 3) {
438 ds.insert(detID);
439 if (ds.size() >= min_ds_hits) {
440 return true;
441 }
442 }
443 }
444 return false;
445}
Definition veto.h:16

◆ FastNoiseFilterScifi_Boards()

bool MCEventBuilder::FastNoiseFilterScifi_Boards ( TClonesArray *  scifiArray,
const std::map< Int_t, std::map< Int_t, std::array< float, 2 > > > &  siPMFibres 
)
private

Definition at line 534 of file MCEventBuilder.cxx.

537 {
538 //A DAQ board reads one Scifi mat and we use the first three digits of the detID_geo, representing the unique combination of StationPlaneMat, to count boards with fired channels.
539 std::set<int> dividedDetIds;
540
541 for (int i = 0; i < scifiArray->GetEntriesFast(); ++i) {
542 auto* p = static_cast<ScifiPoint*>(scifiArray->At(i));
543 int point_detID = p->GetDetectorID();
544 int locFibreID = point_detID % 100000;
545
546 for (const auto& sipmChan : siPMFibresMap.at(locFibreID)) {
547 int channel = sipmChan.first;
548 int detID_geo = (point_detID / 100000) * 100000 + channel;
549 dividedDetIds.insert(detID_geo / 10000);
550
551 if (dividedDetIds.size() >= min_scifi_boards){
552 return true;
553 }
554 }
555 }
556 return false;
557 }

◆ FastNoiseFilterScifi_Hits()

bool MCEventBuilder::FastNoiseFilterScifi_Hits ( TClonesArray *  scifiArray,
const std::map< Int_t, std::map< Int_t, std::array< float, 2 > > > &  siPMFibres 
)
private

Definition at line 511 of file MCEventBuilder.cxx.

514 {
515 std::set<int> detIDs;
516 for (int i = 0; i < scifiArray->GetEntriesFast(); ++i) {
517 auto* p = static_cast<ScifiPoint*>(scifiArray->At(i));
518 int point_detID = p->GetDetectorID();
519 int locFibreID = point_detID % 100000;
520
521 for (const auto& sipmChan : siPMFibresMap.at(locFibreID)) {
522 int channel = sipmChan.first;
523 int detID_geo = (point_detID / 100000) * 100000 + channel;
524 detIDs.insert(detID_geo);
525
526 if (detIDs.size() >= min_scifi_hits){
527 return true;
528 }
529 }
530 }
531 return false;
532 }

◆ FinishTask()

void MCEventBuilder::FinishTask ( )
virtual

Definition at line 583 of file MCEventBuilder.cxx.

583 {
584 // Make the output file structure compatible with FairTasks to facilitate
585 // runing the digitization task.
586 fOutFile->cd();
587 TFolder* folder = new TFolder("cbmroot", "Main Folder");
588 TFolder* fol_Stack = folder->AddFolder("Stack", "Stack Folder");
589 auto mcTrackArr = new TClonesArray("ShipMCTrack");
590 mcTrackArr->SetName("MCTrack"); // avoid default plural names
591 fol_Stack->Add(mcTrackArr);
592 TFolder* fol_veto = folder->AddFolder("veto", "veto Folder");
593 auto vetoArr = new TClonesArray("vetoPoint");
594 vetoArr->SetName("vetoPoint");
595 fol_veto->Add(vetoArr);
596 TFolder* fol_EmulsionDet = folder->AddFolder("EmulsionDet", "EmulsionDetector Folder");
597 auto emuArr = new TClonesArray("EmulsionDetPoint");
598 emuArr->SetName("EmulsionDetPoint");
599 fol_EmulsionDet->Add(emuArr);
600 TFolder* fol_Scifi = folder->AddFolder("Scifi", "Scifi Folder");
601 auto scifiArr = new TClonesArray("ScifiPoint");
602 scifiArr->SetName("ScifiPoint");
603 fol_Scifi->Add(scifiArr);
604 TFolder* fol_MuFilter = folder->AddFolder("MuFilter", "MuFilter Folder");
605 auto mufArr = new TClonesArray("MuFilterPoint");
606 mufArr->SetName("MuFilterPoint");
607 fol_MuFilter->Add(mufArr);
608 TFolder* fol_Event = folder->AddFolder("Event", "Event Folder");
609 folder->Write();
610
611 fOutFile->cd();
612 TList* branch_list = new TList();
613 branch_list->SetName("BranchList");
614 branch_list->Add(new TObjString("MCTrack"));
615 branch_list->Add(new TObjString("vetoPoint"));
616 branch_list->Add(new TObjString("EmulsionDetPoint"));
617 branch_list->Add(new TObjString("ScifiPoint"));
618 branch_list->Add(new TObjString("MuFilterPoint"));
619 branch_list->Add(new TObjString("MCEventHeader."));
620 fOutFile->WriteObject(branch_list, "BranchList");
621 auto* timebased_branch_list = new TList();
622 timebased_branch_list->SetName("TimeBasedBranchList");
623 fOutFile->WriteObject(timebased_branch_list, "TimeBasedBranchList");
624 fOutFile->cd();
625 auto* fileHeader = new FairFileHeader();
626 fileHeader->SetTitle("FileHeader");
627 fileHeader->Write("FileHeader");
628 fOutFile->cd();
629 if (fOutTree) {
630 fOutTree->Write();
631 }
632 LOG(INFO) << "Writing and closing output file: " << fOutputFileName;
633}

◆ Init()

InitStatus MCEventBuilder::Init ( )
virtual

Definition at line 80 of file MCEventBuilder.cxx.

80 {
81 LOG(INFO) << "Initializing MCEventBuilder";
82
83 FairRootManager* ioman = FairRootManager::Instance();
84 if (!ioman) {
85 LOG(FATAL) << "MCEventBuilder::Init: RootManager not instantiated!";
86 }
87
88 // —---------- INPUT BRANCHES —-----------
89 // Only standard ROOT way of reading branches works for the header, otherwise nullptr
90 ioman->GetInTree()->SetBranchAddress("MCEventHeader.",&fInHeader);
91 fInMufiArray = static_cast<TClonesArray*>(ioman->GetObject("MuFilterPoint"));
92 fInSciFiArray = static_cast<TClonesArray*>(ioman->GetObject("ScifiPoint"));
93 fInMCTrackArray = static_cast<TClonesArray*>(ioman->GetObject("MCTrack"));
94
95 if (!fInMufiArray && !fInSciFiArray) {
96 LOG(ERROR) << "No Scifi and no MuFilter MC points array!";
97 return kERROR;
98 }
99
100 // —---------- OUTPUT FILE & TREE —----------
101 fOutFile = new TFile(fOutputFileName.c_str(), "RECREATE");
102 fOutTree = new TTree("cbmsim", "RebuiltEvents");
103 fOutTree->SetTitle("/cbmroot_0");
104
105 fOutHeader = new FairMCEventHeader();
106 fOutMufiArray = new TClonesArray("MuFilterPoint");
107 fOutSciFiArray = new TClonesArray("ScifiPoint");
108 fOutMCTrackArray = new TClonesArray("ShipMCTrack");
109 fOutMufiArray->SetName("MuFilterPoint");
110 fOutSciFiArray->SetName("ScifiPoint");
111 fOutMCTrackArray->SetName("ShipMCTrack");
112
113 fOutTree->Branch("MuFilterPoint", &fOutMufiArray, 32000, 1);
114 fOutTree->Branch("ScifiPoint", &fOutSciFiArray, 32000, 1);
115 fOutTree->Branch("MCTrack", &fOutMCTrackArray, 32000, 1);
116 fOutTree->Branch("MCEventHeader.", &fOutHeader);
117
118 // --------------Global Variables--------------------
119 //Muon Filter
120 MuFilterDet = dynamic_cast<MuFilter*>(gROOT->GetListOfGlobals()->FindObject("MuFilter"));
121 if (!MuFilterDet) {
122 LOG(ERROR) << "MuFilter detector not found in gROOT";
123 return kERROR;
124 }
125 DsPropSpeed = MuFilterDet->GetConfParF("MuFilter/DsPropSpeed");
126 VandUpPropSpeed = MuFilterDet->GetConfParF("MuFilter/VandUpPropSpeed");
127
128 //Scifi
129 ScifiDet = dynamic_cast<Scifi*>(gROOT->GetListOfGlobals()->FindObject("Scifi"));
130 if (!ScifiDet) {
131 LOG(ERROR) << "Scifi detector not found in gROOT";
132 return kERROR;
133 }
134 ScifisignalSpeed = ScifiDet->GetConfParF("Scifi/signalSpeed");
135 ScifiDet->SiPMmapping();
136 siPMFibres = ScifiDet->GetFibresMap();
137
138 LOG(INFO) << "MCEventBuilder initialized successfully.";
139 return kSUCCESS;
140}
Definition Scifi.h:20
Float_t GetConfParF(TString name)
Definition Scifi.h:50

◆ OrderedIds()

std::vector< int > MCEventBuilder::OrderedIds ( const std::vector< double > &  times,
double  firstTime 
) const
private

Definition at line 146 of file MCEventBuilder.cxx.

147 {
148 size_t n = times.size();
149 std::vector<long long> bins(n);
150 for (size_t i = 0; i < n; ++i) {
151 bins[i] = static_cast<long long>((times[i] - firstTime) / timeWindow);
152 }
153
154 std::vector<int> ids(n, 0);
155 for (size_t i = 1; i < n; ++i) {
156 if ((bins[i] - bins[i - 1]) < 1) {
157 ids[i] = ids[i - 1];
158 } else {
159 ids[i] = ids[i - 1] + 1;
160 }
161 }
162 return ids;
163}
real(mps), dimension(0:8) times
cpu time counters
Definition mpmod.f90:134

◆ ProcessEvent()

void MCEventBuilder::ProcessEvent ( )
private

Definition at line 165 of file MCEventBuilder.cxx.

165 {
166
167 fOutHeader->SetRunID(fInHeader->GetRunID());
168 fOutHeader->SetEventID(fInHeader->GetEventID());
169 fOutHeader->SetVertex(fInHeader->GetX(), fInHeader->GetY(), fInHeader->GetZ());
170 fOutHeader->SetTime(fInHeader->GetT());
171 fOutHeader->SetB(fInHeader->GetB());
172 fOutHeader->SetNPrim(fInHeader->GetNPrim());
173 fOutHeader->MarkSet(fInHeader->IsSet());
174 fOutHeader->SetRotX(fInHeader->GetRotX());
175 fOutHeader->SetRotY(fInHeader->GetRotY());
176 fOutHeader->SetRotZ(fInHeader->GetRotZ());
177
178 //---------------------------Muon filter-------------------------------------
179 std::vector<MuFilterPoint*> muFilterPoints;
180 std::vector<double> muArrivalTimes;
181 std::vector<int> muTrackIDs;
182
183 for (auto* p : ROOT::RRangeCast<MuFilterPoint*, false, decltype(*fInMufiArray)>(*fInMufiArray)) {
184 muFilterPoints.push_back(p);
185 muTrackIDs.push_back(p->GetTrackID());
186
187 int detID = p->GetDetectorID();
188
189 float propspeed;
190 if (floor(detID / 10000) == 3)
191 propspeed = DsPropSpeed;
192 else
193 propspeed = VandUpPropSpeed;
194
195 TVector3 vLeft,vRight;
196 TVector3 impact(p->GetX(), p->GetY(), p->GetZ());
197 MuFilterDet->GetPosition(detID, vLeft, vRight);
198 TVector3 vTop = vLeft; //Used for vertical bars
199
200 //Vertical bars with only 1 readout at the top
201 if ( (floor(detID/10000)==3&&detID%1000>59) ||
202 (floor(detID/10000)==1&&int(detID/1000)%10==2) ) {
203 double arrivalTime = p->GetTime() + (vTop - impact).Mag() / propspeed;
204 muArrivalTimes.push_back(arrivalTime);
205 }
206 //Horizontal
207 else{
208 double tLeft = p->GetTime() + (vLeft - impact).Mag() / propspeed;
209 double tRight = p->GetTime() + (vRight - impact).Mag() / propspeed;
210 double arrivalTime = std::min(tLeft, tRight);
211 muArrivalTimes.push_back(arrivalTime);
212 }
213 }
214
215 std::vector<size_t> idxM(muArrivalTimes.size());
216 std::iota(idxM.begin(), idxM.end(), 0);
217 std::sort(idxM.begin(), idxM.end(), [&](size_t a, size_t b) {
218 return muArrivalTimes[a] < muArrivalTimes[b];
219 });
220
221 std::vector<MuFilterPoint*> sortedMuPoints;
222 std::vector<double> sortedMuArrivalTimes;
223 std::vector<int> sortedMuTrackIDs;
224
225 sortedMuPoints.reserve(muFilterPoints.size());
226 sortedMuArrivalTimes.reserve(muArrivalTimes.size());
227 sortedMuTrackIDs.reserve(muTrackIDs.size());
228
229 for (auto i : idxM) {
230 sortedMuPoints.push_back(muFilterPoints[i]);
231 sortedMuArrivalTimes.push_back(muArrivalTimes[i]);
232 sortedMuTrackIDs.push_back(muTrackIDs[i]);
233 }
234
235 //---------------------------Scifi-------------------------------------
236 std::vector<ScifiPoint*> scifiPoints;
237 std::vector<double> scifiArrivalTimes;
238 std::vector<int> scifiTrackIDs;
239
240 float signalSpeed = ScifisignalSpeed;
241
242 for (auto* p : ROOT::RRangeCast<ScifiPoint*, false, decltype(*fInSciFiArray)>(*fInSciFiArray)) {
243 scifiPoints.push_back(p);
244 scifiTrackIDs.push_back(p->GetTrackID());
245
246 TVector3 impact(p->GetX(), p->GetY(), p->GetZ());
247 int point_detID = p->GetDetectorID();
248 int localFiberID = (point_detID)%100000;
249 int a_sipmChan = static_cast<int>(siPMFibres[localFiberID].begin()->first);
250 int detID_geo = int(point_detID/100000)*100000+a_sipmChan;
251
252 TVector3 a, b;
253 ScifiDet->GetSiPMPosition(detID_geo, a, b);
254 bool verticalHit = int(detID_geo / 100000) % 10 == 1;
255 double distance;
256 if (verticalHit) {
257 distance = (b - impact).Mag();
258 } else {
259 distance = (impact - a).Mag();
260 }
261 double arrivalTime = p->GetTime() + distance / signalSpeed;
262 scifiArrivalTimes.push_back(arrivalTime);
263 }
264
265 std::vector<size_t> idxS(scifiArrivalTimes.size());
266 std::iota(idxS.begin(), idxS.end(), 0);
267 std::sort(idxS.begin(), idxS.end(), [&](size_t a, size_t b) {
268 return scifiArrivalTimes[a] < scifiArrivalTimes[b];
269 });
270
271 std::vector<ScifiPoint*> sortedScifiPoints;
272 std::vector<double> sortedScifiArrivalTimes;
273 std::vector<int> sortedScifiTrackIDs;
274
275 sortedScifiPoints.reserve(scifiPoints.size());
276 sortedScifiArrivalTimes.reserve(scifiArrivalTimes.size());
277 sortedScifiTrackIDs.reserve(scifiTrackIDs.size());
278
279 for (auto i : idxS) {
280 sortedScifiPoints.push_back(scifiPoints[i]);
281 sortedScifiArrivalTimes.push_back(scifiArrivalTimes[i]);
282 sortedScifiTrackIDs.push_back(scifiTrackIDs[i]);
283 }
284
285 //----------------------------------Tracks-------------------------------------
286 std::vector<ShipMCTrack*> mcTrackClones;
287 for (auto* t : ROOT::RRangeCast<ShipMCTrack*, false, decltype(*fInMCTrackArray)>(*fInMCTrackArray)) {
288 mcTrackClones.push_back(t);
289 }
290
291 //---------Finding the earliest time between Scifi and MuFilter-----------------
292 double tMufi = sortedMuArrivalTimes.empty() ? -1 : sortedMuArrivalTimes.front();
293 double tScifi = sortedScifiArrivalTimes.empty() ? -1 : sortedScifiArrivalTimes.front();
294 bool hasMCPoints = (tMufi >= 0 || tScifi >= 0);
295 double firstT = hasMCPoints ? (tMufi < 0 ? tScifi : (tScifi < 0 ? tMufi : std::min(tMufi, tScifi))): 0;
296
297 if (!hasMCPoints) {
298 fOutMufiArray->Delete();
299 fOutSciFiArray->Delete();
300 fOutMCTrackArray->Delete();
301 fOutTree->Fill();
302 return;
303 }
304
305 //------------------Preparations before chunking the data---------------------
306 auto idsMufi = OrderedIds(sortedMuArrivalTimes, firstT);
307 auto idsScifi = OrderedIds(sortedScifiArrivalTimes, firstT);
308
309 bool FirstEvent = true;
310
311 std::vector<int> used;
312 int i_mufi = 0, i_scifi = 0, sliceMufi = 0, sliceScifi = 0;
313
314 while (i_mufi < (int)sortedMuArrivalTimes.size() || i_scifi < (int)sortedScifiArrivalTimes.size()) {
315 // To avoid re-copy of points or tracks from chunk N to following chunks, one needs to free the memory
316 // e.g. chunk N has points with track ids(indices) T-10, T-5, T, chunk N+1 holds points with track id T+1.
317 // If memory is not de-allocated, all tracks of smaller index than T+1 will be also written to chunk N+1.
318 // This will be incorrect since there are no points of tracks with id < T in the N+1 chunck.
319 fOutMufiArray->Delete();
320 fOutSciFiArray->Delete();
321 fOutMCTrackArray->Delete();
322
323 std::vector<MuFilterPoint*> muSlicePoints;
324 std::vector<ScifiPoint*> scifiSlicePoints;
325 fOutMCTrackArray->ExpandCreate(mcTrackClones.size());
326
327 //Adding the mother track to the first event
328 if (sliceMufi==0 && sliceScifi==0){
329 ShipMCTrack* newTrack = new ((*fOutMCTrackArray)[0]) ShipMCTrack(*mcTrackClones[0]);
330 // In NC neutrino events, the outgoing neutrino is also saved in the first chunk,
331 // otherwise it will be lost. This is to facilitate NC event tagging in analysis.
332 if (mcTrackClones.size()>1 && std::set({12,14,16}).count(fabs(mcTrackClones[1]->GetPdgCode()))){
333 ShipMCTrack* newTrack = new ((*fOutMCTrackArray)[1]) ShipMCTrack(*mcTrackClones[1]);
334 }
335 }
336
337 //MUON FILTER POINTS CHUNKING
338 while (i_mufi < (int)sortedMuArrivalTimes.size() && idsMufi[i_mufi] == sliceMufi) {
339 MuFilterPoint* origMu = sortedMuPoints[i_mufi];
340 muSlicePoints.push_back(origMu);
341 Int_t trackID = origMu->GetTrackID();
342 Int_t detID = origMu->GetDetectorID();
343 TVector3 pos; origMu->Position(pos);
344 TVector3 mom; origMu->Momentum(mom);
345 Double_t time = origMu->GetTime();
346 Double_t length = origMu->GetLength();
347 Double_t eLoss = origMu->GetEnergyLoss();
348 Int_t pdgCode = origMu->PdgCode();
349
350 new ((*fOutMufiArray)[fOutMufiArray->GetEntriesFast()])
351 MuFilterPoint(trackID, detID, pos, mom, time, length, eLoss, pdgCode);
352
353 int tid = sortedMuTrackIDs[i_mufi++];
354 //MC tracks having ID=-2 are below the pre-set energy threshold (default is 100MeV) and cannot be linked to their corresponding MC point. Such tracks are not saved to the chunked events.
355 if (tid != -2) {
356 ShipMCTrack* newTrack = new ((*fOutMCTrackArray)[tid]) ShipMCTrack(*mcTrackClones[tid]);
357 }
358 }
359 ++sliceMufi;
360
361 //SCIFI POINTS CHUNKING
362 while (i_scifi < (int)sortedScifiArrivalTimes.size() && idsScifi[i_scifi] == sliceScifi) {
363 ScifiPoint* origSci = sortedScifiPoints[i_scifi];
364 scifiSlicePoints.push_back(origSci);
365 Int_t trackID = origSci->GetTrackID();
366 Int_t detID = origSci->GetDetectorID();
367 TVector3 pos; origSci->Position(pos);
368 TVector3 mom; origSci->Momentum(mom);
369 Double_t time = origSci->GetTime();
370 Double_t length = origSci->GetLength();
371 Double_t eLoss = origSci->GetEnergyLoss();
372 Int_t pdgCode = origSci->PdgCode();
373
374 new ((*fOutSciFiArray)[fOutSciFiArray->GetEntriesFast()])
375 ScifiPoint(trackID, detID, pos, mom, time, length, eLoss, pdgCode);
376
377 int tid = sortedScifiTrackIDs[i_scifi++];
378 //MC tracks having ID=-2 are below the pre-set energy threshold (default is 100MeV) and cannot be linked to their corresponding MC point. Such tracks are not saved to the chunked events.
379 if (tid != -2) {
380 ShipMCTrack* newTrack = new ((*fOutMCTrackArray)[tid]) ShipMCTrack(*mcTrackClones[tid]);
381 }
382 }
383 ++sliceScifi;
384
385 //Noise Filters
387 fOutMufiArray->Delete();
388 fOutSciFiArray->Delete();
389 fOutMCTrackArray->Delete();
390 FirstEvent = false;
391 continue;
392 }
394 fOutMufiArray->Delete();
395 fOutSciFiArray->Delete();
396 fOutMCTrackArray->Delete();
397 FirstEvent = false;
398 continue;
399 }
400
402 fOutMufiArray->Delete();
403 fOutSciFiArray->Delete();
404 fOutMCTrackArray->Delete();
405 FirstEvent = false;
406 continue;
407 }
408
409 //--------Save only first 25ns chunks mode-------
410 if (fSaveFirst25nsOnly) {
411 if (FirstEvent) {
412 fOutTree->Fill();
413 FirstEvent = false;
414 break;
415 }
416 }
417 //-----Save all chunks mode-------
418 else {
419 fOutTree->Fill();
420 }
421 }
422}
bool FastNoiseFilterScifi_Hits(TClonesArray *scifiArray, const std::map< Int_t, std::map< Int_t, std::array< float, 2 > > > &siPMFibres)
bool FastNoiseFilterMu_Boards(TClonesArray *muArray)
bool FastNoiseFilterMu_Hits(TClonesArray *muArray)
std::vector< int > OrderedIds(const std::vector< double > &times, double firstTime) const
bool AdvancedNoiseFilterScifi(TClonesArray *scifiArray, const std::map< Int_t, std::map< Int_t, std::array< float, 2 > > > &siPMFibres)
bool FastNoiseFilterScifi_Boards(TClonesArray *scifiArray, const std::map< Int_t, std::map< Int_t, std::array< float, 2 > > > &siPMFibres)
bool AdvancedNoiseFilterMu(TClonesArray *muArray)
Int_t PdgCode() const
Int_t PdgCode() const
Definition ScifiPoint.h:39

Member Data Documentation

◆ fInHeader

FairMCEventHeader* MCEventBuilder::fInHeader
private

Definition at line 55 of file MCEventBuilder.h.

◆ fInMCTrackArray

TClonesArray* MCEventBuilder::fInMCTrackArray
private

Definition at line 58 of file MCEventBuilder.h.

◆ fInMufiArray

TClonesArray* MCEventBuilder::fInMufiArray
private

Definition at line 56 of file MCEventBuilder.h.

◆ fInSciFiArray

TClonesArray* MCEventBuilder::fInSciFiArray
private

Definition at line 57 of file MCEventBuilder.h.

◆ fOutFile

TFile* MCEventBuilder::fOutFile
private

Definition at line 62 of file MCEventBuilder.h.

◆ fOutHeader

FairMCEventHeader* MCEventBuilder::fOutHeader
private

Definition at line 64 of file MCEventBuilder.h.

◆ fOutMCTrackArray

TClonesArray* MCEventBuilder::fOutMCTrackArray
private

Definition at line 67 of file MCEventBuilder.h.

◆ fOutMufiArray

TClonesArray* MCEventBuilder::fOutMufiArray
private

Definition at line 65 of file MCEventBuilder.h.

◆ fOutputFileName

std::string MCEventBuilder::fOutputFileName
private

Definition at line 61 of file MCEventBuilder.h.

◆ fOutSciFiArray

TClonesArray* MCEventBuilder::fOutSciFiArray
private

Definition at line 66 of file MCEventBuilder.h.

◆ fOutTree

TTree* MCEventBuilder::fOutTree
private

Definition at line 63 of file MCEventBuilder.h.

◆ fSaveFirst25nsOnly

bool MCEventBuilder::fSaveFirst25nsOnly
private

Definition at line 54 of file MCEventBuilder.h.


The documentation for this class was generated from the following files: