Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4DNATransformElectronModel.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: G4DNATransformElectronModel.cc 64057 2012-10-30 15:04:49Z gcosmo $
27 //
29 #include "G4SystemOfUnits.hh"
31 #include "G4Electron.hh"
32 #include "G4NistManager.hh"
33 #include "G4DNAChemistryManager.hh"
35 
37  const G4String& nam):
38  G4VEmModel(nam),fIsInitialised(false)
39 {
40  fVerboseLevel = 0 ;
42  SetHighEnergyLimit(0.025*eV);
44  // fNistWater = G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER");
45  fpWaterDensity = 0;
46  fpWaterDensity = 0;
47  fEpsilon = 0.0001*eV;
48 }
49 
50 //______________________________________________________________________
52 {}
53 
54 //______________________________________________________________________
56  const G4DataVector&)
57 {
58 #ifdef G4VERBOSE
59  if (fVerboseLevel)
60  G4cout << "Calling G4DNATransformElectronModel::Initialise()" << G4endl;
61 #endif
62 
63  if (particleDefinition != G4Electron::ElectronDefinition())
64  {
65  G4ExceptionDescription exceptionDescription ;
66  exceptionDescription << "Attempting to calculate cross section for wrong particle";
67  G4Exception("G4DNATransformElectronModel::CrossSectionPerVolume","G4DNATransformElectronModel001",
68  FatalErrorInArgument,exceptionDescription);
69  return;
70  }
71 
72  // Initialize water density pointer
74 
75  if(!fIsInitialised)
76  {
77  fIsInitialised = true;
79  }
80 }
81 
82 //______________________________________________________________________
84  const G4ParticleDefinition*,
85  G4double ekin,
86  G4double,
87  G4double)
88 {
89 #if G4VERBOSE
90  if (fVerboseLevel > 1)
91  G4cout << "Calling CrossSectionPerVolume() of G4DNATransformElectronModel" << G4endl;
92 #endif
93 
94  if(ekin - fEpsilon > HighEnergyLimit())
95  {
96  return 0.0;
97  }
98 
99  G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
100 
101  if(waterDensity!= 0.0)
102  // if (material == nistwater || material->GetBaseMaterial() == nistwater)
103  {
104  if (ekin - fEpsilon <= HighEnergyLimit())
105  {
106  return DBL_MAX;
107  }
108  }
109 
110  return 0.0 ;
111 }
112 
113 //______________________________________________________________________
114 void G4DNATransformElectronModel::SampleSecondaries(std::vector<G4DynamicParticle*>* /*fvect*/,
115  const G4MaterialCutsCouple*,
116  const G4DynamicParticle* particle,
117  G4double,
118  G4double)
119 {
120 #if G4VERBOSE
121  if (fVerboseLevel)
122  G4cout << "Calling SampleSecondaries() of G4DNATransformElectronModel" << G4endl;
123 #endif
124 
125  G4double k = particle->GetKineticEnergy();
126 
127 // if (k - fEpsilon <= HighEnergyLimit())
128 // {
133 // }
134 }
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4double GetKineticEnergy() const
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:592
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:578
void CreateSolvatedElectron(const G4Track *, G4ThreeVector *finalPosition=0)
size_t GetIndex() const
Definition: G4Material.hh:260
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
G4DNATransformElectronModel(const G4ParticleDefinition *p=0, const G4String &nam="DNATransformElectronModel")
string material
Definition: eplot.py:19
G4GLOB_DLL std::ostream G4cout
const std::vector< double > * GetNumMolPerVolTableFor(const G4Material *) const
G4ParticleChangeForGamma * fParticleChangeForGamma
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4DNAChemistryManager * Instance()
static G4DNAMolecularMaterial * Instance()
const G4Track * GetCurrentTrack() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:690
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
#define DBL_MAX
Definition: templates.hh:83
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:121
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)