Geant4  10.02
G4DNAOneStepThermalizationModel.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: G4DNAOneStepThermalizationModel.cc 94218 2015-11-09 08:24: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 #include "G4ITNavigator.hh"
49 #include "G4Navigator.hh"
51 #include "G4ITNavigator.hh"
52 
55  const G4String& nam) :
56  G4VEmModel(nam), fIsInitialised(false)
57 {
58  fVerboseLevel = 0;
61  SetHighEnergyLimit(exStructure.ExcitationEnergy(0));
63  fpWaterDensity = 0;
64  fNavigator = 0;
65 }
66 
67 //------------------------------------------------------------------------------
68 
70 {
71  if(fNavigator)
72  {
73  if(fNavigator->GetNavigatorState())
74  delete fNavigator->GetNavigatorState();
75  delete fNavigator;
76  }
77 }
78 
79 //------------------------------------------------------------------------------
80 
82 Initialise(const G4ParticleDefinition* particleDefinition,
83  const G4DataVector&)
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 
106  G4Navigator* navigator =
108  GetNavigatorForTracking();
109 
110  fNavigator = new G4ITNavigator();
111 
112  fNavigator->SetWorldVolume(navigator->GetWorldVolume());
113  fNavigator->NewNavigatorState();
114 
117  GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
118 }
119 
120 //------------------------------------------------------------------------------
121 
124  const G4ParticleDefinition*,
125  G4double ekin,
126  G4double,
127  G4double)
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 }
151 
152 //------------------------------------------------------------------------------
153 
155 RadialDistributionOfProducts(G4double expectationValue) const
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 }
206 
207 //------------------------------------------------------------------------------
208 
210 SampleSecondaries(std::vector<G4DynamicParticle*>*,
211  const G4MaterialCutsCouple*,
212  const G4DynamicParticle* particle,
213  G4double,
214  G4double)
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 =
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 
273  }
274 }
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:604
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:100
const G4ThreeVector & GetPosition() const
G4DNAOneStepThermalizationModel(const G4ParticleDefinition *p=0, const G4String &nam="DNASancheSolvatationModel")
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
#define G4UniformRand()
Definition: Randomize.hh:97
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 G4DNAChemistryManager * Instance()
static const double eV
Definition: G4SIunits.hh:212
static G4DNAMolecularMaterial * Instance()
const std::vector< G4double > * fpWaterDensity
static const double pi
Definition: G4SIunits.hh:74
const G4VTouchable * GetTouchable() const
const G4Track * GetCurrentTrack() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
#define DBL_MAX
Definition: templates.hh:83
G4VPhysicalVolume * GetWorldVolume() const
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
G4ParticleChangeForGamma * fParticleChangeForGamma
G4ThreeVector RadialDistributionOfProducts(G4double Rrms) const