Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4eCoulombScatteringModel Class Reference

#include <G4eCoulombScatteringModel.hh>

+ Inheritance 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 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

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 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 77 of file G4eCoulombScatteringModel.hh.

Constructor & Destructor Documentation

◆ G4eCoulombScatteringModel()

G4eCoulombScatteringModel::G4eCoulombScatteringModel ( G4bool  combined = true)
explicit

Definition at line 68 of file G4eCoulombScatteringModel.cc.

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

◆ ~G4eCoulombScatteringModel()

G4eCoulombScatteringModel::~G4eCoulombScatteringModel ( )
virtual

Definition at line 97 of file G4eCoulombScatteringModel.cc.

98{
99 delete wokvi;
100}

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 173 of file G4eCoulombScatteringModel.cc.

178{
179 /*
180 G4cout << "### G4eCoulombScatteringModel::ComputeCrossSectionPerAtom for "
181 << p->GetParticleName()<<" Z= "<<Z<<" e(MeV)= "<< kinEnergy/MeV
182 << G4endl;
183 */
184 G4double cross = 0.0;
185 elecRatio = 0.0;
186 if(p != particle) { SetupParticle(p); }
187
188 // cross section is set to zero to avoid problems in sample secondary
189 if(kinEnergy <= 0.0) { return cross; }
191 G4double costmin = wokvi->SetupKinematic(kinEnergy, currentMaterial);
192
193 //G4cout << "cosThetaMax= "<<cosThetaMax<<" costmin= "<<costmin<< G4endl;
194
195 if(cosThetaMax < costmin) {
196 G4int iz = G4lrint(Z);
197 G4double cut = (0.0 < fixedCut) ? fixedCut : cutEnergy;
198 costmin = wokvi->SetupTarget(iz, cut);
199 //G4cout << "SetupTarget: Z= " << iz << " cut= " << cut << " "
200 // << costmin << G4endl;
201 G4double costmax = (1 == iz && particle == theProton && cosThetaMax < 0.0)
202 ? 0.0 : cosThetaMax;
203 if(costmin > costmax) {
204 cross = wokvi->ComputeNuclearCrossSection(costmin, costmax)
205 + wokvi->ComputeElectronCrossSection(costmin, costmax);
206 }
207 /*
208 if(p->GetParticleName() == "e-")
209 G4cout << "Z= " << Z << " e(MeV)= " << kinEnergy/MeV
210 << " cross(b)= " << cross/barn << " 1-costmin= " << 1-costmin
211 << " 1-costmax= " << 1-costmax
212 << " 1-cosThetaMax= " << 1-cosThetaMax
213 << " " << currentMaterial->GetName()
214 << G4endl;
215 */
216 }
217 //G4cout << "====== cross= " << cross << G4endl;
218 return cross;
219}
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:480
G4double SetupTarget(G4int Z, G4double cut)
G4double ComputeElectronCrossSection(G4double CosThetaMin, G4double CosThetaMax)
G4double SetupKinematic(G4double kinEnergy, const G4Material *mat)
G4double ComputeNuclearCrossSection(G4double CosThetaMin, G4double CosThetaMax)
void DefineMaterial(const G4MaterialCutsCouple *)
void SetupParticle(const G4ParticleDefinition *)
int G4lrint(double ad)
Definition: templates.hh:134

◆ DefineMaterial()

void G4eCoulombScatteringModel::DefineMaterial ( const G4MaterialCutsCouple cup)
inlineprotected

Definition at line 166 of file G4eCoulombScatteringModel.hh.

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

Referenced by ComputeCrossSectionPerAtom(), and SampleSecondaries().

◆ GetFixedCut()

G4double G4eCoulombScatteringModel::GetFixedCut ( ) const
inline

Definition at line 204 of file G4eCoulombScatteringModel.hh.

205{
206 return fixedCut;
207}

◆ Initialise()

void G4eCoulombScatteringModel::Initialise ( const G4ParticleDefinition part,
const G4DataVector cuts 
)
overridevirtual

Implements G4VEmModel.

Definition at line 104 of file G4eCoulombScatteringModel.cc.

106{
107 SetupParticle(part);
108 currentCouple = nullptr;
109
110 // defined theta limit between single and multiple scattering
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 pCuts = &cuts;
120 /*
121 G4cout << "G4eCoulombScatteringModel::Initialise for "
122 << part->GetParticleName() << " 1-cos(TetMin)= " << 1.0 - cosThetaMin
123 << " 1-cos(TetMax)= " << 1. - cosThetaMax << G4endl;
124 G4cout << "cut[0]= " << (*pCuts)[0] << G4endl;
125 */
126 if(!fParticleChange) {
127 fParticleChange = GetParticleChangeForGamma();
128 }
129 if(IsMaster() && mass < GeV && part->GetParticleName() != "GenericIon") {
130 InitialiseElementSelectors(part, cuts);
131 }
132}
G4double PolarAngleLimit() const
Definition: G4VEmModel.hh:673
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:133
G4bool IsMaster() const
Definition: G4VEmModel.hh:736
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:148
void Initialise(const G4ParticleDefinition *, G4double CosThetaLim)

◆ InitialiseLocal()

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

Reimplemented from G4VEmModel.

Definition at line 136 of file G4eCoulombScatteringModel.cc.

138{
140}
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:842
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:834

◆ MinPrimaryEnergy()

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

Reimplemented from G4VEmModel.

Definition at line 145 of file G4eCoulombScatteringModel.cc.

148{
149 SetupParticle(part);
150
151 // define cut using cuts for proton
152 G4double cut =
153 std::max(recoilThreshold, (*pCuts)[CurrentCouple()->GetIndex()]);
154
155 // find out lightest element
156 const G4ElementVector* theElementVector = material->GetElementVector();
157 G4int nelm = material->GetNumberOfElements();
158
159 // select lightest element
160 G4int Z = 300;
161 for (G4int j=0; j<nelm; ++j) {
162 Z = std::min(Z,(*theElementVector)[j]->GetZasInt());
163 }
164 G4int A = G4lrint(fNistManager->GetAtomicMassAmu(Z));
166 G4double t = std::max(cut, 0.5*(cut + sqrt(2*cut*targetMass)));
167
168 return t;
169}
double A(double temperature)
std::vector< G4Element * > G4ElementVector
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
G4double GetAtomicMassAmu(const G4String &symb) const
static G4double GetNuclearMass(const G4double A, const G4double Z)

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 223 of file G4eCoulombScatteringModel.cc.

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

◆ SetFixedCut()

void G4eCoulombScatteringModel::SetFixedCut ( G4double  val)
inline

Definition at line 197 of file G4eCoulombScatteringModel.hh.

198{
199 fixedCut = val;
200}

◆ SetLowEnergyThreshold()

void G4eCoulombScatteringModel::SetLowEnergyThreshold ( G4double  val)
inline

◆ SetRecoilThreshold()

void G4eCoulombScatteringModel::SetRecoilThreshold ( G4double  eth)
inline

Definition at line 190 of file G4eCoulombScatteringModel.hh.

191{
192 recoilThreshold = eth;
193}

◆ SetupParticle()

void G4eCoulombScatteringModel::SetupParticle ( const G4ParticleDefinition p)
inlineprotected

Definition at line 178 of file G4eCoulombScatteringModel.hh.

179{
180 // Initialise mass and charge
181 if(p != particle) {
182 particle = p;
183 mass = particle->GetPDGMass();
184 wokvi->SetupParticle(p);
185 }
186}
void SetupParticle(const G4ParticleDefinition *)

Referenced by ComputeCrossSectionPerAtom(), Initialise(), MinPrimaryEnergy(), and SampleSecondaries().


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