Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4PolarizedPEEffectModel Class Reference

#include <G4PolarizedPEEffectModel.hh>

Inheritance diagram for G4PolarizedPEEffectModel:
G4PEEffectFluoModel G4VEmModel

Public Member Functions

 G4PolarizedPEEffectModel (const G4ParticleDefinition *p=0, const G4String &nam="Polarized-PhotoElectric")
 
void Initialise (const G4ParticleDefinition *pd, const G4DataVector &dv)
 
virtual ~G4PolarizedPEEffectModel ()
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
- Public Member Functions inherited from G4PEEffectFluoModel
 G4PEEffectFluoModel (const G4String &nam="PhotoElectric")
 
virtual ~G4PEEffectFluoModel ()
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double, G4double)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
 
- 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 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 *)
 
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)
 
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)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
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
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
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 SetForceBuildTable (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 

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 std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

Detailed Description

Definition at line 64 of file G4PolarizedPEEffectModel.hh.

Constructor & Destructor Documentation

G4PolarizedPEEffectModel::G4PolarizedPEEffectModel ( const G4ParticleDefinition p = 0,
const G4String nam = "Polarized-PhotoElectric" 
)

Definition at line 68 of file G4PolarizedPEEffectModel.cc.

70  : G4PEEffectFluoModel(nam),crossSectionCalculator(0),verboseLevel(0)
71 {
72 }
G4PEEffectFluoModel(const G4String &nam="PhotoElectric")
G4PolarizedPEEffectModel::~G4PolarizedPEEffectModel ( )
virtual

Definition at line 76 of file G4PolarizedPEEffectModel.cc.

77 {
78  if (crossSectionCalculator) delete crossSectionCalculator;
79 }

Member Function Documentation

void G4PolarizedPEEffectModel::Initialise ( const G4ParticleDefinition pd,
const G4DataVector dv 
)
virtual

Reimplemented from G4PEEffectFluoModel.

Definition at line 83 of file G4PolarizedPEEffectModel.cc.

References G4PEEffectFluoModel::Initialise().

85 {
87  if (!crossSectionCalculator)
88  crossSectionCalculator = new G4PolarizedPEEffectCrossSection();
89 }
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
void G4PolarizedPEEffectModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  vdp,
const G4MaterialCutsCouple couple,
const G4DynamicParticle dp,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Reimplemented from G4PEEffectFluoModel.

Definition at line 93 of file G4PolarizedPEEffectModel.cc.

References CLHEP::Hep3Vector::cross(), G4cout, G4VEmModel::GetCurrentElement(), G4PolarizationHelper::GetFrame(), G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMomentumDirection(), G4PolarizedPEEffectCrossSection::GetPol2(), G4DynamicParticle::GetPolarization(), G4PolarizationHelper::GetRandomFrame(), G4PolarizedPEEffectCrossSection::Initialize(), G4StokesVector::InvRotateAz(), CLHEP::Hep3Vector::mag(), G4StokesVector::p1(), G4StokesVector::p2(), G4StokesVector::p3(), G4StokesVector::RotateAz(), G4PEEffectFluoModel::SampleSecondaries(), G4VPolarizedCrossSection::SetMaterial(), G4StokesVector::SetPhoton(), and G4StokesVector::ZERO.

98 {
99  // std::vector<G4DynamicParticle*>* vdp =
100  G4PEEffectFluoModel::SampleSecondaries(vdp,couple, dp, tmin, maxEnergy);
101 
102  if(vdp && vdp->size()>0) {
103  G4double gamEnergy0 = dp->GetKineticEnergy();
104  G4double lepEnergy1 = (*vdp)[0]->GetKineticEnergy();
105  G4double sintheta = dp->GetMomentumDirection().cross((*vdp)[0]->GetMomentumDirection()).mag();
106  if (sintheta>1.) sintheta=1.;
107 
108  G4StokesVector beamPol = dp->GetPolarization();
109  beamPol.SetPhoton();
110  // G4cout<<" beamPol "<<beamPol<<G4endl;
111 
112  // determine interaction plane
113  G4ThreeVector nInteractionFrame =
115  (*vdp)[0]->GetMomentumDirection());
116  // G4cout<<" nInteractionFrame = "<<nInteractionFrame<<G4endl;
117  if (dp->GetMomentumDirection().cross((*vdp)[0]->GetMomentumDirection()).mag()<1.e-10) {
118  // G4cout<<" nInteractionFrame not well defined "<<G4endl;
119  // G4cout<<" choosing random interaction frame "<<G4endl;
120 
122  // G4cout<<"new nInteractionFrame = "<<nInteractionFrame<<G4endl;
123  }
124  /*
125  else {
126  G4ThreeVector mom1=dp->GetMomentumDirection();
127  G4ThreeVector mom2=(*vdp)[0]->GetMomentumDirection();
128  G4cout<<" mom1 = "<<mom1<<G4endl;
129  G4cout<<" mom2 = "<<mom2<<G4endl;
130  G4ThreeVector x=mom1.cross(mom2);
131  G4cout<<" mom1 x mom2 = "<<x<<" "<<x.mag()<<G4endl;
132  G4cout<<" norm = "<<(1./x.mag()*x)<<" "<<G4endl;
133  }
134  */
135 
136  // transform polarization into interaction frame
137  beamPol.InvRotateAz(nInteractionFrame,dp->GetMomentumDirection());
138 
139  // calulcate polarization transfer
140  crossSectionCalculator->SetMaterial(GetCurrentElement()->GetN(), // number of nucleons
141  GetCurrentElement()->GetZ(),
142  GetCurrentElement()->GetfCoulomb());
143  crossSectionCalculator->Initialize(gamEnergy0, lepEnergy1, sintheta,
144  beamPol, G4StokesVector::ZERO);
145 
146  // deterimine final state polarization
147  G4StokesVector lep1Pol = crossSectionCalculator->GetPol2();
148  // G4cout<<" lepPol "<<lep1Pol<<G4endl;
149  lep1Pol.RotateAz(nInteractionFrame,(*vdp)[0]->GetMomentumDirection());
150  (*vdp)[0]->SetPolarization(lep1Pol.p1(),
151  lep1Pol.p2(),
152  lep1Pol.p3());
153 
154  // G4cout<<" lepPol "<<lep1Pol<<G4endl;
155  size_t num = vdp->size();
156  if (num!=1) G4cout<<" WARNING "<<num<<" secondaries in polarized photo electric effect not supported!\n";
157  }
158 }
static G4ThreeVector GetRandomFrame(const G4ThreeVector &)
G4double GetKineticEnergy() const
G4double p2() const
void SetMaterial(G4double A, G4double Z, G4double coul)
G4double p3() const
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
static G4ThreeVector GetFrame(const G4ThreeVector &, const G4ThreeVector &)
G4double p1() const
void InvRotateAz(G4ThreeVector nInteractionFrame, G4ThreeVector particleDirection)
const G4ThreeVector & GetPolarization() const
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
Hep3Vector cross(const Hep3Vector &) const
virtual void Initialize(G4double aGammaE, G4double aLept0E, G4double sintheta, const G4StokesVector &beamPol, const G4StokesVector &, G4int flag=0)
double G4double
Definition: G4Types.hh:76
static const G4StokesVector ZERO
double mag() const
void RotateAz(G4ThreeVector nInteractionFrame, G4ThreeVector particleDirection)
const G4Element * GetCurrentElement() const
Definition: G4VEmModel.hh:440

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