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

#include <G4PenelopeRayleighModel.hh>

+ Inheritance diagram for G4PenelopeRayleighModel:

Public Member Functions

 G4PenelopeRayleighModel (const G4ParticleDefinition *p=0, const G4String &processName="PenRayleigh")
 
virtual ~G4PenelopeRayleighModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
void SetVerbosityLevel (G4int lev)
 
G4int GetVerbosityLevel ()
 
void DumpFormFactorTable (const G4Material *)
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectIsotopeNumber (const G4Element *)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *)
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool LPMFlag () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy)
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void ForceBuildTable (G4bool val)
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 

Protected Attributes

G4ParticleChangeForGammafParticleChange
 
- Protected Attributes inherited from G4VEmModel
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Detailed Description

Definition at line 57 of file G4PenelopeRayleighModel.hh.

Constructor & Destructor Documentation

◆ G4PenelopeRayleighModel()

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

Definition at line 53 of file G4PenelopeRayleighModel.cc.

55 :G4VEmModel(nam),fParticleChange(0),isInitialised(false),logAtomicCrossSection(0),
56 atomicFormFactor(0),logFormFactorTable(0),pMaxTable(0),samplingTable(0)
57{
58 fIntrinsicLowEnergyLimit = 100.0*eV;
59 fIntrinsicHighEnergyLimit = 100.0*GeV;
60 // SetLowEnergyLimit(fIntrinsicLowEnergyLimit);
61 SetHighEnergyLimit(fIntrinsicHighEnergyLimit);
62 //
63 verboseLevel= 0;
64 // Verbosity scale:
65 // 0 = nothing
66 // 1 = warning for energy non-conservation
67 // 2 = details of energy budget
68 // 3 = calculation of cross sections, file openings, sampling of atoms
69 // 4 = entering in methods
70
71 //build the energy grid. It is the same for all materials
72 G4double logenergy = std::log(fIntrinsicLowEnergyLimit/2.);
73 G4double logmaxenergy = std::log(1.5*fIntrinsicHighEnergyLimit);
74 //finer grid below 160 keV
75 G4double logtransitionenergy = std::log(160*keV);
76 G4double logfactor1 = std::log(10.)/250.;
77 G4double logfactor2 = logfactor1*10;
78 logEnergyGridPMax.push_back(logenergy);
79 do{
80 if (logenergy < logtransitionenergy)
81 logenergy += logfactor1;
82 else
83 logenergy += logfactor2;
84 logEnergyGridPMax.push_back(logenergy);
85 }while (logenergy < logmaxenergy);
86}
double G4double
Definition: G4Types.hh:64
G4ParticleChangeForGamma * fParticleChange
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585

◆ ~G4PenelopeRayleighModel()

G4PenelopeRayleighModel::~G4PenelopeRayleighModel ( )
virtual

Definition at line 90 of file G4PenelopeRayleighModel.cc.

91{
92 std::map <G4int,G4PhysicsFreeVector*>::iterator i;
93 if (logAtomicCrossSection)
94 {
95 for (i=logAtomicCrossSection->begin();i != logAtomicCrossSection->end();i++)
96 if (i->second) delete i->second;
97 delete logAtomicCrossSection;
98 }
99
100 if (atomicFormFactor)
101 {
102 for (i=atomicFormFactor->begin();i != atomicFormFactor->end();i++)
103 if (i->second) delete i->second;
104 delete atomicFormFactor;
105 }
106
107 ClearTables();
108}

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 186 of file G4PenelopeRayleighModel.cc.

192{
193 // Cross section of Rayleigh scattering in Penelope v2008 is calculated by the EPDL97
194 // tabulation, Cuellen et al. (1997), with non-relativistic form factors from Hubbel
195 // et al. J. Phys. Chem. Ref. Data 4 (1975) 471; Erratum ibid. 6 (1977) 615.
196
197 if (verboseLevel > 3)
198 G4cout << "Calling CrossSectionPerAtom() of G4PenelopeRayleighModel" << G4endl;
199
200 G4int iZ = (G4int) Z;
201
202 //read data files
203 if (!logAtomicCrossSection->count(iZ))
204 ReadDataFile(iZ);
205 //now it should be ok
206 if (!logAtomicCrossSection->count(iZ))
207 {
208 G4Exception("G4PenelopeRayleighModel::ComputeCrossSectionPerAtom()",
209 "em2040",FatalException,"Unable to load the cross section table");
210 }
211
212 G4double cross = 0;
213
214 G4PhysicsFreeVector* atom = logAtomicCrossSection->find(iZ)->second;
215 if (!atom)
216 {
218 ed << "Unable to find Z=" << iZ << " in the atomic cross section table" << G4endl;
219 G4Exception("G4PenelopeRayleighModel::ComputeCrossSectionPerAtom()",
220 "em2041",FatalException,ed);
221 return 0;
222 }
223 G4double logene = std::log(energy);
224 G4double logXS = atom->Value(logene);
225 cross = std::exp(logXS);
226
227 if (verboseLevel > 2)
228 G4cout << "Rayleigh cross section at " << energy/keV << " keV for Z=" << Z <<
229 " = " << cross/barn << " barn" << G4endl;
230 return cross;
231}
@ FatalException
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4double Value(G4double theEnergy)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76

◆ DumpFormFactorTable()

void G4PenelopeRayleighModel::DumpFormFactorTable ( const G4Material mat)

Definition at line 1160 of file G4PenelopeRayleighModel.cc.

1161{
1162 G4cout << "*****************************************************************" << G4endl;
1163 G4cout << "G4PenelopeRayleighModel: Form Factor Table for " << mat->GetName() << G4endl;
1164 //try to use the same format as Penelope-Fortran, namely Q (/m_e*c) and F
1165 G4cout << "Q/(m_e*c) F(Q) " << G4endl;
1166 G4cout << "*****************************************************************" << G4endl;
1167 if (!logFormFactorTable->count(mat))
1168 BuildFormFactorTable(mat);
1169
1170 G4PhysicsFreeVector* theVec = logFormFactorTable->find(mat)->second;
1171 for (size_t i=0;i<theVec->GetVectorLength();i++)
1172 {
1173 G4double logQ2 = theVec->GetLowEdgeEnergy(i);
1174 G4double Q = std::exp(0.5*logQ2);
1175 G4double logF2 = (*theVec)[i];
1176 G4double F = std::exp(0.5*logF2);
1177 G4cout << Q << " " << F << G4endl;
1178 }
1179 //DONE
1180 return;
1181}
const G4String & GetName() const
Definition: G4Material.hh:177
size_t GetVectorLength() const
virtual G4double GetLowEdgeEnergy(size_t binNumber) const

◆ GetVerbosityLevel()

G4int G4PenelopeRayleighModel::GetVerbosityLevel ( )
inline

Definition at line 82 of file G4PenelopeRayleighModel.hh.

82{return verboseLevel;};

◆ Initialise()

void G4PenelopeRayleighModel::Initialise ( const G4ParticleDefinition ,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 145 of file G4PenelopeRayleighModel.cc.

147{
148 if (verboseLevel > 3)
149 G4cout << "Calling G4PenelopeRayleighModel::Initialise()" << G4endl;
150
151 //clear tables depending on materials, not the atomic ones
152 ClearTables();
153
154 //create new tables
155 //
156 // logAtomicCrossSection and atomicFormFactor are created only once,
157 // since they are never cleared
158 if (!logAtomicCrossSection)
159 logAtomicCrossSection = new std::map<G4int,G4PhysicsFreeVector*>;
160 if (!atomicFormFactor)
161 atomicFormFactor = new std::map<G4int,G4PhysicsFreeVector*>;
162
163 if (!logFormFactorTable)
164 logFormFactorTable = new std::map<const G4Material*,G4PhysicsFreeVector*>;
165 if (!pMaxTable)
166 pMaxTable = new std::map<const G4Material*,G4PhysicsFreeVector*>;
167 if (!samplingTable)
168 samplingTable = new std::map<const G4Material*,G4PenelopeSamplingData*>;
169
170
171 if (verboseLevel > 0) {
172 G4cout << "Penelope Rayleigh model v2008 is initialized " << G4endl
173 << "Energy range: "
174 << LowEnergyLimit() / keV << " keV - "
175 << HighEnergyLimit() / GeV << " GeV"
176 << G4endl;
177 }
178
179 if(isInitialised) return;
181 isInitialised = true;
182}
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:529
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 306 of file G4PenelopeRayleighModel.cc.

311{
312 // Sampling of the Rayleigh final state (namely, scattering angle of the photon)
313 // from the Penelope2008 model. The scattering angle is sampled from the atomic
314 // cross section dOmega/d(cosTheta) from Born ("Atomic Phyisics", 1969), disregarding
315 // anomalous scattering effects. The Form Factor F(Q) function which appears in the
316 // analytical cross section is retrieved via the method GetFSquared(); atomic data
317 // are tabulated for F(Q). Form factor for compounds is calculated according to
318 // the additivity rule. The sampling from the F(Q) is made via a Rational Inverse
319 // Transform with Aliasing (RITA) algorithm; RITA parameters are calculated once
320 // for each material and managed by G4PenelopeSamplingData objects.
321 // The sampling algorithm (rejection method) has efficiency 67% at low energy, and
322 // increases with energy. For E=100 keV the efficiency is 100% and 86% for
323 // hydrogen and uranium, respectively.
324
325 if (verboseLevel > 3)
326 G4cout << "Calling SamplingSecondaries() of G4PenelopeRayleighModel" << G4endl;
327
328 G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy();
329
330 if (photonEnergy0 <= fIntrinsicLowEnergyLimit)
331 {
335 return ;
336 }
337
338 G4ParticleMomentum photonDirection0 = aDynamicGamma->GetMomentumDirection();
339
340 const G4Material* theMat = couple->GetMaterial();
341
342 //1) Verify if tables are ready
343 if (!pMaxTable || !samplingTable)
344 {
345 G4Exception("G4PenelopeRayleighModel::SampleSecondaries()",
346 "em2043",FatalException,"Invalid model initialization");
347 return;
348 }
349
350 //2) retrieve or build the sampling table
351 if (!(samplingTable->count(theMat)))
352 InitializeSamplingAlgorithm(theMat);
353 G4PenelopeSamplingData* theDataTable = samplingTable->find(theMat)->second;
354
355 //3) retrieve or build the pMax data
356 if (!pMaxTable->count(theMat))
357 GetPMaxTable(theMat);
358 G4PhysicsFreeVector* thePMax = pMaxTable->find(theMat)->second;
359
360 G4double cosTheta = 1.0;
361
362 //OK, ready to go!
363 G4double qmax = 2.0*photonEnergy0/electron_mass_c2; //this is non-dimensional now
364
365 if (qmax < 1e-10) //very low momentum transfer
366 {
367 G4bool loopAgain=false;
368 do
369 {
370 loopAgain = false;
371 cosTheta = 1.0-2.0*G4UniformRand();
372 G4double G = 0.5*(1+cosTheta*cosTheta);
373 if (G4UniformRand()>G)
374 loopAgain = true;
375 }while(loopAgain);
376 }
377 else //larger momentum transfer
378 {
379 size_t nData = theDataTable->GetNumberOfStoredPoints();
380 G4double LastQ2inTheTable = theDataTable->GetX(nData-1);
381 G4double q2max = std::min(qmax*qmax,LastQ2inTheTable);
382
383 G4bool loopAgain = false;
384 G4double MaxPValue = thePMax->Value(photonEnergy0);
385 G4double xx=0;
386
387 //Sampling by rejection method. The rejection function is
388 //G = 0.5*(1+cos^2(theta))
389 //
390 do{
391 loopAgain = false;
392 G4double RandomMax = G4UniformRand()*MaxPValue;
393 xx = theDataTable->SampleValue(RandomMax);
394 //xx is a random value of q^2 in (0,q2max),sampled according to
395 //F(Q^2) via the RITA algorithm
396 if (xx > q2max)
397 loopAgain = true;
398 cosTheta = 1.0-2.0*xx/q2max;
399 G4double G = 0.5*(1+cosTheta*cosTheta);
400 if (G4UniformRand()>G)
401 loopAgain = true;
402 }while(loopAgain);
403 }
404
405 G4double sinTheta = std::sqrt(1-cosTheta*cosTheta);
406
407 // Scattered photon angles. ( Z - axis along the parent photon)
408 G4double phi = twopi * G4UniformRand() ;
409 G4double dirX = sinTheta*std::cos(phi);
410 G4double dirY = sinTheta*std::sin(phi);
411 G4double dirZ = cosTheta;
412
413 // Update G4VParticleChange for the scattered photon
414 G4ThreeVector photonDirection1(dirX, dirY, dirZ);
415
416 photonDirection1.rotateUz(photonDirection0);
417 fParticleChange->ProposeMomentumDirection(photonDirection1) ;
419
420 return;
421}
@ fStopAndKill
bool G4bool
Definition: G4Types.hh:67
#define G4UniformRand()
Definition: Randomize.hh:53
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
const G4Material * GetMaterial() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4double GetX(size_t index)
G4double SampleValue(G4double rndm)
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

◆ SetVerbosityLevel()

void G4PenelopeRayleighModel::SetVerbosityLevel ( G4int  lev)
inline

Definition at line 81 of file G4PenelopeRayleighModel.hh.

81{verboseLevel = lev;};

Member Data Documentation

◆ fParticleChange

G4ParticleChangeForGamma* G4PenelopeRayleighModel::fParticleChange
protected

Definition at line 88 of file G4PenelopeRayleighModel.hh.

Referenced by Initialise(), and SampleSecondaries().


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