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

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

maj sur la beta de geant 4.9.3

File size: 6.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// $Id: G4eIonisationCrossSectionHandler.cc,v 1.12 2009/01/29 08:13:34 pandola Exp $
27// GEANT4 tag $Name: geant4-09-03-beta-cand-01 $
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// 28 Jan 2009  L.Pandola    Added public method to make a easier migration of
45//                           G4LowEnergyIonisation to G4LivermoreIonisationModel
46//
47// -------------------------------------------------------------------
48
49#include "G4eIonisationCrossSectionHandler.hh"
50#include "G4VEnergySpectrum.hh"
51#include "G4DataVector.hh"
52#include "G4CompositeEMDataSet.hh"
53#include "G4VDataSetAlgorithm.hh"
54#include "G4LinLogLogInterpolation.hh"
55#include "G4SemiLogInterpolation.hh"
56#include "G4VEMDataSet.hh"
57#include "G4EMDataSet.hh"
58#include "G4Material.hh"
59#include "G4ProductionCutsTable.hh"
60
61
62G4eIonisationCrossSectionHandler::G4eIonisationCrossSectionHandler(
63    const G4VEnergySpectrum* spec, G4VDataSetAlgorithm* alg,
64    G4double emin, G4double emax, G4int nbin)
65 :  G4VCrossSectionHandler(),
66    theParam(spec)
67{
68  G4VCrossSectionHandler::Initialise(alg, emin, emax, nbin);
69  interp = new G4LinLogLogInterpolation();
70}
71
72
73G4eIonisationCrossSectionHandler::~G4eIonisationCrossSectionHandler()
74{
75  delete interp;
76}
77
78
79std::vector<G4VEMDataSet*>* G4eIonisationCrossSectionHandler::BuildCrossSectionsForMaterials(
80                        const G4DataVector& energyVector,
81                        const G4DataVector* energyCuts)
82{
83  G4int verbose = 0;
84  std::vector<G4VEMDataSet*>* set = new std::vector<G4VEMDataSet*>;
85
86  G4DataVector* energies;
87  G4DataVector* cs;
88  G4int nOfBins = energyVector.size();
89
90  const G4ProductionCutsTable* theCoupleTable=
91        G4ProductionCutsTable::GetProductionCutsTable();
92  size_t numOfCouples = theCoupleTable->GetTableSize();
93
94  for (size_t m=0; m<numOfCouples; m++) {
95
96    const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(m);
97    const G4Material* material= couple->GetMaterial();
98    const G4ElementVector* elementVector = material->GetElementVector();
99    const G4double* nAtomsPerVolume = material->GetAtomicNumDensityVector();
100    G4int nElements = material->GetNumberOfElements();
101
102    if(verbose > 0) 
103      {
104        G4cout << "eIonisation CS for " << m << "th material "
105               << material->GetName()
106               << "  eEl= " << nElements << G4endl;
107      }
108
109    G4double tcut  = (*energyCuts)[m];
110
111    G4VDataSetAlgorithm* algo = interp->Clone();
112    G4VEMDataSet* setForMat = new G4CompositeEMDataSet(algo,1.,1.);
113
114    for (G4int i=0; i<nElements; i++) {
115
116      G4int Z = (G4int) (*elementVector)[i]->GetZ();
117      G4int nShells = NumberOfComponents(Z);
118      energies = new G4DataVector;
119      cs       = new G4DataVector;
120      G4double density = nAtomsPerVolume[i];
121
122      for (G4int bin=0; bin<nOfBins; bin++) {
123
124        G4double e = energyVector[bin];
125        energies->push_back(e);
126        G4double value = 0.0;
127
128        if(e > tcut) {
129          for (G4int n=0; n<nShells; n++) {
130            G4double cross = FindValue(Z, e, n);
131            G4double p = theParam->Probability(Z, tcut, e, e, n);
132            value += cross * p * density;
133
134            if(verbose>0 && m == 0 && e>=1. && e<=0.) 
135            {
136              G4cout << "G4eIonCrossSH: e(MeV)= " << e/MeV
137                     << " n= " << n
138                     << " cross= " << cross
139                     << " p= " << p
140                     << " value= " << value
141                     << " tcut(MeV)= " << tcut/MeV
142                     << " rho= " << density
143                     << " Z= " << Z
144                     << G4endl;
145              }
146
147          }
148        }
149        cs->push_back(value);
150      }
151      G4VDataSetAlgorithm* algo = interp->Clone();
152      G4VEMDataSet* elSet = new G4EMDataSet(i,energies,cs,algo,1.,1.);
153      setForMat->AddComponent(elSet);
154    }
155    set->push_back(setForMat);
156  }
157
158  return set;
159}
160
161G4double G4eIonisationCrossSectionHandler::GetCrossSectionAboveThresholdForElement(G4double energy,
162                                                                                   G4double cutEnergy,
163                                                                                   G4int Z)
164{
165  G4int nShells = NumberOfComponents(Z);
166  G4double value = 0.;
167  if(energy > cutEnergy) 
168    {
169      for (G4int n=0; n<nShells; n++) {
170        G4double cross = FindValue(Z, energy, n);
171        G4double p = theParam->Probability(Z, cutEnergy, energy, energy, n);
172        value += cross * p;
173      }
174    }
175  return value;
176}
Note: See TracBrowser for help on using the repository browser.