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

#include <G4eDPWACoulombScatteringModel.hh>

+ Inheritance diagram for G4eDPWACoulombScatteringModel:

Public Member Functions

 G4eDPWACoulombScatteringModel (G4bool ismixed=false, G4bool isscpcor=true, G4double mumin=0.0)
 
 ~G4eDPWACoulombScatteringModel () override
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *) override
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double ekin, G4double Z, G4double A, G4double prodcut, G4double emax) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double) override
 
void SetTheDCS (G4eDPWAElasticDCS *theDCS)
 
G4eDPWAElasticDCSGetTheDCS ()
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerShell (const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double &length, G4double &eloss)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
virtual void FillNumberOfSecondaries (G4int &numberOfTriplets, G4int &numberOfRecoil)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectTargetAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementGetCurrentElement (const G4Material *mat=nullptr) const
 
G4int SelectRandomAtomNumber (const G4Material *) const
 
const G4IsotopeGetCurrentIsotope (const G4Element *elm=nullptr) const
 
G4int SelectIsotopeNumber (const G4Element *) const
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
G4VEmModelGetTripletModel ()
 
void SetTripletModel (G4VEmModel *)
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool LPMFlag () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy) const
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetFluctuationFlag (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4PhysicsTablexSectionTable = nullptr
 
const G4MaterialpBaseMaterial = nullptr
 
const std::vector< G4double > * theDensityFactor = nullptr
 
const std::vector< G4int > * theDensityIdx = nullptr
 
G4double inveplus
 
G4double pFactor = 1.0
 
size_t currentCoupleIndex = 0
 
size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 

Detailed Description

Definition at line 68 of file G4eDPWACoulombScatteringModel.hh.

Constructor & Destructor Documentation

◆ G4eDPWACoulombScatteringModel()

G4eDPWACoulombScatteringModel::G4eDPWACoulombScatteringModel ( G4bool  ismixed = false,
G4bool  isscpcor = true,
G4double  mumin = 0.0 
)

Constructor.

Parameters
[in]ismixedIndicates if the model is for mixed or for pure single Coulomb scattering. Different type of tables are pre- pared for sampling polar angle of Coulomb scattering for mixed and for pure single scattering models: cosine of the polar scattering angle can be sampled in a restriced inteval (see mumin input parameter below).
[in]isscpcorIndicates if scattering power correction should be used. Note, scattering power correction accounts the effects angular deflections due to sub-threshold ionisations when ionisation is described by using condensed history model (should be active only in this case).
[in]muminWhen the model is used for mixed simulation, Coulomb scatterings, resulting in a minimum t_c polar angular deflection, modelled explicitly. Therefore, cross sections are computed, and angular deflections are sampled ina resricted [\theta_c,\pi] interval. The minimum of this interval is determined by the mumin parameter as: \mu_{min} = \mu(\theta_c)=0.5[1-\cos(\theta_c)]

Definition at line 62 of file G4eDPWACoulombScatteringModel.cc.

63: G4VEmModel("eDPWACoulombScattering"),
64 fIsMixedModel(ismixed),
65 fIsScpCorrection(isscpcor),
66 fMuMin(mumin),
67 fTheDCS(nullptr),
68 fParticleChange(nullptr)
69{
70 SetLowEnergyLimit ( 0.0*CLHEP::eV); // ekin = 10 eV is used if (E< 10 eV)
71 SetHighEnergyLimit(100.0*CLHEP::MeV); // ekin = 100 MeV is used if (E>100 MeV)
72}
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:746
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:753

◆ ~G4eDPWACoulombScatteringModel()

G4eDPWACoulombScatteringModel::~G4eDPWACoulombScatteringModel ( )
override

Definition at line 75 of file G4eDPWACoulombScatteringModel.cc.

76{
77 if (IsMaster()) {
78 delete fTheDCS;
79 }
80}
G4bool IsMaster() const
Definition: G4VEmModel.hh:725

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4eDPWACoulombScatteringModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition ,
G4double  ekin,
G4double  Z,
G4double  A,
G4double  prodcut,
G4double  emax 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 125 of file G4eDPWACoulombScatteringModel.cc.

131{
132 // Cross sections are computed by numerical integration of the pre-computed
133 // DCS data between the muMin, muMax limits where mu(theta)=0.5[1-cos(theta)].
134 // In case of single scattering model (i.e. when fMuMin=0): [muMin=0, muMax=1]
135 // In case of mixed simulation model (i.e. when fMuMin>0): [fMuMin , muMax=1]
136 // NOTE: cross sections will be zero if the kinetic enrgy is out of the
137 // [10 eV-100 MeV] range for which DCS data has been computed.
138 //
139 G4double elCS = 0.0; // elastic cross section
140 G4double tr1CS = 0.0; // first transport cross section
141 G4double tr2CS = 0.0; // second transport cross section
142 const G4double muMin = fMuMin;
143 const G4double muMax = 1.0;
144 fTheDCS->ComputeCSPerAtom((G4int)Z, ekin, elCS, tr1CS, tr2CS, muMin, muMax);
145 // scattering power correction: should be only in condensed history ioni!
146 if (fIsScpCorrection && CurrentCouple()) {
147 const G4double theScpCor = fTheDCS->ComputeScatteringPowerCorrection(CurrentCouple(), ekin);
148 elCS *= (theScpCor*(1.0+1.0/Z));
149 }
150 return std::max(0.0, elCS);
151}
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:486
void ComputeCSPerAtom(G4int iz, G4double ekin, G4double &elcs, G4double &tr1cs, G4double &tr2cs, G4double mumin=0.0, G4double mumax=1.0)
G4double ComputeScatteringPowerCorrection(const G4MaterialCutsCouple *matcut, G4double ekin)

◆ GetTheDCS()

G4eDPWAElasticDCS * G4eDPWACoulombScatteringModel::GetTheDCS ( )
inline

Definition at line 122 of file G4eDPWACoulombScatteringModel.hh.

122{ return fTheDCS; }

Referenced by InitialiseLocal().

◆ Initialise()

void G4eDPWACoulombScatteringModel::Initialise ( const G4ParticleDefinition pdef,
const G4DataVector prodcuts 
)
overridevirtual

Implements G4VEmModel.

Definition at line 83 of file G4eDPWACoulombScatteringModel.cc.

85{
86 if(!fParticleChange) {
87 fParticleChange = GetParticleChangeForGamma();
88 }
89 fMuMin = 0.5*(1.0-std::cos(PolarAngleLimit()));
90 fIsMixedModel = (fMuMin > 0.0);
91 if(IsMaster()) {
92 // clean the G4eDPWAElasticDCS object if any
93 delete fTheDCS;
94 fTheDCS = new G4eDPWAElasticDCS(pdef==G4Electron::Electron(), fIsMixedModel);
95 // init only for the elements that are used in the geometry
97 G4int numOfCouples = (G4int)theCpTable->GetTableSize();
98 for(G4int j=0; j<numOfCouples; ++j) {
99 const G4Material* mat = theCpTable->GetMaterialCutsCouple(j)->GetMaterial();
100 const G4ElementVector* elV = mat->GetElementVector();
101 std::size_t numOfElem = mat->GetNumberOfElements();
102 for (std::size_t ie = 0; ie < numOfElem; ++ie) {
103 fTheDCS->InitialiseForZ((*elV)[ie]->GetZasInt());
104 }
105 }
106 // init scattering power correction
107 if (fIsScpCorrection) {
109 }
110 // will make use of the cross sections so the above needs to be done before
111 InitialiseElementSelectors(pdef, prodcuts);
112 }
113}
std::vector< const G4Element * > G4ElementVector
static G4Electron * Electron()
Definition: G4Electron.cc:93
const G4Material * GetMaterial() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:185
size_t GetNumberOfElements() const
Definition: G4Material.hh:181
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
G4double PolarAngleLimit() const
Definition: G4VEmModel.hh:662
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:124
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:139
void InitialiseForZ(std::size_t iz)
void InitSCPCorrection(G4double lowEnergyLimit, G4double highEnergyLimit)

◆ InitialiseLocal()

void G4eDPWACoulombScatteringModel::InitialiseLocal ( const G4ParticleDefinition ,
G4VEmModel masterModel 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 116 of file G4eDPWACoulombScatteringModel.cc.

118{
120 SetTheDCS(static_cast<G4eDPWACoulombScatteringModel*>(masterModel)->GetTheDCS());
121}
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:831
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:823
void SetTheDCS(G4eDPWAElasticDCS *theDCS)

◆ MinPrimaryEnergy()

G4double G4eDPWACoulombScatteringModel::MinPrimaryEnergy ( const G4Material ,
const G4ParticleDefinition ,
G4double   
)
inlineoverridevirtual

Reimplemented from G4VEmModel.

Definition at line 117 of file G4eDPWACoulombScatteringModel.hh.

118 { return 10.0*CLHEP::eV; }

◆ SampleSecondaries()

void G4eDPWACoulombScatteringModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  ,
const G4MaterialCutsCouple cp,
const G4DynamicParticle dp,
G4double  tmin,
G4double  maxEnergy 
)
overridevirtual

Implements G4VEmModel.

Definition at line 155 of file G4eDPWACoulombScatteringModel.cc.

159{
160 const G4double ekin = dp->GetKineticEnergy();
161 const G4double lekin = dp->GetLogKineticEnergy();
162 const G4Element* target = SelectTargetAtom(cp, dp->GetParticleDefinition(), ekin, lekin);
163 const G4int izet = target->GetZasInt();
164 // sample cosine of the polar scattering angle in (hard) elastic insteraction
165 CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
166 G4double cost = 1.0;
167 if (!fIsMixedModel) {
168 G4double rndm[3];
169 rndmEngine->flatArray(3, rndm);
170 cost = fTheDCS->SampleCosineTheta(izet, lekin, rndm[0], rndm[1], rndm[2]);
171 } else {
172 //sample cost between costMax,costMin where costMax = 1-2xfMuMin;
173 const G4double costMax = 1.0-2.0*fMuMin;
174 const G4double costMin = -1.0;
175 G4double rndm[2];
176 rndmEngine->flatArray(2, rndm);
177 cost = fTheDCS->SampleCosineThetaRestricted(izet, lekin, rndm[0], rndm[1], costMin, costMax);
178 }
179 // compute the new direction in the scattering frame
180 const G4double sint = std::sqrt((1.0-cost)*(1.0+cost));
181 const G4double phi = CLHEP::twopi*rndmEngine->flat();
182 G4ThreeVector theNewDirection(sint*std::cos(phi), sint*std::sin(phi), cost);
183 // get original direction in lab frame and rotate new direction to lab frame
184 G4ThreeVector theOrgDirectionLab = dp->GetMomentumDirection();
185 theNewDirection.rotateUz(theOrgDirectionLab);
186 // set new direction
187 fParticleChange->ProposeMomentumDirection(theNewDirection);
188}
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
virtual double flat()=0
virtual void flatArray(const int size, double *vect)=0
const G4ThreeVector & GetMomentumDirection() const
G4double GetLogKineticEnergy() const
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetKineticEnergy() const
G4int GetZasInt() const
Definition: G4Element.hh:132
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
const G4Element * SelectTargetAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:577
G4double SampleCosineTheta(std::size_t iz, G4double lekin, G4double r1, G4double r2, G4double r3)
G4double SampleCosineThetaRestricted(std::size_t iz, G4double lekin, G4double r1, G4double r2, G4double costMax, G4double costMin)

◆ SetTheDCS()

void G4eDPWACoulombScatteringModel::SetTheDCS ( G4eDPWAElasticDCS theDCS)
inline

Definition at line 120 of file G4eDPWACoulombScatteringModel.hh.

120{ fTheDCS = theDCS; }

Referenced by InitialiseLocal().


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