CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
BaseCkfTrajectoryBuilder.cc
Go to the documentation of this file.
2 
7 
8 
15 
19 
21 
22 
29  const TransientTrackingRecHitBuilder* recHitBuilder,
30  const MeasurementTracker* measurementTracker,
31  const TrajectoryFilter* filter,
32  const TrajectoryFilter* inOutFilter):
33  theUpdator(updator),
34  thePropagatorAlong(propagatorAlong),thePropagatorOpposite(propagatorOpposite),
35  theEstimator(estimator),theTTRHBuilder(recHitBuilder),
36  theMeasurementTracker(measurementTracker),
37  theLayerMeasurements(new LayerMeasurements(theMeasurementTracker)),
38  theForwardPropagator(0),theBackwardPropagator(0),
39  theFilter(filter),
40  theInOutFilter(inOutFilter)
41 {}
42 
44  delete theLayerMeasurements;
45 }
46 
47 
48 void
49 BaseCkfTrajectoryBuilder::seedMeasurements(const TrajectorySeed& seed, std::vector<TrajectoryMeasurement> & result) const
50 {
51  TrajectoryStateTransform tsTransform;
52 
53  TrajectorySeed::range hitRange = seed.recHits();
54  for (TrajectorySeed::const_iterator ihit = hitRange.first;
55  ihit != hitRange.second; ihit++) {
57  const GeomDet* hitGeomDet =
58  theMeasurementTracker->geomTracker()->idToDet( ihit->geographicalId());
59 
60  const DetLayer* hitLayer =
61  theMeasurementTracker->geometricSearchTracker()->detLayer(ihit->geographicalId());
62 
63  TSOS invalidState( new BasicSingleTrajectoryState( hitGeomDet->surface()));
64  if (ihit == hitRange.second - 1) {
65  // the seed trajectory state should correspond to this hit
66  PTrajectoryStateOnDet pState( seed.startingState());
67  const GeomDet* gdet = theMeasurementTracker->geomTracker()->idToDet( DetId(pState.detId()));
68  if (&gdet->surface() != &hitGeomDet->surface()) {
69  edm::LogError("CkfPattern") << "CkfTrajectoryBuilder error: the seed state is not on the surface of the detector of the last seed hit";
70  return; // FIXME: should throw exception
71  }
72 
73  TSOS updatedState = tsTransform.transientState( pState, &(gdet->surface()),
75  result.push_back(TM( invalidState, updatedState, recHit, 0, hitLayer));
76  }
77  else {
78  PTrajectoryStateOnDet pState( seed.startingState());
79 
80  TSOS outerState = tsTransform.transientState(pState,
82  (hitRange.second - 1)->geographicalId()))->surface()),
84  TSOS innerState = theBackwardPropagator->propagate(outerState,hitGeomDet->surface());
85  if(innerState.isValid()) {
86  TSOS innerUpdated = theUpdator->update(innerState,*recHit);
87  result.push_back(TM( invalidState, innerUpdated, recHit, 0, hitLayer));
88  }
89  }
90  }
91 
92  // method for debugging
93  fillSeedHistoDebugger(result.begin(),result.end());
94 
95 }
96 
97 
100 {
101  TempTrajectory result( seed, seed.direction());
102  if ( seed.direction() == alongMomentum) {
103  theForwardPropagator = &(*thePropagatorAlong);
104  theBackwardPropagator = &(*thePropagatorOpposite);
105  }
106  else {
107  theForwardPropagator = &(*thePropagatorOpposite);
108  theBackwardPropagator = &(*thePropagatorAlong);
109  }
110 
111  std::vector<TM> seedMeas;
112  seedMeasurements(seed, seedMeas);
113  for (std::vector<TM>::const_iterator i=seedMeas.begin(); i!=seedMeas.end(); i++)
114  result.push(*i);
115 
116  LogDebug("CkfPattern")
117  <<" initial trajectory from the seed: "<<PrintoutHelper::dumpCandidate(result,true);
118 
119  return result;
120 }
121 
122 
124 {
125  if (traj.measurements().size() > 400) {
126  edm::LogError("BaseCkfTrajectoryBuilder_InfiniteLoop");
127  LogTrace("BaseCkfTrajectoryBuilder_InfiniteLoop") <<
128  "Cropping Track After 400 Measurements:\n" <<
129  " Last predicted state: " << traj.lastMeasurement().predictedState() << "\n" <<
130  " Last layer subdetector: " << (traj.lastLayer() ? traj.lastLayer()->subDetector() : -1) << "\n" <<
131  " Found hits: " << traj.foundHits() << ", lost hits: " << traj.lostHits() << "\n\n";
132  return false;
133  }
134  // Called after each new hit is added to the trajectory, to see if it is
135  // worth continuing to build this track candidate.
136  if (inOut) {
137  if (theInOutFilter == 0) edm::LogError("CkfPattern") << "CkfTrajectoryBuilder error: trying to use dedicated filter for in-out tracking phase, when none specified";
138  return theInOutFilter->toBeContinued(traj);
139  } else {
140  return theFilter->toBeContinued(traj);
141  }
142 }
143 
144 
145  bool BaseCkfTrajectoryBuilder::qualityFilter( const TempTrajectory& traj, bool inOut) const
146 {
147  // Called after building a trajectory is completed, to see if it is good enough
148  // to keep.
149  if (inOut) {
150  if (theInOutFilter == 0) edm::LogError("CkfPattern") << "CkfTrajectoryBuilder error: trying to use dedicated filter for in-out tracking phase, when none specified";
151  return theInOutFilter->qualityFilter(traj);
152  } else {
153  return theFilter->qualityFilter(traj);
154  }
155 }
156 
157 
158 void
161  bool inOut) const
162 {
163  // quality check
164  if ( !qualityFilter(tmptraj, inOut) ) return;
165  Trajectory traj = tmptraj.toTrajectory();
166  // discard latest dummy measurements
167  while (!traj.empty() && !traj.lastMeasurement().recHit()->isValid()) traj.pop();
168  LogDebug("CkfPattern")<<inOut<<"=inOut option. pushing a Trajectory with: "<<traj.foundHits()<<" found hits. "<<traj.lostHits()
169  <<" lost hits. Popped :"<<(tmptraj.measurements().size())-(traj.measurements().size())<<" hits.";
170  result.push_back( traj);
171 }
172 void
175  bool inOut) const
176 {
177  // quality check
178  if ( !qualityFilter(tmptraj, inOut) ) return;
179  // discard latest dummy measurements
180  TempTrajectory traj = tmptraj;
181  while (!traj.empty() && !traj.lastMeasurement().recHit()->isValid()) traj.pop();
182  LogDebug("CkfPattern")<<inOut<<"=inOut option. pushing a TempTrajectory with: "<<traj.foundHits()<<" found hits. "<<traj.lostHits()
183  <<" lost hits. Popped :"<<(tmptraj.measurements().size())-(traj.measurements().size())<<" hits.";
184  result.push_back( traj );
185 }
186 
187 
188 
191 {
192  if (traj.empty())
193  {
194  //set the currentState to be the one from the trajectory seed starting point
195  PTrajectoryStateOnDet ptod = traj.seed().startingState();
196  DetId id(ptod.detId());
198  const Surface * surface=&g->surface();
199  TrajectoryStateTransform tsTransform;
200 
201  TSOS currentState = TrajectoryStateOnSurface(tsTransform.transientState(ptod,surface,theForwardPropagator->magneticField()));
203  return StateAndLayers(currentState,lastLayer->nextLayers( *currentState.freeState(), traj.direction()) );
204  }
205  else
206  {
207  TSOS currentState = traj.lastMeasurement().updatedState();
208  return StateAndLayers(currentState,traj.lastLayer()->nextLayers( *currentState.freeState(), traj.direction()) );
209  }
210 }
#define LogDebug(id)
PropagationDirection direction() const
bool empty() const
True if trajectory has no measurements.
Definition: Trajectory.h:200
const TrajectorySeed & seed() const
Access to the seed used to reconstruct the Trajectory.
int foundHits() const
Definition: Trajectory.h:190
int i
Definition: DBlmapReader.cc:9
int lostHits() const
Definition: Trajectory.h:197
const Propagator * theBackwardPropagator
bool empty() const
True if trajectory has no measurements.
virtual TrajectoryStateOnSurface update(const TrajectoryStateOnSurface &, const TransientTrackingRecHit &) const =0
void addToResult(TempTrajectory &traj, TrajectoryContainer &result, bool inOut=false) const
const DataContainer & measurements() const
virtual std::vector< const DetLayer * > nextLayers(NavigationDirection direction) const
Definition: DetLayer.cc:35
virtual SubDetector subDetector() const =0
The type of detector (PixelBarrel, PixelEndcap, TIB, TOB, TID, TEC, CSC, DT, RPCBarrel, RPCEndcap)
const TrajectoryStateUpdator * theUpdator
ConstRecHitPointer recHit() const
int foundHits() const
obsolete name, use measurements() instead.
BaseCkfTrajectoryBuilder(const edm::ParameterSet &conf, const TrajectoryStateUpdator *updator, const Propagator *propagatorAlong, const Propagator *propagatorOpposite, const Chi2MeasurementEstimatorBase *estimator, const TransientTrackingRecHitBuilder *RecHitBuilder, const MeasurementTracker *measurementTracker, const TrajectoryFilter *filter, const TrajectoryFilter *inOutFilter=0)
const TrajectoryFilter * theFilter
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
void seedMeasurements(const TrajectorySeed &seed, std::vector< TrajectoryMeasurement > &result) const
const TrajectoryMeasurement & lastMeasurement() const
bool qualityFilter(const TempTrajectory &traj, bool inOut=false) const
const LayerMeasurements * theLayerMeasurements
PropagationDirection direction() const
DataContainer const & measurements() const
Definition: Trajectory.h:169
FreeTrajectoryState * freeState(bool withErrors=true) const
Trajectory toTrajectory() const
Convert to a standard Trajectory.
recHitContainer::const_iterator const_iterator
virtual RecHitPointer build(const TrackingRecHit *p) const =0
build a tracking rechit from an existing rechit
TrajectoryMeasurement const & lastMeasurement() const
Definition: Trajectory.h:147
const TransientTrackingRecHitBuilder * theTTRHBuilder
tuple result
Definition: query.py:137
const DetLayer * detLayer(const DetId &id) const
obsolete method. Use idToLayer() instead.
std::pair< const_iterator, const_iterator > range
const TrackingGeometry * geomTracker() const
TrajectoryStateOnSurface transientState(const PTrajectoryStateOnDet &ts, const Surface *surface, const MagneticField *field) const
TrajectoryStateOnSurface updatedState() const
virtual bool qualityFilter(const TempTrajectory &) const =0
TrajectoryStateOnSurface predictedState() const
#define LogTrace(id)
static std::string dumpCandidate(const Candidate &candidate, bool showErrors=false)
tuple conf
Definition: dbtoconf.py:185
virtual TrajectoryStateOnSurface propagate(const FreeTrajectoryState &, const Surface &) const
Definition: Propagator.cc:9
Definition: DetId.h:20
const MeasurementTracker * theMeasurementTracker
PTrajectoryStateOnDet const & startingState() const
tuple filter
USE THIS FOR SKIMMED TRACKS process.p = cms.Path(process.hltLevel1GTSeed*process.skimming*process.offlineBeamSpot*process.TrackRefitter2) OTHERWISE USE THIS.
Definition: align_tpl.py:86
virtual const GeomDet * idToDet(DetId) const =0
void pop()
Definition: Trajectory.cc:17
std::vector< TempTrajectory > TempTrajectoryContainer
const unsigned int detId() const
range recHits() const
virtual const MagneticField * magneticField() const =0
bool toBeContinued(TempTrajectory &traj, bool inOut=false) const
std::pair< TSOS, std::vector< const DetLayer * > > StateAndLayers
StateAndLayers findStateAndLayers(const TempTrajectory &traj) const
virtual bool toBeContinued(TempTrajectory &) const =0
size_type size() const
Definition: bqueue.h:143
const TrajectoryFilter * theInOutFilter
const DetLayer * lastLayer() const
Redundant method, returns the layer of lastMeasurement() .
const GeometricSearchTracker * geometricSearchTracker() const
virtual void fillSeedHistoDebugger(std::vector< TrajectoryMeasurement >::iterator begin, std::vector< TrajectoryMeasurement >::iterator end) const
virtual const BoundPlane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
std::vector< Trajectory > TrajectoryContainer
TempTrajectory createStartingTrajectory(const TrajectorySeed &seed) const
int lostHits() const
const Propagator * theForwardPropagator