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 // Calculation of the total, elastic and inelastic cross-sections 27 // based on Barashenkov parametrisations of pion data 28 // 29 // 16.08.06 V.Ivanchenko - first implementation 30 // 22.01.07 V.Ivanchenko - add cross section interfaces with Z and A 31 // 05.03.07 V.Ivanchenko - add IfZAApplicable 32 // 33 34 #ifndef G4UPiNuclearCrossSection_h 35 #define G4UPiNuclearCrossSection_h 36 37 #include "G4VCrossSectionDataSet.hh" 38 #include "G4DynamicParticle.hh" 39 #include "G4ParticleDefinition.hh" 40 #include "globals.hh" 41 42 class G4PhysicsTable; 43 44 class G4UPiNuclearCrossSection : public G4VCrossSectionDataSet 45 { 46 public: 47 48 explicit G4UPiNuclearCrossSection(); 49 50 ~G4UPiNuclearCrossSection() override = default; 51 52 G4bool IsElementApplicable(const G4DynamicParticle* aParticle, 53 G4int Z, const G4Material*) final; 54 55 inline 56 G4double GetElasticCrossSection(const G4DynamicParticle* aParticle, 57 G4int Z, G4int A) const; 58 59 inline 60 G4double GetInelasticCrossSection(const G4DynamicParticle* aParticle, 61 G4int Z, G4int A) const; 62 63 void BuildPhysicsTable(const G4ParticleDefinition&) final; 64 65 void DumpPhysicsTable(const G4ParticleDefinition&) final; 66 67 void CrossSectionDescription(std::ostream&) const final; 68 69 private: 70 71 G4double Interpolate(G4int Z, G4int A, G4double ekin, 72 const G4PhysicsTable*) const; 73 74 void AddDataSet(const G4String& p, const G4double* tot, 75 const G4double* in, const G4double* e, G4int n); 76 77 void LoadData(); 78 79 const G4ParticleDefinition* piPlus; 80 const G4ParticleDefinition* piMinus; 81 82 static const G4int NZ = 16; 83 static G4int theZ[NZ]; 84 static G4int idxZ[93]; 85 86 static G4double theA[NZ]; 87 static G4double APower[93]; 88 89 static G4PhysicsTable* piPlusElastic; 90 static G4PhysicsTable* piPlusInelastic; 91 static G4PhysicsTable* piMinusElastic; 92 static G4PhysicsTable* piMinusInelastic; 93 94 G4double aPower{0.75}; 95 G4double elow; 96 97 G4bool spline{false}; 98 }; 99 100 inline G4double 101 G4UPiNuclearCrossSection::GetElasticCrossSection( 102 const G4DynamicParticle* dp, G4int Z, G4int A) const 103 { 104 const G4PhysicsTable* table = 105 (dp->GetDefinition() == piPlus) ? piPlusElastic : piMinusElastic; 106 return Interpolate(Z, A, dp->GetKineticEnergy(), table); 107 } 108 109 inline G4double 110 G4UPiNuclearCrossSection::GetInelasticCrossSection( 111 const G4DynamicParticle* dp, G4int Z, G4int A) const 112 { 113 const G4PhysicsTable* table = 114 (dp->GetDefinition() == piPlus) ? piPlusInelastic : piMinusInelastic; 115 return Interpolate(Z, A, dp->GetKineticEnergy(), table); 116 } 117 118 #endif 119