00001 #include "DataFormats/MuonReco/interface/MuonCocktails.h"
00002 #include "CommonTools/Statistics/interface/ChiSquaredProbability.h"
00003 #include "DataFormats/TrackReco/interface/Track.h"
00004
00005 #include <TROOT.h>
00006
00007
00008
00009
00010 reco::TrackRef muon::tevOptimized( const reco::TrackRef& combinedTrack,
00011 const reco::TrackRef& trackerTrack,
00012 const reco::TrackToTrackMap tevMap1,
00013 const reco::TrackToTrackMap tevMap2,
00014 const reco::TrackToTrackMap tevMap3 ) {
00015
00016 std::vector<reco::TrackRef> refit(4);
00017 bool ok[4];
00018 ok[0] = true;
00019
00020 reco::TrackToTrackMap::const_iterator gmrTrack = tevMap1.find(combinedTrack);
00021 reco::TrackToTrackMap::const_iterator fmsTrack = tevMap2.find(combinedTrack);
00022 reco::TrackToTrackMap::const_iterator pmrTrack = tevMap3.find(combinedTrack);
00023
00024 ok[1] = gmrTrack != tevMap1.end();
00025 ok[2] = fmsTrack != tevMap2.end();
00026 ok[3] = pmrTrack != tevMap3.end();
00027
00028 double prob[4];
00029 int chosen=3;
00030
00031 if (ok[0]) refit[0] = trackerTrack;
00032 if (ok[1]) refit[1] = (*gmrTrack).val;
00033 if (ok[2]) refit[2] = (*fmsTrack).val;
00034 if (ok[3]) refit[3] = (*pmrTrack).val;
00035
00036 for (unsigned int i=0; i<4; i++)
00037 prob[i] = (ok[i] && refit[i]->numberOfValidHits())
00038 ? trackProbability(refit[i]) : 0.0;
00039
00040
00041
00042 if (!prob[3])
00043 if (prob[2]) chosen=2; else
00044 if (prob[1]) chosen=1; else
00045 if (prob[0]) chosen=0;
00046
00047 if ( prob[0] && prob[3] && ((prob[3]-prob[0]) > 48.) ) chosen=0;
00048 if ( prob[0] && prob[1] && ((prob[1]-prob[0]) < 3.) ) chosen=1;
00049 if ( prob[2] && ((prob[chosen]-prob[2]) > 9.) ) chosen=2;
00050
00051 return refit.at(chosen);
00052 }
00053
00054
00055
00056
00057 reco::TrackRef muon::tevOptimizedOld( const reco::TrackRef& combinedTrack,
00058 const reco::TrackRef& trackerTrack,
00059 const reco::TrackToTrackMap tevMap1,
00060 const reco::TrackToTrackMap tevMap2,
00061 const reco::TrackToTrackMap tevMap3 ) {
00062
00063 std::vector<reco::TrackRef> refit(4);
00064 reco::TrackRef result;
00065 bool ok[4];
00066 ok[0] = true;
00067
00068 reco::TrackToTrackMap::const_iterator gmrTrack = tevMap1.find(combinedTrack);
00069 reco::TrackToTrackMap::const_iterator fmsTrack = tevMap2.find(combinedTrack);
00070 reco::TrackToTrackMap::const_iterator pmrTrack = tevMap3.find(combinedTrack);
00071
00072 ok[1] = gmrTrack != tevMap1.end();
00073 ok[2] = fmsTrack != tevMap2.end();
00074 ok[3] = pmrTrack != tevMap3.end();
00075
00076 double prob[4];
00077
00078 if (ok[0]) refit[0] = trackerTrack;
00079 if (ok[1]) refit[1] = (*gmrTrack).val;
00080 if (ok[2]) refit[2] = (*fmsTrack).val;
00081 if (ok[3]) refit[3] = (*pmrTrack).val;
00082
00083 for (unsigned int i=0; i<4; i++)
00084 prob[i] = (ok[i] && refit[i]->numberOfValidHits())
00085 ? trackProbability(refit[i]) : 0.0;
00086
00087
00088
00089 if (prob[1] ) result = refit[1];
00090 if ((prob[1] == 0) && prob[3]) result = refit[3];
00091
00092 if (prob[1] && prob[3] && ((prob[1] - prob[3]) > 0.05 )) result = refit[3];
00093
00094 if (prob[0] && prob[2] && fabs(prob[2] - prob[0]) > 30.) {
00095 result = refit[0];
00096 return result;
00097 }
00098
00099 if ((prob[1] == 0) && (prob[3] == 0) && prob[2]) result = refit[2];
00100
00101 reco::TrackRef tmin;
00102 double probmin = 0.0;
00103
00104 if (prob[1] && prob[3]) {
00105 probmin = prob[3]; tmin = refit[3];
00106 if ( prob[1] < prob[3] ) { probmin = prob[1]; tmin = refit[1]; }
00107 } else if ((prob[3] == 0) && prob[1]) {
00108 probmin = prob[1]; tmin = refit[1];
00109 } else if ((prob[1] == 0) && prob[3]) {
00110 probmin = prob[3]; tmin = refit[3];
00111 }
00112
00113 if (probmin && prob[2] && ( (probmin - prob[2]) > 3.5 )) {
00114 result = refit[2];
00115 }
00116
00117 return result;
00118 }
00119
00120
00121
00122
00123 reco::TrackRef muon::tevOptimized( const reco::Muon& muon,
00124 const reco::TrackToTrackMap tevMap1,
00125 const reco::TrackToTrackMap tevMap2,
00126 const reco::TrackToTrackMap tevMap3 ) {
00127 return muon::tevOptimized(muon.combinedMuon(), muon.track(),
00128 tevMap1, tevMap2, tevMap3);
00129 }
00130
00131
00132
00133
00134 reco::TrackRef muon::tevOptimizedOld( const reco::Muon& muon,
00135 const reco::TrackToTrackMap tevMap1,
00136 const reco::TrackToTrackMap tevMap2,
00137 const reco::TrackToTrackMap tevMap3 ) {
00138 return muon::tevOptimizedOld(muon.combinedMuon(), muon.track(),
00139 tevMap1, tevMap2, tevMap3);
00140 }
00141
00142
00143
00144
00145 double muon::trackProbability(const reco::TrackRef track) {
00146
00147 int nDOF = (int)track->ndof();
00148 if ( nDOF > 0 && track->chi2()> 0) {
00149 return -LnChiSquaredProbability(track->chi2(), nDOF);
00150 } else {
00151 return 0.0;
00152 }
00153
00154 }
00155
00156