CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CkfTrajectoryBuilder.cc
Go to the documentation of this file.
2 
4 
6 
10 
16 
17 
20 
22 
26 
28 
29 using namespace std;
30 
31 
38  const TransientTrackingRecHitBuilder* recHitBuilder,
39  const MeasurementTracker* measurementTracker,
40  const TrajectoryFilter* filter):
41 
43  updator, propagatorAlong,propagatorOpposite,
44  estimator, recHitBuilder, measurementTracker,filter)
45 {
46  theMaxCand = conf.getParameter<int>("maxCand");
47  theLostHitPenalty = conf.getParameter<double>("lostHitPenalty");
48  theIntermediateCleaning = conf.getParameter<bool>("intermediateCleaning");
49  theAlwaysUseInvalidHits = conf.getParameter<bool>("alwaysUseInvalidHits");
50  /*
51  theSharedSeedCheck = conf.getParameter<bool>("SharedSeedCheck");
52  std::stringstream ss;
53  ss<<"CkfTrajectoryBuilder_"<<conf.getParameter<std::string>("ComponentName")<<"_"<<this;
54  theUniqueName = ss.str();
55  LogDebug("CkfPattern")<<"my unique name is: "<<theUniqueName;
56  */
57 }
58 
60 {
61  theMeasurementTracker->update(event);
62 }
63 
66 {
68  result.reserve(5);
69  trajectories(seed, result);
70  return result;
71 }
72 
73 /*
74  void CkfTrajectoryBuilder::rememberSeedAndTrajectories(const TrajectorySeed& seed,
75  CkfTrajectoryBuilder::TrajectoryContainer &result) const
76  {
77 
78  //result ----> theCachedTrajectories
79  //every first iteration on event. forget about everything that happened before
80  if (edm::Service<UpdaterService>()->checkOnce(theUniqueName))
81  theCachedTrajectories.clear();
82 
83  if (seed.nHits()==0) return;
84 
85  //then remember those trajectories
86  for (TrajectoryContainer::iterator traj=result.begin();
87  traj!=result.end(); ++traj) {
88  theCachedTrajectories.insert(std::make_pair(seed.recHits().first->geographicalId(),*traj));
89  }
90  }
91 
92  bool CkfTrajectoryBuilder::sharedSeed(const TrajectorySeed& s1,const TrajectorySeed& s2) const{
93  //quit right away on nH=0
94  if (s1.nHits()==0 || s2.nHits()==0) return false;
95  //quit right away if not the same number of hits
96  if (s1.nHits()!=s2.nHits()) return false;
97  TrajectorySeed::range r1=s1.recHits();
98  TrajectorySeed::range r2=s2.recHits();
99  TrajectorySeed::const_iterator i1,i2;
100  TrajectorySeed::const_iterator & i1_e=r1.second,&i2_e=r2.second;
101  TrajectorySeed::const_iterator & i1_b=r1.first,&i2_b=r2.first;
102  //quit right away if first detId does not match. front exist because of ==0 ->quit test
103  if(i1_b->geographicalId() != i2_b->geographicalId()) return false;
104  //then check hit by hit if they are the same
105  for (i1=i1_b,i2=i2_b;i1!=i1_e,i2!=i2_e;++i1,++i2){
106  if (!i1->sharesInput(&(*i2),TrackingRecHit::all)) return false;
107  }
108  return true;
109  }
110  bool CkfTrajectoryBuilder::seedAlreadyUsed(const TrajectorySeed& seed,
111  CkfTrajectoryBuilder::TempTrajectoryContainer &candidates) const
112  {
113  //theCachedTrajectories ---> candidates
114  if (seed.nHits()==0) return false;
115  bool answer=false;
116  pair<SharedTrajectory::const_iterator, SharedTrajectory::const_iterator> range =
117  theCachedTrajectories.equal_range(seed.recHits().first->geographicalId());
118  SharedTrajectory::const_iterator trajP;
119  for (trajP = range.first; trajP!=range.second;++trajP){
120  //check whether seeds are identical
121  if (sharedSeed(trajP->second.seed(),seed)){
122  candidates.push_back(trajP->second);
123  answer=true;
124  }//already existing trajectory shares the seed.
125  }//loop already made trajectories
126 
127  return answer;
128  }
129 */
130 
131 void
133 {
134  // analyseSeed( seed);
135  /*
136  if (theSharedSeedCheck){
137  TempTrajectoryContainer candidates;
138  if (seedAlreadyUsed(seed,candidates))
139  {
140  //start with those candidates already made before
141  limitedCandidates(candidates,result);
142  //and quit
143  return;
144  }
145  }
146  */
147 
148  TempTrajectory startingTraj = createStartingTrajectory( seed );
149 
152  limitedCandidates( startingTraj, result);
153 
154  /*
155  //and remember what you just did
156  if (theSharedSeedCheck) rememberSeedAndTrajectories(seed,result);
157  */
158 
159  // analyseResult(result);
160 }
161 
165 {
166  TempTrajectoryContainer candidates;
167  candidates.push_back( startingTraj);
168  limitedCandidates(candidates,result);
169 }
170 
174 {
175  unsigned int nIter=1;
176  // TempTrajectoryContainer candidates; // = TrajectoryContainer();
177  TempTrajectoryContainer newCand; // = TrajectoryContainer();
178  // candidates.push_back( startingTraj);
179 
180  while ( !candidates.empty()) {
181 
182  newCand.clear();
183  for (TempTrajectoryContainer::iterator traj=candidates.begin();
184  traj!=candidates.end(); traj++) {
185  std::vector<TM> meas;
186  findCompatibleMeasurements(*traj, meas);
187 
188  // --- method for debugging
189  if(!analyzeMeasurementsDebugger(*traj,meas,
192  theTTRHBuilder)) return;
193  // ---
194 
195  if ( meas.empty()) {
196  if ( qualityFilter( *traj)) addToResult( *traj, result);
197  }
198  else {
199  std::vector<TM>::const_iterator last;
200  if ( theAlwaysUseInvalidHits) last = meas.end();
201  else {
202  if (meas.front().recHit()->isValid()) {
203  last = find_if( meas.begin(), meas.end(), RecHitIsInvalid());
204  }
205  else last = meas.end();
206  }
207 
208  for( std::vector<TM>::const_iterator itm = meas.begin();
209  itm != last; itm++) {
210  TempTrajectory newTraj = *traj;
211  updateTrajectory( newTraj, *itm);
212 
213  if ( toBeContinued(newTraj)) {
214  newCand.push_back(newTraj);
215  }
216  else {
217  if ( qualityFilter(newTraj)) addToResult( newTraj, result);
219  }
220  }
221  }
222 
223  if ((int)newCand.size() > theMaxCand) {
224  sort( newCand.begin(), newCand.end(), TrajCandLess<TempTrajectory>(theLostHitPenalty));
225  newCand.erase( newCand.begin()+theMaxCand, newCand.end());
226  }
227  }
228 
230 
231  candidates.swap(newCand);
232 
233  LogDebug("CkfPattern") <<result.size()<<" candidates after "<<nIter++<<" CKF iteration: \n"
235  <<"\n "<<candidates.size()<<" running candidates are: \n"
236  <<PrintoutHelper::dumpCandidates(candidates);
237 
238  }
239 }
240 
241 
242 
244  const TM& tm) const
245 {
246  TSOS predictedState = tm.predictedState();
248 
249  if ( hit->isValid()) {
250  TM tmp = TM( predictedState, theUpdator->update( predictedState, *hit),
251  hit, tm.estimate(), tm.layer());
252  traj.push(tmp );
253  }
254  else {
255  traj.push( TM( predictedState, hit, 0, tm.layer()));
256  }
257 }
258 
259 
260 void
262  std::vector<TrajectoryMeasurement> & result) const
263 {
264  int invalidHits = 0;
265  std::pair<TSOS,std::vector<const DetLayer*> > stateAndLayers = findStateAndLayers(traj);
266  if (stateAndLayers.second.empty()) return;
267 
268  vector<const DetLayer*>::iterator layerBegin = stateAndLayers.second.begin();
269  vector<const DetLayer*>::iterator layerEnd = stateAndLayers.second.end();
270  LogDebug("CkfPattern")<<"looping on "<< stateAndLayers.second.size()<<" layers.";
271  for (vector<const DetLayer*>::iterator il = layerBegin;
272  il != layerEnd; il++) {
273 
274  LogDebug("CkfPattern")<<"looping on a layer in findCompatibleMeasurements.\n last layer: "<<traj.lastLayer()<<" current layer: "<<(*il);
275 
276  TSOS stateToUse = stateAndLayers.first;
277  if ((*il)==traj.lastLayer())
278  {
279  LogDebug("CkfPattern")<<" self propagating in findCompatibleMeasurements.\n from: \n"<<stateToUse;
280  //self navigation case
281  // go to a middle point first
283  GlobalPoint center(0,0,0);
284  stateToUse = middle.extrapolate(stateToUse, center, *theForwardPropagator);
285 
286  if (!stateToUse.isValid()) continue;
287  LogDebug("CkfPattern")<<"to: "<<stateToUse;
288  }
289 
290  vector<TrajectoryMeasurement> tmp = theLayerMeasurements->measurements((**il),stateToUse, *theForwardPropagator, *theEstimator);
291 
292  if ( !tmp.empty()) {
293  if ( result.empty()) result = tmp;
294  else {
295  // keep one dummy TM at the end, skip the others
296  result.insert( result.end()-invalidHits, tmp.begin(), tmp.end());
297  }
298  invalidHits++;
299  }
300  }
301 
302  // sort the final result, keep dummy measurements at the end
303  if ( result.size() > 1) {
304  sort( result.begin(), result.end()-invalidHits, TrajMeasLessEstim());
305  }
306 
307  LogDebug("CkfPattern")<<"starting from:\n"
308  <<"x: "<<stateAndLayers.first.globalPosition()<<"\n"
309  <<"p: "<<stateAndLayers.first.globalMomentum()<<"\n"
311 
312 #ifdef DEBUG_INVALID
313  bool afterInvalid = false;
314  for (vector<TM>::const_iterator i=result.begin();
315  i!=result.end(); i++) {
316  if ( ! i->recHit().isValid()) afterInvalid = true;
317  if (afterInvalid && i->recHit().isValid()) {
318  edm::LogError("CkfPattern") << "CkfTrajectoryBuilder error: valid hit after invalid!" ;
319  }
320  }
321 #endif
322 
323  //analyseMeasurements( result, traj);
324 
325 }
326 
#define LogDebug(id)
T getParameter(std::string const &) const
std::vector< TrajectoryMeasurement > measurements(const DetLayer &layer, const TrajectoryStateOnSurface &startingState, const Propagator &prop, const MeasurementEstimator &est) const
int i
Definition: DBlmapReader.cc:9
static std::string dumpCandidates(collection &candidates)
std::vector< Trajectory > TrajectoryContainer
static void clean(TempTrajectoryContainer &tracks)
virtual TrajectoryStateOnSurface update(const TrajectoryStateOnSurface &, const TransientTrackingRecHit &) const =0
void addToResult(TempTrajectory &traj, TrajectoryContainer &result, bool inOut=false) const
const TrajectoryStateUpdator * theUpdator
ConstRecHitPointer recHit() const
static std::string dumpMeasurements(const std::vector< TrajectoryMeasurement > &v)
bool qualityFilter(const TempTrajectory &traj, bool inOut=false) const
const LayerMeasurements * theLayerMeasurements
virtual void setEvent(const edm::Event &event) const
set Event for the internal MeasurementTracker data member
const TransientTrackingRecHitBuilder * theTTRHBuilder
tuple result
Definition: query.py:137
virtual void findCompatibleMeasurements(const TempTrajectory &traj, std::vector< TrajectoryMeasurement > &result) const
CkfTrajectoryBuilder(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 DetLayer * layer() const
TrajectoryStateOnSurface predictedState() const
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
void limitedCandidates(TempTrajectory &startingTraj, TrajectoryContainer &result) const
tuple conf
Definition: dbtoconf.py:185
void updateTrajectory(TempTrajectory &traj, const TM &tm) const
const MeasurementTracker * theMeasurementTracker
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 TrajectoryContainer trajectories(const TrajectorySeed &seed) const
trajectories building starting from a seed
std::vector< TempTrajectory > TempTrajectoryContainer
TrajectoryStateOnSurface extrapolate(const FreeTrajectoryState &fts, const GlobalPoint &vtx) const
extrapolation with default (=geometrical) propagator
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
bool toBeContinued(TempTrajectory &traj, bool inOut=false) const
StateAndLayers findStateAndLayers(const TempTrajectory &traj) const
virtual bool analyzeMeasurementsDebugger(Trajectory &traj, std::vector< TrajectoryMeasurement > meas, const MeasurementTracker *theMeasurementTracker, const Propagator *theForwardPropagator, const Chi2MeasurementEstimatorBase *theEstimator, const TransientTrackingRecHitBuilder *theTTRHBuilder) const
const DetLayer * lastLayer() const
Redundant method, returns the layer of lastMeasurement() .
std::vector< Trajectory > TrajectoryContainer
const Chi2MeasurementEstimatorBase * theEstimator
TempTrajectory createStartingTrajectory(const TrajectorySeed &seed) const
void push(const TrajectoryMeasurement &tm)
const Propagator * theForwardPropagator