00001
00017 #include "RecoMuon/TrackingTools/interface/MuonTrajectoryUpdator.h"
00018 #include "RecoMuon/TrackingTools/interface/MuonPatternRecoDumper.h"
00019
00020 #include "TrackingTools/KalmanUpdators/interface/Chi2MeasurementEstimator.h"
00021 #include "TrackingTools/PatternTools/interface/TrajectoryMeasurement.h"
00022 #include "TrackingTools/PatternTools/interface/Trajectory.h"
00023 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
00024 #include "TrackingTools/GeomPropagators/interface/Propagator.h"
00025 #include "TrackingTools/KalmanUpdators/interface/KFUpdator.h"
00026 #include "TrackingTools/DetLayers/interface/DetLayer.h"
00027
00028 #include "RecoMuon/TransientTrackingRecHit/interface/MuonTransientTrackingRecHit.h"
00029 #include "RecoMuon/TransientTrackingRecHit/interface/MuonTransientTrackingRecHitBreaker.h"
00030
00031 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00032 #include "FWCore/Utilities/interface/Exception.h"
00033 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00034
00035 #include <algorithm>
00036
00037 using namespace edm;
00038 using namespace std;
00039
00041 MuonTrajectoryUpdator::MuonTrajectoryUpdator(const edm::ParameterSet& par,
00042 NavigationDirection fitDirection): theFitDirection(fitDirection){
00043
00044
00045 theMaxChi2 = par.getParameter<double>("MaxChi2");
00046 theEstimator = new Chi2MeasurementEstimator(theMaxChi2);
00047
00048
00049 theUpdator= new KFUpdator();
00050
00051
00052 theGranularity = par.getParameter<int>("Granularity");
00053
00054
00055 theRescaleErrorFlag = par.getParameter<bool>("RescaleError");
00056
00057 if(theRescaleErrorFlag)
00058
00059 theRescaleFactor = par.getParameter<double>("RescaleErrorFactor");
00060
00061
00062 theFirstTSOSFlag = true;
00063 }
00064
00065 MuonTrajectoryUpdator::MuonTrajectoryUpdator( NavigationDirection fitDirection,
00066 double chi2, int granularity): theMaxChi2(chi2),
00067 theGranularity(granularity),
00068 theFitDirection(fitDirection){
00069 theEstimator = new Chi2MeasurementEstimator(theMaxChi2);
00070
00071
00072 theUpdator= new KFUpdator();
00073 }
00074
00076 MuonTrajectoryUpdator::~MuonTrajectoryUpdator(){
00077 delete theEstimator;
00078 delete theUpdator;
00079 }
00080
00081 void MuonTrajectoryUpdator::makeFirstTime(){
00082 theFirstTSOSFlag = true;
00083 }
00084
00085
00086 pair<bool,TrajectoryStateOnSurface>
00087 MuonTrajectoryUpdator::update(const TrajectoryMeasurement* measurement,
00088 Trajectory& trajectory,
00089 const Propagator *propagator){
00090
00091 const std::string metname = "Muon|RecoMuon|MuonTrajectoryUpdator";
00092
00093 MuonPatternRecoDumper muonDumper;
00094
00095
00096 bool updated=false;
00097
00098 if(!measurement) return pair<bool,TrajectoryStateOnSurface>(updated,TrajectoryStateOnSurface() );
00099
00100
00101 const DetLayer* detLayer=measurement->layer();
00102
00103
00104 TransientTrackingRecHit::ConstRecHitPointer muonRecHit = measurement->recHit();
00105
00106
00107 TransientTrackingRecHit::ConstRecHitContainer recHitsForFit =
00108 MuonTransientTrackingRecHitBreaker::breakInSubRecHits(muonRecHit,theGranularity);
00109
00110
00111 sort(recHitsForFit,detLayer);
00112
00113 TrajectoryStateOnSurface lastUpdatedTSOS = measurement->predictedState();
00114
00115 LogTrace(metname)<<"Number of rechits for the fit: "<<recHitsForFit.size()<<endl;
00116
00117 TransientTrackingRecHit::ConstRecHitContainer::iterator recHit;
00118 for(recHit = recHitsForFit.begin(); recHit != recHitsForFit.end(); ++recHit ) {
00119 if ((*recHit)->isValid() ) {
00120
00121
00122 TrajectoryStateOnSurface propagatedTSOS = propagateState(lastUpdatedTSOS, measurement,
00123 *recHit, propagator);
00124
00125 if ( propagatedTSOS.isValid() ) {
00126 pair<bool,double> thisChi2 = estimator()->estimate(propagatedTSOS, *((*recHit).get()));
00127
00128 LogTrace(metname) << "Estimation for Kalman Fit. Chi2: " << thisChi2.second;
00129
00130
00131
00132
00133 if ( thisChi2.first ) {
00134 updated=true;
00135
00136 LogTrace(metname) << endl
00137 << " Kalman Start" << "\n" << "\n";
00138 LogTrace(metname) << " Meas. Position : " << (**recHit).globalPosition() << "\n"
00139 << " Pred. Position : " << propagatedTSOS.globalPosition()
00140 << " Pred Direction : " << propagatedTSOS.globalDirection()<< endl;
00141
00142 if(theRescaleErrorFlag && theFirstTSOSFlag){
00143 propagatedTSOS.rescaleError(theRescaleFactor);
00144 theFirstTSOSFlag = false;
00145 }
00146
00147 lastUpdatedTSOS = measurementUpdator()->update(propagatedTSOS,*((*recHit).get()));
00148
00149 LogTrace(metname) << " Fit Position : " << lastUpdatedTSOS.globalPosition()
00150 << " Fit Direction : " << lastUpdatedTSOS.globalDirection()
00151 << "\n"
00152 << " Fit position radius : "
00153 << lastUpdatedTSOS.globalPosition().perp()
00154 << "filter updated" << endl;
00155
00156 LogTrace(metname) << muonDumper.dumpTSOS(lastUpdatedTSOS);
00157
00158 LogTrace(metname) << "\n\n Kalman End" << "\n" << "\n";
00159
00160 TrajectoryMeasurement updatedMeasurement = updateMeasurement( propagatedTSOS, lastUpdatedTSOS,
00161 *recHit,thisChi2.second,detLayer,
00162 measurement);
00163
00164 trajectory.push(updatedMeasurement, thisChi2.second);
00165 }
00166 }
00167 }
00168 }
00169 recHitsForFit.clear();
00170 return pair<bool,TrajectoryStateOnSurface>(updated,lastUpdatedTSOS);
00171 }
00172
00173 TrajectoryStateOnSurface
00174 MuonTrajectoryUpdator::propagateState(const TrajectoryStateOnSurface& state,
00175 const TrajectoryMeasurement* measurement,
00176 const TransientTrackingRecHit::ConstRecHitPointer & current,
00177 const Propagator *propagator) const{
00178
00179 const TransientTrackingRecHit::ConstRecHitPointer mother = measurement->recHit();
00180
00181 if( current->geographicalId() == mother->geographicalId() )
00182 return measurement->predictedState();
00183
00184 const TrajectoryStateOnSurface tsos =
00185 propagator->propagate(state, current->det()->surface());
00186 return tsos;
00187
00188 }
00189
00190
00191 TrajectoryMeasurement MuonTrajectoryUpdator::updateMeasurement( const TrajectoryStateOnSurface &propagatedTSOS,
00192 const TrajectoryStateOnSurface &lastUpdatedTSOS,
00193 const TransientTrackingRecHit::ConstRecHitPointer &recHit,
00194 const double &chi2, const DetLayer *detLayer,
00195 const TrajectoryMeasurement *initialMeasurement){
00196 return TrajectoryMeasurement(propagatedTSOS, lastUpdatedTSOS,
00197 recHit,chi2,detLayer);
00198
00199
00200
00201
00202
00203
00204
00205
00206
00207
00208
00209
00210
00211
00212 }
00213
00214
00215 void MuonTrajectoryUpdator::sort(TransientTrackingRecHit::ConstRecHitContainer& recHitsForFit,
00216 const DetLayer* detLayer){
00217
00218 if(detLayer->subDetector()==GeomDetEnumerators::DT){
00219 if(fitDirection() == insideOut)
00220 stable_sort(recHitsForFit.begin(),recHitsForFit.end(), RadiusComparatorInOut() );
00221 else if(fitDirection() == outsideIn)
00222 stable_sort(recHitsForFit.begin(),recHitsForFit.end(),RadiusComparatorOutIn() );
00223 else
00224 LogError("Muon|RecoMuon|MuonTrajectoryUpdator") <<"MuonTrajectoryUpdator::sort: Wrong propagation direction!!";
00225 }
00226
00227 else if(detLayer->subDetector()==GeomDetEnumerators::CSC){
00228 if(fitDirection() == insideOut)
00229 stable_sort(recHitsForFit.begin(),recHitsForFit.end(), ZedComparatorInOut() );
00230 else if(fitDirection() == outsideIn)
00231 stable_sort(recHitsForFit.begin(),recHitsForFit.end(), ZedComparatorOutIn() );
00232 else
00233 LogError("Muon|RecoMuon|MuonTrajectoryUpdator") <<"MuonTrajectoryUpdator::sort: Wrong propagation direction!!";
00234 }
00235 }