Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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// $Id$
27//
28
29// Created by Z. Francis
30
32#include "G4SystemOfUnits.hh"
35
36//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
37
38using namespace std;
39
40//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
41
43 const G4String& nam)
44 :G4VEmModel(nam),isInitialised(false)
45{
46// nistwater = G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER");
47 fpWaterDensity = 0;
48
49 lowEnergyLimit = 4 * eV;
50 lowEnergyLimitOfModel = 4 * eV;
51 highEnergyLimit = 13 * eV;
52 SetLowEnergyLimit(lowEnergyLimit);
53 SetHighEnergyLimit(highEnergyLimit);
54
55 verboseLevel= 0;
56 // Verbosity scale:
57 // 0 = nothing
58 // 1 = warning for energy non-conservation
59 // 2 = details of energy budget
60 // 3 = calculation of cross sections, file openings, sampling of atoms
61 // 4 = entering in methods
62
63 if( verboseLevel>0 )
64 {
65 G4cout << "Melton Attachment model is constructed " << G4endl
66 << "Energy range: "
67 << lowEnergyLimit / eV << " eV - "
68 << highEnergyLimit / eV << " eV"
69 << G4endl;
70 }
72 fDissociationFlag = true;
73}
74
75//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
76
78{
79 // For total cross section
80
81 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
82
83 for (pos = tableData.begin(); pos != tableData.end(); ++pos)
84 {
85 G4DNACrossSectionDataSet* table = pos->second;
86 delete table;
87 }
88
89 // For final state
90
91}
92
93//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
94
96 const G4DataVector& /*cuts*/)
97{
98
99 if (verboseLevel > 3)
100 G4cout << "Calling G4DNAMeltonAttachmentModel::Initialise()" << G4endl;
101
102 // Energy limits
103
104 if (LowEnergyLimit() < lowEnergyLimit)
105 {
106 G4cout << "G4DNAMeltonAttachmentModel: low energy limit increased from " <<
107 LowEnergyLimit()/eV << " eV to " << lowEnergyLimit/eV << " eV" << G4endl;
108 SetLowEnergyLimit(lowEnergyLimit);
109 }
110
111 if (HighEnergyLimit() > highEnergyLimit)
112 {
113 G4cout << "G4DNAMeltonAttachmentModel: high energy limit decreased from " <<
114 HighEnergyLimit()/eV << " eV to " << highEnergyLimit/eV << " eV" << G4endl;
115 SetHighEnergyLimit(highEnergyLimit);
116 }
117
118 // Reading of data files
119
120 G4double scaleFactor = 1e-18*cm*cm;
121
122 G4String fileElectron("dna/sigma_attachment_e_melton");
123
125 G4String electron;
126
127 // ELECTRON
128
129 // For total cross section
130
131 electron = electronDef->GetParticleName();
132
133 tableFile[electron] = fileElectron;
134
136 tableE->LoadData(fileElectron);
137 tableData[electron] = tableE;
138
139 //
140
141 if (verboseLevel > 2)
142 G4cout << "Loaded cross section data for Melton Attachment model" << G4endl;
143
144 if( verboseLevel>0 )
145 {
146 G4cout << "Melton Attachment model is initialized " << G4endl
147 << "Energy range: "
148 << LowEnergyLimit() / eV << " eV - "
149 << HighEnergyLimit() / eV << " eV"
150 << G4endl;
151 }
152 // Initialize water density pointer
154
155 if (isInitialised) { return; }
157 isInitialised = true;
158
159}
160
161//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
162
164 const G4ParticleDefinition* particleDefinition,
165 G4double ekin,
166 G4double,
167 G4double)
168{
169 if (verboseLevel > 3)
170 G4cout << "Calling CrossSectionPerVolume() of G4DNAMeltonAttachmentModel" << G4endl;
171
172 // Calculate total cross section for model
173
174 G4double sigma=0;
175
176 G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
177
178 if(waterDensity!= 0.0)
179 // if (material == nistwater || material->GetBaseMaterial() == nistwater)
180 {
181 const G4String& particleName = particleDefinition->GetParticleName();
182
183 if (ekin >= lowEnergyLimit && ekin < highEnergyLimit)
184 {
185
186 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
187 pos = tableData.find(particleName);
188
189 if (pos != tableData.end())
190 {
191 G4DNACrossSectionDataSet* table = pos->second;
192 if (table != 0)
193 {
194 sigma = table->FindValue(ekin);
195 }
196 }
197 else
198 {
199 G4Exception("G4DNAMeltonAttachmentModel::ComputeCrossSectionPerVolume","em0002",
200 FatalException,"Model not applicable to particle type.");
201 }
202 }
203
204 if (verboseLevel > 2)
205 {
206 G4cout << "__________________________________" << G4endl;
207 G4cout << "°°° G4DNAMeltonAttachmentModel - XS INFO START" << G4endl;
208 G4cout << "°°° Kinetic energy(eV)=" << ekin/eV << " particle : " << particleDefinition->GetParticleName() << G4endl;
209 G4cout << "°°° Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
210 G4cout << "°°° Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
211 // G4cout << " - Cross section per water molecule (cm^-1)=" << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
212 G4cout << "°°° G4DNAMeltonAttachmentModel - XS INFO END" << G4endl;
213 }
214
215 } // if water
216
217 return sigma*waterDensity;
218// return sigma*material->GetAtomicNumDensityVector()[1];
219}
220
221//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
222
223void G4DNAMeltonAttachmentModel::SampleSecondaries(std::vector<G4DynamicParticle*>* /*fvect*/,
224 const G4MaterialCutsCouple* /*couple*/,
225 const G4DynamicParticle* aDynamicElectron,
226 G4double,
227 G4double)
228{
229
230 if (verboseLevel > 3)
231 G4cout << "Calling SampleSecondaries() of G4DNAMeltonAttachmentModel" << G4endl;
232
233 // Electron is killed
234
235 G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
239
240 if(fDissociationFlag)
241 {
244 }
245 return ;
246}
@ eDissociativeAttachment
@ FatalException
@ fStopAndKill
double G4double
Definition: G4Types.hh:64
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
static G4DNAChemistryManager * Instance()
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
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()
const std::vector< double > * GetNumMolPerVolTableFor(const G4Material *) const
G4double GetKineticEnergy() const
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
size_t GetIndex() const
Definition: G4Material.hh:261
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:576
const G4Track * GetCurrentTrack() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:529
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:592
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41