|
||||
File indexing completed on 2025-01-18 09:58:54
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 // Author: Luciano Pandola 0028 // 0029 // History: 0030 // ----------- 0031 // 09 Mar 2012 L. Pandola 1st implementation. 0032 // 0033 // ------------------------------------------------------------------- 0034 // 0035 //! Class description: 0036 //! This class is meant to calculate, store and provide the shell-per-shell 0037 //! and the total ionisation cross section calculated by the Penelope model. 0038 //! The information is provided to other physics models, notably: Penelope 0039 //! ionisation model and Penelope "PIXE" model. Cross sections are calculated 0040 //! per material-cut couple and stored as G4PenelopeCrossSection objects. 0041 //! 0042 // ------------------------------------------------------------------- 0043 0044 #ifndef G4PENELOPEIONISATIONXSHANDLER_HH 0045 #define G4PENELOPEIONISATIONXSHANDLER_HH 1 0046 0047 #include "globals.hh" 0048 #include "G4DataVector.hh" 0049 #include "G4Material.hh" 0050 #include <map> 0051 0052 class G4PhysicsFreeVector; 0053 class G4PhysicsLogVector; 0054 class G4ParticleDefinition; 0055 class G4PenelopeOscillatorManager; 0056 class G4PenelopeOscillator; 0057 class G4PenelopeCrossSection; 0058 0059 class G4PenelopeIonisationXSHandler 0060 { 0061 public: 0062 //! Constructor. nBins is the number of intervals in the 0063 //! energy grid. By default the energy grid goes from 100 eV 0064 //! to 100 GeV. 0065 explicit G4PenelopeIonisationXSHandler(size_t nBins=200); 0066 0067 //! Destructor. Clean all tables. 0068 virtual ~G4PenelopeIonisationXSHandler(); 0069 0070 //!Returns the density coeection for the material at the given energy 0071 G4double GetDensityCorrection(const G4Material*,const G4double energy) const; 0072 //! Returns the table of cross sections for the given particle, given 0073 //! material and given cut as a G4PenelopeCrossSection* pointer. 0074 const G4PenelopeCrossSection* GetCrossSectionTableForCouple(const G4ParticleDefinition*, 0075 const G4Material*,const G4double cut) const; 0076 //!Setter for the verbosity level 0077 void SetVerboseLevel(G4int vl){fVerboseLevel = vl;}; 0078 0079 //! This can be inkoved only by the master 0080 void BuildXSTable(const G4Material*,G4double cut, 0081 const G4ParticleDefinition*,G4bool isMaster=true); 0082 0083 G4PenelopeIonisationXSHandler & operator=(const G4PenelopeIonisationXSHandler &right) = delete; 0084 G4PenelopeIonisationXSHandler(const G4PenelopeIonisationXSHandler&) = delete; 0085 0086 private: 0087 void BuildDeltaTable(const G4Material*); 0088 0089 G4DataVector* ComputeShellCrossSectionsElectron(G4PenelopeOscillator* , 0090 G4double energy,G4double cut, 0091 G4double delta); 0092 0093 G4DataVector* ComputeShellCrossSectionsPositron(G4PenelopeOscillator* , 0094 G4double energy,G4double cut, 0095 G4double delta); 0096 0097 //Oscillator manager 0098 G4PenelopeOscillatorManager* fOscManager; 0099 0100 //G4PenelopeCrossSection takes care of the logs 0101 std::map< std::pair<const G4Material*,G4double>, G4PenelopeCrossSection*> *fXSTableElectron; 0102 std::map< std::pair<const G4Material*,G4double>, G4PenelopeCrossSection*> *fXSTablePositron; 0103 0104 //delta vs. log(energy) 0105 std::map<const G4Material*,G4PhysicsFreeVector*> *fDeltaTable; 0106 0107 //energy grid 0108 G4PhysicsLogVector* fEnergyGrid; 0109 0110 G4int fVerboseLevel; 0111 size_t fNBins; 0112 }; 0113 0114 #endif 0115
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |