Geant4 Cross Reference |
1 // 2 // ******************************************************************** 3 // * License and Disclaimer * 4 // * * 5 // * The Geant4 software is copyright of the Copyright Holders of * 6 // * the Geant4 Collaboration. It is provided under the terms and * 7 // * conditions of the Geant4 Software License, included in the file * 8 // * LICENSE and available at http://cern.ch/geant4/license . These * 9 // * include a list of copyright holders. * 10 // * * 11 // * Neither the authors of this software system, nor their employing * 12 // * institutes,nor the agencies providing financial support for this * 13 // * work make any representation or warranty, express or implied, * 14 // * regarding this software system or assume any liability for its * 15 // * use. Please see the license in the file LICENSE and URL above * 16 // * for the full disclaimer and the limitation of liability. * 17 // * * 18 // * This code implementation is the result of the scientific and * 19 // * technical work of the GEANT4 collaboration. * 20 // * By using, copying, modifying or distributing the software (or * 21 // * any work based on the software) you agree to acknowledge its * 22 // * use in resulting scientific publications, and indicate your * 23 // * acceptance of all terms of the Geant4 Software license. * 24 // ******************************************************************** 25 // 26 // ------------------------------------------------------------ 27 // GEANT 4 class header file 28 // 29 // History: 30 // 21-5-98 1 version , V. Grichine 31 // 28-05-01, V.Ivanchenko minor changes to provide ANSI -wall compilation 32 // 23-05-06, H. Burkhardt: Energy spectrum from function rather than table 33 // 34 // ------------------------------------------------------------ 35 36 #ifndef G4SynchrotronRadiation_h 37 #define G4SynchrotronRadiation_h 1 38 39 #include "globals.hh" 40 #include "G4Step.hh" 41 #include "G4ThreeVector.hh" 42 #include "G4Track.hh" 43 #include "G4VDiscreteProcess.hh" 44 45 class G4LossTableManager; 46 class G4ParticleDefinition; 47 class G4PropagatorInField; 48 class G4VEmAngularDistribution; 49 50 class G4SynchrotronRadiation : public G4VDiscreteProcess 51 { 52 public: 53 explicit G4SynchrotronRadiation(const G4String& pName = "SynRad", 54 G4ProcessType type = fElectromagnetic); 55 56 virtual ~G4SynchrotronRadiation(); 57 58 G4SynchrotronRadiation& operator=(const G4SynchrotronRadiation& right) = 59 delete; 60 G4SynchrotronRadiation(const G4SynchrotronRadiation&) = delete; 61 62 virtual G4double GetMeanFreePath(const G4Track& track, 63 G4double previousStepSize, 64 G4ForceCondition* condition) override; 65 66 virtual G4VParticleChange* PostStepDoIt(const G4Track& track, 67 const G4Step& Step) override; 68 69 G4double GetPhotonEnergy(const G4Track& trackData, const G4Step& stepData); 70 71 G4double GetRandomEnergySR(G4double, G4double, G4double); 72 73 G4double InvSynFracInt(G4double x); 74 G4double Chebyshev(G4double a, G4double b, const G4double c[], G4int n, 75 G4double x); 76 77 virtual G4bool IsApplicable(const G4ParticleDefinition&) override; 78 virtual void BuildPhysicsTable(const G4ParticleDefinition&) override; 79 80 void ProcessDescription(std::ostream&) const override; 81 void DumpInfo() const override { ProcessDescription(G4cout); }; 82 83 void SetAngularGenerator(G4VEmAngularDistribution* p); 84 85 private: 86 G4LossTableManager* theManager; 87 G4VEmAngularDistribution* genAngle; 88 G4ParticleDefinition* theGamma; 89 G4PropagatorInField* fFieldPropagator; 90 91 G4bool FirstTime; 92 G4bool FirstTime1; 93 94 G4int secID = -1; // creator modelID 95 }; 96 97 ////////////////////////// INLINE METHODS ///////////////////////////// 98 99 inline G4double G4SynchrotronRadiation::Chebyshev(G4double a, G4double b, 100 const G4double c[], G4int n, 101 G4double x) 102 { 103 G4double y; 104 G4double y2 = 2.0 * (y = (2.0 * x - a - b) / (b - a)); // Change of variable. 105 G4double d = 0., dd = 0.; 106 for(G4int j = n - 1; j >= 1; --j) // Clenshaw's recurrence. 107 { 108 G4double sv = d; 109 d = y2 * d - dd + c[j]; 110 dd = sv; 111 } 112 return y * d - dd + 0.5 * c[0]; 113 } 114 115 #endif // end of G4SynchrotronRadiation.hh 116