Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4eeCrossSections.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//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class header file
30//
31//
32// File name: G4eeCrossSections
33//
34// Author: Vladimir Ivanchenko
35//
36// Creation date: 15.08.2004
37//
38// Modifications:
39//
40
41//
42// Class Description:
43//
44
45// -------------------------------------------------------------------
46//
47
48#ifndef G4eeCrossSections_h
49#define G4eeCrossSections_h 1
50
51#include "G4ThreeVector.hh"
52#include "globals.hh"
53#include <complex>
54
55class G4PhysicsVector;
56
58{
59
60public:
61
63
65
67
69
71
73
75
77
78 std::complex<G4double> DpRho(G4double e);
79
80 // hide assignment operator
83
84private:
85
86 void Initialise();
87
88 G4double Width2p(G4double s_inv, G4double mres, G4double gconst,
89 G4double br, G4double mp);
90
91 G4double Width3p(G4double s_inv, G4double mres, G4double gconst,
92 G4double br);
93
94 G4double PhaseSpace3p(G4double e);
95
96 G4double WidthPg(G4double s_inv, G4double mres, G4double gconst,
97 G4double br, G4double mp);
98
99 G4double WidthRho(G4double e);
100
101 G4double WidthOm(G4double e);
102
103 G4double WidthPhi(G4double e);
104
105 std::complex<G4double> DpOm(G4double e);
106
107 std::complex<G4double> DpPhi(G4double e);
108
109 G4double MsPi, MsPi0, MsEta, MsEtap, MsKs, MsKc, MsRho, MsOm;
110 G4double MsF0, MsA0, MsPhi, MsK892, MsK0892;
111 G4double GRho, GOm, GPhi, GK892, GK0892, PhRho, PhOm, PhPhi, PhRhoPi;
112 G4double BrRhoPiG, BrRhoPi0G, BrRhoEtaG, BrRhoEe, BrOm3Pi;
113 G4double BrOmPi0G, BrOmEtaG, BrOm2Pi, PhOm2Pi, BrOmEe;
114 G4double BrPhi2Kc, BrPhiKsKl, BrPhi3Pi, BrPhiPi0G, BrPhiEtaG;
115 G4double BrPhi2Pi, PhPhi2Pi, BrPhiEe, MeVnb, Alpha;
116 G4double AOmRho, ARhoPRho, cterm, mssig, gsig, brsigpipi;
117 G4double msrho1450, msrho1700, grho1450, grho1700;
118 G4double arhoompi0, arho1450ompi0, arho1700ompi0, phrhoompi0;
119 G4double phrho1450ompi0, phrho1700ompi0, aomrhopi0, phomrhopi0;
120 G4double arhopi0pi0g, aompi0pi0g, phrhopi0pi0g, phompi0pi0g;
121 G4double brrho1450ompi0, brrho1450pipi, brrho1700ompi0;
122 G4double brrho1700pipi, aphirhopi0, phphirhopi0;
123 G4double arhosigg, phrhosigg, aomsigg, phomsigg;
124
125 G4double MsRho3, MsOm3, MsPhi3;
126
127 G4PhysicsVector* ph3p;
128
129};
130
131
132//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
133
134#endif
double G4double
Definition G4Types.hh:83
G4double CrossSection2pi(G4double)
G4double CrossSection2Kcharged(G4double)
G4double CrossSection2Kneutral(G4double)
G4double CrossSectionEtaG(G4double)
G4double CrossSection3pi(G4double)
G4eeCrossSections(const G4eeCrossSections &)=delete
G4eeCrossSections & operator=(const G4eeCrossSections &right)=delete
G4double CrossSectionPi0G(G4double)
std::complex< G4double > DpRho(G4double e)