Geant4 9.6.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 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 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
 

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
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 

Detailed Description

Definition at line 39 of file G4LivermoreRayleighModel.hh.

Constructor & Destructor Documentation

◆ G4LivermoreRayleighModel()

G4LivermoreRayleighModel::G4LivermoreRayleighModel ( )

Definition at line 40 of file G4LivermoreRayleighModel.cc.

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

◆ ~G4LivermoreRayleighModel()

G4LivermoreRayleighModel::~G4LivermoreRayleighModel ( )
virtual

Definition at line 65 of file G4LivermoreRayleighModel.cc.

66{
67 for(G4int i=0; i<=maxZ; ++i) { delete dataCS[i]; }
68}
int G4int
Definition: G4Types.hh:66

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 178 of file G4LivermoreRayleighModel.cc.

183{
184 if (verboseLevel > 1)
185 {
186 G4cout << "Calling ComputeCrossSectionPerAtom() of G4LivermoreRayleighModel"
187 << G4endl;
188 }
189
190 if(GammaEnergy < lowEnergyLimit) { return 0.0; }
191
192 G4double xs = 0.0;
193
194 G4int intZ = G4lrint(Z);
195
196 if(intZ < 1 || intZ > maxZ) { return xs; }
197
198 G4LPhysicsFreeVector* pv = dataCS[intZ];
199
200 // element was not initialised
201 if(!pv)
202 {
203 char* path = getenv("G4LEDATA");
204 ReadData(intZ, path);
205 pv = dataCS[intZ];
206 if(!pv) { return xs; }
207 }
208
209 G4int n = pv->GetVectorLength() - 1;
210 G4double e = GammaEnergy/MeV;
211 if(e >= pv->Energy(n)) {
212 xs = (*pv)[n]/(e*e);
213 } else if(e >= pv->Energy(0)) {
214 xs = pv->Value(e)/(e*e);
215 }
216
217 if(verboseLevel > 0)
218 {
219 G4cout << "****** DEBUG: tcs value for Z=" << Z << " at energy (MeV)=" << e << G4endl;
220 G4cout << " cs (Geant4 internal unit)=" << xs << G4endl;
221 G4cout << " -> first E*E*cs value in CS data file (iu) =" << (*pv)[0] << G4endl;
222 G4cout << " -> last E*E*cs value in CS data file (iu) =" << (*pv)[n] << G4endl;
223 G4cout << "*********************************************************" << G4endl;
224 }
225 return xs;
226}
double G4double
Definition: G4Types.hh:64
G4double Value(G4double theEnergy)
size_t GetVectorLength() const
G4double Energy(size_t index) const
int G4lrint(double ad)
Definition: templates.hh:163

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 72 of file G4LivermoreRayleighModel.cc.

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 230 of file G4LivermoreRayleighModel.cc.

235{
236 if (verboseLevel > 1) {
237 G4cout << "Calling SampleSecondaries() of G4LivermoreRayleighModel"
238 << G4endl;
239 }
240 G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy();
241
242 // absorption of low-energy gamma
243 if (photonEnergy0 <= lowEnergyLimit)
244 {
245 fParticleChange->ProposeTrackStatus(fStopAndKill);
246 fParticleChange->SetProposedKineticEnergy(0.);
247 fParticleChange->ProposeLocalEnergyDeposit(photonEnergy0);
248 return ;
249 }
250
251 // Select randomly one element in the current material
252 const G4ParticleDefinition* particle = aDynamicGamma->GetDefinition();
253 const G4Element* elm = SelectRandomAtom(couple,particle,photonEnergy0);
254 G4int Z = G4lrint(elm->GetZ());
255
256 // Sample the angle of the scattered photon
257
258 G4ThreeVector photonDirection =
259 GetAngularDistribution()->SampleDirection(aDynamicGamma,
260 photonEnergy0, Z, couple->GetMaterial());
261 fParticleChange->ProposeMomentumDirection(photonDirection);
262}
@ fStopAndKill
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetZ() const
Definition: G4Element.hh:131
void SetProposedKineticEnergy(G4double proposedKinEnergy)
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:508
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:459
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

◆ SetLowEnergyThreshold()

void G4LivermoreRayleighModel::SetLowEnergyThreshold ( G4double  val)
inline

Definition at line 85 of file G4LivermoreRayleighModel.hh.

86{
87 lowEnergyLimit = val;
88}

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