source: trunk/source/processes/electromagnetic/lowenergy/src/G4eIonisationCrossSectionHandler.cc @ 924

Last change on this file since 924 was 819, checked in by garnier, 16 years ago

import all except CVS

File size: 5.6 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// $Id: G4eIonisationCrossSectionHandler.cc,v 1.11 2006/06/29 19:42:00 gunter Exp $
27// GEANT4 tag $Name:  $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class file
32//
33//
34// File name:     G4eIonisationCrossSectionHandler
35//
36// Author:        V.Ivanchenko (Vladimir.Ivanchenko@cern.ch)
37//
38// Creation date: 25 Sept 2001
39//
40// Modifications:
41// 10 Oct 2001  M.G. Pia     Revision to improve code quality and consistency with design
42// 19 Jul 2002   VI          Create composite data set for material
43// 21 Jan 2003  V.Ivanchenko Cut per region
44//
45// -------------------------------------------------------------------
46
47#include "G4eIonisationCrossSectionHandler.hh"
48#include "G4VEnergySpectrum.hh"
49#include "G4DataVector.hh"
50#include "G4CompositeEMDataSet.hh"
51#include "G4VDataSetAlgorithm.hh"
52#include "G4LinLogLogInterpolation.hh"
53#include "G4SemiLogInterpolation.hh"
54#include "G4VEMDataSet.hh"
55#include "G4EMDataSet.hh"
56#include "G4Material.hh"
57#include "G4ProductionCutsTable.hh"
58
59
60G4eIonisationCrossSectionHandler::G4eIonisationCrossSectionHandler(
61    const G4VEnergySpectrum* spec, G4VDataSetAlgorithm* alg,
62    G4double emin, G4double emax, G4int nbin)
63 :  G4VCrossSectionHandler(),
64    theParam(spec)
65{
66  G4VCrossSectionHandler::Initialise(alg, emin, emax, nbin);
67  interp = new G4LinLogLogInterpolation();
68}
69
70
71G4eIonisationCrossSectionHandler::~G4eIonisationCrossSectionHandler()
72{
73  delete interp;
74}
75
76
77std::vector<G4VEMDataSet*>* G4eIonisationCrossSectionHandler::BuildCrossSectionsForMaterials(
78                        const G4DataVector& energyVector,
79                        const G4DataVector* energyCuts)
80{
81  G4int verbose = 0;
82  std::vector<G4VEMDataSet*>* set = new std::vector<G4VEMDataSet*>;
83
84  G4DataVector* energies;
85  G4DataVector* cs;
86  G4int nOfBins = energyVector.size();
87
88  const G4ProductionCutsTable* theCoupleTable=
89        G4ProductionCutsTable::GetProductionCutsTable();
90  size_t numOfCouples = theCoupleTable->GetTableSize();
91
92  for (size_t m=0; m<numOfCouples; m++) {
93
94    const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(m);
95    const G4Material* material= couple->GetMaterial();
96    const G4ElementVector* elementVector = material->GetElementVector();
97    const G4double* nAtomsPerVolume = material->GetAtomicNumDensityVector();
98    G4int nElements = material->GetNumberOfElements();
99
100    if(verbose > 0) {
101      G4cout << "eIonisation CS for " << m << "th material "
102             << material->GetName()
103             << "  eEl= " << nElements << G4endl;
104    }
105
106    G4double tcut  = (*energyCuts)[m];
107
108    G4VDataSetAlgorithm* algo = interp->Clone();
109    G4VEMDataSet* setForMat = new G4CompositeEMDataSet(algo,1.,1.);
110
111    for (G4int i=0; i<nElements; i++) {
112
113      G4int Z = (G4int) (*elementVector)[i]->GetZ();
114      G4int nShells = NumberOfComponents(Z);
115      energies = new G4DataVector;
116      cs       = new G4DataVector;
117      G4double density = nAtomsPerVolume[i];
118
119      for (G4int bin=0; bin<nOfBins; bin++) {
120
121        G4double e = energyVector[bin];
122        energies->push_back(e);
123        G4double value = 0.0;
124
125        if(e > tcut) {
126          for (G4int n=0; n<nShells; n++) {
127            G4double cross = FindValue(Z, e, n);
128            G4double p = theParam->Probability(Z, tcut, e, e, n);
129            value += cross * p * density;
130
131            if(verbose>0 && m == 0 && e>=1. && e<=0.) {
132              G4cout << "G4eIonCrossSH: e(MeV)= " << e/MeV
133                     << " n= " << n
134                     << " cross= " << cross
135                     << " p= " << p
136                     << " value= " << value
137                     << " tcut(MeV)= " << tcut/MeV
138                     << " rho= " << density
139                     << " Z= " << Z
140                     << G4endl;
141            }
142
143          }
144        }
145        cs->push_back(value);
146      }
147      G4VDataSetAlgorithm* algo = interp->Clone();
148      G4VEMDataSet* elSet = new G4EMDataSet(i,energies,cs,algo,1.,1.);
149      setForMat->AddComponent(elSet);
150    }
151    set->push_back(setForMat);
152  }
153
154  return set;
155}
156
157
Note: See TracBrowser for help on using the repository browser.