Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PolarizedComptonModel.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// -------------------------------------------------------------------
27//
28// Geant4 Class header file
29//
30// File name: G4PolarizedComptonModel
31//
32// Author: Andreas Schaelicke
33//
34// Class Description:
35// Implementation of polarized gamma Compton scattering on free electron
36
37#ifndef G4PolarizedComptonModel_h
38#define G4PolarizedComptonModel_h 1
39
40#include "globals.hh"
42#include "G4StokesVector.hh"
43#include "G4ThreeVector.hh"
44
49class G4VPolarizedXS;
50
52{
53 public:
54 explicit G4PolarizedComptonModel(const G4ParticleDefinition* p = nullptr,
55 const G4String& nam = "Polarized-Compton");
56
57 virtual ~G4PolarizedComptonModel() override;
58
60 G4double kinEnergy, G4double Z,
61 G4double A, G4double cut,
62 G4double emax) override;
64
65 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
67 const G4DynamicParticle*, G4double tmin,
68 G4double maxEnergy) override;
69
70 // polarized routines
71 inline void SetTargetPolarization(const G4ThreeVector& pTarget);
72 inline void SetBeamPolarization(const G4ThreeVector& pBeam);
73 inline const G4ThreeVector& GetTargetPolarization() const;
74 inline const G4ThreeVector& GetBeamPolarization() const;
75 inline const G4ThreeVector& GetFinalGammaPolarization() const;
76 inline const G4ThreeVector& GetFinalElectronPolarization() const;
77
78 // hide assignment operator
80 delete;
82
83 private:
84 void PrintWarning(const G4DynamicParticle*, G4int, G4double grej,
85 G4double onecos, G4double phi, const G4String) const;
86
87 static constexpr G4int fLoopLim = 10000;
88
89 G4VPolarizedXS* fCrossSectionCalculator;
90 // incoming
91 G4StokesVector fBeamPolarization; // photon
92 G4StokesVector fTargetPolarization; // electron
93 // outgoing
94 G4StokesVector fFinalGammaPolarization;
95 G4StokesVector finalElectronPolarization;
96
97 G4int fVerboseLevel;
98};
99
100//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
101
103 const G4ThreeVector& pTarget)
104{
105 fTargetPolarization = G4StokesVector(pTarget);
106}
108 const G4ThreeVector& pBeam)
109{
110 fBeamPolarization = G4StokesVector(pBeam);
111}
113 const
114{
115 return fTargetPolarization;
116}
118{
119 return fBeamPolarization;
120}
122 const
123{
124 return fFinalGammaPolarization;
125}
126inline const G4ThreeVector&
128{
129 return finalElectronPolarization;
130}
131
132#endif
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
void SetTargetPolarization(const G4ThreeVector &pTarget)
const G4ThreeVector & GetTargetPolarization() const
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
void SetBeamPolarization(const G4ThreeVector &pBeam)
G4double ComputeAsymmetryPerAtom(G4double gammaEnergy, G4double Z)
const G4ThreeVector & GetFinalGammaPolarization() const
virtual ~G4PolarizedComptonModel() override
G4PolarizedComptonModel & operator=(const G4PolarizedComptonModel &right)=delete
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax) override
const G4ThreeVector & GetFinalElectronPolarization() const
G4PolarizedComptonModel(const G4PolarizedComptonModel &)=delete
const G4ThreeVector & GetBeamPolarization() const