File indexing completed on 2025-01-18 09:58:38
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 #ifndef G4LivermorePolarizedGammaConversionModel_h
0031 #define G4LivermorePolarizedGammaConversionModel_h 1
0032
0033 #include "G4VEmModel.hh"
0034 #include "G4PhysicsFreeVector.hh"
0035
0036 class G4ParticleChangeForGamma;
0037
0038 class G4LivermorePolarizedGammaConversionModel : public G4VEmModel
0039 {
0040 public:
0041
0042 explicit G4LivermorePolarizedGammaConversionModel(const G4ParticleDefinition* p = nullptr,
0043 const G4String& nam = "LivermorePolarizedGammaConversion");
0044 virtual ~G4LivermorePolarizedGammaConversionModel();
0045
0046 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
0047
0048 void InitialiseLocal(const G4ParticleDefinition*,
0049 G4VEmModel* masterModel) override;
0050
0051 void InitialiseForElement(const G4ParticleDefinition*, G4int Z) override;
0052
0053 G4double ComputeCrossSectionPerAtom(
0054 const G4ParticleDefinition*,
0055 G4double kinEnergy,
0056 G4double Z,
0057 G4double A=0,
0058 G4double cut=0,
0059 G4double emax=DBL_MAX) override;
0060
0061 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
0062 const G4MaterialCutsCouple*,
0063 const G4DynamicParticle*,
0064 G4double tmin,
0065 G4double maxEnergy) override;
0066
0067 G4double MinPrimaryEnergy(const G4Material*,
0068 const G4ParticleDefinition*,
0069 G4double) override;
0070
0071 G4LivermorePolarizedGammaConversionModel & operator=(const G4LivermorePolarizedGammaConversionModel &right) = delete;
0072 G4LivermorePolarizedGammaConversionModel(const G4LivermorePolarizedGammaConversionModel&) = delete;
0073
0074 private:
0075 void ReadData(size_t Z, const char* path = 0);
0076
0077 G4double ScreenFunction1(G4double screenVariable);
0078 G4double ScreenFunction2(G4double screenVariable);
0079
0080
0081 G4ThreeVector GetRandomPolarization(G4ThreeVector& direction0);
0082 G4ThreeVector GetPerpendicularPolarization(const G4ThreeVector& direction0, const G4ThreeVector& polarization0) const;
0083 G4ThreeVector SetPerpendicularVector(G4ThreeVector& a);
0084 void SystemOfRefChange(G4ThreeVector& direction0, G4ThreeVector& direction1,
0085 G4ThreeVector& polarization0);
0086
0087
0088 G4double SetPhi(G4double);
0089 G4double SetPsi(G4double, G4double);
0090
0091 void SetTheta(G4double*, G4double*, G4double);
0092
0093 G4double Poli(G4double , G4double, G4double, G4double);
0094 G4double Fln(G4double, G4double, G4double);
0095
0096 G4double Encu(G4double*, G4double*, G4double);
0097
0098 G4double Flor(G4double*, G4double);
0099 G4double Glor(G4double*, G4double);
0100
0101 G4double Fdlor(G4double*, G4double);
0102 G4double Fintlor(G4double*, G4double);
0103 G4double Finvlor(G4double*, G4double,G4double);
0104
0105 G4double Ftan(G4double*, G4double);
0106 G4double Fdtan(G4double*, G4double);
0107 G4double Finttan(G4double*, G4double);
0108 G4double Finvtan(G4double*, G4double, G4double);
0109
0110
0111 G4ParticleChangeForGamma* fParticleChange;
0112
0113 static const G4int maxZ = 99;
0114 static G4PhysicsFreeVector* data[100];
0115
0116 G4double lowEnergyLimit;
0117 G4double smallEnergy;
0118 G4double Psi, Phi;
0119 G4int verboseLevel;
0120 G4bool isInitialised;
0121 };
0122
0123
0124
0125 #endif