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
0033
0034
0035
0036
0037
0038
0039
0040
0041
0042
0043
0044
0045
0046
0047
0048
0049
0050
0051
0052
0053
0054
0055
0056
0057
0058
0059
0060
0061
0062
0063
0064
0065
0066
0067
0068 #ifndef G4LowEPPolarizedComptonModel_h
0069 #define G4LowEPPolarizedComptonModel_h 1
0070
0071 #include "G4VEmModel.hh"
0072 #include "G4PhysicsFreeVector.hh"
0073 #include <limits>
0074 #include "G4Electron.hh"
0075 #include "G4ParticleChangeForGamma.hh"
0076 #include "G4LossTableManager.hh"
0077 #include "G4VAtomDeexcitation.hh"
0078 #include "G4AtomicShell.hh"
0079 #include "G4Gamma.hh"
0080 #include "G4ShellData.hh"
0081 #include "G4DopplerProfile.hh"
0082 #include "G4Log.hh"
0083 #include "G4Exp.hh"
0084 #include "G4ForceCondition.hh"
0085
0086 class G4ParticleChangeForGamma;
0087 class G4VAtomDeexcitation;
0088 class G4ShellData;
0089 class G4DopplerProfile;
0090
0091 class G4LowEPPolarizedComptonModel : public G4VEmModel
0092 {
0093 public:
0094 explicit G4LowEPPolarizedComptonModel(const G4ParticleDefinition* p = nullptr,
0095 const G4String& nam = "LowEPComptonModel");
0096 virtual ~G4LowEPPolarizedComptonModel();
0097
0098 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
0099
0100 void InitialiseLocal(const G4ParticleDefinition*,
0101 G4VEmModel* masterModel) override;
0102
0103 void InitialiseForElement(const G4ParticleDefinition*, G4int Z) override;
0104
0105 G4double ComputeCrossSectionPerAtom( const G4ParticleDefinition*,
0106 G4double kinEnergy,
0107 G4double Z,
0108 G4double A=0,
0109 G4double cut=0,
0110 G4double emax=DBL_MAX ) override;
0111
0112 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
0113 const G4MaterialCutsCouple*,
0114 const G4DynamicParticle*,
0115 G4double tmin,
0116 G4double maxEnergy) override;
0117
0118 G4LowEPPolarizedComptonModel & operator=(const G4LowEPPolarizedComptonModel &right) = delete;
0119 G4LowEPPolarizedComptonModel(const G4LowEPPolarizedComptonModel&) = delete;
0120
0121 private:
0122 void ReadData(size_t Z, const char* path = 0);
0123
0124 G4double ComputeScatteringFunction(G4double x, G4int Z);
0125 G4ThreeVector GetRandomPolarization(G4ThreeVector& direction0);
0126 G4ThreeVector GetPerpendicularPolarization(const G4ThreeVector& direction0, const G4ThreeVector& polarization0) const;
0127 G4ThreeVector SetNewPolarization(G4double LowEPPCepsilon, G4double sinT2,
0128 G4double phi, G4double cosTheta);
0129 G4double SetPhi(G4double, G4double);
0130
0131 void SystemOfRefChange(G4ThreeVector& direction0, G4ThreeVector& direction1,
0132 G4ThreeVector& polarization0, G4ThreeVector& polarization1);
0133
0134 void SystemOfRefChangeElect(G4ThreeVector& pdirection, G4ThreeVector& edirection,
0135 G4ThreeVector& ppolarization);
0136
0137 G4ThreeVector SetPerpendicularVector(G4ThreeVector& a);
0138
0139 G4ParticleChangeForGamma* fParticleChange;
0140 G4VAtomDeexcitation* fAtomDeexcitation;
0141
0142 static G4ShellData* shellData;
0143 static G4DopplerProfile* profileData;
0144
0145 static const G4int maxZ = 99;
0146 static G4PhysicsFreeVector* data[100];
0147 static const G4double ScatFuncFitParam[101][9];
0148
0149 G4int verboseLevel;
0150 G4bool isInitialised;
0151 };
0152
0153
0154
0155 #endif