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

#include <G4eSingleCoulombScatteringModel.hh>

+ Inheritance diagram for G4eSingleCoulombScatteringModel:

Public Member Functions

 G4eSingleCoulombScatteringModel (const G4String &nam="eSingleCoulombScattering")
 
virtual ~G4eSingleCoulombScatteringModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
void SetRecoilThreshold (G4double eth)
 
- 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 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 ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., 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 *, G4double &eloss, G4double &niel, G4double length)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
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)
 
G4int SelectIsotopeNumber (const G4Element *)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *)
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
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
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy)
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void ForceBuildTable (G4bool val)
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 

Protected Member Functions

void DefineMaterial (const G4MaterialCutsCouple *)
 
void SetupParticle (const G4ParticleDefinition *)
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Protected Attributes

G4ParticleTabletheParticleTable
 
G4ParticleChangeForGammafParticleChange
 
G4NistManagerfNistManager
 
G4ScreeningMottCrossSectionMottcross
 
const std::vector< G4double > * pCuts
 
const G4MaterialCutsCouplecurrentCouple
 
const G4MaterialcurrentMaterial
 
const G4ElementcurrentElement
 
G4int currentMaterialIndex
 
G4double cosThetaMin
 
G4double recoilThreshold
 
const G4ParticleDefinitionparticle
 
G4double mass
 
G4double lowEnergyLimit
 
- Protected Attributes inherited from G4VEmModel
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 

Detailed Description

Definition at line 66 of file G4eSingleCoulombScatteringModel.hh.

Constructor & Destructor Documentation

◆ G4eSingleCoulombScatteringModel()

G4eSingleCoulombScatteringModel::G4eSingleCoulombScatteringModel ( const G4String nam = "eSingleCoulombScattering")

Definition at line 70 of file G4eSingleCoulombScatteringModel.cc.

71 : G4VEmModel(nam),
72
73 cosThetaMin(1.0),
74 isInitialised(false)
75{
79
80 pCuts=0;
83 currentCouple = 0;
84
85 lowEnergyLimit = 0*eV;
86 recoilThreshold = 0.*eV;
87 particle = 0;
88 mass=0;
90
92
93}
static G4NistManager * Instance()
static G4ParticleTable * GetParticleTable()

◆ ~G4eSingleCoulombScatteringModel()

G4eSingleCoulombScatteringModel::~G4eSingleCoulombScatteringModel ( )
virtual

Definition at line 98 of file G4eSingleCoulombScatteringModel.cc.

99{ delete Mottcross;}

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4eSingleCoulombScatteringModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition p,
G4double  kinEnergy,
G4double  Z,
G4double  A,
G4double  cut,
G4double  emax 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 123 of file G4eSingleCoulombScatteringModel.cc.

130{
131
132 SetupParticle(p);
133
134 G4double cross =0.0;
135 if(kinEnergy < lowEnergyLimit) return cross;
136
138
139 //Total Cross section
140 Mottcross->SetupKinematic(kinEnergy, Z);
142
143 //cout<< "....cross "<<G4BestUnit(cross,"Surface") << " cm2 "<< cross/cm2 <<endl;
144 return cross;
145}
double G4double
Definition: G4Types.hh:64
void SetupKinematic(G4double kinEnergy, G4double Z)
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:377
void SetupParticle(const G4ParticleDefinition *)
void DefineMaterial(const G4MaterialCutsCouple *)

◆ DefineMaterial()

void G4eSingleCoulombScatteringModel::DefineMaterial ( const G4MaterialCutsCouple cup)
inlineprotected

Definition at line 148 of file G4eSingleCoulombScatteringModel.hh.

149{
150 if(cup != currentCouple) {
151 currentCouple = cup;
154
155 }
156}
const G4Material * GetMaterial() const

Referenced by ComputeCrossSectionPerAtom(), and SampleSecondaries().

◆ Initialise()

void G4eSingleCoulombScatteringModel::Initialise ( const G4ParticleDefinition p,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 103 of file G4eSingleCoulombScatteringModel.cc.

105{
106 SetupParticle(p);
107 currentCouple = 0;
111
113
114
115 if(!isInitialised) {
116 isInitialised = true;
118 }
119}
const std::vector< G4double > * GetEnergyCutsVector(size_t pcIdx) const
static G4ProductionCutsTable * GetProductionCutsTable()
void Initialise(const G4ParticleDefinition *, G4double cosThetaLim)
G4double PolarAngleLimit() const
Definition: G4VEmModel.hh:550
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109

◆ SampleSecondaries()

void G4eSingleCoulombScatteringModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple couple,
const G4DynamicParticle dp,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 149 of file G4eSingleCoulombScatteringModel.cc.

155{
156 G4double kinEnergy = dp->GetKineticEnergy();
157 //cout<<"--- kinEnergy "<<kinEnergy<<endl;
158
159
160 if(kinEnergy < lowEnergyLimit) return;
161
162 DefineMaterial(couple);
164
165 // Choose nucleus
167 kinEnergy,cutEnergy,kinEnergy);//last two :cutEnergy= min e kinEnergy=max
168
170 G4int iz = G4int(Z);
172
173 //cout<<"Element "<<currentElement->GetName()<<endl;;
174
176
177
178 if(cross == 0.0)return;
179
180 G4ThreeVector dir = dp->GetMomentumDirection(); //old direction
181 G4ThreeVector newDirection=Mottcross->GetNewDirection();//new direction
182 newDirection.rotateUz(dir);
183
185
186 //Recoil energy
187 G4double trec= Mottcross->GetTrec();
188 //Energy after scattering
189 G4double finalT = kinEnergy - trec;
190
191
192 if(finalT <= lowEnergyLimit) {
193 trec = kinEnergy;
194 finalT = 0.0;
195 }
196
198
200 if(pCuts) { tcut= std::min(tcut,(*pCuts)[currentMaterialIndex]);
201 }
202
203
204 if(trec > tcut) {
205
206 //cout<<"Trec "<<trec/eV<<endl;
207 G4ParticleDefinition* ion = theParticleTable->GetIon(iz, ia, 0.0);
208
209 //incident before scattering
210 G4double ptot=sqrt(Mottcross->GetMom2Lab());
211 //incident after scattering
212 G4double plab = sqrt(finalT*(finalT + 2.0*mass));
213 G4ThreeVector p2 = (ptot*dir - plab*newDirection).unit();
214 //secondary particle
215 G4DynamicParticle* newdp = new G4DynamicParticle(ion, p2, trec);
216 fvect->push_back(newdp);
217 }
218
219 else if(trec > 0.0) {
221 if(trec< tcut) fParticleChange->ProposeLocalEnergyDeposit(trec);
222 }
223
224
225}
int G4int
Definition: G4Types.hh:66
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetZ() const
Definition: G4Element.hh:131
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4ParticleDefinition * GetIon(G4int atomicNumber, G4int atomicMass, G4double excitationEnergy)
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:459
G4int SelectIsotopeNumber(const G4Element *)
Definition: G4VEmModel.hh:478
void ProposeNonIonizingEnergyDeposit(G4double anEnergyPart)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

◆ SetRecoilThreshold()

void G4eSingleCoulombScatteringModel::SetRecoilThreshold ( G4double  eth)
inline

Definition at line 171 of file G4eSingleCoulombScatteringModel.hh.

172{
173 recoilThreshold = eth;
174}

◆ SetupParticle()

void G4eSingleCoulombScatteringModel::SetupParticle ( const G4ParticleDefinition p)
inlineprotected

Definition at line 161 of file G4eSingleCoulombScatteringModel.hh.

162{
163 if(p != particle) {
164 particle = p;
167 }
168}
void SetupParticle(const G4ParticleDefinition *)

Referenced by ComputeCrossSectionPerAtom(), Initialise(), and SampleSecondaries().

Member Data Documentation

◆ cosThetaMin

G4double G4eSingleCoulombScatteringModel::cosThetaMin
protected

Definition at line 131 of file G4eSingleCoulombScatteringModel.hh.

Referenced by Initialise().

◆ currentCouple

const G4MaterialCutsCouple* G4eSingleCoulombScatteringModel::currentCouple
protected

◆ currentElement

const G4Element* G4eSingleCoulombScatteringModel::currentElement
protected

◆ currentMaterial

const G4Material* G4eSingleCoulombScatteringModel::currentMaterial
protected

◆ currentMaterialIndex

G4int G4eSingleCoulombScatteringModel::currentMaterialIndex
protected

◆ fNistManager

G4NistManager* G4eSingleCoulombScatteringModel::fNistManager
protected

Definition at line 121 of file G4eSingleCoulombScatteringModel.hh.

Referenced by G4eSingleCoulombScatteringModel().

◆ fParticleChange

G4ParticleChangeForGamma* G4eSingleCoulombScatteringModel::fParticleChange
protected

◆ lowEnergyLimit

G4double G4eSingleCoulombScatteringModel::lowEnergyLimit
protected

◆ mass

G4double G4eSingleCoulombScatteringModel::mass
protected

◆ Mottcross

◆ particle

const G4ParticleDefinition* G4eSingleCoulombScatteringModel::particle
protected

◆ pCuts

const std::vector<G4double>* G4eSingleCoulombScatteringModel::pCuts
protected

◆ recoilThreshold

G4double G4eSingleCoulombScatteringModel::recoilThreshold
protected

◆ theParticleTable

G4ParticleTable* G4eSingleCoulombScatteringModel::theParticleTable
protected

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