CMS 3D CMS Logo

All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TrackingRegion.h
Go to the documentation of this file.
1 #ifndef TrackingRegion_H
2 #define TrackingRegion_H
3 
11 
18 
19 
23 
24 #include <utility>
25 
26 #include <sstream>
27 
28 #include <vector>
29 #include <string>
30 
32 
33 class DetLayer;
34 class HitRZCompatibility;
35 class BarrelDetLayer;
36 class ForwardDetLayer;
37 
38 namespace edm { class Event; }
39 
41 public:
42 
43 public:
44  virtual ~TrackingRegion(){}
47  typedef std::vector<Hit> Hits;
48 
49 
50 public:
51 
53  const GlobalPoint & originPos,
54  const Range & invPtRange,
55  const float & originRBound,
56  const float & originZBound)
57  : theDirection( direction), theVertexPos( originPos),
58  theInvPtRange( invPtRange),
59  thePtMin(1.f/std::max( std::abs(invPtRange.max()), std::abs(invPtRange.min()) )),
60  theVertexRBound( originRBound),
61  theVertexZBound( originZBound) { }
62 
63 
65  GlobalVector const & direction() const { return theDirection; }
66 
71  GlobalPoint const & origin() const { return theVertexPos; }
72 
74  float originRBound() const { return theVertexRBound; }
75 
77  float originZBound() const { return theVertexZBound; }
78 
80  float ptMin() const { return thePtMin;}
81 
83  Range invPtRange() const { return theInvPtRange; }
84 
85 
88  virtual HitRZCompatibility * checkRZ(const DetLayer* layer,
89  const Hit & outerHit,
90  const edm::EventSetup& iSetup,
91  const DetLayer* outerlayer=0,
92  float lr=0, float gz=0, float dr=0, float dz=0) const = 0;
93 
94 
96  virtual Hits hits(
97  const edm::Event& ev,
98  const edm::EventSetup& es,
99  const ctfseeding::SeedingLayer* layer) const = 0;
100 
103  const float & originRBound, const float & originZBound) const {
104  TrackingRegion* restr = clone();
105  restr->theVertexPos = originPos;
106  restr->theVertexRBound = originRBound;
107  restr->theVertexZBound = originZBound;
108  return restr;
109  }
110 
111  virtual TrackingRegion* clone() const = 0;
112 
113  virtual std::string name() const { return "TrackingRegion"; }
114  virtual std::string print() const {
115  std::ostringstream str;
116  str << name() <<" dir:"<<theDirection<<" vtx:"<<theVertexPos
117  <<" dr:"<<theVertexRBound<<" dz:"<<theVertexZBound<<" pt:"<<1./theInvPtRange.max();
118  return str.str();
119  }
120 
122 
123 private:
124 
128  float thePtMin;
131 
132 };
133 
135 
136 #endif
float originRBound() const
bounds the particle vertex in the transverse plane
virtual Hits hits(const edm::Event &ev, const edm::EventSetup &es, const ctfseeding::SeedingLayer *layer) const =0
get hits from layer compatible with region constraints
GlobalPoint const & origin() const
T max() const
#define abs(x)
Definition: mlp_lapack.h:159
GlobalVector theDirection
TrackingRegion(const GlobalVector &direction, const GlobalPoint &originPos, const Range &invPtRange, const float &originRBound, const float &originZBound)
#define min(a, b)
Definition: mlp_lapack.h:161
void setDirection(const GlobalVector &dir)
virtual ~TrackingRegion()
GlobalVector const & direction() const
the direction around which region is constructed
GlobalPoint theVertexPos
virtual HitRZCompatibility * checkRZ(const DetLayer *layer, const Hit &outerHit, const edm::EventSetup &iSetup, const DetLayer *outerlayer=0, float lr=0, float gz=0, float dr=0, float dz=0) const =0
const T & max(const T &a, const T &b)
PixelRecoRange< float > Range
double f[11][100]
TrackingRegion * restrictedRegion(const GlobalPoint &originPos, const float &originRBound, const float &originZBound) const
clone region with new vertex position
virtual std::string name() const
float originZBound() const
bounds the particle vertex in the longitudinal plane
virtual TrackingRegion * clone() const =0
TransientTrackingRecHit::ConstRecHitPointer Hit
std::vector< Hit > Hits
float ptMin() const
minimal pt of interest
virtual std::string print() const
Range invPtRange() const
inverse pt range
dbl *** dir
Definition: mlp_gen.cc:35