CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch13/src/Validation/TrackingMCTruth/plugins/TrackingTruthValid.cc

Go to the documentation of this file.
00001 #include "Validation/TrackingMCTruth/interface/TrackingTruthValid.h"
00002 #include "FWCore/Framework/interface/EDAnalyzer.h"
00003 #include "FWCore/Framework/interface/Event.h"
00004 #include "DataFormats/Common/interface/Handle.h"
00005 #include "FWCore/Framework/interface/EventSetup.h"
00006 #include "FWCore/Framework/interface/Event.h"
00007 #include "FWCore/Framework/interface/MakerMacros.h"
00008 #include "FWCore/PluginManager/interface/ModuleDef.h"
00009 
00010 #include "DataFormats/Common/interface/EDProduct.h"
00011 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00012 #include "DQMServices/Core/interface/DQMStore.h"
00013 
00014 #include "SimDataFormats/CrossingFrame/interface/CrossingFrame.h"
00015 #include "SimDataFormats/CrossingFrame/interface/MixCollection.h"
00016 #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticleFwd.h"
00017 #include "SimDataFormats/TrackingAnalysis/interface/TrackingVertexContainer.h"
00018 
00019 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00020 
00021 using namespace std;
00022 using namespace ROOT::Math;
00023 using namespace edm;
00024 
00025 typedef edm::RefVector< std::vector<TrackingParticle> > TrackingParticleContainer;
00026 typedef std::vector<TrackingParticle>                   TrackingParticleCollection;
00027 
00028 typedef TrackingParticleRefVector::iterator               tp_iterator;
00029 typedef TrackingParticle::g4t_iterator                   g4t_iterator;
00030 typedef TrackingParticle::genp_iterator                 genp_iterator;
00031 typedef TrackingVertex::genv_iterator                   genv_iterator;
00032 typedef TrackingVertex::g4v_iterator                     g4v_iterator;
00033 
00034 
00035 void TrackingTruthValid::beginJob(const edm::ParameterSet& conf) {}
00036 
00037 TrackingTruthValid::TrackingTruthValid(const edm::ParameterSet& conf) {
00038   
00039   outputFile = conf.getParameter<std::string>("outputFile");
00040   src_ =  conf.getParameter<edm::InputTag>( "src" );
00041   
00042   dbe_  = edm::Service<DQMStore>().operator->();
00043   dbe_->setCurrentFolder("Tracking/TrackingMCTruth/TrackingParticle");
00044   
00045 
00046   meTPMass = dbe_->book1D("TPMass","Tracking Particle Mass",100, -1,+5.);
00047   meTPCharge = dbe_->book1D("TPCharge","Tracking Particle Charge",10, -5, 5);
00048   meTPId = dbe_->book1D("TPId","Tracking Particle Id",500, -5000, 5000);
00049   meTPProc = dbe_->book1D("TPProc","Tracking Particle Proc",20, -0.5, 19.5);
00050   meTPAllHits = dbe_->book1D("TPAllHits", "Tracking Particle All Hits", 200, -0.5, 199.5);
00051   meTPMatchedHits = dbe_->book1D("TPMatchedHits", "Tracking Particle Matched Hits", 100, -0.5, 99.5);
00052   meTPPt = dbe_->book1D("TPPt", "Tracking Particle Pt",100, 0, 100.);
00053   meTPEta = dbe_->book1D("TPEta", "Tracking Particle Eta",100, -7., 7.);
00054   meTPPhi = dbe_->book1D("TPPhi", "Tracking Particle Phi",100, -4., 4);
00055   meTPVtxX = dbe_->book1D("TPVtxX", "Tracking Particle VtxX",100, -100, 100.);
00056   meTPVtxY = dbe_->book1D("TPVtxY", "Tracking Particle VtxY",100, -100, 100.);
00057   meTPVtxZ = dbe_->book1D("TPVtxZ", "Tracking Particle VtxZ",100, -100, 100.);
00058   meTPtip = dbe_->book1D("TPtip", "Tracking Particle tip",100, 0, 1000.);
00059   meTPlip = dbe_->book1D("TPlip", "Tracking Particle lip",100, 0, 100.);
00060   
00061   
00062   // Prepare Axes Labels for Processes
00063   meTPProc->setBinLabel( 1,"Undefined");            // value =  0
00064   meTPProc->setBinLabel( 2,"Unknown");              // value =  1
00065   meTPProc->setBinLabel( 3,"Primary");              // value =  2
00066   meTPProc->setBinLabel( 4,"Hadronic");             // value =  3   
00067   meTPProc->setBinLabel( 5,"Decay");                // value =  4 
00068   meTPProc->setBinLabel( 6,"Compton");              // value =  5
00069   meTPProc->setBinLabel( 7,"Annihilation");         // value =  6
00070   meTPProc->setBinLabel( 8,"EIoni");                // value =  7
00071   meTPProc->setBinLabel( 9,"HIoni");                // value =  8
00072   meTPProc->setBinLabel(10,"MuIoni");               // value =  9
00073   meTPProc->setBinLabel(11,"Photon");               // value = 10
00074   meTPProc->setBinLabel(12,"MuPairProd");           // value = 11
00075   meTPProc->setBinLabel(13,"Conversions");          // value = 12
00076   meTPProc->setBinLabel(14,"EBrem");                // value = 13
00077   meTPProc->setBinLabel(15,"SynchrotronRadiation"); // value = 14
00078   meTPProc->setBinLabel(16,"MuBrem");               // value = 15
00079   meTPProc->setBinLabel(17,"MuNucl");               // value = 16
00080   meTPProc->setBinLabel(18,"");
00081   meTPProc->setBinLabel(19,"");
00082   meTPProc->setBinLabel(20,"");
00083 }
00084 
00085 void TrackingTruthValid::analyze(const edm::Event& event, const edm::EventSetup& c){
00086   using namespace std;
00087 
00088   edm::Handle<TrackingParticleCollection>  TruthTrackContainer ;
00089   //  edm::Handle<TrackingVertexCollection>    TruthVertexContainer;
00090 
00091 
00092   event.getByLabel(src_,TruthTrackContainer );
00093 
00094   //  event.getByLabel(src_,TruthVertexContainer);
00095 
00096   //  std::cout << "Using Collection " << src_ << std::endl;
00097   
00098   const TrackingParticleCollection *tPC   = TruthTrackContainer.product();
00099 
00100   //  const TrackingVertexCollection   *tVC   = TruthVertexContainer.product();
00101 
00102   /*
00103   // Get and print HepMC event for comparison
00104   edm::Handle<edm::HepMCProduct> hepMC;
00105   event.getByLabel("source",hepMC);
00106   const edm::HepMCProduct *mcp = hepMC.product();
00107   //  const HepMC::GenEvent *genEvent = mcp -> GetEvent();
00108   */
00109 
00110   //  cout << "Found " << tPC -> size() << " tracks and " << tVC -> size() << " vertices." <<endl;
00111    
00112 
00113 // Loop over TrackingParticle's
00114 
00115   for (TrackingParticleCollection::const_iterator t = tPC -> begin(); t != tPC -> end(); ++t) {
00116     //if(t -> trackerPSimHit().size() ==0) cout << " Track with 0 SimHit " << endl;
00117 
00118 
00119     meTPMass->Fill(t->mass());
00120 
00121     meTPCharge->Fill(t->charge() );
00122 
00123     meTPId->Fill(t->pdgId());
00124 
00125     meTPPt->Fill(sqrt(t->momentum().perp2()));
00126 
00127     meTPEta->Fill(t->momentum().eta());
00128 
00129     meTPPhi->Fill(t->momentum().Phi());
00130     std::vector<PSimHit> trackerPSimHit( t->trackPSimHit(DetId::Tracker) );
00131     meTPAllHits->Fill(trackerPSimHit.size());
00132     //get the process of the first hit
00133     if(trackerPSimHit.size() !=0) meTPProc->Fill( trackerPSimHit.front().processType());
00134     meTPMatchedHits->Fill(t->matchedHit());
00135     meTPVtxX->Fill(t->vx());
00136     meTPVtxY->Fill(t->vy());
00137     meTPVtxZ->Fill(t->vz());
00138     meTPtip->Fill(sqrt(t->vertex().perp2()));
00139     meTPlip->Fill(t->vz());
00140 
00141       /*
00142    // Compare momenta from sources
00143     cout << "T.P.   Track mass, Momentum, q , ID, & Event # "
00144           << t -> mass()  << " " 
00145           << t -> p4()    << " " << t -> charge() << " "
00146           << t -> pdgId() << " "
00147           << t -> eventId().bunchCrossing() << "." << t -> eventId().event() << endl;
00148 
00149     if(t->mass() < 0) cout << "======= WARNING, this particle has negative mass: " << t->mass()  
00150                            << " and pdgId: " << t->pdgId() << endl;
00151     if(t->pdgId() == 0) cout << "======= WARNING, this particle has pdgId = 0: "     << t->pdgId() << endl;
00152     cout << " Hits for this track: " << t -> trackerPSimHit().size() << endl;
00153     */
00154 
00155     /*
00156       std::cout << std::endl << "### Tracking Particle ###" << std::endl;
00157       std::cout << (*t) << std::endl;
00158       std::cout << "\t Tracker: " << t->trackerPSimHit().size() << std::endl;
00159       std::cout << "\t Muon: "    << t->muonPSimHit().size()    << std::endl;
00160       std::cout << (*t) << std::endl;
00161     */
00162   }  // End loop over TrackingParticle
00163   
00164   // Loop over TrackingVertex's
00165   /*  
00166   cout << "Dumping sample vertex info" << endl;
00167   for (TrackingVertexCollection::const_iterator v = tVC -> begin(); v != tVC -> end(); ++v) {
00168     cout << " Vertex Position & Event #" << v -> position() << " " << v -> eventId().bunchCrossing() << "." << v -> eventId().event() << endl;
00169     cout << "  Associated with " << v -> daughterTracks().size() << " tracks" << endl;
00170     // Get Geant and HepMC positions
00171     for (genv_iterator genV = v -> genVertices_begin(); genV != v -> genVertices_end(); ++genV) {
00172       cout << "  HepMC vertex position " << (*(*genV)).position() << endl;
00173     }
00174     for (g4v_iterator g4V = v -> g4Vertices_begin(); g4V != v -> g4Vertices_end(); ++g4V) {
00175       cout << "  Geant vertex position " << (*g4V).position() << endl;
00176       // Probably empty all the time, currently
00177     }
00178 
00179     // Loop over daughter track(s)
00180     for (tp_iterator iTP = v -> daughterTracks_begin(); iTP != v -> daughterTracks_end(); ++iTP) {
00181       cout << "  Daughter starts:      " << (*(*iTP)).vertex();
00182       for (g4t_iterator g4T  = (*(*iTP)).g4Track_begin(); g4T != (*(*iTP)).g4Track_end(); ++g4T) {
00183         cout << " p " << g4T->momentum();
00184       }
00185       cout << endl;
00186     }
00187 
00188     // Loop over source track(s) (can be multiple since vertices are collapsed)
00189     for (tp_iterator iTP = v -> sourceTracks_begin(); iTP != v -> sourceTracks_end(); ++iTP) {
00190       cout << "  Source   starts: " << (*(*iTP)).vertex();
00191       for (g4t_iterator g4T  = (*iTP)->g4Track_begin(); g4T != (*iTP)->g4Track_end(); ++g4T) {
00192         cout << ", p " <<  g4T ->momentum();
00193       }
00194       cout << endl;
00195     }
00196   }  // End loop over TrackingVertex
00197   */
00198 
00199 
00200 }
00201 
00202 void TrackingTruthValid::endJob(){ 
00203 
00204   if ( outputFile.size() != 0 && dbe_ ) dbe_->save(outputFile);
00205 
00206 }