61 for(
G4int i=0; i<numOfCouples; ++i)
75 Z =
G4lrint((*theElementVector)[0]->GetZ());
78 fkillBelowEnergy_Au = 10. * eV;
79 flowEnergyLimit = 0 * eV;
80 fhighEnergyLimit = 1 * GeV;
87 if(material==fpBaseWater){
88 flowEnergyLimit = 10. * eV;
89 fhighEnergyLimit = 1 * MeV;
99 G4cout <<
"ELSEPA Elastic model is constructed for "
102 << flowEnergyLimit / eV <<
" eV - "
103 << fhighEnergyLimit / MeV <<
" MeV"
109 fpMolDensity =
nullptr;
122 eEdummyVec_Au.clear();
123 eEdummyVec_H2O.clear();
126 fAngleData_Au.clear();
127 fAngleData_H2O.clear();
135 if (verboseLevel > 3)
136 G4cout <<
"Calling G4DNAELSEPAElasticModel::Initialise()" <<
G4endl;
138 if (isInitialised) {
return;}
142 G4Exception(
"G4DNAELSEPAElasticModel::Initialise",
"em0001",
158 for(
G4int i=0; i<numOfCouples; ++i)
176 G4String fileZElectron(
"dna/sigma_elastic_e_elsepa_Z");
177 std::ostringstream oss;
179 oss.clear(stringstream::goodbit);
181 fileZElectron += oss.str()+
"_muffintin";
186 fpData_Au->LoadData(fileZElectron);
188 std::ostringstream eFullFileNameZ;
193 G4Exception(
"G4DNAELSEPAElasticModel::Initialise",
"em0002",
198 eFullFileNameZ.str(
"");
199 eFullFileNameZ.clear(stringstream::goodbit);
203 <<
"/dna/sigmadiff_cumulated_elastic_e_elsepa_Z"
204 << Z <<
"_muffintin.dat";
206 std::ifstream eDiffCrossSectionZ(eFullFileNameZ.str().c_str());
208 if (!eDiffCrossSectionZ)
210 G4Exception(
"G4DNAELSEPAElasticModel::Initialise",
"em0003",
215 eEdummyVec_Au.clear();
217 fAngleData_Au.clear();
219 eEdummyVec_Au.push_back(0.);
224 eDiffCrossSectionZ>>eDummy>>cumDummy;
225 if (eDummy != eEdummyVec_Au.back())
227 eEdummyVec_Au.push_back(eDummy);
228 eCum_Au[eDummy].push_back(0.);
230 eDiffCrossSectionZ>>fAngleData_Au[eDummy][cumDummy];
231 if (cumDummy != eCum_Au[eDummy].back())
233 eCum_Au[eDummy].push_back(cumDummy);
235 }
while(!eDiffCrossSectionZ.eof());
239 if(material == fpBaseWater && !fpData_H2O){
242 G4cout<<
"G4DNAELSEPAElasticModel: low energy limit increased from "
250 G4cout<<
"G4DNAELSEPAElasticModel: high energy limit decreased from "
256 G4String fileZElectron(
"dna/sigma_elastic_e_elsepa_muffin");
261 fpData_H2O->LoadData(fileZElectron);
263 std::ostringstream eFullFileNameZ;
268 G4Exception(
"G4DNAELSEPAElasticModel::Initialise",
"em0004",
273 eFullFileNameZ.str(
"");
274 eFullFileNameZ.clear(stringstream::goodbit);
278 <<
"/dna/sigmadiff_cumulated_elastic_e_elsepa_muffin.dat";
280 std::ifstream eDiffCrossSectionZ(eFullFileNameZ.str().c_str());
282 if (!eDiffCrossSectionZ)
283 G4Exception(
"G4DNAELSEPAElasticModel::Initialise",
"em0005",
285 "Missing data file for cumulated DCS");
287 eEdummyVec_H2O.clear();
289 fAngleData_H2O.clear();
291 eEdummyVec_H2O.push_back(0.);
297 eDiffCrossSectionZ>>eDummy>>cumDummy;
298 if (eDummy != eEdummyVec_H2O.back())
300 eEdummyVec_H2O.push_back(eDummy);
301 eCum_H2O[eDummy].push_back(0.);
303 eDiffCrossSectionZ>>fAngleData_H2O[eDummy][cumDummy];
304 if (cumDummy != eCum_H2O[eDummy].back()){
305 eCum_H2O[eDummy].push_back(cumDummy);
307 }
while(!eDiffCrossSectionZ.eof());
310 if (verboseLevel > 2)
311 G4cout <<
"Loaded cross section files of ELSEPA Elastic model for"
316 G4cout <<
"ELSEPA elastic model is initialized " <<
G4endl
331 isInitialised =
true;
344 if (verboseLevel > 3)
347 "Calling CrossSectionPerVolume() of G4DNAELSEPAElasticModel"
358 if (material->
GetZ()!=79)
return 0.0;
365 if(atomicNDensity!= 0.0)
367 if (ekin < fhighEnergyLimit)
369 if (ekin < fkillBelowEnergy_Au)
return DBL_MAX;
371 if (ekin < 10*eV) sigma = fpData_Au->FindValue(10*eV);
372 else sigma = fpData_Au->FindValue(ekin);
375 if (verboseLevel > 2)
377 G4cout <<
"__________________________________" <<
G4endl;
378 G4cout <<
"=== G4DNAELSEPAElasticModel - XS INFO START" <<
G4endl;
379 G4cout <<
"=== Material is made of one element with Z =" << Z <<
G4endl;
380 G4cout <<
"=== Kinetic energy(eV)=" << ekin/eV <<
" particle : "
381 << particleName <<
G4endl;
382 G4cout <<
"=== Cross section per atom for Z="<<Z<<
" is (cm^2)"
384 G4cout <<
"=== Cross section per atom for Z="<<Z<<
" is (cm^-1)="
385 << sigma*atomicNDensity/(1./cm) <<
G4endl;
386 G4cout <<
"=== G4DNAELSEPAElasticModel - XS INFO END" <<
G4endl;
391 atomicNDensity = (*fpMolDensity)[material->
GetIndex()];
392 if(atomicNDensity!= 0.0)
396 sigma = fpData_H2O->FindValue(ekin);
399 if (verboseLevel > 2)
401 G4cout <<
"__________________________________" <<
G4endl;
402 G4cout <<
"=== G4DNAELSEPAElasticModel - XS INFO START" <<
G4endl;
403 G4cout <<
"=== Kinetic energy(eV)=" << ekin/eV
405 G4cout <<
"=== Cross section per water molecule (cm^2)="
407 G4cout <<
"=== Cross section per water molecule (cm^-1)="
408 << sigma*atomicNDensity/(1./cm) <<
G4endl;
409 G4cout <<
"=== G4DNAELSEPAElasticModel - XS INFO END" <<
G4endl;
413 return sigma*atomicNDensity;
419 std::vector<G4DynamicParticle*>*,
426 if (verboseLevel > 3){
428 "Calling SampleSecondaries() of G4DNAELSEPAElasticModel"
443 if (electronEnergy0 < fkillBelowEnergy_Au)
452 if(electronEnergy0>= fkillBelowEnergy_Au && electronEnergy0 < fhighEnergyLimit)
455 if (electronEnergy0>=10*eV)
457 cosTheta = RandomizeCosTheta(Z,electronEnergy0);
461 cosTheta = RandomizeCosTheta(Z,10*eV);
470 G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
472 xDir *= std::cos(phi);
473 yDir *= std::sin(phi);
475 G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
483 if(material == fpBaseWater)
486 G4double cosTheta = RandomizeCosTheta(0,electronEnergy0);
494 G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
496 xDir *= std::cos(phi);
497 yDir *= std::sin(phi);
499 G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
529 std::vector<G4double>::iterator e2;
531 e2 = std::upper_bound(eEdummyVec_H2O.begin(),
532 eEdummyVec_H2O.end(), k);
534 e2 = std::upper_bound(eEdummyVec_Au.begin(),
535 eEdummyVec_Au.end(), k);
540 std::vector<G4double>::iterator cum12;
542 cum12 = std::upper_bound(eCum_H2O[(*e1)].begin(),
543 eCum_H2O[(*e1)].end(),integrDiff);
545 cum12 = std::upper_bound(eCum_Au[(*e1)].begin(),
546 eCum_Au[(*e1)].end(),integrDiff);
549 auto cum11 = cum12 - 1;
554 std::vector<G4double>::iterator cum22;
556 cum22 = std::upper_bound(eCum_H2O[(*e2)].begin(),
557 eCum_H2O[(*e2)].end(),integrDiff);
559 cum22 = std::upper_bound(eCum_Au[(*e2)].begin(),
560 eCum_Au[(*e2)].end(),integrDiff);
563 auto cum21 = cum22 - 1;
573 a11 = fAngleData_H2O[valueE1][valuecum11];
574 a12 = fAngleData_H2O[valueE1][valuecum12];
575 a21 = fAngleData_H2O[valueE2][valuecum21];
576 a22 = fAngleData_H2O[valueE2][valuecum22];
578 a11 = fAngleData_Au[valueE1][valuecum11];
579 a12 = fAngleData_Au[valueE1][valuecum12];
580 a21 = fAngleData_Au[valueE2][valuecum21];
581 a22 = fAngleData_Au[valueE2][valuecum22];
586 if (a11 == 0 && a12 == 0 && a21 == 0 && a22 == 0)
return (0.);
588 theta = QuadInterpolator(valuecum11, valuecum12, valuecum21, valuecum22,
589 a11, a12,a21, a22, valueE1, valueE2, k, integrDiff);
603 G4double a = std::log10(e) - std::log10(e1);
604 G4double b = std::log10(e2) - std::log10(e);
605 value = xs1 + a/(a+b)*(xs2-xs1);
610 value = (d1 + (d2 - d1) * (e - e1) / (e2 - e1));
626 G4double value = std::pow(10,(d1 + (d2 - d1) * (e - e1) / (e2 - e1)));
640 G4double value = (d1 + (d2 - d1) * (e - e1) / (e2 - e1));
652 G4double a = (std::log10(xs2) - std::log10(xs1))
653 / (std::log10(e2) - std::log10(e1));
654 G4double b = std::log10(xs2) - a * std::log10(e2);
655 G4double sigma = a * std::log10(e) + b;
656 G4double value = (std::pow(10., sigma));
662G4double G4DNAELSEPAElasticModel::QuadInterpolator(
681 interpolatedvalue1 = LinLogInterpolate(cum11, cum12, cum, a11, a12);
684 interpolatedvalue1 = LinLinInterpolate(cum11, cum12, cum, a11, a12);
687 interpolatedvalue2 = LinLogInterpolate(cum21, cum22, cum, a21, a22);
690 interpolatedvalue2 = LinLinInterpolate(cum21, cum22, cum, a21, a22);
693 value = LogLinInterpolate(e1,e2,t,interpolatedvalue1,interpolatedvalue2);
705 integrdiff = uniformRand;
711 cosTheta = std::cos(theta * CLHEP::pi / 180.);
720 fkillBelowEnergy_Au = threshold;
722 if (threshold < 10 * eV)
724 G4cout<<
"*** WARNING : the G4DNAELSEPAElasticModel model is not "
725 "defined below 10 eV !" <<
G4endl;
std::vector< const G4Element * > G4ElementVector
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
CLHEP::Hep3Vector G4ThreeVector
G4GLOB_DLL std::ostream G4cout
Hep3Vector orthogonal() const
Hep3Vector cross(const Hep3Vector &) const
G4DNAELSEPAElasticModel(const G4ParticleDefinition *particle=nullptr, const G4String &nam="DNAELSEPAElasticModel")
G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *particle, G4double ekin, G4double emin, G4double emax) override
void Initialise(const G4ParticleDefinition *particle, const G4DataVector &) override
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
~G4DNAELSEPAElasticModel() override
void SetKillBelowThreshold(G4double threshold)
G4ParticleChangeForGamma * fParticleChangeForGamma
static G4DNAMolecularMaterial * Instance()
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
static G4Electron * ElectronDefinition()
static G4EmParameters * Instance()
const G4String & GetDirLEDATA() const
const G4Material * GetMaterial() const
const G4ElementVector * GetElementVector() const
const G4Material * GetBaseMaterial() const
const G4double * GetAtomicNumDensityVector() const
std::size_t GetIndex() const
std::size_t GetNumberOfElements() const
const G4String & GetName() const
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
const G4String & GetParticleName() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
void SetHighEnergyLimit(G4double)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
G4double LowEnergyLimit() const
G4double HighEnergyLimit() const
void SetLowEnergyLimit(G4double)
G4VEmModel(const G4String &nam)