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

Protected Attributes

G4ParticleChangeForGammafParticleChange
 
- Protected Attributes inherited from G4VEmModel
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 

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 *)
 

Detailed Description

Definition at line 43 of file G4LivermorePolarizedRayleighModel.hh.

Constructor & Destructor Documentation

◆ G4LivermorePolarizedRayleighModel()

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

Definition at line 53 of file G4LivermorePolarizedRayleighModel.cc.

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

◆ ~G4LivermorePolarizedRayleighModel()

G4LivermorePolarizedRayleighModel::~G4LivermorePolarizedRayleighModel ( )
virtual

Definition at line 83 of file G4LivermorePolarizedRayleighModel.cc.

84{
85 if (crossSectionHandler) delete crossSectionHandler;
86 if (formFactorData) delete formFactorData;
87}

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 145 of file G4LivermorePolarizedRayleighModel.cc.

150{
151 if (verboseLevel > 3)
152 G4cout << "Calling CrossSectionPerAtom() of G4LivermorePolarizedRayleighModel" << G4endl;
153
154 if (GammaEnergy < lowEnergyLimit || GammaEnergy > highEnergyLimit) return 0.0;
155
156 G4double cs = crossSectionHandler->FindValue(G4int(Z), GammaEnergy);
157 return cs;
158}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
G4double FindValue(G4int Z, G4double e) const

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 91 of file G4LivermorePolarizedRayleighModel.cc.

93{
94// Rayleigh process: The Quantum Theory of Radiation
95// W. Heitler, Oxford at the Clarendon Press, Oxford (1954)
96// Scattering function: A simple model of photon transport
97// D.E. Cullen, Nucl. Instr. Meth. in Phys. Res. B 101 (1995) 499-510
98// Polarization of the outcoming photon: Beam test of a prototype detector array for the PoGO astronomical hard X-ray/soft gamma-ray polarimeter
99// T. Mizuno et al., Nucl. Instr. Meth. in Phys. Res. A 540 (2005) 158-168
100
101 if (verboseLevel > 3)
102 G4cout << "Calling G4LivermorePolarizedRayleighModel::Initialise()" << G4endl;
103
104 if (crossSectionHandler)
105 {
106 crossSectionHandler->Clear();
107 delete crossSectionHandler;
108 }
109
110 // Read data files for all materials
111
112 crossSectionHandler = new G4CrossSectionHandler;
113 crossSectionHandler->Clear();
114 G4String crossSectionFile = "rayl/re-cs-";
115 crossSectionHandler->LoadData(crossSectionFile);
116
117 G4VDataSetAlgorithm* ffInterpolation = new G4LogLogInterpolation;
118 G4String formFactorFile = "rayl/re-ff-";
119 formFactorData = new G4CompositeEMDataSet(ffInterpolation,1.,1.);
120 formFactorData->LoadData(formFactorFile);
121
122 InitialiseElementSelectors(particle,cuts);
123
124 //
125 if (verboseLevel > 2)
126 G4cout << "Loaded cross section files for Livermore Polarized Rayleigh model" << G4endl;
127
128 InitialiseElementSelectors(particle,cuts);
129
130 if (verboseLevel > 0) {
131 G4cout << "Livermore Polarized Rayleigh model is initialized " << G4endl
132 << "Energy range: "
133 << LowEnergyLimit() / eV << " eV - "
134 << HighEnergyLimit() / GeV << " GeV"
135 << G4endl;
136 }
137
138 if(isInitialised) return;
140 isInitialised = true;
141}
void LoadData(const G4String &dataFile)
virtual G4bool LoadData(const G4String &fileName)=0
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 G4LivermorePolarizedRayleighModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  ,
const G4MaterialCutsCouple couple,
const G4DynamicParticle aDynamicGamma,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 162 of file G4LivermorePolarizedRayleighModel.cc.

167{
168 if (verboseLevel > 3)
169 G4cout << "Calling SampleSecondaries() of G4LivermorePolarizedRayleighModel" << G4endl;
170
171 G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy();
172
173 if (photonEnergy0 <= lowEnergyLimit)
174 {
178 return ;
179 }
180
181 G4ParticleMomentum photonDirection0 = aDynamicGamma->GetMomentumDirection();
182
183 // Select randomly one element in the current material
184 // G4int Z = crossSectionHandler->SelectRandomAtom(couple,photonEnergy0);
185 const G4ParticleDefinition* particle = aDynamicGamma->GetDefinition();
186 const G4Element* elm = SelectRandomAtom(couple,particle,photonEnergy0);
187 G4int Z = (G4int)elm->GetZ();
188
189 G4double outcomingPhotonCosTheta = GenerateCosTheta(photonEnergy0, Z);
190 G4double outcomingPhotonPhi = GeneratePhi(outcomingPhotonCosTheta);
191 G4double beta=GeneratePolarizationAngle();
192
193 // incomingPhoton reference frame:
194 // z = versor parallel to the incomingPhotonDirection
195 // x = versor parallel to the incomingPhotonPolarization
196 // y = defined as z^x
197
198 // outgoingPhoton reference frame:
199 // z' = versor parallel to the outgoingPhotonDirection
200 // x' = defined as x-x*z'z' normalized
201 // y' = defined as z'^x'
202
203 G4ThreeVector z(aDynamicGamma->GetMomentumDirection().unit());
204 G4ThreeVector x(GetPhotonPolarization(*aDynamicGamma));
205 G4ThreeVector y(z.cross(x));
206
207 // z' = std::cos(phi)*std::sin(theta) x + std::sin(phi)*std::sin(theta) y + std::cos(theta) z
208 G4double xDir;
209 G4double yDir;
210 G4double zDir;
211 zDir=outcomingPhotonCosTheta;
212 xDir=std::sqrt(1-outcomingPhotonCosTheta*outcomingPhotonCosTheta);
213 yDir=xDir;
214 xDir*=std::cos(outcomingPhotonPhi);
215 yDir*=std::sin(outcomingPhotonPhi);
216
217 G4ThreeVector zPrime((xDir*x + yDir*y + zDir*z).unit());
218 G4ThreeVector xPrime(x.perpPart(zPrime).unit());
219 G4ThreeVector yPrime(zPrime.cross(xPrime));
220
221 // outgoingPhotonPolarization is directed as x' std::cos(beta) + y' std::sin(beta)
222 G4ThreeVector outcomingPhotonPolarization(xPrime*std::cos(beta) + yPrime*std::sin(beta));
223
225 fParticleChange->ProposePolarization(outcomingPhotonPolarization);
227
228}
@ fStopAndKill
Hep3Vector unit() const
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetZ() const
Definition: G4Element.hh:131
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:459
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

Member Data Documentation

◆ fParticleChange

G4ParticleChangeForGamma* G4LivermorePolarizedRayleighModel::fParticleChange
protected

Definition at line 71 of file G4LivermorePolarizedRayleighModel.hh.

Referenced by Initialise(), and SampleSecondaries().


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