SimG4Core
SaveSimTrackAction
src
SaveSimTrack.cc
Go to the documentation of this file.
1
#include "
SimG4Core/SaveSimTrackAction/interface/SaveSimTrack.h
"
2
3
#include "
SimG4Core/Notification/interface/BeginOfTrack.h
"
4
#include "
SimG4Core/Notification/interface/TrackInformation.h
"
5
6
#include "
FWCore/MessageLogger/interface/MessageLogger.h
"
7
8
#include "G4PhysicalConstants.hh"
9
#include "G4SystemOfUnits.hh"
10
#include "G4Track.hh"
11
#include <algorithm>
12
13
SaveSimTrack::SaveSimTrack
(
edm::ParameterSet
const
&
p
) {
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
23
SaveSimTrack::~SaveSimTrack
() {}
24
25
void
SaveSimTrack::update
(
const
BeginOfTrack
*trk) {
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
}
SaveSimTrack::SaveSimTrack
SaveSimTrack(edm::ParameterSet const &p)
Definition:
SaveSimTrack.cc:13
MessageLogger.h
AlCaHLTBitMon_ParallelJobs.p
p
Definition:
AlCaHLTBitMon_ParallelJobs.py:153
edm::LogInfo
Definition:
MessageLogger.h:254
dileptonTrigSettings_cff.pdg
pdg
Definition:
dileptonTrigSettings_cff.py:6
MeV
const double MeV
edm::ParameterSet::getUntrackedParameter
T getUntrackedParameter(std::string const &, T const &) const
spr::find
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition:
FindCaloHit.cc:19
SaveSimTrack::update
void update(const BeginOfTrack *trk) override
This routine will be called when the appropriate signal arrives.
Definition:
SaveSimTrack.cc:25
BeginOfTrack.h
SaveSimTrack::~SaveSimTrack
~SaveSimTrack() override
Definition:
SaveSimTrack.cc:23
SaveSimTrack::pdgs_
std::vector< int > pdgs_
Definition:
SaveSimTrack.h:18
BeginOfTrack
Definition:
BeginOfTrack.h:6
dqmdumpme.k
k
Definition:
dqmdumpme.py:60
LogDebug
#define LogDebug(id)
Definition:
MessageLogger.h:670
edm::ParameterSet
Definition:
ParameterSet.h:36
TrackInformation
Definition:
TrackInformation.h:8
TrackInformation::storeTrack
bool storeTrack() const
Definition:
TrackInformation.h:14
TrackInformation.h
SaveSimTrack.h
pdg
Definition:
pdg_functions.h:28
Generated for CMSSW Reference Manual by
1.8.16