CMS 3D CMS Logo

SaveSimTrack.cc
Go to the documentation of this file.
2 
5 
7 
8 #include "G4PhysicalConstants.hh"
9 #include "G4SystemOfUnits.hh"
10 #include "G4Track.hh"
11 #include <algorithm>
12 
14  edm::ParameterSet ps = p.getParameter<edm::ParameterSet>("SaveSimTrack");
15  pdgs_ = ps.getUntrackedParameter<std::vector<int>>("PDGCodes");
16 
17  edm::LogInfo("SaveSimTrack") << "SaveSimTrack:: Save Sim Track if PDG code "
18  << "is one from the list of " << pdgs_.size() << " items";
19  for (unsigned int k = 0; k < pdgs_.size(); ++k)
20  edm::LogInfo("SaveSimTrack") << "[" << k << "] " << pdgs_[k];
21 }
22 
24 
26  G4Track *theTrack = (G4Track *)((*trk)());
27  TrackInformation *trkInfo = (TrackInformation *)(theTrack->GetUserInformation());
28  if (trkInfo) {
29  int pdg = theTrack->GetDefinition()->GetPDGEncoding();
30  if (std::find(pdgs_.begin(), pdgs_.end(), pdg) != pdgs_.end()) {
31  trkInfo->storeTrack(true);
32  LogDebug("SaveSimTrack") << "Save SimTrack the Track " << theTrack->GetTrackID() << " Type "
33  << theTrack->GetDefinition()->GetParticleName() << " Momentum "
34  << theTrack->GetMomentum() / MeV << " MeV/c";
35  }
36  }
37 }
#define LogDebug(id)
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
bool storeTrack() const
~SaveSimTrack() override
Definition: SaveSimTrack.cc:23
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
void update(const BeginOfTrack *trk) override
This routine will be called when the appropriate signal arrives.
Definition: SaveSimTrack.cc:25
const double MeV
std::vector< int > pdgs_
Definition: SaveSimTrack.h:18
SaveSimTrack(edm::ParameterSet const &p)
Definition: SaveSimTrack.cc:13
int k[5][pyjets_maxn]