source: trunk/source/processes/electromagnetic/lowenergy/src/G4PaulKCrossSection.cc @ 1347

Last change on this file since 1347 was 1347, checked in by garnier, 13 years ago

geant4 tag 9.4

File size: 3.8 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// History:
28// -----------
29//  21 Apr 2008   H. Abdelohauwed - 1st implementation
30//  29 Apr 2009   ALF  Major Design Revision
31//
32// -------------------------------------------------------------------
33
34// Class description:
35// Low Energy Electromagnetic Physics, Cross section, p ionisation, K shell
36// Further documentation available from http://www.ge.infn.it/geant4/lowE
37
38// -------------------------------------------------------------------
39
40#include "globals.hh"
41#include "G4ios.hh"
42#include <fstream>
43#include <iomanip>
44//#include "G4CompositeEMDataSet.hh"
45//#include "G4ShellEMDataSet.hh"
46#include "G4EMDataSet.hh"
47//#include "G4VEMDataSet.hh"
48//#include "G4VDataSetAlgorithm.hh"
49#include "G4LogLogInterpolation.hh"
50#include "G4PaulKCrossSection.hh"
51#include "G4Proton.hh"
52#include "G4Alpha.hh"
53
54
55G4PaulKCrossSection::G4PaulKCrossSection()
56{ 
57
58 
59  interpolation = new G4LogLogInterpolation();
60
61  /*
62    G4String path = getenv("G4LEDATA");
63 
64    if (!path)
65    G4Exception("G4paulKCrossSection::G4paulKCrossSection: G4LEDATA environment variable not set");
66    G4cout << path + "/kcsPaul/kcs-" << G4endl;
67  */
68
69
70    for (G4int i=4; i<93; i++) {
71      protonDataSetMap[i] = new G4EMDataSet(i,interpolation);
72      protonDataSetMap[i]->LoadData("pixe/kpcsPaul/kcs-");
73    }
74    for (G4int i=6; i<93; i++) {
75      alphaDataSetMap[i] = new G4EMDataSet(i,interpolation);
76      alphaDataSetMap[i]->LoadData("pixe/kacsPaul/kacs-");
77    }
78
79
80
81
82}
83
84G4PaulKCrossSection::~G4PaulKCrossSection()
85{ 
86
87  protonDataSetMap.clear();
88  alphaDataSetMap.clear();
89
90}
91
92G4double G4PaulKCrossSection::CalculateKCrossSection(G4int zTarget,G4double massIncident, G4double energyIncident)
93{
94 
95  G4Proton* aProtone = G4Proton::Proton();
96  G4Alpha* aAlpha = G4Alpha::Alpha();
97 
98  G4double sigma = 0;
99
100  if (massIncident == aProtone->GetPDGMass() )
101    {
102     
103      sigma = protonDataSetMap[zTarget]->FindValue(energyIncident/MeV); 
104     
105    }
106  else
107    {
108      if (massIncident == aAlpha->GetPDGMass())
109        {
110         
111          sigma = alphaDataSetMap[zTarget]->FindValue(energyIncident/MeV); 
112         
113        }
114      else
115        { 
116          G4cout << "we can treat only Proton or Alpha incident particles " << G4endl;
117          sigma = 0.;
118        }
119    }
120 
121 
122  // sigma is in internal units (mm^2)
123  return sigma;
124}
125
126
127
128
129
130
131
132
133
134
135
136
Note: See TracBrowser for help on using the repository browser.