Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ChargeExchangeProcess.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//
29// Geant4 Hadron Charge Exchange Process -- source file
30//
31// Created 21 April 2006 V.Ivanchenko
32//
33// Modified:
34// 24-Apr-06 V.Ivanchenko add neutron scattering on hydrogen from CHIPS
35// 07-Jun-06 V.Ivanchenko fix problem of rotation of final state
36// 25-Jul-06 V.Ivanchenko add 19 MeV low energy for CHIPS
37// 23-Jan-07 V.Ivanchenko add cross section interfaces with Z and A
38// and do not use CHIPS for cross sections
39// 14-Sep-12 M.Kelsey -- Pass subType code to base ctor
40// 06-Aug-15 A.Ribon migrating to G4Pow
41
43#include "globals.hh"
44#include "G4SystemOfUnits.hh"
48#include "G4Element.hh"
49#include "G4ElementVector.hh"
50#include "G4IsotopeVector.hh"
51#include "G4Neutron.hh"
52#include "G4Proton.hh"
54
55#include "G4Pow.hh"
56
57
59 : G4HadronicProcess(procName,fChargeExchange), first(true)
60{
61 thEnergy = 20.*MeV;
62 pPDG = 0;
63 verboseLevel= 1;
65 theProton = G4Proton::Proton();
66 theNeutron = G4Neutron::Neutron();
67 theAProton = G4AntiProton::AntiProton();
68 theANeutron = G4AntiNeutron::AntiNeutron();
69 thePiPlus = G4PionPlus::PionPlus();
70 thePiMinus = G4PionMinus::PionMinus();
71 thePiZero = G4PionZero::PionZero();
72 theKPlus = G4KaonPlus::KaonPlus();
73 theKMinus = G4KaonMinus::KaonMinus();
76 theL = G4Lambda::Lambda();
77 theAntiL = G4AntiLambda::AntiLambda();
78 theSPlus = G4SigmaPlus::SigmaPlus();
80 theSMinus = G4SigmaMinus::SigmaMinus();
82 theS0 = G4SigmaZero::SigmaZero();
84 theXiMinus = G4XiMinus::XiMinus();
85 theXi0 = G4XiZero::XiZero();
86 theAXiMinus = G4AntiXiMinus::AntiXiMinus();
87 theAXi0 = G4AntiXiZero::AntiXiZero();
88 theOmega = G4OmegaMinus::OmegaMinus();
90 theD = G4Deuteron::Deuteron();
91 theT = G4Triton::Triton();
92 theA = G4Alpha::Alpha();
93 theHe3 = G4He3::He3();
94}
95
97{
98 if (factors) delete factors;
99}
100
102BuildPhysicsTable(const G4ParticleDefinition& aParticleType)
103{
104 if(first) {
105 first = false;
106 theParticle = &aParticleType;
107 pPDG = theParticle->GetPDGEncoding();
108
110
111 const size_t n = 10;
112 if(theParticle == thePiPlus || theParticle == thePiMinus ||
113 theParticle == theKPlus || theParticle == theKMinus ||
114 theParticle == theK0S || theParticle == theK0L) {
115
116 G4double F[n] = {0.33,0.27,0.29,0.31,0.27,0.18,0.13,0.1,0.09,0.07};
117 factors = new G4PhysicsLinearVector(0.0,2.0*GeV,n);
118 for(size_t i=0; i<n; i++) {factors->PutValue(i,F[i]);}
119
120 } else {
121
122 G4double F[n] = {0.50,0.45,0.40,0.35,0.30,0.25,0.06,0.04,0.005,0.0};
123 factors = new G4PhysicsLinearVector(0.0,4.0*GeV,n);
124 for(size_t i=0; i<n; i++) {factors->PutValue(i,F[i]);}
125 }
126
127 if(verboseLevel>1)
128 G4cout << "G4ChargeExchangeProcess for "
129 << theParticle->GetParticleName()
130 << G4endl;
131 }
133}
134
136 const G4DynamicParticle* dp,
137 const G4Element* elm,
138 const G4Material* mat)
139{
140 // gives the microscopic cross section in GEANT4 internal units
141 G4double Z = elm->GetZ();
142 G4int iz = G4int(Z);
143 G4double x = 0.0;
144
145 // The process is effective only above the threshold
146 if(iz == 1 || dp->GetKineticEnergy() < thEnergy) return x;
147
148 if(verboseLevel>1)
149 G4cout << "G4ChargeExchangeProcess compute GHAD CS for element "
150 << elm->GetName()
151 << G4endl;
152 x = store->GetCrossSection(dp, elm, mat);
153
154 if(verboseLevel>1)
155 G4cout << "G4ChargeExchangeProcess cross(mb)= " << x/millibarn
156 << " E(MeV)= " << dp->GetKineticEnergy()
157 << " " << theParticle->GetParticleName()
158 << " in Z= " << iz
159 << G4endl;
160 G4bool b;
161 G4double A = elm->GetN();
162 G4double ptot = dp->GetTotalMomentum();
163 x *= factors->GetValue(ptot, b)/G4Pow::GetInstance()->powA(A, 0.42);
164 if(theParticle == thePiPlus || theParticle == theProton ||
165 theParticle == theKPlus || theParticle == theANeutron)
166 { x *= (1.0 - Z/A); }
167
168 else if(theParticle == thePiMinus || theParticle == theNeutron ||
169 theParticle == theKMinus || theParticle == theAProton)
170 { x *= Z/A; }
171
172 if(theParticle->GetPDGMass() < GeV) {
173 if(ptot > 2.*GeV) x *= 4.0*GeV*GeV/(ptot*ptot);
174 }
175
176 if(verboseLevel>1)
177 G4cout << "Corrected cross(mb)= " << x/millibarn << G4endl;
178
179 return x;
180}
181
183IsApplicable(const G4ParticleDefinition& aParticleType)
184{
185 const G4ParticleDefinition* p = &aParticleType;
186 return (p == thePiPlus || p == thePiMinus ||
187 p == theProton || p == theNeutron ||
188 p == theAProton|| p == theANeutron||
189 p == theKPlus || p == theKMinus ||
190 p == theK0S || p == theK0L ||
191 p == theL);
192}
193
195DumpPhysicsTable(const G4ParticleDefinition& aParticleType)
196{
197 store->DumpPhysicsTable(aParticleType);
198}
@ fChargeExchange
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
const G4double A[17]
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
static G4Alpha * Alpha()
Definition G4Alpha.cc:83
static G4AntiLambda * AntiLambda()
static G4AntiNeutron * AntiNeutron()
static G4AntiOmegaMinus * AntiOmegaMinus()
static G4AntiProton * AntiProton()
static G4AntiSigmaMinus * AntiSigmaMinus()
static G4AntiSigmaPlus * AntiSigmaPlus()
static G4AntiSigmaZero * AntiSigmaZero()
static G4AntiXiMinus * AntiXiMinus()
static G4AntiXiZero * AntiXiZero()
virtual G4double GetElementCrossSection(const G4DynamicParticle *aParticle, const G4Element *anElement, const G4Material *mat=0)
virtual void DumpPhysicsTable(const G4ParticleDefinition &aParticleType)
virtual G4bool IsApplicable(const G4ParticleDefinition &aParticleType)
virtual void BuildPhysicsTable(const G4ParticleDefinition &aParticleType)
G4ChargeExchangeProcess(const G4String &procName="chargeExchange")
void DumpPhysicsTable(const G4ParticleDefinition &)
G4double GetCrossSection(const G4DynamicParticle *, const G4Material *)
static G4Deuteron * Deuteron()
Definition G4Deuteron.cc:90
G4double GetKineticEnergy() const
G4double GetTotalMomentum() const
G4double GetZ() const
Definition G4Element.hh:119
const G4String & GetName() const
Definition G4Element.hh:115
G4double GetN() const
Definition G4Element.hh:123
void AddDataSet(G4VCrossSectionDataSet *aDataSet)
void BuildPhysicsTable(const G4ParticleDefinition &) override
G4CrossSectionDataStore * GetCrossSectionDataStore()
static G4He3 * He3()
Definition G4He3.cc:90
static G4KaonMinus * KaonMinus()
static G4KaonPlus * KaonPlus()
static G4KaonZeroLong * KaonZeroLong()
static G4KaonZeroShort * KaonZeroShort()
static G4Lambda * Lambda()
Definition G4Lambda.cc:105
static G4Neutron * Neutron()
Definition G4Neutron.cc:101
static G4OmegaMinus * OmegaMinus()
const G4String & GetParticleName() const
void PutValue(const std::size_t index, const G4double value)
G4double GetValue(const G4double energy, G4bool &isOutRange) const
static G4PionMinus * PionMinus()
static G4PionPlus * PionPlus()
Definition G4PionPlus.cc:93
static G4PionZero * PionZero()
static G4Pow * GetInstance()
Definition G4Pow.cc:41
G4double powA(G4double A, G4double y) const
Definition G4Pow.hh:230
static G4Proton * Proton()
Definition G4Proton.cc:90
static G4SigmaMinus * SigmaMinus()
static G4SigmaPlus * SigmaPlus()
static G4SigmaZero * SigmaZero()
static G4Triton * Triton()
Definition G4Triton.cc:90
G4int verboseLevel
static G4XiMinus * XiMinus()
Definition G4XiMinus.cc:103
static G4XiZero * XiZero()
Definition G4XiZero.cc:103