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

#include <G4IonCoulombScatteringModel.hh>

+ Inheritance diagram for G4IonCoulombScatteringModel:

Public Member Functions

 G4IonCoulombScatteringModel (const G4String &nam="IonCoulombScattering")
 
virtual ~G4IonCoulombScatteringModel ()
 
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)
 
void SetHeavyIonCorr (G4int b)
 
G4int GetHeavyIonCorr ()
 
- 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
 
G4IonCoulombCrossSectionioncross
 
const std::vector< G4double > * pCuts
 
const G4MaterialCutsCouplecurrentCouple
 
const G4MaterialcurrentMaterial
 
const G4ElementcurrentElement
 
G4int currentMaterialIndex
 
G4int heavycorr
 
G4double cosThetaMin
 
G4double recoilThreshold
 
const G4ParticleDefinitionparticle
 
const G4ParticleDefinitiontheProton
 
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 71 of file G4IonCoulombScatteringModel.hh.

Constructor & Destructor Documentation

◆ G4IonCoulombScatteringModel()

G4IonCoulombScatteringModel::G4IonCoulombScatteringModel ( const G4String nam = "IonCoulombScattering")

Definition at line 73 of file G4IonCoulombScatteringModel.cc.

74 : G4VEmModel(nam),
75
76 cosThetaMin(1.0),
77 isInitialised(false)
78{
82
83 pCuts=0;
86 currentCouple = 0;
87
88 lowEnergyLimit = 100*eV;
89 recoilThreshold = 0.*eV;
90 heavycorr =0;
91 particle = 0;
92 mass=0;
94
96
97}
const G4ParticleDefinition * particle
const std::vector< G4double > * pCuts
const G4ParticleDefinition * theProton
const G4MaterialCutsCouple * currentCouple
static G4NistManager * Instance()
static G4ParticleTable * GetParticleTable()
static G4Proton * Proton()
Definition: G4Proton.cc:93

◆ ~G4IonCoulombScatteringModel()

G4IonCoulombScatteringModel::~G4IonCoulombScatteringModel ( )
virtual

Definition at line 102 of file G4IonCoulombScatteringModel.cc.

103{ delete ioncross;}

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 127 of file G4IonCoulombScatteringModel.cc.

134{
135
136 SetupParticle(p);
137
138 G4double cross =0.0;
139 if(kinEnergy < lowEnergyLimit) return cross;
140
142
143 G4int iz = G4int(Z);
144
145 //from lab to pCM & mu_rel of effective particle
146 ioncross->SetupKinematic(kinEnergy, cutEnergy,iz);
147
148
149 ioncross->SetupTarget(Z, kinEnergy, heavycorr);
150
151 cross = ioncross->NuclearCrossSection();
152
153//cout<< "..........cross "<<G4BestUnit(cross,"Surface") <<endl;
154 return cross;
155}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
void SetupTarget(G4double Z, G4double kinEnergy, G4int heavycorr)
void SetupKinematic(G4double kinEnergy, G4double cut, G4int iz)
void DefineMaterial(const G4MaterialCutsCouple *)
void SetupParticle(const G4ParticleDefinition *)
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:377

Referenced by SampleSecondaries().

◆ DefineMaterial()

void G4IonCoulombScatteringModel::DefineMaterial ( const G4MaterialCutsCouple cup)
inlineprotected

Definition at line 157 of file G4IonCoulombScatteringModel.hh.

158{
159 if(cup != currentCouple) {
160 currentCouple = cup;
163
164 }
165}
const G4Material * GetMaterial() const

Referenced by ComputeCrossSectionPerAtom(), and SampleSecondaries().

◆ GetHeavyIonCorr()

G4int G4IonCoulombScatteringModel::GetHeavyIonCorr ( )
inline

Definition at line 100 of file G4IonCoulombScatteringModel.hh.

100{return heavycorr; };

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 107 of file G4IonCoulombScatteringModel.cc.

109{
110 SetupParticle(p);
111 currentCouple = 0;
115
117
118
119 if(!isInitialised) {
120 isInitialised = true;
122 }
123}
void Initialise(const G4ParticleDefinition *, G4double cosThetaLim)
G4ParticleChangeForGamma * fParticleChange
const std::vector< G4double > * GetEnergyCutsVector(size_t pcIdx) const
static G4ProductionCutsTable * GetProductionCutsTable()
G4double PolarAngleLimit() const
Definition: G4VEmModel.hh:550
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 159 of file G4IonCoulombScatteringModel.cc.

165{
166 G4double kinEnergy = dp->GetKineticEnergy();
167
168 if(kinEnergy < lowEnergyLimit) return;
169
170 DefineMaterial(couple);
171
173
174 // Choose nucleus
176 kinEnergy,cutEnergy,kinEnergy);
177
179 G4int iz = G4int(Z);
182
183
184
185 G4double cross= ComputeCrossSectionPerAtom(particle,kinEnergy, Z,
186 kinEnergy, cutEnergy, kinEnergy) ;
187 if(cross == 0.0) { return; }
188
189 //scattering angle, z1 == (1-cost)
191 if(z1 > 2.0) { z1 = 2.0; }
192 else if(z1 < 0.0) { z1 = 0.0; }
193
194 G4double cost = 1.0 - z1;
195 G4double sint = sqrt(z1*(1.0 + cost));
196 G4double phi = twopi * G4UniformRand();
197
198
199 // kinematics in the Lab system
200 G4double etot = kinEnergy + mass;
201 G4double mom2= kinEnergy*(kinEnergy+2.0*mass);
202 G4double ptot = sqrt(mom2);
203
204 //CM particle 1
205 G4double bet = ptot/(etot + mass2);
206 G4double gam = 1.0/sqrt((1.0 - bet)*(1.0 + bet));
207
208 //CM
209 G4double momCM2= ioncross->GetMomentum2();
210 G4double momCM =std::sqrt(momCM2);
211 //energy & momentum after scattering of incident particle
212 G4double pxCM = momCM*sint*cos(phi);
213 G4double pyCM = momCM*sint*sin(phi);
214 G4double pzCM = momCM*cost;
215 G4double eCM = sqrt(momCM2 + mass*mass);
216
217 //CM--->Lab
218 G4ThreeVector v1(pxCM , pyCM, gam*(pzCM + bet*eCM));
220
221 G4ThreeVector newDirection = v1.unit();
222 newDirection.rotateUz(dir);
223
225
226 // V.Ivanchenko fix of final energies after scattering
227 // recoil.......................................
228 //G4double trec =(1.0 - cost)* mass2*(etot*etot - mass*mass )/
229 // (mass*mass + mass2*mass2+ 2.*mass2*etot);
230 //G4double finalT = kinEnergy - trec;
231
232 // new computation
233 G4double finalT = gam*(eCM + bet*pzCM) - mass;
234 G4double trec = kinEnergy - finalT;
235
236 if(finalT <= lowEnergyLimit) {
237 trec = kinEnergy;
238 finalT = 0.0;
239 } else if(trec < 0.0) {
240 trec = 0.0;
241 finalT = kinEnergy;
242 }
243
245
247 if(pCuts) { tcut= std::max(tcut,(*pCuts)[currentMaterialIndex]);
248
249 //G4cout<<" tcut eV "<<tcut/eV<<endl;
250 }
251
252 if(trec > tcut) {
253 G4ParticleDefinition* ion = theParticleTable->GetIon(iz, ia, 0.0);
254 G4double plab = sqrt(finalT*(finalT + 2.0*mass));
255 G4ThreeVector p2 = (ptot*dir - plab*newDirection).unit();
256 G4DynamicParticle* newdp = new G4DynamicParticle(ion, p2, trec);
257 fvect->push_back(newdp);
258 } else if(trec > 0.0) {
261 }
262
263
264}
#define G4UniformRand()
Definition: Randomize.hh:53
Hep3Vector unit() const
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
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax)
static G4double GetNuclearMass(const G4double A, const G4double Z)
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)

◆ SetHeavyIonCorr()

void G4IonCoulombScatteringModel::SetHeavyIonCorr ( G4int  b)
inline

Definition at line 99 of file G4IonCoulombScatteringModel.hh.

99{heavycorr=b; };

◆ SetRecoilThreshold()

void G4IonCoulombScatteringModel::SetRecoilThreshold ( G4double  eth)
inline

Definition at line 180 of file G4IonCoulombScatteringModel.hh.

181{
182 recoilThreshold = eth;
183}

◆ SetupParticle()

void G4IonCoulombScatteringModel::SetupParticle ( const G4ParticleDefinition p)
inlineprotected

Definition at line 170 of file G4IonCoulombScatteringModel.hh.

171{
172 if(p != particle) {
173 particle = p;
176 }
177}
void SetupParticle(const G4ParticleDefinition *)

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

Member Data Documentation

◆ cosThetaMin

G4double G4IonCoulombScatteringModel::cosThetaMin
protected

Definition at line 139 of file G4IonCoulombScatteringModel.hh.

Referenced by Initialise().

◆ currentCouple

const G4MaterialCutsCouple* G4IonCoulombScatteringModel::currentCouple
protected

◆ currentElement

const G4Element* G4IonCoulombScatteringModel::currentElement
protected

◆ currentMaterial

const G4Material* G4IonCoulombScatteringModel::currentMaterial
protected

Definition at line 132 of file G4IonCoulombScatteringModel.hh.

Referenced by DefineMaterial(), and G4IonCoulombScatteringModel().

◆ currentMaterialIndex

G4int G4IonCoulombScatteringModel::currentMaterialIndex
protected

◆ fNistManager

G4NistManager* G4IonCoulombScatteringModel::fNistManager
protected

Definition at line 127 of file G4IonCoulombScatteringModel.hh.

Referenced by G4IonCoulombScatteringModel().

◆ fParticleChange

G4ParticleChangeForGamma* G4IonCoulombScatteringModel::fParticleChange
protected

Definition at line 126 of file G4IonCoulombScatteringModel.hh.

Referenced by Initialise(), and SampleSecondaries().

◆ heavycorr

G4int G4IonCoulombScatteringModel::heavycorr
protected

◆ ioncross

◆ lowEnergyLimit

G4double G4IonCoulombScatteringModel::lowEnergyLimit
protected

◆ mass

G4double G4IonCoulombScatteringModel::mass
protected

◆ particle

const G4ParticleDefinition* G4IonCoulombScatteringModel::particle
protected

◆ pCuts

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

◆ recoilThreshold

G4double G4IonCoulombScatteringModel::recoilThreshold
protected

◆ theParticleTable

G4ParticleTable* G4IonCoulombScatteringModel::theParticleTable
protected

◆ theProton

const G4ParticleDefinition* G4IonCoulombScatteringModel::theProton
protected

Definition at line 145 of file G4IonCoulombScatteringModel.hh.

Referenced by G4IonCoulombScatteringModel().


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