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

#include <G4PolarizedComptonModel.hh>

+ Inheritance diagram for G4PolarizedComptonModel:

Public Member Functions

 G4PolarizedComptonModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="Polarized-Compton")
 
virtual ~G4PolarizedComptonModel () override
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax) override
 
G4double ComputeAsymmetryPerAtom (G4double gammaEnergy, G4double Z)
 
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 G4ThreeVectorGetTargetPolarization () const
 
const G4ThreeVectorGetBeamPolarization () const
 
const G4ThreeVectorGetFinalGammaPolarization () const
 
const G4ThreeVectorGetFinalElectronPolarization () const
 
G4PolarizedComptonModeloperator= (const G4PolarizedComptonModel &right)=delete
 
 G4PolarizedComptonModel (const G4PolarizedComptonModel &)=delete
 
- Public Member Functions inherited from G4KleinNishinaCompton
 G4KleinNishinaCompton (const G4ParticleDefinition *p=nullptr, const G4String &nam="Klein-Nishina")
 
 ~G4KleinNishinaCompton () override
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel) override
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
G4KleinNishinaComptonoperator= (const G4KleinNishinaCompton &right)=delete
 
 G4KleinNishinaCompton (const G4KleinNishinaCompton &)=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 G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 
- Protected Attributes inherited from G4KleinNishinaCompton
G4ParticleDefinitiontheGamma
 
G4ParticleDefinitiontheElectron
 
G4ParticleChangeForGammafParticleChange
 
G4double lowestSecondaryEnergy
 
- 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 51 of file G4PolarizedComptonModel.hh.

Constructor & Destructor Documentation

◆ G4PolarizedComptonModel() [1/2]

G4PolarizedComptonModel::G4PolarizedComptonModel ( const G4ParticleDefinition p = nullptr,
const G4String nam = "Polarized-Compton" 
)
explicit

Definition at line 47 of file G4PolarizedComptonModel.cc.

49 : G4KleinNishinaCompton(nullptr, nam)
50 , fVerboseLevel(0)
51{
52 fCrossSectionCalculator = new G4PolarizedComptonXS();
53 fBeamPolarization = G4StokesVector::ZERO;
54 fTargetPolarization = G4StokesVector::ZERO;
55}
static const G4StokesVector ZERO

◆ ~G4PolarizedComptonModel()

G4PolarizedComptonModel::~G4PolarizedComptonModel ( )
overridevirtual

Definition at line 58 of file G4PolarizedComptonModel.cc.

59{
60 delete fCrossSectionCalculator;
61}

◆ G4PolarizedComptonModel() [2/2]

G4PolarizedComptonModel::G4PolarizedComptonModel ( const G4PolarizedComptonModel )
delete

Member Function Documentation

◆ ComputeAsymmetryPerAtom()

G4double G4PolarizedComptonModel::ComputeAsymmetryPerAtom ( G4double  gammaEnergy,
G4double  Z 
)

Definition at line 64 of file G4PolarizedComptonModel.cc.

66{
67 G4double asymmetry = 0.0;
68
69 G4double k0 = gammaEnergy / electron_mass_c2;
70 G4double k1 = 1. + 2. * k0;
71
72 asymmetry = -k0;
73 asymmetry *=
74 (k0 + 1.) * sqr(k1) * G4Log(k1) - 2. * k0 * (5. * sqr(k0) + 4. * k0 + 1.);
75 asymmetry /= ((k0 - 2.) * k0 - 2.) * sqr(k1) * G4Log(k1) +
76 2. * k0 * (k0 * (k0 + 1.) * (k0 + 8.) + 2.);
77
78 if(asymmetry > 1.)
79 {
81 ed << "ERROR in G4PolarizedComptonModel::ComputeAsymmetryPerAtom.\n"
82 << " asymmetry = " << asymmetry << "\n";
83 G4Exception("G4PolarizedComptonModel::ComputeAsymmetryPerAtom", "pol035",
84 JustWarning, ed);
85 }
86
87 return asymmetry;
88}
@ 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
G4double G4Log(G4double x)
Definition: G4Log.hh:227
double G4double
Definition: G4Types.hh:83
T sqr(const T &x)
Definition: templates.hh:128

Referenced by ComputeCrossSectionPerAtom().

◆ ComputeCrossSectionPerAtom()

G4double G4PolarizedComptonModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition pd,
G4double  kinEnergy,
G4double  Z,
G4double  A,
G4double  cut,
G4double  emax 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 91 of file G4PolarizedComptonModel.cc.

94{
96 pd, kinEnergy, Z, A, cut, emax);
97 G4double polzz = fBeamPolarization.p3() * fTargetPolarization.z();
98 if(polzz > 0.0)
99 {
100 G4double asym = ComputeAsymmetryPerAtom(kinEnergy, Z);
101 xs *= (1. + polzz * asym);
102 }
103 return xs;
104}
const G4int Z[17]
const G4double A[17]
double z() const
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax) override
G4double ComputeAsymmetryPerAtom(G4double gammaEnergy, G4double Z)
G4double p3() const

◆ GetBeamPolarization()

const G4ThreeVector & G4PolarizedComptonModel::GetBeamPolarization ( ) const
inline

Definition at line 117 of file G4PolarizedComptonModel.hh.

118{
119 return fBeamPolarization;
120}

◆ GetFinalElectronPolarization()

const G4ThreeVector & G4PolarizedComptonModel::GetFinalElectronPolarization ( ) const
inline

Definition at line 127 of file G4PolarizedComptonModel.hh.

128{
129 return finalElectronPolarization;
130}

◆ GetFinalGammaPolarization()

const G4ThreeVector & G4PolarizedComptonModel::GetFinalGammaPolarization ( ) const
inline

Definition at line 121 of file G4PolarizedComptonModel.hh.

123{
124 return fFinalGammaPolarization;
125}

◆ GetTargetPolarization()

const G4ThreeVector & G4PolarizedComptonModel::GetTargetPolarization ( ) const
inline

Definition at line 112 of file G4PolarizedComptonModel.hh.

114{
115 return fTargetPolarization;
116}

◆ operator=()

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

◆ SampleSecondaries()

void G4PolarizedComptonModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple ,
const G4DynamicParticle aDynamicGamma,
G4double  tmin,
G4double  maxEnergy 
)
overridevirtual

Implements G4VEmModel.

Definition at line 107 of file G4PolarizedComptonModel.cc.

110{
111 // do nothing below the threshold
112 if(aDynamicGamma->GetKineticEnergy() <= LowEnergyLimit())
113 {
114 return;
115 }
116
117 const G4Track* aTrack = fParticleChange->GetCurrentTrack();
118 G4VPhysicalVolume* aPVolume = aTrack->GetVolume();
119 G4LogicalVolume* aLVolume = aPVolume->GetLogicalVolume();
120
121 if(fVerboseLevel >= 1)
122 {
123 G4cout << "G4PolarizedComptonModel::SampleSecondaries in "
124 << aLVolume->GetName() << G4endl;
125 }
126 G4PolarizationManager* polarizationManager =
128
129 // obtain polarization of the beam
130 fBeamPolarization = G4StokesVector(aDynamicGamma->GetPolarization());
131 fBeamPolarization.SetPhoton();
132
133 // obtain polarization of the media
134 G4bool targetIsPolarized = polarizationManager->IsPolarized(aLVolume);
135 fTargetPolarization = polarizationManager->GetVolumePolarization(aLVolume);
136
137 // if beam is linear polarized or target is transversely polarized
138 // determine the angle to x-axis
139 // (assumes same PRF as in the polarization definition)
140 G4ThreeVector gamDirection0 = aDynamicGamma->GetMomentumDirection();
141
142 // transfer fTargetPolarization
143 // into the gamma frame (problem electron is at rest)
144 if(targetIsPolarized)
145 {
146 fTargetPolarization.rotateUz(gamDirection0);
147 }
148 // The scattered gamma energy is sampled according to
149 // Klein - Nishina formula.
150 // The random number techniques of Butcher & Messel are used
151 // (Nuc Phys 20(1960),15).
152 // Note : Effects due to binding of atomic electrons are neglected.
153
154 G4double gamEnergy0 = aDynamicGamma->GetKineticEnergy();
155 G4double E0_m = gamEnergy0 / electron_mass_c2;
156
157 // sample the energy rate of the scattered gamma
158 G4double epsilon, sint2;
159 G4double onecost = 0.0;
160 G4double Phi = 0.0;
161 G4double greject = 1.0;
162 G4double cosTeta = 1.0;
163 G4double sinTeta = 0.0;
164
165 G4double eps0 = 1. / (1. + 2. * E0_m);
166 G4double epsilon0sq = eps0 * eps0;
167 G4double alpha1 = -G4Log(eps0);
168 G4double alpha2 = alpha1 + 0.5 * (1. - epsilon0sq);
169
170 G4double polarization = fBeamPolarization.p3() * fTargetPolarization.p3();
171
172 CLHEP::HepRandomEngine* rndmEngineMod = G4Random::getTheEngine();
173 G4int nloop = 0;
174 G4bool end = false;
175
176 G4double rndm[3];
177
178 do
179 {
180 do
181 {
182 ++nloop;
183 // false interaction if too many iterations
184 if(nloop > fLoopLim)
185 {
186 PrintWarning(aDynamicGamma, nloop, greject, onecost, Phi,
187 "too many iterations");
188 return;
189 }
190
191 // 3 random numbers to sample scattering
192 rndmEngineMod->flatArray(3, rndm);
193
194 if(alpha1 > alpha2 * rndm[0])
195 {
196 epsilon = G4Exp(-alpha1 * rndm[1]);
197 }
198 else
199 {
200 epsilon = std::sqrt(epsilon0sq + (1. - epsilon0sq) * rndm[1]);
201 }
202
203 onecost = (1. - epsilon) / (epsilon * E0_m);
204 sint2 = onecost * (2. - onecost);
205
206 G4double gdiced = 2. * (1. / epsilon + epsilon);
207 G4double gdist = 1. / epsilon + epsilon - sint2 -
208 polarization * (1. / epsilon - epsilon) * (1. - onecost);
209
210 greject = gdist / gdiced;
211
212 if(greject > 1.0)
213 {
214 PrintWarning(aDynamicGamma, nloop, greject, onecost, Phi,
215 "theta majoranta wrong");
216 }
217 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
218 } while(greject < rndm[2]);
219
220 // assuming phi loop successful
221 end = true;
222
223 // scattered gamma angles. ( Z - axis along the parent gamma)
224 cosTeta = 1. - onecost;
225 sinTeta = std::sqrt(sint2);
226 do
227 {
228 ++nloop;
229
230 // 2 random numbers to sample scattering
231 rndmEngineMod->flatArray(2, rndm);
232
233 // false interaction if too many iterations
234 Phi = twopi * rndm[0];
235 if(nloop > fLoopLim)
236 {
237 PrintWarning(aDynamicGamma, nloop, greject, onecost, Phi,
238 "too many iterations");
239 return;
240 }
241
242 G4double gdiced = 1. / epsilon + epsilon - sint2 +
243 std::abs(fBeamPolarization.p3()) *
244 (std::abs((1. / epsilon - epsilon) * cosTeta *
245 fTargetPolarization.p3()) +
246 (1. - epsilon) * sinTeta *
247 (std::sqrt(sqr(fTargetPolarization.p1()) +
248 sqr(fTargetPolarization.p2())))) +
249 sint2 * (std::sqrt(sqr(fBeamPolarization.p1()) +
250 sqr(fBeamPolarization.p2())));
251
252 G4double gdist =
253 1. / epsilon + epsilon - sint2 +
254 fBeamPolarization.p3() *
255 ((1. / epsilon - epsilon) * cosTeta * fTargetPolarization.p3() +
256 (1. - epsilon) * sinTeta *
257 (std::cos(Phi) * fTargetPolarization.p1() +
258 std::sin(Phi) * fTargetPolarization.p2())) -
259 sint2 * (std::cos(2. * Phi) * fBeamPolarization.p1() +
260 std::sin(2. * Phi) * fBeamPolarization.p2());
261 greject = gdist / gdiced;
262
263 if(greject > 1.0)
264 {
265 PrintWarning(aDynamicGamma, nloop, greject, onecost, Phi,
266 "phi majoranta wrong");
267 }
268
269 if(greject < 1.e-3)
270 {
271 PrintWarning(aDynamicGamma, nloop, greject, onecost, Phi,
272 "phi loop ineffective");
273 // restart theta loop
274 end = false;
275 break;
276 }
277
278 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
279 } while(greject < rndm[1]);
280 } while(!end);
281 G4double dirx = sinTeta * std::cos(Phi);
282 G4double diry = sinTeta * std::sin(Phi);
283 G4double dirz = cosTeta;
284
285 // update G4VParticleChange for the scattered gamma
286 G4ThreeVector gamDirection1(dirx, diry, dirz);
287 gamDirection1.rotateUz(gamDirection0);
288 G4double gamEnergy1 = epsilon * gamEnergy0;
289
290 G4double edep = 0.0;
291 if(gamEnergy1 > lowestSecondaryEnergy)
292 {
295 }
296 else
297 {
300 edep = gamEnergy1;
301 }
302
303 // calculate Stokes vector of final state photon and electron
304 G4ThreeVector nInteractionFrame =
305 G4PolarizationHelper::GetFrame(gamDirection1, gamDirection0);
306
307 // transfer fBeamPolarization and fTargetPolarization
308 // into the interaction frame (note electron is in gamma frame)
309 if(fVerboseLevel >= 1)
310 {
311 G4cout << "========================================" << G4endl;
312 G4cout << " nInteractionFrame = " << nInteractionFrame << G4endl;
313 G4cout << " GammaDirection0 = " << gamDirection0 << G4endl;
314 G4cout << " gammaPolarization = " << fBeamPolarization << G4endl;
315 G4cout << " electronPolarization = " << fTargetPolarization << G4endl;
316 }
317
318 fBeamPolarization.InvRotateAz(nInteractionFrame, gamDirection0);
319 fTargetPolarization.InvRotateAz(nInteractionFrame, gamDirection0);
320
321 if(fVerboseLevel >= 1)
322 {
323 G4cout << "----------------------------------------" << G4endl;
324 G4cout << " gammaPolarization = " << fBeamPolarization << G4endl;
325 G4cout << " electronPolarization = " << fTargetPolarization << G4endl;
326 G4cout << "----------------------------------------" << G4endl;
327 }
328
329 // initialize the polarization transfer matrix
330 fCrossSectionCalculator->Initialize(epsilon, E0_m, 0., fBeamPolarization,
331 fTargetPolarization, 2);
332
333 if(gamEnergy1 > lowestSecondaryEnergy)
334 {
335 // in interaction frame
336 // calculate polarization transfer to the photon (in interaction plane)
337 fFinalGammaPolarization = fCrossSectionCalculator->GetPol2();
338 if(fVerboseLevel >= 1)
339 {
340 G4cout << " gammaPolarization1 = " << fFinalGammaPolarization << G4endl;
341 }
342 fFinalGammaPolarization.SetPhoton();
343
344 // translate polarization into particle reference frame
345 fFinalGammaPolarization.RotateAz(nInteractionFrame, gamDirection1);
346 if(fFinalGammaPolarization.mag() > 1. + 1.e-8)
347 {
349 ed << "ERROR in Polarizaed Compton Scattering !\n";
350 ed << "Polarization of final photon more than 100%.\n";
351 ed << fFinalGammaPolarization
352 << " mag = " << fFinalGammaPolarization.mag() << "\n";
353 G4Exception("G4PolarizedComptonModel::SampleSecondaries", "pol033",
354 FatalException, ed);
355 }
356 // store polarization vector
357 fParticleChange->ProposePolarization(fFinalGammaPolarization);
358 if(fVerboseLevel >= 1)
359 {
360 G4cout << " gammaPolarization1 = " << fFinalGammaPolarization << G4endl;
361 G4cout << " GammaDirection1 = " << gamDirection1 << G4endl;
362 }
363 }
364
365 // kinematic of the scattered electron
366 G4double eKinEnergy = gamEnergy0 - gamEnergy1;
367
368 if(eKinEnergy > lowestSecondaryEnergy)
369 {
370 G4ThreeVector eDirection =
371 gamEnergy0 * gamDirection0 - gamEnergy1 * gamDirection1;
372 eDirection = eDirection.unit();
373
374 finalElectronPolarization = fCrossSectionCalculator->GetPol3();
375 if(fVerboseLevel >= 1)
376 {
377 G4cout << " electronPolarization1 = " << finalElectronPolarization
378 << G4endl;
379 }
380 // transfer into particle reference frame
381 finalElectronPolarization.RotateAz(nInteractionFrame, eDirection);
382 if(fVerboseLevel >= 1)
383 {
384 G4cout << " electronPolarization1 = " << finalElectronPolarization
385 << G4endl << " ElecDirection = " << eDirection << G4endl;
386 }
387
388 // create G4DynamicParticle object for the electron.
389 G4DynamicParticle* aElectron =
390 new G4DynamicParticle(theElectron, eDirection, eKinEnergy);
391 // store polarization vector
392 if(finalElectronPolarization.mag() > 1. + 1.e-8)
393 {
395 ed << "ERROR in Polarized Compton Scattering !\n";
396 ed << "Polarization of final electron more than 100%.\n";
397 ed << finalElectronPolarization
398 << " mag = " << finalElectronPolarization.mag() << G4endl;
399 G4Exception("G4PolarizedComptonModel::SampleSecondaries", "pol034",
400 FatalException, ed);
401 }
402 aElectron->SetPolarization(finalElectronPolarization.p1(),
403 finalElectronPolarization.p2(),
404 finalElectronPolarization.p3());
405 fvect->push_back(aElectron);
406 }
407 else
408 {
409 edep += eKinEnergy;
410 }
411 // energy balance
412 if(edep > 0.0)
413 {
415 }
416}
G4double epsilon(G4double density, G4double temperature)
@ FatalException
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
@ fStopAndKill
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4double alpha2
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
Hep3Vector unit() const
double mag() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
virtual void flatArray(const int size, double *vect)=0
void SetPolarization(const G4ThreeVector &)
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
const G4ThreeVector & GetPolarization() const
G4ParticleChangeForGamma * fParticleChange
G4ParticleDefinition * theElectron
const G4String & GetName() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposePolarization(const G4ThreeVector &dir)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
static G4ThreeVector GetFrame(const G4ThreeVector &, const G4ThreeVector &)
bool IsPolarized(G4LogicalVolume *lVol) const
const G4StokesVector GetVolumePolarization(G4LogicalVolume *lVol) const
static G4PolarizationManager * GetInstance()
G4double p1() const
void InvRotateAz(G4ThreeVector nInteractionFrame, G4ThreeVector particleDirection)
G4double p2() const
void RotateAz(G4ThreeVector nInteractionFrame, G4ThreeVector particleDirection)
G4VPhysicalVolume * GetVolume() const
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
void ProposeTrackStatus(G4TrackStatus status)
const G4Track * GetCurrentTrack() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4LogicalVolume * GetLogicalVolume() const
virtual void Initialize(G4double, G4double, G4double, const G4StokesVector &p0, const G4StokesVector &p1, G4int flag=0)=0
virtual G4StokesVector GetPol2()
virtual G4StokesVector GetPol3()

◆ SetBeamPolarization()

void G4PolarizedComptonModel::SetBeamPolarization ( const G4ThreeVector pBeam)
inline

Definition at line 107 of file G4PolarizedComptonModel.hh.

109{
110 fBeamPolarization = G4StokesVector(pBeam);
111}

◆ SetTargetPolarization()

void G4PolarizedComptonModel::SetTargetPolarization ( const G4ThreeVector pTarget)
inline

Definition at line 102 of file G4PolarizedComptonModel.hh.

104{
105 fTargetPolarization = G4StokesVector(pTarget);
106}

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