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

#include <G4BraggModel.hh>

+ Inheritance diagram for G4BraggModel:

Public Member Functions

 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
 
G4BraggModeloperator= (const G4BraggModel &right)=delete
 
 G4BraggModel (const G4BraggModel &)=delete
 
- 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
 

Protected Member Functions

G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kinEnergy) final
 
G4double GetChargeSquareRatio () const
 
void SetChargeSquareRatio (G4double val)
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Additional Inherited Members

- 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 71 of file G4BraggModel.hh.

Constructor & Destructor Documentation

◆ G4BraggModel() [1/2]

G4BraggModel::G4BraggModel ( const G4ParticleDefinition p = nullptr,
const G4String nam = "Bragg" 
)
explicit

Definition at line 84 of file G4BraggModel.cc.

85 : G4VEmModel(nam),
86 protonMassAMU(1.007276)
87{
88 SetHighEnergyLimit(2.0*CLHEP::MeV);
89
90 lowestKinEnergy = 0.25*CLHEP::keV;
91 theZieglerFactor = CLHEP::eV*CLHEP::cm2*1.0e-15;
92 theElectron = G4Electron::Electron();
93 expStopPower125 = 0.0;
94
96 if(nullptr != p) { SetParticle(p); }
97 else { SetParticle(theElectron); }
98}
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4LossTableManager * Instance()
G4EmCorrections * EmCorrections()
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:746

◆ ~G4BraggModel()

G4BraggModel::~G4BraggModel ( )
override

Definition at line 102 of file G4BraggModel.cc.

103{
104 if(IsMaster()) {
105 delete fPSTAR;
106 fPSTAR = nullptr;
107 }
108}
G4bool IsMaster() const
Definition: G4VEmModel.hh:725

◆ G4BraggModel() [2/2]

G4BraggModel::G4BraggModel ( const G4BraggModel )
delete

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 207 of file G4BraggModel.cc.

212{
213 return
214 Z*ComputeCrossSectionPerElectron(p,kineticEnergy,cutEnergy,maxEnergy);
215}
const G4int Z[17]
G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)

◆ ComputeCrossSectionPerElectron()

G4double G4BraggModel::ComputeCrossSectionPerElectron ( const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy,
G4double  maxEnergy 
)

Definition at line 177 of file G4BraggModel.cc.

182{
183 G4double cross = 0.0;
184 const G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
185 const G4double maxEnergy = std::min(tmax, maxKinEnergy);
186 const G4double cutEnergy = std::max(cut, lowestKinEnergy*massRate);
187 if(cutEnergy < maxEnergy) {
188
189 const G4double energy = kineticEnergy + mass;
190 const G4double energy2 = energy*energy;
191 const G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
192 cross = (maxEnergy - cutEnergy)/(cutEnergy*maxEnergy)
193 - beta2*G4Log(maxEnergy/cutEnergy)/tmax;
194
195 if( 0.0 < spin ) { cross += 0.5*(maxEnergy - cutEnergy)/energy2; }
196
197 cross *= CLHEP::twopi_mc2_rcl2*chargeSquare/beta2;
198 }
199 // G4cout << "BR: e= " << kineticEnergy << " tmin= " << cutEnergy
200 // << " tmax= " << tmax << " cross= " << cross << G4endl;
201 return cross;
202}
G4double G4Log(G4double x)
Definition: G4Log.hh:227
double G4double
Definition: G4Types.hh:83
G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) final
G4double energy(const ThreeVector &p, const G4double m)

Referenced by ComputeCrossSectionPerAtom(), and CrossSectionPerVolume().

◆ ComputeDEDXPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 231 of file G4BraggModel.cc.

235{
236 const G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
237 const G4double tkin = kineticEnergy/massRate;
238 const G4double cutEnergy = std::max(cut, lowestKinEnergy*massRate);
239 G4double dedx = 0.0;
240
241 if(tkin < lowestKinEnergy) {
242 dedx = DEDX(material, lowestKinEnergy)*std::sqrt(tkin/lowestKinEnergy);
243 } else {
244 dedx = DEDX(material, tkin);
245
246 if (cutEnergy < tmax) {
247 const G4double tau = kineticEnergy/mass;
248 const G4double x = cutEnergy/tmax;
249
250 dedx += (G4Log(x)*(tau + 1.)*(tau + 1.)/(tau * (tau + 2.0)) + 1.0 - x) *
251 CLHEP::twopi_mc2_rcl2 * material->GetElectronDensity();
252 }
253 }
254 dedx = std::max(dedx, 0.0) * chargeSquare;
255
256 //G4cout << "E(MeV)= " << tkin/MeV << " dedx= " << dedx
257 // << " " << material->GetName() << G4endl;
258 return dedx;
259}
G4double GetElectronDensity() const
Definition: G4Material.hh:212

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 219 of file G4BraggModel.cc.

224{
225 return material->GetElectronDensity()
226 *ComputeCrossSectionPerElectron(p,kineticEnergy,cutEnergy,maxEnergy);
227}

◆ GetChargeSquareRatio() [1/2]

G4double G4BraggModel::GetChargeSquareRatio ( ) const
inlineprotected

Definition at line 196 of file G4BraggModel.hh.

197{
198 return chargeSquare;
199}

◆ GetChargeSquareRatio() [2/2]

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

Reimplemented from G4VEmModel.

Definition at line 147 of file G4BraggModel.cc.

150{
151 // this method is called only for ions
152 G4double q2 = corr->EffectiveChargeSquareRatio(p,mat,kineticEnergy);
154 return q2*corr->EffectiveChargeCorrection(p,mat,kineticEnergy);
155}
G4double EffectiveChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double EffectiveChargeCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
virtual void SetParticleAndCharge(const G4ParticleDefinition *, G4double q2)
G4VEmFluctuationModel * GetModelOfFluctuations()
Definition: G4VEmModel.hh:593

◆ GetParticleCharge()

G4double G4BraggModel::GetParticleCharge ( const G4ParticleDefinition p,
const G4Material mat,
G4double  kineticEnergy 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 167 of file G4BraggModel.cc.

170{
171 // this method is called only for ions, so no check if it is an ion
172 return corr->GetParticleCharge(p,mat,kineticEnergy);
173}
G4double GetParticleCharge(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)

◆ Initialise()

void G4BraggModel::Initialise ( const G4ParticleDefinition p,
const G4DataVector  
)
overridevirtual

Implements G4VEmModel.

Definition at line 112 of file G4BraggModel.cc.

114{
115 if(p != particle) { SetParticle(p); }
116
117 // always false before the run
118 SetDeexcitationFlag(false);
119
120 if(IsMaster()) {
121 if(nullptr == fPSTAR) { fPSTAR = new G4PSTARStopping(); }
122 if(particle->GetPDGMass() < CLHEP::GeV) { fPSTAR->Initialise(); }
123 if(G4EmParameters::Instance()->UseICRU90Data()) {
124 if(!fICRU90) {
126 } else if(particle->GetPDGMass() < CLHEP::GeV) { fICRU90->Initialise(); }
127 }
128 }
129
130 if(nullptr == fParticleChange) {
131
134 }
135 G4String pname = particle->GetParticleName();
136 if(particle->GetParticleType() == "nucleus" &&
137 pname != "deuteron" && pname != "triton" &&
138 pname != "alpha+" && pname != "helium" &&
139 pname != "hydrogen") { isIon = true; }
140
141 fParticleChange = GetParticleChangeForLoss();
142 }
143}
static G4EmParameters * Instance()
G4ICRU90StoppingData * GetICRU90StoppingData()
static G4NistManager * Instance()
const G4String & GetParticleType() const
const G4String & GetParticleName() const
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:600
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:802
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:607
G4bool UseAngularGeneratorFlag() const
Definition: G4VEmModel.hh:697
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:109

◆ MaxSecondaryEnergy()

G4double G4BraggModel::MaxSecondaryEnergy ( const G4ParticleDefinition pd,
G4double  kinEnergy 
)
finalprotectedvirtual

Reimplemented from G4VEmModel.

Definition at line 341 of file G4BraggModel.cc.

343{
344 if(pd != particle) { SetParticle(pd); }
345 G4double tau = kinEnergy/mass;
346 G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.) /
347 (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
348 return tmax;
349}

Referenced by ComputeCrossSectionPerElectron(), and ComputeDEDXPerVolume().

◆ MinEnergyCut()

G4double G4BraggModel::MinEnergyCut ( const G4ParticleDefinition ,
const G4MaterialCutsCouple couple 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 159 of file G4BraggModel.cc.

161{
163}
G4double GetMeanExcitationEnergy() const
const G4Material * GetMaterial() const
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:221

◆ operator=()

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

◆ SampleSecondaries()

void G4BraggModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  vdp,
const G4MaterialCutsCouple couple,
const G4DynamicParticle dp,
G4double  tmin,
G4double  maxEnergy 
)
overridevirtual

Implements G4VEmModel.

Definition at line 263 of file G4BraggModel.cc.

268{
269 const G4double tmax = MaxSecondaryKinEnergy(dp);
270 const G4double xmax = std::min(tmax, maxEnergy);
271 const G4double xmin = std::max(lowestKinEnergy*massRate, minEnergy);
272 if(xmin >= xmax) { return; }
273
274 G4double kineticEnergy = dp->GetKineticEnergy();
275 const G4double energy = kineticEnergy + mass;
276 const G4double energy2 = energy*energy;
277 const G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
278 const G4double grej = 1.0;
279 G4double deltaKinEnergy, f;
280
281 CLHEP::HepRandomEngine* rndmEngineMod = G4Random::getTheEngine();
282 G4double rndm[2];
283
284 // sampling follows ...
285 do {
286 rndmEngineMod->flatArray(2, rndm);
287 deltaKinEnergy = xmin*xmax/(xmin*(1.0 - rndm[0]) + xmax*rndm[0]);
288
289 f = 1.0 - beta2*deltaKinEnergy/tmax;
290
291 if(f > grej) {
292 G4cout << "G4BraggModel::SampleSecondary Warning! "
293 << "Majorant " << grej << " < "
294 << f << " for e= " << deltaKinEnergy
295 << G4endl;
296 }
297
298 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
299 } while( grej*rndm[1] >= f );
300
301 G4ThreeVector deltaDirection;
302
304 const G4Material* mat = couple->GetMaterial();
306
307 deltaDirection =
308 GetAngularDistribution()->SampleDirection(dp, deltaKinEnergy, Z, mat);
309
310 } else {
311
312 G4double deltaMomentum =
313 std::sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
314 G4double cost = deltaKinEnergy * (energy + electron_mass_c2) /
315 (deltaMomentum * dp->GetTotalMomentum());
316 if(cost > 1.0) { cost = 1.0; }
317 G4double sint = std::sqrt((1.0 - cost)*(1.0 + cost));
318
319 G4double phi = twopi*rndmEngineMod->flat();
320
321 deltaDirection.set(sint*std::cos(phi),sint*std::sin(phi), cost) ;
322 deltaDirection.rotateUz(dp->GetMomentumDirection());
323 }
324
325 // create G4DynamicParticle object for delta ray
326 auto delta = new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
327
328 // Change kinematics of primary particle
329 kineticEnergy -= deltaKinEnergy;
330 G4ThreeVector finalP = dp->GetMomentum() - delta->GetMomentum();
331 finalP = finalP.unit();
332
333 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
334 fParticleChange->SetProposedMomentumDirection(finalP);
335
336 vdp->push_back(delta);
337}
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
Hep3Vector unit() const
void set(double x, double y, double z)
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 GetKineticEnergy() const
G4ThreeVector GetMomentum() const
G4double GetTotalMomentum() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
G4int SelectRandomAtomNumber(const G4Material *) const
Definition: G4VEmModel.cc:253
G4double MaxSecondaryKinEnergy(const G4DynamicParticle *dynParticle)
Definition: G4VEmModel.hh:501

◆ SetChargeSquareRatio()

void G4BraggModel::SetChargeSquareRatio ( G4double  val)
inlineprotected

Definition at line 201 of file G4BraggModel.hh.

202{
203 chargeSquare = val;
204}

Referenced by G4BraggIonGasModel::ChargeSquareRatio().


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