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

#include <G4DNAIonElasticModel.hh>

Inheritance diagram for G4DNAIonElasticModel:
Collaboration diagram for G4DNAIonElasticModel:

Public Member Functions

 G4DNAIonElasticModel (const G4ParticleDefinition *p=0, const G4String &nam="DNAIonElasticModel")
 
virtual ~G4DNAIonElasticModel ()
 
virtual void Initialise (const G4ParticleDefinition *particuleDefinition, 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 SetKillBelowThreshold (G4double threshold)
 
G4double GetKillBelowThreshold ()
 
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 46 of file G4DNAIonElasticModel.hh.

Constructor & Destructor Documentation

G4DNAIonElasticModel::G4DNAIonElasticModel ( const G4ParticleDefinition p = 0,
const G4String nam = "DNAIonElasticModel" 
)

Definition at line 46 of file G4DNAIonElasticModel.cc.

47  :
48  G4VEmModel(nam), isInitialised(false)
49 {
50  //nistwater = G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER");
51 
52  killBelowEnergy = 100 * eV;
53  lowEnergyLimit = 0 * eV;
54  highEnergyLimit = 1 * MeV;
55  SetLowEnergyLimit(lowEnergyLimit);
56  SetHighEnergyLimit(highEnergyLimit);
57 
58  verboseLevel = 0;
59  // Verbosity scale:
60  // 0 = nothing
61  // 1 = warning for energy non-conservation
62  // 2 = details of energy budget
63  // 3 = calculation of cross sections, file openings, sampling of atoms
64  // 4 = entering in methods
65 
66  if(verboseLevel > 0)
67  {
68  G4cout << "Ion elastic model is constructed " << G4endl<< "Energy range: "
69  << lowEnergyLimit / eV << " eV - "
70  << highEnergyLimit / MeV << " MeV"
71  << G4endl;
72  }
74  fpMolWaterDensity = 0;
75  fpTableData = 0;
76  fParticle_Mass = -1;
77 
78  // Selection of stationary mode
79 
80  statCode = false;
81 }
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
G4GLOB_DLL std::ostream G4cout
G4ParticleChangeForGamma * fParticleChangeForGamma
static constexpr double eV
Definition: G4SIunits.hh:215
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:739

Here is the call graph for this function:

G4DNAIonElasticModel::~G4DNAIonElasticModel ( )
virtual

Definition at line 85 of file G4DNAIonElasticModel.cc.

86 {
87  // For total cross section
88  if(fpTableData) delete fpTableData;
89 }

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 248 of file G4DNAIonElasticModel.cc.

251 {
252  if(verboseLevel > 3)
253  {
254  G4cout << "Calling CrossSectionPerVolume() of G4DNAIonElasticModel"
255  << G4endl;
256  }
257 
258  // Calculate total cross section for model
259 
260  G4double sigma=0;
261 
262  G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
263 
264  if(waterDensity!= 0.0)
265  {
266  const G4String& particleName = p->GetParticleName();
267 
268  if (ekin < highEnergyLimit)
269  {
270  //SI : XS must not be zero otherwise sampling of secondaries method ignored
271  if (ekin < killBelowEnergy) return DBL_MAX;
272  //
273 
274  if (fpTableData != 0) // MK: is this check necessary?
275  {
276  sigma = fpTableData->FindValue(ekin);
277  }
278  else
279  {
280  G4Exception("G4DNAIonElasticModel::ComputeCrossSectionPerVolume","em0002",
281  FatalException,"Model not applicable to particle type.");
282  }
283  }
284 
285  if (verboseLevel > 2)
286  {
287  G4cout << "__________________________________" << G4endl;
288  G4cout << "G4DNAIonElasticModel - XS INFO START" << G4endl;
289  G4cout << "Kinetic energy(eV)=" << ekin/eV << " particle : " << particleName << G4endl;
290  G4cout << "Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
291  G4cout << "Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
292  G4cout << "G4DNAIonElasticModel - XS INFO END" << G4endl;
293  }
294 
295  }
296 
297  return sigma*waterDensity;
298 }
size_t GetIndex() const
Definition: G4Material.hh:262
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
static constexpr double cm
Definition: G4SIunits.hh:119
static constexpr double eV
Definition: G4SIunits.hh:215
virtual G4double FindValue(G4double e, G4int componentId=0) const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83

Here is the call graph for this function:

G4double G4DNAIonElasticModel::GetKillBelowThreshold ( )
inline

Definition at line 75 of file G4DNAIonElasticModel.hh.

76  {
77  return killBelowEnergy;
78  }
void G4DNAIonElasticModel::Initialise ( const G4ParticleDefinition particuleDefinition,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 94 of file G4DNAIonElasticModel.cc.

97 {
98 
99  if(verboseLevel > 3)
100  {
101  G4cout << "Calling G4DNAIonElasticModel::Initialise()" << G4endl;
102  }
103 
104  // Energy limits
105 
106  if (LowEnergyLimit() < lowEnergyLimit)
107  {
108  G4cout << "G4DNAIonElasticModel: low energy limit increased from " <<
109  LowEnergyLimit()/eV << " eV to " << lowEnergyLimit/eV << " eV" << G4endl;
110  SetLowEnergyLimit(lowEnergyLimit);
111  }
112 
113  if (HighEnergyLimit() > highEnergyLimit)
114  {
115  G4cout << "G4DNAIonElasticModel: high energy limit decreased from " <<
116  HighEnergyLimit()/MeV << " MeV to " << highEnergyLimit/MeV << " MeV" << G4endl;
117  SetHighEnergyLimit(highEnergyLimit);
118  }
119 
120  // Reading of data files
121 
122  G4double scaleFactor = 1e-16*cm*cm;
123 
124  char *path = getenv("G4LEDATA");
125 
126  if (!path)
127  {
128  G4Exception("G4IonElasticModel::Initialise","em0006",
129  FatalException,"G4LEDATA environment variable not set.");
130  return;
131  }
132 
133  G4String totalXSFile;
134  std::ostringstream fullFileName;
135 
138  G4ParticleDefinition* protonDef =
140  G4ParticleDefinition* hydrogenDef = instance->GetIon("hydrogen");
141  G4ParticleDefinition* heliumDef = instance->GetIon("helium");
142  G4ParticleDefinition* alphaplusDef = instance->GetIon("alpha+");
143  G4ParticleDefinition* alphaplusplusDef = instance->GetIon("alpha++");
144  G4String proton, hydrogen, helium, alphaplus, alphaplusplus;
145 
146  if (
147  (particleDefinition == protonDef && protonDef != 0)
148  ||
149  (particleDefinition == hydrogenDef && hydrogenDef != 0)
150  )
151  {
152  // For total cross section of p,h
153  fParticle_Mass = 1.;
154  totalXSFile = "dna/sigma_elastic_proton_HTran";
155 
156  // For final state
157  fullFileName << path << "/dna/sigmadiff_cumulated_elastic_proton_HTran.dat";
158  }
159 
160  if (
161  (particleDefinition == instance->GetIon("helium") && heliumDef)
162  ||
163  (particleDefinition == instance->GetIon("alpha+") && alphaplusDef)
164  ||
165  (particleDefinition == instance->GetIon("alpha++") && alphaplusplusDef)
166  )
167  {
168  fParticle_Mass = 4.;
169  // For total cross section of he,he+,he++
170  totalXSFile = "dna/sigma_elastic_alpha_HTran";
171 
172  // For final state
173  fullFileName << path << "/dna/sigmadiff_cumulated_elastic_alpha_HTran.dat";
174  }
175 
176  fpTableData = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, eV,scaleFactor );
177  fpTableData->LoadData(totalXSFile);
178  std::ifstream diffCrossSection(fullFileName.str().c_str());
179 
180  if (!diffCrossSection)
181  {
182  G4ExceptionDescription description;
183  description << "Missing data file:"
184  <<fullFileName.str().c_str()<< G4endl;
185  G4Exception("G4IonElasticModel::Initialise","em0003",
187  description);
188  }
189 
190  // Added clear for MT
191 
192  eTdummyVec.clear();
193  eVecm.clear();
194  fDiffCrossSectionData.clear();
195 
196  //
197 
198  eTdummyVec.push_back(0.);
199 
200  while(!diffCrossSection.eof())
201  {
202  double tDummy;
203  double eDummy;
204  diffCrossSection>>tDummy>>eDummy;
205 
206  // SI : mandatory eVecm initialization
207 
208  if (tDummy != eTdummyVec.back())
209  {
210  eTdummyVec.push_back(tDummy);
211  eVecm[tDummy].push_back(0.);
212  }
213 
214  diffCrossSection>>fDiffCrossSectionData[tDummy][eDummy];
215 
216  if (eDummy != eVecm[tDummy].back()) eVecm[tDummy].push_back(eDummy);
217 
218  }
219 
220  // End final state
221  if( verboseLevel>0 )
222  {
223  if (verboseLevel > 2)
224  {
225  G4cout << "Loaded cross section files for ion elastic model" << G4endl;
226  }
227  G4cout << "Ion elastic model is initialized " << G4endl
228  << "Energy range: "
229  << LowEnergyLimit() / eV << " eV - "
230  << HighEnergyLimit() / MeV << " MeV"
231  << G4endl;
232  }
233 
234  // Initialize water density pointer
236  fpMolWaterDensity = G4DNAMolecularMaterial::Instance()->
237  GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
238 
239  if (isInitialised)
240  { return;}
242  isInitialised = true;
243 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:643
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:636
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:602
virtual G4bool LoadData(const G4String &argFileName)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
G4GLOB_DLL std::ostream G4cout
G4ParticleChangeForGamma * fParticleChangeForGamma
static constexpr double cm
Definition: G4SIunits.hh:119
static constexpr double eV
Definition: G4SIunits.hh:215
static G4DNAGenericIonsManager * Instance(void)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4DNAMolecularMaterial * Instance()
static G4ParticleTable * GetParticleTable()
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
static MCTruthManager * instance
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:739
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:132
G4ParticleDefinition * GetIon(const G4String &name)

Here is the call graph for this function:

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

Implements G4VEmModel.

Definition at line 303 of file G4DNAIonElasticModel.cc.

307 {
308 
309  if(verboseLevel > 3)
310  {
311  G4cout << "Calling SampleSecondaries() of G4DNAIonElasticModel" << G4endl;
312  }
313 
314  G4double particleEnergy0 = aDynamicParticle->GetKineticEnergy();
315 
316  if (particleEnergy0 < killBelowEnergy)
317  {
321  return;
322  }
323 
324  if (particleEnergy0>= killBelowEnergy && particleEnergy0 < highEnergyLimit)
325  {
326  G4double water_mass = 18.;
327 
328  G4double thetaCM = RandomizeThetaCM(particleEnergy0, aDynamicParticle->GetDefinition());
329 
330  //HT:convert to laboratory system
331 
332  G4double theta = std::atan(std::sin(thetaCM*CLHEP::pi/180)
333  /(fParticle_Mass/water_mass+std::cos(thetaCM*CLHEP::pi/180)));
334 
335  G4double cosTheta= std::cos(theta);
336 
337  //
338 
339  G4double phi = 2. * CLHEP::pi * G4UniformRand();
340 
341  G4ThreeVector zVers = aDynamicParticle->GetMomentumDirection();
342  G4ThreeVector xVers = zVers.orthogonal();
343  G4ThreeVector yVers = zVers.cross(xVers);
344 
345  G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
346  G4double yDir = xDir;
347  xDir *= std::cos(phi);
348  yDir *= std::sin(phi);
349 
350  G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
351 
353 
354  G4double depositEnergyCM = 0;
355 
356  //HT: deposited energy
357  depositEnergyCM = 4. * particleEnergy0 * fParticle_Mass * water_mass *
358  (1-std::cos(thetaCM*CLHEP::pi/180))
359  / (2 * std::pow((fParticle_Mass+water_mass),2));
360 
361  //SI: added protection particleEnergy0 >= depositEnergyCM
362  if (!statCode && (particleEnergy0 >= depositEnergyCM) )
363 
364  fParticleChangeForGamma->SetProposedKineticEnergy(particleEnergy0 - depositEnergyCM);
365 
366  else fParticleChangeForGamma->SetProposedKineticEnergy(particleEnergy0);
367 
369  }
370 }
G4double GetKineticEnergy() const
G4ParticleDefinition * GetDefinition() const
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
G4ParticleChangeForGamma * fParticleChangeForGamma
const G4ThreeVector & GetMomentumDirection() const
Hep3Vector orthogonal() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
Hep3Vector cross(const Hep3Vector &) const
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
static constexpr double pi
Definition: SystemOfUnits.h:54

Here is the call graph for this function:

void G4DNAIonElasticModel::SelectStationary ( G4bool  input)
inline

Definition at line 148 of file G4DNAIonElasticModel.hh.

149 {
150  statCode = input;
151 }
void G4DNAIonElasticModel::SetKillBelowThreshold ( G4double  threshold)

Definition at line 508 of file G4DNAIonElasticModel.cc.

509 {
510  killBelowEnergy = threshold;
511 
512  if(killBelowEnergy < 100 * eV)
513  {
514  G4cout << "*** WARNING : the G4DNAIonElasticModel class is not "
515  "activated below 100 eV !"
516  << G4endl;
517  }
518 }
G4GLOB_DLL std::ostream G4cout
static constexpr double eV
Definition: G4SIunits.hh:215
#define G4endl
Definition: G4ios.hh:61

Member Data Documentation

G4ParticleChangeForGamma* G4DNAIonElasticModel::fParticleChangeForGamma
protected

Definition at line 84 of file G4DNAIonElasticModel.hh.


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