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

#include <G4DNADingfelderChargeIncreaseModel.hh>

Inheritance diagram for G4DNADingfelderChargeIncreaseModel:
Collaboration diagram for G4DNADingfelderChargeIncreaseModel:

Public Member Functions

 G4DNADingfelderChargeIncreaseModel (const G4ParticleDefinition *p=0, const G4String &nam="DNADingfelderChargeIncreaseModel")
 
virtual ~G4DNADingfelderChargeIncreaseModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
 
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 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)
 
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 Attributes

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

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 41 of file G4DNADingfelderChargeIncreaseModel.hh.

Constructor & Destructor Documentation

G4DNADingfelderChargeIncreaseModel::G4DNADingfelderChargeIncreaseModel ( const G4ParticleDefinition p = 0,
const G4String nam = "DNADingfelderChargeIncreaseModel" 
)

Definition at line 40 of file G4DNADingfelderChargeIncreaseModel.cc.

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

Definition at line 71 of file G4DNADingfelderChargeIncreaseModel.cc.

72 {}

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 198 of file G4DNADingfelderChargeIncreaseModel.cc.

203 {
204  if (verboseLevel > 3)
205  {
206  G4cout
207  << "Calling CrossSectionPerVolume() of G4DNADingfelderChargeIncreaseModel"
208  << G4endl;
209  }
210 
211  // Calculate total cross section for model
212 
213  G4DNAGenericIonsManager *instance;
215 
216  if (
217  particleDefinition != instance->GetIon("hydrogen")
218  &&
219  particleDefinition != instance->GetIon("alpha+")
220  &&
221  particleDefinition != instance->GetIon("helium")
222  )
223 
224  return 0;
225 
226  G4double lowLim = 0;
227  G4double highLim = 0;
228  G4double totalCrossSection = 0.;
229 
230  G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
231 
232  if(waterDensity!= 0.0)
233  // if (material == nistwater || material->GetBaseMaterial() == nistwater)
234  {
235  const G4String& particleName = particleDefinition->GetParticleName();
236 
237  std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
238  pos1 = lowEnergyLimit.find(particleName);
239 
240  if (pos1 != lowEnergyLimit.end())
241  {
242  lowLim = pos1->second;
243  }
244 
245  std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
246  pos2 = highEnergyLimit.find(particleName);
247 
248  if (pos2 != highEnergyLimit.end())
249  {
250  highLim = pos2->second;
251  }
252 
253  if (k >= lowLim && k < highLim)
254  {
255  //HYDROGEN
256  if (particleDefinition == instance->GetIon("hydrogen"))
257  {
258  const G4double aa = 2.835;
259  const G4double bb = 0.310;
260  const G4double cc = 2.100;
261  const G4double dd = 0.760;
262  const G4double fac = 1.0e-18;
263  const G4double rr = 13.606 * eV;
264 
266  G4double x = t / rr;
267  G4double temp = 4.0 * pi * Bohr_radius/nm * Bohr_radius/nm * fac;
268  G4double sigmal = temp * cc * (std::pow(x,dd));
269  G4double sigmah = temp * (aa * std::log(1.0 + x) + bb) / x;
270  totalCrossSection = 1.0/(1.0/sigmal + 1.0/sigmah) *m*m;
271  }
272  else
273  {
274  totalCrossSection = Sum(k,particleDefinition);
275  }
276  }
277 
278  if (verboseLevel > 2)
279  {
280  G4cout << "__________________________________" << G4endl;
281  G4cout << "G4DNADingfelderChargeIncreaseModel - XS INFO START" << G4endl;
282  G4cout << "Kinetic energy(eV)=" << k/eV << " particle : " << particleName << G4endl;
283  G4cout << "Cross section per water molecule (cm^2)=" << totalCrossSection/cm/cm << G4endl;
284  G4cout << "Cross section per water molecule (cm^-1)=" << totalCrossSection*waterDensity/(1./cm) << G4endl;
285  // G4cout << " - Cross section per water molecule (cm^-1)=" << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
286  G4cout << "G4DNADingfelderChargeIncreaseModel - XS INFO END" << G4endl;
287  }
288 
289  }
290 
291  return totalCrossSection*waterDensity;
292 // return totalCrossSection*material->GetAtomicNumDensityVector()[1];
293 
294 }
static constexpr double proton_mass_c2
static constexpr double Bohr_radius
size_t GetIndex() const
Definition: G4Material.hh:262
static constexpr double electron_mass_c2
static ulg bb
Definition: csz_inflate.cc:344
G4GLOB_DLL std::ostream G4cout
static constexpr double m
Definition: G4SIunits.hh:129
static constexpr double cm
Definition: G4SIunits.hh:119
static constexpr double eV
Definition: G4SIunits.hh:215
static G4DNAGenericIonsManager * Instance(void)
static constexpr double nm
Definition: G4SIunits.hh:112
static const G4double fac
#define G4endl
Definition: G4ios.hh:61
static constexpr double pi
Definition: G4SIunits.hh:75
double G4double
Definition: G4Types.hh:76
G4ParticleDefinition * GetIon(const G4String &name)

Here is the call graph for this function:

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

Implements G4VEmModel.

Definition at line 76 of file G4DNADingfelderChargeIncreaseModel.cc.

78 {
79 
80  if (verboseLevel > 3)
81  {
82  G4cout << "Calling G4DNADingfelderChargeIncreaseModel::Initialise()"
83  << G4endl;
84  }
85 
86  // Energy limits
87 
88  G4DNAGenericIonsManager *instance;
90  G4ParticleDefinition* hydrogenDef = instance->GetIon("hydrogen");
91  G4ParticleDefinition* alphaPlusDef = instance->GetIon("alpha+");
92  G4ParticleDefinition* heliumDef = instance->GetIon("helium");
93 
94  G4String hydrogen;
95  G4String alphaPlus;
96  G4String helium;
97 
98  // LIMITS
99 
100  hydrogen = hydrogenDef->GetParticleName();
101  lowEnergyLimit[hydrogen] = 100. * eV;
102  highEnergyLimit[hydrogen] = 100. * MeV;
103 
104  alphaPlus = alphaPlusDef->GetParticleName();
105  lowEnergyLimit[alphaPlus] = 1. * keV;
106  highEnergyLimit[alphaPlus] = 400. * MeV;
107 
108  helium = heliumDef->GetParticleName();
109  lowEnergyLimit[helium] = 1. * keV;
110  highEnergyLimit[helium] = 400. * MeV;
111 
112  //
113 
114  if (particle==hydrogenDef)
115  {
116  SetLowEnergyLimit(lowEnergyLimit[hydrogen]);
117  SetHighEnergyLimit(highEnergyLimit[hydrogen]);
118  }
119 
120  if (particle==alphaPlusDef)
121  {
122  SetLowEnergyLimit(lowEnergyLimit[alphaPlus]);
123  SetHighEnergyLimit(highEnergyLimit[alphaPlus]);
124  }
125 
126  if (particle==heliumDef)
127  {
128  SetLowEnergyLimit(lowEnergyLimit[helium]);
129  SetHighEnergyLimit(highEnergyLimit[helium]);
130  }
131 
132  // Final state
133 
134  //ALPHA+
135 
136  f0[0][0]=1.;
137  a0[0][0]=2.25;
138  a1[0][0]=-0.75;
139  b0[0][0]=-32.10;
140  c0[0][0]=0.600;
141  d0[0][0]=2.40;
142  x0[0][0]=4.60;
143 
144  x1[0][0]=-1.;
145  b1[0][0]=-1.;
146 
147  numberOfPartialCrossSections[0]=1;
148 
149  //HELIUM
150 
151  f0[0][1]=1.;
152  a0[0][1]=2.25;
153  a1[0][1]=-0.75;
154  b0[0][1]=-30.93;
155  c0[0][1]=0.590;
156  d0[0][1]=2.35;
157  x0[0][1]=4.29;
158 
159  f0[1][1]=1.;
160  a0[1][1]=2.25;
161  a1[1][1]=-0.75;
162  b0[1][1]=-32.61;
163  c0[1][1]=0.435;
164  d0[1][1]=2.70;
165  x0[1][1]=4.45;
166 
167  x1[0][1]=-1.;
168  b1[0][1]=-1.;
169 
170  x1[1][1]=-1.;
171  b1[1][1]=-1.;
172 
173  numberOfPartialCrossSections[1]=2;
174 
175  //
176 
177  if( verboseLevel>0 )
178  {
179  G4cout << "Dingfelder charge increase model is initialized " << G4endl
180  << "Energy range: "
181  << LowEnergyLimit() / keV << " keV - "
182  << HighEnergyLimit() / MeV << " MeV for "
183  << particle->GetParticleName()
184  << G4endl;
185  }
186 
187  // Initialize water density pointer
189 
190  if (isInitialised)
191  { return;}
193  isInitialised = true;
194 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:640
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:633
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:602
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:724
G4GLOB_DLL std::ostream G4cout
const std::vector< double > * GetNumMolPerVolTableFor(const G4Material *) const
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
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:731
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 G4DNADingfelderChargeIncreaseModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple ,
const G4DynamicParticle aDynamicParticle,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 298 of file G4DNADingfelderChargeIncreaseModel.cc.

304 {
305  if (verboseLevel > 3)
306  {
307  G4cout
308  << "Calling SampleSecondaries() of G4DNADingfelderChargeIncreaseModel"
309  << G4endl;
310  }
311 
313 
314  G4ParticleDefinition* definition = aDynamicParticle->GetDefinition();
315 
316  G4double particleMass = definition->GetPDGMass();
317 
318  G4double inK = aDynamicParticle->GetKineticEnergy();
319 
320  G4int finalStateIndex = RandomSelect(inK,definition);
321 
322  G4int n = NumberOfFinalStates(definition,finalStateIndex);
323 
324  G4double outK = 0.;
325 
326  if (!statCode) outK = inK - IncomingParticleBindingEnergyConstant(definition,finalStateIndex);
327 
328  else outK = inK;
329 
330  if (statCode)
331  fParticleChangeForGamma->ProposeLocalEnergyDeposit(IncomingParticleBindingEnergyConstant(definition,finalStateIndex));
332 
334 
335  G4DNAGenericIonsManager* instance;
337 
338  G4double electronK;
339  if (definition == instance->GetIon("hydrogen")) electronK = inK*electron_mass_c2/proton_mass_c2;
340  else electronK = inK*electron_mass_c2/(particleMass);
341 
342  if (outK<0)
343  {
344  G4Exception("G4DNADingfelderChargeIncreaseModel::SampleSecondaries","em0004",
345  FatalException,"Final kinetic energy is negative.");
346  }
347 
348  G4DynamicParticle* dp = new G4DynamicParticle(OutgoingParticleDefinition(definition,finalStateIndex),
349  aDynamicParticle->GetMomentumDirection(),
350  outK);
351 
352  fvect->push_back(dp);
353 
354  n = n - 1;
355 
356  while (n>0)
357  {
358  n--;
359  fvect->push_back(new G4DynamicParticle
360  (G4Electron::Electron(), aDynamicParticle->GetMomentumDirection(), electronK) );
361  }
362 }
static constexpr double proton_mass_c2
G4double GetKineticEnergy() const
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
static constexpr double electron_mass_c2
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
static G4DNAGenericIonsManager * Instance(void)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4double GetPDGMass() const
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
G4ParticleDefinition * GetIon(const G4String &name)

Here is the call graph for this function:

void G4DNADingfelderChargeIncreaseModel::SelectStationary ( G4bool  input)
inline

Definition at line 122 of file G4DNADingfelderChargeIncreaseModel.hh.

123 {
124  statCode = input;
125 }

Member Data Documentation

G4ParticleChangeForGamma* G4DNADingfelderChargeIncreaseModel::fParticleChangeForGamma
protected

Definition at line 69 of file G4DNADingfelderChargeIncreaseModel.hh.


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