ATLAS Offline Software
Loading...
Searching...
No Matches
FastCaloSimCaloTransportation.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
3*/
4
5/* Athena includes */
7#include "GaudiKernel/ToolHandle.h"
8#include "GaudiKernel/IPartPropSvc.h"
9#include "GaudiKernel/ThreadLocalContext.h"
10
12/* Header include */
14
15/* ISF includes */
19
20/* Tracking includes */
22
23/* Geometry primitives */
25
26/* Particle data */
27#include "HepPDT/ParticleDataTable.hh"
28
29/* Transport steps will be return as G4FieldTracks*/
30#include "G4FieldTrack.hh"
32
33
34FastCaloSimCaloTransportation::FastCaloSimCaloTransportation(const std::string& t, const std::string& n, const IInterface* p)
35 : base_class(t,n,p)
36{
37
38}
39
41{
42 ATH_MSG_INFO( "Initializing FastCaloSimCaloTransportation" );
43
44 // Get TimedExtrapolator
45 if(!m_extrapolator.empty()){
46 ATH_CHECK(m_extrapolator.retrieve());
47 ATH_MSG_INFO("Extrapolator retrieved "<< m_extrapolator);
48 }
49
50 return StatusCode::SUCCESS;
51
52}
53
54
56 ATH_MSG_INFO( "Finalizing FastCaloSimCaloTransportation" );
57 return StatusCode::SUCCESS;
58}
59
60
61std::vector<G4FieldTrack> FastCaloSimCaloTransportation::transport(const TFCSTruthState* truth, bool forceNeutral) const{
62 if (auto* eventInfo = AtlasG4EventUserInfo::GetEventUserInfo()) {
63 Gaudi::Hive::setCurrentContext(eventInfo->GetEventContext());
64 }
65
66 // Start calo extrapolation
67 ATH_MSG_DEBUG ("[ fastCaloSim transport ] processing particle "<<truth->pdgid() );
68
69 auto hitVector = std::make_unique<std::vector<Trk::HitInfo>>();
70
71 int pdgId = truth->pdgid();
72 double charge = HepPDT::ParticleID(pdgId).charge();
73 if (forceNeutral) charge = 0.;
74
75 // particle Hypothesis for the extrapolation
76 Trk::ParticleHypothesis pHypothesis = m_pdgToParticleHypothesis.convert(pdgId, charge);
77
78 ATH_MSG_DEBUG ("particle hypothesis "<< pHypothesis);
79
80 // geantinos not handled by PdgToParticleHypothesis - fix there
81 if (pdgId == 999) pHypothesis = Trk::geantino;
82
83 Amg::Vector3D pos = Amg::Vector3D(truth->vertex().X(), truth->vertex().Y(), truth->vertex().Z());
84
85 Amg::Vector3D mom(truth->X(), truth->Y(), truth->Z());
86
87 ATH_MSG_DEBUG( "[ fastCaloSim transport ] x from position eta="<<pos.eta()<<" phi="<<pos.phi()<<" d="<<pos.mag()<<" pT="<<mom.perp() );
88
89 // input parameters : curvilinear parameters
90 Trk::CurvilinearParameters inputPar(pos,mom,charge);
91
92 double freepath = -1.;
93 double tDec = freepath > 0. ? freepath : -1.;
94 int decayProc = 0;
95
96 Trk::TimeLimit timeLim(tDec, 0., decayProc);
97 Trk::PathLimit pathLim(-1., 0);
99
100 // first extrapolation to reach the ID boundary
101 ATH_MSG_DEBUG( "[ fastCaloSim transport ] before calo entrance ");
102
103 // get CaloEntrance if not done already
104 if(!m_caloEntrance.get()){
105 m_caloEntrance.set(m_extrapolator->trackingGeometry()->trackingVolume(m_caloEntranceName));
106 if(!m_caloEntrance.get())
107 ATH_MSG_WARNING("CaloEntrance not found");
108 else
109 ATH_MSG_DEBUG("CaloEntrance found");
110 }
111
112 ATH_MSG_DEBUG( "[ fastCaloSim transport ] after calo entrance ");
113
114 std::unique_ptr<const Trk::TrackParameters> caloEntry = nullptr;
115
116 if(m_caloEntrance.get() && m_caloEntrance.get()->inside(pos, 0.001) && !m_extrapolator->trackingGeometry()->atVolumeBoundary(pos,m_caloEntrance.get(), 0.001)){
117 std::vector<Trk::HitInfo>* dummyHitVector = nullptr;
118 if (charge == 0){
119 caloEntry =
120 m_extrapolator->transportNeutralsWithPathLimit(inputPar,
121 pathLim,
122 timeLim,
124 pHypothesis,
125 dummyHitVector,
126 nextGeoID,
127 m_caloEntrance.get());
128 }else{
129 caloEntry = m_extrapolator->extrapolateWithPathLimit(inputPar,
130 pathLim,
131 timeLim,
133 pHypothesis,
134 dummyHitVector,
135 nextGeoID,
136 m_caloEntrance.get());
137 }
138 } else
139 caloEntry = inputPar.uniqueClone();
140
141 ATH_MSG_DEBUG( "[ fastCaloSim transport ] after calo caloEntry ");
142
143 if(caloEntry){
144 std::unique_ptr<const Trk::TrackParameters> eParameters = nullptr;
145
146 // save Calo entry hit (fallback info)
147 hitVector->push_back(Trk::HitInfo(caloEntry->uniqueClone(), timeLim.time, nextGeoID, 0.));
148
150 "[ fastCaloSim transport ] starting Calo transport from position eta="
151 << caloEntry->position().eta() << " phi=" << caloEntry->position().phi()
152 << " d=" << caloEntry->position().mag());
153
154 std::vector<Trk::HitInfo>* rawHitVector = hitVector.get();
155 if (charge == 0){
156 eParameters =
157 m_extrapolator->transportNeutralsWithPathLimit(*caloEntry,
158 pathLim,
159 timeLim,
161 pHypothesis,
162 rawHitVector,
163 nextGeoID);
164 }else{
165 eParameters =
166 m_extrapolator->extrapolateWithPathLimit(*caloEntry,
167 pathLim,
168 timeLim,
170 pHypothesis,
171 rawHitVector,
172 nextGeoID);
173 }
174 // save Calo exit hit (fallback info)
175 if(eParameters) hitVector->push_back(Trk::HitInfo(std::move(eParameters), timeLim.time, nextGeoID, 0.));
176 } //if caloEntry
177
178 //used to identify ID-Calo boundary in tracking tools
179 int IDCaloBoundary = 3000;
180
181 if(msgLvl(MSG::DEBUG)){
182 std::vector<Trk::HitInfo>::iterator it = hitVector->begin();
183 while (it < hitVector->end()){
184 int sample=(*it).detID;
185 Amg::Vector3D hitPos = (*it).trackParms->position();
186 ATH_MSG_DEBUG(" HIT: layer="<<sample<<" sample="<<sample-IDCaloBoundary<<" eta="<<hitPos.eta()<<" phi="<<hitPos.phi()<<" d="<<hitPos.mag());
187 ++it;
188 }
189 }
190
191 std::vector<Trk::HitInfo>::iterator it2 = hitVector->begin();
192 while(it2 < hitVector->end()){
193 int sample=(*it2).detID;
194 Amg::Vector3D hitPos = (*it2).trackParms->position();
195 ATH_MSG_DEBUG(" HIT: layer="<<sample<<" sample="<<sample-IDCaloBoundary<<" eta="<<hitPos.eta()<<" phi="<<hitPos.phi()<<" r="<<hitPos.perp()<<" z="<<hitPos[Amg::z]);
196 ++it2;
197 }
198
199 // Extrapolation may fail for very low pT charged particles. Enforce charge 0 to prevent this
200 if (!forceNeutral && hitVector->empty()){
201 ATH_MSG_DEBUG("forcing neutral charge in FastCaloSimCaloTransportation::caloHits");
202 transport(truth, true);
203 }
204 // Don't expect this ever to happen. Nevertheless, error handling should be improved.
205 // This may require changes in periphery (adjustments after setting function type to StatusCode)
206 else if(hitVector->empty()) ATH_MSG_ERROR("Empty hitVector even after forcing neutral charge. This may cause a segfault soon.");
207
208
209 std::vector<G4FieldTrack> caloSteps = convertToFieldTrack(*hitVector);
210
211 return caloSteps;
212
213}
214
215
216std::vector<G4FieldTrack> FastCaloSimCaloTransportation::convertToFieldTrack(const std::vector<Trk::HitInfo>& vec) const{
217
218 std::vector<G4FieldTrack> caloSteps;
219 for (auto& step : vec){
220 G4FieldTrack track = G4FieldTrack(' ');
221 track.SetPosition(Amg::EigenToHep3Vector(step.trackParms->position()));
222 track.SetMomentum(Amg::EigenToHep3Vector(step.trackParms->momentum()));
223 track.SetChargeAndMoments(step.trackParms->charge());
224 caloSteps.push_back(track);
225 }
226
227 return caloSteps;
228
229}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
double charge(const T &p)
Definition AtlasPID.h:997
std::vector< size_t > vec
std::vector< FPGATrackSimHit > hitVector
static AtlasG4EventUserInfo * GetEventUserInfo()
virtual std::vector< G4FieldTrack > transport(const TFCSTruthState *truth, bool forceNeutral=false) const override final
Trk::PdgToParticleHypothesis m_pdgToParticleHypothesis
virtual StatusCode initialize() override final
CxxUtils::CachedPointer< const Trk::TrackingVolume > m_caloEntrance
FastCaloSimCaloTransportation(const std::string &t, const std::string &n, const IInterface *p)
virtual StatusCode finalize() override final
PublicToolHandle< Trk::ITimedExtrapolator > m_extrapolator
std::vector< G4FieldTrack > convertToFieldTrack(const std::vector< Trk::HitInfo > &vec) const
int pdgid() const
const TLorentzVector & vertex() const
std::unique_ptr< ParametersBase< DIM, T > > uniqueClone() const
clone method for polymorphic deep copy returning unique_ptr; it is not overriden, but uses the existi...
CLHEP::Hep3Vector EigenToHep3Vector(const Amg::Vector3D &eigenvector)
Converts an Eigen-based Amg::Vector3D into a CLHEP-based CLHEP::Hep3Vector.
Eigen::Matrix< double, 3, 1 > Vector3D
@ alongMomentum
CurvilinearParametersT< TrackParametersDim, Charged, PlaneSurface > CurvilinearParameters
ParticleHypothesis
Enumeration for Particle hypothesis respecting the interaction with material.