CMS 3D CMS Logo

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