Geant4  10.01.p01
G4DNAOneStepSolvatationModel.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // $Id: G4DNAOneStepSolvatationModel.cc 87137 2014-11-25 09:12:48Z gcosmo $
27 //
28 // Author: Mathieu Karamitros (kara (AT) cenbg . in2p3 . fr)
29 //
30 // WARNING : This class is released as a prototype.
31 // It might strongly evolve or even disapear in the next releases.
32 //
33 // History:
34 // -----------
35 // 10 Oct 2011 M.Karamitros created
36 //
37 // -------------------------------------------------------------------
38 
40 #include "G4PhysicalConstants.hh"
41 #include "G4SystemOfUnits.hh"
44 #include "G4Electron.hh"
45 #include "G4NistManager.hh"
46 #include "G4DNAChemistryManager.hh"
48 
50  const G4String& nam) :
51  G4VEmModel(nam), fIsInitialised(false)
52 {
53  fVerboseLevel = 0;
56  SetHighEnergyLimit(exStructure.ExcitationEnergy(0));
58  // fNistWater = G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER");
59  fpWaterDensity = 0;
60 }
61 
62 //______________________________________________________________________
64 {
65 }
66 
67 //______________________________________________________________________
69  const G4DataVector&)
70 {
71 #ifdef G4VERBOSE
72  if (fVerboseLevel)
73  G4cout << "Calling G4SancheSolvatationModel::Initialise()" << G4endl;
74 #endif
75  if (particleDefinition != G4Electron::ElectronDefinition())
76  {
77  G4ExceptionDescription exceptionDescription;
78  exceptionDescription << "Attempting to calculate cross section for wrong particle";
79  G4Exception("G4DNAOneStepSolvatationModel::CrossSectionPerVolume","G4DNAOneStepSolvatationModel001",
80  FatalErrorInArgument,exceptionDescription);
81  return;
82  }
83 
84  if(!fIsInitialised)
85  {
86  fIsInitialised = true;
88  }
89 
90  // Initialize water density pointer
92 }
93 
94 //______________________________________________________________________
96  const G4ParticleDefinition*,
97  G4double ekin,
98  G4double,
99  G4double)
100 {
101 #ifdef G4VERBOSE
102  if (fVerboseLevel > 1)
103  G4cout << "Calling CrossSectionPerVolume() of G4SancheSolvatationModel"
104  << G4endl;
105 #endif
106 
107  if(ekin > HighEnergyLimit())
108  {
109  return 0.0;
110  }
111 
112  G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
113 
114  if(waterDensity!= 0.0)
115  // if (material == nistwater || material->GetBaseMaterial() == nistwater)
116  {
117  if (ekin <= HighEnergyLimit())
118  {
119  return DBL_MAX;
120  }
121  }
122  return 0.;
123 }
124 
125 //______________________________________________________________________
127 {
128  G4double sigma = std::sqrt(1.57) / 2 * expectationValue;
129 
130  G4double XValueForfMax = std::sqrt(2. * sigma * sigma);
131  G4double fMaxValue = std::sqrt(2. / 3.14) * 1. / (sigma * sigma * sigma)
132  * (XValueForfMax * XValueForfMax)
133  * std::exp(-1. / 2. * (XValueForfMax * XValueForfMax) / (sigma * sigma));
134 
135  G4double R;
136 
137  do
138  {
139  G4double aRandomfValue = fMaxValue * G4UniformRand();
140 
141  G4double sign;
142  if(G4UniformRand() > 0.5)
143  {
144  sign = +1.;
145  }
146  else
147  {
148  sign = -1;
149  }
150 
151  R = expectationValue + sign*3.*sigma* G4UniformRand();
152  G4double f = std::sqrt(2./3.14) * 1/std::pow(sigma, 3) * R*R * std::exp(-1./2. * R*R/(sigma*sigma));
153 
154  if(aRandomfValue < f)
155  {
156  break;
157  }
158  }
159  while(1);
160 
161  G4double costheta = (2. * G4UniformRand()-1.);
162  G4double theta = std::acos(costheta);
163  G4double phi = 2. * pi * G4UniformRand();
164 
165  G4double xDirection = R * std::cos(phi) * std::sin(theta);
166  G4double yDirection = R * std::sin(theta) * std::sin(phi);
167  G4double zDirection = R * costheta;
168  G4ThreeVector RandDirection = G4ThreeVector(xDirection, yDirection,
169  zDirection);
170 
171  return RandDirection;
172 }
173 
174 //______________________________________________________________________
177  const G4MaterialCutsCouple*,
178  const G4DynamicParticle* particle,
179  G4double,
180  G4double)
181 {
182 #ifdef G4VERBOSE
183  if (fVerboseLevel)
184  G4cout << "Calling SampleSecondaries() of G4SancheSolvatationModel" << G4endl;
185 #endif
186  G4double k = particle->GetKineticEnergy();
187 
188  if (k <= HighEnergyLimit())
189  {
190  G4double r_mean =
191  (-0.003*std::pow(k/eV,6) + 0.0749*std::pow(k/eV,5) - 0.7197*std::pow(k/eV,4)
192  + 3.1384*std::pow(k/eV,3) - 5.6926*std::pow(k/eV,2) + 5.6237*k/eV - 0.7883)*nanometer;
193 
194  G4ThreeVector displacement = RadialDistributionOfProducts (r_mean);
195 
196  //______________________________________________________________
197  const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
198  G4ThreeVector finalPosition(theIncomingTrack->GetPosition()+displacement);
199 
200  G4DNAChemistryManager::Instance()->CreateSolvatedElectron(theIncomingTrack,&finalPosition);
201 
205  }
206 }
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:615
G4DNAOneStepSolvatationModel(const G4ParticleDefinition *p=0, const G4String &nam="DNASancheSolvatationModel")
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:603
void CreateSolvatedElectron(const G4Track *, G4ThreeVector *finalPosition=0)
On the same idea as the previous method but for solvated electron.
size_t GetIndex() const
Definition: G4Material.hh:262
static const double nanometer
Definition: G4SIunits.hh:91
const G4double pi
const G4ThreeVector & GetPosition() const
const std::vector< G4double > * fpWaterDensity
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:706
#define G4UniformRand()
Definition: Randomize.hh:95
G4GLOB_DLL std::ostream G4cout
const std::vector< double > * GetNumMolPerVolTableFor(const G4Material *) const
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4DNAChemistryManager * Instance()
static const double eV
Definition: G4SIunits.hh:194
static G4DNAMolecularMaterial * Instance()
G4ThreeVector RadialDistributionOfProducts(G4double Rrms) const
G4ParticleChangeForGamma * fParticleChangeForGamma
const G4Track * GetCurrentTrack() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:713
#define DBL_MAX
Definition: templates.hh:83
G4int sign(const T t)
A simple sign function that allows us to port fortran code to c++ more easily.
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:134