Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4MiscLHEPBuilder.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// GEANT4 tag $Name: not supported by cvs2svn $
27//
28//---------------------------------------------------------------------------
29//
30// ClassName: G4MiscLHEPBuilder
31//
32// Author: 2002 J.P. Wellisch
33//
34// Modified:
35// 16.11.2005 G.Folger: don't keep processes as data members, but new these
36// 13.06.2006 G.Folger: (re)move elastic scatterring
37//
38//----------------------------------------------------------------------------
39//
40#include "G4MiscLHEPBuilder.hh"
41
42#include "G4SystemOfUnits.hh"
44#include "G4ParticleTable.hh"
45#include "G4ProcessManager.hh"
46
48 theAntiProtonInelastic(0), theLEAntiProtonModel(0),
49 theHEAntiProtonModel(0),
50 theAntiNeutronInelastic(0), theLEAntiNeutronModel(0),
51 theHEAntiNeutronModel(0),
52 theLambdaInelastic(0), theLELambdaModel(0), theHELambdaModel(0),
53 theAntiLambdaInelastic(0), theLEAntiLambdaModel(0), theHEAntiLambdaModel(0),
54 theSigmaMinusInelastic(0), theLESigmaMinusModel(0), theHESigmaMinusModel(0),
55 theAntiSigmaMinusInelastic(0), theLEAntiSigmaMinusModel(0), theHEAntiSigmaMinusModel(0),
56 theSigmaPlusInelastic(0), theLESigmaPlusModel(0), theHESigmaPlusModel(0),
57 theAntiSigmaPlusInelastic(0), theLEAntiSigmaPlusModel(0), theHEAntiSigmaPlusModel(0),
58 theXiZeroInelastic(0), theLEXiZeroModel(0), theHEXiZeroModel(0),
59 theAntiXiZeroInelastic(0), theLEAntiXiZeroModel(0), theHEAntiXiZeroModel(0),
60 theXiMinusInelastic(0), theLEXiMinusModel(0), theHEXiMinusModel(0),
61 theAntiXiMinusInelastic(0), theLEAntiXiMinusModel(0), theHEAntiXiMinusModel(0),
62 theOmegaMinusInelastic(0), theLEOmegaMinusModel(0), theHEOmegaMinusModel(0),
63 theAntiOmegaMinusInelastic(0), theLEAntiOmegaMinusModel(0), theHEAntiOmegaMinusModel(0),
64 wasActivated(false)
65
66{}
67
68
70{
71}
72
74{
75 G4ProcessManager * aProcMan = 0;
76 wasActivated = true;
77
78 // anti-Proton
79 theAntiProtonInelastic = new G4AntiProtonInelasticProcess();
81 theLEAntiProtonModel = new G4LEAntiProtonInelastic();
82 theHEAntiProtonModel = new G4HEAntiProtonInelastic();
83 theHEAntiProtonModel->SetMaxEnergy(100*TeV);
84 theAntiProtonInelastic->RegisterMe(theLEAntiProtonModel);
85 theAntiProtonInelastic->RegisterMe(theHEAntiProtonModel);
86 aProcMan->AddDiscreteProcess(theAntiProtonInelastic);
87
88 // AntiNeutron
89 theAntiNeutronInelastic = new G4AntiNeutronInelasticProcess();
91 theLEAntiNeutronModel = new G4LEAntiNeutronInelastic();
92 theHEAntiNeutronModel = new G4HEAntiNeutronInelastic();
93 theHEAntiNeutronModel->SetMaxEnergy(100*TeV);
94 theAntiNeutronInelastic->RegisterMe(theLEAntiNeutronModel);
95 theAntiNeutronInelastic->RegisterMe(theHEAntiNeutronModel);
96 aProcMan->AddDiscreteProcess(theAntiNeutronInelastic);
97
98 // Lambda
99 theLambdaInelastic = new G4LambdaInelasticProcess();
100 aProcMan = G4Lambda::Lambda()->GetProcessManager();
101 theLELambdaModel = new G4LELambdaInelastic();
102 theHELambdaModel = new G4HELambdaInelastic();
103 theHELambdaModel->SetMaxEnergy(100*TeV);
104 theLambdaInelastic->RegisterMe(theLELambdaModel);
105 theLambdaInelastic->RegisterMe(theHELambdaModel);
106 aProcMan->AddDiscreteProcess(theLambdaInelastic);
107
108 // AntiLambda
109 theAntiLambdaInelastic = new G4AntiLambdaInelasticProcess();
111 theLEAntiLambdaModel = new G4LEAntiLambdaInelastic();
112 theHEAntiLambdaModel = new G4HEAntiLambdaInelastic();
113 theHEAntiLambdaModel->SetMaxEnergy(100*TeV);
114 theAntiLambdaInelastic->RegisterMe(theLEAntiLambdaModel);
115 theAntiLambdaInelastic->RegisterMe(theHEAntiLambdaModel);
116 aProcMan->AddDiscreteProcess(theAntiLambdaInelastic);
117
118 // SigmaMinus
119 theSigmaMinusInelastic = new G4SigmaMinusInelasticProcess();
121 theLESigmaMinusModel = new G4LESigmaMinusInelastic();
122 theHESigmaMinusModel = new G4HESigmaMinusInelastic();
123 theHESigmaMinusModel->SetMaxEnergy(100*TeV);
124 theSigmaMinusInelastic->RegisterMe(theLESigmaMinusModel);
125 theSigmaMinusInelastic->RegisterMe(theHESigmaMinusModel);
126 aProcMan->AddDiscreteProcess(theSigmaMinusInelastic);
127
128 // anti-SigmaMinus
129 theAntiSigmaMinusInelastic = new G4AntiSigmaMinusInelasticProcess();
131 theLEAntiSigmaMinusModel = new G4LEAntiSigmaMinusInelastic();
132 theHEAntiSigmaMinusModel = new G4HEAntiSigmaMinusInelastic();
133 theHEAntiSigmaMinusModel->SetMaxEnergy(100*TeV);
134 theAntiSigmaMinusInelastic->RegisterMe(theLEAntiSigmaMinusModel);
135 theAntiSigmaMinusInelastic->RegisterMe(theHEAntiSigmaMinusModel);
136 aProcMan->AddDiscreteProcess(theAntiSigmaMinusInelastic);
137
138 // SigmaPlus
139 theSigmaPlusInelastic = new G4SigmaPlusInelasticProcess();
141 theLESigmaPlusModel = new G4LESigmaPlusInelastic();
142 theHESigmaPlusModel = new G4HESigmaPlusInelastic();
143 theHESigmaPlusModel->SetMaxEnergy(100*TeV);
144 theSigmaPlusInelastic->RegisterMe(theLESigmaPlusModel);
145 theSigmaPlusInelastic->RegisterMe(theHESigmaPlusModel);
146 aProcMan->AddDiscreteProcess(theSigmaPlusInelastic);
147
148 // anti-SigmaPlus
149 theAntiSigmaPlusInelastic = new G4AntiSigmaPlusInelasticProcess();
151 theLEAntiSigmaPlusModel = new G4LEAntiSigmaPlusInelastic();
152 theHEAntiSigmaPlusModel = new G4HEAntiSigmaPlusInelastic();
153 theHEAntiSigmaPlusModel->SetMaxEnergy(100*TeV);
154 theAntiSigmaPlusInelastic->RegisterMe(theLEAntiSigmaPlusModel);
155 theAntiSigmaPlusInelastic->RegisterMe(theHEAntiSigmaPlusModel);
156 aProcMan->AddDiscreteProcess(theAntiSigmaPlusInelastic);
157
158 // XiMinus
159 theXiMinusInelastic = new G4XiMinusInelasticProcess();
161 theLEXiMinusModel = new G4LEXiMinusInelastic();
162 theHEXiMinusModel = new G4HEXiMinusInelastic();
163 theHEXiMinusModel->SetMaxEnergy(100*TeV);
164 theXiMinusInelastic->RegisterMe(theLEXiMinusModel);
165 theXiMinusInelastic->RegisterMe(theHEXiMinusModel);
166 aProcMan->AddDiscreteProcess(theXiMinusInelastic);
167
168 // anti-XiMinus
169 theAntiXiMinusInelastic = new G4AntiXiMinusInelasticProcess();
171 theLEAntiXiMinusModel = new G4LEAntiXiMinusInelastic();
172 theHEAntiXiMinusModel = new G4HEAntiXiMinusInelastic();
173 theHEAntiXiMinusModel->SetMaxEnergy(100*TeV);
174 theAntiXiMinusInelastic->RegisterMe(theLEAntiXiMinusModel);
175 theAntiXiMinusInelastic->RegisterMe(theHEAntiXiMinusModel);
176 aProcMan->AddDiscreteProcess(theAntiXiMinusInelastic);
177
178 // XiZero
179 theXiZeroInelastic = new G4XiZeroInelasticProcess();
180 aProcMan = G4XiZero::XiZero()->GetProcessManager();
181 theLEXiZeroModel = new G4LEXiZeroInelastic();
182 theHEXiZeroModel = new G4HEXiZeroInelastic();
183 theHEXiZeroModel->SetMaxEnergy(100*TeV);
184 theXiZeroInelastic->RegisterMe(theLEXiZeroModel);
185 theXiZeroInelastic->RegisterMe(theHEXiZeroModel);
186 aProcMan->AddDiscreteProcess(theXiZeroInelastic);
187
188 // anti-XiZero
189 theAntiXiZeroInelastic = new G4AntiXiZeroInelasticProcess();
191 theLEAntiXiZeroModel = new G4LEAntiXiZeroInelastic();
192 theHEAntiXiZeroModel = new G4HEAntiXiZeroInelastic();
193 theHEAntiXiZeroModel->SetMaxEnergy(100*TeV);
194 theAntiXiZeroInelastic->RegisterMe(theLEAntiXiZeroModel);
195 theAntiXiZeroInelastic->RegisterMe(theHEAntiXiZeroModel);
196 aProcMan->AddDiscreteProcess(theAntiXiZeroInelastic);
197
198 // OmegaMinus
199 theOmegaMinusInelastic = new G4OmegaMinusInelasticProcess();
201 theLEOmegaMinusModel = new G4LEOmegaMinusInelastic();
202 theHEOmegaMinusModel = new G4HEOmegaMinusInelastic();
203 theHEOmegaMinusModel->SetMaxEnergy(100*TeV);
204 theOmegaMinusInelastic->RegisterMe(theLEOmegaMinusModel);
205 theOmegaMinusInelastic->RegisterMe(theHEOmegaMinusModel);
206 aProcMan->AddDiscreteProcess(theOmegaMinusInelastic);
207
208 // anti-OmegaMinus
209 theAntiOmegaMinusInelastic = new G4AntiOmegaMinusInelasticProcess();
211 theLEAntiOmegaMinusModel = new G4LEAntiOmegaMinusInelastic();
212 theHEAntiOmegaMinusModel = new G4HEAntiOmegaMinusInelastic();
213 theHEAntiOmegaMinusModel->SetMaxEnergy(100*TeV);
214 theAntiOmegaMinusInelastic->RegisterMe(theLEAntiOmegaMinusModel);
215 theAntiOmegaMinusInelastic->RegisterMe(theHEAntiOmegaMinusModel);
216 aProcMan->AddDiscreteProcess(theAntiOmegaMinusInelastic);
217}
218
219// 2002 by J.P. Wellisch
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()
void SetMaxEnergy(const G4double anEnergy)
void RegisterMe(G4HadronicInteraction *a)
static G4Lambda * Lambda()
Definition: G4Lambda.cc:108
virtual ~G4MiscLHEPBuilder()
static G4OmegaMinus * OmegaMinus()
G4ProcessManager * GetProcessManager() const
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
static G4SigmaMinus * SigmaMinus()
static G4SigmaPlus * SigmaPlus()
Definition: G4SigmaPlus.cc:108
static G4XiMinus * XiMinus()
Definition: G4XiMinus.cc:106
static G4XiZero * XiZero()
Definition: G4XiZero.cc:106