Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4NuclearStopping.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 // $Id: G4NuclearStopping.cc 105121 2017-07-13 13:38:25Z gcosmo $
27 //
28 // -----------------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 // File name: G4NuclearStopping
33 //
34 // Author: Vladimir Ivanchenko
35 //
36 // Creation date: 20 July 2009
37 //
38 // Modified:
39 //
40 // Warning: this class should be instantiated after G4ionIonisation
41 //
42 // -----------------------------------------------------------------------------
43 //
44 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
45 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
46 
47 #include "G4NuclearStopping.hh"
49 #include "G4SystemOfUnits.hh"
50 #include "G4PhysicalConstants.hh"
51 
52 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
53 
54 using namespace std;
55 
57  : G4VEmProcess(processName)
58 {
59  isInitialized = false;
61  SetBuildTableFlag(false);
62  enableAlongStepDoIt = true;
63  enablePostStepDoIt = false;
64 }
65 
66 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
67 
69 {}
70 
71 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
72 
74 {
75  return (p.GetPDGCharge() != 0.0 && !p.IsShortLived());
76 }
77 
78 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
79 
81 {
82  if(!isInitialized) {
83  isInitialized = true;
84 
85  if(!EmModel(1)) {
87  }
88  AddEmModel(1, EmModel());
89  EmModel()->SetParticleChange(&nParticleChange);
90  }
91 }
92 
93 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
94 
95 G4double
97  const G4Track&, G4double, G4double, G4double&, G4GPILSelection* selection)
98 {
99  *selection = NotCandidateForSelection;
100  return DBL_MAX;
101 }
102 
103 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
104 
106  const G4Step& step)
107 {
108  nParticleChange.InitializeForAlongStep(track);
109 
110  // this line only valid if nuclear stopping
111  // is computed after G4ionIonisation process
112  nParticleChange.SetProposedCharge(step.GetPostStepPoint()->GetCharge());
113 
115 
116  const G4ParticleDefinition* part = track.GetParticleDefinition();
117  G4double Z = std::abs(part->GetPDGCharge()/eplus);
118 
119  if(T2 > 0.0 && T2*proton_mass_c2 < Z*Z*MeV*part->GetPDGMass()) {
120 
121  G4double length = step.GetStepLength();
122  if(length > 0.0) {
123 
124  // primary
126  G4double T = 0.5*(T1 + T2);
127  const G4MaterialCutsCouple* couple = track.GetMaterialCutsCouple();
128  G4VEmModel* mod = SelectModel(T, couple->GetIndex());
129 
130  // sample stopping
131  G4double nloss =
132  length*mod->ComputeDEDXPerVolume(couple->GetMaterial(), part, T);
133  if(nloss > T1) { nloss = T1; }
134  nParticleChange.SetProposedKineticEnergy(T1 - nloss);
135  nParticleChange.ProposeLocalEnergyDeposit(nloss);
136  nParticleChange.ProposeNonIonizingEnergyDeposit(nloss);
137  }
138  }
139  return &nParticleChange;
140 }
141 
142 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
143 
145 {}
146 
147 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
148 
G4NuclearStopping(const G4String &processName="nuclearStopping")
G4VEmModel * SelectModel(G4double &kinEnergy, size_t index)
G4double GetStepLength() const
void SetBuildTableFlag(G4bool val)
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
Definition: G4VEmModel.cc:418
G4VEmModel * EmModel(G4int index=1) const
const char * p
Definition: xmltok.h:285
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
void InitializeForAlongStep(const G4Track &)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4StepPoint * GetPreStepPoint() const
virtual void PrintInfo() final
void SetEmModel(G4VEmModel *, G4int index=1)
void ProposeNonIonizingEnergyDeposit(G4double anEnergyPart)
bool G4bool
Definition: G4Types.hh:79
static constexpr double eplus
Definition: G4SIunits.hh:199
G4double GetCharge() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
const G4ParticleDefinition * GetParticleDefinition() const
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
Definition: G4Step.hh:76
G4bool enablePostStepDoIt
Definition: G4VProcess.hh:352
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)
G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection) final
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
Definition: G4VEmModel.cc:248
G4bool IsApplicable(const G4ParticleDefinition &p) final
G4StepPoint * GetPostStepPoint() const
void SetProposedCharge(G4double theCharge)
G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &step) final
G4bool enableAlongStepDoIt
Definition: G4VProcess.hh:351
G4double GetKineticEnergy() const
double G4double
Definition: G4Types.hh:76
void InitialiseProcess(const G4ParticleDefinition *) final
G4double GetPDGCharge() const
#define DBL_MAX
Definition: templates.hh:83
G4GPILSelection
const G4Material * GetMaterial() const
virtual ~G4NuclearStopping()