source: trunk/source/processes/cuts/src/G4RToEConvForGamma.cc @ 961

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

update processes

File size: 6.4 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: G4RToEConvForGamma.cc,v 1.4 2006/06/29 19:30:24 gunter Exp $
28// GEANT4 tag $Name: geant4-09-02-ref-02 $
29//
30//
31// --------------------------------------------------------------
32//      GEANT 4 class implementation file/  History:
33//    5 Oct. 2002, H.Kuirashige : Structure created based on object model
34// --------------------------------------------------------------
35
36#include "G4RToEConvForGamma.hh"
37#include "G4ParticleDefinition.hh"
38#include "G4ParticleTable.hh"
39#include "G4Material.hh"
40#include "G4PhysicsLogVector.hh"
41
42#include "G4ios.hh"
43
44G4RToEConvForGamma::G4RToEConvForGamma() : G4VRangeToEnergyConverter()
45{   
46  theParticle =  G4ParticleTable::GetParticleTable()->FindParticle("gamma");
47  if (theParticle ==0) {
48#ifdef G4VERBOSE
49    if (GetVerboseLevel()>0) {
50      G4cout << " G4RToEConvForGamma::G4RToEConvForGamma() ";
51      G4cout << " Gamma is not defined !!" << G4endl;
52    }
53#endif
54  } 
55  TotBin = 100;
56}
57
58G4RToEConvForGamma::~G4RToEConvForGamma()
59{ 
60}
61
62
63// ***********************************************************************
64// ******************* BuildAbsorptionLengthVector ***********************
65// ***********************************************************************
66void G4RToEConvForGamma::BuildAbsorptionLengthVector(
67                            const G4Material* aMaterial,
68                            G4double       ,     
69                            G4double       ,
70                            G4RangeVector* absorptionLengthVector )
71{
72  // fill the absorption length vector for this material
73  // absorption length is defined here as
74  //
75  //    absorption length = 5./ macroscopic absorption cross section
76  //
77  const G4CrossSectionTable* aCrossSectionTable = (G4CrossSectionTable*)(theLossTable);
78  const G4ElementVector* elementVector = aMaterial->GetElementVector();
79  const G4double* atomicNumDensityVector = aMaterial->GetAtomicNumDensityVector();
80
81  //  fill absorption length vector
82  G4int NumEl = aMaterial->GetNumberOfElements();
83  G4double absorptionLengthMax = 0.0;
84  for (size_t ibin=0; ibin<size_t(TotBin); ibin++) {
85    G4double lowEdgeEnergy = absorptionLengthVector->GetLowEdgeEnergy(ibin);
86   
87    G4double SIGMA = 0. ;
88   
89    for (size_t iel=0; iel<size_t(NumEl); iel++) {
90      G4bool isOut;
91      G4int IndEl = (*elementVector)[iel]->GetIndex();
92      SIGMA +=  atomicNumDensityVector[iel]*
93                   (*aCrossSectionTable)[IndEl]->GetValue(lowEdgeEnergy,isOut);
94    }
95    //  absorption length=5./SIGMA
96    absorptionLengthVector->PutValue(ibin, 5./SIGMA);
97    if (absorptionLengthMax < 5./SIGMA ) absorptionLengthMax = 5./SIGMA;
98  }
99}
100
101
102
103// ***********************************************************************
104// ********************** ComputeCrossSection ****************************
105// ***********************************************************************
106G4double G4RToEConvForGamma::ComputeCrossSection(G4double AtomicNumber,
107                                                 G4double KineticEnergy) const
108{
109  //  Compute the "absorption" cross section of the photon "absorption"
110  //  cross section means here the sum of the cross sections of the
111  //  pair production, Compton scattering and photoelectric processes
112  static G4double Z; 
113  const  G4double t1keV = 1.*keV;
114  const  G4double t200keV = 200.*keV;
115  const  G4double t100MeV = 100.*MeV;
116
117  static G4double s200keV, s1keV;
118  static G4double tmin, tlow; 
119  static G4double smin, slow;
120  static G4double cmin, clow, chigh;
121  //  compute Z dependent quantities in the case of a new AtomicNumber
122  if(std::abs(AtomicNumber-Z)>0.1)  {
123    Z = AtomicNumber;
124    G4double Zsquare = Z*Z;
125    G4double Zlog = std::log(Z);
126    G4double Zlogsquare = Zlog*Zlog;
127
128    s200keV = (0.2651-0.1501*Zlog+0.02283*Zlogsquare)*Zsquare;
129    tmin = (0.552+218.5/Z+557.17/Zsquare)*MeV;
130    smin = (0.01239+0.005585*Zlog-0.000923*Zlogsquare)*std::exp(1.5*Zlog);
131    cmin=std::log(s200keV/smin)/(std::log(tmin/t200keV)*std::log(tmin/t200keV));
132    tlow = 0.2*std::exp(-7.355/std::sqrt(Z))*MeV;
133    slow = s200keV*std::exp(0.042*Z*std::log(t200keV/tlow)*std::log(t200keV/tlow));
134    s1keV = 300.*Zsquare;
135    clow =std::log(s1keV/slow)/std::log(tlow/t1keV);
136
137    chigh=(7.55e-5-0.0542e-5*Z)*Zsquare*Z/std::log(t100MeV/tmin);
138  }
139
140  //  calculate the cross section (using an approximate empirical formula)
141  G4double s;
142  if ( KineticEnergy<tlow ) {
143    if(KineticEnergy<t1keV) s = slow*std::exp(clow*std::log(tlow/t1keV));
144    else                    s = slow*std::exp(clow*std::log(tlow/KineticEnergy));
145
146  } else if ( KineticEnergy<t200keV ) {
147    s = s200keV
148         * std::exp(0.042*Z*std::log(t200keV/KineticEnergy)*std::log(t200keV/KineticEnergy));
149
150  } else if( KineticEnergy<tmin ){
151    s = smin
152         * std::exp(cmin*std::log(tmin/KineticEnergy)*std::log(tmin/KineticEnergy));
153
154  } else {
155    s = smin + chigh*std::log(KineticEnergy/tmin);
156
157  }
158  return s * barn;
159}
160
Note: See TracBrowser for help on using the repository browser.