Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4AdjointForcedInteractionForGamma.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
28
29#include "G4AdjointCSManager.hh"
30#include "G4AdjointGamma.hh"
32#include "G4ParticleChange.hh"
33#include "G4SystemOfUnits.hh"
34#include "G4VEmAdjointModel.hh"
35
37 G4String process_name)
38 : G4VContinuousDiscreteProcess(process_name)
39 , fAdjointComptonModel(nullptr)
40 , fAdjointBremModel(nullptr)
41{
43 fParticleChange = new G4ParticleChange();
44}
45
46//////////////////////////////////////////////////////////////////////////////
48{
49 if(fParticleChange)
50 delete fParticleChange;
51}
52
53//////////////////////////////////////////////////////////////////////////////
55 std::ostream& out) const
56{
57 out << "Forced interaction for gamma.\n";
58}
59
60//////////////////////////////////////////////////////////////////////////////
63{
64 fCSManager->BuildCrossSectionMatrices(); // it will be done just once
65 fCSManager->BuildTotalSigmaTables();
66}
67
68// Note on weight correction for forced interaction.
69// For the forced interaction applied here we use a truncated exponential law
70// for the probability of survival over a fixed total length. This is done by
71// using a linear transformation of the non-biased probability survival. In
72// math this is written P'(x)=C1P(x)+C2 , with P(x)=exp(-sum(sigma_ixi)) . x and
73// L can cross different volumes with different cross section sigma. For forced
74// interaction, we get the limit conditions:
75// P'(L)=0 and P'(0)=1 (L can be used over different volumes)
76// From simple solving of linear equations we
77// get C1=1/(1-P(L)) and C2=-P(L)/(1-P(L))
78// P'(x)=(P(x)-P(L))/(1-P(L))
79// For the probability over a step x1 to x2, P'(x1->x2)=P'(x2)/P'(x1).
80// The effective cross
81// section is defined -d(P'(x))/dx/P'(x).
82// We get therefore
83// sigma_eff = C1sigmaP(x)/(C1P(x)+C2) = sigmaP(x)/(P(x)+C2/C1)
84// = sigmaP(x)/(P(x)-P(L)) = sigma/(1-P(L)/P(x))
85//////////////////////////////////////////////////////////////////////////////
87 const G4Track& track, const G4Step&)
88{
89 fParticleChange->Initialize(track);
90 // For the free flight gamma no interaction occurs but a gamma with same
91 // properties is produced for further forced interaction. It is done at the
92 // very beginning of the track so that the weight can be the same
93 if(fCopyGammaForForced)
94 {
95 G4ThreeVector theGammaMomentum = track.GetMomentum();
96 fParticleChange->AddSecondary(
97 new G4DynamicParticle(G4AdjointGamma::AdjointGamma(), theGammaMomentum));
98 fParticleChange->SetParentWeightByProcess(false);
99 fParticleChange->SetSecondaryWeightByProcess(false);
100 }
101 else
102 { // Occurrence of forced interaction
103 // Selection of the model to be called
104 G4VEmAdjointModel* theSelectedModel = nullptr;
105 G4bool is_scat_proj_to_proj_case = false;
106 G4double factor=1.;
107 if(!fAdjointComptonModel && !fAdjointBremModel)
108 return fParticleChange;
109 if(!fAdjointComptonModel)
110 {
111 theSelectedModel = fAdjointBremModel;
112 is_scat_proj_to_proj_case = false;
113 // This is needed because the results of it will be used in the post step
114 // do it weight correction inside the model
115 fAdjointBremModel->AdjointCrossSection(track.GetMaterialCutsCouple(),
116 track.GetKineticEnergy(), false);
117 }
118 else if(!fAdjointBremModel)
119 {
120 theSelectedModel = fAdjointComptonModel;
121 is_scat_proj_to_proj_case = true;
122 }
123 else
124 { // Choose the model according to a 50-50 % probability
125 G4double bremAdjCS = fAdjointBremModel->AdjointCrossSection(
126 track.GetMaterialCutsCouple(), track.GetKineticEnergy(), false);
127 if(G4UniformRand() < 0.5)
128 {
129 theSelectedModel = fAdjointBremModel;
130 is_scat_proj_to_proj_case = false;
131 factor=bremAdjCS/fLastAdjCS/0.5;
132 }
133 else
134 {
135 theSelectedModel = fAdjointComptonModel;
136 is_scat_proj_to_proj_case = true;
137 factor=(fLastAdjCS-bremAdjCS)/fLastAdjCS/0.5;
138 }
139 }
140
141 // Compute the weight correction factor
142 G4double invEffectiveAdjointCS =
143 (1. - std::exp(fNbAdjIntLength - fTotNbAdjIntLength)) / fLastAdjCS/fCSBias;
144
145 // Call the selected model without correction of the weight in the model
146 theSelectedModel->SetCorrectWeightForPostStepInModel(false);
147 theSelectedModel
149 factor*fLastAdjCS * invEffectiveAdjointCS);
150 theSelectedModel->SampleSecondaries(track, is_scat_proj_to_proj_case,
151 fParticleChange);
152 theSelectedModel->SetCorrectWeightForPostStepInModel(true);
153
154 fContinueGammaAsNewFreeFlight = true;
155 }
156 return fParticleChange;
157}
158
159//////////////////////////////////////////////////////////////////////////////
161 const G4Track& track, const G4Step&)
162{
163 fParticleChange->Initialize(track);
164 // Compute nb of interactions length over step length
166 G4double stepLength = track.GetStep()->GetStepLength();
167 G4double ekin = track.GetKineticEnergy();
168 fLastAdjCS = fCSManager->GetTotalAdjointCS(track.GetDefinition(), ekin,
169 track.GetMaterialCutsCouple());
170 G4double nb_fwd_interaction_length_over_step =
171 stepLength * fCSManager->GetTotalForwardCS(G4AdjointGamma::AdjointGamma(),
172 ekin,
173 track.GetMaterialCutsCouple());
174
175 G4double nb_adj_interaction_length_over_step = stepLength * fLastAdjCS;
176 G4double fwd_survival_probability =
177 std::exp(-nb_fwd_interaction_length_over_step);
178 G4double mc_induced_survival_probability = 1.;
179
180 if(fFreeFlightGamma)
181 { // for free_flight survival probability stays 1
182 // Accumulate the number of interaction lengths during free flight of gamma
183 fTotNbAdjIntLength += nb_adj_interaction_length_over_step;
184 fAccTrackLength += stepLength;
185 }
186 else
187 {
188 G4double previous_acc_nb_adj_interaction_length = fNbAdjIntLength;
189 fNbAdjIntLength += fCSBias*nb_adj_interaction_length_over_step;
190 theNumberOfInteractionLengthLeft -= fCSBias*nb_adj_interaction_length_over_step;
191
192 // protection against rare race condition
193 if(std::abs(fTotNbAdjIntLength - previous_acc_nb_adj_interaction_length) <=
194 1.e-15)
195 {
196 mc_induced_survival_probability = 1.e50;
197 }
198 else
199 {
200 mc_induced_survival_probability =
201 std::exp(-fNbAdjIntLength) - std::exp(-fTotNbAdjIntLength);
202 mc_induced_survival_probability /=
203 (std::exp(-previous_acc_nb_adj_interaction_length) -
204 std::exp(-fTotNbAdjIntLength));
205 }
206 }
207 G4double weight_correction =
208 fwd_survival_probability / mc_induced_survival_probability;
209
210 // Caution!!!
211 // It is important to select the weight of the post_step_point as the
212 // current weight and not the weight of the track, as the weight of the track
213 // is changed after having applied all the along_step_do_it.
214 G4double new_weight =
215 weight_correction * track.GetStep()->GetPostStepPoint()->GetWeight();
216
217 fParticleChange->SetParentWeightByProcess(false);
218 fParticleChange->SetSecondaryWeightByProcess(false);
219 fParticleChange->ProposeParentWeight(new_weight);
220
221 return fParticleChange;
222}
223
224//////////////////////////////////////////////////////////////////////////////
228{
229 static G4int lastFreeFlightTrackId = 1000;
230 G4int step_id = track.GetCurrentStepNumber();
232 fCopyGammaForForced = false;
233 G4int track_id = track.GetTrackID();
234 fFreeFlightGamma =
235 (track_id != lastFreeFlightTrackId + 1 || fContinueGammaAsNewFreeFlight);
236 if(fFreeFlightGamma)
237 {
238 if(step_id == 1 || fContinueGammaAsNewFreeFlight)
239 {
240 *condition = Forced;
241 // A gamma with same conditions will be generate at next post_step do it
242 // for the forced interaction
243 fCopyGammaForForced = true;
244 lastFreeFlightTrackId = track_id;
245 fAccTrackLength = 0.;
246 fTotNbAdjIntLength = 0.;
247 fContinueGammaAsNewFreeFlight = false;
248 return 1.e-90;
249 }
250 else
251 {
252 return DBL_MAX;
253 }
254 }
255 else
256 { // compute the interaction length for forced interaction
257 if(step_id == 1)
258 {
259 fCSBias=0.000001/fTotNbAdjIntLength;
260 fTotNbAdjIntLength*=fCSBias;
261 G4double min_val = std::exp(-fTotNbAdjIntLength);
263 -std::log(min_val + G4UniformRand() * (1. - min_val));
265 fNbAdjIntLength = 0.;
266 }
267 G4VPhysicalVolume* thePostPhysVolume =
269 G4double ekin = track.GetKineticEnergy();
270 G4double postCS = 0.;
271 if(thePostPhysVolume)
272 {
273 postCS = fCSManager->GetTotalAdjointCS(
275 thePostPhysVolume->GetLogicalVolume()->GetMaterialCutsCouple());
276 }
277 if(postCS > 0.)
278 return theNumberOfInteractionLengthLeft / postCS /fCSBias;
279 else
280 return DBL_MAX;
281 }
282}
283
284////////////////////////////////////////////////////////////////////////////////
287{
288 return DBL_MAX;
289}
290
291////////////////////////////////////////////////////////////////////////////////
292// Not used in this process but should be implemented as virtual method
294 G4double,
296{
297 return 0.;
298}
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
@ NotForced
@ Forced
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
G4double GetTotalForwardCS(G4ParticleDefinition *aPartDef, G4double Ekin, const G4MaterialCutsCouple *aCouple)
static G4AdjointCSManager * GetAdjointCSManager()
G4double GetTotalAdjointCS(G4ParticleDefinition *aPartDef, G4double Ekin, const G4MaterialCutsCouple *aCouple)
G4double GetContinuousStepLimit(const G4Track &aTrack, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety) override
G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &step) override
void ProcessDescription(std::ostream &) const override
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
void BuildPhysicsTable(const G4ParticleDefinition &) override
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &) override
static G4AdjointGamma * AdjointGamma()
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
void AddSecondary(G4Track *aSecondary)
void Initialize(const G4Track &) override
G4VPhysicalVolume * GetPhysicalVolume() const
G4double GetWeight() const
Definition: G4Step.hh:62
G4StepPoint * GetPreStepPoint() const
G4double GetStepLength() const
G4StepPoint * GetPostStepPoint() const
G4int GetTrackID() const
const G4ThreeVector & GetPosition() const
G4int GetCurrentStepNumber() const
G4ThreeVector GetMomentum() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
const G4Step * GetStep() const
virtual void SampleSecondaries(const G4Track &aTrack, G4bool isScatProjToProj, G4ParticleChange *fParticleChange)=0
void SetAdditionalWeightCorrectionFactorForPostStepOutsideModel(G4double factor)
virtual G4double AdjointCrossSection(const G4MaterialCutsCouple *aCouple, G4double primEnergy, G4bool isScatProjToProj)
void SetCorrectWeightForPostStepInModel(G4bool aBool)
void SetSecondaryWeightByProcess(G4bool)
void SetParentWeightByProcess(G4bool)
void ProposeParentWeight(G4double finalWeight)
G4LogicalVolume * GetLogicalVolume() const
G4double theInitialNumberOfInteractionLength
Definition: G4VProcess.hh:342
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:335
#define DBL_MAX
Definition: templates.hh:62