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