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

#include <G4LivermoreRayleighModel.hh>

+ Inheritance diagram for G4LivermoreRayleighModel:

Public Member Functions

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

Additional Inherited Members

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

Detailed Description

Definition at line 39 of file G4LivermoreRayleighModel.hh.

Constructor & Destructor Documentation

◆ G4LivermoreRayleighModel()

G4LivermoreRayleighModel::G4LivermoreRayleighModel ( )

Definition at line 43 of file G4LivermoreRayleighModel.cc.

44 :G4VEmModel("LivermoreRayleigh"),isInitialised(false)
45{
46 fParticleChange = 0;
47 lowEnergyLimit = 10 * eV;
48
50
51 verboseLevel= 0;
52 // Verbosity scale for debugging purposes:
53 // 0 = nothing
54 // 1 = calculation of cross sections, file openings...
55 // 2 = entering in methods
56
57 if(verboseLevel > 0)
58 {
59 G4cout << "G4LivermoreRayleighModel is constructed " << G4endl;
60 }
61
62}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:618

◆ ~G4LivermoreRayleighModel()

G4LivermoreRayleighModel::~G4LivermoreRayleighModel ( )
virtual

Definition at line 66 of file G4LivermoreRayleighModel.cc.

67{
68 if(IsMaster()) {
69 for(G4int i=0; i<maxZ; ++i) {
70 if(dataCS[i]) {
71 delete dataCS[i];
72 dataCS[i] = 0;
73 }
74 }
75 }
76}
int G4int
Definition: G4Types.hh:85
G4bool IsMaster() const
Definition: G4VEmModel.hh:736

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 193 of file G4LivermoreRayleighModel.cc.

198{
199 if (verboseLevel > 1)
200 {
201 G4cout << "G4LivermoreRayleighModel::ComputeCrossSectionPerAtom()"
202 << G4endl;
203 }
204
205 if(GammaEnergy < lowEnergyLimit) { return 0.0; }
206
207 G4double xs = 0.0;
208
209 G4int intZ = G4lrint(Z);
210
211 if(intZ < 1 || intZ > maxZ) { return xs; }
212
213 G4LPhysicsFreeVector* pv = dataCS[intZ];
214
215 // if element was not initialised
216 // do initialisation safely for MT mode
217 if(!pv) {
218 InitialiseForElement(0, intZ);
219 pv = dataCS[intZ];
220 if(!pv) { return xs; }
221 }
222
223 G4int n = pv->GetVectorLength() - 1;
224 G4double e = GammaEnergy/MeV;
225 if(e >= pv->Energy(n)) {
226 xs = (*pv)[n]/(e*e);
227 } else if(e >= pv->Energy(0)) {
228 xs = pv->Value(e)/(e*e);
229 }
230
231 if(verboseLevel > 0)
232 {
233 G4cout << "****** DEBUG: tcs value for Z=" << Z << " at energy (MeV)="
234 << e << G4endl;
235 G4cout << " cs (Geant4 internal unit)=" << xs << G4endl;
236 G4cout << " -> first E*E*cs value in CS data file (iu) =" << (*pv)[0]
237 << G4endl;
238 G4cout << " -> last E*E*cs value in CS data file (iu) =" << (*pv)[n]
239 << G4endl;
240 G4cout << "*********************************************************"
241 << G4endl;
242 }
243 return xs;
244}
double G4double
Definition: G4Types.hh:83
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
G4double Energy(std::size_t index) const
G4double Value(G4double theEnergy, std::size_t &lastidx) const
std::size_t GetVectorLength() const
int G4lrint(double ad)
Definition: templates.hh:134

◆ Initialise()

void G4LivermoreRayleighModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector cuts 
)
virtual

Implements G4VEmModel.

Definition at line 80 of file G4LivermoreRayleighModel.cc.

82{
83 if (verboseLevel > 1)
84 {
85 G4cout << "Calling Initialise() of G4LivermoreRayleighModel." << G4endl
86 << "Energy range: "
87 << LowEnergyLimit() / eV << " eV - "
88 << HighEnergyLimit() / GeV << " GeV"
89 << G4endl;
90 }
91
92 if(IsMaster()) {
93
94 // Initialise element selector
95 InitialiseElementSelectors(particle, cuts);
96
97 // Access to elements
98 char* path = std::getenv("G4LEDATA");
99 G4ProductionCutsTable* theCoupleTable =
101 G4int numOfCouples = theCoupleTable->GetTableSize();
102
103 for(G4int i=0; i<numOfCouples; ++i)
104 {
105 const G4MaterialCutsCouple* couple =
106 theCoupleTable->GetMaterialCutsCouple(i);
107 const G4Material* material = couple->GetMaterial();
108 const G4ElementVector* theElementVector = material->GetElementVector();
109 G4int nelm = material->GetNumberOfElements();
110
111 for (G4int j=0; j<nelm; ++j)
112 {
113 G4int Z = (*theElementVector)[j]->GetZasInt();
114 if(Z < 1) { Z = 1; }
115 else if(Z > maxZ) { Z = maxZ; }
116 if( (!dataCS[Z]) ) { ReadData(Z, path); }
117 }
118 }
119 }
120
121 if(isInitialised) { return; }
122 fParticleChange = GetParticleChangeForGamma();
123 isInitialised = true;
124
125}
std::vector< G4Element * > G4ElementVector
const G4Material * GetMaterial() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:133
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:652
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:148

◆ InitialiseForElement()

void G4LivermoreRayleighModel::InitialiseForElement ( const G4ParticleDefinition ,
G4int  Z 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 290 of file G4LivermoreRayleighModel.cc.

292{
293 G4AutoLock l(&LivermoreRayleighModelMutex);
294 // G4cout << "G4LivermoreRayleighModel::InitialiseForElement Z= "
295 // << Z << G4endl;
296 if(!dataCS[Z]) { ReadData(Z); }
297 l.unlock();
298}

Referenced by ComputeCrossSectionPerAtom().

◆ InitialiseLocal()

void G4LivermoreRayleighModel::InitialiseLocal ( const G4ParticleDefinition ,
G4VEmModel masterModel 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 129 of file G4LivermoreRayleighModel.cc.

131{
133}
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:842
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:834

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 248 of file G4LivermoreRayleighModel.cc.

253{
254 if (verboseLevel > 1) {
255 G4cout << "Calling SampleSecondaries() of G4LivermoreRayleighModel"
256 << G4endl;
257 }
258 G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy();
259
260 // absorption of low-energy gamma
261 /*
262 if (photonEnergy0 <= lowEnergyLimit)
263 {
264 fParticleChange->ProposeTrackStatus(fStopAndKill);
265 fParticleChange->SetProposedKineticEnergy(0.);
266 fParticleChange->ProposeLocalEnergyDeposit(photonEnergy0);
267 return ;
268 }
269 */
270 // Select randomly one element in the current material
271 const G4ParticleDefinition* particle = aDynamicGamma->GetDefinition();
272 const G4Element* elm = SelectRandomAtom(couple,particle,photonEnergy0);
273 G4int Z = G4lrint(elm->GetZ());
274
275 // Sample the angle of the scattered photon
276
277 G4ThreeVector photonDirection =
278 GetAngularDistribution()->SampleDirection(aDynamicGamma,
279 photonEnergy0,
280 Z, couple->GetMaterial());
281 fParticleChange->ProposeMomentumDirection(photonDirection);
282}
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetZ() const
Definition: G4Element.hh:130
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:611
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:570

◆ SetLowEnergyThreshold()

void G4LivermoreRayleighModel::SetLowEnergyThreshold ( G4double  val)
inline

Definition at line 90 of file G4LivermoreRayleighModel.hh.

91{
92 lowEnergyLimit = val;
93}

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