Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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//
27// -----------------------------------------------------------------------------
28//
29// GEANT4 Class file
30//
31// File name: G4NuclearStopping
32//
33// Author: Vladimir Ivanchenko
34//
35// Creation date: 20 July 2009
36//
37// Modified:
38//
39// Warning: this class should be instantiated after G4ionIonisation
40//
41// -----------------------------------------------------------------------------
42//
43//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
44//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
45
46#include "G4NuclearStopping.hh"
48#include "G4SystemOfUnits.hh"
50
51//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
52
53using namespace std;
54
56 : G4VEmProcess(processName)
57{
58 isInitialized = false;
60 SetBuildTableFlag(false);
61 enableAlongStepDoIt = true;
62 enablePostStepDoIt = false;
63}
64
65//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
66
68{}
69
70//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
71
73{
74 return (p.GetPDGCharge() != 0.0 && !p.IsShortLived());
75}
76
77//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
78
80{
81 if(!isInitialized) {
82 isInitialized = true;
83
85 AddEmModel(1, EmModel());
87 EmModel()->SetParticleChange(&nParticleChange);
88 }
89}
90
91//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
92
95 const G4Track&, G4double, G4double, G4double&, G4GPILSelection* selection)
96{
97 *selection = NotCandidateForSelection;
98 return DBL_MAX;
99}
100
101//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
102
104 const G4Step& step)
105{
106 nParticleChange.InitializeForAlongStep(track);
107
108 // this line only valid if nuclear stopping
109 // is computed after G4ionIonisation process
110 nParticleChange.SetProposedCharge(step.GetPostStepPoint()->GetCharge());
111
113
114 const G4ParticleDefinition* part = track.GetParticleDefinition();
115 G4double Z = std::abs(part->GetPDGCharge()/eplus);
116 G4double massR = proton_mass_c2/part->GetPDGMass();
117
118 if(T2 > 0.0 && T2*massR < Z*Z*MeV) {
119
120 G4double length = step.GetStepLength();
121 if(length > 0.0) {
122
123 // primary
125 G4double T = 0.5*(T1 + T2);
126 const G4MaterialCutsCouple* couple = track.GetMaterialCutsCouple();
127 G4double Tscaled = T*massR;
128 G4VEmModel* mod = SelectModel(Tscaled, couple->GetIndex());
129
130 // sample stopping
131 if(mod->IsActive(Tscaled)) {
132 G4double nloss =
133 length*mod->ComputeDEDXPerVolume(couple->GetMaterial(), part, T);
134 nloss = std::min(nloss, T1);
135 nParticleChange.SetProposedKineticEnergy(T1 - nloss);
136 nParticleChange.ProposeLocalEnergyDeposit(nloss);
137 nParticleChange.ProposeNonIonizingEnergyDeposit(nloss);
138 }
139 }
140 }
141 return &nParticleChange;
142}
143
144//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
145
147{}
148
149//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
150
151void G4NuclearStopping::ProcessDescription(std::ostream& out) const
152{
153 out << " Nuclear stopping";
155}
156
157//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
@ fNuclearStopping
G4GPILSelection
@ NotCandidateForSelection
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
const G4Material * GetMaterial() const
G4bool IsApplicable(const G4ParticleDefinition &p) final
void InitialiseProcess(const G4ParticleDefinition *) final
virtual void PrintInfo() final
G4NuclearStopping(const G4String &processName="nuclearStopping")
G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection) final
G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &step) final
virtual void ProcessDescription(std::ostream &) const override
virtual ~G4NuclearStopping()
void InitializeForAlongStep(const G4Track &)
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedCharge(G4double theCharge)
G4double GetPDGCharge() const
G4double GetCharge() const
G4double GetKineticEnergy() const
Definition: G4Step.hh:62
G4StepPoint * GetPreStepPoint() const
G4double GetStepLength() const
G4StepPoint * GetPostStepPoint() const
const G4ParticleDefinition * GetParticleDefinition() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
Definition: G4VEmModel.cc:456
void SetActivationHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:771
G4bool IsActive(G4double kinEnergy) const
Definition: G4VEmModel.hh:785
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
Definition: G4VEmModel.cc:245
G4VEmModel * EmModel(size_t index=0) const
void SetBuildTableFlag(G4bool val)
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)
void SetEmModel(G4VEmModel *, G4int index=0)
virtual void ProcessDescription(std::ostream &outFile) const override
G4VEmModel * SelectModel(G4double kinEnergy, size_t index)
void ProposeNonIonizingEnergyDeposit(G4double anEnergyPart)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4bool enableAlongStepDoIt
Definition: G4VProcess.hh:360
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406
G4bool enablePostStepDoIt
Definition: G4VProcess.hh:361
#define DBL_MAX
Definition: templates.hh:62