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