Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ChargeExchange Class Reference

#include <G4ChargeExchange.hh>

+ Inheritance diagram for G4ChargeExchange:

Public Member Functions

 G4ChargeExchange (G4ChargeExchangeXS *)
 
 ~G4ChargeExchange () override
 
 G4ChargeExchange (const G4ChargeExchange &right)=delete
 
const G4ChargeExchangeoperator= (const G4ChargeExchange &right)=delete
 
G4bool operator== (const G4ChargeExchange &right) const =delete
 
G4bool operator!= (const G4ChargeExchange &right) const =delete
 
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus) override
 
G4double SampleT (const G4ParticleDefinition *theSec, const G4int A, const G4double tmax) const
 
- Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
 
virtual ~G4HadronicInteraction ()
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
virtual G4bool IsApplicable (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
G4double GetMinEnergy () const
 
G4double GetMinEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMinEnergy (G4double anEnergy)
 
void SetMinEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMinEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4double GetMaxEnergy () const
 
G4double GetMaxEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMaxEnergy (const G4double anEnergy)
 
void SetMaxEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMaxEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4int GetVerboseLevel () const
 
void SetVerboseLevel (G4int value)
 
const G4StringGetModelName () const
 
void DeActivateFor (const G4Material *aMaterial)
 
void ActivateFor (const G4Material *aMaterial)
 
void DeActivateFor (const G4Element *anElement)
 
void ActivateFor (const G4Element *anElement)
 
G4bool IsBlocked (const G4Material *aMaterial) const
 
G4bool IsBlocked (const G4Element *anElement) const
 
void SetRecoilEnergyThreshold (G4double val)
 
G4double GetRecoilEnergyThreshold () const
 
virtual const std::pair< G4double, G4doubleGetFatalEnergyCheckLevels () const
 
virtual std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
virtual void ModelDescription (std::ostream &outFile) const
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void InitialiseModel ()
 
 G4HadronicInteraction (const G4HadronicInteraction &right)=delete
 
const G4HadronicInteractionoperator= (const G4HadronicInteraction &right)=delete
 
G4bool operator== (const G4HadronicInteraction &right) const =delete
 
G4bool operator!= (const G4HadronicInteraction &right) const =delete
 

Additional Inherited Members

- Protected Member Functions inherited from G4HadronicInteraction
void SetModelName (const G4String &nam)
 
G4bool IsBlocked () const
 
void Block ()
 
- Protected Attributes inherited from G4HadronicInteraction
G4HadFinalState theParticleChange
 
G4int verboseLevel
 
G4double theMinEnergy
 
G4double theMaxEnergy
 
G4bool isBlocked
 

Detailed Description

Definition at line 51 of file G4ChargeExchange.hh.

Constructor & Destructor Documentation

◆ G4ChargeExchange() [1/2]

G4ChargeExchange::G4ChargeExchange ( G4ChargeExchangeXS * ptr)
explicit

Definition at line 66 of file G4ChargeExchange.cc.

67 : G4HadronicInteraction("ChargeExchange"),
68 fXSection(ptr), fXSWeightFactor(1.0)
69{
70 lowEnergyLimit = 1.*CLHEP::MeV;
71 secID = G4PhysicsModelCatalog::GetModelID( "model_ChargeExchange" );
73 fHandler = new G4ExcitationHandler();
74 if (nullptr != fXSection) {
75 fXSWeightFactor = 1.0/fXSection->GetCrossSectionFactor();
76 }
77}
G4HadronicInteraction(const G4String &modelName="HadronicModel")
static G4NistManager * Instance()
static G4int GetModelID(const G4int modelIndex)

Referenced by G4ChargeExchange(), operator!=(), operator=(), and operator==().

◆ ~G4ChargeExchange()

G4ChargeExchange::~G4ChargeExchange ( )
override

Definition at line 79 of file G4ChargeExchange.cc.

80{
81 delete fHandler;
82}

◆ G4ChargeExchange() [2/2]

G4ChargeExchange::G4ChargeExchange ( const G4ChargeExchange & right)
delete

Member Function Documentation

◆ ApplyYourself()

G4HadFinalState * G4ChargeExchange::ApplyYourself ( const G4HadProjectile & aTrack,
G4Nucleus & targetNucleus )
overridevirtual

Reimplemented from G4HadronicInteraction.

Definition at line 84 of file G4ChargeExchange.cc.

86{
87 theParticleChange.Clear();
88 auto part = aTrack.GetDefinition();
89 G4double ekin = aTrack.GetKineticEnergy();
90
91 G4int A = targetNucleus.GetA_asInt();
92 G4int Z = targetNucleus.GetZ_asInt();
93
94 if (ekin <= lowEnergyLimit) {
95 return &theParticleChange;
96 }
97 theParticleChange.SetWeightChange(fXSWeightFactor);
98
99 G4int projPDG = part->GetPDGEncoding();
100
101 // for hydrogen targets and positive projectile change exchange
102 // is not possible on proton, only on deuteron
103 if (1 == Z && (211 == projPDG || 321 == projPDG)) { A = 2; }
104
105 if (verboseLevel > 1)
106 G4cout << "G4ChargeExchange for " << part->GetParticleName()
107 << " PDGcode= " << projPDG << " on nucleus Z= " << Z
108 << " A= " << A << " N= " << A - Z
109 << G4endl;
110
112 G4LorentzVector lv0 = aTrack.Get4Momentum();
113 G4double etot = mass1 + lv0.e();
114
115 // select final state
116 const G4ParticleDefinition* theSecondary =
117 fXSection->SampleSecondaryType(part, Z, A);
118 G4int pdg = theSecondary->GetPDGEncoding();
119
120 // omega(782) and f2(1270)
121 G4bool isShortLived = (pdg == 223 || pdg == 225);
122
123 // atomic number of the recoil nucleus
124 if (projPDG == -211) { --Z; }
125 else if (projPDG == 211) { ++Z; }
126 else if (projPDG == -321) { --Z; }
127 else if (projPDG == 321) { ++Z; }
128 else if (projPDG == 130) {
129 if (theSecondary->GetPDGCharge() > 0.0) { --Z; }
130 else { ++Z; }
131 } else {
132 // not ready for other projectile
133 return &theParticleChange;
134 }
135
136 // recoil nucleus
137 const G4ParticleDefinition* theRecoil = nullptr;
138 if (Z == 0 && A == 1) { theRecoil = G4Neutron::Neutron(); }
139 else if (Z == 1 && A == 1) { theRecoil = G4Proton::Proton(); }
140 else if (Z == 1 && A == 2) { theRecoil = G4Deuteron::Deuteron(); }
141 else if (Z == 1 && A == 3) { theRecoil = G4Triton::Triton(); }
142 else if (Z == 2 && A == 3) { theRecoil = G4He3::He3(); }
143 else if (Z == 2 && A == 4) { theRecoil = G4Alpha::Alpha(); }
144 else if (nist->GetIsotopeAbundance(Z, A) > 0.0) {
146 ->GetIonTable()->GetIon(Z, A, 0.0);
147 }
148
149 // check if there is enough energy for the final state
150 // and sample mass of produced state
151 const G4double mass0 = theSecondary->GetPDGMass();
152 G4double mass3 = (nullptr == theRecoil) ?
154 G4double mass2 = mass0;
155 if (isShortLived &&
156 !SampleMass(mass2, theSecondary->GetPDGWidth(), etot - mass3)) {
157 return &theParticleChange;
158 }
159
160 // not possible kinematically
161 if (etot <= mass2 + mass3) {
162 return &theParticleChange;
163 }
164
165 // sample kinematics
166 G4LorentzVector lv1(0.0, 0.0, 0.0, mass1);
167 G4LorentzVector lv = lv0 + lv1;
168 G4ThreeVector bst = lv.boostVector();
169 G4double ss = lv.mag2();
170
171 // tmax = 4*momCMS^2
172 G4double e2 = ss + mass2*mass2 - mass3*mass3;
173 G4double tmax = e2*e2/ss - 4*mass2*mass2;
174
175 G4double t = SampleT(theSecondary, A, tmax);
176
177 G4double phi = G4UniformRand()*CLHEP::twopi;
178 G4double cost = 1. - 2.0*t/tmax;
179
180 if (cost > 1.0) { cost = 1.0; }
181 else if(cost < -1.0) { cost = -1.0; }
182
183 G4double sint = std::sqrt((1.0-cost)*(1.0+cost));
184
185 if (verboseLevel>1) {
186 G4cout << " t= " << t << " tmax(GeV^2)= " << tmax/(GeV*GeV)
187 << " cos(t)=" << cost << " sin(t)=" << sint << G4endl;
188 }
189 G4double momentumCMS = 0.5*std::sqrt(tmax);
190 G4LorentzVector lv2(momentumCMS*sint*std::cos(phi),
191 momentumCMS*sint*std::sin(phi),
192 momentumCMS*cost,
193 std::sqrt(momentumCMS*momentumCMS + mass2*mass2));
194
195 // kinematics in the final state, may be a warning should be added if
196 lv2.boost(bst);
197 if (lv2.e() < mass2) {
198 lv2.setE(mass2);
199 }
200 lv -= lv2;
201 if (lv.e() < mass3) {
202 lv.setE(mass3);
203 }
204
205 // prepare secondary particles
206 theParticleChange.SetStatusChange(stopAndKill);
207 theParticleChange.SetEnergyChange(0.0);
208
209 if (!isShortLived) {
210 auto aSec = new G4DynamicParticle(theSecondary, lv2);
211 theParticleChange.AddSecondary(aSec, secID);
212 } else {
213 auto channel = theSecondary->GetDecayTable()->SelectADecayChannel(mass2);
214 auto products = channel->DecayIt(mass2);
215 G4ThreeVector bst1 = lv2.boostVector();
216 G4int N = products->entries();
217 for (G4int i=0; i<N; ++i) {
218 auto p = (*products)[i];
219 auto lvp = p->Get4Momentum();
220 lvp.boost(bst1);
221 p->Set4Momentum(lvp);
222 theParticleChange.AddSecondary(p, secID);
223 }
224 delete products;
225 }
226
227 // recoil is a stable isotope
228 if (nullptr != theRecoil) {
229 auto aRec = new G4DynamicParticle(theRecoil, lv);
230 theParticleChange.AddSecondary(aRec, secID);
231 } else {
232 // recoil is an unstable fragment
233 G4Fragment frag(A, Z, lv);
234 auto products = fHandler->BreakItUp(frag);
235 for (auto & prod : *products) {
236 auto dp = new G4DynamicParticle(prod->GetDefinition(), prod->GetMomentum());
237 theParticleChange.AddSecondary(dp, secID);
238 delete prod;
239 }
240 delete products;
241 }
242 return &theParticleChange;
243}
@ stopAndKill
CLHEP::HepLorentzVector G4LorentzVector
CLHEP::Hep3Vector G4ThreeVector
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
#define G4UniformRand()
Definition Randomize.hh:52
Hep3Vector boostVector() const
static G4Alpha * Alpha()
Definition G4Alpha.cc:83
G4double SampleT(const G4ParticleDefinition *theSec, const G4int A, const G4double tmax) const
G4VDecayChannel * SelectADecayChannel(G4double parentMass=-1.)
static G4Deuteron * Deuteron()
Definition G4Deuteron.cc:90
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
static G4He3 * He3()
Definition G4He3.cc:90
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
static G4Neutron * Neutron()
Definition G4Neutron.cc:101
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4int GetA_asInt() const
Definition G4Nucleus.hh:99
G4int GetZ_asInt() const
Definition G4Nucleus.hh:105
G4DecayTable * GetDecayTable() const
G4IonTable * GetIonTable() const
static G4ParticleTable * GetParticleTable()
static G4Proton * Proton()
Definition G4Proton.cc:90
static G4Triton * Triton()
Definition G4Triton.cc:90
virtual G4DecayProducts * DecayIt(G4double parentMass=-1.0)=0
#define N
Definition crc32.c:57

◆ operator!=()

G4bool G4ChargeExchange::operator!= ( const G4ChargeExchange & right) const
delete

◆ operator=()

const G4ChargeExchange & G4ChargeExchange::operator= ( const G4ChargeExchange & right)
delete

◆ operator==()

G4bool G4ChargeExchange::operator== ( const G4ChargeExchange & right) const
delete

◆ SampleT()

G4double G4ChargeExchange::SampleT ( const G4ParticleDefinition * theSec,
const G4int A,
const G4double tmax ) const

Definition at line 245 of file G4ChargeExchange.cc.

247{
248 G4double aa, bb, cc, dd;
249 G4Pow* g4pow = G4Pow::GetInstance();
250 if (A <= 62.) {
251 aa = g4pow->powZ(A, 1.63);
252 bb = 14.5*g4pow->powZ(A, 0.66);
253 cc = 1.4*g4pow->powZ(A, 0.33);
254 dd = 10.;
255 } else {
256 aa = g4pow->powZ(A, 1.33);
257 bb = 60.*g4pow->powZ(A, 0.33);
258 cc = 0.4*g4pow->powZ(A, 0.40);
259 dd = 10.;
260 }
261 G4double x1 = (1.0 - G4Exp(-tmax*bb))*aa/bb;
262 G4double x2 = (1.0 - G4Exp(-tmax*dd))*cc/dd;
263
264 G4double t;
265 G4double y = bb;
266 if(G4UniformRand()*(x1 + x2) < x2) y = dd;
267
268 for (G4int i=0; i<maxN; ++i) {
269 t = -G4Log(G4UniformRand())/y;
270 if (t <= tmax) { return t; }
271 }
272 return 0.0;
273}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition G4Exp.hh:180
G4double G4Log(G4double x)
Definition G4Log.hh:227
static G4Pow * GetInstance()
Definition G4Pow.cc:41
G4double powZ(G4int Z, G4double y) const
Definition G4Pow.hh:225

Referenced by ApplyYourself().


The documentation for this class was generated from the following files: