Geant4 11.2.2
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 ()
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel) override
 
void InitialiseForElement (const G4ParticleDefinition *, G4int Z) override
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
void SetLowEnergyThreshold (G4double val)
 
void SetDebugVerbosity (G4int val)
 
G4JAEAElasticScatteringModeloperator= (const G4JAEAElasticScatteringModel &right)=delete
 
 G4JAEAElasticScatteringModel (const G4JAEAElasticScatteringModel &)=delete
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
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 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 *, const G4double &length, G4double &eloss)
 
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 FillNumberOfSecondaries (G4int &numberOfTriplets, G4int &numberOfRecoil)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
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)
 
const G4ElementGetCurrentElement (const G4Material *mat=nullptr) const
 
G4int SelectRandomAtomNumber (const G4Material *) const
 
const G4IsotopeGetCurrentIsotope (const G4Element *elm=nullptr) const
 
G4int SelectIsotopeNumber (const G4Element *) const
 
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 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 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 *)
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
void SetLPMFlag (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 = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4PhysicsTablexSectionTable = nullptr
 
const G4MaterialpBaseMaterial = nullptr
 
const std::vector< G4double > * theDensityFactor = nullptr
 
const std::vector< G4int > * theDensityIdx = nullptr
 
G4double inveplus
 
G4double pFactor = 1.0
 
std::size_t currentCoupleIndex = 0
 
std::size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 

Detailed Description

Definition at line 51 of file G4JAEAElasticScatteringModel.hh.

Constructor & Destructor Documentation

◆ G4JAEAElasticScatteringModel() [1/2]

G4JAEAElasticScatteringModel::G4JAEAElasticScatteringModel ( )
explicit

Definition at line 59 of file G4JAEAElasticScatteringModel.cc.

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

◆ ~G4JAEAElasticScatteringModel()

G4JAEAElasticScatteringModel::~G4JAEAElasticScatteringModel ( )
virtual

Definition at line 80 of file G4JAEAElasticScatteringModel.cc.

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

◆ G4JAEAElasticScatteringModel() [2/2]

G4JAEAElasticScatteringModel::G4JAEAElasticScatteringModel ( const G4JAEAElasticScatteringModel & )
delete

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 229 of file G4JAEAElasticScatteringModel.cc.

234{
235 if (verboseLevel > 2)
236 {
237 G4cout << "G4JAEAElasticScatteringModel::ComputeCrossSectionPerAtom()"
238 << G4endl;
239 }
240
241 if(GammaEnergy < lowEnergyLimit) { return 0.0; }
242
243 G4double xs = 0.0;
244
245 G4int intZ = G4lrint(Z);
246
247 if(intZ < 1 || intZ > maxZ) { return xs; }
248
249 G4PhysicsFreeVector* pv = dataCS[intZ];
250
251 // if element was not initialised
252 // do initialisation safely for MT mode
253 if(!pv) {
254 InitialiseForElement(0, intZ);
255 pv = dataCS[intZ];
256 if(!pv) { return xs; }
257 }
258
259 G4int n = G4int(pv->GetVectorLength() - 1);
260
261 G4double e = GammaEnergy;
262 if(e >= pv->Energy(n)) {
263 xs = (*pv)[n];
264 } else if(e >= pv->Energy(0)) {
265 xs = pv->Value(e);
266 }
267
268 if(verboseLevel > 0)
269 {
270 G4cout << "****** DEBUG: tcs value for Z=" << Z << " at energy (MeV)="
271 << e << G4endl;
272 G4cout << " cs (Geant4 internal unit)=" << xs << G4endl;
273 G4cout << " -> first E*E*cs value in CS data file (iu) =" << (*pv)[0]
274 << G4endl;
275 G4cout << " -> last E*E*cs value in CS data file (iu) =" << (*pv)[n]
276 << G4endl;
277 G4cout << "*********************************************************"
278 << G4endl;
279 }
280
281 return (xs);
282}
double G4double
Definition G4Types.hh:83
void InitialiseForElement(const G4ParticleDefinition *, G4int Z) override
G4double Energy(const std::size_t index) const
G4double Value(const G4double energy, 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 )
overridevirtual

Implements G4VEmModel.

Definition at line 98 of file G4JAEAElasticScatteringModel.cc.

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

◆ InitialiseForElement()

void G4JAEAElasticScatteringModel::InitialiseForElement ( const G4ParticleDefinition * ,
G4int Z )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 385 of file G4JAEAElasticScatteringModel.cc.

387{
388 G4AutoLock l(&G4JAEAElasticScatteringModelMutex);
389 if(!dataCS[Z]) { ReadData(Z); }
390 l.unlock();
391}

Referenced by ComputeCrossSectionPerAtom().

◆ InitialiseLocal()

void G4JAEAElasticScatteringModel::InitialiseLocal ( const G4ParticleDefinition * ,
G4VEmModel * masterModel )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 146 of file G4JAEAElasticScatteringModel.cc.

148{
150}
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
std::vector< G4EmElementSelector * > * GetElementSelectors()

◆ operator=()

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 286 of file G4JAEAElasticScatteringModel.cc.

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

◆ SetDebugVerbosity()

void G4JAEAElasticScatteringModel::SetDebugVerbosity ( G4int val)
inline

Definition at line 79 of file G4JAEAElasticScatteringModel.hh.

79{verboseLevel = val;};

◆ SetLowEnergyThreshold()

void G4JAEAElasticScatteringModel::SetLowEnergyThreshold ( G4double val)
inline

Definition at line 78 of file G4JAEAElasticScatteringModel.hh.

78{lowEnergyLimit = val;};

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