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

#include <G4PolarizedMollerBhabhaModel.hh>

Inheritance diagram for G4PolarizedMollerBhabhaModel:
G4MollerBhabhaModel G4VEmModel

Public Member Functions

 G4PolarizedMollerBhabhaModel (const G4ParticleDefinition *p=0, const G4String &nam="PolarizedMollerBhabha")
 
virtual ~G4PolarizedMollerBhabhaModel ()
 
virtual G4double ComputeCrossSectionPerElectron (const G4ParticleDefinition *, G4double kinEnergy, G4double cut, G4double emax)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
void SetTargetPolarization (const G4ThreeVector &pTarget)
 
void SetBeamPolarization (const G4ThreeVector &pBeam)
 
const G4StokesVectorGetTargetPolarization ()
 
const G4StokesVectorGetBeamPolarization ()
 
const G4StokesVectorGetFinalElectronPolarization ()
 
const G4StokesVectorGetFinalPositronPolarization ()
 
- Public Member Functions inherited from G4MollerBhabhaModel
 G4MollerBhabhaModel (const G4ParticleDefinition *p=0, const G4String &nam="MollerBhabha")
 
virtual ~G4MollerBhabhaModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
 
- 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 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 G4MollerBhabhaModel
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kinEnergy)
 
void SetParticle (const G4ParticleDefinition *p)
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 
- Protected Attributes inherited from G4MollerBhabhaModel
const G4ParticleDefinitionparticle
 
G4ParticleDefinitiontheElectron
 
G4ParticleChangeForLossfParticleChange
 
G4bool isElectron
 
G4double twoln10
 
G4double lowLimit
 
- 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 59 of file G4PolarizedMollerBhabhaModel.hh.

Constructor & Destructor Documentation

G4PolarizedMollerBhabhaModel::G4PolarizedMollerBhabhaModel ( const G4ParticleDefinition p = 0,
const G4String nam = "PolarizedMollerBhabha" 
)

Definition at line 65 of file G4PolarizedMollerBhabhaModel.cc.

References G4cout, G4endl, G4MollerBhabhaModel::isElectron, and G4MollerBhabhaModel::theElectron.

67  : G4MollerBhabhaModel(p,nam)
68 {
69 
70  // G4cout<<" particle==electron "<<(p==theElectron)<<G4endl;
71  isElectron=(p==theElectron); // necessary due to wrong order in G4MollerBhabhaModel constructor!
72 
73  if (p==0) {
74 
75  }
76  if (!isElectron) {
77  G4cout<<" buildBhabha cross section "<<isElectron<<G4endl;
78  crossSectionCalculator = new G4PolarizedBhabhaCrossSection();
79  } else {
80  G4cout<<" buildMoller cross section "<<isElectron<<G4endl;
81  crossSectionCalculator = new G4PolarizedMollerCrossSection();
82  }
83 }
G4MollerBhabhaModel(const G4ParticleDefinition *p=0, const G4String &nam="MollerBhabha")
G4GLOB_DLL std::ostream G4cout
G4ParticleDefinition * theElectron
#define G4endl
Definition: G4ios.hh:61
G4PolarizedMollerBhabhaModel::~G4PolarizedMollerBhabhaModel ( )
virtual

Definition at line 87 of file G4PolarizedMollerBhabhaModel.cc.

88 {
89  if (crossSectionCalculator) {
90  delete crossSectionCalculator;
91  }
92 }

Member Function Documentation

G4double G4PolarizedMollerBhabhaModel::ComputeCrossSectionPerElectron ( const G4ParticleDefinition pd,
G4double  kinEnergy,
G4double  cut,
G4double  emax 
)
virtual

Reimplemented from G4MollerBhabhaModel.

Definition at line 97 of file G4PolarizedMollerBhabhaModel.cc.

References G4MollerBhabhaModel::ComputeCrossSectionPerElectron(), python.hepunit::electron_mass_c2, G4InuclParticleNames::gam, G4MollerBhabhaModel::MaxSecondaryEnergy(), G4INCL::Math::min(), and G4StokesVector::ZERO.

102 {
103  G4double xs =
105  cut,emax);
106 // G4cout<<"calc eIoni xsec "<<xs<<G4endl;
107 // G4cout<<" "<<kinEnergy<<" "<<cut<<" "<<emax<<G4endl;
108  G4double factor=1.;
109  if (xs!=0) {
110  // G4cout<<"calc asym"<<G4endl;
111  G4double tmax = MaxSecondaryEnergy(pd, kinEnergy);
112  tmax = std::min(emax, tmax);
113 
114  if (std::fabs(cut/emax-1.)<1.e-10) return xs;
115 
116  if(cut < tmax) {
117 
118  G4double xmin = cut/kinEnergy;
119  G4double xmax = tmax/kinEnergy;
120 // G4cout<<"calc asym "<<xmin<<","<<xmax<<G4endl;
121  G4double gam = kinEnergy/electron_mass_c2 + 1.0;
122 
123  G4double crossPol=crossSectionCalculator->
124  TotalXSection(xmin,xmax,gam,
125  theBeamPolarization,
126  theTargetPolarization);
127  G4double crossUnpol=crossSectionCalculator->
128  TotalXSection(xmin,xmax,gam,
131  if (crossUnpol>0.) factor=crossPol/crossUnpol;
132  // G4cout<<" factor="<<factor<<G4endl;
133  }
134  }
135  return xs*factor;
136 }
float electron_mass_c2
Definition: hepunit.py:274
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
double G4double
Definition: G4Types.hh:76
static const G4StokesVector ZERO
const G4StokesVector& G4PolarizedMollerBhabhaModel::GetBeamPolarization ( )
inline

Definition at line 96 of file G4PolarizedMollerBhabhaModel.hh.

97  {
98  return theBeamPolarization;
99  }
const G4StokesVector& G4PolarizedMollerBhabhaModel::GetFinalElectronPolarization ( )
inline

Definition at line 100 of file G4PolarizedMollerBhabhaModel.hh.

101  {
102  return fElectronPolarization;
103  }
const G4StokesVector& G4PolarizedMollerBhabhaModel::GetFinalPositronPolarization ( )
inline

Definition at line 104 of file G4PolarizedMollerBhabhaModel.hh.

105  {
106  return fPositronPolarization;
107  }
const G4StokesVector& G4PolarizedMollerBhabhaModel::GetTargetPolarization ( )
inline

Definition at line 92 of file G4PolarizedMollerBhabhaModel.hh.

93  {
94  return theTargetPolarization;
95  }
void G4PolarizedMollerBhabhaModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  vdp,
const G4MaterialCutsCouple ,
const G4DynamicParticle dp,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Reimplemented from G4MollerBhabhaModel.

Definition at line 140 of file G4PolarizedMollerBhabhaModel.cc.

References python.hepunit::classic_electr_radius, DBL_MIN, python.hepunit::electron_mass_c2, energy(), G4MollerBhabhaModel::fParticleChange, G4cout, G4endl, G4UniformRand, G4InuclParticleNames::gam, G4ParticleChangeForLoss::GetCurrentTrack(), G4PolarizationHelper::GetFrame(), G4PolarizationManager::GetInstance(), G4DynamicParticle::GetKineticEnergy(), G4VPhysicalVolume::GetLogicalVolume(), G4DynamicParticle::GetMomentumDirection(), G4VPolarizedCrossSection::GetPol2(), G4VPolarizedCrossSection::GetPol3(), G4DynamicParticle::GetPolarization(), G4Track::GetVolume(), G4PolarizationManager::GetVolumePolarization(), G4VPolarizedCrossSection::Initialize(), G4StokesVector::InvRotateAz(), G4MollerBhabhaModel::isElectron, G4PolarizationManager::IsPolarized(), G4VEmModel::MaxSecondaryKinEnergy(), G4INCL::Math::min(), G4ParticleChangeForLoss::ProposePolarization(), G4StokesVector::RotateAz(), CLHEP::Hep3Vector::rotateUz(), G4DynamicParticle::SetPolarization(), G4ParticleChangeForLoss::SetProposedKineticEnergy(), G4ParticleChangeForLoss::SetProposedMomentumDirection(), sqr(), G4MollerBhabhaModel::theElectron, python.hepunit::twopi, CLHEP::Hep3Vector::unit(), test::x, CLHEP::Hep3Vector::x(), G4VPolarizedCrossSection::XSection(), CLHEP::Hep3Vector::y(), z, CLHEP::Hep3Vector::z(), and G4StokesVector::ZERO.

145 {
146  // *** obtain and save target and beam polarization ***
148 
149  const G4Track * aTrack = fParticleChange->GetCurrentTrack();
150 
151  // obtain polarization of the beam
152  theBeamPolarization = dp->GetPolarization();
153 
154  // obtain polarization of the media
155  G4VPhysicalVolume* aPVolume = aTrack->GetVolume();
156  G4LogicalVolume* aLVolume = aPVolume->GetLogicalVolume();
157  const G4bool targetIsPolarized = polarizationManger->IsPolarized(aLVolume);
158  theTargetPolarization = polarizationManger->GetVolumePolarization(aLVolume);
159 
160  // transfer target polarization in interaction frame
161  if (targetIsPolarized)
162  theTargetPolarization.rotateUz(dp->GetMomentumDirection());
163 
164 
165 
166 
167  G4double tmax = std::min(maxEnergy, MaxSecondaryKinEnergy(dp));
168  if(tmin >= tmax) return;
169  // if(tmin > tmax) tmin = tmax;
170 
171  G4double polL = theBeamPolarization.z()*theTargetPolarization.z();
172  polL=std::fabs(polL);
173  G4double polT = theBeamPolarization.x()*theTargetPolarization.x() +
174  theBeamPolarization.y()*theTargetPolarization.y();
175  polT=std::fabs(polT);
176 
177  G4double kineticEnergy = dp->GetKineticEnergy();
178  G4double energy = kineticEnergy + electron_mass_c2;
179  G4double totalMomentum = std::sqrt(kineticEnergy*(energy + electron_mass_c2));
180  G4double xmin = tmin/kineticEnergy;
181  G4double xmax = tmax/kineticEnergy;
182  G4double gam = energy/electron_mass_c2;
183  G4double gamma2 = gam*gam;
184  G4double gmo = gam - 1.;
185  G4double gmo2 = gmo*gmo;
186  G4double gmo3 = gmo2*gmo;
187  G4double gpo = gam + 1.;
188  G4double gpo2 = gpo*gpo;
189  G4double gpo3 = gpo2*gpo;
190  G4double x, y, q, grej, grej2;
191  G4double z = 0.;
192  G4double xs = 0., phi =0.;
193  G4ThreeVector direction = dp->GetMomentumDirection();
194 
195  //(Polarized) Moller (e-e-) scattering
196  if (isElectron) {
197  // *** dice according to polarized cross section
198  G4double G = ((2.0*gam - 1.0)/gamma2)*(1. - polT - polL*gam);
199  G4double H = (sqr(gam - 1.0)/gamma2)*(1. + polT + polL*((gam + 3.)/(gam - 1.)));
200 
201  y = 1.0 - xmax;
202  grej = 1.0 - G*xmax + xmax*xmax*(H + (1.0 - G*y)/(y*y));
203  grej2 = 1.0 - G*xmin + xmin*xmin*(H + (1.0 - G*y)/(y*y));
204  if (grej2 > grej) grej = grej2;
205  G4double prefM = gamma2*classic_electr_radius*classic_electr_radius/(gmo2*(gam + 1.0));
206  grej *= prefM;
207  do {
208  q = G4UniformRand();
209  x = xmin*xmax/(xmin*(1.0 - q) + xmax*q);
210  if (crossSectionCalculator) {
211  crossSectionCalculator->Initialize(x,gam,phi,theBeamPolarization,
212  theTargetPolarization,1);
213  xs=crossSectionCalculator->XSection(G4StokesVector::ZERO,
215  z=xs*sqr(x)*4.;
216  if (grej < z) {
217  G4cout<<"WARNING : error in Moller rejection routine! \n"
218  <<" z = "<<z<<" grej="<<grej<<"\n";
219  }
220  } else {
221  G4cout<<"No calculator in Moller scattering"<<G4endl;
222  }
223  } while(grej * G4UniformRand() > z);
224  //Bhabha (e+e-) scattering
225  } else {
226  // *** dice according to polarized cross section
227  y = xmax*xmax;
228  grej = 0.;
229  grej += y*y*gmo3*(1. + (polL + polT)*(gam + 3.)/gmo);
230  grej += -2.*xmin*xmin*xmin*gam*gmo2*(1. - (polL + polT)*(gam + 3.)/gmo);
231  grej += y*y*gmo*(3.*gamma2 + 6.*gam + 4.)*(1. + (polL*(3.*gam + 1.)*(gamma2 + gam + 1.) + polT*((gam + 2.)*gamma2 + 1.))/(gmo*(3.*gam*(gam + 2.) + 4.)));
232  grej /= gpo3;
233  grej += -xmin*(2.*gamma2 + 4.*gam + 1.)*(1. - gam*(polL*(2.*gam + 1.) + polT)/(2.*gam*(gam + 2.) + 1.))/gpo2;
234  grej += gamma2/(gamma2 - 1.);
236  grej *= prefB;
237 
238  do {
239  q = G4UniformRand();
240  x = xmin*xmax/(xmin*(1.0 - q) + xmax*q);
241  if (crossSectionCalculator) {
242  crossSectionCalculator->Initialize(x,gam,phi,theBeamPolarization,
243  theTargetPolarization,1);
244  xs=crossSectionCalculator->XSection(G4StokesVector::ZERO,
246  z=xs*sqr(x)*4.;
247  } else {
248  G4cout<<"No calculator in Bhabha scattering"<<G4endl;
249  }
250 
251  if(z > grej) {
252  G4cout<<"&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&"<<G4endl;
253  G4cout << "G4PolarizedMollerBhabhaModel::SampleSecondaries Warning! "<<G4endl
254  << "Majorant " << grej << " < "
255  << z << " for x= " << x<<G4endl
256  << " e+e- (Bhabha) scattering"<<" at KinEnergy "<<kineticEnergy<<G4endl;
257  G4cout<<"&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&"<<G4endl;
258  }
259  } while(grej * G4UniformRand() > z);
260  }
261  //
262  //
263  // polar asymmetries (due to transverse polarizations)
264  //
265  //
266  if (crossSectionCalculator) {
267  // grej*=1./(sqr(x)*sqr(gamma2-1))*sqr(gam*(1+gam));
268  grej=xs*2.;
269  do {
270  phi = twopi * G4UniformRand() ;
271  crossSectionCalculator->Initialize(x,gam,phi,theBeamPolarization,
272  theTargetPolarization,1);
273  xs=crossSectionCalculator->XSection(G4StokesVector::ZERO,
275  if(xs > grej) {
276  if (isElectron){
277  G4cout << "G4PolarizedMollerBhabhaModel::SampleSecondaries Warning! "<<G4endl
278  << "Majorant " << grej << " < "
279  << xs << " for phi= " << phi<<G4endl
280  << " e-e- (Moller) scattering"<< G4endl
281  <<"PHI DICING"<<G4endl;
282  } else {
283  G4cout << "G4PolarizedMollerBhabhaModel::SampleSecondaries Warning! "<<G4endl
284  << "Majorant " << grej << " < "
285  << xs << " for phi= " << phi<<G4endl
286  << " e+e- (Bhabha) scattering"<< G4endl
287  <<"PHI DICING"<<G4endl;
288  }
289  }
290  } while(grej * G4UniformRand() > xs);
291  }
292 
293  // fix kinematics of delta electron
294  G4double deltaKinEnergy = x * kineticEnergy;
295  G4double deltaMomentum =
296  std::sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
297  G4double cost = deltaKinEnergy * (energy + electron_mass_c2) /
298  (deltaMomentum * totalMomentum);
299  G4double sint = 1.0 - cost*cost;
300  if(sint > 0.0) sint = std::sqrt(sint);
301 
302 
303  G4ThreeVector deltaDirection(-sint*std::cos(phi),-sint*std::sin(phi), cost) ;
304  deltaDirection.rotateUz(direction);
305 
306  // primary change
307  kineticEnergy -= deltaKinEnergy;
309 
310  if(kineticEnergy > DBL_MIN) {
311  G4ThreeVector dir = totalMomentum*direction - deltaMomentum*deltaDirection;
312  direction = dir.unit();
314  }
315 
316  // create G4DynamicParticle object for delta ray
317  G4DynamicParticle* delta = new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
318  vdp->push_back(delta);
319 
320  // get interaction frame
321  G4ThreeVector nInteractionFrame =
322  G4PolarizationHelper::GetFrame(direction,deltaDirection);
323 
324  if (crossSectionCalculator) {
325  // calculate mean final state polarizations
326 
327  theBeamPolarization.InvRotateAz(nInteractionFrame,direction);
328  theTargetPolarization.InvRotateAz(nInteractionFrame,direction);
329  crossSectionCalculator->Initialize(x,gam,phi,theBeamPolarization,
330  theTargetPolarization,2);
331 
332  // electron/positron
333  fPositronPolarization=crossSectionCalculator->GetPol2();
334  fPositronPolarization.RotateAz(nInteractionFrame,direction);
335 
336  fParticleChange->ProposePolarization(fPositronPolarization);
337 
338  // electron
339  fElectronPolarization=crossSectionCalculator->GetPol3();
340  fElectronPolarization.RotateAz(nInteractionFrame,deltaDirection);
341  delta->SetPolarization(fElectronPolarization.x(),
342  fElectronPolarization.y(),
343  fElectronPolarization.z());
344  }
345  else {
346  fPositronPolarization=G4ThreeVector();
347  fElectronPolarization=G4ThreeVector();
348  }
349 }
void ProposePolarization(const G4ThreeVector &dir)
G4double MaxSecondaryKinEnergy(const G4DynamicParticle *dynParticle)
Definition: G4VEmModel.hh:448
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
double x() const
G4double z
Definition: TRTMaterials.hh:39
static G4PolarizationManager * GetInstance()
virtual G4double XSection(const G4StokesVector &pol2, const G4StokesVector &pol3)=0
double z() const
virtual G4StokesVector GetPol2()
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
G4ParticleChangeForLoss * fParticleChange
bool G4bool
Definition: G4Types.hh:79
const G4ThreeVector & GetMomentumDirection() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
void SetProposedKineticEnergy(G4double proposedKinEnergy)
static G4ThreeVector GetFrame(const G4ThreeVector &, const G4ThreeVector &)
float electron_mass_c2
Definition: hepunit.py:274
void SetPolarization(G4double polX, G4double polY, G4double polZ)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
void InvRotateAz(G4ThreeVector nInteractionFrame, G4ThreeVector particleDirection)
G4ParticleDefinition * theElectron
G4LogicalVolume * GetLogicalVolume() const
Hep3Vector unit() const
double y() const
bool IsPolarized(G4LogicalVolume *lVol) const
const G4ThreeVector & GetPolarization() const
#define DBL_MIN
Definition: templates.hh:75
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4VPhysicalVolume * GetVolume() const
virtual G4StokesVector GetPol3()
#define G4endl
Definition: G4ios.hh:61
const G4Track * GetCurrentTrack() const
T sqr(const T &x)
Definition: templates.hh:145
double G4double
Definition: G4Types.hh:76
static const G4StokesVector ZERO
const G4ThreeVector & GetVolumePolarization(G4LogicalVolume *lVol) const
void RotateAz(G4ThreeVector nInteractionFrame, G4ThreeVector particleDirection)
virtual void Initialize(G4double, G4double, G4double, const G4StokesVector &p0, const G4StokesVector &p1, G4int flag=0)
void G4PolarizedMollerBhabhaModel::SetBeamPolarization ( const G4ThreeVector pBeam)
inline

Definition at line 88 of file G4PolarizedMollerBhabhaModel.hh.

Referenced by G4ePolarizedIonisation::ComputeAsymmetry().

89  {
90  theBeamPolarization = pBeam;
91  }
void G4PolarizedMollerBhabhaModel::SetTargetPolarization ( const G4ThreeVector pTarget)
inline

Definition at line 84 of file G4PolarizedMollerBhabhaModel.hh.

Referenced by G4ePolarizedIonisation::ComputeAsymmetry().

85  {
86  theTargetPolarization = pTarget;
87  }

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