ATLAS Offline Software
Loading...
Searching...
No Matches
DiscLayer.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// DiscLayer.cxx, (c) ATLAS Detector software
8
9// Trk
11
18// CLHEP
20
22 std::shared_ptr<const Trk::DiscBounds> dbounds,
23 const Trk::LayerMaterialProperties& laymatprop,
24 double thickness,
25 std::unique_ptr<Trk::OverlapDescriptor> olap,
26 int laytyp)
27 : DiscSurface(transform, std::move(dbounds))
28 , Layer(laymatprop, thickness, std::move(olap), laytyp)
29 , m_approachDescriptor(nullptr)
30{
32}
33
35 const Trk::LayerMaterialProperties& laymatprop,
36 double thickness,
37 std::unique_ptr<Trk::OverlapDescriptor> olap,
38 int laytyp)
39 : DiscSurface(*disc)
40 , Layer(laymatprop, thickness, std::move(olap), laytyp)
41 , m_approachDescriptor(nullptr)
42{
44}
45
47 std::shared_ptr<const Trk::DiscBounds> dbounds,
48 std::unique_ptr<Trk::SurfaceArray> surfaceArray,
49 double thickness,
50 std::unique_ptr<Trk::OverlapDescriptor> olap,
52 int laytyp)
53 : DiscSurface(transform, std::move(dbounds))
54 , Layer(std::move(surfaceArray), thickness, std::move(olap), laytyp)
56{
59 // register the layer
60 if (ades) m_approachDescriptor->registerLayer(*this);
61}
62
64 std::shared_ptr<const Trk::DiscBounds> dbounds,
65 std::unique_ptr<Trk::SurfaceArray> surfaceArray,
66 const Trk::LayerMaterialProperties& laymatprop,
67 double thickness,
68 std::unique_ptr<Trk::OverlapDescriptor> olap,
70 int laytyp)
71 : DiscSurface(transform, std::move(dbounds))
72 , Layer(std::move(surfaceArray), laymatprop, thickness, std::move(olap), laytyp)
74{
76 if (!ades && m_surfaceArray)
78 // register the layer
79 if (ades)
80 m_approachDescriptor->registerLayer(*this);
81}
82
88
95
97 if (this != &dlay) {
98 // call the assignments of the base classes
101 m_approachDescriptor.reset();
102 if (m_surfaceArray) {
104 }
106 }
107 return (*this);
108}
109
110const Trk::DiscSurface&
112{
113 return (*this);
114}
117{
118 return (*this);
119}
120
122 Trk::PropDirection dir) const {
124 // const Amg::Vector3D& parmPos = parm.position();
125 if (Trk::DiscSurface::normal().dot(dir * parm.momentum().normalized()) > 0.)
126 return Trk::Layer::m_layerMaterialProperties->alongPreFactor();
127 return Trk::Layer::m_layerMaterialProperties->oppositePreFactor();
128}
129
131 const Trk::TrackParameters& parm, Trk::PropDirection dir) const {
133 // const Amg::Vector3D& parmPos = parm.position();
134 if (Trk::DiscSurface::normal().dot(dir * parm.momentum().normalized()) > 0.)
135 return Trk::Layer::m_layerMaterialProperties->alongPostFactor();
136 return Trk::Layer::m_layerMaterialProperties->oppositePostFactor();
137}
138
140 Amg::Transform3D transf = shift * m_transforms->transform;
141 m_transforms = std::make_unique<Transforms>(transf, transf.translation(),
142 transf.rotation().col(2));
143 // rebuild that - deletes the current one
144 if (m_approachDescriptor && m_approachDescriptor->rebuild()) {
146 }
147}
148
149void Trk::DiscLayer::resizeLayer(const VolumeBounds& bounds, double envelope) {
150 // only do this if the volume bounds a CylinderVolumeBounds
151 const Trk::CylinderVolumeBounds* cvb =
152 dynamic_cast<const Trk::CylinderVolumeBounds*>(&bounds);
153 if (cvb) {
154 // get the dimensions
155 double rInner = cvb->innerRadius();
156 double rOuter = cvb->outerRadius();
157 // (0) first, resize the layer itself
158 Trk::DiscBounds* rDiscBounds =
159 new Trk::DiscBounds(rInner + envelope, rOuter - envelope);
161 std::shared_ptr<const Trk::SurfaceBounds>(rDiscBounds);
162 // (1) resize the material properties by updating the BinUtility, assuming
163 // r/phi binning
165 const BinUtility* layerMaterialBU =
167 if (layerMaterialBU && layerMaterialBU->binningValue(0) == Trk::binR) {
168 size_t binsR = layerMaterialBU->max(0) + 1;
169 // create a new binning with the new dimensions
170 Trk::BinUtility* rBinUtility = new Trk::BinUtility(
171 binsR, rInner + envelope, rOuter - envelope, Trk::open, Trk::binR);
172 // copy the second dimension over if exist
173 if (layerMaterialBU->dimensions() > 1) {
174 (*rBinUtility) += Trk::BinUtility(layerMaterialBU->max(1) + 1, -M_PI,
176 }
177 Trk::Layer::m_layerMaterialProperties->updateBinning(rBinUtility);
178 }
179 }
180 }
181
182 if (m_approachDescriptor && m_approachDescriptor->rebuild()) {
183 // rebuild the approach descriptor - delete the current approach descriptor
185 }
186}
187
191 const Amg::Vector3D& pos, const Amg::Vector3D& dir,
192 const Trk::BoundaryCheck& bcheck) const {
194 // get the test surfaces from the approach Descriptor
195 const Trk::ApproachSurfaces* surfacesOnApproach =
196 m_approachDescriptor->approachSurfaces(pos, dir);
197 if (surfacesOnApproach) {
198 // test the intersections and go
199 std::vector<Trk::Intersection> sfIntersections;
200 const Trk::Surface* aSurface = nullptr;
201 double aPathLength = 10e10;
203 for (const auto& sfIter : (*surfacesOnApproach)) {
204 // get the intersection with the surface
205 Trk::Intersection sIntersection =
206 sfIter->straightLineIntersection(pos, dir, true, bcheck);
207 // validation
208 if (sIntersection.valid && sIntersection.pathLength < aPathLength) {
209 aPathLength = sIntersection.pathLength;
210 aSurface = sfIter;
211 }
212 }
213 if (aSurface) return (*aSurface);
214 }
215 }
216 return surfaceRepresentation();
217}
218
222 const Amg::Vector3D& pos, const Amg::Vector3D& mom, Trk::PropDirection pDir,
223 const Trk::BoundaryCheck& bcheck, bool resolveSubSurfaces,
224 const Trk::ICompatibilityEstimator*) const {
225 // resolve the surfaces
226 if (m_approachDescriptor && resolveSubSurfaces) {
227 // resolve based on straight line intersection
228 return approachSurface(pos, double(pDir) * mom.unit(), bcheck);
229 }
230 return surfaceRepresentation();
231}
232
235 // create the surface container
236 auto aSurfaces = std::make_unique<Trk::ApproachSurfaces>();
237 // get the center
238 const auto halfThickness = 0.5 * thickness() * normal();
239 const Amg::Vector3D aspPosition(center() + halfThickness);
240 const Amg::Vector3D asnPosition(center() - halfThickness);
241 // create the new surfaces
242 const Trk::DiscBounds* db =
243 dynamic_cast<const Trk::DiscBounds*>(m_bounds.get());
244 if (db) {
245 // create new surfaces
246 Amg::Transform3D asnTransform =
248 Amg::Transform3D aspTransform =
250 aSurfaces->push_back(new Trk::DiscSurface(aspTransform, std::make_shared<Trk::DiscBounds>(*db)));
251 aSurfaces->push_back(new Trk::DiscSurface(asnTransform, std::make_shared<Trk::DiscBounds>(*db)));
252 // set the layer and make TGOwn
253 for (auto& sIter : (*aSurfaces)) {
254 sIter->associateLayer(*this);
255 sIter->setOwner(Trk::TGOwn);
256 }
257 }
258 // The approach descriptor takes owneship
260 std::make_unique<Trk::ApproachDescriptor>(std::move(aSurfaces));
261}
262
264 const Amg::Vector3D& vCenter,
265 double envelope) {
266 // resize first of all
267 resizeLayer(vBounds, envelope);
268 // now reposition to the potentially center if necessary, do not change layers
269 // with no transform
270 const Trk::CylinderVolumeBounds* cvb =
271 dynamic_cast<const Trk::CylinderVolumeBounds*>(&vBounds);
272 if (cvb) {
273 // get the halflength
274 double hLengthZ = cvb->halflengthZ();
275 Amg::Vector3D nDiscCenter =
276 center().z() < 0.
277 ? Amg::Vector3D(vCenter -
278 Amg::Vector3D(0., 0., hLengthZ - 0.5 * thickness()))
280 vCenter +
281 Amg::Vector3D(0., 0., hLengthZ - 0.5 * thickness()));
282 if (center().isApprox(nDiscCenter)) return;
283 // else set to the new volume center
284 Amg::Transform3D transf((Amg::Translation3D(nDiscCenter)));
285 Amg::Vector3D center(nDiscCenter);
286 m_transforms = std::make_unique<Transforms>(transf, center);
287 }
288 // rebuild the approaching layer
291}
#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
BinningValue binningValue(size_t ba=0) const
The type/value of the binning.
Definition BinUtility.h:230
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 Volume, the decomposeToSurfaces method creates a vector of up to 6 surfaces:
double innerRadius() const
This method returns the inner radius.
double halflengthZ() const
This method returns the halflengthZ.
double outerRadius() const
This method returns the outer radius.
Class to describe the bounds for a planar DiscSurface.
Definition DiscBounds.h:44
Class to describe a disc-like detector layer for tracking, it inhertis from both, Layer base class an...
Definition DiscLayer.h:45
void buildApproachDescriptor()
build approach surfaces
virtual const Surface & surfaceOnApproach(const Amg::Vector3D &pos, const Amg::Vector3D &mom, 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()
virtual void resizeAndRepositionLayer(const VolumeBounds &vBounds, const Amg::Vector3D &cCenter, double envelop) override final
Resize the layer to the tracking volume - not implemented.
DiscLayer()=default
Default Constructor.
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()
virtual void resizeLayer(const VolumeBounds &vBounds, double envelope) override final
Resize the layer to the tracking volume - only works for CylinderVolumeBouns.
DiscLayer & operator=(const DiscLayer &)
Assignment operator for DiscLayers.
Definition DiscLayer.cxx:96
virtual void moveLayer(Amg::Transform3D &shift) override final
move the Layer non-const
std::unique_ptr< IApproachDescriptor > m_approachDescriptor
< surface for approaching
Definition DiscLayer.h:146
virtual double preUpdateMaterialFactor(const Trk::TrackParameters &par, Trk::PropDirection dir) const override final
getting the MaterialProperties back - for pre-update
virtual double postUpdateMaterialFactor(const Trk::TrackParameters &par, Trk::PropDirection dir) const override final
getting the MaterialProperties back - for post-update
virtual const DiscSurface & surfaceRepresentation() const override final
Transforms the layer into a Surface representation for extrapolation.
Class for a DiscSurface in the ATLAS detector.
Definition DiscSurface.h:54
DiscSurface()
Default Constructor.
DiscSurface & operator=(const DiscSurface &dsf)
Assignement operator.
std::shared_ptr< const SurfaceBounds > m_bounds
reference Point on the Surface
const SurfaceBounds & bounds() const override final
This method returns the bounds by reference.
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.
Abstract Base Class for tracking surfaces.
void associateLayer(const Layer &lay)
method to associate a Trk::Layer.
virtual const Amg::Vector3D & normal() const
Returns the normal vector of the Surface (i.e.
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.
Eigen::Affine3d Transform3D
Eigen::Matrix< double, 3, 1 > Vector3D
Eigen::Translation< double, 3 > Translation3D
PropDirection
PropDirection, enum for direction of the propagation.
@ open
Definition BinningType.h:40
@ closed
Definition BinningType.h:41
ParametersBase< TrackParametersDim, Charged > TrackParameters
@ binR
Definition BinningType.h:50
@ binPhi
Definition BinningType.h:51
Definition dot.py:1
STL namespace.