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

#include <G4DNAMillerGreenExcitationModel.hh>

+ Inheritance diagram for G4DNAMillerGreenExcitationModel:

Public Member Functions

 G4DNAMillerGreenExcitationModel (const G4ParticleDefinition *p=0, const G4String &nam="DNAMillerGreenExcitationModel")
 
virtual ~G4DNAMillerGreenExcitationModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int, const G4ParticleDefinition *, G4double)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
void SelectStationary (G4bool input)
 
- 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 Attributes

G4ParticleChangeForGammafParticleChangeForGamma
 
- 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
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Detailed Description

Definition at line 43 of file G4DNAMillerGreenExcitationModel.hh.

Constructor & Destructor Documentation

◆ G4DNAMillerGreenExcitationModel()

G4DNAMillerGreenExcitationModel::G4DNAMillerGreenExcitationModel ( const G4ParticleDefinition p = 0,
const G4String nam = "DNAMillerGreenExcitationModel" 
)

Definition at line 40 of file G4DNAMillerGreenExcitationModel.cc.

42:G4VEmModel(nam),isInitialised(false)
43{
44 fpMolWaterDensity = 0;
45
46 nLevels=0;
47 kineticEnergyCorrection[0]=0.;
48 kineticEnergyCorrection[1]=0.;
49 kineticEnergyCorrection[2]=0.;
50 kineticEnergyCorrection[3]=0.;
51
52 verboseLevel= 0;
53 // Verbosity scale:
54 // 0 = nothing
55 // 1 = warning for energy non-conservation
56 // 2 = details of energy budget
57 // 3 = calculation of cross sections, file openings, sampling of atoms
58 // 4 = entering in methods
59
60 if( verboseLevel>0 )
61 {
62 G4cout << "Miller & Green excitation model is constructed " << G4endl;
63 }
65
66 // Selection of stationary mode
67
68 statCode = false;
69}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4ParticleChangeForGamma * fParticleChangeForGamma

◆ ~G4DNAMillerGreenExcitationModel()

G4DNAMillerGreenExcitationModel::~G4DNAMillerGreenExcitationModel ( )
virtual

Definition at line 73 of file G4DNAMillerGreenExcitationModel.cc.

74{}

Member Function Documentation

◆ CrossSectionPerVolume()

G4double G4DNAMillerGreenExcitationModel::CrossSectionPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  ekin,
G4double  emin,
G4double  emax 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 224 of file G4DNAMillerGreenExcitationModel.cc.

229{
230 if (verboseLevel > 3)
231 G4cout << "Calling CrossSectionPerVolume() of G4DNAMillerGreenExcitationModel" << G4endl;
232
233 // Calculate total cross section for model
234
235 G4DNAGenericIonsManager *instance;
237
238 if (
239 particleDefinition != G4Proton::ProtonDefinition()
240 &&
241 particleDefinition != instance->GetIon("hydrogen")
242 &&
243 particleDefinition != instance->GetIon("alpha++")
244 &&
245 particleDefinition != instance->GetIon("alpha+")
246 &&
247 particleDefinition != instance->GetIon("helium")
248 )
249
250 return 0;
251
252 G4double lowLim = 0;
253 G4double highLim = 0;
254 G4double crossSection = 0.;
255
256 G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
257
258 const G4String& particleName = particleDefinition->GetParticleName();
259
260 std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
261 pos1 = lowEnergyLimit.find(particleName);
262
263 if (pos1 != lowEnergyLimit.end())
264 {
265 lowLim = pos1->second;
266 }
267
268 std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
269 pos2 = highEnergyLimit.find(particleName);
270
271 if (pos2 != highEnergyLimit.end())
272 {
273 highLim = pos2->second;
274 }
275
276 if (k >= lowLim && k <= highLim)
277 {
278 crossSection = Sum(k,particleDefinition);
279
280 // add ONE or TWO electron-water excitation for alpha+ and helium
281 /*
282 if ( particleDefinition == instance->GetIon("alpha+")
283 ||
284 particleDefinition == instance->GetIon("helium")
285 )
286 {
287
288 G4DNAEmfietzoglouExcitationModel * excitationXS = new G4DNAEmfietzoglouExcitationModel();
289 excitationXS->Initialise(G4Electron::ElectronDefinition());
290
291 G4double sigmaExcitation=0;
292 G4double tmp =0.;
293
294 if (k*0.511/3728 > 8.23*eV && k*0.511/3728 < 10*MeV ) sigmaExcitation =
295 excitationXS->CrossSectionPerVolume(material,G4Electron::ElectronDefinition(),k*0.511/3728,tmp,tmp)
296 /material->GetAtomicNumDensityVector()[1];
297
298 if ( particleDefinition == instance->GetIon("alpha+") )
299 crossSection = crossSection + sigmaExcitation ;
300
301 if ( particleDefinition == instance->GetIon("helium") )
302 crossSection = crossSection + 2*sigmaExcitation ;
303
304 delete excitationXS;
305
306 // Alternative excitation model
307
308 G4DNABornExcitationModel * excitationXS = new G4DNABornExcitationModel();
309 excitationXS->Initialise(G4Electron::ElectronDefinition());
310
311 G4double sigmaExcitation=0;
312 G4double tmp=0;
313
314 if (k*0.511/3728 > 9*eV && k*0.511/3728 < 1*MeV ) sigmaExcitation =
315 excitationXS->CrossSectionPerVolume(material,G4Electron::ElectronDefinition(),k*0.511/3728,tmp,tmp)
316 /material->GetAtomicNumDensityVector()[1];
317
318 if ( particleDefinition == instance->GetIon("alpha+") )
319 crossSection = crossSection + sigmaExcitation ;
320
321 if ( particleDefinition == instance->GetIon("helium") )
322 crossSection = crossSection + 2*sigmaExcitation ;
323
324 delete excitationXS;
325
326 }
327 */
328
329 }
330
331 if (verboseLevel > 2)
332 {
333 G4cout << "__________________________________" << G4endl;
334 G4cout << "G4DNAMillerGreenExcitationModel - XS INFO START" << G4endl;
335 G4cout << "Kinetic energy(eV)=" << k/eV << " particle : " << particleDefinition->GetParticleName() << G4endl;
336 G4cout << "Cross section per water molecule (cm^2)=" << crossSection/cm/cm << G4endl;
337 G4cout << "Cross section per water molecule (cm^-1)=" << crossSection*waterDensity/(1./cm) << G4endl;
338 // G4cout << " - Cross section per water molecule (cm^-1)=" << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
339 G4cout << "G4DNAMillerGreenExcitationModel - XS INFO END" << G4endl;
340 }
341
342 return crossSection*waterDensity;
343}
double G4double
Definition: G4Types.hh:83
static G4DNAGenericIonsManager * Instance(void)
G4ParticleDefinition * GetIon(const G4String &name)
size_t GetIndex() const
Definition: G4Material.hh:258
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:87

◆ GetPartialCrossSection()

G4double G4DNAMillerGreenExcitationModel::GetPartialCrossSection ( const G4Material ,
G4int  level,
const G4ParticleDefinition particleDefinition,
G4double  kineticEnergy 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 387 of file G4DNAMillerGreenExcitationModel.cc.

391{
392 return PartialCrossSection(kineticEnergy, level, particleDefinition);
393}

◆ Initialise()

void G4DNAMillerGreenExcitationModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 78 of file G4DNAMillerGreenExcitationModel.cc.

80{
81
82 if (verboseLevel > 3)
83 G4cout << "Calling G4DNAMillerGreenExcitationModel::Initialise()" << G4endl;
84
85 // Energy limits
86
90 G4ParticleDefinition* hydrogenDef = instance->GetIon("hydrogen");
91 G4ParticleDefinition* alphaPlusPlusDef = instance->GetIon("alpha++");
92 G4ParticleDefinition* alphaPlusDef = instance->GetIon("alpha+");
93 G4ParticleDefinition* heliumDef = instance->GetIon("helium");
94
96 G4String hydrogen;
97 G4String alphaPlusPlus;
98 G4String alphaPlus;
99 G4String helium;
100
101 // LIMITS AND CONSTANTS
102
103 proton = protonDef->GetParticleName();
104 lowEnergyLimit[proton] = 10. * eV;
105 highEnergyLimit[proton] = 500. * keV;
106
107 kineticEnergyCorrection[0] = 1.;
108 slaterEffectiveCharge[0][0] = 0.;
109 slaterEffectiveCharge[1][0] = 0.;
110 slaterEffectiveCharge[2][0] = 0.;
111 sCoefficient[0][0] = 0.;
112 sCoefficient[1][0] = 0.;
113 sCoefficient[2][0] = 0.;
114
115 hydrogen = hydrogenDef->GetParticleName();
116 lowEnergyLimit[hydrogen] = 10. * eV;
117 highEnergyLimit[hydrogen] = 500. * keV;
118
119 kineticEnergyCorrection[0] = 1.;
120 slaterEffectiveCharge[0][0] = 0.;
121 slaterEffectiveCharge[1][0] = 0.;
122 slaterEffectiveCharge[2][0] = 0.;
123 sCoefficient[0][0] = 0.;
124 sCoefficient[1][0] = 0.;
125 sCoefficient[2][0] = 0.;
126
127 alphaPlusPlus = alphaPlusPlusDef->GetParticleName();
128 lowEnergyLimit[alphaPlusPlus] = 1. * keV;
129 highEnergyLimit[alphaPlusPlus] = 400. * MeV;
130
131 kineticEnergyCorrection[1] = 0.9382723/3.727417;
132 slaterEffectiveCharge[0][1]=0.;
133 slaterEffectiveCharge[1][1]=0.;
134 slaterEffectiveCharge[2][1]=0.;
135 sCoefficient[0][1]=0.;
136 sCoefficient[1][1]=0.;
137 sCoefficient[2][1]=0.;
138
139 alphaPlus = alphaPlusDef->GetParticleName();
140 lowEnergyLimit[alphaPlus] = 1. * keV;
141 highEnergyLimit[alphaPlus] = 400. * MeV;
142
143 kineticEnergyCorrection[2] = 0.9382723/3.727417;
144 slaterEffectiveCharge[0][2]=2.0;
145
146 // Following values provided by M. Dingfelder
147 slaterEffectiveCharge[1][2]=2.00;
148 slaterEffectiveCharge[2][2]=2.00;
149 //
150 sCoefficient[0][2]=0.7;
151 sCoefficient[1][2]=0.15;
152 sCoefficient[2][2]=0.15;
153
154 helium = heliumDef->GetParticleName();
155 lowEnergyLimit[helium] = 1. * keV;
156 highEnergyLimit[helium] = 400. * MeV;
157
158 kineticEnergyCorrection[3] = 0.9382723/3.727417;
159 slaterEffectiveCharge[0][3]=1.7;
160 slaterEffectiveCharge[1][3]=1.15;
161 slaterEffectiveCharge[2][3]=1.15;
162 sCoefficient[0][3]=0.5;
163 sCoefficient[1][3]=0.25;
164 sCoefficient[2][3]=0.25;
165
166 //
167
168 if (particle==protonDef)
169 {
170 SetLowEnergyLimit(lowEnergyLimit[proton]);
171 SetHighEnergyLimit(highEnergyLimit[proton]);
172 }
173
174 if (particle==hydrogenDef)
175 {
176 SetLowEnergyLimit(lowEnergyLimit[hydrogen]);
177 SetHighEnergyLimit(highEnergyLimit[hydrogen]);
178 }
179
180 if (particle==alphaPlusPlusDef)
181 {
182 SetLowEnergyLimit(lowEnergyLimit[alphaPlusPlus]);
183 SetHighEnergyLimit(highEnergyLimit[alphaPlusPlus]);
184 }
185
186 if (particle==alphaPlusDef)
187 {
188 SetLowEnergyLimit(lowEnergyLimit[alphaPlus]);
189 SetHighEnergyLimit(highEnergyLimit[alphaPlus]);
190 }
191
192 if (particle==heliumDef)
193 {
194 SetLowEnergyLimit(lowEnergyLimit[helium]);
195 SetHighEnergyLimit(highEnergyLimit[helium]);
196 }
197
198 //
199
200 nLevels = waterExcitation.NumberOfLevels();
201
202 //
203 if( verboseLevel>0 )
204 {
205 G4cout << "Miller & Green excitation model is initialized " << G4endl
206 << "Energy range: "
207 << LowEnergyLimit() / eV << " eV - "
208 << HighEnergyLimit() / keV << " keV for "
209 << particle->GetParticleName()
210 << G4endl;
211 }
212
213 // Initialize water density pointer
215
216 if (isInitialised) { return; }
218 isInitialised = true;
219
220}
const std::vector< G4double > * GetNumMolPerVolTableFor(const G4Material *) const
Retrieve a table of molecular densities (number of molecules per unit volume) in the G4 unit system f...
static G4DNAMolecularMaterial * Instance()
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:651
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:757
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:133
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:652
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:764

◆ SampleSecondaries()

void G4DNAMillerGreenExcitationModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  ,
const G4MaterialCutsCouple ,
const G4DynamicParticle aDynamicParticle,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 347 of file G4DNAMillerGreenExcitationModel.cc.

352{
353
354 if (verboseLevel > 3)
355 G4cout << "Calling SampleSecondaries() of G4DNAMillerGreenExcitationModel" << G4endl;
356
357 G4double particleEnergy0 = aDynamicParticle->GetKineticEnergy();
358
359 G4int level = RandomSelect(particleEnergy0,aDynamicParticle->GetDefinition());
360
361 // Dingfelder's excitation levels
362 const G4double excitation[]={ 8.17*eV, 10.13*eV, 11.31*eV, 12.91*eV, 14.50*eV};
363 G4double excitationEnergy = excitation[level];
364
365 G4double newEnergy = 0.;
366
367 if (!statCode) newEnergy = particleEnergy0 - excitationEnergy;
368
369 else newEnergy = particleEnergy0;
370
371 if (newEnergy>0)
372 {
376
377 const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
379 level, theIncomingTrack);
380
381 }
382
383}
@ eExcitedMolecule
int G4int
Definition: G4Types.hh:85
static G4DNAChemistryManager * Instance()
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4Track * GetCurrentTrack() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

◆ SelectStationary()

void G4DNAMillerGreenExcitationModel::SelectStationary ( G4bool  input)
inline

Definition at line 138 of file G4DNAMillerGreenExcitationModel.hh.

139{
140 statCode = input;
141}

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNAMillerGreenExcitationModel::fParticleChangeForGamma
protected

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