CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/SimTracker/TrackAssociation/src/TrackAssociatorByChi2.cc

Go to the documentation of this file.
00001 #include "SimTracker/TrackAssociation/interface/TrackAssociatorByChi2.h"
00002 #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticle.h"
00003 
00004 #include "DataFormats/GeometrySurface/interface/Line.h"
00005 #include "DataFormats/GeometryVector/interface/Pi.h"
00006 #include "TrackingTools/TrajectoryState/interface/FreeTrajectoryState.h"
00007 #include "TrackingTools/PatternTools/interface/TSCBLBuilderNoMaterial.h"
00008 using namespace edm;
00009 using namespace reco;
00010 using namespace std;
00011 
00012 double TrackAssociatorByChi2::compareTracksParam ( TrackCollection::const_iterator rt, 
00013                                                    SimTrackContainer::const_iterator st, 
00014                                                    const math::XYZTLorentzVectorD vertexPosition, 
00015                                                    GlobalVector magField,
00016                                                    TrackBase::CovarianceMatrix invertedCovariance,
00017                                                    const reco::BeamSpot& bs) const{
00018   
00019   Basic3DVector<double> momAtVtx(st->momentum().x(),st->momentum().y(),st->momentum().z());
00020   Basic3DVector<double> vert = (Basic3DVector<double>) vertexPosition;
00021 
00022   std::pair<bool,reco::TrackBase::ParameterVector> params = parametersAtClosestApproach(vert, momAtVtx, st->charge(), bs);
00023   if (params.first){
00024     TrackBase::ParameterVector sParameters = params.second;
00025     TrackBase::ParameterVector rParameters = rt->parameters();
00026 
00027     TrackBase::ParameterVector diffParameters = rParameters - sParameters;
00028     double chi2 = ROOT::Math::Dot(diffParameters * invertedCovariance, diffParameters);
00029     
00030     return chi2;
00031   } else {
00032     return 10000000000.;
00033   }
00034 }
00035 
00036 
00037 TrackAssociatorByChi2::RecoToSimPairAssociation 
00038 TrackAssociatorByChi2::compareTracksParam(const TrackCollection& rtColl,
00039                                           const SimTrackContainer& stColl,
00040                                           const SimVertexContainer& svColl,
00041                                           const reco::BeamSpot& bs) const{
00042   
00043   RecoToSimPairAssociation outputVec;
00044 
00045   for (TrackCollection::const_iterator track=rtColl.begin(); track!=rtColl.end(); track++){
00046      Chi2SimMap outMap;
00047 
00048     TrackBase::ParameterVector rParameters = track->parameters();
00049 
00050     TrackBase::CovarianceMatrix recoTrackCovMatrix = track->covariance();
00051     if (onlyDiagonal){
00052       for (unsigned int i=0;i<5;i++){
00053         for (unsigned int j=0;j<5;j++){
00054           if (i!=j) recoTrackCovMatrix(i,j)=0;
00055         }
00056       }
00057     }
00058     recoTrackCovMatrix.Invert();
00059 
00060     for (SimTrackContainer::const_iterator st=stColl.begin(); st!=stColl.end(); st++){
00061 
00062       Basic3DVector<double> momAtVtx(st->momentum().x(),st->momentum().y(),st->momentum().z());
00063       Basic3DVector<double> vert = (Basic3DVector<double>)  svColl[st->vertIndex()].position();
00064 
00065       std::pair<bool,reco::TrackBase::ParameterVector> params = parametersAtClosestApproach(vert, momAtVtx, st->charge(), bs);
00066       if (params.first){
00067         TrackBase::ParameterVector sParameters = params.second;
00068       
00069         TrackBase::ParameterVector diffParameters = rParameters - sParameters;
00070         double chi2 = ROOT::Math::Dot(diffParameters * recoTrackCovMatrix, diffParameters);
00071         chi2/=5;
00072         if (chi2<chi2cut) outMap[chi2]=*st;
00073       }
00074     }
00075     outputVec.push_back(RecoToSimPair(*track,outMap));
00076   }
00077   return outputVec;
00078 }
00079 
00080 
00081 double TrackAssociatorByChi2::associateRecoToSim( TrackCollection::const_iterator rt, 
00082                                                   TrackingParticleCollection::const_iterator tp, 
00083                                                   const reco::BeamSpot& bs) const{
00084   
00085   double chi2;
00086   
00087   TrackBase::ParameterVector rParameters = rt->parameters();
00088 
00089   TrackBase::CovarianceMatrix recoTrackCovMatrix = rt->covariance();
00090       if (onlyDiagonal) {
00091         for (unsigned int i=0;i<5;i++){
00092           for (unsigned int j=0;j<5;j++){
00093             if (i!=j) recoTrackCovMatrix(i,j)=0;
00094           }
00095         }
00096       }
00097   recoTrackCovMatrix.Invert();
00098   
00099   Basic3DVector<double> momAtVtx(tp->momentum().x(),tp->momentum().y(),tp->momentum().z());
00100   Basic3DVector<double> vert(tp->vertex().x(),tp->vertex().y(),tp->vertex().z());
00101 
00102   std::pair<bool,reco::TrackBase::ParameterVector> params = parametersAtClosestApproach(vert, momAtVtx, tp->charge(), bs);
00103   if (params.first){
00104     TrackBase::ParameterVector sParameters=params.second;
00105     
00106     TrackBase::ParameterVector diffParameters = rParameters - sParameters;
00107     chi2 = ROOT::Math::Dot(diffParameters * recoTrackCovMatrix, diffParameters);
00108     chi2 /= 5;
00109     
00110     LogDebug("TrackAssociator") << "====NEW RECO TRACK WITH PT=" << rt->pt() << "====\n" 
00111                                 << "qoverp sim: " << sParameters[0] << "\n" 
00112                                 << "lambda sim: " << sParameters[1] << "\n" 
00113                                 << "phi    sim: " << sParameters[2] << "\n" 
00114                                 << "dxy    sim: " << sParameters[3] << "\n" 
00115                                 << "dsz    sim: " << sParameters[4] << "\n" 
00116                                 << ": " /*<< */ << "\n" 
00117                                 << "qoverp rec: " << rt->qoverp()/*rParameters[0]*/ << "\n" 
00118                                 << "lambda rec: " << rt->lambda()/*rParameters[1]*/ << "\n" 
00119                                 << "phi    rec: " << rt->phi()/*rParameters[2]*/ << "\n" 
00120                                 << "dxy    rec: " << rt->dxy(bs.position())/*rParameters[3]*/ << "\n" 
00121                                 << "dsz    rec: " << rt->dsz(bs.position())/*rParameters[4]*/ << "\n" 
00122                                 << ": " /*<< */ << "\n" 
00123                                 << "chi2: " << chi2 << "\n";
00124     
00125     return chi2;  
00126   } else {
00127     return 10000000000.;
00128   }
00129 }
00130 
00131 pair<bool,TrackBase::ParameterVector> 
00132 TrackAssociatorByChi2::parametersAtClosestApproach(Basic3DVector<double> vertex,
00133                                                    Basic3DVector<double> momAtVtx,
00134                                                    float charge,
00135                                                    const BeamSpot& bs) const{
00136   
00137   TrackBase::ParameterVector sParameters;
00138   try {
00139     FreeTrajectoryState ftsAtProduction(GlobalPoint(vertex.x(),vertex.y(),vertex.z()),
00140                                         GlobalVector(momAtVtx.x(),momAtVtx.y(),momAtVtx.z()),
00141                                         TrackCharge(charge),
00142                                         theMF.product());
00143     TSCBLBuilderNoMaterial tscblBuilder;
00144     TrajectoryStateClosestToBeamLine tsAtClosestApproach = tscblBuilder(ftsAtProduction,bs);//as in TrackProducerAlgorithm
00145     
00146     GlobalPoint v = tsAtClosestApproach.trackStateAtPCA().position();
00147     GlobalVector p = tsAtClosestApproach.trackStateAtPCA().momentum();
00148     sParameters[0] = tsAtClosestApproach.trackStateAtPCA().charge()/p.mag();
00149     sParameters[1] = Geom::halfPi() - p.theta();
00150     sParameters[2] = p.phi();
00151     sParameters[3] = (-v.x()*sin(p.phi())+v.y()*cos(p.phi()));
00152     sParameters[4] = v.z()*p.perp()/p.mag() - (v.x()*p.x()+v.y()*p.y())/p.perp() * p.z()/p.mag();
00153     
00154     return pair<bool,TrackBase::ParameterVector>(true,sParameters);
00155   } catch ( ... ) {
00156     return pair<bool,TrackBase::ParameterVector>(false,sParameters);
00157   }
00158 }
00159 
00160 RecoToSimCollection TrackAssociatorByChi2::associateRecoToSim(const edm::RefToBaseVector<reco::Track>& tC, 
00161                                                               const edm::RefVector<TrackingParticleCollection>& tPCH,
00162                                                               const edm::Event * e,
00163                                                               const edm::EventSetup *setup ) const{
00164   edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
00165   e->getByLabel(bsSrc,recoBeamSpotHandle);
00166   reco::BeamSpot bs = *recoBeamSpotHandle;      
00167 
00168   RecoToSimCollection  outputCollection;
00169   double chi2;
00170 
00171   TrackingParticleCollection tPC;
00172   if (tPCH.size()!=0)  tPC = *tPCH.product();
00173 
00174   int tindex=0;
00175   for (RefToBaseVector<reco::Track>::const_iterator rt=tC.begin(); rt!=tC.end(); rt++, tindex++){
00176 
00177     LogDebug("TrackAssociator") << "=========LOOKING FOR ASSOCIATION===========" << "\n"
00178                                 << "rec::Track #"<<tindex<<" with pt=" << (*rt)->pt() <<  "\n"
00179                                 << "===========================================" << "\n";
00180  
00181     TrackBase::ParameterVector rParameters = (*rt)->parameters();
00182 
00183     TrackBase::CovarianceMatrix recoTrackCovMatrix = (*rt)->covariance();
00184     if (onlyDiagonal){
00185       for (unsigned int i=0;i<5;i++){
00186         for (unsigned int j=0;j<5;j++){
00187           if (i!=j) recoTrackCovMatrix(i,j)=0;
00188         }
00189       }
00190     } 
00191 
00192     recoTrackCovMatrix.Invert();
00193 
00194     int tpindex =0;
00195     for (TrackingParticleCollection::const_iterator tp=tPC.begin(); tp!=tPC.end(); tp++, ++tpindex){
00196         
00197       //skip tps with a very small pt
00198       //if (sqrt(tp->momentum().perp2())<0.5) continue;
00199       Basic3DVector<double> momAtVtx(tp->momentum().x(),tp->momentum().y(),tp->momentum().z());
00200       Basic3DVector<double> vert=(Basic3DVector<double>) tp->vertex();
00201 
00202       std::pair<bool,reco::TrackBase::ParameterVector> params = parametersAtClosestApproach(vert, momAtVtx, tp->charge(), bs);
00203       if (params.first){
00204         TrackBase::ParameterVector sParameters=params.second;
00205         
00206         TrackBase::ParameterVector diffParameters = rParameters - sParameters;
00207         
00208         chi2 = ROOT::Math::Similarity(diffParameters, recoTrackCovMatrix);
00209         chi2 /= 5;
00210         
00211         LogTrace("TrackAssociator") << "====TP index=" << tpindex << "  RECO TRACK index="<<tindex<<" WITH PT=" << (*rt)->pt() << "====\n" 
00212                                     << "qoverp sim: " << sParameters[0] << "\n" 
00213                                     << "lambda sim: " << sParameters[1] << "\n" 
00214                                     << "phi    sim: " << sParameters[2] << "\n" 
00215                                     << "dxy    sim: " << sParameters[3] << "\n" 
00216                                     << "dsz    sim: " << sParameters[4] << "\n" 
00217                                     << ": " /*<< */ << "\n" 
00218                                     << "qoverp rec: " << (*rt)->qoverp() << " err: " << (*rt)->error(0) << "\n"
00219                                     << "lambda rec: " << (*rt)->lambda() << " err: " << (*rt)->error(1) << "\n"
00220                                     << "phi    rec: " << (*rt)->phi() << " err: " << (*rt)->error(2) << "\n"
00221                                     << "dxy    rec: " << (*rt)->dxy(bs.position()) << " err: " << (*rt)->error(3) << "\n"
00222                                     << "dsz    rec: " << (*rt)->dsz(bs.position()) << " err: " << (*rt)->error(4) << "\n"
00223                                     << ": " /*<< */ << "\n" 
00224                                     << "chi2: " << chi2 << "\n";
00225         
00226         if (chi2<chi2cut) {
00227           outputCollection.insert(tC[tindex], 
00228                                   std::make_pair(edm::Ref<TrackingParticleCollection>(tPCH, tpindex),
00229                                                  -chi2));//-chi2 because the Association Map is ordered using std::greater
00230         }
00231       }
00232     }
00233   }
00234   outputCollection.post_insert();
00235   return outputCollection;
00236 }
00237 
00238 
00239 
00240 SimToRecoCollection TrackAssociatorByChi2::associateSimToReco(const edm::RefToBaseVector<reco::Track>& tC, 
00241                                                               const edm::RefVector<TrackingParticleCollection>& tPCH,
00242                                                               const edm::Event * e,
00243                                                               const edm::EventSetup *setup ) const {
00244   edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
00245   e->getByLabel(bsSrc,recoBeamSpotHandle);
00246   reco::BeamSpot bs = *recoBeamSpotHandle;      
00247 
00248   SimToRecoCollection  outputCollection;
00249   double chi2;
00250 
00251   TrackingParticleCollection tPC;
00252   if (tPCH.size()!=0)  tPC = *tPCH.product();
00253 
00254   int tpindex =0;
00255   for (TrackingParticleCollection::const_iterator tp=tPC.begin(); tp!=tPC.end(); tp++, ++tpindex){
00256     
00257     //skip tps with a very small pt
00258     //if (sqrt(tp->momentum().perp2())<0.5) continue;
00259     
00260     LogDebug("TrackAssociator") << "=========LOOKING FOR ASSOCIATION===========" << "\n"
00261                                 << "TrackingParticle #"<<tpindex<<" with pt=" << sqrt(tp->momentum().perp2()) << "\n"
00262                                 << "===========================================" << "\n";
00263     
00264     Basic3DVector<double> momAtVtx(tp->momentum().x(),tp->momentum().y(),tp->momentum().z());
00265     Basic3DVector<double> vert(tp->vertex().x(),tp->vertex().y(),tp->vertex().z());
00266     
00267     std::pair<bool,reco::TrackBase::ParameterVector> params = parametersAtClosestApproach(vert, momAtVtx, tp->charge(), bs);
00268     if (params.first){
00269       TrackBase::ParameterVector sParameters=params.second;
00270       
00271       int tindex=0;
00272       for (RefToBaseVector<reco::Track>::const_iterator rt=tC.begin(); rt!=tC.end(); rt++, tindex++){
00273         
00274         TrackBase::ParameterVector rParameters = (*rt)->parameters();
00275 
00276         TrackBase::CovarianceMatrix recoTrackCovMatrix = (*rt)->covariance();
00277         if (onlyDiagonal) {
00278           for (unsigned int i=0;i<5;i++){
00279             for (unsigned int j=0;j<5;j++){
00280               if (i!=j) recoTrackCovMatrix(i,j)=0;
00281             }
00282           }
00283         }
00284         
00285         recoTrackCovMatrix.Invert();
00286         
00287         TrackBase::ParameterVector diffParameters = rParameters - sParameters;
00288         
00289         chi2 = ROOT::Math::Similarity(recoTrackCovMatrix, diffParameters);
00290         chi2 /= 5;
00291         
00292         LogTrace("TrackAssociator") << "====TP index=" << tpindex << "  RECO TRACK index="<<tindex<<" WITH PT=" << (*rt)->pt() << "====\n" 
00293                                     << "qoverp sim: " << sParameters[0] << "\n" 
00294                                     << "lambda sim: " << sParameters[1] << "\n" 
00295                                     << "phi    sim: " << sParameters[2] << "\n" 
00296                                     << "dxy    sim: " << sParameters[3] << "\n" 
00297                                     << "dsz    sim: " << sParameters[4] << "\n" 
00298                                     << ": " /*<< */ << "\n" 
00299                                     << "qoverp rec: " << (*rt)->qoverp() << " err: " << (*rt)->error(0) << "\n"
00300                                     << "lambda rec: " << (*rt)->lambda() << " err: " << (*rt)->error(1) << "\n"
00301                                     << "phi    rec: " << (*rt)->phi() << " err: " << (*rt)->error(2) << "\n"
00302                                     << "dxy    rec: " << (*rt)->dxy(bs.position()) << " err: " << (*rt)->error(3) << "\n"
00303                                     << "dsz    rec: " << (*rt)->dsz(bs.position()) << " err: " << (*rt)->error(4) << "\n"
00304                                     << ": " /*<< */ << "\n" 
00305                                     << "chi2: " << chi2 << "\n";
00306         
00307         if (chi2<chi2cut) {
00308           outputCollection.insert(edm::Ref<TrackingParticleCollection>(tPCH, tpindex),
00309                                   std::make_pair(tC[tindex],
00310                                                  -chi2));//-chi2 because the Association Map is ordered using std::greater
00311         }
00312       }
00313     }
00314   }
00315   outputCollection.post_insert();
00316   return outputCollection;
00317 }