CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
KillSecondariesTrackAction.cc
Go to the documentation of this file.
2 
5 
7 
8 #include "G4Track.hh"
9 
11 
12  killHeavy = p.getParameter<bool>("KillHeavy");
13  kmaxIon = p.getParameter<double>("IonThreshold")*MeV;
14  kmaxProton = p.getParameter<double>("ProtonThreshold")*MeV;
15  kmaxNeutron = p.getParameter<double>("NeutronThreshold")*MeV;
16 
17  edm::LogInfo("KillSecondaries") << "KillSecondariesTrackAction:: Killing"
18  << " Flag " << killHeavy << " protons below "
19  << kmaxProton << " MeV, neutrons below "
20  << kmaxNeutron << " MeV and ions below "
21  << kmaxIon << " MeV\n";
22 }
23 
25 
27 
28  if (killHeavy) {
29  G4Track* theTrack = (G4Track*)((*trk)());
30  TrackInformation * trkInfo = (TrackInformation *)(theTrack->GetUserInformation());
31  if (trkInfo) {
32  int pdg = theTrack->GetDefinition()->GetPDGEncoding();
33  if (!(trkInfo->isPrimary())) { // Only secondary particles
34  double ke = theTrack->GetKineticEnergy()/MeV;
35  if ((((pdg/1000000000 == 1 && ((pdg/10000)%100) > 0 &&
36  ((pdg/10)%100) > 0)) && (ke<kmaxIon)) ||
37  ((pdg == 2212) && (ke < kmaxProton)) ||
38  ((pdg == 2112) && (ke < kmaxNeutron))) {
39  theTrack->SetTrackStatus(fStopAndKill);
40  edm::LogInfo("KillSecondaries") << "Kill Track " << theTrack->GetTrackID()
41  << " Type " << theTrack->GetDefinition()->GetParticleName()
42  << " Kinetic Energy " << ke <<" MeV";
43  }
44  }
45  }
46  }
47 }
48 
T getParameter(std::string const &) const
KillSecondariesTrackAction(edm::ParameterSet const &p)
void update(const BeginOfTrack *trk)
This routine will be called when the appropriate signal arrives.
int ke
bool isPrimary() const