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

#include <G4BraggIonModel.hh>

Inheritance diagram for G4BraggIonModel:

Public Member Functions

 G4BraggIonModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="BraggIon")
 ~G4BraggIonModel () override
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy) override
void CorrectionsAlongStep (const G4Material *, const G4ParticleDefinition *, const G4double kinEnergy, const G4double cutEnergy, const G4double &length, G4double &eloss) override
G4BraggIonModeloperator= (const G4BraggIonModel &right)=delete
 G4BraggIonModel (const G4BraggIonModel &)=delete
Public Member Functions inherited from G4BraggModel
 G4BraggModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="Bragg")
 ~G4BraggModel () override
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *couple) override
G4double ComputeCrossSectionPerElectron (const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy) override
G4double GetParticleCharge (const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy) override
void SetChargeSquareRatio (G4double val)
G4BraggModeloperator= (const G4BraggModel &right)=delete
 G4BraggModel (const G4BraggModel &)=delete
Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
virtual ~G4VEmModel ()
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
virtual G4double GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
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 void StartTracking (G4Track *)
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
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 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 SetDeexcitationFlag (G4bool val)
void SetForceBuildTable (G4bool val)
void SetFluctuationFlag (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)
void SetLPMFlag (G4bool)
void SetMasterThread (G4bool)
G4VEmModeloperator= (const G4VEmModel &right)=delete
 G4VEmModel (const G4VEmModel &)=delete

Additional Inherited Members

Protected Member Functions inherited from G4BraggModel
void SetParticle (const G4ParticleDefinition *p)
G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kinEnergy) final
Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
G4ParticleChangeForGammaGetParticleChangeForGamma ()
const G4MaterialCutsCoupleCurrentCouple () const
void SetCurrentElement (const G4Element *)
Protected Attributes inherited from G4BraggModel
const G4ParticleDefinitionparticle = nullptr
G4ParticleDefinitiontheElectron = nullptr
G4ParticleChangeForLossfParticleChange = nullptr
const G4MaterialcurrentMaterial = nullptr
const G4MaterialbaseMaterial = nullptr
G4EmCorrectionscorr = nullptr
G4double mass = 0.0
G4double spin = 0.0
G4double chargeSquareRatio = 1.0
G4double massRate = 1.0
G4double ratio = 1.0
G4double protonMassAMU = 1.007276
G4double lowestKinEnergy
G4double theZieglerFactor
G4double expStopPower125
G4int iMolecula = -1
G4int iPSTAR = -1
G4int iICRU90 = -1
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
std::size_t currentCoupleIndex = 0
std::size_t basedCoupleIndex = 0
G4bool lossFlucFlag = true
Static Protected Attributes inherited from G4BraggModel
static G4ICRU90StoppingDatafICRU90 = nullptr
static G4PSTARStoppingfPSTAR = nullptr

Detailed Description

Definition at line 63 of file G4BraggIonModel.hh.

Constructor & Destructor Documentation

◆ G4BraggIonModel() [1/2]

G4BraggIonModel::G4BraggIonModel ( const G4ParticleDefinition * p = nullptr,
const G4String & nam = "BraggIon" )
explicit

Definition at line 83 of file G4BraggIonModel.cc.

85 : G4BraggModel(p, nam)
86{
87 HeMass = 3.727417*CLHEP::GeV;
88 massFactor = 1000.*CLHEP::amu_c2/HeMass;
89}
G4BraggModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="Bragg")

Referenced by G4BraggIonModel(), G4BraggNoDeltaModel::G4BraggNoDeltaModel(), and operator=().

◆ ~G4BraggIonModel()

G4BraggIonModel::~G4BraggIonModel ( )
override

Definition at line 93 of file G4BraggIonModel.cc.

94{
95 if(isFirstAlpha) {
96 delete fASTAR;
97 fASTAR = nullptr;
98 }
99}

◆ G4BraggIonModel() [2/2]

G4BraggIonModel::G4BraggIonModel ( const G4BraggIonModel & )
delete

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4BraggIonModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition * p,
G4double kineticEnergy,
G4double Z,
G4double A,
G4double cutEnergy,
G4double maxEnergy )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 136 of file G4BraggIonModel.cc.

142{
143 G4double sigma =
144 Z*ComputeCrossSectionPerElectron(p,kinEnergy,cutEnergy,maxEnergy);
145 if(isAlpha) {
146 sigma *= (HeEffChargeSquare(Z, kinEnergy/CLHEP::MeV)/chargeSquareRatio);
147 }
148 return sigma;
149}
double G4double
Definition G4Types.hh:83
G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
G4double chargeSquareRatio

◆ ComputeDEDXPerVolume()

G4double G4BraggIonModel::ComputeDEDXPerVolume ( const G4Material * material,
const G4ParticleDefinition * p,
G4double kineticEnergy,
G4double cutEnergy )
overridevirtual

Reimplemented from G4VEmModel.

Reimplemented in G4BraggNoDeltaModel.

Definition at line 172 of file G4BraggIonModel.cc.

176{
177 const G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
178 const G4double tlim = lowestKinEnergy*massRate;
179 const G4double tmin = std::max(std::min(cut, tmax), tlim);
180 G4double dedx = 0.0;
181
182 if(kineticEnergy < tlim) {
183 dedx = HeDEDX(material, tlim)*std::sqrt(kineticEnergy/tlim);
184 } else {
185 dedx = HeDEDX(material, kineticEnergy);
186
187 if (tmin < tmax) {
188 const G4double tau = kineticEnergy/mass;
189 const G4double x = tmin/tmax;
190
191 G4double del =
192 (G4Log(x)*(tau + 1.)*(tau + 1.)/(tau * (tau + 2.0)) + 1.0 - x) *
193 CLHEP::twopi_mc2_rcl2*material->GetElectronDensity();
194 if(isAlpha) {
195 const G4double zeff = material->GetTotNbOfElectPerVolume()/
196 material->GetTotNbOfAtomsPerVolume();
197 heChargeSquare = HeEffChargeSquare(zeff, kineticEnergy/CLHEP::MeV);
198 del *= heChargeSquare;
199 }
200 dedx += del;
201 }
202 }
203 dedx = std::max(dedx, 0.0);
204 /*
205 G4cout << "BraggIon: " << material->GetName()
206 << " E(MeV)=" << kineticEnergy/MeV
207 << " Tmin(MeV)=" << tmin << " dedx(MeV*cm^2/g)="
208 << dedx*gram/(MeV*cm2*material->GetDensity())
209 << " q2=" << chargeSquare << G4endl;
210 */
211 return dedx;
212}
G4double G4Log(G4double x)
Definition G4Log.hh:169
G4double massRate
G4double lowestKinEnergy
G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) final
G4double GetTotNbOfAtomsPerVolume() const
G4double GetTotNbOfElectPerVolume() const
G4double GetElectronDensity() const

Referenced by G4BraggNoDeltaModel::ComputeDEDXPerVolume().

◆ CorrectionsAlongStep()

void G4BraggIonModel::CorrectionsAlongStep ( const G4Material * mat,
const G4ParticleDefinition * p,
const G4double kinEnergy,
const G4double cutEnergy,
const G4double & length,
G4double & eloss )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 216 of file G4BraggIonModel.cc.

222{
223 // no correction for alpha
224 if (isAlpha) { return; }
225
226 // no correction at a small step at the last step
227 if(eloss >= preKinEnergy || eloss < preKinEnergy*0.05) { return; }
228
229 // corrections only for ions
230 if(p != particle) { SetParticle(p); }
231
232 // effective energy and charge at a step
233 const G4double e = std::max(preKinEnergy - eloss*0.5, preKinEnergy*0.5);
234 const G4double q20 = corr->EffectiveChargeSquareRatio(p, mat, preKinEnergy);
235 const G4double q2 = corr->EffectiveChargeSquareRatio(p, mat, e);
236 const G4double qfactor = q2/q20;
237 /*
238 G4cout << "G4BraggIonModel::CorrectionsAlongStep: Epre(MeV)="
239 << preKinEnergy << " Eeff(MeV)=" << e
240 << " eloss=" << eloss << " elossnew=" << eloss*qfactor
241 << " qfactor=" << qfactor << " Qpre=" << q20
242 << p->GetParticleName() <<G4endl;
243 */
244 eloss *= qfactor;
245}
void SetParticle(const G4ParticleDefinition *p)
G4EmCorrections * corr
const G4ParticleDefinition * particle

◆ CrossSectionPerVolume()

G4double G4BraggIonModel::CrossSectionPerVolume ( const G4Material * material,
const G4ParticleDefinition * p,
G4double kineticEnergy,
G4double cutEnergy,
G4double maxEnergy )
overridevirtual

Reimplemented from G4VEmModel.

Reimplemented in G4BraggNoDeltaModel.

Definition at line 153 of file G4BraggIonModel.cc.

159{
160 G4double sigma = material->GetElectronDensity()*
161 ComputeCrossSectionPerElectron(p,kinEnergy,cutEnergy,maxEnergy);
162 if(isAlpha) {
163 const G4double zeff = material->GetTotNbOfElectPerVolume()/
164 material->GetTotNbOfAtomsPerVolume();
165 sigma *= (HeEffChargeSquare(zeff, kinEnergy/CLHEP::MeV)/chargeSquareRatio);
166 }
167 return sigma;
168}

◆ GetChargeSquareRatio()

G4double G4BraggIonModel::GetChargeSquareRatio ( const G4ParticleDefinition * p,
const G4Material * mat,
G4double kineticEnergy )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 125 of file G4BraggIonModel.cc.

128{
129 // this method is called only for ions, so no check if it is an ion
130 if(isAlpha) { return 1.0; }
131 return G4BraggModel::GetChargeSquareRatio(p, mat, kinEnergy);
132}
G4double GetChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy) override

◆ Initialise()

void G4BraggIonModel::Initialise ( const G4ParticleDefinition * p,
const G4DataVector & ref )
overridevirtual

Implements G4VEmModel.

Definition at line 103 of file G4BraggIonModel.cc.

105{
107 const G4String& pname = particle->GetParticleName();
108 if(pname == "alpha") { isAlpha = true; }
109 if(isAlpha && fASTAR == nullptr) {
110 G4AutoLock l(&alphaMutex);
111 if(fASTAR == nullptr) {
112 isFirstAlpha = true;
113 fASTAR = new G4ASTARStopping();
114 }
115 l.unlock();
116 }
117 if(isFirstAlpha) {
118 fASTAR->Initialise();
119 }
120}
G4TemplateAutoLock< G4Mutex > G4AutoLock
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override

◆ operator=()

G4BraggIonModel & G4BraggIonModel::operator= ( const G4BraggIonModel & right)
delete

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