Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 // $Id: G4DNAOneStepThermalizationModel.cc 96841 2016-05-12 13:40:37Z matkara $
27 //
28 // Author: Mathieu Karamitros
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 // 13 Nov 2016 M.Karamitros created
36 //
37 // -------------------------------------------------------------------
38 
39 #include "G4PhysicalConstants.hh"
40 #include "G4SystemOfUnits.hh"
43 #include "G4NistManager.hh"
44 #include "G4DNAChemistryManager.hh"
47 #include "G4ITNavigator.hh"
48 #include "G4Navigator.hh"
49 
50 //#define MODEL_VERBOSE
51 
52 //------------------------------------------------------------------------------
53 
54 template<typename MODEL>
57  const G4String& nam) :
58 G4VEmModel(nam), fIsInitialised(false)
59 {
60  fVerboseLevel = 0;
63  SetHighEnergyLimit(exStructure.ExcitationEnergy(0));
65  fpWaterDensity = 0;
66  fNavigator = 0;
67 }
68 
69 //------------------------------------------------------------------------------
70 
71 template<typename MODEL>
73 {
74  if(fNavigator)
75  {
76  // if(fNavigator->GetNavigatorState())
77  // delete fNavigator->GetNavigatorState();
78  delete fNavigator;
79  }
80 }
81 
82 //------------------------------------------------------------------------------
83 template<typename MODEL>
85 Initialise(const G4ParticleDefinition* particleDefinition,
86  const G4DataVector&)
87 {
88 #ifdef MODEL_VERBOSE
89  if(fVerboseLevel)
90  G4cout << "Calling G4DNAOneStepThermalizationModel::Initialise()"
91  << G4endl;
92 #endif
93  if (particleDefinition->GetParticleName() != "e-")
94  {
96  errMsg << "G4DNAOneStepThermalizationModel can only be applied "
97  "to electrons";
98  G4Exception("G4DNAOneStepThermalizationModel::CrossSectionPerVolume",
99  "G4DNAOneStepThermalizationModel001",
100  FatalErrorInArgument,errMsg);
101  return;
102  }
103 
104  if(!fIsInitialised)
105  {
106  fIsInitialised = true;
107  fParticleChangeForGamma = GetParticleChangeForGamma();
108  }
109 
110  G4Navigator* navigator =
112  GetNavigatorForTracking();
113 
114  fNavigator = new G4Navigator();
115 
116  if(navigator){ // add these checks for testing mode
117  auto world=navigator->GetWorldVolume();
118  if(world){
119  fNavigator->SetWorldVolume(world);
120  //fNavigator->NewNavigatorState();
121  }
122  }
123 
124  fpWaterDensity =
126  GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
127 }
128 
129 //------------------------------------------------------------------------------
130 template<typename MODEL>
133  const G4ParticleDefinition*,
134  G4double ekin,
135  G4double,
136  G4double)
137 {
138 #ifdef MODEL_VERBOSE
139  if(fVerboseLevel > 1)
140  G4cout << "Calling CrossSectionPerVolume() of G4DNAOneStepThermalizationModel"
141  << G4endl;
142 #endif
143 
144  if(ekin > HighEnergyLimit()){
145  return 0.0;
146  }
147 
148  G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
149 
150  if(waterDensity!= 0.0){
151  return DBL_MAX;
152  }
153  return 0.;
154 }
155 
156 //------------------------------------------------------------------------------
157 template<typename MODEL>
159  return MODEL::GetRmean(k);
160 }
161 
162 
163 //------------------------------------------------------------------------------
164 
165 template<typename MODEL>
168 {
169  return MODEL::GetPenetration(k, displacement);
170 }
171 
172 //------------------------------------------------------------------------------
173 template<typename MODEL>
175 SampleSecondaries(std::vector<G4DynamicParticle*>*,
176  const G4MaterialCutsCouple*,
177  const G4DynamicParticle* particle,
178  G4double,
179  G4double)
180 {
181 #ifdef MODEL_VERBOSE
182  if(fVerboseLevel)
183  G4cout << "Calling SampleSecondaries() of G4DNAOneStepThermalizationModel"
184  << G4endl;
185 #endif
186 
187  G4double k = particle->GetKineticEnergy();
188 
189  if (k <= HighEnergyLimit())
190  {
191  fParticleChangeForGamma->ProposeTrackStatus(fStopAndKill);
192  fParticleChangeForGamma->ProposeLocalEnergyDeposit(k);
193 
195  {
196  G4ThreeVector displacement(0,0,0);
197  GetPenetration(k, displacement);
198 
199  //______________________________________________________________
200  const G4Track * theIncomingTrack =
201  fParticleChangeForGamma->GetCurrentTrack();
202  G4ThreeVector finalPosition(theIncomingTrack->GetPosition()+displacement);
203 
204  fNavigator->SetWorldVolume(theIncomingTrack->GetTouchable()->
205  GetVolume(theIncomingTrack->GetTouchable()->
206  GetHistoryDepth()));
207 
208  double displacementMag = displacement.mag();
209  double safety = DBL_MAX;
210  G4ThreeVector direction = displacement/displacementMag;
211 
212  //--
213  // 6/09/16 - recupere de molecular dissocation
214  double mag_displacement = displacement.mag();
215  G4ThreeVector displacement_direction = displacement/mag_displacement;
216 
217  // double step = DBL_MAX;
218  // step = fNavigator->CheckNextStep(theIncomingTrack->GetPosition(),
219  // displacement_direction,
220  // mag_displacement,
221  // safety);
222  //
223  //
224  // if(safety < mag_displacement)
225  // {
227  // finalPosition = theIncomingTrack->GetPosition()
228  // + (displacement/displacementMag)*safety*0.80;
229  // }
230  //--
231 
232  fNavigator->ResetHierarchyAndLocate(theIncomingTrack->GetPosition(),
233  direction,
234  *((G4TouchableHistory*)
235  theIncomingTrack->GetTouchable()));
236 
237  fNavigator->ComputeStep(theIncomingTrack->GetPosition(),
238  displacement/displacementMag,
239  displacementMag,
240  safety);
241 
242  if(safety <= displacementMag)
243  {
244  finalPosition = theIncomingTrack->GetPosition()
245  + (displacement/displacementMag)*safety*0.80;
246  }
247 
249  &finalPosition);
250 
251  fParticleChangeForGamma->SetProposedKineticEnergy(25.e-3*eV);
252  }
253  }
254 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4double GetKineticEnergy() const
const std::vector< G4double > * fpWaterDensity
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:602
void CreateSolvatedElectron(const G4Track *, G4ThreeVector *finalPosition=0)
size_t GetIndex() const
Definition: G4Material.hh:262
const G4ThreeVector & GetPosition() const
G4TDNAOneStepThermalizationModel(const G4ParticleDefinition *p=0, const G4String &nam="DNAOneStepThermalizationModel")
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:724
G4GLOB_DLL std::ostream G4cout
void GetPenetration(G4double energy, G4ThreeVector &displacement)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
static constexpr double eV
Definition: G4SIunits.hh:215
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 G4TransportationManager * GetTransportationManager()
static G4DNAChemistryManager * Instance()
static G4DNAMolecularMaterial * Instance()
const G4VTouchable * GetTouchable() const
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:731
double mag() const
#define DBL_MAX
Definition: templates.hh:83
G4VPhysicalVolume * GetWorldVolume() const