Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4MicroElecInelasticModel_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// G4MicroElecInelasticModel_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......
75
76
77#include "globals.hh"
80#include "G4SystemOfUnits.hh"
81#include "G4ios.hh"
82#include "G4UnitsTable.hh"
84#include "G4LossTableManager.hh"
87#include "G4DeltaAngle.hh"
88
89#include <sstream>
90
91//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
92
93using namespace std;
94
95//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
96
98 const G4ParticleDefinition*, const G4String& nam)
99 :G4VEmModel(nam),isInitialised(false)
100{
101
102 verboseLevel= 0;
103 // Verbosity scale:
104 // 0 = nothing
105 // 1 = warning for energy non-conservation
106 // 2 = details of energy budget
107 // 3 = calculation of cross sections, file openings, sampling of atoms
108 // 4 = entering in methods
109
110 if( verboseLevel>0 )
111 {
112 G4cout << "MicroElec inelastic model is constructed " << G4endl;
113 }
114
115 //Mark this model as "applicable" for atomic deexcitation
117 fAtomDeexcitation = nullptr;
118 fParticleChangeForGamma = nullptr;
119
120 // default generator
122
123 // Selection of computation method
124 fasterCode = true;
125 SEFromFermiLevel = false;
126}
127
128//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
129
131{
132 // Cross section
133 // (0)
134 TCSMap::iterator pos2;
135 for (pos2 = tableTCS.begin(); pos2 != tableTCS.end(); ++pos2) {
136 MapData* tableData = pos2->second;
137 for (auto pos = tableData->begin(); pos != tableData->end(); ++pos)
138 {
139 G4MicroElecCrossSectionDataSet_new* table = pos->second;
140 delete table;
141 }
142 delete tableData;
143 }
144 tableTCS.clear();
145
146 dataDiffCSMap::iterator iterator_proba;
147 // (1)
148 for (iterator_proba = eNrjTransStorage.begin(); iterator_proba != eNrjTransStorage.end(); ++iterator_proba) {
149 vector<TriDimensionMap>* eNrjTransfData = iterator_proba->second;
150 eNrjTransfData->clear();
151 delete eNrjTransfData;
152 }
153 eNrjTransStorage.clear();
154
155 for (iterator_proba = pNrjTransStorage.begin(); iterator_proba != pNrjTransStorage.end(); ++iterator_proba) {
156 vector<TriDimensionMap>* pNrjTransfData = iterator_proba->second;
157 pNrjTransfData->clear();
158 delete pNrjTransfData;
159 }
160 pNrjTransStorage.clear();
161
162 // (2)
163 for (iterator_proba = eDiffDatatable.begin(); iterator_proba != eDiffDatatable.end(); ++iterator_proba) {
164 vector<TriDimensionMap>* eDiffCrossSectionData = iterator_proba->second;
165 eDiffCrossSectionData->clear();
166 delete eDiffCrossSectionData;
167 }
168 eDiffDatatable.clear();
169
170 for (iterator_proba = pDiffDatatable.begin(); iterator_proba != pDiffDatatable.end(); ++iterator_proba) {
171 vector<TriDimensionMap>* pDiffCrossSectionData = iterator_proba->second;
172 pDiffCrossSectionData->clear();
173 delete pDiffCrossSectionData;
174 }
175 pDiffDatatable.clear();
176
177 // (3)
178 dataProbaShellMap::iterator iterator_probaShell;
179
180 for (iterator_probaShell = eProbaShellStorage.begin(); iterator_probaShell != eProbaShellStorage.end(); ++iterator_probaShell) {
181 vector<VecMap>* eProbaShellMap = iterator_probaShell->second;
182 eProbaShellMap->clear();
183 delete eProbaShellMap;
184 }
185 eProbaShellStorage.clear();
186
187 for (iterator_probaShell = pProbaShellStorage.begin(); iterator_probaShell != pProbaShellStorage.end(); ++iterator_probaShell) {
188 vector<VecMap>* pProbaShellMap = iterator_probaShell->second;
189 pProbaShellMap->clear();
190 delete pProbaShellMap;
191 }
192 pProbaShellStorage.clear();
193
194 // (4)
195 TranfEnergyMap::iterator iterator_nrjtransf;
196 for (iterator_nrjtransf = eVecmStorage.begin(); iterator_nrjtransf != eVecmStorage.end(); ++iterator_nrjtransf) {
197 VecMap* eVecm = iterator_nrjtransf->second;
198 eVecm->clear();
199 delete eVecm;
200 }
201 eVecmStorage.clear();
202 for (iterator_nrjtransf = pVecmStorage.begin(); iterator_nrjtransf != pVecmStorage.end(); ++iterator_nrjtransf) {
203 VecMap* pVecm = iterator_nrjtransf->second;
204 pVecm->clear();
205 delete pVecm;
206 }
207 pVecmStorage.clear();
208
209 // (5)
210 incidentEnergyMap::iterator iterator_energy;
211 for (iterator_energy = eIncidentEnergyStorage.begin(); iterator_energy != eIncidentEnergyStorage.end(); ++iterator_energy) {
212 std::vector<G4double>* eTdummyVec = iterator_energy->second;
213 eTdummyVec->clear();
214 delete eTdummyVec;
215 }
216 eIncidentEnergyStorage.clear();
217
218 for (iterator_energy = pIncidentEnergyStorage.begin(); iterator_energy != pIncidentEnergyStorage.end(); ++iterator_energy) {
219 std::vector<G4double>* pTdummyVec = iterator_energy->second;
220 pTdummyVec->clear();
221 delete pTdummyVec;
222 }
223 pIncidentEnergyStorage.clear();
224
225 // (6)
226 MapStructure::iterator iterator_matStructure;
227 for (iterator_matStructure = tableMaterialsStructures.begin();
228 iterator_matStructure != tableMaterialsStructures.end(); ++iterator_matStructure) {
229 currentMaterialStructure = iterator_matStructure->second;
230 delete currentMaterialStructure;
231 }
232 tableMaterialsStructures.clear();
233 currentMaterialStructure = nullptr;
234}
235
236//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
237
239 const G4DataVector& /*cuts*/)
240{
241 if (isInitialised) { return; }
242
243 if (verboseLevel > 3)
244 G4cout << "Calling G4MicroElecInelasticModel_new::Initialise()" << G4endl;
245
246 const char* path = G4FindDataDir("G4LEDATA");
247 if (!path)
248 {
249 G4Exception("G4MicroElecElasticModel_new::Initialise","em0006",FatalException,"G4LEDATA environment variable not set.");
250 return;
251 }
252
253 G4String modelName = "mermin";
254 G4cout << "****************************" << G4endl;
255 G4cout << modelName << " model loaded !" << G4endl;
256
257 // Energy limits
260 G4String electron = electronDef->GetParticleName();
261 G4String proton = protonDef->GetParticleName();
262
263 G4double scaleFactor = 1.0;
264
265 // *** ELECTRON
266 lowEnergyLimit[electron] = 2 * eV;
267 highEnergyLimit[electron] = 10.0 * MeV;
268
269 // Cross section
271 G4int numOfCouples = (G4int)theCoupleTable->GetTableSize();
272
273 for (G4int i = 0; i < numOfCouples; ++i) {
274 const G4Material* material = theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
275 G4cout << "Material " << i + 1 << " / " << numOfCouples << " : " << material->GetName() << G4endl;
276 if (material->GetName() == "Vacuum") continue;
277 G4String mat = material->GetName().substr(3, material->GetName().size());
278 MapData* tableData = new MapData;
279 currentMaterialStructure = new G4MicroElecMaterialStructure(mat);
280
281 tableMaterialsStructures[mat] = currentMaterialStructure;
282 if (particle == electronDef) {
283 //TCS
284 G4String fileElectron("Inelastic/" + modelName + "_sigma_inelastic_e-_" + mat);
285 G4cout << fileElectron << G4endl;
287 tableE->LoadData(fileElectron);
288 tableData->insert(make_pair(electron, tableE));
289
290 // DCS
291 std::ostringstream eFullFileName;
292 if (fasterCode) {
293 eFullFileName << path << "/microelec/Inelastic/cumulated_" + modelName + "_sigmadiff_inelastic_e-_" + mat + ".dat";
294 G4cout << "Faster code = true" << G4endl;
295 G4cout << "Inelastic/cumulated_" + modelName + "_sigmadiff_inelastic_e-_" + mat + ".dat" << G4endl;
296 }
297 else {
298 eFullFileName << path << "/microelec/Inelastic/" + modelName + "_sigmadiff_inelastic_e-_" + mat + ".dat";
299 G4cout << "Faster code = false" << G4endl;
300 G4cout << "Inelastic/" + modelName + "_sigmadiff_inelastic_e-_" + mat + ".dat" << G4endl;
301 }
302
303 std::ifstream eDiffCrossSection(eFullFileName.str().c_str());
304 if (!eDiffCrossSection)
305 {
306 std::stringstream ss;
307 ss << "Missing data " << eFullFileName.str().c_str();
308 std::string sortieString = ss.str();
309
310 if (fasterCode) G4Exception("G4MicroElecInelasticModel_new::Initialise", "em0003",
311 FatalException, sortieString.c_str());
312
313 else {
314 G4Exception("G4MicroElecInelasticModel_new::Initialise", "em0003",
315 FatalException, "Missing data file:/microelec/sigmadiff_inelastic_e_Si.dat");
316 }
317 }
318
319 // Clear the arrays for re-initialization case (MT mode)
320 // Octobre 22nd, 2014 - Melanie Raine
321 //Creating vectors of maps for DCS and Cumulated DCS for the current material.
322 //Each vector is storing one map for each shell.
323 vector<TriDimensionMap>* eDiffCrossSectionData =
324 new vector<TriDimensionMap>; //Storage of [IncidentEnergy, TransfEnergy, DCS values], used in slower code
325 vector<TriDimensionMap>* eNrjTransfData =
326 new vector<TriDimensionMap>; //Storage of possible transfer energies by shell
327 vector<VecMap>* eProbaShellMap = new vector<VecMap>; //Storage of the vectors containing all cumulated DCS values for an initial energy, by shell
328 vector<G4double>* eTdummyVec = new vector<G4double>; //Storage of incident energies for interpolation
329 VecMap* eVecm = new VecMap; //Transfered energy map for slower code
330
331 for (int j = 0; j < currentMaterialStructure->NumberOfLevels(); j++) //Filling the map vectors with an empty map for each shell
332 {
333 eDiffCrossSectionData->push_back(TriDimensionMap());
334 eNrjTransfData->push_back(TriDimensionMap());
335 eProbaShellMap->push_back(VecMap());
336 }
337
338 eTdummyVec->push_back(0.);
339 while (!eDiffCrossSection.eof())
340 {
341 G4double tDummy; //incident energy
342 G4double eDummy; //transfered energy
343 eDiffCrossSection >> tDummy >> eDummy;
344 if (tDummy != eTdummyVec->back()) eTdummyVec->push_back(tDummy);
345
346 G4double tmp; //probability
347 for (int j = 0; j < currentMaterialStructure->NumberOfLevels(); j++)
348 {
349 eDiffCrossSection >> tmp;
350 (*eDiffCrossSectionData)[j][tDummy][eDummy] = tmp;
351
352 if (fasterCode)
353 {
354 (*eNrjTransfData)[j][tDummy][(*eDiffCrossSectionData)[j][tDummy][eDummy]] = eDummy;
355 (*eProbaShellMap)[j][tDummy].push_back((*eDiffCrossSectionData)[j][tDummy][eDummy]);
356 }
357 else { // SI - only if eof is not reached !
358 if (!eDiffCrossSection.eof()) (*eDiffCrossSectionData)[j][tDummy][eDummy] *= scaleFactor;
359 (*eVecm)[tDummy].push_back(eDummy);
360 }
361 }
362 }
363 //
364 G4cout << "add to material vector" << G4endl;
365
366 //Filing maps for the current material into the master maps
367 if (fasterCode) {
368 eNrjTransStorage[mat] = eNrjTransfData;
369 eProbaShellStorage[mat] = eProbaShellMap;
370 }
371 else {
372 eDiffDatatable[mat] = eDiffCrossSectionData;
373 eVecmStorage[mat] = eVecm;
374 }
375 eIncidentEnergyStorage[mat] = eTdummyVec;
376
377 //Cleanup support vectors
378 // delete eProbaShellMap;
379 // delete eDiffCrossSectionData;
380 // delete eNrjTransfData;
381 }
382
383 // *** PROTON
384 if (particle == protonDef)
385 {
386 // Cross section
387 G4String fileProton("Inelastic/" + modelName + "_sigma_inelastic_p_" + mat); G4cout << fileProton << G4endl;
389 tableP->LoadData(fileProton);
390 tableData->insert(make_pair(proton, tableP));
391
392 // DCS
393 std::ostringstream pFullFileName;
394 if (fasterCode) {
395 pFullFileName << path << "/microelec/Inelastic/cumulated_" + modelName + "_sigmadiff_inelastic_p_" + mat + ".dat";
396 G4cout << "Faster code = true" << G4endl;
397 G4cout << "Inelastic/cumulated_" + modelName + "_sigmadiff_inelastic_p_" + mat + ".dat" << G4endl;
398 }
399 else {
400 pFullFileName << path << "/microelec/Inelastic/" + modelName + "_sigmadiff_inelastic_p_" + mat + ".dat";
401 G4cout << "Faster code = false" << G4endl;
402 G4cout << "Inelastic/" + modelName + "_sigmadiff_inelastic_e-_" + mat + ".dat" << G4endl;
403 }
404
405 std::ifstream pDiffCrossSection(pFullFileName.str().c_str());
406 if (!pDiffCrossSection)
407 {
408 if (fasterCode) G4Exception("G4MicroElecInelasticModel_new::Initialise", "em0003",
409 FatalException, "Missing data file:/microelec/sigmadiff_cumulated_inelastic_p_Si.dat");
410 else {
411 G4Exception("G4MicroElecInelasticModel_new::Initialise", "em0003",
412 FatalException, "Missing data file:/microelec/sigmadiff_inelastic_p_Si.dat");
413 }
414 }
415
416 //
417 // Clear the arrays for re-initialization case (MT mode)
418 // Octobre 22nd, 2014 - Melanie Raine
419 //Creating vectors of maps for DCS and Cumulated DCS for the current material.
420 //Each vector is storing one map for each shell.
421
422 vector<TriDimensionMap>* pDiffCrossSectionData =
423 new vector<TriDimensionMap>; //Storage of [IncidentEnergy, TransfEnergy, DCS values], used in slower code
424 vector<TriDimensionMap>* pNrjTransfData =
425 new vector<TriDimensionMap>; //Storage of possible transfer energies by shell
426 vector<VecMap>* pProbaShellMap =
427 new vector<VecMap>; //Storage of the vectors containing all cumulated DCS values for an initial energy, by shell
428 vector<G4double>* pTdummyVec =
429 new vector<G4double>; //Storage of incident energies for interpolation
430 VecMap* eVecm = new VecMap; //Transfered energy map for slower code
431
432 for (int j = 0; j < currentMaterialStructure->NumberOfLevels(); ++j)
433 //Filling the map vectors with an empty map for each shell
434 {
435 pDiffCrossSectionData->push_back(TriDimensionMap());
436 pNrjTransfData->push_back(TriDimensionMap());
437 pProbaShellMap->push_back(VecMap());
438 }
439
440 pTdummyVec->push_back(0.);
441 while (!pDiffCrossSection.eof())
442 {
443 G4double tDummy; //incident energy
444 G4double eDummy; //transfered energy
445 pDiffCrossSection >> tDummy >> eDummy;
446 if (tDummy != pTdummyVec->back()) pTdummyVec->push_back(tDummy);
447
448 G4double tmp; //probability
449 for (int j = 0; j < currentMaterialStructure->NumberOfLevels(); j++)
450 {
451 pDiffCrossSection >> tmp;
452 (*pDiffCrossSectionData)[j][tDummy][eDummy] = tmp;
453 // ArrayofMaps[j] -> fill with 3DMap(incidentEnergy,
454 // 2Dmap (transferedEnergy,proba=tmp) ) -> fill map for shell j
455 // with proba for transfered energy eDummy
456
457 if (fasterCode)
458 {
459 (*pNrjTransfData)[j][tDummy][(*pDiffCrossSectionData)[j][tDummy][eDummy]] = eDummy;
460 (*pProbaShellMap)[j][tDummy].push_back((*pDiffCrossSectionData)[j][tDummy][eDummy]);
461 }
462 else { // SI - only if eof is not reached !
463 if (!pDiffCrossSection.eof()) (*pDiffCrossSectionData)[j][tDummy][eDummy] *= scaleFactor;
464 (*eVecm)[tDummy].push_back(eDummy);
465 }
466 }
467 }
468
469 //Filing maps for the current material into the master maps
470 if (fasterCode) {
471 pNrjTransStorage[mat] = pNrjTransfData;
472 pProbaShellStorage[mat] = pProbaShellMap;
473 }
474 else {
475 pDiffDatatable[mat] = pDiffCrossSectionData;
476 pVecmStorage[mat] = eVecm;
477 }
478 pIncidentEnergyStorage[mat] = pTdummyVec;
479
480 //Cleanup support vectors
481 // delete pNrjTransfData;
482 // delete eVecm;
483 // delete pDiffCrossSectionData;
484 // delete pProbaShellMap;
485 }
486 tableTCS[mat] = tableData;
487}
488 if (particle==electronDef)
489 {
490 SetLowEnergyLimit(lowEnergyLimit[electron]);
491 SetHighEnergyLimit(highEnergyLimit[electron]);
492 }
493
494 if (particle==protonDef)
495 {
496 SetLowEnergyLimit(100*eV);
497 SetHighEnergyLimit(10*MeV);
498 }
499
500 if( verboseLevel>1 )
501 {
502 G4cout << "MicroElec Inelastic model is initialized " << G4endl
503 << "Energy range: "
504 << LowEnergyLimit() / keV << " keV - "
505 << HighEnergyLimit() / MeV << " MeV for "
506 << particle->GetParticleName()
507 << " with mass (amu) " << particle->GetPDGMass()/proton_mass_c2
508 << " and charge " << particle->GetPDGCharge()
509 << G4endl << G4endl ;
510 }
511
512 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
513
515 isInitialised = true;
516}
517
518//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
519
521 const G4ParticleDefinition* particleDefinition,
522 G4double ekin,
523 G4double,
524 G4double)
525{
526 if (verboseLevel > 3) G4cout << "Calling CrossSectionPerVolume() of G4MicroElecInelasticModel" << G4endl;
527
528 G4double density = material->GetTotNbOfAtomsPerVolume();
529 currentMaterial = material->GetName().substr(3, material->GetName().size());
530
531 MapStructure::iterator structPos;
532 structPos = tableMaterialsStructures.find(currentMaterial);
533
534 // Calculate total cross section for model
535 TCSMap::iterator tablepos;
536 tablepos = tableTCS.find(currentMaterial);
537
538 if (tablepos == tableTCS.end() )
539 {
540 G4String str = "Material ";
541 str += currentMaterial + " TCS Table not found!";
542 G4Exception("G4MicroElecInelasticModel_new::ComputeCrossSectionPerVolume", "em0002", FatalException, str);
543 return 0;
544 }
545 else if(structPos == tableMaterialsStructures.end())
546 {
547 G4String str = "Material ";
548 str += currentMaterial + " Structure not found!";
549 G4Exception("G4MicroElecInelasticModel_new::ComputeCrossSectionPerVolume", "em0002", FatalException, str);
550 return 0;
551 }
552 else {
553 MapData* tableData = tablepos->second;
554 currentMaterialStructure = structPos->second;
555
556 G4double sigma = 0;
557
558 const G4String& particleName = particleDefinition->GetParticleName();
559 G4String nameLocal = particleName;
560 G4int pdg = particleDefinition->GetPDGEncoding();
561 G4int Z = particleDefinition->GetAtomicNumber();
562
563 G4double Zeff = 1.0, Zeff2 = Zeff*Zeff;
564 G4double Mion_c2 = particleDefinition->GetPDGMass();
565
566 if (Mion_c2 > proton_mass_c2)
567 {
568 ekin *= proton_mass_c2 / Mion_c2;
569 nameLocal = "proton";
570 }
571
572 G4double lowLim = currentMaterialStructure->GetInelasticModelLowLimit(pdg);
573 G4double highLim = currentMaterialStructure->GetInelasticModelHighLimit(pdg);
574
575 if (ekin >= lowLim && ekin < highLim)
576 {
577 std::map< G4String, G4MicroElecCrossSectionDataSet_new*, std::less<G4String> >::iterator pos;
578 pos = tableData->find(nameLocal); //find particle type
579
580 if (pos != tableData->end())
581 {
582 G4MicroElecCrossSectionDataSet_new* table = pos->second;
583 if (table != 0)
584 {
585 sigma = table->FindValue(ekin);
586
587 if (Mion_c2 > proton_mass_c2) {
588 sigma = 0.;
589 for (G4int i = 0; i < currentMaterialStructure->NumberOfLevels(); i++) {
590 Zeff = BKZ(ekin / (proton_mass_c2 / Mion_c2), Mion_c2 / c_squared, Z, currentMaterialStructure->Energy(i)); // il faut garder le vrai ekin car le calcul à l'interieur de la methode convertie l'énergie en vitesse
591 Zeff2 = Zeff*Zeff;
592 sigma += Zeff2*table->FindShellValue(ekin, i);
593// il faut utiliser le ekin mis à l'echelle pour chercher la bonne
594// valeur dans les tables proton
595
596 }
597 }
598 else {
599 sigma = table->FindValue(ekin);
600 }
601 }
602 }
603 else
604 {
605 G4Exception("G4MicroElecInelasticModel_new::CrossSectionPerVolume",
606 "em0002", FatalException,
607 "Model not applicable to particle type.");
608 }
609 }
610 else
611 {
612 return 1 / DBL_MAX;
613 }
614
615 if (verboseLevel > 3)
616 {
617 G4cout << "---> Kinetic energy (eV)=" << ekin / eV << G4endl;
618 G4cout << " - Cross section per Si atom (cm^2)=" << sigma / cm2 << G4endl;
619 G4cout << " - Cross section per Si atom (cm^-1)=" << sigma*density / (1. / cm) << G4endl;
620 }
621
622 return (sigma)*density;}
623
624}
625
626//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
627
628void G4MicroElecInelasticModel_new::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
629 const G4MaterialCutsCouple* couple,
630 const G4DynamicParticle* particle,
631 G4double,
632 G4double)
633{
634
635 if (verboseLevel > 3)
636 G4cout << "Calling SampleSecondaries() of G4MicroElecInelasticModel" << G4endl;
637
638 G4int pdg = particle->GetParticleDefinition()->GetPDGEncoding();
639 G4double lowLim = currentMaterialStructure->GetInelasticModelLowLimit(pdg);
640 G4double highLim = currentMaterialStructure->GetInelasticModelHighLimit(pdg);
641
642 G4double ekin = particle->GetKineticEnergy();
643 G4double k = ekin ;
644
645 G4ParticleDefinition* PartDef = particle->GetDefinition();
646 const G4String& particleName = PartDef->GetParticleName();
647 G4String nameLocal2 = particleName ;
648 G4double particleMass = particle->GetDefinition()->GetPDGMass();
649 G4double originalMass = particleMass; // a passer en argument dans samplesecondaryenergy pour évaluer correctement Qmax
650 G4int originalZ = particle->GetDefinition()->GetAtomicNumber();
651
652 if (particleMass > proton_mass_c2)
653 {
654 k *= proton_mass_c2/particleMass ;
655 PartDef = G4Proton::ProtonDefinition();
656 nameLocal2 = "proton" ;
657 }
658
659 if (k >= lowLim && k < highLim)
660 {
661 G4ParticleMomentum primaryDirection = particle->GetMomentumDirection();
662 G4double totalEnergy = ekin + particleMass;
663 G4double pSquare = ekin * (totalEnergy + particleMass);
664 G4double totalMomentum = std::sqrt(pSquare);
665
666 G4int Shell = 1;
667
668 Shell = RandomSelect(k,nameLocal2,originalMass, originalZ);
669
670 G4double bindingEnergy = currentMaterialStructure->Energy(Shell);
671 G4double limitEnergy = currentMaterialStructure->GetLimitEnergy(Shell);
672
673 if (verboseLevel > 3)
674 {
675 G4cout << "---> Kinetic energy (eV)=" << k/eV << G4endl ;
676 G4cout << "Shell: " << Shell << ", energy: " << bindingEnergy/eV << G4endl;
677 }
678
679 // sample deexcitation
680
681 std::size_t secNumberInit = 0; // need to know at a certain point the energy of secondaries
682 std::size_t secNumberFinal = 0; // So I'll make the difference and then sum the energies
683
684 //SI: additional protection if tcs interpolation method is modified
685 //if (k<bindingEnergy) return;
686 if (k<limitEnergy) return;
687 // G4cout << currentMaterial << G4endl;
688 G4int Z = currentMaterialStructure->GetZ(Shell);
689 G4int shellEnum = currentMaterialStructure->GetEADL_Enumerator(Shell);
690 if (currentMaterialStructure->IsShellWeaklyBound(Shell)) { shellEnum = -1; }
691
692 if(fAtomDeexcitation && shellEnum >=0)
693 {
694 // G4cout << "enter if deex and shell 0" << G4endl;
696 const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);
697 secNumberInit = fvect->size();
698 fAtomDeexcitation->GenerateParticles(fvect, shell, Z, 0, 0);
699 secNumberFinal = fvect->size();
700 }
701
702 G4double secondaryKinetic=-1000*eV;
703 SEFromFermiLevel = false;
704 if (!fasterCode)
705 {
706 secondaryKinetic = RandomizeEjectedElectronEnergy(PartDef, k, Shell, originalMass, originalZ);
707 }
708 else
709 {
710 secondaryKinetic = RandomizeEjectedElectronEnergyFromCumulatedDcs(PartDef, k, Shell) ;
711 }
712
713 if (verboseLevel > 3)
714 {
715 G4cout << "Ionisation process" << G4endl;
716 G4cout << "Shell: " << Shell << " Kin. energy (eV)=" << k/eV
717 << " Sec. energy (eV)=" << secondaryKinetic/eV << G4endl;
718 }
719 G4ThreeVector deltaDirection =
720 GetAngularDistribution()->SampleDirectionForShell(particle, secondaryKinetic,
721 Z, Shell,
722 couple->GetMaterial());
723
725 {
726 G4double deltaTotalMomentum = std::sqrt(secondaryKinetic*(secondaryKinetic + 2.*electron_mass_c2 ));
727
728 G4double finalPx = totalMomentum*primaryDirection.x() - deltaTotalMomentum*deltaDirection.x();
729 G4double finalPy = totalMomentum*primaryDirection.y() - deltaTotalMomentum*deltaDirection.y();
730 G4double finalPz = totalMomentum*primaryDirection.z() - deltaTotalMomentum*deltaDirection.z();
731 G4double finalMomentum = std::sqrt(finalPx*finalPx + finalPy*finalPy + finalPz*finalPz);
732 finalPx /= finalMomentum;
733 finalPy /= finalMomentum;
734 finalPz /= finalMomentum;
735
736 G4ThreeVector direction;
737 direction.set(finalPx,finalPy,finalPz);
738
740 }
742
743 // note that secondaryKinetic is the energy of the delta ray, not of all secondaries.
744 G4double deexSecEnergy = 0;
745 for (std::size_t j=secNumberInit; j < secNumberFinal; ++j) {
746 deexSecEnergy = deexSecEnergy + (*fvect)[j]->GetKineticEnergy();
747 }
748 if (SEFromFermiLevel) limitEnergy = currentMaterialStructure->GetEnergyGap();
749 fParticleChangeForGamma->SetProposedKineticEnergy(ekin - secondaryKinetic - limitEnergy); //Ef = Ei-(Q-El)-El = Ei-Q
750 fParticleChangeForGamma->ProposeLocalEnergyDeposit(limitEnergy - deexSecEnergy);
751
752 if (secondaryKinetic>0)
753 {
754 G4DynamicParticle* dp = new G4DynamicParticle(G4Electron::Electron(), deltaDirection, secondaryKinetic); //Esec = Q-El
755 fvect->push_back(dp);
756 }
757 }
758}
759
760//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
761
762G4double G4MicroElecInelasticModel_new::RandomizeEjectedElectronEnergy(
763 const G4ParticleDefinition* particleDefinition,
764 G4double k, G4int shell, G4double originalMass, G4int)
765{
766 G4double secondaryElectronKineticEnergy=0.;
767 if (particleDefinition == G4Electron::ElectronDefinition())
768 {
769 G4double maximumEnergyTransfer=k;
770 G4double crossSectionMaximum = 0.;
771 G4double minEnergy = currentMaterialStructure->GetLimitEnergy(shell);
772 G4double maxEnergy = maximumEnergyTransfer;
773 G4int nEnergySteps = 100;
774
775 G4double value(minEnergy);
776 G4double stpEnergy(std::pow(maxEnergy/value, 1./static_cast<G4double>(nEnergySteps-1)));
777 G4int step(nEnergySteps);
778 while (step>0)
779 {
780 --step;
781 G4double differentialCrossSection =
782 DifferentialCrossSection(particleDefinition, k, value, shell);
783 crossSectionMaximum = std::max(crossSectionMaximum, differentialCrossSection);
784 value*=stpEnergy;
785 }
786
787 do
788 {
789 secondaryElectronKineticEnergy = G4UniformRand() *
790 (maximumEnergyTransfer-currentMaterialStructure->GetLimitEnergy(shell));
791 } while(G4UniformRand()*crossSectionMaximum >
792 DifferentialCrossSection(particleDefinition, k,
793 (secondaryElectronKineticEnergy+currentMaterialStructure->GetLimitEnergy(shell)),shell));
794 }
795 else if (particleDefinition == G4Proton::ProtonDefinition())
796 {
797 G4double maximumEnergyTransfer =
798 ComputeElasticQmax(k/(proton_mass_c2/originalMass),
799 currentMaterialStructure->Energy(shell),
800 originalMass/c_squared, electron_mass_c2/c_squared);
801
802 G4double crossSectionMaximum = 0.;
803
804 G4double minEnergy = currentMaterialStructure->GetLimitEnergy(shell);
805 G4double maxEnergy = maximumEnergyTransfer;
806 G4int nEnergySteps = 100;
807
808 G4double value(minEnergy);
809 G4double stpEnergy(std::pow(maxEnergy/value, 1./static_cast<G4double>(nEnergySteps-1)));
810 G4int step(nEnergySteps);
811
812 while (step>0)
813 {
814 --step;
815 G4double differentialCrossSection =
816 DifferentialCrossSection(particleDefinition, k, value, shell);
817 crossSectionMaximum = std::max(crossSectionMaximum, differentialCrossSection);
818 value*=stpEnergy;
819 }
820
821 G4double energyTransfer = 0.;
822 do
823 {
824 energyTransfer = G4UniformRand() * maximumEnergyTransfer;
825 } while(G4UniformRand()*crossSectionMaximum >
826 DifferentialCrossSection(particleDefinition, k,energyTransfer,shell));
827
828 secondaryElectronKineticEnergy =
829 energyTransfer-currentMaterialStructure->GetLimitEnergy(shell);
830
831 }
832 return std::max(secondaryElectronKineticEnergy, 0.0);
833}
834
835//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
836
837G4double G4MicroElecInelasticModel_new::RandomizeEjectedElectronEnergyFromCumulatedDcs(
838 const G4ParticleDefinition* particleDefinition, G4double k, G4int shell)
839{
840 G4double secondaryElectronKineticEnergy = 0.;
841 G4double random = G4UniformRand();
842
843 secondaryElectronKineticEnergy = TransferedEnergy(particleDefinition, k, shell, random)
844 - currentMaterialStructure->GetLimitEnergy(shell) ;
845
846 if (isnan(secondaryElectronKineticEnergy)) { secondaryElectronKineticEnergy = k - currentMaterialStructure->GetLimitEnergy(shell); }
847
848 if (secondaryElectronKineticEnergy < 0.) {
849 secondaryElectronKineticEnergy = k - currentMaterialStructure->GetEnergyGap();
850 SEFromFermiLevel = true;
851 }
852 return secondaryElectronKineticEnergy;
853}
854
855//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
856
857G4double G4MicroElecInelasticModel_new::TransferedEnergy(
858 const G4ParticleDefinition* particleDefinition,
859 G4double k,
860 G4int ionizationLevelIndex,
861 G4double random)
862{
863 G4double nrj = 0.;
864 G4double valueK1 = 0;
865 G4double valueK2 = 0;
866 G4double valuePROB21 = 0;
867 G4double valuePROB22 = 0;
868 G4double valuePROB12 = 0;
869 G4double valuePROB11 = 0;
870 G4double nrjTransf11 = 0;
871 G4double nrjTransf12 = 0;
872 G4double nrjTransf21 = 0;
873 G4double nrjTransf22 = 0;
874
875 G4double maximumEnergyTransfer1 = 0;
876 G4double maximumEnergyTransfer2 = 0;
877 G4double maximumEnergyTransferP = 4.* (electron_mass_c2 / proton_mass_c2) * k;
878 G4double bindingEnergy = currentMaterialStructure->GetLimitEnergy(ionizationLevelIndex);
879
880 if (particleDefinition == G4Electron::ElectronDefinition())
881 {
882 dataDiffCSMap::iterator iterator_Nrj;
883 iterator_Nrj = eNrjTransStorage.find(currentMaterial);
884
885 dataProbaShellMap::iterator iterator_Proba;
886 iterator_Proba = eProbaShellStorage.find(currentMaterial);
887
888 incidentEnergyMap::iterator iterator_Tdummy;
889 iterator_Tdummy = eIncidentEnergyStorage.find(currentMaterial);
890
891 if(iterator_Nrj == eNrjTransStorage.end() || iterator_Proba == eProbaShellStorage.end() ||
892 iterator_Tdummy == eIncidentEnergyStorage.end())
893 {
894 G4String str = "Material ";
895 str += currentMaterial + " not found!";
896 G4Exception("G4MicroElecInelasticModel_new::TransferedEnergy", "em0002",
897 FatalException, str);
898 }
899 else {
900 vector<TriDimensionMap>* eNrjTransfData = iterator_Nrj->second; //Storage of possible transfer energies
901 vector<VecMap>* eProbaShellMap = iterator_Proba->second; //Storage of probabilities for energy transfer
902 vector<G4double>* eTdummyVec = iterator_Tdummy->second; //Incident energies for interpolation
903
904 // k should be in eV
905 auto k2 = std::upper_bound(eTdummyVec->begin(),
906 eTdummyVec->end(),
907 k);
908 auto k1 = k2 - 1;
909
910 // SI : the following condition avoids situations where random >last vector element
911 if (random <= (*eProbaShellMap)[ionizationLevelIndex][(*k1)].back()
912 && random <= (*eProbaShellMap)[ionizationLevelIndex][(*k2)].back())
913 {
914 auto prob12 =
915 std::upper_bound((*eProbaShellMap)[ionizationLevelIndex][(*k1)].begin(),
916 (*eProbaShellMap)[ionizationLevelIndex][(*k1)].end(),
917 random);
918
919 auto prob11 = prob12 - 1;
920
921 auto prob22 =
922 std::upper_bound((*eProbaShellMap)[ionizationLevelIndex][(*k2)].begin(),
923 (*eProbaShellMap)[ionizationLevelIndex][(*k2)].end(),
924 random);
925
926 auto prob21 = prob22 - 1;
927
928 valueK1 = *k1;
929 valueK2 = *k2;
930 valuePROB21 = *prob21;
931 valuePROB22 = *prob22;
932 valuePROB12 = *prob12;
933 valuePROB11 = *prob11;
934
935 // The following condition avoid getting transfered energy < binding energy and forces cumxs = 1 for maximum energy transfer.
936 if (valuePROB11 == 0) nrjTransf11 = bindingEnergy;
937 else nrjTransf11 = (*eNrjTransfData)[ionizationLevelIndex][valueK1][valuePROB11];
938 if (valuePROB12 == 1)
939 {
940 if ((valueK1 + bindingEnergy) / 2. > valueK1)
941 maximumEnergyTransfer1 = valueK1;
942 else
943 maximumEnergyTransfer1 = (valueK1 + bindingEnergy) / 2.;
944
945 nrjTransf12 = maximumEnergyTransfer1;
946 }
947 else
948 nrjTransf12 = (*eNrjTransfData)[ionizationLevelIndex][valueK1][valuePROB12];
949
950 if (valuePROB21 == 0) nrjTransf21 = bindingEnergy;
951 else nrjTransf21 = (*eNrjTransfData)[ionizationLevelIndex][valueK2][valuePROB21];
952 if (valuePROB22 == 1)
953 {
954 if ((valueK2 + bindingEnergy) / 2. > valueK2) maximumEnergyTransfer2 = valueK2;
955 else maximumEnergyTransfer2 = (valueK2 + bindingEnergy) / 2.;
956
957 nrjTransf22 = maximumEnergyTransfer2;
958 }
959 else nrjTransf22 = (*eNrjTransfData)[ionizationLevelIndex][valueK2][valuePROB22];
960
961 }
962 // Avoids cases where cum xs is zero for k1 and is not for k2 (with always k1<k2)
963 if (random > (*eProbaShellMap)[ionizationLevelIndex][(*k1)].back())
964 {
965 auto prob22 =
966 std::upper_bound((*eProbaShellMap)[ionizationLevelIndex][(*k2)].begin(),
967 (*eProbaShellMap)[ionizationLevelIndex][(*k2)].end(),
968 random);
969 auto prob21 = prob22 - 1;
970
971 valueK1 = *k1;
972 valueK2 = *k2;
973 valuePROB21 = *prob21;
974 valuePROB22 = *prob22;
975
976 nrjTransf21 = (*eNrjTransfData)[ionizationLevelIndex][valueK2][valuePROB21];
977 nrjTransf22 = (*eNrjTransfData)[ionizationLevelIndex][valueK2][valuePROB22];
978
979 G4double interpolatedvalue2 = Interpolate(valuePROB21,
980 valuePROB22,
981 random,
982 nrjTransf21,
983 nrjTransf22);
984
985 // zeros are explicitly set
986 G4double value = Interpolate(valueK1, valueK2, k, 0., interpolatedvalue2);
987
988 return value;
989 }
990 }
991 }
992 else if (particleDefinition == G4Proton::ProtonDefinition())
993 {
994 // k should be in eV
995 dataDiffCSMap::iterator iterator_Nrj;
996 iterator_Nrj = pNrjTransStorage.find(currentMaterial);
997
998 dataProbaShellMap::iterator iterator_Proba;
999 iterator_Proba = pProbaShellStorage.find(currentMaterial);
1000
1001 incidentEnergyMap::iterator iterator_Tdummy;
1002 iterator_Tdummy = pIncidentEnergyStorage.find(currentMaterial);
1003
1004 if (iterator_Nrj == pNrjTransStorage.end() || iterator_Proba == pProbaShellStorage.end() ||
1005 iterator_Tdummy == pIncidentEnergyStorage.end())
1006 {
1007 G4String str = "Material ";
1008 str += currentMaterial + " not found!";
1009 G4Exception("G4MicroElecInelasticModel_new::TransferedEnergy", "em0002",
1010 FatalException, str);
1011 }
1012 else
1013 {
1014 vector<TriDimensionMap>* pNrjTransfData = iterator_Nrj->second; //Storage of possible transfer energies
1015 vector<VecMap>* pProbaShellMap = iterator_Proba->second; //Storage of probabilities for energy transfer
1016 vector<G4double>* pTdummyVec = iterator_Tdummy->second; //Incident energies for interpolation
1017
1018 auto k2 = std::upper_bound(pTdummyVec->begin(),
1019 pTdummyVec->end(),
1020 k);
1021
1022 auto k1 = k2 - 1;
1023
1024 // SI : the following condition avoids situations where random > last vector element,
1025 // for eg. when the last element is zero
1026 if (random <= (*pProbaShellMap)[ionizationLevelIndex][(*k1)].back()
1027 && random <= (*pProbaShellMap)[ionizationLevelIndex][(*k2)].back())
1028 {
1029 auto prob12 =
1030 std::upper_bound((*pProbaShellMap)[ionizationLevelIndex][(*k1)].begin(),
1031 (*pProbaShellMap)[ionizationLevelIndex][(*k1)].end(),
1032 random);
1033 auto prob11 = prob12 - 1;
1034 auto prob22 =
1035 std::upper_bound((*pProbaShellMap)[ionizationLevelIndex][(*k2)].begin(),
1036 (*pProbaShellMap)[ionizationLevelIndex][(*k2)].end(),
1037 random);
1038 auto prob21 = prob22 - 1;
1039
1040 valueK1 = *k1;
1041 valueK2 = *k2;
1042 valuePROB21 = *prob21;
1043 valuePROB22 = *prob22;
1044 valuePROB12 = *prob12;
1045 valuePROB11 = *prob11;
1046
1047 // The following condition avoid getting transfered energy < binding energy
1048 // and forces cumxs = 1 for maximum energy transfer.
1049 if (valuePROB11 == 0) nrjTransf11 = bindingEnergy;
1050 else nrjTransf11 = (*pNrjTransfData)[ionizationLevelIndex][valueK1][valuePROB11];
1051
1052 if (valuePROB12 == 1) nrjTransf12 = maximumEnergyTransferP;
1053 else nrjTransf12 = (*pNrjTransfData)[ionizationLevelIndex][valueK1][valuePROB12];
1054
1055 if (valuePROB21 == 0) nrjTransf21 = bindingEnergy;
1056 else nrjTransf21 = (*pNrjTransfData)[ionizationLevelIndex][valueK2][valuePROB21];
1057
1058 if (valuePROB22 == 1) nrjTransf22 = maximumEnergyTransferP;
1059 else nrjTransf22 = (*pNrjTransfData)[ionizationLevelIndex][valueK2][valuePROB22];
1060
1061 }
1062
1063 // Avoids cases where cum xs is zero for k1 and is not for k2 (with always k1<k2)
1064 if (random > (*pProbaShellMap)[ionizationLevelIndex][(*k1)].back())
1065 {
1066 auto prob22 =
1067 std::upper_bound((*pProbaShellMap)[ionizationLevelIndex][(*k2)].begin(),
1068 (*pProbaShellMap)[ionizationLevelIndex][(*k2)].end(),
1069 random);
1070
1071 auto prob21 = prob22 - 1;
1072
1073 valueK1 = *k1;
1074 valueK2 = *k2;
1075 valuePROB21 = *prob21;
1076 valuePROB22 = *prob22;
1077
1078 nrjTransf21 = (*pNrjTransfData)[ionizationLevelIndex][valueK2][valuePROB21];
1079 nrjTransf22 = (*pNrjTransfData)[ionizationLevelIndex][valueK2][valuePROB22];
1080
1081 G4double interpolatedvalue2 = Interpolate(valuePROB21,
1082 valuePROB22,
1083 random,
1084 nrjTransf21,
1085 nrjTransf22);
1086
1087 // zeros are explicitly set
1088 G4double value = Interpolate(valueK1, valueK2, k, 0., interpolatedvalue2);
1089 return value;
1090 }
1091 }
1092 }
1093 // End electron and proton cases
1094
1095 G4double nrjTransfProduct = nrjTransf11 * nrjTransf12 * nrjTransf21 * nrjTransf22;
1096
1097 if (nrjTransfProduct != 0.)
1098 {
1099 nrj = QuadInterpolator(valuePROB11,
1100 valuePROB12,
1101 valuePROB21,
1102 valuePROB22,
1103 nrjTransf11,
1104 nrjTransf12,
1105 nrjTransf21,
1106 nrjTransf22,
1107 valueK1,
1108 valueK2,
1109 k,
1110 random);
1111 }
1112
1113 return nrj;
1114}
1115
1116//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1117
1119 const G4ParticleDefinition * particleDefinition,
1120 G4double k,
1121 G4double energyTransfer,
1122 G4int LevelIndex)
1123{
1124 G4double sigma = 0.;
1125
1126 if (energyTransfer >= currentMaterialStructure->GetLimitEnergy(LevelIndex))
1127 {
1128 G4double valueT1 = 0;
1129 G4double valueT2 = 0;
1130 G4double valueE21 = 0;
1131 G4double valueE22 = 0;
1132 G4double valueE12 = 0;
1133 G4double valueE11 = 0;
1134
1135 G4double xs11 = 0;
1136 G4double xs12 = 0;
1137 G4double xs21 = 0;
1138 G4double xs22 = 0;
1139
1140 if (particleDefinition == G4Electron::ElectronDefinition())
1141 {
1142
1143 dataDiffCSMap::iterator iterator_Proba;
1144 iterator_Proba = eDiffDatatable.find(currentMaterial);
1145
1146 incidentEnergyMap::iterator iterator_Nrj;
1147 iterator_Nrj = eIncidentEnergyStorage.find(currentMaterial);
1148
1149 TranfEnergyMap::iterator iterator_TransfNrj;
1150 iterator_TransfNrj = eVecmStorage.find(currentMaterial);
1151
1152 if (iterator_Proba != eDiffDatatable.end() && iterator_Nrj != eIncidentEnergyStorage.end()
1153 && iterator_TransfNrj!= eVecmStorage.end())
1154 {
1155 vector<TriDimensionMap>* eDiffCrossSectionData = (iterator_Proba->second);
1156 vector<G4double>* eTdummyVec = iterator_Nrj->second; //Incident energies for interpolation
1157 VecMap* eVecm = iterator_TransfNrj->second;
1158
1159 // k should be in eV and energy transfer eV also
1160 auto t2 = std::upper_bound(eTdummyVec->begin(), eTdummyVec->end(), k);
1161 auto t1 = t2 - 1;
1162 // SI : the following condition avoids situations where energyTransfer >last vector element
1163 if (energyTransfer <= (*eVecm)[(*t1)].back() && energyTransfer <= (*eVecm)[(*t2)].back())
1164 {
1165 auto e12 = std::upper_bound((*eVecm)[(*t1)].begin(), (*eVecm)[(*t1)].end(), energyTransfer);
1166 auto e11 = e12 - 1;
1167 auto e22 = std::upper_bound((*eVecm)[(*t2)].begin(), (*eVecm)[(*t2)].end(), energyTransfer);
1168 auto e21 = e22 - 1;
1169
1170 valueT1 = *t1;
1171 valueT2 = *t2;
1172 valueE21 = *e21;
1173 valueE22 = *e22;
1174 valueE12 = *e12;
1175 valueE11 = *e11;
1176
1177 xs11 = (*eDiffCrossSectionData)[LevelIndex][valueT1][valueE11];
1178 xs12 = (*eDiffCrossSectionData)[LevelIndex][valueT1][valueE12];
1179 xs21 = (*eDiffCrossSectionData)[LevelIndex][valueT2][valueE21];
1180 xs22 = (*eDiffCrossSectionData)[LevelIndex][valueT2][valueE22];
1181 }
1182 }
1183 else {
1184 G4String str = "Material ";
1185 str += currentMaterial + " not found!";
1186 G4Exception("G4MicroElecDielectricModels::DifferentialCrossSection", "em0002", FatalException, str);
1187 }
1188 }
1189
1190 if (particleDefinition == G4Proton::ProtonDefinition())
1191 {
1192 dataDiffCSMap::iterator iterator_Proba;
1193 iterator_Proba = pDiffDatatable.find(currentMaterial);
1194
1195 incidentEnergyMap::iterator iterator_Nrj;
1196 iterator_Nrj = pIncidentEnergyStorage.find(currentMaterial);
1197
1198 TranfEnergyMap::iterator iterator_TransfNrj;
1199 iterator_TransfNrj = pVecmStorage.find(currentMaterial);
1200
1201 if (iterator_Proba != pDiffDatatable.end() && iterator_Nrj != pIncidentEnergyStorage.end()
1202 && iterator_TransfNrj != pVecmStorage.end())
1203 {
1204 vector<TriDimensionMap>* pDiffCrossSectionData = (iterator_Proba->second);
1205 vector<G4double>* pTdummyVec = iterator_Nrj->second; //Incident energies for interpolation
1206 VecMap* pVecm = iterator_TransfNrj->second;
1207
1208 // k should be in eV and energy transfer eV also
1209 auto t2 =
1210 std::upper_bound(pTdummyVec->begin(), pTdummyVec->end(), k);
1211 auto t1 = t2 - 1;
1212 if (energyTransfer <= (*pVecm)[(*t1)].back() && energyTransfer <= (*pVecm)[(*t2)].back())
1213 {
1214 auto e12 = std::upper_bound((*pVecm)[(*t1)].begin(), (*pVecm)[(*t1)].end(), energyTransfer);
1215 auto e11 = e12 - 1;
1216 auto e22 = std::upper_bound((*pVecm)[(*t2)].begin(), (*pVecm)[(*t2)].end(), energyTransfer);
1217 auto e21 = e22 - 1;
1218
1219 valueT1 = *t1;
1220 valueT2 = *t2;
1221 valueE21 = *e21;
1222 valueE22 = *e22;
1223 valueE12 = *e12;
1224 valueE11 = *e11;
1225
1226 xs11 = (*pDiffCrossSectionData)[LevelIndex][valueT1][valueE11];
1227 xs12 = (*pDiffCrossSectionData)[LevelIndex][valueT1][valueE12];
1228 xs21 = (*pDiffCrossSectionData)[LevelIndex][valueT2][valueE21];
1229 xs22 = (*pDiffCrossSectionData)[LevelIndex][valueT2][valueE22];
1230 }
1231 }
1232 else {
1233 G4String str = "Material ";
1234 str += currentMaterial + " not found!";
1235 G4Exception("G4MicroElecDielectricModels::DifferentialCrossSection", "em0002", FatalException, str);
1236 }
1237 }
1238
1239 G4double xsProduct = xs11 * xs12 * xs21 * xs22;
1240 if (xsProduct != 0.)
1241 {
1242 sigma = QuadInterpolator( valueE11, valueE12,
1243 valueE21, valueE22,
1244 xs11, xs12,
1245 xs21, xs22,
1246 valueT1, valueT2,
1247 k, energyTransfer);
1248 }
1249
1250 }
1251
1252 return sigma;
1253}
1254
1255//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1256
1257
1258G4double G4MicroElecInelasticModel_new::Interpolate(G4double e1,
1259 G4double e2,
1260 G4double e,
1261 G4double xs1,
1262 G4double xs2)
1263{
1264 G4double value = 0.;
1265
1266 // Log-log interpolation by default
1267 if (e1 != 0 && e2 != 0 && (e2-e1) != 0 && !fasterCode)
1268 {
1269 G4double a = std::log(xs2/xs1)/ std::log(e2/e1);
1270 G4double b = std::log(xs2) - a * std::log(e2);
1271 G4double sigma = a * std::log(e) + b;
1272 value = (std::exp(sigma));
1273 }
1274
1275 // Switch to log-lin interpolation for faster code
1276 if ((e2 - e1) != 0 && xs1 != 0 && xs2 != 0 && fasterCode)
1277 {
1278 G4double d1 = std::log(xs1);
1279 G4double d2 = std::log(xs2);
1280 value = std::exp((d1 + (d2 - d1) * (e - e1) / (e2 - e1)));
1281 }
1282
1283 // Switch to lin-lin interpolation for faster code
1284 // in case one of xs1 or xs2 (=cum proba) value is zero
1285 if ((e2 - e1) != 0 && (xs1 == 0 || xs2 == 0) && fasterCode)
1286 {
1287 G4double d1 = xs1;
1288 G4double d2 = xs2;
1289 value = (d1 + (d2 - d1) * (e - e1) / (e2 - e1));
1290 }
1291
1292 return value;
1293}
1294
1295//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1296
1297G4double G4MicroElecInelasticModel_new::QuadInterpolator(G4double e11, G4double e12,
1298 G4double e21, G4double e22,
1299 G4double xs11, G4double xs12,
1300 G4double xs21, G4double xs22,
1301 G4double t1, G4double t2,
1302 G4double t, G4double e)
1303{
1304 G4double interpolatedvalue1 = Interpolate(e11, e12, e, xs11, xs12);
1305 G4double interpolatedvalue2 = Interpolate(e21, e22, e, xs21, xs22);
1306 G4double value = Interpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2);
1307 return value;
1308}
1309
1310//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1311
1312G4int G4MicroElecInelasticModel_new::RandomSelect(G4double k, const G4String& particle, G4double originalMass, G4int originalZ )
1313{
1314 G4int level = 0;
1315
1316 TCSMap::iterator tablepos;
1317 tablepos = tableTCS.find(currentMaterial);
1318 MapData* tableData = tablepos->second;
1319
1320 std::map< G4String,G4MicroElecCrossSectionDataSet_new*,std::less<G4String> >::iterator pos;
1321 pos = tableData->find(particle);
1322
1323 std::vector<G4double> Zeff(currentMaterialStructure->NumberOfLevels(), 1.0);
1324 if(originalMass>proton_mass_c2) {
1325 for(G4int nl=0;nl<currentMaterialStructure->NumberOfLevels();nl++) {
1326 Zeff[nl] = BKZ(k/(proton_mass_c2/originalMass), originalMass/c_squared, originalZ, currentMaterialStructure->Energy(nl));
1327 }
1328 }
1329
1330 if (pos != tableData->end())
1331 {
1332 G4MicroElecCrossSectionDataSet_new* table = pos->second;
1333
1334 if (table != 0)
1335 {
1336 G4double* valuesBuffer = new G4double[table->NumberOfComponents()];
1337 const G4int n = (G4int)table->NumberOfComponents();
1338 G4int i = (G4int)n;
1339 G4double value = 0.;
1340
1341 while (i>0)
1342 {
1343 --i;
1344 valuesBuffer[i] = table->GetComponent(i)->FindValue(k)*Zeff[i]*Zeff[i];
1345 value += valuesBuffer[i];
1346 }
1347 value *= G4UniformRand();
1348
1349 i = n;
1350
1351 while (i > 0)
1352 {
1353 --i;
1354
1355 if (valuesBuffer[i] > value)
1356 {
1357 delete[] valuesBuffer;
1358 return i;
1359 }
1360 value -= valuesBuffer[i];
1361 }
1362
1363 if (valuesBuffer) delete[] valuesBuffer;
1364
1365 }
1366 }
1367 else
1368 {
1369 G4Exception("G4MicroElecInelasticModel_new::RandomSelect","em0002",FatalException,"Model not applicable to particle type.");
1370 }
1371
1372 return level;
1373}
1374
1375//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1376
1378 G4double x = E/mass;
1379 return c_light*std::sqrt(x*(x + 2.0))/(x + 1.0);
1380}
1381
1382//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1383
1385 G4double v1i = ComputeRelativistVelocity(T1i, M1);
1386 G4double v2i = ComputeRelativistVelocity(T2i, M2);
1387
1388 G4double v2f = 2*M1/(M1+M2)*v1i + (M2-M1)/(M1+M2)*-1*v2i;
1389 G4double vtransfer2a = v2f*v2f-v2i*v2i;
1390
1391 v2f = 2*M1/(M1+M2)*v1i + (M2-M1)/(M1+M2)*v2i;
1392 G4double vtransfer2b = v2f*v2f-v2i*v2i;
1393
1394 G4double vtransfer2 = std::max(vtransfer2a, vtransfer2b);
1395 return 0.5*M2*vtransfer2;
1396}
1397
1398//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1399
1401 return (x < 0.) ? 1.0 : 0.0;
1402}
1403
1404//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1405
1407{
1408 G4double r = vF*( std::pow(v/vF+1., 3.) - fabs(std::pow(v/vF-1., 3.))
1409 + 4.*(v/vF)*(v/vF) ) + stepFunc(v/vF-1.) * (3./2.*v/vF -
1410 4.*(v/vF)*(v/vF) + 3.*std::pow(v/vF, 3.)
1411 - 0.5*std::pow(v/vF, 5.));
1412 return r/(10.*v/vF);
1413}
1414
1415//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1416
1418{
1419 // need atomic unit conversion
1420 G4double hbar = hbar_Planck, hbar2 = hbar*hbar, me = electron_mass_c2/c_squared, Ry = me*elm_coupling*elm_coupling/(2*hbar2);
1421 G4double hartree = 2*Ry, a0 = Bohr_radius, velocity = a0*hartree/hbar;
1423
1424 vp /= velocity;
1425
1426 G4double wp = Eplasmon/hartree;
1427 G4double a = std::pow(4./9./CLHEP::pi, 1./3.);
1428 G4double vF = std::pow(wp*wp/(3.*a*a*a), 1./3.);
1429 G4double c = 0.9;
1430 G4double vr = vrkreussler(vp /*in u.a*/, vF /*in u.a*/);
1431 G4double yr = vr/std::pow(Zp, 2./3.);
1432 G4double q = 0.;
1433 if(Zp==2) q = 1-exp(-c*vr/(Zp-5./16.));
1434 else q = 1.-exp(-c*(yr-0.07));
1435 G4double Neq = Zp*(1.-q);
1436 G4double l0 = 0.;
1437 if(Neq<=2) l0 = 3./(Zp-0.3*(Neq-1))/2.;
1438 else l0 = 0.48*std::pow(Neq, 2./3.)/(Zp-Neq/7.);
1439 if(Zp==2) c = 1.0;
1440 else c = 3./2.;
1441 return Zp*(q + c*(1.-q)/vF/vF/2.0 * log(1.+std::pow(2.*l0*vF,2.)));
1442}
1443
1444//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
G4AtomicShellEnumerator
const char * G4FindDataDir(const char *)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
const G4double a0
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
double z() const
Hep3Vector unit() const
double x() const
double y() const
void set(double x, double y, double z)
const G4ThreeVector & GetMomentumDirection() const
const G4ParticleDefinition * GetParticleDefinition() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:88
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
const G4Material * GetMaterial() const
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:204
const G4String & GetName() const
Definition: G4Material.hh:172
G4double FindShellValue(G4double argEnergy, G4int shell) const
const G4VEMDataSet * GetComponent(G4int componentId) const override
G4bool LoadData(const G4String &argFileName) override
G4double FindValue(G4double e, G4int componentId=0) const override
G4double DifferentialCrossSection(const G4ParticleDefinition *aParticleDefinition, G4double k, G4double energyTransfer, G4int shell)
G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
G4MicroElecInelasticModel_new(const G4ParticleDefinition *p=nullptr, const G4String &nam="MicroElecInelasticModel")
G4double ComputeElasticQmax(G4double T1i, G4double T2i, G4double m1, G4double m2)
G4ParticleChangeForGamma * fParticleChangeForGamma
G4double BKZ(G4double Ep, G4double mp, G4int Zp, G4double EF)
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4double ComputeRelativistVelocity(G4double E, G4double mass)
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4double vrkreussler(G4double v, G4double vF)
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
G4int GetAtomicNumber() const
G4double GetPDGCharge() const
const G4String & GetParticleName() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:87
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
virtual G4double FindValue(G4double x, G4int componentId=0) const =0
virtual G4ThreeVector & SampleDirectionForShell(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, G4int shellID, const G4Material *)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:746
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:600
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 SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:802
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:607
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double bindingEnergy(G4int A, G4int Z)
#define DBL_MAX
Definition: templates.hh:62