ATLAS Offline Software
Loading...
Searching...
No Matches
CylinderLayer.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
6// CylinderLayer.cxx, (c) ATLAS Detector software
8
9// Trk
11
19// Amg
21
24 std::shared_ptr<const Trk::CylinderBounds> cbounds,
25 const Trk::LayerMaterialProperties& laymatprop,
26 double thickness,
27 std::unique_ptr<Trk::OverlapDescriptor> olap,
28 int laytyp)
29 : CylinderSurface(transform, std::move(cbounds))
30 , Layer(laymatprop, thickness, std::move(olap), laytyp)
31 , m_approachDescriptor(nullptr)
32{
34}
35
38 const Trk::LayerMaterialProperties& laymatprop,
39 double thickness,
40 std::unique_ptr<Trk::OverlapDescriptor> olap,
41 int laytyp)
42 : CylinderSurface(*cyl)
43 , Layer(laymatprop, thickness, std::move(olap), laytyp)
44 , m_approachDescriptor(nullptr)
45{
47}
48
51 std::shared_ptr<const Trk::CylinderBounds> cbounds,
52 std::unique_ptr<Trk::SurfaceArray> surfaceArray,
53 double thickness,
54 std::unique_ptr<Trk::OverlapDescriptor> olap,
56 int laytyp)
57 : CylinderSurface(transform, std::move(cbounds))
58 , Layer(std::move(surfaceArray), thickness, std::move(olap), laytyp)
60{
62 if (!ades && m_surfaceArray)
64 // register the layer
65 if (ades)
66 m_approachDescriptor->registerLayer(*this);
67}
68
71 std::shared_ptr<const Trk::CylinderBounds> cbounds,
72 std::unique_ptr<Trk::SurfaceArray> surfaceArray,
73 const Trk::LayerMaterialProperties& laymatprop,
74 double thickness,
75 std::unique_ptr<Trk::OverlapDescriptor> olap,
77 int laytyp)
78 : CylinderSurface(transform, std::move(cbounds))
79 , Layer(std::move(surfaceArray), laymatprop, thickness, std::move(olap), laytyp)
81{
83 if (!ades && m_surfaceArray)
85 // register the layer
86 if (ades)
87 m_approachDescriptor->registerLayer(*this);
88}
89
91 std::shared_ptr<const Trk::CylinderBounds> cbounds,
92 const Trk::LayerMaterialProperties& laymatprop,
93 double thickness,
94 std::unique_ptr<Trk::OverlapDescriptor> olap,
95 int laytyp)
96 : CylinderSurface(std::move(cbounds)),
97 Layer(laymatprop, thickness, std::move(olap), laytyp),
98 m_approachDescriptor(nullptr) {
100}
101
103 std::shared_ptr<const Trk::CylinderBounds> cbounds,
104 std::unique_ptr<Trk::SurfaceArray> surfaceArray,
105 double thickness,
106 std::unique_ptr<Trk::OverlapDescriptor> olap,
107 Trk::IApproachDescriptor* ades, int laytyp)
108 : CylinderSurface(std::move(cbounds)),
109 Layer(std::move(surfaceArray), thickness, std::move(olap), laytyp),
113 // register the layer
114 if (ades) m_approachDescriptor->registerLayer(*this);
115}
116
118 std::shared_ptr<const Trk::CylinderBounds> cbounds,
119 std::unique_ptr<Trk::SurfaceArray> surfaceArray,
120 const Trk::LayerMaterialProperties& laymatprop,
121 double thickness,
122 std::unique_ptr<Trk::OverlapDescriptor> olap,
123 Trk::IApproachDescriptor* ades, int laytyp)
124 : CylinderSurface(std::move(cbounds)),
125 Layer(std::move(surfaceArray), laymatprop, thickness, std::move(olap), laytyp),
129 // register the layer
130 if (ades) m_approachDescriptor->registerLayer(*this);
131}
132
142
144 const Amg::Transform3D& transf)
145 : CylinderSurface(clay, transf),
146 Layer(clay),
147 m_approachDescriptor(nullptr) {
148 if (m_surfaceArray) {
150 }
151}
152
154 if (this != &clay) {
155 // call the assignments of the base classes
160 }
161 return (*this);
162}
163
166{
167 return (*this);
168}
169
172{
173 return (*this);
174}
175
177 const Trk::TrackParameters& parm, Trk::PropDirection dir) const {
179 // calculate the direction to the normal
180 const Amg::Vector3D& parmPos = parm.position();
181 Amg::Vector3D pastStep(parmPos + dir * parm.momentum().normalized());
182 if (pastStep.perp() > parm.position().perp())
183 return Trk::Layer::m_layerMaterialProperties->alongPreFactor();
184 return Trk::Layer::m_layerMaterialProperties->oppositePreFactor();
185}
186
188 const Trk::TrackParameters& parm, Trk::PropDirection dir) const {
190 const Amg::Vector3D& parmPos = parm.position();
191 Amg::Vector3D pastStep(parmPos + dir * parm.momentum().normalized());
192 if (pastStep.perp() > parm.position().perp())
193 return Trk::Layer::m_layerMaterialProperties->alongPostFactor();
194 return Trk::Layer::m_layerMaterialProperties->oppositePostFactor();
195}
196
199 std::make_unique<Transforms>(shift * (m_transforms->transform));
200
201 if (m_approachDescriptor && m_approachDescriptor->rebuild()) {
202 // build the new approach descriptor - deletes the current one
204 }
205}
206
208 double envelope) {
209 // only do this if the volume bounds a CylinderVolumeBounds
210 const Trk::CylinderVolumeBounds* cvb =
211 dynamic_cast<const Trk::CylinderVolumeBounds*>(&bounds);
212 if (cvb) {
213 // get the dimensions
214 double hLengthZ = cvb->halflengthZ();
215 double r = surfaceRepresentation().bounds().r();
216 // (0) first, resize the layer itself
217 Trk::CylinderBounds* rCylinderBounds =
218 new Trk::CylinderBounds(r, hLengthZ - envelope);
220 std::shared_ptr<const Trk::CylinderBounds>(rCylinderBounds);
221 // (1) resize the material properties by updating the BinUtility, assuming
222 // rphi/z binning
224 const BinUtility* layerMaterialBU =
226 if (layerMaterialBU && layerMaterialBU->dimensions() > 1) {
227 size_t binsRPhi = layerMaterialBU->max(0) + 1;
228 size_t binsZ = layerMaterialBU->max(1) + 1;
229 // create a new binning with the new dimensions
230 Trk::BinUtility* rBinUtility = new Trk::BinUtility(
231 binsRPhi, -r * M_PI, r * M_PI, Trk::closed, Trk::binRPhi);
232 (*rBinUtility) +=
233 Trk::BinUtility(binsZ, -hLengthZ + envelope, hLengthZ - envelope,
235 Trk::Layer::m_layerMaterialProperties->updateBinning(rBinUtility);
236 }
237 }
238 }
239
240 if (m_approachDescriptor && m_approachDescriptor->rebuild()) {
241 // build the approach descriptor - delete the current approach descriptor
243 }
244}
245
249 const Amg::Vector3D& pos, const Amg::Vector3D& dir,
250 const Trk::BoundaryCheck& bcheck) const {
252 // get the test surfaces from the approach Descriptor
253 const Trk::ApproachSurfaces* surfacesOnApproach =
254 m_approachDescriptor->approachSurfaces(pos, dir);
255 if (surfacesOnApproach) {
256 // test the intersections and go
257 std::vector<Trk::Intersection> sfIntersections;
258 const Trk::Surface* aSurface = nullptr;
259 double aPathLength = 10e10;
260 // get the surfaces
261 for (const auto& sfIter : (*surfacesOnApproach)) {
262 // get the intersection with the surface
263 Trk::Intersection sIntersection =
264 sfIter->straightLineIntersection(pos, dir, true, bcheck);
265 // validation
266 if (sIntersection.valid && sIntersection.pathLength < aPathLength) {
267 aPathLength = sIntersection.pathLength;
268 aSurface = sfIter;
269 }
270 }
271 if (aSurface) return (*aSurface);
272 }
273 }
274 return surfaceRepresentation();
275}
276
280 const Amg::Vector3D& pos, const Amg::Vector3D& mom, Trk::PropDirection pDir,
281 const Trk::BoundaryCheck& bcheck, bool resolveSubSurfaces,
282 const Trk::ICompatibilityEstimator*) const {
283 // resolve the surfaces
284 if (m_approachDescriptor && resolveSubSurfaces) {
285 // resolve based on straight line intersection
286 return approachSurface(pos, double(pDir) * mom.unit(), bcheck);
287 }
288 return surfaceRepresentation();
289}
290
293 // delete it
294 // delete the surfaces
295 auto aSurfaces = std::make_unique<Trk::ApproachSurfaces>();
296 // create the new surfaces
297 if (m_transforms) {
298 Amg::Transform3D asTransform = m_transforms->transform;
299 aSurfaces->push_back(new Trk::CylinderSurface(
300 asTransform, m_bounds->r() - 0.5 * thickness(), m_bounds->halflengthZ()));
301 aSurfaces->push_back(new Trk::CylinderSurface(
302 asTransform, m_bounds->r() + 0.5 * thickness(), m_bounds->halflengthZ()));
303 } else {
304 aSurfaces->push_back(new Trk::CylinderSurface(
305 m_bounds->r() - 0.5 * thickness(), m_bounds->halflengthZ()));
306 aSurfaces->push_back(new Trk::CylinderSurface(
307 m_bounds->r() + 0.5 * thickness(), m_bounds->halflengthZ()));
308 }
309 // set the layer and make TGOwn
310 for (auto& sIter : (*aSurfaces)) {
311 sIter->associateLayer(*this);
312 sIter->setOwner(Trk::TGOwn);
313 }
315 std::make_unique<Trk::ApproachDescriptor>(std::move(aSurfaces));
316}
317
319 const Amg::Vector3D& vCenter,
320 double envelope) {
321 // resize first of all
322 resizeLayer(vBounds, envelope);
323 // now reposition to the potentially center if necessary, do not change layers
324 // with no transform
325 if (Trk::CylinderSurface::m_transforms || center().isApprox(vCenter)) {
326 return;
327 }
328 Amg::Transform3D transf(vCenter);
330 std::make_unique<Transforms>(Amg::Transform3D(vCenter), vCenter);
331 // rebuild approaching layers if needed
334}
335
#define M_PI
just implement the delete on the objects
A generic symmetric BinUtility, for fully symmetric binning in terms of binning grid and binning type...
Definition BinUtility.h:39
size_t max(size_t ba=0) const
First bin maximal value.
Definition BinUtility.h:212
size_t dimensions() const
First bin maximal value.
Definition BinUtility.h:209
The BoundaryCheck class allows to steer the way surface boundaries are used for inside/outside checks...
Bounds for a cylindrical Surface.
Class to describe a cylindrical detector layer for tracking, it inhertis from both,...
virtual double postUpdateMaterialFactor(const Trk::TrackParameters &par, Trk::PropDirection dir) const override final
getting the MaterialProperties back - for post-update
virtual void resizeAndRepositionLayer(const VolumeBounds &vBounds, const Amg::Vector3D &cCenter, double envelope) override final
Resize the layer to the tracking volume.
virtual const Surface & surfaceOnApproach(const Amg::Vector3D &pos, const Amg::Vector3D &dir, PropDirection pdir, const BoundaryCheck &bcheck, bool resolveSubSurfaces=0, const ICompatibilityEstimator *ice=nullptr) const override final
Surface seen on approach - if not defined differently, it is the surfaceRepresentation()
CylinderLayer()
Default Constructor.
CylinderLayer & operator=(const CylinderLayer &)
Assignment operator for CylinderLayers.
const Surface & approachSurface(const Amg::Vector3D &pos, const Amg::Vector3D &dir, const BoundaryCheck &bcheck) const
Surface seen on approach - if not defined differently, it is the surfaceRepresentation()
void buildApproachDescriptor()
build approach surfaces
virtual void moveLayer(Amg::Transform3D &shift) override final
move the Layer
virtual double preUpdateMaterialFactor(const Trk::TrackParameters &par, Trk::PropDirection dir) const override final
getting the MaterialProperties back - for pre-update
virtual void resizeLayer(const VolumeBounds &vBounds, double envelope) override final
Resize the layer to the tracking volume - only works for CylinderVolumeBouns.
std::unique_ptr< IApproachDescriptor > m_approachDescriptor
surfaces on approach to the layer
virtual const CylinderSurface & surfaceRepresentation() const override final
Transforms the layer into a Surface representation for extrapolation.
Class for a CylinderSurface in the ATLAS detector.
CylinderSurface()
Default Constructor.
std::shared_ptr< const CylinderBounds > m_bounds
The global reference point (== a point on the surface)
CylinderSurface & operator=(const CylinderSurface &csf)
Assignment operator.
virtual const CylinderBounds & bounds() const override final
This method returns the CylinderBounds by reference (NoBounds is not possible for cylinder)
Bounds for a cylindrical Volume, the decomposeToSurfaces method creates a vector of up to 6 surfaces:
double halflengthZ() const
This method returns the halflengthZ.
CVirtual class to decide and return which approaching surface to be taken.
This virtual base class encapsulates the logics to build pre/post/full update material for Layer stru...
Layer()=default
Default Constructor.
const SurfaceArray * surfaceArray() const
Return the entire SurfaceArray, returns nullptr if no SurfaceArray.
double thickness() const
Return the Thickness of the Layer.
Layer & operator=(const Layer &lay)
Assignment operator for Derived classes.
Definition Layer.cxx:84
std::unique_ptr< LayerMaterialProperties > m_layerMaterialProperties
thickness of the Layer
Definition Layer.h:288
std::unique_ptr< SurfaceArray > m_surfaceArray
MaterialPoperties of this layer Surface.
Definition Layer.h:286
const Amg::Vector3D & momentum() const
Access method for the momentum.
const Amg::Vector3D & position() const
Access method for the position.
Abstract Base Class for tracking surfaces.
void associateLayer(const Layer &lay)
method to associate a Trk::Layer.
const Amg::Transform3D & transform() const
Returns HepGeom::Transform3D by reference.
std::unique_ptr< Transforms > m_transforms
Unique Pointer to the Transforms struct.
const Amg::Vector3D & center() const
Returns the center position of the Surface.
Pure Absract Base Class for Volume bounds.
int r
Definition globals.cxx:22
Eigen::Affine3d Transform3D
Eigen::Matrix< double, 3, 1 > Vector3D
PropDirection
PropDirection, enum for direction of the propagation.
@ open
Definition BinningType.h:40
@ closed
Definition BinningType.h:41
ParametersBase< TrackParametersDim, Charged > TrackParameters
@ binRPhi
Definition BinningType.h:52
@ binZ
Definition BinningType.h:49
STL namespace.