ATLAS Offline Software
Loading...
Searching...
No Matches
Trk::CylinderVolumeBounds Class Referencefinal

Bounds for a cylindrical Volume, the decomposeToSurfaces method creates a vector of up to 6 surfaces: More...

#include <CylinderVolumeBounds.h>

Inheritance diagram for Trk::CylinderVolumeBounds:
Collaboration diagram for Trk::CylinderVolumeBounds:

Public Member Functions

 CylinderVolumeBounds ()
 Default Constructor.
 CylinderVolumeBounds (double radius, double halez)
 Constructor - full cylinder.
 CylinderVolumeBounds (double rinner, double router, double halez)
 Constructor - extruded cylinder.
 CylinderVolumeBounds (double rinner, double router, double halfPhiSector, double halez)
 Constructor - extruded cylinder segment.
 CylinderVolumeBounds (const CylinderVolumeBounds &cylbo)
 Copy Constructor.
virtual ~CylinderVolumeBounds ()
 Destructor.
CylinderVolumeBoundsoperator= (const CylinderVolumeBounds &cylbo)
 Assignment operator.
CylinderVolumeBoundsclone () const override
 Virtual constructor.
bool inside (const Amg::Vector3D &, double tol=0.) const override final
 This method checks if position in the 3D volume frame is inside the cylinder.
virtual std::vector< std::unique_ptr< Trk::Surface > > decomposeToSurfaces (const Amg::Transform3D &transform) override final
 Method to decompose the Bounds into boundarySurfaces.
ObjectAccessor boundarySurfaceAccessor (const Amg::Vector3D &gp, const Amg::Vector3D &dir, bool forceInside=false) const override final
 Provide accessor for BoundarySurfaces.
double innerRadius () const
 This method returns the inner radius.
double outerRadius () const
 This method returns the outer radius.
double mediumRadius () const
 This method returns the medium radius.
double deltaRadius () const
 This method returns the delta radius.
double halfPhiSector () const
 This method returns the halfPhiSector angle.
double halflengthZ () const
 This method returns the halflengthZ.
MsgStream & dump (MsgStream &sl) const override final
 Output Method for MsgStream.
std::ostream & dump (std::ostream &sl) const override final
 Output Method for std::ostream.

Private Member Functions

std::shared_ptr< CylinderBoundsinnerCylinderBounds () const
 This method returns the associated CylinderBounds of the inner CylinderSurfaces.
std::shared_ptr< CylinderBoundsouterCylinderBounds () const
 This method returns the associated CylinderBounds of the outer CylinderSurfaces.
std::shared_ptr< DiscBoundsbottomDiscBounds () const
 This method returns the associated DiscBounds for the bottom/top DiscSurface.
std::shared_ptr< DiscBoundstopDiscBounds () const
 This method returns the associated DiscBounds for the bottom/top DiscSurface.
std::shared_ptr< RectangleBoundssectorPlaneBounds () const
 This method returns the associated PlaneBounds limiting a sectoral CylinderVolume.
void createBoundarySurfaceAccessors ()
 Private method to construct the accessors.

Private Attributes

double m_innerRadius
double m_outerRadius
double m_halfPhiSector
double m_halfZ
CylinderVolumeBoundaryAccessors m_boundaryAccessors
 Accessors for Boundary surface access - static is not possible due to mismatched delete() / free () with TrkMagFieldUtils.

Static Private Attributes

static const double s_numericalStable
 numerical stability

Detailed Description

Bounds for a cylindrical Volume, the decomposeToSurfaces method creates a vector of up to 6 surfaces:

case A) 3 Surfaces (full cylindrical tube): BoundarySurfaceFace [index]:

case B) 4 Surfaces (tube with inner and outer radius): BoundarySurfaceFace [index]:

case C) 6 Surfaces (sectoral tube with inner and outer radius): BoundarySurfaceFace [index]:

  • negativeFaceXY [0] : Trk::DiscSurface with \( r_{inner}>0 \) and \( \phi < \pi \), parallel to \( xy \) plane at negative \( z \)
  • positiveFaceXY [1] : Trk::DiscSurface with \( r_{inner}>0 \) and \( \phi < \pi \), parallel to \( xy \) plane at positive \( z \)
  • tubeSectorOuterCover [2] : Trk::CylinderSurface with \( r = r_{outer} \)
  • tubeSectorInnerCover [3] : Trk::CylinderSurface with \( r = r_{inner} \)
  • tubeSectorNegativePhi [4] : Rectangular Trk::PlaneSurface attached to [0] and [1] at negative \( \phi \)
  • tubeSectorNegativePhi [5] : Rectangular Trk::PlaneSurface attached to [0] and [1] at positive \( \phi \)
Author
Andre.nosp@m.as.S.nosp@m.alzbu.nosp@m.rger.nosp@m.@cern.nosp@m..ch
Christos Anastopoulos (Athena MT modifications)

Definition at line 70 of file CylinderVolumeBounds.h.

Constructor & Destructor Documentation

◆ CylinderVolumeBounds() [1/5]

Trk::CylinderVolumeBounds::CylinderVolumeBounds ( )

Default Constructor.

Definition at line 30 of file CylinderVolumeBounds.cxx.

31 : VolumeBounds()
32 , m_innerRadius(0.)
33 , m_outerRadius(0.)
34 , m_halfPhiSector(0.)
35 , m_halfZ(0.)
37{}
CylinderVolumeBoundaryAccessors m_boundaryAccessors
Accessors for Boundary surface access - static is not possible due to mismatched delete() / free () w...
VolumeBounds()
Default Constructor.

◆ CylinderVolumeBounds() [2/5]

Trk::CylinderVolumeBounds::CylinderVolumeBounds ( double radius,
double halez )

Constructor - full cylinder.

Definition at line 39 of file CylinderVolumeBounds.cxx.

40 : VolumeBounds()
41 , m_innerRadius(0.)
42 , m_outerRadius(std::abs(radius))
44 , m_halfZ(std::abs(halez))
46{}
#define M_PI

◆ CylinderVolumeBounds() [3/5]

Trk::CylinderVolumeBounds::CylinderVolumeBounds ( double rinner,
double router,
double halez )

Constructor - extruded cylinder.

Definition at line 48 of file CylinderVolumeBounds.cxx.

52 : VolumeBounds()
53 , m_innerRadius(std::abs(rinner))
54 , m_outerRadius(std::abs(router))
56 , m_halfZ(std::abs(halez))
58{}

◆ CylinderVolumeBounds() [4/5]

Trk::CylinderVolumeBounds::CylinderVolumeBounds ( double rinner,
double router,
double halfPhiSector,
double halez )

Constructor - extruded cylinder segment.

Definition at line 60 of file CylinderVolumeBounds.cxx.

65 : VolumeBounds()
66 , m_innerRadius(std::abs(rinner))
67 , m_outerRadius(std::abs(router))
68 , m_halfPhiSector(haphi)
69 , m_halfZ(std::abs(halez))
71{}

◆ CylinderVolumeBounds() [5/5]

Trk::CylinderVolumeBounds::CylinderVolumeBounds ( const CylinderVolumeBounds & cylbo)

Copy Constructor.

Definition at line 73 of file CylinderVolumeBounds.cxx.

75 : VolumeBounds()
76 , m_innerRadius(cylbo.m_innerRadius)
77 , m_outerRadius(cylbo.m_outerRadius)
78 , m_halfPhiSector(cylbo.m_halfPhiSector)
79 , m_halfZ(cylbo.m_halfZ)
81{}

◆ ~CylinderVolumeBounds()

Trk::CylinderVolumeBounds::~CylinderVolumeBounds ( )
virtualdefault

Destructor.

Member Function Documentation

◆ bottomDiscBounds()

std::shared_ptr< Trk::DiscBounds > Trk::CylinderVolumeBounds::bottomDiscBounds ( ) const
private

This method returns the associated DiscBounds for the bottom/top DiscSurface.

Definition at line 364 of file CylinderVolumeBounds.cxx.

365{
366 return std::make_shared<Trk::DiscBounds>(m_innerRadius, m_outerRadius, m_halfPhiSector);
367}

◆ boundarySurfaceAccessor()

Trk::ObjectAccessor Trk::CylinderVolumeBounds::boundarySurfaceAccessor ( const Amg::Vector3D & gp,
const Amg::Vector3D & dir,
bool forceInside = false ) const
finaloverridevirtual

Provide accessor for BoundarySurfaces.

Implements Trk::VolumeBounds.

Definition at line 163 of file CylinderVolumeBounds.cxx.

167{
168 // the tube case - most likely
169 if (
171 std::abs(m_halfPhiSector - M_PI) < s_numericalStable) {
172 // prepare the data
173 double posZ = gp.z();
174 double posR = gp.perp();
175 // difference Vector
176 Amg::Vector3D diff(gp + dir);
177 // differences
178 double deltaZ = diff.z() - posZ;
179 double deltaR = diff.perp() - posR;
180
181 // the isOnSurface cases + switchers (that change with the cases)
182 bool isOnFace = false;
183 bool intersectionIndicator = (deltaR > 0.);
184 bool choiceIndicator = false;
185
186 // on surface or look-a-likes
187 // on zMin or slightly outside
188 if (
189 std::abs(posZ + m_halfZ) < s_numericalStable ||
190 (posZ < (-m_halfZ) && deltaZ > 0.)) {
191 isOnFace = true;
192 } else if (
193 std::abs(posZ - m_halfZ) < s_numericalStable ||
194 (posZ > m_halfZ && deltaZ < 0.)) {
195 // on zMax or slighly outside
196 isOnFace = true;
197 } else if (
198 std::abs(posR - m_innerRadius) < s_numericalStable ||
199 (posR < m_innerRadius && deltaR > 0.)) {
200 // on innerRadius or slightly outside
201 isOnFace = true;
202 choiceIndicator =
203 (deltaZ > 0.); // the side choice indicator <========== WRONG ==========
204 } else if (
205 std::abs(posR - m_outerRadius) < s_numericalStable ||
206 (posR > m_outerRadius && deltaR < 0.)) {
207 // on outRadius or slightly outside
208 isOnFace = true;
209 choiceIndicator =
210 (deltaZ > 0.); // the side choice indicator <========== WRONG ==========
211 }
212
213 // the onSurface case
214 // =================================================================================
215 if (isOnFace) {
216 if (intersectionIndicator) {
217 // intersect the Rmax / with floatint point exception check
218 Trk::CylinderIntersector intersectRmax(
219 posR, posZ, deltaR != 0 ? deltaZ / deltaR : 0, m_outerRadius);
220 double zOfIntersect = intersectRmax.yOfX;
221 // now check if the intersect is inside m_halfZ
222 if (std::abs(zOfIntersect) <= m_halfZ)
223 return {(choiceIndicator || zOfIntersect > 0.)
224 ? m_boundaryAccessors.tubeAccessor(
226 : m_boundaryAccessors.tubeAccessor(
228 // if the intersect is outside
229 return {
230 (choiceIndicator || zOfIntersect > 0.)
232 : m_boundaryAccessors.tubeAccessor(
234 }
235 // intersect the Rmin
236 Trk::CylinderIntersector intersectRmin(
237 posR, posZ, deltaR != 0 ? deltaZ / deltaR : 0, m_innerRadius);
238 double zOfIntersect = intersectRmin.yOfX;
239 // now check if the intersect is inside m_halfZ
240 if (std::abs(zOfIntersect) <= m_halfZ)
241 return {
242 (choiceIndicator || zOfIntersect > 0.)
244 : m_boundaryAccessors.tubeAccessor(
246 // if the intersect is outside
247 return {
248 (choiceIndicator || zOfIntersect > 0.)
251 }
252 // =================================================================================================
253
254 // ======================= the inside/outside part remains
255 // =========================================
256 // (a) outside cases
257 if (posR < m_innerRadius && deltaR < 0.)
259 if (posR > m_outerRadius && deltaR > 0.)
261 if (posZ < -m_halfZ && deltaZ < 0.)
263 if (posZ > m_halfZ && deltaZ > 0.)
265 // (b) inside cases
266 // the increase R case
267 if (deltaR > 0.) {
268 // solve the linear equation for the outer Radius - with floating point
269 // exception check
270 Trk::CylinderIntersector intersectRmax(
271 posR, posZ, deltaR != 0 ? deltaZ / deltaR : 0, m_outerRadius);
272 double zOfIntersect = intersectRmax.yOfX;
273
274 if (std::abs(zOfIntersect) <= m_halfZ && zOfIntersect > 0.)
276 if (std::abs(zOfIntersect) <= m_halfZ && zOfIntersect < 0.)
278 if (std::abs(zOfIntersect) > m_halfZ && zOfIntersect < 0.)
280 if (std::abs(zOfIntersect) > m_halfZ && zOfIntersect > 0.)
282
283 } else {
284 // solve the linear equation for the inner Radius with floating point
285 // exceptin check
286 Trk::CylinderIntersector intersectRmin(
287 posR, posZ, deltaR != 0 ? deltaZ / deltaR : 0, m_innerRadius);
288 double zOfIntersect = intersectRmin.yOfX;
289
290 if (std::abs(zOfIntersect) <= m_halfZ && zOfIntersect > 0.)
292 if (std::abs(zOfIntersect) <= m_halfZ && zOfIntersect < 0.)
294 if (std::abs(zOfIntersect) > m_halfZ && zOfIntersect > 0.)
297 }
298 }
299 // the cylinder case
300 if (
302 std::abs(m_halfPhiSector - M_PI) < s_numericalStable) {
303
304 // ----------------------------- z / r information
305 double posZ = gp.z();
306 double posR = gp.perp();
307
308 // is on R case
309 bool isOnCylinder = (std::abs(posR - m_outerRadius) < s_numericalStable);
310
311 // difference Vector
312 Amg::Vector3D diff(gp + dir);
313 // differences
314 double deltaZ = diff.z() - posZ;
315 double deltaR = diff.perp() - posR;
316 // take the sign of the intersect
317 int radiusSign = deltaR > 0. ? 1 : -1;
318
319 // solve the linear equation for the cylinder Radius / with floating point
320 // exceptin check
321 Trk::CylinderIntersector intersectR(
322 posR,
323 posZ,
324 deltaR != 0 ? deltaZ / deltaR : 0,
325 radiusSign * m_outerRadius);
326
327 // the intersection point - a little trick for the OnFace cases
328 double zOfIntersect = intersectR.yOfX;
329 // check if intersection exists
330 bool intersectsCylinder =
331 !isOnCylinder ? (zOfIntersect * zOfIntersect <= m_halfZ * m_halfZ)
332 : false;
333
334 // return the cases for going through the cylinder
335 if (intersectsCylinder && zOfIntersect > 0.)
336 return {m_boundaryAccessors.cylinderAccessor(Trk::CylinderZincrease)};
337 if (intersectsCylinder && zOfIntersect <= 0.)
338 return {m_boundaryAccessors.cylinderAccessor(Trk::CylinderZdecrease)};
339 if (deltaZ > 0.)
340 return {m_boundaryAccessors.cylinderAccessor(Trk::CylinderPositiveFace)};
341 return {m_boundaryAccessors.cylinderAccessor(Trk::CylinderNegativeFace)};
342 }
343 // the sectoral cylinder case
344 if (m_innerRadius != 0. && std::abs(m_halfPhiSector - M_PI) > 10e-3)
345 return {m_boundaryAccessors.sectoralCylinderAccessor(
347 // it remains the sectoral tube case
348 return {m_boundaryAccessors.sectoralTubeAccessor(Trk::StandardSectoralTube)};
349}
Scalar deltaR(const MatrixBase< Derived > &vec) const
void diff(const Jet &rJet1, const Jet &rJet2, std::map< std::string, double > varDiff)
Difference between jets - Non-Class function required by trigger.
Definition Jet.cxx:631
static const double s_numericalStable
numerical stability
Eigen::Matrix< double, 3, 1 > Vector3D
@ CylinderZincrease
Cylinder hit, then pos.
@ CylinderZdecrease
Cylinder hit, the neg.
@ TubeOutsideZminZdecrease
Accessor type [ 0,3,2,1 ] - inverse case.
@ TubeRdecreaseZdecrease
Accessor type [ 3,0,1,2 ].
@ TubeOutsideZmaxZincrease
Accessor type [ 1,3,2,0 ] - inverse case.
@ TubeZincreaseRincrease
Accessor type [ 1,2,0,3 ].
@ TubeRincreaseZdecrease
Accessor type [ 2,0,1,3 ].
@ TubeOutsideRminRdecrease
Accessor type [ 3,1,0,2] - inverse case.
@ TubeZdecreaseRdecrease
Accessor type [ 0,3,1,2 ].
@ TubeRincreaseZincrease
Accessor type [ 2,1,0,3 ].
@ TubeZincreaseRdecrease
Accessor type [ 1,3,0,2 ].
@ TubeRdecreaseZincrease
Accessor type [ 3,1,0,2 ].
@ TubeOutsideRmaxRincrease
Accessor type [ 2,1,0,3 ] - inverse case.
@ TubeZdecreaseRincrease
Accessor type [ 0,2,1,3 ].

◆ clone()

CylinderVolumeBounds * Trk::CylinderVolumeBounds::clone ( ) const
inlineoverridevirtual

Virtual constructor.

Implements Trk::VolumeBounds.

Definition at line 171 of file CylinderVolumeBounds.h.

171 {
172 return new CylinderVolumeBounds(*this);
173}
CylinderVolumeBounds()
Default Constructor.

◆ createBoundarySurfaceAccessors()

void Trk::CylinderVolumeBounds::createBoundarySurfaceAccessors ( )
private

Private method to construct the accessors.

◆ decomposeToSurfaces()

std::vector< std::unique_ptr< Trk::Surface > > Trk::CylinderVolumeBounds::decomposeToSurfaces ( const Amg::Transform3D & transform)
finaloverridevirtual

Method to decompose the Bounds into boundarySurfaces.

Implements Trk::VolumeBounds.

Definition at line 98 of file CylinderVolumeBounds.cxx.

100{
101 auto retsf = std::vector<std::unique_ptr<Trk::Surface>>();
102 // memory optimisation --- reserve the maximum
103 retsf.reserve(6);
104
105 // check if the transform is approximatively the identity
106 bool isConcentric = transform.isApprox(Amg::Transform3D::Identity());
107
108 Amg::RotationMatrix3D discRot = transform.linear();
109 Amg::Vector3D cylCenter(transform.translation());
110
111 // bottom Disc (negative z)
112 Amg::RotationMatrix3D bottomDiscRot;
113 bottomDiscRot.col(0) = discRot.col(1);
114 bottomDiscRot.col(1) = discRot.col(0);
115 bottomDiscRot.col(2) = -discRot.col(2);
116 retsf.push_back(std::make_unique<Trk::DiscSurface>(
118 transform * Amg::AngleAxis3D(M_PI, Amg::Vector3D(1., 0., 0.)) *
121 // top Disc (positive z)
122 retsf.push_back(std::make_unique<Trk::DiscSurface>(
124 topDiscBounds()));
125 // outer Cylinder
126 if (!isConcentric)
127 retsf.push_back(std::make_unique<Trk::CylinderSurface>(
129 else
130 retsf.push_back(std::make_unique<Trk::CylinderSurface>(outerCylinderBounds()));
131
132 // innermost Cylinder
134 if (!isConcentric)
135 retsf.push_back(std::make_unique<Trk::CylinderSurface>(
137 else
138 retsf.push_back(std::make_unique<Trk::CylinderSurface>(innerCylinderBounds()));
139 }
140
141 if (std::abs(halfPhiSector() - M_PI) > s_numericalStable) {
142 // sectorPlane 1 (negative phi)
143 retsf.push_back(std::make_unique<Trk::PlaneSurface>(
145 transform *
148 Amg::AngleAxis3D(M_PI / 2, Amg::Vector3D(1., 0., 0.))),
150 // sectorPlane 2 (positive phi)
151 retsf.push_back(std::make_unique<Trk::PlaneSurface>(
153 transform *
156 Amg::AngleAxis3D(-M_PI / 2, Amg::Vector3D(1., 0., 0.))),
158 }
159 return retsf;
160}
double mediumRadius() const
This method returns the medium radius.
double innerRadius() const
This method returns the inner radius.
double halflengthZ() const
This method returns the halflengthZ.
std::shared_ptr< CylinderBounds > outerCylinderBounds() const
This method returns the associated CylinderBounds of the outer CylinderSurfaces.
std::shared_ptr< DiscBounds > bottomDiscBounds() const
This method returns the associated DiscBounds for the bottom/top DiscSurface.
std::shared_ptr< DiscBounds > topDiscBounds() const
This method returns the associated DiscBounds for the bottom/top DiscSurface.
double halfPhiSector() const
This method returns the halfPhiSector angle.
std::shared_ptr< CylinderBounds > innerCylinderBounds() const
This method returns the associated CylinderBounds of the inner CylinderSurfaces.
std::shared_ptr< RectangleBounds > sectorPlaneBounds() const
This method returns the associated PlaneBounds limiting a sectoral CylinderVolume.
Eigen::AngleAxisd AngleAxis3D
Eigen::Matrix< double, 3, 3 > RotationMatrix3D
Eigen::Affine3d Transform3D
Amg::Vector3D transform(Amg::Vector3D &v, Amg::Transform3D &tr)
Transform a point from a Trasformation3D.
Eigen::Translation< double, 3 > Translation3D

◆ deltaRadius()

double Trk::CylinderVolumeBounds::deltaRadius ( ) const
inline

This method returns the delta radius.

Definition at line 199 of file CylinderVolumeBounds.h.

199 {
200 return (m_outerRadius - m_innerRadius);
201}

◆ dump() [1/2]

MsgStream & Trk::CylinderVolumeBounds::dump ( MsgStream & sl) const
finaloverridevirtual

Output Method for MsgStream.

Implements Trk::VolumeBounds.

Definition at line 379 of file CylinderVolumeBounds.cxx.

380{
381 std::stringstream temp_sl;
382 temp_sl << std::setiosflags(std::ios::fixed);
383 temp_sl << std::setprecision(2);
384 temp_sl << "Trk::CylinderVolumeBounds: (rMin, rMax, halfPhi, halfZ) = ";
385 temp_sl << m_innerRadius << ", " << m_outerRadius << ", " << m_halfPhiSector
386 << ", " << m_halfZ;
387 sl << temp_sl.str();
388 return sl;
389}

◆ dump() [2/2]

std::ostream & Trk::CylinderVolumeBounds::dump ( std::ostream & sl) const
finaloverridevirtual

Output Method for std::ostream.

Implements Trk::VolumeBounds.

Definition at line 392 of file CylinderVolumeBounds.cxx.

393{
394 std::stringstream temp_sl;
395 temp_sl << std::setiosflags(std::ios::fixed);
396 temp_sl << std::setprecision(2);
397 temp_sl << "Trk::CylinderVolumeBounds: (rMin, rMax, halfPhi, halfZ) = ";
398 temp_sl << m_innerRadius << ", " << m_outerRadius << ", " << m_halfPhiSector
399 << ", " << m_halfZ;
400 sl << temp_sl.str();
401 return sl;
402}

◆ halflengthZ()

double Trk::CylinderVolumeBounds::halflengthZ ( ) const
inline

This method returns the halflengthZ.

Definition at line 207 of file CylinderVolumeBounds.h.

207{ return m_halfZ; }

◆ halfPhiSector()

double Trk::CylinderVolumeBounds::halfPhiSector ( ) const
inline

This method returns the halfPhiSector angle.

Definition at line 203 of file CylinderVolumeBounds.h.

203 {
204 return m_halfPhiSector;
205}

◆ innerCylinderBounds()

std::shared_ptr< Trk::CylinderBounds > Trk::CylinderVolumeBounds::innerCylinderBounds ( ) const
private

This method returns the associated CylinderBounds of the inner CylinderSurfaces.

Definition at line 352 of file CylinderVolumeBounds.cxx.

353{
354 return std::make_shared<Trk::CylinderBounds>(m_innerRadius, m_halfPhiSector, m_halfZ);
355}

◆ innerRadius()

double Trk::CylinderVolumeBounds::innerRadius ( ) const
inline

This method returns the inner radius.

Definition at line 187 of file CylinderVolumeBounds.h.

187 {
188 return m_innerRadius;
189}

◆ inside()

bool Trk::CylinderVolumeBounds::inside ( const Amg::Vector3D & pos,
double tol = 0. ) const
inlinefinaloverridevirtual

This method checks if position in the 3D volume frame is inside the cylinder.

Implements Trk::VolumeBounds.

Definition at line 175 of file CylinderVolumeBounds.h.

176 {
177 double ros = pos.perp();
178 // bool insidePhi = fabs(pos.phi()) <= m_halfPhiSector + tol;
179 bool insidePhi = cos(pos.phi()) >= cos(m_halfPhiSector) - tol;
180 bool insideR =
181 insidePhi ? ((ros >= m_innerRadius - tol) && (ros <= m_outerRadius + tol))
182 : false;
183 bool insideZ = insideR ? (fabs(pos.z()) <= m_halfZ + tol) : false;
184 return (insideZ && insideR && insidePhi);
185}

◆ mediumRadius()

double Trk::CylinderVolumeBounds::mediumRadius ( ) const
inline

This method returns the medium radius.

Definition at line 195 of file CylinderVolumeBounds.h.

195 {
196 return 0.5 * (m_innerRadius + m_outerRadius);
197}

◆ operator=()

Trk::CylinderVolumeBounds & Trk::CylinderVolumeBounds::operator= ( const CylinderVolumeBounds & cylbo)

Assignment operator.

Definition at line 86 of file CylinderVolumeBounds.cxx.

87{
88 if (this != &cylbo) {
89 m_innerRadius = cylbo.m_innerRadius;
90 m_outerRadius = cylbo.m_outerRadius;
91 m_halfPhiSector = cylbo.m_halfPhiSector;
92 m_halfZ = cylbo.m_halfZ;
93 }
94 return *this;
95}

◆ outerCylinderBounds()

std::shared_ptr< Trk::CylinderBounds > Trk::CylinderVolumeBounds::outerCylinderBounds ( ) const
private

This method returns the associated CylinderBounds of the outer CylinderSurfaces.

Definition at line 358 of file CylinderVolumeBounds.cxx.

359{
360 return std::make_shared<Trk::CylinderBounds>(m_outerRadius, m_halfPhiSector, m_halfZ);
361}

◆ outerRadius()

double Trk::CylinderVolumeBounds::outerRadius ( ) const
inline

This method returns the outer radius.

Definition at line 191 of file CylinderVolumeBounds.h.

191 {
192 return m_outerRadius;
193}

◆ sectorPlaneBounds()

std::shared_ptr< Trk::RectangleBounds > Trk::CylinderVolumeBounds::sectorPlaneBounds ( ) const
private

This method returns the associated PlaneBounds limiting a sectoral CylinderVolume.

Definition at line 370 of file CylinderVolumeBounds.cxx.

371{
372 return std::make_shared<Trk::RectangleBounds>(
374}

◆ topDiscBounds()

std::shared_ptr< DiscBounds > Trk::CylinderVolumeBounds::topDiscBounds ( ) const
inlineprivate

This method returns the associated DiscBounds for the bottom/top DiscSurface.

Definition at line 209 of file CylinderVolumeBounds.h.

209 {
210 return this->bottomDiscBounds();
211}

Member Data Documentation

◆ m_boundaryAccessors

CylinderVolumeBoundaryAccessors Trk::CylinderVolumeBounds::m_boundaryAccessors
private

Accessors for Boundary surface access - static is not possible due to mismatched delete() / free () with TrkMagFieldUtils.

Definition at line 165 of file CylinderVolumeBounds.h.

◆ m_halfPhiSector

double Trk::CylinderVolumeBounds::m_halfPhiSector
private

Definition at line 159 of file CylinderVolumeBounds.h.

◆ m_halfZ

double Trk::CylinderVolumeBounds::m_halfZ
private

Definition at line 160 of file CylinderVolumeBounds.h.

◆ m_innerRadius

double Trk::CylinderVolumeBounds::m_innerRadius
private

Definition at line 157 of file CylinderVolumeBounds.h.

◆ m_outerRadius

double Trk::CylinderVolumeBounds::m_outerRadius
private

Definition at line 158 of file CylinderVolumeBounds.h.

◆ s_numericalStable

const double Trk::CylinderVolumeBounds::s_numericalStable
staticprivate
Initial value:
=
10e-2 * Gaudi::Units::mm

numerical stability

Definition at line 168 of file CylinderVolumeBounds.h.


The documentation for this class was generated from the following files: