112 G4cout <<
"MicroElec inelastic model is constructed " <<
G4endl;
117 fAtomDeexcitation =
nullptr;
118 fParticleChangeForGamma =
nullptr;
125 SEFromFermiLevel =
false;
134 for (
auto const& p : tableTCS) {
135 MapData* tableData = p.second;
136 for (
auto const& pos : *tableData) {
delete pos.second; }
142 for (
auto const & obj : eDiffDatatable) {
143 auto ptr = obj.second;
148 for (
auto const & obj : pDiffDatatable) {
149 auto ptr = obj.second;
155 for (
auto const & obj : eNrjTransStorage) {
158 for (
auto const & obj : pNrjTransStorage) {
163 for (
auto const& p : eProbaShellStorage) {
167 for (
auto const& p : pProbaShellStorage) {
172 for (
auto const& p : eIncidentEnergyStorage) {
176 for (
auto const& p : pIncidentEnergyStorage) {
181 for (
auto const& p : eVecmStorage) {
185 for (
auto const& p : pVecmStorage) {
190 for (
auto const& p : tableMaterialsStructures) {
200 if (isInitialised) {
return; }
202 if (verboseLevel > 3)
203 G4cout <<
"Calling G4MicroElecInelasticModel_new::Initialise()" <<
G4endl;
225 lowEnergyLimit[electron] = 2 * eV;
226 highEnergyLimit[electron] = 10.0 * MeV;
232 for (
G4int i = 0; i < numOfCouples; ++i) {
234 G4cout <<
"Material " << i + 1 <<
" / " << numOfCouples <<
" : " << material->
GetName() <<
G4endl;
235 if (material->
GetName() ==
"Vacuum")
continue;
237 MapData* tableData =
new MapData;
240 tableMaterialsStructures[mat] = currentMaterialStructure;
241 if (particle == electronDef) {
243 G4String fileElectron(
"Inelastic/" + modelName +
"_sigma_inelastic_e-_" + mat);
247 tableData->insert(make_pair(electron, tableE));
250 std::ostringstream eFullFileName;
252 eFullFileName << path <<
"/microelec/Inelastic/cumulated_" + modelName +
"_sigmadiff_inelastic_e-_" + mat +
".dat";
254 G4cout <<
"Inelastic/cumulated_" + modelName +
"_sigmadiff_inelastic_e-_" + mat +
".dat" <<
G4endl;
257 eFullFileName << path <<
"/microelec/Inelastic/" + modelName +
"_sigmadiff_inelastic_e-_" + mat +
".dat";
259 G4cout <<
"Inelastic/" + modelName +
"_sigmadiff_inelastic_e-_" + mat +
".dat" <<
G4endl;
262 std::ifstream eDiffCrossSection(eFullFileName.str().c_str());
263 if (!eDiffCrossSection)
265 std::stringstream ss;
266 ss <<
"Missing data " << eFullFileName.str().c_str();
267 std::string sortieString = ss.str();
269 if (fasterCode)
G4Exception(
"G4MicroElecInelasticModel_new::Initialise",
"em0003",
272 G4Exception(
"G4MicroElecInelasticModel_new::Initialise",
"em0003",
273 FatalException,
"Missing data file:/microelec/sigmadiff_inelastic_e_Si.dat");
282 vector<TriDimensionMap>* eDiffCrossSectionData =
283 new vector<TriDimensionMap>;
284 vector<TriDimensionMap>* eNrjTransfData =
285 new vector<TriDimensionMap>;
286 vector<VecMap>* eProbaShellMap =
new vector<VecMap>;
287 vector<G4double>* eTdummyVec =
new vector<G4double>;
288 VecMap* eVecm =
new VecMap;
292 eDiffCrossSectionData->push_back(TriDimensionMap());
293 eNrjTransfData->push_back(TriDimensionMap());
294 eProbaShellMap->push_back(VecMap());
297 eTdummyVec->push_back(0.);
298 while (!eDiffCrossSection.eof())
302 eDiffCrossSection >> tDummy >> eDummy;
303 if (tDummy != eTdummyVec->back()) eTdummyVec->push_back(tDummy);
308 eDiffCrossSection >> tmp;
309 (*eDiffCrossSectionData)[j][tDummy][eDummy] = tmp;
313 (*eNrjTransfData)[j][tDummy][(*eDiffCrossSectionData)[j][tDummy][eDummy]] = eDummy;
314 (*eProbaShellMap)[j][tDummy].push_back((*eDiffCrossSectionData)[j][tDummy][eDummy]);
317 (*eDiffCrossSectionData)[j][tDummy][eDummy] *= scaleFactor;
318 (*eVecm)[tDummy].push_back(eDummy);
328 eNrjTransStorage[mat] = eNrjTransfData;
329 eProbaShellStorage[mat] = eProbaShellMap;
332 eDiffDatatable[mat] = eDiffCrossSectionData;
333 eVecmStorage[mat] = eVecm;
335 eIncidentEnergyStorage[mat] = eTdummyVec;
339 delete eProbaShellMap;
340 delete eNrjTransfData;
342 delete eDiffCrossSectionData;
348 if (particle == protonDef)
351 G4String fileProton(
"Inelastic/" + modelName +
"_sigma_inelastic_p_" + mat);
G4cout << fileProton <<
G4endl;
354 tableData->insert(make_pair(proton, tableP));
357 std::ostringstream pFullFileName;
359 pFullFileName << path <<
"/microelec/Inelastic/cumulated_" + modelName +
"_sigmadiff_inelastic_p_" + mat +
".dat";
361 G4cout <<
"Inelastic/cumulated_" + modelName +
"_sigmadiff_inelastic_p_" + mat +
".dat" <<
G4endl;
364 pFullFileName << path <<
"/microelec/Inelastic/" + modelName +
"_sigmadiff_inelastic_p_" + mat +
".dat";
366 G4cout <<
"Inelastic/" + modelName +
"_sigmadiff_inelastic_e-_" + mat +
".dat" <<
G4endl;
369 std::ifstream pDiffCrossSection(pFullFileName.str().c_str());
370 if (!pDiffCrossSection)
372 if (fasterCode)
G4Exception(
"G4MicroElecInelasticModel_new::Initialise",
"em0003",
373 FatalException,
"Missing data file:/microelec/sigmadiff_cumulated_inelastic_p_Si.dat");
375 G4Exception(
"G4MicroElecInelasticModel_new::Initialise",
"em0003",
376 FatalException,
"Missing data file:/microelec/sigmadiff_inelastic_p_Si.dat");
387 vector<TriDimensionMap>* pDiffCrossSectionData =
388 new vector<TriDimensionMap>;
389 vector<TriDimensionMap>* pNrjTransfData =
390 new vector<TriDimensionMap>;
391 vector<VecMap>* pProbaShellMap =
393 vector<G4double>* pTdummyVec =
394 new vector<G4double>;
395 VecMap* pVecm =
new VecMap;
400 pDiffCrossSectionData->push_back(TriDimensionMap());
401 pNrjTransfData->push_back(TriDimensionMap());
402 pProbaShellMap->push_back(VecMap());
405 pTdummyVec->push_back(0.);
406 while (!pDiffCrossSection.eof())
410 pDiffCrossSection >> tDummy >> eDummy;
411 if (tDummy != pTdummyVec->back()) pTdummyVec->push_back(tDummy);
416 pDiffCrossSection >> tmp;
417 (*pDiffCrossSectionData)[j][tDummy][eDummy] = tmp;
424 (*pNrjTransfData)[j][tDummy][(*pDiffCrossSectionData)[j][tDummy][eDummy]] = eDummy;
425 (*pProbaShellMap)[j][tDummy].push_back((*pDiffCrossSectionData)[j][tDummy][eDummy]);
428 (*pDiffCrossSectionData)[j][tDummy][eDummy] *= scaleFactor;
429 (*pVecm)[tDummy].push_back(eDummy);
437 pNrjTransStorage[mat] = pNrjTransfData;
438 pProbaShellStorage[mat] = pProbaShellMap;
441 pDiffDatatable[mat] = pDiffCrossSectionData;
442 pVecmStorage[mat] = pVecm;
444 pIncidentEnergyStorage[mat] = pTdummyVec;
448 delete pProbaShellMap;
449 delete pNrjTransfData;
451 delete pDiffCrossSectionData;
455 tableTCS[mat] = tableData;
457 if (particle==electronDef)
463 if (particle==protonDef)
471 G4cout <<
"MicroElec Inelastic model is initialized " <<
G4endl
476 <<
" with mass (amu) " << particle->
GetPDGMass()/proton_mass_c2
484 isInitialised =
true;
495 if (verboseLevel > 3)
G4cout <<
"Calling CrossSectionPerVolume() of G4MicroElecInelasticModel" <<
G4endl;
498 currentMaterial = material->
GetName().substr(3, material->
GetName().size());
500 MapStructure::iterator structPos;
501 structPos = tableMaterialsStructures.find(currentMaterial);
504 TCSMap::iterator tablepos;
505 tablepos = tableTCS.find(currentMaterial);
507 if (tablepos == tableTCS.end() )
510 str += currentMaterial +
" TCS Table not found!";
514 else if(structPos == tableMaterialsStructures.end())
517 str += currentMaterial +
" Structure not found!";
522 MapData* tableData = tablepos->second;
523 currentMaterialStructure = structPos->second;
532 G4double Zeff = 1.0, Zeff2 = Zeff*Zeff;
535 if (Mion_c2 > proton_mass_c2)
537 ekin *= proton_mass_c2 / Mion_c2;
538 nameLocal =
"proton";
544 if (ekin >= lowLim && ekin < highLim)
546 std::map< G4String, G4MicroElecCrossSectionDataSet_new*, std::less<G4String> >::iterator pos;
547 pos = tableData->find(nameLocal);
549 if (pos != tableData->end())
556 if (Mion_c2 > proton_mass_c2) {
559 Zeff =
BKZ(ekin / (proton_mass_c2 / Mion_c2), Mion_c2 / c_squared, Z, currentMaterialStructure->
Energy(i));
574 G4Exception(
"G4MicroElecInelasticModel_new::CrossSectionPerVolume",
576 "Model not applicable to particle type.");
584 if (verboseLevel > 3)
586 G4cout <<
"---> Kinetic energy (eV)=" << ekin / eV <<
G4endl;
587 G4cout <<
" - Cross section per Si atom (cm^2)=" << sigma / cm2 <<
G4endl;
588 G4cout <<
" - Cross section per Si atom (cm^-1)=" << sigma*density / (1. / cm) <<
G4endl;
591 return (sigma)*density;
604 if (verboseLevel > 3)
605 G4cout <<
"Calling SampleSecondaries() of G4MicroElecInelasticModel" <<
G4endl;
616 G4String nameLocal2 = particleName ;
618 G4double originalMass = particleMass;
621 if (particleMass > proton_mass_c2)
623 k *= proton_mass_c2/particleMass ;
625 nameLocal2 =
"proton" ;
628 if (k >= lowLim && k < highLim)
631 G4double totalEnergy = ekin + particleMass;
632 G4double pSquare = ekin * (totalEnergy + particleMass);
633 G4double totalMomentum = std::sqrt(pSquare);
637 Shell = RandomSelect(k,nameLocal2,originalMass, originalZ);
644 if (verboseLevel > 3)
646 G4cout <<
"---> Kinetic energy (eV)=" << k/eV <<
G4endl ;
647 G4cout <<
"Shell: " << Shell <<
", energy: " << bindingEnergy/eV <<
G4endl;
652 if (weaklyBound && k > currentMaterialStructure->
GetEnergyGap()) {
659 std::size_t secNumberInit = 0;
660 std::size_t secNumberFinal = 0;
663 G4int Z = currentMaterialStructure->
GetZ(Shell);
667 if(fAtomDeexcitation && shellEnum >=0)
672 secNumberInit = fvect->size();
674 secNumberFinal = fvect->size();
678 SEFromFermiLevel =
false;
681 secondaryKinetic = RandomizeEjectedElectronEnergy(PartDef, k, Shell, originalMass, originalZ);
685 secondaryKinetic = RandomizeEjectedElectronEnergyFromCumulatedDcs(PartDef, k, Shell) ;
688 if (verboseLevel > 3)
691 G4cout <<
"Shell: " << Shell <<
" Kin. energy (eV)=" << k/eV
692 <<
" Sec. energy (eV)=" << secondaryKinetic/eV <<
G4endl;
701 G4double deltaTotalMomentum = std::sqrt(secondaryKinetic*(secondaryKinetic + 2.*electron_mass_c2 ));
703 G4double finalPx = totalMomentum*primaryDirection.
x() - deltaTotalMomentum*deltaDirection.
x();
704 G4double finalPy = totalMomentum*primaryDirection.
y() - deltaTotalMomentum*deltaDirection.
y();
705 G4double finalPz = totalMomentum*primaryDirection.
z() - deltaTotalMomentum*deltaDirection.
z();
706 G4double finalMomentum = std::sqrt(finalPx*finalPx + finalPy*finalPy + finalPz*finalPz);
707 finalPx /= finalMomentum;
708 finalPy /= finalMomentum;
709 finalPz /= finalMomentum;
712 direction.
set(finalPx,finalPy,finalPz);
720 for (std::size_t j=secNumberInit; j < secNumberFinal; ++j) {
721 deexSecEnergy = deexSecEnergy + (*fvect)[j]->GetKineticEnergy();
734 if (secondaryKinetic>0)
737 fvect->push_back(dp);
744G4double G4MicroElecInelasticModel_new::RandomizeEjectedElectronEnergy(
748 G4double secondaryElectronKineticEnergy=0.;
754 G4double maxEnergy = maximumEnergyTransfer;
755 G4int nEnergySteps = 100;
758 G4double stpEnergy(std::pow(maxEnergy/value, 1./
static_cast<G4double>(nEnergySteps-1)));
759 G4int step(nEnergySteps);
765 crossSectionMaximum = std::max(crossSectionMaximum, differentialCrossSection);
772 (maximumEnergyTransfer-currentMaterialStructure->
GetLimitEnergy(shell));
775 (secondaryElectronKineticEnergy+currentMaterialStructure->
GetLimitEnergy(shell)),shell));
777 return secondaryElectronKineticEnergy;
783 currentMaterialStructure->
Energy(shell),
784 originalMass/c_squared, electron_mass_c2/c_squared);
789 G4double maxEnergy = maximumEnergyTransfer;
790 G4int nEnergySteps = 100;
793 G4double stpEnergy(std::pow(maxEnergy/value, 1./
static_cast<G4double>(nEnergySteps-1)));
794 G4int step(nEnergySteps);
801 crossSectionMaximum = std::max(crossSectionMaximum, differentialCrossSection);
812 secondaryElectronKineticEnergy =
816 return std::max(secondaryElectronKineticEnergy, 0.0);
821G4double G4MicroElecInelasticModel_new::RandomizeEjectedElectronEnergyFromCumulatedDcs(
824 G4double secondaryElectronKineticEnergy = 0.;
827 G4double transf = TransferedEnergy(particleDefinition, k, shell, random);
829 secondaryElectronKineticEnergy = transf - currentMaterialStructure->
GetLimitEnergy(shell);
830 if(secondaryElectronKineticEnergy <= 0.) {
831 secondaryElectronKineticEnergy = 0.0;
835 secondaryElectronKineticEnergy = transf - currentMaterialStructure->
GetLimitEnergy(shell);
837 if (secondaryElectronKineticEnergy <= 0.) {
838 secondaryElectronKineticEnergy = 0.0;
839 SEFromFermiLevel =
true;
843 return secondaryElectronKineticEnergy;
848G4double G4MicroElecInelasticModel_new::TransferedEnergy(
851 G4int ionizationLevelIndex,
866 G4double maximumEnergyTransfer1 = 0;
867 G4double maximumEnergyTransfer2 = 0;
868 G4double maximumEnergyTransferP = 4.* (electron_mass_c2 / proton_mass_c2) * k;
873 dataDiffCSMap::iterator iterator_Nrj;
874 iterator_Nrj = eNrjTransStorage.find(currentMaterial);
876 dataProbaShellMap::iterator iterator_Proba;
877 iterator_Proba = eProbaShellStorage.find(currentMaterial);
879 incidentEnergyMap::iterator iterator_Tdummy;
880 iterator_Tdummy = eIncidentEnergyStorage.find(currentMaterial);
882 if(iterator_Nrj == eNrjTransStorage.end() || iterator_Proba == eProbaShellStorage.end() ||
883 iterator_Tdummy == eIncidentEnergyStorage.end())
886 str += currentMaterial +
" not found!";
887 G4Exception(
"G4MicroElecInelasticModel_new::TransferedEnergy",
"em0002",
891 vector<TriDimensionMap>* eNrjTransfData = iterator_Nrj->second;
892 vector<VecMap>* eProbaShellMap = iterator_Proba->second;
893 vector<G4double>* eTdummyVec = iterator_Tdummy->second;
896 auto k2 = std::upper_bound(eTdummyVec->begin(),
902 if (random <= (*eProbaShellMap)[ionizationLevelIndex][(*k1)].back()
903 && random <= (*eProbaShellMap)[ionizationLevelIndex][(*k2)].back())
906 std::upper_bound((*eProbaShellMap)[ionizationLevelIndex][(*k1)].begin(),
907 (*eProbaShellMap)[ionizationLevelIndex][(*k1)].end(),
910 auto prob11 = prob12 - 1;
913 std::upper_bound((*eProbaShellMap)[ionizationLevelIndex][(*k2)].begin(),
914 (*eProbaShellMap)[ionizationLevelIndex][(*k2)].end(),
917 auto prob21 = prob22 - 1;
921 valuePROB21 = *prob21;
922 valuePROB22 = *prob22;
923 valuePROB12 = *prob12;
924 valuePROB11 = *prob11;
928 else nrjTransf11 = (*eNrjTransfData)[ionizationLevelIndex][valueK1][valuePROB11];
929 if (valuePROB12 == 1)
931 if ((valueK1 + bindingEnergy) / 2. > valueK1)
932 maximumEnergyTransfer1 = valueK1;
936 nrjTransf12 = maximumEnergyTransfer1;
939 nrjTransf12 = (*eNrjTransfData)[ionizationLevelIndex][valueK1][valuePROB12];
942 else nrjTransf21 = (*eNrjTransfData)[ionizationLevelIndex][valueK2][valuePROB21];
943 if (valuePROB22 == 1)
945 if ((valueK2 + bindingEnergy) / 2. > valueK2) maximumEnergyTransfer2 = valueK2;
946 else maximumEnergyTransfer2 = (valueK2 +
bindingEnergy) / 2.;
948 nrjTransf22 = maximumEnergyTransfer2;
950 else nrjTransf22 = (*eNrjTransfData)[ionizationLevelIndex][valueK2][valuePROB22];
954 if (random > (*eProbaShellMap)[ionizationLevelIndex][(*k1)].back())
957 std::upper_bound((*eProbaShellMap)[ionizationLevelIndex][(*k2)].begin(),
958 (*eProbaShellMap)[ionizationLevelIndex][(*k2)].end(),
960 auto prob21 = prob22 - 1;
964 valuePROB21 = *prob21;
965 valuePROB22 = *prob22;
967 nrjTransf21 = (*eNrjTransfData)[ionizationLevelIndex][valueK2][valuePROB21];
968 nrjTransf22 = (*eNrjTransfData)[ionizationLevelIndex][valueK2][valuePROB22];
970 G4double interpolatedvalue2 = Interpolate(valuePROB21,
977 G4double value = Interpolate(valueK1, valueK2, k, 0., interpolatedvalue2);
986 dataDiffCSMap::iterator iterator_Nrj;
987 iterator_Nrj = pNrjTransStorage.find(currentMaterial);
989 dataProbaShellMap::iterator iterator_Proba;
990 iterator_Proba = pProbaShellStorage.find(currentMaterial);
992 incidentEnergyMap::iterator iterator_Tdummy;
993 iterator_Tdummy = pIncidentEnergyStorage.find(currentMaterial);
995 if (iterator_Nrj == pNrjTransStorage.end() || iterator_Proba == pProbaShellStorage.end() ||
996 iterator_Tdummy == pIncidentEnergyStorage.end())
999 str += currentMaterial +
" not found!";
1000 G4Exception(
"G4MicroElecInelasticModel_new::TransferedEnergy",
"em0002",
1005 vector<TriDimensionMap>* pNrjTransfData = iterator_Nrj->second;
1006 vector<VecMap>* pProbaShellMap = iterator_Proba->second;
1007 vector<G4double>* pTdummyVec = iterator_Tdummy->second;
1009 auto k2 = std::upper_bound(pTdummyVec->begin(),
1017 if (random <= (*pProbaShellMap)[ionizationLevelIndex][(*k1)].back()
1018 && random <= (*pProbaShellMap)[ionizationLevelIndex][(*k2)].back())
1021 std::upper_bound((*pProbaShellMap)[ionizationLevelIndex][(*k1)].begin(),
1022 (*pProbaShellMap)[ionizationLevelIndex][(*k1)].end(),
1024 auto prob11 = prob12 - 1;
1026 std::upper_bound((*pProbaShellMap)[ionizationLevelIndex][(*k2)].begin(),
1027 (*pProbaShellMap)[ionizationLevelIndex][(*k2)].end(),
1029 auto prob21 = prob22 - 1;
1033 valuePROB21 = *prob21;
1034 valuePROB22 = *prob22;
1035 valuePROB12 = *prob12;
1036 valuePROB11 = *prob11;
1041 else nrjTransf11 = (*pNrjTransfData)[ionizationLevelIndex][valueK1][valuePROB11];
1043 if (valuePROB12 == 1) nrjTransf12 = maximumEnergyTransferP;
1044 else nrjTransf12 = (*pNrjTransfData)[ionizationLevelIndex][valueK1][valuePROB12];
1047 else nrjTransf21 = (*pNrjTransfData)[ionizationLevelIndex][valueK2][valuePROB21];
1049 if (valuePROB22 == 1) nrjTransf22 = maximumEnergyTransferP;
1050 else nrjTransf22 = (*pNrjTransfData)[ionizationLevelIndex][valueK2][valuePROB22];
1055 if (random > (*pProbaShellMap)[ionizationLevelIndex][(*k1)].back())
1058 std::upper_bound((*pProbaShellMap)[ionizationLevelIndex][(*k2)].begin(),
1059 (*pProbaShellMap)[ionizationLevelIndex][(*k2)].end(),
1062 auto prob21 = prob22 - 1;
1066 valuePROB21 = *prob21;
1067 valuePROB22 = *prob22;
1069 nrjTransf21 = (*pNrjTransfData)[ionizationLevelIndex][valueK2][valuePROB21];
1070 nrjTransf22 = (*pNrjTransfData)[ionizationLevelIndex][valueK2][valuePROB22];
1072 G4double interpolatedvalue2 = Interpolate(valuePROB21,
1079 G4double value = Interpolate(valueK1, valueK2, k, 0., interpolatedvalue2);
1086 G4double nrjTransfProduct = nrjTransf11 * nrjTransf12 * nrjTransf21 * nrjTransf22;
1088 if (nrjTransfProduct != 0.)
1090 nrj = QuadInterpolator(valuePROB11,
1117 if (energyTransfer >= currentMaterialStructure->
GetLimitEnergy(LevelIndex))
1134 dataDiffCSMap::iterator iterator_Proba;
1135 iterator_Proba = eDiffDatatable.find(currentMaterial);
1137 incidentEnergyMap::iterator iterator_Nrj;
1138 iterator_Nrj = eIncidentEnergyStorage.find(currentMaterial);
1140 TranfEnergyMap::iterator iterator_TransfNrj;
1141 iterator_TransfNrj = eVecmStorage.find(currentMaterial);
1143 if (iterator_Proba != eDiffDatatable.end() && iterator_Nrj != eIncidentEnergyStorage.end()
1144 && iterator_TransfNrj!= eVecmStorage.end())
1146 vector<TriDimensionMap>* eDiffCrossSectionData = (iterator_Proba->second);
1147 vector<G4double>* eTdummyVec = iterator_Nrj->second;
1148 VecMap* eVecm = iterator_TransfNrj->second;
1151 auto t2 = std::upper_bound(eTdummyVec->begin(), eTdummyVec->end(), k);
1154 if (energyTransfer <= (*eVecm)[(*t1)].back() && energyTransfer <= (*eVecm)[(*t2)].back())
1156 auto e12 = std::upper_bound((*eVecm)[(*t1)].begin(), (*eVecm)[(*t1)].end(), energyTransfer);
1158 auto e22 = std::upper_bound((*eVecm)[(*t2)].begin(), (*eVecm)[(*t2)].end(), energyTransfer);
1168 xs11 = (*eDiffCrossSectionData)[LevelIndex][valueT1][valueE11];
1169 xs12 = (*eDiffCrossSectionData)[LevelIndex][valueT1][valueE12];
1170 xs21 = (*eDiffCrossSectionData)[LevelIndex][valueT2][valueE21];
1171 xs22 = (*eDiffCrossSectionData)[LevelIndex][valueT2][valueE22];
1176 str += currentMaterial +
" not found!";
1183 dataDiffCSMap::iterator iterator_Proba;
1184 iterator_Proba = pDiffDatatable.find(currentMaterial);
1186 incidentEnergyMap::iterator iterator_Nrj;
1187 iterator_Nrj = pIncidentEnergyStorage.find(currentMaterial);
1189 TranfEnergyMap::iterator iterator_TransfNrj;
1190 iterator_TransfNrj = pVecmStorage.find(currentMaterial);
1192 if (iterator_Proba != pDiffDatatable.end() && iterator_Nrj != pIncidentEnergyStorage.end()
1193 && iterator_TransfNrj != pVecmStorage.end())
1195 vector<TriDimensionMap>* pDiffCrossSectionData = (iterator_Proba->second);
1196 vector<G4double>* pTdummyVec = iterator_Nrj->second;
1197 VecMap* pVecm = iterator_TransfNrj->second;
1201 std::upper_bound(pTdummyVec->begin(), pTdummyVec->end(), k);
1203 if (energyTransfer <= (*pVecm)[(*t1)].back() && energyTransfer <= (*pVecm)[(*t2)].back())
1205 auto e12 = std::upper_bound((*pVecm)[(*t1)].begin(), (*pVecm)[(*t1)].end(), energyTransfer);
1207 auto e22 = std::upper_bound((*pVecm)[(*t2)].begin(), (*pVecm)[(*t2)].end(), energyTransfer);
1217 xs11 = (*pDiffCrossSectionData)[LevelIndex][valueT1][valueE11];
1218 xs12 = (*pDiffCrossSectionData)[LevelIndex][valueT1][valueE12];
1219 xs21 = (*pDiffCrossSectionData)[LevelIndex][valueT2][valueE21];
1220 xs22 = (*pDiffCrossSectionData)[LevelIndex][valueT2][valueE22];
1225 str += currentMaterial +
" not found!";
1230 G4double xsProduct = xs11 * xs12 * xs21 * xs22;
1231 if (xsProduct != 0.)
1233 sigma = QuadInterpolator( valueE11, valueE12,
1256 if (e == e1 || e1 == e2 ) {
return xs1; }
1257 if (e == e2) {
return xs2; }
1260 if (e1 > 0. && e2 > 0. && xs1 > 0. && xs2 > 0. && !fasterCode)
1262 G4double a = std::log(xs2/xs1)/ std::log(e2/e1);
1263 G4double b = std::log(xs2) - a * std::log(e2);
1264 G4double sigma = a * std::log(e) + b;
1265 value = (std::exp(sigma));
1269 else if (xs1 > 0. && xs2 > 0. && fasterCode)
1273 value = std::exp((d1 + (d2 - d1) * (e - e1) / (e2 - e1)));
1282 value = (d1 + (d2 - d1) * (e - e1) / (e2 - e1));
1297 G4double interpolatedvalue1 = Interpolate(e11, e12, e, xs11, xs12);
1298 G4double interpolatedvalue2 = Interpolate(e21, e22, e, xs21, xs22);
1299 G4double value = Interpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2);
1309 TCSMap::iterator tablepos;
1310 tablepos = tableTCS.find(currentMaterial);
1311 MapData* tableData = tablepos->second;
1313 std::map< G4String,G4MicroElecCrossSectionDataSet_new*,std::less<G4String> >::iterator
pos;
1314 pos = tableData->find(particle);
1316 std::vector<G4double> Zeff(currentMaterialStructure->
NumberOfLevels(), 1.0);
1317 if(originalMass>proton_mass_c2) {
1319 Zeff[nl] =
BKZ(k/(proton_mass_c2/originalMass), originalMass/c_squared, originalZ, currentMaterialStructure->
Energy(nl));
1323 if (pos != tableData->end())
1338 value += valuesBuffer[i];
1348 if (valuesBuffer[i] > value)
1350 delete[] valuesBuffer;
1353 value -= valuesBuffer[i];
1356 if (valuesBuffer)
delete[] valuesBuffer;
1362 G4Exception(
"G4MicroElecInelasticModel_new::RandomSelect",
"em0002",
FatalException,
"Model not applicable to particle type.");
1372 return c_light*std::sqrt(x*(x + 2.0))/(x + 1.0);
1381 G4double v2f = 2*M1/(M1+M2)*v1i + (M2-M1)/(M1+M2)*-1*v2i;
1382 G4double vtransfer2a = v2f*v2f-v2i*v2i;
1384 v2f = 2*M1/(M1+M2)*v1i + (M2-M1)/(M1+M2)*v2i;
1385 G4double vtransfer2b = v2f*v2f-v2i*v2i;
1387 G4double vtransfer2 = std::max(vtransfer2a, vtransfer2b);
1388 return 0.5*M2*vtransfer2;
1394 return (x < 0.) ? 1.0 : 0.0;
1401 G4double r = vF*( std::pow(v/vF+1., 3.) - fabs(std::pow(v/vF-1., 3.))
1402 + 4.*(v/vF)*(v/vF) ) +
stepFunc(v/vF-1.) * (3./2.*v/vF -
1403 4.*(v/vF)*(v/vF) + 3.*std::pow(v/vF, 3.)
1404 - 0.5*std::pow(v/vF, 5.));
1405 return r/(10.*v/vF);
1413 G4double hbar = hbar_Planck, hbar2 = hbar*hbar, me = electron_mass_c2/c_squared, Ry = me*elm_coupling*elm_coupling/(2*hbar2);
1414 G4double hartree = 2*Ry,
a0 = Bohr_radius, velocity =
a0*hartree/hbar;
1420 G4double a = std::pow(4./9./CLHEP::pi, 1./3.);
1421 G4double vF = std::pow(wp*wp/(3.*a*a*a), 1./3.);
1424 G4double yr = vr/std::pow(Zp, 2./3.);
1426 if(Zp==2) q = 1-exp(-c*vr/(Zp-5./16.));
1427 else q = 1.-exp(-c*(yr-0.07));
1430 if(Neq<=2) l0 = 3./(Zp-0.3*(Neq-1))/2.;
1431 else l0 = 0.48*std::pow(Neq, 2./3.)/(Zp-Neq/7.);
1434 return Zp*(q + c*(1.-q)/vF/vF/2.0 * log(1.+std::pow(2.*l0*vF,2.)));
const char * G4FindDataDir(const char *)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
G4GLOB_DLL std::ostream G4cout
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()
static G4Electron * Electron()
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
const G4Material * GetMaterial() const
G4double GetTotNbOfAtomsPerVolume() const
const G4String & GetName() const
G4double FindShellValue(G4double argEnergy, G4int shell) const
const G4VEMDataSet * GetComponent(G4int componentId) const override
size_t NumberOfComponents(void) 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 stepFunc(G4double x)
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")
~G4MicroElecInelasticModel_new() override
G4double ComputeElasticQmax(G4double T1i, G4double T2i, G4double m1, G4double m2)
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)
G4double GetInelasticModelHighLimit(G4int pdg)
G4bool IsShellWeaklyBound(G4int level)
G4double GetLimitEnergy(G4int level)
G4int GetEADL_Enumerator(G4int shell)
G4double GetZ(G4int Shell)
G4double Energy(G4int level)
G4double GetInelasticModelLowLimit(G4int pdg)
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
G4int GetAtomicNumber() const
G4double GetPDGMass() const
G4int GetPDGEncoding() 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()
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)
G4VEmAngularDistribution * GetAngularDistribution()
G4ParticleChangeForGamma * GetParticleChangeForGamma()
G4double LowEnergyLimit() const
G4double HighEnergyLimit() const
void SetLowEnergyLimit(G4double)
void SetDeexcitationFlag(G4bool val)
void SetAngularDistribution(G4VEmAngularDistribution *)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double bindingEnergy(G4int A, G4int Z)