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

#include <G4PolarizedPhotoElectricModel.hh>

+ Inheritance diagram for G4PolarizedPhotoElectricModel:

Public Member Functions

 G4PolarizedPhotoElectricModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="Polarized-PhotoElectric")
 
virtual ~G4PolarizedPhotoElectricModel () override
 
void Initialise (const G4ParticleDefinition *pd, const G4DataVector &dv) override
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
G4PolarizedPhotoElectricModeloperator= (const G4PolarizedPhotoElectricModel &right)=delete
 
 G4PolarizedPhotoElectricModel (const G4PolarizedPhotoElectricModel &)=delete
 
- Public Member Functions inherited from G4PEEffectFluoModel
 G4PEEffectFluoModel (const G4String &nam="PhotoElectric")
 
 ~G4PEEffectFluoModel () override
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double, G4double) override
 
G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
 
G4PEEffectFluoModeloperator= (const G4PEEffectFluoModel &right)=delete
 
 G4PEEffectFluoModel (const G4PEEffectFluoModel &)=delete
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
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 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 47 of file G4PolarizedPhotoElectricModel.hh.

Constructor & Destructor Documentation

◆ G4PolarizedPhotoElectricModel() [1/2]

G4PolarizedPhotoElectricModel::G4PolarizedPhotoElectricModel ( const G4ParticleDefinition * p = nullptr,
const G4String & nam = "Polarized-PhotoElectric" )
explicit

Definition at line 46 of file G4PolarizedPhotoElectricModel.cc.

49 , fCrossSectionCalculator(nullptr)
50 , fVerboseLevel(0)
51{}
G4PEEffectFluoModel(const G4String &nam="PhotoElectric")

◆ ~G4PolarizedPhotoElectricModel()

G4PolarizedPhotoElectricModel::~G4PolarizedPhotoElectricModel ( )
overridevirtual

Definition at line 54 of file G4PolarizedPhotoElectricModel.cc.

55{
56 if(fCrossSectionCalculator)
57 delete fCrossSectionCalculator;
58}

◆ G4PolarizedPhotoElectricModel() [2/2]

G4PolarizedPhotoElectricModel::G4PolarizedPhotoElectricModel ( const G4PolarizedPhotoElectricModel & )
delete

Member Function Documentation

◆ Initialise()

void G4PolarizedPhotoElectricModel::Initialise ( const G4ParticleDefinition * pd,
const G4DataVector & dv )
overridevirtual

Reimplemented from G4PEEffectFluoModel.

Definition at line 61 of file G4PolarizedPhotoElectricModel.cc.

63{
65 if(!fCrossSectionCalculator)
66 fCrossSectionCalculator = new G4PolarizedPhotoElectricXS();
67}
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override

◆ operator=()

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

◆ SampleSecondaries()

void G4PolarizedPhotoElectricModel::SampleSecondaries ( std::vector< G4DynamicParticle * > * vdp,
const G4MaterialCutsCouple * couple,
const G4DynamicParticle * dp,
G4double tmin,
G4double maxEnergy )
overridevirtual

Reimplemented from G4PEEffectFluoModel.

Definition at line 70 of file G4PolarizedPhotoElectricModel.cc.

73{
74 G4PEEffectFluoModel::SampleSecondaries(vdp, couple, dp, tmin, maxEnergy);
75
76 if(fVerboseLevel >= 1)
77 {
78 G4cout << "G4PolarizedPhotoElectricModel::SampleSecondaries" << G4endl;
79 }
80
81 if(vdp && !vdp->empty())
82 {
83 G4double gamEnergy0 = dp->GetKineticEnergy();
84 G4double lepEnergy1 = (*vdp)[0]->GetKineticEnergy();
85 G4double sintheta =
86 dp->GetMomentumDirection().cross((*vdp)[0]->GetMomentumDirection()).mag();
87 if(sintheta > 1.)
88 sintheta = 1.;
89
91 beamPol.SetPhoton();
92
93 // determine interaction plane
95 dp->GetMomentumDirection(), (*vdp)[0]->GetMomentumDirection());
97 .cross((*vdp)[0]->GetMomentumDirection())
98 .mag() < 1.e-10)
99 {
100 nInteractionFrame =
102 }
103
104 // transform polarization into interaction frame
105 beamPol.InvRotateAz(nInteractionFrame, dp->GetMomentumDirection());
106
107 // calulcate polarization transfer
108 fCrossSectionCalculator->SetMaterial(
109 GetCurrentElement()->GetN(), // number of nucleons
110 GetCurrentElement()->GetZ(), GetCurrentElement()->GetfCoulomb());
111 fCrossSectionCalculator->Initialize(gamEnergy0, lepEnergy1, sintheta,
112 beamPol, G4StokesVector::ZERO);
113
114 // determine final state polarization
115 G4StokesVector lep1Pol = fCrossSectionCalculator->GetPol2();
116 lep1Pol.RotateAz(nInteractionFrame, (*vdp)[0]->GetMomentumDirection());
117 (*vdp)[0]->SetPolarization(lep1Pol.p1(), lep1Pol.p2(), lep1Pol.p3());
118
119 if(vdp->size() != 1)
120 {
122 ed << " WARNING " << vdp->size()
123 << " secondaries in polarized photo electric effect not supported!\n";
124 G4Exception("G4PolarizedPhotoElectricModel::SampleSecondaries", "pol024",
125 JustWarning, ed);
126 }
127 }
128}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
double G4double
Definition G4Types.hh:83
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
Hep3Vector cross(const Hep3Vector &) const
double mag() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
const G4ThreeVector & GetPolarization() const
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
static G4ThreeVector GetFrame(const G4ThreeVector &, const G4ThreeVector &)
static G4ThreeVector GetRandomFrame(const G4ThreeVector &)
void Initialize(G4double aGammaE, G4double aLept0E, G4double sintheta, const G4StokesVector &beamPol, const G4StokesVector &, G4int flag=0) override
G4double p3() const
G4double p1() const
static const G4StokesVector ZERO
void InvRotateAz(G4ThreeVector nInteractionFrame, G4ThreeVector particleDirection)
G4double p2() const
void RotateAz(G4ThreeVector nInteractionFrame, G4ThreeVector particleDirection)
const G4Element * GetCurrentElement(const G4Material *mat=nullptr) const
void SetMaterial(G4double A, G4double Z, G4double coul)

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