Geant4 9.6.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)
 

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 53 of file G4AdjointeIonisationModel.hh.

Constructor & Destructor Documentation

◆ G4AdjointeIonisationModel()

G4AdjointeIonisationModel::G4AdjointeIonisationModel ( )

Definition at line 42 of file G4AdjointeIonisationModel.cc.

42 :
43 G4VEmAdjointModel("Inv_eIon_model")
44
45{
46
47 UseMatrix =true;
49 ApplyCutInRange = true;
52 WithRapidSampling = false;
53
58}
static G4AdjointElectron * AdjointElectron()
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4ParticleDefinition * theDirectPrimaryPartDef
G4bool UseOnlyOneMatrixForAllElements
G4ParticleDefinition * theAdjEquivOfDirectSecondPartDef
G4ParticleDefinition * theAdjEquivOfDirectPrimPartDef

◆ ~G4AdjointeIonisationModel()

G4AdjointeIonisationModel::~G4AdjointeIonisationModel ( )
virtual

Definition at line 61 of file G4AdjointeIonisationModel.cc.

62{;}

Member Function Documentation

◆ DiffCrossSectionPerAtomPrimToSecond()

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

Reimplemented from G4VEmAdjointModel.

Definition at line 176 of file G4AdjointeIonisationModel.cc.

181{
182 G4double dSigmadEprod=0;
183 G4double Emax_proj = GetSecondAdjEnergyMaxForProdToProjCase(kinEnergyProd);
184 G4double Emin_proj = GetSecondAdjEnergyMinForProdToProjCase(kinEnergyProd);
185
186
187 if (kinEnergyProj>Emin_proj && kinEnergyProj<=Emax_proj){ //the produced particle should have a kinetic energy smaller than the projectile
188 dSigmadEprod=Z*DiffCrossSectionMoller(kinEnergyProj,kinEnergyProd);
189 }
190 return dSigmadEprod;
191
192
193
194}
double G4double
Definition: G4Types.hh:64
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 65 of file G4AdjointeIonisationModel.cc.

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