Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4AntiNuclElastic.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// Geant4 Header : G4AntiNuclElastic
29//
30// Author : A. Galoyan
31//
32// Class Description
33// Model for AntiNuclear Nuclear Elastic Scattering;
34// Class Description - End
35
36#ifndef G4AntiNuclElastic_h
37#define G4AntiNuclElastic_h 1
38
39#include "G4HadronElastic.hh"
40#include "globals.hh"
41#include "G4Nucleus.hh"
43
45
47{
48public:
49
50 explicit G4AntiNuclElastic();
51
52 ~G4AntiNuclElastic() override;
53
55 G4double plab, G4int Z, G4int A) override;
56
58 G4int Z, G4int A);
59
61 G4double plab, G4int Z, G4int A);
62
63
65 G4double momentum );
66
68
70
72
74
76
78
80
82
83private:
84
85 // Assignment operator and copy constructor
86 G4AntiNuclElastic & operator=(const G4AntiNuclElastic &right);
88
89 G4ComponentAntiNuclNuclearXS* cs; //cross section of antiA-A interaction
90 const G4ParticleDefinition* fParticle;
91
92 G4double fTetaCMS; // sampled Theta in CMS
93 G4double fThetaLab; //sampled Theta in Lab system
94 G4double fWaveVector;
95 G4double fBeta; // velosity of projectile
96 G4double fZommerfeld; // parameter of Zommerfeld for calculation of Coulomb cross-section
97 G4double fAm; // parameter for calculation of Coulomb cross-section
98 G4double fRa; // Radius of target
99 G4double fRef; // Effective radiuse for Calculation of hadron cross-section
100 G4double fceff; // Effective diffuse parameter
101
102 G4ThreeVector fbst; // boost vector
103 G4double fptot; // momentum of projectile in CMS system
104 G4double fTmax;
105
106 G4ParticleDefinition* theAProton;
107 G4ParticleDefinition* theANeutron;
108 G4ParticleDefinition* theADeuteron;
109 G4ParticleDefinition* theATriton;
110 G4ParticleDefinition* theAAlpha;
111 G4ParticleDefinition* theAHe3;
112
113 G4ParticleDefinition* theProton;
114 G4ParticleDefinition* theNeutron;
115 G4ParticleDefinition* theDeuteron;
116 G4ParticleDefinition* theAlpha;
117
118};
119
125
126#endif
127
128
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
const G4double A[17]
G4double SampleInvariantT(const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A) override
~G4AntiNuclElastic() override
G4double BesselOneByArg(G4double z)
G4double CalculateParticleBeta(const G4ParticleDefinition *particle, G4double momentum)
G4double SampleThetaLab(const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
G4double CalculateZommerfeld(G4double beta, G4double Z1, G4double Z2)
G4double DampFactor(G4double z)
G4double BesselJone(G4double z)
G4double GetcosTeta1(G4double plab, G4int A)
G4double BesselJzero(G4double z)
G4double SampleThetaCMS(const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
G4ComponentAntiNuclNuclearXS * GetComponentCrossSection()
G4double CalculateAm(G4double momentum, G4double n, G4double Z)
G4HadronElastic(const G4String &name="hElasticLHEP")