Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PenelopeRayleighModelMI Class Reference

#include <G4PenelopeRayleighModelMI.hh>

+ Inheritance diagram for G4PenelopeRayleighModelMI:

Public Member Functions

 G4PenelopeRayleighModelMI (const G4ParticleDefinition *p=nullptr, const G4String &processName="PenRayleighMI")
 
virtual ~G4PenelopeRayleighModelMI ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel) override
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX) override
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0., G4double maxEnergy=DBL_MAX) override
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
void SetVerbosityLevel (G4int lev)
 
G4int GetVerbosityLevel ()
 
void DumpFormFactorTable (const G4Material *)
 
void SetMIActive (G4bool val)
 
G4bool IsMIActive ()
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerShell (const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
virtual G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectTargetAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
G4VEmModelGetTripletModel ()
 
void SetTripletModel (G4VEmModel *)
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool LPMFlag () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy) const
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetFluctuationFlag (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const G4MaterialpBaseMaterial
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 
G4double inveplus
 
G4double pFactor
 

Detailed Description

Definition at line 69 of file G4PenelopeRayleighModelMI.hh.

Constructor & Destructor Documentation

◆ G4PenelopeRayleighModelMI()

G4PenelopeRayleighModelMI::G4PenelopeRayleighModelMI ( const G4ParticleDefinition p = nullptr,
const G4String processName = "PenRayleighMI" 
)

Definition at line 68 of file G4PenelopeRayleighModelMI.cc.

69 :
70 G4VEmModel(nam),
71 fParticleChange(nullptr),fParticle(nullptr),isInitialised(false),logAtomicCrossSection(nullptr),
72 atomicFormFactor(nullptr),MolInterferenceData(nullptr),logFormFactorTable(nullptr),
73 pMaxTable(nullptr),samplingTable(nullptr),fLocalTable(false),fIsMIActive(true),
74 angularFunction(nullptr), fKnownMaterials(nullptr)
75{
76 fIntrinsicLowEnergyLimit = 100.0*eV;
77 fIntrinsicHighEnergyLimit = 100.0*GeV;
78 //SetLowEnergyLimit(fIntrinsicLowEnergyLimit);
79 SetHighEnergyLimit(fIntrinsicHighEnergyLimit);
80
81 if (part) SetParticle(part);
82
83 verboseLevel = 0;
84 // Verbosity scale:
85 // 0 = nothing
86 // 1 = warning for energy non-conservation
87 // 2 = details of energy budget
88 // 3 = calculation of FF and CS, file openings, sampling of atoms
89 // 4 = entering in methods
90
91 //build the energy grid. It is the same for all materials
92 G4double logenergy = G4Log(fIntrinsicLowEnergyLimit/2.);
93 G4double logmaxenergy = G4Log(1.5*fIntrinsicHighEnergyLimit);
94 //finer grid below 160 keV
95 G4double logtransitionenergy = G4Log(160*keV);
96 G4double logfactor1 = G4Log(10.)/250.;
97 G4double logfactor2 = logfactor1*10;
98 logEnergyGridPMax.push_back(logenergy);
99 do {
100 if (logenergy < logtransitionenergy)
101 logenergy += logfactor1;
102 else
103 logenergy += logfactor2;
104 logEnergyGridPMax.push_back(logenergy);
105 } while (logenergy < logmaxenergy);
106}
G4double G4Log(G4double x)
Definition: G4Log.hh:226
double G4double
Definition: G4Types.hh:83
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:757

◆ ~G4PenelopeRayleighModelMI()

G4PenelopeRayleighModelMI::~G4PenelopeRayleighModelMI ( )
virtual

Definition at line 110 of file G4PenelopeRayleighModelMI.cc.

111{
112 if (IsMaster() || fLocalTable) {
113
114 if (logAtomicCrossSection) {
115 for (auto& item : (*logAtomicCrossSection))
116 if (item.second) delete item.second;
117 delete logAtomicCrossSection;
118 logAtomicCrossSection = nullptr;
119 }
120
121 if (atomicFormFactor) {
122 for (auto& item : (*atomicFormFactor))
123 if (item.second) delete item.second;
124 delete atomicFormFactor;
125 atomicFormFactor = nullptr;
126 }
127
128 if (MolInterferenceData) {
129 for (auto& item : (*MolInterferenceData))
130 if (item.second) delete item.second;
131 delete MolInterferenceData;
132 MolInterferenceData = nullptr;
133 }
134
135 if (fKnownMaterials)
136 {
137 delete fKnownMaterials;
138 fKnownMaterials = nullptr;
139 }
140
141 if (angularFunction)
142 {
143 delete angularFunction;
144 angularFunction = nullptr;
145 }
146
147 ClearTables();
148
149 }
150}
G4bool IsMaster() const
Definition: G4VEmModel.hh:736

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4PenelopeRayleighModelMI::ComputeCrossSectionPerAtom ( const G4ParticleDefinition ,
G4double  kinEnergy,
G4double  Z,
G4double  A = 0,
G4double  cut = 0,
G4double  emax = DBL_MAX 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 343 of file G4PenelopeRayleighModelMI.cc.

349{
350 //Cross section of Rayleigh scattering in Penelope v2008 is calculated by the EPDL97
351 //tabulation, Cuellen et al. (1997), with non-relativistic form factors from Hubbel
352 //et al. J. Phys. Chem. Ref. Data 4 (1975) 471; Erratum ibid. 6 (1977) 615.
353
354 if (verboseLevel > 3)
355 G4cout << "Calling CrossSectionPerAtom() of G4PenelopeRayleighModelMI" << G4endl;
356
357 G4int iZ = (G4int) Z;
358
359 //Either Initialize() was not called, or we are in a slave and InitializeLocal() was
360 //not invoked
361 if (!logAtomicCrossSection) {
362 //create a **thread-local** version of the table. Used only for G4EmCalculator and
363 //Unit Tests
364 fLocalTable = true;
365 logAtomicCrossSection = new std::map<G4int,G4PhysicsFreeVector*>;
366 }
367
368 //now it should be ok
369 if (!logAtomicCrossSection->count(iZ)) {
370 //If we are here, it means that Initialize() was inkoved, but the MaterialTable was
371 //not filled up. This can happen in a UnitTest or via G4EmCalculator
372 if (verboseLevel > 0) {
373 //Issue a G4Exception (warning) only in verbose mode
375 ed << "Unable to retrieve the cross section table for Z=" << iZ << G4endl;
376 ed << "This can happen only in Unit Tests or via G4EmCalculator" << G4endl;
377 G4Exception("G4PenelopeRayleighModelMI::ComputeCrossSectionPerAtom()",
378 "em2040",JustWarning,ed);
379 }
380
381 //protect file reading via autolock
382 G4AutoLock lock(&PenelopeRayleighModelMutex);
383 ReadDataFile(iZ);
384 lock.unlock();
385
386 }
387
388 G4double cross = 0;
389
390 G4PhysicsFreeVector* atom = logAtomicCrossSection->find(iZ)->second;
391 if (!atom) {
393 ed << "Unable to find Z=" << iZ << " in the atomic cross section table" << G4endl;
394 G4Exception("G4PenelopeRayleighModelMI::ComputeCrossSectionPerAtom()",
395 "em2041",FatalException,ed);
396 return 0;
397 }
398
399 G4double logene = G4Log(energy);
400 G4double logXS = atom->Value(logene);
401 cross = G4Exp(logXS);
402
403 if (verboseLevel > 2) {
404 G4cout << "Rayleigh cross section at " << energy/keV << " keV for Z="
405 << Z << " = " << cross/barn << " barn" << G4endl;
406 }
407
408 return cross;
409}
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4double Value(G4double theEnergy, std::size_t &lastidx) const
G4double energy(const ThreeVector &p, const G4double m)

◆ CrossSectionPerVolume()

G4double G4PenelopeRayleighModelMI::CrossSectionPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy = 0.,
G4double  maxEnergy = DBL_MAX 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 450 of file G4PenelopeRayleighModelMI.cc.

455{
456 //check if we are in a Unit Test (only for the first time)
457 static G4bool amInAUnitTest = false;
458 if (G4ProductionCutsTable::GetProductionCutsTable()->GetTableSize() == 0 && !amInAUnitTest)
459 {
460 amInAUnitTest = true;
462 ed << "The ProductionCuts table is empty " << G4endl;
463 ed << "This should happen only in Unit Tests" << G4endl;
464 G4Exception("G4PenelopeRayleighModelMI::CrossSectionPerVolume()",
465 "em2019",JustWarning,ed);
466 }
467 //If the material does not have a MIFF, continue with the old-style calculation
468 G4String matname = material->GetName();
469 if (amInAUnitTest)
470 {
471 //No need to lock, as this is always called in a master
472 const G4ElementVector* theElementVector = material->GetElementVector();
473 //protect file reading via autolock
474 for (size_t j=0;j<material->GetNumberOfElements();j++) {
475 G4int iZ = theElementVector->at(j)->GetZasInt();
476 if (!logAtomicCrossSection->count(iZ)) {
477 ReadDataFile(iZ);
478 }
479 }
480 if (fIsMIActive)
481 ReadMolInterferenceData(matname);
482 if (!logFormFactorTable->count(material))
483 BuildFormFactorTable(material);
484 if (!(samplingTable->count(material)))
485 InitializeSamplingAlgorithm(material);
486 if (!pMaxTable->count(material))
487 GetPMaxTable(material);
488 }
489
490 G4bool useMIFF = fIsMIActive && (MolInterferenceData->count(matname) || matname.find("MedMat") != std::string::npos);
491 if (!useMIFF)
492 {
493 if (verboseLevel > 2)
494 G4cout << "Rayleigh CS of: " << matname << " calculated through CSperAtom!" << G4endl;
495 return G4VEmModel::CrossSectionPerVolume(material,p,energy);
496 }
497
498 // If we are here, it means that we have to integrate the cross section
499 if (verboseLevel > 2)
500 G4cout << "Rayleigh CS of: " << matname
501 << " calculated through integration of the DCS" << G4endl;
502
503
504 G4double cs = 0;
505
506 //force null cross-section if below the low-energy edge of the table
507 if (energy < LowEnergyLimit()) return cs;
508
509
510
511
512 //if the material is a CRYSTAL, forbid this process
513 if (material->IsExtended() && material->GetName() != "CustomMat") {
514 G4ExtendedMaterial* extendedmaterial = (G4ExtendedMaterial*)material;
515 G4CrystalExtension* crystalExtension = (G4CrystalExtension*)extendedmaterial->RetrieveExtension("crystal");
516 if (crystalExtension != 0) {
517 G4cout << "The material has a crystalline structure, a dedicated diffraction model is used!" << G4endl;
518 return 0;
519 }
520 }
521
522
523 //GET MATERIAL INFORMATION
524
525 G4double atomDensity = material->GetTotNbOfAtomsPerVolume();
526 G4int nElements = material->GetNumberOfElements();
527 const G4ElementVector* elementVector = material->GetElementVector();
528 const G4double* fractionVector = material->GetFractionVector();
529 //const G4double* theAtomNumDensityVector = material-> GetVecNbOfAtomsPerVolume();
530
531 //Stoichiometric factors
532 std::vector<G4double> *StoichiometricFactors = new std::vector<G4double>;
533 for (G4int i=0;i<nElements;i++) {
534 G4double fraction = fractionVector[i];
535 G4double atomicWeigth = (*elementVector)[i]->GetA()/(g/mole);
536 StoichiometricFactors->push_back(fraction/atomicWeigth);
537 }
538 G4double MaxStoichiometricFactor = 0.;
539 for (G4int i=0;i<nElements;i++) {
540 if ((*StoichiometricFactors)[i] > MaxStoichiometricFactor)
541 MaxStoichiometricFactor = (*StoichiometricFactors)[i];
542 }
543 for (G4int i=0;i<nElements;i++) {
544 (*StoichiometricFactors)[i] /= MaxStoichiometricFactor;
545 //G4cout << "StoichiometricFactors[" << i+1 << "]: " << (*StoichiometricFactors)[i] << G4endl;
546 }
547
548 //Equivalent atoms per molecule
549 G4double atPerMol = 0;
550 for (G4int i=0;i<nElements;i++)
551 atPerMol += (*StoichiometricFactors)[i];
552 G4double moleculeDensity = 0.;
553 if (atPerMol) moleculeDensity = atomDensity/atPerMol;
554
555 if (verboseLevel > 2)
556 G4cout << "Material " << material->GetName() << " has " << atPerMol << " atoms "
557 << "per molecule and " << moleculeDensity/(cm*cm*cm) << " molecule/cm3" << G4endl;
558
559 //Equivalent molecular weight (dimensionless)
560 G4double MolWeight = 0.;
561 for (G4int i=0;i<nElements;i++)
562 MolWeight += (*StoichiometricFactors)[i]*(*elementVector)[i]->GetA()/(g/mole);
563
564 if (verboseLevel > 2)
565 G4cout << "Molecular weight of " << matname << ": " << MolWeight << " g/mol" << G4endl;
566
567 G4double IntegrandFun[Ntheta];
568 for (G4int k=0; k<Ntheta; k++) {
569 G4double theta = angularFunction->Energy(k); //the x-value is called "Energy", but is an angle
570 G4double F2 = GetFSquared(material,CalculateQSquared(theta,energy));
571 IntegrandFun[k] = (*angularFunction)[k]*F2;
572 }
573
574 G4double constant = pi*classic_electr_radius*classic_electr_radius;
575 cs = constant*IntegrateFun(IntegrandFun,Ntheta,fDTheta);
576
577 //Now cs is the cross section per molecule, let's calculate the cross section per volume
578 G4double csvolume = cs*moleculeDensity;
579
580
581 //print CS and mfp
582 if (verboseLevel > 2)
583 G4cout << "Rayleigh CS of " << matname << " at " << energy/keV
584 << " keV: " << cs/barn << " barn"
585 << ", mean free path: " << 1./csvolume/mm << " mm" << G4endl;
586
587
588 delete StoichiometricFactors;
589
590 //return CS
591 return csvolume;
592}
std::vector< G4Element * > G4ElementVector
bool G4bool
Definition: G4Types.hh:86
G4VMaterialExtension * RetrieveExtension(const G4String &name)
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:207
virtual G4bool IsExtended() const
Definition: G4Material.cc:795
const G4double * GetFractionVector() const
Definition: G4Material.hh:192
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
const G4String & GetName() const
Definition: G4Material.hh:175
G4double Energy(std::size_t index) const
static G4ProductionCutsTable * GetProductionCutsTable()
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:652
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:254
const G4double pi

◆ DumpFormFactorTable()

void G4PenelopeRayleighModelMI::DumpFormFactorTable ( const G4Material mat)

Definition at line 1784 of file G4PenelopeRayleighModelMI.cc.

1785{
1786 G4cout << "*****************************************************************" << G4endl;
1787 G4cout << "G4PenelopeRayleighModelMI: Form Factor Table for " << mat->GetName() << G4endl;
1788 //try to use the same format as Penelope-Fortran, namely Q (/m_e*c) and F
1789 G4cout << "Q/(m_e*c) F(Q) " << G4endl;
1790 G4cout << "*****************************************************************" << G4endl;
1791 if (!logFormFactorTable->count(mat))
1792 BuildFormFactorTable(mat);
1793
1794 G4PhysicsFreeVector* theVec = logFormFactorTable->find(mat)->second;
1795 for (size_t i=0;i<theVec->GetVectorLength();i++)
1796 {
1797 G4double logQ2 = theVec->GetLowEdgeEnergy(i);
1798 G4double Q = G4Exp(0.5*logQ2);
1799 G4double logF2 = (*theVec)[i];
1800 G4double F = G4Exp(0.5*logF2);
1801 G4cout << Q << " " << F << G4endl;
1802 }
1803 //DONE
1804 return;
1805}
G4double GetLowEdgeEnergy(std::size_t binNumber) const
std::size_t GetVectorLength() const

◆ GetVerbosityLevel()

G4int G4PenelopeRayleighModelMI::GetVerbosityLevel ( )
inline

Definition at line 103 of file G4PenelopeRayleighModelMI.hh.

103{return verboseLevel;};

◆ Initialise()

void G4PenelopeRayleighModelMI::Initialise ( const G4ParticleDefinition part,
const G4DataVector  
)
overridevirtual

Implements G4VEmModel.

Definition at line 189 of file G4PenelopeRayleighModelMI.cc.

191{
192 if (verboseLevel > 3)
193 G4cout << "Calling G4PenelopeRayleighModelMI::Initialise()" << G4endl;
194
195 SetParticle(part);
196
197
198 if (verboseLevel)
199 G4cout << "# Molecular Interference is " << (fIsMIActive ? "ON" : "OFF") << " #" << G4endl;
200
201 //Only the master model creates/fills/destroys the tables
202 if (IsMaster() && part == fParticle) {
203 //clear tables depending on materials, not the atomic ones
204 ClearTables();
205
206 //Use here the highest verbosity, from G4EmParameter or local
207 G4int globVerb = G4EmParameters::Instance()->Verbose();
208 if (globVerb > verboseLevel)
209 {
210 verboseLevel = globVerb;
211 if (verboseLevel)
212 G4cout << "Verbosity level of G4PenelopeRayleighModelMI set to " << verboseLevel <<
213 " from G4EmParameters()" << G4endl;
214 }
215
216 if (verboseLevel > 3)
217 G4cout << "Calling G4PenelopeRayleighModelMI::Initialise() [master]" << G4endl;
218
219 //Load the list of known materials and the DCS integration grid
220 if (fIsMIActive)
221 {
222 if (!fKnownMaterials)
223 fKnownMaterials = new std::map<G4String,G4String>;
224 if (!fKnownMaterials->size())
225 LoadKnownMIFFMaterials();
226 if (!angularFunction)
227 {
228 //Create and fill once
229 angularFunction = new G4PhysicsFreeVector(Ntheta);
230 CalculateThetaAndAngFun(); //angular funtion for DCS integration
231 }
232 }
233
234 //create new tables
235 //logAtomicCrossSection and atomicFormFactor are created only once,
236 //since they are never cleared
237 if (!logAtomicCrossSection)
238 logAtomicCrossSection = new std::map<G4int,G4PhysicsFreeVector*>;
239
240 if (!atomicFormFactor)
241 atomicFormFactor = new std::map<G4int,G4PhysicsFreeVector*>;
242
243 if (fIsMIActive && !MolInterferenceData)
244 MolInterferenceData = new std::map<G4String,G4PhysicsFreeVector*>; //G. PaternĂ²
245
246 if (!logFormFactorTable)
247 logFormFactorTable = new std::map<const G4Material*,G4PhysicsFreeVector*>;
248
249 if (!pMaxTable)
250 pMaxTable = new std::map<const G4Material*,G4PhysicsFreeVector*>;
251
252 if (!samplingTable)
253 samplingTable = new std::map<const G4Material*,G4PenelopeSamplingData*>;
254
255 //loop on the used materials
257
258 for (size_t i=0;i<theCoupleTable->GetTableSize();i++) {
259 const G4Material* material =
260 theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
261 const G4ElementVector* theElementVector = material->GetElementVector();
262
263 for (size_t j=0;j<material->GetNumberOfElements();j++) {
264 G4int iZ = theElementVector->at(j)->GetZasInt();
265 //read data files only in the master
266 if (!logAtomicCrossSection->count(iZ))
267 ReadDataFile(iZ);
268 }
269
270 //1) Read MI form factors
271 if (fIsMIActive && !MolInterferenceData->count(material->GetName()))
272 ReadMolInterferenceData(material->GetName());
273
274 //2) If the table has not been built for the material, do it!
275 if (!logFormFactorTable->count(material))
276 BuildFormFactorTable(material);
277
278 //3) retrieve or build the sampling table
279 if (!(samplingTable->count(material)))
280 InitializeSamplingAlgorithm(material);
281
282 //4) retrieve or build the pMax data
283 if (!pMaxTable->count(material))
284 GetPMaxTable(material);
285 }
286
287 if (verboseLevel > 1) {
288 G4cout << G4endl << "Penelope Rayleigh model v2008 is initialized" << G4endl
289 << "Energy range: "
290 << LowEnergyLimit() / keV << " keV - "
291 << HighEnergyLimit() / GeV << " GeV"
292 << G4endl;
293 }
294
295 }
296
297 if (isInitialised)
298 return;
299 fParticleChange = GetParticleChangeForGamma();
300 isInitialised = true;
301}
static G4EmParameters * Instance()
G4int Verbose() const
const G4Material * GetMaterial() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:133
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645

◆ InitialiseLocal()

void G4PenelopeRayleighModelMI::InitialiseLocal ( const G4ParticleDefinition part,
G4VEmModel masterModel 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 305 of file G4PenelopeRayleighModelMI.cc.

307{
308 if (verboseLevel > 3)
309 G4cout << "Calling G4PenelopeRayleighModelMI::InitialiseLocal()" << G4endl;
310
311 //Check that particle matches: one might have multiple master models
312 //(e.g. for e+ and e-)
313 if (part == fParticle) {
314
315 //Get the const table pointers from the master to the workers
316 const G4PenelopeRayleighModelMI* theModel =
317 static_cast<G4PenelopeRayleighModelMI*> (masterModel);
318
319 //Copy pointers to the data tables
320 logAtomicCrossSection = theModel->logAtomicCrossSection;
321 atomicFormFactor = theModel->atomicFormFactor;
322 MolInterferenceData = theModel->MolInterferenceData;
323 logFormFactorTable = theModel->logFormFactorTable;
324 pMaxTable = theModel->pMaxTable;
325 samplingTable = theModel->samplingTable;
326 fKnownMaterials = theModel->fKnownMaterials;
327 angularFunction = theModel->angularFunction;
328
329 //Copy the G4DataVector with the grid
330 logQSquareGrid = theModel->logQSquareGrid;
331
332 //Same verbosity for all workers, as the master
333 verboseLevel = theModel->verboseLevel;
334
335 }
336
337 return;
338}

◆ IsMIActive()

G4bool G4PenelopeRayleighModelMI::IsMIActive ( )
inline

Definition at line 110 of file G4PenelopeRayleighModelMI.hh.

110{return fIsMIActive;};

◆ SampleSecondaries()

void G4PenelopeRayleighModelMI::SampleSecondaries ( std::vector< G4DynamicParticle * > *  ,
const G4MaterialCutsCouple couple,
const G4DynamicParticle aDynamicGamma,
G4double  tmin,
G4double  maxEnergy 
)
overridevirtual

Implements G4VEmModel.

Definition at line 775 of file G4PenelopeRayleighModelMI.cc.

780{
781 // Sampling of the Rayleigh final state (namely, scattering angle of the photon)
782 // from the Penelope2008 model. The scattering angle is sampled from the atomic
783 // cross section dOmega/d(cosTheta) from Born ("Atomic Phyisics", 1969), disregarding
784 // anomalous scattering effects. The Form Factor F(Q) function which appears in the
785 // analytical cross section is retrieved via the method GetFSquared(); atomic data
786 // are tabulated for F(Q). Form factor for compounds is calculated according to
787 // the additivity rule. The sampling from the F(Q) is made via a Rational Inverse
788 // Transform with Aliasing (RITA) algorithm; RITA parameters are calculated once
789 // for each material and managed by G4PenelopeSamplingData objects.
790 // The sampling algorithm (rejection method) has efficiency 67% at low energy, and
791 // increases with energy. For E=100 keV the efficiency is 100% and 86% for
792 // hydrogen and uranium, respectively.
793
794 if (verboseLevel > 3)
795 G4cout << "Calling SamplingSecondaries() of G4PenelopeRayleighModelMI" << G4endl;
796
797 G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy();
798
799 if (photonEnergy0 <= fIntrinsicLowEnergyLimit) {
800 fParticleChange->ProposeTrackStatus(fStopAndKill);
801 fParticleChange->SetProposedKineticEnergy(0.);
802 fParticleChange->ProposeLocalEnergyDeposit(photonEnergy0);
803 return;
804 }
805
806 G4ParticleMomentum photonDirection0 = aDynamicGamma->GetMomentumDirection();
807
808 const G4Material* theMat = couple->GetMaterial();
809
810 //1) Verify if tables are ready
811 //Either Initialize() was not called, or we are in a slave and InitializeLocal() was
812 //not invoked
813 if (!pMaxTable || !samplingTable || !logAtomicCrossSection || !atomicFormFactor ||
814 !logFormFactorTable) {
815 //create a **thread-local** version of the table. Used only for G4EmCalculator and
816 //Unit Tests
817 fLocalTable = true;
818 if (!logAtomicCrossSection)
819 logAtomicCrossSection = new std::map<G4int,G4PhysicsFreeVector*>;
820 if (!atomicFormFactor)
821 atomicFormFactor = new std::map<G4int,G4PhysicsFreeVector*>;
822 if (!logFormFactorTable)
823 logFormFactorTable = new std::map<const G4Material*,G4PhysicsFreeVector*>;
824 if (!pMaxTable)
825 pMaxTable = new std::map<const G4Material*,G4PhysicsFreeVector*>;
826 if (!samplingTable)
827 samplingTable = new std::map<const G4Material*,G4PenelopeSamplingData*>;
828 if (fIsMIActive && !MolInterferenceData)
829 MolInterferenceData = new std::map<G4String,G4PhysicsFreeVector*>;
830 }
831
832 if (!samplingTable->count(theMat)) {
833 //If we are here, it means that Initialize() was inkoved, but the MaterialTable was
834 //not filled up. This can happen in a UnitTest
835 if (verboseLevel > 0) {
836 //Issue a G4Exception (warning) only in verbose mode
838 ed << "Unable to find the samplingTable data for " <<
839 theMat->GetName() << G4endl;
840 ed << "This can happen only in Unit Tests" << G4endl;
841 G4Exception("G4PenelopeRayleighModelMI::SampleSecondaries()",
842 "em2019",JustWarning,ed);
843 }
844 const G4ElementVector* theElementVector = theMat->GetElementVector();
845 //protect file reading via autolock
846 G4AutoLock lock(&PenelopeRayleighModelMutex);
847 for (size_t j=0;j<theMat->GetNumberOfElements();j++) {
848 G4int iZ = theElementVector->at(j)->GetZasInt();
849 if (!logAtomicCrossSection->count(iZ)) {
850 lock.lock();
851 ReadDataFile(iZ);
852 lock.unlock();
853 }
854 }
855 lock.lock();
856
857 //1) If the table has not been built for the material, do it!
858 if (!logFormFactorTable->count(theMat))
859 BuildFormFactorTable(theMat);
860
861 //2) retrieve or build the sampling table
862 if (!(samplingTable->count(theMat)))
863 InitializeSamplingAlgorithm(theMat);
864
865 //3) retrieve or build the pMax data
866 if (!pMaxTable->count(theMat))
867 GetPMaxTable(theMat);
868
869 lock.unlock();
870 }
871
872 //Ok, restart the job
873 G4PenelopeSamplingData* theDataTable = samplingTable->find(theMat)->second;
874 G4PhysicsFreeVector* thePMax = pMaxTable->find(theMat)->second;
875
876 G4double cosTheta = 1.0;
877
878 //OK, ready to go!
879 G4double qmax = 2.0*photonEnergy0/electron_mass_c2; //this is non-dimensional now
880
881 if (qmax < 1e-10) { //very low momentum transfer
882 G4bool loopAgain=false;
883 do {
884 loopAgain = false;
885 cosTheta = 1.0-2.0*G4UniformRand();
886 G4double G = 0.5*(1+cosTheta*cosTheta);
887 if (G4UniformRand()>G)
888 loopAgain = true;
889 } while(loopAgain);
890 }
891
892 else { //larger momentum transfer
893
894 size_t nData = theDataTable->GetNumberOfStoredPoints();
895 G4double LastQ2inTheTable = theDataTable->GetX(nData-1);
896 G4double q2max = std::min(qmax*qmax,LastQ2inTheTable);
897
898 G4bool loopAgain = false;
899 G4double MaxPValue = thePMax->Value(photonEnergy0);
900 G4double xx=0;
901
902 //Sampling by rejection method. The rejection function is
903 //G = 0.5*(1+cos^2(theta))
904 do {
905 loopAgain = false;
906 G4double RandomMax = G4UniformRand()*MaxPValue;
907 xx = theDataTable->SampleValue(RandomMax);
908
909 //xx is a random value of q^2 in (0,q2max),sampled according to
910 //F(Q^2) via the RITA algorithm
911 if (xx > q2max)
912 loopAgain = true;
913 cosTheta = 1.0-2.0*xx/q2max;
914 G4double G = 0.5*(1+cosTheta*cosTheta);
915 if (G4UniformRand()>G)
916 loopAgain = true;
917 } while(loopAgain);
918 }
919
920 G4double sinTheta = std::sqrt(1-cosTheta*cosTheta);
921
922 //Scattered photon angles. ( Z - axis along the parent photon)
923 G4double phi = twopi * G4UniformRand() ;
924 G4double dirX = sinTheta*std::cos(phi);
925 G4double dirY = sinTheta*std::sin(phi);
926 G4double dirZ = cosTheta;
927
928 //Update G4VParticleChange for the scattered photon
929 G4ThreeVector photonDirection1(dirX, dirY, dirZ);
930
931 photonDirection1.rotateUz(photonDirection0);
932 fParticleChange->ProposeMomentumDirection(photonDirection1) ;
933 fParticleChange->SetProposedKineticEnergy(photonEnergy0) ;
934
935 return;
936}
@ fStopAndKill
#define G4UniformRand()
Definition: Randomize.hh:52
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4double GetX(size_t index)
G4double SampleValue(G4double rndm)
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

◆ SetMIActive()

void G4PenelopeRayleighModelMI::SetMIActive ( G4bool  val)
inline

Definition at line 109 of file G4PenelopeRayleighModelMI.hh.

109{fIsMIActive = val;};

◆ SetVerbosityLevel()

void G4PenelopeRayleighModelMI::SetVerbosityLevel ( G4int  lev)
inline

Definition at line 102 of file G4PenelopeRayleighModelMI.hh.

102{verboseLevel = lev;};

The documentation for this class was generated from the following files: