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

#include <G4LivermorePolarizedRayleighModel.hh>

+ Inheritance diagram for G4LivermorePolarizedRayleighModel:

Public Member Functions

 G4LivermorePolarizedRayleighModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="LivermorePolarizedRayleigh")
 
virtual ~G4LivermorePolarizedRayleighModel ()
 
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
 
G4LivermorePolarizedRayleighModeloperator= (const G4LivermorePolarizedRayleighModel &right)=delete
 
 G4LivermorePolarizedRayleighModel (const G4LivermorePolarizedRayleighModel &)=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 40 of file G4LivermorePolarizedRayleighModel.hh.

Constructor & Destructor Documentation

◆ G4LivermorePolarizedRayleighModel() [1/2]

G4LivermorePolarizedRayleighModel::G4LivermorePolarizedRayleighModel ( const G4ParticleDefinition * p = nullptr,
const G4String & nam = "LivermorePolarizedRayleigh" )
explicit

Definition at line 59 of file G4LivermorePolarizedRayleighModel.cc.

61 :G4VEmModel(nam),fParticleChange(nullptr),isInitialised(false)
62{
63 lowEnergyLimit = 250 * CLHEP::eV;
64 //
65 verboseLevel= 0;
66 // Verbosity scale:
67 // 0 = nothing
68 // 1 = warning for energy non-conservation
69 // 2 = details of energy budget
70 // 3 = calculation of cross sections, file openings, sampling of atoms
71 // 4 = entering in methods
72
73 if(verboseLevel > 0) {
74 G4cout << "Livermore Polarized Rayleigh is constructed " << G4endl
75 << "Energy range: " << LowEnergyLimit() / CLHEP::eV << " eV - "
76 << HighEnergyLimit() / CLHEP::GeV << " GeV"
77 << G4endl;
78 }
79}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4double LowEnergyLimit() const
G4double HighEnergyLimit() const
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67

◆ ~G4LivermorePolarizedRayleighModel()

G4LivermorePolarizedRayleighModel::~G4LivermorePolarizedRayleighModel ( )
virtual

Definition at line 83 of file G4LivermorePolarizedRayleighModel.cc.

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

◆ G4LivermorePolarizedRayleighModel() [2/2]

G4LivermorePolarizedRayleighModel::G4LivermorePolarizedRayleighModel ( const G4LivermorePolarizedRayleighModel & )
delete

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 210 of file G4LivermorePolarizedRayleighModel.cc.

215{
216 if (verboseLevel > 1) {
217 G4cout << "G4LivermoreRayleighModel::ComputeCrossSectionPerAtom()"
218 << G4endl;
219 }
220
221 if(GammaEnergy < lowEnergyLimit) { return 0.0; }
222
223 G4double xs = 0.0;
224
225 G4int intZ = G4lrint(Z);
226 if(intZ < 1 || intZ > maxZ) { return xs; }
227
228 G4PhysicsFreeVector* pv = dataCS[intZ];
229
230 // if element was not initialised
231 // do initialisation safely for MT mode
232 if(nullptr == pv) {
233 InitialiseForElement(0, intZ);
234 pv = dataCS[intZ];
235 if(nullptr == pv) { return xs; }
236 }
237
238 G4int n = G4int(pv->GetVectorLength() - 1);
239 G4double e = GammaEnergy/MeV;
240 if(e >= pv->Energy(n)) {
241 xs = (*pv)[n]/(e*e);
242 } else if(e >= pv->Energy(0)) {
243 xs = pv->Value(e)/(e*e);
244 }
245 return xs;
246}
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 G4LivermorePolarizedRayleighModel::Initialise ( const G4ParticleDefinition * particle,
const G4DataVector & cuts )
overridevirtual

Implements G4VEmModel.

Definition at line 99 of file G4LivermorePolarizedRayleighModel.cc.

101{
102 // Rayleigh process: The Quantum Theory of Radiation
103 // W. Heitler, Oxford at the Clarendon Press, Oxford (1954)
104 // Scattering function: A simple model of photon transport
105 // D.E. Cullen, Nucl. Instr. Meth. in Phys. Res. B 101 (1995) 499-510
106 // Polarization of the outcoming photon: Beam test of a prototype detector array for the PoGO astronomical hard
107 // X-ray/soft gamma-ray polarimeter
108 // T. Mizuno et al., Nucl. Instr. Meth. in Phys. Res. A 540 (2005) 158-168
109 if (verboseLevel > 3)
110 G4cout << "Calling G4LivermorePolarizedRayleighModel::Initialise()" << G4endl;
111
112 if(IsMaster()) {
113
114 // Initialise element selector
115 InitialiseElementSelectors(particle, cuts);
116
117 // Access to elements
118 const char* path = G4FindDataDir("G4LEDATA");
119 auto elmTable = G4Element::GetElementTable();
120 for (auto const & elm : *elmTable) {
121 G4int Z = std::min(elm->GetZasInt(), maxZ);
122 if( nullptr == dataCS[Z] ) { ReadData(Z, path); }
123 }
124 }
125
126 if(isInitialised) { return; }
127 fParticleChange = GetParticleChangeForGamma();
128 isInitialised = true;
129}
const char * G4FindDataDir(const char *)
static G4ElementTable * GetElementTable()
Definition G4Element.cc:389
G4ParticleChangeForGamma * GetParticleChangeForGamma()
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)

◆ InitialiseForElement()

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

Reimplemented from G4VEmModel.

Definition at line 450 of file G4LivermorePolarizedRayleighModel.cc.

452{
453 G4AutoLock l(&LivermorePolarizedRayleighModelMutex);
454 if(nullptr == dataCS[Z]) { ReadData(Z); }
455 l.unlock();
456}

Referenced by ComputeCrossSectionPerAtom().

◆ InitialiseLocal()

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

Reimplemented from G4VEmModel.

Definition at line 134 of file G4LivermorePolarizedRayleighModel.cc.

136{
138}
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
std::vector< G4EmElementSelector * > * GetElementSelectors()

◆ operator=()

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 250 of file G4LivermorePolarizedRayleighModel.cc.

255{
256 if (verboseLevel > 3)
257 G4cout << "Calling SampleSecondaries() of G4LivermorePolarizedRayleighModel" << G4endl;
258
259 G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy();
260
261 if (photonEnergy0 <= lowEnergyLimit)
262 {
263 fParticleChange->ProposeTrackStatus(fStopAndKill);
264 fParticleChange->SetProposedKineticEnergy(0.);
265 fParticleChange->ProposeLocalEnergyDeposit(photonEnergy0);
266 return;
267 }
268
269 G4ParticleMomentum photonDirection0 = aDynamicGamma->GetMomentumDirection();
270
271 // Select randomly one element in the current material
272 const G4ParticleDefinition* particle = aDynamicGamma->GetDefinition();
273 const G4Element* elm = SelectRandomAtom(couple,particle,photonEnergy0);
274 G4int Z = elm->GetZasInt();
275
276 G4double outcomingPhotonCosTheta = GenerateCosTheta(photonEnergy0, Z);
277 G4double outcomingPhotonPhi = GeneratePhi(outcomingPhotonCosTheta);
278 G4double beta = GeneratePolarizationAngle();
279
280 // incomingPhoton reference frame:
281 // z = versor parallel to the incomingPhotonDirection
282 // x = versor parallel to the incomingPhotonPolarization
283 // y = defined as z^x
284
285 // outgoingPhoton reference frame:
286 // z' = versor parallel to the outgoingPhotonDirection
287 // x' = defined as x-x*z'z' normalized
288 // y' = defined as z'^x'
289 G4ThreeVector z(aDynamicGamma->GetMomentumDirection().unit());
290 G4ThreeVector x(GetPhotonPolarization(*aDynamicGamma));
291 G4ThreeVector y(z.cross(x));
292
293 // z' = std::cos(phi)*std::sin(theta)
294 // x + std::sin(phi)*std::sin(theta) y + std::cos(theta) z
295 G4double xDir;
296 G4double yDir;
297 G4double zDir;
298 zDir=outcomingPhotonCosTheta;
299 xDir=std::sqrt(1-outcomingPhotonCosTheta*outcomingPhotonCosTheta);
300 yDir=xDir;
301 xDir*=std::cos(outcomingPhotonPhi);
302 yDir*=std::sin(outcomingPhotonPhi);
303
304 G4ThreeVector zPrime((xDir*x + yDir*y + zDir*z).unit());
305 G4ThreeVector xPrime(x.perpPart(zPrime).unit());
306 G4ThreeVector yPrime(zPrime.cross(xPrime));
307
308 // outgoingPhotonPolarization is directed as
309 // x' std::cos(beta) + y' std::sin(beta)
310 G4ThreeVector outcomingPhotonPolarization(xPrime*std::cos(beta) + yPrime*std::sin(beta));
311
312 fParticleChange->ProposeMomentumDirection(zPrime);
313 fParticleChange->ProposePolarization(outcomingPhotonPolarization);
314 fParticleChange->SetProposedKineticEnergy(photonEnergy0);
315}
@ fStopAndKill
Hep3Vector unit() const
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4int GetZasInt() const
Definition G4Element.hh:120
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposePolarization(const G4ThreeVector &dir)
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)

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