Geant4-11
G4ChipsProtonInelasticXS.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// GEANT4 physics class: G4ChipsProtonInelasticXS -- header file
27// Created: M.V. Kossov, CERN/ITEP(Moscow), 20-Dec-01
28// The last update: M.V. Kossov, CERN/ITEP (Moscow) 17-May-02
29//
30// ****************************************************************************************
31// Short description: Cross-sections extracted (by W.Pokorski) from the CHIPS package for
32// proton-nuclear interactions. Original author: M. Kossov
33// -------------------------------------------------------------------------------------
34//
35
36#ifndef G4ChipsProtonInelasticXS_h
37#define G4ChipsProtonInelasticXS_h 1
38
39#include "G4ParticleTable.hh"
40#include "G4NucleiProperties.hh"
41#include <vector>
43
45{
46
47public:
48
50
52
53 static const char* Default_Name() {return "ChipsProtonInelasticXS";}
54
55 virtual void CrossSectionDescription(std::ostream&) const;
56
58 const G4Element* elm,
59 const G4Material* mat );
60
61 // At present momentum (pMom) in MeV/c, CS in mb (@@ Units)
63 const G4Isotope* iso = 0,
64 const G4Element* elm = 0,
65 const G4Material* mat = 0);
66
67 virtual G4double GetChipsCrossSection(G4double momentum, G4int Z, G4int N, G4int pdg);
68
69private:
71 G4int N, G4double Momentum);
72
76 G4double ThresholdMomentum(G4int targZ, G4int targN); // Threshold of pA reaction (MeV/c)
78// Body
79private:
80 G4double* lastLEN; // Pointer to the last array of LowEnergy cross sections
81 G4double* lastHEN; // Pointer to the last array of HighEnergy cross sections
82 G4int lastN; // The last N of calculated nucleus
83 G4int lastZ; // The last Z of calculated nucleus
84 G4double lastP; // Last used in the cross section Momentum
85 G4double lastTH; // Last value of the Momentum Threshold
86 G4double lastCS; // Last value of the Cross Section
87 G4int lastI; // The last position in the DAMDB
88 std::vector<G4double*>* LEN; // Vector of pointers to LowEnProtonCrossSection
89 std::vector<G4double*>* HEN; // Vector of pointers to HighEnProtonCrossSection
90
91 G4int j=0; // A#0f Z/N-records already tested in AMDB
92 std::vector <G4int> colN; // Vector of N for calculated nuclei (isotops)
93 std::vector <G4int> colZ; // Vector of Z for calculated nuclei (isotops)
94 std::vector <G4double> colP; // Vector of last momenta for the reaction
95 std::vector <G4double> colTH; // Vector of energy thresholds for the reaction
96 std::vector <G4double> colCS; // Vector of last cross sections for the reaction
97
98};
99
100#endif
G4double Y(G4double density)
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
G4double CalculateCrossSection(G4int F, G4int I, G4int PDG, G4int Z, G4int N, G4double Momentum)
std::vector< G4double > colTH
std::vector< G4double * > * LEN
virtual void CrossSectionDescription(std::ostream &) const
G4double CrossSectionFormula(G4int targZ, G4int targN, G4double P, G4double lP)
G4double EquLinearFit(G4double X, G4int N, G4double X0, G4double DX, G4double *Y)
static const char * Default_Name()
std::vector< G4double * > * HEN
G4double ThresholdMomentum(G4int targZ, G4int targN)
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int tgZ, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
G4double CrossSectionLog(G4int targZ, G4int targN, G4double lP)
G4double CrossSectionLin(G4int targZ, G4int targN, G4double P)
std::vector< G4double > colP
virtual G4bool IsIsoApplicable(const G4DynamicParticle *Pt, G4int Z, G4int A, const G4Element *elm, const G4Material *mat)
virtual G4double GetChipsCrossSection(G4double momentum, G4int Z, G4int N, G4int pdg)
std::vector< G4double > colCS
static double P[]