Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4HadProcesses.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//---------------------------------------------------------------------------
27//
28// ClassName: G4HadProcesses
29//
30// Author: 8 July 2020 V.Ivanchenko
31//
32// Modified:
33//
34//----------------------------------------------------------------------------
35//
36
37#include "G4HadProcesses.hh"
38#include "G4SystemOfUnits.hh"
39#include "G4PhysListUtil.hh"
40#include "G4ParticleTable.hh"
41#include "G4Neutron.hh"
52#include "G4NeutronElasticXS.hh"
53#include "G4NeutronCaptureXS.hh"
54
56{
58}
59
61{
63}
64
66{
67 return FindInelasticProcess( FindParticle(pname) );
68}
69
71{
73}
74
76{
77 return FindElasticProcess( FindParticle(pname) );
78}
79
81{
83}
84
86{
88}
89
91{
92 G4CrossSectionInelastic* xs = nullptr;
94 if( comp != nullptr ) {
95 xs = new G4CrossSectionInelastic(comp);
96 } else if( "Glauber-Gribov" == compName ) {
98 } else if( "Glauber-Gribov Nucl-nucl" == compName ) {
100 } else if( "AntiAGlauber" == compName ) {
102 }
103 return xs;
104}
105
107{
108 G4CrossSectionElastic* xs = nullptr;
110 if( comp != nullptr ) {
111 xs = new G4CrossSectionElastic(comp);
112 } else if( "Glauber-Gribov" == compName ) {
114 } else if( "Glauber-Gribov Nucl-nucl" == compName ) {
116 } else if( "AntiAGlauber" == compName ) {
118 }
119 return xs;
120}
121
124{
125 G4bool isOK(false);
126 if( ptr != nullptr ) {
128 if( had != nullptr ) {
129 isOK = true;
130 had->AddDataSet( xs );
131 }
132 }
133 return isOK;
134}
135
137{
138 return AddInelasticCrossSection( FindParticle(pname), xs );
139}
140
142{
143 G4bool isOK(false);
144 if( ptr != nullptr ) {
146 if( had != nullptr ) {
147 isOK = true;
148 had->AddDataSet( xs );
149 }
150 }
151 return isOK;
152}
153
155{
156 return AddElasticCrossSection( FindParticle(pname), xs );
157}
158
160{
161 G4bool isOK(false);
163 if( had != nullptr ) {
164 isOK = true;
165 had->AddDataSet( xs );
166 }
167 return isOK;
168}
169
171{
172 G4bool isOK(false);
174 if( had != nullptr ) {
175 isOK = true;
176 had->AddDataSet( xs );
177 }
178 return isOK;
179}
180
182{
184 G4bool useNeutronGeneral = param->EnableNeutronGeneralProcess();
185
186 G4HadronicProcess* nCap = new G4NeutronCaptureProcess("nCapture");
187 nCap->RegisterMe(new G4NeutronRadCapture());
188
189 if(useNeutronGeneral) {
191 nGen->SetInelasticProcess(nInel);
192 nGen->SetCaptureProcess(nCap);
193 } else {
194 auto neutron = G4Neutron::Neutron();
196 nInel->AddDataSet(new G4NeutronInelasticXS());
197 ph->RegisterProcess(nInel, neutron);
198 ph->RegisterProcess(nCap, neutron);
199 if( param->ApplyFactorXS() ) {
201 }
202 }
203}
204
206{
208 G4bool useNeutronGeneral = param->EnableNeutronGeneralProcess();
209
210 if(useNeutronGeneral) {
212 nGen->SetElasticProcess(nEl);
213 } else {
214 auto neutron = G4Neutron::Neutron();
215 nEl->AddDataSet(new G4NeutronElasticXS());
217 ph->RegisterProcess(nEl, neutron);
218 if( param->ApplyFactorXS() ) {
220 }
221 }
222}
bool G4bool
Definition: G4Types.hh:86
G4VComponentCrossSection * GetComponentCrossSection(const G4String &name)
static G4CrossSectionDataSetRegistry * Instance()
static G4bool AddElasticCrossSection(const G4ParticleDefinition *, G4VCrossSectionDataSet *)
static G4CrossSectionElastic * ElasticXS(const G4String &componentName)
static G4HadronicProcess * FindInelasticProcess(const G4ParticleDefinition *)
static void BuildNeutronInelasticAndCapture(G4HadronicProcess *)
static const G4ParticleDefinition * FindParticle(const G4String &)
static G4bool AddFissionCrossSection(G4VCrossSectionDataSet *)
static G4HadronicProcess * FindCaptureProcess()
static G4bool AddInelasticCrossSection(const G4ParticleDefinition *, G4VCrossSectionDataSet *)
static G4bool AddCaptureCrossSection(G4VCrossSectionDataSet *)
static G4HadronicProcess * FindFissionProcess()
static G4CrossSectionInelastic * InelasticXS(const G4String &componentName)
static G4HadronicProcess * FindElasticProcess(const G4ParticleDefinition *)
static void BuildNeutronElastic(G4HadronicProcess *)
static G4HadronicParameters * Instance()
G4double XSFactorNucleonElastic() const
G4bool EnableNeutronGeneralProcess() const
G4double XSFactorNucleonInelastic() const
void AddDataSet(G4VCrossSectionDataSet *aDataSet)
void MultiplyCrossSectionBy(G4double factor)
void RegisterMe(G4HadronicInteraction *a)
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
static G4NeutronGeneralProcess * FindNeutronGeneralProcess()
static G4HadronicProcess * FindElasticProcess(const G4ParticleDefinition *)
static G4HadronicProcess * FindInelasticProcess(const G4ParticleDefinition *)
static G4HadronicProcess * FindCaptureProcess(const G4ParticleDefinition *)
static G4HadronicProcess * FindFissionProcess(const G4ParticleDefinition *)
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static G4PhysicsListHelper * GetPhysicsListHelper()