Geant4
9.6.p02
Main Page
Related Pages
Modules
Namespaces
Classes
Files
File List
File Members
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Groups
Pages
geant4_9_6_p02
source
processes
electromagnetic
dna
models
src
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
//
28
#include "
G4DNATransformElectronModel.hh
"
29
#include "
G4SystemOfUnits.hh
"
30
#include "
G4ParticleChangeForGamma.hh
"
31
#include "
G4Electron.hh
"
32
#include "
G4NistManager.hh
"
33
#include "
G4DNAChemistryManager.hh
"
34
#include "
G4DNAMolecularMaterial.hh
"
35
36
G4DNATransformElectronModel::G4DNATransformElectronModel
(
const
G4ParticleDefinition
*,
37
const
G4String
& nam):
38
G4VEmModel
(nam),fIsInitialised(false)
39
{
40
fVerboseLevel = 0 ;
41
SetLowEnergyLimit
(0.*
eV
);
42
SetHighEnergyLimit
(0.025*
eV
);
43
fParticleChangeForGamma
= 0;
44
// fNistWater = G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER");
45
fpWaterDensity = 0;
46
fpWaterDensity = 0;
47
fEpsilon = 0.0001*
eV
;
48
}
49
50
//______________________________________________________________________
51
G4DNATransformElectronModel::~G4DNATransformElectronModel
()
52
{}
53
54
//______________________________________________________________________
55
void
G4DNATransformElectronModel::Initialise
(
const
G4ParticleDefinition
* particleDefinition,
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
73
fpWaterDensity =
G4DNAMolecularMaterial::Instance
()->
GetNumMolPerVolTableFor
(
G4Material::GetMaterial
(
"G4_WATER"
));
74
75
if
(!fIsInitialised)
76
{
77
fIsInitialised =
true
;
78
fParticleChangeForGamma
=
GetParticleChangeForGamma
();
79
}
80
}
81
82
//______________________________________________________________________
83
G4double
G4DNATransformElectronModel::CrossSectionPerVolume
(
const
G4Material
*
material
,
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
// {
129
const
G4Track
* track =
fParticleChangeForGamma
->
GetCurrentTrack
();
130
G4DNAChemistryManager::Instance
()->
CreateSolvatedElectron
(track);
131
fParticleChangeForGamma
->
ProposeTrackStatus
(
fStopAndKill
);
132
fParticleChangeForGamma
->
ProposeLocalEnergyDeposit
(k);
133
// }
134
}
Generated on Sat May 25 2013 14:33:31 for Geant4 by
1.8.4