CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
KalmanAlignmentTrackRefitter.cc
Go to the documentation of this file.
1 
3 
5 
7 
10 
13 
15 
17 
18 #include "CLHEP/GenericFunctions/CumulativeChiSquare.hh"
19 
20 #include <iostream>
21 
22 using namespace std;
23 using namespace reco;
24 using namespace Genfun;
25 
26 
28  theRefitterAlgo( config ),
29  theNavigator( navigator ),
30  theDebugFlag( config.getUntrackedParameter<bool>( "debug", true ) )
31 {
33  // --- GPetrucc: I can't understand where anything is read from the event, and who's the consumer.
34  // If there is one anywhere, it should do the consumes<T> calls and pass that to the setSrc
35  //TrackProducerBase< reco::Track >::setSrc( consumes<TrackCandidateCollection>(iConfig.getParameter<edm::InputTag>( "src" )),
36  // consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>( "beamSpot" )));
37 }
38 
39 
41 
42 
45  const AlignmentSetupCollection& algoSetups,
47  const reco::BeamSpot* beamSpot )
48 {
49  // Retrieve what we need from the EventSetup
51  edm::ESHandle< MagneticField > aMagneticField;
52  edm::ESHandle< TrajectoryFitter > aTrajectoryFitter;
53  edm::ESHandle< Propagator > aPropagator;
56 
57  getFromES( setup,
58  aGeometry,
59  aMagneticField,
60  aTrajectoryFitter,
61  aPropagator,
62  theMeasTk,
63  aRecHitBuilder );
64 
66  TrackCollection fullTracks;
67 
68  ConstTrajTrackPairCollection refittedFullTracks;
69  ConstTrajTrackPairCollection::const_iterator itTrack;
70 
71  for( itTrack = tracks.begin(); itTrack != tracks.end(); ++itTrack )
72  {
73  TransientTrack fullTrack( *(*itTrack).second, aMagneticField.product() );
74 
75  AlignmentSetupCollection::const_iterator itSetup;
76  for ( itSetup = algoSetups.begin(); itSetup != algoSetups.end(); ++itSetup )
77  {
78  RecHitContainer trackingRecHits;
79  RecHitContainer externalTrackingRecHits;
80 
82  RecHitContainer zPlusRecHits;
83  RecHitContainer zMinusRecHits;
84 
85  // Extract collection with TrackingRecHits
86  Trajectory::ConstRecHitContainer hits = (*itTrack).first->recHits();
87  Trajectory::ConstRecHitContainer::iterator itHits;
88 
89  for ( itHits = hits.begin(); itHits != hits.end(); ++itHits )
90  {
91  if ( !(*itHits)->isValid() ) continue;
92 
93  try
94  {
95  //if ( !theNavigator->alignableFromDetId( (*itHits)->geographicalId() )->alignmentParameters() ) continue;
96  theNavigator->alignableFromDetId( (*itHits)->geographicalId() );
97  } catch(...) { continue; }
98 
99  if ( (*itSetup)->useForTracking( *itHits ) )
100  {
101  trackingRecHits.push_back( (*itHits)->hit()->clone() );
102 
103  ( (*itHits)->det()->position().z() > 0. ) ?
104  zPlusRecHits.push_back( (*itHits)->hit()->clone() ) :
105  zMinusRecHits.push_back( (*itHits)->hit()->clone() );
106 
109  }
110  else if ( (*itSetup)->useForExternalTracking( *itHits ) )
111  {
112  externalTrackingRecHits.push_back( (*itHits)->hit()->clone() );
113  }
114  }
115 
116  //edm::LogInfo( "KalmanAlignmentTrackRefitter" ) << "Hits for tracking/external: " << trackingRecHits.size() << "/" << externalTrackingRecHits.size();
117 
118  //if ( !zPlusRecHits.size() || !zMinusRecHits.size() ) continue;
120 
121  if ( trackingRecHits.empty() ) continue;
122 
123  if ( externalTrackingRecHits.empty() )
124  {
125  if ( ( (*itSetup)->getExternalTrackingSubDetIds().size() == 0 ) && // O.K., no external hits expected,
126  ( trackingRecHits.size() >= (*itSetup)->minTrackingHits() ) )
127  {
128  TrajTrackPairCollection refitted = refitSingleTracklet( aGeometry.product(), aMagneticField.product(),
129  (*itSetup)->fitter(), (*itSetup)->propagator(),
130  aRecHitBuilder.product(), fullTrack,
131  trackingRecHits, beamSpot,
132  (*itSetup)->sortingDirection(), false, true );
133 
134  // The refitting did not work ... Try next!
135  if ( refitted.empty() ) continue;
136  if ( rejectTrack( refitted.front().second ) ) continue;
137 
138  if ( theDebugFlag )
139  {
140  debugTrackData( (*itSetup)->id(), refitted.front().first, refitted.front().second, beamSpot );
141  debugTrackData( "OrigFullTrack", (*itTrack).first, (*itTrack).second, beamSpot );
142  }
143 
144 
145  TrackletPtr trackletPtr( new KalmanAlignmentTracklet( refitted.front(), *itSetup ) );
146  result.push_back( trackletPtr );
147  }
148  else { continue; } // Expected external hits but found none or not enough hits.
149  }
150  else if ( ( trackingRecHits.size() >= (*itSetup)->minTrackingHits() ) &&
151  ( externalTrackingRecHits.size() >= (*itSetup)->minExternalHits() ) )
152  {
153  // Create an instance of KalmanAlignmentTracklet with an external prediction.
154 
155  TrajTrackPairCollection external = refitSingleTracklet( aGeometry.product(), aMagneticField.product(),
156  (*itSetup)->externalFitter(), (*itSetup)->externalPropagator(),
157  aRecHitBuilder.product(), fullTrack,
158  externalTrackingRecHits, beamSpot,
159  (*itSetup)->externalSortingDirection(),
160  false, true );
161  //if ( external.empty() || rejectTrack( external.front().second ) ) { continue; }
162  if ( external.empty() ) { continue; }
163 
164  TransientTrack externalTrack( *external.front().second, aMagneticField.product() );
165 
166  TrajTrackPairCollection refitted = refitSingleTracklet( aGeometry.product(), aMagneticField.product(),
167  (*itSetup)->fitter(), (*itSetup)->propagator(),
168  aRecHitBuilder.product(), externalTrack,
169  trackingRecHits, beamSpot,
170  (*itSetup)->sortingDirection(),
171  false, true, (*itSetup)->id() );
172 
173  if ( refitted.empty() ) { continue; }
174  if ( rejectTrack( refitted.front().second ) ) continue;
175 
176  //const Surface& surface = refitted.front().first->firstMeasurement().updatedState().surface();
177  const Surface& surface = refitted.front().first->lastMeasurement().updatedState().surface();
178  TrajectoryStateOnSurface externalTsos = externalTrack.impactPointState();
179  AnalyticalPropagator externalPredictionPropagator( aMagneticField.product(), anyDirection );
180  TrajectoryStateOnSurface externalPrediction = externalPredictionPropagator.propagate( externalTsos, surface );
181  if ( !externalPrediction.isValid() ) continue;
182 
183  if ( theDebugFlag )
184  {
185  debugTrackData( string("External") + (*itSetup)->id(), external.front().first, external.front().second, beamSpot );
186  debugTrackData( (*itSetup)->id(), refitted.front().first, refitted.front().second, beamSpot );
187  debugTrackData( "OrigFullTrack", (*itTrack).first, (*itTrack).second, beamSpot );
188  }
189 
190  TrackletPtr trackletPtr( new KalmanAlignmentTracklet( refitted.front(), externalPrediction, *itSetup ) );
191  result.push_back( trackletPtr );
192 
193  delete external.front().first;
194  delete external.front().second;
195  }
196  }
197  }
198 
199  return result;
200 }
201 
202 
205  const MagneticField* magneticField,
206  const TrajectoryFitter* fitter,
207  const Propagator* propagator,
208  const TransientTrackingRecHitBuilder* recHitBuilder,
209  const TransientTrack& fullTrack,
210  RecHitContainer& recHits,
211  const reco::BeamSpot* beamSpot,
212  const SortingDirection& sortingDir,
213  bool useExternalEstimate,
214  bool reuseMomentumEstimate,
215  const string identifier )
216 {
217 
219 
220  if ( recHits.size() < 2 ) return result;
221 
222  sortRecHits( recHits, recHitBuilder, sortingDir );
223 
225  RecHitContainer::iterator itRecHit;
226  for ( itRecHit = recHits.begin(); itRecHit != recHits.end(); ++itRecHit )
227  hits.push_back( recHitBuilder->build( &(*itRecHit) ) );
228 
229  TransientTrackingRecHit::ConstRecHitPointer firstHit = hits.front();
230 
231  AnalyticalPropagator firstStatePropagator( magneticField, anyDirection );
232  TrajectoryStateOnSurface firstState = firstStatePropagator.propagate( fullTrack.impactPointState(), firstHit->det()->surface() );
233 
234  KalmanAlignmentDataCollector::fillHistogram( identifier + string("_IPPt"),
235  1e-2*fullTrack.impactPointState().globalParameters().momentum().perp() );
236 
237  if ( !firstState.isValid() ) return result;
238 
239 // LocalTrajectoryError startError;
240 
241 // const double startErrorValue = 100;
242 // const unsigned int nTrajParam = 5;
243 
244 // if ( useExternalEstimate ) {
245 // startError = firstState.localError();
246 // } else {
247 // if ( reuseMomentumEstimate )
248 // {
249 // AlgebraicSymMatrix firstStateError( asHepMatrix( firstState.localError().matrix() ) );
250 // AlgebraicSymMatrix startErrorMatrix( nTrajParam, 0 );
251 // startErrorMatrix[0][0] = 1e-10;
252 // //startErrorMatrix[0][0] = firstStateError[0][0];
253 // startErrorMatrix[1][1] = startErrorValue;//firstStateError[1][1];
254 // startErrorMatrix[2][2] = startErrorValue;//firstStateError[2][2];
255 // startErrorMatrix[3][3] = startErrorValue;
256 // startErrorMatrix[4][4] = startErrorValue;
257 // startError = LocalTrajectoryError( startErrorMatrix );
258 // } else {
259 // AlgebraicSymMatrix startErrorMatrix( startErrorValue*AlgebraicSymMatrix( nTrajParam, 1 ) );
260 // startError = LocalTrajectoryError( startErrorMatrix );
261 // }
262 
263 // }
264 
265 // // MOMENTUM ESTIMATE FOR COSMICS. P = 1.5 GeV
266 // LocalTrajectoryParameters firstStateParameters = firstState.localParameters();
267 // AlgebraicVector firstStateParamVec = asHepVector( firstStateParameters.mixedFormatVector() );
268 // firstStateParamVec[0] = 1./1.5;
269 // LocalTrajectoryParameters cosmicsStateParameters( firstStateParamVec, firstStateParameters.pzSign(), true );
270 // TrajectoryStateOnSurface tsos( cosmicsStateParameters, startError, firstState.surface(), magneticField );
271 
272  KalmanAlignmentDataCollector::fillHistogram( identifier + string("_FSPt"),
273  1e-2*firstState.globalParameters().momentum().perp() );
274 
275  firstState.rescaleError( 100 );
276  TrajectoryStateOnSurface tsos( firstState.localParameters(), firstState.localError(),
277  firstState.surface(), magneticField );
278 
279  // Generate a trajectory seed.
280  TrajectorySeed seed( PTrajectoryStateOnDet(), recHits, propagator->propagationDirection() );
281 
282  // Generate track candidate.
283 
284  PTrajectoryStateOnDet state = trajectoryStateTransform::persistentState( tsos, firstHit->det()->geographicalId().rawId() );
285  TrackCandidate candidate( recHits, seed, state );
286 
287  AlgoProductCollection algoResult;
288 
289  int charge = static_cast<int>( tsos.charge() );
290  double momentum = firstState.localParameters().momentum().mag();
292  TRecHit1DMomConstraint::build( charge, momentum, 1e-10, &tsos.surface() );
293 
294  //no insert in OwnVector...
296  tmpHits.push_back(testhit);
297  for (TransientTrackingRecHit::RecHitContainer::const_iterator i=hits.begin(); i!=hits.end(); i++){
298  tmpHits.push_back(*i);
299  }
300  hits.swap(tmpHits);
301 
302  theRefitterAlgo.buildTrack( fitter, propagator, algoResult, hits, tsos, seed, 0, *beamSpot, candidate.seedRef());
303 
304  for ( AlgoProductCollection::iterator it = algoResult.begin(); it != algoResult.end(); ++it )
305  result.push_back( make_pair( (*it).first, (*it).second.first ) );
306 
307  return result;
308 }
309 
310 
312  const TransientTrackingRecHitBuilder* builder,
313  const SortingDirection& sortingDir ) const
314 {
315  // Don't start sorting if there is only 1 or even 0 elements.
316  if ( hits.size() < 2 ) return;
317 
318  TransientTrackingRecHit::RecHitPointer firstHit = builder->build( &hits.front() );
319  double firstRadius = firstHit->det()->surface().toGlobal( firstHit->localPosition() ).mag();
320  double firstY = firstHit->det()->surface().toGlobal( firstHit->localPosition() ).y();
321 
322  TransientTrackingRecHit::RecHitPointer lastHit = builder->build( &hits.back() );
323  double lastRadius = lastHit->det()->surface().toGlobal( lastHit->localPosition() ).mag();
324  double lastY = lastHit->det()->surface().toGlobal( lastHit->localPosition() ).y();
325 
326  bool insideOut = firstRadius < lastRadius;
327  bool upsideDown = lastY < firstY;
328 
329  if ( ( insideOut && ( sortingDir == KalmanAlignmentSetup::sortInsideOut ) ) ||
330  ( !insideOut && ( sortingDir == KalmanAlignmentSetup::sortOutsideIn ) ) ||
331  ( upsideDown && ( sortingDir == KalmanAlignmentSetup::sortUpsideDown ) ) ||
332  ( !upsideDown && ( sortingDir == KalmanAlignmentSetup::sortDownsideUp ) ) ) return;
333 
334  // Fill temporary container with reversed hits.
336  RecHitContainer::iterator itHit = hits.end();
337  do { --itHit; tmp.push_back( ( *itHit ).clone() ); } while ( itHit != hits.begin() );
338 
339  // Swap the content of the temporary and the input container.
340  hits.swap( tmp );
341 
342  return;
343 }
344 
345 
347 {
348  double trackChi2 = track->chi2();
349  unsigned int ndof = static_cast<unsigned int>( track->ndof() );
350  if ( trackChi2 <= 0. || ndof <= 0 ) return false;
351 
352  //FIXME: should be configurable (via KalmanAlignmentSetup)
353  double minChi2Prob = 0;//1e-6;
354  double maxChi2Prob = 1.0;
355 
356  GENFUNCTION cumulativeChi2 = Genfun::CumulativeChiSquare( ndof );
357  double chi2Prob = 1. - cumulativeChi2( trackChi2 );
358  return ( chi2Prob < minChi2Prob ) || ( chi2Prob > maxChi2Prob );
359 }
360 
361 
362 void KalmanAlignmentTrackRefitter::debugTrackData( const string identifier,
363  const Trajectory* traj,
364  const Track* track,
365  const reco::BeamSpot* bs )
366 {
367  unsigned int ndof = static_cast<unsigned int>( track->ndof() );
368  double trackChi2 = track->chi2();
369  if ( ( trackChi2 > 0. ) && ( ndof > 0 ) )
370  {
371  GENFUNCTION cumulativeChi2 = Genfun::CumulativeChiSquare( ndof );
372  KalmanAlignmentDataCollector::fillHistogram( identifier + string("_CumChi2"), 1. - cumulativeChi2( trackChi2 ) );
373  } else if ( ndof == 0 ) {
374  KalmanAlignmentDataCollector::fillHistogram( identifier + string("_CumChi2"), -1. );
375  } else {
376  KalmanAlignmentDataCollector::fillHistogram( identifier + string("_CumChi2"), -2. );
377  }
378 
379  KalmanAlignmentDataCollector::fillHistogram( identifier + string("_NHits"), traj->foundHits() );
380  KalmanAlignmentDataCollector::fillHistogram( identifier + string("_Pt"), 1e-2*track->pt() );
381  KalmanAlignmentDataCollector::fillHistogram( identifier + string("_Eta"), track->eta() );
382  KalmanAlignmentDataCollector::fillHistogram( identifier + string("_Phi"), track->phi() );
383  KalmanAlignmentDataCollector::fillHistogram( identifier + string("_NormChi2"), track->normalizedChi2() );
384  KalmanAlignmentDataCollector::fillHistogram( identifier + string("_DZ"), track->dz() );
385 
386  KalmanAlignmentDataCollector::fillHistogram( identifier + string("_DXY_BS"), fabs( track->dxy( bs->position() ) ) );
387  KalmanAlignmentDataCollector::fillHistogram( identifier + string("_DXY"), fabs( track->dxy() ) );
388  //KalmanAlignmentDataCollector::fillHistogram( identifier + string("_D0"), fabs( track->d0() ) );
389 }
int foundHits() const
Definition: Trajectory.h:234
int i
Definition: DBlmapReader.cc:9
TrackletCollection refitTracks(const edm::EventSetup &eventSetup, const AlignmentSetupCollection &algoSetups, const ConstTrajTrackPairCollection &tracks, const reco::BeamSpot *beamSpot)
reference back()
Definition: OwnVector.h:327
T perp() const
Definition: PV3DBase.h:72
static void fillHistogram(std::string histo_name, float data)
virtual PropagationDirection propagationDirection() const
Definition: Propagator.h:155
bool rejectTrack(const reco::Track *track) const
const LocalTrajectoryParameters & localParameters() const
double normalizedChi2() const
chi-squared divided by n.d.o.f. (or chi-squared * 1e6 if n.d.o.f. is zero)
Definition: TrackBase.h:537
std::vector< reco::Track > TrackCollection
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
size_type size() const
Definition: OwnVector.h:254
PTrajectoryStateOnDet persistentState(const TrajectoryStateOnSurface &ts, unsigned int detid)
std::vector< ConstRecHitPointer > RecHitContainer
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:621
std::vector< KalmanAlignmentSetup * > AlignmentSetupCollection
void setConf(const edm::ParameterSet &conf)
Set parameter set.
bool buildTrack(const TrajectoryFitter *, const Propagator *, AlgoProductCollection &, TransientTrackingRecHit::RecHitContainer &, TrajectoryStateOnSurface &, const TrajectorySeed &, float, const reco::BeamSpot &, SeedRef seedRef=SeedRef(), int qualityMask=0, signed char nLoops=0)
Construct Tracks to be put in the event.
iterator begin()
Definition: OwnVector.h:234
KalmanAlignmentTracklet::TrajTrackPairCollection TrajTrackPairCollection
void push_back(D *&d)
Definition: OwnVector.h:280
const SurfaceType & surface() const
T mag() const
Definition: PV3DBase.h:67
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:627
double chi2() const
chi-squared of the fit
Definition: TrackBase.h:525
std::shared_ptr< TrackingRecHit const > ConstRecHitPointer
double ndof() const
number of degrees of freedom of the fit
Definition: TrackBase.h:531
virtual RecHitPointer build(const TrackingRecHit *p) const =0
build a tracking rechit from an existing rechit
bool empty() const
Definition: OwnVector.h:259
double pt() const
track transverse momentum
Definition: TrackBase.h:597
tuple result
Definition: query.py:137
LocalVector momentum() const
Momentum vector in the local frame.
const LocalTrajectoryError & localError() const
std::vector< AlgoProduct > AlgoProductCollection
TrackingRecHit::ConstRecHitContainer ConstRecHitContainer
Definition: Trajectory.h:47
std::shared_ptr< TrackingRecHit const > RecHitPointer
void sortRecHits(RecHitContainer &hits, const TransientTrackingRecHitBuilder *builder, const SortingDirection &sortingDir) const
TrajTrackPairCollection refitSingleTracklet(const TrackingGeometry *geometry, const MagneticField *magneticField, const TrajectoryFitter *fitter, const Propagator *propagator, const TransientTrackingRecHitBuilder *recHitBuilder, const reco::TransientTrack &originalTrack, RecHitContainer &recHits, const reco::BeamSpot *beamSpot, const SortingDirection &sortingDir, bool useExternalEstimate, bool reuseMomentumEstimate, const std::string identifier=std::string("RefitSingle_"))
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:585
void debugTrackData(const std::string identifier, const Trajectory *traj, const reco::Track *track, const reco::BeamSpot *bs)
iterator end()
Definition: OwnVector.h:239
const GlobalTrajectoryParameters & globalParameters() const
tuple tracks
Definition: testEve_cfg.py:39
AlignmentAlgorithmBase::ConstTrajTrackPairCollection ConstTrajTrackPairCollection
T const * product() const
Definition: ESHandle.h:86
std::vector< TrackletPtr > TrackletCollection
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
ESHandle< TrackerGeometry > geometry
const Point & position() const
position
Definition: BeamSpot.h:62
TrajectoryStateOnSurface impactPointState() const
static RecHitPointer build(const int charge, const double mom, const double err, const Surface *surface)
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:567
void swap(OwnVector< T, P > &other)
Definition: OwnVector.h:407
AlignableDetOrUnitPtr alignableFromDetId(const DetId &detid)
Returns AlignableDetOrUnitPtr corresponding to given DetId.
reference front()
Definition: OwnVector.h:355
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
KalmanAlignmentTrackRefitter(const edm::ParameterSet &config, AlignableNavigator *navigator)
Constructor.
TrackProducerAlgorithm< reco::Track > theRefitterAlgo
virtual void getFromES(const edm::EventSetup &, edm::ESHandle< TrackerGeometry > &, edm::ESHandle< MagneticField > &, edm::ESHandle< TrajectoryFitter > &, edm::ESHandle< Propagator > &, edm::ESHandle< MeasurementTracker > &, edm::ESHandle< TransientTrackingRecHitBuilder > &)
Get needed services from the Event Setup.