CbmRoot
Loading...
Searching...
No Matches
CbmCaIdealHitProducerDet.cxx
Go to the documentation of this file.
1/* Copyright (C) 2023 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Sergei Zharko [committer] */
4
9
11
12//templateClassImp(cbm::ca::IdealHitProducerDet);
13
14using namespace cbm::ca;
15
16// ------------------------------------------------------------------------------------------------------------------
17//
18template<ca::EDetectorID DetID>
23
24// ------------------------------------------------------------------------------------------------------------------
25//
26template<ca::EDetectorID DetID>
28{
29 if (fpBrHitsTmp) {
30 fpBrHitsTmp->Delete();
31 delete fpBrHitsTmp;
32 fpBrHitsTmp = nullptr;
33 }
34
36 fpBrHitMatchesTmp->Delete();
37 delete fpBrHitMatchesTmp;
38 fpBrHitMatchesTmp = nullptr;
39 }
40}
41
42// ------------------------------------------------------------------------------------------------------------------
43//
44template<ca::EDetectorID DetID>
46try {
47 auto CheckBranch = [](auto pBranch, const char* brName) {
48 if (!pBranch) {
49 throw std::logic_error(Form("Branch %s is not available", brName));
50 }
51 };
52
53 // ------ Input branch initialization
54 auto* pFairManager = FairRootManager::Instance();
55 auto* pMcManager = dynamic_cast<CbmMCDataManager*>(pFairManager->GetObject("MCDataManager"));
56 CheckBranch(pMcManager, "MCDataManager");
57
58 fpTimeSlice = dynamic_cast<CbmTimeSlice*>(pFairManager->GetObject("TimeSlice."));
59 fpRecoEvents = dynamic_cast<TClonesArray*>(pFairManager->GetObject("CbmEvent"));
60 fpMCEventList = dynamic_cast<CbmMCEventList*>(pFairManager->GetObject("MCEventList."));
61
62 CheckBranch(fpTimeSlice, "TimeSlice.");
63 CheckBranch(fpMCEventList, "MCEventList.");
64 fpMCEventHeader = pMcManager->GetObject("MCEventHeader.");
65 fpBrMCTracks = pMcManager->InitBranch("MCTrack");
66 CheckBranch(fpBrMCTracks, "MCTrack");
67
68 const auto* pRecoSetupManager = cbm::RecoSetupManager::Instance();
69 if (!pRecoSetupManager->IsInitialized()) {
70 return kFATAL;
71 }
72
73 fpDetInterface = pRecoSetupManager->GetSetup().Get<ToCbmModuleId(DetID)>();
74 switch (DetID) {
76 fpBrPoints = pMcManager->InitBranch("MvdPoint");
77 fpBrHits = dynamic_cast<TClonesArray*>(pFairManager->GetObject("MvdHit"));
78 fpBrHitMatches = dynamic_cast<TClonesArray*>(pFairManager->GetObject("MvdHitMatch"));
80 CheckBranch(fpBrPoints, "MvdPoint");
81 CheckBranch(fpBrHits, "MvdHit");
82 CheckBranch(fpBrHitMatches, "MvdHitMatch");
83 break;
85 fpBrPoints = pMcManager->InitBranch("StsPoint");
86 fpBrHits = dynamic_cast<TClonesArray*>(pFairManager->GetObject("StsHit"));
87 fpBrHitMatches = dynamic_cast<TClonesArray*>(pFairManager->GetObject("StsHitMatch"));
89 CheckBranch(fpBrPoints, "StsPoint");
90 CheckBranch(fpBrHits, "StsHit");
91 CheckBranch(fpBrHitMatches, "StsHitMatch");
92 break;
94 fpBrPoints = pMcManager->InitBranch("MuchPoint");
95 fpBrHits = dynamic_cast<TClonesArray*>(pFairManager->GetObject("MuchPixelHit"));
96 fpBrHitMatches = dynamic_cast<TClonesArray*>(pFairManager->GetObject("MuchPixelHitMatch"));
98 CheckBranch(fpBrPoints, "MuchPoint");
99 CheckBranch(fpBrHits, "MuchPixelHit");
100 CheckBranch(fpBrHitMatches, "MuchPixelHitMatch");
101 break;
103 fpBrPoints = pMcManager->InitBranch("TrdPoint");
104 fpBrHits = dynamic_cast<TClonesArray*>(pFairManager->GetObject("TrdHit"));
105 fpBrHitMatches = dynamic_cast<TClonesArray*>(pFairManager->GetObject("TrdHitMatch"));
107 CheckBranch(fpBrPoints, "TrdPoint");
108 CheckBranch(fpBrHits, "TrdHit");
109 CheckBranch(fpBrHitMatches, "TrdHitMatch");
110 break;
112 fpBrPoints = pMcManager->InitBranch("TofPoint");
113 fpBrHits = dynamic_cast<TClonesArray*>(pFairManager->GetObject("TofHit"));
114 fpBrHitMatches = dynamic_cast<TClonesArray*>(pFairManager->GetObject("TofHitMatch"));
116 CheckBranch(fpBrPoints, "TofPoint");
117 CheckBranch(fpBrHits, "TofHit");
118 CheckBranch(fpBrHitMatches, "TofHitMatch");
119 break;
121 LOG(fatal) << "CA ideal hit producer for " << kDetName[DetID]
122 << ": Wrong template argument is used: ca::EDetectorID::END";
123 break;
124 }
125
126 if constexpr (DetID == ca::EDetectorID::kTof) {
127 if (auto pGeoNodeMap = pRecoSetupManager->GetGeoNodeMap()) {
128 fTofRpcZpos = std::move(
129 pGeoNodeMap->GetVolumeProperties([](const auto& vol) { return vol.GetCenterZ(); }, ECbmModuleId::kTof));
130 }
131 else {
132 throw std::runtime_error("geo node map was not built in cbm::RecoSetupManager. Please request it with "
133 "cbm::RecoSetupManager::BuildGeoNodeMaps(), when adding the manager as a task "
134 "to the FairRun");
135 }
136 }
137
138 // ------ Checks of the settings of hit creation/adjustment mode
139 if (fsConfigName.size() != 0) {
140 this->ParseConfig();
141 // Skip event/time slice processing, if there are no requests for hit modification of the detector
142 // (fbRunTheRoutine = true, if there is at least one station with fMode > 0)
143 fbRunTheRoutine =
144 std::any_of(fvStationPars.begin(), fvStationPars.end(), [](const auto& p) { return p.fMode > 0; });
145 }
146 else {
147 fbRunTheRoutine = false;
148 }
149
150 std::stringstream msg;
151 msg << "\033[1;31mATTENTION! IDEAL HIT PRODUCER IS USED FOR " << kDetName[DetID] << "\033[0m";
152 LOG_IF(info, fbRunTheRoutine) << msg.str();
153
154 return kSUCCESS;
155}
156catch (const std::logic_error& error) {
157 LOG(error) << "CA ideal hit producer for " << kDetName[DetID] << ": initialization failed. Reason: " << error.what();
158 return kERROR;
159}
160
161// -------------------------------------------------------------------------------------------------------------------
162//
163template<ca::EDetectorID DetID>
165{
166 int nStations = fpDetInterface->GetNofTrackingStations();
167 fvStationPars.clear();
168 fvStationPars.resize(nStations);
169
170 // Read file
171 YAML::Node config;
172 try {
173 config = YAML::LoadFile(fsConfigName)["ideal_hit_producer"];
174
175 std::string detBranchName = "";
176 if constexpr (ca::EDetectorID::kMvd == DetID) {
177 detBranchName = "mvd";
178 }
179 else if constexpr (ca::EDetectorID::kSts == DetID) {
180 detBranchName = "sts";
181 }
182 else if constexpr (ca::EDetectorID::kMuch == DetID) {
183 detBranchName = "much";
184 }
185 else if constexpr (ca::EDetectorID::kTrd == DetID) {
186 detBranchName = "trd";
187 }
188 else if constexpr (ca::EDetectorID::kTof == DetID) {
189 detBranchName = "tof";
190 }
191
192 const auto& parNode = config["parameters"][detBranchName];
193 const auto& flagNode = config["flags"][detBranchName];
194
195 for (int iSt = 0; iSt < nStations; ++iSt) {
196 auto& par = fvStationPars[iSt];
197 auto& node = parNode[iSt];
198 short mode = flagNode[iSt].as<int>(0);
199 if (mode > 0) {
200 par.fPhiU = node["du"]["angle"].as<double>() * 180. / M_PI;
201 par.fDu = node["du"]["value"].as<double>();
202 switch (node["du"]["pdf"].as<char>()) {
203 case 'g': par.fPdfU = 0; break;
204 case 'u': par.fPdfU = 1; break;
205 default: par.fPdfU = -1; break;
206 }
207 par.fPhiV = node["dv"]["angle"].as<double>() * 180. / M_PI;
208 par.fDv = node["dv"]["value"].as<double>();
209 switch (node["dv"]["pdf"].as<char>()) {
210 case 'g': par.fPdfV = 0; break;
211 case 'u': par.fPdfV = 1; break;
212 default: par.fPdfV = -1; break;
213 }
214 par.fDt = node["dt"]["value"].as<double>();
215 switch (node["dt"]["pdf"].as<char>()) {
216 case 'g': par.fPdfT = 0; break;
217 case 'u': par.fPdfT = 1; break;
218 default: par.fPdfT = -1; break;
219 }
220 switch (mode) {
221 case 1:
222 par.fMode = 1;
223 par.fbSmear = false;
224 break;
225 case 2:
226 par.fMode = 1;
227 par.fbSmear = true;
228 break;
229 case 3:
230 par.fMode = 2;
231 par.fbSmear = false;
232 break;
233 case 4:
234 par.fMode = 2;
235 par.fbSmear = true;
236 break;
237 }
238 }
239 }
240 }
241 catch (const YAML::BadFile& err) {
242 std::stringstream msg;
243 msg << "YAML configuration file " << fsConfigName << " does not exist";
244 throw std::runtime_error(msg.str());
245 }
246 catch (const YAML::ParserException& err) {
247 std::stringstream msg;
248 msg << "YAML configuration file " << fsConfigName << " is formatted improperly";
249 throw std::runtime_error(msg.str());
250 }
251}
252
253// -------------------------------------------------------------------------------------------------------------------
254//
255template<ca::EDetectorID DetID>
257{
258 // ------ Check, if there are any requirement to run the routine for this detector
259 if (!fbRunTheRoutine) {
260 return;
261 }
262
263 // ------ Copy hit and match branches to the temporary array and clear the main arrays
264 if (fpBrHits) {
265 fpBrHitsTmp = static_cast<TClonesArray*>(fpBrHits->Clone("tmpHit"));
266 fpBrHits->Delete();
267 }
268 if (fpBrHitMatches) {
269 fpBrHitMatchesTmp = static_cast<TClonesArray*>(fpBrHitMatches->Clone("tmpHitMatch"));
270 fpBrHitMatches->Delete();
271 }
272
273 assert(fpBrHits);
274 assert(fpBrHitMatches);
275 assert(fpBrHitsTmp);
276 assert(fpBrHitMatchesTmp);
277
278 std::vector<double> vHitMcEventTime;
279
280 // ------ Fill main hit array
281 fHitCounter = 0;
282 int iCluster = 0;
283 for (int iH = 0; iH < fpBrHitsTmp->GetEntriesFast(); ++iH) {
284 auto* pHit = static_cast<Hit_t*>(fpBrHitsTmp->At(iH));
285 double tMC = pHit->GetTime();
286
287 // Get setting
288 int iSt = fpDetInterface->GetTrackingStationId(pHit->GetAddress());
289 if (iSt < 0) {
290 continue;
291 }
292 const auto& setting = fvStationPars[iSt];
293
294 // Skip hits from stations, where new hits are to be created from points
295 if (setting.fMode == 2) {
296 continue;
297 }
298
299 // ------ Adjust hit to matched MC point
300 if (setting.fMode == 1) {
301 // Access matched MC point
302 auto oLinkToPoint = GetMatchedPointLink(iH);
303 if (oLinkToPoint == std::nullopt) {
304 continue;
305 } // Hit had no links, so it's fake
306 auto& link = oLinkToPoint.value();
307 auto* pPoint = static_cast<Point_t*>(fpBrPoints->Get(link));
308
309 // Get point position and time
310 auto pos = TVector3{};
311 pPoint->Position(pos);
312
313 // set the hit position to the middle of the point entry and exit positions
314 //if constexpr (ca::EDetectorID::kTof != DetID) {
315 // auto posOut = TVector3{};
316 // pPoint->PositionOut(posOut);
317 // pos = 0.5 * (pos + posOut);
318 //}
319
320 double x = pos.X();
321 double y = pos.Y();
322 double z = pos.Z();
323 double t = pPoint->GetTime() + fpMCEventList->GetEventTime(link);
324 tMC = fpMCEventList->GetEventTime(link);
325
326 if (setting.fbSmear) {
327 double dx2 = pHit->GetDx() * pHit->GetDx();
328 double dxy = pHit->GetDxy();
329 double dy2 = pHit->GetDy() * pHit->GetDy();
330
331 CaUvConverter conv(setting.fPhiU, dx2, dxy, dy2);
332
333 auto [u, v] = conv.ConvertXYtoUV(x, y);
334 auto [du2, duv, dv2] = conv.ConvertCovMatrixXYtoUV(dx2, dxy, dy2);
335
336 double dt = pHit->GetTimeError();
337 this->SmearValue(u, sqrt(fabs(du2)), setting.fPdfU);
338 this->SmearValue(v, sqrt(fabs(dv2)), setting.fPdfV);
339 this->SmearValue(t, dt, setting.fPdfT);
340
341 std::tie(x, y) = conv.ConvertUVtoXY(u, v);
342 }
343 pHit->SetX(x);
344 pHit->SetY(y);
345 pHit->SetZ(z);
346 pHit->SetTime(t);
347 } // IF setting.fMode == 1
348
349 // Check out the index of cluster
350 if constexpr (ca::EDetectorID::kSts == DetID) {
351 iCluster = std::max(pHit->GetFrontClusterId(), iCluster);
352 iCluster = std::max(pHit->GetBackClusterId(), iCluster);
353 }
354
355 PushBackHit(pHit);
356 vHitMcEventTime.push_back(tMC);
357 auto* pHitMatchNew = new ((*fpBrHitMatches)[fHitCounter]) CbmMatch();
358 pHitMatchNew->AddLinks(*static_cast<CbmMatch*>(fpBrHitMatchesTmp->At(iH)));
359 ++fHitCounter;
360 }
361
362 // ------ Create hits from points
363 int nEvents = fpMCEventList->GetNofEvents();
364 for (int iE = 0; iE < nEvents; ++iE) {
365 int fileId = fpMCEventList->GetFileIdByIndex(iE);
366 int eventId = fpMCEventList->GetEventIdByIndex(iE);
367 int nPoints = fpBrPoints->Size(fileId, eventId);
368
369 // TODO: Write point selector for TOF
370 // NOTE: does not implemented for MVD
371 this->FillPointIsLegit(iE);
372 for (int iP = 0; iP < nPoints; ++iP) {
373 if (!fvbPointIsLegit[iP]) {
374 continue;
375 }
376 auto* pPoint = static_cast<Point_t*>(fpBrPoints->Get(fileId, eventId, iP));
377 int iSt = fpDetInterface->GetTrackingStationId(pPoint->GetDetectorID());
378 if (iSt < 0) {
379 continue;
380 }
381
382 const auto& setting = fvStationPars[iSt];
383 if (setting.fMode != 2) {
384 continue;
385 }
386 double du = setting.fDu;
387 double dv = setting.fDv;
388 double dt = setting.fDt;
389 double dz = 0.;
390
391 CaUvConverter conv(setting.fPhiU, setting.fPhiV);
392 auto [dx2, dxy, dy2] = conv.ConvertCovMatrixUVtoXY(du * du, 0., dv * dv);
393
394 // Get point position and time
395 auto pos = TVector3{};
396 pPoint->Position(pos);
397
398 // set the hit position to the middle of the point entry and exit positions
399 //if constexpr (ca::EDetectorID::kTof != DetID) {
400 // auto posOut = TVector3{};
401 // pPoint->PositionOut(posOut);
402 // pos = 0.5 * (pos + posOut);
403 //}
404
405 double x = pos.X();
406 double y = pos.Y();
407 double z = pos.Z();
408 double t = pPoint->GetTime() + fpMCEventList->GetEventTime(eventId, fileId);
409
410 double tMC = fpMCEventList->GetEventTime(eventId, fileId);
411
412 if (setting.fbSmear) { // TODO: Provide more realistic profiles for TRD
413 auto [u, v] = conv.ConvertXYtoUV(x, y);
414 this->SmearValue(u, du, setting.fPdfU);
415 this->SmearValue(v, dv, setting.fPdfV);
416 this->SmearValue(t, dt, setting.fPdfT);
417 std::tie(x, y) = conv.ConvertUVtoXY(u, v);
418 }
419 pos.SetX(x);
420 pos.SetY(y);
421 pos.SetZ(z);
422 auto dpos = TVector3{std::sqrt(dx2), std::sqrt(dy2), dz};
423
424 // Push back a new artificial hit
425 if constexpr (ca::EDetectorID::kMvd == DetID) {
426 new ((*fpBrHits)[fHitCounter]) CbmMvdHit();
427 }
428 else {
429 auto address = pPoint->GetDetectorID();
430 if constexpr (ca::EDetectorID::kSts == DetID) {
431 int32_t iClustF = iCluster;
432 int32_t iClustB = iCluster + 1;
433 new ((*fpBrHits)[fHitCounter]) CbmStsHit(address, pos, dpos, dxy, iClustF, iClustB, t, dt, du, dv);
434 iCluster += 2;
435 }
436 else if constexpr (ca::EDetectorID::kMuch == DetID) {
437 // TODO: What is planeID in MuCh?
438 int32_t refId = -1;
439 int32_t planeId = -1;
440 new ((*fpBrHits)[fHitCounter]) CbmMuchPixelHit(address, pos, dpos, dxy, refId, planeId, t, dt);
441 }
442 else if constexpr (ca::EDetectorID::kTrd == DetID) {
443 int32_t refId = -1;
444 auto eLoss = pPoint->GetEnergyLoss();
445 new ((*fpBrHits)[fHitCounter]) CbmTrdHit(address, pos, dpos, dxy, refId, eLoss, t, dt);
446 }
447 else if constexpr (ca::EDetectorID::kTof == DetID) {
448 int32_t refId = 0;
449 int32_t flag = 0;
450 int32_t channel = 0;
451 new ((*fpBrHits)[fHitCounter]) CbmTofHit(address, pos, dpos, refId, t, dt, flag, channel);
452 }
453 }
454
455 // Update hit match
456 auto* pHitMatchNew = new ((*fpBrHitMatches)[fHitCounter]) CbmMatch();
457 pHitMatchNew->AddLink(1., iP, eventId, fileId);
458 vHitMcEventTime.push_back(tMC);
459 ++fHitCounter;
460 assert(fpBrHitMatches->GetEntriesFast() == fHitCounter);
461 } // iP
462 }
463
464 // ------ Clear temporary hits array
465 fpBrHitsTmp->Delete();
466 delete fpBrHitsTmp;
467 fpBrHitsTmp = nullptr;
468
469 fpBrHitMatchesTmp->Delete();
470 delete fpBrHitMatchesTmp;
471 fpBrHitMatchesTmp = nullptr;
472
473 // --- Set new hit data to the reconstructed events
474 if (fpRecoEvents) {
475 for (Int_t iEvent = 0; iEvent < fpRecoEvents->GetEntriesFast(); iEvent++) {
476 CbmEvent* event = static_cast<CbmEvent*>(fpRecoEvents->At(iEvent));
477 event->ClearData(fHitDataType);
478 double tStart = event->GetStartTime();
479 double tEnd = event->GetEndTime();
480 for (UInt_t iH = 0; iH < vHitMcEventTime.size(); iH++) {
481 double t = vHitMcEventTime[iH];
482 if ((t >= tStart && t <= tEnd) || tStart < 0 || tEnd < 0) {
483 event->AddData(fHitDataType, iH);
484 }
485 }
486 }
487 }
488}
489
490// -------------------------------------------------------------------------------------------------------------------
491//
492template<ca::EDetectorID DetID>
494{
495 int iFile = fpMCEventList->GetFileIdByIndex(iEvGlob);
496 int iEvent = fpMCEventList->GetEventIdByIndex(iEvGlob);
497 int nPoints = fpBrPoints->Size(iFile, iEvent);
498
499 fvbPointIsLegit.clear();
500 fvbPointIsLegit.resize(nPoints);
501 if constexpr (ca::EDetectorID::kTof == DetID) {
502 constexpr int kNofBitsRpcAddress = 11;
503 std::map<std::pair<int, int>, int> mMatchedPointId; // map (iTr, addr) -> is matched
504 for (int iH = 0; iH < fpBrHitMatches->GetEntriesFast(); ++iH) {
505 auto* pHitMatch = dynamic_cast<CbmMatch*>(fpBrHitMatches->At(iH));
506 LOG_IF(fatal, !pHitMatch) << "CA MC Module: hit match was not found for TOF hit " << iH;
507 double bestWeight = 0;
508 for (int iLink = 0; iLink < pHitMatch->GetNofLinks(); ++iLink) {
509 const auto& link = pHitMatch->GetLink(iLink);
510 if (link.GetFile() != iFile || link.GetEntry() != iEvent) {
511 continue;
512 }
513 int iPext = link.GetIndex();
514 if (iPext < 0) {
515 continue;
516 }
517 auto* pExtPoint = static_cast<CbmTofPoint*>(fpBrPoints->Get(link));
518 int trkId = pExtPoint->GetTrackID();
519 int rpcAddr = pExtPoint->GetDetectorID() << kNofBitsRpcAddress; // FIXME:
520 auto key = std::make_pair(trkId, rpcAddr);
521 auto prev = mMatchedPointId.find(key);
522 if (prev == mMatchedPointId.end()) {
523 bestWeight = link.GetWeight();
524 mMatchedPointId[key] = iPext;
525 }
526 else { // If we find two links for the same interaction, we select the one with the largest weight
527 if (bestWeight < link.GetWeight()) {
528 mMatchedPointId[key] = iPext;
529 }
530 }
531 }
532 } // iH
533
534 // In TOF each RPC contains up to 10 points, and one have to select only
535 // one of them
536 {
537 int iPointSelected = -1;
538 int iTmcCurr = -1;
539 int rpcAddrCurr = -1;
540 bool bTrkHasHits = false;
541 double zDist = std::numeric_limits<double>::max();
542 double zCell = std::numeric_limits<double>::signaling_NaN();
543 for (int iP = 0; iP < nPoints; ++iP) {
544 auto* pExtPoint = dynamic_cast<CbmTofPoint*>(fpBrPoints->Get(iFile, iEvent, iP));
545 LOG_IF(fatal, !pExtPoint) << "NO POINT: TOF: file, event, index = " << iFile << ' ' << iEvent << ' ' << iP;
546
547 int iTmc = pExtPoint->GetTrackID();
548 int rpcAddr = pExtPoint->GetDetectorID() << kNofBitsRpcAddress;
549 // New point
550 if (rpcAddrCurr != rpcAddr || iTmcCurr != iTmc) { // The new interaction of the MC track with the TOF RPC
551 if (iPointSelected != -1) {
552 fvbPointIsLegit[iPointSelected] = true;
553 }
554 iTmcCurr = iTmc;
555 rpcAddrCurr = rpcAddr;
556 auto key = std::make_pair(iTmc, rpcAddr);
557 auto found = mMatchedPointId.find(key);
558 bTrkHasHits = found != mMatchedPointId.end();
559 if (bTrkHasHits) {
560 iPointSelected = found->second;
561 }
562 else {
563 // First iteration
564 zCell = fTofRpcZpos.at(cbm::GeoNodeMap::ApplyMask(pExtPoint->GetDetectorID(), ECbmModuleId::kTof));
565 zDist = std::fabs(zCell - pExtPoint->GetZ());
566 iPointSelected = iP;
567 }
568 }
569 else {
570 if (!bTrkHasHits) {
571 auto newDist = std::fabs(pExtPoint->GetZ() - zCell);
572 if (newDist < zDist) {
573 zDist = newDist;
574 iPointSelected = iP;
575 }
576 }
577 }
578 }
579 // Add the last point
580 if (iPointSelected != -1) {
581 fvbPointIsLegit[iPointSelected] = true;
582 }
583 }
584 }
585 else {
586 // In MVD, STS, MuCh and TRD every MC point is legit
587 std::fill(fvbPointIsLegit.begin(), fvbPointIsLegit.end(), true);
588 }
589}
590
@ kTof
Time-of-flight Detector.
Definition CbmDefs.h:52
TVector3 dpos
friend fvec sqrt(const fvec &a)
fscal v[fmask::Size]
Definition KfSimdPseudo.h:4
int Int_t
Class characterising one event by a collection of links (indices) to data objects,...
Definition CbmEvent.h:34
void ClearData(ECbmDataType type)
Definition CbmEvent.cxx:37
Task class creating and managing CbmMCDataArray objects.
Container class for MC events with number, file and start time.
void AddLinks(const CbmMatch &match)
Definition CbmMatch.cxx:39
void AddLink(const CbmLink &newLink)
Definition CbmMatch.cxx:47
data class for a reconstructed 3-d hit in the STS
Definition CbmStsHit.h:35
Bookkeeping of time-slice content.
Geometric intersection of a MC track with a TOFb detector.
Definition CbmTofPoint.h:44
data class for a reconstructed Energy-4D measurement in the TRD
Definition CbmTrdHit.h:40
static uint32_t ApplyMask(uint32_t address, ECbmModuleId moduleId)
Applies mask.
void BuildGeoNodeMaps()
Requests geo-node map building.
static RecoSetupManager * Instance()
Instance access.
A class to convert XY coordinates to UV coordinates.
std::tuple< double, double, double > ConvertCovMatrixUVtoXY(double du2, double duv, double dv2) const
Conversion function (du2, duv, dv2) -> (dx2, dxy, dy2)
std::pair< double, double > ConvertUVtoXY(double u, double v) const
Conversion function (x,y) -> (u,v)
std::tuple< double, double, double > ConvertCovMatrixXYtoUV(double dx2, double dxy, double dy2) const
Conversion function (dx2, dxy, dy2) -> (du2, duv, dv2)
std::pair< double, double > ConvertXYtoUV(double x, double y) const
Conversion function (x,y) -> (u,v)
TClonesArray * fpRecoEvents
Array of reconstructed events.
std::vector< unsigned char > fvbPointIsLegit
Map of used point index.
void Exec(Option_t *option)
Execution of the task.
TClonesArray * fpBrHitsTmp
Temporary array of hits.
PointTypes_t::at< DetID > Point_t
CbmMCEventList * fpMCEventList
MC event list.
TClonesArray * fpBrHitMatchesTmp
Temporary array of hit matches.
std::vector< HitParameters > fvStationPars
Parameters, stored for each station.
const cbm::algo::RecoSetupUnit_t< ToCbmModuleId(DetID)> * fpDetInterface
CbmMCDataObject * fpMCEventHeader
MC event header.
int fHitCounter
Hit counter in a new branch.
CbmTimeSlice * fpTimeSlice
Current time slice.
void PushBackHit(const Hit_t *pHit)
Pushes back a hit into the hits branch.
std::optional< CbmLink > GetMatchedPointLink(int iH)
Gets pointer to matched MC point.
bool fbRunTheRoutine
Management flag, which does run the routine if there was a request.
void FillPointIsLegit(int iEvGlob)
Fills map of legit points (inside MC event!)
std::string fsConfigName
Name of configuration file.
InitStatus Init()
Initialization of the task.
CbmMCDataArray * fpBrMCTracks
Branch: array of MC tracks.
void ParseConfig()
Parses the YAML configuration file.
void SmearValue(double &value, double sigma, int opt)
Smears the value by a given standard deviation with a selected function.
constexpr DetIdArr_t< const char * > kDetName
Names of detector subsystems.