Geant4 10.7.0
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=0, const G4String &nam="LivermorePolarizedRayleigh")
 
virtual ~G4LivermorePolarizedRayleighModel ()
 
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)
 
- 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 49 of file G4LivermorePolarizedRayleighModel.hh.

Constructor & Destructor Documentation

◆ G4LivermorePolarizedRayleighModel()

G4LivermorePolarizedRayleighModel::G4LivermorePolarizedRayleighModel ( const G4ParticleDefinition p = 0,
const G4String nam = "LivermorePolarizedRayleigh" 
)

Definition at line 58 of file G4LivermorePolarizedRayleighModel.cc.

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

◆ ~G4LivermorePolarizedRayleighModel()

G4LivermorePolarizedRayleighModel::~G4LivermorePolarizedRayleighModel ( )
virtual

Definition at line 86 of file G4LivermorePolarizedRayleighModel.cc.

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

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 226 of file G4LivermorePolarizedRayleighModel.cc.

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

Implements G4VEmModel.

Definition at line 103 of file G4LivermorePolarizedRayleighModel.cc.

105{
106// Rayleigh process: The Quantum Theory of Radiation
107// W. Heitler, Oxford at the Clarendon Press, Oxford (1954)
108// Scattering function: A simple model of photon transport
109// D.E. Cullen, Nucl. Instr. Meth. in Phys. Res. B 101 (1995) 499-510
110// Polarization of the outcoming photon: Beam test of a prototype detector array for the PoGO astronomical hard X-ray/soft gamma-ray polarimeter
111// T. Mizuno et al., Nucl. Instr. Meth. in Phys. Res. A 540 (2005) 158-168
112
113 if (verboseLevel > 3)
114 G4cout << "Calling G4LivermorePolarizedRayleighModel::Initialise()" << G4endl;
115
116
117 if(IsMaster()) {
118
119 // Form Factor
120
121 G4VDataSetAlgorithm* ffInterpolation = new G4LogLogInterpolation;
122 G4String formFactorFile = "rayl/re-ff-";
123 formFactorData = new G4CompositeEMDataSet(ffInterpolation,1.,1.);
124 formFactorData->LoadData(formFactorFile);
125
126 // Initialise element selector
127 InitialiseElementSelectors(particle, cuts);
128
129 // Access to elements
130 char* path = std::getenv("G4LEDATA");
131 G4ProductionCutsTable* theCoupleTable =
133 G4int numOfCouples = theCoupleTable->GetTableSize();
134
135 for(G4int i=0; i<numOfCouples; ++i)
136 {
137 const G4MaterialCutsCouple* couple =
138 theCoupleTable->GetMaterialCutsCouple(i);
139 const G4Material* material = couple->GetMaterial();
140 const G4ElementVector* theElementVector = material->GetElementVector();
141 G4int nelm = material->GetNumberOfElements();
142
143 for (G4int j=0; j<nelm; ++j)
144 {
145 G4int Z = G4lrint((*theElementVector)[j]->GetZ());
146 if(Z < 1) { Z = 1; }
147 else if(Z > maxZ) { Z = maxZ; }
148 if( (!dataCS[Z]) ) { ReadData(Z, path); }
149 }
150 }
151 }
152
153 if(isInitialised) { return; }
154 fParticleChange = GetParticleChangeForGamma();
155 isInitialised = true;
156
157}
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()
virtual G4bool LoadData(const G4String &fileName)=0
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:133
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:148

◆ InitialiseForElement()

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

Reimplemented from G4VEmModel.

Definition at line 493 of file G4LivermorePolarizedRayleighModel.cc.

495 {
496 G4AutoLock l(&LivermorePolarizedRayleighModelMutex);
497 // G4cout << "G4LivermoreRayleighModel::InitialiseForElement Z= "
498 // << Z << G4endl;
499 if(!dataCS[Z]) { ReadData(Z); }
500 l.unlock();
501 }

Referenced by ComputeCrossSectionPerAtom().

◆ InitialiseLocal()

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

Reimplemented from G4VEmModel.

Definition at line 162 of file G4LivermorePolarizedRayleighModel.cc.

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 283 of file G4LivermorePolarizedRayleighModel.cc.

288{
289 if (verboseLevel > 3)
290 G4cout << "Calling SampleSecondaries() of G4LivermorePolarizedRayleighModel" << G4endl;
291
292 G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy();
293
294 if (photonEnergy0 <= lowEnergyLimit)
295 {
296 fParticleChange->ProposeTrackStatus(fStopAndKill);
297 fParticleChange->SetProposedKineticEnergy(0.);
298 fParticleChange->ProposeLocalEnergyDeposit(photonEnergy0);
299 return ;
300 }
301
302 G4ParticleMomentum photonDirection0 = aDynamicGamma->GetMomentumDirection();
303
304 // Select randomly one element in the current material
305 // G4int Z = crossSectionHandler->SelectRandomAtom(couple,photonEnergy0);
306 const G4ParticleDefinition* particle = aDynamicGamma->GetDefinition();
307 const G4Element* elm = SelectRandomAtom(couple,particle,photonEnergy0);
308 G4int Z = (G4int)elm->GetZ();
309
310 G4double outcomingPhotonCosTheta = GenerateCosTheta(photonEnergy0, Z);
311 G4double outcomingPhotonPhi = GeneratePhi(outcomingPhotonCosTheta);
312 G4double beta=GeneratePolarizationAngle();
313
314 // incomingPhoton reference frame:
315 // z = versor parallel to the incomingPhotonDirection
316 // x = versor parallel to the incomingPhotonPolarization
317 // y = defined as z^x
318
319 // outgoingPhoton reference frame:
320 // z' = versor parallel to the outgoingPhotonDirection
321 // x' = defined as x-x*z'z' normalized
322 // y' = defined as z'^x'
323
324 G4ThreeVector z(aDynamicGamma->GetMomentumDirection().unit());
325 G4ThreeVector x(GetPhotonPolarization(*aDynamicGamma));
326 G4ThreeVector y(z.cross(x));
327
328 // z' = std::cos(phi)*std::sin(theta) x + std::sin(phi)*std::sin(theta) y + std::cos(theta) z
329 G4double xDir;
330 G4double yDir;
331 G4double zDir;
332 zDir=outcomingPhotonCosTheta;
333 xDir=std::sqrt(1-outcomingPhotonCosTheta*outcomingPhotonCosTheta);
334 yDir=xDir;
335 xDir*=std::cos(outcomingPhotonPhi);
336 yDir*=std::sin(outcomingPhotonPhi);
337
338 G4ThreeVector zPrime((xDir*x + yDir*y + zDir*z).unit());
339 G4ThreeVector xPrime(x.perpPart(zPrime).unit());
340 G4ThreeVector yPrime(zPrime.cross(xPrime));
341
342 // outgoingPhotonPolarization is directed as x' std::cos(beta) + y' std::sin(beta)
343 G4ThreeVector outcomingPhotonPolarization(xPrime*std::cos(beta) + yPrime*std::sin(beta));
344
345 fParticleChange->ProposeMomentumDirection(zPrime);
346 fParticleChange->ProposePolarization(outcomingPhotonPolarization);
347 fParticleChange->SetProposedKineticEnergy(photonEnergy0);
348
349}
@ fStopAndKill
Hep3Vector unit() const
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetZ() const
Definition: G4Element.hh:130
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposePolarization(const G4ThreeVector &dir)
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)

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