Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Member Functions | Protected Attributes
G4DNAOneStepSolvatationModel Class Reference

#include <G4DNAOneStepSolvatationModel.hh>

Inheritance diagram for G4DNAOneStepSolvatationModel:
G4VEmModel

Public Member Functions

 G4DNAOneStepSolvatationModel (const G4ParticleDefinition *p=0, const G4String &nam="DNASancheSolvatationModel")
 
virtual ~G4DNAOneStepSolvatationModel ()
 
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 ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., 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 *)
 
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
 

Protected Member Functions

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

Protected Attributes

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

Detailed Description

When an electron reaches the highest energy domain of G4DNAOneStepSolvatationModel, 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 54 of file G4DNAOneStepSolvatationModel.hh.

Constructor & Destructor Documentation

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

Definition at line 49 of file G4DNAOneStepSolvatationModel.cc.

References G4DNAWaterExcitationStructure::ExcitationEnergy(), fParticleChangeForGamma, fpWaterDensity, fVerboseLevel, G4VEmModel::SetHighEnergyLimit(), and G4VEmModel::SetLowEnergyLimit().

50  :
51  G4VEmModel(nam),fIsInitialised(false)
52 {
53  fVerboseLevel = 0 ;
55  G4DNAWaterExcitationStructure exStructure ;
56  SetHighEnergyLimit(exStructure.ExcitationEnergy(0));
58  // fNistWater = G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER");
59  fpWaterDensity = 0;
60 }
const std::vector< G4double > * fpWaterDensity
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:65
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
G4ParticleChangeForGamma * fParticleChangeForGamma
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:690
G4DNAOneStepSolvatationModel::~G4DNAOneStepSolvatationModel ( )
virtual

Definition at line 63 of file G4DNAOneStepSolvatationModel.cc.

64 {}

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 94 of file G4DNAOneStepSolvatationModel.cc.

References DBL_MAX, fVerboseLevel, G4cout, G4endl, G4Material::GetIndex(), and G4VEmModel::HighEnergyLimit().

99 {
100 #ifdef G4VERBOSE
101  if (fVerboseLevel > 1)
102  G4cout << "Calling CrossSectionPerVolume() of G4SancheSolvatationModel" << G4endl;
103 #endif
104 
105  if(ekin > HighEnergyLimit())
106  {
107  return 0.0;
108  }
109 
110  G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
111 
112  if(waterDensity!= 0.0)
113  // if (material == nistwater || material->GetBaseMaterial() == nistwater)
114  {
115  if (ekin <= HighEnergyLimit())
116  {
117  return DBL_MAX;
118  }
119  }
120  return 0. ;
121 }
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:592
size_t GetIndex() const
Definition: G4Material.hh:260
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
void G4DNAOneStepSolvatationModel::Initialise ( const G4ParticleDefinition particleDefinition,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 67 of file G4DNAOneStepSolvatationModel.cc.

References G4Electron::ElectronDefinition(), FatalErrorInArgument, fIsInitialised, fParticleChangeForGamma, fpWaterDensity, fVerboseLevel, G4cout, G4endl, G4Exception(), G4Material::GetMaterial(), G4DNAMolecularMaterial::GetNumMolPerVolTableFor(), G4VEmModel::GetParticleChangeForGamma(), and G4DNAMolecularMaterial::Instance().

69 {
70 #ifdef G4VERBOSE
71  if (fVerboseLevel)
72  G4cout << "Calling G4SancheSolvatationModel::Initialise()" << G4endl;
73 #endif
74  if (particleDefinition != G4Electron::ElectronDefinition())
75  {
76  G4ExceptionDescription exceptionDescription ;
77  exceptionDescription << "Attempting to calculate cross section for wrong particle";
78  G4Exception("G4DNAOneStepSolvatationModel::CrossSectionPerVolume","G4DNAOneStepSolvatationModel001",
79  FatalErrorInArgument,exceptionDescription);
80  return;
81  }
82 
83  if(!fIsInitialised)
84  {
85  fIsInitialised = true;
87  }
88 
89  // Initialize water density pointer
91 }
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:578
const std::vector< G4double > * fpWaterDensity
G4GLOB_DLL std::ostream G4cout
const std::vector< double > * GetNumMolPerVolTableFor(const G4Material *) const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4DNAMolecularMaterial * Instance()
G4ParticleChangeForGamma * fParticleChangeForGamma
#define G4endl
Definition: G4ios.hh:61
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:121
G4ThreeVector G4DNAOneStepSolvatationModel::RadialDistributionOfProducts ( G4double  Rrms) const
protected

Definition at line 124 of file G4DNAOneStepSolvatationModel.cc.

References G4UniformRand, python.hepunit::pi, and G4INCL::Math::sign().

Referenced by SampleSecondaries().

125 {
126  G4double sigma = std::sqrt(1.57)/2*expectationValue;
127 
128  G4double XValueForfMax = std::sqrt(2.*sigma*sigma);
129  G4double fMaxValue = std::sqrt(2./3.14) * 1./(sigma*sigma*sigma) *
130  (XValueForfMax*XValueForfMax)*
131  std::exp(-1./2. * (XValueForfMax*XValueForfMax)/(sigma*sigma));
132 
133  G4double R;
134 
135  do
136  {
137  G4double aRandomfValue = fMaxValue*G4UniformRand();
138 
139  G4double sign;
140  if(G4UniformRand() > 0.5)
141  {
142  sign = +1.;
143  }
144  else
145  {
146  sign = -1;
147  }
148 
149  R = expectationValue + sign*3.*sigma* G4UniformRand();
150  G4double f = std::sqrt(2./3.14) * 1/std::pow(sigma, 3) * R*R * std::exp(-1./2. * R*R/(sigma*sigma));
151 
152  if(aRandomfValue < f)
153  {
154  break;
155  }
156  }
157  while(1);
158 
159  G4double costheta = (2.*G4UniformRand()-1.);
160  G4double theta = std::acos (costheta);
161  G4double phi = 2.*pi*G4UniformRand();
162 
163  G4double xDirection = R*std::cos(phi)* std::sin(theta);
164  G4double yDirection = R*std::sin(theta)*std::sin(phi);
165  G4double zDirection = R*costheta;
166  G4ThreeVector RandDirection = G4ThreeVector(xDirection, yDirection, zDirection);
167 
168  return RandDirection;
169 }
CLHEP::Hep3Vector G4ThreeVector
#define G4UniformRand()
Definition: Randomize.hh:87
double G4double
Definition: G4Types.hh:76
G4int sign(const T t)
void G4DNAOneStepSolvatationModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  ,
const G4MaterialCutsCouple ,
const G4DynamicParticle particle,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 172 of file G4DNAOneStepSolvatationModel.cc.

References G4DNAChemistryManager::CreateSolvatedElectron(), python.hepunit::eV, fParticleChangeForGamma, fStopAndKill, fVerboseLevel, G4cout, G4endl, G4ParticleChangeForGamma::GetCurrentTrack(), G4DynamicParticle::GetKineticEnergy(), G4Track::GetPosition(), G4VEmModel::HighEnergyLimit(), G4DNAChemistryManager::Instance(), python.hepunit::nanometer, G4VParticleChange::ProposeLocalEnergyDeposit(), G4VParticleChange::ProposeTrackStatus(), RadialDistributionOfProducts(), and G4ParticleChangeForGamma::SetProposedKineticEnergy().

177 {
178 #ifdef G4VERBOSE
179  if (fVerboseLevel)
180  G4cout << "Calling SampleSecondaries() of G4SancheSolvatationModel" << G4endl;
181 #endif
182  G4double k = particle->GetKineticEnergy();
183 
184  if (k <= HighEnergyLimit())
185  {
186  G4double r_mean =
187  (-0.003*std::pow(k/eV,6) + 0.0749*std::pow(k/eV,5) - 0.7197*std::pow(k/eV,4)
188  + 3.1384*std::pow(k/eV,3) - 5.6926*std::pow(k/eV,2) + 5.6237*k/eV - 0.7883)*nanometer;
189 
190  G4ThreeVector displacement = RadialDistributionOfProducts (r_mean);
191 
192  //______________________________________________________________
193  const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
194  G4ThreeVector finalPosition(theIncomingTrack->GetPosition()+displacement);
195 
196  G4DNAChemistryManager::Instance()->CreateSolvatedElectron(theIncomingTrack,&finalPosition);
197 
201  }
202 }
G4double GetKineticEnergy() const
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:592
void CreateSolvatedElectron(const G4Track *, G4ThreeVector *finalPosition=0)
const G4ThreeVector & GetPosition() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
int nanometer
Definition: hepunit.py:35
G4GLOB_DLL std::ostream G4cout
static G4DNAChemistryManager * Instance()
G4ThreeVector RadialDistributionOfProducts(G4double Rrms) const
G4ParticleChangeForGamma * fParticleChangeForGamma
const G4Track * GetCurrentTrack() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
void G4DNAOneStepSolvatationModel::SetVerbose ( int  flag)
inline

Definition at line 92 of file G4DNAOneStepSolvatationModel.hh.

References fVerboseLevel.

Field Documentation

G4bool G4DNAOneStepSolvatationModel::fIsInitialised
protected

Definition at line 84 of file G4DNAOneStepSolvatationModel.hh.

Referenced by Initialise().

G4ParticleChangeForGamma* G4DNAOneStepSolvatationModel::fParticleChangeForGamma
protected
const std::vector<G4double>* G4DNAOneStepSolvatationModel::fpWaterDensity
protected

Definition at line 79 of file G4DNAOneStepSolvatationModel.hh.

Referenced by G4DNAOneStepSolvatationModel(), and Initialise().

G4int G4DNAOneStepSolvatationModel::fVerboseLevel
protected

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