CMS 3D CMS Logo

TrackingTruthValid Class Reference

#include <Validation/TrackingMCTruth/interface/TrackingTruthValid.h>

Inheritance diagram for TrackingTruthValid:

edm::EDAnalyzer

List of all members.

Public Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
void beginJob (const edm::ParameterSet &conf)
void endJob ()
 TrackingTruthValid (const edm::ParameterSet &conf)
 ~TrackingTruthValid ()

Private Attributes

edm::ParameterSet conf_
DQMStoredbe_
MonitorElementmeTPAllHits
MonitorElementmeTPCharge
MonitorElementmeTPEta
MonitorElementmeTPId
MonitorElementmeTPlip
MonitorElementmeTPMass
MonitorElementmeTPMatchedHits
MonitorElementmeTPPhi
MonitorElementmeTPProc
MonitorElementmeTPPt
MonitorElementmeTPtip
MonitorElementmeTPVtxX
MonitorElementmeTPVtxY
MonitorElementmeTPVtxZ
std::string outputFile
edm::InputTag src_


Detailed Description

Definition at line 20 of file TrackingTruthValid.h.


Constructor & Destructor Documentation

TrackingTruthValid::TrackingTruthValid ( const edm::ParameterSet conf  )  [explicit]

Definition at line 37 of file TrackingTruthValid.cc.

References dbe_, edm::ParameterSet::getParameter(), meTPAllHits, meTPCharge, meTPEta, meTPId, meTPlip, meTPMass, meTPMatchedHits, meTPPhi, meTPProc, meTPPt, meTPtip, meTPVtxX, meTPVtxY, meTPVtxZ, outputFile, MonitorElement::setBinLabel(), and src_.

00037                                                                   {
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("TrackingMCTruthV/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 }

TrackingTruthValid::~TrackingTruthValid (  )  [inline]

Definition at line 25 of file TrackingTruthValid.h.

00025 {} ;


Member Function Documentation

void TrackingTruthValid::analyze ( const edm::Event event,
const edm::EventSetup c 
) [virtual]

Implements edm::EDAnalyzer.

Definition at line 85 of file TrackingTruthValid.cc.

References begin, end, MonitorElement::Fill(), prof2calltree::front, meTPAllHits, meTPCharge, meTPEta, meTPId, meTPlip, meTPMass, meTPMatchedHits, meTPPhi, meTPProc, meTPPt, meTPtip, meTPVtxX, meTPVtxY, meTPVtxZ, edm::Handle< T >::product(), funct::sqrt(), src_, std, and t.

00085                                                                              {
00086   using namespace std;
00087 
00088   edm::Handle<TrackingParticleCollection>  TruthTrackContainer ;
00089   edm::Handle<TrackingVertexCollection>    TruthVertexContainer;
00090 
00091   event.getByLabel(src_,TruthTrackContainer );
00092   event.getByLabel(src_,TruthVertexContainer);
00093 
00094   //  std::cout << "Using Collection " << src_ << std::endl;
00095   
00096   TrackingParticleCollection *tPC   = const_cast<TrackingParticleCollection*>(TruthTrackContainer.product());
00097   const TrackingVertexCollection   *tVC   = TruthVertexContainer.product();
00098   
00099   /*
00100   // Get and print HepMC event for comparison
00101   edm::Handle<edm::HepMCProduct> hepMC;
00102   event.getByLabel("source",hepMC);
00103   const edm::HepMCProduct *mcp = hepMC.product();
00104   //  const HepMC::GenEvent *genEvent = mcp -> GetEvent();
00105   */
00106 
00107   //  cout << "Found " << tPC -> size() << " tracks and " << tVC -> size() << " vertices." <<endl;
00108    
00109 
00110 // Loop over TrackingParticle's
00111 
00112   for (TrackingParticleCollection::iterator t = tPC -> begin(); t != tPC -> end(); ++t) {
00113     //if(t -> trackerPSimHit().size() ==0) cout << " Track with 0 SimHit " << endl;
00114 
00115 
00116     meTPMass->Fill(t->mass());
00117     meTPCharge->Fill(t->charge() );
00118     meTPId->Fill(t->pdgId());
00119     meTPPt->Fill(sqrt(t->momentum().perp2()));
00120     meTPEta->Fill(t->momentum().eta());
00121     meTPPhi->Fill(t->momentum().Phi());
00122     meTPAllHits->Fill(t->trackerPSimHit().size());
00123     //get the process of the first hit
00124     if(t -> trackerPSimHit().size() !=0) meTPProc->Fill( t -> trackerPSimHit().front().processType());
00125     meTPMatchedHits->Fill(t->matchedHit());
00126     meTPVtxX->Fill(sqrt(t->vertex().x()));
00127     meTPVtxY->Fill(sqrt(t->vertex().y()));
00128     meTPVtxZ->Fill(sqrt(t->vertex().z()));
00129     meTPtip->Fill(sqrt(t->vertex().perp2()));
00130     meTPlip->Fill(sqrt(t->vertex().z()));
00131 
00132       /*
00133    // Compare momenta from sources
00134     cout << "T.P.   Track mass, Momentum, q , ID, & Event # "
00135           << t -> mass()  << " " 
00136           << t -> p4()    << " " << t -> charge() << " "
00137           << t -> pdgId() << " "
00138           << t -> eventId().bunchCrossing() << "." << t -> eventId().event() << endl;
00139 
00140     if(t->mass() < 0) cout << "======= WARNING, this particle has negative mass: " << t->mass()  
00141                            << " and pdgId: " << t->pdgId() << endl;
00142     if(t->pdgId() == 0) cout << "======= WARNING, this particle has pdgId = 0: "     << t->pdgId() << endl;
00143     cout << " Hits for this track: " << t -> trackerPSimHit().size() << endl;
00144     */
00145 
00146     /*
00147       std::cout << std::endl << "### Tracking Particle ###" << std::endl;
00148       std::cout << (*t) << std::endl;
00149       std::cout << "\t Tracker: " << t->trackerPSimHit().size() << std::endl;
00150       std::cout << "\t Muon: "    << t->muonPSimHit().size()    << std::endl;
00151       std::cout << (*t) << std::endl;
00152     */
00153   }  // End loop over TrackingParticle
00154   
00155   // Loop over TrackingVertex's
00156   /*  
00157   cout << "Dumping sample vertex info" << endl;
00158   for (TrackingVertexCollection::const_iterator v = tVC -> begin(); v != tVC -> end(); ++v) {
00159     cout << " Vertex Position & Event #" << v -> position() << " " << v -> eventId().bunchCrossing() << "." << v -> eventId().event() << endl;
00160     cout << "  Associated with " << v -> daughterTracks().size() << " tracks" << endl;
00161     // Get Geant and HepMC positions
00162     for (genv_iterator genV = v -> genVertices_begin(); genV != v -> genVertices_end(); ++genV) {
00163       cout << "  HepMC vertex position " << (*(*genV)).position() << endl;
00164     }
00165     for (g4v_iterator g4V = v -> g4Vertices_begin(); g4V != v -> g4Vertices_end(); ++g4V) {
00166       cout << "  Geant vertex position " << (*g4V).position() << endl;
00167       // Probably empty all the time, currently
00168     }
00169 
00170     // Loop over daughter track(s)
00171     for (tp_iterator iTP = v -> daughterTracks_begin(); iTP != v -> daughterTracks_end(); ++iTP) {
00172       cout << "  Daughter starts:      " << (*(*iTP)).vertex();
00173       for (g4t_iterator g4T  = (*(*iTP)).g4Track_begin(); g4T != (*(*iTP)).g4Track_end(); ++g4T) {
00174         cout << " p " << g4T->momentum();
00175       }
00176       cout << endl;
00177     }
00178 
00179     // Loop over source track(s) (can be multiple since vertices are collapsed)
00180     for (tp_iterator iTP = v -> sourceTracks_begin(); iTP != v -> sourceTracks_end(); ++iTP) {
00181       cout << "  Source   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   }  // End loop over TrackingVertex
00188   */
00189 
00190 
00191 }

void TrackingTruthValid::beginJob ( const edm::ParameterSet conf  ) 

Definition at line 35 of file TrackingTruthValid.cc.

00035 {}

void TrackingTruthValid::endJob ( void   )  [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 193 of file TrackingTruthValid.cc.

References dbe_, outputFile, and DQMStore::save().

00193                                { 
00194 
00195   if ( outputFile.size() != 0 && dbe_ ) dbe_->save(outputFile);
00196 
00197 } 


Member Data Documentation

edm::ParameterSet TrackingTruthValid::conf_ [private]

Definition at line 34 of file TrackingTruthValid.h.

DQMStore* TrackingTruthValid::dbe_ [private]

Definition at line 33 of file TrackingTruthValid.h.

Referenced by endJob(), and TrackingTruthValid().

MonitorElement* TrackingTruthValid::meTPAllHits [private]

Definition at line 42 of file TrackingTruthValid.h.

Referenced by analyze(), and TrackingTruthValid().

MonitorElement* TrackingTruthValid::meTPCharge [private]

Definition at line 39 of file TrackingTruthValid.h.

Referenced by analyze(), and TrackingTruthValid().

MonitorElement* TrackingTruthValid::meTPEta [private]

Definition at line 45 of file TrackingTruthValid.h.

Referenced by analyze(), and TrackingTruthValid().

MonitorElement* TrackingTruthValid::meTPId [private]

Definition at line 40 of file TrackingTruthValid.h.

Referenced by analyze(), and TrackingTruthValid().

MonitorElement* TrackingTruthValid::meTPlip [private]

Definition at line 51 of file TrackingTruthValid.h.

Referenced by analyze(), and TrackingTruthValid().

MonitorElement* TrackingTruthValid::meTPMass [private]

Definition at line 38 of file TrackingTruthValid.h.

Referenced by analyze(), and TrackingTruthValid().

MonitorElement* TrackingTruthValid::meTPMatchedHits [private]

Definition at line 43 of file TrackingTruthValid.h.

Referenced by analyze(), and TrackingTruthValid().

MonitorElement* TrackingTruthValid::meTPPhi [private]

Definition at line 46 of file TrackingTruthValid.h.

Referenced by analyze(), and TrackingTruthValid().

MonitorElement* TrackingTruthValid::meTPProc [private]

Definition at line 41 of file TrackingTruthValid.h.

Referenced by analyze(), and TrackingTruthValid().

MonitorElement* TrackingTruthValid::meTPPt [private]

Definition at line 44 of file TrackingTruthValid.h.

Referenced by analyze(), and TrackingTruthValid().

MonitorElement* TrackingTruthValid::meTPtip [private]

Definition at line 50 of file TrackingTruthValid.h.

Referenced by analyze(), and TrackingTruthValid().

MonitorElement* TrackingTruthValid::meTPVtxX [private]

Definition at line 47 of file TrackingTruthValid.h.

Referenced by analyze(), and TrackingTruthValid().

MonitorElement* TrackingTruthValid::meTPVtxY [private]

Definition at line 48 of file TrackingTruthValid.h.

Referenced by analyze(), and TrackingTruthValid().

MonitorElement* TrackingTruthValid::meTPVtxZ [private]

Definition at line 49 of file TrackingTruthValid.h.

Referenced by analyze(), and TrackingTruthValid().

std::string TrackingTruthValid::outputFile [private]

Definition at line 35 of file TrackingTruthValid.h.

Referenced by endJob(), and TrackingTruthValid().

edm::InputTag TrackingTruthValid::src_ [private]

Definition at line 36 of file TrackingTruthValid.h.

Referenced by analyze(), and TrackingTruthValid().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:34:14 2009 for CMSSW by  doxygen 1.5.4