Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4BraggModel Class Reference

#include <G4BraggModel.hh>

Inheritance diagram for G4BraggModel:
Collaboration 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 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 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)
 
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)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
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
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
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 SetForceBuildTable (G4bool val)
 
void SetFluctuationFlag (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () 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)
 

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 ()
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Additional Inherited Members

- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 71 of file G4BraggModel.hh.

Constructor & Destructor Documentation

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

Definition at line 84 of file G4BraggModel.cc.

85  : G4VEmModel(nam),
86  particle(0),
87  currentMaterial(0),
88  protonMassAMU(1.007276),
89  iMolecula(-1),
90  iPSTAR(-1),
91  isIon(false)
92 {
93  fParticleChange = nullptr;
95 
96  lowestKinEnergy = 1.0*keV;
97  theZieglerFactor = eV*cm2*1.0e-15;
98  theElectron = G4Electron::Electron();
99  expStopPower125 = 0.0;
100 
102  if(p) { SetParticle(p); }
103  else { SetParticle(theElectron); }
104 }
static G4LossTableManager * Instance()
static constexpr double cm2
Definition: G4SIunits.hh:120
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
G4EmCorrections * EmCorrections()
static constexpr double eV
Definition: G4SIunits.hh:215
static G4Electron * Electron()
Definition: G4Electron.cc:94
static constexpr double MeV
Definition: G4SIunits.hh:214
static constexpr double keV
Definition: G4SIunits.hh:216

Here is the call graph for this function:

G4BraggModel::~G4BraggModel ( )
virtual

Definition at line 108 of file G4BraggModel.cc.

109 {
110  if(IsMaster()) { delete fPSTAR; fPSTAR = nullptr; }
111 }
G4bool IsMaster() const
Definition: G4VEmModel.hh:725

Here is the call graph for this function:

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 196 of file G4BraggModel.cc.

202 {
204  (p,kineticEnergy,cutEnergy,maxEnergy);
205  return cross;
206 }
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

Definition at line 167 of file G4BraggModel.cc.

172 {
173  G4double cross = 0.0;
174  G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
175  G4double maxEnergy = std::min(tmax,maxKinEnergy);
176  if(cutEnergy < maxEnergy) {
177 
178  G4double energy = kineticEnergy + mass;
179  G4double energy2 = energy*energy;
180  G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
181  cross = (maxEnergy - cutEnergy)/(cutEnergy*maxEnergy)
182  - beta2*G4Log(maxEnergy/cutEnergy)/tmax;
183 
184  if( 0.0 < spin ) { cross += 0.5*(maxEnergy - cutEnergy)/energy2; }
185 
186  cross *= twopi_mc2_rcl2*chargeSquare/beta2;
187  }
188  // G4cout << "BR: e= " << kineticEnergy << " tmin= " << cutEnergy
189  // << " tmax= " << tmax << " cross= " << cross << G4endl;
190 
191  return cross;
192 }
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) final
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double energy(const ThreeVector &p, const G4double m)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 225 of file G4BraggModel.cc.

229 {
230  G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
231  G4double tkin = kineticEnergy/massRate;
232  G4double dedx = 0.0;
233 
234  if(tkin < lowestKinEnergy) {
235  dedx = DEDX(material, lowestKinEnergy)*sqrt(tkin/lowestKinEnergy);
236  } else {
237  dedx = DEDX(material, tkin);
238  }
239 
240  if (cutEnergy < tmax) {
241 
242  G4double tau = kineticEnergy/mass;
243  G4double gam = tau + 1.0;
244  G4double bg2 = tau * (tau+2.0);
245  G4double beta2 = bg2/(gam*gam);
246  G4double x = cutEnergy/tmax;
247 
248  dedx += (G4Log(x) + (1.0 - x)*beta2) * twopi_mc2_rcl2
249  * (material->GetElectronDensity())/beta2;
250  }
251 
252  // now compute the total ionization loss
253 
254  dedx = std::max(dedx, 0.0) * chargeSquare;
255 
256  //G4cout << "E(MeV)= " << tkin/MeV << " dedx= " << dedx
257  // << " " << material->GetName() << G4endl;
258 
259  return dedx;
260 }
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) final
tuple x
Definition: test.py:50
G4double GetElectronDensity() const
Definition: G4Material.hh:217
G4double G4Log(G4double x)
Definition: G4Log.hh:230
T max(const T t1, const T t2)
brief Return the largest of the two arguments
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 210 of file G4BraggModel.cc.

216 {
217  G4double eDensity = material->GetElectronDensity();
219  (p,kineticEnergy,cutEnergy,maxEnergy);
220  return cross;
221 }
G4double GetElectronDensity() const
Definition: G4Material.hh:217
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 145 of file G4BraggModel.cc.

148 {
149  // this method is called only for ions
150  G4double q2 = corr->EffectiveChargeSquareRatio(p,mat,kineticEnergy);
152  return q2*corr->EffectiveChargeCorrection(p,mat,kineticEnergy);
153 }
G4double EffectiveChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double EffectiveChargeCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4VEmFluctuationModel * GetModelOfFluctuations()
Definition: G4VEmModel.hh:612
double G4double
Definition: G4Types.hh:76
virtual void SetParticleAndCharge(const G4ParticleDefinition *, G4double q2)

Here is the call graph for this function:

G4double G4BraggModel::GetChargeSquareRatio ( ) const
inlineprotected

Definition at line 191 of file G4BraggModel.hh.

192 {
193  return chargeSquare;
194 }
G4double G4BraggModel::GetParticleCharge ( const G4ParticleDefinition p,
const G4Material mat,
G4double  kineticEnergy 
)
overridevirtual

Reimplemented from G4VEmModel.

Reimplemented in G4BraggIonGasModel.

Definition at line 157 of file G4BraggModel.cc.

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

Here is the call graph for this function:

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

Implements G4VEmModel.

Definition at line 115 of file G4BraggModel.cc.

117 {
118  if(p != particle) { SetParticle(p); }
119 
120  // always false before the run
121  SetDeexcitationFlag(false);
122 
123  if(IsMaster()) {
124  if(nullptr == fPSTAR) { fPSTAR = new G4PSTARStopping(); }
125  if(particle->GetPDGMass() < GeV) { fPSTAR->Initialise(); }
126  }
127 
128  if(nullptr == fParticleChange) {
129 
132  }
133  G4String pname = particle->GetParticleName();
134  if(particle->GetParticleType() == "nucleus" &&
135  pname != "deuteron" && pname != "triton" &&
136  pname != "alpha+" && pname != "helium" &&
137  pname != "hydrogen") { isIon = true; }
138 
139  fParticleChange = GetParticleChangeForLoss();
140  }
141 }
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:118
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:619
G4bool IsMaster() const
Definition: G4VEmModel.hh:725
const G4String & GetParticleName() const
G4bool UseAngularGeneratorFlag() const
Definition: G4VEmModel.hh:699
const G4String & GetParticleType() const
string pname
Definition: eplot.py:33
G4double GetPDGMass() const
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:626
static constexpr double GeV
Definition: G4SIunits.hh:217
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:788

Here is the call graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 342 of file G4BraggModel.cc.

344 {
345  if(pd != particle) { SetParticle(pd); }
346  G4double tau = kinEnergy/mass;
347  G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.) /
348  (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
349  return tmax;
350 }
float electron_mass_c2
Definition: hepunit.py:274
double G4double
Definition: G4Types.hh:76

Here is the caller graph for this function:

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

Implements G4VEmModel.

Definition at line 264 of file G4BraggModel.cc.

269 {
270  G4double tmax = MaxSecondaryKinEnergy(dp);
271  G4double xmax = std::min(tmax, maxEnergy);
272  if(xmin >= xmax) { return; }
273 
274  G4double kineticEnergy = dp->GetKineticEnergy();
275  G4double energy = kineticEnergy + mass;
276  G4double energy2 = energy*energy;
277  G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
278  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  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 = sqrt((1.0 - cost)*(1.0 + cost));
318 
319  G4double phi = twopi*rndmEngineMod->flat();
320 
321  deltaDirection.set(sint*cos(phi),sint*sin(phi), cost) ;
322  deltaDirection.rotateUz(dp->GetMomentumDirection());
323  }
324 
325  // create G4DynamicParticle object for delta ray
326  G4DynamicParticle* delta =
327  new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
328 
329  // Change kinematics of primary particle
330  kineticEnergy -= deltaKinEnergy;
331  G4ThreeVector finalP = dp->GetMomentum() - delta->GetMomentum();
332  finalP = finalP.unit();
333 
334  fParticleChange->SetProposedKineticEnergy(kineticEnergy);
335  fParticleChange->SetProposedMomentumDirection(finalP);
336 
337  vdp->push_back(delta);
338 }
void set(double x, double y, double z)
G4double MaxSecondaryKinEnergy(const G4DynamicParticle *dynParticle)
Definition: G4VEmModel.hh:484
G4double GetKineticEnergy() const
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:619
virtual double flat()=0
int G4int
Definition: G4Types.hh:78
static constexpr double twopi
Definition: G4SIunits.hh:76
G4double GetTotalMomentum() const
G4GLOB_DLL std::ostream G4cout
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
G4bool UseAngularGeneratorFlag() const
Definition: G4VEmModel.hh:699
const G4ThreeVector & GetMomentumDirection() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
void SetProposedKineticEnergy(G4double proposedKinEnergy)
float electron_mass_c2
Definition: hepunit.py:274
void SetProposedMomentumDirection(const G4ThreeVector &dir)
G4double energy(const ThreeVector &p, const G4double m)
Hep3Vector unit() const
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
virtual void flatArray(const int size, double *vect)=0
G4ThreeVector GetMomentum() const
const G4Material * GetMaterial() const
G4int SelectRandomAtomNumber(const G4Material *)
Definition: G4VEmModel.hh:564

Here is the call graph for this function:

void G4BraggModel::SetChargeSquareRatio ( G4double  val)
inlineprotected

Definition at line 196 of file G4BraggModel.hh.

197 {
198  chargeSquare = val;
199 }

Here is the caller graph for this function:


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