Geant4 10.7.0
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 ()
 
virtual ~G4AdjointeIonisationModel ()
 
virtual void SampleSecondaries (const G4Track &aTrack, G4bool IsScatProjToProjCase, G4ParticleChange *fParticleChange)
 
virtual G4double DiffCrossSectionPerAtomPrimToSecond (G4double kinEnergyProj, G4double kinEnergyProd, G4double Z, G4double A=0.)
 
- 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)
 
void SetCorrectWeightForPostStepInModel (G4bool aBool)
 
void SetAdditionalWeightCorrectionFactorForPostStepOutsideModel (G4double factor)
 

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
 
G4bool correct_weight_for_post_step_in_model
 
G4double additional_weight_correction_factor_for_post_step_outside_model
 

Detailed Description

Definition at line 51 of file G4AdjointeIonisationModel.hh.

Constructor & Destructor Documentation

◆ G4AdjointeIonisationModel()

G4AdjointeIonisationModel::G4AdjointeIonisationModel ( )

Definition at line 41 of file G4AdjointeIonisationModel.cc.

41 :
42 G4VEmAdjointModel("Inv_eIon_model")
43
44{
45
46 UseMatrix =true;
48 ApplyCutInRange = true;
51 WithRapidSampling = false;
52
57}
static G4AdjointElectron * AdjointElectron()
static G4Electron * Electron()
Definition: G4Electron.cc:93
G4ParticleDefinition * theDirectPrimaryPartDef
G4bool UseOnlyOneMatrixForAllElements
G4ParticleDefinition * theAdjEquivOfDirectSecondPartDef
G4ParticleDefinition * theAdjEquivOfDirectPrimPartDef

◆ ~G4AdjointeIonisationModel()

G4AdjointeIonisationModel::~G4AdjointeIonisationModel ( )
virtual

Definition at line 60 of file G4AdjointeIonisationModel.cc.

61{;}

Member Function Documentation

◆ DiffCrossSectionPerAtomPrimToSecond()

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

Reimplemented from G4VEmAdjointModel.

Definition at line 175 of file G4AdjointeIonisationModel.cc.

180{
181 G4double dSigmadEprod=0;
182 G4double Emax_proj = GetSecondAdjEnergyMaxForProdToProjCase(kinEnergyProd);
183 G4double Emin_proj = GetSecondAdjEnergyMinForProdToProjCase(kinEnergyProd);
184
185
186 if (kinEnergyProj>Emin_proj && kinEnergyProj<=Emax_proj){ //the produced particle should have a kinetic energy smaller than the projectile
187 dSigmadEprod=Z*DiffCrossSectionMoller(kinEnergyProj,kinEnergyProd);
188 }
189 return dSigmadEprod;
190
191
192
193}
double G4double
Definition: G4Types.hh:83
virtual G4double GetSecondAdjEnergyMaxForProdToProjCase(G4double PrimAdjEnergy)
virtual G4double GetSecondAdjEnergyMinForProdToProjCase(G4double PrimAdjEnergy)

◆ SampleSecondaries()

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

Implements G4VEmAdjointModel.

Definition at line 64 of file G4AdjointeIonisationModel.cc.

67{
68
69
70 const G4DynamicParticle* theAdjointPrimary =aTrack.GetDynamicParticle();
71
72 //Elastic inverse scattering
73 //---------------------------------------------------------
74 G4double adjointPrimKinEnergy = theAdjointPrimary->GetKineticEnergy();
75 G4double adjointPrimP =theAdjointPrimary->GetTotalMomentum();
76
77 if (adjointPrimKinEnergy>HighEnergyLimit*0.999){
78 return;
79 }
80
81 //Sample secondary energy
82 //-----------------------
83 G4double projectileKinEnergy;
84 if (!WithRapidSampling ) { //used by default
85 projectileKinEnergy = SampleAdjSecEnergyFromCSMatrix(adjointPrimKinEnergy, IsScatProjToProjCase);
86
87 CorrectPostStepWeight(fParticleChange,
88 aTrack.GetWeight(),
89 adjointPrimKinEnergy,
90 projectileKinEnergy,
91 IsScatProjToProjCase); //Caution !!!this weight correction should be always applied
92 }
93 else { //only for test at the moment
94
95 G4double Emin,Emax;
96 if (IsScatProjToProjCase) {
98 Emax=GetSecondAdjEnergyMaxForScatProjToProjCase(adjointPrimKinEnergy);
99 }
100 else {
101 Emin=GetSecondAdjEnergyMinForProdToProjCase(adjointPrimKinEnergy);
102 Emax=GetSecondAdjEnergyMaxForProdToProjCase(adjointPrimKinEnergy);
103 }
104 projectileKinEnergy = Emin*std::pow(Emax/Emin,G4UniformRand());
105
106
107
109 if ( !IsScatProjToProjCase) lastCS=lastAdjointCSForProdToProjCase;
110
111 G4double new_weight=aTrack.GetWeight();
112 G4double used_diffCS=lastCS*std::log(Emax/Emin)/projectileKinEnergy;
113 G4double needed_diffCS=adjointPrimKinEnergy/projectileKinEnergy;
114 if (!IsScatProjToProjCase) needed_diffCS *=DiffCrossSectionPerVolumePrimToSecond(currentMaterial,projectileKinEnergy,adjointPrimKinEnergy);
115 else needed_diffCS *=DiffCrossSectionPerVolumePrimToScatPrim(currentMaterial,projectileKinEnergy,adjointPrimKinEnergy);
116 new_weight*=needed_diffCS/used_diffCS;
117 fParticleChange->SetParentWeightByProcess(false);
118 fParticleChange->SetSecondaryWeightByProcess(false);
119 fParticleChange->ProposeParentWeight(new_weight);
120
121
122 }
123
124
125
126 //Kinematic:
127 //we consider a two body elastic scattering for the forward processes where the projectile knock on an e- at rest and gives
128 // him part of its energy
129 //----------------------------------------------------------------------------------------
130
132 G4double projectileTotalEnergy = projectileM0+projectileKinEnergy;
133 G4double projectileP2 = projectileTotalEnergy*projectileTotalEnergy - projectileM0*projectileM0;
134
135
136
137 //Companion
138 //-----------
140 if (IsScatProjToProjCase) {
142 }
143 G4double companionTotalEnergy =companionM0+ projectileKinEnergy-adjointPrimKinEnergy;
144 G4double companionP2 = companionTotalEnergy*companionTotalEnergy - companionM0*companionM0;
145
146
147 //Projectile momentum
148 //--------------------
149 G4double P_parallel = (adjointPrimP*adjointPrimP + projectileP2 - companionP2)/(2.*adjointPrimP);
150 G4double P_perp = std::sqrt( projectileP2 - P_parallel*P_parallel);
151 G4ThreeVector dir_parallel=theAdjointPrimary->GetMomentumDirection();
152 G4double phi =G4UniformRand()*2.*3.1415926;
153 G4ThreeVector projectileMomentum = G4ThreeVector(P_perp*std::cos(phi),P_perp*std::sin(phi),P_parallel);
154 projectileMomentum.rotateUz(dir_parallel);
155
156
157
158 if (!IsScatProjToProjCase ){ //kill the primary and add a secondary
159 fParticleChange->ProposeTrackStatus(fStopAndKill);
160 fParticleChange->AddSecondary(new G4DynamicParticle(theAdjEquivOfDirectPrimPartDef,projectileMomentum));
161 //G4cout<<"projectileMomentum "<<projectileMomentum<<G4endl;
162 }
163 else {
164 fParticleChange->ProposeEnergy(projectileKinEnergy);
165 fParticleChange->ProposeMomentumDirection(projectileMomentum.unit());
166 }
167
168
169
170
171}
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
virtual G4double DiffCrossSectionPerVolumePrimToScatPrim(const G4Material *aMaterial, G4double kinEnergyProj, G4double kinEnergyScatProj)
G4double lastAdjointCSForProdToProjCase
G4double lastAdjointCSForScatProjToProjCase
virtual void CorrectPostStepWeight(G4ParticleChange *fParticleChange, G4double old_weight, G4double adjointPrimKinEnergy, G4double projectileKinEnergy, G4bool IsScatProjToProjCase)
virtual G4double DiffCrossSectionPerVolumePrimToSecond(const G4Material *aMaterial, G4double kinEnergyProj, G4double kinEnergyProd)
G4double SampleAdjSecEnergyFromCSMatrix(size_t MatrixIndex, G4double prim_energy, G4bool IsScatProjToProjCase)
virtual G4double GetSecondAdjEnergyMaxForScatProjToProjCase(G4double PrimAdjEnergy)
G4Material * currentMaterial
G4double currentTcutForDirectSecond
virtual G4double GetSecondAdjEnergyMinForScatProjToProjCase(G4double PrimAdjEnergy, G4double Tcut=0)
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: