Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4MicroElecElasticModel_new.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// G4MicroElecElasticModel_new.cc, 2011/08/29 A.Valentin, M. Raine are with CEA [a]
28// 2020/05/20 P. Caron, C. Inguimbert are with ONERA [b]
29// Q. Gibaru is with CEA [a], ONERA [b] and CNES [c]
30// M. Raine and D. Lambert are with CEA [a]
31//
32// A part of this work has been funded by the French space agency(CNES[c])
33// [a] CEA, DAM, DIF - 91297 ARPAJON, France
34// [b] ONERA - DPHY, 2 avenue E.Belin, 31055 Toulouse, France
35// [c] CNES, 18 av.E.Belin, 31401 Toulouse CEDEX, France
36//
37// Based on the following publications
38// - A.Valentin, M. Raine,
39// Inelastic cross-sections of low energy electrons in silicon
40// for the simulation of heavy ion tracks with the Geant4-DNA toolkit,
41// NSS Conf. Record 2010, pp. 80-85
42// https://doi.org/10.1109/NSSMIC.2010.5873720
43//
44// - A.Valentin, M. Raine, M.Gaillardin, P.Paillet
45// Geant4 physics processes for microdosimetry simulation:
46// very low energy electromagnetic models for electrons in Silicon,
47// https://doi.org/10.1016/j.nimb.2012.06.007
48// NIM B, vol. 288, pp. 66-73, 2012, part A
49// heavy ions in Si, NIM B, vol. 287, pp. 124-129, 2012, part B
50// https://doi.org/10.1016/j.nimb.2012.07.028
51//
52// - M. Raine, M. Gaillardin, P. Paillet
53// Geant4 physics processes for silicon microdosimetry simulation:
54// Improvements and extension of the energy-range validity up to 10 GeV/nucleon
55// NIM B, vol. 325, pp. 97-100, 2014
56// https://doi.org/10.1016/j.nimb.2014.01.014
57//
58// - J. Pierron, C. Inguimbert, M. Belhaj, T. Gineste, J. Puech, M. Raine
59// Electron emission yield for low energy electrons:
60// Monte Carlo simulation and experimental comparison for Al, Ag, and Si
61// Journal of Applied Physics 121 (2017) 215107.
62// https://doi.org/10.1063/1.4984761
63//
64// - P. Caron,
65// Study of Electron-Induced Single-Event Upset in Integrated Memory Devices
66// PHD, 16th October 2019
67//
68// - Q.Gibaru, C.Inguimbert, P.Caron, M.Raine, D.Lambert, J.Puech,
69// Geant4 physics processes for microdosimetry and secondary electron emission simulation :
70// Extension of MicroElec to very low energies and new materials
71// NIM B, 2020, in review.
72//
73//
74//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
77#include "G4SystemOfUnits.hh"
78#include "G4Exp.hh"
79#include "G4Material.hh"
80#include "G4String.hh"
81
82//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
83
84using namespace std;
85
86//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
87
89 const G4String& nam)
90 :G4VEmModel(nam), isInitialised(false)
91{
92 killBelowEnergy = 0.1*eV; // Minimum e- energy for energy loss by excitation
93 lowEnergyLimit = 0.1 * eV;
94 lowEnergyLimitOfModel = 10 * eV; // The model lower energy is 10 eV
95 highEnergyLimit = 500. * keV;
96 SetLowEnergyLimit(lowEnergyLimit);
97 SetHighEnergyLimit(highEnergyLimit);
98
99 verboseLevel= 0;
100 // Verbosity scale:
101 // 0 = nothing
102 // 1 = warning for energy non-conservation
103 // 2 = details of energy budget
104 // 3 = calculation of cross sections, file openings, sampling of atoms
105 // 4 = entering in methods
106
107 if( verboseLevel>0 )
108 {
109 G4cout << "MicroElec Elastic model is constructed " << G4endl
110 << "Energy range: "
111 << lowEnergyLimit / eV << " eV - "
112 << highEnergyLimit / MeV << " MeV"
113 << G4endl;
114 }
116
117 killElectron = false;
118 acousticModelEnabled = false;
119 currentMaterialName = "";
120 isOkToBeInitialised = false;
121}
122
123//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
124
126{
127 // For total cross section
128 TCSMap::iterator pos2;
129 for (pos2 = tableTCS.begin(); pos2 != tableTCS.end(); ++pos2) {
130 MapData* tableData = pos2->second;
131 std::map< G4String, G4MicroElecCrossSectionDataSet_new*, std::less<G4String> >::iterator pos;
132 for (pos = tableData->begin(); pos != tableData->end(); ++pos)
133 {
134 G4MicroElecCrossSectionDataSet_new* table = pos->second;
135 delete table;
136 }
137 delete tableData;
138 }
139
140 //Clearing DCS maps
141
142 ThetaMap::iterator iterator_angle;
143 for (iterator_angle = thetaDataStorage.begin(); iterator_angle != thetaDataStorage.end(); ++iterator_angle) {
144 TriDimensionMap* eDiffCrossSectionData = iterator_angle->second;
145 eDiffCrossSectionData->clear();
146 delete eDiffCrossSectionData;
147 }
148
149 energyMap::iterator iterator_energy;
150 for (iterator_energy = eIncidentEnergyStorage.begin(); iterator_energy != eIncidentEnergyStorage.end(); ++iterator_energy) {
151 std::vector<G4double>* eTdummyVec = iterator_energy->second;
152 eTdummyVec->clear();
153 delete eTdummyVec;
154 }
155
156 ProbaMap::iterator iterator_proba;
157 for (iterator_proba = eProbaStorage.begin(); iterator_proba != eProbaStorage.end(); ++iterator_proba) {
158 VecMap* eVecm = iterator_proba->second;
159 eVecm->clear();
160 delete eVecm;
161 }
162
163}
164
165//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
166
168 const G4DataVector& /*cuts*/)
169{
170 if (isOkToBeInitialised == true && isInitialised == false) {
171
172 if (verboseLevel > -1)
173 G4cout << "Calling G4MicroElecElasticModel_new::Initialise()" << G4endl;
174 // Energy limits
175 // Reading of data files
176
177 G4double scaleFactor = 1e-18 * cm * cm;
178
179 G4ProductionCutsTable* theCoupleTable =
181 G4int numOfCouples = (G4int)theCoupleTable->GetTableSize();
182
183 for (G4int i = 0; i < numOfCouples; ++i) {
184 const G4Material* material =
185 theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
186
187 //theCoupleTable->GetMaterialCutsCouple(i)->;
188
189 G4cout << "MicroElasticModel, Material " << i + 1 << " / " << numOfCouples << " : " << material->GetName() << G4endl;
190 if (material->GetName() == "Vacuum") continue;
191
192 G4String matName = material->GetName().substr(3, material->GetName().size());
193 G4cout<< matName<< G4endl;
194
195 currentMaterialStructure = new G4MicroElecMaterialStructure(matName);
196 lowEnergyLimitTable[matName]=currentMaterialStructure->GetElasticModelLowLimit();
197 highEnergyLimitTable[matName]=currentMaterialStructure->GetElasticModelHighLimit();
198 workFunctionTable[matName] = currentMaterialStructure->GetWorkFunction();
199
200 delete currentMaterialStructure;
201
202 G4cout << "Reading TCS file" << G4endl;
203 G4String fileElectron = "Elastic/elsepa_elastic_cross_e_" + matName;
204 G4cout << "Elastic Total Cross file : " << fileElectron << G4endl;
205
207 G4String electron = electronDef->GetParticleName();
208
209 // For total cross section
210 MapData* tableData = new MapData();
211
213 tableE->LoadData(fileElectron);
214 tableData->insert(make_pair(electron, tableE));
215 tableTCS[matName] = tableData; //Storage of TCS
216
217 // For final state
218 const char* path = G4FindDataDir("G4LEDATA");
219 if (!path)
220 {
221 G4Exception("G4MicroElecElasticModel_new::Initialise","em0006",FatalException,"G4LEDATA environment variable not set.");
222 return;
223 }
224
225 //Reading DCS file
226 std::ostringstream eFullFileName;
227 eFullFileName << path << "/microelec/Elastic/elsepa_elastic_cumulated_diffcross_e_" + matName + ".dat";
228 G4cout << "Elastic Cumulated Diff Cross : " << eFullFileName.str().c_str() << G4endl;
229 std::ifstream eDiffCrossSection(eFullFileName.str().c_str());
230
231 if (!eDiffCrossSection)
232 G4Exception("G4MicroElecElasticModel_new::Initialise", "em0003", FatalException, "Missing data file: /microelec/sigmadiff_cumulated_elastic_e_Si.dat");
233
234 // October 21th, 2014 - Melanie Raine
235 // Added clear for MT
236 // Diff Cross Sections in cumulated mode
237 TriDimensionMap* eDiffCrossSectionData = new TriDimensionMap(); //Angles
238 std::vector<G4double>* eTdummyVec = new std::vector<G4double>; //Incident energy vector
239 VecMap* eProbVec = new VecMap; //Probabilities
240
241 eTdummyVec->push_back(0.);
242
243 while (!eDiffCrossSection.eof())
244 {
245 G4double tDummy; //incident energy
246 G4double eProb; //Proba
247 eDiffCrossSection >> tDummy >> eProb;
248
249 // SI : mandatory eVecm initialization
250 if (tDummy != eTdummyVec->back())
251 {
252 eTdummyVec->push_back(tDummy); //adding values for incident energy points
253 (*eProbVec)[tDummy].push_back(0.); //adding probability for the first angle, equal to 0
254 }
255
256 eDiffCrossSection >> (*eDiffCrossSectionData)[tDummy][eProb]; //adding Angle Value to map
257
258 if (eProb != (*eProbVec)[tDummy].back()) {
259 (*eProbVec)[tDummy].push_back(eProb); //Adding cumulated proba to map
260 }
261
262 }
263
264 //Filling maps for the material
265 thetaDataStorage[matName] = eDiffCrossSectionData;
266 eIncidentEnergyStorage[matName] = eTdummyVec;
267 eProbaStorage[matName] = eProbVec;
268 }
269 // End final state
270
271 if (verboseLevel > 2)
272 G4cout << "Loaded cross section files for MicroElec Elastic model" << G4endl;
273
274 if (verboseLevel > 0)
275 {
276 G4cout << "MicroElec Elastic model is initialized " << G4endl
277 << "Energy range: "
278 << LowEnergyLimit() / eV << " eV - "
279 << HighEnergyLimit() / MeV << " MeV"
280 << G4endl; // system("pause"); linux doesn't like
281 }
282
283 if (isInitialised) { return; }
285 isInitialised = true;
286 }
287}
288
289//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
290
292 const G4ParticleDefinition* p,
293 G4double ekin,
294 G4double,
295 G4double)
296{
297 if (verboseLevel > 3)
298 G4cout << "Calling CrossSectionPerVolume() of G4MicroElecElasticModel" << G4endl;
299
300 isOkToBeInitialised = true;
301 currentMaterialName = material->GetName().substr(3, material->GetName().size());
302 const G4DataVector cuts;
303 Initialise(p, cuts);
304 // Calculate total cross section for model
305 MapEnergy::iterator lowEPos;
306 lowEPos = lowEnergyLimitTable.find(currentMaterialName);
307
308 MapEnergy::iterator highEPos;
309 highEPos = highEnergyLimitTable.find(currentMaterialName);
310
311 MapEnergy::iterator killEPos;
312 killEPos = workFunctionTable.find(currentMaterialName);
313
314 if (lowEPos == lowEnergyLimitTable.end() || highEPos == highEnergyLimitTable.end() || killEPos == workFunctionTable.end())
315 {
316 G4String str = "Material ";
317 str += currentMaterialName + " not found!";
318 G4Exception("G4MicroElecElasticModel_new::EnergyLimits", "em0002", FatalException, str);
319 return 0;
320 }
321 else {
322 // G4cout << "normal elastic " << G4endl;
323 lowEnergyLimit = lowEPos->second;
324 highEnergyLimit = highEPos->second;
325 killBelowEnergy = killEPos->second;
326
327 }
328
329 if (ekin < killBelowEnergy) {
330
331 return DBL_MAX; }
332
333 G4double sigma=0;
334
335 //Phonon for SiO2
336 if (currentMaterialName == "SILICON_DIOXIDE" && ekin < 100 * eV) {
337 acousticModelEnabled = true;
338
339 //Values for SiO2
340 G4double kbz = 11.54e9,
341 rho = 2.2 * 1000, // [g/cm3] * 1000
342 cs = 3560, //Sound speed
343 Ebz = 5.1 * 1.6e-19,
344 Aac = 17 * Ebz, //A screening parameter
345 Eac = 3.5 * 1.6e-19, //C deformation potential
346 prefactor = 2.2;// Facteur pour modifier les MFP
347
348 return AcousticCrossSectionPerVolume(ekin, kbz, rho, cs, Aac, Eac, prefactor);
349 }
350
351 //Elastic
352 else {
353 acousticModelEnabled = false;
354
355 G4double density = material->GetTotNbOfAtomsPerVolume();
356 const G4String& particleName = p->GetParticleName();
357
358 TCSMap::iterator tablepos;
359 tablepos = tableTCS.find(currentMaterialName);
360
361 if (tablepos != tableTCS.end())
362 {
363 MapData* tableData = tablepos->second;
364
365 if (ekin >= lowEnergyLimit && ekin < highEnergyLimit)
366 {
367 std::map< G4String, G4MicroElecCrossSectionDataSet_new*, std::less<G4String> >::iterator pos;
368 pos = tableData->find(particleName);
369
370 if (pos != tableData->end())
371 {
372 G4MicroElecCrossSectionDataSet_new* table = pos->second;
373 if (table != 0)
374 {
375 sigma = table->FindValue(ekin);
376 }
377 }
378 else
379 {
380 G4Exception("G4MicroElecElasticModel_new::ComputeCrossSectionPerVolume", "em0002", FatalException, "Model not applicable to particle type.");
381 }
382 }
383 else return 1 / DBL_MAX;
384 }
385 else
386 {
387 G4String str = "Material ";
388 str += currentMaterialName + " TCS table not found!";
389 G4Exception("G4MicroElecElasticModel_new::ComputeCrossSectionPerVolume", "em0002", FatalException, str);
390 }
391
392 if (verboseLevel > 3)
393 {
394 G4cout << "---> Kinetic energy(eV)=" << ekin / eV << G4endl;
395 G4cout << " - Cross section per Si atom (cm^2)=" << sigma / cm / cm << G4endl;
396 G4cout << " - Cross section per Si atom (cm^-1)=" << sigma*density / (1. / cm) << G4endl;
397 }
398 return sigma*density;
399 }
400}
401
402//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
403
405 G4double kbz,
406 G4double rho,
407 G4double cs,
408 G4double Aac,
409 G4double Eac,
410 G4double prefactor)
411{
412
413 G4double e = 1.6e-19,
414 m0 = 9.10938356e-31,
415 h = 1.0546e-34,
416 kb = 1.38e-23;
417
418 G4double E = (ekin / eV) * e;
419 G4double D = (2 / (std::sqrt(2) * std::pow(pi, 2) * std::pow(h, 3))) * (1 + 2 * E) * std::pow(m0, 1.5) * std::sqrt(E);
420
421 // Parametres SiO2
422 G4double T = 300,
423 Ebz = (std::pow(h, 2) * std::pow(kbz, 2)) / (2 * m0),
424 hwbz = cs * kbz * h,
425 nbz = 1.0 / (exp(hwbz / (kb * T)) - 1),
426 Pac;
427
428 if (E < Ebz / 4.0)
429 {
430 Pac = ((pi * kb * T) / (h * std::pow(cs, 2) * rho)) * (std::pow(Eac, 2) * D) / (1 + (E / Aac));
431 }
432
433 else if (E > Ebz) //Screened relationship
434 {
435 Pac = ((2 * pi * m0 * (2 * nbz + 1)) / (h * rho * hwbz)) * std::pow(Eac, 2) * D * E * 2 * std::pow((Aac / E), 2) * (((-E / Aac) / (1 + (E / Aac))) + log(1 + (E / Aac)));
436 }
437 else //Linear interpolation
438 {
439 G4double fEbz = ((2 * pi * m0 * (2 * nbz + 1)) / (h * rho * hwbz)) * std::pow(Eac, 2) * D * Ebz * 2 * std::pow((Aac / Ebz), 2) * (((-Ebz / Aac) / (1 + (Ebz / Aac))) + log(1 + (Ebz / Aac)));
440 G4double fEbz4 = ((pi * kb * T) / (h * std::pow(cs, 2) * rho)) * (std::pow(Eac, 2) * D) / (1 + ((Ebz / 4) / Aac));
441 G4double alpha = ((fEbz - fEbz4) / (Ebz - (Ebz / 4)));
442 Pac = alpha * E + (fEbz - alpha * Ebz);
443 }
444
445 G4double MFP = (std::sqrt(2 * E / m0) / (prefactor * Pac)) * m;
446
447 return 1 / MFP;
448}
449
450
451//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
452
453void G4MicroElecElasticModel_new::SampleSecondaries(std::vector<G4DynamicParticle*>* /*fvect*/,
454 const G4MaterialCutsCouple* /*couple*/,
455 const G4DynamicParticle* aDynamicElectron,
456 G4double,
457 G4double)
458{
459
460 if (verboseLevel > 3)
461 G4cout << "Calling SampleSecondaries() of G4MicroElecElasticModel" << G4endl;
462
463 G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
464
465 if (electronEnergy0 < killBelowEnergy)
466 {
470 return;
471 }
472
473 if (electronEnergy0 < highEnergyLimit)
474 {
475 G4double cosTheta = 0;
476 if (acousticModelEnabled)
477 {
478 cosTheta = 1 - 2 * G4UniformRand(); //Isotrope
479 }
480 else if (electronEnergy0 >= lowEnergyLimit)
481 {
482 cosTheta = RandomizeCosTheta(electronEnergy0);
483 }
484
485 G4double phi = 2. * pi * G4UniformRand();
486
487 G4ThreeVector zVers = aDynamicElectron->GetMomentumDirection();
488 G4ThreeVector xVers = zVers.orthogonal();
489 G4ThreeVector yVers = zVers.cross(xVers);
490
491 G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
492 G4double yDir = xDir;
493 xDir *= std::cos(phi);
494 yDir *= std::sin(phi);
495
496 G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
497
500 }
501}
502
503//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
504
506{
507 //.................. T in eV!!!!!!!!!!!!!
508 G4double Z2= Z;
509 G4double M2= A;
510 G4double k_d;
511 G4double epsilon_d;
512 G4double g_epsilon_d;
513 G4double E_nu;
514
515 k_d=0.1334*std::pow(Z2,(2./3.))*std::pow(M2,(-1./2.));
516 epsilon_d=0.01014*std::pow(Z2,(-7./3.))*(T/eV);
517 g_epsilon_d= epsilon_d+0.40244*std::pow(epsilon_d,(3./4.))+3.4008*std::pow(epsilon_d,(1./6.));
518
519 E_nu=1./(1.+ k_d*g_epsilon_d);
520
521 return E_nu;
522}
523
524//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
525
526G4double G4MicroElecElasticModel_new::Theta
527 (G4ParticleDefinition * particleDefinition, G4double k, G4double integrDiff)
528{
529
530 G4double theta = 0.;
531 G4double valueT1 = 0;
532 G4double valueT2 = 0;
533 G4double valueE21 = 0;
534 G4double valueE22 = 0;
535 G4double valueE12 = 0;
536 G4double valueE11 = 0;
537 G4double xs11 = 0;
538 G4double xs12 = 0;
539 G4double xs21 = 0;
540 G4double xs22 = 0;
541
542 if (particleDefinition == G4Electron::ElectronDefinition())
543 {
544 ThetaMap::iterator iterator_angle;
545 iterator_angle = thetaDataStorage.find(currentMaterialName);
546
547 energyMap::iterator iterator_energy;
548 iterator_energy = eIncidentEnergyStorage.find(currentMaterialName);
549
550 ProbaMap::iterator iterator_proba;
551 iterator_proba = eProbaStorage.find(currentMaterialName);
552
553 if (iterator_angle != thetaDataStorage.end() && iterator_energy != eIncidentEnergyStorage.end() && iterator_proba != eProbaStorage.end())
554 {
555 TriDimensionMap* eDiffCrossSectionData = iterator_angle->second; //Theta points
556 std::vector<G4double>* eTdummyVec = iterator_energy->second;
557 VecMap* eVecm = iterator_proba->second;
558
559 auto t2 = std::upper_bound(eTdummyVec->begin(), eTdummyVec->end(), k);
560 auto t1 = t2 - 1;
561 auto e12 = std::upper_bound((*eVecm)[(*t1)].begin(), (*eVecm)[(*t1)].end(), integrDiff);
562 auto e11 = e12 - 1;
563 auto e22 = std::upper_bound((*eVecm)[(*t2)].begin(), (*eVecm)[(*t2)].end(), integrDiff);
564 auto e21 = e22 - 1;
565
566 valueT1 = *t1;
567 valueT2 = *t2;
568 valueE21 = *e21;
569 valueE22 = *e22;
570 valueE12 = *e12;
571 valueE11 = *e11;
572
573 xs11 = (*eDiffCrossSectionData)[valueT1][valueE11];
574 xs12 = (*eDiffCrossSectionData)[valueT1][valueE12];
575 xs21 = (*eDiffCrossSectionData)[valueT2][valueE21];
576 xs22 = (*eDiffCrossSectionData)[valueT2][valueE22];
577 }
578 else
579 {
580 G4String str = "Material ";
581 str += currentMaterialName + " not found!";
582 G4Exception("G4MicroElecElasticModel_new::ComputeCrossSectionPerVolume", "em0002", FatalException, str);
583 }
584
585 }
586
587 if (xs11==0 || xs12==0 ||xs21==0 ||xs22==0) return (0.);
588
589 theta = QuadInterpolator( valueE11, valueE12,
590 valueE21, valueE22,
591 xs11, xs12,
592 xs21, xs22,
593 valueT1, valueT2,
594 k, integrDiff );
595
596 return theta;
597}
598
599//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
600
601G4double G4MicroElecElasticModel_new::LinLogInterpolate(G4double e1,
602 G4double e2,
603 G4double e,
604 G4double xs1,
605 G4double xs2)
606{
607 G4double d1 = std::log(xs1);
608 G4double d2 = std::log(xs2);
609 G4double value = G4Exp(d1 + (d2 - d1)*(e - e1)/ (e2 - e1));
610 return value;
611}
612
613//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
614
615G4double G4MicroElecElasticModel_new::LinLinInterpolate(G4double e1,
616 G4double e2,
617 G4double e,
618 G4double xs1,
619 G4double xs2)
620{
621 G4double d1 = xs1;
622 G4double d2 = xs2;
623 G4double value = (d1 + (d2 - d1)*(e - e1)/ (e2 - e1));
624 return value;
625}
626
627//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
628
629G4double G4MicroElecElasticModel_new::LogLogInterpolate(G4double e1,
630 G4double e2,
631 G4double e,
632 G4double xs1,
633 G4double xs2)
634{
635 G4double a = (std::log10(xs2)-std::log10(xs1)) / (std::log10(e2)-std::log10(e1));
636 G4double b = std::log10(xs2) - a*std::log10(e2);
637 G4double sigma = a*std::log10(e) + b;
638 G4double value = (std::pow(10.,sigma));
639 return value;
640}
641
642//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
643
644G4double G4MicroElecElasticModel_new::QuadInterpolator(G4double e11, G4double e12,
645 G4double e21, G4double e22,
646 G4double xs11, G4double xs12,
647 G4double xs21, G4double xs22,
648 G4double t1, G4double t2,
649 G4double t, G4double e)
650{
651
652
653 // Lin-Lin
654 G4double interpolatedvalue1 = LinLinInterpolate(e11, e12, e, xs11, xs12);
655 G4double interpolatedvalue2 = LinLinInterpolate(e21, e22, e, xs21, xs22);
656 G4double value = LinLinInterpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2);
657
658 return value;
659}
660
661//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
662
663G4double G4MicroElecElasticModel_new::RandomizeCosTheta(G4double k)
664{
665 G4double integrdiff=0;
666 G4double uniformRand=G4UniformRand();
667 integrdiff = uniformRand;
668
669 G4double theta=0.;
670 G4double cosTheta=0.;
671 theta = Theta(G4Electron::ElectronDefinition(),k/eV,integrdiff);
672
673 cosTheta= std::cos(theta*pi/180.);
674
675 return cosTheta;
676}
677
678//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
679
680
682{
683 killBelowEnergy = threshold;
684
685 if (threshold < 5*CLHEP::eV)
686 {
687 G4Exception ("*** WARNING : the G4MicroElecElasticModel class is not validated below 5 eV !","",JustWarning,"") ;
688 threshold = 5*CLHEP::eV;
689 }
690}
691
692//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
G4double D(G4double temp)
const char * G4FindDataDir(const char *)
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector unit() const
Hep3Vector orthogonal() const
Hep3Vector cross(const Hep3Vector &) const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:88
const G4Material * GetMaterial() const
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:204
const G4String & GetName() const
Definition: G4Material.hh:172
G4bool LoadData(const G4String &argFileName) override
G4double FindValue(G4double e, G4int componentId=0) const override
G4double AcousticCrossSectionPerVolume(G4double ekin, G4double kbz, G4double rho, G4double cs, G4double Aac, G4double Eac, G4double prefactor)
void SetKillBelowThreshold(G4double threshold)
G4MicroElecElasticModel_new(const G4ParticleDefinition *p=0, const G4String &nam="MicroElecElasticModel")
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4double DamageEnergy(G4double T, G4double A, G4double Z)
G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
G4ParticleChangeForGamma * fParticleChangeForGamma
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
const G4String & GetParticleName() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
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)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
#define DBL_MAX
Definition: templates.hh:62