BmnRoot
Loading...
Searching...
No Matches
BmnFHCal.cxx
Go to the documentation of this file.
1#include "BmnFHCal.h"
2
3#include "CbmStack.h"
4#include "FairLinkManager.h"
5#include "FairRootManager.h"
6#include "FairRun.h"
7#include "FairRuntimeDb.h"
8#include "FairVolume.h"
9#include "TGeoManager.h"
10#include "TLorentzVector.h"
11#include "TParticle.h"
12#include "TVirtualMC.h"
13
14#include <cassert>
15#include <iostream>
16
18 : FairDetector("BmnFHCal", true)
19 , fCollection(new TClonesArray("BmnFHCalPoint"))
20 , fGeoHandler(std::make_unique<BmnFHCalGeo>())
21 , fPoint(nullptr)
22{
23 LOG(debug3) << "BmnFHCal default constructor called";
24}
25
26BmnFHCal::BmnFHCal(const char* name, Bool_t active)
27 : FairDetector(name, active)
28 , fCollection(new TClonesArray("BmnFHCalPoint"))
29 , fGeoHandler(std::make_unique<BmnFHCalGeo>())
30 , fPoint(nullptr)
31{
32 LOG(debug3) << "BmnFHCal constructor called with name: " << name << " and active: " << active;
33}
34
36{
37 LOG(debug3) << "BmnFHCal destructor called";
38 if (fCollection) {
39 fCollection->Delete();
40 delete fCollection;
41 fCollection = nullptr;
42 }
43
44 fMultiLinkMap.clear();
45 fGeoHandler.reset();
46 fPoint.reset();
47}
48
50{
51 LOG(debug3) << "BmnFHCal Initialize method called";
52 FairDetector::Initialize();
53}
54
56{
57 LOG(debug3) << "BmnFHCal BeginEvent method called";
58 fMultiLinkMap.clear();
59}
60
62{
63 LOG(debug4) << "BmnFHCal GetHit method called";
64 Int_t nEntries = fCollection->GetEntriesFast();
65 for (Int_t i = nEntries - 1; i >= 0; i--) {
66 auto* hit = static_cast<BmnFHCalPoint*>(fCollection->At(i));
67 if (hit->GetAddress() == address) {
68 return hit;
69 }
70 }
71 return nullptr;
72}
73
74static double CalculateEnergyLoss(double step, double edep)
75{
76 constexpr double BirkConst = 12.6; // Static constant.
77 return (step > 0) ? edep / (1.0 + (BirkConst / step) * edep) : edep;
78}
79
80Bool_t BmnFHCal::ProcessHits(FairVolume* vol)
81{
82 LOG(debug4) << "BmnFHCal ProcessHits method called for event " << gMC->CurrentEvent();
83
84 // Skip non-sensitive volumes or invalid energies
85 if (!CheckIfSensitive(gMC->CurrentVolName()) || gMC->Edep() < 0) {
86 return kFALSE;
87 }
88
89 LOG(debug4) << "BmnFHCal::ProcessHits. trackID = " << gMC->GetStack()->GetCurrentTrackNumber() << " "
90 << BmnFHCalAddress::GetInfoString(fGeoHandler->GetAddressFromPath(gMC->CurrentVolPath()))
91 << ". Edep = " << gMC->Edep() << " Time = " << gMC->TrackTime() * 1.0e09;
92
93 if (gMC->IsTrackEntering() || !fPoint) {
94 fPoint.reset(new BmnFHCalPoint());
95 fPoint->SetEventID(gMC->CurrentEvent());
96 fPoint->SetAddress(fGeoHandler->GetAddressFromPath(gMC->CurrentVolPath()));
97 TLorentzVector tPos, tMom;
98 gMC->TrackPosition(tPos);
99 gMC->TrackMomentum(tMom);
100 fPoint->SetPosition(tPos.Vect());
101 fPoint->SetMomentum(tMom.Vect());
102 fPoint->SetEnergyLoss(CalculateEnergyLoss(gMC->TrackStep(), gMC->Edep()));
103 fPoint->SetTime(gMC->TrackTime() * 1.0e09);
104 fPoint->SetLength(gMC->TrackLength());
105 }
106
107 assert(fPoint);
108 assert(fPoint->GetAddress() == fGeoHandler->GetAddressFromPath(gMC->CurrentVolPath()));
109 double dEdep = CalculateEnergyLoss(gMC->TrackStep(), gMC->Edep());
110 double weighted_time = fPoint->GetEnergyLoss() * fPoint->GetTime() + dEdep * gMC->TrackTime() * 1.0e09;
111 fPoint->SetEnergyLoss(fPoint->GetEnergyLoss() + dEdep);
112 if (fPoint->GetEnergyLoss() > 0.)
113 fPoint->SetTime(weighted_time / fPoint->GetEnergyLoss());
114
115 if (gMC->IsTrackExiting() || gMC->IsTrackStop() || gMC->IsTrackDisappeared()) {
116 if (fPoint->GetEnergyLoss() > std::numeric_limits<double>::epsilon()) {
117 int surface_track = GetSurfaceMCTrack(gMC->GetStack()->GetCurrentTrackNumber());
118 if (surface_track >= 0) {
119 fPoint->SetTrackID(surface_track);
120 fMultiLinkMap[fPoint->GetAddress()].AddLink(
121 FairLink("MCTrack", surface_track, fPoint->GetEnergyLoss()));
122 ((CbmStack*)gMC->GetStack())->AddPoint(kFHCAL, surface_track);
123 }
124 auto* existing = GetHit(fPoint->GetAddress());
125 if (existing)
126 UpdateHit(*existing, *fPoint);
127 else
128 AddHit(fPoint.get());
129 }
130 fPoint.reset();
131 }
132
133 return kTRUE;
134}
135
137{
138 LOG(debug3) << "BmnFHCal EndOfEvent method called";
139 Print("");
140 Reset();
141}
142
144{
145 LOG(debug3) << "BmnFHCal Register method called";
146 FairRootManager::Instance()->Register("FHCalPoint", "FHCal", fCollection, kTRUE);
147}
148
149TClonesArray* BmnFHCal::GetCollection(Int_t iColl) const
150{
151 LOG(debug3) << "BmnFHCal GetCollection method called with iColl: " << iColl;
152 for (Int_t i = 0; i < fCollection->GetEntriesFast(); i++) {
153 auto* hit = static_cast<BmnFHCalPoint*>(fCollection->At(i));
154 hit->SetEntryNr(FairLink("FHCalPoint", i, hit->GetEnergyLoss()));
155 uint32_t address = hit->GetAddress();
156 if (fMultiLinkMap.find(address) != fMultiLinkMap.end()) {
157 hit->SetLinks(fMultiLinkMap.at(address)); // Overwrite
158 LOG(debug3) << "BmnFHCal::GetCollection Links: " << hit->GetNLinks();
159 if (FairLogger::GetLogger()->IsLogNeeded(fair::Severity::debug3) && hit->GetNLinks()) {
160 hit->PrintLinkInfo();
161 printf("\n");
162 }
163 }
164 }
165 return (iColl == 0) ? fCollection : nullptr;
166}
167
168void BmnFHCal::Print(Option_t*) const
169{
170 LOG(debug3) << "BmnFHCal Print method called";
171 Int_t nHits = fCollection->GetEntriesFast();
172 std::cout << "-I- BmnFHCal: " << nHits << " points registered in event " << gMC->CurrentEvent() << std::endl;
173
174 if (FairLogger::GetLogger()->IsLogNeeded(fair::Severity::debug4)) {
175 for (Int_t i = 0; i < nHits; ++i) {
176 auto hit = static_cast<BmnFHCalPoint*>((*fCollection)[i]);
177 hit->Print();
178 }
179 }
180}
181
183{
184 LOG(debug3) << "BmnFHCal Reset method called";
185 fCollection->Delete();
186}
187
189{
190 TString fileName = GetGeometryFileName();
191 if (fileName.EndsWith(".root")) {
192 LOG(info) << "Constructing FHCal geometry from ROOT file %s" << fileName.Data();
193 ConstructRootGeometry();
194 }
195
196 FairGeoLoader* geoLoad = FairGeoLoader::Instance();
197 FairGeoInterface* geoFace = geoLoad->getGeoInterface();
198 fGeoHandler->setGeomFile(GetGeometryFileName());
199 geoFace->addGeoModule(fGeoHandler.get());
200
201 Bool_t rc = geoFace->readSet(fGeoHandler.get());
202 if (rc)
203 fGeoHandler->create(geoLoad->getGeoBuilder());
204 TList* volList = fGeoHandler->getListOfVolumes();
205
206 // store geo parameter
207 FairRun* fRun = FairRun::Instance();
208 FairRuntimeDb* rtdb = FairRun::Instance()->GetRuntimeDb();
209 BmnFHCalGeoPar* par = (BmnFHCalGeoPar*)(rtdb->getContainer("BmnFHCalGeoPar"));
210 TObjArray* fSensNodes = par->GetGeoSensitiveNodes();
211 TObjArray* fPassNodes = par->GetGeoPassiveNodes();
212
213 TListIter iter(volList);
214 FairGeoNode* node = NULL;
215 FairGeoVolume* aVol = NULL;
216 while ((node = (FairGeoNode*)iter.Next())) {
217 aVol = dynamic_cast<FairGeoVolume*>(node);
218 if (node->isSensitive()) {
219 fSensNodes->AddLast(aVol);
220 } else {
221 fPassNodes->AddLast(aVol);
222 }
223 }
224 par->setChanged();
225 par->setInputVersion(fRun->GetRunId(), 1);
226
227 ProcessNodes(volList);
228}
229
230Bool_t BmnFHCal::CheckIfSensitive(std::string name)
231{
232 LOG(debug4) << "BmnFHCal CheckIfSensitive method called with name: " << name;
233 TString tsname = name.c_str();
234 return tsname.Contains(BmnFHCalGeoPar::SensitiveVolume_name);
235}
236
238{
239 TClonesArray& clref = *fCollection;
240 Int_t size = clref.GetEntriesFast();
241 return new (clref[size]) BmnFHCalPoint(*point);
242}
243
244void BmnFHCal::UpdateHit(BmnFHCalPoint& existing, const BmnFHCalPoint& update)
245{
246 LOG(debug4) << "UpdateHit: Updating hit at address " << existing.GetAddress()
247 << " with current E_loss = " << existing.GetEnergyLoss()
248 << " update E_loss += " << update.GetEnergyLoss();
249
250 const double eloss_existing = existing.GetEnergyLoss();
251 const double eloss_update = update.GetEnergyLoss();
252 const double eloss_total = eloss_existing + eloss_update;
253
254 if (eloss_total == 0)
255 return;
256
257 const double weight_existing = eloss_existing / eloss_total;
258 const double weight_update = eloss_update / eloss_total;
259
260 TVector3 pos_existing, mom_existing;
261 existing.Position(pos_existing);
262 existing.Momentum(mom_existing);
263
264 TVector3 pos_update, mom_update;
265 update.Position(pos_update);
266 update.Momentum(mom_update);
267
268 existing.SetPosition(weight_existing * pos_existing + weight_update * pos_update);
269 existing.SetMomentum(weight_existing * mom_existing + weight_update * mom_update);
270 existing.SetTime(weight_existing * existing.GetTime() + weight_update * update.GetTime());
271 existing.SetLength(weight_existing * existing.GetLength() + weight_update * update.GetLength());
272 if (weight_update > weight_existing)
273 existing.SetTrackID(update.GetTrackID());
274 existing.SetEnergyLoss(eloss_total);
275}
276
277int BmnFHCal::GetSurfaceMCTrack(int start_track_id)
278{
279 int current_id = start_track_id;
280 TParticle* current = ((CbmStack*)gMC->GetStack())->GetParticle(current_id);
281 while (current) {
282 if (current->IsPrimary()) {
283 LOG(debug4) << Form("BmnFHCal::GetSurfaceMCTrack. Current particle %d is primary", current_id);
284 return current_id;
285 }
286 int mother_id = current->GetFirstMother();
287 TParticle* mother = ((CbmStack*)gMC->GetStack())->GetParticle(mother_id);
288 if (!mother) {
289 LOG(error) << Form("BmnFHCal::GetSurfaceMCTrack. Current track with id %d has mother but it is "
290 "not stored in MC stack",
291 current_id);
292 return -1;
293 }
294 if (mother_id == current_id) {
295 LOG(error) << "BmnFHCal::GetSurfaceMCTrack. Infinite loop detected for track ID: " << current_id;
296 return -1;
297 }
298
299 TVector3 current_vtx(current->Vx(), current->Vy(), current->Vz());
300 TVector3 mother_vtx(mother->Vx(), mother->Vy(), mother->Vz());
301 if (fGeoHandler->IsPointInside(current_vtx) && !fGeoHandler->IsPointInside(mother_vtx))
302 { // daughter vtx is inside detector, mother vtx is not
303 LOG(debug4) << Form("BmnFHCal::GetSurfaceMCTrack. Current track with id %d has vtx (%.2f %.2f "
304 "%.2f) inside detector. Its mother %d vtx (%.2f %.2f %.2f) is outside",
305 current_id, current_vtx.X(), current_vtx.Y(), current_vtx.Z(), mother_id,
306 mother_vtx.X(), mother_vtx.Y(), mother_vtx.Z());
307 return mother_id;
308 }
309 current = mother;
310 current_id = mother_id;
311 }
312 return -1;
313}
int i
Definition P4_F32vec4.h:22
@ kFHCAL
static std::string GetInfoString(uint32_t address)
Return a formatted string with all address components.
TObjArray * GetGeoSensitiveNodes()
static const TString SensitiveVolume_name
TObjArray * GetGeoPassiveNodes()
virtual void Print(const Option_t *opt="") const override
uint32_t GetAddress() const
int GetSurfaceMCTrack(int start_track_id)
Definition BmnFHCal.cxx:277
BmnFHCalPoint * AddHit(BmnFHCalPoint *point)
Definition BmnFHCal.cxx:237
BmnFHCalPoint * GetHit(uint32_t address)
Definition BmnFHCal.cxx:61
virtual ~BmnFHCal()
Definition BmnFHCal.cxx:35
virtual void Reset() override
Definition BmnFHCal.cxx:182
virtual void Print(Option_t *) const override
Definition BmnFHCal.cxx:168
virtual void ConstructGeometry() override
Definition BmnFHCal.cxx:188
void UpdateHit(BmnFHCalPoint &existing, const BmnFHCalPoint &update)
Definition BmnFHCal.cxx:244
virtual void EndOfEvent() override
Definition BmnFHCal.cxx:136
virtual void Initialize() override
Definition BmnFHCal.cxx:49
virtual void Register() override
Definition BmnFHCal.cxx:143
virtual void BeginEvent() override
Definition BmnFHCal.cxx:55
virtual TClonesArray * GetCollection(Int_t iColl) const override
Definition BmnFHCal.cxx:149
virtual Bool_t ProcessHits(FairVolume *vol=nullptr) override
Definition BmnFHCal.cxx:80
virtual Bool_t CheckIfSensitive(std::string name) override
Definition BmnFHCal.cxx:230
STL namespace.