Geant4-11
G4DNAOneStepThermalizationModel.hpp
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//
27// Author: Mathieu Karamitros
28//
29// WARNING : This class is released as a prototype.
30// It might strongly evolve or even disappear in the next releases.
31//
32// History:
33// -----------
34// 13 Nov 2016 M.Karamitros created
35//
36// -------------------------------------------------------------------
37
39#include "G4SystemOfUnits.hh"
42#include "G4NistManager.hh"
46#include "G4ITNavigator.hh"
47#include "G4Navigator.hh"
48
49//#define MODEL_VERBOSE
50
51//------------------------------------------------------------------------------
52
53template<typename MODEL>
56 const G4String& nam) :
57G4VEmModel(nam), fIsInitialised(false)
58{
59 fVerboseLevel = 0;
65}
66
67//------------------------------------------------------------------------------
68
69template<typename MODEL>
71{
72 // if(fpNavigator && fpNavigator->GetNavigatorState())
73 // delete fpNavigator->GetNavigatorState();
74}
75
76//------------------------------------------------------------------------------
77template<typename MODEL>
79Initialise(const G4ParticleDefinition* particleDefinition,
80 const G4DataVector&)
81{
82#ifdef MODEL_VERBOSE
83 if(fVerboseLevel)
84 G4cout << "Calling G4DNAOneStepThermalizationModel::Initialise()"
85 << G4endl;
86#endif
87 if (particleDefinition->GetParticleName() != "e-")
88 {
90 errMsg << "G4DNAOneStepThermalizationModel can only be applied "
91 "to electrons";
92 G4Exception("G4DNAOneStepThermalizationModel::CrossSectionPerVolume",
93 "G4DNAOneStepThermalizationModel001",
95 return;
96 }
97
98 if(!fIsInitialised)
99 {
100 fIsInitialised = true;
101 fpParticleChangeForGamma = GetParticleChangeForGamma();
102 }
103
106 GetNavigatorForTracking();
107
108 fpNavigator.reset(new G4Navigator());
109
110 if(navigator){ // add these checks for testing mode
111 auto world=navigator->GetWorldVolume();
112 if(world){
113 fpNavigator->SetWorldVolume(world);
114 //fNavigator->NewNavigatorState();
115 }
116 }
117
118 fpWaterDensity =
120 GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
121}
122
123//------------------------------------------------------------------------------
124template<typename MODEL>
128 G4double ekin,
129 G4double,
130 G4double)
131{
132#ifdef MODEL_VERBOSE
133 if(fVerboseLevel > 1)
134 G4cout << "Calling CrossSectionPerVolume() of G4DNAOneStepThermalizationModel"
135 << G4endl;
136#endif
137
138 if(ekin > HighEnergyLimit()){
139 return 0.0;
140 }
141
142 G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
143
144 if(waterDensity!= 0.0){
145 return DBL_MAX;
146 }
147 return 0.;
148}
149
150//------------------------------------------------------------------------------
151template<typename MODEL>
153 return MODEL::GetRmean(k);
154}
155
156
157//------------------------------------------------------------------------------
158
159template<typename MODEL>
161GetPenetration(G4double k, G4ThreeVector& displacement)
162{
163 return MODEL::GetPenetration(k, displacement);
164}
165
166//------------------------------------------------------------------------------
167template<typename MODEL>
169SampleSecondaries(std::vector<G4DynamicParticle*>*,
171 const G4DynamicParticle* particle,
172 G4double,
173 G4double)
174{
175#ifdef MODEL_VERBOSE
176 if(fVerboseLevel)
177 G4cout << "Calling SampleSecondaries() of G4DNAOneStepThermalizationModel"
178 << G4endl;
179#endif
180
181 G4double k = particle->GetKineticEnergy();
182
183 if (k <= HighEnergyLimit())
184 {
185 fpParticleChangeForGamma->ProposeTrackStatus(fStopAndKill);
186 fpParticleChangeForGamma->ProposeLocalEnergyDeposit(k);
187
189 {
190 G4ThreeVector displacement(0,0,0);
191 GetPenetration(k, displacement);
192
193 //______________________________________________________________
194 const G4Track * theIncomingTrack =
195 fpParticleChangeForGamma->GetCurrentTrack();
196 G4ThreeVector finalPosition(theIncomingTrack->GetPosition()+displacement);
197
198 fpNavigator->SetWorldVolume(theIncomingTrack->GetTouchable()->
199 GetVolume(theIncomingTrack->GetTouchable()->
200 GetHistoryDepth()));
201
202 double displacementMag = displacement.mag();
203 double safety = DBL_MAX;
204 G4ThreeVector direction = displacement/displacementMag;
205
206 //--
207 // 6/09/16 - recupere de molecular dissocation
208 double mag_displacement = displacement.mag();
209 G4ThreeVector displacement_direction = displacement/mag_displacement;
210
211 // double step = DBL_MAX;
212 // step = fNavigator->CheckNextStep(theIncomingTrack->GetPosition(),
213 // displacement_direction,
214 // mag_displacement,
215 // safety);
216 //
217 //
218 // if(safety < mag_displacement)
219 // {
221 // finalPosition = theIncomingTrack->GetPosition()
222 // + (displacement/displacementMag)*safety*0.80;
223 // }
224 //--
225
226 fpNavigator->ResetHierarchyAndLocate(theIncomingTrack->GetPosition(),
227 direction,
229 theIncomingTrack->GetTouchable()));
230
231 fpNavigator->ComputeStep(theIncomingTrack->GetPosition(),
232 displacement/displacementMag,
233 displacementMag,
234 safety);
235
236 if(safety <= displacementMag)
237 {
238 finalPosition = theIncomingTrack->GetPosition()
239 + (displacement/displacementMag)*safety*0.80;
240 }
241
243 &finalPosition);
244
245 fpParticleChangeForGamma->SetProposedKineticEnergy(25.e-3*eV);
246 }
247 }
248}
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
static constexpr double eV
Definition: G4SIunits.hh:201
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double mag() const
void CreateSolvatedElectron(const G4Track *, G4ThreeVector *pFinalPosition=nullptr)
static G4DNAChemistryManager * Instance()
static G4DNAMolecularMaterial * Instance()
G4double GetKineticEnergy() const
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:686
const G4String & GetParticleName() const
const std::vector< G4double > * fpWaterDensity
void GetPenetration(G4double energy, G4ThreeVector &displacement)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
G4TDNAOneStepThermalizationModel(const G4ParticleDefinition *p=0, const G4String &nam="DNAOneStepThermalizationModel")
const G4ThreeVector & GetPosition() const
const G4VTouchable * GetTouchable() const
static G4TransportationManager * GetTransportationManager()
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:767
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:774
string material
Definition: eplot.py:19
#define DBL_MAX
Definition: templates.hh:62