Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/lowenergy/src/G4CrossSectionHandler.cc

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 // 1  Aug 2001   MGP        Created
 33 // 19 Jul 2002   VI         Create composite data set for material
 34 // 24 Apr 2003   VI         Cut per region mfpt
 35 //
 36 // 15 Jul 2009   Nicolas A. Karakatsanis
 37 //
 38 //                           - BuildCrossSectionForMaterials method was revised in order to calculate the 
 39 //                             logarithmic values of the loaded data. 
 40 //                             It retrieves the data values from the G4EMLOW data files but, then, calculates the
 41 //                             respective log values and loads them to seperate data structures.
 42 //                             The EM data sets, initialized this way, contain both non-log and log values.
 43 //                             These initialized data sets can enhance the computing performance of data interpolation
 44 //                             operations
 45 //
 46 // -------------------------------------------------------------------
 47 
 48 #include "G4CrossSectionHandler.hh"
 49 #include "G4VDataSetAlgorithm.hh"
 50 #include "G4VEMDataSet.hh"
 51 #include "G4EMDataSet.hh"
 52 #include "G4CompositeEMDataSet.hh"
 53 #include "G4ShellEMDataSet.hh"
 54 #include "G4ProductionCutsTable.hh"
 55 #include "G4Material.hh"
 56 #include "G4Element.hh"
 57 #include "Randomize.hh"
 58 #include <map>
 59 #include <vector>
 60 
 61 #include "G4LogLogInterpolation.hh"
 62 
 63 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 64 
 65 G4CrossSectionHandler::G4CrossSectionHandler()
 66 { }
 67 
 68 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 69 
 70 G4CrossSectionHandler::~G4CrossSectionHandler()
 71 { }
 72 
 73 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 74 
 75 std::vector<G4VEMDataSet*>*
 76 G4CrossSectionHandler::BuildCrossSectionsForMaterials(const G4DataVector& energyVector,
 77                   const G4DataVector*)
 78 {
 79   G4DataVector* energies;
 80   G4DataVector* data;
 81 
 82   G4DataVector* log_energies;
 83   G4DataVector* log_data;
 84 
 85   std::vector<G4VEMDataSet*>* matCrossSections = new std::vector<G4VEMDataSet*>;
 86 
 87   const G4ProductionCutsTable* theCoupleTable=
 88         G4ProductionCutsTable::GetProductionCutsTable();
 89   G4int numOfCouples = (G4int)theCoupleTable->GetTableSize();
 90 
 91   std::size_t nOfBins = energyVector.size();
 92   const G4VDataSetAlgorithm* interpolationAlgo = CreateInterpolation();
 93 
 94   for (G4int mLocal=0; mLocal<numOfCouples; ++mLocal)
 95     {
 96       const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(mLocal);
 97       const G4Material* material= couple->GetMaterial();
 98       G4int nElements = (G4int)material->GetNumberOfElements();
 99       const G4ElementVector* elementVector = material->GetElementVector();
100       const G4double* nAtomsPerVolume = material->GetAtomicNumDensityVector();
101 
102       G4VDataSetAlgorithm* algo = interpolationAlgo->Clone();
103 
104       G4VEMDataSet* setForMat = new G4CompositeEMDataSet(algo,1.,1.);
105 
106       for (G4int i=0; i<nElements; ++i) {
107  
108         G4int Z = (G4int) (*elementVector)[i]->GetZ();
109         G4double density = nAtomsPerVolume[i];
110 
111         energies = new G4DataVector;
112         data = new G4DataVector;
113 
114         log_energies = new G4DataVector;
115         log_data = new G4DataVector;
116 
117         for (std::size_t bin=0; bin<nOfBins; ++bin)
118     {
119       G4double e = energyVector[bin];
120       energies->push_back(e);
121             if (e==0.) e=1e-300;
122             log_energies->push_back(std::log10(e));
123       G4double cross = density*FindValue(Z,e);
124       data->push_back(cross);
125             if (cross==0.) cross=1e-300;
126             log_data->push_back(std::log10(cross));
127     }
128 
129         G4VDataSetAlgorithm* algo1 = interpolationAlgo->Clone();
130         G4VEMDataSet* elSet = new G4EMDataSet(i,energies,data,log_energies,log_data,algo1,1.,1.);
131 
132         setForMat->AddComponent(elSet);
133       }
134 
135       matCrossSections->push_back(setForMat);
136     }
137   delete interpolationAlgo;
138   return matCrossSections;
139 }
140 
141