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

#include <G4VPolarizedXS.hh>

+ Inheritance diagram for G4VPolarizedXS:

Public Member Functions

 G4VPolarizedXS ()
 
virtual ~G4VPolarizedXS ()
 
virtual void Initialize (G4double, G4double, G4double, const G4StokesVector &p0, const G4StokesVector &p1, G4int flag=0)=0
 
virtual G4double XSection (const G4StokesVector &pol2, const G4StokesVector &pol3)=0
 
virtual G4double TotalXSection (G4double xmin, G4double xmax, G4double y, const G4StokesVector &pol0, const G4StokesVector &pol1)
 
virtual G4StokesVector GetPol2 ()
 
virtual G4StokesVector GetPol3 ()
 
G4double GetYmin ()
 
virtual G4double GetXmin (G4double y)
 
virtual G4double GetXmax (G4double y)
 
void SetMaterial (G4double A, G4double Z, G4double coul)
 
G4VPolarizedXSoperator= (const G4VPolarizedXS &right)=delete
 
 G4VPolarizedXS (const G4VPolarizedXS &)=delete
 

Protected Member Functions

void SetXmin (G4double xmin)
 
void SetXmax (G4double xmax)
 
void SetYmin (G4double ymin)
 

Protected Attributes

G4double fXmin
 
G4double fXmax
 
G4double fYmin
 
G4double fA
 
G4double fZ
 
G4double fCoul
 

Detailed Description

Definition at line 42 of file G4VPolarizedXS.hh.

Constructor & Destructor Documentation

◆ G4VPolarizedXS() [1/2]

G4VPolarizedXS::G4VPolarizedXS ( )

Definition at line 40 of file G4VPolarizedXS.cc.

41 : fXmin(0)
42 , fXmax(1.)
43 , fYmin(1.)
44 , fA(1.)
45 , fZ(1.)
46 , fCoul(0.)
47{}

◆ ~G4VPolarizedXS()

G4VPolarizedXS::~G4VPolarizedXS ( )
virtual

Definition at line 49 of file G4VPolarizedXS.cc.

49{}

◆ G4VPolarizedXS() [2/2]

G4VPolarizedXS::G4VPolarizedXS ( const G4VPolarizedXS & )
delete

Member Function Documentation

◆ GetPol2()

G4StokesVector G4VPolarizedXS::GetPol2 ( )
virtual

Reimplemented in G4PolarizedAnnihilationXS, G4PolarizedBremsstrahlungXS, G4PolarizedComptonXS, G4PolarizedGammaConversionXS, G4PolarizedIonisationBhabhaXS, G4PolarizedIonisationMollerXS, and G4PolarizedPhotoElectricXS.

Definition at line 56 of file G4VPolarizedXS.cc.

57{
58 // neglects correlation effects!
59 G4double invXsecTotal =
65 invXsecTotal * xsPol1, invXsecTotal * xsPol2, invXsecTotal * xsPol3));
66}
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
static const G4StokesVector P3
static const G4StokesVector ZERO
static const G4StokesVector P2
static const G4StokesVector P1
virtual G4double XSection(const G4StokesVector &pol2, const G4StokesVector &pol3)=0

Referenced by G4PolarizedBremsstrahlungModel::SampleSecondaries(), G4PolarizedComptonModel::SampleSecondaries(), and G4PolarizedIonisationModel::SampleSecondaries().

◆ GetPol3()

◆ GetXmax()

G4double G4VPolarizedXS::GetXmax ( G4double y)
virtual

Reimplemented in G4PolarizedAnnihilationXS.

Definition at line 84 of file G4VPolarizedXS.cc.

84{ return fXmax; }

◆ GetXmin()

G4double G4VPolarizedXS::GetXmin ( G4double y)
virtual

Reimplemented in G4PolarizedAnnihilationXS.

Definition at line 81 of file G4VPolarizedXS.cc.

81{ return fXmin; }

◆ GetYmin()

G4double G4VPolarizedXS::GetYmin ( )
inline

Definition at line 65 of file G4VPolarizedXS.hh.

65{ return fYmin; }

◆ Initialize()

◆ operator=()

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

◆ SetMaterial()

void G4VPolarizedXS::SetMaterial ( G4double A,
G4double Z,
G4double coul )
inline

◆ SetXmax()

void G4VPolarizedXS::SetXmax ( G4double xmax)
inlineprotected

◆ SetXmin()

void G4VPolarizedXS::SetXmin ( G4double xmin)
inlineprotected

Definition at line 83 of file G4VPolarizedXS.hh.

83{ fXmin = xmin; }

◆ SetYmin()

void G4VPolarizedXS::SetYmin ( G4double ymin)
inlineprotected

Definition at line 85 of file G4VPolarizedXS.hh.

85{ fYmin = ymin; }

Referenced by G4PolarizedComptonXS::G4PolarizedComptonXS().

◆ TotalXSection()

G4double G4VPolarizedXS::TotalXSection ( G4double xmin,
G4double xmax,
G4double y,
const G4StokesVector & pol0,
const G4StokesVector & pol1 )
virtual

Reimplemented in G4PolarizedAnnihilationXS, G4PolarizedComptonXS, G4PolarizedIonisationBhabhaXS, and G4PolarizedIonisationMollerXS.

Definition at line 86 of file G4VPolarizedXS.cc.

89{
91 ed << "WARNING virtual function G4VPolarizedXS::TotalXSection() "
92 "called.\n";
93 G4Exception("G4VPolarizedXS::TotalXSection", "pol032", FatalException, ed);
94 return 0.;
95}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription

Referenced by G4PolarizedIonisationModel::ComputeCrossSectionPerElectron().

◆ XSection()

Member Data Documentation

◆ fA

G4double G4VPolarizedXS::fA
protected

Definition at line 90 of file G4VPolarizedXS.hh.

Referenced by SetMaterial().

◆ fCoul

G4double G4VPolarizedXS::fCoul
protected

◆ fXmax

G4double G4VPolarizedXS::fXmax
protected

Definition at line 88 of file G4VPolarizedXS.hh.

Referenced by GetXmax(), and SetXmax().

◆ fXmin

G4double G4VPolarizedXS::fXmin
protected

Definition at line 88 of file G4VPolarizedXS.hh.

Referenced by GetXmin(), and SetXmin().

◆ fYmin

G4double G4VPolarizedXS::fYmin
protected

Definition at line 88 of file G4VPolarizedXS.hh.

Referenced by GetYmin(), and SetYmin().

◆ fZ


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