Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
HadronPhysicsQGSP_BERT_95.cc
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// $Id$
27//
28//---------------------------------------------------------------------------
29//
30// Modified:
31//
32//----------------------------------------------------------------------------
33//
34#include <iomanip>
35
37
38#include "globals.hh"
39#include "G4ios.hh"
40#include "G4SystemOfUnits.hh"
42#include "G4ParticleTable.hh"
43
44#include "G4MesonConstructor.hh"
47
54
55#include "G4PhysListUtil.hh"
56
57// factory
59//
61
63 : G4VPhysicsConstructor("hInelastic QGSP_BERT_95")
64 , theNeutrons(0)
65 , theLEPNeutron(0)
66 , theQGSPNeutron(0)
67 , theBertiniNeutron(0)
68 , thePiK(0)
69 , theLEPPiK(0)
70 , theQGSPPiK(0)
71 , theBertiniPiK(0)
72 , thePro(0)
73 , theLEPPro(0)
74 , theQGSPPro(0)
75 , theBertiniPro(0)
76 , theMiscLHEP(0)
77 , QuasiElastic(true)
78 , ProjectileDiffraction(false)
79 , xsBarashenkovGGPion(0)
80 , xsGeisha(0)
81 , xsAxenWellischGGProton(0)
82 , xsLaidlawWellischGGNeutron(0)
83{
84}
85
88 , theNeutrons(0)
89 , theLEPNeutron(0)
90 , theQGSPNeutron(0)
91 , theBertiniNeutron(0)
92 , thePiK(0)
93 , theLEPPiK(0)
94 , theQGSPPiK(0)
95 , theBertiniPiK(0)
96 , thePro(0)
97 , theLEPPro(0)
98 , theQGSPPro(0)
99 , theBertiniPro(0)
100 , theMiscLHEP(0)
101 , QuasiElastic(quasiElastic)
102 , ProjectileDiffraction(false)
103 , xsBarashenkovGGPion(0)
104 , xsGeisha(0)
105 , xsAxenWellischGGProton(0)
106 , xsLaidlawWellischGGNeutron(0)
107{
108}
109
110void HadronPhysicsQGSP_BERT_95::CreateModels()
111{
112 theNeutrons=new G4NeutronBuilder;
113 theNeutrons->RegisterMe(theQGSPNeutron=new G4QGSPNeutronBuilder(QuasiElastic, ProjectileDiffraction));
114 theNeutrons->RegisterMe(theLEPNeutron=new G4LEPNeutronBuilder);
115 theLEPNeutron->SetMinInelasticEnergy(9.5*GeV);
116 theLEPNeutron->SetMaxInelasticEnergy(25*GeV);
117
118 theNeutrons->RegisterMe(theBertiniNeutron=new G4BertiniNeutronBuilder);
119 theBertiniNeutron->SetMinEnergy(0.0*GeV);
120 theBertiniNeutron->SetMaxEnergy(9.9*GeV);
121
122 thePro=new G4ProtonBuilder;
123 thePro->RegisterMe(theQGSPPro=new G4QGSPProtonBuilder(QuasiElastic, ProjectileDiffraction));
124 thePro->RegisterMe(theLEPPro=new G4LEPProtonBuilder);
125 theLEPPro->SetMinEnergy(9.5*GeV);
126 theLEPPro->SetMaxEnergy(25*GeV);
127
128 thePro->RegisterMe(theBertiniPro=new G4BertiniProtonBuilder);
129 theBertiniPro->SetMaxEnergy(9.9*GeV);
130
131 thePiK=new G4PiKBuilder;
132 thePiK->RegisterMe(theQGSPPiK=new G4QGSPPiKBuilder(QuasiElastic));
133 thePiK->RegisterMe(theLEPPiK=new G4LEPPiKBuilder);
134 theLEPPiK->SetMaxEnergy(25*GeV);
135 theLEPPiK->SetMinEnergy(9.5*GeV);
136
137 thePiK->RegisterMe(theBertiniPiK=new G4BertiniPiKBuilder);
138 theBertiniPiK->SetMaxEnergy(9.9*GeV);
139
140 theMiscLHEP=new G4MiscLHEPBuilder;
141}
142
144{
145 delete theMiscLHEP;
146 delete theQGSPNeutron;
147 delete theLEPNeutron;
148 delete theNeutrons;
149 delete theBertiniNeutron;
150 delete theQGSPPro;
151 delete theLEPPro;
152 delete thePro;
153 delete theBertiniPro;
154 delete theQGSPPiK;
155 delete theLEPPiK;
156 delete theBertiniPiK;
157 delete thePiK;
158
159 delete xsBarashenkovGGPion;
160 delete xsGeisha;
161 delete xsAxenWellischGGProton;
162 delete xsLaidlawWellischGGNeutron;
163}
164
166{
167 G4MesonConstructor pMesonConstructor;
168 pMesonConstructor.ConstructParticle();
169
170 G4BaryonConstructor pBaryonConstructor;
171 pBaryonConstructor.ConstructParticle();
172
173 G4ShortLivedConstructor pShortLivedConstructor;
174 pShortLivedConstructor.ConstructParticle();
175}
176
177#include "G4ProcessManager.hh"
179{
180 CreateModels();
181 theNeutrons->Build();
182 thePro->Build();
183 thePiK->Build();
184 theMiscLHEP->Build();
185
186 // Inelastic cross sections
187
188 // --- Pions ---
189 // Use Barashenkov inelastic pion cross section up to 91 GeV,
190 // and Glauber-Gribov above
191 xsBarashenkovGGPion = new G4CrossSectionPairGG(new G4PiNuclearCrossSection(), 91*GeV);
194
195 // --- Kaons ---
196 // Use Geisha inelastic cross sections
197 xsGeisha = new G4HadronInelasticDataSet();
202
203 // --- Protons ---
204 // Use Axen-Wellisch inelastic proton cross section up to 91 GeV,
205 // and Glauber-Gribov above
206 xsAxenWellischGGProton = new G4CrossSectionPairGG(new G4ProtonInelasticCrossSection(), 91*GeV);
208
209 // --- Neutrons ---
210 // Use Laidlaw-Wellisch inelastic neutron cross section up to 91 GeV,
211 // and Glauber-Gribov above
212 xsLaidlawWellischGGNeutron = new G4CrossSectionPairGG(new G4NeutronInelasticCrossSection(), 91*GeV);
214
215 // --- Hyperons ---
216 // Use Geisha inelastic cross sections
229
230 // --- AntiBaryons ---
231 // Use Geisha inelastic cross sections
234
235}
#define G4_DECLARE_PHYSCONSTR_FACTORY(physics_constructor)
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
static G4AntiLambda * AntiLambda()
static G4AntiNeutron * AntiNeutron()
static G4AntiOmegaMinus * AntiOmegaMinus()
static G4AntiProton * AntiProton()
Definition: G4AntiProton.cc:93
static G4AntiSigmaMinus * AntiSigmaMinus()
static G4AntiSigmaPlus * AntiSigmaPlus()
static G4AntiXiMinus * AntiXiMinus()
static G4AntiXiZero * AntiXiZero()
static void ConstructParticle()
void SetMaxEnergy(G4double aM)
void AddDataSet(G4VCrossSectionDataSet *aDataSet)
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:113
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:113
static G4KaonZeroLong * KaonZeroLong()
static G4KaonZeroShort * KaonZeroShort()
void SetMaxInelasticEnergy(G4double aM)
void SetMinInelasticEnergy(G4double aM)
void SetMaxEnergy(G4double aM)
void SetMinEnergy(G4double aM)
void SetMinEnergy(G4double aM)
void SetMaxEnergy(G4double aM)
static G4Lambda * Lambda()
Definition: G4Lambda.cc:108
static void ConstructParticle()
void RegisterMe(G4VNeutronBuilder *aB)
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4OmegaMinus * OmegaMinus()
static G4HadronicProcess * FindInelasticProcess(const G4ParticleDefinition *)
void RegisterMe(G4VPiKBuilder *aB)
Definition: G4PiKBuilder.hh:58
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
void RegisterMe(G4VProtonBuilder *aB)
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4SigmaMinus * SigmaMinus()
static G4SigmaPlus * SigmaPlus()
Definition: G4SigmaPlus.cc:108
static G4XiMinus * XiMinus()
Definition: G4XiMinus.cc:106
static G4XiZero * XiZero()
Definition: G4XiZero.cc:106