source: trunk/source/processes/electromagnetic/lowenergy/src/G4DNAMeltonAttachmentModel.cc

Last change on this file was 1350, checked in by garnier, 15 years ago

update to last version 4.9.4

File size: 7.4 KB
Line 
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: G4DNAMeltonAttachmentModel.cc,v 1.2 2010/09/15 05:47:33 sincerti Exp $
27// GEANT4 tag $Name: geant4-09-04-ref-00 $
28//
29
30// Created by Z. Francis
31
32#include "G4DNAMeltonAttachmentModel.hh"
33
34//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
35
36using namespace std;
37
38//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
39
40G4DNAMeltonAttachmentModel::G4DNAMeltonAttachmentModel(const G4ParticleDefinition*,
41 const G4String& nam)
42:G4VEmModel(nam),isInitialised(false)
43{
44
45 lowEnergyLimit = 4 * eV;
46 lowEnergyLimitOfModel = 4 * eV;
47 highEnergyLimit = 13 * eV;
48 SetLowEnergyLimit(lowEnergyLimit);
49 SetHighEnergyLimit(highEnergyLimit);
50
51 verboseLevel= 0;
52 // Verbosity scale:
53 // 0 = nothing
54 // 1 = warning for energy non-conservation
55 // 2 = details of energy budget
56 // 3 = calculation of cross sections, file openings, sampling of atoms
57 // 4 = entering in methods
58
59 if( verboseLevel>0 )
60 {
61 G4cout << "Melton Attachment model is constructed " << G4endl
62 << "Energy range: "
63 << lowEnergyLimit / eV << " eV - "
64 << highEnergyLimit / eV << " eV"
65 << G4endl;
66 }
67}
68
69//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
70
71G4DNAMeltonAttachmentModel::~G4DNAMeltonAttachmentModel()
72{
73 // For total cross section
74
75 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
76
77 for (pos = tableData.begin(); pos != tableData.end(); ++pos)
78 {
79 G4DNACrossSectionDataSet* table = pos->second;
80 delete table;
81 }
82
83 // For final state
84
85}
86
87//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
88
89void G4DNAMeltonAttachmentModel::Initialise(const G4ParticleDefinition* /*particle*/,
90 const G4DataVector& /*cuts*/)
91{
92
93 if (verboseLevel > 3)
94 G4cout << "Calling G4DNAMeltonAttachmentModel::Initialise()" << G4endl;
95
96 // Energy limits
97
98 if (LowEnergyLimit() < lowEnergyLimit)
99 {
100 G4cout << "G4DNAMeltonAttachmentModel: low energy limit increased from " <<
101 LowEnergyLimit()/eV << " eV to " << lowEnergyLimit/eV << " eV" << G4endl;
102 SetLowEnergyLimit(lowEnergyLimit);
103 }
104
105 if (HighEnergyLimit() > highEnergyLimit)
106 {
107 G4cout << "G4DNAMeltonAttachmentModel: high energy limit decreased from " <<
108 HighEnergyLimit()/eV << " eV to " << highEnergyLimit/eV << " eV" << G4endl;
109 SetHighEnergyLimit(highEnergyLimit);
110 }
111
112 // Reading of data files
113
114 G4double scaleFactor = 1e-18*cm*cm;
115
116 G4String fileElectron("dna/sigma_attachment_e_melton");
117
118 G4ParticleDefinition* electronDef = G4Electron::ElectronDefinition();
119 G4String electron;
120
121 if (electronDef != 0)
122 {
123 // For total cross section
124
125 electron = electronDef->GetParticleName();
126
127 tableFile[electron] = fileElectron;
128
129 G4DNACrossSectionDataSet* tableE = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, eV,scaleFactor );
130 tableE->LoadData(fileElectron);
131 tableData[electron] = tableE;
132
133 }
134 else G4Exception("G4DNAMeltonAttachmentModel::Initialise: electron is not defined");
135
136 if (verboseLevel > 2)
137 G4cout << "Loaded cross section data for Melton Attachment model" << G4endl;
138
139 if( verboseLevel>0 )
140 {
141 G4cout << "Melton Attachment model is initialized " << G4endl
142 << "Energy range: "
143 << LowEnergyLimit() / eV << " eV - "
144 << HighEnergyLimit() / eV << " eV"
145 << G4endl;
146 }
147
148 if(!isInitialised)
149 {
150 isInitialised = true;
151
152 if(pParticleChange)
153 fParticleChangeForGamma = reinterpret_cast<G4ParticleChangeForGamma*>(pParticleChange);
154 else
155 fParticleChangeForGamma = new G4ParticleChangeForGamma();
156 }
157
158}
159
160//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
161
162G4double G4DNAMeltonAttachmentModel::CrossSectionPerVolume(const G4Material* material,
163 const G4ParticleDefinition* p,
164 G4double ekin,
165 G4double,
166 G4double)
167{
168 if (verboseLevel > 3)
169 G4cout << "Calling CrossSectionPerVolume() of G4DNAMeltonAttachmentModel" << G4endl;
170
171 // Calculate total cross section for model
172
173 G4double sigma=0;
174
175 if (material->GetName() == "G4_WATER")
176 {
177 const G4String& particleName = p->GetParticleName();
178
179 if (ekin >= lowEnergyLimit && ekin < highEnergyLimit)
180 {
181
182 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
183 pos = tableData.find(particleName);
184
185 if (pos != tableData.end())
186 {
187 G4DNACrossSectionDataSet* table = pos->second;
188 if (table != 0)
189 {
190 sigma = table->FindValue(ekin);
191 }
192 }
193 else
194 {
195 G4Exception("G4DNAMeltonAttachmentModel::ComputeCrossSectionPerVolume: attempting to calculate cross section for wrong particle");
196 }
197 }
198
199 if (verboseLevel > 3)
200 {
201 G4cout << "---> Kinetic energy(eV)=" << ekin/eV << G4endl;
202 G4cout << " - Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
203 G4cout << " - Cross section per water molecule (cm^-1)=" << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
204 }
205
206 } // if water
207
208 return sigma*material->GetAtomicNumDensityVector()[1];
209}
210
211//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
212
213void G4DNAMeltonAttachmentModel::SampleSecondaries(std::vector<G4DynamicParticle*>* /*fvect*/,
214 const G4MaterialCutsCouple* /*couple*/,
215 const G4DynamicParticle* aDynamicElectron,
216 G4double,
217 G4double)
218{
219
220 if (verboseLevel > 3)
221 G4cout << "Calling SampleSecondaries() of G4DNAMeltonAttachmentModel" << G4endl;
222
223 // Electron is killed
224
225 G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
226 fParticleChangeForGamma->ProposeTrackStatus(fStopAndKill);
227 fParticleChangeForGamma->ProposeLocalEnergyDeposit(electronEnergy0);
228
229 return ;
230}
231
232
233
234
235
236
237
238
Note: See TracBrowser for help on using the repository browser.