CMS 3D CMS Logo

MuonTrajectoryUpdator.cc
Go to the documentation of this file.
1 
21 
29 
32 
36 
37 #include <algorithm>
38 
39 using namespace edm;
40 using namespace std;
41 
44  NavigationDirection fitDirection): theFitDirection(fitDirection){
45 
46  // The max allowed chi2 to accept a rechit in the fit
47  theMaxChi2 = par.getParameter<double>("MaxChi2");
49 
50  // The KF updator
51  theUpdator= new KFUpdator();
52 
53  // The granularity
54  theGranularity = par.getParameter<int>("Granularity");
55 
56  // Rescale the error of the first state?
57  theRescaleErrorFlag = par.getParameter<bool>("RescaleError");
58 
59  if(theRescaleErrorFlag)
60  // The rescale factor
61  theRescaleFactor = par.getParameter<double>("RescaleErrorFactor");
62 
63  // Use the invalid hits?
64  useInvalidHits = par.getParameter<bool>("UseInvalidHits");
65 
66  // Flag needed for the rescaling
67  theFirstTSOSFlag = true;
68 
69  // Exlude RPC from the fit?
70  theRPCExFlag = par.getParameter<bool>("ExcludeRPCFromFit");
71 }
72 
74  double chi2, int granularity): theMaxChi2(chi2),
75  theGranularity(granularity),
76  theFitDirection(fitDirection){
78 
79  // The KF updator
80  theUpdator= new KFUpdator();
81 }
82 
85  delete theEstimator;
86  delete theUpdator;
87 }
88 
90  theFirstTSOSFlag = true;
91 }
92 
93 
94 pair<bool,TrajectoryStateOnSurface>
96  Trajectory& trajectory,
97  const Propagator *propagator){
98 
99  const std::string metname = "Muon|RecoMuon|MuonTrajectoryUpdator";
100 
101  MuonPatternRecoDumper muonDumper;
102 
103  // Status of the updating
104  bool updated=false;
105 
106  if(!measurement) return pair<bool,TrajectoryStateOnSurface>(updated,TrajectoryStateOnSurface() );
107 
108  // measurement layer
109  const DetLayer* detLayer=measurement->layer();
110 
111  // these are the 4D segment for the CSC/DT and a point for the RPC
112  TransientTrackingRecHit::ConstRecHitPointer muonRecHit = measurement->recHit();
113 
114  // The KFUpdator takes TransientTrackingRecHits as arg.
117 
118  // sort the container in agreement with the porpagation direction
119  sort(recHitsForFit,detLayer);
120 
121  TrajectoryStateOnSurface lastUpdatedTSOS = measurement->predictedState();
122 
123  LogTrace(metname)<<"Number of rechits for the fit: "<<recHitsForFit.size()<<endl;
124 
125  TransientTrackingRecHit::ConstRecHitContainer::iterator recHit;
126  for(recHit = recHitsForFit.begin(); recHit != recHitsForFit.end(); ++recHit ) {
127  if ((*recHit)->isValid() ) {
128 
129  // propagate the TSOS onto the rechit plane
130  TrajectoryStateOnSurface propagatedTSOS = propagateState(lastUpdatedTSOS, measurement,
131  *recHit, propagator);
132 
133  if ( propagatedTSOS.isValid() ) {
134  pair<bool,double> thisChi2 = estimator()->estimate(propagatedTSOS, *((*recHit).get()));
135 
136  LogTrace(metname) << "Estimation for Kalman Fit. Chi2: " << thisChi2.second;
137 
138  // Is an RPC hit? Prepare the variable to possibly exluding it from the fit
139  bool wantIncludeThisHit = true;
140  if (theRPCExFlag &&
141  (*recHit)->geographicalId().det() == DetId::Muon &&
142  (*recHit)->geographicalId().subdetId() == MuonSubdetId::RPC){
143  wantIncludeThisHit = false;
144  LogTrace(metname) << "This is an RPC hit and the present configuration is such that it will be excluded from the fit";
145  }
146 
147 
148  // The Chi2 cut was already applied in the estimator, which
149  // returns 0 if the chi2 is bigger than the cut defined in its
150  // constructor
151  if (thisChi2.first) {
152  updated=true;
153  if (wantIncludeThisHit) { // This split is a trick to have the RPC hits counted as updatable (in used chamber counting), while are not actually included in the fit when the proper obtion is activated.
154 
155  LogTrace(metname) << endl
156  << " Kalman Start" << "\n" << "\n";
157  LogTrace(metname) << " Meas. Position : " << (**recHit).globalPosition() << "\n"
158  << " Pred. Position : " << propagatedTSOS.globalPosition()
159  << " Pred Direction : " << propagatedTSOS.globalDirection()<< endl;
160 
162  propagatedTSOS.rescaleError(theRescaleFactor);
163  theFirstTSOSFlag = false;
164  }
165 
166  lastUpdatedTSOS = measurementUpdator()->update(propagatedTSOS,*((*recHit).get()));
167 
168  LogTrace(metname) << " Fit Position : " << lastUpdatedTSOS.globalPosition()
169  << " Fit Direction : " << lastUpdatedTSOS.globalDirection()
170  << "\n"
171  << " Fit position radius : "
172  << lastUpdatedTSOS.globalPosition().perp()
173  << "filter updated" << endl;
174 
175  LogTrace(metname) << muonDumper.dumpTSOS(lastUpdatedTSOS);
176 
177  LogTrace(metname) << "\n\n Kalman End" << "\n" << "\n";
178 
179  TrajectoryMeasurement && updatedMeasurement = updateMeasurement( propagatedTSOS, lastUpdatedTSOS,
180  *recHit, thisChi2.second, detLayer,
181  measurement);
182  // FIXME: check!
183  trajectory.push(std::move(updatedMeasurement), thisChi2.second);
184  }
185  else {
186  LogTrace(metname) << " Compatible RecHit with good chi2 but made with RPC when it was decided to not include it in the fit"
187  << " --> trajectory NOT updated, invalid RecHit added." << endl;
188 
189  MuonTransientTrackingRecHit::MuonRecHitPointer invalidRhPtr = MuonTransientTrackingRecHit::specificBuild( (*recHit)->det(), (*recHit)->hit() );
190  invalidRhPtr->invalidateHit();
191  TrajectoryMeasurement invalidRhMeasurement(propagatedTSOS, propagatedTSOS, invalidRhPtr, thisChi2.second, detLayer);
192  trajectory.push(std::move(invalidRhMeasurement), thisChi2.second);
193  }
194  }
195  else {
196  if(useInvalidHits) {
197  LogTrace(metname) << " Compatible RecHit with too large chi2"
198  << " --> trajectory NOT updated, invalid RecHit added." << endl;
199 
200  MuonTransientTrackingRecHit::MuonRecHitPointer invalidRhPtr = MuonTransientTrackingRecHit::specificBuild( (*recHit)->det(), (*recHit)->hit() );
201  invalidRhPtr->invalidateHit();
202  TrajectoryMeasurement invalidRhMeasurement(propagatedTSOS, propagatedTSOS, invalidRhPtr, thisChi2.second, detLayer);
203  trajectory.push(std::move(invalidRhMeasurement), thisChi2.second);
204  }
205  }
206  }
207  }
208  }
209  recHitsForFit.clear();
210  return pair<bool,TrajectoryStateOnSurface>(updated,lastUpdatedTSOS);
211 }
212 
215  const TrajectoryMeasurement* measurement,
217  const Propagator *propagator) const{
218 
219  const TransientTrackingRecHit::ConstRecHitPointer& mother = measurement->recHit();
220 
221  if( current->geographicalId() == mother->geographicalId() )
222  return measurement->predictedState();
223 
224  const TrajectoryStateOnSurface tsos =
225  propagator->propagate(state, current->det()->surface());
226  return tsos;
227 
228 }
229 
230 // FIXME: would I a different threatment for the two prop dirrections??
232  const TrajectoryStateOnSurface &lastUpdatedTSOS,
234  const double &chi2, const DetLayer *detLayer,
235  const TrajectoryMeasurement *initialMeasurement){
236  return TrajectoryMeasurement(propagatedTSOS, lastUpdatedTSOS,
237  recHit,chi2,detLayer);
238 
239  // // FIXME: put a better check! One could fit in first out-in and then in - out
240  // if(propagator()->propagationDirection() == alongMomentum)
241  // return TrajectoryMeasurement(propagatedTSOS, lastUpdatedTSOS,
242  // recHit,thisChi2.second,detLayer);
243 
244  // // FIXME: Check this carefully!!
245  // else if(propagator()->propagationDirection() == oppositeToMomentum)
246  // return TrajectoryMeasurement(initialMeasurement->forwardPredictedState(),
247  // propagatedTSOS, lastUpdatedTSOS,
248  // recHit,thisChi2.second,detLayer);
249  // else{
250  // LogError("MuonTrajectoryUpdator::updateMeasurement") <<"Wrong propagation direction!!";
251  // }
252 }
253 
254 
256  const DetLayer* detLayer){
257 
258  if(detLayer->subDetector()==GeomDetEnumerators::DT){
259  if(fitDirection() == insideOut)
260  stable_sort(recHitsForFit.begin(),recHitsForFit.end(), RadiusComparatorInOut() );
261  else if(fitDirection() == outsideIn)
262  stable_sort(recHitsForFit.begin(),recHitsForFit.end(),RadiusComparatorOutIn() );
263  else
264  LogError("Muon|RecoMuon|MuonTrajectoryUpdator") <<"MuonTrajectoryUpdator::sort: Wrong propagation direction!!";
265  }
266 
267  else if(detLayer->subDetector()==GeomDetEnumerators::CSC){
268  if(fitDirection() == insideOut)
269  stable_sort(recHitsForFit.begin(),recHitsForFit.end(), ZedComparatorInOut() );
270  else if(fitDirection() == outsideIn)
271  stable_sort(recHitsForFit.begin(),recHitsForFit.end(), ZedComparatorOutIn() );
272  else
273  LogError("Muon|RecoMuon|MuonTrajectoryUpdator") <<"MuonTrajectoryUpdator::sort: Wrong propagation direction!!";
274  }
275 
276  else if(detLayer->subDetector()==GeomDetEnumerators::GEM){
277  if(fitDirection() == insideOut)
278  stable_sort(recHitsForFit.begin(),recHitsForFit.end(), ZedComparatorInOut() );
279  else if(fitDirection() == outsideIn)
280  stable_sort(recHitsForFit.begin(),recHitsForFit.end(), ZedComparatorOutIn() );
281  else
282  LogError("Muon|RecoMuon|MuonTrajectoryUpdator") <<"MuonTrajectoryUpdator::sort: Wrong propagation direction!!";
283  }
284  else if(detLayer->subDetector()==GeomDetEnumerators::ME0){
285  if(fitDirection() == insideOut)
286  stable_sort(recHitsForFit.begin(),recHitsForFit.end(), ZedComparatorInOut() );
287  else if(fitDirection() == outsideIn)
288  stable_sort(recHitsForFit.begin(),recHitsForFit.end(), ZedComparatorOutIn() );
289  else
290  LogError("Muon|RecoMuon|MuonTrajectoryUpdator") <<"MuonTrajectoryUpdator::sort: Wrong propagation direction!!";
291  }
292 }
T getParameter(std::string const &) const
Ordering along decreasing radius (for DT rechits)
TrajectoryStateOnSurface const & predictedState() const
const MeasurementEstimator * estimator() const
accasso at the propagator
T perp() const
Definition: PV3DBase.h:72
ConstRecHitPointer const & recHit() const
const std::string metname
void makeFirstTime()
reset the theFirstTSOSFlag
Ordering along increasing zed (for CSC rechits)
virtual SubDetector subDetector() const =0
The type of detector (PixelBarrel, PixelEndcap, TIB, TOB, TID, TEC, CSC, DT, RPCBarrel, RPCEndcap)
NavigationDirection fitDirection()
get the fit direction
GlobalPoint globalPosition() const
virtual ~MuonTrajectoryUpdator()
Destructor.
TrajectoryStateOnSurface propagateState(const TrajectoryStateOnSurface &state, const TrajectoryMeasurement *measurement, const TransientTrackingRecHit::ConstRecHitPointer &current, const Propagator *propagator) const
MuonTrajectoryUpdator(const edm::ParameterSet &par, NavigationDirection fitDirection)
Constructor from Propagator and Parameter set.
virtual TrajectoryStateOnSurface update(const TrajectoryStateOnSurface &, const TrackingRecHit &) const =0
static TransientTrackingRecHit::ConstRecHitContainer breakInSubRecHits(TransientTrackingRecHit::ConstRecHitPointer, int granularity)
takes a muon rechit and returns its sub-rechits given a certain granularity
void sort(TransientTrackingRecHit::ConstRecHitContainer &, const DetLayer *)
std::string dumpTSOS(const TrajectoryStateOnSurface &tsos) const
const SurfaceType & surface() const
std::shared_ptr< MuonTransientTrackingRecHit > MuonRecHitPointer
std::shared_ptr< TrackingRecHit const > ConstRecHitPointer
TrajectoryMeasurement updateMeasurement(const TrajectoryStateOnSurface &propagatedTSOS, const TrajectoryStateOnSurface &lastUpdatedTSOS, const TransientTrackingRecHit::ConstRecHitPointer &recHit, const double &chi2, const DetLayer *detLayer, const TrajectoryMeasurement *initialMeasurement)
Return the trajectory measurement. It handles both the fw and the bw propagation. ...
double theMaxChi2
the max chi2 allowed
const DetLayer * layer() const
#define LogTrace(id)
std::vector< ConstRecHitPointer > ConstRecHitContainer
Ordering along increasing radius (for DT rechits)
const TrajectoryStateUpdator * measurementUpdator() const
TrajectoryStateUpdator * theUpdator
TrajectoryStateOnSurface propagate(STA const &state, SUR const &surface) const
Definition: Propagator.h:53
HLT enums.
static const int RPC
Definition: MuonSubdetId.h:14
virtual std::pair< bool, TrajectoryStateOnSurface > update(const TrajectoryMeasurement *measurement, Trajectory &trajectory, const Propagator *propagator)
update the Trajectory with the TrajectoryMeasurement
virtual HitReturnType estimate(const TrajectoryStateOnSurface &ts, const TrackingRecHit &hit) const =0
NavigationDirection theFitDirection
static MuonRecHitPointer specificBuild(const GeomDet *geom, const TrackingRecHit *rh)
Ordering along decreasing zed (for CSC rechits)
void push(const TrajectoryMeasurement &tm)
Definition: Trajectory.cc:50
def move(src, dest)
Definition: eostools.py:510
GlobalVector globalDirection() const
MeasurementEstimator * theEstimator