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
0031
0032 #ifndef G4LivermorePolarizedRayleighModel_h
0033 #define G4LivermorePolarizedRayleighModel_h 1
0034
0035 #include "G4VEmModel.hh"
0036 #include "G4ParticleChangeForGamma.hh"
0037 #include "G4PhysicsFreeVector.hh"
0038 #include "G4ProductionCutsTable.hh"
0039
0040 class G4LivermorePolarizedRayleighModel : public G4VEmModel
0041 {
0042 public:
0043
0044 explicit G4LivermorePolarizedRayleighModel(
0045 const G4ParticleDefinition* p = nullptr,
0046 const G4String& nam = "LivermorePolarizedRayleigh");
0047
0048 virtual ~G4LivermorePolarizedRayleighModel();
0049
0050 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
0051 void InitialiseLocal(const G4ParticleDefinition*,
0052 G4VEmModel* masterModel) override;
0053 void InitialiseForElement(const G4ParticleDefinition*, G4int Z) override;
0054
0055 G4double ComputeCrossSectionPerAtom(
0056 const G4ParticleDefinition*,
0057 G4double kinEnergy,
0058 G4double Z,
0059 G4double A=0,
0060 G4double cut=0,
0061 G4double emax=DBL_MAX) override;
0062
0063 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
0064 const G4MaterialCutsCouple*,
0065 const G4DynamicParticle*,
0066 G4double tmin,
0067 G4double maxEnergy) override;
0068
0069 G4LivermorePolarizedRayleighModel & operator=
0070 (const G4LivermorePolarizedRayleighModel &right) = delete;
0071 G4LivermorePolarizedRayleighModel(const G4LivermorePolarizedRayleighModel&) = delete;
0072
0073 private:
0074
0075 void ReadData(size_t Z, const char* path = nullptr);
0076
0077
0078
0079
0080
0081
0082 G4double GenerateCosTheta(G4double incomingPhotonEnergy, G4int zAtom) const;
0083
0084
0085
0086
0087 G4double GeneratePhi(G4double cosTheta) const;
0088
0089
0090
0091 G4double GeneratePolarizationAngle(void) const;
0092 G4ThreeVector GetPhotonPolarization(const G4DynamicParticle& photon);
0093
0094 G4ParticleChangeForGamma* fParticleChange;
0095
0096 const G4int maxZ = 100;
0097 static G4PhysicsFreeVector* dataCS[101];
0098 static G4PhysicsFreeVector* formFactorData[101];
0099
0100 G4double lowEnergyLimit;
0101 G4int verboseLevel;
0102 G4bool isInitialised;
0103 };
0104
0105 #endif