diff --git a/core/data/CbmEvent.cxx b/core/data/CbmEvent.cxx index 37857fa7e8df52adda2c1ed44a9fe8d7ac1c1b83..5339587459d0af4342c274715e812578d88ed0f3 100644 --- a/core/data/CbmEvent.cxx +++ b/core/data/CbmEvent.cxx @@ -33,6 +33,10 @@ CbmEvent::CbmEvent(const CbmEvent& rhs) void CbmEvent::AddData(ECbmDataType type, uint32_t index) { fIndexMap[type].push_back(index); } // ------------------------------------------------------------------------- +// ----- Clear a specific data branch --------------------------------------------- +void CbmEvent::ClearData(ECbmDataType type) { fIndexMap[type].clear(); } +// ------------------------------------------------------------------------- + // ----- Get a data index ---------------------------------------------- uint32_t CbmEvent::GetIndex(ECbmDataType type, uint32_t iData) diff --git a/core/data/CbmEvent.h b/core/data/CbmEvent.h index 01f416f26bcfa089c00cd92ae23d3c4d78f7f25f..e51c70fcc19a8bda518f0db7342997abaf2fe7c7 100644 --- a/core/data/CbmEvent.h +++ b/core/data/CbmEvent.h @@ -66,6 +66,10 @@ public: /** Overload TObject Clear to clear the map! **/ void Clear(Option_t*) { fIndexMap.clear(); } + /** Clear a specific data branch in the index map + ** @param DataType Type of data (for values see CbmDetectorList.h) + */ + void ClearData(ECbmDataType type); /** Add a data object to the index map ** @param DataType Type of data (for values see CbmDetectorList.h) diff --git a/reco/L1/utils/CbmCaIdealHitProducerDet.h b/reco/L1/utils/CbmCaIdealHitProducerDet.h index e7c5475ae3f8e6598a4f1af87eb13857c4fdc498..6393f354c3b6dc416aad60ff3b59e27ae4b7b9cc 100644 --- a/reco/L1/utils/CbmCaIdealHitProducerDet.h +++ b/reco/L1/utils/CbmCaIdealHitProducerDet.h @@ -10,6 +10,7 @@ #ifndef CbmCaIdealHitProducerDet_h #define CbmCaIdealHitProducerDet_h 1 +#include "CbmEvent.h" #include "CbmL1DetectorID.h" #include "CbmMCDataArray.h" #include "CbmMCDataManager.h" @@ -158,6 +159,7 @@ namespace cbm::ca // ----- Input branches: CbmTimeSlice* fpTimeSlice = nullptr; ///< Current time slice + TClonesArray* fpRecoEvents = nullptr; ///< Array of reconstructed events CbmMCEventList* fpMCEventList = nullptr; ///< MC event list CbmMCDataObject* fpMCEventHeader = nullptr; ///< MC event header CbmMCDataArray* fpBrPoints = nullptr; ///< Branch: array of MC points @@ -168,6 +170,8 @@ namespace cbm::ca TClonesArray* fpBrHitsTmp = nullptr; ///< Temporary array of hits TClonesArray* fpBrHitMatchesTmp = nullptr; ///< Temporary array of hit matches + ECbmDataType fHitDataType = ECbmDataType::kUnknown; ///< Hit data type + std::string fsConfigName = ""; ///< Name of configuration file std::vector<HitParameters> fvStationPars; ///< Parameters, stored for each station @@ -219,6 +223,7 @@ namespace cbm::ca CheckBranch(pMcManager, "MCDataManager"); fpTimeSlice = dynamic_cast<CbmTimeSlice*>(pFairManager->GetObject("TimeSlice.")); + fpRecoEvents = dynamic_cast<TClonesArray*>(pFairManager->GetObject("CbmEvent")); fpMCEventList = dynamic_cast<CbmMCEventList*>(pFairManager->GetObject("MCEventList.")); CheckBranch(fpTimeSlice, "TimeSlice."); CheckBranch(fpMCEventList, "MCEventList."); @@ -230,6 +235,7 @@ namespace cbm::ca fpBrPoints = pMcManager->InitBranch("MvdPoint"); fpBrHits = dynamic_cast<TClonesArray*>(pFairManager->GetObject("MvdHit")); fpBrHitMatches = dynamic_cast<TClonesArray*>(pFairManager->GetObject("MvdHitMatch")); + fHitDataType = ECbmDataType::kMvdHit; CheckBranch(fpBrPoints, "MvdPoint"); CheckBranch(fpBrHits, "MvdHit"); CheckBranch(fpBrHitMatches, "MvdHitMatch"); @@ -239,6 +245,7 @@ namespace cbm::ca fpBrPoints = pMcManager->InitBranch("StsPoint"); fpBrHits = dynamic_cast<TClonesArray*>(pFairManager->GetObject("StsHit")); fpBrHitMatches = dynamic_cast<TClonesArray*>(pFairManager->GetObject("StsHitMatch")); + fHitDataType = ECbmDataType::kStsHit; CheckBranch(fpBrPoints, "StsPoint"); CheckBranch(fpBrHits, "StsHit"); CheckBranch(fpBrHitMatches, "StsHitMatch"); @@ -248,6 +255,7 @@ namespace cbm::ca fpBrPoints = pMcManager->InitBranch("MuchPoint"); fpBrHits = dynamic_cast<TClonesArray*>(pFairManager->GetObject("MuchPixelHit")); fpBrHitMatches = dynamic_cast<TClonesArray*>(pFairManager->GetObject("MuchPixelHitMatch")); + fHitDataType = ECbmDataType::kMuchPixelHit; CheckBranch(fpBrPoints, "MuchPoint"); CheckBranch(fpBrHits, "MuchPixelHit"); CheckBranch(fpBrHitMatches, "MuchPixelHitMatch"); @@ -257,6 +265,7 @@ namespace cbm::ca fpBrPoints = pMcManager->InitBranch("TrdPoint"); fpBrHits = dynamic_cast<TClonesArray*>(pFairManager->GetObject("TrdHit")); fpBrHitMatches = dynamic_cast<TClonesArray*>(pFairManager->GetObject("TrdHitMatch")); + fHitDataType = ECbmDataType::kTrdHit; CheckBranch(fpBrPoints, "TrdPoint"); CheckBranch(fpBrHits, "TrdHit"); CheckBranch(fpBrHitMatches, "TrdHitMatch"); @@ -266,6 +275,7 @@ namespace cbm::ca fpBrPoints = pMcManager->InitBranch("TofPoint"); fpBrHits = dynamic_cast<TClonesArray*>(pFairManager->GetObject("TofHit")); fpBrHitMatches = dynamic_cast<TClonesArray*>(pFairManager->GetObject("TofHitMatch")); + fHitDataType = ECbmDataType::kTofHit; CheckBranch(fpBrPoints, "TofPoint"); CheckBranch(fpBrHits, "TofHit"); CheckBranch(fpBrHitMatches, "TofHitMatch"); @@ -499,11 +509,14 @@ namespace cbm::ca assert(fpBrHitsTmp); assert(fpBrHitMatchesTmp); - // ------ Fill main hits array + std::vector<double> vHitMcEventTime; + + // ------ Fill main hit array fHitCounter = 0; int iCluster = 0; ///< Current cluster index for (int iH = 0; iH < fpBrHitsTmp->GetEntriesFast(); ++iH) { auto* pHit = static_cast<Hit_t*>(fpBrHitsTmp->At(iH)); + double tMC = pHit->GetTime(); // Get setting int iSt = fpDetInterface->GetTrackingStationIndex(pHit); @@ -541,6 +554,7 @@ namespace cbm::ca double y = pos.Y(); double z = pos.Z(); double t = pPoint->GetTime() + fpMCEventList->GetEventTime(link); + tMC = fpMCEventList->GetEventTime(link); if (setting.fbSmear) { double dx2 = pHit->GetDx() * pHit->GetDx(); @@ -572,6 +586,7 @@ namespace cbm::ca } PushBackHit(pHit); + vHitMcEventTime.push_back(tMC); auto* pHitMatchNew = new ((*fpBrHitMatches)[fHitCounter]) CbmMatch(); pHitMatchNew->AddLinks(*static_cast<CbmMatch*>(fpBrHitMatchesTmp->At(iH))); ++fHitCounter; @@ -617,6 +632,8 @@ namespace cbm::ca double z = pos.Z(); double t = pPoint->GetTime() + fpMCEventList->GetEventTime(eventId, fileId); + double tMC = fpMCEventList->GetEventTime(eventId, fileId); + if (setting.fbSmear) { // TODO: Provide more realistic profiles for TRD auto [u, v] = conv.ConvertXYtoUV(x, y); this->SmearValue(u, du, setting.fPdfU); @@ -661,7 +678,9 @@ namespace cbm::ca // Update hit match auto* pHitMatchNew = new ((*fpBrHitMatches)[fHitCounter]) CbmMatch(); pHitMatchNew->AddLink(1., iP, eventId, fileId); + vHitMcEventTime.push_back(tMC); ++fHitCounter; + assert(fpBrHitMatches->GetEntriesFast() == fHitCounter); } // iP } @@ -673,6 +692,20 @@ namespace cbm::ca fpBrHitMatchesTmp->Delete(); delete fpBrHitMatchesTmp; fpBrHitMatchesTmp = nullptr; + + // --- Set new hit data to the reconstructed events + if (fpRecoEvents) { + for (Int_t iEvent = 0; iEvent < fpRecoEvents->GetEntriesFast(); iEvent++) { + CbmEvent* event = static_cast<CbmEvent*>(fpRecoEvents->At(iEvent)); + event->ClearData(fHitDataType); + double tStart = event->GetStartTime(); + double tEnd = event->GetEndTime(); + for (UInt_t iH = 0; iH < vHitMcEventTime.size(); iH++) { + double t = vHitMcEventTime[iH]; + if ((t >= tStart && t <= tEnd) || tStart < 0 || tEnd < 0) { event->AddData(fHitDataType, iH); } + } + } + } } } // namespace cbm::ca