Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4DNAMillerGreenExcitationModel Class Reference

#include <G4DNAMillerGreenExcitationModel.hh>

Inheritance diagram for G4DNAMillerGreenExcitationModel:
Collaboration 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 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 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)
 
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 SetFluctuationFlag (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 Attributes

G4ParticleChangeForGammafParticleChangeForGamma
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 

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 *)
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 44 of file G4DNAMillerGreenExcitationModel.hh.

Constructor & Destructor Documentation

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

Definition at line 42 of file G4DNAMillerGreenExcitationModel.cc.

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

Definition at line 76 of file G4DNAMillerGreenExcitationModel.cc.

77 {}

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 227 of file G4DNAMillerGreenExcitationModel.cc.

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

Here is the call graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 408 of file G4DNAMillerGreenExcitationModel.cc.

412 {
413  return PartialCrossSection(kineticEnergy, level, particleDefinition);
414 }
void G4DNAMillerGreenExcitationModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 81 of file G4DNAMillerGreenExcitationModel.cc.

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

Here is the call graph for this function:

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

Implements G4VEmModel.

Definition at line 365 of file G4DNAMillerGreenExcitationModel.cc.

370 {
371 
372  if (verboseLevel > 3)
373  G4cout << "Calling SampleSecondaries() of G4DNAMillerGreenExcitationModel" << G4endl;
374 
375  G4double particleEnergy0 = aDynamicParticle->GetKineticEnergy();
376 
377  G4int level = RandomSelect(particleEnergy0,aDynamicParticle->GetDefinition());
378 
379  // G4double excitationEnergy = waterExcitation.ExcitationEnergy(level);
380 
381  // Dingfelder's excitation levels
382  const G4double excitation[]={ 8.17*eV, 10.13*eV, 11.31*eV, 12.91*eV, 14.50*eV};
383  G4double excitationEnergy = excitation[level];
384 
385  G4double newEnergy = 0.;
386 
387  if (!statCode) newEnergy = particleEnergy0 - excitationEnergy;
388 
389  else newEnergy = particleEnergy0;
390 
391  if (newEnergy>0)
392  {
396 
397  const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
399  level,
400  theIncomingTrack);
401 
402  }
403 
404 }
G4double GetKineticEnergy() const
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4GLOB_DLL std::ostream G4cout
G4ParticleChangeForGamma * fParticleChangeForGamma
const G4ThreeVector & GetMomentumDirection() const
static constexpr double eV
Definition: G4SIunits.hh:215
static G4DNAChemistryManager * Instance()
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
const G4Track * GetCurrentTrack() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

void G4DNAMillerGreenExcitationModel::SelectStationary ( G4bool  input)
inline

Definition at line 139 of file G4DNAMillerGreenExcitationModel.hh.

140 {
141  statCode = input;
142 }

Member Data Documentation

G4ParticleChangeForGamma* G4DNAMillerGreenExcitationModel::fParticleChangeForGamma
protected

Definition at line 77 of file G4DNAMillerGreenExcitationModel.hh.


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