Geant4  10.02.p03
G4DNAOneStepThermalizationModel Class Reference

#include <G4DNAOneStepThermalizationModel.hh>

Inheritance diagram for G4DNAOneStepThermalizationModel:
Collaboration diagram for G4DNAOneStepThermalizationModel:

Public Member Functions

 G4DNAOneStepThermalizationModel (const G4ParticleDefinition *p=0, const G4String &nam="DNASancheSolvatationModel")
 
virtual ~G4DNAOneStepThermalizationModel ()
 
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 SetVerbose (int)
 
- 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, const G4ParticleDefinition *, G4double)
 
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 *> *)
 
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=0)
 
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 Member Functions

G4ThreeVector RadialDistributionOfProducts (G4double Rrms) const
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLoss * GetParticleChangeForLoss ()
 
G4ParticleChangeForGamma * GetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Protected Attributes

const std::vector< G4double > * fpWaterDensity
 
G4ParticleChangeForGamma * fParticleChangeForGamma
 
G4bool fIsInitialised
 
G4int fVerboseLevel
 
G4ITNavigator * fNavigator
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChange * pParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

Private Member Functions

G4DNAOneStepThermalizationModeloperator= (const G4DNAOneStepThermalizationModel &right)
 
 G4DNAOneStepThermalizationModel (const G4DNAOneStepThermalizationModel &)
 

Additional Inherited Members

- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

When an electron reaches the highest energy domain of G4DNAOneStepThermalizationModel, it is then automatically converted into a solvated electron and displace from its original position using a published thermalization statistic.

Article: Jintana Meesungnoen, Jean-Paul Jay-Gerin, Abdelali Filali-Mouhim, and Samlee Mankhetkorn (2002) Low-Energy Electron Penetration Range in Liquid Water. Radiation Research: November 2002, Vol. 158, No. 5, pp. 657-660.

Definition at line 65 of file G4DNAOneStepThermalizationModel.hh.

Constructor & Destructor Documentation

◆ G4DNAOneStepThermalizationModel() [1/2]

G4DNAOneStepThermalizationModel::G4DNAOneStepThermalizationModel ( const G4ParticleDefinition p = 0,
const G4String nam = "DNASancheSolvatationModel" 
)

Definition at line 54 of file G4DNAOneStepThermalizationModel.cc.

55  :
56  G4VEmModel(nam), fIsInitialised(false)
57 {
58  fVerboseLevel = 0;
60  G4DNAWaterExcitationStructure exStructure;
61  SetHighEnergyLimit(exStructure.ExcitationEnergy(0));
63  fpWaterDensity = 0;
64  fNavigator = 0;
65 }
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:69
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
const std::vector< G4double > * fpWaterDensity
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
G4ParticleChangeForGamma * fParticleChangeForGamma
Here is the call graph for this function:

◆ ~G4DNAOneStepThermalizationModel()

G4DNAOneStepThermalizationModel::~G4DNAOneStepThermalizationModel ( )
virtual

Definition at line 69 of file G4DNAOneStepThermalizationModel.cc.

70 {
71  if(fNavigator)
72  {
73  if(fNavigator->GetNavigatorState())
74  delete fNavigator->GetNavigatorState();
75  delete fNavigator;
76  }
77 }
Here is the call graph for this function:

◆ G4DNAOneStepThermalizationModel() [2/2]

G4DNAOneStepThermalizationModel::G4DNAOneStepThermalizationModel ( const G4DNAOneStepThermalizationModel )
private

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 123 of file G4DNAOneStepThermalizationModel.cc.

128 {
129 #ifdef G4VERBOSE
130  if(fVerboseLevel > 1)
131  G4cout << "Calling CrossSectionPerVolume() of G4DNAOneStepThermalizationModel"
132  << G4endl;
133 #endif
134 
135  if(ekin > HighEnergyLimit())
136  {
137  return 0.0;
138  }
139 
140  G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
141 
142  if(waterDensity!= 0.0)
143  {
144 // if (ekin <= HighEnergyLimit()) // already tested
145  {
146  return DBL_MAX;
147  }
148  }
149  return 0.;
150 }
size_t GetIndex() const
Definition: G4Material.hh:262
G4GLOB_DLL std::ostream G4cout
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
#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:
Here is the caller graph for this function:

◆ Initialise()

void G4DNAOneStepThermalizationModel::Initialise ( const G4ParticleDefinition particleDefinition,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 82 of file G4DNAOneStepThermalizationModel.cc.

84 {
85 #ifdef G4VERBOSE
86  if(fVerboseLevel)
87  G4cout << "Calling G4DNAOneStepThermalizationModel::Initialise()" << G4endl;
88 #endif
89  if (particleDefinition != G4Electron::ElectronDefinition())
90  {
91  G4ExceptionDescription exceptionDescription;
92  exceptionDescription << "G4DNAOneStepThermalizationModel can only be applied "
93  "to electrons";
94  G4Exception("G4DNAOneStepThermalizationModel::CrossSectionPerVolume",
95  "G4DNAOneStepThermalizationModel001",
96  FatalErrorInArgument,exceptionDescription);
97  return;
98  }
99 
100  if(!fIsInitialised)
101  {
102  fIsInitialised = true;
104  }
105 
108  GetNavigatorForTracking();
109 
110  fNavigator = new G4ITNavigator();
111 
112  fNavigator->SetWorldVolume(navigator->GetWorldVolume());
113  fNavigator->NewNavigatorState();
114 
117  GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
118 }
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:604
G4GLOB_DLL std::ostream G4cout
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4TransportationManager * GetTransportationManager()
static G4DNAMolecularMaterial * Instance()
const std::vector< G4double > * fpWaterDensity
#define G4endl
Definition: G4ios.hh:61
G4VPhysicalVolume * GetWorldVolume() const
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:134
G4ParticleChangeForGamma * fParticleChangeForGamma
Here is the call graph for this function:
Here is the caller graph for this function:

◆ operator=()

G4DNAOneStepThermalizationModel& G4DNAOneStepThermalizationModel::operator= ( const G4DNAOneStepThermalizationModel right)
private

◆ RadialDistributionOfProducts()

G4ThreeVector G4DNAOneStepThermalizationModel::RadialDistributionOfProducts ( G4double  Rrms) const
protected

Definition at line 155 of file G4DNAOneStepThermalizationModel.cc.

156 {
157  G4double sigma = std::sqrt(1.57) / 2 * expectationValue;
158 
159  G4double XValueForfMax = std::sqrt(2. * sigma * sigma);
160  G4double fMaxValue = std::sqrt(2. / 3.14)
161  * 1. / (sigma * sigma * sigma)
162  * (XValueForfMax * XValueForfMax)
163  * std::exp(-1. / 2. * (XValueForfMax * XValueForfMax)
164  / (sigma * sigma));
165 
166  G4double R;
167 
168  do
169  {
170  G4double aRandomfValue = fMaxValue * G4UniformRand();
171 
172  G4double sign;
173  if(G4UniformRand() > 0.5)
174  {
175  sign = +1.;
176  }
177  else
178  {
179  sign = -1;
180  }
181 
182  R = expectationValue + sign*3.*sigma* G4UniformRand();
183  G4double f = std::sqrt(2./3.14) * 1/std::pow(sigma, 3)
184  * R*R * std::exp(-1./2. * R*R/(sigma*sigma));
185 
186  if(aRandomfValue < f)
187  {
188  break;
189  }
190  }
191  while(1);
192 
193  G4double costheta = (2. * G4UniformRand()-1.);
194  G4double theta = std::acos(costheta);
195  G4double phi = 2. * pi * G4UniformRand();
196 
197  G4double xDirection = R * std::cos(phi) * std::sin(theta);
198  G4double yDirection = R * std::sin(theta) * std::sin(phi);
199  G4double zDirection = R * costheta;
200  G4ThreeVector RandDirection = G4ThreeVector(xDirection,
201  yDirection,
202  zDirection);
203 
204  return RandDirection;
205 }
CLHEP::Hep3Vector G4ThreeVector
#define G4UniformRand()
Definition: Randomize.hh:97
static const double pi
Definition: G4SIunits.hh:74
double G4double
Definition: G4Types.hh:76
G4int sign(const T t)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 210 of file G4DNAOneStepThermalizationModel.cc.

215 {
216 #ifdef G4VERBOSE
217  if(fVerboseLevel)
218  G4cout << "Calling SampleSecondaries() of G4DNAOneStepThermalizationModel"
219  << G4endl;
220 #endif
221  G4double k = particle->GetKineticEnergy();
222 
223  if (k <= HighEnergyLimit())
224  {
225  G4double k_eV = k/eV;
226 
227  G4double r_mean =
228  (-0.003*std::pow(k_eV,6)
229  + 0.0749*std::pow(k_eV,5)
230  - 0.7197*std::pow(k_eV,4)
231  + 3.1384*std::pow(k_eV,3)
232  - 5.6926*std::pow(k_eV,2)
233  + 5.6237*k_eV
234  - 0.7883)*nanometer;
235 
236  G4ThreeVector displacement = RadialDistributionOfProducts (r_mean);
237 
238  //______________________________________________________________
239  const G4Track * theIncomingTrack =
240  fParticleChangeForGamma->GetCurrentTrack();
241  G4ThreeVector finalPosition(theIncomingTrack->GetPosition()+displacement);
242 
243  fNavigator->SetWorldVolume(theIncomingTrack->GetTouchable()->
244  GetVolume(theIncomingTrack->GetTouchable()->
245  GetHistoryDepth()));
246 
247  double displacementMag = displacement.mag();
248  double safety = DBL_MAX;
249  G4ThreeVector direction = displacement/displacementMag;
250 
251  fNavigator->ResetHierarchyAndLocate(theIncomingTrack->GetPosition(),
252  direction,
253  *((G4TouchableHistory*)
254  theIncomingTrack->GetTouchable()));
255 
256  fNavigator->ComputeStep(theIncomingTrack->GetPosition(),
257  displacement/displacementMag,
258  displacementMag,
259  safety);
260 
261  if(safety <= displacementMag)
262  {
263  finalPosition = theIncomingTrack->GetPosition()
264  + (displacement/displacementMag)*safety*0.80;
265  }
266 
268  &finalPosition);
269 
270  fParticleChangeForGamma->SetProposedKineticEnergy(25.e-3*eV);
271  fParticleChangeForGamma->ProposeTrackStatus(fStopAndKill);
272  fParticleChangeForGamma->ProposeLocalEnergyDeposit(k);
273  }
274 }
G4ThreeVector RadialDistributionOfProducts(G4double Rrms) const
void CreateSolvatedElectron(const G4Track *, G4ThreeVector *finalPosition=0)
static const double nanometer
Definition: G4SIunits.hh:100
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
double mag() const
static G4DNAChemistryManager * Instance()
static const double eV
Definition: G4SIunits.hh:212
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
G4ParticleChangeForGamma * fParticleChangeForGamma
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetVerbose()

void G4DNAOneStepThermalizationModel::SetVerbose ( int  flag)
inline

Definition at line 105 of file G4DNAOneStepThermalizationModel.hh.

Member Data Documentation

◆ fIsInitialised

G4bool G4DNAOneStepThermalizationModel::fIsInitialised
protected

Definition at line 95 of file G4DNAOneStepThermalizationModel.hh.

◆ fNavigator

G4ITNavigator* G4DNAOneStepThermalizationModel::fNavigator
protected

Definition at line 97 of file G4DNAOneStepThermalizationModel.hh.

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNAOneStepThermalizationModel::fParticleChangeForGamma
protected

Definition at line 93 of file G4DNAOneStepThermalizationModel.hh.

◆ fpWaterDensity

const std::vector<G4double>* G4DNAOneStepThermalizationModel::fpWaterDensity
protected

Definition at line 90 of file G4DNAOneStepThermalizationModel.hh.

◆ fVerboseLevel

G4int G4DNAOneStepThermalizationModel::fVerboseLevel
protected

Definition at line 96 of file G4DNAOneStepThermalizationModel.hh.


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