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

#include <G4PolarizedIonisationModel.hh>

+ Inheritance diagram for G4PolarizedIonisationModel:

Public Member Functions

 G4PolarizedIonisationModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="PolarizedMollerBhabha")
 
virtual ~G4PolarizedIonisationModel () override
 
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
 
 G4PolarizedIonisationModel (G4PolarizedIonisationModel &)=delete
 
G4PolarizedIonisationModeloperator= (const G4PolarizedIonisationModel &right)=delete
 
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")
 
 ~G4MollerBhabhaModel () override
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
virtual G4double ComputeCrossSectionPerElectron (const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
 
G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
 
G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
G4MollerBhabhaModeloperator= (const G4MollerBhabhaModel &right)=delete
 
 G4MollerBhabhaModel (const G4MollerBhabhaModel &)=delete
 
- 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 *, const G4double &length, G4double &eloss)
 
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 FillNumberOfSecondaries (G4int &numberOfTriplets, G4int &numberOfRecoil)
 
virtual void ModelDescription (std::ostream &outFile) const
 
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)
 
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)
 
const G4ElementGetCurrentElement (const G4Material *mat=nullptr) const
 
G4int SelectRandomAtomNumber (const G4Material *) const
 
const G4IsotopeGetCurrentIsotope (const G4Element *elm=nullptr) const
 
G4int SelectIsotopeNumber (const G4Element *) const
 
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 *)
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Additional Inherited Members

- Protected Member Functions inherited from G4MollerBhabhaModel
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 = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4PhysicsTablexSectionTable = nullptr
 
const G4MaterialpBaseMaterial = nullptr
 
const std::vector< G4double > * theDensityFactor = nullptr
 
const std::vector< G4int > * theDensityIdx = nullptr
 
G4double inveplus
 
G4double pFactor = 1.0
 
size_t currentCoupleIndex = 0
 
size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 

Detailed Description

Definition at line 52 of file G4PolarizedIonisationModel.hh.

Constructor & Destructor Documentation

◆ G4PolarizedIonisationModel() [1/2]

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

Definition at line 49 of file G4PolarizedIonisationModel.cc.

51 : G4MollerBhabhaModel(p, nam)
52 , fCrossSectionCalculator(nullptr)
53{
54 fBeamPolarization = G4StokesVector::ZERO;
55 fTargetPolarization = G4StokesVector::ZERO;
56
57 fPositronPolarization = G4StokesVector::ZERO;
58 fElectronPolarization = G4StokesVector::ZERO;
59
60 isElectron = (p == theElectron); // necessary due to wrong order in
61 // G4MollerBhabhaModel constructor!
62
63 if(!isElectron)
64 {
65 G4cout << " buildBhabha cross section " << isElectron << G4endl;
66 fCrossSectionCalculator = new G4PolarizedIonisationBhabhaXS();
67 }
68 else
69 {
70 G4cout << " buildMoller cross section " << isElectron << G4endl;
71 fCrossSectionCalculator = new G4PolarizedIonisationMollerXS();
72 }
73}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4ParticleDefinition * theElectron
static const G4StokesVector ZERO

◆ ~G4PolarizedIonisationModel()

G4PolarizedIonisationModel::~G4PolarizedIonisationModel ( )
overridevirtual

Definition at line 76 of file G4PolarizedIonisationModel.cc.

77{
78 if(fCrossSectionCalculator)
79 {
80 delete fCrossSectionCalculator;
81 }
82}

◆ G4PolarizedIonisationModel() [2/2]

G4PolarizedIonisationModel::G4PolarizedIonisationModel ( G4PolarizedIonisationModel )
delete

Member Function Documentation

◆ ComputeCrossSectionPerElectron()

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

Reimplemented from G4MollerBhabhaModel.

Definition at line 85 of file G4PolarizedIonisationModel.cc.

88{
90 pd, kinEnergy, cut, emax);
91 G4double factor = 1.;
92 if(xs != 0.)
93 {
94 G4double tmax = MaxSecondaryEnergy(pd, kinEnergy);
95 tmax = std::min(emax, tmax);
96
97 if(std::fabs(cut / emax - 1.) < 1.e-10)
98 return xs;
99
100 if(cut < tmax)
101 {
102 G4double xmin = cut / kinEnergy;
103 G4double xmax = tmax / kinEnergy;
104 G4double gam = kinEnergy / electron_mass_c2 + 1.0;
105 G4double crossPol = fCrossSectionCalculator->TotalXSection(
106 xmin, xmax, gam, fBeamPolarization, fTargetPolarization);
107 G4double crossUnpol = fCrossSectionCalculator->TotalXSection(
109 if(crossUnpol > 0.)
110 factor = crossPol / crossUnpol;
111 }
112 }
113 return xs * factor;
114}
double G4double
Definition: G4Types.hh:83
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) final
virtual G4double TotalXSection(G4double xmin, G4double xmax, G4double y, const G4StokesVector &pol0, const G4StokesVector &pol1)

◆ GetBeamPolarization()

const G4StokesVector & G4PolarizedIonisationModel::GetBeamPolarization ( )
inline

Definition at line 84 of file G4PolarizedIonisationModel.hh.

84{ return fBeamPolarization; }

◆ GetFinalElectronPolarization()

const G4StokesVector & G4PolarizedIonisationModel::GetFinalElectronPolarization ( )
inline

Definition at line 85 of file G4PolarizedIonisationModel.hh.

86 {
87 return fElectronPolarization;
88 }

◆ GetFinalPositronPolarization()

const G4StokesVector & G4PolarizedIonisationModel::GetFinalPositronPolarization ( )
inline

Definition at line 89 of file G4PolarizedIonisationModel.hh.

90 {
91 return fPositronPolarization;
92 }

◆ GetTargetPolarization()

const G4StokesVector & G4PolarizedIonisationModel::GetTargetPolarization ( )
inline

Definition at line 83 of file G4PolarizedIonisationModel.hh.

83{ return fTargetPolarization; }

◆ operator=()

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 117 of file G4PolarizedIonisationModel.cc.

120{
121 // *** obtain and save target and beam polarization ***
122 G4PolarizationManager* polarizationManager =
124
125 const G4Track* aTrack = fParticleChange->GetCurrentTrack();
126
127 // obtain polarization of the beam
128 fBeamPolarization = G4StokesVector(dp->GetPolarization());
129
130 // obtain polarization of the media
131 G4VPhysicalVolume* aPVolume = aTrack->GetVolume();
132 G4LogicalVolume* aLVolume = aPVolume->GetLogicalVolume();
133 const G4bool targetIsPolarized = polarizationManager->IsPolarized(aLVolume);
134 fTargetPolarization = polarizationManager->GetVolumePolarization(aLVolume);
135
136 // transfer target polarization in interaction frame
137 if(targetIsPolarized)
138 fTargetPolarization.rotateUz(dp->GetMomentumDirection());
139
140 G4double tmax = std::min(maxEnergy, MaxSecondaryKinEnergy(dp));
141 if(tmin >= tmax)
142 return;
143
144 G4double polL = fBeamPolarization.z() * fTargetPolarization.z();
145 polL = std::fabs(polL);
146 G4double polT = fBeamPolarization.x() * fTargetPolarization.x() +
147 fBeamPolarization.y() * fTargetPolarization.y();
148 polT = std::fabs(polT);
149
150 G4double kineticEnergy = dp->GetKineticEnergy();
151 G4double energy = kineticEnergy + electron_mass_c2;
152 G4double totalMomentum =
153 std::sqrt(kineticEnergy * (energy + electron_mass_c2));
154 G4double xmin = tmin / kineticEnergy;
155 G4double xmax = tmax / kineticEnergy;
156 G4double gam = energy / electron_mass_c2;
157 G4double gamma2 = gam * gam;
158 G4double gmo = gam - 1.;
159 G4double gmo2 = gmo * gmo;
160 G4double gmo3 = gmo2 * gmo;
161 G4double gpo = gam + 1.;
162 G4double gpo2 = gpo * gpo;
163 G4double gpo3 = gpo2 * gpo;
164 G4double x, y, q, grej, grej2;
165 G4double z = 0.;
166 G4double xs = 0., phi = 0.;
167 G4ThreeVector direction = dp->GetMomentumDirection();
168
169 //(Polarized) Moller (e-e-) scattering
170 if(isElectron)
171 {
172 // *** dice according to polarized cross section
173 G4double G = ((2.0 * gam - 1.0) / gamma2) * (1. - polT - polL * gam);
174 G4double H = (sqr(gam - 1.0) / gamma2) *
175 (1. + polT + polL * ((gam + 3.) / (gam - 1.)));
176
177 y = 1.0 - xmax;
178 grej = 1.0 - G * xmax + xmax * xmax * (H + (1.0 - G * y) / (y * y));
179 grej2 = 1.0 - G * xmin + xmin * xmin * (H + (1.0 - G * y) / (y * y));
180 if(grej2 > grej)
181 grej = grej2;
182 G4double prefM = gamma2 * classic_electr_radius * classic_electr_radius /
183 (gmo2 * (gam + 1.0));
184 grej *= prefM;
185 do
186 {
187 q = G4UniformRand();
188 x = xmin * xmax / (xmin * (1.0 - q) + xmax * q);
189 if(fCrossSectionCalculator)
190 {
191 fCrossSectionCalculator->Initialize(x, gam, phi, fBeamPolarization,
192 fTargetPolarization, 1);
193 xs = fCrossSectionCalculator->XSection(G4StokesVector::ZERO,
195 z = xs * sqr(x) * 4.;
196 if(grej < z)
197 {
199 ed << "WARNING : error in Moller rejection routine! \n"
200 << " z = " << z << " grej=" << grej << "\n";
201 G4Exception("G4PolarizedIonisationModel::SampleSecondaries", "pol019",
202 JustWarning, ed);
203 }
204 }
205 else
206 {
207 G4cout << "No calculator in Moller scattering" << G4endl;
208 }
209 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
210 } while(grej * G4UniformRand() > z);
211 // Bhabha (e+e-) scattering
212 }
213 else
214 {
215 // *** dice according to polarized cross section
216 y = xmax * xmax;
217 grej = 0.;
218 grej += y * y * gmo3 * (1. + (polL + polT) * (gam + 3.) / gmo);
219 grej += -2. * xmin * xmin * xmin * gam * gmo2 *
220 (1. - (polL + polT) * (gam + 3.) / gmo);
221 grej += y * y * gmo * (3. * gamma2 + 6. * gam + 4.) *
222 (1. + (polL * (3. * gam + 1.) * (gamma2 + gam + 1.) +
223 polT * ((gam + 2.) * gamma2 + 1.)) /
224 (gmo * (3. * gam * (gam + 2.) + 4.)));
225 grej /= gpo3;
226 grej += -xmin * (2. * gamma2 + 4. * gam + 1.) *
227 (1. - gam * (polL * (2. * gam + 1.) + polT) /
228 (2. * gam * (gam + 2.) + 1.)) /
229 gpo2;
230 grej += gamma2 / (gamma2 - 1.);
231 G4double prefB =
232 classic_electr_radius * classic_electr_radius / (gam - 1.0);
233 grej *= prefB;
234
235 do
236 {
237 q = G4UniformRand();
238 x = xmin * xmax / (xmin * (1.0 - q) + xmax * q);
239 if(fCrossSectionCalculator)
240 {
241 fCrossSectionCalculator->Initialize(x, gam, phi, fBeamPolarization,
242 fTargetPolarization, 1);
243 xs = fCrossSectionCalculator->XSection(G4StokesVector::ZERO,
245 z = xs * sqr(x) * 4.;
246 }
247 else
248 {
249 G4cout << "No calculator in Bhabha scattering" << G4endl;
250 }
251
252 if(z > grej)
253 {
255 ed << "G4PolarizedIonisationModel::SampleSecondaries Warning!\n "
256 << "Majorant " << grej << " < " << z << " for x= " << x << G4endl
257 << " e+e- (Bhabha) scattering"
258 << " at KinEnergy " << kineticEnergy << G4endl;
259 G4Exception("G4PolarizedIonisationModel::SampleSecondaries", "pol020",
260 JustWarning, ed);
261 }
262 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
263 } while(grej * G4UniformRand() > z);
264 }
265
266 // polar asymmetries (due to transverse polarizations)
267 if(fCrossSectionCalculator)
268 {
269 grej = xs * 2.;
270 do
271 {
272 phi = twopi * G4UniformRand();
273 fCrossSectionCalculator->Initialize(x, gam, phi, fBeamPolarization,
274 fTargetPolarization, 1);
275 xs = fCrossSectionCalculator->XSection(G4StokesVector::ZERO,
277 if(xs > grej)
278 {
279 if(isElectron)
280 {
282 ed << "Majorant " << grej << " < " << xs << " for phi= " << phi
283 << "\n"
284 << " e-e- (Moller) scattering\n"
285 << "PHI DICING\n";
286 G4Exception("G4PolarizedIonisationModel::SampleSecondaries", "pol021",
287 JustWarning, ed);
288 }
289 else
290 {
292 ed << "Majorant " << grej << " < " << xs << " for phi= " << phi
293 << "\n"
294 << " e+e- (Bhabha) scattering\n"
295 << "PHI DICING\n";
296 G4Exception("G4PolarizedIonisationModel::SampleSecondaries", "pol022",
297 JustWarning, ed);
298 }
299 }
300 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
301 } while(grej * G4UniformRand() > xs);
302 }
303
304 // fix kinematics of delta electron
305 G4double deltaKinEnergy = x * kineticEnergy;
306 G4double deltaMomentum =
307 std::sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0 * electron_mass_c2));
308 G4double cost = deltaKinEnergy * (energy + electron_mass_c2) /
309 (deltaMomentum * totalMomentum);
310 G4double sint = 1.0 - cost * cost;
311 if(sint > 0.0)
312 sint = std::sqrt(sint);
313
314 G4ThreeVector deltaDirection(-sint * std::cos(phi), -sint * std::sin(phi),
315 cost);
316 deltaDirection.rotateUz(direction);
317
318 // primary change
319 kineticEnergy -= deltaKinEnergy;
321
322 if(kineticEnergy > DBL_MIN)
323 {
324 G4ThreeVector dir =
325 totalMomentum * direction - deltaMomentum * deltaDirection;
326 direction = dir.unit();
328 }
329
330 // create G4DynamicParticle object for delta ray
331 G4DynamicParticle* delta =
332 new G4DynamicParticle(theElectron, deltaDirection, deltaKinEnergy);
333 vdp->push_back(delta);
334
335 // get interaction frame
336 G4ThreeVector nInteractionFrame =
337 G4PolarizationHelper::GetFrame(direction, deltaDirection);
338
339 if(fCrossSectionCalculator)
340 {
341 // calculate mean final state polarizations
342 fBeamPolarization.InvRotateAz(nInteractionFrame, direction);
343 fTargetPolarization.InvRotateAz(nInteractionFrame, direction);
344 fCrossSectionCalculator->Initialize(x, gam, phi, fBeamPolarization,
345 fTargetPolarization, 2);
346
347 // electron/positron
348 fPositronPolarization = fCrossSectionCalculator->GetPol2();
349 fPositronPolarization.RotateAz(nInteractionFrame, direction);
350
351 fParticleChange->ProposePolarization(fPositronPolarization);
352
353 // electron
354 fElectronPolarization = fCrossSectionCalculator->GetPol3();
355 fElectronPolarization.RotateAz(nInteractionFrame, deltaDirection);
356 delta->SetPolarization(fElectronPolarization.x(), fElectronPolarization.y(),
357 fElectronPolarization.z());
358 }
359 else
360 {
361 fPositronPolarization = G4StokesVector::ZERO;
362 fElectronPolarization = G4StokesVector::ZERO;
363 }
364}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
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)
void ProposePolarization(const G4ThreeVector &dir)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
static G4ThreeVector GetFrame(const G4ThreeVector &, const G4ThreeVector &)
bool IsPolarized(G4LogicalVolume *lVol) const
const G4StokesVector GetVolumePolarization(G4LogicalVolume *lVol) const
static G4PolarizationManager * GetInstance()
void InvRotateAz(G4ThreeVector nInteractionFrame, G4ThreeVector particleDirection)
void RotateAz(G4ThreeVector nInteractionFrame, G4ThreeVector particleDirection)
G4VPhysicalVolume * GetVolume() const
G4double MaxSecondaryKinEnergy(const G4DynamicParticle *dynParticle)
Definition: G4VEmModel.hh:501
const G4Track * GetCurrentTrack() const
G4LogicalVolume * GetLogicalVolume() const
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 G4StokesVector GetPol2()
virtual G4StokesVector GetPol3()
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 G4PolarizedIonisationModel::SetBeamPolarization ( const G4ThreeVector pBeam)
inline

Definition at line 79 of file G4PolarizedIonisationModel.hh.

80 {
81 fBeamPolarization = G4StokesVector(pBeam);
82 }

◆ SetTargetPolarization()

void G4PolarizedIonisationModel::SetTargetPolarization ( const G4ThreeVector pTarget)
inline

Definition at line 75 of file G4PolarizedIonisationModel.hh.

76 {
77 fTargetPolarization = G4StokesVector(pTarget);
78 }

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