CMS 3D CMS Logo

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