Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4DNATransformElectronModel.cc
이 파일의 문서화 페이지로 가기
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 96861 2016-05-13 13:43:04Z gcosmo $
27 //
29 #include "G4SystemOfUnits.hh"
31 #include "G4NistManager.hh"
32 #include "G4DNAChemistryManager.hh"
34 
35 //#define MODEL_VERBOSE
36 
39  const G4String& nam) :
40  G4VEmModel(nam),fIsInitialised(false)
41 {
42  fVerboseLevel = 0;
43  SetLowEnergyLimit(0. * eV);
44  SetHighEnergyLimit(0.025 * eV);
46  fpWaterDensity = 0;
47  fpWaterDensity = 0;
48  fEpsilon = 0.0001 * eV;
49 }
50 
51 //______________________________________________________________________
53 {}
54 
55 //______________________________________________________________________
57 Initialise(const G4ParticleDefinition* particleDefinition,
58  const G4DataVector&)
59 {
60 #ifdef MODEL_VERBOSE
61  if (fVerboseLevel)
62  G4cout << "Calling G4DNATransformElectronModel::Initialise()" << G4endl;
63 #endif
64 
65  if(particleDefinition->GetParticleName() != "e-")
66  {
68  errMsg << "Attempting to calculate cross section for wrong particle";
69  G4Exception("G4DNATransformElectronModel::CrossSectionPerVolume",
70  "G4DNATransformElectronModel001", FatalErrorInArgument, errMsg);
71  return;
72  }
73 
74  // Initialize water density pointer
76  GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
77 
78  if(!fIsInitialised)
79  {
80  fIsInitialised = true;
82  }
83 }
84 
85 //______________________________________________________________________
88  const G4ParticleDefinition*,
89  G4double ekin,
90  G4double,
91  G4double)
92 {
93 #if MODEL_VERBOSE
94  if (fVerboseLevel > 1)
95  G4cout << "Calling CrossSectionPerVolume() of G4DNATransformElectronModel" << G4endl;
96 #endif
97 
98  if(ekin - fEpsilon > HighEnergyLimit())
99  {
100  return 0.0;
101  }
102 
103  G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
104 
105  if(waterDensity != 0.0)
106  {
107  return DBL_MAX;
108  }
109 
110  return 0.0;
111 }
112 
113 //______________________________________________________________________
114 void
115 G4DNATransformElectronModel::SampleSecondaries(std::vector<G4DynamicParticle*>* /*fvect*/,
116  const G4MaterialCutsCouple*,
117  const G4DynamicParticle* particle,
118  G4double,
119  G4double)
120 {
121 #if MODEL_VERBOSE
122  if (fVerboseLevel)
123  G4cout << "Calling SampleSecondaries() of G4DNATransformElectronModel" << G4endl;
124 #endif
125 
126  G4double k = particle->GetKineticEnergy();
127 
128 // if (k - fEpsilon <= HighEnergyLimit()) // should be already checked
129 // {
134 // }
135 }
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:45
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:707
G4DNATransformElectronModel(const G4ParticleDefinition *p=0, const G4String &nam="DNATransformElectronModel")
size_t GetIndex() const
Definition: G4Material.hh:262
#define G4endl
Definition: G4ios.hh:61
const G4String & GetParticleName() const
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:608
void CreateSolvatedElectron(const G4Track *, G4ThreeVector *finalPosition=0)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
double G4double
Definition: G4Types.hh:76
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
static G4DNAMolecularMaterial * Instance()
static constexpr double eV
Definition: G4SIunits.hh:215
const G4Track * GetCurrentTrack() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.hh:65
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:714
G4ParticleChangeForGamma * fParticleChangeForGamma
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:131
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:609
#define DBL_MAX
Definition: templates.hh:83
const std::vector< G4double > * fpWaterDensity
void ProposeTrackStatus(G4TrackStatus status)
static G4DNAChemistryManager * Instance()