Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4ePolarizedBremsstrahlungModel Class Reference

#include <G4ePolarizedBremsstrahlungModel.hh>

Inheritance diagram for G4ePolarizedBremsstrahlungModel:
Collaboration diagram for G4ePolarizedBremsstrahlungModel:

Public Member Functions

 G4ePolarizedBremsstrahlungModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="PolBrem")
 
virtual ~G4ePolarizedBremsstrahlungModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
const G4ElementSelectedAtom ()
 
- Public Member Functions inherited from G4SeltzerBergerModel
 G4SeltzerBergerModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="eBremSB")
 
virtual ~G4SeltzerBergerModel ()
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z) override
 
void SetBicubicInterpolationFlag (G4bool)
 
- Public Member Functions inherited from G4eBremsstrahlungRelModel
 G4eBremsstrahlungRelModel (const G4ParticleDefinition *p=0, const G4String &nam="eBremLPM")
 
virtual ~G4eBremsstrahlungRelModel ()
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel) override
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double tkin, G4double Z, G4double, G4double cutEnergy, G4double maxEnergy=DBL_MAX) override
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double) override
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut) override
 
void SetLPMconstant (G4double val)
 
G4double LPMconstant () const
 
void SetLowestKinEnergy (G4double)
 
G4double LowestKinEnergy () const
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
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 *, G4double &eloss, G4double &niel, G4double length)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
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)
 
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=nullptr)
 
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
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 

Protected Attributes

G4VPolarizedCrossSectioncrossSectionCalculator
 
- Protected Attributes inherited from G4eBremsstrahlungRelModel
G4NistManagernist
 
const G4ParticleDefinitionparticle
 
G4ParticleDefinitiontheGamma
 
G4ParticleChangeForLossfParticleChange
 
G4double bremFactor
 
G4double particleMass
 
G4double kinEnergy
 
G4double totalEnergy
 
G4double densityFactor
 
G4double densityCorr
 
G4int currentZ
 
G4bool isElectron
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

Additional Inherited Members

- Protected Member Functions inherited from G4SeltzerBergerModel
virtual G4double ComputeDXSectionPerAtom (G4double gammaEnergy) override
 
virtual G4String DirectoryPath () const
 
- Protected Member Functions inherited from G4eBremsstrahlungRelModel
void SetCurrentElement (G4int)
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 58 of file G4ePolarizedBremsstrahlungModel.hh.

Constructor & Destructor Documentation

G4ePolarizedBremsstrahlungModel::G4ePolarizedBremsstrahlungModel ( const G4ParticleDefinition p = nullptr,
const G4String nam = "PolBrem" 
)
explicit

Definition at line 56 of file G4ePolarizedBremsstrahlungModel.cc.

58  : G4SeltzerBergerModel(p,nam),
59  crossSectionCalculator(nullptr)
60 {
61 }
G4SeltzerBergerModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="eBremSB")
G4ePolarizedBremsstrahlungModel::~G4ePolarizedBremsstrahlungModel ( )
virtual

Definition at line 65 of file G4ePolarizedBremsstrahlungModel.cc.

66 {
68 }

Member Function Documentation

void G4ePolarizedBremsstrahlungModel::Initialise ( const G4ParticleDefinition p,
const G4DataVector d 
)
overridevirtual

Reimplemented from G4SeltzerBergerModel.

Definition at line 72 of file G4ePolarizedBremsstrahlungModel.cc.

Here is the call graph for this function:

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

Reimplemented from G4SeltzerBergerModel.

Definition at line 84 of file G4ePolarizedBremsstrahlungModel.cc.

89 {
90  G4SeltzerBergerModel::SampleSecondaries(vdp,couple,dp,tmin,maxEnergy);
91  G4int num = vdp->size();
92 
93  if(num > 0) {
94  G4double lepEnergy0 = dp->GetKineticEnergy();
95  G4double gamEnergy1 = (*vdp)[0]->GetKineticEnergy();
96  G4double sintheta = dp->GetMomentumDirection().cross((*vdp)[0]->GetMomentumDirection()).mag();
97  if (sintheta>1.) sintheta=1.;
98 
99 
100  G4StokesVector beamPol = dp->GetPolarization();
101 
102  // determine interaction plane
103  G4ThreeVector nInteractionFrame =
106 
107  // transform polarization into interaction frame
108  beamPol.InvRotateAz(nInteractionFrame,dp->GetMomentumDirection());
109 
110  // calulcate polarization transfer
111  crossSectionCalculator->SetMaterial(GetCurrentElement()->GetN(), // number of nucleons
112  GetCurrentElement()->GetZ(),
113  GetCurrentElement()->GetfCoulomb());
114  crossSectionCalculator->Initialize(lepEnergy0, gamEnergy1, sintheta,
115  beamPol, G4StokesVector::ZERO);
116 
117  // deterimine final state polarization
119  newBeamPol.RotateAz(nInteractionFrame,
122 
123  if (num!=1) G4cout<<" WARNING "<<num<<" secondaries in polarized bremsstrahlung not supported!\n";
124  for (G4int i=0; i<num; i++) {
126  photonPol.SetPhoton();
127  photonPol.RotateAz(nInteractionFrame,(*vdp)[i]->GetMomentumDirection());
128  (*vdp)[i]->SetPolarization(photonPol.p1(),
129  photonPol.p2(),
130  photonPol.p3());
131  }
132  }
133  return;
134 }
void ProposePolarization(const G4ThreeVector &dir)
G4double GetKineticEnergy() const
G4double p2() const
void SetMaterial(G4double A, G4double Z, G4double coul)
const G4ThreeVector & GetProposedMomentumDirection() const
int G4int
Definition: G4Types.hh:78
G4double p3() const
virtual G4StokesVector GetPol2()
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
static G4ThreeVector GetFrame(const G4ThreeVector &, const G4ThreeVector &)
G4double p1() const
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double cutEnergy, G4double maxEnergy) override
void InvRotateAz(G4ThreeVector nInteractionFrame, G4ThreeVector particleDirection)
const G4ThreeVector & GetPolarization() const
virtual G4StokesVector GetPol3()
G4ParticleChangeForLoss * fParticleChange
Hep3Vector cross(const Hep3Vector &) const
double G4double
Definition: G4Types.hh:76
static const G4StokesVector ZERO
double mag() const
void RotateAz(G4ThreeVector nInteractionFrame, G4ThreeVector particleDirection)
virtual void Initialize(G4double, G4double, G4double, const G4StokesVector &p0, const G4StokesVector &p1, G4int flag=0)
const G4Element * GetCurrentElement() const
Definition: G4VEmModel.hh:466

Here is the call graph for this function:

const G4Element * G4ePolarizedBremsstrahlungModel::SelectedAtom ( )
inline

Definition at line 87 of file G4ePolarizedBremsstrahlungModel.hh.

88 {
89  return GetCurrentElement();
90 }
const G4Element * GetCurrentElement() const
Definition: G4VEmModel.hh:466

Here is the call graph for this function:

Member Data Documentation

G4VPolarizedCrossSection* G4ePolarizedBremsstrahlungModel::crossSectionCalculator
protected

Definition at line 80 of file G4ePolarizedBremsstrahlungModel.hh.


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