Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4AdjointPhotoElectricModel Class Reference

#include <G4AdjointPhotoElectricModel.hh>

+ Inheritance diagram for G4AdjointPhotoElectricModel:

Public Member Functions

 G4AdjointPhotoElectricModel ()
 
 ~G4AdjointPhotoElectricModel ()
 
virtual void SampleSecondaries (const G4Track &aTrack, G4bool IsScatProjToProjCase, G4ParticleChange *fParticleChange)
 
virtual G4double AdjointCrossSection (const G4MaterialCutsCouple *aCouple, G4double primEnergy, G4bool IsScatProjToProjCase)
 
virtual G4double GetAdjointCrossSection (const G4MaterialCutsCouple *aCouple, G4double primEnergy, G4bool IsScatProjToProjCase)
 
G4double AdjointCrossSectionPerAtom (const G4Element *anElement, G4double electronEnergy)
 
void SetTheDirectPEEffectModel (G4PEEffectModel *aModel)
 
virtual void CorrectPostStepWeight (G4ParticleChange *fParticleChange, G4double old_weight, G4double adjointPrimKinEnergy, G4double projectileKinEnergy, G4bool IsScatProjToProjCase)
 
- Public Member Functions inherited from G4VEmAdjointModel
 G4VEmAdjointModel (const G4String &nam)
 
virtual ~G4VEmAdjointModel ()
 
virtual void SampleSecondaries (const G4Track &aTrack, G4bool IsScatProjToProjCase, G4ParticleChange *fParticleChange)=0
 
virtual G4double AdjointCrossSection (const G4MaterialCutsCouple *aCouple, G4double primEnergy, G4bool IsScatProjToProjCase)
 
virtual G4double GetAdjointCrossSection (const G4MaterialCutsCouple *aCouple, G4double primEnergy, G4bool IsScatProjToProjCase)
 
virtual G4double DiffCrossSectionPerAtomPrimToSecond (G4double kinEnergyProj, G4double kinEnergyProd, G4double Z, G4double A=0.)
 
virtual G4double DiffCrossSectionPerAtomPrimToScatPrim (G4double kinEnergyProj, G4double kinEnergyScatProj, G4double Z, G4double A=0.)
 
virtual G4double DiffCrossSectionPerVolumePrimToSecond (const G4Material *aMaterial, G4double kinEnergyProj, G4double kinEnergyProd)
 
virtual G4double DiffCrossSectionPerVolumePrimToScatPrim (const G4Material *aMaterial, G4double kinEnergyProj, G4double kinEnergyScatProj)
 
virtual G4double GetSecondAdjEnergyMaxForScatProjToProjCase (G4double PrimAdjEnergy)
 
virtual G4double GetSecondAdjEnergyMinForScatProjToProjCase (G4double PrimAdjEnergy, G4double Tcut=0)
 
virtual G4double GetSecondAdjEnergyMaxForProdToProjCase (G4double PrimAdjEnergy)
 
virtual G4double GetSecondAdjEnergyMinForProdToProjCase (G4double PrimAdjEnergy)
 
void DefineCurrentMaterial (const G4MaterialCutsCouple *couple)
 
std::vector< std::vector< double > * > ComputeAdjointCrossSectionVectorPerAtomForSecond (G4double kinEnergyProd, G4double Z, G4double A=0., G4int nbin_pro_decade=10)
 
std::vector< std::vector< double > * > ComputeAdjointCrossSectionVectorPerAtomForScatProj (G4double kinEnergyProd, G4double Z, G4double A=0., G4int nbin_pro_decade=10)
 
std::vector< std::vector< double > * > ComputeAdjointCrossSectionVectorPerVolumeForSecond (G4Material *aMaterial, G4double kinEnergyProd, G4int nbin_pro_decade=10)
 
std::vector< std::vector< double > * > ComputeAdjointCrossSectionVectorPerVolumeForScatProj (G4Material *aMaterial, G4double kinEnergyProd, G4int nbin_pro_decade=10)
 
void SetCSMatrices (std::vector< G4AdjointCSMatrix * > *Vec1CSMatrix, std::vector< G4AdjointCSMatrix * > *Vec2CSMatrix)
 
G4ParticleDefinitionGetAdjointEquivalentOfDirectPrimaryParticleDefinition ()
 
G4ParticleDefinitionGetAdjointEquivalentOfDirectSecondaryParticleDefinition ()
 
G4double GetHighEnergyLimit ()
 
G4double GetLowEnergyLimit ()
 
void SetHighEnergyLimit (G4double aVal)
 
void SetLowEnergyLimit (G4double aVal)
 
void DefineDirectEMModel (G4VEmModel *aModel)
 
void SetAdjointEquivalentOfDirectPrimaryParticleDefinition (G4ParticleDefinition *aPart)
 
void SetAdjointEquivalentOfDirectSecondaryParticleDefinition (G4ParticleDefinition *aPart)
 
void SetSecondPartOfSameType (G4bool aBool)
 
G4bool GetSecondPartOfSameType ()
 
void SetUseMatrix (G4bool aBool)
 
void SetUseMatrixPerElement (G4bool aBool)
 
void SetUseOnlyOneMatrixForAllElements (G4bool aBool)
 
void SetApplyCutInRange (G4bool aBool)
 
G4bool GetUseMatrix ()
 
G4bool GetUseMatrixPerElement ()
 
G4bool GetUseOnlyOneMatrixForAllElements ()
 
G4bool GetApplyCutInRange ()
 
G4String GetName ()
 
virtual void SetCSBiasingFactor (G4double aVal)
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmAdjointModel
G4double DiffCrossSectionFunction1 (G4double kinEnergyProj)
 
G4double DiffCrossSectionFunction2 (G4double kinEnergyProj)
 
G4double DiffCrossSectionPerVolumeFunctionForIntegrationOverEkinProj (G4double EkinProd)
 
G4double SampleAdjSecEnergyFromCSMatrix (size_t MatrixIndex, G4double prim_energy, G4bool IsScatProjToProjCase)
 
G4double SampleAdjSecEnergyFromCSMatrix (G4double prim_energy, G4bool IsScatProjToProjCase)
 
void SelectCSMatrix (G4bool IsScatProjToProjCase)
 
virtual G4double SampleAdjSecEnergyFromDiffCrossSectionPerAtom (G4double prim_energy, G4bool IsScatProjToProjCase)
 
virtual void CorrectPostStepWeight (G4ParticleChange *fParticleChange, G4double old_weight, G4double adjointPrimKinEnergy, G4double projectileKinEnergy, G4bool IsScatProjToProjCase)
 
- Protected Attributes inherited from G4VEmAdjointModel
G4VEmModeltheDirectEMModel
 
G4VParticleChangepParticleChange
 
const G4String name
 
G4int ASelectedNucleus
 
G4int ZSelectedNucleus
 
G4MaterialSelectedMaterial
 
G4double kinEnergyProdForIntegration
 
G4double kinEnergyScatProjForIntegration
 
G4double kinEnergyProjForIntegration
 
std::vector< G4AdjointCSMatrix * > * pOnCSMatrixForProdToProjBackwardScattering
 
std::vector< G4AdjointCSMatrix * > * pOnCSMatrixForScatProjToProjBackwardScattering
 
std::vector< G4doubleCS_Vs_ElementForScatProjToProjCase
 
std::vector< G4doubleCS_Vs_ElementForProdToProjCase
 
G4double lastCS
 
G4double lastAdjointCSForScatProjToProjCase
 
G4double lastAdjointCSForProdToProjCase
 
G4ParticleDefinitiontheAdjEquivOfDirectPrimPartDef
 
G4ParticleDefinitiontheAdjEquivOfDirectSecondPartDef
 
G4ParticleDefinitiontheDirectPrimaryPartDef
 
G4bool second_part_of_same_type
 
G4double preStepEnergy
 
G4MaterialcurrentMaterial
 
G4MaterialCutsCouplecurrentCouple
 
size_t currentMaterialIndex
 
size_t currentCoupleIndex
 
G4double currentTcutForDirectPrim
 
G4double currentTcutForDirectSecond
 
G4bool ApplyCutInRange
 
G4double mass_ratio_product
 
G4double mass_ratio_projectile
 
G4double HighEnergyLimit
 
G4double LowEnergyLimit
 
G4double CS_biasing_factor
 
G4bool UseMatrix
 
G4bool UseMatrixPerElement
 
G4bool UseOnlyOneMatrixForAllElements
 
size_t indexOfUsedCrossSectionMatrix
 
size_t model_index
 

Detailed Description

Definition at line 62 of file G4AdjointPhotoElectricModel.hh.

Constructor & Destructor Documentation

◆ G4AdjointPhotoElectricModel()

G4AdjointPhotoElectricModel::G4AdjointPhotoElectricModel ( )

Definition at line 42 of file G4AdjointPhotoElectricModel.cc.

42 :
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}
static G4AdjointElectron * AdjointElectron()
static G4AdjointGamma * AdjointGamma()
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
void SetUseMatrix(G4bool aBool)
G4ParticleDefinition * theDirectPrimaryPartDef
G4ParticleDefinition * theAdjEquivOfDirectSecondPartDef
void SetApplyCutInRange(G4bool aBool)
G4ParticleDefinition * theAdjEquivOfDirectPrimPartDef

◆ ~G4AdjointPhotoElectricModel()

G4AdjointPhotoElectricModel::~G4AdjointPhotoElectricModel ( )

Definition at line 66 of file G4AdjointPhotoElectricModel.cc.

67{;}

Member Function Documentation

◆ AdjointCrossSection()

G4double G4AdjointPhotoElectricModel::AdjointCrossSection ( const G4MaterialCutsCouple aCouple,
G4double  primEnergy,
G4bool  IsScatProjToProjCase 
)
virtual

Reimplemented from G4VEmAdjointModel.

Definition at line 188 of file G4AdjointPhotoElectricModel.cc.

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}
std::vector< G4Element * > G4ElementVector
G4double AdjointCrossSectionPerAtom(const G4Element *anElement, G4double electronEnergy)
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:189
size_t GetNumberOfElements() const
Definition: G4Material.hh:185
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:205
G4Material * currentMaterial
G4MaterialCutsCouple * currentCouple

Referenced by GetAdjointCrossSection(), and SampleSecondaries().

◆ AdjointCrossSectionPerAtom()

G4double G4AdjointPhotoElectricModel::AdjointCrossSectionPerAtom ( const G4Element anElement,
G4double  electronEnergy 
)

Definition at line 231 of file G4AdjointPhotoElectricModel.cc.

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}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
G4double GetZ() const
Definition: G4Element.hh:131
G4int GetNbOfAtomicShells() const
Definition: G4Element.hh:146
G4double GetAtomicShell(G4int index) const
Definition: G4Element.cc:367
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double, G4double)

Referenced by AdjointCrossSection().

◆ CorrectPostStepWeight()

void G4AdjointPhotoElectricModel::CorrectPostStepWeight ( G4ParticleChange fParticleChange,
G4double  old_weight,
G4double  adjointPrimKinEnergy,
G4double  projectileKinEnergy,
G4bool  IsScatProjToProjCase 
)
virtual

Reimplemented from G4VEmAdjointModel.

Definition at line 168 of file G4AdjointPhotoElectricModel.cc.

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}
G4double GetPostStepWeightCorrection()
static G4AdjointCSManager * GetAdjointCSManager()
void SetSecondaryWeightByProcess(G4bool)
void SetParentWeightByProcess(G4bool)
void ProposeParentWeight(G4double finalWeight)

Referenced by SampleSecondaries().

◆ GetAdjointCrossSection()

G4double G4AdjointPhotoElectricModel::GetAdjointCrossSection ( const G4MaterialCutsCouple aCouple,
G4double  primEnergy,
G4bool  IsScatProjToProjCase 
)
virtual

Reimplemented from G4VEmAdjointModel.

Definition at line 223 of file G4AdjointPhotoElectricModel.cc.

226{ return AdjointCrossSection(aCouple,electronEnergy,IsScatProjToProjCase);
227}
virtual G4double AdjointCrossSection(const G4MaterialCutsCouple *aCouple, G4double primEnergy, G4bool IsScatProjToProjCase)

◆ SampleSecondaries()

void G4AdjointPhotoElectricModel::SampleSecondaries ( const G4Track aTrack,
G4bool  IsScatProjToProjCase,
G4ParticleChange fParticleChange 
)
virtual

Implements G4VEmAdjointModel.

Definition at line 71 of file G4AdjointPhotoElectricModel.cc.

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}
@ fStopAndKill
#define G4UniformRand()
Definition: Randomize.hh:53
virtual void CorrectPostStepWeight(G4ParticleChange *fParticleChange, G4double old_weight, G4double adjointPrimKinEnergy, G4double projectileKinEnergy, G4bool IsScatProjToProjCase)
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
void AddSecondary(G4Track *aSecondary)
G4double GetWeight() const
const G4DynamicParticle * GetDynamicParticle() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
void ProposeTrackStatus(G4TrackStatus status)

◆ SetTheDirectPEEffectModel()

void G4AdjointPhotoElectricModel::SetTheDirectPEEffectModel ( G4PEEffectModel aModel)
inline

Definition at line 86 of file G4AdjointPhotoElectricModel.hh.

86 {theDirectPEEffectModel = aModel;
87 DefineDirectEMModel(aModel);}
void DefineDirectEMModel(G4VEmModel *aModel)

The documentation for this class was generated from the following files: