Geant4 9.6.0
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=0, const G4String &nam="Polarized-Compton")
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax)
 
G4double ComputeAsymmetryPerAtom (G4double gammaEnergy, G4double Z)
 
virtual ~G4PolarizedComptonModel ()
 
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 G4ThreeVectorGetTargetPolarization () const
 
const G4ThreeVectorGetBeamPolarization () const
 
const G4ThreeVectorGetFinalGammaPolarization () const
 
const G4ThreeVectorGetFinalElectronPolarization () const
 
- Public Member Functions inherited from G4KleinNishinaCompton
 G4KleinNishinaCompton (const G4ParticleDefinition *p=0, const G4String &nam="Klein-Nishina")
 
virtual ~G4KleinNishinaCompton ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
- 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 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 ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., 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 *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
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)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *)
 
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
 
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 ForceBuildTable (G4bool val)
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 

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 lowestGammaEnergy
 
- Protected Attributes inherited from G4VEmModel
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 

Detailed Description

Definition at line 61 of file G4PolarizedComptonModel.hh.

Constructor & Destructor Documentation

◆ G4PolarizedComptonModel()

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

Definition at line 67 of file G4PolarizedComptonModel.cc.

◆ ~G4PolarizedComptonModel()

G4PolarizedComptonModel::~G4PolarizedComptonModel ( )
virtual

Definition at line 77 of file G4PolarizedComptonModel.cc.

78{
79 if (crossSectionCalculator) delete crossSectionCalculator;
80}

Member Function Documentation

◆ ComputeAsymmetryPerAtom()

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

Definition at line 84 of file G4PolarizedComptonModel.cc.

87{
88 G4double asymmetry = 0.0 ;
89
90 G4double k0 = gammaEnergy / electron_mass_c2 ;
91 G4double k1 = 1 + 2*k0 ;
92
93 asymmetry = -k0;
94 asymmetry *= (k0 + 1.)*sqr(k1)*std::log(k1) - 2.*k0*(5.*sqr(k0) + 4.*k0 + 1.);
95 asymmetry /= ((k0 - 2.)*k0 -2.)*sqr(k1)*std::log(k1) + 2.*k0*(k0*(k0 + 1.)*(k0 + 8.) + 2.);
96
97 // G4cout<<"energy = "<<GammaEnergy<<" asymmetry = "<<asymmetry<<"\t\t GAM = "<<k0<<G4endl;
98 if (asymmetry>1.) G4cout<<"ERROR in G4PolarizedComptonModel::ComputeAsymmetryPerAtom"<<G4endl;
99
100 return asymmetry;
101}
double G4double
Definition: G4Types.hh:64
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
T sqr(const T &x)
Definition: templates.hh:145

Referenced by ComputeCrossSectionPerAtom().

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4KleinNishinaCompton.

Definition at line 104 of file G4PolarizedComptonModel.cc.

111{
112 double xs =
114 Z,A,cut,emax);
115 G4double polzz = theBeamPolarization.p3()*theTargetPolarization.z();
116 if (polzz!=0) {
117 G4double asym=ComputeAsymmetryPerAtom(kinEnergy, Z);
118 xs*=(1.+polzz*asym);
119 }
120 return xs;
121}
double z() const
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax)
G4double ComputeAsymmetryPerAtom(G4double gammaEnergy, G4double Z)
G4double p3() const

◆ GetBeamPolarization()

const G4ThreeVector & G4PolarizedComptonModel::GetBeamPolarization ( ) const
inline

Definition at line 123 of file G4PolarizedComptonModel.hh.

124{
125 return theBeamPolarization;
126}

◆ GetFinalElectronPolarization()

const G4ThreeVector & G4PolarizedComptonModel::GetFinalElectronPolarization ( ) const
inline

Definition at line 131 of file G4PolarizedComptonModel.hh.

132{
133 return finalElectronPolarization;
134}

◆ GetFinalGammaPolarization()

const G4ThreeVector & G4PolarizedComptonModel::GetFinalGammaPolarization ( ) const
inline

Definition at line 127 of file G4PolarizedComptonModel.hh.

128{
129 return finalGammaPolarization;
130}

◆ GetTargetPolarization()

const G4ThreeVector & G4PolarizedComptonModel::GetTargetPolarization ( ) const
inline

Definition at line 119 of file G4PolarizedComptonModel.hh.

120{
121 return theTargetPolarization;
122}

◆ SampleSecondaries()

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

Reimplemented from G4KleinNishinaCompton.

Definition at line 126 of file G4PolarizedComptonModel.cc.

131{
132 const G4Track * aTrack = fParticleChange->GetCurrentTrack();
133 G4VPhysicalVolume* aPVolume = aTrack->GetVolume();
134 G4LogicalVolume* aLVolume = aPVolume->GetLogicalVolume();
135
136 if (verboseLevel>=1)
137 G4cout<<"G4PolarizedComptonModel::SampleSecondaries in "
138 << aLVolume->GetName() <<G4endl;
139
141
142 // obtain polarization of the beam
143 theBeamPolarization = aDynamicGamma->GetPolarization();
144 theBeamPolarization.SetPhoton();
145
146 // obtain polarization of the media
147 const G4bool targetIsPolarized = polarizationManager->IsPolarized(aLVolume);
148 theTargetPolarization = polarizationManager->GetVolumePolarization(aLVolume);
149
150 // if beam is linear polarized or target is transversely polarized
151 // determine the angle to x-axis
152 // (assumes same PRF as in the polarization definition)
153
154 G4ThreeVector gamDirection0 = aDynamicGamma->GetMomentumDirection();
155
156 // transfere theTargetPolarization
157 // into the gamma frame (problem electron is at rest)
158 if (targetIsPolarized)
159 theTargetPolarization.rotateUz(gamDirection0);
160
161 // The scattered gamma energy is sampled according to Klein - Nishina formula.
162 // The random number techniques of Butcher & Messel are used
163 // (Nuc Phys 20(1960),15).
164 // Note : Effects due to binding of atomic electrons are negliged.
165
166 G4double gamEnergy0 = aDynamicGamma->GetKineticEnergy();
167 G4double E0_m = gamEnergy0 / electron_mass_c2 ;
168
169
170 //
171 // sample the energy rate of the scattered gamma
172 //
173
174 G4double epsilon, epsilonsq, onecost, sint2, greject ;
175
176 G4double eps0 = 1./(1. + 2.*E0_m);
177 G4double epsilon0sq = eps0*eps0;
178 G4double alpha1 = - std::log(eps0);
179 G4double alpha2 = 0.5*(1.- epsilon0sq);
180
181 G4double polarization = theBeamPolarization.p3()*theTargetPolarization.p3();
182 do {
183 if ( alpha1/(alpha1+alpha2) > G4UniformRand() ) {
184 epsilon = std::exp(-alpha1*G4UniformRand()); // epsilon0**r
185 epsilonsq = epsilon*epsilon;
186
187 } else {
188 epsilonsq = epsilon0sq + (1.- epsilon0sq)*G4UniformRand();
189 epsilon = std::sqrt(epsilonsq);
190 };
191
192 onecost = (1.- epsilon)/(epsilon*E0_m);
193 sint2 = onecost*(2.-onecost);
194
195
196 G4double gdiced = 2.*(1./epsilon+epsilon);
197 G4double gdist = 1./epsilon + epsilon - sint2
198 - polarization*(1./epsilon-epsilon)*(1.-onecost);
199
200 greject = gdist/gdiced;
201
202 if (greject>1) G4cout<<"ERROR in PolarizedComptonScattering::PostStepDoIt\n"
203 <<" costh rejection does not work properly: "<<greject<<G4endl;
204
205 } while (greject < G4UniformRand());
206
207 //
208 // scattered gamma angles. ( Z - axis along the parent gamma)
209 //
210
211 G4double cosTeta = 1. - onecost;
212 G4double sinTeta = std::sqrt (sint2);
213 G4double Phi;
214 do {
215 Phi = twopi * G4UniformRand();
216 G4double gdiced = 1./epsilon + epsilon - sint2
217 + std::abs(theBeamPolarization.p3())*
218 ( std::abs((1./epsilon-epsilon)*cosTeta*theTargetPolarization.p3())
219 +(1.-epsilon)*sinTeta*(std::sqrt(sqr(theTargetPolarization.p1())
220 + sqr(theTargetPolarization.p2()))))
221 +sint2*(std::sqrt(sqr(theBeamPolarization.p1()) + sqr(theBeamPolarization.p2())));
222
223 G4double gdist = 1./epsilon + epsilon - sint2
224 + theBeamPolarization.p3()*
225 ((1./epsilon-epsilon)*cosTeta*theTargetPolarization.p3()
226 +(1.-epsilon)*sinTeta*(std::cos(Phi)*theTargetPolarization.p1()+
227 std::sin(Phi)*theTargetPolarization.p2()))
228 -sint2*(std::cos(2.*Phi)*theBeamPolarization.p1()
229 +std::sin(2.*Phi)*theBeamPolarization.p2());
230 greject = gdist/gdiced;
231
232 if (greject>1.+1.e-10 || greject<0) G4cout<<"ERROR in PolarizedComptonScattering::PostStepDoIt\n"
233 <<" phi rejection does not work properly: "<<greject<<G4endl;
234
235 if (greject<1.e-3) {
236 G4cout<<"ERROR in PolarizedComptonScattering::PostStepDoIt\n"
237 <<" phi rejection does not work properly: "<<greject<<"\n";
238 G4cout<<" greject="<<greject<<" phi="<<Phi<<" cost="<<cosTeta<<"\n";
239 G4cout<<" gdiced="<<gdiced<<" gdist="<<gdist<<"\n";
240 G4cout<<" eps="<<epsilon<<" 1/eps="<<1./epsilon<<"\n";
241 }
242
243 } while (greject < G4UniformRand());
244 G4double dirx = sinTeta*std::cos(Phi), diry = sinTeta*std::sin(Phi), dirz = cosTeta;
245
246 //
247 // update G4VParticleChange for the scattered gamma
248 //
249
250 G4ThreeVector gamDirection1 ( dirx,diry,dirz );
251 gamDirection1.rotateUz(gamDirection0);
252 G4double gamEnergy1 = epsilon*gamEnergy0;
254
255
256 if(gamEnergy1 > lowestGammaEnergy) {
258 } else {
260 gamEnergy1 += fParticleChange->GetLocalEnergyDeposit();
262 }
263
264 //
265 // kinematic of the scattered electron
266 //
267
268 G4double eKinEnergy = gamEnergy0 - gamEnergy1;
269 G4ThreeVector eDirection = gamEnergy0*gamDirection0 - gamEnergy1*gamDirection1;
270 eDirection = eDirection.unit();
271
272 //
273 // calculate Stokesvector of final state photon and electron
274 //
275 G4ThreeVector nInteractionFrame;
276 if((gamEnergy1 > lowestGammaEnergy) ||
277 (eKinEnergy > DBL_MIN)) {
278
279 // determine interaction plane
280// nInteractionFrame =
281// G4PolarizationHelper::GetFrame(gamDirection1,eDirection);
282 if (gamEnergy1 > lowestGammaEnergy)
283 nInteractionFrame = G4PolarizationHelper::GetFrame(gamDirection1,gamDirection0);
284 else
285 nInteractionFrame = G4PolarizationHelper::GetFrame(gamDirection0, eDirection);
286
287 // transfere theBeamPolarization and theTargetPolarization
288 // into the interaction frame (note electron is in gamma frame)
289 if (verboseLevel>=1) {
290 G4cout << "========================================\n";
291 G4cout << " nInteractionFrame = " <<nInteractionFrame<<"\n";
292 G4cout << " GammaDirection0 = " <<gamDirection0<<"\n";
293 G4cout << " gammaPolarization = " <<theBeamPolarization<<"\n";
294 G4cout << " electronPolarization = " <<theTargetPolarization<<"\n";
295 }
296
297 theBeamPolarization.InvRotateAz(nInteractionFrame,gamDirection0);
298 theTargetPolarization.InvRotateAz(nInteractionFrame,gamDirection0);
299
300 if (verboseLevel>=1) {
301 G4cout << "----------------------------------------\n";
302 G4cout << " gammaPolarization = " <<theBeamPolarization<<"\n";
303 G4cout << " electronPolarization = " <<theTargetPolarization<<"\n";
304 G4cout << "----------------------------------------\n";
305 }
306
307 // initialize the polarization transfer matrix
308 crossSectionCalculator->Initialize(epsilon,E0_m,0.,
309 theBeamPolarization,
310 theTargetPolarization,2);
311 }
312
313 // if(eKinEnergy > DBL_MIN)
314 {
315 // in interaction frame
316 // calculate polarization transfer to the photon (in interaction plane)
317 finalGammaPolarization = crossSectionCalculator->GetPol2();
318 if (verboseLevel>=1) G4cout << " gammaPolarization1 = " <<finalGammaPolarization<<"\n";
319 finalGammaPolarization.SetPhoton();
320
321 // translate polarization into particle reference frame
322 finalGammaPolarization.RotateAz(nInteractionFrame,gamDirection1);
323 //store polarization vector
324 fParticleChange->ProposePolarization(finalGammaPolarization);
325 if (finalGammaPolarization.mag() > 1.+1.e-8){
326 G4cout<<"ERROR in Polarizaed Compton Scattering !"<<G4endl;
327 G4cout<<"Polarization of final photon more than 100%"<<G4endl;
328 G4cout<<finalGammaPolarization<<" mag = "<<finalGammaPolarization.mag()<<G4endl;
329 }
330 if (verboseLevel>=1) {
331 G4cout << " gammaPolarization1 = " <<finalGammaPolarization<<"\n";
332 G4cout << " GammaDirection1 = " <<gamDirection1<<"\n";
333 }
334 }
335
336 // if (ElecKineEnergy > fminimalEnergy) {
337 {
338 finalElectronPolarization = crossSectionCalculator->GetPol3();
339 if (verboseLevel>=1)
340 G4cout << " electronPolarization1 = " <<finalElectronPolarization<<"\n";
341
342 // transfer into particle reference frame
343 finalElectronPolarization.RotateAz(nInteractionFrame,eDirection);
344 if (verboseLevel>=1) {
345 G4cout << " electronPolarization1 = " <<finalElectronPolarization<<"\n";
346 G4cout << " ElecDirection = " <<eDirection<<"\n";
347 }
348 }
349 if (verboseLevel>=1)
350 G4cout << "========================================\n";
351
352
353 if(eKinEnergy > DBL_MIN) {
354
355 // create G4DynamicParticle object for the electron.
356 G4DynamicParticle* aElectron = new G4DynamicParticle(theElectron,eDirection,eKinEnergy);
357 //store polarization vector
358 if (finalElectronPolarization.mag() > 1.+1.e-8){
359 G4cout<<"ERROR in Polarizaed Compton Scattering !"<<G4endl;
360 G4cout<<"Polarization of final electron more than 100%"<<G4endl;
361 G4cout<<finalElectronPolarization<<" mag = "<<finalElectronPolarization.mag()<<G4endl;
362 }
363 aElectron->SetPolarization(finalElectronPolarization.p1(),
364 finalElectronPolarization.p2(),
365 finalElectronPolarization.p3());
366 fvect->push_back(aElectron);
367 }
368}
@ fStopAndKill
bool G4bool
Definition: G4Types.hh:67
#define G4UniformRand()
Definition: Randomize.hh:53
Hep3Vector unit() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
const G4ThreeVector & GetMomentumDirection() const
void SetPolarization(G4double polX, G4double polY, G4double polZ)
G4double GetKineticEnergy() const
const G4ThreeVector & GetPolarization() const
G4ParticleChangeForGamma * fParticleChange
G4ParticleDefinition * theElectron
G4String GetName() const
const G4Track * GetCurrentTrack() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposePolarization(const G4ThreeVector &dir)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
static G4ThreeVector GetFrame(const G4ThreeVector &, const G4ThreeVector &)
bool IsPolarized(G4LogicalVolume *lVol) const
static G4PolarizationManager * GetInstance()
const G4ThreeVector & GetVolumePolarization(G4LogicalVolume *lVol) const
virtual void Initialize(G4double eps, G4double X, G4double phi, const G4StokesVector &p0, const G4StokesVector &p1, G4int flag=0)
G4double p1() const
void InvRotateAz(G4ThreeVector nInteractionFrame, G4ThreeVector particleDirection)
G4double p2() const
void RotateAz(G4ThreeVector nInteractionFrame, G4ThreeVector particleDirection)
G4VPhysicalVolume * GetVolume() const
void ProposeTrackStatus(G4TrackStatus status)
G4double GetLocalEnergyDeposit() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4LogicalVolume * GetLogicalVolume() const
#define DBL_MIN
Definition: templates.hh:75

◆ SetBeamPolarization()

void G4PolarizedComptonModel::SetBeamPolarization ( const G4ThreeVector pBeam)
inline

Definition at line 115 of file G4PolarizedComptonModel.hh.

116{
117 theBeamPolarization = pBeam;
118}

Referenced by G4PolarizedCompton::ComputeAsymmetry().

◆ SetTargetPolarization()

void G4PolarizedComptonModel::SetTargetPolarization ( const G4ThreeVector pTarget)
inline

Definition at line 111 of file G4PolarizedComptonModel.hh.

112{
113 theTargetPolarization = pTarget;
114}

Referenced by G4PolarizedCompton::ComputeAsymmetry().


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