CMS 3D CMS Logo

DiskSectorBounds.h
Go to the documentation of this file.
1 #ifndef RecoTracker_TkDetLayers_DiskSectorBounds_h
2 #define RecoTracker_TkDetLayers_DiskSectorBounds_h
3 
7 #include <algorithm>
8 #include <cmath>
9 #include <cassert>
10 
11 #pragma GCC visibility push(hidden)
12 class DiskSectorBounds final : public Bounds {
13 public:
14  DiskSectorBounds(float rmin, float rmax, float zmin, float zmax, float phiExt)
15  : theRmin(rmin), theRmax(rmax), theZmin(zmin), theZmax(zmax), thePhiExtH(0.5f * phiExt) {
16  assert(thePhiExtH > 0);
17  if (theRmin > theRmax)
19  if (theZmin > theZmax)
21  theOffset = theRmin + 0.5f * (theRmax - theRmin);
22  }
23 
24  float length() const override { return theRmax - theRmin * std::cos(thePhiExtH); }
25  float width() const override { return 2.f * theRmax * std::sin(thePhiExtH); }
26  float thickness() const override { return theZmax - theZmin; }
27 
28  bool inside(const Local3DPoint& p) const override;
29 
30  bool inside(const Local3DPoint& p, const LocalError& err, float scale) const override;
31 
32  Bounds* clone() const override { return new DiskSectorBounds(*this); }
33 
34  float innerRadius() const { return theRmin; }
35  float outerRadius() const { return theRmax; }
36  float phiHalfExtension() const { return thePhiExtH; }
37 
38 private:
39  float theRmin;
40  float theRmax;
41  float theZmin;
42  float theZmax;
43  float thePhiExtH;
44  float theOffset;
45 };
46 
47 #pragma GCC visibility pop
48 #endif
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
float innerRadius() const
bool inside(const Local3DPoint &p) const override
Determine if the point is inside the bounds.
void swap(edm::DataFrameContainer &lhs, edm::DataFrameContainer &rhs)
float length() const override
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
float thickness() const override
double f[11][100]
Bounds * clone() const override
float phiHalfExtension() const
float width() const override
Definition: Bounds.h:20
float outerRadius() const
DiskSectorBounds(float rmin, float rmax, float zmin, float zmax, float phiExt)