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

Class to describe the bounds for a planar DiscSurface. More...

#include <DiscTrapezoidalBounds.h>

Inheritance diagram for Trk::DiscTrapezoidalBounds:
Collaboration diagram for Trk::DiscTrapezoidalBounds:

Public Types

enum  BoundValues {
  bv_rMin = 0 , bv_rMax = 1 , bv_minHalfX = 2 , bv_maxHalfX = 3 ,
  bv_halfY = 4 , bv_halfPhiSector = 5 , bv_averagePhi = 6 , bv_rCenter = 7 ,
  bv_stereo = 8 , bv_length = 9
}
 enumeration for readability More...
enum  BoundsType {
  Cone = 0 , Cylinder = 1 , Diamond = 2 , Disc = 3 ,
  Ellipse = 5 , Rectangle = 6 , RotatedTrapezoid = 7 , Trapezoid = 8 ,
  Triangle = 9 , DiscTrapezoidal = 10 , Annulus = 11 , Other = 12
}
 This enumerator simplifies the persistency, by saving a dynamic_cast to happen. More...

Public Member Functions

 DiscTrapezoidalBounds ()
 Default Constructor.
 DiscTrapezoidalBounds (const DiscTrapezoidalBounds &disctrbo)
 Copy constructor.
DiscTrapezoidalBoundsoperator= (const DiscTrapezoidalBounds &disctrbo)
 Assignment operator.
 DiscTrapezoidalBounds (DiscTrapezoidalBounds &&disctrbo) noexcept=default
 Move constructor.
DiscTrapezoidalBoundsoperator= (DiscTrapezoidalBounds &&disctrbo) noexcept=default
 Move Assignment operator.
virtual ~DiscTrapezoidalBounds ()=default
 Destructor.
 DiscTrapezoidalBounds (double minhalfx, double maxhalfx, double rMin, double rMax, double avephi, double stereo=0.)
 Constructor for a symmetric Trapezoid giving min X lenght, max X lenght, Rmin and R max.
virtual bool operator== (const SurfaceBounds &sbo) const override
 Equality operator.
virtual DiscTrapezoidalBoundsclone () const override
 Virtual constructor.
virtual SurfaceBounds::BoundsType type () const override final
 Return the type - mainly for persistency.
virtual bool inside (const Amg::Vector2D &locpo, double tol1=0., double tol2=0.) const override final
 This method cheks if the radius given in the LocalPosition is inside [rMin,rMax] if only tol1 is given and additional in the phi sector is tol2 is given.
virtual bool inside (const Amg::Vector2D &locpo, const BoundaryCheck &bchk) const override final
virtual bool insideLoc1 (const Amg::Vector2D &locpo, double tol1=0.) const override final
 This method checks inside bounds in loc1.
virtual bool insideLoc2 (const Amg::Vector2D &locpo, double tol2=0.) const override final
 This method checks inside bounds in loc2.
virtual double minDistance (const Amg::Vector2D &pos) const override final
 Minimal distance to boundary ( > 0 if outside and <=0 if inside)
double rMin () const
 This method returns inner radius.
double rMax () const
 This method returns outer radius.
virtual double r () const override
 This method returns the maximum expansion on the plane (=rMax)
double averagePhi () const
 This method returns the average phi.
double rCenter () const
 This method returns the center radius.
double stereo () const
 This method returns the stereo angle.
double halfPhiSector () const
 This method returns the halfPhiSector which is covered by the disc.
double minHalflengthX () const
 This method returns the minimal halflength in X.
double maxHalflengthX () const
 This method returns the maximal halflength in X.
double halflengthY () const
 This method returns the halflength in Y (this is Rmax -Rmin)
virtual MsgStream & dump (MsgStream &sl) const override
 Output Method for MsgStream.
virtual std::ostream & dump (std::ostream &sl) const override
 Output Method for std::ostream.
virtual bool operator!= (const SurfaceBounds &sb) const
 Non-Equality operator.

Protected Member Functions

void swap (double &b1, double &b2)
 Swap method to be called from DiscBounds or TrapezoidalBounds.
virtual void initCache ()
 virtual initCache method for object persistency

Private Attributes

std::vector< TDD_real_tm_boundValues
 Internal members of the bounds (float/double)

Detailed Description

Class to describe the bounds for a planar DiscSurface.

By providing an argument for hphisec, the bounds can be restricted to a phirange around the center position.

Author
Noemi.nosp@m..Cal.nosp@m.ace@c.nosp@m.ern..nosp@m.ch , Andre.nosp@m.as.S.nosp@m.alzbu.nosp@m.rger.nosp@m.@cern.nosp@m..ch

Definition at line 41 of file DiscTrapezoidalBounds.h.

Member Enumeration Documentation

◆ BoundsType

This enumerator simplifies the persistency, by saving a dynamic_cast to happen.

Other is reserved for the GeometrySurfaces implementation.

Enumerator
Cone 
Cylinder 
Diamond 
Disc 
Ellipse 
Rectangle 
RotatedTrapezoid 
Trapezoid 
Triangle 
DiscTrapezoidal 
Annulus 
Other 

Definition at line 58 of file SurfaceBounds.h.

◆ BoundValues

enumeration for readability

Enumerator
bv_rMin 
bv_rMax 
bv_minHalfX 
bv_maxHalfX 
bv_halfY 
bv_halfPhiSector 
bv_averagePhi 
bv_rCenter 
bv_stereo 
bv_length 

Definition at line 46 of file DiscTrapezoidalBounds.h.

Constructor & Destructor Documentation

◆ DiscTrapezoidalBounds() [1/4]

Trk::DiscTrapezoidalBounds::DiscTrapezoidalBounds ( )

Default Constructor.

Definition at line 17 of file DiscTrapezoidalBounds.cxx.

19{}
std::vector< TDD_real_t > m_boundValues
Internal members of the bounds (float/double)

◆ DiscTrapezoidalBounds() [2/4]

Trk::DiscTrapezoidalBounds::DiscTrapezoidalBounds ( const DiscTrapezoidalBounds & disctrbo)

Copy constructor.

Definition at line 57 of file DiscTrapezoidalBounds.cxx.

58 : Trk::SurfaceBounds()
59 , m_boundValues(disctrbo.m_boundValues)
60{}

◆ DiscTrapezoidalBounds() [3/4]

Trk::DiscTrapezoidalBounds::DiscTrapezoidalBounds ( DiscTrapezoidalBounds && disctrbo)
defaultnoexcept

Move constructor.

◆ ~DiscTrapezoidalBounds()

virtual Trk::DiscTrapezoidalBounds::~DiscTrapezoidalBounds ( )
virtualdefault

Destructor.

◆ DiscTrapezoidalBounds() [4/4]

Trk::DiscTrapezoidalBounds::DiscTrapezoidalBounds ( double minhalfx,
double maxhalfx,
double rMin,
double rMax,
double avephi,
double stereo = 0. )

Constructor for a symmetric Trapezoid giving min X lenght, max X lenght, Rmin and R max.

Definition at line 21 of file DiscTrapezoidalBounds.cxx.

28{
31
36
41
44
45 double hmax =
48
49 double hmin =
52
55}
double stereo() const
This method returns the stereo angle.
void swap(double &b1, double &b2)
Swap method to be called from DiscBounds or TrapezoidalBounds.
double hmax(TH1 *&h)
Definition listroot.cxx:115

Member Function Documentation

◆ averagePhi()

double Trk::DiscTrapezoidalBounds::averagePhi ( ) const

This method returns the average phi.

◆ clone()

virtual DiscTrapezoidalBounds * Trk::DiscTrapezoidalBounds::clone ( ) const
overridevirtual

Virtual constructor.

Implements Trk::SurfaceBounds.

◆ dump() [1/2]

MsgStream & Trk::DiscTrapezoidalBounds::dump ( MsgStream & sl) const
overridevirtual

Output Method for MsgStream.

Implements Trk::SurfaceBounds.

Definition at line 369 of file DiscTrapezoidalBounds.cxx.

370{
371 sl << std::setiosflags(std::ios::fixed);
372 sl << std::setprecision(7);
373 sl << "Trk::DiscTrapezoidalBounds: (innerRadius, outerRadius, hMinX, hMaxX, hlengthY, hPhiSector, averagePhi, "
374 "rCenter, stereo) = ";
375 sl << "(" << this->rMin() << ", " << this->rMax() << ", " << this->minHalflengthX() << ", " << this->maxHalflengthX()
376 << ", " << this->halflengthY() << ", " << this->halfPhiSector() << ", " << this->averagePhi() << ", "
377 << this->rCenter() << ", " << this->stereo() << ")";
378 sl << std::setprecision(-1);
379 return sl;
380}
double rMin() const
This method returns inner radius.
double halflengthY() const
This method returns the halflength in Y (this is Rmax -Rmin)
double averagePhi() const
This method returns the average phi.
double maxHalflengthX() const
This method returns the maximal halflength in X.
double halfPhiSector() const
This method returns the halfPhiSector which is covered by the disc.
double minHalflengthX() const
This method returns the minimal halflength in X.
double rMax() const
This method returns outer radius.
double rCenter() const
This method returns the center radius.

◆ dump() [2/2]

std::ostream & Trk::DiscTrapezoidalBounds::dump ( std::ostream & sl) const
overridevirtual

Output Method for std::ostream.

Implements Trk::SurfaceBounds.

Definition at line 383 of file DiscTrapezoidalBounds.cxx.

384{
385 sl << std::setiosflags(std::ios::fixed);
386 sl << std::setprecision(7);
387 sl << "Trk::DiscTrapezoidalBounds: (innerRadius, outerRadius, hMinX, hMaxX, hlengthY, hPhiSector, averagePhi, "
388 "rCenter, stereo) = ";
389 sl << "(" << this->rMin() << ", " << this->rMax() << ", " << this->minHalflengthX() << ", " << this->maxHalflengthX()
390 << ", " << this->halflengthY() << ", " << this->halfPhiSector() << ", " << this->averagePhi() << ", "
391 << this->rCenter() << ", " << this->stereo() << ")";
392 sl << std::setprecision(-1);
393 return sl;
394}

◆ halflengthY()

double Trk::DiscTrapezoidalBounds::halflengthY ( ) const

This method returns the halflength in Y (this is Rmax -Rmin)

◆ halfPhiSector()

double Trk::DiscTrapezoidalBounds::halfPhiSector ( ) const

This method returns the halfPhiSector which is covered by the disc.

◆ initCache()

virtual void Trk::SurfaceBounds::initCache ( )
inlineprotectedvirtualinherited

virtual initCache method for object persistency

Reimplemented in Trk::ConeBounds, Trk::DiamondBounds, Trk::RotatedDiamondBounds, and Trk::RotatedTrapezoidBounds.

Definition at line 129 of file SurfaceBounds.h.

129{}

◆ inside() [1/2]

bool Trk::DiscTrapezoidalBounds::inside ( const Amg::Vector2D & locpo,
const BoundaryCheck & bchk ) const
finaloverridevirtual

Implements Trk::SurfaceBounds.

Definition at line 125 of file DiscTrapezoidalBounds.cxx.

127{
128 if (bchk.bcType == 0 || bchk.nSigmas == 0 ||
131 locpo, bchk.toleranceLoc1, bchk.toleranceLoc2);
132 double alpha =
134 if (alpha > M_PI)
135 alpha = 2 * M_PI - alpha;
136 // a fast FALSE
137 sincosCache scResult = bchk.FastSinCos(locpo(1, 0));
138 double dx = bchk.nSigmas * sqrt(bchk.lCovariance(0, 0));
139 double dy =
140 bchk.nSigmas * sqrt(scResult.sinC * scResult.sinC * bchk.lCovariance(0, 0) +
141 locpo(0, 0) * locpo(0, 0) * scResult.cosC *
142 scResult.cosC * bchk.lCovariance(1, 1) +
143 2 * scResult.cosC * scResult.sinC * locpo(0, 0) *
144 bchk.lCovariance(0, 1));
145 double max_ell = dx > dy ? dx : dy;
146 if (locpo(0, 0) >
149 cos(alpha) +
150 max_ell))
151 return false;
152 // a fast TRUE
153 double min_ell = dx < dy ? dx : dy;
154 if (locpo(0, 0) <
157 cos(alpha) +
158 min_ell))
159 return true;
160
161 // we are not using the KDOP approach here but rather a highly optimized one
162 class EllipseCollisionTest
163 {
164 private:
165 int m_maxIterations;
166 bool iterate(double x,
167 double y,
168 double c0x,
169 double c0y,
170 double c2x,
171 double c2y,
172 double rr) const
173 {
174 std::vector<double> innerPolygonCoef(m_maxIterations + 1);
175 std::vector<double> outerPolygonCoef(m_maxIterations + 1);
176 /*
177 t2______t4
178 --_ \
179 --_ \ /¨¨¨ ¨¨\
180 t1 = (0, 0) ( t )
181 | \ \__ _ /
182 | \
183 | t3
184 | /
185 | /
186 t0
187 */
188 for (int t = 1; t <= m_maxIterations; t++) {
189 int numNodes = 4 << t;
190 innerPolygonCoef[t] = 0.5 / cos(2.0f*M_PI / numNodes);
191 double c1x = (c0x + c2x) * innerPolygonCoef[t];
192 double c1y = (c0y + c2y) * innerPolygonCoef[t];
193 double tx = x - c1x; // t indicates a translated coordinate
194 double ty = y - c1y;
195 if (tx * tx + ty * ty <= rr) {
196 return true; // collision with t1
197 }
198 double t2x = c2x - c1x;
199 double t2y = c2y - c1y;
200 if (tx * t2x + ty * t2y >= 0 &&
201 tx * t2x + ty * t2y <= t2x * t2x + t2y * t2y &&
202 (ty * t2x - tx * t2y >= 0 ||
203 rr * (t2x * t2x + t2y * t2y) >=
204 (ty * t2x - tx * t2y) * (ty * t2x - tx * t2y))) {
205 return true; // collision with t1---t2
206 }
207 double t0x = c0x - c1x;
208 double t0y = c0y - c1y;
209 if (tx * t0x + ty * t0y >= 0 &&
210 tx * t0x + ty * t0y <= t0x * t0x + t0y * t0y &&
211 (ty * t0x - tx * t0y <= 0 ||
212 rr * (t0x * t0x + t0y * t0y) >=
213 (ty * t0x - tx * t0y) * (ty * t0x - tx * t0y))) {
214 return true; // collision with t1---t0
215 }
216 outerPolygonCoef[t] =
217 0.5 / (std::cos(M_PI / numNodes) * std::cos(M_PI / numNodes));
218 double c3x = (c0x + c1x) * outerPolygonCoef[t];
219 double c3y = (c0y + c1y) * outerPolygonCoef[t];
220 if ((c3x - x) * (c3x - x) + (c3y - y) * (c3y - y) < rr) {
221 c2x = c1x;
222 c2y = c1y;
223 continue; // t3 is inside circle
224 }
225 double c4x = c1x - c3x + c1x;
226 double c4y = c1y - c3y + c1y;
227 if ((c4x - x) * (c4x - x) + (c4y - y) * (c4y - y) < rr) {
228 c0x = c1x;
229 c0y = c1y;
230 continue; // t4 is inside circle
231 }
232 double t3x = c3x - c1x;
233 double t3y = c3y - c1y;
234 if (ty * t3x - tx * t3y <= 0 ||
235 rr * (t3x * t3x + t3y * t3y) >
236 (ty * t3x - tx * t3y) * (ty * t3x - tx * t3y)) {
237 if (tx * t3x + ty * t3y > 0) {
238 if (std::abs(tx * t3x + ty * t3y) <= t3x * t3x + t3y * t3y ||
239 (x - c3x) * (c0x - c3x) + (y - c3y) * (c0y - c3y) >= 0) {
240 c2x = c1x;
241 c2y = c1y;
242 continue; // circle center is inside t0---t1---t3
243 }
244 } else if (-(tx * t3x + ty * t3y) <= t3x * t3x + t3y * t3y ||
245 (x - c4x) * (c2x - c4x) + (y - c4y) * (c2y - c4y) >= 0) {
246 c0x = c1x;
247 c0y = c1y;
248 continue; // circle center is inside t1---t2---t4
249 }
250 }
251 return false; // no collision possible
252 }
253 return false; // out of iterations so it is unsure if there was a
254 // collision. But have to return something.
255 }
256
257 public:
258 // test for collision between an ellipse of horizontal radius w and vertical
259 // radius h at (x0, y0) and a circle of radius r at (x1, y1)
260 bool collide(double x0,
261 double y0,
262 double w,
263 double h,
264 double x1,
265 double y1,
266 double r) const
267 {
268 double x = std::abs(x1 - x0);
269 double y = std::abs(y1 - y0);
270 if (x * x + (h - y) * (h - y) <= r * r ||
271 (w - x) * (w - x) + y * y <= r * r ||
272 x * h + y * w <= w * h // collision with (0, h)
273 || ((x * h + y * w - w * h) * (x * h + y * w - w * h) <=
274 r * r * (w * w + h * h) &&
275 x * w - y * h >= -h * h &&
276 x * w - y * h <= w * w)) { // collision with (0, h)---(w, 0)
277 return true;
278 } else {
279 if ((x - w) * (x - w) + (y - h) * (y - h) <= r * r ||
280 (x <= w && y - r <= h) || (y <= h && x - r <= w)) {
281 return iterate(
282 x, y, w, 0, 0, h, r * r); // collision within triangle (0, h) (w, h)
283 // (0, 0) is possible
284 }
285 return false;
286 }
287 }
288 explicit EllipseCollisionTest(int maxIterations) : m_maxIterations(maxIterations)
289 {
290
291 }
292 };
293
294 EllipseCollisionTest test(4);
295 // convert to cartesian coordinates
296 AmgMatrix(2, 2) covRotMatrix;
297 // cppcheck-suppress constStatement
298 covRotMatrix << scResult.cosC, -locpo(0, 0) * scResult.sinC, scResult.sinC,
299 locpo(0, 0) * scResult.cosC;
300 AmgMatrix(2, 2) lCovarianceCar =
301 covRotMatrix * bchk.lCovariance * covRotMatrix.transpose();
302 Amg::Vector2D locpoCar(covRotMatrix(1, 1), -covRotMatrix(0, 1));
303
304 // ellipse is always at (0,0), surface is moved to ellipse position and then
305 // rotated
306 double w = bchk.nSigmas * sqrt(lCovarianceCar(0, 0));
307 double h = bchk.nSigmas * sqrt(lCovarianceCar(1, 1));
308 double x0 = 0;
309 double y0 = 0;
310 float theta =
311 (lCovarianceCar(1, 0) != 0 &&
312 (lCovarianceCar(1, 1) - lCovarianceCar(0, 0)) != 0)
313 ? .5 * bchk.FastArcTan(2 * lCovarianceCar(1, 0) /
314 (lCovarianceCar(1, 1) - lCovarianceCar(0, 0)))
315 : 0.;
316 scResult = bchk.FastSinCos(theta);
317 AmgMatrix(2, 2) rotMatrix;
318 rotMatrix << scResult.cosC, scResult.sinC, -scResult.sinC, scResult.cosC;
319 Amg::Vector2D tmp = rotMatrix * (-locpoCar);
320 double x1 = tmp(0, 0);
321 double y1 = tmp(1, 0);
323 // check if ellipse and circle overlap and return result
324 return test.collide(x0, y0, w, h, x1, y1, r);
325}
const boost::regex rr(r_r)
#define M_PI
#define AmgMatrix(rows, cols)
#define y
#define x
virtual double r() const override
This method returns the maximum expansion on the plane (=rMax)
DiscTrapezoidalBounds()
Default Constructor.
virtual bool inside(const Amg::Vector2D &locpo, double tol1=0., double tol2=0.) const override final
This method cheks if the radius given in the LocalPosition is inside [rMin,rMax] if only tol1 is give...
int r
Definition globals.cxx:22
iterate(ROOT.TDirectory thisdir, ROOT.TDirectory targetdir, str prefix, typing.Pattern regex, bool excludeTrees)
@ theta
Definition ParamDefs.h:66
@ locPhi
local polar
Definition ParamDefs.h:45

◆ inside() [2/2]

bool Trk::DiscTrapezoidalBounds::inside ( const Amg::Vector2D & locpo,
double tol1 = 0.,
double tol2 = 0. ) const
finaloverridevirtual

This method cheks if the radius given in the LocalPosition is inside [rMin,rMax] if only tol1 is given and additional in the phi sector is tol2 is given.

Implements Trk::SurfaceBounds.

Definition at line 82 of file DiscTrapezoidalBounds.cxx.

83{
86
87 Amg::Vector2D cartCenter(rMedium * cos(phi), rMedium * sin(phi));
88
89 Amg::Vector2D cartPos(locpo[Trk::locR] * cos(locpo[Trk::locPhi]),
90 locpo[Trk::locR] * sin(locpo[Trk::locPhi]));
91
92 Amg::Vector2D Pos = cartPos - cartCenter;
93
94 Amg::Vector2D DeltaPos(Pos[Trk::locX] * sin(phi) - Pos[Trk::locY] * cos(phi),
95 Pos[Trk::locY] * sin(phi) + Pos[Trk::locX] * cos(phi));
96
97 bool insideX = (std::abs(DeltaPos[Trk::locX]) <
99 bool insideY = (std::abs(DeltaPos[Trk::locY]) <
101
102 if (!insideX || !insideY)
103 return false;
104
105 if (std::abs(DeltaPos[Trk::locX]) <
107 return true;
108
112
118
119 bool inside = (DeltaPos[Trk::locY] > (m * std::abs(DeltaPos[Trk::locX]) + q));
120
121 return inside;
122}
Eigen::Matrix< double, 2, 1 > Vector2D
@ locY
local cartesian
Definition ParamDefs.h:38
@ locX
Definition ParamDefs.h:37
@ locR
Definition ParamDefs.h:44
@ phi
Definition ParamDefs.h:75

◆ insideLoc1()

virtual bool Trk::DiscTrapezoidalBounds::insideLoc1 ( const Amg::Vector2D & locpo,
double tol1 = 0. ) const
finaloverridevirtual

This method checks inside bounds in loc1.

  • loc1/loc2 correspond to the natural coordinates of the surface

Implements Trk::SurfaceBounds.

◆ insideLoc2()

virtual bool Trk::DiscTrapezoidalBounds::insideLoc2 ( const Amg::Vector2D & locpo,
double tol2 = 0. ) const
finaloverridevirtual

This method checks inside bounds in loc2.

  • loc1/loc2 correspond to the natural coordinates of the surface

Implements Trk::SurfaceBounds.

◆ maxHalflengthX()

double Trk::DiscTrapezoidalBounds::maxHalflengthX ( ) const

This method returns the maximal halflength in X.

◆ minDistance()

double Trk::DiscTrapezoidalBounds::minDistance ( const Amg::Vector2D & pos) const
finaloverridevirtual

Minimal distance to boundary ( > 0 if outside and <=0 if inside)

Implements Trk::SurfaceBounds.

Definition at line 329 of file DiscTrapezoidalBounds.cxx.

330{
331 const double pi2 = 2. * M_PI;
332 double alpha = std::abs(pos[locPhi]);
333 if (alpha > M_PI)
334 alpha = pi2 - alpha;
335
336 double r = pos[locR];
337 if (r == 0.)
339 cos(alpha);
340
341 // check if it is external in R
344 r;
345 if (sr0 > 0.)
346 return sr0;
349 if (sr1 > 0.)
350 return sr1;
351
352 // check if it is external in phi
354 return r * std::abs(sin(alpha - m_boundValues[DiscTrapezoidalBounds::bv_halfPhiSector]));
355
356 // if here it is inside:
357 // Evaluate the distance from the 4 segments
358 double dist[4] = { sr0,
359 sr1,
362
363 return *std::max_element(dist, dist + 4);
364}

◆ minHalflengthX()

double Trk::DiscTrapezoidalBounds::minHalflengthX ( ) const

This method returns the minimal halflength in X.

◆ operator!=()

bool Trk::SurfaceBounds::operator!= ( const SurfaceBounds & sb) const
inlinevirtualinherited

Non-Equality operator.

Reimplemented in Trk::InvalidBounds.

Definition at line 141 of file SurfaceBounds.h.

142{
143 return !((*this) == sb);
144}

◆ operator=() [1/2]

Trk::DiscTrapezoidalBounds & Trk::DiscTrapezoidalBounds::operator= ( const DiscTrapezoidalBounds & disctrbo)

Assignment operator.

Definition at line 64 of file DiscTrapezoidalBounds.cxx.

65{
66 if (this != &disctrbo)
67 m_boundValues = disctrbo.m_boundValues;
68 return *this;
69}

◆ operator=() [2/2]

DiscTrapezoidalBounds & Trk::DiscTrapezoidalBounds::operator= ( DiscTrapezoidalBounds && disctrbo)
defaultnoexcept

Move Assignment operator.

◆ operator==()

bool Trk::DiscTrapezoidalBounds::operator== ( const SurfaceBounds & sbo) const
overridevirtual

Equality operator.

Implements Trk::SurfaceBounds.

Definition at line 72 of file DiscTrapezoidalBounds.cxx.

73{
74 // check the type first not to compare apples with oranges
75 const Trk::DiscTrapezoidalBounds* disctrbo = dynamic_cast<const Trk::DiscTrapezoidalBounds*>(&sbo);
76 if (!disctrbo)
77 return false;
78 return (m_boundValues == disctrbo->m_boundValues);
79}

◆ r()

virtual double Trk::DiscTrapezoidalBounds::r ( ) const
overridevirtual

This method returns the maximum expansion on the plane (=rMax)

Implements Trk::SurfaceBounds.

◆ rCenter()

double Trk::DiscTrapezoidalBounds::rCenter ( ) const

This method returns the center radius.

◆ rMax()

double Trk::DiscTrapezoidalBounds::rMax ( ) const

This method returns outer radius.

◆ rMin()

double Trk::DiscTrapezoidalBounds::rMin ( ) const

This method returns inner radius.

◆ stereo()

double Trk::DiscTrapezoidalBounds::stereo ( ) const

This method returns the stereo angle.

◆ swap()

void Trk::SurfaceBounds::swap ( double & b1,
double & b2 )
inlineprotectedinherited

Swap method to be called from DiscBounds or TrapezoidalBounds.

Definition at line 133 of file SurfaceBounds.h.

134{
135 double tmp = b1;
136 b1 = b2;
137 b2 = tmp;
138}

◆ type()

virtual SurfaceBounds::BoundsType Trk::DiscTrapezoidalBounds::type ( ) const
inlinefinaloverridevirtual

Return the type - mainly for persistency.

Implements Trk::SurfaceBounds.

Definition at line 88 of file DiscTrapezoidalBounds.h.

Member Data Documentation

◆ m_boundValues

std::vector<TDD_real_t> Trk::DiscTrapezoidalBounds::m_boundValues
private

Internal members of the bounds (float/double)

Definition at line 144 of file DiscTrapezoidalBounds.h.


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