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 // GEANT4 physics class: G4PhotoNuclearCrossSection -- header file 28 // Created: M.V. Kossov, CERN/ITEP(Moscow), 10-OCT-01 29 // The last update: M.V. Kossov, CERN/ITEP (Moscow) 17-May-02 30 31 #ifndef G4PhotoNuclearCrossSection_h 32 #define G4PhotoNuclearCrossSection_h 1 33 34 #include "G4VCrossSectionDataSet.hh" 35 #include "G4DynamicParticle.hh" 36 #include "G4Element.hh" 37 #include "G4ParticleTable.hh" 38 #include "G4NucleiProperties.hh" 39 #include "G4NistManager.hh" 40 #include <vector> 41 42 class G4PhotoNuclearCrossSection : public G4VCrossSectionDataSet 43 { 44 public: 45 46 G4PhotoNuclearCrossSection(); 47 ~G4PhotoNuclearCrossSection() override; 48 49 static const char* Default_Name() {return "PhotoNuclearXS";} 50 51 void CrossSectionDescription(std::ostream&) const override; 52 53 G4bool IsIsoApplicable(const G4DynamicParticle*, G4int Z, G4int A, 54 const G4Element* elm = nullptr, 55 const G4Material* mat = nullptr) override; 56 57 G4bool IsElementApplicable(const G4DynamicParticle*, G4int Z, 58 const G4Material* mat = nullptr) override; 59 60 G4double GetIsoCrossSection(const G4DynamicParticle*, 61 G4int Z, G4int A, 62 const G4Isotope* iso = nullptr, 63 const G4Element* elm = nullptr, 64 const G4Material* mat = nullptr) override; 65 66 G4double GetElementCrossSection(const G4DynamicParticle*, G4int Z, 67 const G4Material*) override; 68 69 G4double ComputeElementXSection(G4double energy, G4int Z); 70 71 G4double ComputeIsoXSection(G4double energy, G4int Z, G4int A); 72 73 G4PhotoNuclearCrossSection& operator= 74 (const G4PhotoNuclearCrossSection& right) = delete; 75 G4PhotoNuclearCrossSection(const G4PhotoNuclearCrossSection&) = delete; 76 77 private: 78 79 G4int GetFunctions(G4double a, G4double* y, G4double* z); 80 G4double EquLinearFit(G4double X, G4int N, const G4double X0, 81 const G4double XD, const G4double* Y); 82 G4double ThresholdEnergy(G4int Z, G4int N); 83 84 G4int lastZ = 0; // The last Z of calculated nucleus 85 G4double lastSig = 0.0; // Last value of the Cross Section 86 G4double* lastGDR = nullptr; // Pointer to the last array of GDR cross sections 87 G4double* lastHEN = nullptr; // Pointer to the last array of HEn cross sections 88 G4double lastE = 0.0; // Last used in the cross section Energy 89 G4double lastTH = 0.0; // Last value of the Energy Threshold (A-dependent) 90 G4double lastSP = 0.0; // Last value of the ShadowingPomeron (A-dependent) 91 92 // Vector of pointers to the GDRPhotonuclearCrossSection 93 std::vector <G4double*> GDR; 94 95 // store deuteron, triton, He3 XS 96 G4double* deuteron_GDR = nullptr; 97 G4double* deuteron_HR = nullptr; 98 G4double deuteron_TH = 0.0; 99 G4double deuteron_SP = 0.0; 100 G4double* triton_GDR = nullptr; 101 G4double* triton_HR = nullptr; 102 G4double triton_TH = 0.0; 103 G4double triton_SP = 0.0; 104 G4double* he3_GDR = nullptr; 105 G4double* he3_HR = nullptr; 106 G4double he3_TH = 0.0; 107 G4double he3_SP = 0.0; 108 109 // Vector of pointers to the HighEnPhotonuclearCrossSect 110 std::vector <G4double*> HEN; 111 112 std::vector <G4double> spA; // shadowing coefficients (A-dependent) 113 std::vector <G4double> eTH; // energy threshold (A-dependent) 114 115 G4NistManager* nistmngr; 116 117 G4double mNeut; 118 G4double mProt; 119 }; 120 121 #endif 122