File indexing completed on 2025-01-18 09:58:56
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035
0036
0037
0038
0039 #ifndef G4PolarizedAnnihilationModel_h
0040 #define G4PolarizedAnnihilationModel_h 1
0041
0042 #include "globals.hh"
0043 #include "G4eeToTwoGammaModel.hh"
0044 #include "G4ThreeVector.hh"
0045 #include "G4StokesVector.hh"
0046
0047 class G4DynamicParticle;
0048 class G4MaterialCutsCouple;
0049 class G4ParticleChangeForGamma;
0050 class G4ParticleDefinition;
0051 class G4PolarizedAnnihilationXS;
0052
0053 class G4PolarizedAnnihilationModel : public G4eeToTwoGammaModel
0054 {
0055 public:
0056 explicit G4PolarizedAnnihilationModel(
0057 const G4ParticleDefinition* p = nullptr,
0058 const G4String& nam = "Polarized-Annihilation");
0059
0060 virtual ~G4PolarizedAnnihilationModel() override;
0061
0062 virtual void Initialise(const G4ParticleDefinition*,
0063 const G4DataVector&) override;
0064
0065 virtual G4double ComputeCrossSectionPerElectron(G4double kinEnergy) override;
0066
0067 void ComputeAsymmetriesPerElectron(G4double gammaEnergy, G4double& valueX,
0068 G4double& valueA, G4double& valueT);
0069
0070 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
0071 const G4MaterialCutsCouple*,
0072 const G4DynamicParticle*, G4double tmin,
0073 G4double maxEnergy) final;
0074
0075
0076 inline void SetTargetPolarization(const G4ThreeVector& pTarget);
0077 inline void SetBeamPolarization(const G4ThreeVector& pBeam);
0078 inline const G4ThreeVector& GetTargetPolarization() const;
0079 inline const G4ThreeVector& GetBeamPolarization() const;
0080 inline const G4ThreeVector& GetFinalGamma1Polarization() const;
0081 inline const G4ThreeVector& GetFinalGamma2Polarization() const;
0082
0083 G4PolarizedAnnihilationModel& operator =(
0084 const G4PolarizedAnnihilationModel& right) = delete;
0085 G4PolarizedAnnihilationModel(const G4PolarizedAnnihilationModel&) = delete;
0086
0087 private:
0088 G4PolarizedAnnihilationXS* fCrossSectionCalculator;
0089 G4ParticleChangeForGamma* fParticleChange;
0090
0091
0092 G4StokesVector fBeamPolarization;
0093 G4StokesVector fTargetPolarization;
0094
0095 G4StokesVector fFinalGamma1Polarization;
0096 G4StokesVector fFinalGamma2Polarization;
0097
0098 G4int fVerboseLevel;
0099 };
0100
0101
0102
0103 inline void G4PolarizedAnnihilationModel::SetTargetPolarization(
0104 const G4ThreeVector& pTarget)
0105 {
0106 fTargetPolarization = G4StokesVector(pTarget);
0107 }
0108 inline void G4PolarizedAnnihilationModel::SetBeamPolarization(
0109 const G4ThreeVector& pBeam)
0110 {
0111 fBeamPolarization = G4StokesVector(pBeam);
0112 }
0113 inline const G4ThreeVector&
0114 G4PolarizedAnnihilationModel::GetTargetPolarization() const
0115 {
0116 return fTargetPolarization;
0117 }
0118 inline const G4ThreeVector& G4PolarizedAnnihilationModel::GetBeamPolarization()
0119 const
0120 {
0121 return fBeamPolarization;
0122 }
0123 inline const G4ThreeVector&
0124 G4PolarizedAnnihilationModel::GetFinalGamma1Polarization() const
0125 {
0126 return fFinalGamma1Polarization;
0127 }
0128 inline const G4ThreeVector&
0129 G4PolarizedAnnihilationModel::GetFinalGamma2Polarization() const
0130 {
0131 return fFinalGamma2Polarization;
0132 }
0133
0134 #endif