Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/hadronic/cross_sections/include/G4NeutronElasticXS.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 //
 27 // -------------------------------------------------------------------
 28 //
 29 // GEANT4 Class header file
 30 //
 31 //
 32 // File name:    G4NeutronElasticXS
 33 //
 34 // Author  Ivantchenko, Geant4, 3-AUG-09
 35 //
 36  
 37 // Class Description:
 38 // This is a base class for neutron elastic hadronic cross section based on
 39 // data files from G4PARTICLEXSDATA data set 
 40 // Class Description - End
 41 
 42 #ifndef G4NeutronElasticXS_h
 43 #define G4NeutronElasticXS_h 1
 44 
 45 #include "G4VCrossSectionDataSet.hh"
 46 #include "globals.hh"
 47 #include "G4PhysicsVector.hh"
 48 #include <vector>
 49 
 50 class G4DynamicParticle;
 51 class G4ParticleDefinition;
 52 class G4Element;
 53 class G4VComponentCrossSection;
 54 
 55 class G4NeutronElasticXS final : public G4VCrossSectionDataSet
 56 {
 57 public: 
 58 
 59   G4NeutronElasticXS();
 60 
 61   ~G4NeutronElasticXS() final;
 62     
 63   static const char* Default_Name() {return "G4NeutronElasticXS";}
 64 
 65   G4bool IsElementApplicable(const G4DynamicParticle*, 
 66            G4int Z, const G4Material*) final;
 67 
 68   G4bool IsIsoApplicable(const G4DynamicParticle*, G4int Z, G4int A,
 69        const G4Element*, const G4Material*) final;
 70 
 71   G4double GetElementCrossSection(const G4DynamicParticle*, 
 72                 G4int Z, const G4Material*) final; 
 73 
 74   G4double GetIsoCrossSection(const G4DynamicParticle*, G4int Z, G4int A,
 75                               const G4Isotope* iso,
 76                               const G4Element* elm,
 77                               const G4Material* mat) final;
 78 
 79   G4double ComputeCrossSectionPerElement(G4double kinEnergy, G4double loge,
 80                                          const G4ParticleDefinition*,
 81                                          const G4Element*,
 82                                          const G4Material*) final;
 83   
 84   G4double ComputeIsoCrossSection(G4double kinEnergy, G4double loge,
 85                                   const G4ParticleDefinition*,
 86                                   G4int Z, G4int A,
 87                                   const G4Isotope* iso,
 88                                   const G4Element* elm,
 89                                   const G4Material* mat) final;
 90 
 91   const G4Isotope* SelectIsotope(const G4Element*, 
 92                                  G4double kinEnergy, G4double logE) final;
 93 
 94   void BuildPhysicsTable(const G4ParticleDefinition&) final;
 95 
 96   void CrossSectionDescription(std::ostream&) const final;
 97 
 98   G4double ElementCrossSection(G4double kinEnergy, G4double loge, G4int Z);
 99 
100   G4NeutronElasticXS & operator=(const G4NeutronElasticXS &right) = delete;
101   G4NeutronElasticXS(const G4NeutronElasticXS&) = delete;
102   
103 private: 
104 
105   void Initialise(G4int Z);
106 
107   void InitialiseOnFly(G4int Z);
108 
109   const G4String& FindDirectoryPath();
110 
111   inline G4PhysicsVector* GetPhysicsVector(G4int Z);
112 
113   G4VComponentCrossSection* ggXsection = nullptr;
114   const G4ParticleDefinition* neutron;
115 
116   G4bool isFirst = false;
117 
118   static const G4int MAXZEL = 93;
119   static G4PhysicsVector* data[MAXZEL];
120   static G4double coeff[MAXZEL];
121   static G4String gDataDirectory;
122   static G4bool fLock;
123 };
124 
125 inline
126 G4PhysicsVector* G4NeutronElasticXS::GetPhysicsVector(G4int Z)
127 {
128   if(nullptr == data[Z]) { InitialiseOnFly(Z); }
129   return data[Z];
130 }
131 
132 #endif
133