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 #ifndef G4PolarizedComptonModel_h
0038 #define G4PolarizedComptonModel_h 1
0039
0040 #include "globals.hh"
0041 #include "G4KleinNishinaCompton.hh"
0042 #include "G4StokesVector.hh"
0043 #include "G4ThreeVector.hh"
0044
0045 class G4DynamicParticle;
0046 class G4MaterialCutsCouple;
0047 class G4ParticleDefinition;
0048 class G4ParticleChangeForGamma;
0049 class G4VPolarizedXS;
0050
0051 class G4PolarizedComptonModel : public G4KleinNishinaCompton
0052 {
0053 public:
0054 explicit G4PolarizedComptonModel(const G4ParticleDefinition* p = nullptr,
0055 const G4String& nam = "Polarized-Compton");
0056
0057 virtual ~G4PolarizedComptonModel() override;
0058
0059 virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition*,
0060 G4double kinEnergy, G4double Z,
0061 G4double A, G4double cut,
0062 G4double emax) override;
0063 G4double ComputeAsymmetryPerAtom(G4double gammaEnergy, G4double Z);
0064
0065 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
0066 const G4MaterialCutsCouple*,
0067 const G4DynamicParticle*, G4double tmin,
0068 G4double maxEnergy) override;
0069
0070
0071 inline void SetTargetPolarization(const G4ThreeVector& pTarget);
0072 inline void SetBeamPolarization(const G4ThreeVector& pBeam);
0073 inline const G4ThreeVector& GetTargetPolarization() const;
0074 inline const G4ThreeVector& GetBeamPolarization() const;
0075 inline const G4ThreeVector& GetFinalGammaPolarization() const;
0076 inline const G4ThreeVector& GetFinalElectronPolarization() const;
0077
0078
0079 G4PolarizedComptonModel& operator=(const G4PolarizedComptonModel& right) =
0080 delete;
0081 G4PolarizedComptonModel(const G4PolarizedComptonModel&) = delete;
0082
0083 private:
0084 void PrintWarning(const G4DynamicParticle*, G4int, G4double grej,
0085 G4double onecos, G4double phi, const G4String) const;
0086
0087 static constexpr G4int fLoopLim = 10000;
0088
0089 G4VPolarizedXS* fCrossSectionCalculator;
0090
0091 G4StokesVector fBeamPolarization;
0092 G4StokesVector fTargetPolarization;
0093
0094 G4StokesVector fFinalGammaPolarization;
0095 G4StokesVector finalElectronPolarization;
0096
0097 G4int fVerboseLevel;
0098 };
0099
0100
0101
0102 inline void G4PolarizedComptonModel::SetTargetPolarization(
0103 const G4ThreeVector& pTarget)
0104 {
0105 fTargetPolarization = G4StokesVector(pTarget);
0106 }
0107 inline void G4PolarizedComptonModel::SetBeamPolarization(
0108 const G4ThreeVector& pBeam)
0109 {
0110 fBeamPolarization = G4StokesVector(pBeam);
0111 }
0112 inline const G4ThreeVector& G4PolarizedComptonModel::GetTargetPolarization()
0113 const
0114 {
0115 return fTargetPolarization;
0116 }
0117 inline const G4ThreeVector& G4PolarizedComptonModel::GetBeamPolarization() const
0118 {
0119 return fBeamPolarization;
0120 }
0121 inline const G4ThreeVector& G4PolarizedComptonModel::GetFinalGammaPolarization()
0122 const
0123 {
0124 return fFinalGammaPolarization;
0125 }
0126 inline const G4ThreeVector&
0127 G4PolarizedComptonModel::GetFinalElectronPolarization() const
0128 {
0129 return finalElectronPolarization;
0130 }
0131
0132 #endif