Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4AdjointPhotoElectricModel.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$
27//
29#include "G4AdjointCSManager.hh"
30
32#include "G4Integrator.hh"
33#include "G4TrackStatus.hh"
34#include "G4ParticleChange.hh"
35#include "G4AdjointElectron.hh"
36#include "G4Gamma.hh"
37#include "G4AdjointGamma.hh"
38
39
40////////////////////////////////////////////////////////////////////////////////
41//
43 G4VEmAdjointModel("AdjointPEEffect")
44
45{ SetUseMatrix(false);
46 SetApplyCutInRange(false);
47
48 //Initialization
49 current_eEnergy =0.;
50 totAdjointCS=0.;
51 factorCSBiasing =1.;
52 post_step_AdjointCS =0.;
53 pre_step_AdjointCS =0.;
54 totBiasedAdjointCS =0.;
55
56 index_element=0;
57
62 theDirectPEEffectModel = new G4PEEffectModel();
63}
64////////////////////////////////////////////////////////////////////////////////
65//
67{;}
68
69////////////////////////////////////////////////////////////////////////////////
70//
72 G4bool IsScatProjToProjCase,
73 G4ParticleChange* fParticleChange)
74{ if (IsScatProjToProjCase) return ;
75
76 //Compute the totAdjointCS vectors if not already done for the current couple and electron energy
77 //-----------------------------------------------------------------------------------------------
78 const G4MaterialCutsCouple* aCouple = aTrack.GetMaterialCutsCouple();
79 const G4DynamicParticle* aDynPart = aTrack.GetDynamicParticle() ;
80 G4double electronEnergy = aDynPart->GetKineticEnergy();
81 G4ThreeVector electronDirection= aDynPart->GetMomentumDirection() ;
82 pre_step_AdjointCS = totAdjointCS; //The last computed CS was at pre step point
83 post_step_AdjointCS = AdjointCrossSection(aCouple, electronEnergy,IsScatProjToProjCase);
84 post_step_AdjointCS = totAdjointCS;
85
86
87
88
89 //Sample element
90 //-------------
91 const G4ElementVector* theElementVector = currentMaterial->GetElementVector();
92 size_t nelm = currentMaterial->GetNumberOfElements();
93 G4double rand_CS= G4UniformRand()*xsec[nelm-1];
94 for (index_element=0; index_element<nelm-1; index_element++){
95 if (rand_CS<xsec[index_element]) break;
96 }
97
98 //Sample shell and binding energy
99 //-------------
100 G4int nShells = (*theElementVector)[index_element]->GetNbOfAtomicShells();
101 rand_CS= shell_prob[index_element][nShells-1]*G4UniformRand();
102 G4int i = 0;
103 for (i=0; i<nShells-1; i++){
104 if (rand_CS<shell_prob[index_element][i]) break;
105 }
106 G4double gammaEnergy= electronEnergy+(*theElementVector)[index_element]->GetAtomicShell(i);
107
108 //Sample cos theta
109 //Copy of the G4PEEffectModel cos theta sampling method ElecCosThetaDistribution.
110 //This method cannot be used directly from G4PEEffectModel because it is a friend method. I should ask Vladimir to change that
111 //------------------------------------------------------------------------------------------------
112 //G4double cos_theta = theDirectPEEffectModel->ElecCosThetaDistribution(electronEnergy);
113
114 G4double cos_theta = 1.;
115 G4double gamma = 1. + electronEnergy/electron_mass_c2;
116 if (gamma <= 5.) {
117 G4double beta = std::sqrt(gamma*gamma-1.)/gamma;
118 G4double b = 0.5*gamma*(gamma-1.)*(gamma-2);
119
120 G4double rndm,term,greject,grejsup;
121 if (gamma < 2.) grejsup = gamma*gamma*(1.+b-beta*b);
122 else grejsup = gamma*gamma*(1.+b+beta*b);
123
124 do { rndm = 1.-2*G4UniformRand();
125 cos_theta = (rndm+beta)/(rndm*beta+1.);
126 term = 1.-beta*cos_theta;
127 greject = (1.-cos_theta*cos_theta)*(1.+b*term)/(term*term);
128 } while(greject < G4UniformRand()*grejsup);
129 }
130
131 // direction of the adjoint gamma electron
132 //---------------------------------------
133
134
135 G4double sin_theta = std::sqrt(1.-cos_theta*cos_theta);
136 G4double Phi = twopi * G4UniformRand();
137 G4double dirx = sin_theta*std::cos(Phi),diry = sin_theta*std::sin(Phi),dirz = cos_theta;
138 G4ThreeVector adjoint_gammaDirection(dirx,diry,dirz);
139 adjoint_gammaDirection.rotateUz(electronDirection);
140
141
142
143 //Weight correction
144 //-----------------------
145 CorrectPostStepWeight(fParticleChange, aTrack.GetWeight(), electronEnergy,gammaEnergy,IsScatProjToProjCase);
146
147
148
149 //Create secondary and modify fParticleChange
150 //--------------------------------------------
151 G4DynamicParticle* anAdjointGamma = new G4DynamicParticle (
152 G4AdjointGamma::AdjointGamma(),adjoint_gammaDirection, gammaEnergy);
153
154
155
156
157
158 fParticleChange->ProposeTrackStatus(fStopAndKill);
159 fParticleChange->AddSecondary(anAdjointGamma);
160
161
162
163
164}
165
166////////////////////////////////////////////////////////////////////////////////
167//
169 G4double old_weight,
170 G4double adjointPrimKinEnergy,
171 G4double projectileKinEnergy ,
172 G4bool )
173{
174 G4double new_weight=old_weight;
175
177 w_corr*=post_step_AdjointCS/pre_step_AdjointCS;
178 new_weight*=w_corr;
179 new_weight*=projectileKinEnergy/adjointPrimKinEnergy;
180 fParticleChange->SetParentWeightByProcess(false);
181 fParticleChange->SetSecondaryWeightByProcess(false);
182 fParticleChange->ProposeParentWeight(new_weight);
183}
184
185////////////////////////////////////////////////////////////////////////////////
186//
187
189 G4double electronEnergy,
190 G4bool IsScatProjToProjCase)
191{
192
193
194 if (IsScatProjToProjCase) return 0.;
195
196
197 if (aCouple !=currentCouple || current_eEnergy !=electronEnergy) {
198 totAdjointCS = 0.;
199 DefineCurrentMaterialAndElectronEnergy(aCouple, electronEnergy);
200 const G4ElementVector* theElementVector = currentMaterial->GetElementVector();
201 const double* theAtomNumDensityVector = currentMaterial->GetVecNbOfAtomsPerVolume();
202 size_t nelm = currentMaterial->GetNumberOfElements();
203 for (index_element=0;index_element<nelm;index_element++){
204
205 totAdjointCS +=AdjointCrossSectionPerAtom((*theElementVector)[index_element],electronEnergy)*theAtomNumDensityVector[index_element];
206 xsec[index_element] = totAdjointCS;
207 }
208
209 totBiasedAdjointCS=std::min(totAdjointCS,0.01);
210// totBiasedAdjointCS=totAdjointCS;
211 factorCSBiasing = totBiasedAdjointCS/totAdjointCS;
212 lastCS=totBiasedAdjointCS;
213
214
215 }
216 return totBiasedAdjointCS;
217
218
219}
220////////////////////////////////////////////////////////////////////////////////
221//
222
224 G4double electronEnergy,
225 G4bool IsScatProjToProjCase)
226{ return AdjointCrossSection(aCouple,electronEnergy,IsScatProjToProjCase);
227}
228////////////////////////////////////////////////////////////////////////////////
229//
230
232{
233 G4int nShells = anElement->GetNbOfAtomicShells();
234 G4double Z= anElement->GetZ();
235 G4int i = 0;
236 G4double B0=anElement->GetAtomicShell(0);
237 G4double gammaEnergy = electronEnergy+B0;
238 G4double CS= theDirectPEEffectModel->ComputeCrossSectionPerAtom(G4Gamma::Gamma(),gammaEnergy,Z,0.,0.,0.);
239 G4double adjointCS =0.;
240 if (CS >0) adjointCS += CS/gammaEnergy;
241 shell_prob[index_element][0] = adjointCS;
242 for (i=1;i<nShells;i++){
243 //G4cout<<i<<G4endl;
244 G4double Bi_= anElement->GetAtomicShell(i-1);
245 G4double Bi = anElement->GetAtomicShell(i);
246 //G4cout<<Bi_<<'\t'<<Bi<<G4endl;
247 if (electronEnergy <Bi_-Bi) {
248 gammaEnergy = electronEnergy+Bi;
249
250 CS=theDirectPEEffectModel->ComputeCrossSectionPerAtom(G4Gamma::Gamma(),gammaEnergy,Z,0.,0.,0.);
251 if (CS>0) adjointCS +=CS/gammaEnergy;
252 }
253 shell_prob[index_element][i] = adjointCS;
254
255 }
256 adjointCS*=electronEnergy;
257 return adjointCS;
258
259}
260////////////////////////////////////////////////////////////////////////////////
261//
262
263void G4AdjointPhotoElectricModel::DefineCurrentMaterialAndElectronEnergy(const G4MaterialCutsCouple* couple, G4double anEnergy)
264{ currentCouple = const_cast<G4MaterialCutsCouple*> (couple);
265 currentMaterial = const_cast<G4Material*> (couple->GetMaterial());
266 currentCoupleIndex = couple->GetIndex();
268 current_eEnergy = anEnergy;
269}
std::vector< G4Element * > G4ElementVector
@ fStopAndKill
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4UniformRand()
Definition: Randomize.hh:53
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
G4double GetPostStepWeightCorrection()
static G4AdjointCSManager * GetAdjointCSManager()
static G4AdjointElectron * AdjointElectron()
static G4AdjointGamma * AdjointGamma()
virtual void SampleSecondaries(const G4Track &aTrack, G4bool IsScatProjToProjCase, G4ParticleChange *fParticleChange)
virtual void CorrectPostStepWeight(G4ParticleChange *fParticleChange, G4double old_weight, G4double adjointPrimKinEnergy, G4double projectileKinEnergy, G4bool IsScatProjToProjCase)
virtual G4double AdjointCrossSection(const G4MaterialCutsCouple *aCouple, G4double primEnergy, G4bool IsScatProjToProjCase)
G4double AdjointCrossSectionPerAtom(const G4Element *anElement, G4double electronEnergy)
virtual G4double GetAdjointCrossSection(const G4MaterialCutsCouple *aCouple, G4double primEnergy, G4bool IsScatProjToProjCase)
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4double GetZ() const
Definition: G4Element.hh:131
G4int GetNbOfAtomicShells() const
Definition: G4Element.hh:146
G4double GetAtomicShell(G4int index) const
Definition: G4Element.cc:367
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
const G4Material * GetMaterial() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:189
size_t GetNumberOfElements() const
Definition: G4Material.hh:185
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:205
size_t GetIndex() const
Definition: G4Material.hh:261
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double, G4double)
void AddSecondary(G4Track *aSecondary)
G4double GetWeight() const
const G4DynamicParticle * GetDynamicParticle() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
void SetUseMatrix(G4bool aBool)
G4ParticleDefinition * theDirectPrimaryPartDef
G4Material * currentMaterial
G4MaterialCutsCouple * currentCouple
G4ParticleDefinition * theAdjEquivOfDirectSecondPartDef
void SetApplyCutInRange(G4bool aBool)
G4ParticleDefinition * theAdjEquivOfDirectPrimPartDef
void ProposeTrackStatus(G4TrackStatus status)
void SetSecondaryWeightByProcess(G4bool)
void SetParentWeightByProcess(G4bool)
void ProposeParentWeight(G4double finalWeight)