Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAQuinnPlasmonExcitationModel.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// Created on 2016/04/08
27//
28// Authors: D. Sakata, S. Incerti
29//
30// This class perform transmission term of volume plasmon excitation,
31// based on Quinn Model, see Phys. Rev. vol 126, number 4 (1962)
32
34#include "G4SystemOfUnits.hh"
35#include "G4RandomDirection.hh"
36
37//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
38
39using namespace std;
40
41//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
42
45 const G4String& nam):
46G4VEmModel(nam), isInitialised(false)
47{
48 fpMaterialDensity = 0;
49 fLowEnergyLimit = 10 * eV;
50 fHighEnergyLimit = 1.0 * GeV;
51
52 for(G4int i=0;i<100;i++) nValenceElectron[i]=0;
53
54 verboseLevel = 0;
55
56 if (verboseLevel > 0)
57 {
58 G4cout << "Quinn plasmon excitation model is constructed " << G4endl;
59 }
61 statCode = false;
62}
63
64//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
65
67{
68}
69
70//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
71
73 (const G4ParticleDefinition* particle,
74 const G4DataVector& /*cuts*/)
75{
76 for(G4int i=0;i<100;i++) nValenceElectron[i]=0;
77
78 if (verboseLevel > 3)
79 {
80 G4cout <<
81 "Calling G4DNAQuinnPlasmonExcitationModel::Initialise()"
82 << G4endl;
83 }
84
85 if(particle == G4Electron::ElectronDefinition())
86 {
87 fLowEnergyLimit = 10 * eV;
88 fHighEnergyLimit = 1.0 * GeV;
89 }
90 else
91 {
92 G4Exception("G4DNAQuinnPlasmonExcitationModel::Initialise","em0001",
93 FatalException,"Not defined for other particles than electrons.");
94 return;
95 }
96
97 // Get Number of valence electrons
98 G4ProductionCutsTable* theCoupleTable =
100
101 G4int numOfCouples = (G4int)theCoupleTable->GetTableSize();
102
103 for(G4int i=0;i<numOfCouples;i++){
104
105 const G4MaterialCutsCouple* couple =
106 theCoupleTable->GetMaterialCutsCouple(i);
107
108 const G4Material* material = couple->GetMaterial();
109
110 const G4ElementVector* theElementVector =material->GetElementVector();
111
112 std::size_t nelm = material->GetNumberOfElements();
113 if (nelm==1) // Protection: only for single element
114 {
115 G4int z = G4lrint((*theElementVector)[0]->GetZ());
116 if(z<=100)
117 {
118 nValenceElectron[z] = GetNValenceElectron(z);
119 }
120 else
121 {
122 G4Exception("G4DNAQuinnPlasmonExcitationModel::Initialise","em0002",
123 FatalException,"The model is not applied for z>100");
124 }
125 }
126 //for(G4int j=0;j<nelm;j++){
127 // G4int z=G4lrint((*theElementVector)[j]->GetZ());
128 // if(z<=100){nValenceElectron[z] = GetNValenceElectron(z);}
129 //}
130 }
131
132 if( verboseLevel>0 )
133 {
134 G4cout << "Quinn plasmon excitation model is initialized " << G4endl
135 << "Energy range: "
136 << LowEnergyLimit() / eV << " eV - "
137 << HighEnergyLimit() / keV << " keV for "
138 << particle->GetParticleName()
139 << G4endl;
140 }
141
142 if (isInitialised){return;}
144 isInitialised = true;
145}
146
147//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
148
150 (const G4Material* material,
151 const G4ParticleDefinition* particleDefinition,
152 G4double ekin,
153 G4double,
154 G4double)
155{
156 if (verboseLevel > 3)
157 {
158 G4cout <<
159 "Calling CrossSectionPerVolume() of G4DNAQuinnPlasmonExcitationModel"
160 << G4endl;
161 }
162
163 // Protection: only for single element
164 if(material->GetNumberOfElements()>1) return 0.;
165 G4double z = material->GetZ();
166
167 // Protection: only for Gold
168 if (z!=79){return 0.;}
169
170
171 G4double sigma = 0;
172 G4double atomicNDensity = material->GetAtomicNumDensityVector()[0];
173
174 if(atomicNDensity!= 0.0)
175 {
176 if (ekin >= fLowEnergyLimit && ekin < fHighEnergyLimit)
177 {
178 sigma = GetCrossSection(material,particleDefinition,ekin);
179 }
180
181 if (verboseLevel > 2)
182 {
183 G4cout<<"__________________________________" << G4endl;
184 G4cout<<"=== G4DNAQuinnPlasmonExcitationModel - XS INFO START"<<G4endl;
185 G4cout<<"=== Kinetic energy (eV)=" << ekin/eV << " particle : "
186 <<particleDefinition->GetParticleName() << G4endl;
187 G4cout<<"=== Cross section per atom for Z="<<z<<" is (cm^2)"
188 <<sigma/cm/cm << G4endl;
189 G4cout<<"=== Cross section per atom for Z="<<z<<" is (cm^-1)="
190 <<sigma*atomicNDensity/(1./cm) << G4endl;
191 G4cout<<"=== G4DNAQuinnPlasmonExcitationModel - XS INFO END" << G4endl;
192 }
193 }
194
195 return sigma*atomicNDensity;
196}
197
198//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
199
201 (std::vector<G4DynamicParticle*>* /*fvect*/,
202 const G4MaterialCutsCouple* couple,
203 const G4DynamicParticle* aDynamicParticle,
205{
206
207 if (verboseLevel > 3)
208 {
209 G4cout <<
210 "Calling SampleSecondaries() of G4DNAQuinnPlasmonExcitationModel"
211 << G4endl;
212 }
213
214 const G4Material *material = couple->GetMaterial();
215
216 G4ParticleDefinition* particle = aDynamicParticle->GetDefinition();
217
218 G4double k = aDynamicParticle->GetKineticEnergy();
219
220 if(particle == G4Electron::ElectronDefinition())
221 {
222 G4double e = 1.;
223 G4int z = material->GetZ();
224 G4int Nve = 0;
225
226 //TODO: have to be change to realistic!!
227 if(z<100) Nve = nValenceElectron[z];
228
229 G4double A = material->GetA()/g/mole;
230 G4double Dens = material->GetDensity()/g*cm*cm*cm;
231 G4double veDens = Dens*CLHEP::Avogadro*Nve/A;
232
233 G4double omega_p = std::sqrt(veDens*std::pow(e,2)/
234 (CLHEP::epsilon0/(1./cm)*CLHEP::electron_mass_c2
235 /(CLHEP::c_squared/cm/cm)));
236
237 G4double excitationEnergy = CLHEP::hbar_Planck*omega_p;
238 G4double newEnergy = k - excitationEnergy;
239
240
241 if (newEnergy > 0)
242 {
244 ProposeMomentumDirection(aDynamicParticle->GetMomentumDirection());
245
247
248 if(!statCode)
249 {
251 }
252 else
253 {
255
256 }
257 }
258 }
259}
260
261//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
262
264 (const G4Material* material,
265 const G4ParticleDefinition* particle,
266 G4double kineticEnergy)
267{
268 G4double value=0;
269
270 if(particle == G4Electron::ElectronDefinition())
271 {
272 G4double e = 1.;
273 G4int z = material->GetZ();
274 G4int Nve = 0;
275 if(z<100) Nve = nValenceElectron[z];
276 G4double A = material->GetA()/g/mole;
277 G4double Dens = material->GetDensity()/g*cm*cm*cm;
278 G4double veDens = Dens*CLHEP::Avogadro*Nve/A;
279
280 G4double omega_p = std::sqrt(veDens*std::pow(e,2)
281 /(CLHEP::epsilon0/(1./cm)*CLHEP::electron_mass_c2/
282 (CLHEP::c_squared/cm/cm)));
283
284 G4double fEnergy = std::pow(CLHEP::h_Planck,2)/(8*CLHEP::electron_mass_c2)*
285 std::pow(3*veDens/CLHEP::pi,2./3.)/e
286 *(CLHEP::c_squared/cm/cm);
287
288 G4double p0 = sqrt(2*CLHEP::electron_mass_c2
289 /(CLHEP::c_squared/cm/cm)*fEnergy);
290
291 G4double p = sqrt(2*CLHEP::electron_mass_c2
292 /(CLHEP::c_squared/cm/cm)*kineticEnergy);
293
294 G4double mfp = 2*CLHEP::Bohr_radius/cm*kineticEnergy
295 /(CLHEP::hbar_Planck*omega_p)/
296 (G4Log((std::pow(std::pow(p0,2)
297 +2*CLHEP::electron_mass_c2/
298 (CLHEP::c_squared/cm/cm)*omega_p
299 *CLHEP::hbar_Planck,1./2.)-p0)
300 /(p-std::pow(std::pow(p,2)-2*CLHEP::electron_mass_c2/
301 (CLHEP::c_squared/cm/cm)*omega_p
302 *CLHEP::hbar_Planck,1./2.))));
303
304 G4double excitationEnergy = CLHEP::hbar_Planck*omega_p;
305
306 if((0<mfp)&&(0<veDens)&&(excitationEnergy<kineticEnergy)){
307 value = 1./(veDens*mfp);
308 }
309 }
310 return value*cm*cm;
311}
312
313//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
314
315G4int G4DNAQuinnPlasmonExcitationModel::GetNValenceElectron(G4int z)
316{
317
318 G4int Nve=0;
319
320 // Current limitation to gold
321 if (z!=79){return 0.;}
322
323 if (verboseLevel > 3)
324 {
325 G4cout <<
326 "Calling GetNValenceElectron() of G4DNAQuinnPlasmonExcitationModel"
327 << G4endl;
328 }
329
330 const char *datadir=0;
331
332 if(!datadir)
333 {
334 datadir = G4FindDataDir("G4LEDATA");
335 if(!datadir)
336 {
337 G4Exception("G4DNAQuinnPlasmonExcitationModel::GetNValenceElectron()"
338 ,"em0002",FatalException,
339 "Enviroment variable G4LEDATA not defined");
340 return 0;
341 }
342 }
343
344 std::ostringstream targetfile;
345 targetfile.str("");
346 targetfile.clear(stringstream::goodbit);
347 targetfile << datadir <<"/dna/atomicstate_Z"<< z <<".dat";
348 std::ifstream fin(targetfile.str().c_str());
349
350 if(!fin)
351 {
352 G4cout<< " Error : "<< targetfile.str() <<" is not found "<<endl;
353 G4Exception("G4DNAQuinnPlasmonExcitationModel::GetNValenceElectron()"
354 ,"em0003",FatalException,
355 "There is no target file");
356 return 0;
357 }
358
359 string buff0,buff1,buff2,buff3,buff4,buff5,buff6;
360 fin >> buff0 >>buff1>>buff2>>buff3>>buff4>>buff5>>buff6;
361
362 while(true){
363 fin >> buff0 >>buff1>>buff2>>buff3>>buff4>>buff5>>buff6;
364 if(!fin.eof())
365 {
366 Nve = stoi(buff3);
367 }
368 else
369 {
370 break;
371 }
372 }
373 return Nve;
374}
375
std::vector< const G4Element * > G4ElementVector
const char * G4FindDataDir(const char *)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
G4double G4Log(G4double x)
Definition: G4Log.hh:227
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4double A[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &= *(new G4DataVector()))
virtual G4double GetCrossSection(const G4Material *material, const G4ParticleDefinition *, G4double kineticEnergy)
G4DNAQuinnPlasmonExcitationModel(const G4ParticleDefinition *p=0, const G4String &nam="DNAQuinnPlasmonExcitationModel")
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:88
const G4Material * GetMaterial() const
G4double GetDensity() const
Definition: G4Material.hh:175
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:185
G4double GetZ() const
Definition: G4Material.cc:745
size_t GetNumberOfElements() const
Definition: G4Material.hh:181
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:211
G4double GetA() const
Definition: G4Material.cc:759
void SetProposedKineticEnergy(G4double proposedKinEnergy)
const G4String & GetParticleName() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:124
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
int G4lrint(double ad)
Definition: templates.hh:134