Geant4  10.00.p03
G4ElectronCapture.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // This example is provided by the Geant4-DNA collaboration
27 // Any report or published results obtained using the Geant4-DNA software
28 // shall cite the following Geant4-DNA collaboration publication:
29 // Med. Phys. 37 (2010) 4692-4708
30 // The Geant4-DNA web site is available at http://geant4-dna.org
31 //
32 //---------------------------------------------------------------------------
33 //
34 // ClassName: G4ElectronCapture
35 //
36 // Description: The process to kill particles to save CPU
37 //
38 // Author: V.Ivanchenko 31 August 2010
39 //
40 //----------------------------------------------------------------------------
41 //
42 
43 #include "G4ElectronCapture.hh"
44 #include "G4SystemOfUnits.hh"
45 #include "G4ParticleDefinition.hh"
46 #include "G4Step.hh"
47 #include "G4Track.hh"
48 #include "G4Region.hh"
49 #include "G4RegionStore.hh"
50 #include "G4Electron.hh"
51 
52 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
53 
55  : G4VDiscreteProcess("eCapture", fElectromagnetic), kinEnergyThreshold(ekinlim),
56  regionName(regName), region(0)
57 {
58  if(regName == "" || regName == "world") {
59  regionName = "DefaultRegionForTheWorld";
60  }
62 }
63 
64 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
65 
67 {}
68 
69 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
70 
72 {
73  kinEnergyThreshold = val;
74  if(verboseLevel > 0) {
75  G4cout << "### G4ElectronCapture: Tracking cut E(MeV) = "
77  }
78 }
79 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
80 
82 {
84  if(region && verboseLevel > 0) {
85  G4cout << "### G4ElectronCapture: Tracking cut E(MeV) = "
86  << kinEnergyThreshold/MeV << " is assigned to " << regionName
87  << G4endl;
88  }
89 }
90 
91 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
92 
94 {
95  return true;
96 }
97 
98 G4double
100  G4double,
102 {
103  // condition is set to "Not Forced"
104  *condition = NotForced;
105 
106  G4double limit = DBL_MAX;
107  if(region) {
108  if(aTrack.GetVolume()->GetLogicalVolume()->GetRegion() == region &&
109  aTrack.GetKineticEnergy() < kinEnergyThreshold) { limit = 0.0; }
110  }
111  return limit;
112 }
113 
114 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
115 
117  const G4Step&)
118 {
119  pParticleChange->Initialize(aTrack);
123  return pParticleChange;
124 }
125 
126 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
127 
130 {
131  return DBL_MAX;
132 }
133 
134 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
135 
136 
G4double condition(const G4ErrorSymMatrix &m)
virtual void Initialize(const G4Track &)
static const double MeV
Definition: G4SIunits.hh:193
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
G4int verboseLevel
Definition: G4VProcess.hh:368
virtual ~G4ElectronCapture()
G4Region * GetRegion() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
static G4RegionStore * GetInstance()
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
virtual G4double GetMeanFreePath(const G4Track &, G4double, G4ForceCondition *)
bool G4bool
Definition: G4Types.hh:79
G4ElectronCapture(const G4String &regName, G4double ekinlimit)
Definition: G4Step.hh:76
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
G4LogicalVolume * GetLogicalVolume() const
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
G4VPhysicalVolume * GetVolume() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
G4ParticleChangeForGamma fParticleChange
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
G4ForceCondition
#define DBL_MAX
Definition: templates.hh:83
virtual G4bool IsApplicable(const G4ParticleDefinition &)
void SetKinEnergyLimit(G4double)