File indexing completed on 2025-01-18 09:59:08
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 #ifndef G4SynchrotronRadiation_h
0037 #define G4SynchrotronRadiation_h 1
0038
0039 #include "globals.hh"
0040 #include "G4Step.hh"
0041 #include "G4ThreeVector.hh"
0042 #include "G4Track.hh"
0043 #include "G4VDiscreteProcess.hh"
0044
0045 class G4LossTableManager;
0046 class G4ParticleDefinition;
0047 class G4PropagatorInField;
0048 class G4VEmAngularDistribution;
0049
0050 class G4SynchrotronRadiation : public G4VDiscreteProcess
0051 {
0052 public:
0053 explicit G4SynchrotronRadiation(const G4String& pName = "SynRad",
0054 G4ProcessType type = fElectromagnetic);
0055
0056 virtual ~G4SynchrotronRadiation();
0057
0058 G4SynchrotronRadiation& operator=(const G4SynchrotronRadiation& right) =
0059 delete;
0060 G4SynchrotronRadiation(const G4SynchrotronRadiation&) = delete;
0061
0062 virtual G4double GetMeanFreePath(const G4Track& track,
0063 G4double previousStepSize,
0064 G4ForceCondition* condition) override;
0065
0066 virtual G4VParticleChange* PostStepDoIt(const G4Track& track,
0067 const G4Step& Step) override;
0068
0069 G4double GetPhotonEnergy(const G4Track& trackData, const G4Step& stepData);
0070
0071 G4double GetRandomEnergySR(G4double, G4double, G4double);
0072
0073 G4double InvSynFracInt(G4double x);
0074 G4double Chebyshev(G4double a, G4double b, const G4double c[], G4int n,
0075 G4double x);
0076
0077 virtual G4bool IsApplicable(const G4ParticleDefinition&) override;
0078 virtual void BuildPhysicsTable(const G4ParticleDefinition&) override;
0079
0080 void ProcessDescription(std::ostream&) const override;
0081 void DumpInfo() const override { ProcessDescription(G4cout); };
0082
0083 void SetAngularGenerator(G4VEmAngularDistribution* p);
0084
0085 private:
0086 G4LossTableManager* theManager;
0087 G4VEmAngularDistribution* genAngle;
0088 G4ParticleDefinition* theGamma;
0089 G4PropagatorInField* fFieldPropagator;
0090
0091 G4bool FirstTime;
0092 G4bool FirstTime1;
0093
0094 G4int secID = -1;
0095 };
0096
0097
0098
0099 inline G4double G4SynchrotronRadiation::Chebyshev(G4double a, G4double b,
0100 const G4double c[], G4int n,
0101 G4double x)
0102 {
0103 G4double y;
0104 G4double y2 = 2.0 * (y = (2.0 * x - a - b) / (b - a));
0105 G4double d = 0., dd = 0.;
0106 for(G4int j = n - 1; j >= 1; --j)
0107 {
0108 G4double sv = d;
0109 d = y2 * d - dd + c[j];
0110 dd = sv;
0111 }
0112 return y * d - dd + 0.5 * c[0];
0113 }
0114
0115 #endif