Geant4 10.7.0
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")
 
virtual ~G4BraggModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
virtual G4double ComputeCrossSectionPerElectron (const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy) override
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy) override
 
- 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 *, G4double &eloss, G4double &niel, G4double length)
 
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 ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
virtual 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)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
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 *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Member Functions

virtual 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
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const G4MaterialpBaseMaterial
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 
G4double inveplus
 
G4double pFactor
 

Detailed Description

Definition at line 71 of file G4BraggModel.hh.

Constructor & Destructor Documentation

◆ G4BraggModel()

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

Definition at line 85 of file G4BraggModel.cc.

86 : G4VEmModel(nam),
87 particle(nullptr),
88 fICRU90(nullptr),
89 currentMaterial(nullptr),
90 baseMaterial(nullptr),
91 protonMassAMU(1.007276),
92 iMolecula(-1),
93 iPSTAR(-1),
94 iICRU90(-1),
95 isIon(false)
96{
97 fParticleChange = nullptr;
98 SetHighEnergyLimit(2.0*MeV);
99
100 lowestKinEnergy = 1.0*keV;
101 theZieglerFactor = eV*cm2*1.0e-15;
102 theElectron = G4Electron::Electron();
103 expStopPower125 = 0.0;
104
106 if(p) { SetParticle(p); }
107 else { SetParticle(theElectron); }
108}
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4LossTableManager * Instance()
G4EmCorrections * EmCorrections()
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:757

◆ ~G4BraggModel()

G4BraggModel::~G4BraggModel ( )
virtual

Definition at line 112 of file G4BraggModel.cc.

113{
114 if(IsMaster()) { delete fPSTAR; fPSTAR = nullptr; }
115}
G4bool IsMaster() const
Definition: G4VEmModel.hh:736

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 206 of file G4BraggModel.cc.

211{
212 return
213 Z*ComputeCrossSectionPerElectron(p,kineticEnergy,cutEnergy,maxEnergy);
214}
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)

◆ ComputeCrossSectionPerElectron()

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

Definition at line 176 of file G4BraggModel.cc.

181{
182 G4double cross = 0.0;
183 G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
184 G4double maxEnergy = std::min(tmax,maxKinEnergy);
185 if(cutEnergy < maxEnergy) {
186
187 G4double energy = kineticEnergy + mass;
188 G4double energy2 = energy*energy;
189 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
190 cross = (maxEnergy - cutEnergy)/(cutEnergy*maxEnergy)
191 - beta2*G4Log(maxEnergy/cutEnergy)/tmax;
192
193 if( 0.0 < spin ) { cross += 0.5*(maxEnergy - cutEnergy)/energy2; }
194
195 cross *= twopi_mc2_rcl2*chargeSquare/beta2;
196 }
197 // G4cout << "BR: e= " << kineticEnergy << " tmin= " << cutEnergy
198 // << " tmax= " << tmax << " cross= " << cross << G4endl;
199
200 return cross;
201}
G4double G4Log(G4double x)
Definition: G4Log.hh:226
double G4double
Definition: G4Types.hh:83
virtual 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 230 of file G4BraggModel.cc.

234{
235 G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
236 G4double tkin = kineticEnergy/massRate;
237 G4double dedx = 0.0;
238
239 if(tkin < lowestKinEnergy) {
240 dedx = DEDX(material, lowestKinEnergy)*sqrt(tkin/lowestKinEnergy);
241 } else {
242 dedx = DEDX(material, tkin);
243 }
244
245 if (cutEnergy < tmax) {
246
247 G4double tau = kineticEnergy/mass;
248 G4double gam = tau + 1.0;
249 G4double bg2 = tau * (tau+2.0);
250 G4double beta2 = bg2/(gam*gam);
251 G4double x = cutEnergy/tmax;
252
253 dedx += (G4Log(x) + (1.0 - x)*beta2) * twopi_mc2_rcl2
254 * (material->GetElectronDensity())/beta2;
255 }
256
257 // now compute the total ionization loss
258
259 dedx = std::max(dedx, 0.0) * chargeSquare;
260
261 //G4cout << "E(MeV)= " << tkin/MeV << " dedx= " << dedx
262 // << " " << material->GetName() << G4endl;
263
264 return dedx;
265}
G4double GetElectronDensity() const
Definition: G4Material.hh:215

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 218 of file G4BraggModel.cc.

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

◆ GetChargeSquareRatio() [1/2]

G4double G4BraggModel::GetChargeSquareRatio ( ) const
inlineprotected

Definition at line 194 of file G4BraggModel.hh.

195{
196 return chargeSquare;
197}

◆ GetChargeSquareRatio() [2/2]

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

Reimplemented from G4VEmModel.

Definition at line 154 of file G4BraggModel.cc.

157{
158 // this method is called only for ions
159 G4double q2 = corr->EffectiveChargeSquareRatio(p,mat,kineticEnergy);
161 return q2*corr->EffectiveChargeCorrection(p,mat,kineticEnergy);
162}
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:604

◆ GetParticleCharge()

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

Reimplemented from G4VEmModel.

Reimplemented in G4BraggIonGasModel.

Definition at line 166 of file G4BraggModel.cc.

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

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 119 of file G4BraggModel.cc.

121{
122 if(p != particle) { SetParticle(p); }
123
124 // always false before the run
125 SetDeexcitationFlag(false);
126
127 if(IsMaster()) {
128 if(nullptr == fPSTAR) { fPSTAR = new G4PSTARStopping(); }
129 if(particle->GetPDGMass() < GeV) { fPSTAR->Initialise(); }
130 if(G4EmParameters::Instance()->UseICRU90Data()) {
131 if(!fICRU90) {
133 } else if(particle->GetPDGMass() < GeV) { fICRU90->Initialise(); }
134 }
135 }
136
137 if(nullptr == fParticleChange) {
138
141 }
142 G4String pname = particle->GetParticleName();
143 if(particle->GetParticleType() == "nucleus" &&
144 pname != "deuteron" && pname != "triton" &&
145 pname != "alpha+" && pname != "helium" &&
146 pname != "hydrogen") { isIon = true; }
147
148 fParticleChange = GetParticleChangeForLoss();
149 }
150}
static G4EmParameters * Instance()
G4ICRU90StoppingData * GetICRU90StoppingData()
static G4NistManager * Instance()
const G4String & GetParticleType() const
const G4String & GetParticleName() const
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:611
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:813
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:618
G4bool UseAngularGeneratorFlag() const
Definition: G4VEmModel.hh:708
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:118

◆ MaxSecondaryEnergy()

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

Reimplemented from G4VEmModel.

Definition at line 347 of file G4BraggModel.cc.

349{
350 if(pd != particle) { SetParticle(pd); }
351 G4double tau = kinEnergy/mass;
352 G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.) /
353 (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
354 return tmax;
355}

Referenced by ComputeCrossSectionPerElectron(), and ComputeDEDXPerVolume().

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 269 of file G4BraggModel.cc.

274{
276 G4double xmax = std::min(tmax, maxEnergy);
277 if(xmin >= xmax) { return; }
278
279 G4double kineticEnergy = dp->GetKineticEnergy();
280 G4double energy = kineticEnergy + mass;
281 G4double energy2 = energy*energy;
282 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
283 G4double grej = 1.0;
284 G4double deltaKinEnergy, f;
285
286 CLHEP::HepRandomEngine* rndmEngineMod = G4Random::getTheEngine();
287 G4double rndm[2];
288
289 // sampling follows ...
290 do {
291 rndmEngineMod->flatArray(2, rndm);
292 deltaKinEnergy = xmin*xmax/(xmin*(1.0 - rndm[0]) + xmax*rndm[0]);
293
294 f = 1.0 - beta2*deltaKinEnergy/tmax;
295
296 if(f > grej) {
297 G4cout << "G4BraggModel::SampleSecondary Warning! "
298 << "Majorant " << grej << " < "
299 << f << " for e= " << deltaKinEnergy
300 << G4endl;
301 }
302
303 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
304 } while( grej*rndm[1] >= f );
305
306 G4ThreeVector deltaDirection;
307
309 const G4Material* mat = couple->GetMaterial();
311
312 deltaDirection =
313 GetAngularDistribution()->SampleDirection(dp, deltaKinEnergy, Z, mat);
314
315 } else {
316
317 G4double deltaMomentum =
318 sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
319 G4double cost = deltaKinEnergy * (energy + electron_mass_c2) /
320 (deltaMomentum * dp->GetTotalMomentum());
321 if(cost > 1.0) { cost = 1.0; }
322 G4double sint = sqrt((1.0 - cost)*(1.0 + cost));
323
324 G4double phi = twopi*rndmEngineMod->flat();
325
326 deltaDirection.set(sint*cos(phi),sint*sin(phi), cost) ;
327 deltaDirection.rotateUz(dp->GetMomentumDirection());
328 }
329
330 // create G4DynamicParticle object for delta ray
331 G4DynamicParticle* delta =
332 new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
333
334 // Change kinematics of primary particle
335 kineticEnergy -= deltaKinEnergy;
336 G4ThreeVector finalP = dp->GetMomentum() - delta->GetMomentum();
337 finalP = finalP.unit();
338
339 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
340 fParticleChange->SetProposedMomentumDirection(finalP);
341
342 vdp->push_back(delta);
343}
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
const G4Material * GetMaterial() 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 *)
Definition: G4VEmModel.cc:315
G4double MaxSecondaryKinEnergy(const G4DynamicParticle *dynParticle)
Definition: G4VEmModel.hh:510

◆ SetChargeSquareRatio()

void G4BraggModel::SetChargeSquareRatio ( G4double  val)
inlineprotected

Definition at line 199 of file G4BraggModel.hh.

200{
201 chargeSquare = val;
202}

Referenced by G4BraggIonGasModel::ChargeSquareRatio().


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