Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/hadronic/cross_sections/include/G4NucleonNuclearCrossSection.hh

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 // author: Vladimir.Grichine@cern.ch
 27 //
 28 // Implements data from: Barashenkov V.S., Nucleon-Nucleus Cross Section,
 29 // Preprint JINR P2-89-770, p. 12, Dubna 1989 (scanned version from KEK)
 30 // Based on G. Folger version of G4PiNuclearCrossSection class
 31 //
 32 // Modified:
 33 // 05.03.2024 V.Ivanchenko removed obsolete methods and calls
 34 //
 35 
 36 #ifndef G4NucleonNuclearCrossSection_h
 37 #define G4NucleonNuclearCrossSection_h
 38 
 39 #include "G4VCrossSectionDataSet.hh"
 40 #include "G4ParticleDefinition.hh"
 41 #include "globals.hh"
 42 
 43 class G4ComponentBarNucleonNucleusXsc;
 44 class G4ParticleDefinition;
 45 
 46 class G4NucleonNuclearCrossSection : public G4VCrossSectionDataSet
 47 {
 48 public:
 49   
 50   G4NucleonNuclearCrossSection();
 51   ~G4NucleonNuclearCrossSection() override = default;
 52     
 53   static const char* Default_Name() { return "BarashenkovNucleonXS"; }
 54 
 55   G4bool IsElementApplicable(const G4DynamicParticle* aParticle,
 56            G4int Z, const G4Material* mat) final;
 57 
 58   // return inelastic x-section
 59   G4double GetElementCrossSection(const G4DynamicParticle* aParticle, 
 60                                   G4int Z, const G4Material* mat=nullptr) final;
 61 
 62   void CrossSectionDescription(std::ostream&) const final;
 63 
 64   // return elastic x-section
 65   inline G4double GetElasticCrossSection(const G4DynamicParticle* aParticle,
 66            G4int Z);
 67 
 68   // access methods should be called after ComputeCrossSection(...)
 69   inline G4double GetTotalXsc()     { return fTotalXsc; };
 70   inline G4double GetInelasticXsc() { return fInelasticXsc; };
 71   inline G4double GetElasticXsc()   { return fElasticXsc; };
 72 
 73   G4NucleonNuclearCrossSection& operator=
 74   (const G4NucleonNuclearCrossSection &right) = delete;
 75   G4NucleonNuclearCrossSection(const G4NucleonNuclearCrossSection&) = delete;
 76 
 77 private:
 78 
 79   void ComputeCrossSections(const G4ParticleDefinition*, 
 80                             G4double kinEnergy, G4int Z);
 81 
 82   G4ComponentBarNucleonNucleusXsc* fBarash;
 83 
 84   G4double fTotalXsc{0.0};
 85   G4double fInelasticXsc{0.0};
 86   G4double fElasticXsc{0.0};
 87 };
 88 
 89 inline
 90 G4double G4NucleonNuclearCrossSection::GetElasticCrossSection(
 91          const G4DynamicParticle* dp, G4int Z)
 92 {
 93   ComputeCrossSections(dp->GetDefinition(), dp->GetKineticEnergy(), Z);
 94   return fElasticXsc;
 95 }
 96 
 97 #endif
 98