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

#include <G4DNADingfelderChargeDecreaseModel.hh>

Inheritance diagram for G4DNADingfelderChargeDecreaseModel:
Collaboration diagram for G4DNADingfelderChargeDecreaseModel:

Public Member Functions

 G4DNADingfelderChargeDecreaseModel (const G4ParticleDefinition *p=0, const G4String &nam="DNADingfelderChargeDecreaseModel")
 
virtual ~G4DNADingfelderChargeDecreaseModel ()
 
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 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 40 of file G4DNADingfelderChargeDecreaseModel.hh.

Constructor & Destructor Documentation

G4DNADingfelderChargeDecreaseModel::G4DNADingfelderChargeDecreaseModel ( const G4ParticleDefinition p = 0,
const G4String nam = "DNADingfelderChargeDecreaseModel" 
)

Definition at line 41 of file G4DNADingfelderChargeDecreaseModel.cc.

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

Definition at line 72 of file G4DNADingfelderChargeDecreaseModel.cc.

73 {
74 }

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 211 of file G4DNADingfelderChargeDecreaseModel.cc.

216 {
217  if (verboseLevel > 3)
218  {
219  G4cout
220  << "Calling CrossSectionPerVolume() of G4DNADingfelderChargeDecreaseModel"
221  << G4endl;
222  }
223 
224  // Calculate total cross section for model
225 
228 
229  if (
230  particleDefinition != G4Proton::ProtonDefinition()
231  &&
232  particleDefinition != instance->GetIon("alpha++")
233  &&
234  particleDefinition != instance->GetIon("alpha+")
235  )
236 
237  return 0;
238 
239  G4double lowLim = 0;
240  G4double highLim = 0;
241  G4double crossSection = 0.;
242 
243  G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
244 
245  if(waterDensity!= 0.0)
246  // if (material == nistwater || material->GetBaseMaterial() == nistwater)
247  {
248  const G4String& particleName = particleDefinition->GetParticleName();
249 
250  std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
251  pos1 = lowEnergyLimit.find(particleName);
252 
253  if (pos1 != lowEnergyLimit.end())
254  {
255  lowLim = pos1->second;
256  }
257 
258  std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
259  pos2 = highEnergyLimit.find(particleName);
260 
261  if (pos2 != highEnergyLimit.end())
262  {
263  highLim = pos2->second;
264  }
265 
266  if (k >= lowLim && k < highLim)
267  {
268  crossSection = Sum(k,particleDefinition);
269  }
270 
271  if (verboseLevel > 2)
272  {
273  G4cout << "_______________________________________" << G4endl;
274  G4cout << "G4DNADingfelderChargeDecreaeModel" << G4endl;
275  G4cout << "Kinetic energy(eV)=" << k/eV << "particle :" << particleName << G4endl;
276  G4cout << "Cross section per water molecule (cm^2)=" << crossSection/cm/cm << G4endl;
277  G4cout << "Cross section per water molecule (cm^-1)=" << crossSection*
278  waterDensity/(1./cm) << G4endl;
279 // material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
280  }
281  }
282 
283  return crossSection*waterDensity;
284 // return crossSection*material->GetAtomicNumDensityVector()[1];
285 
286 }
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:

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

Implements G4VEmModel.

Definition at line 78 of file G4DNADingfelderChargeDecreaseModel.cc.

80 {
81 
82  if (verboseLevel > 3)
83  {
84  G4cout << "Calling G4DNADingfelderChargeDecreaseModel::Initialise()"
85  << G4endl;
86  }
87 
88  // Energy limits
89 
93  G4ParticleDefinition* alphaPlusPlusDef = instance->GetIon("alpha++");
94  G4ParticleDefinition* alphaPlusDef = instance->GetIon("alpha+");
95 
97  G4String alphaPlusPlus;
98  G4String alphaPlus;
99 
100  // LIMITS
101 
102  proton = protonDef->GetParticleName();
103  lowEnergyLimit[proton] = 100. * eV;
104  highEnergyLimit[proton] = 100. * MeV;
105 
106  alphaPlusPlus = alphaPlusPlusDef->GetParticleName();
107  lowEnergyLimit[alphaPlusPlus] = 1. * keV;
108  highEnergyLimit[alphaPlusPlus] = 400. * MeV;
109 
110  alphaPlus = alphaPlusDef->GetParticleName();
111  lowEnergyLimit[alphaPlus] = 1. * keV;
112  highEnergyLimit[alphaPlus] = 400. * MeV;
113 
114  //
115 
116  if (particle==protonDef)
117  {
118  SetLowEnergyLimit(lowEnergyLimit[proton]);
119  SetHighEnergyLimit(highEnergyLimit[proton]);
120  }
121 
122  if (particle==alphaPlusPlusDef)
123  {
124  SetLowEnergyLimit(lowEnergyLimit[alphaPlusPlus]);
125  SetHighEnergyLimit(highEnergyLimit[alphaPlusPlus]);
126  }
127 
128  if (particle==alphaPlusDef)
129  {
130  SetLowEnergyLimit(lowEnergyLimit[alphaPlus]);
131  SetHighEnergyLimit(highEnergyLimit[alphaPlus]);
132  }
133 
134  // Final state
135 
136  //PROTON
137  f0[0][0]=1.;
138  a0[0][0]=-0.180;
139  a1[0][0]=-3.600;
140  b0[0][0]=-18.22;
141  b1[0][0]=-1.997;
142  c0[0][0]=0.215;
143  d0[0][0]=3.550;
144  x0[0][0]=3.450;
145  x1[0][0]=5.251;
146 
147  numberOfPartialCrossSections[0] = 1;
148 
149  //ALPHA++
150  f0[0][1]=1.; a0[0][1]=0.95;
151  a1[0][1]=-2.75;
152  b0[0][1]=-23.00;
153  c0[0][1]=0.215;
154  d0[0][1]=2.95;
155  x0[0][1]=3.50;
156 
157  f0[1][1]=1.;
158  a0[1][1]=0.95;
159  a1[1][1]=-2.75;
160  b0[1][1]=-23.73;
161  c0[1][1]=0.250;
162  d0[1][1]=3.55;
163  x0[1][1]=3.72;
164 
165  x1[0][1]=-1.;
166  b1[0][1]=-1.;
167 
168  x1[1][1]=-1.;
169  b1[1][1]=-1.;
170 
171  numberOfPartialCrossSections[1] = 2;
172 
173  // ALPHA+
174  f0[0][2]=1.;
175  a0[0][2]=0.65;
176  a1[0][2]=-2.75;
177  b0[0][2]=-21.81;
178  c0[0][2]=0.232;
179  d0[0][2]=2.95;
180  x0[0][2]=3.53;
181 
182  x1[0][2]=-1.;
183  b1[0][2]=-1.;
184 
185  numberOfPartialCrossSections[2] = 1;
186 
187  //
188 
189  if( verboseLevel>0 )
190  {
191  G4cout << "Dingfelder charge decrease model is initialized " << G4endl
192  << "Energy range: "
193  << LowEnergyLimit() / keV << " keV - "
194  << HighEnergyLimit() / MeV << " MeV for "
195  << particle->GetParticleName()
196  << G4endl;
197  }
198 
199  // Initialize water density pointer
201 
202  if (isInitialised)
203  { return;}
205  isInitialised = true;
206 
207 }
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
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 G4DNADingfelderChargeDecreaseModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple ,
const G4DynamicParticle aDynamicParticle,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 290 of file G4DNADingfelderChargeDecreaseModel.cc.

296 {
297  if (verboseLevel > 3)
298  {
299  G4cout
300  << "Calling SampleSecondaries() of G4DNADingfelderChargeDecreaseModel"
301  << G4endl;
302  }
303 
304  G4double inK = aDynamicParticle->GetKineticEnergy();
305 
306  G4ParticleDefinition* definition = aDynamicParticle->GetDefinition();
307 
308  G4double particleMass = definition->GetPDGMass();
309 
310  G4int finalStateIndex = RandomSelect(inK,definition);
311 
312  G4int n = NumberOfFinalStates(definition, finalStateIndex);
313  G4double waterBindingEnergy = WaterBindingEnergyConstant(definition, finalStateIndex);
314  G4double outgoingParticleBindingEnergy = OutgoingParticleBindingEnergyConstant(definition, finalStateIndex);
315 
316  G4double outK = 0.;
317 
318  if (definition==G4Proton::Proton())
319  {
320  if (!statCode) outK = inK - n*(inK*electron_mass_c2/proton_mass_c2)
321  - waterBindingEnergy + outgoingParticleBindingEnergy;
322  else outK = inK;
323  }
324 
325  else
326  {
327  if (!statCode) outK = inK - n*(inK*electron_mass_c2/particleMass)
328  - waterBindingEnergy + outgoingParticleBindingEnergy;
329  else outK = inK;
330  }
331 
332  if (outK<0)
333  {
334  G4Exception("G4DNADingfelderChargeDecreaseModel::SampleSecondaries","em0004",
335  FatalException,"Final kinetic energy is negative.");
336  }
337 
339 
340  if (!statCode) fParticleChangeForGamma->ProposeLocalEnergyDeposit(waterBindingEnergy);
341 
342  else
343 
344  {
345  if (definition==G4Proton::Proton())
347  + waterBindingEnergy - outgoingParticleBindingEnergy);
348  else
350  + waterBindingEnergy - outgoingParticleBindingEnergy);
351  }
352 
353  G4DynamicParticle* dp = new G4DynamicParticle (OutgoingParticleDefinition(definition, finalStateIndex),
354  aDynamicParticle->GetMomentumDirection(),
355  outK);
356  fvect->push_back(dp);
357 
358  const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
360  1,
361  theIncomingTrack);
362 
363 }
G4double GetKineticEnergy() const
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
static G4Proton * Proton()
Definition: G4Proton.cc:93
float proton_mass_c2
Definition: hepunit.py:275
float electron_mass_c2
Definition: hepunit.py:274
const G4int n
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4DNAChemistryManager * Instance()
G4double GetPDGMass() const
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
const G4Track * GetCurrentTrack() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)

Here is the call graph for this function:

void G4DNADingfelderChargeDecreaseModel::SelectStationary ( G4bool  input)
inline

Definition at line 122 of file G4DNADingfelderChargeDecreaseModel.hh.

123 {
124  statCode = input;
125 }

Member Data Documentation

G4ParticleChangeForGamma* G4DNADingfelderChargeDecreaseModel::fParticleChangeForGamma
protected

Definition at line 68 of file G4DNADingfelderChargeDecreaseModel.hh.


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