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

#include <G4PolarizedMollerBhabhaModel.hh>

+ Inheritance diagram for G4PolarizedMollerBhabhaModel:

Public Member Functions

 G4PolarizedMollerBhabhaModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="PolarizedMollerBhabha")
 
virtual ~G4PolarizedMollerBhabhaModel ()
 
virtual G4double ComputeCrossSectionPerElectron (const G4ParticleDefinition *, G4double kinEnergy, G4double cut, G4double emax) override
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
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=nullptr, const G4String &nam="MollerBhabha")
 
virtual ~G4MollerBhabhaModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
virtual G4double ComputeCrossSectionPerElectron (const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
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 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 ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
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 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 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)
 
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)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
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 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) const
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
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 *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Additional Inherited Members

- Protected Member Functions inherited from G4MollerBhabhaModel
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kinEnergy) final
 
void SetParticle (const G4ParticleDefinition *p)
 
- 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 G4MollerBhabhaModel
const G4ParticleDefinitionparticle
 
G4ParticleDefinitiontheElectron
 
G4ParticleChangeForLossfParticleChange
 
G4bool isElectron
 
G4double twoln10
 
G4double lowLimit
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const G4MaterialpBaseMaterial
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 
G4double inveplus
 
G4double pFactor
 

Detailed Description

Definition at line 58 of file G4PolarizedMollerBhabhaModel.hh.

Constructor & Destructor Documentation

◆ G4PolarizedMollerBhabhaModel()

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

Definition at line 64 of file G4PolarizedMollerBhabhaModel.cc.

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

◆ ~G4PolarizedMollerBhabhaModel()

G4PolarizedMollerBhabhaModel::~G4PolarizedMollerBhabhaModel ( )
virtual

Definition at line 86 of file G4PolarizedMollerBhabhaModel.cc.

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

Member Function Documentation

◆ ComputeCrossSectionPerElectron()

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

Reimplemented from G4MollerBhabhaModel.

Definition at line 96 of file G4PolarizedMollerBhabhaModel.cc.

101{
102 G4double xs =
104 cut,emax);
105// G4cout<<"calc eIoni xsec "<<xs<<G4endl;
106// G4cout<<" "<<kinEnergy<<" "<<cut<<" "<<emax<<G4endl;
107 G4double factor=1.;
108 if (xs!=0.) {
109 // G4cout<<"calc asym"<<G4endl;
110 G4double tmax = MaxSecondaryEnergy(pd, kinEnergy);
111 tmax = std::min(emax, tmax);
112
113 if (std::fabs(cut/emax-1.)<1.e-10) return xs;
114
115 if(cut < tmax) {
116
117 G4double xmin = cut/kinEnergy;
118 G4double xmax = tmax/kinEnergy;
119// G4cout<<"calc asym "<<xmin<<","<<xmax<<G4endl;
120 G4double gam = kinEnergy/electron_mass_c2 + 1.0;
121
122 G4double crossPol=crossSectionCalculator->
123 TotalXSection(xmin,xmax,gam,
124 theBeamPolarization,
125 theTargetPolarization);
126 G4double crossUnpol=crossSectionCalculator->
127 TotalXSection(xmin,xmax,gam,
130 if (crossUnpol>0.) factor=crossPol/crossUnpol;
131 // G4cout<<" factor="<<factor<<G4endl;
132 }
133 }
134 return xs*factor;
135}
double G4double
Definition: G4Types.hh:83
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) final
static const G4StokesVector ZERO

◆ GetBeamPolarization()

const G4StokesVector & G4PolarizedMollerBhabhaModel::GetBeamPolarization ( )
inline

Definition at line 95 of file G4PolarizedMollerBhabhaModel.hh.

96 {
97 return theBeamPolarization;
98 }

◆ GetFinalElectronPolarization()

const G4StokesVector & G4PolarizedMollerBhabhaModel::GetFinalElectronPolarization ( )
inline

Definition at line 99 of file G4PolarizedMollerBhabhaModel.hh.

100 {
101 return fElectronPolarization;
102 }

◆ GetFinalPositronPolarization()

const G4StokesVector & G4PolarizedMollerBhabhaModel::GetFinalPositronPolarization ( )
inline

Definition at line 103 of file G4PolarizedMollerBhabhaModel.hh.

104 {
105 return fPositronPolarization;
106 }

◆ GetTargetPolarization()

const G4StokesVector & G4PolarizedMollerBhabhaModel::GetTargetPolarization ( )
inline

Definition at line 91 of file G4PolarizedMollerBhabhaModel.hh.

92 {
93 return theTargetPolarization;
94 }

◆ SampleSecondaries()

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

Reimplemented from G4MollerBhabhaModel.

Definition at line 139 of file G4PolarizedMollerBhabhaModel.cc.

144{
145 // *** obtain and save target and beam polarization ***
147
148 const G4Track * aTrack = fParticleChange->GetCurrentTrack();
149
150 // obtain polarization of the beam
151 theBeamPolarization = dp->GetPolarization();
152
153 // obtain polarization of the media
154 G4VPhysicalVolume* aPVolume = aTrack->GetVolume();
155 G4LogicalVolume* aLVolume = aPVolume->GetLogicalVolume();
156 const G4bool targetIsPolarized = polarizationManger->IsPolarized(aLVolume);
157 theTargetPolarization = polarizationManger->GetVolumePolarization(aLVolume);
158
159 // transfer target polarization in interaction frame
160 if (targetIsPolarized)
161 theTargetPolarization.rotateUz(dp->GetMomentumDirection());
162
163
164
165
166 G4double tmax = std::min(maxEnergy, MaxSecondaryKinEnergy(dp));
167 if(tmin >= tmax) return;
168 // if(tmin > tmax) tmin = tmax;
169
170 G4double polL = theBeamPolarization.z()*theTargetPolarization.z();
171 polL=std::fabs(polL);
172 G4double polT = theBeamPolarization.x()*theTargetPolarization.x() +
173 theBeamPolarization.y()*theTargetPolarization.y();
174 polT=std::fabs(polT);
175
176 G4double kineticEnergy = dp->GetKineticEnergy();
177 G4double energy = kineticEnergy + electron_mass_c2;
178 G4double totalMomentum = std::sqrt(kineticEnergy*(energy + electron_mass_c2));
179 G4double xmin = tmin/kineticEnergy;
180 G4double xmax = tmax/kineticEnergy;
181 G4double gam = energy/electron_mass_c2;
182 G4double gamma2 = gam*gam;
183 G4double gmo = gam - 1.;
184 G4double gmo2 = gmo*gmo;
185 G4double gmo3 = gmo2*gmo;
186 G4double gpo = gam + 1.;
187 G4double gpo2 = gpo*gpo;
188 G4double gpo3 = gpo2*gpo;
189 G4double x, y, q, grej, grej2;
190 G4double z = 0.;
191 G4double xs = 0., phi =0.;
192 G4ThreeVector direction = dp->GetMomentumDirection();
193
194 //(Polarized) Moller (e-e-) scattering
195 if (isElectron) {
196 // *** dice according to polarized cross section
197 G4double G = ((2.0*gam - 1.0)/gamma2)*(1. - polT - polL*gam);
198 G4double H = (sqr(gam - 1.0)/gamma2)*(1. + polT + polL*((gam + 3.)/(gam - 1.)));
199
200 y = 1.0 - xmax;
201 grej = 1.0 - G*xmax + xmax*xmax*(H + (1.0 - G*y)/(y*y));
202 grej2 = 1.0 - G*xmin + xmin*xmin*(H + (1.0 - G*y)/(y*y));
203 if (grej2 > grej) grej = grej2;
204 G4double prefM = gamma2*classic_electr_radius*classic_electr_radius/(gmo2*(gam + 1.0));
205 grej *= prefM;
206 do {
207 q = G4UniformRand();
208 x = xmin*xmax/(xmin*(1.0 - q) + xmax*q);
209 if (crossSectionCalculator) {
210 crossSectionCalculator->Initialize(x,gam,phi,theBeamPolarization,
211 theTargetPolarization,1);
212 xs=crossSectionCalculator->XSection(G4StokesVector::ZERO,
214 z=xs*sqr(x)*4.;
215 if (grej < z) {
216 G4cout<<"WARNING : error in Moller rejection routine! \n"
217 <<" z = "<<z<<" grej="<<grej<<"\n";
218 }
219 } else {
220 G4cout<<"No calculator in Moller scattering"<<G4endl;
221 }
222 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
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.);
235 G4double prefB = classic_electr_radius*classic_electr_radius/(gam - 1.0);
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 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
260 } while(grej * G4UniformRand() > z);
261 }
262 //
263 //
264 // polar asymmetries (due to transverse polarizations)
265 //
266 //
267 if (crossSectionCalculator) {
268 // grej*=1./(sqr(x)*sqr(gamma2-1))*sqr(gam*(1+gam));
269 grej=xs*2.;
270 do {
271 phi = twopi * G4UniformRand() ;
272 crossSectionCalculator->Initialize(x,gam,phi,theBeamPolarization,
273 theTargetPolarization,1);
274 xs=crossSectionCalculator->XSection(G4StokesVector::ZERO,
276 if(xs > grej) {
277 if (isElectron){
278 G4cout << "G4PolarizedMollerBhabhaModel::SampleSecondaries Warning! "<<G4endl
279 << "Majorant " << grej << " < "
280 << xs << " for phi= " << phi<<G4endl
281 << " e-e- (Moller) scattering"<< G4endl
282 <<"PHI DICING"<<G4endl;
283 } else {
284 G4cout << "G4PolarizedMollerBhabhaModel::SampleSecondaries Warning! "<<G4endl
285 << "Majorant " << grej << " < "
286 << xs << " for phi= " << phi<<G4endl
287 << " e+e- (Bhabha) scattering"<< G4endl
288 <<"PHI DICING"<<G4endl;
289 }
290 }
291 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
292 } while(grej * G4UniformRand() > xs);
293 }
294
295 // fix kinematics of delta electron
296 G4double deltaKinEnergy = x * kineticEnergy;
297 G4double deltaMomentum =
298 std::sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
299 G4double cost = deltaKinEnergy * (energy + electron_mass_c2) /
300 (deltaMomentum * totalMomentum);
301 G4double sint = 1.0 - cost*cost;
302 if(sint > 0.0) sint = std::sqrt(sint);
303
304
305 G4ThreeVector deltaDirection(-sint*std::cos(phi),-sint*std::sin(phi), cost) ;
306 deltaDirection.rotateUz(direction);
307
308 // primary change
309 kineticEnergy -= deltaKinEnergy;
311
312 if(kineticEnergy > DBL_MIN) {
313 G4ThreeVector dir = totalMomentum*direction - deltaMomentum*deltaDirection;
314 direction = dir.unit();
316 }
317
318 // create G4DynamicParticle object for delta ray
319 G4DynamicParticle* delta = new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
320 vdp->push_back(delta);
321
322 // get interaction frame
323 G4ThreeVector nInteractionFrame =
324 G4PolarizationHelper::GetFrame(direction,deltaDirection);
325
326 if (crossSectionCalculator) {
327 // calculate mean final state polarizations
328
329 theBeamPolarization.InvRotateAz(nInteractionFrame,direction);
330 theTargetPolarization.InvRotateAz(nInteractionFrame,direction);
331 crossSectionCalculator->Initialize(x,gam,phi,theBeamPolarization,
332 theTargetPolarization,2);
333
334 // electron/positron
335 fPositronPolarization=crossSectionCalculator->GetPol2();
336 fPositronPolarization.RotateAz(nInteractionFrame,direction);
337
338 fParticleChange->ProposePolarization(fPositronPolarization);
339
340 // electron
341 fElectronPolarization=crossSectionCalculator->GetPol3();
342 fElectronPolarization.RotateAz(nInteractionFrame,deltaDirection);
343 delta->SetPolarization(fElectronPolarization.x(),
344 fElectronPolarization.y(),
345 fElectronPolarization.z());
346 }
347 else {
348 fPositronPolarization=G4ThreeVector();
349 fElectronPolarization=G4ThreeVector();
350 }
351}
CLHEP::Hep3Vector G4ThreeVector
bool G4bool
Definition: G4Types.hh:86
#define G4UniformRand()
Definition: Randomize.hh:52
double z() const
Hep3Vector unit() const
double x() const
double y() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
void SetPolarization(const G4ThreeVector &)
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
const G4ThreeVector & GetPolarization() const
G4ParticleChangeForLoss * fParticleChange
void SetProposedKineticEnergy(G4double proposedKinEnergy)
const G4Track * GetCurrentTrack() const
void ProposePolarization(const G4ThreeVector &dir)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
static G4ThreeVector GetFrame(const G4ThreeVector &, const G4ThreeVector &)
bool IsPolarized(G4LogicalVolume *lVol) const
static G4PolarizationManager * GetInstance()
const G4ThreeVector & GetVolumePolarization(G4LogicalVolume *lVol) const
void InvRotateAz(G4ThreeVector nInteractionFrame, G4ThreeVector particleDirection)
void RotateAz(G4ThreeVector nInteractionFrame, G4ThreeVector particleDirection)
G4VPhysicalVolume * GetVolume() const
G4double MaxSecondaryKinEnergy(const G4DynamicParticle *dynParticle)
Definition: G4VEmModel.hh:510
G4LogicalVolume * GetLogicalVolume() const
virtual G4double XSection(const G4StokesVector &pol2, const G4StokesVector &pol3)=0
virtual G4StokesVector GetPol3()
virtual G4StokesVector GetPol2()
virtual void Initialize(G4double, G4double, G4double, const G4StokesVector &p0, const G4StokesVector &p1, G4int flag=0)
G4double energy(const ThreeVector &p, const G4double m)
T sqr(const T &x)
Definition: templates.hh:128
#define DBL_MIN
Definition: templates.hh:54

◆ SetBeamPolarization()

void G4PolarizedMollerBhabhaModel::SetBeamPolarization ( const G4ThreeVector pBeam)
inline

Definition at line 87 of file G4PolarizedMollerBhabhaModel.hh.

88 {
89 theBeamPolarization = pBeam;
90 }

◆ SetTargetPolarization()

void G4PolarizedMollerBhabhaModel::SetTargetPolarization ( const G4ThreeVector pTarget)
inline

Definition at line 83 of file G4PolarizedMollerBhabhaModel.hh.

84 {
85 theTargetPolarization = pTarget;
86 }

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