Geant4-11
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes | Static Private Attributes
G4PenelopeRayleighModel Class Reference

#include <G4PenelopeRayleighModel.hh>

Inheritance diagram for G4PenelopeRayleighModel:
G4VEmModel

Public Member Functions

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

Protected Member Functions

const G4MaterialCutsCoupleCurrentCouple () const
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
void SetCurrentElement (const G4Element *)
 

Protected Attributes

size_t basedCoupleIndex = 0
 
size_t currentCoupleIndex = 0
 
G4ElementDatafElementData = nullptr
 
const G4ParticleDefinitionfParticle
 
G4ParticleChangeForGammafParticleChange
 
G4double inveplus
 
G4bool lossFlucFlag = true
 
const G4MaterialpBaseMaterial = nullptr
 
G4double pFactor = 1.0
 
G4VParticleChangepParticleChange = nullptr
 
const std::vector< G4double > * theDensityFactor = nullptr
 
const std::vector< G4int > * theDensityIdx = nullptr
 
G4PhysicsTablexSectionTable = nullptr
 

Private Member Functions

void BuildFormFactorTable (const G4Material *)
 
void ClearTables ()
 
G4double GetFSquared (const G4Material *, const G4double)
 
void GetPMaxTable (const G4Material *)
 
void InitializeSamplingAlgorithm (const G4Material *)
 
void ReadDataFile (G4int)
 
void SetParticle (const G4ParticleDefinition *)
 

Private Attributes

G4VEmAngularDistributionanglModel = nullptr
 
std::vector< G4EmElementSelector * > * elmSelectors = nullptr
 
G4double eMaxActive = DBL_MAX
 
G4double eMinActive = 0.0
 
const G4MaterialCutsCouplefCurrentCouple = nullptr
 
const G4ElementfCurrentElement = nullptr
 
const G4IsotopefCurrentIsotope = nullptr
 
G4LossTableManagerfEmManager
 
G4double fIntrinsicHighEnergyLimit
 
G4double fIntrinsicLowEnergyLimit
 
G4bool fIsInitialised
 
G4bool flagDeexcitation = false
 
G4bool flagForceBuildTable = false
 
G4bool fLocalTable
 
G4DataVector fLogEnergyGridPMax
 
std::map< const G4Material *, G4PhysicsFreeVector * > * fLogFormFactorTable
 
G4DataVector fLogQSquareGrid
 
G4VEmFluctuationModelflucModel = nullptr
 
std::map< const G4Material *, G4PhysicsFreeVector * > * fPMaxTable
 
std::map< const G4Material *, G4PenelopeSamplingData * > * fSamplingTable
 
G4VEmModelfTripletModel = nullptr
 
G4int fVerboseLevel
 
G4double highLimit
 
G4bool isLocked = false
 
G4bool isMaster = true
 
G4bool localElmSelectors = true
 
G4bool localTable = true
 
G4double lowLimit
 
const G4String name
 
G4int nsec = 5
 
G4int nSelectors = 0
 
G4double polarAngleLimit
 
G4double secondaryThreshold = DBL_MAX
 
G4bool theLPMflag = false
 
G4bool useAngularGenerator = false
 
G4bool useBaseMaterials = false
 
std::vector< G4doublexsec
 

Static Private Attributes

static G4PhysicsFreeVectorfAtomicFormFactor [fMaxZ+1] = {nullptr}
 
static G4PhysicsFreeVectorfLogAtomicCrossSection [fMaxZ+1] = {nullptr}
 
static const G4int fMaxZ =99
 

Detailed Description

Definition at line 58 of file G4PenelopeRayleighModel.hh.

Constructor & Destructor Documentation

◆ G4PenelopeRayleighModel() [1/2]

G4PenelopeRayleighModel::G4PenelopeRayleighModel ( const G4ParticleDefinition p = nullptr,
const G4String processName = "PenRayleigh" 
)
explicit

Definition at line 59 of file G4PenelopeRayleighModel.cc.

61 :G4VEmModel(nam),fParticleChange(nullptr),fParticle(nullptr),
62 fLogFormFactorTable(nullptr),fPMaxTable(nullptr),fSamplingTable(nullptr),
63 fIsInitialised(false),fLocalTable(false)
64{
67 // SetLowEnergyLimit(fIntrinsicLowEnergyLimit);
69
70 if (part)
71 SetParticle(part);
72
73 //
75 // Verbosity scale:
76 // 0 = nothing
77 // 1 = warning for energy non-conservation
78 // 2 = details of energy budget
79 // 3 = calculation of cross sections, file openings, sampling of atoms
80 // 4 = entering in methods
81
82 //build the energy grid. It is the same for all materials
84 G4double logmaxenergy = G4Log(1.5*fIntrinsicHighEnergyLimit);
85 //finer grid below 160 keV
86 G4double logtransitionenergy = G4Log(160*keV);
87 G4double logfactor1 = G4Log(10.)/250.;
88 G4double logfactor2 = logfactor1*10;
89 fLogEnergyGridPMax.push_back(logenergy);
90 do{
91 if (logenergy < logtransitionenergy)
92 logenergy += logfactor1;
93 else
94 logenergy += logfactor2;
95 fLogEnergyGridPMax.push_back(logenergy);
96 }while (logenergy < logmaxenergy);
97}
G4double G4Log(G4double x)
Definition: G4Log.hh:226
static constexpr double keV
Definition: G4SIunits.hh:202
static constexpr double eV
Definition: G4SIunits.hh:201
static constexpr double GeV
Definition: G4SIunits.hh:203
double G4double
Definition: G4Types.hh:83
std::map< const G4Material *, G4PenelopeSamplingData * > * fSamplingTable
const G4ParticleDefinition * fParticle
G4ParticleChangeForGamma * fParticleChange
std::map< const G4Material *, G4PhysicsFreeVector * > * fLogFormFactorTable
void SetParticle(const G4ParticleDefinition *)
std::map< const G4Material *, G4PhysicsFreeVector * > * fPMaxTable
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:767
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:66

References eV, fIntrinsicHighEnergyLimit, fIntrinsicLowEnergyLimit, fLogEnergyGridPMax, fVerboseLevel, G4Log(), GeV, keV, G4VEmModel::SetHighEnergyLimit(), and SetParticle().

◆ ~G4PenelopeRayleighModel()

G4PenelopeRayleighModel::~G4PenelopeRayleighModel ( )
virtual

Definition at line 101 of file G4PenelopeRayleighModel.cc.

102{
103 if (IsMaster() || fLocalTable)
104 {
105
106 for(G4int i=0; i<=fMaxZ; ++i)
107 {
109 {
110 delete fLogAtomicCrossSection[i];
111 fLogAtomicCrossSection[i] = nullptr;
112 }
113 if(fAtomicFormFactor[i])
114 {
115 delete fAtomicFormFactor[i];
116 fAtomicFormFactor[i] = nullptr;
117 }
118 }
119 ClearTables();
120 }
121}
int G4int
Definition: G4Types.hh:85
static G4PhysicsFreeVector * fLogAtomicCrossSection[fMaxZ+1]
static G4PhysicsFreeVector * fAtomicFormFactor[fMaxZ+1]
G4bool IsMaster() const
Definition: G4VEmModel.hh:746

References ClearTables(), fAtomicFormFactor, fLocalTable, fLogAtomicCrossSection, fMaxZ, and G4VEmModel::IsMaster().

◆ G4PenelopeRayleighModel() [2/2]

G4PenelopeRayleighModel::G4PenelopeRayleighModel ( const G4PenelopeRayleighModel )
delete

Member Function Documentation

◆ BuildFormFactorTable()

void G4PenelopeRayleighModel::BuildFormFactorTable ( const G4Material material)
private

Definition at line 320 of file G4PenelopeRayleighModel.cc.

321{
322 /*
323 1) get composition and equivalent molecular density
324 */
325 G4int nElements = material->GetNumberOfElements();
326 const G4ElementVector* elementVector = material->GetElementVector();
327 const G4double* fractionVector = material->GetFractionVector();
328
329 std::vector<G4double> *StechiometricFactors = new std::vector<G4double>;
330 for (G4int i=0;i<nElements;i++)
331 {
332 G4double fraction = fractionVector[i];
333 G4double atomicWeigth = (*elementVector)[i]->GetA()/(g/mole);
334 StechiometricFactors->push_back(fraction/atomicWeigth);
335 }
336 //Find max
337 G4double MaxStechiometricFactor = 0.;
338 for (G4int i=0;i<nElements;i++)
339 {
340 if ((*StechiometricFactors)[i] > MaxStechiometricFactor)
341 MaxStechiometricFactor = (*StechiometricFactors)[i];
342 }
343 if (MaxStechiometricFactor<1e-16)
344 {
346 ed << "Inconsistent data of atomic composition for " <<
347 material->GetName() << G4endl;
348 G4Exception("G4PenelopeRayleighModel::BuildFormFactorTable()",
349 "em2042",FatalException,ed);
350 }
351 //Normalize
352 for (G4int i=0;i<nElements;i++)
353 (*StechiometricFactors)[i] /= MaxStechiometricFactor;
354
355 /*
356 CREATE THE FORM FACTOR TABLE
357 */
358 G4PhysicsFreeVector* theFFVec = new G4PhysicsFreeVector(fLogQSquareGrid.size(),/*spline=*/true);
359
360 for (size_t k=0;k<fLogQSquareGrid.size();k++)
361 {
362 G4double ff2 = 0; //squared form factor
363 for (G4int i=0;i<nElements;i++)
364 {
365 G4int iZ = (*elementVector)[i]->GetZasInt();
366 G4PhysicsFreeVector* theAtomVec = fAtomicFormFactor[iZ];
367 G4double f = (*theAtomVec)[k]; //the q-grid is always the same
368 ff2 += f*f*(*StechiometricFactors)[i];
369 }
370 if (ff2)
371 theFFVec->PutValue(k,fLogQSquareGrid[k],G4Log(ff2)); //NOTICE: THIS IS log(Q^2) vs. log(F^2)
372 }
373 theFFVec->FillSecondDerivatives(); //vector is filled!
374 fLogFormFactorTable->insert(std::make_pair(material,theFFVec));
375
376 delete StechiometricFactors;
377 return;
378}
std::vector< const G4Element * > G4ElementVector
@ 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
static constexpr double mole
Definition: G4SIunits.hh:279
static constexpr double g
Definition: G4SIunits.hh:168
#define G4endl
Definition: G4ios.hh:57
void PutValue(const std::size_t index, const G4double e, const G4double value)
void FillSecondDerivatives(const G4SplineType=G4SplineType::Base, const G4double dir1=0.0, const G4double dir2=0.0)
string material
Definition: eplot.py:19

References FatalException, fAtomicFormFactor, G4PhysicsVector::FillSecondDerivatives(), fLogFormFactorTable, fLogQSquareGrid, g, G4endl, G4Exception(), G4Log(), eplot::material, mole, and G4PhysicsFreeVector::PutValue().

Referenced by DumpFormFactorTable(), Initialise(), and SampleSecondaries().

◆ ChargeSquareRatio()

G4double G4VEmModel::ChargeSquareRatio ( const G4Track track)
virtualinherited

Reimplemented in G4BraggIonGasModel, and G4BetheBlochIonGasModel.

Definition at line 374 of file G4VEmModel.cc.

375{
377 track.GetMaterial(), track.GetKineticEnergy());
378}
const G4ParticleDefinition * GetParticleDefinition() const
G4Material * GetMaterial() const
G4double GetKineticEnergy() const
virtual G4double GetChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:382

References G4VEmModel::GetChargeSquareRatio(), G4Track::GetKineticEnergy(), G4Track::GetMaterial(), and G4Track::GetParticleDefinition().

Referenced by G4VEnergyLossProcess::PostStepGetPhysicalInteractionLength().

◆ ClearTables()

void G4PenelopeRayleighModel::ClearTables ( )
private

Definition at line 124 of file G4PenelopeRayleighModel.cc.

125{
127 {
128 for (auto& item : (*fLogFormFactorTable))
129 if (item.second) delete item.second;
130 delete fLogFormFactorTable;
131 fLogFormFactorTable = nullptr; //zero explicitly
132 }
133 if (fPMaxTable)
134 {
135 for (auto& item : (*fPMaxTable))
136 if (item.second) delete item.second;
137 delete fPMaxTable;
138 fPMaxTable = nullptr; //zero explicitly
139 }
140 if (fSamplingTable)
141 {
142 for (auto& item : (*fSamplingTable))
143 if (item.second) delete item.second;
144 delete fSamplingTable;
145 fSamplingTable = nullptr; //zero explicitly
146 }
147 return;
148}

References fLogFormFactorTable, fPMaxTable, and fSamplingTable.

Referenced by Initialise(), and ~G4PenelopeRayleighModel().

◆ ComputeCrossSectionPerAtom() [1/2]

G4double G4VEmModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition part,
const G4Element elm,
G4double  kinEnergy,
G4double  cutEnergy = 0.0,
G4double  maxEnergy = DBL_MAX 
)
inlineinherited

Definition at line 566 of file G4VEmModel.hh.

571{
573 return ComputeCrossSectionPerAtom(part,kinEnergy,elm->GetZ(),elm->GetN(),
574 cutEnergy,maxEnergy);
575}
G4double GetZ() const
Definition: G4Element.hh:131
G4double GetN() const
Definition: G4Element.hh:135
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:341
void SetCurrentElement(const G4Element *)
Definition: G4VEmModel.hh:497

References G4VEmModel::ComputeCrossSectionPerAtom(), G4Element::GetN(), G4Element::GetZ(), and G4VEmModel::SetCurrentElement().

◆ ComputeCrossSectionPerAtom() [2/2]

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

Reimplemented from G4VEmModel.

Definition at line 261 of file G4PenelopeRayleighModel.cc.

267{
268 // Cross section of Rayleigh scattering in Penelope v2008 is calculated by the EPDL97
269 // tabulation, Cuellen et al. (1997), with non-relativistic form factors from Hubbel
270 // et al. J. Phys. Chem. Ref. Data 4 (1975) 471; Erratum ibid. 6 (1977) 615.
271
272 if (fVerboseLevel > 3)
273 G4cout << "Calling CrossSectionPerAtom() of G4PenelopeRayleighModel" << G4endl;
274
275 G4int iZ = G4int(Z);
276
277 if (!fLogAtomicCrossSection[iZ])
278 {
279 //If we are here, it means that Initialize() was inkoved, but the MaterialTable was
280 //not filled up. This can happen in a UnitTest or via G4EmCalculator
281 if (fVerboseLevel > 0)
282 {
283 //Issue a G4Exception (warning) only in verbose mode
285 ed << "Unable to retrieve the cross section table for Z=" << iZ << G4endl;
286 ed << "This can happen only in Unit Tests or via G4EmCalculator" << G4endl;
287 G4Exception("G4PenelopeRayleighModel::ComputeCrossSectionPerAtom()",
288 "em2040",JustWarning,ed);
289 }
290 //protect file reading via autolock
292 ReadDataFile(iZ);
293 lock.unlock();
294 }
295
296 G4double cross = 0;
298 if (!atom)
299 {
301 ed << "Unable to find Z=" << iZ << " in the atomic cross section table" << G4endl;
302 G4Exception("G4PenelopeRayleighModel::ComputeCrossSectionPerAtom()",
303 "em2041",FatalException,ed);
304 return 0;
305 }
306 G4double logene = G4Log(energy);
307 G4double logXS = atom->Value(logene);
308 cross = G4Exp(logXS);
309
310 if (fVerboseLevel > 2)
311 {
312 G4cout << "Rayleigh cross section at " << energy/keV << " keV for Z=" << Z <<
313 " = " << cross/barn << " barn" << G4endl;
314 }
315 return cross;
316}
@ JustWarning
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
static constexpr double barn
Definition: G4SIunits.hh:85
const G4int Z[17]
G4GLOB_DLL std::ostream G4cout
G4double Value(const G4double energy, std::size_t &lastidx) const
G4double energy(const ThreeVector &p, const G4double m)

References barn, G4INCL::KinematicsUtils::energy(), FatalException, fLogAtomicCrossSection, fVerboseLevel, G4cout, G4endl, G4Exception(), G4Exp(), G4Log(), JustWarning, keV, anonymous_namespace{G4PenelopeRayleighModel.cc}::PenelopeRayleighModelMutex, ReadDataFile(), G4TemplateAutoLock< _Mutex_t >::unlock(), G4PhysicsVector::Value(), and Z.

◆ ComputeCrossSectionPerShell()

G4double G4VEmModel::ComputeCrossSectionPerShell ( const G4ParticleDefinition ,
G4int  Z,
G4int  shellIdx,
G4double  kinEnergy,
G4double  cutEnergy = 0.0,
G4double  maxEnergy = DBL_MAX 
)
virtualinherited

Definition at line 351 of file G4VEmModel.cc.

354{
355 return 0.0;
356}

Referenced by G4EmCalculator::ComputeCrossSectionPerShell().

◆ ComputeDEDX()

G4double G4VEmModel::ComputeDEDX ( const G4MaterialCutsCouple couple,
const G4ParticleDefinition part,
G4double  kineticEnergy,
G4double  cutEnergy = DBL_MAX 
)
inlineinherited

Definition at line 528 of file G4VEmModel.hh.

532{
533 SetCurrentCouple(couple);
534 return pFactor*ComputeDEDXPerVolume(pBaseMaterial,part,kinEnergy,cutEnergy);
535}
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:472
const G4Material * pBaseMaterial
Definition: G4VEmModel.hh:427
G4double pFactor
Definition: G4VEmModel.hh:432
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
Definition: G4VEmModel.cc:228

References G4VEmModel::ComputeDEDXPerVolume(), G4VEmModel::pBaseMaterial, G4VEmModel::pFactor, and G4VEmModel::SetCurrentCouple().

◆ ComputeDEDXPerVolume()

G4double G4VEmModel::ComputeDEDXPerVolume ( const G4Material ,
const G4ParticleDefinition ,
G4double  kineticEnergy,
G4double  cutEnergy = DBL_MAX 
)
virtualinherited

◆ ComputeMeanFreePath()

G4double G4VEmModel::ComputeMeanFreePath ( const G4ParticleDefinition part,
G4double  kineticEnergy,
const G4Material material,
G4double  cutEnergy = 0.0,
G4double  maxEnergy = DBL_MAX 
)
inlineinherited

Definition at line 553 of file G4VEmModel.hh.

558{
559 G4double cross = CrossSectionPerVolume(material,part,ekin,emin,emax);
560 return (cross > 0.0) ? 1./cross : DBL_MAX;
561}
static const G4double emax
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:237
#define DBL_MAX
Definition: templates.hh:62

References G4VEmModel::CrossSectionPerVolume(), DBL_MAX, emax, and eplot::material.

◆ CorrectionsAlongStep()

void G4VEmModel::CorrectionsAlongStep ( const G4MaterialCutsCouple ,
const G4DynamicParticle ,
const G4double length,
G4double eloss 
)
virtualinherited

◆ CrossSection()

G4double G4VEmModel::CrossSection ( const G4MaterialCutsCouple couple,
const G4ParticleDefinition part,
G4double  kineticEnergy,
G4double  cutEnergy = 0.0,
G4double  maxEnergy = DBL_MAX 
)
inlineinherited

◆ CrossSectionPerVolume()

G4double G4VEmModel::CrossSectionPerVolume ( const G4Material mat,
const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy = 0.0,
G4double  maxEnergy = DBL_MAX 
)
virtualinherited

Reimplemented in G4LivermorePhotoElectricModel, G4PAIModel, G4PAIPhotModel, G4BetheBlochNoDeltaModel, G4BraggNoDeltaModel, G4eeToHadronsModel, G4eeToHadronsMultiModel, G4ICRU73NoDeltaModel, G4MuBetheBlochModel, G4AtimaEnergyLossModel, G4BetheBlochModel, G4BraggIonModel, G4BraggModel, G4ICRU73QOModel, G4LindhardSorensenIonModel, G4MollerBhabhaModel, G4PEEffectFluoModel, G4PenelopeRayleighModelMI, G4LEPTSAttachmentModel, G4LEPTSDissociationModel, G4LEPTSElasticModel, G4LEPTSExcitationModel, G4LEPTSIonisationModel, G4LEPTSPositroniumModel, G4LEPTSRotExcitationModel, G4LEPTSVibExcitationModel, G4eplusTo2GammaOKVIModel, G4eplusTo3GammaOKVIModel, G4PenelopeComptonModel, G4eeToTwoGammaModel, G4GoudsmitSaundersonMscModel, G4IonParametrisedLossModel, G4DNABornExcitationModel1, G4DNABornExcitationModel2, G4DNABornIonisationModel1, G4DNABornIonisationModel2, G4DNAChampionElasticModel, G4DNACPA100ElasticModel, G4DNACPA100ExcitationModel, G4DNACPA100IonisationModel, G4DNADingfelderChargeDecreaseModel, G4DNADingfelderChargeIncreaseModel, G4DNADiracRMatrixExcitationModel, G4DNAEmfietzoglouExcitationModel, G4DNAEmfietzoglouIonisationModel, G4DNAIonElasticModel, G4DNAMeltonAttachmentModel, G4DNAMillerGreenExcitationModel, G4DNAModelInterface, G4TDNAOneStepThermalizationModel< MODEL >, G4DNAQuinnPlasmonExcitationModel, G4DNARelativisticIonisationModel, G4DNARuddIonisationExtendedModel, G4DNARuddIonisationModel, G4DNASancheExcitationModel, G4DNAScreenedRutherfordElasticModel, G4DNATransformElectronModel, G4DNAUeharaScreenedRutherfordElasticModel, G4MicroElecElasticModel, G4MicroElecElasticModel_new, G4MicroElecInelasticModel, G4MicroElecInelasticModel_new, G4MicroElecLOPhononModel, G4DNAELSEPAElasticModel, G4PenelopeBremsstrahlungModel, and G4PenelopeIonisationModel.

Definition at line 237 of file G4VEmModel.cc.

242{
243 SetupForMaterial(p, mat, ekin);
244 const G4double* theAtomNumDensityVector = mat->GetVecNbOfAtomsPerVolume();
245 G4int nelm = mat->GetNumberOfElements();
246 if(nelm > nsec) {
247 xsec.resize(nelm);
248 nsec = nelm;
249 }
250 G4double cross = 0.0;
251 for (G4int i=0; i<nelm; ++i) {
252 cross += theAtomNumDensityVector[i]*
253 ComputeCrossSectionPerAtom(p,mat->GetElement(i),ekin,emin,emax);
254 xsec[i] = cross;
255 }
256 return cross;
257}
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:198
size_t GetNumberOfElements() const
Definition: G4Material.hh:182
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:202
std::vector< G4double > xsec
Definition: G4VEmModel.hh:466
G4int nsec
Definition: G4VEmModel.hh:444
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:440

References G4VEmModel::ComputeCrossSectionPerAtom(), emax, G4Material::GetElement(), G4Material::GetNumberOfElements(), G4Material::GetVecNbOfAtomsPerVolume(), G4VEmModel::nsec, G4VEmModel::SetupForMaterial(), and G4VEmModel::xsec.

Referenced by G4AdjointBremsstrahlungModel::AdjointCrossSection(), G4EmCalculator::ComputeCrossSectionPerVolume(), G4VEmProcess::ComputeCurrentLambda(), G4VEmModel::ComputeMeanFreePath(), G4TablesForExtrapolator::ComputeTrasportXS(), G4VEmModel::CrossSection(), G4LivermorePhotoElectricModel::CrossSectionPerVolume(), G4PenelopeRayleighModelMI::CrossSectionPerVolume(), G4DNADummyModel::CrossSectionPerVolume(), G4VEnergyLossProcess::CrossSectionPerVolume(), G4VEmAdjointModel::DiffCrossSectionPerVolumePrimToSecond(), G4VMscModel::GetTransportMeanFreePath(), G4VEmModel::SelectRandomAtom(), and G4VEmModel::Value().

◆ CurrentCouple()

const G4MaterialCutsCouple * G4VEmModel::CurrentCouple ( ) const
inlineprotectedinherited

◆ DeexcitationFlag()

G4bool G4VEmModel::DeexcitationFlag ( ) const
inlineinherited

Definition at line 704 of file G4VEmModel.hh.

705{
706 return flagDeexcitation;
707}
G4bool flagDeexcitation
Definition: G4VEmModel.hh:455

References G4VEmModel::flagDeexcitation.

Referenced by G4EmModelManager::DumpModelList().

◆ DefineForRegion()

void G4VEmModel::DefineForRegion ( const G4Region )
virtualinherited

Reimplemented in G4PAIModel, and G4PAIPhotModel.

Definition at line 360 of file G4VEmModel.cc.

361{}

Referenced by G4EmModelManager::AddEmModel().

◆ DumpFormFactorTable()

void G4PenelopeRayleighModel::DumpFormFactorTable ( const G4Material mat)

Definition at line 1264 of file G4PenelopeRayleighModel.cc.

1265{
1266 G4cout << "*****************************************************************" << G4endl;
1267 G4cout << "G4PenelopeRayleighModel: Form Factor Table for " << mat->GetName() << G4endl;
1268 //try to use the same format as Penelope-Fortran, namely Q (/m_e*c) and F
1269 G4cout << "Q/(m_e*c) F(Q) " << G4endl;
1270 G4cout << "*****************************************************************" << G4endl;
1271 if (!fLogFormFactorTable->count(mat))
1273
1274 G4PhysicsFreeVector* theVec = fLogFormFactorTable->find(mat)->second;
1275 for (size_t i=0;i<theVec->GetVectorLength();i++)
1276 {
1277 G4double logQ2 = theVec->GetLowEdgeEnergy(i);
1278 G4double Q = G4Exp(0.5*logQ2);
1279 G4double logF2 = (*theVec)[i];
1280 G4double F = G4Exp(0.5*logF2);
1281 G4cout << Q << " " << F << G4endl;
1282 }
1283 //DONE
1284 return;
1285}
const G4String & GetName() const
Definition: G4Material.hh:173
void BuildFormFactorTable(const G4Material *)
G4double GetLowEdgeEnergy(const std::size_t index) const
std::size_t GetVectorLength() const
static double Q[]

References BuildFormFactorTable(), fLogFormFactorTable, G4cout, G4endl, G4Exp(), G4PhysicsVector::GetLowEdgeEnergy(), G4Material::GetName(), G4PhysicsVector::GetVectorLength(), and Q.

◆ FillNumberOfSecondaries()

void G4VEmModel::FillNumberOfSecondaries ( G4int numberOfTriplets,
G4int numberOfRecoil 
)
virtualinherited

Definition at line 365 of file G4VEmModel.cc.

367{
368 numberOfTriplets = 0;
369 numberOfRecoil = 0;
370}

Referenced by G4VEmProcess::PostStepDoIt(), and G4VEnergyLossProcess::PostStepDoIt().

◆ ForceBuildTableFlag()

G4bool G4VEmModel::ForceBuildTableFlag ( ) const
inlineinherited

Definition at line 711 of file G4VEmModel.hh.

712{
713 return flagForceBuildTable;
714}
G4bool flagForceBuildTable
Definition: G4VEmModel.hh:456

References G4VEmModel::flagForceBuildTable.

Referenced by G4VMscModel::GetParticleChangeForMSC().

◆ GetAngularDistribution()

G4VEmAngularDistribution * G4VEmModel::GetAngularDistribution ( )
inlineinherited

Definition at line 621 of file G4VEmModel.hh.

622{
623 return anglModel;
624}
G4VEmAngularDistribution * anglModel
Definition: G4VEmModel.hh:414

References G4VEmModel::anglModel.

Referenced by G4EmModelManager::DumpModelList(), G4AtimaEnergyLossModel::Initialise(), G4BetheBlochModel::Initialise(), G4BraggIonModel::Initialise(), G4BraggModel::Initialise(), G4ICRU73QOModel::Initialise(), G4LindhardSorensenIonModel::Initialise(), G4MollerBhabhaModel::Initialise(), G4AdjointBremsstrahlungModel::RapidSampleSecondaries(), G4LivermoreBremsstrahlungModel::SampleSecondaries(), G4eBremParametrizedModel::SampleSecondaries(), G4eBremsstrahlungRelModel::SampleSecondaries(), G4SeltzerBergerModel::SampleSecondaries(), G4DNABornIonisationModel1::SampleSecondaries(), G4DNABornIonisationModel2::SampleSecondaries(), G4DNAEmfietzoglouIonisationModel::SampleSecondaries(), G4DNARuddIonisationExtendedModel::SampleSecondaries(), G4DNARuddIonisationModel::SampleSecondaries(), G4PAIModel::SampleSecondaries(), G4PAIPhotModel::SampleSecondaries(), G4LivermoreIonisationModel::SampleSecondaries(), G4LivermorePhotoElectricModel::SampleSecondaries(), G4LivermoreRayleighModel::SampleSecondaries(), G4MicroElecInelasticModel::SampleSecondaries(), G4MicroElecInelasticModel_new::SampleSecondaries(), G4MuBremsstrahlungModel::SampleSecondaries(), G4MuPairProductionModel::SampleSecondaries(), G4AtimaEnergyLossModel::SampleSecondaries(), G4BetheBlochModel::SampleSecondaries(), G4BetheHeitlerModel::SampleSecondaries(), G4BraggIonModel::SampleSecondaries(), G4BraggModel::SampleSecondaries(), G4ICRU73QOModel::SampleSecondaries(), G4LindhardSorensenIonModel::SampleSecondaries(), G4MollerBhabhaModel::SampleSecondaries(), G4PairProductionRelModel::SampleSecondaries(), G4PEEffectFluoModel::SampleSecondaries(), and G4IonParametrisedLossModel::SampleSecondaries().

◆ GetChargeSquareRatio()

G4double G4VEmModel::GetChargeSquareRatio ( const G4ParticleDefinition p,
const G4Material ,
G4double  kineticEnergy 
)
virtualinherited

◆ GetCrossSectionTable()

G4PhysicsTable * G4VEmModel::GetCrossSectionTable ( )
inlineinherited

◆ GetCurrentElement()

const G4Element * G4VEmModel::GetCurrentElement ( ) const
inlineinherited

◆ GetCurrentIsotope()

const G4Isotope * G4VEmModel::GetCurrentIsotope ( ) const
inlineinherited

Definition at line 512 of file G4VEmModel.hh.

513{
514 return fCurrentIsotope;
515}
const G4Isotope * fCurrentIsotope
Definition: G4VEmModel.hh:418

References G4VEmModel::fCurrentIsotope.

Referenced by G4VEmProcess::GetTargetIsotope().

◆ GetElementData()

G4ElementData * G4VEmModel::GetElementData ( )
inlineinherited

◆ GetElementSelectors()

std::vector< G4EmElementSelector * > * G4VEmModel::GetElementSelectors ( )
inlineinherited

◆ GetFSquared()

G4double G4PenelopeRayleighModel::GetFSquared ( const G4Material mat,
const G4double  QSquared 
)
private

Definition at line 665 of file G4PenelopeRayleighModel.cc.

666{
667 G4double f2 = 0;
668 //Input value QSquared could be zero: protect the log() below against
669 //the FPE exception
670 //If Q<1e-10, set Q to 1e-10
671 G4double logQSquared = (QSquared>1e-10) ? G4Log(QSquared) : -23.;
672 //last value of the table
673 G4double maxlogQ2 = fLogQSquareGrid[fLogQSquareGrid.size()-1];
674
675 //now it should be all right
676 G4PhysicsFreeVector* theVec = fLogFormFactorTable->find(mat)->second;
677
678 if (!theVec)
679 {
681 ed << "Unable to retrieve F squared table for " << mat->GetName() << G4endl;
682 G4Exception("G4PenelopeRayleighModel::GetFSquared()",
683 "em2046",FatalException,ed);
684 return 0;
685 }
686 if (logQSquared < -20) // Q < 1e-9
687 {
688 G4double logf2 = (*theVec)[0]; //first value of the table
689 f2 = G4Exp(logf2);
690 }
691 else if (logQSquared > maxlogQ2)
692 f2 =0;
693 else
694 {
695 //log(Q^2) vs. log(F^2)
696 G4double logf2 = theVec->Value(logQSquared);
697 f2 = G4Exp(logf2);
698
699 }
700 if (fVerboseLevel > 3)
701 {
702 G4cout << "G4PenelopeRayleighModel::GetFSquared() in " << mat->GetName() << G4endl;
703 G4cout << "Q^2 = " << QSquared << " (units of 1/(m_e*c); F^2 = " << f2 << G4endl;
704 }
705 return f2;
706}

References FatalException, fLogFormFactorTable, fLogQSquareGrid, fVerboseLevel, G4cout, G4endl, G4Exception(), G4Exp(), G4Log(), G4Material::GetName(), and G4PhysicsVector::Value().

Referenced by InitializeSamplingAlgorithm().

◆ GetModelOfFluctuations()

G4VEmFluctuationModel * G4VEmModel::GetModelOfFluctuations ( )
inlineinherited

◆ GetName()

const G4String & G4VEmModel::GetName ( ) const
inlineinherited

◆ GetPartialCrossSection()

G4double G4VEmModel::GetPartialCrossSection ( const G4Material ,
G4int  level,
const G4ParticleDefinition ,
G4double  kineticEnergy 
)
virtualinherited

◆ GetParticleChangeForGamma()

G4ParticleChangeForGamma * G4VEmModel::GetParticleChangeForGamma ( )
protectedinherited

Definition at line 123 of file G4VEmModel.cc.

124{
125 G4ParticleChangeForGamma* p = nullptr;
126 if (pParticleChange != nullptr) {
127 p = static_cast<G4ParticleChangeForGamma*>(pParticleChange);
128 } else {
129 p = new G4ParticleChangeForGamma();
130 pParticleChange = p;
131 }
132 if(fTripletModel != nullptr) { fTripletModel->SetParticleChange(p); }
133 return p;
134}
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
Definition: G4VEmModel.cc:447
G4VParticleChange * pParticleChange
Definition: G4VEmModel.hh:425
G4VEmModel * fTripletModel
Definition: G4VEmModel.hh:415

References G4VEmModel::fTripletModel, G4VEmModel::pParticleChange, and G4VEmModel::SetParticleChange().

Referenced by G4MicroElecLOPhononModel::G4MicroElecLOPhononModel(), G4DNAChampionElasticModel::Initialise(), G4DNACPA100ElasticModel::Initialise(), G4DNADingfelderChargeDecreaseModel::Initialise(), G4DNADingfelderChargeIncreaseModel::Initialise(), G4DNAMeltonAttachmentModel::Initialise(), G4DNAMillerGreenExcitationModel::Initialise(), G4DNARuddIonisationExtendedModel::Initialise(), G4DNARuddIonisationModel::Initialise(), G4DNASancheExcitationModel::Initialise(), G4DNAScreenedRutherfordElasticModel::Initialise(), G4DNATransformElectronModel::Initialise(), G4DNAUeharaScreenedRutherfordElasticModel::Initialise(), G4LEPTSAttachmentModel::Initialise(), G4LEPTSDissociationModel::Initialise(), G4LEPTSElasticModel::Initialise(), G4LEPTSExcitationModel::Initialise(), G4LEPTSIonisationModel::Initialise(), G4LEPTSPositroniumModel::Initialise(), G4LEPTSRotExcitationModel::Initialise(), G4LEPTSVibExcitationModel::Initialise(), G4BoldyshevTripletModel::Initialise(), G4eplusTo3GammaOKVIModel::Initialise(), G4eSingleCoulombScatteringModel::Initialise(), G4IonCoulombScatteringModel::Initialise(), G4eeToHadronsMultiModel::Initialise(), G4JAEAElasticScatteringModel::Initialise(), G4JAEAPolarizedElasticScatteringModel::Initialise(), G4LivermoreComptonModel::Initialise(), G4LivermoreGammaConversion5DModel::Initialise(), G4LivermoreGammaConversionModel::Initialise(), G4LivermoreNuclearGammaConversionModel::Initialise(), G4LivermorePhotoElectricModel::Initialise(), G4LivermorePolarizedComptonModel::Initialise(), G4LivermorePolarizedGammaConversionModel::Initialise(), G4LivermorePolarizedRayleighModel::Initialise(), G4LivermoreRayleighModel::Initialise(), G4LowEPComptonModel::Initialise(), G4LowEPPolarizedComptonModel::Initialise(), G4MicroElecElasticModel::Initialise(), G4MicroElecElasticModel_new::Initialise(), G4MicroElecInelasticModel::Initialise(), G4MicroElecInelasticModel_new::Initialise(), G4MicroElecLOPhononModel::Initialise(), G4PenelopeAnnihilationModel::Initialise(), G4PenelopeComptonModel::Initialise(), G4PenelopeGammaConversionModel::Initialise(), G4PenelopePhotoElectricModel::Initialise(), Initialise(), G4PenelopeRayleighModelMI::Initialise(), G4PolarizedAnnihilationModel::Initialise(), G4BetheHeitlerModel::Initialise(), G4eCoulombScatteringModel::Initialise(), G4eDPWACoulombScatteringModel::Initialise(), G4eeToTwoGammaModel::Initialise(), G4eplusTo2GammaOKVIModel::Initialise(), G4hCoulombScatteringModel::Initialise(), G4KleinNishinaCompton::Initialise(), G4KleinNishinaModel::Initialise(), G4PairProductionRelModel::Initialise(), G4PEEffectFluoModel::Initialise(), G4XrayRayleighModel::Initialise(), G4DNABornExcitationModel1::Initialise(), G4DNABornExcitationModel2::Initialise(), G4DNABornIonisationModel1::Initialise(), G4DNABornIonisationModel2::Initialise(), G4DNACPA100ExcitationModel::Initialise(), G4DNACPA100IonisationModel::Initialise(), G4DNADiracRMatrixExcitationModel::Initialise(), G4DNAEmfietzoglouExcitationModel::Initialise(), G4DNAEmfietzoglouIonisationModel::Initialise(), G4DNAQuinnPlasmonExcitationModel::Initialise(), G4DNARelativisticIonisationModel::Initialise(), G4DNAELSEPAElasticModel::Initialise(), G4DNAModelInterface::Initialise(), and G4DNAIonElasticModel::Initialise().

◆ GetParticleChangeForLoss()

G4ParticleChangeForLoss * G4VEmModel::GetParticleChangeForLoss ( )
protectedinherited

◆ GetParticleCharge()

G4double G4VEmModel::GetParticleCharge ( const G4ParticleDefinition p,
const G4Material ,
G4double  kineticEnergy 
)
virtualinherited

◆ GetPMaxTable()

void G4PenelopeRayleighModel::GetPMaxTable ( const G4Material mat)
private

Definition at line 1131 of file G4PenelopeRayleighModel.cc.

1132{
1133 if (!fPMaxTable)
1134 {
1135 G4cout << "G4PenelopeRayleighModel::BuildPMaxTable" << G4endl;
1136 G4cout << "Going to instanziate the fPMaxTable !" << G4endl;
1137 G4cout << "That should _not_ be here! " << G4endl;
1138 fPMaxTable = new std::map<const G4Material*,G4PhysicsFreeVector*>;
1139 }
1140 //check if the table is already there
1141 if (fPMaxTable->count(mat))
1142 return;
1143
1144 //otherwise build it
1145 if (!fSamplingTable)
1146 {
1147 G4Exception("G4PenelopeRayleighModel::GetPMaxTable()",
1148 "em2052",FatalException,
1149 "SamplingTable is not properly instantiated");
1150 return;
1151 }
1152
1153 //This should not be: the sampling table is built before the p-table
1154 if (!fSamplingTable->count(mat))
1155 {
1157 ed << "Sampling table for material " << mat->GetName() << " not found";
1158 G4Exception("G4PenelopeRayleighModel::GetPMaxTable()",
1159 "em2052",FatalException,
1160 ed);
1161 return;
1162 }
1163
1164 G4PenelopeSamplingData *theTable = fSamplingTable->find(mat)->second;
1165 size_t tablePoints = theTable->GetNumberOfStoredPoints();
1166
1167 size_t nOfEnergyPoints = fLogEnergyGridPMax.size();
1168 G4PhysicsFreeVector* theVec = new G4PhysicsFreeVector(nOfEnergyPoints);
1169
1170 const size_t nip = 51; //hard-coded in Penelope
1171
1172 for (size_t ie=0;ie<fLogEnergyGridPMax.size();ie++)
1173 {
1175 G4double Qm = 2.0*energy/electron_mass_c2; //this is non-dimensional now
1176 G4double Qm2 = Qm*Qm;
1177 G4double firstQ2 = theTable->GetX(0);
1178 G4double lastQ2 = theTable->GetX(tablePoints-1);
1179 G4double thePMax = 0;
1180
1181 if (Qm2 > firstQ2)
1182 {
1183 if (Qm2 < lastQ2)
1184 {
1185 //bisection to look for the index of Qm
1186 size_t lowerBound = 0;
1187 size_t upperBound = tablePoints-1;
1188 while (lowerBound <= upperBound)
1189 {
1190 size_t midBin = (lowerBound + upperBound)/2;
1191 if( Qm2 < theTable->GetX(midBin))
1192 { upperBound = midBin-1; }
1193 else
1194 { lowerBound = midBin+1; }
1195 }
1196 //upperBound is the output (but also lowerBounf --> should be the same!)
1197 G4double Q1 = theTable->GetX(upperBound);
1198 G4double Q2 = Qm2;
1199 G4double DQ = (Q2-Q1)/((G4double)(nip-1));
1200 G4double theA = theTable->GetA(upperBound);
1201 G4double theB = theTable->GetB(upperBound);
1202 G4double thePAC = theTable->GetPAC(upperBound);
1203 G4DataVector* fun = new G4DataVector();
1204 for (size_t k=0;k<nip;k++)
1205 {
1206 G4double qi = Q1 + k*DQ;
1207 G4double tau = (qi-Q1)/
1208 (theTable->GetX(upperBound+1)-Q1);
1209 G4double con1 = 2.0*theB*tau;
1210 G4double ci = 1.0+theA+theB;
1211 G4double con2 = ci-theA*tau;
1212 G4double etap = 0;
1213 if (std::fabs(con1) > 1.0e-16*std::fabs(con2))
1214 etap = con2*(1.0-std::sqrt(1.0-2.0*tau*con1/(con2*con2)))/con1;
1215 else
1216 etap = tau/con2;
1217 G4double theFun = (theTable->GetPAC(upperBound+1)-thePAC)*
1218 (1.0+(theA+theB*etap)*etap)*(1.0+(theA+theB*etap)*etap)/
1219 ((1.0-theB*etap*etap)*ci*(theTable->GetX(upperBound+1)-Q1));
1220 fun->push_back(theFun);
1221 }
1222 //Now intergrate numerically the fun Cavalieri-Simpson's method
1223 G4DataVector* sum = new G4DataVector;
1224 G4double CONS = DQ*(1./12.);
1225 G4double HCONS = 0.5*CONS;
1226 sum->push_back(0.);
1227 G4double secondPoint = (*sum)[0] +
1228 (5.0*(*fun)[0]+8.0*(*fun)[1]-(*fun)[2])*CONS;
1229 sum->push_back(secondPoint);
1230 for (size_t hh=2;hh<nip-1;hh++)
1231 {
1232 G4double previous = (*sum)[hh-1];
1233 G4double next = previous+(13.0*((*fun)[hh-1]+(*fun)[hh])-
1234 (*fun)[hh+1]-(*fun)[hh-2])*HCONS;
1235 sum->push_back(next);
1236 }
1237 G4double last = (*sum)[nip-2]+(5.0*(*fun)[nip-1]+8.0*(*fun)[nip-2]-
1238 (*fun)[nip-3])*CONS;
1239 sum->push_back(last);
1240 thePMax = thePAC + (*sum)[sum->size()-1]; //last point
1241 delete fun;
1242 delete sum;
1243 }
1244 else
1245 {
1246 thePMax = 1.0;
1247 }
1248 }
1249 else
1250 {
1251 thePMax = theTable->GetPAC(0);
1252 }
1253
1254 //Write number in the table
1255 theVec->PutValue(ie,energy,thePMax);
1256 }
1257
1258 fPMaxTable->insert(std::make_pair(mat,theVec));
1259 return;
1260}
G4double GetA(size_t index)
G4double GetPAC(size_t index)
G4double GetX(size_t index)
G4double GetB(size_t index)
float electron_mass_c2
Definition: hepunit.py:273

References source.hepunit::electron_mass_c2, G4INCL::KinematicsUtils::energy(), FatalException, fLogEnergyGridPMax, fPMaxTable, fSamplingTable, G4cout, G4endl, G4Exception(), G4Exp(), G4PenelopeSamplingData::GetA(), G4PenelopeSamplingData::GetB(), G4Material::GetName(), G4PenelopeSamplingData::GetNumberOfStoredPoints(), G4PenelopeSamplingData::GetPAC(), G4PenelopeSamplingData::GetX(), and G4PhysicsFreeVector::PutValue().

Referenced by Initialise(), and SampleSecondaries().

◆ GetTripletModel()

G4VEmModel * G4VEmModel::GetTripletModel ( )
inlineinherited

◆ GetVerbosityLevel()

G4int G4PenelopeRayleighModel::GetVerbosityLevel ( )
inline

Definition at line 85 of file G4PenelopeRayleighModel.hh.

85{return fVerboseLevel;};

References fVerboseLevel.

◆ HighEnergyActivationLimit()

G4double G4VEmModel::HighEnergyActivationLimit ( ) const
inlineinherited

◆ HighEnergyLimit()

G4double G4VEmModel::HighEnergyLimit ( ) const
inlineinherited

Definition at line 655 of file G4VEmModel.hh.

656{
657 return highLimit;
658}
G4double highLimit
Definition: G4VEmModel.hh:437

References G4VEmModel::highLimit.

Referenced by G4DNAChampionElasticModel::CrossSectionPerVolume(), G4DNACPA100ElasticModel::CrossSectionPerVolume(), G4DNACPA100ExcitationModel::CrossSectionPerVolume(), G4DNACPA100IonisationModel::CrossSectionPerVolume(), G4DNAEmfietzoglouExcitationModel::CrossSectionPerVolume(), G4DNAEmfietzoglouIonisationModel::CrossSectionPerVolume(), G4DNAMeltonAttachmentModel::CrossSectionPerVolume(), G4DNASancheExcitationModel::CrossSectionPerVolume(), G4DNAScreenedRutherfordElasticModel::CrossSectionPerVolume(), G4DNATransformElectronModel::CrossSectionPerVolume(), G4DNAELSEPAElasticModel::CrossSectionPerVolume(), G4DNAChampionElasticModel::G4DNAChampionElasticModel(), G4DNACPA100ElasticModel::G4DNACPA100ElasticModel(), G4DNAMeltonAttachmentModel::G4DNAMeltonAttachmentModel(), G4DNASancheExcitationModel::G4DNASancheExcitationModel(), G4DNAUeharaScreenedRutherfordElasticModel::G4DNAUeharaScreenedRutherfordElasticModel(), G4eeToHadronsModel::G4eeToHadronsModel(), G4IonParametrisedLossModel::G4IonParametrisedLossModel(), G4LivermorePolarizedRayleighModel::G4LivermorePolarizedRayleighModel(), G4ContinuousGainOfEnergy::GetContinuousStepLimit(), G4PenelopeBremsstrahlungModel::GetCrossSectionTableForCouple(), G4VMscModel::GetParticleChangeForMSC(), G4DNAChampionElasticModel::Initialise(), G4DNACPA100ElasticModel::Initialise(), G4DNADingfelderChargeDecreaseModel::Initialise(), G4DNADingfelderChargeIncreaseModel::Initialise(), G4DNAMeltonAttachmentModel::Initialise(), G4DNAMillerGreenExcitationModel::Initialise(), G4DNARuddIonisationExtendedModel::Initialise(), G4DNARuddIonisationModel::Initialise(), G4DNASancheExcitationModel::Initialise(), G4DNAScreenedRutherfordElasticModel::Initialise(), G4DNAUeharaScreenedRutherfordElasticModel::Initialise(), G4BoldyshevTripletModel::Initialise(), G4PAIModel::Initialise(), G4PAIPhotModel::Initialise(), G4JAEAElasticScatteringModel::Initialise(), G4JAEAPolarizedElasticScatteringModel::Initialise(), G4LivermoreComptonModel::Initialise(), G4LivermoreGammaConversion5DModel::Initialise(), G4LivermoreGammaConversionModel::Initialise(), G4LivermoreIonisationModel::Initialise(), G4LivermoreNuclearGammaConversionModel::Initialise(), G4LivermorePolarizedComptonModel::Initialise(), G4LivermorePolarizedGammaConversionModel::Initialise(), G4LivermoreRayleighModel::Initialise(), G4LowEPComptonModel::Initialise(), G4LowEPPolarizedComptonModel::Initialise(), G4MicroElecElasticModel::Initialise(), G4MicroElecElasticModel_new::Initialise(), G4MicroElecInelasticModel::Initialise(), G4MicroElecInelasticModel_new::Initialise(), G4PenelopeAnnihilationModel::Initialise(), G4PenelopeBremsstrahlungModel::Initialise(), G4PenelopeComptonModel::Initialise(), G4PenelopeGammaConversionModel::Initialise(), G4PenelopeIonisationModel::Initialise(), G4PenelopePhotoElectricModel::Initialise(), Initialise(), G4PenelopeRayleighModelMI::Initialise(), G4MuBremsstrahlungModel::Initialise(), G4MuPairProductionModel::Initialise(), G4eBremsstrahlungRelModel::Initialise(), G4eDPWACoulombScatteringModel::Initialise(), G4GoudsmitSaundersonMscModel::Initialise(), G4PairProductionRelModel::Initialise(), G4SeltzerBergerModel::Initialise(), G4WentzelVIModel::Initialise(), G4DNABornExcitationModel1::Initialise(), G4DNABornExcitationModel2::Initialise(), G4DNABornIonisationModel1::Initialise(), G4DNABornIonisationModel2::Initialise(), G4DNACPA100ExcitationModel::Initialise(), G4DNACPA100IonisationModel::Initialise(), G4DNADiracRMatrixExcitationModel::Initialise(), G4DNAEmfietzoglouExcitationModel::Initialise(), G4DNAEmfietzoglouIonisationModel::Initialise(), G4DNAQuinnPlasmonExcitationModel::Initialise(), G4DNARelativisticIonisationModel::Initialise(), G4EmModelManager::Initialise(), G4DNAELSEPAElasticModel::Initialise(), G4DNADummyModel::Initialise(), G4DNAIonElasticModel::Initialise(), G4mplIonisation::InitialiseEnergyLossProcess(), G4ionIonisation::InitialiseEnergyLossProcess(), G4MuBremsstrahlungModel::InitialiseLocal(), G4MuPairProductionModel::InitialiseLocal(), G4eBremsstrahlungRelModel::InitialiseLocal(), G4PairProductionRelModel::InitialiseLocal(), G4CoulombScattering::InitialiseProcess(), G4VEmProcess::PostStepDoIt(), G4VEmProcess::PreparePhysicsTable(), G4VEnergyLossProcess::PreparePhysicsTable(), G4VMultipleScattering::PreparePhysicsTable(), G4DNACPA100IonisationModel::SampleSecondaries(), G4DNAEmfietzoglouIonisationModel::SampleSecondaries(), G4DNASancheExcitationModel::SampleSecondaries(), G4EmConfigurator::SetExtraEmModel(), G4mplIonisationModel::SetParticle(), G4mplIonisationWithDeltaModel::SetParticle(), G4eBremsstrahlung::StreamProcessInfo(), and G4EmConfigurator::UpdateModelEnergyRange().

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 152 of file G4PenelopeRayleighModel.cc.

154{
155 if (fVerboseLevel > 3)
156 G4cout << "Calling G4PenelopeRayleighModel::Initialise()" << G4endl;
157
158 SetParticle(part);
159
160 //Only the master model creates/fills/destroys the tables
161 if (IsMaster() && part == fParticle)
162 {
163 //clear tables depending on materials, not the atomic ones
164 ClearTables();
165
166 if (fVerboseLevel > 3)
167 G4cout << "Calling G4PenelopeRayleighModel::Initialise() [master]" << G4endl;
168
169 //create new tables
171 fLogFormFactorTable = new std::map<const G4Material*,G4PhysicsFreeVector*>;
172 if (!fPMaxTable)
173 fPMaxTable = new std::map<const G4Material*,G4PhysicsFreeVector*>;
174 if (!fSamplingTable)
175 fSamplingTable = new std::map<const G4Material*,G4PenelopeSamplingData*>;
176
177 G4ProductionCutsTable* theCoupleTable =
179
180 for (size_t i=0;i<theCoupleTable->GetTableSize();i++)
181 {
182 const G4Material* material =
183 theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
184 const G4ElementVector* theElementVector = material->GetElementVector();
185
186 for (size_t j=0;j<material->GetNumberOfElements();j++)
187 {
188 G4int iZ = theElementVector->at(j)->GetZasInt();
189 //read data files only in the master
190 if (!fLogAtomicCrossSection[iZ])
191 ReadDataFile(iZ);
192 }
193
194 //1) If the table has not been built for the material, do it!
195 if (!fLogFormFactorTable->count(material))
197
198 //2) retrieve or build the sampling table
199 if (!(fSamplingTable->count(material)))
201
202 //3) retrieve or build the pMax data
203 if (!fPMaxTable->count(material))
205 }
206
207 if (fVerboseLevel > 1) {
208 G4cout << "Penelope Rayleigh model v2008 is initialized " << G4endl
209 << "Energy range: "
210 << LowEnergyLimit() / keV << " keV - "
211 << HighEnergyLimit() / GeV << " GeV"
212 << G4endl;
213 }
214 }
215
216 if(fIsInitialised) return;
218 fIsInitialised = true;
219}
const G4Material * GetMaterial() const
void GetPMaxTable(const G4Material *)
void InitializeSamplingAlgorithm(const G4Material *)
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:123
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:662
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:655

References BuildFormFactorTable(), ClearTables(), fIsInitialised, fLogAtomicCrossSection, fLogFormFactorTable, fParticle, fParticleChange, fPMaxTable, fSamplingTable, fVerboseLevel, G4cout, G4endl, G4MaterialCutsCouple::GetMaterial(), G4ProductionCutsTable::GetMaterialCutsCouple(), G4VEmModel::GetParticleChangeForGamma(), GetPMaxTable(), G4ProductionCutsTable::GetProductionCutsTable(), G4ProductionCutsTable::GetTableSize(), GeV, G4VEmModel::HighEnergyLimit(), InitializeSamplingAlgorithm(), G4VEmModel::IsMaster(), keV, G4VEmModel::LowEnergyLimit(), eplot::material, ReadDataFile(), and SetParticle().

◆ InitialiseElementSelectors()

void G4VEmModel::InitialiseElementSelectors ( const G4ParticleDefinition part,
const G4DataVector cuts 
)
inherited

Definition at line 138 of file G4VEmModel.cc.

140{
141 // using spline for element selectors should be investigated in details
142 // because small number of points may provide biased results
143 // large number of points requires significant increase of memory
144 G4bool spline = false;
145
146 //G4cout << "IES: for " << GetName() << " Emin(MeV)= " << lowLimit/MeV
147 // << " Emax(MeV)= " << highLimit/MeV << G4endl;
148
149 // two times less bins because probability functon is normalized
150 // so correspondingly is more smooth
151 if(highLimit <= lowLimit) { return; }
152
154
155 G4ProductionCutsTable* theCoupleTable=
157 G4int numOfCouples = theCoupleTable->GetTableSize();
158
159 // prepare vector
160 if(!elmSelectors) {
161 elmSelectors = new std::vector<G4EmElementSelector*>;
162 }
163 if(numOfCouples > nSelectors) {
164 for(G4int i=nSelectors; i<numOfCouples; ++i) {
165 elmSelectors->push_back(nullptr);
166 }
167 nSelectors = numOfCouples;
168 }
169
170 // initialise vector
171 for(G4int i=0; i<numOfCouples; ++i) {
172
173 // no need in element selectors for infinite cuts
174 if(cuts[i] == DBL_MAX) { continue; }
175
176 auto couple = theCoupleTable->GetMaterialCutsCouple(i);
177 auto material = couple->GetMaterial();
178 SetCurrentCouple(couple);
179
180 // selector already exist then delete
181 delete (*elmSelectors)[i];
182
183 G4double emin = std::max(lowLimit, MinPrimaryEnergy(material, part, cuts[i]));
184 G4double emax = std::max(highLimit, 10*emin);
185 static const G4double invlog106 = 1.0/(6*G4Log(10.));
186 G4int nbins = (G4int)(nbinsPerDec*G4Log(emax/emin)*invlog106);
187 nbins = std::max(nbins, 3);
188
189 (*elmSelectors)[i] = new G4EmElementSelector(this,material,nbins,
190 emin,emax,spline);
191 ((*elmSelectors)[i])->Initialise(part, cuts[i]);
192 /*
193 G4cout << "G4VEmModel::InitialiseElmSelectors i= " << i
194 << " " << part->GetParticleName()
195 << " for " << GetName() << " cut= " << cuts[i]
196 << " " << (*elmSelectors)[i] << G4endl;
197 ((*elmSelectors)[i])->Dump(part);
198 */
199 }
200}
bool G4bool
Definition: G4Types.hh:86
static G4EmParameters * Instance()
G4int NumberOfBinsPerDecade() const
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
Definition: G4VEmModel.cc:415
G4double lowLimit
Definition: G4VEmModel.hh:436
G4int nSelectors
Definition: G4VEmModel.hh:443
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)=0
T max(const T t1, const T t2)
brief Return the largest of the two arguments

References DBL_MAX, G4VEmModel::elmSelectors, emax, G4Log(), G4ProductionCutsTable::GetMaterialCutsCouple(), G4ProductionCutsTable::GetProductionCutsTable(), G4ProductionCutsTable::GetTableSize(), G4VEmModel::highLimit, G4VEmModel::Initialise(), G4EmParameters::Instance(), G4VEmModel::lowLimit, eplot::material, G4INCL::Math::max(), G4VEmModel::MinPrimaryEnergy(), G4VEmModel::nSelectors, G4EmParameters::NumberOfBinsPerDecade(), and G4VEmModel::SetCurrentCouple().

Referenced by G4eSingleCoulombScatteringModel::Initialise(), G4PAIModel::Initialise(), G4PAIPhotModel::Initialise(), G4JAEAElasticScatteringModel::Initialise(), G4JAEAPolarizedElasticScatteringModel::Initialise(), G4LivermoreComptonModel::Initialise(), G4LivermoreGammaConversion5DModel::Initialise(), G4LivermoreGammaConversionModel::Initialise(), G4LivermoreNuclearGammaConversionModel::Initialise(), G4LivermorePolarizedComptonModel::Initialise(), G4LivermorePolarizedGammaConversionModel::Initialise(), G4LivermorePolarizedRayleighModel::Initialise(), G4LivermoreRayleighModel::Initialise(), G4LowEPComptonModel::Initialise(), G4LowEPPolarizedComptonModel::Initialise(), G4PenelopePhotoElectricModel::Initialise(), G4MuBremsstrahlungModel::Initialise(), G4MuPairProductionModel::Initialise(), G4BetheHeitlerModel::Initialise(), G4eBremParametrizedModel::Initialise(), G4eBremsstrahlungRelModel::Initialise(), G4eCoulombScatteringModel::Initialise(), G4eDPWACoulombScatteringModel::Initialise(), G4hCoulombScatteringModel::Initialise(), G4KleinNishinaCompton::Initialise(), G4KleinNishinaModel::Initialise(), G4PairProductionRelModel::Initialise(), G4SeltzerBergerModel::Initialise(), and G4XrayRayleighModel::Initialise().

◆ InitialiseForElement()

void G4VEmModel::InitialiseForElement ( const G4ParticleDefinition ,
G4int  Z 
)
virtualinherited

◆ InitialiseForMaterial()

void G4VEmModel::InitialiseForMaterial ( const G4ParticleDefinition part,
const G4Material material 
)
virtualinherited

Definition at line 209 of file G4VEmModel.cc.

211{
212 if(material != nullptr) {
213 size_t n = material->GetNumberOfElements();
214 for(size_t i=0; i<n; ++i) {
215 G4int Z = material->GetElement(i)->GetZasInt();
216 InitialiseForElement(part, Z);
217 }
218 }
219}
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
Definition: G4VEmModel.cc:223

References G4VEmModel::InitialiseForElement(), eplot::material, CLHEP::detail::n, and Z.

Referenced by G4EmCalculator::FindEmModel().

◆ InitialiseLocal()

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

Reimplemented from G4VEmModel.

Definition at line 223 of file G4PenelopeRayleighModel.cc.

225{
226 if (fVerboseLevel > 3)
227 G4cout << "Calling G4PenelopeRayleighModel::InitialiseLocal()" << G4endl;
228 //
229 //Check that particle matches: one might have multiple master models (e.g.
230 //for e+ and e-).
231 //
232 if (part == fParticle)
233 {
234 //Get the const table pointers from the master to the workers
235 const G4PenelopeRayleighModel* theModel =
236 static_cast<G4PenelopeRayleighModel*> (masterModel);
237
238 //Copy pointers to the data tables
239 for(G4int i=0; i<=fMaxZ; ++i)
240 {
242 fAtomicFormFactor[i] = theModel->fAtomicFormFactor[i];
243 }
245 fPMaxTable = theModel->fPMaxTable;
246 fSamplingTable = theModel->fSamplingTable;
247
248 //copy the G4DataVector with the grid
250
251 //Same verbosity for all workers, as the master
252 fVerboseLevel = theModel->fVerboseLevel;
253 }
254
255 return;
256}

References fAtomicFormFactor, fLogAtomicCrossSection, fLogFormFactorTable, fLogQSquareGrid, fMaxZ, fParticle, fPMaxTable, fSamplingTable, fVerboseLevel, G4cout, and G4endl.

◆ InitializeSamplingAlgorithm()

void G4PenelopeRayleighModel::InitializeSamplingAlgorithm ( const G4Material mat)
private

Definition at line 710 of file G4PenelopeRayleighModel.cc.

711{
712 G4double q2min = 0;
713 G4double q2max = 0;
714 const size_t np = 150; //hard-coded in Penelope
715 //G4cout << "Init N= " << fLogQSquareGrid.size() << G4endl;
716 for (size_t i=1;i<fLogQSquareGrid.size();i++)
717 {
719 if (GetFSquared(mat,Q2) > 1e-35)
720 {
721 q2max = G4Exp(fLogQSquareGrid[i-1]);
722 }
723 //G4cout << "Q2= " << Q2 << " q2max= " << q2max << G4endl;
724 }
725
726 size_t nReducedPoints = np/4;
727
728 //check for errors
729 if (np < 16)
730 {
731 G4Exception("G4PenelopeRayleighModel::InitializeSamplingAlgorithm()",
732 "em2047",FatalException,
733 "Too few points to initialize the sampling algorithm");
734 }
735 if (q2min > (q2max-1e-10))
736 {
737 G4cout << "q2min= " << q2min << " q2max= " << q2max << G4endl;
738 G4Exception("G4PenelopeRayleighModel::InitializeSamplingAlgorithm()",
739 "em2048",FatalException,
740 "Too narrow grid to initialize the sampling algorithm");
741 }
742
743 //This is subroutine RITAI0 of Penelope
744 //Create an object of type G4PenelopeRayleighSamplingData --> store in a map::Material*
745
746 //temporary vectors --> Then everything is stored in G4PenelopeSamplingData
747 G4DataVector* x = new G4DataVector();
748
749 /*******************************************************************************
750 Start with a grid of NUNIF points uniformly spaced in the interval q2min,q2max
751 ********************************************************************************/
752 size_t NUNIF = std::min(std::max(((size_t)8),nReducedPoints),np/2);
753 const G4int nip = 51; //hard-coded in Penelope
754
755 G4double dx = (q2max-q2min)/((G4double) NUNIF-1);
756 x->push_back(q2min);
757 for (size_t i=1;i<NUNIF-1;i++)
758 {
759 G4double app = q2min + i*dx;
760 x->push_back(app); //increase
761 }
762 x->push_back(q2max);
763
764 if (fVerboseLevel> 3)
765 G4cout << "Vector x has " << x->size() << " points, while NUNIF = " << NUNIF << G4endl;
766
767 //Allocate temporary storage vectors
768 G4DataVector* area = new G4DataVector();
769 G4DataVector* a = new G4DataVector();
770 G4DataVector* b = new G4DataVector();
771 G4DataVector* c = new G4DataVector();
772 G4DataVector* err = new G4DataVector();
773
774 for (size_t i=0;i<NUNIF-1;i++) //build all points but the last
775 {
776 //Temporary vectors for this loop
777 G4DataVector* pdfi = new G4DataVector();
778 G4DataVector* pdfih = new G4DataVector();
779 G4DataVector* sumi = new G4DataVector();
780
781 G4double dxi = ((*x)[i+1]-(*x)[i])/(G4double (nip-1));
782 G4double pdfmax = 0;
783 for (G4int k=0;k<nip;k++)
784 {
785 G4double xik = (*x)[i]+k*dxi;
786 G4double pdfk = std::max(GetFSquared(mat,xik),0.);
787 pdfi->push_back(pdfk);
788 pdfmax = std::max(pdfmax,pdfk);
789 if (k < (nip-1))
790 {
791 G4double xih = xik + 0.5*dxi;
792 G4double pdfIK = std::max(GetFSquared(mat,xih),0.);
793 pdfih->push_back(pdfIK);
794 pdfmax = std::max(pdfmax,pdfIK);
795 }
796 }
797
798 //Simpson's integration
799 G4double cons = dxi*0.5*(1./3.);
800 sumi->push_back(0.);
801 for (G4int k=1;k<nip;k++)
802 {
803 G4double previous = (*sumi)[k-1];
804 G4double next = previous + cons*((*pdfi)[k-1]+4.0*(*pdfih)[k-1]+(*pdfi)[k]);
805 sumi->push_back(next);
806 }
807
808 G4double lastIntegral = (*sumi)[sumi->size()-1];
809 area->push_back(lastIntegral);
810 //Normalize cumulative function
811 G4double factor = 1.0/lastIntegral;
812 for (size_t k=0;k<sumi->size();k++)
813 (*sumi)[k] *= factor;
814
815 //When the PDF vanishes at one of the interval end points, its value is modified
816 if ((*pdfi)[0] < 1e-35)
817 (*pdfi)[0] = 1e-5*pdfmax;
818 if ((*pdfi)[pdfi->size()-1] < 1e-35)
819 (*pdfi)[pdfi->size()-1] = 1e-5*pdfmax;
820
821 G4double pli = (*pdfi)[0]*factor;
822 G4double pui = (*pdfi)[pdfi->size()-1]*factor;
823 G4double B_temp = 1.0-1.0/(pli*pui*dx*dx);
824 G4double A_temp = (1.0/(pli*dx))-1.0-B_temp;
825 G4double C_temp = 1.0+A_temp+B_temp;
826 if (C_temp < 1e-35)
827 {
828 a->push_back(0.);
829 b->push_back(0.);
830 c->push_back(1.);
831 }
832 else
833 {
834 a->push_back(A_temp);
835 b->push_back(B_temp);
836 c->push_back(C_temp);
837 }
838
839 //OK, now get ERR(I), the integral of the absolute difference between the rational interpolation
840 //and the true pdf, extended over the interval (X(I),X(I+1))
841 G4int icase = 1; //loop code
842 G4bool reLoop = false;
843 err->push_back(0.);
844 do
845 {
846 reLoop = false;
847 (*err)[i] = 0.; //zero variable
848 for (G4int k=0;k<nip;k++)
849 {
850 G4double rr = (*sumi)[k];
851 G4double pap = (*area)[i]*(1.0+((*a)[i]+(*b)[i]*rr)*rr)*(1.0+((*a)[i]+(*b)[i]*rr)*rr)/
852 ((1.0-(*b)[i]*rr*rr)*(*c)[i]*((*x)[i+1]-(*x)[i]));
853 if (k == 0 || k == nip-1)
854 (*err)[i] += 0.5*std::fabs(pap-(*pdfi)[k]);
855 else
856 (*err)[i] += std::fabs(pap-(*pdfi)[k]);
857 }
858 (*err)[i] *= dxi;
859
860 //If err(I) is too large, the pdf is approximated by a uniform distribution
861 if ((*err)[i] > 0.1*(*area)[i] && icase == 1)
862 {
863 (*b)[i] = 0;
864 (*a)[i] = 0;
865 (*c)[i] = 1.;
866 icase = 2;
867 reLoop = true;
868 }
869 }while(reLoop);
870 delete pdfi;
871 delete pdfih;
872 delete sumi;
873 } //end of first loop over i
874
875 //Now assign last point
876 (*x)[x->size()-1] = q2max;
877 a->push_back(0.);
878 b->push_back(0.);
879 c->push_back(0.);
880 err->push_back(0.);
881 area->push_back(0.);
882
883 if (x->size() != NUNIF || a->size() != NUNIF ||
884 err->size() != NUNIF || area->size() != NUNIF)
885 {
887 ed << "Problem in building the Table for Sampling: array dimensions do not match" << G4endl;
888 G4Exception("G4PenelopeRayleighModel::InitializeSamplingAlgorithm()",
889 "em2049",FatalException,ed);
890 }
891
892 /*******************************************************************************
893 New grid points are added by halving the sub-intervals with the largest absolute error
894 This is done up to np=150 points in the grid
895 ********************************************************************************/
896 do
897 {
898 G4double maxError = 0.0;
899 size_t iErrMax = 0;
900 for (size_t i=0;i<err->size()-2;i++)
901 {
902 //maxError is the lagest of the interval errors err[i]
903 if ((*err)[i] > maxError)
904 {
905 maxError = (*err)[i];
906 iErrMax = i;
907 }
908 }
909
910 //OK, now I have to insert one new point in the position iErrMax
911 G4double newx = 0.5*((*x)[iErrMax]+(*x)[iErrMax+1]);
912
913 x->insert(x->begin()+iErrMax+1,newx);
914 //Add place-holders in the other vectors
915 area->insert(area->begin()+iErrMax+1,0.);
916 a->insert(a->begin()+iErrMax+1,0.);
917 b->insert(b->begin()+iErrMax+1,0.);
918 c->insert(c->begin()+iErrMax+1,0.);
919 err->insert(err->begin()+iErrMax+1,0.);
920
921 //Now calculate the other parameters
922 for (size_t i=iErrMax;i<=iErrMax+1;i++)
923 {
924 //Temporary vectors for this loop
925 G4DataVector* pdfi = new G4DataVector();
926 G4DataVector* pdfih = new G4DataVector();
927 G4DataVector* sumi = new G4DataVector();
928
929 G4double dxLocal = (*x)[i+1]-(*x)[i];
930 G4double dxi = ((*x)[i+1]-(*x)[i])/(G4double (nip-1));
931 G4double pdfmax = 0;
932 for (G4int k=0;k<nip;k++)
933 {
934 G4double xik = (*x)[i]+k*dxi;
935 G4double pdfk = std::max(GetFSquared(mat,xik),0.);
936 pdfi->push_back(pdfk);
937 pdfmax = std::max(pdfmax,pdfk);
938 if (k < (nip-1))
939 {
940 G4double xih = xik + 0.5*dxi;
941 G4double pdfIK = std::max(GetFSquared(mat,xih),0.);
942 pdfih->push_back(pdfIK);
943 pdfmax = std::max(pdfmax,pdfIK);
944 }
945 }
946
947 //Simpson's integration
948 G4double cons = dxi*0.5*(1./3.);
949 sumi->push_back(0.);
950 for (G4int k=1;k<nip;k++)
951 {
952 G4double previous = (*sumi)[k-1];
953 G4double next = previous + cons*((*pdfi)[k-1]+4.0*(*pdfih)[k-1]+(*pdfi)[k]);
954 sumi->push_back(next);
955 }
956 G4double lastIntegral = (*sumi)[sumi->size()-1];
957 (*area)[i] = lastIntegral;
958
959 //Normalize cumulative function
960 G4double factor = 1.0/lastIntegral;
961 for (size_t k=0;k<sumi->size();k++)
962 (*sumi)[k] *= factor;
963
964 //When the PDF vanishes at one of the interval end points, its value is modified
965 if ((*pdfi)[0] < 1e-35)
966 (*pdfi)[0] = 1e-5*pdfmax;
967 if ((*pdfi)[pdfi->size()-1] < 1e-35)
968 (*pdfi)[pdfi->size()-1] = 1e-5*pdfmax;
969
970 G4double pli = (*pdfi)[0]*factor;
971 G4double pui = (*pdfi)[pdfi->size()-1]*factor;
972 G4double B_temp = 1.0-1.0/(pli*pui*dxLocal*dxLocal);
973 G4double A_temp = (1.0/(pli*dxLocal))-1.0-B_temp;
974 G4double C_temp = 1.0+A_temp+B_temp;
975 if (C_temp < 1e-35)
976 {
977 (*a)[i]= 0.;
978 (*b)[i] = 0.;
979 (*c)[i] = 1;
980 }
981 else
982 {
983 (*a)[i]= A_temp;
984 (*b)[i] = B_temp;
985 (*c)[i] = C_temp;
986 }
987 //OK, now get ERR(I), the integral of the absolute difference between the rational interpolation
988 //and the true pdf, extended over the interval (X(I),X(I+1))
989 G4int icase = 1; //loop code
990 G4bool reLoop = false;
991 do
992 {
993 reLoop = false;
994 (*err)[i] = 0.; //zero variable
995 for (G4int k=0;k<nip;k++)
996 {
997 G4double rr = (*sumi)[k];
998 G4double pap = (*area)[i]*(1.0+((*a)[i]+(*b)[i]*rr)*rr)*(1.0+((*a)[i]+(*b)[i]*rr)*rr)/
999 ((1.0-(*b)[i]*rr*rr)*(*c)[i]*((*x)[i+1]-(*x)[i]));
1000 if (k == 0 || k == nip-1)
1001 (*err)[i] += 0.5*std::fabs(pap-(*pdfi)[k]);
1002 else
1003 (*err)[i] += std::fabs(pap-(*pdfi)[k]);
1004 }
1005 (*err)[i] *= dxi;
1006
1007 //If err(I) is too large, the pdf is approximated by a uniform distribution
1008 if ((*err)[i] > 0.1*(*area)[i] && icase == 1)
1009 {
1010 (*b)[i] = 0;
1011 (*a)[i] = 0;
1012 (*c)[i] = 1.;
1013 icase = 2;
1014 reLoop = true;
1015 }
1016 }while(reLoop);
1017 delete pdfi;
1018 delete pdfih;
1019 delete sumi;
1020 }
1021 }while(x->size() < np);
1022
1023 if (x->size() != np || a->size() != np ||
1024 err->size() != np || area->size() != np)
1025 {
1026 G4Exception("G4PenelopeRayleighModel::InitializeSamplingAlgorithm()",
1027 "em2050",FatalException,
1028 "Problem in building the extended Table for Sampling: array dimensions do not match ");
1029 }
1030
1031 /*******************************************************************************
1032 Renormalization
1033 ********************************************************************************/
1034 G4double ws = 0;
1035 for (size_t i=0;i<np-1;i++)
1036 ws += (*area)[i];
1037 ws = 1.0/ws;
1038 G4double errMax = 0;
1039 for (size_t i=0;i<np-1;i++)
1040 {
1041 (*area)[i] *= ws;
1042 (*err)[i] *= ws;
1043 errMax = std::max(errMax,(*err)[i]);
1044 }
1045
1046 //Vector with the normalized cumulative distribution
1047 G4DataVector* PAC = new G4DataVector();
1048 PAC->push_back(0.);
1049 for (size_t i=0;i<np-1;i++)
1050 {
1051 G4double previous = (*PAC)[i];
1052 PAC->push_back(previous+(*area)[i]);
1053 }
1054 (*PAC)[PAC->size()-1] = 1.;
1055
1056 /*******************************************************************************
1057 Pre-calculated limits for the initial binary search for subsequent sampling
1058 ********************************************************************************/
1059 std::vector<size_t> *ITTL = new std::vector<size_t>;
1060 std::vector<size_t> *ITTU = new std::vector<size_t>;
1061
1062 //Just create place-holders
1063 for (size_t i=0;i<np;i++)
1064 {
1065 ITTL->push_back(0);
1066 ITTU->push_back(0);
1067 }
1068
1069 G4double bin = 1.0/(np-1);
1070 (*ITTL)[0]=0;
1071 for (size_t i=1;i<(np-1);i++)
1072 {
1073 G4double ptst = i*bin;
1074 G4bool found = false;
1075 for (size_t j=(*ITTL)[i-1];j<np && !found;j++)
1076 {
1077 if ((*PAC)[j] > ptst)
1078 {
1079 (*ITTL)[i] = j-1;
1080 (*ITTU)[i-1] = j;
1081 found = true;
1082 }
1083 }
1084 }
1085 (*ITTU)[ITTU->size()-2] = ITTU->size()-1;
1086 (*ITTU)[ITTU->size()-1] = ITTU->size()-1;
1087 (*ITTL)[ITTL->size()-1] = ITTU->size()-2;
1088
1089 if (ITTU->size() != np || ITTU->size() != np)
1090 {
1091 G4Exception("G4PenelopeRayleighModel::InitializeSamplingAlgorithm()",
1092 "em2051",FatalException,
1093 "Problem in building the Limit Tables for Sampling: array dimensions do not match");
1094 }
1095
1096 /********************************************************************************
1097 Copy tables
1098 ********************************************************************************/
1100 for (size_t i=0;i<np;i++)
1101 {
1102 theTable->AddPoint((*x)[i],(*PAC)[i],(*a)[i],(*b)[i],(*ITTL)[i],(*ITTU)[i]);
1103 }
1104
1105 if (fVerboseLevel > 2)
1106 {
1107 G4cout << "*************************************************************************" <<
1108 G4endl;
1109 G4cout << "Sampling table for Penelope Rayleigh scattering in " << mat->GetName() << G4endl;
1110 theTable->DumpTable();
1111 }
1112 fSamplingTable->insert(std::make_pair(mat,theTable));
1113
1114 //Clean up temporary vectors
1115 delete x;
1116 delete a;
1117 delete b;
1118 delete c;
1119 delete err;
1120 delete area;
1121 delete PAC;
1122 delete ITTL;
1123 delete ITTU;
1124
1125 //DONE!
1126 return;
1127}
G4double GetFSquared(const G4Material *, const G4double)
void AddPoint(G4double x0, G4double pac0, G4double a0, G4double b0, size_t ITTL0, size_t ITTU0)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
app
Definition: demo.py:189

References G4PenelopeSamplingData::AddPoint(), demo::app, G4PenelopeSamplingData::DumpTable(), FatalException, fLogQSquareGrid, fSamplingTable, fVerboseLevel, G4cout, G4endl, G4Exception(), G4Exp(), GetFSquared(), G4Material::GetName(), G4INCL::Math::max(), and G4INCL::Math::min().

Referenced by Initialise(), and SampleSecondaries().

◆ IsActive()

G4bool G4VEmModel::IsActive ( G4double  kinEnergy) const
inlineinherited

◆ IsLocked()

G4bool G4VEmModel::IsLocked ( ) const
inlineinherited

◆ IsMaster()

G4bool G4VEmModel::IsMaster ( ) const
inlineinherited

Definition at line 746 of file G4VEmModel.hh.

747{
748 return isMaster;
749}
G4bool isMaster
Definition: G4VEmModel.hh:457

References G4VEmModel::isMaster.

Referenced by G4PenelopeBremsstrahlungModel::BuildXSTable(), G4PenelopeBremsstrahlungModel::ClearTables(), G4MuPairProductionModel::DataCorrupted(), G4PenelopePhotoElectricModel::GetNumberOfShellXS(), G4VMscModel::GetParticleChangeForMSC(), G4BoldyshevTripletModel::Initialise(), G4eSingleCoulombScatteringModel::Initialise(), G4PAIModel::Initialise(), G4PAIPhotModel::Initialise(), G4EmMultiModel::Initialise(), G4mplIonisationModel::Initialise(), G4mplIonisationWithDeltaModel::Initialise(), G4JAEAElasticScatteringModel::Initialise(), G4JAEAPolarizedElasticScatteringModel::Initialise(), G4LivermoreBremsstrahlungModel::Initialise(), G4LivermoreComptonModel::Initialise(), G4LivermoreGammaConversion5DModel::Initialise(), G4LivermoreGammaConversionModel::Initialise(), G4LivermoreNuclearGammaConversionModel::Initialise(), G4LivermorePhotoElectricModel::Initialise(), G4LivermorePolarizedComptonModel::Initialise(), G4LivermorePolarizedGammaConversionModel::Initialise(), G4LivermorePolarizedRayleighModel::Initialise(), G4LivermoreRayleighModel::Initialise(), G4LowEPComptonModel::Initialise(), G4LowEPPolarizedComptonModel::Initialise(), G4PenelopeAnnihilationModel::Initialise(), G4PenelopeBremsstrahlungModel::Initialise(), G4PenelopeComptonModel::Initialise(), G4PenelopeGammaConversionModel::Initialise(), G4PenelopeIonisationModel::Initialise(), G4PenelopePhotoElectricModel::Initialise(), Initialise(), G4PenelopeRayleighModelMI::Initialise(), G4MuBremsstrahlungModel::Initialise(), G4MuPairProductionModel::Initialise(), G4BetheBlochModel::Initialise(), G4BetheHeitlerModel::Initialise(), G4BraggIonModel::Initialise(), G4BraggModel::Initialise(), G4eBremParametrizedModel::Initialise(), G4eBremsstrahlungRelModel::Initialise(), G4eCoulombScatteringModel::Initialise(), G4eDPWACoulombScatteringModel::Initialise(), G4eeToTwoGammaModel::Initialise(), G4eplusTo2GammaOKVIModel::Initialise(), G4GoudsmitSaundersonMscModel::Initialise(), G4hCoulombScatteringModel::Initialise(), G4KleinNishinaCompton::Initialise(), G4KleinNishinaModel::Initialise(), G4LindhardSorensenIonModel::Initialise(), G4PairProductionRelModel::Initialise(), G4SeltzerBergerModel::Initialise(), G4WentzelVIModel::Initialise(), G4PenelopeBremsstrahlungModel::InitialiseLocal(), G4PenelopeGammaConversionModel::ReadDataFile(), G4PenelopePhotoElectricModel::ReadDataFile(), G4BetheHeitlerModel::~G4BetheHeitlerModel(), G4BoldyshevTripletModel::~G4BoldyshevTripletModel(), G4BraggIonModel::~G4BraggIonModel(), G4BraggModel::~G4BraggModel(), G4eBremsstrahlungRelModel::~G4eBremsstrahlungRelModel(), G4eDPWACoulombScatteringModel::~G4eDPWACoulombScatteringModel(), G4GoudsmitSaundersonMscModel::~G4GoudsmitSaundersonMscModel(), G4JAEAElasticScatteringModel::~G4JAEAElasticScatteringModel(), G4JAEAPolarizedElasticScatteringModel::~G4JAEAPolarizedElasticScatteringModel(), G4LivermoreBremsstrahlungModel::~G4LivermoreBremsstrahlungModel(), G4LivermoreComptonModel::~G4LivermoreComptonModel(), G4LivermoreGammaConversion5DModel::~G4LivermoreGammaConversion5DModel(), G4LivermoreGammaConversionModel::~G4LivermoreGammaConversionModel(), G4LivermoreNuclearGammaConversionModel::~G4LivermoreNuclearGammaConversionModel(), G4LivermorePhotoElectricModel::~G4LivermorePhotoElectricModel(), G4LivermorePolarizedComptonModel::~G4LivermorePolarizedComptonModel(), G4LivermorePolarizedGammaConversionModel::~G4LivermorePolarizedGammaConversionModel(), G4LivermorePolarizedRayleighModel::~G4LivermorePolarizedRayleighModel(), G4LivermoreRayleighModel::~G4LivermoreRayleighModel(), G4LowEPComptonModel::~G4LowEPComptonModel(), G4LowEPPolarizedComptonModel::~G4LowEPPolarizedComptonModel(), G4mplIonisationModel::~G4mplIonisationModel(), G4mplIonisationWithDeltaModel::~G4mplIonisationWithDeltaModel(), G4PAIModel::~G4PAIModel(), G4PAIPhotModel::~G4PAIPhotModel(), G4PairProductionRelModel::~G4PairProductionRelModel(), G4PenelopeBremsstrahlungModel::~G4PenelopeBremsstrahlungModel(), G4PenelopeGammaConversionModel::~G4PenelopeGammaConversionModel(), G4PenelopeIonisationModel::~G4PenelopeIonisationModel(), G4PenelopePhotoElectricModel::~G4PenelopePhotoElectricModel(), ~G4PenelopeRayleighModel(), G4PenelopeRayleighModelMI::~G4PenelopeRayleighModelMI(), G4SeltzerBergerModel::~G4SeltzerBergerModel(), and G4WentzelVIModel::~G4WentzelVIModel().

◆ LowEnergyActivationLimit()

G4double G4VEmModel::LowEnergyActivationLimit ( ) const
inlineinherited

◆ LowEnergyLimit()

G4double G4VEmModel::LowEnergyLimit ( ) const
inlineinherited

Definition at line 662 of file G4VEmModel.hh.

663{
664 return lowLimit;
665}

References G4VEmModel::lowLimit.

Referenced by G4eBremsstrahlungRelModel::ComputeCrossSectionPerAtom(), G4KleinNishinaCompton::ComputeCrossSectionPerAtom(), G4KleinNishinaModel::ComputeCrossSectionPerAtom(), G4LivermoreComptonModel::ComputeCrossSectionPerAtom(), G4LivermorePolarizedComptonModel::ComputeCrossSectionPerAtom(), G4LowEPComptonModel::ComputeCrossSectionPerAtom(), G4LowEPPolarizedComptonModel::ComputeCrossSectionPerAtom(), G4mplIonisationWithDeltaModel::ComputeCrossSectionPerElectron(), G4EmCalculator::ComputeDEDX(), G4eBremsstrahlungRelModel::ComputeDEDXPerVolume(), G4mplIonisationWithDeltaModel::ComputeDEDXPerVolume(), G4IonParametrisedLossModel::ComputeDEDXPerVolume(), G4IonParametrisedLossModel::CorrectionsAlongStep(), G4PenelopeRayleighModelMI::CrossSectionPerVolume(), G4PenelopeComptonModel::CrossSectionPerVolume(), G4DNAChampionElasticModel::CrossSectionPerVolume(), G4DNACPA100ElasticModel::CrossSectionPerVolume(), G4DNACPA100ExcitationModel::CrossSectionPerVolume(), G4DNACPA100IonisationModel::CrossSectionPerVolume(), G4DNAEmfietzoglouExcitationModel::CrossSectionPerVolume(), G4DNAEmfietzoglouIonisationModel::CrossSectionPerVolume(), G4DNAMeltonAttachmentModel::CrossSectionPerVolume(), G4DNASancheExcitationModel::CrossSectionPerVolume(), G4DNAScreenedRutherfordElasticModel::CrossSectionPerVolume(), G4DNAELSEPAElasticModel::CrossSectionPerVolume(), G4EmCalculator::FindEmModel(), G4DNAChampionElasticModel::G4DNAChampionElasticModel(), G4DNACPA100ElasticModel::G4DNACPA100ElasticModel(), G4DNAMeltonAttachmentModel::G4DNAMeltonAttachmentModel(), G4DNASancheExcitationModel::G4DNASancheExcitationModel(), G4DNAUeharaScreenedRutherfordElasticModel::G4DNAUeharaScreenedRutherfordElasticModel(), G4eeToHadronsModel::G4eeToHadronsModel(), G4LivermorePolarizedRayleighModel::G4LivermorePolarizedRayleighModel(), G4PenelopeBremsstrahlungModel::GetCrossSectionTableForCouple(), G4VMscModel::GetParticleChangeForMSC(), G4DNAChampionElasticModel::Initialise(), G4DNACPA100ElasticModel::Initialise(), G4DNADingfelderChargeDecreaseModel::Initialise(), G4DNADingfelderChargeIncreaseModel::Initialise(), G4DNAMeltonAttachmentModel::Initialise(), G4DNAMillerGreenExcitationModel::Initialise(), G4DNARuddIonisationExtendedModel::Initialise(), G4DNARuddIonisationModel::Initialise(), G4DNASancheExcitationModel::Initialise(), G4DNAScreenedRutherfordElasticModel::Initialise(), G4DNAUeharaScreenedRutherfordElasticModel::Initialise(), G4BoldyshevTripletModel::Initialise(), G4PAIModel::Initialise(), G4PAIPhotModel::Initialise(), G4JAEAElasticScatteringModel::Initialise(), G4JAEAPolarizedElasticScatteringModel::Initialise(), G4LivermoreComptonModel::Initialise(), G4LivermoreGammaConversion5DModel::Initialise(), G4LivermoreGammaConversionModel::Initialise(), G4LivermoreIonisationModel::Initialise(), G4LivermoreNuclearGammaConversionModel::Initialise(), G4LivermorePolarizedComptonModel::Initialise(), G4LivermorePolarizedGammaConversionModel::Initialise(), G4LivermoreRayleighModel::Initialise(), G4LowEPComptonModel::Initialise(), G4LowEPPolarizedComptonModel::Initialise(), G4MicroElecElasticModel::Initialise(), G4MicroElecElasticModel_new::Initialise(), G4MicroElecInelasticModel::Initialise(), G4MicroElecInelasticModel_new::Initialise(), G4PenelopeAnnihilationModel::Initialise(), G4PenelopeBremsstrahlungModel::Initialise(), G4PenelopeComptonModel::Initialise(), G4PenelopeGammaConversionModel::Initialise(), G4PenelopeIonisationModel::Initialise(), G4PenelopePhotoElectricModel::Initialise(), Initialise(), G4PenelopeRayleighModelMI::Initialise(), G4MuPairProductionModel::Initialise(), G4eBremParametrizedModel::Initialise(), G4eBremsstrahlungRelModel::Initialise(), G4eDPWACoulombScatteringModel::Initialise(), G4GoudsmitSaundersonMscModel::Initialise(), G4PairProductionRelModel::Initialise(), G4SeltzerBergerModel::Initialise(), G4WentzelVIModel::Initialise(), G4DNABornExcitationModel1::Initialise(), G4DNABornExcitationModel2::Initialise(), G4DNABornIonisationModel1::Initialise(), G4DNABornIonisationModel2::Initialise(), G4DNACPA100ExcitationModel::Initialise(), G4DNACPA100IonisationModel::Initialise(), G4DNADiracRMatrixExcitationModel::Initialise(), G4DNAEmfietzoglouExcitationModel::Initialise(), G4DNAEmfietzoglouIonisationModel::Initialise(), G4DNAQuinnPlasmonExcitationModel::Initialise(), G4DNARelativisticIonisationModel::Initialise(), G4EmModelManager::Initialise(), G4DNAELSEPAElasticModel::Initialise(), G4DNADummyModel::Initialise(), G4DNAIonElasticModel::Initialise(), G4mplIonisation::InitialiseEnergyLossProcess(), G4eBremsstrahlungRelModel::InitialiseLocal(), G4PairProductionRelModel::InitialiseLocal(), G4CoulombScattering::InitialiseProcess(), G4VEmProcess::PostStepDoIt(), G4eBremsstrahlungRelModel::SampleSecondaries(), G4DNACPA100IonisationModel::SampleSecondaries(), G4DNAEmfietzoglouIonisationModel::SampleSecondaries(), G4LivermoreComptonModel::SampleSecondaries(), G4LivermorePolarizedComptonModel::SampleSecondaries(), G4LowEPComptonModel::SampleSecondaries(), G4LowEPPolarizedComptonModel::SampleSecondaries(), G4PenelopeComptonModel::SampleSecondaries(), G4PolarizedComptonModel::SampleSecondaries(), G4KleinNishinaCompton::SampleSecondaries(), G4KleinNishinaModel::SampleSecondaries(), G4EmConfigurator::SetExtraEmModel(), G4mplIonisationModel::SetParticle(), G4mplIonisationWithDeltaModel::SetParticle(), and G4EmConfigurator::UpdateModelEnergyRange().

◆ LPMFlag()

G4bool G4VEmModel::LPMFlag ( ) const
inlineinherited

◆ MaxSecondaryEnergy()

G4double G4VEmModel::MaxSecondaryEnergy ( const G4ParticleDefinition ,
G4double  kineticEnergy 
)
protectedvirtualinherited

◆ MaxSecondaryKinEnergy()

G4double G4VEmModel::MaxSecondaryKinEnergy ( const G4DynamicParticle dynParticle)
inlineinherited

◆ MinEnergyCut()

G4double G4VEmModel::MinEnergyCut ( const G4ParticleDefinition ,
const G4MaterialCutsCouple  
)
virtualinherited

◆ MinPrimaryEnergy()

G4double G4VEmModel::MinPrimaryEnergy ( const G4Material ,
const G4ParticleDefinition ,
G4double  cut = 0.0 
)
virtualinherited

◆ ModelDescription()

void G4VEmModel::ModelDescription ( std::ostream &  outFile) const
virtualinherited

Reimplemented in G4eeToHadronsMultiModel.

Definition at line 469 of file G4VEmModel.cc.

470{
471 outFile << "The description for this model has not been written yet.\n";
472}

◆ operator=()

G4PenelopeRayleighModel & G4PenelopeRayleighModel::operator= ( const G4PenelopeRayleighModel right)
delete

◆ PolarAngleLimit()

G4double G4VEmModel::PolarAngleLimit ( ) const
inlineinherited

◆ ReadDataFile()

void G4PenelopeRayleighModel::ReadDataFile ( G4int  Z)
private

Definition at line 546 of file G4PenelopeRayleighModel.cc.

547{
548 if (fVerboseLevel > 2)
549 {
550 G4cout << "G4PenelopeRayleighModel::ReadDataFile()" << G4endl;
551 G4cout << "Going to read Rayleigh data files for Z=" << Z << G4endl;
552 }
553 char* path = std::getenv("G4LEDATA");
554 if (!path)
555 {
556 G4String excep = "G4LEDATA environment variable not set!";
557 G4Exception("G4PenelopeRayleighModel::ReadDataFile()",
558 "em0006",FatalException,excep);
559 return;
560 }
561
562 /*
563 Read first the cross section file
564 */
565 std::ostringstream ost;
566 if (Z>9)
567 ost << path << "/penelope/rayleigh/pdgra" << Z << ".p08";
568 else
569 ost << path << "/penelope/rayleigh/pdgra0" << Z << ".p08";
570 std::ifstream file(ost.str().c_str());
571 if (!file.is_open())
572 {
573 G4String excep = "Data file " + G4String(ost.str()) + " not found!";
574 G4Exception("G4PenelopeRayleighModel::ReadDataFile()",
575 "em0003",FatalException,excep);
576 }
577 G4int readZ =0;
578 size_t nPoints= 0;
579 file >> readZ >> nPoints;
580 //check the right file is opened.
581 if (readZ != Z || nPoints <= 0 || nPoints >= 5000)
582 {
584 ed << "Corrupted data file for Z=" << Z << G4endl;
585 G4Exception("G4PenelopeRayleighModel::ReadDataFile()",
586 "em0005",FatalException,ed);
587 return;
588 }
589
590 fLogAtomicCrossSection[Z] = new G4PhysicsFreeVector((size_t)nPoints);
591 G4double ene=0,f1=0,f2=0,xs=0;
592 for (size_t i=0;i<nPoints;i++)
593 {
594 file >> ene >> f1 >> f2 >> xs;
595 //dimensional quantities
596 ene *= eV;
597 xs *= cm2;
599 if (file.eof() && i != (nPoints-1)) //file ended too early
600 {
602 ed << "Corrupted data file for Z=" << Z << G4endl;
603 ed << "Found less than " << nPoints << "entries " <<G4endl;
604 G4Exception("G4PenelopeRayleighModel::ReadDataFile()",
605 "em0005",FatalException,ed);
606 }
607 }
608 file.close();
609
610 /*
611 Then read the form factor file
612 */
613 std::ostringstream ost2;
614 if (Z>9)
615 ost2 << path << "/penelope/rayleigh/pdaff" << Z << ".p08";
616 else
617 ost2 << path << "/penelope/rayleigh/pdaff0" << Z << ".p08";
618 file.open(ost2.str().c_str());
619 if (!file.is_open())
620 {
621 G4String excep = "Data file " + G4String(ost2.str()) + " not found!";
622 G4Exception("G4PenelopeRayleighModel::ReadDataFile()",
623 "em0003",FatalException,excep);
624 }
625 file >> readZ >> nPoints;
626 //check the right file is opened.
627 if (readZ != Z || nPoints <= 0 || nPoints >= 5000)
628 {
630 ed << "Corrupted data file for Z=" << Z << G4endl;
631 G4Exception("G4PenelopeRayleighModel::ReadDataFile()",
632 "em0005",FatalException,ed);
633 return;
634 }
635 fAtomicFormFactor[Z] = new G4PhysicsFreeVector((size_t)nPoints);
636 G4double q=0,ff=0,incoh=0;
637 G4bool fillQGrid = false;
638 //fill this vector only the first time.
639 if (!fLogQSquareGrid.size())
640 fillQGrid = true;
641 for (size_t i=0;i<nPoints;i++)
642 {
643 file >> q >> ff >> incoh;
644 //q and ff are dimensionless (q is in units of (m_e*c)
645 fAtomicFormFactor[Z]->PutValue(i,q,ff);
646 if (fillQGrid)
647 {
648 fLogQSquareGrid.push_back(2.0*G4Log(q));
649 }
650 if (file.eof() && i != (nPoints-1)) //file ended too early
651 {
653 ed << "Corrupted data file for Z=" << Z << G4endl;
654 ed << "Found less than " << nPoints << "entries " <<G4endl;
655 G4Exception("G4PenelopeRayleighModel::ReadDataFile()",
656 "em0005",FatalException,ed);
657 }
658 }
659 file.close();
660 return;
661}
static constexpr double cm2
Definition: G4SIunits.hh:100

References cm2, eV, FatalException, fAtomicFormFactor, geant4_check_module_cycles::file, fLogAtomicCrossSection, fLogQSquareGrid, fVerboseLevel, G4cout, G4endl, G4Exception(), G4Log(), G4PhysicsFreeVector::PutValue(), and Z.

Referenced by ComputeCrossSectionPerAtom(), Initialise(), and SampleSecondaries().

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 382 of file G4PenelopeRayleighModel.cc.

387{
388 // Sampling of the Rayleigh final state (namely, scattering angle of the photon)
389 // from the Penelope2008 model. The scattering angle is sampled from the atomic
390 // cross section dOmega/d(cosTheta) from Born ("Atomic Phyisics", 1969), disregarding
391 // anomalous scattering effects. The Form Factor F(Q) function which appears in the
392 // analytical cross section is retrieved via the method GetFSquared(); atomic data
393 // are tabulated for F(Q). Form factor for compounds is calculated according to
394 // the additivity rule. The sampling from the F(Q) is made via a Rational Inverse
395 // Transform with Aliasing (RITA) algorithm; RITA parameters are calculated once
396 // for each material and managed by G4PenelopeSamplingData objects.
397 // The sampling algorithm (rejection method) has efficiency 67% at low energy, and
398 // increases with energy. For E=100 keV the efficiency is 100% and 86% for
399 // hydrogen and uranium, respectively.
400
401 if (fVerboseLevel > 3)
402 G4cout << "Calling SamplingSecondaries() of G4PenelopeRayleighModel" << G4endl;
403
404 G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy();
405
406 if (photonEnergy0 <= fIntrinsicLowEnergyLimit)
407 {
411 return ;
412 }
413
414 G4ParticleMomentum photonDirection0 = aDynamicGamma->GetMomentumDirection();
415
416 const G4Material* theMat = couple->GetMaterial();
417
418 //1) Verify if tables are ready
419 //Either Initialize() was not called, or we are in a slave and InitializeLocal() was
420 //not invoked
422 {
423 //create a **thread-local** version of the table. Used only for G4EmCalculator and
424 //Unit Tests
425 fLocalTable = true;
427 fLogFormFactorTable = new std::map<const G4Material*,G4PhysicsFreeVector*>;
428 if (!fPMaxTable)
429 fPMaxTable = new std::map<const G4Material*,G4PhysicsFreeVector*>;
430 if (!fSamplingTable)
431 fSamplingTable = new std::map<const G4Material*,G4PenelopeSamplingData*>;
432 }
433
434 if (!fSamplingTable->count(theMat))
435 {
436 //If we are here, it means that Initialize() was inkoved, but the MaterialTable was
437 //not filled up. This can happen in a UnitTest
438 if (fVerboseLevel > 0)
439 {
440 //Issue a G4Exception (warning) only in verbose mode
442 ed << "Unable to find the fSamplingTable data for " <<
443 theMat->GetName() << G4endl;
444 ed << "This can happen only in Unit Tests" << G4endl;
445 G4Exception("G4PenelopeRayleighModel::SampleSecondaries()",
446 "em2019",JustWarning,ed);
447 }
448 const G4ElementVector* theElementVector = theMat->GetElementVector();
449 //protect file reading via autolock
451 for (size_t j=0;j<theMat->GetNumberOfElements();j++)
452 {
453 G4int iZ = theElementVector->at(j)->GetZasInt();
454 if (!fLogAtomicCrossSection[iZ])
455 {
456 lock.lock();
457 ReadDataFile(iZ);
458 lock.unlock();
459 }
460 }
461 lock.lock();
462 //1) If the table has not been built for the material, do it!
463 if (!fLogFormFactorTable->count(theMat))
464 BuildFormFactorTable(theMat);
465
466 //2) retrieve or build the sampling table
467 if (!(fSamplingTable->count(theMat)))
469
470 //3) retrieve or build the pMax data
471 if (!fPMaxTable->count(theMat))
472 GetPMaxTable(theMat);
473 lock.unlock();
474 }
475
476 //Ok, restart the job
477 G4PenelopeSamplingData* theDataTable = fSamplingTable->find(theMat)->second;
478 G4PhysicsFreeVector* thePMax = fPMaxTable->find(theMat)->second;
479
480 G4double cosTheta = 1.0;
481
482 //OK, ready to go!
483 G4double qmax = 2.0*photonEnergy0/electron_mass_c2; //this is non-dimensional now
484
485 if (qmax < 1e-10) //very low momentum transfer
486 {
487 G4bool loopAgain=false;
488 do
489 {
490 loopAgain = false;
491 cosTheta = 1.0-2.0*G4UniformRand();
492 G4double G = 0.5*(1+cosTheta*cosTheta);
493 if (G4UniformRand()>G)
494 loopAgain = true;
495 }while(loopAgain);
496 }
497 else //larger momentum transfer
498 {
499 size_t nData = theDataTable->GetNumberOfStoredPoints();
500 G4double LastQ2inTheTable = theDataTable->GetX(nData-1);
501 G4double q2max = std::min(qmax*qmax,LastQ2inTheTable);
502
503 G4bool loopAgain = false;
504 G4double MaxPValue = thePMax->Value(photonEnergy0);
505 G4double xx=0;
506
507 //Sampling by rejection method. The rejection function is
508 //G = 0.5*(1+cos^2(theta))
509 //
510 do{
511 loopAgain = false;
512 G4double RandomMax = G4UniformRand()*MaxPValue;
513 xx = theDataTable->SampleValue(RandomMax);
514 //xx is a random value of q^2 in (0,q2max),sampled according to
515 //F(Q^2) via the RITA algorithm
516 if (xx > q2max)
517 loopAgain = true;
518 cosTheta = 1.0-2.0*xx/q2max;
519 G4double G = 0.5*(1+cosTheta*cosTheta);
520 if (G4UniformRand()>G)
521 loopAgain = true;
522 }while(loopAgain);
523 }
524
525 G4double sinTheta = std::sqrt(1-cosTheta*cosTheta);
526
527 // Scattered photon angles. ( Z - axis along the parent photon)
528 G4double phi = twopi * G4UniformRand() ;
529 G4double dirX = sinTheta*std::cos(phi);
530 G4double dirY = sinTheta*std::sin(phi);
531 G4double dirZ = cosTheta;
532
533 // Update G4VParticleChange for the scattered photon
534 G4ThreeVector photonDirection1(dirX, dirY, dirZ);
535
536 photonDirection1.rotateUz(photonDirection0);
537 fParticleChange->ProposeMomentumDirection(photonDirection1) ;
539
540 return;
541}
static constexpr double twopi
Definition: G4SIunits.hh:56
@ fStopAndKill
#define G4UniformRand()
Definition: Randomize.hh:52
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:186
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4double SampleValue(G4double rndm)
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

References BuildFormFactorTable(), source.hepunit::electron_mass_c2, fIntrinsicLowEnergyLimit, fLocalTable, fLogAtomicCrossSection, fLogFormFactorTable, fParticleChange, fPMaxTable, fSamplingTable, fStopAndKill, fVerboseLevel, G4cout, G4endl, G4Exception(), G4UniformRand, G4Material::GetElementVector(), G4DynamicParticle::GetKineticEnergy(), G4MaterialCutsCouple::GetMaterial(), G4DynamicParticle::GetMomentumDirection(), G4Material::GetName(), G4Material::GetNumberOfElements(), G4PenelopeSamplingData::GetNumberOfStoredPoints(), GetPMaxTable(), G4PenelopeSamplingData::GetX(), InitializeSamplingAlgorithm(), JustWarning, G4TemplateAutoLock< _Mutex_t >::lock(), G4INCL::Math::min(), anonymous_namespace{G4PenelopeRayleighModel.cc}::PenelopeRayleighModelMutex, G4VParticleChange::ProposeLocalEnergyDeposit(), G4ParticleChangeForGamma::ProposeMomentumDirection(), G4VParticleChange::ProposeTrackStatus(), ReadDataFile(), CLHEP::Hep3Vector::rotateUz(), G4PenelopeSamplingData::SampleValue(), G4ParticleChangeForGamma::SetProposedKineticEnergy(), twopi, G4TemplateAutoLock< _Mutex_t >::unlock(), and G4PhysicsVector::Value().

◆ SecondaryThreshold()

G4double G4VEmModel::SecondaryThreshold ( ) const
inlineinherited

◆ SelectIsotopeNumber()

G4int G4VEmModel::SelectIsotopeNumber ( const G4Element elm)
inherited

Definition at line 319 of file G4VEmModel.cc.

320{
322 const size_t ni = elm->GetNumberOfIsotopes();
323 fCurrentIsotope = elm->GetIsotope(0);
324 size_t idx = 0;
325 if(ni > 1) {
326 const G4double* ab = elm->GetRelativeAbundanceVector();
328 for(; idx<ni; ++idx) {
329 x -= ab[idx];
330 if (x <= 0.0) {
331 fCurrentIsotope = elm->GetIsotope(idx);
332 break;
333 }
334 }
335 }
336 return fCurrentIsotope->GetN();
337}
static const G4double ab
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:167
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:170
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:159
G4int GetN() const
Definition: G4Isotope.hh:93

References ab, G4VEmModel::fCurrentIsotope, G4UniformRand, G4Element::GetIsotope(), G4Isotope::GetN(), G4Element::GetNumberOfIsotopes(), G4Element::GetRelativeAbundanceVector(), and G4VEmModel::SetCurrentElement().

Referenced by G4eSingleCoulombScatteringModel::SampleSecondaries(), G4IonCoulombScatteringModel::SampleSecondaries(), G4eCoulombScatteringModel::SampleSecondaries(), G4hCoulombScatteringModel::SampleSecondaries(), and G4BetheHeitler5DModel::SampleSecondaries().

◆ SelectRandomAtom() [1/2]

const G4Element * G4VEmModel::SelectRandomAtom ( const G4Material mat,
const G4ParticleDefinition pd,
G4double  kineticEnergy,
G4double  cutEnergy = 0.0,
G4double  maxEnergy = DBL_MAX 
)
inherited

Definition at line 275 of file G4VEmModel.cc.

280{
281 size_t n = mat->GetNumberOfElements();
282 fCurrentElement = mat->GetElement(0);
283 if (n > 1) {
284 const G4double x = G4UniformRand()*
285 G4VEmModel::CrossSectionPerVolume(mat,pd,kinEnergy,tcut,tmax);
286 for(size_t i=0; i<n; ++i) {
287 if (x <= xsec[i]) {
288 fCurrentElement = mat->GetElement(i);
289 break;
290 }
291 }
292 }
293 return fCurrentElement;
294}

References G4VEmModel::CrossSectionPerVolume(), G4VEmModel::fCurrentElement, G4UniformRand, G4Material::GetElement(), G4Material::GetNumberOfElements(), CLHEP::detail::n, and G4VEmModel::xsec.

◆ SelectRandomAtom() [2/2]

const G4Element * G4VEmModel::SelectRandomAtom ( const G4MaterialCutsCouple couple,
const G4ParticleDefinition part,
G4double  kineticEnergy,
G4double  cutEnergy = 0.0,
G4double  maxEnergy = DBL_MAX 
)
inlineinherited

Definition at line 580 of file G4VEmModel.hh.

585{
586 SetCurrentCouple(couple);
588 ((*elmSelectors)[couple->GetIndex()])->SelectRandomAtom(kinEnergy) :
589 SelectRandomAtom(pBaseMaterial,part,kinEnergy,cutEnergy,maxEnergy);
590 fCurrentIsotope = nullptr;
591 return fCurrentElement;
592}
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:580

References G4VEmModel::elmSelectors, G4VEmModel::fCurrentElement, G4VEmModel::fCurrentIsotope, G4MaterialCutsCouple::GetIndex(), G4VEmModel::nSelectors, G4VEmModel::pBaseMaterial, G4VEmModel::SelectRandomAtom(), and G4VEmModel::SetCurrentCouple().

Referenced by G4AdjointBremsstrahlungModel::RapidSampleSecondaries(), G4LivermoreBremsstrahlungModel::SampleSecondaries(), G4JAEAElasticScatteringModel::SampleSecondaries(), G4JAEAPolarizedElasticScatteringModel::SampleSecondaries(), G4LivermoreComptonModel::SampleSecondaries(), G4LivermoreNuclearGammaConversionModel::SampleSecondaries(), G4LivermorePhotoElectricModel::SampleSecondaries(), G4LivermorePolarizedComptonModel::SampleSecondaries(), G4LivermorePolarizedGammaConversionModel::SampleSecondaries(), G4LivermorePolarizedRayleighModel::SampleSecondaries(), G4LivermoreRayleighModel::SampleSecondaries(), G4LowEPComptonModel::SampleSecondaries(), G4LowEPPolarizedComptonModel::SampleSecondaries(), G4PenelopePhotoElectricModel::SampleSecondaries(), G4MuBremsstrahlungModel::SampleSecondaries(), G4MuPairProductionModel::SampleSecondaries(), G4hCoulombScatteringModel::SampleSecondaries(), G4KleinNishinaModel::SampleSecondaries(), G4PEEffectFluoModel::SampleSecondaries(), G4VEmModel::SelectRandomAtom(), and G4VEmModel::SelectTargetAtom().

◆ SelectRandomAtomNumber()

G4int G4VEmModel::SelectRandomAtomNumber ( const G4Material mat)
inherited

◆ SelectTargetAtom()

const G4Element * G4VEmModel::SelectTargetAtom ( const G4MaterialCutsCouple couple,
const G4ParticleDefinition part,
G4double  kineticEnergy,
G4double  logKineticEnergy,
G4double  cutEnergy = 0.0,
G4double  maxEnergy = DBL_MAX 
)
inlineinherited

◆ SetActivationHighEnergyLimit()

void G4VEmModel::SetActivationHighEnergyLimit ( G4double  val)
inlineinherited

◆ SetActivationLowEnergyLimit()

void G4VEmModel::SetActivationLowEnergyLimit ( G4double  val)
inlineinherited

◆ SetAngularDistribution()

void G4VEmModel::SetAngularDistribution ( G4VEmAngularDistribution p)
inlineinherited

Definition at line 628 of file G4VEmModel.hh.

629{
630 if(p != anglModel) {
631 delete anglModel;
632 anglModel = p;
633 }
634}

References G4VEmModel::anglModel.

Referenced by G4EmLivermorePhysics::ConstructProcess(), G4EmLowEPPhysics::ConstructProcess(), G4EmStandardPhysics::ConstructProcess(), G4EmStandardPhysics_option3::ConstructProcess(), G4EmStandardPhysics_option4::ConstructProcess(), G4EmStandardPhysicsSS::ConstructProcess(), G4BetheHeitlerModel::G4BetheHeitlerModel(), G4DNABornIonisationModel1::G4DNABornIonisationModel1(), G4DNABornIonisationModel2::G4DNABornIonisationModel2(), G4DNAEmfietzoglouIonisationModel::G4DNAEmfietzoglouIonisationModel(), G4DNARuddIonisationExtendedModel::G4DNARuddIonisationExtendedModel(), G4DNARuddIonisationModel::G4DNARuddIonisationModel(), G4eBremParametrizedModel::G4eBremParametrizedModel(), G4eBremsstrahlungRelModel::G4eBremsstrahlungRelModel(), G4IonParametrisedLossModel::G4IonParametrisedLossModel(), G4LivermoreBremsstrahlungModel::G4LivermoreBremsstrahlungModel(), G4LivermoreIonisationModel::G4LivermoreIonisationModel(), G4LivermorePhotoElectricModel::G4LivermorePhotoElectricModel(), G4LivermoreRayleighModel::G4LivermoreRayleighModel(), G4MicroElecInelasticModel::G4MicroElecInelasticModel(), G4MicroElecInelasticModel_new::G4MicroElecInelasticModel_new(), G4MuBremsstrahlungModel::G4MuBremsstrahlungModel(), G4MuPairProductionModel::G4MuPairProductionModel(), G4PAIModel::G4PAIModel(), G4PAIPhotModel::G4PAIPhotModel(), G4PairProductionRelModel::G4PairProductionRelModel(), G4PEEffectFluoModel::G4PEEffectFluoModel(), G4SeltzerBergerModel::G4SeltzerBergerModel(), G4AtimaEnergyLossModel::Initialise(), G4BetheBlochModel::Initialise(), G4BraggIonModel::Initialise(), G4BraggModel::Initialise(), G4ICRU73QOModel::Initialise(), G4LindhardSorensenIonModel::Initialise(), and G4MollerBhabhaModel::Initialise().

◆ SetAngularGeneratorFlag()

void G4VEmModel::SetAngularGeneratorFlag ( G4bool  val)
inlineinherited

Definition at line 725 of file G4VEmModel.hh.

726{
728}
G4bool useAngularGenerator
Definition: G4VEmModel.hh:461

References G4VEmModel::useAngularGenerator.

Referenced by G4VEnergyLossProcess::PreparePhysicsTable().

◆ SetCrossSectionTable()

void G4VEmModel::SetCrossSectionTable ( G4PhysicsTable p,
G4bool  isLocal 
)
inherited

Definition at line 455 of file G4VEmModel.cc.

456{
457 if(p != xSectionTable) {
458 if(xSectionTable != nullptr && localTable) {
460 delete xSectionTable;
461 }
462 xSectionTable = p;
463 }
464 localTable = isLocal;
465}
void clearAndDestroy()
G4bool localTable
Definition: G4VEmModel.hh:459

References G4PhysicsTable::clearAndDestroy(), G4VEmModel::localTable, and G4VEmModel::xSectionTable.

Referenced by G4VMultipleScattering::BuildPhysicsTable().

◆ SetCurrentCouple()

void G4VEmModel::SetCurrentCouple ( const G4MaterialCutsCouple ptr)
inlineinherited

Definition at line 472 of file G4VEmModel.hh.

473{
474 if(fCurrentCouple != ptr) {
475 fCurrentCouple = ptr;
477 pBaseMaterial = ptr->GetMaterial();
478 pFactor = 1.0;
479 if(useBaseMaterials) {
480 basedCoupleIndex = (*theDensityIdx)[currentCoupleIndex];
481 if(nullptr != pBaseMaterial->GetBaseMaterial())
483 pFactor = (*theDensityFactor)[currentCoupleIndex];
484 }
485 }
486}
const G4Material * GetBaseMaterial() const
Definition: G4Material.hh:229
G4bool useBaseMaterials
Definition: G4VEmModel.hh:462
size_t currentCoupleIndex
Definition: G4VEmModel.hh:448
size_t basedCoupleIndex
Definition: G4VEmModel.hh:449

References G4VEmModel::basedCoupleIndex, G4VEmModel::currentCoupleIndex, G4VEmModel::fCurrentCouple, G4Material::GetBaseMaterial(), G4MaterialCutsCouple::GetIndex(), G4MaterialCutsCouple::GetMaterial(), G4VEmModel::pBaseMaterial, G4VEmModel::pFactor, and G4VEmModel::useBaseMaterials.

Referenced by G4VMultipleScattering::AlongStepGetPhysicalInteractionLength(), G4EmMultiModel::ComputeCrossSectionPerAtom(), G4VEmModel::ComputeDEDX(), G4TablesForExtrapolator::ComputeTrasportXS(), G4UrbanAdjointMscModel::ComputeTruePathLengthLimit(), G4GoudsmitSaundersonMscModel::ComputeTruePathLengthLimit(), G4UrbanMscModel::ComputeTruePathLengthLimit(), G4VEmModel::CrossSection(), G4AdjointPhotoElectricModel::DefineCurrentMaterialAndElectronEnergy(), G4WentzelVIModel::DefineMaterial(), G4WentzelVIRelModel::DefineMaterial(), G4EmCalculator::GetCrossSectionPerVolume(), G4LivermoreGammaConversion5DModel::Initialise(), G4VEmModel::InitialiseElementSelectors(), G4PEEffectFluoModel::SampleSecondaries(), G4EmMultiModel::SampleSecondaries(), G4VEnergyLossProcess::SelectModel(), G4VEmProcess::SelectModel(), G4VEmModel::SelectRandomAtom(), G4VEmModel::SelectTargetAtom(), and G4VEmModel::Value().

◆ SetCurrentElement()

void G4VEmModel::SetCurrentElement ( const G4Element elm)
inlineprotectedinherited

◆ SetDeexcitationFlag()

void G4VEmModel::SetDeexcitationFlag ( G4bool  val)
inlineinherited

Definition at line 823 of file G4VEmModel.hh.

824{
825 flagDeexcitation = val;
826}

References G4VEmModel::flagDeexcitation.

Referenced by G4DNABornIonisationModel1::G4DNABornIonisationModel1(), G4DNABornIonisationModel2::G4DNABornIonisationModel2(), G4DNACPA100IonisationModel::G4DNACPA100IonisationModel(), G4DNAEmfietzoglouIonisationModel::G4DNAEmfietzoglouIonisationModel(), G4DNARelativisticIonisationModel::G4DNARelativisticIonisationModel(), G4DNARuddIonisationExtendedModel::G4DNARuddIonisationExtendedModel(), G4DNARuddIonisationModel::G4DNARuddIonisationModel(), G4KleinNishinaModel::G4KleinNishinaModel(), G4LEPTSIonisationModel::G4LEPTSIonisationModel(), G4LivermoreComptonModel::G4LivermoreComptonModel(), G4LivermorePhotoElectricModel::G4LivermorePhotoElectricModel(), G4LivermorePolarizedComptonModel::G4LivermorePolarizedComptonModel(), G4LowEPComptonModel::G4LowEPComptonModel(), G4LowEPPolarizedComptonModel::G4LowEPPolarizedComptonModel(), G4MicroElecInelasticModel::G4MicroElecInelasticModel(), G4MicroElecInelasticModel_new::G4MicroElecInelasticModel_new(), G4PEEffectFluoModel::G4PEEffectFluoModel(), G4PenelopeBremsstrahlungModel::G4PenelopeBremsstrahlungModel(), G4PenelopeComptonModel::G4PenelopeComptonModel(), G4PenelopeIonisationModel::G4PenelopeIonisationModel(), G4PenelopePhotoElectricModel::G4PenelopePhotoElectricModel(), G4AtimaEnergyLossModel::Initialise(), G4BetheBlochModel::Initialise(), G4BraggIonModel::Initialise(), G4BraggModel::Initialise(), G4ICRU73QOModel::Initialise(), and G4LindhardSorensenIonModel::Initialise().

◆ SetElementSelectors()

void G4VEmModel::SetElementSelectors ( std::vector< G4EmElementSelector * > *  p)
inlineinherited

Definition at line 852 of file G4VEmModel.hh.

853{
854 if(p != elmSelectors) {
855 elmSelectors = p;
856 nSelectors = (nullptr != elmSelectors) ? G4int(elmSelectors->size()) : 0;
857 localElmSelectors = false;
858 }
859}
G4bool localElmSelectors
Definition: G4VEmModel.hh:460

References G4VEmModel::elmSelectors, G4VEmModel::localElmSelectors, and G4VEmModel::nSelectors.

Referenced by G4eDPWACoulombScatteringModel::InitialiseLocal(), G4eSingleCoulombScatteringModel::InitialiseLocal(), G4PAIModel::InitialiseLocal(), G4PAIPhotModel::InitialiseLocal(), G4JAEAElasticScatteringModel::InitialiseLocal(), G4JAEAPolarizedElasticScatteringModel::InitialiseLocal(), G4LivermoreComptonModel::InitialiseLocal(), G4LivermoreNuclearGammaConversionModel::InitialiseLocal(), G4LivermorePolarizedComptonModel::InitialiseLocal(), G4LivermorePolarizedGammaConversionModel::InitialiseLocal(), G4LivermorePolarizedRayleighModel::InitialiseLocal(), G4LivermoreRayleighModel::InitialiseLocal(), G4LowEPComptonModel::InitialiseLocal(), G4LowEPPolarizedComptonModel::InitialiseLocal(), G4PenelopePhotoElectricModel::InitialiseLocal(), G4MuBremsstrahlungModel::InitialiseLocal(), G4MuPairProductionModel::InitialiseLocal(), G4BetheHeitlerModel::InitialiseLocal(), G4eBremParametrizedModel::InitialiseLocal(), G4eBremsstrahlungRelModel::InitialiseLocal(), G4eCoulombScatteringModel::InitialiseLocal(), G4hCoulombScatteringModel::InitialiseLocal(), G4KleinNishinaCompton::InitialiseLocal(), G4KleinNishinaModel::InitialiseLocal(), and G4PairProductionRelModel::InitialiseLocal().

◆ SetFluctuationFlag()

void G4VEmModel::SetFluctuationFlag ( G4bool  val)
inlineinherited

Definition at line 732 of file G4VEmModel.hh.

733{
734 lossFlucFlag = val;
735}
G4bool lossFlucFlag
Definition: G4VEmModel.hh:450

References G4VEmModel::lossFlucFlag.

Referenced by G4EmCalculator::ComputeNuclearDEDX().

◆ SetForceBuildTable()

void G4VEmModel::SetForceBuildTable ( G4bool  val)
inlineinherited

Definition at line 830 of file G4VEmModel.hh.

831{
833}

References G4VEmModel::flagForceBuildTable.

◆ SetHighEnergyLimit()

void G4VEmModel::SetHighEnergyLimit ( G4double  val)
inlineinherited

Definition at line 767 of file G4VEmModel.hh.

768{
769 highLimit = val;
770}

References G4VEmModel::highLimit.

Referenced by G4EmModelActivator::ActivateEmOptions(), G4EmModelActivator::ActivatePAI(), LBE::ConstructEM(), G4EmDNAPhysics_option2::ConstructProcess(), G4EmDNAPhysics_option5::ConstructProcess(), G4EmDNAPhysics_option7::ConstructProcess(), G4EmDNAPhysics_stationary::ConstructProcess(), G4EmLivermorePhysics::ConstructProcess(), G4EmLowEPPhysics::ConstructProcess(), G4EmPenelopePhysics::ConstructProcess(), G4EmStandardPhysics::ConstructProcess(), G4EmStandardPhysics_option1::ConstructProcess(), G4EmStandardPhysics_option2::ConstructProcess(), G4EmStandardPhysics_option4::ConstructProcess(), G4EmStandardPhysicsGS::ConstructProcess(), G4EmStandardPhysicsWVI::ConstructProcess(), G4BraggIonModel::G4BraggIonModel(), G4BraggModel::G4BraggModel(), G4DNAChampionElasticModel::G4DNAChampionElasticModel(), G4DNACPA100ElasticModel::G4DNACPA100ElasticModel(), G4DNACPA100ExcitationModel::G4DNACPA100ExcitationModel(), G4DNACPA100IonisationModel::G4DNACPA100IonisationModel(), G4DNAELSEPAElasticModel::G4DNAELSEPAElasticModel(), G4DNAEmfietzoglouExcitationModel::G4DNAEmfietzoglouExcitationModel(), G4DNAEmfietzoglouIonisationModel::G4DNAEmfietzoglouIonisationModel(), G4DNAIonElasticModel::G4DNAIonElasticModel(), G4DNAMeltonAttachmentModel::G4DNAMeltonAttachmentModel(), G4DNASancheExcitationModel::G4DNASancheExcitationModel(), G4DNAScreenedRutherfordElasticModel::G4DNAScreenedRutherfordElasticModel(), G4DNATransformElectronModel::G4DNATransformElectronModel(), G4DNAUeharaScreenedRutherfordElasticModel::G4DNAUeharaScreenedRutherfordElasticModel(), G4eDPWACoulombScatteringModel::G4eDPWACoulombScatteringModel(), G4ICRU73QOModel::G4ICRU73QOModel(), G4MicroElecElasticModel::G4MicroElecElasticModel(), G4MicroElecElasticModel_new::G4MicroElecElasticModel_new(), G4PenelopeAnnihilationModel::G4PenelopeAnnihilationModel(), G4PenelopeBremsstrahlungModel::G4PenelopeBremsstrahlungModel(), G4PenelopeComptonModel::G4PenelopeComptonModel(), G4PenelopeGammaConversionModel::G4PenelopeGammaConversionModel(), G4PenelopeIonisationModel::G4PenelopeIonisationModel(), G4PenelopePhotoElectricModel::G4PenelopePhotoElectricModel(), G4PenelopeRayleighModel(), G4PenelopeRayleighModelMI::G4PenelopeRayleighModelMI(), G4TDNAOneStepThermalizationModel< MODEL >::G4TDNAOneStepThermalizationModel(), G4XrayRayleighModel::G4XrayRayleighModel(), G4VLEPTSModel::Init(), G4DNAChampionElasticModel::Initialise(), G4DNACPA100ElasticModel::Initialise(), G4DNADingfelderChargeDecreaseModel::Initialise(), G4DNADingfelderChargeIncreaseModel::Initialise(), G4DNAMeltonAttachmentModel::Initialise(), G4DNAMillerGreenExcitationModel::Initialise(), G4DNARuddIonisationExtendedModel::Initialise(), G4DNARuddIonisationModel::Initialise(), G4DNASancheExcitationModel::Initialise(), G4MicroElecElasticModel::Initialise(), G4MicroElecInelasticModel::Initialise(), G4MicroElecInelasticModel_new::Initialise(), G4DNABornExcitationModel1::Initialise(), G4DNABornExcitationModel2::Initialise(), G4DNABornIonisationModel1::Initialise(), G4DNABornIonisationModel2::Initialise(), G4DNAELSEPAElasticModel::Initialise(), G4DNAModelInterface::Initialise(), G4DNAIonElasticModel::Initialise(), G4hhIonisation::InitialiseEnergyLossProcess(), G4mplIonisation::InitialiseEnergyLossProcess(), G4ePairProduction::InitialiseEnergyLossProcess(), G4MuBremsstrahlung::InitialiseEnergyLossProcess(), G4MuIonisation::InitialiseEnergyLossProcess(), G4MuPairProduction::InitialiseEnergyLossProcess(), G4PolarizedBremsstrahlung::InitialiseEnergyLossProcess(), G4PolarizedIonisation::InitialiseEnergyLossProcess(), G4eBremsstrahlung::InitialiseEnergyLossProcess(), G4eIonisation::InitialiseEnergyLossProcess(), G4hIonisation::InitialiseEnergyLossProcess(), G4ionIonisation::InitialiseEnergyLossProcess(), G4DNAAttachment::InitialiseProcess(), G4DNAChargeDecrease::InitialiseProcess(), G4DNAChargeIncrease::InitialiseProcess(), G4DNADissociation::InitialiseProcess(), G4DNAElastic::InitialiseProcess(), G4DNAExcitation::InitialiseProcess(), G4DNAIonisation::InitialiseProcess(), G4DNAPlasmonExcitation::InitialiseProcess(), G4DNAPositronium::InitialiseProcess(), G4DNARotExcitation::InitialiseProcess(), G4DNAVibExcitation::InitialiseProcess(), G4PolarizedCompton::InitialiseProcess(), G4PolarizedGammaConversion::InitialiseProcess(), G4PolarizedPhotoElectric::InitialiseProcess(), G4ComptonScattering::InitialiseProcess(), G4CoulombScattering::InitialiseProcess(), G4eplusAnnihilation::InitialiseProcess(), G4GammaConversion::InitialiseProcess(), G4PhotoElectricEffect::InitialiseProcess(), G4VEmProcess::PreparePhysicsTable(), G4VEnergyLossProcess::PreparePhysicsTable(), G4VMultipleScattering::PreparePhysicsTable(), G4DNAUeharaScreenedRutherfordElasticModel::SelectHighEnergyLimit(), G4VEmAdjointModel::SetHighEnergyLimit(), G4DNAELSEPAElasticModel::SetMaximumEnergy(), G4mplIonisationModel::SetParticle(), G4mplIonisationWithDeltaModel::SetParticle(), and G4EmConfigurator::UpdateModelEnergyRange().

◆ SetLocked()

void G4VEmModel::SetLocked ( G4bool  val)
inlineinherited

◆ SetLowEnergyLimit()

void G4VEmModel::SetLowEnergyLimit ( G4double  val)
inlineinherited

Definition at line 774 of file G4VEmModel.hh.

775{
776 lowLimit = val;
777}

References G4VEmModel::lowLimit.

Referenced by G4EmModelActivator::ActivatePAI(), G4EmDNAPhysics_option2::ConstructProcess(), G4EmDNAPhysics_option5::ConstructProcess(), G4EmDNAPhysics_stationary::ConstructProcess(), G4EmLivermorePhysics::ConstructProcess(), G4EmLowEPPhysics::ConstructProcess(), G4EmPenelopePhysics::ConstructProcess(), G4EmStandardPhysics::ConstructProcess(), G4EmStandardPhysics_option1::ConstructProcess(), G4EmStandardPhysics_option2::ConstructProcess(), G4EmStandardPhysics_option3::ConstructProcess(), G4EmStandardPhysics_option4::ConstructProcess(), G4EmStandardPhysicsGS::ConstructProcess(), G4EmStandardPhysicsWVI::ConstructProcess(), G4DNASancheExcitationModel::ExtendLowEnergyLimit(), G4AtimaEnergyLossModel::G4AtimaEnergyLossModel(), G4BetheBlochModel::G4BetheBlochModel(), G4BetheHeitler5DModel::G4BetheHeitler5DModel(), G4DNAChampionElasticModel::G4DNAChampionElasticModel(), G4DNACPA100ElasticModel::G4DNACPA100ElasticModel(), G4DNACPA100ExcitationModel::G4DNACPA100ExcitationModel(), G4DNACPA100IonisationModel::G4DNACPA100IonisationModel(), G4DNAELSEPAElasticModel::G4DNAELSEPAElasticModel(), G4DNAEmfietzoglouExcitationModel::G4DNAEmfietzoglouExcitationModel(), G4DNAEmfietzoglouIonisationModel::G4DNAEmfietzoglouIonisationModel(), G4DNAIonElasticModel::G4DNAIonElasticModel(), G4DNAMeltonAttachmentModel::G4DNAMeltonAttachmentModel(), G4DNASancheExcitationModel::G4DNASancheExcitationModel(), G4DNAScreenedRutherfordElasticModel::G4DNAScreenedRutherfordElasticModel(), G4DNATransformElectronModel::G4DNATransformElectronModel(), G4DNAUeharaScreenedRutherfordElasticModel::G4DNAUeharaScreenedRutherfordElasticModel(), G4DummyModel::G4DummyModel(), G4eBremParametrizedModel::G4eBremParametrizedModel(), G4eBremsstrahlungRelModel::G4eBremsstrahlungRelModel(), G4eDPWACoulombScatteringModel::G4eDPWACoulombScatteringModel(), G4LivermoreBremsstrahlungModel::G4LivermoreBremsstrahlungModel(), G4MicroElecElasticModel::G4MicroElecElasticModel(), G4MicroElecElasticModel_new::G4MicroElecElasticModel_new(), G4SeltzerBergerModel::G4SeltzerBergerModel(), G4TDNAOneStepThermalizationModel< MODEL >::G4TDNAOneStepThermalizationModel(), G4VLEPTSModel::Init(), G4DNAChampionElasticModel::Initialise(), G4DNACPA100ElasticModel::Initialise(), G4DNADingfelderChargeDecreaseModel::Initialise(), G4DNADingfelderChargeIncreaseModel::Initialise(), G4DNAMeltonAttachmentModel::Initialise(), G4DNAMillerGreenExcitationModel::Initialise(), G4DNARuddIonisationExtendedModel::Initialise(), G4DNARuddIonisationModel::Initialise(), G4MicroElecElasticModel::Initialise(), G4MicroElecInelasticModel::Initialise(), G4MicroElecInelasticModel_new::Initialise(), G4BetheHeitler5DModel::Initialise(), G4DNABornExcitationModel1::Initialise(), G4DNABornExcitationModel2::Initialise(), G4DNABornIonisationModel1::Initialise(), G4DNABornIonisationModel2::Initialise(), G4DNAELSEPAElasticModel::Initialise(), G4DNAModelInterface::Initialise(), G4DNAIonElasticModel::Initialise(), G4hhIonisation::InitialiseEnergyLossProcess(), G4mplIonisation::InitialiseEnergyLossProcess(), G4ePairProduction::InitialiseEnergyLossProcess(), G4MuBremsstrahlung::InitialiseEnergyLossProcess(), G4MuIonisation::InitialiseEnergyLossProcess(), G4MuPairProduction::InitialiseEnergyLossProcess(), G4PolarizedBremsstrahlung::InitialiseEnergyLossProcess(), G4PolarizedIonisation::InitialiseEnergyLossProcess(), G4eBremsstrahlung::InitialiseEnergyLossProcess(), G4eIonisation::InitialiseEnergyLossProcess(), G4hIonisation::InitialiseEnergyLossProcess(), G4ionIonisation::InitialiseEnergyLossProcess(), G4DNAAttachment::InitialiseProcess(), G4DNAChargeDecrease::InitialiseProcess(), G4DNAChargeIncrease::InitialiseProcess(), G4DNADissociation::InitialiseProcess(), G4DNAElastic::InitialiseProcess(), G4DNAExcitation::InitialiseProcess(), G4DNAIonisation::InitialiseProcess(), G4DNAPlasmonExcitation::InitialiseProcess(), G4DNAPositronium::InitialiseProcess(), G4DNARotExcitation::InitialiseProcess(), G4DNAVibExcitation::InitialiseProcess(), G4PolarizedCompton::InitialiseProcess(), G4PolarizedGammaConversion::InitialiseProcess(), G4PolarizedPhotoElectric::InitialiseProcess(), G4ComptonScattering::InitialiseProcess(), G4CoulombScattering::InitialiseProcess(), G4eplusAnnihilation::InitialiseProcess(), G4GammaConversion::InitialiseProcess(), G4PhotoElectricEffect::InitialiseProcess(), G4VEmAdjointModel::SetLowEnergyLimit(), G4mplIonisationModel::SetParticle(), G4mplIonisationWithDeltaModel::SetParticle(), and G4EmConfigurator::UpdateModelEnergyRange().

◆ SetLPMFlag()

void G4VEmModel::SetLPMFlag ( G4bool  val)
inlineinherited

◆ SetMasterThread()

void G4VEmModel::SetMasterThread ( G4bool  val)
inlineinherited

◆ SetParticle()

void G4PenelopeRayleighModel::SetParticle ( const G4ParticleDefinition p)
private

Definition at line 1289 of file G4PenelopeRayleighModel.cc.

1290{
1291 if(!fParticle) {
1292 fParticle = p;
1293 }
1294}

References fParticle.

Referenced by G4PenelopeRayleighModel(), and Initialise().

◆ SetParticleChange()

void G4VEmModel::SetParticleChange ( G4VParticleChange p,
G4VEmFluctuationModel f = nullptr 
)
inherited

◆ SetPolarAngleLimit()

void G4VEmModel::SetPolarAngleLimit ( G4double  val)
inlineinherited

◆ SetSecondaryThreshold()

void G4VEmModel::SetSecondaryThreshold ( G4double  val)
inlineinherited

◆ SetTripletModel()

void G4VEmModel::SetTripletModel ( G4VEmModel p)
inlineinherited

Definition at line 645 of file G4VEmModel.hh.

646{
647 if(p != fTripletModel) {
648 delete fTripletModel;
649 fTripletModel = p;
650 }
651}

References G4VEmModel::fTripletModel.

Referenced by G4eplusTo2GammaOKVIModel::G4eplusTo2GammaOKVIModel().

◆ SetupForMaterial()

void G4VEmModel::SetupForMaterial ( const G4ParticleDefinition ,
const G4Material ,
G4double  kineticEnergy 
)
virtualinherited

◆ SetUseBaseMaterials()

void G4VEmModel::SetUseBaseMaterials ( G4bool  val)
inlineinherited

◆ SetVerbosityLevel()

void G4PenelopeRayleighModel::SetVerbosityLevel ( G4int  lev)
inline

Definition at line 84 of file G4PenelopeRayleighModel.hh.

84{fVerboseLevel = lev;};

References fVerboseLevel.

◆ StartTracking()

void G4VEmModel::StartTracking ( G4Track )
virtualinherited

◆ UseAngularGeneratorFlag()

G4bool G4VEmModel::UseAngularGeneratorFlag ( ) const
inlineinherited

◆ UseBaseMaterials()

G4bool G4VEmModel::UseBaseMaterials ( ) const
inlineinherited

Definition at line 760 of file G4VEmModel.hh.

761{
762 return useBaseMaterials;
763}

References G4VEmModel::useBaseMaterials.

◆ Value()

G4double G4VEmModel::Value ( const G4MaterialCutsCouple couple,
const G4ParticleDefinition p,
G4double  kineticEnergy 
)
virtualinherited

Field Documentation

◆ anglModel

G4VEmAngularDistribution* G4VEmModel::anglModel = nullptr
privateinherited

◆ basedCoupleIndex

size_t G4VEmModel::basedCoupleIndex = 0
protectedinherited

◆ currentCoupleIndex

size_t G4VEmModel::currentCoupleIndex = 0
protectedinherited

Definition at line 448 of file G4VEmModel.hh.

Referenced by G4VEmModel::SetCurrentCouple().

◆ elmSelectors

std::vector<G4EmElementSelector*>* G4VEmModel::elmSelectors = nullptr
privateinherited

◆ eMaxActive

G4double G4VEmModel::eMaxActive = DBL_MAX
privateinherited

◆ eMinActive

G4double G4VEmModel::eMinActive = 0.0
privateinherited

◆ fAtomicFormFactor

G4PhysicsFreeVector * G4PenelopeRayleighModel::fAtomicFormFactor = {nullptr}
staticprivate

◆ fCurrentCouple

const G4MaterialCutsCouple* G4VEmModel::fCurrentCouple = nullptr
privateinherited

Definition at line 416 of file G4VEmModel.hh.

Referenced by G4VEmModel::CurrentCouple(), and G4VEmModel::SetCurrentCouple().

◆ fCurrentElement

const G4Element* G4VEmModel::fCurrentElement = nullptr
privateinherited

◆ fCurrentIsotope

const G4Isotope* G4VEmModel::fCurrentIsotope = nullptr
privateinherited

◆ fElementData

G4ElementData* G4VEmModel::fElementData = nullptr
protectedinherited

◆ fEmManager

G4LossTableManager* G4VEmModel::fEmManager
privateinherited

Definition at line 420 of file G4VEmModel.hh.

Referenced by G4VEmModel::G4VEmModel(), and G4VEmModel::~G4VEmModel().

◆ fIntrinsicHighEnergyLimit

G4double G4PenelopeRayleighModel::fIntrinsicHighEnergyLimit
private

Definition at line 123 of file G4PenelopeRayleighModel.hh.

Referenced by G4PenelopeRayleighModel().

◆ fIntrinsicLowEnergyLimit

G4double G4PenelopeRayleighModel::fIntrinsicLowEnergyLimit
private

Definition at line 122 of file G4PenelopeRayleighModel.hh.

Referenced by G4PenelopeRayleighModel(), and SampleSecondaries().

◆ fIsInitialised

G4bool G4PenelopeRayleighModel::fIsInitialised
private

Definition at line 125 of file G4PenelopeRayleighModel.hh.

Referenced by Initialise().

◆ flagDeexcitation

G4bool G4VEmModel::flagDeexcitation = false
privateinherited

Definition at line 455 of file G4VEmModel.hh.

Referenced by G4VEmModel::DeexcitationFlag(), and G4VEmModel::SetDeexcitationFlag().

◆ flagForceBuildTable

G4bool G4VEmModel::flagForceBuildTable = false
privateinherited

◆ fLocalTable

G4bool G4PenelopeRayleighModel::fLocalTable
private

Definition at line 127 of file G4PenelopeRayleighModel.hh.

Referenced by SampleSecondaries(), and ~G4PenelopeRayleighModel().

◆ fLogAtomicCrossSection

G4PhysicsFreeVector * G4PenelopeRayleighModel::fLogAtomicCrossSection = {nullptr}
staticprivate

◆ fLogEnergyGridPMax

G4DataVector G4PenelopeRayleighModel::fLogEnergyGridPMax
private

Definition at line 111 of file G4PenelopeRayleighModel.hh.

Referenced by G4PenelopeRayleighModel(), and GetPMaxTable().

◆ fLogFormFactorTable

std::map<const G4Material*,G4PhysicsFreeVector*>* G4PenelopeRayleighModel::fLogFormFactorTable
private

◆ fLogQSquareGrid

G4DataVector G4PenelopeRayleighModel::fLogQSquareGrid
private

◆ flucModel

G4VEmFluctuationModel* G4VEmModel::flucModel = nullptr
privateinherited

◆ fMaxZ

const G4int G4PenelopeRayleighModel::fMaxZ =99
staticprivate

Definition at line 116 of file G4PenelopeRayleighModel.hh.

Referenced by InitialiseLocal(), and ~G4PenelopeRayleighModel().

◆ fParticle

const G4ParticleDefinition* G4PenelopeRayleighModel::fParticle
protected

Definition at line 95 of file G4PenelopeRayleighModel.hh.

Referenced by Initialise(), InitialiseLocal(), and SetParticle().

◆ fParticleChange

G4ParticleChangeForGamma* G4PenelopeRayleighModel::fParticleChange
protected

Definition at line 94 of file G4PenelopeRayleighModel.hh.

Referenced by Initialise(), and SampleSecondaries().

◆ fPMaxTable

std::map<const G4Material*,G4PhysicsFreeVector*>* G4PenelopeRayleighModel::fPMaxTable
private

◆ fSamplingTable

std::map<const G4Material*,G4PenelopeSamplingData*>* G4PenelopeRayleighModel::fSamplingTable
private

◆ fTripletModel

G4VEmModel* G4VEmModel::fTripletModel = nullptr
privateinherited

◆ fVerboseLevel

G4int G4PenelopeRayleighModel::fVerboseLevel
private

◆ highLimit

G4double G4VEmModel::highLimit
privateinherited

◆ inveplus

G4double G4VEmModel::inveplus
protectedinherited

◆ isLocked

G4bool G4VEmModel::isLocked = false
privateinherited

◆ isMaster

G4bool G4VEmModel::isMaster = true
privateinherited

◆ localElmSelectors

G4bool G4VEmModel::localElmSelectors = true
privateinherited

Definition at line 460 of file G4VEmModel.hh.

Referenced by G4VEmModel::SetElementSelectors(), and G4VEmModel::~G4VEmModel().

◆ localTable

G4bool G4VEmModel::localTable = true
privateinherited

Definition at line 459 of file G4VEmModel.hh.

Referenced by G4VEmModel::SetCrossSectionTable(), and G4VEmModel::~G4VEmModel().

◆ lossFlucFlag

G4bool G4VEmModel::lossFlucFlag = true
protectedinherited

◆ lowLimit

G4double G4VEmModel::lowLimit
privateinherited

◆ name

const G4String G4VEmModel::name
privateinherited

◆ nsec

G4int G4VEmModel::nsec = 5
privateinherited

Definition at line 444 of file G4VEmModel.hh.

Referenced by G4VEmModel::CrossSectionPerVolume(), and G4VEmModel::G4VEmModel().

◆ nSelectors

G4int G4VEmModel::nSelectors = 0
privateinherited

◆ pBaseMaterial

const G4Material* G4VEmModel::pBaseMaterial = nullptr
protectedinherited

◆ pFactor

G4double G4VEmModel::pFactor = 1.0
protectedinherited

◆ polarAngleLimit

G4double G4VEmModel::polarAngleLimit
privateinherited

Definition at line 441 of file G4VEmModel.hh.

Referenced by G4VEmModel::PolarAngleLimit(), and G4VEmModel::SetPolarAngleLimit().

◆ pParticleChange

G4VParticleChange* G4VEmModel::pParticleChange = nullptr
protectedinherited

◆ secondaryThreshold

G4double G4VEmModel::secondaryThreshold = DBL_MAX
privateinherited

◆ theDensityFactor

const std::vector<G4double>* G4VEmModel::theDensityFactor = nullptr
protectedinherited

Definition at line 428 of file G4VEmModel.hh.

Referenced by G4VEmModel::G4VEmModel().

◆ theDensityIdx

const std::vector<G4int>* G4VEmModel::theDensityIdx = nullptr
protectedinherited

Definition at line 429 of file G4VEmModel.hh.

Referenced by G4VEmModel::G4VEmModel().

◆ theLPMflag

G4bool G4VEmModel::theLPMflag = false
privateinherited

Definition at line 454 of file G4VEmModel.hh.

Referenced by G4VEmModel::LPMFlag(), and G4VEmModel::SetLPMFlag().

◆ useAngularGenerator

G4bool G4VEmModel::useAngularGenerator = false
privateinherited

◆ useBaseMaterials

G4bool G4VEmModel::useBaseMaterials = false
privateinherited

◆ xsec

std::vector<G4double> G4VEmModel::xsec
privateinherited

◆ xSectionTable

G4PhysicsTable* G4VEmModel::xSectionTable = nullptr
protectedinherited

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