source: trunk/source/processes/electromagnetic/lowenergy/src/G4CrossSectionHandler.cc @ 991

Last change on this file since 991 was 991, checked in by garnier, 15 years ago

update

File size: 4.2 KB
Line 
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// $Id: G4CrossSectionHandler.cc,v 1.18 2006/06/29 19:38:48 gunter Exp $
28// GEANT4 tag $Name: geant4-09-02 $
29//
30// Author: Maria Grazia Pia (Maria.Grazia.Pia@cern.ch)
31//
32// History:
33// -----------
34// 1  Aug 2001   MGP        Created
35// 19 Jul 2002   VI         Create composite data set for material
36// 24 Apr 2003   VI         Cut per region mfpt
37//
38// -------------------------------------------------------------------
39
40#include "G4CrossSectionHandler.hh"
41#include "G4VDataSetAlgorithm.hh"
42#include "G4VEMDataSet.hh"
43#include "G4EMDataSet.hh"
44#include "G4CompositeEMDataSet.hh"
45#include "G4ShellEMDataSet.hh"
46#include "G4ProductionCutsTable.hh"
47#include "G4Material.hh"
48#include "G4Element.hh"
49#include "Randomize.hh"
50#include <map>
51#include <vector>
52
53#include "G4LogLogInterpolation.hh"
54
55G4CrossSectionHandler::G4CrossSectionHandler()
56{ }
57
58G4CrossSectionHandler::~G4CrossSectionHandler()
59{ }
60
61std::vector<G4VEMDataSet*>*
62G4CrossSectionHandler::BuildCrossSectionsForMaterials(const G4DataVector& energyVector,
63                                                      const G4DataVector*)
64{
65  G4DataVector* energies;
66  G4DataVector* data;
67
68  std::vector<G4VEMDataSet*>* matCrossSections = new std::vector<G4VEMDataSet*>;
69
70  const G4ProductionCutsTable* theCoupleTable=
71        G4ProductionCutsTable::GetProductionCutsTable();
72  size_t numOfCouples = theCoupleTable->GetTableSize();
73
74  size_t nOfBins = energyVector.size();
75  const G4VDataSetAlgorithm* interpolationAlgo = CreateInterpolation();
76
77  for (size_t m=0; m<numOfCouples; m++)
78    {
79      const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(m);
80      const G4Material* material= couple->GetMaterial();
81      G4int nElements = material->GetNumberOfElements();
82      const G4ElementVector* elementVector = material->GetElementVector();
83      const G4double* nAtomsPerVolume = material->GetAtomicNumDensityVector();
84
85      G4VDataSetAlgorithm* algo = interpolationAlgo->Clone();
86
87      G4VEMDataSet* setForMat = new G4CompositeEMDataSet(algo,1.,1.);
88
89      for (G4int i=0; i<nElements; i++) {
90 
91        G4int Z = (G4int) (*elementVector)[i]->GetZ();
92        G4double density = nAtomsPerVolume[i];
93
94        energies = new G4DataVector;
95        data = new G4DataVector;
96
97
98        for (size_t bin=0; bin<nOfBins; bin++)
99          {
100            G4double e = energyVector[bin];
101            energies->push_back(e);
102            G4double cross = density*FindValue(Z,e);
103            data->push_back(cross);
104          }
105
106        G4VDataSetAlgorithm* algo1 = interpolationAlgo->Clone();
107        G4VEMDataSet* elSet = new G4EMDataSet(i,energies,data,algo1,1.,1.);
108        setForMat->AddComponent(elSet);
109      }
110
111      matCrossSections->push_back(setForMat);
112    }
113  return matCrossSections;
114}
115
Note: See TracBrowser for help on using the repository browser.