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

#include <G4JAEAElasticScatteringModel.hh>

+ Inheritance diagram for G4JAEAElasticScatteringModel:

Public Member Functions

 G4JAEAElasticScatteringModel ()
 
virtual ~G4JAEAElasticScatteringModel ()
 
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 val)
 
void SetDebugVerbosity (G4int val)
 
- 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 53 of file G4JAEAElasticScatteringModel.hh.

Constructor & Destructor Documentation

◆ G4JAEAElasticScatteringModel()

G4JAEAElasticScatteringModel::G4JAEAElasticScatteringModel ( )

Definition at line 60 of file G4JAEAElasticScatteringModel.cc.

61 :G4VEmModel("G4JAEAElasticScatteringModel"),isInitialised(false)
62{
63 fParticleChange = 0;
64//Low energy limit for G4JAEAElasticScatteringModel process.
65 lowEnergyLimit = 10 * keV;
66
67 verboseLevel= 0;
68 // Verbosity scale for debugging purposes:
69 // 0 = nothing
70 // 1 = calculation of cross sections, file openings...
71 // 2 = entering in methods
72
73 if(verboseLevel > 0)
74 {
75 G4cout << "G4JAEAElasticScatteringModel is constructed " << G4endl;
76 }
77
78}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout

◆ ~G4JAEAElasticScatteringModel()

G4JAEAElasticScatteringModel::~G4JAEAElasticScatteringModel ( )
virtual

Definition at line 82 of file G4JAEAElasticScatteringModel.cc.

83{
84 if(IsMaster()) {
85 for(G4int i=0; i<=maxZ; ++i) {
86 if(dataCS[i]) {
87 delete dataCS[i];
88 dataCS[i] = nullptr;
89 }
90 if (ES_Data[i]){
91 delete ES_Data[i];
92 ES_Data[i] = nullptr;
93 }
94 }
95 }
96}
int G4int
Definition: G4Types.hh:85
G4bool IsMaster() const
Definition: G4VEmModel.hh:736

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 239 of file G4JAEAElasticScatteringModel.cc.

244{
245
246 if (verboseLevel > 2)
247 {
248 G4cout << "G4JAEAElasticScatteringModel::ComputeCrossSectionPerAtom()"
249 << G4endl;
250 }
251
252 if(GammaEnergy < lowEnergyLimit) { return 0.0; }
253
254 G4double xs = 0.0;
255
256 G4int intZ = G4lrint(Z);
257
258 if(intZ < 1 || intZ > maxZ) { return xs; }
259
260 G4LPhysicsFreeVector* pv = dataCS[intZ];
261
262 // if element was not initialised
263 // do initialisation safely for MT mode
264 if(!pv) {
265 InitialiseForElement(0, intZ);
266 pv = dataCS[intZ];
267 if(!pv) { return xs; }
268 }
269
270 G4int n = pv->GetVectorLength() - 1;
271
272 G4double e = GammaEnergy;
273 if(e >= pv->Energy(n)) {
274 xs = (*pv)[n];
275 } else if(e >= pv->Energy(0)) {
276 xs = pv->Value(e);
277 }
278
279 if(verboseLevel > 0)
280 {
281 G4cout << "****** DEBUG: tcs value for Z=" << Z << " at energy (MeV)="
282 << e << G4endl;
283 G4cout << " cs (Geant4 internal unit)=" << xs << G4endl;
284 G4cout << " -> first E*E*cs value in CS data file (iu) =" << (*pv)[0]
285 << G4endl;
286 G4cout << " -> last E*E*cs value in CS data file (iu) =" << (*pv)[n]
287 << G4endl;
288 G4cout << "*********************************************************"
289 << G4endl;
290 }
291
292 return (xs);
293}
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 G4JAEAElasticScatteringModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector cuts 
)
virtual

Implements G4VEmModel.

Definition at line 100 of file G4JAEAElasticScatteringModel.cc.

102{
103 if (verboseLevel > 1)
104 {
105 G4cout << "Calling Initialise() of G4JAEAElasticScatteringModel." << G4endl
106 << "Energy range: "
107 << LowEnergyLimit() / eV << " eV - "
108 << HighEnergyLimit() / GeV << " GeV"
109 << G4endl;
110 }
111
112 if(IsMaster()) {
113
114 // Initialise element selector
115 InitialiseElementSelectors(particle, cuts);
116
117 // Access to elements
118 char* path = std::getenv("G4LEDATA");
119 G4ProductionCutsTable* theCoupleTable =
121 G4int numOfCouples = theCoupleTable->GetTableSize();
122
123 for(G4int i=0; i<numOfCouples; ++i)
124 {
125 const G4MaterialCutsCouple* couple =
126 theCoupleTable->GetMaterialCutsCouple(i);
127 const G4Material* material = couple->GetMaterial();
128 const G4ElementVector* theElementVector = material->GetElementVector();
129 G4int nelm = material->GetNumberOfElements();
130
131 for (G4int j=0; j<nelm; ++j)
132 {
133 G4int Z = G4lrint((*theElementVector)[j]->GetZ());
134 if(Z < 1) { Z = 1; }
135 else if(Z > maxZ) { Z = maxZ; }
136 if( (!dataCS[Z]) ) { ReadData(Z, path); }
137 }
138 }
139 }
140
141 if(isInitialised) { return; }
142 fParticleChange = GetParticleChangeForGamma();
143 isInitialised = true;
144
145}
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 G4JAEAElasticScatteringModel::InitialiseForElement ( const G4ParticleDefinition ,
G4int  Z 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 402 of file G4JAEAElasticScatteringModel.cc.

404{
405 G4AutoLock l(&G4JAEAElasticScatteringModelMutex);
406 if(!dataCS[Z]) { ReadData(Z); }
407 l.unlock();
408}

Referenced by ComputeCrossSectionPerAtom().

◆ InitialiseLocal()

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

Reimplemented from G4VEmModel.

Definition at line 149 of file G4JAEAElasticScatteringModel.cc.

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 297 of file G4JAEAElasticScatteringModel.cc.

302{
303 if (verboseLevel > 2) {
304 G4cout << "Calling SampleSecondaries() of G4JAEAElasticScatteringModel."
305 << G4endl;
306 }
307
308 G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy();
309
310 // Absorption of low-energy gamma
311 if (photonEnergy0 <= lowEnergyLimit)
312 {
313 fParticleChange->ProposeTrackStatus(fStopAndKill);
314 fParticleChange->SetProposedKineticEnergy(0.);
315 fParticleChange->ProposeLocalEnergyDeposit(photonEnergy0);
316 return;
317 }
318
319
320//Warning if the incoming photon has polarization
321
322G4double Xi1=0, Xi2=0, Xi3=0;
323 G4ThreeVector gammaPolarization0 = aDynamicGamma->GetPolarization();
324 Xi1=gammaPolarization0.x();
325 Xi2=gammaPolarization0.y();
326 Xi3=gammaPolarization0.z();
327
328 G4double polarization_magnitude=Xi1*Xi1+Xi2*Xi2+Xi3*Xi3;
329 if ((polarization_magnitude)>0 || (Xi1*Xi1>0) || (Xi2*Xi2>0) || (Xi3*Xi3>0))
330 {
331 G4cout<<"WARNING: G4JAEAElasticScatteringModel is only compatible with non-polarized photons."<<G4endl;
332 G4cout<<"The event is ignored."<<G4endl;
333 return;
334 }
335
336 // Select randomly one element in the current material
337 const G4ParticleDefinition* particle = aDynamicGamma->GetDefinition();
338 const G4Element* elm = SelectRandomAtom(couple,particle,photonEnergy0);
339 G4int Z = G4lrint(elm->GetZ());
340
341
342G4int energyindex=round(100*photonEnergy0)-1;
343/*
344Getting the normalized probablity distrbution function and
345normalization factor to create the probability distribution function
346*/
347G4double a1=0, a2=0, a3=0,a4=0;
348G4double normdist=0;
349for (G4int i=0;i<=180;++i)
350 {
351 a1=ES_Data[Z]->at(4*i+300+181*4*(energyindex));
352 a2=ES_Data[Z]->at(4*i+1+300+181*4*(energyindex));
353 a3=ES_Data[Z]->at(4*i+2+300+181*4*(energyindex));
354 a4=ES_Data[Z]->at(4*i+3+300+181*4*(energyindex));
355 distribution[i]=a1*a1+a2*a2+a3*a3+a4*a4;
356 normdist += distribution[i];
357 }
358
359
360//Create the cummulative distribution function (cdf)
361 for (G4int i =0;i<=180;++i)
362 pdf[i]=distribution[i]/normdist;
363 cdf[0]=0;
364 G4double cdfsum =0;
365 for (G4int i=0; i<=180;++i)
366 {
367 cdfsum=cdfsum+pdf[i];
368 cdf[i]=cdfsum;
369 }
370 //Sampling the polar angle by inverse transform uing cdf.
372 G4double *cdfptr=lower_bound(cdf,cdf+181,r);
373 G4int cdfindex = (G4int)(cdfptr-cdf-1);
374 G4double cdfinv = (r-cdf[cdfindex])/(cdf[cdfindex+1]-cdf[cdfindex]);
375 G4double theta = (cdfindex+cdfinv)/180.;
376//polar is now ready
377 theta = theta*CLHEP::pi;
378
379
380/* Alternative sampling using CLHEP functions
381 CLHEP::RandGeneral GenDistTheta(distribution,181);
382 G4double theta = CLHEP::pi*GenDistTheta.shoot();
383 theta =theta*CLHEP::pi; //polar is now ready
384*/
385
386//Azimuth is uniformally distributed
387G4double phi = CLHEP::twopi*G4UniformRand();
388
389G4ThreeVector finaldirection(sin(theta)*cos(phi), sin(theta)*sin(phi), cos(theta));
390finaldirection.rotateUz(aDynamicGamma->GetMomentumDirection());
391//Sampling the Final State
392 fParticleChange->ProposeMomentumDirection(finaldirection);
393 fParticleChange->SetProposedKineticEnergy(photonEnergy0);
394}
@ fStopAndKill
#define G4UniformRand()
Definition: Randomize.hh:52
double z() const
double x() const
double y() const
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4ThreeVector & GetPolarization() const
G4double GetZ() const
Definition: G4Element.hh:130
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:570
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

◆ SetDebugVerbosity()

void G4JAEAElasticScatteringModel::SetDebugVerbosity ( G4int  val)
inline

Definition at line 84 of file G4JAEAElasticScatteringModel.hh.

84{verboseLevel = val;};

◆ SetLowEnergyThreshold()

void G4JAEAElasticScatteringModel::SetLowEnergyThreshold ( G4double  val)
inline

Definition at line 83 of file G4JAEAElasticScatteringModel.hh.

83{lowEnergyLimit = val;};

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