Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ANSTOecpssrMixsModel.hh
Go to the documentation of this file.
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// History:
27// -----------
28// 10 Nov 2021 S. Guatelli & S. Bakr - 1st implementation
29//
30// Class description
31// ----------------
32// Computation of K, L & M shell ECPSSR ionisation cross sections for protons and alphas
33// Based on the work of
34// - S. Bakr et al. (2021) NIM B, 507:11-19.
35// - S. Bakr et al (2018), NIMB B, 436: 285-291.
36// ---------------------------------------------------------------------------------------
37
38#ifndef G4ANSTOecpssrMixsModel_HH
39#define G4ANSTOecpssrMixsModel_HH 1
40
41#include "G4VecpssrMiModel.hh"
42#include "globals.hh"
43#include <map>
44
46class G4VEMDataSet;
47
49
50{
51public:
52
54
56
57 G4double CalculateM1CrossSection (G4int zTarget, G4double massIncident, G4double energyIncident);
58 G4double CalculateM2CrossSection (G4int zTarget, G4double massIncident, G4double energyIncident);
59 G4double CalculateM3CrossSection (G4int zTarget, G4double massIncident, G4double energyIncident);
60 G4double CalculateM4CrossSection (G4int zTarget, G4double massIncident, G4double energyIncident);
61 G4double CalculateM5CrossSection (G4int zTarget, G4double massIncident, G4double energyIncident);
62private:
63
64 G4double CalculateMiCrossSection (G4int zTarget, G4double massIncident, G4double energyIncident, G4int mShellId);
65
67 G4ANSTOecpssrMixsModel & operator = (const G4ANSTOecpssrMixsModel &right);
68
69 G4VDataSetAlgorithm* interpolation;
70
71 std::vector< std::map<G4int,G4VEMDataSet*> > protonMiXsVector;
72 std::vector< std::map<G4int,G4VEMDataSet*> > alphaMiXsVector;
73 std::vector< std::map<G4int,G4VEMDataSet*> > carbonMiXsVector;
74
75 std::map< G4int , G4VEMDataSet* > protonM1DataSetMap;
76 std::map< G4int , G4VEMDataSet* > protonM2DataSetMap;
77 std::map< G4int , G4VEMDataSet* > protonM3DataSetMap;
78 std::map< G4int , G4VEMDataSet* > protonM4DataSetMap;
79 std::map< G4int , G4VEMDataSet* > protonM5DataSetMap;
80
81 std::map< G4int , G4VEMDataSet* > alphaM1DataSetMap;
82 std::map< G4int , G4VEMDataSet* > alphaM2DataSetMap;
83 std::map< G4int , G4VEMDataSet* > alphaM3DataSetMap;
84 std::map< G4int , G4VEMDataSet* > alphaM4DataSetMap;
85 std::map< G4int , G4VEMDataSet* > alphaM5DataSetMap;
86
87 std::map< G4int , G4VEMDataSet* > carbonM1DataSetMap;
88 std::map< G4int , G4VEMDataSet* > carbonM2DataSetMap;
89 std::map< G4int , G4VEMDataSet* > carbonM3DataSetMap;
90 std::map< G4int , G4VEMDataSet* > carbonM4DataSetMap;
91 std::map< G4int , G4VEMDataSet* > carbonM5DataSetMap;
92};
93
94#endif
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
G4double CalculateM2CrossSection(G4int zTarget, G4double massIncident, G4double energyIncident)
G4double CalculateM5CrossSection(G4int zTarget, G4double massIncident, G4double energyIncident)
G4double CalculateM4CrossSection(G4int zTarget, G4double massIncident, G4double energyIncident)
G4double CalculateM1CrossSection(G4int zTarget, G4double massIncident, G4double energyIncident)
G4double CalculateM3CrossSection(G4int zTarget, G4double massIncident, G4double energyIncident)