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 
24 
25 #include <utility>
26 
27 #include <sstream>
28 
29 #include <vector>
30 #include <string>
31 
33 
34 class DetLayer;
35 class HitRZCompatibility;
36 class BarrelDetLayer;
37 class ForwardDetLayer;
38 
39 namespace edm { class Event; }
40 
42 public:
43 
44 public:
45  virtual ~TrackingRegion(){}
48  typedef std::vector<Hit> Hits;
49 
50 
51 public:
52 
54  const GlobalPoint & originPos,
55  const Range & invPtRange,
56  const float & originRBound,
57  const float & originZBound)
58  : theDirection( direction), theUnitDirection(direction.unit()), theVertexPos( originPos),
59  theInvPtRange( invPtRange), thePhi(direction.barePhi()),
60  thePtMin(1.f/std::max( std::abs(invPtRange.max()), std::abs(invPtRange.min()) )),
61  theVertexRBound( originRBound),
62  theVertexZBound( originZBound) { }
63 
64 
66  GlobalVector const & direction() const { return theDirection; }
67  GlobalVector const & unitDirection() const { return theUnitDirection; }
68 
69  float phiDirection() const { return thePhi;}
70 
75  GlobalPoint const & origin() const { return theVertexPos; }
76 
78  float originRBound() const { return theVertexRBound; }
79 
81  float originZBound() const { return theVertexZBound; }
82 
84  float ptMin() const { return thePtMin;}
85 
87  Range invPtRange() const { return theInvPtRange; }
88 
89 
92  virtual HitRZCompatibility * checkRZ(const DetLayer* layer,
93  const Hit & outerHit,
94  const edm::EventSetup& iSetup,
95  const DetLayer* outerlayer=0,
96  float lr=0, float gz=0, float dr=0, float dz=0) const = 0;
97 
98 
100  virtual Hits hits(
101  const edm::Event& ev,
102  const edm::EventSetup& es,
103  const ctfseeding::SeedingLayer* layer) const = 0;
104 
105  virtual Hits hits(
106  const edm::Event& ev,
107  const edm::EventSetup& es,
108  const SeedingLayerSetsHits::SeedingLayer& layer) const = 0;
109 
112  const float & originRBound, const float & originZBound) const {
113  TrackingRegion* restr = clone();
114  restr->theVertexPos = originPos;
115  restr->theVertexRBound = originRBound;
116  restr->theVertexZBound = originZBound;
117  return restr;
118  }
119 
120  virtual TrackingRegion* clone() const = 0;
121 
122  virtual std::string name() const { return "TrackingRegion"; }
123  virtual std::string print() const {
124  std::ostringstream str;
125  str << name() <<" dir:"<<theDirection<<" vtx:"<<theVertexPos
126  <<" dr:"<<theVertexRBound<<" dz:"<<theVertexZBound<<" pt:"<<1./theInvPtRange.max();
127  return str.str();
128  }
129 
130  // void setDirection(const GlobalVector & dir ) { theDirection = dir; }
131 
132 private:
133 
134  GlobalVector theDirection; // this is not direction is momentum...
135  GlobalVector theUnitDirection; // the real direction
138  float thePhi;
139  float thePtMin;
142 
143 };
144 
146 
147 #endif
float originRBound() const
bounds the particle vertex in the transverse plane
T barePhi() const
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
GlobalVector theUnitDirection
GlobalVector theDirection
TrackingRegion(const GlobalVector &direction, const GlobalPoint &originPos, const Range &invPtRange, const float &originRBound, const float &originZBound)
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
string unit
Definition: csvLumiCalc.py:46
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
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 phiDirection() const
float originZBound() const
bounds the particle vertex in the longitudinal plane
virtual TrackingRegion * clone() const =0
GlobalVector const & unitDirection() const
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