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

#include <G4AdjointeIonisationModel.hh>

+ Inheritance diagram for G4AdjointeIonisationModel:

Public Member Functions

 G4AdjointeIonisationModel ()
 
 ~G4AdjointeIonisationModel () override
 
void SampleSecondaries (const G4Track &aTrack, G4bool isScatProjToProj, G4ParticleChange *fParticleChange) override
 
G4double DiffCrossSectionPerAtomPrimToSecond (G4double kinEnergyProj, G4double kinEnergyProd, G4double Z, G4double A=0.) override
 
 G4AdjointeIonisationModel (G4AdjointeIonisationModel &)=delete
 
G4AdjointeIonisationModeloperator= (const G4AdjointeIonisationModel &right)=delete
 
- Public Member Functions inherited from G4VEmAdjointModel
 G4VEmAdjointModel (const G4String &nam)
 
virtual ~G4VEmAdjointModel ()
 
virtual void SampleSecondaries (const G4Track &aTrack, G4bool isScatProjToProj, G4ParticleChange *fParticleChange)=0
 
virtual G4double AdjointCrossSection (const G4MaterialCutsCouple *aCouple, G4double primEnergy, G4bool isScatProjToProj)
 
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 GetSecondAdjEnergyMaxForScatProjToProj (G4double primAdjEnergy)
 
virtual G4double GetSecondAdjEnergyMinForScatProjToProj (G4double primAdjEnergy, G4double tcut=0.)
 
virtual G4double GetSecondAdjEnergyMaxForProdToProj (G4double primAdjEnergy)
 
virtual G4double GetSecondAdjEnergyMinForProdToProj (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)
 
void SetCorrectWeightForPostStepInModel (G4bool aBool)
 
void SetAdditionalWeightCorrectionFactorForPostStepOutsideModel (G4double factor)
 
 G4VEmAdjointModel (G4VEmAdjointModel &)=delete
 
G4VEmAdjointModeloperator= (const G4VEmAdjointModel &right)=delete
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmAdjointModel
G4double DiffCrossSectionFunction1 (G4double kinEnergyProj)
 
G4double DiffCrossSectionFunction2 (G4double kinEnergyProj)
 
G4double SampleAdjSecEnergyFromCSMatrix (std::size_t MatrixIndex, G4double prim_energy, G4bool isScatProjToProj)
 
G4double SampleAdjSecEnergyFromCSMatrix (G4double prim_energy, G4bool isScatProjToProj)
 
void SelectCSMatrix (G4bool isScatProjToProj)
 
virtual G4double SampleAdjSecEnergyFromDiffCrossSectionPerAtom (G4double prim_energy, G4bool isScatProjToProj)
 
virtual void CorrectPostStepWeight (G4ParticleChange *fParticleChange, G4double old_weight, G4double adjointPrimKinEnergy, G4double projectileKinEnergy, G4bool isScatProjToProj)
 
- Protected Attributes inherited from G4VEmAdjointModel
G4AdjointCSManagerfCSManager
 
G4VEmModelfDirectModel = nullptr
 
const G4String fName
 
G4MaterialfSelectedMaterial = nullptr
 
G4MaterialfCurrentMaterial = nullptr
 
G4MaterialCutsCouplefCurrentCouple = nullptr
 
G4ParticleDefinitionfAdjEquivDirectPrimPart = nullptr
 
G4ParticleDefinitionfAdjEquivDirectSecondPart = nullptr
 
G4ParticleDefinitionfDirectPrimaryPart = nullptr
 
std::vector< G4AdjointCSMatrix * > * fCSMatrixProdToProjBackScat = nullptr
 
std::vector< G4AdjointCSMatrix * > * fCSMatrixProjToProjBackScat = nullptr
 
std::vector< G4doublefElementCSScatProjToProj
 
std::vector< G4doublefElementCSProdToProj
 
G4double fKinEnergyProdForIntegration = 0.
 
G4double fKinEnergyScatProjForIntegration = 0.
 
G4double fLastCS = 0.
 
G4double fLastAdjointCSForScatProjToProj = 0.
 
G4double fLastAdjointCSForProdToProj = 0.
 
G4double fPreStepEnergy = 0.
 
G4double fTcutPrim = 0.
 
G4double fTcutSecond = 0.
 
G4double fHighEnergyLimit = 0.
 
G4double fLowEnergyLimit = 0.
 
G4double fCsBiasingFactor = 1.
 
G4double fOutsideWeightFactor = 1.
 
G4int fASelectedNucleus = 0
 
G4int fZSelectedNucleus = 0
 
std::size_t fCSMatrixUsed = 0
 
G4bool fSecondPartSameType = false
 
G4bool fInModelWeightCorr
 
G4bool fApplyCutInRange = true
 
G4bool fUseMatrix = false
 
G4bool fUseMatrixPerElement = false
 
G4bool fOneMatrixForAllElements = false
 

Detailed Description

Definition at line 41 of file G4AdjointeIonisationModel.hh.

Constructor & Destructor Documentation

◆ G4AdjointeIonisationModel() [1/2]

G4AdjointeIonisationModel::G4AdjointeIonisationModel ( )

Definition at line 36 of file G4AdjointeIonisationModel.cc.

37 : G4VEmAdjointModel("Inv_eIon_model")
38
39{
40 fUseMatrix = true;
42 fApplyCutInRange = true;
44
49}
static G4AdjointElectron * AdjointElectron()
static G4Electron * Electron()
Definition: G4Electron.cc:93
G4ParticleDefinition * fAdjEquivDirectSecondPart
G4ParticleDefinition * fAdjEquivDirectPrimPart
G4ParticleDefinition * fDirectPrimaryPart

◆ ~G4AdjointeIonisationModel()

G4AdjointeIonisationModel::~G4AdjointeIonisationModel ( )
override

Definition at line 52 of file G4AdjointeIonisationModel.cc.

52{}

◆ G4AdjointeIonisationModel() [2/2]

G4AdjointeIonisationModel::G4AdjointeIonisationModel ( G4AdjointeIonisationModel )
delete

Member Function Documentation

◆ DiffCrossSectionPerAtomPrimToSecond()

G4double G4AdjointeIonisationModel::DiffCrossSectionPerAtomPrimToSecond ( G4double  kinEnergyProj,
G4double  kinEnergyProd,
G4double  Z,
G4double  A = 0. 
)
overridevirtual

Reimplemented from G4VEmAdjointModel.

Definition at line 164 of file G4AdjointeIonisationModel.cc.

166{
167 G4double dSigmadEprod = 0.;
168 G4double Emax_proj = GetSecondAdjEnergyMaxForProdToProj(kinEnergyProd);
169 G4double Emin_proj = GetSecondAdjEnergyMinForProdToProj(kinEnergyProd);
170
171 // the produced particle should have a kinetic energy smaller than the
172 // projectile
173 if(kinEnergyProj > Emin_proj && kinEnergyProj <= Emax_proj)
174 {
175 dSigmadEprod = Z * DiffCrossSectionMoller(kinEnergyProj, kinEnergyProd);
176 }
177 return dSigmadEprod;
178}
double G4double
Definition: G4Types.hh:83
const G4int Z[17]
virtual G4double GetSecondAdjEnergyMaxForProdToProj(G4double primAdjEnergy)
virtual G4double GetSecondAdjEnergyMinForProdToProj(G4double primAdjEnergy)

◆ operator=()

G4AdjointeIonisationModel & G4AdjointeIonisationModel::operator= ( const G4AdjointeIonisationModel right)
delete

◆ SampleSecondaries()

void G4AdjointeIonisationModel::SampleSecondaries ( const G4Track aTrack,
G4bool  isScatProjToProj,
G4ParticleChange fParticleChange 
)
overridevirtual

Implements G4VEmAdjointModel.

Definition at line 55 of file G4AdjointeIonisationModel.cc.

58{
59 const G4DynamicParticle* theAdjointPrimary = aTrack.GetDynamicParticle();
60
61 // Elastic inverse scattering
62 G4double adjointPrimKinEnergy = theAdjointPrimary->GetKineticEnergy();
63 G4double adjointPrimP = theAdjointPrimary->GetTotalMomentum();
64
65 if(adjointPrimKinEnergy > GetHighEnergyLimit() * 0.999)
66 {
67 return;
68 }
69
70 // Sample secondary energy
71 G4double projectileKinEnergy;
72 if(!fWithRapidSampling)
73 { // used by default
74 projectileKinEnergy =
75 SampleAdjSecEnergyFromCSMatrix(adjointPrimKinEnergy, IsScatProjToProj);
76
77 // Caution!!! this weight correction should be always applied
78 CorrectPostStepWeight(fParticleChange, aTrack.GetWeight(),
79 adjointPrimKinEnergy, projectileKinEnergy,
80 IsScatProjToProj);
81 }
82 else
83 { // only for testing
84 G4double Emin, Emax;
85 if(IsScatProjToProj)
86 {
87 Emin = GetSecondAdjEnergyMinForScatProjToProj(adjointPrimKinEnergy,
89 Emax = GetSecondAdjEnergyMaxForScatProjToProj(adjointPrimKinEnergy);
90 }
91 else
92 {
93 Emin = GetSecondAdjEnergyMinForProdToProj(adjointPrimKinEnergy);
94 Emax = GetSecondAdjEnergyMaxForProdToProj(adjointPrimKinEnergy);
95 }
96 projectileKinEnergy = Emin * std::pow(Emax / Emin, G4UniformRand());
97
99 if(!IsScatProjToProj)
101
102 G4double new_weight = aTrack.GetWeight();
103 G4double used_diffCS =
104 fLastCS * std::log(Emax / Emin) / projectileKinEnergy;
105 G4double needed_diffCS = adjointPrimKinEnergy / projectileKinEnergy;
106 if(!IsScatProjToProj)
108 fCurrentMaterial, projectileKinEnergy, adjointPrimKinEnergy);
109 else
111 fCurrentMaterial, projectileKinEnergy, adjointPrimKinEnergy);
112 new_weight *= needed_diffCS / used_diffCS;
113 fParticleChange->SetParentWeightByProcess(false);
114 fParticleChange->SetSecondaryWeightByProcess(false);
115 fParticleChange->ProposeParentWeight(new_weight);
116 }
117
118 // Kinematic:
119 // we consider a two body elastic scattering for the forward processes where
120 // the projectile knock on an e- at rest and gives it part of its energy
122 G4double projectileTotalEnergy = projectileM0 + projectileKinEnergy;
123 G4double projectileP2 =
124 projectileTotalEnergy * projectileTotalEnergy - projectileM0 * projectileM0;
125
126 // Companion
128 if(IsScatProjToProj)
129 {
130 companionM0 = fAdjEquivDirectSecondPart->GetPDGMass();
131 }
132 G4double companionTotalEnergy =
133 companionM0 + projectileKinEnergy - adjointPrimKinEnergy;
134 G4double companionP2 =
135 companionTotalEnergy * companionTotalEnergy - companionM0 * companionM0;
136
137 // Projectile momentum
138 G4double P_parallel =
139 (adjointPrimP * adjointPrimP + projectileP2 - companionP2) /
140 (2. * adjointPrimP);
141 G4double P_perp = std::sqrt(projectileP2 - P_parallel * P_parallel);
142 G4ThreeVector dir_parallel = theAdjointPrimary->GetMomentumDirection();
143 G4double phi = G4UniformRand() * twopi;
144 G4ThreeVector projectileMomentum =
145 G4ThreeVector(P_perp * std::cos(phi), P_perp * std::sin(phi), P_parallel);
146 projectileMomentum.rotateUz(dir_parallel);
147
148 if(!IsScatProjToProj)
149 { // kill the primary and add a secondary
150 fParticleChange->ProposeTrackStatus(fStopAndKill);
151 fParticleChange->AddSecondary(
152 new G4DynamicParticle(fAdjEquivDirectPrimPart, projectileMomentum));
153 }
154 else
155 {
156 fParticleChange->ProposeEnergy(projectileKinEnergy);
157 fParticleChange->ProposeMomentumDirection(projectileMomentum.unit());
158 }
159}
CLHEP::Hep3Vector G4ThreeVector
@ fStopAndKill
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector unit() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4double GetTotalMomentum() const
void AddSecondary(G4Track *aSecondary)
void ProposeEnergy(G4double finalEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4double GetWeight() const
const G4DynamicParticle * GetDynamicParticle() const
G4double fLastAdjointCSForScatProjToProj
virtual G4double DiffCrossSectionPerVolumePrimToScatPrim(const G4Material *aMaterial, G4double kinEnergyProj, G4double kinEnergyScatProj)
virtual G4double GetSecondAdjEnergyMinForScatProjToProj(G4double primAdjEnergy, G4double tcut=0.)
virtual G4double GetSecondAdjEnergyMaxForScatProjToProj(G4double primAdjEnergy)
G4Material * fCurrentMaterial
virtual G4double DiffCrossSectionPerVolumePrimToSecond(const G4Material *aMaterial, G4double kinEnergyProj, G4double kinEnergyProd)
virtual void CorrectPostStepWeight(G4ParticleChange *fParticleChange, G4double old_weight, G4double adjointPrimKinEnergy, G4double projectileKinEnergy, G4bool isScatProjToProj)
G4double SampleAdjSecEnergyFromCSMatrix(std::size_t MatrixIndex, G4double prim_energy, G4bool isScatProjToProj)
G4double GetHighEnergyLimit()
G4double fLastAdjointCSForProdToProj
void ProposeTrackStatus(G4TrackStatus status)
void SetSecondaryWeightByProcess(G4bool)
void SetParentWeightByProcess(G4bool)
void ProposeParentWeight(G4double finalWeight)

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