ATLAS Offline Software
Loading...
Searching...
No Matches
PlaneSurface.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
6// PlaneSurface.h, (c) ATLAS Detector software
8
9#ifndef TRKSURFACES_PLANESURFACE_H
10#define TRKSURFACES_PLANESURFACE_H
11
12// Trk
13#include <memory>
16#include "TrkSurfaces/Surface.h"
18// Amg
21
22class MsgStream;
23class Identifier;
24template<class SURFACE, class BOUNDS_CNV>
26template<class SURFACE, class BOUNDS_CNV>
28
29namespace Trk {
30
31class LocalDirection;
32class LocalParameters;
34class RectangleBounds;
35class TriangleBounds;
36class AnnulusBounds;
37class TrapezoidBounds;
39class DiamondBounds;
40class EllipseBounds;
41class CurvilinearUVT;
42template<int DIM, class T, class S>
43class ParametersT;
44
62
63class PlaneSurface : public Surface
64{
65public:
68
71
73 PlaneSurface(const PlaneSurface& psf) = default;
74
76 PlaneSurface& operator=(const PlaneSurface& psf) = default;
77
79 PlaneSurface(PlaneSurface&& psf) noexcept = default;
80
82 PlaneSurface& operator=(PlaneSurface&& psf) noexcept = default;
83
85 virtual ~PlaneSurface() = default;
86
88 PlaneSurface(const PlaneSurface& psf, const Amg::Transform3D& transf);
89
91 PlaneSurface(const Amg::Vector3D& position, const CurvilinearUVT& curvUVT);
92
94 PlaneSurface(const TrkDetElementBase& detelement,
95 const Amg::Transform3D& transf);
96
98 PlaneSurface(const TrkDetElementBase& detelement);
99
102 PlaneSurface(const TrkDetElementBase& detelement,
103 const Identifier& id,
104 const Amg::Transform3D & transf);
105
108 PlaneSurface(const TrkDetElementBase& detelement,
109 const Identifier& id);
110
112 PlaneSurface(const Amg::Transform3D& htrans);
113
114
116 PlaneSurface(const Amg::Transform3D & htrans, double halephi, double haleta);
117
119 PlaneSurface(const Amg::Transform3D & htrans,
120 double minhalephi,
121 double maxhalephi,
122 double haleta);
123
124
126 PlaneSurface(const Amg::Transform3D& htrans,
127 std::shared_ptr<const Trk::SurfaceBounds> sbounds);
128
130 virtual bool operator==(const Surface& sf) const override;
131 bool operator==(const PlaneSurface& cf) const;
132
134 virtual PlaneSurface* clone() const override;
135
137 constexpr virtual SurfaceType type() const override final;
138
142 double l1,
143 double l2,
144 double phi,
145 double theta,
146 double qop,
147 std::optional<AmgSymMatrix(5)> cov = std::nullopt) const override final;
148
152 const Amg::Vector3D& position,
153 const Amg::Vector3D& momentum,
154 double charge,
155 std::optional<AmgSymMatrix(5)> cov = std::nullopt) const override final;
156
160 double l1,
161 double l2,
162 double phi,
163 double theta,
164 double oop,
165 std::optional<AmgSymMatrix(5)> cov = std::nullopt) const override final;
166
170 const Amg::Vector3D& position,
171 const Amg::Vector3D& momentum,
172 double charge = 0.,
173 std::optional<AmgSymMatrix(5)> cov = std::nullopt) const override final;
174
176 template<int DIM, class T>
178 double l1,
179 double l2,
180 double phi,
181 double theta,
182 double qop,
183 std::optional<AmgSymMatrix(DIM)> cov = std::nullopt) const;
184
186 template<int DIM, class T>
188 const Amg::Vector3D& position,
189 const Amg::Vector3D& momentum,
190 double charge,
191 std::optional<AmgSymMatrix(DIM)> cov = std::nullopt) const;
192
195 virtual const SurfaceBounds& bounds() const override final;
196
198 virtual bool insideBounds(const Amg::Vector2D& locpos,
199 double tol1 = 0.,
200 double tol2 = 0.) const override;
201
202 virtual bool insideBoundsCheck(
203 const Amg::Vector2D& locpos,
204 const BoundaryCheck& bchk) const override final;
205
209 virtual bool isOnSurface(const Amg::Vector3D& glopo,
210 const BoundaryCheck& bchk = true,
211 double tol1 = 0.,
212 double tol2 = 0.) const override final;
213
216 virtual void localToGlobal(const Amg::Vector2D& locp,
217 const Amg::Vector3D& mom,
218 Amg::Vector3D& glob) const override final;
219
223 virtual bool globalToLocal(const Amg::Vector3D& glob,
224 const Amg::Vector3D& mom,
225 Amg::Vector2D& loc) const override final;
226
230 Amg::Vector3D& globdir) const;
231
235 Trk::LocalDirection& locdir) const;
236
259 const Amg::Vector3D& pos,
260 const Amg::Vector3D& dir,
261 bool forceDir,
262 Trk::BoundaryCheck bchk) const override final;
263
266 const Amg::Vector3D& pos,
267 const Amg::Vector3D& dir) const override final;
268
271 const Amg::Vector3D& pos,
272 const Amg::Vector3D& dir,
273 bool Bound) const override final;
274
276 virtual std::string name() const override;
277
278protected:
279 template<class SURFACE, class BOUNDS_CNV>
280 friend class ::BoundSurfaceCnv_p1;
281 template<class SURFACE, class BOUNDS_CNV>
282 friend class ::BoundSurfaceCnv_p2;
283
287};
288
289} // end of namespace
290
291#include "TrkSurfaces/PlaneSurface.icc"
292#endif // TRKSURFACES_PLANESURFACE_H
double charge(const T &p)
Definition AtlasPID.h:997
#define protected
Eigen::Matrix< double, 3, 1 > Vector3D
PlaneSurface()
Default Constructor - needed for persistency.
Bounds for a annulus-like, planar Surface.
The BoundaryCheck class allows to steer the way surface boundaries are used for inside/outside checks...
simple class that constructs the curvilinear vectors curvU and curvV from a given momentum direction ...
Bounds for a double trapezoidal ("diamond"), planar Surface.
Access to distance solutions.
Class to describe the bounds for a planar EllipseSurface, i.e.
represents the three-dimensional global direction with respect to a planar surface frame.
Bounds object for a boundless surface (...)
Definition NoBounds.h:30
Dummy class used to allow special convertors to be called for surfaces owned by a detector element.
Definition ParametersT.h:49
void localToGlobalDirection(const Trk::LocalDirection &locdir, Amg::Vector3D &globdir) const
This method transforms a local direction wrt the plane to a global direction.
virtual bool insideBounds(const Amg::Vector2D &locpos, double tol1=0., double tol2=0.) const override
This method calls the inside() method of the Bounds.
PlaneSurface(const PlaneSurface &psf)=default
Copy Constructor.
PlaneSurface & operator=(const PlaneSurface &psf)=default
Assignment operator.
static constexpr SurfaceType staticType
The surface type static constexpr.
bool operator==(const PlaneSurface &cf) const
virtual DistanceSolution straightLineDistanceEstimate(const Amg::Vector3D &pos, const Amg::Vector3D &dir) const override final
fast straight line distance evaluation to Surface
virtual NeutralTrackParametersUniquePtr createUniqueNeutralParameters(double l1, double l2, double phi, double theta, double oop, std::optional< AmgSymMatrix(5)> cov=std::nullopt) const override final
Use the Surface as a ParametersBase constructor, from local parameters - neutral.
std::unique_ptr< ParametersT< DIM, T, PlaneSurface > > createUniqueParameters(double l1, double l2, double phi, double theta, double qop, std::optional< AmgSymMatrix(DIM)> cov=std::nullopt) const
Use the Surface as a ParametersBase constructor, from local parameters.
virtual bool isOnSurface(const Amg::Vector3D &glopo, const BoundaryCheck &bchk=true, double tol1=0., double tol2=0.) const override final
This method returns true if the GlobalPosition is on the Surface for both, within or without check of...
virtual Surface::ChargedTrackParametersUniquePtr createUniqueTrackParameters(double l1, double l2, double phi, double theta, double qop, std::optional< AmgSymMatrix(5)> cov=std::nullopt) const override final
Use the Surface as a ParametersBase constructor, from local parameters - charged.
virtual constexpr SurfaceType type() const override final
Return the surface type.
virtual void localToGlobal(const Amg::Vector2D &locp, const Amg::Vector3D &mom, Amg::Vector3D &glob) const override final
Specified for PlaneSurface: LocalToGlobal method without dynamic memory allocation.
static const NoBounds s_boundless
virtual bool insideBoundsCheck(const Amg::Vector2D &locpos, const BoundaryCheck &bchk) const override final
friend class ::BoundSurfaceCnv_p1
< data members
PlaneSurface()
Default Constructor - needed for persistency.
PlaneSurface & operator=(PlaneSurface &&psf) noexcept=default
Move assignment operator.
virtual ~PlaneSurface()=default
Destructor.
virtual Intersection straightLineIntersection(const Amg::Vector3D &pos, const Amg::Vector3D &dir, bool forceDir, Trk::BoundaryCheck bchk) const override final
fast straight line intersection schema - standard: provides closest intersection and (signed) path le...
virtual bool globalToLocal(const Amg::Vector3D &glob, const Amg::Vector3D &mom, Amg::Vector2D &loc) const override final
Specified for PlaneSurface: GlobalToLocal method without dynamic memory allocation - boolean checks i...
virtual PlaneSurface * clone() const override
Virtual constructor.
virtual bool operator==(const Surface &sf) const override
Equality operator.
void globalToLocalDirection(const Amg::Vector3D &glodir, Trk::LocalDirection &locdir) const
This method transforms the global direction to a local direction wrt the plane.
friend class ::BoundSurfaceCnv_p2
PlaneSurface(PlaneSurface &&psf) noexcept=default
Move Constructor.
std::shared_ptr< const SurfaceBounds > m_bounds
bounds (shared)
virtual const SurfaceBounds & bounds() const override final
This method returns the bounds by reference, static NoBounds in case of no boundaries.
virtual std::string name() const override
Return properly formatted class name for screen output.
Bounds for a rectangular, planar surface.
Bounds for a rotated trapezoidal, planar Surface.
Abstract base class for surface bounds to be specified.
std::unique_ptr< ParametersBase< 5, Trk::Charged > > ChargedTrackParametersUniquePtr
Unique ptr types.
Surface()
Default Constructor for inheriting classes.
std::unique_ptr< ParametersBase< 5, Trk::Neutral > > NeutralTrackParametersUniquePtr
Bounds for a trapezoidal, planar Surface.
Bounds for a triangular, planar surface.
This is the base class for all tracking detector elements with read-out relevant information.
STL class.
Definition of ATLAS Math & Geometry primitives (Amg)
Eigen::Affine3d Transform3D
Eigen::Matrix< double, 3, 1 > Vector3D
Ensure that the ATLAS eigen extensions are properly loaded.
SurfaceType
This enumerator simplifies the persistency & calculations,.
@ theta
Definition ParamDefs.h:66
@ phi
Definition ParamDefs.h:75
AmgSymMatrix(5) &GXFTrackState
STL namespace.