CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CRackTrajectoryBuilder.h
Go to the documentation of this file.
1 #ifndef CRackTrajectoryBuilder_h
2 #define CRackTrajectoryBuilder_h
3 
4 //
5 // Package: RecoTracker/SingleTrackPattern
6 // Class: CRackTrajectoryBuilder
7 // Original Author: Michele Pioppi-INFN perugia
8 
9 #include <string>
10 
13 
24 
26 
40 
42 
43 
44 //to sort hits by the det position
46  public:
48  bool operator()( const TrackingRecHit *rh1,
49  const TrackingRecHit *rh2)
50  {
51  const GeomDet* detPos1 = _tracker.idToDet(rh1->geographicalId());
52  const GeomDet* detPos2 = _tracker.idToDet(rh2->geographicalId());
53 
54  GlobalPoint gp1 = detPos1->position();
55  GlobalPoint gp2 = detPos2->position();
56 
57  if (gp1.y()>gp2.y())
58  return true;
59  if (gp1.y()<gp2.y())
60  return false;
61  // if (gp1.y()== gp2.y())
62  //
63  return (rh1->geographicalId() < rh2->geographicalId());
64  };
65  private:
66  // edm::ESHandle<TrackerGeometry> _tracker;
68  };
69 
71  public:
73  bool operator()( const TrackingRecHit *rh1,
74  const TrackingRecHit *rh2)
75  {
76  const GeomDet* detPos1 = _tracker.idToDet(rh1->geographicalId());
77  const GeomDet* detPos2 = _tracker.idToDet(rh2->geographicalId());
78 
79  GlobalPoint gp1 = detPos1->position();
80  GlobalPoint gp2 = detPos2->position();
81 
82  if (gp1.y()<gp2.y())
83  return true;
84  if (gp1.y()>gp2.y())
85  return false;
86  // if (gp1.y()== gp2.y())
87  //
88  return (rh1->geographicalId() < rh2->geographicalId());
89  };
90  private:
91  // edm::ESHandle<TrackerGeometry> _tracker;
93  };
94 
95 #ifndef TrajectoryBuilder_CompareHitY
96 #define TrajectoryBuilder_CompareHitY
97 
98 class CompareHitY {
99  public:
101  bool operator()( const TrackingRecHit *rh1,
102  const TrackingRecHit *rh2)
103  {
104  GlobalPoint gp1=_tracker.idToDet(rh1->geographicalId())->surface().toGlobal(rh1->localPosition());
105  GlobalPoint gp2=_tracker.idToDet(rh2->geographicalId())->surface().toGlobal(rh2->localPosition());
106  return gp1.y()<gp2.y();};
107  private:
108  // edm::ESHandle<TrackerGeometry> _tracker;
109  const TrackerGeometry& _tracker;
110 };
111 
112  class CompareHitY_plus {
113  public:
115  bool operator()( const TrackingRecHit *rh1,
116  const TrackingRecHit *rh2)
117  {
118  GlobalPoint gp1=_tracker.idToDet(rh1->geographicalId())->surface().toGlobal(rh1->localPosition());
119  GlobalPoint gp2=_tracker.idToDet(rh2->geographicalId())->surface().toGlobal(rh2->localPosition());
120  return gp1.y()>gp2.y();};
121  private:
122  // edm::ESHandle<TrackerGeometry> _tracker;
123  const TrackerGeometry& _tracker;
124  };
125 
126 #endif
127 
129 {
130 // using namespace std;
131 
134 
135  typedef std::vector<const TrackingRecHit*>::iterator TrackingRecHitIterator;
136 
137  typedef std::pair<TrackingRecHitIterator, TrackingRecHitIterator> TrackingRecHitRange;
138  typedef std::vector<TrackingRecHitRange>::iterator TrackingRecHitRangeIterator;
139 
140  // typedef std::pair<TrackingRecHitIterator, TSOS> PairTrackingRecHitTsos;
141  typedef std::pair<TrackingRecHitRangeIterator, TSOS> PairTrackingRecHitTsos;
142 
143  public:
145  friend class CompareDetByTraj;
146 
148  public:
149  CompareDetByTraj(const TSOS& tSos ):_tSos(tSos)
150  {};
151  bool operator()( const std::pair<TrackingRecHitRangeIterator, TSOS> rh1,
152  const std::pair<TrackingRecHitRangeIterator, TSOS> rh2)
153  {
154  GlobalPoint gp1 = rh1.second.globalPosition();
155  GlobalPoint gp2 = rh2.second.globalPosition();
156 
158  GlobalVector gpDiff1 = gp1-gpT;
159  GlobalVector gpDiff2 = gp2-gpT;
160 
161  //this might have a better performance ...
162  // float dist1 = ( gp1.x()-gpT.x() ) * ( gp1.x()-gpT.x() ) + ( gp1.y()-gpT.y() ) * ( gp1.y()-gpT.y() ) + ( gp1.z()-gpT.z() ) * ( gp1.z()-gpT.z() );
163  // float dist2 = ( gp2.x()-gpT.x() ) * ( gp2.x()-gpT.x() ) + ( gp2.y()-gpT.y() ) * ( gp2.y()-gpT.y() ) + ( gp2.z()-gpT.z() ) * ( gp2.z()-gpT.z() );
164  //if ( dist1<dist2 )
165 
166  // if ( gpDiff1.mag2() < gpDiff2.mag2() )
167 
168 
169  float dist1 = gpDiff1 * _tSos.globalDirection();
170  float dist2 = gpDiff2 * _tSos.globalDirection();
171 
172  if (dist1 < 0)
173  return false;
174  if ( dist1<dist2 )
175  return true;
176 
177  return false;
178  };
179  private:
181  };
182 
183 
184 
185  public:
186 
189 
191 
192  void run(const TrajectorySeedCollection &collseed,
193  const SiStripRecHit2DCollection &collstereo,
194  const SiStripRecHit2DCollection &collrphi ,
195  const SiStripMatchedRecHit2DCollection &collmatched,
196  const SiPixelRecHitCollection &collpixel,
197  const edm::EventSetup& es,
198  edm::Event& e,
199  std::vector<Trajectory> &trajoutput);
200 
201  void init(const edm::EventSetup& es,bool);
203  private:
204  std::vector<TrajectoryMeasurement> seedMeasurements(const TrajectorySeed& seed) const;
205 
206 
207  std::vector<const TrackingRecHit*> SortHits(const SiStripRecHit2DCollection &collstereo,
208  const SiStripRecHit2DCollection &collrphi ,
209  const SiStripMatchedRecHit2DCollection &collmatched,
210  const SiPixelRecHitCollection &collpixel,
211  const TrajectorySeed &seed,
212  const bool bAddSeedHits
213  );
214 
215 
216  // std::vector<TrackingRecHitRange> SortByTrajectory (const std::vector<TrackingRecHitRange>& inputHits);
217 
218 
219  TSOS startingTSOS(const TrajectorySeed& seed)const;
220  void updateTrajectory( Trajectory& traj,
221  const TM& tm,
222  const TransientTrackingRecHit& hit) const;
223 
224  void AddHit(Trajectory &traj,
225  std::vector<const TrackingRecHit*>Hits,
226  Propagator *currPropagator
227  );
228  // edm::OwnVector<TransientTrackingRecHit> hits);
229  bool qualityFilter(Trajectory traj);
230 
231  bool isDifferentStripReHit2D (const SiStripRecHit2D& hitA, const SiStripRecHit2D& hitB );
232 
233  std::pair<TrajectoryStateOnSurface, const GeomDet*>
234  innerState( const Trajectory& traj) const;
235 
236  private:
240 
243 
244 // AnalyticalPropagator *thePropagator;
245 // AnalyticalPropagator *thePropagatorOp;
246 
252 // const KFTrajectoryFitter * theFitterOp;
253 
257 
259  double chi2cut;
260  std::vector<Trajectory> trajFit;
261  //RC edm::OwnVector<const TransientTrackingRecHit> hits;
263  bool seed_plus;
264  std::string geometry;
265 // TransientInitialStateEstimator* theInitialState;
266 };
267 
268 #endif
edm::ESHandle< MagneticField > magfield
std::vector< TrackingRecHitRange >::iterator TrackingRecHitRangeIterator
Trajectory createStartingTrajectory(const TrajectorySeed &seed) const
bool operator()(const TrackingRecHit *rh1, const TrackingRecHit *rh2)
bool qualityFilter(Trajectory traj)
CompareHitY_plus(const TrackerGeometry &tracker)
Chi2MeasurementEstimator * theEstimator
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:47
T y() const
Definition: PV3DBase.h:62
const TrackerGeometry & _tracker
GlobalPoint globalPosition() const
const TrackerGeometry & _tracker
const TransientTrackingRecHitBuilder * RHBuilder
TransientTrackingRecHit::RecHitContainer hits
std::vector< ConstRecHitPointer > RecHitContainer
std::vector< const TrackingRecHit * >::iterator TrackingRecHitIterator
CompareDetY_minus(const TrackerGeometry &tracker)
const KFTrajectorySmoother * theSmoother
std::pair< TrajectoryStateOnSurface, const GeomDet * > innerState(const Trajectory &traj) const
std::vector< Trajectory > trajFit
std::pair< TrackingRecHitRangeIterator, TSOS > PairTrackingRecHitTsos
const Surface::PositionType & position() const
The position (origin of the R.F.)
Definition: GeomDet.h:41
bool operator()(const TrackingRecHit *rh1, const TrackingRecHit *rh2)
std::vector< TrajectorySeed > TrajectorySeedCollection
void run(const TrajectorySeedCollection &collseed, const SiStripRecHit2DCollection &collstereo, const SiStripRecHit2DCollection &collrphi, const SiStripMatchedRecHit2DCollection &collmatched, const SiPixelRecHitCollection &collpixel, const edm::EventSetup &es, edm::Event &e, std::vector< Trajectory > &trajoutput)
Runs the algorithm.
CompareDetY_plus(const TrackerGeometry &tracker)
edm::ESHandle< TrackerGeometry > tracker
CompareHitY(const TrackerGeometry &tracker)
PropagatorWithMaterial * thePropagatorOp
std::vector< TrajectoryMeasurement > seedMeasurements(const TrajectorySeed &seed) const
virtual const GeomDet * idToDet(DetId) const
TSOS startingTSOS(const TrajectorySeed &seed) const
tuple conf
Definition: dbtoconf.py:185
const TrajectoryStateOnSurface & _tSos
bool operator()(const TrackingRecHit *rh1, const TrackingRecHit *rh2)
bool operator()(const std::pair< TrackingRecHitRangeIterator, TSOS > rh1, const std::pair< TrackingRecHitRangeIterator, TSOS > rh2)
void AddHit(Trajectory &traj, std::vector< const TrackingRecHit * >Hits, Propagator *currPropagator)
CRackTrajectoryBuilder(const edm::ParameterSet &conf)
bool isDifferentStripReHit2D(const SiStripRecHit2D &hitA, const SiStripRecHit2D &hitB)
std::pair< TrackingRecHitIterator, TrackingRecHitIterator > TrackingRecHitRange
void init(const edm::EventSetup &es, bool)
TrajectoryStateOnSurface TSOS
PropagatorWithMaterial * thePropagator
const KFTrajectoryFitter * theFitter
void updateTrajectory(Trajectory &traj, const TM &tm, const TransientTrackingRecHit &hit) const
std::vector< const TrackingRecHit * > SortHits(const SiStripRecHit2DCollection &collstereo, const SiStripRecHit2DCollection &collrphi, const SiStripMatchedRecHit2DCollection &collmatched, const SiPixelRecHitCollection &collpixel, const TrajectorySeed &seed, const bool bAddSeedHits)
DetId geographicalId() const
virtual LocalPoint localPosition() const =0
const TrackerGeometry & _tracker
const TrackerGeometry & _tracker
bool operator()(const TrackingRecHit *rh1, const TrackingRecHit *rh2)
GlobalVector globalDirection() const
TrajectoryMeasurement TM