Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:29:43

0001 #include "SimDataFormats/TrackingAnalysis/interface/TrackingVertex.h"
0002 #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticle.h"
0003 
0004 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
0005 
0006 #include <FWCore/MessageLogger/interface/MessageLogger.h>
0007 
0008 const unsigned int TrackingParticle::longLivedTag = 65536;
0009 
0010 TrackingParticle::TrackingParticle() {
0011   // No operation
0012 }
0013 
0014 TrackingParticle::TrackingParticle(const SimTrack& simtrk, const TrackingVertexRef& parentVertex) {
0015   addG4Track(simtrk);
0016   setParentVertex(parentVertex);
0017 }
0018 
0019 TrackingParticle::~TrackingParticle() {}
0020 
0021 void TrackingParticle::addGenParticle(const reco::GenParticleRef& ref) { genParticles_.push_back(ref); }
0022 
0023 void TrackingParticle::addG4Track(const SimTrack& t) { g4Tracks_.push_back(t); }
0024 
0025 TrackingParticle::genp_iterator TrackingParticle::genParticle_begin() const { return genParticles_.begin(); }
0026 
0027 TrackingParticle::genp_iterator TrackingParticle::genParticle_end() const { return genParticles_.end(); }
0028 
0029 TrackingParticle::g4t_iterator TrackingParticle::g4Track_begin() const { return g4Tracks_.begin(); }
0030 
0031 TrackingParticle::g4t_iterator TrackingParticle::g4Track_end() const { return g4Tracks_.end(); }
0032 
0033 void TrackingParticle::setParentVertex(const TrackingVertexRef& ref) { parentVertex_ = ref; }
0034 
0035 void TrackingParticle::addDecayVertex(const TrackingVertexRef& ref) { decayVertices_.push_back(ref); }
0036 
0037 void TrackingParticle::clearParentVertex() { parentVertex_ = TrackingVertexRef(); }
0038 
0039 void TrackingParticle::clearDecayVertices() { decayVertices_.clear(); }
0040 
0041 int TrackingParticle::matchedHit() const {
0042   edm::LogWarning("TrackingParticle")
0043       << "The method matchedHit() has been deprecated. Use numberOfTrackerLayers() instead.";
0044   return numberOfTrackerLayers_;
0045 }
0046 
0047 void TrackingParticle::setNumberOfHits(int numberOfHits) { numberOfHits_ = numberOfHits; }
0048 
0049 void TrackingParticle::setNumberOfTrackerHits(int numberOfTrackerHits) { numberOfTrackerHits_ = numberOfTrackerHits; }
0050 
0051 void TrackingParticle::setNumberOfTrackerLayers(const int numberOfTrackerLayers) {
0052   numberOfTrackerLayers_ = numberOfTrackerLayers;
0053 }
0054 
0055 std::ostream& operator<<(std::ostream& s, TrackingParticle const& tp) {
0056   s << "TP momentum, q, ID, & Event #: " << tp.p4() << " " << tp.charge() << " " << tp.pdgId() << " "
0057     << tp.eventId().bunchCrossing() << "." << tp.eventId().event() << std::endl;
0058 
0059   for (TrackingParticle::genp_iterator hepT = tp.genParticle_begin(); hepT != tp.genParticle_end(); ++hepT) {
0060     s << " HepMC Track Momentum " << (*hepT)->momentum().rho() << std::endl;
0061   }
0062 
0063   for (TrackingParticle::g4t_iterator g4T = tp.g4Track_begin(); g4T != tp.g4Track_end(); ++g4T) {
0064     s << " Geant Track Momentum  " << g4T->momentum() << std::endl;
0065     s << " Geant Track ID & type " << g4T->trackId() << " " << g4T->type() << std::endl;
0066     if (g4T->type() != tp.pdgId()) {
0067       s << " Mismatch b/t TrackingParticle and Geant types" << std::endl;
0068     }
0069   }
0070   // Loop over decay vertices
0071   s << " TP Vertex " << tp.vertex() << std::endl;
0072   s << " Source vertex: " << tp.parentVertex()->position() << std::endl;
0073   s << " " << tp.decayVertices().size() << " Decay vertices" << std::endl;
0074   for (tv_iterator iTV = tp.decayVertices_begin(); iTV != tp.decayVertices_end(); ++iTV) {
0075     s << " Decay vertices:      " << (**iTV).position() << std::endl;
0076   }
0077 
0078   return s;
0079 }