ATLAS Offline Software
Loading...
Searching...
No Matches
SiHitCollectionCnv_p3.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
10
11#include <cmath>
12
13//CLHEP
14#include "CLHEP/Geometry/Point3D.h"
15// Gaudi
16#include "GaudiKernel/MsgStream.h"
17#include "GaudiKernel/ThreadLocalContext.h"
18
19// Athena
23
24// * * * stolen from eflowRec * * * //
25inline double phicorr(double a)
26{
27 if (a <= -M_PI)
28 {
29 return a+(2*M_PI*floor(-(a-M_PI)/(2*M_PI)));
30 }
31 else if (a > M_PI)
32 {
33 return a-(2*M_PI*floor((a+M_PI)/(2*M_PI)));
34 }
35 else
36 {
37 return a;
38 }
39}
40
41// * * * stolen from eflowRec * * * //
42inline double cycle(double a, double b)
43{
44 double del = b-a;
45 if (del > M_PI)
46 {
47 return a+2.0*M_PI;
48 }
49 else if (del < -M_PI)
50 {
51 return a-2.0*M_PI;
52 }
53 else
54 {
55 return a;
56 }
57}
58
59
60const double SiHitCollectionCnv_p3::m_persEneUnit = 1.0e-5;
61const double SiHitCollectionCnv_p3::m_persLenUnit = 1.0e-5;
62const double SiHitCollectionCnv_p3::m_persAngUnit = 1.0e-5;
63const double SiHitCollectionCnv_p3::m_2bHalfMaximum = pow(2.0, 15.0);
64const int SiHitCollectionCnv_p3::m_2bMaximum = (unsigned short)(-1);
65
66
67void SiHitCollectionCnv_p3::transToPers(const SiHitCollection* transCont, SiHitCollection_p3* persCont, MsgStream &log)
68{
69 // Finds hits belonging to a "string" (in which the end point of one hit is the same as the start point of the next) and
70 // persistifies the end point of each hit plus the start point of the first hit in each string.
71 //
72 // Further compression is achieved by optimising the storage of the position vectors:- start (x,y,z) and (theta,phi) of
73 // first hit are stored as floats, (delta_theta,delta_phi) relative to the fisrst hit are stored as 2 byte numbers and
74 // used to specify the hit direction. All hit lengths are stored as 2 byte numbers.
75 //
76 // Additional savings are achieved by storing the energy loss for each hit as a 2 byte number and only storing the mean
77 // time of the first hit per string.
78 //
79 // See http://indico.cern.ch/getFile.py/access?contribId=11&resId=2&materialId=slides&confId=30893 for more info.
80
81 static const double dRcut = 1.0e-7;
82 static const double dTcut = 1.0;
83
84 const EventContext& ctx = Gaudi::Hive::currentContext();
85 const IProxyDict* proxy = Atlas::getExtendedEventContext(ctx).proxy();
86 int lastIndex{-1};
87 int lastBarcode{-1};
88 int lastId = -1;
89 double stringFirstTheta = 0.0;
90 double stringFirstPhi = 0.0;
91 double lastT = 0.0;
92 double persSumE = 0.0;
93 double transSumE = 0.0;
94 unsigned int idx = 0;
95 unsigned int endBC = 0;
96 unsigned int endId = 0;
97 unsigned int endHit = 0;
98 HepGeom::Point3D<double> lastTransEnd(0.0, 0.0, 0.0);
99 HepGeom::Point3D<double> lastPersEnd(0.0, 0.0, 0.0);
100
101 for (SiHitCollection::const_iterator it = transCont->begin(); it != transCont->end(); ++it) {
102
104 const HepMcParticleLink * currentLink = &(siHit->particleLink());
105 const int barcode = currentLink->barcode();
106 int index{0};
108 proxy).at(0) != 0) {
109 index = currentLink->eventIndex();
110 }
111
112 if ( lastBarcode != barcode || lastIndex != index || (idx-endBC)==USHRT_MAX ) {
113 lastBarcode = barcode;
114 lastIndex = index;
115 const unsigned short persIndex = static_cast<unsigned short>(index);
116 if (static_cast<HepMcParticleLink::index_type>(lastIndex) != static_cast<HepMcParticleLink::index_type>(persIndex)) {
117 log << MSG::WARNING << "Attempting to persistify an eventIndex larger than max unsigned short!" << endmsg;
118 }
119 // store barcode, eventIndex and McEventCollection once for set
120 // of consecutive hits with the same barcode and eventIndex
121 persCont->m_barcode.push_back(static_cast<unsigned long>(lastBarcode));
122 persCont->m_mcEvtIndex.push_back(persIndex);
123 persCont->m_evtColl.push_back('a'); // Hard-coding as this only ever had a single value in production
124
125 if (idx > 0) {
126 persCont->m_nBC.push_back(idx - endBC);
127 endBC = idx;
128 }
129 }
130
131 if ( ( (int)siHit->identify() != lastId ) || (idx-endId)==USHRT_MAX) {
132
133 // store SiHitIdentifier once for set of consecutive hits with
134 // the same SiHitIdentifier
135
136 lastId = siHit->identify();
137 persCont->m_id.push_back(lastId);
138
139 if (idx > 0) {
140 persCont->m_nId.push_back(idx - endId);
141 endId = idx;
142 }
143 }
144
145 HepGeom::Point3D<double> st = siHit->localStartPosition();
146 HepGeom::Point3D<double> en = siHit->localEndPosition();
147
148 const double dx = st.x() - lastTransEnd.x();
149 const double dy = st.y() - lastTransEnd.y();
150 const double dz = st.z() - lastTransEnd.z();
151 const double t = siHit->meanTime();
152
153 const double dRLast = sqrt(dx * dx + dy * dy + dz * dz); // dR between end of previous hit and start of current one
154 const double dTLast = fabs(t - lastT);
155
156 CLHEP::Hep3Vector direction(0.0, 0.0, 0.0);
157 double theta = 0.0;
158 double phi = 0.0;
159 bool startNewString = (dRLast >= dRcut || dTLast >= dTcut || (idx - endHit) == USHRT_MAX);
160
161 if (!startNewString) {
162
163 // hit is part of existing string
164
165 direction = CLHEP::Hep3Vector( en.x() - lastPersEnd.x(), en.y() - lastPersEnd.y(), en.z() - lastPersEnd.z() );
166
167 theta = direction.theta();
168 phi = phicorr( direction.phi() );
169
170 const int dTheta_2b = (int)( (theta - stringFirstTheta) / m_persAngUnit + m_2bHalfMaximum + 0.5 );
171 const int dPhi_2b = (int)( (cycle(phi, stringFirstPhi) - stringFirstPhi) / m_persAngUnit + m_2bHalfMaximum + 0.5 );
172
173 if ( dTheta_2b < m_2bMaximum && dTheta_2b >= 0 && dPhi_2b < m_2bMaximum && dPhi_2b >= 0) {
174 persCont->m_dTheta.push_back(dTheta_2b);
175 persCont->m_dPhi.push_back(dPhi_2b);
176 theta = stringFirstTheta + ( (double)dTheta_2b - m_2bHalfMaximum ) * m_persAngUnit;
177 phi = stringFirstPhi + ( (double)dPhi_2b - m_2bHalfMaximum ) * m_persAngUnit;
178 phi = phicorr(phi);
179 }
180 else {
181 startNewString = true;
182 }
183 }
184
185 if (startNewString) {
186
187 // begin new hit string
188
189 direction = CLHEP::Hep3Vector( en.x() - st.x(), en.y() - st.y(), en.z() - st.z() );
190
191 theta = direction.theta();
192 phi = phicorr( direction.phi() );
193
194 persCont->m_hit1_meanTime.push_back(t);
195 persCont->m_hit1_x0.push_back(st.x());
196 persCont->m_hit1_y0.push_back(st.y());
197 persCont->m_hit1_z0.push_back(st.z());
198 persCont->m_hit1_theta.push_back(theta);
199 persCont->m_hit1_phi.push_back(phi);
200
201 lastPersEnd = HepGeom::Point3D<double>(st.x(), st.y(), st.z());
202
203 stringFirstTheta = theta;
204 stringFirstPhi = phi;
205
206 if (idx > 0) {
207 persCont->m_nHits.push_back(idx - endHit);
208 endHit = idx;
209 }
210 }
211
212 lastTransEnd = HepGeom::Point3D<double>(en.x(), en.y(), en.z());
213 transSumE += siHit->energyLoss();
214
215 const int eneLoss_2b = (int)((transSumE - persSumE) / m_persEneUnit + 0.5); // calculated to allow recovery sum over
216 // whole hit string to chosen precision
217
218 const int hitLength_2b = (int)(direction.mag() / m_persLenUnit + 0.5); // calculated to give the correct position to
219 // the chosen precision, NOT the length of the
220 // hit (small difference in practice).
221 double eneLoss = 0.0;
222
223 if (eneLoss_2b >= m_2bMaximum) {
224 eneLoss = siHit->energyLoss();
225 persCont->m_hitEne_2b.push_back(m_2bMaximum);
226 persCont->m_hitEne_4b.push_back(eneLoss);
227 }
228 else {
229 eneLoss = eneLoss_2b * m_persEneUnit;
230 persCont->m_hitEne_2b.push_back(eneLoss_2b);
231 }
232
233 double length = 0.0;
234
235 if (hitLength_2b >= m_2bMaximum) {
236 length = direction.mag();
237 persCont->m_hitLength_2b.push_back(m_2bMaximum);
238 persCont->m_hitLength_4b.push_back(direction.mag());
239 }
240 else {
241 length = hitLength_2b * m_persLenUnit;
242 persCont->m_hitLength_2b.push_back(hitLength_2b);
243 }
244
245 CLHEP::Hep3Vector persDir(length, 0.0, 0.0);
246 persDir.setTheta(theta);
247 persDir.setPhi(phi);
248
249 lastPersEnd = (CLHEP::Hep3Vector)lastPersEnd + persDir;
250 persSumE += eneLoss;
251 lastT = t;
252
253 ++idx;
254 }
255
256 persCont->m_nBC.push_back(idx - endBC);
257 persCont->m_nId.push_back(idx - endId);
258 persCont->m_nHits.push_back(idx - endHit);
259#ifdef ENABLE_SANITY_CHECKS
260 // Sanity check
261 const unsigned int init(0);
262 const unsigned int transContSize = transCont->size();
263 if (std::accumulate(persCont->m_nBC.begin(), persCont->m_nBC.end(), init)!=transContSize) {
264 log << MSG::ERROR << "transToPers: sum of entries of persCont->m_nBC (" << std::accumulate(persCont->m_nBC.begin(), persCont->m_nBC.end(), init) << ") does not match transient container size = " << transContSize << endmsg;
265 }
266 if (std::accumulate(persCont->m_nId.begin(), persCont->m_nId.end(), init)!=transContSize) {
267 log << MSG::ERROR << "transToPers: sum of entries of persCont->m_nId (" << std::accumulate(persCont->m_nId.begin(), persCont->m_nId.end(), init) << ") does not match transient container size = " << transContSize << endmsg;
268 }
269 if (std::accumulate(persCont->m_nHits.begin(), persCont->m_nHits.end(), init)!=transContSize) {
270 log << MSG::ERROR << "transToPers: sum of entries of persCont->m_nHits (" << std::accumulate(persCont->m_nHits.begin(), persCont->m_nHits.end(), init) << ") does not match transient container size = " << transContSize << endmsg;
271 }
272#endif
273}
274
275
277 std::unique_ptr<SiHitCollection> trans(std::make_unique<SiHitCollection>("DefaultCollectionName",persObj->m_nHits.size()));
278 persToTrans(persObj, trans.get(), log);
279 return(trans.release());
280}
281
282
283#ifdef ENABLE_SANITY_CHECKS
284void SiHitCollectionCnv_p3::persToTrans(const SiHitCollection_p3* persCont, SiHitCollection* transCont, MsgStream &log)
285#else
286void SiHitCollectionCnv_p3::persToTrans(const SiHitCollection_p3* persCont, SiHitCollection* transCont, MsgStream &/*log*/)
287#endif
288{
289 const EventContext& ctx = Gaudi::Hive::currentContext();
290#ifdef ENABLE_SANITY_CHECKS
291 // Sanity check
292 const unsigned int transContSize = persCont->m_hitEne_2b.size(); // this vector has one entry per transient SiHit, so its size can be used as a proxy for the transient Container size
293 const unsigned int init(0);
294 if (std::accumulate(persCont->m_nBC.begin(), persCont->m_nBC.end(), init)!=transContSize) {
295 log << MSG::ERROR << "persToTrans: sum of entries of persCont->m_nBC (" << std::accumulate(persCont->m_nBC.begin(), persCont->m_nBC.end(), init) << ") does not match transient container size = " << transContSize << endmsg;
296 }
297 if (std::accumulate(persCont->m_nId.begin(), persCont->m_nId.end(), init)!=transContSize) {
298 log << MSG::ERROR << "persToTrans: sum of entries of persCont->m_nId (" << std::accumulate(persCont->m_nId.begin(), persCont->m_nId.end(), init) << ") does not match transient container size = " << transContSize << endmsg;
299 }
300 if (std::accumulate(persCont->m_nHits.begin(), persCont->m_nHits.end(), init)!=transContSize) {
301 log << MSG::ERROR << "persToTrans: sum of entries of persCont->m_nHits (" << std::accumulate(persCont->m_nHits.begin(), persCont->m_nHits.end(), init) << ") does not match transient container size = " << transContSize << endmsg;
302 }
303#endif
304
305 unsigned int hitCount = 0;
306 unsigned int angleCount = 0;
307 unsigned int idxBC = 0;
308 unsigned int idxId = 0;
309 unsigned int idxEne4b = 0;
310 unsigned int idxLen4b = 0;
311 unsigned int endHit = 0;
312 unsigned int endBC = 0;
313 unsigned int endId = 0;
314
315 for (unsigned int i = 0; i < persCont->m_nHits.size(); i++) {
316
317 if (persCont->m_nHits[i]) {
318
319 const unsigned int start = endHit;
320 endHit += persCont->m_nHits[i];
321
322 const double t0 = persCont->m_hit1_meanTime[i];
323 const double theta0 = persCont->m_hit1_theta[i];
324 const double phi0 = persCont->m_hit1_phi[i];
325 HepGeom::Point3D<double> endLast(persCont->m_hit1_x0[i], persCont->m_hit1_y0[i], persCont->m_hit1_z0[i]);
326
327 for (unsigned int j = start; j < endHit; j++) {
328
329 if (j >= endBC + persCont->m_nBC[idxBC])
330 endBC += persCont->m_nBC[idxBC++];
331
332 if (j >= endId + persCont->m_nId[idxId])
333 endId += persCont->m_nId[idxId++];
334
335 const double eneLoss_2b = persCont->m_hitEne_2b[hitCount];
336 const double hitLength_2b = persCont->m_hitLength_2b[hitCount];
337
338 const double eneLoss = (eneLoss_2b < m_2bMaximum) ? eneLoss_2b * m_persEneUnit : persCont->m_hitEne_4b[idxEne4b++];
339 const double length = (hitLength_2b < m_2bMaximum) ? hitLength_2b * m_persLenUnit : persCont->m_hitLength_4b[idxLen4b++];
340
341 const double dTheta = (j > start) ? ((double)persCont->m_dTheta[angleCount] - m_2bHalfMaximum) * m_persAngUnit : 0.0;
342 const double dPhi = (j > start) ? ((double)persCont->m_dPhi[angleCount] - m_2bHalfMaximum) * m_persAngUnit : 0.0;
343
344 const double meanTime = t0;
345 const double theta = theta0 + dTheta;
346 const double phi = phicorr(phi0 + dPhi);
347
348 CLHEP::Hep3Vector r(length, 0.0, 0.0);
349 r.setTheta(theta);
350 r.setPhi(phi);
351
352 HepGeom::Point3D<double> endThis( endLast + r );
353
355 if (persCont->m_mcEvtIndex[idxBC] == 0) {
357 }
358 HepMcParticleLink partLink( persCont->m_barcode[idxBC], persCont->m_mcEvtIndex[idxBC], flag, HepMcParticleLink::IS_BARCODE, ctx );
359 if ( HepMC::BarcodeBased::is_truth_suppressed_pileup(static_cast<int>(persCont->m_barcode[idxBC])) ) {
361 }
362 transCont->Emplace( endLast, endThis, eneLoss, meanTime, partLink, persCont->m_id[idxId]);
363
364 endLast = std::move(endThis);
365
366 ++hitCount;
367 if (j > start) ++angleCount;
368 }
369 }
370 }
371}
#define M_PI
Scalar phi() const
phi method
Scalar theta() const
theta method
#define endmsg
double length(const pvec &v)
static Double_t a
static Double_t t0
@ EBC_PU_SUPPRESSED
double phicorr(double a)
double cycle(double a, double b)
double phicorr(double a)
double cycle(double a, double b)
AtlasHitsVector< SiHit > SiHitCollection
constexpr int pow(int base, int exp) noexcept
CONT::const_iterator const_iterator
const_iterator begin() const
void Emplace(Args &&... args)
size_type size() const
const_iterator end() const
static const double m_2bHalfMaximum
static const double m_persEneUnit
static const double m_persLenUnit
static const double m_persAngUnit
virtual void transToPers(const SiHitCollection *transCont, SiHitCollection_p3 *persCont, MsgStream &log)
virtual SiHitCollection * createTransient(const SiHitCollection_p3 *persObj, MsgStream &log)
virtual void persToTrans(const SiHitCollection_p3 *persCont, SiHitCollection *transCont, MsgStream &log)
std::vector< unsigned short > m_nBC
std::vector< float > m_hit1_theta
std::vector< float > m_hit1_meanTime
std::vector< float > m_hit1_y0
std::vector< float > m_hit1_x0
std::vector< char > m_evtColl
std::vector< unsigned long > m_barcode
std::vector< float > m_hitEne_4b
std::vector< unsigned short > m_mcEvtIndex
std::vector< unsigned short > m_nId
std::vector< unsigned short > m_dPhi
std::vector< float > m_hitLength_4b
std::vector< unsigned short > m_nHits
std::vector< unsigned short > m_dTheta
std::vector< unsigned short > m_hitLength_2b
std::vector< float > m_hit1_z0
std::vector< unsigned short > m_hitEne_2b
std::vector< float > m_hit1_phi
std::vector< unsigned long > m_id
int r
Definition globals.cxx:22
const ExtendedEventContext & getExtendedEventContext(const EventContext &ctx)
Retrieve an extended context from a context object.
bool is_truth_suppressed_pileup(const T &p)
Method to establish if a particle (or barcode) corresponds to truth-suppressed pile-up.
Definition index.py:1