File indexing completed on 2025-01-18 09:59:28
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 #ifndef G4VXTRenergyLoss_h
0040 #define G4VXTRenergyLoss_h 1
0041
0042 #include "globals.hh"
0043 #include "G4Gamma.hh"
0044 #include "G4LogicalVolume.hh"
0045 #include "G4Material.hh"
0046 #include "G4ParticleChange.hh"
0047 #include "G4PhysicsTable.hh"
0048 #include "G4Step.hh"
0049 #include "G4Track.hh"
0050 #include "G4VDiscreteProcess.hh"
0051
0052 class G4SandiaTable;
0053 class G4VParticleChange;
0054 class G4PhysicsFreeVector;
0055 class G4PhysicsLinearVector;
0056 class G4PhysicsLogVector;
0057
0058 class G4VXTRenergyLoss : public G4VDiscreteProcess
0059 {
0060 public:
0061 explicit G4VXTRenergyLoss(G4LogicalVolume* anEnvelope, G4Material*,
0062 G4Material*, G4double, G4double, G4int,
0063 const G4String& processName = "XTRenergyLoss",
0064 G4ProcessType type = fElectromagnetic);
0065 virtual ~G4VXTRenergyLoss();
0066
0067 virtual void ProcessDescription(std::ostream&) const override;
0068 virtual void DumpInfo() const override { ProcessDescription(G4cout); };
0069
0070 G4VXTRenergyLoss(G4VXTRenergyLoss&) = delete;
0071 G4VXTRenergyLoss& operator=(const G4VXTRenergyLoss& right) = delete;
0072
0073
0074 virtual G4double GetStackFactor(G4double energy, G4double gamma,
0075 G4double varAngle);
0076
0077 virtual G4bool IsApplicable(const G4ParticleDefinition&) override;
0078
0079 virtual G4VParticleChange* PostStepDoIt(const G4Track& aTrack,
0080 const G4Step& aStep) override;
0081
0082 virtual G4double GetMeanFreePath(const G4Track& aTrack,
0083 G4double previousStepSize,
0084 G4ForceCondition* condition) override;
0085
0086 virtual void BuildPhysicsTable(const G4ParticleDefinition&) override;
0087 void BuildEnergyTable();
0088 void BuildAngleForEnergyBank();
0089
0090 void BuildTable(){};
0091 void BuildAngleTable();
0092 void BuildGlobalAngleTable();
0093
0094 G4complex OneInterfaceXTRdEdx(G4double energy, G4double gamma,
0095 G4double varAngle);
0096
0097 G4double SpectralAngleXTRdEdx(G4double varAngle);
0098
0099 virtual G4double SpectralXTRdEdx(G4double energy);
0100
0101 G4double AngleSpectralXTRdEdx(G4double energy);
0102
0103 G4double AngleXTRdEdx(G4double varAngle);
0104
0105 G4double OneBoundaryXTRNdensity(G4double energy, G4double gamma,
0106 G4double varAngle) const;
0107
0108
0109 G4double XTRNSpectralAngleDensity(G4double varAngle);
0110 G4double XTRNSpectralDensity(G4double energy);
0111
0112
0113 G4double XTRNAngleSpectralDensity(G4double energy);
0114 G4double XTRNAngleDensity(G4double varAngle);
0115
0116 void GetNumberOfPhotons();
0117
0118
0119 G4double GetPlateFormationZone(G4double, G4double, G4double);
0120 G4complex GetPlateComplexFZ(G4double, G4double, G4double);
0121 void ComputePlatePhotoAbsCof();
0122 G4double GetPlateLinearPhotoAbs(G4double);
0123 void GetPlateZmuProduct();
0124 G4double GetPlateZmuProduct(G4double, G4double, G4double);
0125
0126 G4double GetGasFormationZone(G4double, G4double, G4double);
0127 G4complex GetGasComplexFZ(G4double, G4double, G4double);
0128 void ComputeGasPhotoAbsCof();
0129 G4double GetGasLinearPhotoAbs(G4double);
0130 void GetGasZmuProduct();
0131 G4double GetGasZmuProduct(G4double, G4double, G4double);
0132
0133 G4double GetPlateCompton(G4double);
0134 G4double GetGasCompton(G4double);
0135 G4double GetComptonPerAtom(G4double, G4double);
0136
0137 G4double GetXTRrandomEnergy(G4double scaledTkin, G4int iTkin);
0138 G4double GetXTRenergy(G4int iPlace, G4double position, G4int iTransfer);
0139
0140 G4double GetRandomAngle(G4double energyXTR, G4int iTkin);
0141 G4double GetAngleXTR(G4int iTR, G4double position, G4int iAngle);
0142
0143
0144
0145 void SetGamma(G4double gamma) { fGamma = gamma; };
0146 G4double GetGamma() { return fGamma; };
0147 void SetEnergy(G4double energy) { fEnergy = energy; };
0148 G4double GetEnergy() { return fEnergy; };
0149 void SetVarAngle(G4double varAngle) { fVarAngle = varAngle; };
0150 G4double GetVarAngle() { return fVarAngle; };
0151 void SetCompton(G4bool pC) { fCompton = pC; };
0152 G4bool GetCompton() { return fCompton; };
0153
0154 G4int GetKrange(){ return fKrange;};
0155 void SetKrange( G4int kk ){ fKrange = kk;};
0156
0157
0158 void SetAlphaGas(G4double ag){ fAlphaGas = ag;};
0159 G4double GetAlphaGas() { return fAlphaGas; };
0160 void SetAlphaPlate(G4double ap){ fAlphaPlate = ap;};
0161 G4double GetAlphaPlate() { return fAlphaPlate; };
0162
0163 void SetTheMinEnergyTR(G4double minetr){ fTheMinEnergyTR = minetr;};
0164 G4double GetTheMinEnergyTR() { return fTheMinEnergyTR; };
0165 void SetTheMaxEnergyTR(G4double maxetr){ fTheMaxEnergyTR = maxetr;};
0166 G4double GetTheMaxEnergyTR() { return fTheMaxEnergyTR; };
0167
0168 void SetMinEnergyTR(G4double minetr){ fMinEnergyTR = minetr;};
0169 G4double GetMinEnergyTR() { return fMinEnergyTR; };
0170 void SetMaxEnergyTR(G4double maxetr){ fMaxEnergyTR = maxetr;};
0171 G4double GetMaxEnergyTR() { return fMaxEnergyTR; };
0172
0173 void SetTheMinAngle(G4double minang){ fTheMinAngle = minang;};
0174 G4double GetTheMinAngle() { return fTheMinAngle; };
0175 void SetTheMaxAngle(G4double maxang){ fTheMaxAngle = maxang;};
0176 G4double GetTheMaxAngle() { return fTheMaxAngle; };
0177
0178 void SetMinThetaTR(G4double minatr){ fMinThetaTR = minatr;};
0179 G4double GetMinThetaTR() { return fMinThetaTR; };
0180 void SetMaxThetaTR(G4double maxatr){ fMaxThetaTR = maxatr;};
0181 G4double GetMaxThetaTR() { return fMaxThetaTR; };
0182
0183
0184
0185 void SetFastAngle(G4bool fatr){ fFastAngle = fatr;};
0186 G4bool GetFastAngle() { return fFastAngle; };
0187 void SetAngleRadDistr(G4bool fatr){ fAngleRadDistr = fatr;};
0188 G4bool GetAngleRadDistr() { return fAngleRadDistr; };
0189
0190
0191
0192
0193 G4PhysicsLogVector* GetProtonVector() { return fProtonEnergyVector; };
0194 G4int GetTotBin() { return fTotBin; };
0195 G4PhysicsFreeVector* GetAngleVector(G4double energy, G4int n);
0196
0197 protected:
0198
0199 G4double fTheMinEnergyTR;
0200
0201 G4double fTheMaxEnergyTR;
0202 G4double fTheMinAngle;
0203 G4double fTheMaxAngle;
0204
0205
0206
0207
0208 static constexpr G4double fMinProtonTkin = 100. * CLHEP::GeV;
0209
0210 static constexpr G4double fMaxProtonTkin = 100. * CLHEP::TeV;
0211
0212 static constexpr G4double fPlasmaCof =
0213 4. * CLHEP::pi * CLHEP::fine_structure_const * CLHEP::hbarc * CLHEP::hbarc *
0214 CLHEP::hbarc / CLHEP::electron_mass_c2;
0215 static constexpr G4double fCofTR = CLHEP::fine_structure_const / CLHEP::pi;
0216
0217 G4int fTotBin;
0218 G4int fBinTR;
0219 G4int fKrange;
0220 G4ParticleDefinition* fPtrGamma;
0221
0222 G4double* fGammaCutInKineticEnergy;
0223 G4LogicalVolume* fEnvelope;
0224 G4PhysicsTable* fAngleDistrTable;
0225 G4PhysicsTable* fEnergyDistrTable;
0226 G4PhysicsTable* fAngleForEnergyTable;
0227 G4PhysicsLogVector* fProtonEnergyVector;
0228 G4PhysicsLogVector* fXTREnergyVector;
0229 G4SandiaTable* fPlatePhotoAbsCof;
0230 G4SandiaTable* fGasPhotoAbsCof;
0231
0232 G4ParticleChange fParticleChange;
0233 std::vector<G4PhysicsTable*> fAngleBank;
0234
0235 G4double fGammaTkinCut;
0236 G4double fMinEnergyTR;
0237 G4double fMaxEnergyTR;
0238 G4double fMinThetaTR, fMaxThetaTR;
0239 G4double fTotalDist;
0240 G4double fPlateThick;
0241 G4double fGasThick;
0242 G4double fAlphaPlate;
0243 G4double fAlphaGas;
0244 G4double fGamma;
0245 G4double fEnergy;
0246 G4double fVarAngle;
0247 G4double fLambda;
0248 G4double fSigma1;
0249 G4double fSigma2;
0250
0251 G4int fMatIndex1;
0252 G4int fMatIndex2;
0253 G4int fPlateNumber;
0254
0255 G4bool fExitFlux;
0256 G4bool fFastAngle, fAngleRadDistr;
0257 G4bool fCompton;
0258
0259 G4int secID = -1;
0260 };
0261
0262 #endif