Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4BraggIonModel Class Reference

#include <G4BraggIonModel.hh>

Inheritance diagram for G4BraggIonModel:
Collaboration diagram for G4BraggIonModel:

Public Member Functions

 G4BraggIonModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="BraggIon")
 
virtual ~G4BraggIonModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *couple) 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
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length) 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 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 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 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
 
- 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
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 68 of file G4BraggIonModel.hh.

Constructor & Destructor Documentation

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

Definition at line 78 of file G4BraggIonModel.cc.

80  : G4VEmModel(nam),
81  corr(nullptr),
82  particle(nullptr),
83  fParticleChange(nullptr),
84  currentMaterial(nullptr),
85  iMolecula(-1),
86  iASTAR(-1),
87  isIon(false)
88 {
90 
91  HeMass = 3.727417*GeV;
92  rateMassHe2p = HeMass/proton_mass_c2;
93  lowestKinEnergy = 1.0*keV/rateMassHe2p;
94  massFactor = 1000.*amu_c2/HeMass;
95  theZieglerFactor = eV*cm2*1.0e-15;
96  theElectron = G4Electron::Electron();
97  corrFactor = 1.0;
98  if(p) { SetParticle(p); }
99  else { SetParticle(theElectron); }
100 }
static constexpr double proton_mass_c2
static constexpr double cm2
Definition: G4SIunits.hh:120
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:724
static constexpr double eV
Definition: G4SIunits.hh:215
static constexpr double amu_c2
static constexpr double GeV
Definition: G4SIunits.hh:217
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:

G4BraggIonModel::~G4BraggIonModel ( )
virtual

Definition at line 104 of file G4BraggIonModel.cc.

105 {
106  if(IsMaster()) { delete fASTAR; fASTAR = nullptr; }
107 }
G4bool IsMaster() const
Definition: G4VEmModel.hh:717

Here is the call graph for this function:

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 207 of file G4BraggIonModel.cc.

213 {
215  (p,kineticEnergy,cutEnergy,maxEnergy);
216  return cross;
217 }
double G4double
Definition: G4Types.hh:76
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)

Here is the call graph for this function:

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

Definition at line 178 of file G4BraggIonModel.cc.

183 {
184  G4double cross = 0.0;
185  G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
186  G4double maxEnergy = std::min(tmax,maxKinEnergy);
187  if(cutEnergy < tmax) {
188 
189  G4double energy = kineticEnergy + mass;
190  G4double energy2 = energy*energy;
191  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 *= twopi_mc2_rcl2*chargeSquare/beta2;
198  }
199  // G4cout << "BR: e= " << kineticEnergy << " tmin= " << cutEnergy
200  // << " tmax= " << tmax << " cross= " << cross << G4endl;
201 
202  return cross;
203 }
static constexpr double twopi_mc2_rcl2
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 G4BraggIonModel::ComputeDEDXPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy 
)
overridevirtual

Reimplemented from G4VEmModel.

Reimplemented in G4BraggNoDeltaModel.

Definition at line 236 of file G4BraggIonModel.cc.

240 {
241  G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
242  G4double tmin = min(cutEnergy, tmax);
243  G4double tkin = kineticEnergy/massRate;
244  G4double dedx = 0.0;
245 
246  if(tkin < lowestKinEnergy) {
247  dedx = DEDX(material, lowestKinEnergy)*sqrt(tkin/lowestKinEnergy);
248  } else {
249  dedx = DEDX(material, tkin);
250  }
251 
252  if (cutEnergy < tmax) {
253 
254  G4double tau = kineticEnergy/mass;
255  G4double gam = tau + 1.0;
256  G4double bg2 = tau * (tau+2.0);
257  G4double beta2 = bg2/(gam*gam);
258  G4double x = tmin/tmax;
259 
260  dedx += (G4Log(x) + (1.0 - x)*beta2) * twopi_mc2_rcl2
261  * (material->GetElectronDensity())/beta2;
262  }
263 
264  // now compute the total ionization loss
265 
266  dedx = std::max(dedx, 0.0);
267  dedx *= chargeSquare;
268  /*
269  G4cout << "Bragg: tkin(MeV) = " << tkin/MeV << " dedx(MeVxcm^2/g) = "
270  << dedx*gram/(MeV*cm2*material->GetDensity())
271  << " q2 = " << chargeSquare << G4endl;
272  */
273  return dedx;
274 }
static constexpr double twopi_mc2_rcl2
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) final
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
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:

void G4BraggIonModel::CorrectionsAlongStep ( const G4MaterialCutsCouple couple,
const G4DynamicParticle dp,
G4double eloss,
G4double niel,
G4double  length 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 278 of file G4BraggIonModel.cc.

283 {
284  // this method is called only for ions
285  const G4ParticleDefinition* p = dp->GetDefinition();
286  const G4Material* mat = couple->GetMaterial();
287  G4double preKinEnergy = dp->GetKineticEnergy();
288  G4double e = preKinEnergy - eloss*0.5;
289  if(e < 0.0) { e = preKinEnergy*0.5; }
290 
291  G4double q2 = corr->EffectiveChargeSquareRatio(p,mat,e);
293  G4double qfactor = q2*corr->EffectiveChargeCorrection(p,mat,e)/corrFactor;
294  eloss *= qfactor;
295 
296  //G4cout << "G4BraggIonModel::CorrectionsAlongStep e= " << e
297  // << " qfactor= " << qfactor << " " << p->GetParticleName() <<G4endl;
298 }
G4double GetKineticEnergy() const
G4double EffectiveChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
const char * p
Definition: xmltok.h:285
G4double EffectiveChargeCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4ParticleDefinition * GetDefinition() const
G4VEmFluctuationModel * GetModelOfFluctuations()
Definition: G4VEmModel.hh:609
double G4double
Definition: G4Types.hh:76
virtual void SetParticleAndCharge(const G4ParticleDefinition *, G4double q2)
const G4Material * GetMaterial() const

Here is the call graph for this function:

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 221 of file G4BraggIonModel.cc.

227 {
228  G4double eDensity = material->GetElectronDensity();
230  (p,kineticEnergy,cutEnergy,maxEnergy);
231  return cross;
232 }
G4double GetElectronDensity() const
Definition: G4Material.hh:217
double G4double
Definition: G4Types.hh:76
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)

Here is the call graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 153 of file G4BraggIonModel.cc.

156 {
157  //G4cout<<"G4BraggIonModel::GetChargeSquareRatio e= "<<kineticEnergy<<G4endl;
158  // this method is called only for ions
159  G4double q2 = corr->EffectiveChargeSquareRatio(p,mat,kineticEnergy);
160  corrFactor = q2*corr->EffectiveChargeCorrection(p,mat,kineticEnergy);
161  return corrFactor;
162 }
G4double EffectiveChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double EffectiveChargeCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 166 of file G4BraggIonModel.cc.

169 {
170  //G4cout<<"G4BraggIonModel::GetParticleCharge e= "<<kineticEnergy <<
171  // " q= " << corr->GetParticleCharge(p,mat,kineticEnergy) <<G4endl;
172  // this method is called only for ions
173  return corr->GetParticleCharge(p,mat,kineticEnergy);
174 }
G4double GetParticleCharge(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)

Here is the call graph for this function:

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

Implements G4VEmModel.

Definition at line 111 of file G4BraggIonModel.cc.

113 {
114  if(p != particle) { SetParticle(p); }
115 
116  corrFactor = chargeSquare;
117 
118  // always false before the run
119  SetDeexcitationFlag(false);
120 
121  if(IsMaster()) {
122  if(nullptr == fASTAR) { fASTAR = new G4ASTARStopping(); }
123  if(particle->GetPDGMass() < GeV) { fASTAR->Initialise(); }
124  }
125 
126  if(nullptr == fParticleChange) {
127 
130  }
131  G4String pname = particle->GetParticleName();
132  if(particle->GetParticleType() == "nucleus" &&
133  pname != "deuteron" && pname != "triton" &&
134  pname != "alpha+" && pname != "helium" &&
135  pname != "hydrogen") { isIon = true; }
136 
138 
139  fParticleChange = GetParticleChangeForLoss();
140  }
141 }
static G4LossTableManager * Instance()
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:118
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:616
G4bool IsMaster() const
Definition: G4VEmModel.hh:717
const G4String & GetParticleName() const
G4bool UseAngularGeneratorFlag() const
Definition: G4VEmModel.hh:696
G4EmCorrections * EmCorrections()
const G4String & GetParticleType() const
G4double GetPDGMass() const
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:623
static constexpr double GeV
Definition: G4SIunits.hh:217
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:780

Here is the call graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 380 of file G4BraggIonModel.cc.

382 {
383  if(pd != particle) { SetParticle(pd); }
384  G4double tau = kinEnergy/mass;
385  G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.) /
386  (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
387  return tmax;
388 }
static constexpr double electron_mass_c2
double G4double
Definition: G4Types.hh:76

Here is the caller graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 145 of file G4BraggIonModel.cc.

147 {
148  return couple->GetMaterial()->GetIonisation()->GetMeanExcitationEnergy();
149 }
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:226
G4double GetMeanExcitationEnergy() const
const G4Material * GetMaterial() const

Here is the call graph for this function:

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

Implements G4VEmModel.

Definition at line 302 of file G4BraggIonModel.cc.

307 {
308  G4double tmax = MaxSecondaryKinEnergy(dp);
309  G4double xmax = std::min(tmax, maxEnergy);
310  if(xmin >= xmax) { return; }
311 
312  G4double kineticEnergy = dp->GetKineticEnergy();
313  G4double energy = kineticEnergy + mass;
314  G4double energy2 = energy*energy;
315  G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
316  G4double grej = 1.0;
317  G4double deltaKinEnergy, f;
318 
319  CLHEP::HepRandomEngine* rndmEngineMod = G4Random::getTheEngine();
320  G4double rndm[2];
321 
322  // sampling follows ...
323  do {
324  rndmEngineMod->flatArray(2, rndm);
325  deltaKinEnergy = xmin*xmax/(xmin*(1.0 - rndm[0]) + xmax*rndm[0]);
326 
327  f = 1.0 - beta2*deltaKinEnergy/tmax;
328 
329  if(f > grej) {
330  G4cout << "G4BraggIonModel::SampleSecondary Warning! "
331  << "Majorant " << grej << " < "
332  << f << " for e= " << deltaKinEnergy
333  << G4endl;
334  }
335 
336  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
337  } while( grej*rndm[1] >= f );
338 
339  G4ThreeVector deltaDirection;
340 
342  const G4Material* mat = couple->GetMaterial();
344 
345  deltaDirection =
346  GetAngularDistribution()->SampleDirection(dp, deltaKinEnergy, Z, mat);
347 
348  } else {
349 
350  G4double deltaMomentum =
351  sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
352  G4double cost = deltaKinEnergy * (energy + electron_mass_c2) /
353  (deltaMomentum * dp->GetTotalMomentum());
354  if(cost > 1.0) { cost = 1.0; }
355  G4double sint = sqrt((1.0 - cost)*(1.0 + cost));
356 
357  G4double phi = twopi*rndmEngineMod->flat();
358 
359  deltaDirection.set(sint*cos(phi),sint*sin(phi), cost) ;
360  deltaDirection.rotateUz(dp->GetMomentumDirection());
361  }
362 
363  // create G4DynamicParticle object for delta ray
364  G4DynamicParticle* delta =
365  new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
366 
367  vdp->push_back(delta);
368 
369  // Change kinematics of primary particle
370  kineticEnergy -= deltaKinEnergy;
371  G4ThreeVector finalP = dp->GetMomentum() - delta->GetMomentum();
372  finalP = finalP.unit();
373 
374  fParticleChange->SetProposedKineticEnergy(kineticEnergy);
375  fParticleChange->SetProposedMomentumDirection(finalP);
376 }
void set(double x, double y, double z)
G4double MaxSecondaryKinEnergy(const G4DynamicParticle *dynParticle)
Definition: G4VEmModel.hh:481
G4double GetKineticEnergy() const
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:616
virtual double flat()=0
int G4int
Definition: G4Types.hh:78
static constexpr double twopi
Definition: G4SIunits.hh:76
static constexpr double electron_mass_c2
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:696
const G4ThreeVector & GetMomentumDirection() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
void SetProposedKineticEnergy(G4double proposedKinEnergy)
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:561

Here is the call graph for this function:


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