Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:58:04

0001 //
0002 // ********************************************************************
0003 // * License and Disclaimer                                           *
0004 // *                                                                  *
0005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
0006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
0007 // * conditions of the Geant4 Software License,  included in the file *
0008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
0009 // * include a list of copyright holders.                             *
0010 // *                                                                  *
0011 // * Neither the authors of this software system, nor their employing *
0012 // * institutes,nor the agencies providing financial support for this *
0013 // * work  make  any representation or  warranty, express or implied, *
0014 // * regarding  this  software system or assume any liability for its *
0015 // * use.  Please see the license in the file  LICENSE  and URL above *
0016 // * for the full disclaimer and the limitation of liability.         *
0017 // *                                                                  *
0018 // * This  code  implementation is the result of  the  scientific and *
0019 // * technical work of the GEANT4 collaboration.                      *
0020 // * By using,  copying,  modifying or  distributing the software (or *
0021 // * any work based  on the software)  you  agree  to acknowledge its *
0022 // * use  in  resulting  scientific  publications,  and indicate your *
0023 // * acceptance of all terms of the Geant4 Software license.          *
0024 // ********************************************************************
0025 //
0026 //
0027 // -------------------------------------------------------------------
0028 // File name:     G4CrossSectionDataStore
0029 //
0030 // Modifications:
0031 // 23.01.2009 V.Ivanchenko move constructor and destructor to source,
0032 //                         use STL vector instead of C-array
0033 //
0034 // August 2011  Re-designed
0035 //              by G. Folger, V. Ivantchenko, T. Koi and D.H. Wright
0036 
0037 // Class Description
0038 // This is the class to which cross section data sets may be registered. 
0039 // An instance of it is contained in each hadronic process, allowing
0040 // the use of the AddDataSet() method to tailor the cross sections to
0041 // your application.
0042 // Class Description - End
0043 
0044 #ifndef G4CrossSectionDataStore_h
0045 #define G4CrossSectionDataStore_h 1
0046 
0047 #include "globals.hh"
0048 #include "G4VCrossSectionDataSet.hh"
0049 #include "G4DynamicParticle.hh"
0050 #include "G4PhysicsVector.hh"
0051 #include <vector>
0052 #include <iostream>
0053 
0054 class G4Nucleus;
0055 class G4ParticleDefinition;
0056 class G4Isotope;
0057 class G4Element;
0058 class G4Material;
0059 class G4NistManager;
0060 
0061 class G4CrossSectionDataStore
0062 {
0063 public:
0064 
0065   G4CrossSectionDataStore();
0066 
0067   ~G4CrossSectionDataStore() = default;
0068 
0069   // Run time cross section per unit volume
0070   inline G4double GetCrossSection(const G4DynamicParticle*, const G4Material*);
0071   G4double ComputeCrossSection(const G4DynamicParticle*, const G4Material*);
0072 
0073   // Cross section per element
0074   G4double GetCrossSection(const G4DynamicParticle*, 
0075                const G4Element*, const G4Material*);
0076 
0077   // Cross section per isotope 
0078   G4double GetCrossSection(const G4DynamicParticle*, G4int Z, G4int A,
0079                            const G4Isotope*,
0080                const G4Element*, const G4Material*);
0081 
0082   // Sample Z and A of a target nucleus and upload into G4Nucleus
0083   const G4Element* SampleZandA(const G4DynamicParticle*, const G4Material*,
0084                    G4Nucleus& target);
0085 
0086   // Initialisation before run
0087   void BuildPhysicsTable(const G4ParticleDefinition&);
0088 
0089   // Dump store to G4cout
0090   void DumpPhysicsTable(const G4ParticleDefinition&);
0091 
0092   // Dump store as html
0093   void DumpHtml(const G4ParticleDefinition&, std::ofstream&) const;
0094   void PrintCrossSectionHtml(const G4VCrossSectionDataSet *cs,
0095                              const G4String&, const G4String&) const;
0096   
0097   void AddDataSet(G4VCrossSectionDataSet*);
0098   void AddDataSet(G4VCrossSectionDataSet*, std::size_t);
0099   inline const std::vector<G4VCrossSectionDataSet*>& GetDataSetList() const;
0100 
0101   inline void SetVerboseLevel(G4int value);
0102 
0103   // may be used by special processes
0104   inline void SetForcedElement(const G4Element*);
0105 
0106   G4CrossSectionDataStore & operator=
0107   (const G4CrossSectionDataStore &right) = delete;
0108   G4CrossSectionDataStore(const G4CrossSectionDataStore&) = delete;
0109 
0110 private:
0111 
0112   G4double GetIsoCrossSection(const G4DynamicParticle*, G4int Z, G4int A,
0113                   const G4Isotope*, const G4Element*,
0114                               const G4Material*, const G4int index);
0115 
0116   G4String HtmlFileName(const G4String & in) const;
0117 
0118   G4NistManager* nist;
0119   const G4Material* currentMaterial = nullptr;
0120   const G4ParticleDefinition* matParticle = nullptr;
0121   const G4Element* forcedElement = nullptr;
0122   G4double matKinEnergy = 0.0;
0123   G4double matCrossSection = 0.0;
0124 
0125   G4int nDataSetList = 0;
0126   G4int verboseLevel = 1;
0127 
0128   std::vector<G4VCrossSectionDataSet*> dataSetList;
0129   std::vector<G4double> xsecelm;
0130   std::vector<G4double> xseciso;
0131 };
0132 
0133 inline void G4CrossSectionDataStore::SetVerboseLevel(G4int value)
0134 {
0135   verboseLevel = value;
0136 }
0137 
0138 inline void G4CrossSectionDataStore::SetForcedElement(const G4Element* ptr)
0139 {
0140   forcedElement = ptr;
0141 }
0142 
0143 inline const std::vector<G4VCrossSectionDataSet*>&
0144 G4CrossSectionDataStore::GetDataSetList() const
0145 {
0146   return dataSetList;
0147 }
0148 
0149 inline G4double 
0150 G4CrossSectionDataStore::GetCrossSection(const G4DynamicParticle* dp,
0151                                          const G4Material* mat)
0152 {
0153   if(dp->GetKineticEnergy() != matKinEnergy || mat != currentMaterial ||
0154      dp->GetDefinition() != matParticle) {
0155     ComputeCrossSection(dp, mat);
0156   }
0157   return matCrossSection;
0158 }
0159 
0160 #endif