CMS 3D CMS Logo

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