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

#include <G4eCoulombScatteringModel.hh>

Inheritance diagram for G4eCoulombScatteringModel:
Collaboration diagram for G4eCoulombScatteringModel:

Public Member Functions

 G4eCoulombScatteringModel (G4bool combined=true)
 
virtual ~G4eCoulombScatteringModel ()
 
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 78 of file G4eCoulombScatteringModel.hh.

Constructor & Destructor Documentation

G4eCoulombScatteringModel::G4eCoulombScatteringModel ( G4bool  combined = true)
explicit

Definition at line 80 of file G4eCoulombScatteringModel.cc.

81  : G4VEmModel("eCoulombScattering"),
82  cosThetaMin(1.0),
83  cosThetaMax(-1.0),
84  isCombined(combined)
85 {
86  fParticleChange = nullptr;
87  fNistManager = G4NistManager::Instance();
89  theProton = G4Proton::Proton();
90  currentMaterial = 0;
91  fixedCut = -1.0;
92 
93  pCuts = nullptr;
94 
95  recoilThreshold = 0.*keV; // by default does not work
96 
97  particle = nullptr;
98  currentCouple = nullptr;
99  wokvi = new G4WentzelOKandVIxSection(combined);
100 
101  currentMaterialIndex = 0;
102  mass = proton_mass_c2;
103  elecRatio = 0.0;
104 }
static constexpr double proton_mass_c2
static G4NistManager * Instance()
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
G4IonTable * GetIonTable() const
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4ParticleTable * GetParticleTable()
static constexpr double keV
Definition: G4SIunits.hh:216

Here is the call graph for this function:

G4eCoulombScatteringModel::~G4eCoulombScatteringModel ( )
virtual

Definition at line 108 of file G4eCoulombScatteringModel.cc.

109 {
110  delete wokvi;
111 }

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 190 of file G4eCoulombScatteringModel.cc.

195 {
196  //G4cout << "### G4eCoulombScatteringModel::ComputeCrossSectionPerAtom for "
197  //<< p->GetParticleName()<<" Z= "<<Z<<" e(MeV)= "<< kinEnergy/MeV << G4endl;
198  G4double cross = 0.0;
199  elecRatio = 0.0;
200  if(p != particle) { SetupParticle(p); }
201 
202  // cross section is set to zero to avoid problems in sample secondary
203  if(kinEnergy <= 0.0) { return cross; }
205  G4double costmin = wokvi->SetupKinematic(kinEnergy, currentMaterial);
206  if(cosThetaMax < costmin) {
207  G4int iz = G4lrint(Z);
208  G4double cut = cutEnergy;
209  if(fixedCut > 0.0) { cut = fixedCut; }
210  costmin = wokvi->SetupTarget(iz, cut);
211  G4double costmax = cosThetaMax;
212  if(iz == 1 && costmax < 0.0 && particle == theProton) {
213  costmax = 0.0;
214  }
215  if(costmin > costmax) {
216  cross = wokvi->ComputeNuclearCrossSection(costmin, costmax)
217  + wokvi->ComputeElectronCrossSection(costmin, costmax);
218  }
219  /*
220  if(p->GetParticleName() == "mu+")
221  G4cout << "e(MeV)= " << kinEnergy/MeV << " cross(b)= " << cross/barn
222  << " 1-costmin= " << 1-costmin
223  << " 1-costmax= " << 1-costmax
224  << " 1-cosThetaMax= " << 1-cosThetaMax
225  << G4endl;
226  */
227  }
228  return cross;
229 }
int G4int
Definition: G4Types.hh:78
G4double SetupTarget(G4int Z, G4double cut=DBL_MAX)
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:451
G4double SetupKinematic(G4double kinEnergy, const G4Material *mat)
void DefineMaterial(const G4MaterialCutsCouple *)
G4double ComputeElectronCrossSection(G4double CosThetaMin, G4double CosThetaMax)
int G4lrint(double ad)
Definition: templates.hh:163
G4double ComputeNuclearCrossSection(G4double CosThetaMin, G4double CosThetaMax)
double G4double
Definition: G4Types.hh:76
void SetupParticle(const G4ParticleDefinition *)

Here is the call graph for this function:

void G4eCoulombScatteringModel::DefineMaterial ( const G4MaterialCutsCouple cup)
inlineprotected

Definition at line 167 of file G4eCoulombScatteringModel.hh.

168 {
169  if(cup != currentCouple) {
170  currentCouple = cup;
171  currentMaterial = cup->GetMaterial();
172  currentMaterialIndex = currentCouple->GetIndex();
173  }
174 }
const G4Material * GetMaterial() const

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4eCoulombScatteringModel::GetFixedCut ( ) const
inline

Definition at line 205 of file G4eCoulombScatteringModel.hh.

206 {
207  return fixedCut;
208 }
void G4eCoulombScatteringModel::Initialise ( const G4ParticleDefinition part,
const G4DataVector cuts 
)
overridevirtual

Implements G4VEmModel.

Definition at line 115 of file G4eCoulombScatteringModel.cc.

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

Here is the call graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 153 of file G4eCoulombScatteringModel.cc.

155 {
157 }
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 G4eCoulombScatteringModel::MinPrimaryEnergy ( const G4Material material,
const G4ParticleDefinition part,
G4double   
)
finalvirtual

Reimplemented from G4VEmModel.

Definition at line 162 of file G4eCoulombScatteringModel.cc.

165 {
166  SetupParticle(part);
167 
168  // define cut using cuts for proton
169  G4double cut =
170  std::max(recoilThreshold, (*pCuts)[CurrentCouple()->GetIndex()]);
171 
172  // find out lightest element
173  const G4ElementVector* theElementVector = material->GetElementVector();
174  G4int nelm = material->GetNumberOfElements();
175 
176  // select lightest element
177  G4int Z = 300;
178  for (G4int j=0; j<nelm; ++j) {
179  Z = std::min(Z,(*theElementVector)[j]->GetZasInt());
180  }
181  G4int A = G4lrint(fNistManager->GetAtomicMassAmu(Z));
182  G4double targetMass = G4NucleiProperties::GetNuclearMass(A, Z);
183  G4double t = std::max(cut, 0.5*(cut + sqrt(2*cut*targetMass)));
184 
185  return t;
186 }
static G4double GetNuclearMass(const G4double A, const G4double Z)
std::vector< G4Element * > G4ElementVector
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
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
double G4double
Definition: G4Types.hh:76
void SetupParticle(const G4ParticleDefinition *)

Here is the call graph for this function:

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

Implements G4VEmModel.

Definition at line 233 of file G4eCoulombScatteringModel.cc.

239 {
240  G4double kinEnergy = dp->GetKineticEnergy();
242  DefineMaterial(couple);
243  /*
244  G4cout << "G4eCoulombScatteringModel::SampleSecondaries e(MeV)= "
245  << kinEnergy << " " << particle->GetParticleName()
246  << " cut= " << cutEnergy<< G4endl;
247  */
248  // Choose nucleus
249  G4double cut = cutEnergy;
250  if(fixedCut > 0.0) { cut = fixedCut; }
251 
252  wokvi->SetupKinematic(kinEnergy, currentMaterial);
253 
254  const G4Element* currentElement =
255  SelectRandomAtom(couple,particle,kinEnergy,cut,kinEnergy);
256 
257  G4int iz = currentElement->GetZasInt();
258 
259  G4double costmin = wokvi->SetupTarget(iz, cut);
260  G4double costmax = cosThetaMax;
261  if(iz == 1 && costmax < 0.0 && particle == theProton) {
262  costmax = 0.0;
263  }
264 
265  if(costmin > costmax) {
266  G4double cross = wokvi->ComputeNuclearCrossSection(costmin, costmax);
267  G4double ecross = wokvi->ComputeElectronCrossSection(costmin, costmax);
268  G4double ratio = ecross/(cross + ecross);
269 
270  G4int ia = SelectIsotopeNumber(currentElement);
271  G4double targetMass = G4NucleiProperties::GetNuclearMass(ia, iz);
272  wokvi->SetTargetMass(targetMass);
273 
274  G4ThreeVector newDirection =
275  wokvi->SampleSingleScattering(costmin, costmax, ratio);
276  G4double cost = newDirection.z();
277  /*
278  G4cout << "SampleSec: e(MeV)= " << kinEnergy/MeV
279  << " 1-costmin= " << 1-costmin
280  << " 1-costmax= " << 1-costmax
281  << " 1-cost= " << 1-cost
282  << " ratio= " << ratio
283  << G4endl;
284  */
285  G4ThreeVector direction = dp->GetMomentumDirection();
286  newDirection.rotateUz(direction);
287 
288  fParticleChange->ProposeMomentumDirection(newDirection);
289 
290  // recoil sampling assuming a small recoil
291  // and first order correction to primary 4-momentum
292  G4double mom2 = wokvi->GetMomentumSquare();
293  G4double trec = mom2*(1.0 - cost)
294  /(targetMass + (mass + kinEnergy)*(1.0 - cost));
295 
296  // the check likely not needed
297  if(trec > kinEnergy) { trec = kinEnergy; }
298  G4double finalT = kinEnergy - trec;
299  G4double edep = 0.0;
300  /*
301  G4cout<<"G4eCoulombScatteringModel: finalT= "<<finalT<<" Trec= "
302  <<trec << " Z= " << iz << " A= " << ia
303  << " tcut(keV)= " << (*pCuts)[currentMaterialIndex]/keV << G4endl;
304  */
305  G4double tcut = recoilThreshold;
306  if(pCuts) { tcut= std::max(tcut,(*pCuts)[currentMaterialIndex]); }
307 
308  if(trec > tcut) {
309  G4ParticleDefinition* ion = theIonTable->GetIon(iz, ia, 0);
310  G4ThreeVector dir = (direction*sqrt(mom2) -
311  newDirection*sqrt(finalT*(2*mass + finalT))).unit();
312  G4DynamicParticle* newdp = new G4DynamicParticle(ion, dir, trec);
313  fvect->push_back(newdp);
314  } else {
315  edep = trec;
316  fParticleChange->ProposeNonIonizingEnergyDeposit(edep);
317  }
318 
319  // finelize primary energy and energy balance
320  // this threshold may be applied only because for low-enegry
321  // e+e- msc model is applied
322  if(finalT < 0.0) {
323  edep += finalT;
324  finalT = 0.0;
325  if(edep < 0.0) { edep = 0.0; }
326  }
327  fParticleChange->SetProposedKineticEnergy(finalT);
328  fParticleChange->ProposeLocalEnergyDeposit(edep);
329  }
330 }
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4double GetKineticEnergy() const
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:503
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4double SetupTarget(G4int Z, G4double cut=DBL_MAX)
double z() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void ProposeNonIonizingEnergyDeposit(G4double anEnergyPart)
const G4ThreeVector & GetMomentumDirection() const
G4double SetupKinematic(G4double kinEnergy, const G4Material *mat)
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
void DefineMaterial(const G4MaterialCutsCouple *)
G4int SelectIsotopeNumber(const G4Element *)
Definition: G4VEmModel.hh:584
G4double ComputeElectronCrossSection(G4double CosThetaMin, G4double CosThetaMax)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4double ComputeNuclearCrossSection(G4double CosThetaMin, G4double CosThetaMax)
G4int GetZasInt() const
Definition: G4Element.hh:132
void SetProposedKineticEnergy(G4double proposedKinEnergy)
double G4double
Definition: G4Types.hh:76
void SetupParticle(const G4ParticleDefinition *)
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:541
G4ThreeVector & SampleSingleScattering(G4double CosThetaMin, G4double CosThetaMax, G4double elecRatio=0.0)

Here is the call graph for this function:

void G4eCoulombScatteringModel::SetFixedCut ( G4double  val)
inline

Definition at line 198 of file G4eCoulombScatteringModel.hh.

199 {
200  fixedCut = val;
201 }
void G4eCoulombScatteringModel::SetLowEnergyThreshold ( G4double  val)
inline
void G4eCoulombScatteringModel::SetRecoilThreshold ( G4double  eth)
inline

Definition at line 191 of file G4eCoulombScatteringModel.hh.

192 {
193  recoilThreshold = eth;
194 }
void G4eCoulombScatteringModel::SetupParticle ( const G4ParticleDefinition p)
inlineprotected

Definition at line 179 of file G4eCoulombScatteringModel.hh.

180 {
181  // Initialise mass and charge
182  if(p != particle) {
183  particle = p;
184  mass = particle->GetPDGMass();
185  wokvi->SetupParticle(p);
186  }
187 }
void SetupParticle(const G4ParticleDefinition *)
const char * p
Definition: xmltok.h:285
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: