Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/pii/include/G4PixeCrossSectionHandler.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 // Author: Maria Grazia Pia (Maria.Grazia.Pia@cern.ch)
 29 //
 30 // History:
 31 // -----------
 32 // 16 Jun 2008   MGP           Created on the basis of G4CrossSectionHandler
 33 
 34 // -------------------------------------------------------------------
 35 // Class description:
 36 // Cross section manager for hadron impact ionization
 37 // Documented in:
 38 // M.G. Pia et al., PIXE Simulation With Geant4,
 39 // IEEE Trans. Nucl. Sci., vol. 56, no. 6, pp. 3614-3649, Dec. 2009.
 40 
 41 // -------------------------------------------------------------------
 42 
 43 #ifndef G4PIXECROSSSECTIONHANDLER_HH
 44 #define G4PIXECROSSSECTIONHANDLER_HH 1
 45 
 46 #include <map>
 47 #include <vector>
 48 #include <CLHEP/Units/SystemOfUnits.h>
 49 
 50 #include "globals.hh"
 51 #include "G4DataVector.hh"
 52 #include "G4MaterialCutsCouple.hh"
 53 
 54 class G4IInterpolator;
 55 class G4IDataSet;
 56 class G4Material;
 57 class G4Element;
 58 
 59 class G4PixeCrossSectionHandler {
 60 
 61 public:
 62 
 63   G4PixeCrossSectionHandler();
 64 
 65   G4PixeCrossSectionHandler(G4IInterpolator* interpolation,
 66           const G4String& modelK="ecpssr",
 67           const G4String& modelL="ecpssr",
 68           const G4String& modelM="ecpssr",
 69           G4double minE = 1*CLHEP::keV,
 70                             G4double maxE = 0.1*CLHEP::GeV,
 71           G4int nBins = 200,
 72           G4double unitE = CLHEP::MeV,
 73                             G4double unitData = CLHEP::barn,
 74           G4int minZ = 6, G4int maxZ = 92);
 75   
 76   virtual ~G4PixeCrossSectionHandler();
 77 
 78   void Initialise(G4IInterpolator* interpolation,
 79       const G4String& modelK="ecpssr",
 80       const G4String& modelL="ecpssr",
 81       const G4String& modelM="ecpssr",
 82       G4double minE = 1*CLHEP::keV,
 83                   G4double maxE = 0.1*CLHEP::GeV,
 84       G4int nBins = 200,
 85       G4double unitE = CLHEP::MeV,
 86                   G4double unitData = CLHEP::barn,
 87       G4int minZ = 6, G4int maxZ = 92);
 88 
 89   G4int SelectRandomAtom(const G4Material* material, G4double e) const;
 90 
 91   G4int SelectRandomShell(G4int Z, G4double e) const;
 92 
 93   G4double FindValue(G4int Z, G4double e) const;
 94 
 95   G4double FindValue(G4int Z, G4double e, G4int shellIndex) const;
 96 
 97   G4double ValueForMaterial(const G4Material* material, G4double e) const;
 98 
 99   // void LoadData(const G4String& dataFile);
100 
101   void LoadShellData(const G4String& dataFile);
102 
103   // Ionisation cross section as in Geant4 Physics Reference Manual
104   G4double MicroscopicCrossSection(const G4ParticleDefinition* particleDef,
105            G4double kineticEnergy,
106            G4double Z,
107            G4double deltaCut) const;
108 
109   void PrintData() const;
110 
111   void Clear();
112 
113 private:
114 
115  // Hide copy constructor and assignment operator
116   G4PixeCrossSectionHandler(const G4PixeCrossSectionHandler&);
117   G4PixeCrossSectionHandler & operator=(const G4PixeCrossSectionHandler &right);
118 
119   G4int NumberOfComponents(G4int Z) const;
120 
121   void ActiveElements();
122 
123   void BuildForMaterials();
124   // Factory method
125   std::vector<G4IDataSet*>* BuildCrossSectionsForMaterials(const G4DataVector& energyVector);
126 
127   // Factory method
128   G4IInterpolator* CreateInterpolation();
129 
130   const G4IInterpolator* GetInterpolation() const { return interpolation; }
131  
132   G4IInterpolator* interpolation;
133 
134   G4double eMin;
135   G4double eMax;
136   G4int nBins;
137 
138   G4double unit1;
139   G4double unit2;
140 
141   G4int zMin;
142   G4int zMax;
143 
144   G4DataVector activeZ;
145 
146   // Map of PixeShellDataSets with the shell cross sections for each element
147   std::map<G4int,G4IDataSet*,std::less<G4int> > dataMap;
148 
149   // Vector of composite cross sections for each material in the MaterialTable 
150   // The composite cross section is composed of cross sections for each element in material
151   std::vector<G4IDataSet*>* crossSections;
152 
153   std::vector<G4String> crossModel;
154 
155 };
156 
157 #endif
158 
159 
160 
161 
162 
163 
164 
165 
166 
167 
168 
169