CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10/src/MuonAnalysis/MuonAssociators/src/MatcherUsingTracksAlgorithm.cc

Go to the documentation of this file.
00001 #include "MuonAnalysis/MuonAssociators/interface/MatcherUsingTracksAlgorithm.h"
00002 
00003 #include "DataFormats/Math/interface/deltaR.h"
00004 #include "DataFormats/Math/interface/deltaPhi.h"
00005 
00006 #include "DataFormats/RecoCandidate/interface/RecoCandidate.h"
00007 #include "DataFormats/TrackReco/interface/Track.h"
00008 
00009 #include "TrackingTools/Records/interface/TrackingComponentsRecord.h"
00010 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00011 #include "Geometry/Records/interface/GlobalTrackingGeometryRecord.h"
00012 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
00013 #include "TrackingTools/PatternTools/interface/TSCPBuilderNoMaterial.h"
00014 #include "TrackingTools/TrajectoryState/interface/PerigeeConversions.h"
00015 
00016 #include "FWCore/Utilities/interface/Exception.h"
00017 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00018 
00019 // a few template-related workarounds
00020 namespace reco {
00021 template<>
00022 inline double deltaR<GlobalVector>(const GlobalVector &v1, const GlobalVector &v2) {
00023     return deltaR<float>(v1.eta(),v1.phi(),v2.eta(),v2.phi());
00024 }
00025 template<>
00026 inline double deltaR2<GlobalVector>(const GlobalVector &v1, const GlobalVector &v2) {
00027     return deltaR2<float>(v1.eta(),v1.phi(),v2.eta(),v2.phi());
00028 }
00029 }
00030 
00031 MatcherUsingTracksAlgorithm::MatcherUsingTracksAlgorithm(const edm::ParameterSet & iConfig) :
00032     whichTrack1_(None),     whichTrack2_(None), 
00033     whichState1_(AtVertex), whichState2_(AtVertex),
00034     srcCut_(iConfig.existsAs<std::string>("srcPreselection") ? iConfig.getParameter<std::string>("srcPreselection") : ""),
00035     matchedCut_(iConfig.existsAs<std::string>("matchedPreselection") ? iConfig.getParameter<std::string>("matchedPreselection") : ""),
00036     requireSameCharge_(iConfig.existsAs<bool>("requireSameCharge") ? iConfig.getParameter<bool>("requireSameCharge") : false)
00037 {
00038     std::string algo = iConfig.getParameter<std::string>("algorithm");
00039     if      (algo == "byTrackRef")           { algo_ = ByTrackRef; }
00040     else if (algo == "byPropagatingSrc")     { algo_ = ByPropagatingSrc; }
00041     else if (algo == "byPropagatingMatched") { algo_ = ByPropagatingMatched; }
00042     else if (algo == "byDirectComparison")   { algo_ = ByDirectComparison; }
00043     else throw cms::Exception("Configuration") << "Value '" << algo << "' for algorithm not yet implemented.\n";
00044 
00045     getConf(iConfig, "src",     whichTrack1_, whichState1_);
00046     getConf(iConfig, "matched", whichTrack2_, whichState2_);
00047 
00048     if (algo_ == ByTrackRef) {
00049         // validate the config
00050         if (whichTrack1_ == None || whichTrack2_ == None) throw cms::Exception("Configuration") << "Algorithm 'byTrackRef' needs tracks not to be 'none'.\n";
00051     } else if (algo_ == ByPropagatingSrc || algo_ == ByPropagatingMatched || algo_ == ByDirectComparison) {
00052         // read matching cuts
00053         maxLocalPosDiff_      = iConfig.getParameter<double>("maxDeltaLocalPos");
00054         maxGlobalMomDeltaR_   = iConfig.getParameter<double>("maxDeltaR");
00055         maxGlobalMomDeltaEta_ = iConfig.existsAs<double>("maxDeltaEta") ? iConfig.getParameter<double>("maxDeltaEta") : maxGlobalMomDeltaR_;
00056         maxGlobalMomDeltaPhi_ = iConfig.existsAs<double>("maxDeltaPhi") ? iConfig.getParameter<double>("maxDeltaPhi") : maxGlobalMomDeltaR_;
00057         maxGlobalDPtRel_      = iConfig.getParameter<double>("maxDeltaPtRel");
00058 
00059         // choice of sorting variable
00060         std::string sortBy = iConfig.getParameter<std::string>("sortBy");
00061         if      (sortBy == "deltaLocalPos") sortBy_ = LocalPosDiff;
00062         else if (sortBy == "deltaPtRel")    sortBy_ = GlobalDPtRel;
00063         else if (sortBy == "deltaR")        sortBy_ = GlobalMomDeltaR;
00064         else if (sortBy == "deltaEta")      sortBy_ = GlobalMomDeltaEta;
00065         else if (sortBy == "deltaPhi")      sortBy_ = GlobalMomDeltaPhi;
00066         else if (sortBy == "chi2")          sortBy_ = Chi2;
00067         else throw cms::Exception("Configuration") << "Parameter 'sortBy' must be one of: deltaLocalPos, deltaPtRel, deltaR, chi2.\n";
00068         // validate the config
00069         if (algo_ == ByPropagatingSrc) {
00070             if (whichTrack2_ == None || whichState2_ == AtVertex) {
00071                 algo_ = ByPropagatingSrcTSCP;
00072                 //throw cms::Exception("Configuration") << "Destination track must be non-null, and state must not be 'AtVertex' (not yet).\n";
00073             }
00074         } else if (algo_ == ByPropagatingMatched) {
00075             if (whichTrack1_ == None || whichState1_ == AtVertex) {
00076                 algo_ = ByPropagatingMatchedTSCP;
00077                 //throw cms::Exception("Configuration") << "Destination track must be non-null, and state must not be 'AtVertex' (not yet).\n";
00078             }
00079         } else if (algo_ == ByDirectComparison) {
00080             bool firstAtVertex = (whichTrack1_ == None || whichState1_ == AtVertex);
00081             bool   secAtVertex = (whichTrack2_ == None || whichState2_ == AtVertex);
00082             if (firstAtVertex) {
00083                 if (!secAtVertex)  throw cms::Exception("Configuration") << "When using 'byDirectComparison' with 'src' at vertex (or None), 'matched' must be at vertex too.\n";
00084             } else {
00085                 if (secAtVertex)  throw cms::Exception("Configuration") << "When using 'byDirectComparison' with 'src' not at vertex, 'matched' can't be at vertex or None.\n";
00086                 if (whichState1_ != whichState2_) throw cms::Exception("Configuration") << "You can't use 'byDirectComparison' with non-matching states.\n";
00087             }
00088         }
00089 
00090         useChi2_ = iConfig.existsAs<bool>("computeChi2") ? iConfig.getParameter<bool>("computeChi2") : false;
00091         if (useChi2_) {
00092            if (whichTrack1_ == None && whichTrack2_ == None) throw cms::Exception("Configuration") << "Can't compute chi2s if both tracks are set to 'none'.\n";
00093            maxChi2_ = iConfig.getParameter<double>("maxChi2");
00094            chi2DiagonalOnly_ = iConfig.getParameter<bool>("chi2DiagonalOnly");
00095            if (algo_ == ByPropagatingSrc || algo_ == ByPropagatingMatched) {
00096                chi2UseVertex_ = iConfig.getParameter<bool>("chi2UsePosition");
00097            } else {
00098                chi2UseVertex_ = iConfig.getParameter<bool>("chi2UseVertex");
00099                if (algo_ == ByDirectComparison) {
00100                     std::string choice = iConfig.getParameter<std::string>("chi2MomentumForDxy");
00101                     if (choice == "src") chi2FirstMomentum_ = true;
00102                     else if (choice != "matched") throw cms::Exception("Configuration") << "chi2MomentumForDxy must be 'src' or 'matched'\n";
00103                }
00104            }
00105         } else maxChi2_ = 1;
00106 
00107         if (sortBy_ == Chi2 && !useChi2_) throw cms::Exception("Configuration") << "Can't sort by chi2s if 'computeChi2s' is not set to true.\n";
00108     }
00109 }
00110 
00111 void 
00112 MatcherUsingTracksAlgorithm::getConf(const edm::ParameterSet & iConfig, const std::string &whatFor, WhichTrack &whichTrack, WhichState &whichState) {
00113     std::string s_whichTrack = iConfig.getParameter<std::string>(whatFor+"Track");
00114     if      (s_whichTrack == "none")    { whichTrack = None; }
00115     else if (s_whichTrack == "tracker") { whichTrack = TrackerTk; }
00116     else if (s_whichTrack == "muon")    { whichTrack = MuonTk; }
00117     else if (s_whichTrack == "global")  { whichTrack = GlobalTk; }
00118     else throw cms::Exception("Configuration") << "Parameter 'useTrack' must be 'none', 'tracker', 'muon', 'global'\n";
00119     if ((whichTrack != None) && (algo_ != ByTrackRef)) {
00120         std::string s_whichState = iConfig.getParameter<std::string>(whatFor+"State");
00121         if      (s_whichState == "atVertex")  { whichState = AtVertex; }
00122         else if (s_whichState == "innermost") { whichState = Innermost; }
00123         else if (s_whichState == "outermost") { whichState = Outermost; }
00124         else throw cms::Exception("Configuration") << "Parameter 'useState' must be 'atVertex', 'innermost', 'outermost'\n";
00125     }
00126 }
00127 
00130 bool 
00131 MatcherUsingTracksAlgorithm::match(const reco::Candidate &c1, const reco::Candidate &c2, float &deltR, float &deltEta, float &deltPhi, float &deltaLocalPos, float &deltaPtRel, float &chi2) const {
00132     if (!(srcCut_(c1) && matchedCut_(c2))) return false;
00133     if (requireSameCharge_ && (c1.charge() != c2.charge())) return false;
00134     switch (algo_) {
00135         case ByTrackRef: { 
00136             reco::TrackRef t1 = getTrack(c1, whichTrack1_); 
00137             reco::TrackRef t2 = getTrack(c2, whichTrack2_); 
00138             if (t1.isNonnull()) {
00139                 if (t1 == t2) return true; 
00140                 if (t1.id() != t2.id()) { edm::LogWarning("MatcherUsingTracksAlgorithm") << "Trying to match by reference tracks coming from different collections.\n"; }
00141             } 
00142             }
00143         case ByPropagatingSrc: {
00144             FreeTrajectoryState start = startingState(c1, whichTrack1_, whichState1_);
00145             TrajectoryStateOnSurface target = targetState(c2, whichTrack2_, whichState2_);
00146             return matchWithPropagation(start, target, deltR, deltEta, deltPhi, deltaLocalPos, deltaPtRel, chi2);
00147             }
00148         case ByPropagatingMatched: { 
00149             FreeTrajectoryState start = startingState(c2, whichTrack2_, whichState2_);
00150             TrajectoryStateOnSurface target = targetState(c1, whichTrack1_, whichState1_);
00151             return matchWithPropagation(start, target, deltR, deltEta, deltPhi, deltaLocalPos, deltaPtRel, chi2);
00152             }
00153         case ByPropagatingSrcTSCP: {
00154             FreeTrajectoryState start = startingState(c1, whichTrack1_, whichState1_);
00155             FreeTrajectoryState target = startingState(c2, whichTrack2_, whichState2_);
00156             return matchWithPropagation(start, target, deltR, deltEta, deltPhi, deltaLocalPos, deltaPtRel, chi2);
00157             }
00158         case ByPropagatingMatchedTSCP: { 
00159             FreeTrajectoryState start = startingState(c2, whichTrack2_, whichState2_);
00160             FreeTrajectoryState target = startingState(c1, whichTrack1_, whichState1_);
00161             return matchWithPropagation(start, target, deltR, deltEta, deltPhi, deltaLocalPos, deltaPtRel, chi2);
00162             }
00163 
00164         case ByDirectComparison: {
00165             FreeTrajectoryState start = startingState(c1, whichTrack1_, whichState1_);
00166             FreeTrajectoryState otherstart = startingState(c2, whichTrack2_, whichState2_);
00167             return matchByDirectComparison(start, otherstart, deltR, deltEta, deltPhi, deltaLocalPos, deltaPtRel, chi2);
00168             }
00169     }
00170     return false;
00171 }
00172 
00176 int 
00177 MatcherUsingTracksAlgorithm::match(const reco::Candidate &c1, const edm::View<reco::Candidate> &c2s, float &deltR, float &deltEta, float &deltPhi, float &deltaLocalPos, float &deltaPtRel, float &chi2) const {
00178     if (!srcCut_(c1)) return -1;
00179     
00180     // working and output variables
00181     FreeTrajectoryState start; TrajectoryStateOnSurface target;
00182     int match = -1;
00183 
00184     // pre-fetch some states if needed
00185     if (algo_ == ByPropagatingSrc || algo_ == ByPropagatingSrcTSCP || algo_ == ByPropagatingMatchedTSCP || algo_ == ByDirectComparison) {
00186         start = startingState(c1, whichTrack1_, whichState1_);
00187     } else if (algo_ == ByPropagatingMatched) target = targetState(c1, whichTrack1_, whichState1_);
00188 
00189     // loop on the  collection
00190     edm::View<reco::Candidate>::const_iterator it, ed; int i;
00191     for (it = c2s.begin(), ed = c2s.end(), i = 0; it != ed; ++it, ++i) {
00192         if (!matchedCut_(*it)) continue;
00193         if (requireSameCharge_ && (c1.charge() != it->charge()) ) continue;
00194         bool exit = false;
00195         switch (algo_) {
00196             case ByTrackRef: {
00197                      reco::TrackRef t1 = getTrack(c1, whichTrack1_); 
00198                      reco::TrackRef t2 = getTrack(*it, whichTrack2_); 
00199                      if (t1.isNonnull()) {
00200                          if (t1 == t2) { match = i; exit = true; }
00201                          if (t1.id() != t2.id()) { edm::LogWarning("MatcherUsingTracksAlgorithm") << "Trying to match by reference tracks coming from different collections.\n"; }
00202                      } 
00203                  } break;
00204             case ByPropagatingSrc: 
00205             case ByPropagatingMatched: {
00206                      if (algo_ == ByPropagatingMatched)  start = startingState(*it, whichTrack2_, whichState2_);
00207                      else if (algo_ == ByPropagatingSrc) target = targetState(*it, whichTrack2_, whichState2_);
00208                      if (matchWithPropagation(start, target, deltR, deltEta, deltPhi, deltaLocalPos, deltaPtRel, chi2)) { 
00209                          match = i; 
00210                      }
00211                  } break;
00212             case ByDirectComparison: {
00213                      FreeTrajectoryState otherstart = startingState(*it, whichTrack2_, whichState2_);
00214                      if (matchByDirectComparison(start, otherstart, deltR, deltEta, deltPhi, deltaLocalPos, deltaPtRel, chi2)) {
00215                          match = i;
00216                      }
00217                  } break;
00218             case ByPropagatingSrcTSCP: {
00219                     FreeTrajectoryState otherstart = startingState(*it, whichTrack2_, whichState2_);
00220                     if (matchWithPropagation(start, otherstart, deltR, deltEta, deltPhi, deltaLocalPos, deltaPtRel, chi2)) {
00221                          match = i;
00222                     }
00223                 } break;
00224             case ByPropagatingMatchedTSCP: {
00225                     FreeTrajectoryState otherstart = startingState(*it, whichTrack2_, whichState2_);
00226                     if (matchWithPropagation(otherstart, start, deltR, deltEta, deltPhi, deltaLocalPos, deltaPtRel, chi2)) {
00227                          match = i;
00228                     }
00229                 } break;
00230         }
00231         if (exit) break;
00232     }
00233 
00234     return match;
00235 }
00236 
00237 void 
00238 MatcherUsingTracksAlgorithm::init(const edm::EventSetup & iSetup) {
00239     iSetup.get<IdealMagneticFieldRecord>().get(magfield_);
00240     iSetup.get<TrackingComponentsRecord>().get("SteppingHelixPropagatorAny", propagator_);
00241     iSetup.get<GlobalTrackingGeometryRecord>().get(geometry_);
00242 }
00243 
00244 
00245 reco::TrackRef
00246 MatcherUsingTracksAlgorithm::getTrack(const reco::Candidate &reco, WhichTrack whichTrack) const {
00247     reco::TrackRef tk;
00248     const reco::RecoCandidate *rc = dynamic_cast<const reco::RecoCandidate *>(&reco);
00249     if (rc == 0) throw cms::Exception("Invalid Data") << "Input object is not a RecoCandidate.\n";
00250     switch (whichTrack) {
00251         case TrackerTk: tk = rc->track();          break; 
00252         case MuonTk   : tk = rc->standAloneMuon(); break;
00253         case GlobalTk : tk = rc->combinedMuon();   break;
00254         default: break; // just to make gcc happy
00255     }
00256     return tk;
00257 }
00258 
00259 FreeTrajectoryState 
00260 MatcherUsingTracksAlgorithm::startingState(const reco::Candidate &reco, WhichTrack whichTrack, WhichState whichState) const {
00261     FreeTrajectoryState ret;
00262     if (whichTrack != None) {
00263         reco::TrackRef tk = getTrack(reco, whichTrack);
00264         if (tk.isNull()) {
00265             ret = FreeTrajectoryState();
00266         } else {
00267             switch (whichState) {
00268                 case AtVertex:  ret = trajectoryStateTransform::initialFreeState(*tk, magfield_.product()); break;
00269                 case Innermost: ret = trajectoryStateTransform::innerFreeState(  *tk, magfield_.product()); break;
00270                 case Outermost: ret = trajectoryStateTransform::outerFreeState(  *tk, magfield_.product()); break;
00271             }
00272         }
00273     } else {
00274         ret = FreeTrajectoryState(  GlobalPoint( reco.vx(), reco.vy(), reco.vz()),
00275                                     GlobalVector(reco.px(), reco.py(), reco.pz()),
00276                                     reco.charge(),
00277                                     magfield_.product());
00278     }
00279     return ret;
00280 }
00281 
00282 TrajectoryStateOnSurface 
00283 MatcherUsingTracksAlgorithm::targetState(const reco::Candidate &reco, WhichTrack whichTrack, WhichState whichState) const {
00284     TrajectoryStateOnSurface ret;
00285     reco::TrackRef tk = getTrack(reco, whichTrack);
00286     if (tk.isNonnull()) {
00287         switch (whichState) {
00288             case Innermost: ret = trajectoryStateTransform::innerStateOnSurface(  *tk, *geometry_, magfield_.product()); break;
00289             case Outermost: ret = trajectoryStateTransform::outerStateOnSurface(  *tk, *geometry_, magfield_.product()); break;
00290             default: break; // just to make gcc happy
00291         }
00292     }
00293     return ret;
00294 }
00295 
00296 bool
00297 MatcherUsingTracksAlgorithm::matchWithPropagation(const FreeTrajectoryState &start, 
00298                                               const TrajectoryStateOnSurface &target, 
00299                                               float &lastDeltaR,
00300                                               float &lastDeltaEta,
00301                                               float &lastDeltaPhi, 
00302                                               float &lastDeltaLocalPos,
00303                                               float &lastGlobalDPtRel,
00304                                               float &lastChi2) const 
00305 {
00306     if ((start.momentum().mag() == 0) || !target.isValid()) return false;
00307 
00308     TrajectoryStateOnSurface tsos = propagator_->propagate(start, target.surface());
00309 
00310     bool isBest = false;
00311     if (tsos.isValid()) {
00312         float thisLocalPosDiff = (tsos.localPosition()-target.localPosition()).mag();
00313         float thisGlobalMomDeltaR = deltaR(tsos.globalMomentum(), target.globalMomentum());
00314         float thisGlobalMomDeltaPhi = fabs(deltaPhi(tsos.globalMomentum().phi(), target.globalMomentum().phi()));
00315         float thisGlobalMomDeltaEta = fabs(tsos.globalMomentum().eta() - target.globalMomentum().eta());
00316         float thisGlobalDPtRel = (tsos.globalMomentum().perp() - target.globalMomentum().perp())/target.globalMomentum().perp();
00317 
00318         if ((thisLocalPosDiff       < maxLocalPosDiff_) &&
00319             (thisGlobalMomDeltaR    < maxGlobalMomDeltaR_) &&
00320             (thisGlobalMomDeltaEta  < maxGlobalMomDeltaEta_) &&
00321             (thisGlobalMomDeltaPhi  < maxGlobalMomDeltaPhi_) &&
00322             (fabs(thisGlobalDPtRel) < maxGlobalDPtRel_)) {
00323             float thisChi2         = useChi2_ ? getChi2(target, tsos, chi2DiagonalOnly_, chi2UseVertex_) : 0;
00324             if (thisChi2 >= maxChi2_) return false;
00325             switch (sortBy_) {
00326                 case LocalPosDiff:      isBest = (thisLocalPosDiff      < lastDeltaLocalPos); break;
00327                 case GlobalMomDeltaR:   isBest = (thisGlobalMomDeltaR   < lastDeltaR);        break;
00328                 case GlobalMomDeltaEta: isBest = (thisGlobalMomDeltaEta < lastDeltaEta);      break;
00329                 case GlobalMomDeltaPhi: isBest = (thisGlobalMomDeltaPhi < lastDeltaPhi);      break;
00330                 case GlobalDPtRel:      isBest = (thisGlobalDPtRel      < lastGlobalDPtRel);  break;
00331                 case Chi2:              isBest = (thisChi2              < lastChi2);          break;
00332             }
00333             if (isBest) {
00334                 lastDeltaLocalPos = thisLocalPosDiff;
00335                 lastDeltaR        = thisGlobalMomDeltaR;
00336                 lastDeltaEta      = thisGlobalMomDeltaEta;
00337                 lastDeltaPhi      = thisGlobalMomDeltaPhi;
00338                 lastGlobalDPtRel  = thisGlobalDPtRel;
00339                 lastChi2          = thisChi2;
00340             } 
00341         }  // if match
00342     }
00343 
00344     return isBest;
00345 }
00346 
00347 bool
00348 MatcherUsingTracksAlgorithm::matchWithPropagation(const FreeTrajectoryState &start, 
00349                                               const FreeTrajectoryState &target, 
00350                                               float &lastDeltaR,
00351                                               float &lastDeltaEta,
00352                                               float &lastDeltaPhi, 
00353                                               float &lastDeltaLocalPos,
00354                                               float &lastGlobalDPtRel,
00355                                               float &lastChi2) const 
00356 {
00357     if ((start.momentum().mag() == 0) || (target.momentum().mag() == 0)) return false;
00358     TSCPBuilderNoMaterial propagator;
00359     /*2.2.X*/ try { 
00360     TrajectoryStateClosestToPoint tscp = propagator(start, target.position());
00361     // if (!tscp.isValid()) return false;  // in 3.1.X
00362     
00363     bool isBest = false;
00364     float thisLocalPosDiff = (tscp.position()-target.position()).mag();
00365     float thisGlobalMomDeltaR   = deltaR(tscp.momentum(), target.momentum());
00366     float thisGlobalMomDeltaPhi = fabs(deltaPhi(tscp.momentum().phi(), target.momentum().phi()));
00367     float thisGlobalMomDeltaEta = fabs(tscp.momentum().eta() - target.momentum().eta());
00368     float thisGlobalDPtRel = (tscp.momentum().perp() - target.momentum().perp())/target.momentum().perp();
00369 
00370     if ((thisLocalPosDiff       < maxLocalPosDiff_) &&
00371             (thisGlobalMomDeltaR    < maxGlobalMomDeltaR_) &&
00372             (thisGlobalMomDeltaEta  < maxGlobalMomDeltaEta_) &&
00373             (thisGlobalMomDeltaPhi  < maxGlobalMomDeltaPhi_) &&
00374             (fabs(thisGlobalDPtRel) < maxGlobalDPtRel_)) {
00375         float thisChi2         = useChi2_ ? getChi2(target, tscp, chi2DiagonalOnly_, chi2UseVertex_) : 0;
00376         if (thisChi2 >= maxChi2_) return false;
00377         switch (sortBy_) {
00378             case LocalPosDiff:      isBest = (thisLocalPosDiff      < lastDeltaLocalPos); break;
00379             case GlobalMomDeltaR:   isBest = (thisGlobalMomDeltaR   < lastDeltaR);        break;
00380             case GlobalMomDeltaEta: isBest = (thisGlobalMomDeltaEta < lastDeltaEta);      break;
00381             case GlobalMomDeltaPhi: isBest = (thisGlobalMomDeltaPhi < lastDeltaPhi);      break;
00382             case GlobalDPtRel:      isBest = (thisGlobalDPtRel      < lastGlobalDPtRel);  break;
00383             case Chi2:              isBest = (thisChi2              < lastChi2);          break;
00384         }
00385         if (isBest) {
00386             lastDeltaLocalPos = thisLocalPosDiff;
00387             lastDeltaR        = thisGlobalMomDeltaR;
00388             lastDeltaEta      = thisGlobalMomDeltaEta;
00389             lastDeltaPhi      = thisGlobalMomDeltaPhi;
00390             lastGlobalDPtRel  = thisGlobalDPtRel;
00391             lastChi2          = thisChi2;
00392         } 
00393     }  // if match
00394 
00395     return isBest;
00396     /*2.2.X*/ } catch (const TrajectoryStateException &err) { return false; } 
00397 }
00398 
00399 
00400 bool
00401 MatcherUsingTracksAlgorithm::matchByDirectComparison(const FreeTrajectoryState &start, 
00402                                                  const FreeTrajectoryState &target, 
00403                                                  float &lastDeltaR,
00404                                                  float &lastDeltaEta,
00405                                                  float &lastDeltaPhi, 
00406                                                  float &lastDeltaLocalPos,
00407                                                  float &lastGlobalDPtRel,
00408                                                  float &lastChi2) const 
00409 {
00410     if ((start.momentum().mag() == 0) || target.momentum().mag() == 0) return false;
00411 
00412     bool isBest = false;
00413     float thisLocalPosDiff = (start.position()-target.position()).mag();
00414     float thisGlobalMomDeltaR   = deltaR(start.momentum(), target.momentum());
00415     float thisGlobalMomDeltaPhi = fabs(deltaPhi(start.momentum().phi(), target.momentum().phi()));
00416     float thisGlobalMomDeltaEta = fabs(start.momentum().eta() - target.momentum().eta());
00417     float thisGlobalDPtRel = (start.momentum().perp() - target.momentum().perp())/target.momentum().perp();
00418 
00419     if ((thisLocalPosDiff       < maxLocalPosDiff_) &&
00420             (thisGlobalMomDeltaR    < maxGlobalMomDeltaR_) &&
00421             (thisGlobalMomDeltaEta  < maxGlobalMomDeltaEta_) &&
00422             (thisGlobalMomDeltaPhi  < maxGlobalMomDeltaPhi_) &&
00423             (fabs(thisGlobalDPtRel) < maxGlobalDPtRel_)) {
00424         float thisChi2 = useChi2_ ? getChi2(start, target, chi2DiagonalOnly_, chi2UseVertex_, chi2FirstMomentum_) : 0;
00425         if (thisChi2 >= maxChi2_) return false;
00426         switch (sortBy_) {
00427             case LocalPosDiff:      isBest = (thisLocalPosDiff      < lastDeltaLocalPos); break;
00428             case GlobalMomDeltaR:   isBest = (thisGlobalMomDeltaR   < lastDeltaR);        break;
00429             case GlobalMomDeltaEta: isBest = (thisGlobalMomDeltaEta < lastDeltaEta);      break;
00430             case GlobalMomDeltaPhi: isBest = (thisGlobalMomDeltaPhi < lastDeltaPhi);      break;
00431             case GlobalDPtRel:      isBest = (thisGlobalDPtRel      < lastGlobalDPtRel);  break;
00432             case Chi2:              isBest = (thisChi2              < lastChi2);          break;
00433         }
00434         if (isBest) {
00435             lastDeltaLocalPos = thisLocalPosDiff;
00436             lastDeltaR        = thisGlobalMomDeltaR;
00437             lastDeltaEta      = thisGlobalMomDeltaEta;
00438             lastDeltaPhi      = thisGlobalMomDeltaPhi;
00439             lastGlobalDPtRel  = thisGlobalDPtRel;
00440             lastChi2          = thisChi2;
00441         } 
00442     }  // if match
00443 
00444     return isBest;
00445 }
00446 
00447 double
00448 MatcherUsingTracksAlgorithm::getChi2(const FreeTrajectoryState &start, const FreeTrajectoryState &other, bool diagonalOnly, bool useVertex, bool useFirstMomentum) {
00449     if (!start.hasError() && !other.hasError()) throw cms::Exception("LogicError") << "At least one of the two states must have errors to make chi2s.\n";
00450     AlgebraicSymMatrix55 cov;
00451     if (start.hasError()) cov += start.curvilinearError().matrix();
00452     if (other.hasError()) cov += other.curvilinearError().matrix();
00453     cropAndInvert(cov, diagonalOnly, !useVertex);
00454     GlobalVector p1 = start.momentum(), p2 = other.momentum();
00455     GlobalPoint  x1 = start.position(), x2 = other.position();
00456     GlobalVector p = useFirstMomentum ? p1 : p2; double pt = p.perp(), pm = p.mag();
00457     double dsz = (x1.z()-x2.z())*pt/pm - ( (x1.x()-x2.x())*p.x() + (x1.y()-x2.y())*p.y() ) / pt * p.z()/pm;
00458     double dxy = ( -(x1.x()-x2.x())*p.y() + (x1.y()-x2.y())*p.x() )/pt;
00459     AlgebraicVector5 diff(start.charge()/p1.mag()-other.charge()/p2.mag(),
00460                           p1.theta()-p2.theta(),
00461                           (p1.phi()-p2.phi()).value(), 
00462                           dxy,
00463                           dsz);
00464     return ROOT::Math::Similarity(diff, cov);
00465 }
00466 
00467 double
00468 MatcherUsingTracksAlgorithm::getChi2(const FreeTrajectoryState &start, const TrajectoryStateClosestToPoint &other, bool diagonalOnly, bool useVertex) {
00469     if (!start.hasError() && !other.hasError()) throw cms::Exception("LogicError") << "At least one of the two states must have errors to make chi2s.\n";
00470     PerigeeConversions pgconvert; double pt; // needed by pgconvert
00471     AlgebraicSymMatrix55 cov;
00472     if (start.hasError()) cov += pgconvert.ftsToPerigeeError(start).covarianceMatrix();
00473     if (other.hasError()) cov += other.perigeeError().covarianceMatrix();
00474     cropAndInvert(cov, diagonalOnly, !useVertex);
00475     AlgebraicVector5 pgpar1 = pgconvert.ftsToPerigeeParameters(start,other.referencePoint(),pt).vector();
00476     AlgebraicVector5 pgpar2 = other.perigeeParameters().vector();
00477     AlgebraicVector5 diff(pgpar1-pgpar2);
00478     return ROOT::Math::Similarity(diff, cov);
00479 }
00480 
00481 double
00482 MatcherUsingTracksAlgorithm::getChi2(const TrajectoryStateOnSurface &start, const TrajectoryStateOnSurface &other, bool diagonalOnly, bool usePosition) {
00483     if (!start.hasError() && !other.hasError()) throw cms::Exception("LogicError") << "At least one of the two states must have errors to make chi2s.\n";
00484     AlgebraicSymMatrix55 cov;
00485     if (start.hasError()) cov += start.localError().matrix();
00486     if (other.hasError()) cov += other.localError().matrix();
00487     cropAndInvert(cov, diagonalOnly, !usePosition);
00488     AlgebraicVector5 diff(start.localParameters().mixedFormatVector() - other.localParameters().mixedFormatVector());
00489     return ROOT::Math::Similarity(diff, cov);
00490 }
00491 
00492 void
00493 MatcherUsingTracksAlgorithm::cropAndInvert(AlgebraicSymMatrix55 &cov, bool diagonalOnly, bool top3by3only) {
00494     if (!top3by3only) {
00495         if (diagonalOnly) {
00496             for (size_t i = 0; i < 5; ++i) { for (size_t j = i+1; j < 5; ++j) {
00497                 cov(i,j) = 0;
00498             } }
00499         }
00500         cov.Invert();
00501     } else {
00502         // get 3x3 covariance
00503         AlgebraicSymMatrix33 momCov = cov.Sub<AlgebraicSymMatrix33>(0,0); // get 3x3 matrix
00504         if (diagonalOnly) { momCov(0,1) = 0; momCov(0,2) = 0; momCov(1,2) = 0; }
00505         // invert        
00506         momCov.Invert();
00507         // place it
00508         cov.Place_at(momCov,0,0);
00509         // zero the rest
00510         for (size_t i = 3; i < 5; ++i) { for (size_t j = i; j < 5; ++j) { 
00511             cov(i,j) = 0;
00512         } }
00513     }
00514 }
00515 
00516