CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/RecoTracker/DebugTools/plugins/TrackAlgoCompareUtil.cc

Go to the documentation of this file.
00001 #include "RecoTracker/DebugTools/plugins/TrackAlgoCompareUtil.h"
00002 
00003 using namespace std;
00004 using namespace edm;
00005 
00006 
00007 // constructors and destructor
00008 TrackAlgoCompareUtil::TrackAlgoCompareUtil(const edm::ParameterSet& iConfig)
00009 {
00010     //now do what ever other initialization is needed
00011     trackLabel_algoA = iConfig.getParameter<edm::InputTag>("trackLabel_algoA");
00012     trackLabel_algoB = iConfig.getParameter<edm::InputTag>("trackLabel_algoB");
00013     trackingParticleLabel_fakes = iConfig.getParameter<edm::InputTag>("trackingParticleLabel_fakes");
00014     trackingParticleLabel_effic = iConfig.getParameter<edm::InputTag>("trackingParticleLabel_effic");
00015     vertexLabel_algoA = iConfig.getParameter<edm::InputTag>("vertexLabel_algoA");
00016     vertexLabel_algoB = iConfig.getParameter<edm::InputTag>("vertexLabel_algoB");
00017     beamSpotLabel = iConfig.getParameter<edm::InputTag>("beamSpotLabel");
00018     assocLabel_algoA = iConfig.getUntrackedParameter<std::string>("assocLabel_algoA", "TrackAssociatorByHits");
00019     assocLabel_algoB = iConfig.getUntrackedParameter<std::string>("assocLabel_algoB", "TrackAssociatorByHits");
00020     associatormap_algoA = iConfig.getParameter< edm::InputTag >("associatormap_algoA");
00021     associatormap_algoB = iConfig.getParameter< edm::InputTag >("associatormap_algoB");
00022     UseAssociators = iConfig.getParameter< bool >("UseAssociators");
00023     UseVertex = iConfig.getParameter< bool >("UseVertex");
00024   
00025     produces<RecoTracktoTPCollection>("AlgoA");
00026     produces<RecoTracktoTPCollection>("AlgoB");
00027     produces<TPtoRecoTrackCollection>("TP");
00028 }
00029 
00030 
00031 TrackAlgoCompareUtil::~TrackAlgoCompareUtil()
00032 {
00033 }
00034 
00035 
00036 // ------------ method called once each job just before starting event loop  ------------
00037 void TrackAlgoCompareUtil::beginJob()
00038 {
00039 }
00040 
00041 
00042 // ------------ method called to produce the data  ------------
00043 void
00044 TrackAlgoCompareUtil::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00045 {
00046      // create output collection instance
00047     std::auto_ptr<RecoTracktoTPCollection> outputAlgoA(new RecoTracktoTPCollection());
00048     std::auto_ptr<RecoTracktoTPCollection> outputAlgoB(new RecoTracktoTPCollection());
00049     std::auto_ptr<TPtoRecoTrackCollection> outputTP(new TPtoRecoTrackCollection());
00050   
00051     // Get Inputs
00052     edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
00053     iEvent.getByLabel(beamSpotLabel, recoBeamSpotHandle);
00054     reco::BeamSpot beamSpot = *recoBeamSpotHandle; 
00055   
00056     edm::Handle<View<reco::Track> > trackCollAlgoA;
00057     iEvent.getByLabel(trackLabel_algoA, trackCollAlgoA);
00058   
00059     edm::Handle< View<reco::Track> > trackCollAlgoB;
00060     iEvent.getByLabel(trackLabel_algoB, trackCollAlgoB);
00061   
00062     edm::Handle<TrackingParticleCollection> trackingParticleCollFakes;
00063     iEvent.getByLabel(trackingParticleLabel_fakes, trackingParticleCollFakes);
00064   
00065     edm::Handle<TrackingParticleCollection> trackingParticleCollEffic;
00066     iEvent.getByLabel(trackingParticleLabel_effic, trackingParticleCollEffic);
00067   
00068     edm::Handle<reco::VertexCollection> vertexCollAlgoA;
00069     edm::Handle<reco::VertexCollection> vertexCollAlgoB;
00070     if(UseVertex) 
00071     {
00072         iEvent.getByLabel(vertexLabel_algoA, vertexCollAlgoA);
00073         iEvent.getByLabel(vertexLabel_algoB, vertexCollAlgoB);
00074     }
00075   
00076     // call the associator functions:
00077     reco::RecoToSimCollection recSimColl_AlgoA; 
00078     reco::RecoToSimCollection recSimColl_AlgoB;
00079   
00080     reco::SimToRecoCollection simRecColl_AlgoA;
00081     reco::SimToRecoCollection simRecColl_AlgoB;
00082 
00083     if(UseAssociators)
00084     {
00085         edm::ESHandle<TrackAssociatorBase> theAssociator_algoA;
00086         iSetup.get<TrackAssociatorRecord>().get(assocLabel_algoA, theAssociator_algoA);
00087   
00088         edm::ESHandle<TrackAssociatorBase> theAssociator_algoB;
00089         iSetup.get<TrackAssociatorRecord>().get(assocLabel_algoB, theAssociator_algoB);
00090   
00091         recSimColl_AlgoA = theAssociator_algoA->associateRecoToSim(trackCollAlgoA, trackingParticleCollFakes, &iEvent);
00092         recSimColl_AlgoB = theAssociator_algoB->associateRecoToSim(trackCollAlgoB, trackingParticleCollFakes, &iEvent);
00093 
00094         simRecColl_AlgoA = theAssociator_algoA->associateSimToReco(trackCollAlgoA, trackingParticleCollEffic, &iEvent);
00095         simRecColl_AlgoB = theAssociator_algoB->associateSimToReco(trackCollAlgoB, trackingParticleCollEffic, &iEvent);
00096     }
00097     else
00098     {
00099         Handle<reco::RecoToSimCollection > recotosimCollectionH_AlgoA;
00100         iEvent.getByLabel(associatormap_algoA,recotosimCollectionH_AlgoA);
00101         recSimColl_AlgoA  = *(recotosimCollectionH_AlgoA.product());
00102         
00103         Handle<reco::RecoToSimCollection > recotosimCollectionH_AlgoB;
00104         iEvent.getByLabel(associatormap_algoB,recotosimCollectionH_AlgoB);
00105         recSimColl_AlgoB  = *(recotosimCollectionH_AlgoB.product());
00106         
00107         Handle<reco::SimToRecoCollection > simtorecoCollectionH_AlgoA;
00108         iEvent.getByLabel(associatormap_algoA, simtorecoCollectionH_AlgoA);
00109         simRecColl_AlgoA = *(simtorecoCollectionH_AlgoA.product());
00110 
00111         Handle<reco::SimToRecoCollection > simtorecoCollectionH_AlgoB;
00112         iEvent.getByLabel(associatormap_algoB, simtorecoCollectionH_AlgoB);
00113         simRecColl_AlgoB = *(simtorecoCollectionH_AlgoB.product());
00114     }
00115     
00116     // define the vector of references to trackingParticleColl associated with a given reco::Track
00117     std::vector<std::pair<TrackingParticleRef, double> > associatedTrackingParticles;
00118   
00119     // define the vector of references to trackColl associated with a given TrackingParticle
00120     std::vector<std::pair<reco::TrackBaseRef, double> > associatedRecoTracks;
00121   
00122     // Get the magnetic field data from the event (used to calculate the point of closest TrackingParticle)
00123     edm::ESHandle<MagneticField> theMagneticField;
00124     iSetup.get<IdealMagneticFieldRecord>().get(theMagneticField);
00125     const MagneticField *magneticField = theMagneticField.product();
00126     
00127     // fill collection algoA
00128     for(View<reco::Track>::size_type i = 0; i < trackCollAlgoA->size(); ++i)
00129     {
00130         // get recoTrack algo A
00131         reco::TrackBaseRef recoTrack(trackCollAlgoA, i);
00132         RecoTracktoTP  recoTracktoTP;
00133         recoTracktoTP.SetRecoTrack(recoTrack);
00134         recoTracktoTP.SetBeamSpot(beamSpot.position());
00135         
00136         // get the associated trackingParticle
00137         if(recSimColl_AlgoA.find(recoTrack) != recSimColl_AlgoA.end())
00138         {
00139             associatedTrackingParticles = recSimColl_AlgoA[recoTrack];
00140             recoTracktoTP.SetTrackingParticle( associatedTrackingParticles.begin()->first );
00141             recoTracktoTP.SetShared( associatedTrackingParticles.begin()->second );
00142             SetTrackingParticleD0Dz(associatedTrackingParticles.begin()->first, beamSpot, magneticField, recoTracktoTP);
00143         }
00144         else
00145         {           
00146             recoTracktoTP.SetTrackingParticle(TrackingParticleRef());
00147             recoTracktoTP.SetShared(-1.0);
00148         }       
00149         
00150         // get the reco primary vertex info
00151         if(UseVertex && vertexCollAlgoA->size())  
00152         {
00153             recoTracktoTP.SetRecoVertex( reco::VertexRef(vertexCollAlgoA, 0) );
00154         }
00155         else
00156         {
00157             recoTracktoTP.SetRecoVertex( reco::VertexRef() );
00158         }   
00159        
00160         outputAlgoA->push_back(recoTracktoTP);
00161     }
00162   
00163   
00164     // fill collection algoB
00165     for(reco::TrackCollection::size_type i = 0; i < trackCollAlgoB->size(); ++i)
00166     {
00167         // get recoTrack algo B
00168         reco::TrackBaseRef recoTrack(trackCollAlgoB, i);
00169         RecoTracktoTP  recoTracktoTP;
00170         recoTracktoTP.SetRecoTrack(recoTrack);
00171         recoTracktoTP.SetBeamSpot(beamSpot.position());
00172       
00173         // get the associated trackingParticle
00174         if(recSimColl_AlgoB.find(recoTrack) != recSimColl_AlgoB.end())
00175         {
00176             associatedTrackingParticles = recSimColl_AlgoB[recoTrack];
00177             recoTracktoTP.SetTrackingParticle( associatedTrackingParticles.begin()->first );
00178             recoTracktoTP.SetShared( associatedTrackingParticles.begin()->second );
00179             SetTrackingParticleD0Dz(associatedTrackingParticles.begin()->first, beamSpot, magneticField, recoTracktoTP);
00180         }
00181         else
00182         {           
00183             recoTracktoTP.SetTrackingParticle(TrackingParticleRef());
00184             recoTracktoTP.SetShared(-1.0);
00185         }       
00186         
00187         // get the reco primary vertex info
00188         if(UseVertex && vertexCollAlgoB->size())  
00189         {
00190             recoTracktoTP.SetRecoVertex( reco::VertexRef(vertexCollAlgoB, 0) );
00191         }
00192         else
00193         {
00194             recoTracktoTP.SetRecoVertex( reco::VertexRef() );
00195         }
00196         
00197         outputAlgoB->push_back(recoTracktoTP);
00198     }
00199   
00200   
00201     for(TrackingParticleCollection::size_type i = 0; i < trackingParticleCollEffic->size(); ++i)
00202     {
00203         // initialize the trackingParticle (sim) info
00204         TrackingParticleRef tparticle(trackingParticleCollEffic, i);
00205         TPtoRecoTrack  tptoRecoTrack;
00206         tptoRecoTrack.SetBeamSpot(beamSpot.position());
00207         tptoRecoTrack.SetTrackingParticle(tparticle);
00208         SetTrackingParticleD0Dz(tparticle, beamSpot,  magneticField, tptoRecoTrack);
00209       
00210         // get the assocated recoTrack algoA
00211         if(simRecColl_AlgoA.find(tparticle) != simRecColl_AlgoA.end())
00212         {
00213             associatedRecoTracks = simRecColl_AlgoA[tparticle];
00214             tptoRecoTrack.SetRecoTrack_AlgoA(associatedRecoTracks.begin()->first );
00215             tptoRecoTrack.SetShared_AlgoA(associatedRecoTracks.begin()->second );
00216         }
00217         else
00218         {
00219             tptoRecoTrack.SetRecoTrack_AlgoA(reco::TrackBaseRef());
00220             tptoRecoTrack.SetShared_AlgoA(-1.0);
00221         }
00222         
00223         // get the recoVertex algo A
00224         if(UseVertex && vertexCollAlgoA->size())
00225         {
00226             tptoRecoTrack.SetRecoVertex_AlgoA( reco::VertexRef(vertexCollAlgoA, 0) );
00227         }
00228         else
00229         {
00230             tptoRecoTrack.SetRecoVertex_AlgoA( reco::VertexRef() );
00231         }
00232       
00233         // get the assocated recoTrack algoB
00234         if(simRecColl_AlgoB.find(tparticle) != simRecColl_AlgoB.end())
00235         {
00236             associatedRecoTracks = simRecColl_AlgoB[tparticle];
00237             tptoRecoTrack.SetRecoTrack_AlgoB(associatedRecoTracks.begin()->first );
00238             tptoRecoTrack.SetShared_AlgoB(associatedRecoTracks.begin()->second );
00239         }
00240         else
00241         {
00242             tptoRecoTrack.SetRecoTrack_AlgoB(reco::TrackBaseRef());
00243             tptoRecoTrack.SetShared_AlgoB(-1.0);
00244         }
00245         // get the recoVertex algo B
00246         if(UseVertex && vertexCollAlgoB->size())
00247         {
00248             tptoRecoTrack.SetRecoVertex_AlgoB( reco::VertexRef(vertexCollAlgoB, 0) );
00249         }
00250         else
00251         {
00252             tptoRecoTrack.SetRecoVertex_AlgoB( reco::VertexRef() );
00253         }
00254       
00255         outputTP->push_back(tptoRecoTrack);
00256     }
00257 
00258 
00259     // put the collection in the event record
00260     iEvent.put(outputAlgoA, "AlgoA");
00261     iEvent.put(outputAlgoB, "AlgoB");
00262     iEvent.put(outputTP, "TP");
00263 }
00264 
00265 // ------------ method called once each job just after ending the event loop  ------------
00266 void TrackAlgoCompareUtil::endJob() 
00267 {
00268 }
00269 
00270 
00271 // ------------ Producer Specific Meber Fucntions ----------------------------------------
00272 void TrackAlgoCompareUtil::SetTrackingParticleD0Dz(TrackingParticleRef tp, const reco::BeamSpot &bs, const MagneticField *bf, TPtoRecoTrack& TPRT)
00273 {
00274     GlobalPoint trackingParticleVertex( tp->vertex().x(), tp->vertex().y(), tp->vertex().z() );
00275     GlobalVector trackingParticleP3(tp->g4Track_begin()->momentum().x(),
00276                                     tp->g4Track_begin()->momentum().y(),
00277                                     tp->g4Track_begin()->momentum().z() );
00278     TrackCharge trackingParticleCharge(tp->charge());
00279   
00280     FreeTrajectoryState ftsAtProduction( trackingParticleVertex, trackingParticleP3, trackingParticleCharge, bf );
00281     TSCBLBuilderNoMaterial tscblBuilder;
00282     TrajectoryStateClosestToBeamLine tsAtClosestApproach = tscblBuilder(ftsAtProduction, bs);  //as in TrackProducerAlgorithm
00283   
00284     if(tsAtClosestApproach.isValid())
00285     {
00286         GlobalPoint v1 = tsAtClosestApproach.trackStateAtPCA().position();
00287         GlobalVector p = tsAtClosestApproach.trackStateAtPCA().momentum();
00288 
00289         TPRT.SetTrackingParticleMomentumPCA(p);
00290         TPRT.SetTrackingParticlePCA(v1);
00291     }
00292     else
00293     {
00294         TPRT.SetTrackingParticleMomentumPCA(GlobalVector(-9999.0, -9999.0, -9999.0));
00295         TPRT.SetTrackingParticlePCA(GlobalPoint(-9999.0, -9999.0, -9999.0));
00296     }
00297 }
00298 
00299 
00300 void TrackAlgoCompareUtil::SetTrackingParticleD0Dz(TrackingParticleRef tp, const reco::BeamSpot &bs, const MagneticField *bf, RecoTracktoTP& RTTP)
00301 {
00302     GlobalPoint trackingParticleVertex( tp->vertex().x(), tp->vertex().y(), tp->vertex().z() );
00303     GlobalVector trackingParticleP3(tp->g4Track_begin()->momentum().x(),
00304                                     tp->g4Track_begin()->momentum().y(),
00305                                     tp->g4Track_begin()->momentum().z() );
00306     TrackCharge trackingParticleCharge(tp->charge());
00307 
00308     FreeTrajectoryState ftsAtProduction( trackingParticleVertex, trackingParticleP3, trackingParticleCharge, bf );
00309     TSCBLBuilderNoMaterial tscblBuilder;
00310     TrajectoryStateClosestToBeamLine tsAtClosestApproach = tscblBuilder(ftsAtProduction, bs);  //as in TrackProducerAlgorithm
00311 
00312     if(tsAtClosestApproach.isValid())
00313     {
00314         GlobalPoint v1 = tsAtClosestApproach.trackStateAtPCA().position();
00315         GlobalVector p = tsAtClosestApproach.trackStateAtPCA().momentum();
00316 
00317         RTTP.SetTrackingParticleMomentumPCA(p);
00318         RTTP.SetTrackingParticlePCA(v1);
00319     }
00320     else
00321     {
00322         RTTP.SetTrackingParticleMomentumPCA(GlobalVector(-9999.0, -9999.0, -9999.0));
00323         RTTP.SetTrackingParticlePCA(GlobalPoint(-9999.0, -9999.0, -9999.0));
00324     }
00325 }
00326