Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4DNAMeltonAttachmentModel.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//
27
28// Created by Z. Francis
29
31#include "G4SystemOfUnits.hh"
34
35//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
36
37using namespace std;
38
39//#define MELTON_VERBOSE // prevent checking conditions at run time
40
41//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
42
44 const G4String& nam) :
45G4VEmModel(nam), isInitialised(false)
46{
47 fpWaterDensity = 0;
48
49 SetLowEnergyLimit(4.*eV);
50 SetHighEnergyLimit(13.*eV);
51
52 verboseLevel = 0;
53 // Verbosity scale:
54 // 0 = nothing
55 // 1 = warning for energy non-conservation
56 // 2 = details of energy budget
57 // 3 = calculation of cross sections, file openings, sampling of atoms
58 // 4 = entering in methods
59
60#ifdef MELTON_VERBOSE
61 if (verboseLevel > 0)
62 {
63 G4cout << "Melton Attachment model is constructed "
64 << G4endl
65 << "Energy range: "
66 << LowEnergyLimit() / eV << " eV - "
67 << HighEnergyLimit() / eV << " eV"
68 << G4endl;
69 }
70#endif
71
73 fDissociationFlag = true;
74 fData = 0;
75
76 // Selection of stationary mode
77
78 statCode = false;
79}
80
81//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
82
84{
85 if(fData) delete fData;
86}
87
88//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
89
91 const G4DataVector& /*cuts*/)
92{
93#ifdef MELTON_VERBOSE
94 if (verboseLevel > 3)
95 G4cout
96 << "Calling G4DNAMeltonAttachmentModel::Initialise()" << G4endl;
97#endif
98
99 // Only electron
100
101 if(particle->GetParticleName() != "e-")
102 {
103 G4Exception("G4DNAMeltonAttachmentModel::Initialise",
104 "em0002",
106 "Model not applicable to particle type.");
107 }
108
109 // Energy limits
110
111 if (LowEnergyLimit() < 4.*eV)
112 {
114 errMsg << "G4DNAMeltonAttachmentModel: low energy limit increased from " <<
115 LowEnergyLimit()/eV << " eV to " << 4. << " eV" << G4endl;
116
117 G4Exception("G4DNAMeltonAttachmentModel::Initialise",
118 "Melton_LowerEBoundary",
120 errMsg);
121
122 SetLowEnergyLimit(4*eV);
123 }
124
125 if (HighEnergyLimit() > 13.*eV)
126 {
128 errMsg << "G4DNAMeltonAttachmentModel: high energy limit decreased from " <<
129 HighEnergyLimit()/eV << " eV to " << 13. << " eV" << G4endl;
130
131 G4Exception("G4DNAMeltonAttachmentModel::Initialise",
132 "Melton_HigherEBoundary",
134 errMsg);
135
136 SetHighEnergyLimit(13.*eV);
137 }
138
139 // Reading of data files
140
141 G4double scaleFactor = 1e-18*cm2;
142
143 // For total cross section
144 G4String fileElectron("dna/sigma_attachment_e_melton");
145
147 eV, scaleFactor);
148 fData->LoadData(fileElectron);
149
150
151#ifdef MELTON_VERBOSE
152 if( verboseLevel >0)
153 {
154 if (verboseLevel > 2)
155 {
156 G4cout << "Loaded cross section data for Melton Attachment model" << G4endl;
157 }
158
159 G4cout << "Melton Attachment model is initialized " << G4endl
160 << "Energy range: "
161 << LowEnergyLimit() / eV << " eV - "
162 << HighEnergyLimit() / eV << " eV"
163 << G4endl;
164 }
165#endif
166
167 // Initialize water density pointer
168 fpWaterDensity = G4DNAMolecularMaterial::Instance()->
169 GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
170
171 if (isInitialised) return;
172
174 isInitialised = true;
175}
176
177//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
178
182 G4double ekin,
183 G4double,
184 G4double)
185{
186#ifdef MELTON_VERBOSE
187 if (verboseLevel > 3)
188 G4cout
189 << "Calling CrossSectionPerVolume() of G4DNAMeltonAttachmentModel"
190 << G4endl;
191#endif
192
193 // Calculate total cross section for model
194
195 G4double sigma = 0.;
196
197 G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
198
199 if (ekin >= LowEnergyLimit() && ekin <= HighEnergyLimit())
200 sigma = fData->FindValue(ekin);
201
202#ifdef MELTON_VERBOSE
203 if (verboseLevel > 2)
204 {
205 G4cout << "__________________________________" << G4endl;
206 G4cout << "=== G4DNAMeltonAttachmentModel - XS INFO START" << G4endl;
207 G4cout << "--- Kinetic energy(eV)=" << ekin/eV
208 << " particle : " << particleDefinition->GetParticleName()
209 << G4endl;
210 G4cout << "--- Cross section per water molecule (cm^2)="
211 << sigma/cm/cm << G4endl;
212 G4cout << "--- Cross section per water molecule (cm^-1)="
213 << sigma*waterDensity/(1./cm) << G4endl;
214 G4cout << "--- G4DNAMeltonAttachmentModel - XS INFO END" << G4endl;
215 }
216#endif
217
218 return sigma*waterDensity;
219}
220
221//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
222
223void
225SampleSecondaries(std::vector<G4DynamicParticle*>* /*fvect*/,
226 const G4MaterialCutsCouple* /*couple*/,
227 const G4DynamicParticle* aDynamicElectron,
228 G4double,
229 G4double)
230{
231
232#ifdef MELTON_VERBOSE
233 if (verboseLevel > 3)
234 G4cout
235 << "Calling SampleSecondaries() of G4DNAMeltonAttachmentModel" << G4endl;
236#endif
237
238 // Electron is killed
239
240 G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
241
242 if (!statCode)
243 {
247 }
248
249 else
250 {
253 }
254
255 if(fDissociationFlag)
256 {
258 CreateWaterMolecule(eDissociativeAttachment,
259 -1,
261 }
262 return;
263}
@ eDissociativeAttachment
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4DNAChemistryManager * Instance()
virtual G4double FindValue(G4double e, G4int componentId=0) const
virtual G4bool LoadData(const G4String &argFileName)
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
G4DNAMeltonAttachmentModel(const G4ParticleDefinition *p=0, const G4String &nam="DNAMeltonAttachmentModel")
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4ParticleChangeForGamma * fParticleChangeForGamma
static G4DNAMolecularMaterial * Instance()
G4double GetKineticEnergy() const
size_t GetIndex() const
Definition: G4Material.hh:255
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:691
void SetProposedKineticEnergy(G4double proposedKinEnergy)
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:746
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:124
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:753
void ProposeTrackStatus(G4TrackStatus status)
const G4Track * GetCurrentTrack() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)