File indexing completed on 2025-01-18 09:58:54
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 #ifndef G4PENELOPEBREMSSTRAHLUNGMODEL_HH
0043 #define G4PENELOPEBREMSSTRAHLUNGMODEL_HH 1
0044
0045 #include "globals.hh"
0046 #include "G4VEmModel.hh"
0047 #include "G4DataVector.hh"
0048 #include "G4ParticleChangeForLoss.hh"
0049
0050
0051 class G4PhysicsFreeVector;
0052 class G4PhysicsLogVector;
0053 class G4ParticleDefinition;
0054 class G4DynamicParticle;
0055 class G4MaterialCutsCouple;
0056 class G4Material;
0057 class G4PenelopeOscillatorManager;
0058 class G4PenelopeCrossSection;
0059 class G4PenelopeBremsstrahlungFS;
0060 class G4PenelopeBremsstrahlungAngular;
0061
0062 class G4PenelopeBremsstrahlungModel : public G4VEmModel
0063 {
0064 public:
0065 explicit G4PenelopeBremsstrahlungModel(const G4ParticleDefinition* p=nullptr,
0066 const G4String& processName ="PenBrem");
0067 virtual ~G4PenelopeBremsstrahlungModel();
0068
0069 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
0070 void InitialiseLocal(const G4ParticleDefinition*,
0071 G4VEmModel*) override;
0072
0073
0074 G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition* theParticle,
0075 G4double kinEnergy,
0076 G4double Z,
0077 G4double A=0,
0078 G4double cut=0,
0079 G4double emax=DBL_MAX) override;
0080
0081
0082 G4double CrossSectionPerVolume(const G4Material* material,
0083 const G4ParticleDefinition* theParticle,
0084 G4double kineticEnergy,
0085 G4double cutEnergy,
0086 G4double maxEnergy = DBL_MAX) override;
0087
0088 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
0089 const G4MaterialCutsCouple*,
0090 const G4DynamicParticle*,
0091 G4double tmin,
0092 G4double maxEnergy) override;
0093
0094 G4double ComputeDEDXPerVolume(const G4Material*,
0095 const G4ParticleDefinition*,
0096 G4double kineticEnergy,
0097 G4double cutEnergy) override;
0098
0099
0100 G4double MinEnergyCut(const G4ParticleDefinition*,
0101 const G4MaterialCutsCouple*) override;
0102
0103 void SetVerbosityLevel(G4int lev){fVerboseLevel = lev;};
0104 G4int GetVerbosityLevel(){return fVerboseLevel;};
0105
0106 G4PenelopeBremsstrahlungModel & operator=
0107 (const G4PenelopeBremsstrahlungModel &right) = delete;
0108 G4PenelopeBremsstrahlungModel(const G4PenelopeBremsstrahlungModel&) = delete;
0109
0110 protected:
0111 G4ParticleChangeForLoss* fParticleChange;
0112 const G4ParticleDefinition* fParticle;
0113
0114 private:
0115 void ClearTables();
0116 G4PenelopeCrossSection* GetCrossSectionTableForCouple(const G4ParticleDefinition*,
0117 const G4Material*,G4double cut);
0118 void SetParticle(const G4ParticleDefinition*);
0119
0120
0121 void BuildXSTable(const G4Material* material,G4double cut);
0122 G4double GetPositronXSCorrection(const G4Material*,G4double energy);
0123
0124
0125 G4PenelopeOscillatorManager* fOscManager;
0126 G4PenelopeBremsstrahlungFS* fPenelopeFSHelper;
0127 G4PenelopeBremsstrahlungAngular* fPenelopeAngular;
0128
0129
0130 G4PhysicsLogVector* fEnergyGrid;
0131 size_t nBins;
0132
0133 std::map< std::pair<const G4Material*,G4double>, G4PenelopeCrossSection*> *fXSTableElectron;
0134 std::map< std::pair<const G4Material*,G4double>, G4PenelopeCrossSection*> *fXSTablePositron;
0135
0136
0137 G4double fIntrinsicLowEnergyLimit;
0138 G4double fIntrinsicHighEnergyLimit;
0139
0140 G4int fVerboseLevel;
0141 G4bool fIsInitialised;
0142
0143 G4bool fLocalTable;
0144 };
0145
0146 #endif
0147