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

#include <G4hCoulombScatteringModel.hh>

Inheritance diagram for G4hCoulombScatteringModel:
Collaboration diagram for G4hCoulombScatteringModel:

Public Member Functions

 G4hCoulombScatteringModel (G4bool combined=true)
 
virtual ~G4hCoulombScatteringModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel) override
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax) override
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double) final
 
void SetLowEnergyThreshold (G4double val)
 
void SetRecoilThreshold (G4double eth)
 
void SetFixedCut (G4double)
 
G4double GetFixedCut () const
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
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 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 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 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

void DefineMaterial (const G4MaterialCutsCouple *)
 
void SetupParticle (const G4ParticleDefinition *)
 
- 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 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 67 of file G4hCoulombScatteringModel.hh.

Constructor & Destructor Documentation

G4hCoulombScatteringModel::G4hCoulombScatteringModel ( G4bool  combined = true)
explicit

Definition at line 70 of file G4hCoulombScatteringModel.cc.

71  : G4VEmModel("hCoulombScattering"),
72  cosThetaMin(1.0),
73  cosThetaMax(-1.0),
74  isCombined(combined)
75 {
76  fParticleChange = nullptr;
77  fNistManager = G4NistManager::Instance();
79  theProton = G4Proton::Proton();
80  currentMaterial = nullptr;
81  fixedCut = -1.0;
82 
83  pCuts = 0;
84 
85  recoilThreshold = 0.*CLHEP::keV; // by default does not work
86 
87  particle = nullptr;
88  currentCouple = nullptr;
89  wokvi = new G4WentzelVIRelXSection(combined);
90 
91  currentMaterialIndex = 0;
92  mass = CLHEP::proton_mass_c2;
93  elecRatio = 0.0;
94 }
static constexpr double proton_mass_c2
static constexpr double keV
static G4NistManager * Instance()
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
G4IonTable * GetIonTable() const
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4ParticleTable * GetParticleTable()

Here is the call graph for this function:

G4hCoulombScatteringModel::~G4hCoulombScatteringModel ( )
virtual

Definition at line 98 of file G4hCoulombScatteringModel.cc.

99 {
100  delete wokvi;
101 }

Member Function Documentation

G4double G4hCoulombScatteringModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition p,
G4double  kinEnergy,
G4double  Z,
G4double  A,
G4double  cut,
G4double  emax 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 181 of file G4hCoulombScatteringModel.cc.

186 {
187  //G4cout << "### G4hCoulombScatteringModel::ComputeCrossSectionPerAtom for "
188  //<< p->GetParticleName()<<" Z= "<<Z<<" e(MeV)= "<< kinEnergy/MeV << G4endl;
189  G4double cross = 0.0;
190  elecRatio = 0.0;
191  if(p != particle) { SetupParticle(p); }
192 
193  // cross section is set to zero to avoid problems in sample secondary
194  if(kinEnergy <= 0.0) { return cross; }
196 
197  G4int iz = G4int(Z);
198  G4double tmass = proton_mass_c2;
199  if(1 < iz) {
200  tmass = fNistManager->GetAtomicMassAmu(iz)*amu_c2;
201  }
202  G4double costmin =
203  wokvi->SetupKinematic(kinEnergy, currentMaterial, cutEnergy, tmass);
204  if(cosThetaMax < costmin) {
205  G4double cut = cutEnergy;
206  if(fixedCut > 0.0) { cut = fixedCut; }
207  costmin = wokvi->SetupTarget(iz, cut);
208  G4double costmax = cosThetaMax;
209  if(iz == 1 && costmax < 0.0 && particle == theProton) {
210  costmax = 0.0;
211  }
212  if(costmin > costmax) {
213  cross = wokvi->ComputeNuclearCrossSection(costmin, costmax)
214  + wokvi->ComputeElectronCrossSection(costmin, costmax);
215  }
216  /*
217  if(p->GetParticleName() == "mu+")
218  G4cout << "e(MeV)= " << kinEnergy/MeV << " cross(b)= " << cross/barn
219  << " 1-costmin= " << 1-costmin
220  << " 1-costmax= " << 1-costmax
221  << " 1-cosThetaMax= " << 1-cosThetaMax
222  << G4endl;
223  */
224  }
225  return cross;
226 }
G4double ComputeNuclearCrossSection(G4double CosThetaMin, G4double CosThetaMax)
static constexpr double proton_mass_c2
void DefineMaterial(const G4MaterialCutsCouple *)
void SetupParticle(const G4ParticleDefinition *)
int G4int
Definition: G4Types.hh:78
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:451
G4double SetupKinematic(G4double kinEnergy, const G4Material *mat, G4double cut, G4double tmass)
G4double SetupTarget(G4int Z, G4double cut)
static constexpr double amu_c2
G4double GetAtomicMassAmu(const G4String &symb) const
G4double ComputeElectronCrossSection(G4double CosThetaMin, G4double CosThetaMax)
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

void G4hCoulombScatteringModel::DefineMaterial ( const G4MaterialCutsCouple cup)
inlineprotected

Definition at line 154 of file G4hCoulombScatteringModel.hh.

155 {
156  if(cup != currentCouple) {
157  currentCouple = cup;
158  currentMaterial = cup->GetMaterial();
159  currentMaterialIndex = currentCouple->GetIndex();
160  }
161 }
const G4Material * GetMaterial() const

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4hCoulombScatteringModel::GetFixedCut ( ) const
inline

Definition at line 192 of file G4hCoulombScatteringModel.hh.

193 {
194  return fixedCut;
195 }
void G4hCoulombScatteringModel::Initialise ( const G4ParticleDefinition part,
const G4DataVector cuts 
)
overridevirtual

Implements G4VEmModel.

Definition at line 105 of file G4hCoulombScatteringModel.cc.

107 {
108  SetupParticle(part);
109  currentCouple = 0;
110 
111  if(isCombined) {
112  cosThetaMin = 1.0;
114  if(tet >= pi) { cosThetaMin = -1.0; }
115  else if(tet > 0.0) { cosThetaMin = cos(tet); }
116  }
117 
118  wokvi->Initialise(part, cosThetaMin);
119  /*
120  G4cout << "G4hCoulombScatteringModel: " << particle->GetParticleName()
121  << " 1-cos(ThetaLimit)= " << 1 - cosThetaMin
122  << " cos(thetaMax)= " << cosThetaMax
123  << G4endl;
124  */
125  pCuts = &cuts;
126  //G4ProductionCutsTable::GetProductionCutsTable()->GetEnergyCutsVector(3);
127  /*
128  G4cout << "!!! G4hCoulombScatteringModel::Initialise for "
129  << part->GetParticleName() << " cos(TetMin)= " << cosThetaMin
130  << " cos(TetMax)= " << cosThetaMax <<G4endl;
131  G4cout << "cut= " << (*pCuts)[0] << " cut1= " << (*pCuts)[1] << G4endl;
132  */
133  if(!fParticleChange) {
134  fParticleChange = GetParticleChangeForGamma();
135  }
136  if(IsMaster() && mass < GeV && part->GetParticleName() != "GenericIon") {
137  InitialiseElementSelectors(part, cuts);
138  }
139 }
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:146
G4bool IsMaster() const
Definition: G4VEmModel.hh:717
void SetupParticle(const G4ParticleDefinition *)
static G4double tet[DIM]
G4double PolarAngleLimit() const
Definition: G4VEmModel.hh:661
static constexpr double pi
Definition: G4SIunits.hh:75
void Initialise(const G4ParticleDefinition *, G4double CosThetaLim)
double G4double
Definition: G4Types.hh:76
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:132

Here is the call graph for this function:

void G4hCoulombScatteringModel::InitialiseLocal ( const G4ParticleDefinition ,
G4VEmModel masterModel 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 143 of file G4hCoulombScatteringModel.cc.

145 {
147 }
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:801
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:809

Here is the call graph for this function:

G4double G4hCoulombScatteringModel::MinPrimaryEnergy ( const G4Material material,
const G4ParticleDefinition part,
G4double   
)
finalvirtual

Reimplemented from G4VEmModel.

Definition at line 152 of file G4hCoulombScatteringModel.cc.

155 {
156  SetupParticle(part);
157 
158  // define cut using cuts for proton
159  G4double cut =
160  std::max(recoilThreshold, (*pCuts)[CurrentCouple()->GetIndex()]);
161 
162  // find out lightest element
163  const G4ElementVector* theElementVector = material->GetElementVector();
164  G4int nelm = material->GetNumberOfElements();
165 
166  // select lightest element
167  G4int Z = 300;
168  for (G4int j=0; j<nelm; ++j) {
169  G4int iz = G4lrint((*theElementVector)[j]->GetZ());
170  if(iz < Z) { Z = iz; }
171  }
172  G4int A = G4lrint(fNistManager->GetAtomicMassAmu(Z));
173  G4double targetMass = G4NucleiProperties::GetNuclearMass(A, Z);
174  G4double t = std::max(cut, 0.5*(cut + sqrt(2*cut*targetMass)));
175 
176  return t;
177 }
static G4double GetNuclearMass(const G4double A, const G4double Z)
std::vector< G4Element * > G4ElementVector
void SetupParticle(const G4ParticleDefinition *)
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
int G4int
Definition: G4Types.hh:78
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:451
double A(double temperature)
int G4lrint(double ad)
Definition: templates.hh:163
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4double GetAtomicMassAmu(const G4String &symb) const
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

Implements G4VEmModel.

Definition at line 230 of file G4hCoulombScatteringModel.cc.

236 {
237  G4double kinEnergy = dp->GetKineticEnergy();
239  DefineMaterial(couple);
240 
241  // Choose nucleus
242  const G4Element* elm = SelectRandomAtom(couple,particle,
243  kinEnergy,cutEnergy,kinEnergy);
244 
245  G4int iz = elm->GetZasInt();
246  G4int ia = SelectIsotopeNumber(elm);
248 
249  wokvi->SetupKinematic(kinEnergy, currentMaterial, cutEnergy, mass2);
250  G4double costmin = wokvi->SetupTarget(iz, cutEnergy);
251  G4double costmax = cosThetaMax;
252  if(iz == 1 && costmax < 0.0 && particle == theProton) {
253  costmax = 0.0;
254  }
255 
256  if(costmin <= costmax) { return; }
257  G4double cross = wokvi->ComputeNuclearCrossSection(costmin, costmax);
258  G4double ecross = wokvi->ComputeElectronCrossSection(costmin, costmax);
259  G4double ratio = ecross/(cross + ecross);
260 
261  G4ThreeVector newDirection =
262  wokvi->SampleSingleScattering(costmin, costmax, ratio);
263 
264  // kinematics in the Lab system
265  G4double ptot = dp->GetTotalMomentum();
266  G4double e1 = dp->GetTotalEnergy();
267 
268  // Lab. system kinematics along projectile direction
269  G4LorentzVector v0 = G4LorentzVector(0, 0, ptot, e1 + mass2);
270  G4double bet = ptot/v0.e();
271  G4double gam = 1.0/sqrt((1.0 - bet)*(1.0 + bet));
272 
273  // CM projectile
274  G4double momCM = gam*(ptot - bet*e1);
275  G4double eCM = gam*(e1 - bet*ptot);
276  // energy & momentum after scattering of incident particle
277  G4double pxCM = momCM*newDirection.x();
278  G4double pyCM = momCM*newDirection.y();
279  G4double pzCM = momCM*newDirection.z();
280 
281  // CM--->Lab
282  G4LorentzVector v1(pxCM , pyCM, gam*(pzCM + bet*eCM), gam*(eCM + bet*pzCM));
283 
284  G4ThreeVector dir = dp->GetMomentumDirection();
285  newDirection = v1.vect().unit();
286  newDirection.rotateUz(dir);
287 
288  fParticleChange->ProposeMomentumDirection(newDirection);
289 
290  // recoil
291  v0 -= v1;
292  G4double trec = v0.e() - mass2;
293  G4double edep = 0.0;
294 
295  G4double tcut = recoilThreshold;
296  if(pCuts) {
297  tcut= std::max(tcut,(*pCuts)[currentMaterialIndex]);
298  //G4cout<<" tcut eV "<<tcut/eV<<endl;
299  }
300 
301  if(trec > tcut) {
302  G4ParticleDefinition* ion = theIonTable->GetIon(iz, ia, 0);
303  newDirection = v0.vect().unit();
304  newDirection.rotateUz(dir);
305  G4DynamicParticle* newdp = new G4DynamicParticle(ion, newDirection, trec);
306  fvect->push_back(newdp);
307  } else if(trec > 0.0) {
308  edep = trec;
309  fParticleChange->ProposeNonIonizingEnergyDeposit(edep);
310  }
311 
312  // finelize primary energy and energy balance
313  G4double finalT = v1.e() - mass;
314  if(finalT < 0.0) {
315  edep += finalT;
316  finalT = 0.0;
317  }
318  edep = std::max(edep, 0.0);
319  fParticleChange->SetProposedKineticEnergy(finalT);
320  fParticleChange->ProposeLocalEnergyDeposit(edep);
321 }
G4double ComputeNuclearCrossSection(G4double CosThetaMin, G4double CosThetaMax)
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4double GetKineticEnergy() const
G4double GetTotalEnergy() const
double x() const
void DefineMaterial(const G4MaterialCutsCouple *)
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:503
G4ParticleDefinition * GetDefinition() const
void SetupParticle(const G4ParticleDefinition *)
int G4int
Definition: G4Types.hh:78
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
double z() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double GetTotalMomentum() const
Hep3Vector vect() const
void ProposeNonIonizingEnergyDeposit(G4double anEnergyPart)
G4double SetupKinematic(G4double kinEnergy, const G4Material *mat, G4double cut, G4double tmass)
G4double SetupTarget(G4int Z, G4double cut)
const G4ThreeVector & GetMomentumDirection() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
G4int SelectIsotopeNumber(const G4Element *)
Definition: G4VEmModel.hh:584
T max(const T t1, const T t2)
brief Return the largest of the two arguments
Hep3Vector unit() const
G4int GetZasInt() const
Definition: G4Element.hh:132
double y() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
G4double ComputeElectronCrossSection(G4double CosThetaMin, G4double CosThetaMax)
double G4double
Definition: G4Types.hh:76
G4ThreeVector & SampleSingleScattering(G4double CosThetaMin, G4double CosThetaMax, G4double elecRatio)
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:541
CLHEP::HepLorentzVector G4LorentzVector

Here is the call graph for this function:

void G4hCoulombScatteringModel::SetFixedCut ( G4double  val)
inline

Definition at line 185 of file G4hCoulombScatteringModel.hh.

186 {
187  fixedCut = val;
188 }
void G4hCoulombScatteringModel::SetLowEnergyThreshold ( G4double  val)
inline
void G4hCoulombScatteringModel::SetRecoilThreshold ( G4double  eth)
inline

Definition at line 178 of file G4hCoulombScatteringModel.hh.

179 {
180  recoilThreshold = eth;
181 }
void G4hCoulombScatteringModel::SetupParticle ( const G4ParticleDefinition p)
inlineprotected

Definition at line 166 of file G4hCoulombScatteringModel.hh.

167 {
168  // Initialise mass and charge
169  if(p != particle) {
170  particle = p;
171  mass = particle->GetPDGMass();
172  wokvi->SetupParticle(p);
173  }
174 }
const char * p
Definition: xmltok.h:285
void SetupParticle(const G4ParticleDefinition *)
G4double GetPDGMass() const

Here is the call graph for this function:

Here is the caller graph for this function:


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