Geant4 11.1.1
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 ()
 
 ~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
 
void SetLowestEnergyLimit (G4double value)
 
G4double SampleT (G4double p, G4int A)
 
- Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
 
virtual ~G4HadronicInteraction ()
 
virtual G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
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 ( )
explicit

Definition at line 56 of file G4ChargeExchange.cc.

56 : G4HadronicInteraction("Charge Exchange"), secID(-1)
57{
58 SetMinEnergy( 0.0*GeV );
60
61 lowestEnergyLimit = 1.*MeV;
62
63 theProton = G4Proton::Proton();
64 theNeutron = G4Neutron::Neutron();
65 theAProton = G4AntiProton::AntiProton();
66 theANeutron = G4AntiNeutron::AntiNeutron();
67 thePiPlus = G4PionPlus::PionPlus();
68 thePiMinus = G4PionMinus::PionMinus();
69 thePiZero = G4PionZero::PionZero();
70 theKPlus = G4KaonPlus::KaonPlus();
71 theKMinus = G4KaonMinus::KaonMinus();
74 theL = G4Lambda::Lambda();
75 theAntiL = G4AntiLambda::AntiLambda();
76 theSPlus = G4SigmaPlus::SigmaPlus();
78 theSMinus = G4SigmaMinus::SigmaMinus();
80 theS0 = G4SigmaZero::SigmaZero();
82 theXiMinus = G4XiMinus::XiMinus();
83 theXi0 = G4XiZero::XiZero();
84 theAXiMinus = G4AntiXiMinus::AntiXiMinus();
85 theAXi0 = G4AntiXiZero::AntiXiZero();
86 theOmega = G4OmegaMinus::OmegaMinus();
88 theD = G4Deuteron::Deuteron();
89 theT = G4Triton::Triton();
90 theA = G4Alpha::Alpha();
91 theHe3 = G4He3::He3();
92
93 secID = G4PhysicsModelCatalog::GetModelID( "model_ChargeExchange" );
94}
static G4Alpha * Alpha()
Definition: G4Alpha.cc:88
static G4AntiLambda * AntiLambda()
static G4AntiNeutron * AntiNeutron()
static G4AntiOmegaMinus * AntiOmegaMinus()
static G4AntiProton * AntiProton()
Definition: G4AntiProton.cc:92
static G4AntiSigmaMinus * AntiSigmaMinus()
static G4AntiSigmaPlus * AntiSigmaPlus()
static G4AntiSigmaZero * AntiSigmaZero()
static G4AntiXiMinus * AntiXiMinus()
static G4AntiXiZero * AntiXiZero()
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:93
void SetMinEnergy(G4double anEnergy)
void SetMaxEnergy(const G4double anEnergy)
static G4HadronicParameters * Instance()
static G4He3 * He3()
Definition: G4He3.cc:93
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:112
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:112
static G4KaonZeroLong * KaonZeroLong()
static G4KaonZeroShort * KaonZeroShort()
static G4Lambda * Lambda()
Definition: G4Lambda.cc:107
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
static G4OmegaMinus * OmegaMinus()
static G4int GetModelID(const G4int modelIndex)
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:97
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:97
static G4PionZero * PionZero()
Definition: G4PionZero.cc:107
static G4Proton * Proton()
Definition: G4Proton.cc:92
static G4SigmaMinus * SigmaMinus()
static G4SigmaPlus * SigmaPlus()
Definition: G4SigmaPlus.cc:107
static G4SigmaZero * SigmaZero()
Definition: G4SigmaZero.cc:101
static G4Triton * Triton()
Definition: G4Triton.cc:93
static G4XiMinus * XiMinus()
Definition: G4XiMinus.cc:105
static G4XiZero * XiZero()
Definition: G4XiZero.cc:105

◆ ~G4ChargeExchange()

G4ChargeExchange::~G4ChargeExchange ( )
override

Definition at line 96 of file G4ChargeExchange.cc.

97{}

◆ 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 99 of file G4ChargeExchange.cc.

101{
103 const G4HadProjectile* aParticle = &aTrack;
104 G4double ekin = aParticle->GetKineticEnergy();
105
106 G4int A = targetNucleus.GetA_asInt();
107 G4int Z = targetNucleus.GetZ_asInt();
108
109 if(ekin <= lowestEnergyLimit || A < 3) {
112 return &theParticleChange;
113 }
114
115 G4double plab = aParticle->GetTotalMomentum();
116
117 if (verboseLevel > 1)
118 G4cout << "G4ChargeExchange::DoIt: Incident particle plab="
119 << plab/GeV << " GeV/c "
120 << " ekin(MeV) = " << ekin/MeV << " "
121 << aParticle->GetDefinition()->GetParticleName() << G4endl;
122
123 // Scattered particle referred to axis of incident particle
124 const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
125
126 G4int N = A - Z;
127 G4int projPDG = theParticle->GetPDGEncoding();
128 if (verboseLevel > 1)
129 G4cout << "G4ChargeExchange for " << theParticle->GetParticleName()
130 << " PDGcode= " << projPDG << " on nucleus Z= " << Z
131 << " A= " << A << " N= " << N
132 << G4endl;
133
134 const G4ParticleDefinition* theDef = nullptr;
135
137 G4LorentzVector lv1 = aParticle->Get4Momentum();
138 G4LorentzVector lv0(0.0,0.0,0.0,mass2);
139
140 G4LorentzVector lv = lv0 + lv1;
141 G4ThreeVector bst = lv.boostVector();
142 lv1.boost(-bst);
143 lv0.boost(-bst);
144
145 // Sample final particles
146 G4bool theHyperon = false;
147 const G4ParticleDefinition* theRecoil = nullptr;
148 const G4ParticleDefinition* theSecondary = nullptr;
149
150 if(theParticle == theProton) {
151 theSecondary = theNeutron;
152 Z++;
153 } else if(theParticle == theNeutron) {
154 theSecondary = theProton;
155 Z--;
156 } else if(theParticle == thePiPlus) {
157 theSecondary = thePiZero;
158 Z++;
159 } else if(theParticle == thePiMinus) {
160 theSecondary = thePiZero;
161 Z--;
162 } else if(theParticle == theKPlus) {
163 if(G4UniformRand()<0.5) theSecondary = theK0S;
164 else theSecondary = theK0L;
165 Z++;
166 } else if(theParticle == theKMinus) {
167 if(G4UniformRand()<0.5) theSecondary = theK0S;
168 else theSecondary = theK0L;
169 Z--;
170 } else if(theParticle == theK0S || theParticle == theK0L) {
171 if(G4UniformRand()*A < G4double(Z)) {
172 theSecondary = theKPlus;
173 Z--;
174 } else {
175 theSecondary = theKMinus;
176 Z++;
177 }
178 } else if(theParticle == theANeutron) {
179 theSecondary = theAProton;
180 Z++;
181 } else if(theParticle == theAProton) {
182 theSecondary = theANeutron;
183 Z--;
184 } else if(theParticle == theL) {
186 if(G4UniformRand()*A < G4double(Z)) {
187 if(x < 0.2) {
188 theSecondary = theS0;
189 } else if (x < 0.4) {
190 theSecondary = theSPlus;
191 Z--;
192 } else if (x < 0.6) {
193 theSecondary = theProton;
194 theRecoil = theL;
195 theHyperon = true;
196 A--;
197 } else if (x < 0.8) {
198 theSecondary = theProton;
199 theRecoil = theS0;
200 theHyperon = true;
201 A--;
202 } else {
203 theSecondary = theNeutron;
204 theRecoil = theSPlus;
205 theHyperon = true;
206 A--;
207 }
208 } else {
209 if(x < 0.2) {
210 theSecondary = theS0;
211 } else if (x < 0.4) {
212 theSecondary = theSMinus;
213 Z++;
214 } else if (x < 0.6) {
215 theSecondary = theNeutron;
216 theRecoil = theL;
217 A--;
218 theHyperon = true;
219 } else if (x < 0.8) {
220 theSecondary = theNeutron;
221 theRecoil = theS0;
222 theHyperon = true;
223 A--;
224 } else {
225 theSecondary = theProton;
226 theRecoil = theSMinus;
227 theHyperon = true;
228 A--;
229 }
230 }
231 }
232
233 if (Z == 1 && A == 2) theDef = theD;
234 else if (Z == 1 && A == 3) theDef = theT;
235 else if (Z == 2 && A == 3) theDef = theHe3;
236 else if (Z == 2 && A == 4) theDef = theA;
237 else {
238 theDef =
240 }
241 if(!theSecondary) { return &theParticleChange; }
242
243 G4double m11 = theSecondary->GetPDGMass();
244 G4double m21 = theDef->GetPDGMass();
245 if(theRecoil) { m21 += theRecoil->GetPDGMass(); }
246 else { theRecoil = theDef; }
247
248 G4double etot = lv0.e() + lv1.e();
249
250 // kinematiacally impossible
251 if(etot < m11 + m21) {
254 return &theParticleChange;
255 }
256
257 G4ThreeVector p1 = lv1.vect();
258 G4double e1 = 0.5*etot*(1.0 - (m21*m21 - m11*m11)/(etot*etot));
259 // G4double e2 = etot - e1;
260 G4double ptot = std::sqrt(e1*e1 - m11*m11);
261
262 G4double tmax = 4.0*ptot*ptot;
263 G4double g2 = GeV*GeV;
264
265 G4double t = g2*SampleT(tmax/g2, A);
266
267 if(verboseLevel>1) {
268 G4cout <<"## G4ChargeExchange t= " << t << " tmax= " << tmax
269 << " ptot= " << ptot << G4endl;
270 }
271 // Sampling in CM system
272 G4double phi = G4UniformRand()*twopi;
273 G4double cost = 1. - 2.0*t/tmax;
274 if(std::abs(cost) > 1.0) cost = 1.0;
275 G4double sint = std::sqrt((1.0-cost)*(1.0+cost));
276
277 //if (verboseLevel > 1)
278 // G4cout << "cos(t)=" << cost << " std::sin(t)=" << sint << G4endl;
279
280 G4ThreeVector v1(sint*std::cos(phi),sint*std::sin(phi),cost);
281 v1 *= ptot;
282 G4LorentzVector nlv1(v1.x(),v1.y(),v1.z(),e1);
283 G4LorentzVector nlv0 = lv0 + lv1 - nlv1;
284
285 nlv0.boost(bst);
286 nlv1.boost(bst);
287
290 G4DynamicParticle * aSec = new G4DynamicParticle(theSecondary, nlv1);
291 theParticleChange.AddSecondary(aSec, secID);
292
293 G4double erec = std::max(nlv0.e() - m21, 0.0);
294
295 //G4cout << "erec= " <<erec << " Esec= " << aSec->GetKineticEnergy() << G4endl;
296
297 if(theHyperon) {
299 aSec = new G4DynamicParticle();
300 aSec->SetDefinition(theRecoil);
301 aSec->SetKineticEnergy(0.0);
302 } else if(erec > GetRecoilEnergyThreshold()) {
303 aSec = new G4DynamicParticle(theRecoil, nlv0);
304 theParticleChange.AddSecondary(aSec, secID);
305 } else {
307 }
308 return &theParticleChange;
309}
@ stopAndKill
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
G4double SampleT(G4double p, G4int A)
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
void SetKineticEnergy(G4double aEnergy)
void SetStatusChange(G4HadFinalStateStatus aS)
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
void SetEnergyChange(G4double anEnergy)
void SetMomentumChange(const G4ThreeVector &aV)
void SetLocalEnergyDeposit(G4double aE)
G4double GetTotalMomentum() const
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
G4double GetRecoilEnergyThreshold() const
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:522
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
const G4String & GetParticleName() const
G4IonTable * GetIonTable() const
static G4ParticleTable * GetParticleTable()
#define N
Definition: crc32.c:56

◆ 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 ( G4double  p,
G4int  A 
)

Definition at line 311 of file G4ChargeExchange.cc.

312{
313 G4double aa, bb, cc, dd;
314 G4Pow* g4pow = G4Pow::GetInstance();
315 if (A <= 62.) {
316 aa = g4pow->powZ(A, 1.63);
317 bb = 14.5*g4pow->powZ(A, 0.66);
318 cc = 1.4*g4pow->powZ(A, 0.33);
319 dd = 10.;
320 } else {
321 aa = g4pow->powZ(A, 1.33);
322 bb = 60.*g4pow->powZ(A, 0.33);
323 cc = 0.4*g4pow->powZ(A, 0.40);
324 dd = 10.;
325 }
326 G4double x1 = (1.0 - G4Exp(-tmax*bb))*aa/bb;
327 G4double x2 = (1.0 - G4Exp(-tmax*dd))*cc/dd;
328
329 G4double t;
330 G4double y = bb;
331 if(G4UniformRand()*(x1 + x2) < x2) y = dd;
332
333 const G4int maxNumberOfLoops = 10000;
334 G4int loopCounter = 0;
335 do {
336 t = -G4Log(G4UniformRand())/y;
337 } while ( (t > tmax) &&
338 ++loopCounter < maxNumberOfLoops ); /* Loop checking, 10.08.2015, A.Ribon */
339 if ( loopCounter >= maxNumberOfLoops ) {
340 t = 0.0;
341 }
342 return t;
343}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
G4double G4Log(G4double x)
Definition: G4Log.hh:227
Definition: G4Pow.hh:49
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double powZ(G4int Z, G4double y) const
Definition: G4Pow.hh:225

Referenced by ApplyYourself().

◆ SetLowestEnergyLimit()

void G4ChargeExchange::SetLowestEnergyLimit ( G4double  value)
inline

Definition at line 107 of file G4ChargeExchange.hh.

108{
109 lowestEnergyLimit = value;
110}

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