Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:31:01

0001 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0002 #include "FWCore/Framework/interface/Event.h"
0003 #include "FWCore/Framework/interface/EventSetup.h"
0004 #include "FWCore/Utilities/interface/InputTag.h"
0005 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0006 #include "MagneticField/Engine/interface/MagneticField.h"
0007 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
0008 
0009 #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticle.h"
0010 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0011 #include "SimDataFormats/Associations/interface/TrackToTrackingParticleAssociator.h"
0012 
0013 #include <iostream>
0014 #include <map>
0015 #include <memory>
0016 #include <set>
0017 #include <string>
0018 
0019 // class TrackAssociator;
0020 class TrackAssociatorByHits;
0021 class TrackerHitAssociator;
0022 
0023 namespace reco {
0024   class TrackToTrackingParticleAssociator;
0025 }
0026 
0027 class testTrackAssociator : public edm::one::EDAnalyzer<> {
0028 public:
0029   testTrackAssociator(const edm::ParameterSet &conf);
0030   ~testTrackAssociator() override = default;
0031   void beginJob() override {}
0032   void analyze(const edm::Event &, const edm::EventSetup &) override;
0033 
0034 private:
0035   reco::TrackToTrackingParticleAssociator const *associatorByChi2;
0036   reco::TrackToTrackingParticleAssociator const *associatorByHits;
0037   edm::InputTag tracksTag, tpTag, simtracksTag, simvtxTag;
0038   edm::ESGetToken<MagneticField, IdealMagneticFieldRecord> tokenMF_;
0039 };
0040 
0041 using namespace reco;
0042 using namespace std;
0043 using namespace edm;
0044 
0045 testTrackAssociator::testTrackAssociator(edm::ParameterSet const &conf) {
0046   tracksTag = conf.getParameter<edm::InputTag>("tracksTag");
0047   tpTag = conf.getParameter<edm::InputTag>("tpTag");
0048   simtracksTag = conf.getParameter<edm::InputTag>("simtracksTag");
0049   simvtxTag = conf.getParameter<edm::InputTag>("simvtxTag");
0050   tokenMF_ = esConsumes<MagneticField, IdealMagneticFieldRecord>();
0051 }
0052 
0053 void testTrackAssociator::analyze(const edm::Event &event, const edm::EventSetup &setup) {
0054   //const auto &theMF = setup.getHandle(tokenMF_);
0055   edm::Handle<reco::TrackToTrackingParticleAssociator> theChiAssociator;
0056   event.getByLabel("trackAssociatorByChi2", theChiAssociator);
0057   associatorByChi2 = theChiAssociator.product();
0058   edm::Handle<reco::TrackToTrackingParticleAssociator> theHitsAssociator;
0059   event.getByLabel("trackAssociatorByHits", theHitsAssociator);
0060   associatorByHits = theHitsAssociator.product();
0061 
0062   Handle<View<Track>> trackCollectionH;
0063   event.getByLabel(tracksTag, trackCollectionH);
0064   const View<Track> &tC = *(trackCollectionH.product());
0065 
0066   Handle<SimTrackContainer> simTrackCollection;
0067   event.getByLabel(simtracksTag, simTrackCollection);
0068   const SimTrackContainer &simTC = *(simTrackCollection.product());
0069 
0070   Handle<SimVertexContainer> simVertexCollection;
0071   event.getByLabel(simvtxTag, simVertexCollection);
0072 
0073   edm::Handle<TrackingParticleCollection> TPCollectionH;
0074   event.getByLabel(tpTag, TPCollectionH);
0075 
0076   cout << "\nEvent ID = " << event.id() << endl;
0077 
0078   // RECOTOSIM
0079   cout << "                      ****************** Reco To Sim "
0080           "****************** "
0081        << endl;
0082   cout << "-- Associator by hits --" << endl;
0083   reco::RecoToSimCollection p = associatorByHits->associateRecoToSim(trackCollectionH, TPCollectionH);
0084   for (View<Track>::size_type i = 0; i < tC.size(); ++i) {
0085     RefToBase<Track> track(trackCollectionH, i);
0086     try {
0087       std::vector<std::pair<TrackingParticleRef, double>> tp = p[track];
0088       cout << "Reco Track pT: " << setw(6) << track->pt() << " matched to " << tp.size() << " MC Tracks" << std::endl;
0089       for (std::vector<std::pair<TrackingParticleRef, double>>::const_iterator it = tp.begin(); it != tp.end(); ++it) {
0090         TrackingParticleRef tpr = it->first;
0091         double assocChi2 = it->second;
0092         cout << "\t\tMCTrack " << setw(2) << tpr.index() << " pT: " << setw(6) << tpr->pt() << " NShared: " << assocChi2
0093              << endl;
0094       }
0095     } catch (Exception const &) {
0096       cout << "->   Track pT: " << setprecision(2) << setw(6) << track->pt() << " matched to 0  MC Tracks" << endl;
0097     }
0098   }
0099   cout << "-- Associator by chi2 --" << endl;
0100   p = associatorByChi2->associateRecoToSim(trackCollectionH, TPCollectionH);
0101   for (View<Track>::size_type i = 0; i < tC.size(); ++i) {
0102     RefToBase<Track> track(trackCollectionH, i);
0103     try {
0104       std::vector<std::pair<TrackingParticleRef, double>> tp = p[track];
0105       cout << "Reco Track pT: " << setw(6) << track->pt() << " matched to " << tp.size() << " MC Tracks" << std::endl;
0106       for (std::vector<std::pair<TrackingParticleRef, double>>::const_iterator it = tp.begin(); it != tp.end(); ++it) {
0107         TrackingParticleRef tpr = it->first;
0108         double assocChi2 = it->second;
0109         cout << "\t\tMCTrack " << setw(2) << tpr.index() << " pT: " << setw(6) << tpr->pt() << " chi2: " << assocChi2
0110              << endl;
0111       }
0112     } catch (Exception const &) {
0113       cout << "->   Track pT: " << setprecision(2) << setw(6) << track->pt() << " matched to 0  MC Tracks" << endl;
0114     }
0115   }
0116   // SIMTORECO
0117   cout << "                      ****************** Sim To Reco "
0118           "****************** "
0119        << endl;
0120   cout << "-- Associator by hits --" << endl;
0121   reco::SimToRecoCollection q = associatorByHits->associateSimToReco(trackCollectionH, TPCollectionH);
0122   for (SimTrackContainer::size_type i = 0; i < simTC.size(); ++i) {
0123     TrackingParticleRef tp(TPCollectionH, i);
0124     try {
0125       std::vector<std::pair<RefToBase<Track>, double>> trackV = q[tp];
0126       cout << "Sim Track " << setw(2) << tp.index() << " pT: " << setw(6) << tp->pt() << " matched to " << trackV.size()
0127            << " reco::Tracks" << std::endl;
0128       for (std::vector<std::pair<RefToBase<Track>, double>>::const_iterator it = trackV.begin(); it != trackV.end();
0129            ++it) {
0130         RefToBase<Track> tr = it->first;
0131         double assocChi2 = it->second;
0132         cout << "\t\treco::Track pT: " << setw(6) << tr->pt() << " NShared: " << assocChi2 << endl;
0133       }
0134     } catch (Exception const &) {
0135       cout << "->   TrackingParticle " << setw(2) << tp.index() << " pT: " << setprecision(2) << setw(6) << tp->pt()
0136            << " matched to 0  reco::Tracks" << endl;
0137     }
0138   }
0139   cout << "-- Associator by chi2 --" << endl;
0140   q = associatorByChi2->associateSimToReco(trackCollectionH, TPCollectionH);
0141   for (SimTrackContainer::size_type i = 0; i < simTC.size(); ++i) {
0142     TrackingParticleRef tp(TPCollectionH, i);
0143     try {
0144       std::vector<std::pair<RefToBase<Track>, double>> trackV = q[tp];
0145       cout << "Sim Track " << setw(2) << tp.index() << " pT: " << setw(6) << tp->pt() << " matched to " << trackV.size()
0146            << " reco::Tracks" << std::endl;
0147       for (std::vector<std::pair<RefToBase<Track>, double>>::const_iterator it = trackV.begin(); it != trackV.end();
0148            ++it) {
0149         RefToBase<Track> tr = it->first;
0150         double assocChi2 = it->second;
0151         cout << "\t\treco::Track pT: " << setw(6) << tr->pt() << " chi2: " << assocChi2 << endl;
0152       }
0153     } catch (Exception const &) {
0154       cout << "->   TrackingParticle " << setw(2) << tp.index() << " pT: " << setprecision(2) << setw(6) << tp->pt()
0155            << " matched to 0  reco::Tracks" << endl;
0156     }
0157   }
0158 }
0159 
0160 #include "FWCore/Framework/interface/MakerMacros.h"
0161 #include "FWCore/PluginManager/interface/ModuleDef.h"
0162 
0163 DEFINE_FWK_MODULE(testTrackAssociator);