#include <KillSecondariesTrackAction.h>
Public Member Functions | |
KillSecondariesTrackAction (edm::ParameterSet const &p) | |
void | update (const BeginOfTrack *trk) |
This routine will be called when the appropriate signal arrives. | |
~KillSecondariesTrackAction () | |
Private Attributes | |
bool | killHeavy |
double | kmaxIon |
double | kmaxNeutron |
double | kmaxProton |
Definition at line 10 of file KillSecondariesTrackAction.h.
KillSecondariesTrackAction::KillSecondariesTrackAction | ( | edm::ParameterSet const & | p | ) |
Definition at line 10 of file KillSecondariesTrackAction.cc.
References edm::ParameterSet::getParameter(), killHeavy, kmaxIon, kmaxNeutron, and kmaxProton.
{ killHeavy = p.getParameter<bool>("KillHeavy"); kmaxIon = p.getParameter<double>("IonThreshold")*MeV; kmaxProton = p.getParameter<double>("ProtonThreshold")*MeV; kmaxNeutron = p.getParameter<double>("NeutronThreshold")*MeV; edm::LogInfo("KillSecondaries") << "KillSecondariesTrackAction:: Killing" << " Flag " << killHeavy << " protons below " << kmaxProton << " MeV, neutrons below " << kmaxNeutron << " MeV and ions below " << kmaxIon << " MeV\n"; }
KillSecondariesTrackAction::~KillSecondariesTrackAction | ( | ) |
Definition at line 24 of file KillSecondariesTrackAction.cc.
{}
void KillSecondariesTrackAction::update | ( | const BeginOfTrack * | ) | [virtual] |
This routine will be called when the appropriate signal arrives.
Implements Observer< const BeginOfTrack * >.
Definition at line 26 of file KillSecondariesTrackAction.cc.
References TrackInformation::isPrimary(), ke, killHeavy, kmaxIon, kmaxNeutron, and kmaxProton.
{ if (killHeavy) { G4Track* theTrack = (G4Track*)((*trk)()); TrackInformation * trkInfo = (TrackInformation *)(theTrack->GetUserInformation()); if (trkInfo) { int pdg = theTrack->GetDefinition()->GetPDGEncoding(); if (!(trkInfo->isPrimary())) { // Only secondary particles double ke = theTrack->GetKineticEnergy()/MeV; if ((((pdg/1000000000 == 1 && ((pdg/10000)%100) > 0 && ((pdg/10)%100) > 0)) && (ke<kmaxIon)) || ((pdg == 2212) && (ke < kmaxProton)) || ((pdg == 2112) && (ke < kmaxNeutron))) { theTrack->SetTrackStatus(fStopAndKill); edm::LogInfo("KillSecondaries") << "Kill Track " << theTrack->GetTrackID() << " Type " << theTrack->GetDefinition()->GetParticleName() << " Kinetic Energy " << ke <<" MeV"; } } } } }
bool KillSecondariesTrackAction::killHeavy [private] |
Definition at line 19 of file KillSecondariesTrackAction.h.
Referenced by KillSecondariesTrackAction(), and update().
double KillSecondariesTrackAction::kmaxIon [private] |
Definition at line 20 of file KillSecondariesTrackAction.h.
Referenced by KillSecondariesTrackAction(), and update().
double KillSecondariesTrackAction::kmaxNeutron [private] |
Definition at line 20 of file KillSecondariesTrackAction.h.
Referenced by KillSecondariesTrackAction(), and update().
double KillSecondariesTrackAction::kmaxProton [private] |
Definition at line 20 of file KillSecondariesTrackAction.h.
Referenced by KillSecondariesTrackAction(), and update().