Geant4 11.2.2
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=default
 
 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 49 of file G4ChargeExchange.hh.

Constructor & Destructor Documentation

◆ G4ChargeExchange() [1/2]

G4ChargeExchange::G4ChargeExchange ( G4ChargeExchangeXS * ptr)
explicit

Definition at line 54 of file G4ChargeExchange.cc.

55 : G4HadronicInteraction("ChargeExchange"),
56 fXSection(ptr)
57{
58 lowEnergyLimit = 1.*CLHEP::MeV;
59 secID = G4PhysicsModelCatalog::GetModelID( "model_ChargeExchange" );
60}
G4HadronicInteraction(const G4String &modelName="HadronicModel")
static G4int GetModelID(const G4int modelIndex)

◆ ~G4ChargeExchange()

G4ChargeExchange::~G4ChargeExchange ( )
overridedefault

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

64{
66 auto part = aTrack.GetDefinition();
67 G4int pdg = part->GetPDGEncoding();
68 G4double ekin = aTrack.GetKineticEnergy();
69
70 G4int A = targetNucleus.GetA_asInt();
71 G4int Z = targetNucleus.GetZ_asInt();
72
73 if (ekin <= lowEnergyLimit) {
74 return &theParticleChange;
75 }
76
77 G4int projPDG = part->GetPDGEncoding();
78 if (verboseLevel > 1)
79 G4cout << "G4ChargeExchange for " << part->GetParticleName()
80 << " PDGcode= " << projPDG << " on nucleus Z= " << Z
81 << " A= " << A << " N= " << A - Z
82 << G4endl;
83
85 G4LorentzVector lv0 = aTrack.Get4Momentum();
86 G4double etot = mass1 + lv0.e();
87
88 // select final state
89 const G4ParticleDefinition* theRecoil = nullptr;
90 const G4ParticleDefinition* theSecondary =
91 fXSection->SampleSecondaryType(part, Z, A);
92
93 // atomic number of the recoil nucleus
94 if (pdg == -211) { --Z; }
95 else if (pdg == 211) { ++Z; }
96 else if (pdg == -321) { --Z; }
97 else if (pdg == 321) { ++Z; }
98 else if (pdg == 130) {
99 if (theSecondary->GetPDGCharge() > 0.0) { --Z; }
100 else { ++Z; }
101 } else {
102 // not ready for other projectile
103 return &theParticleChange;
104 }
105
106 if (Z == 0 && A == 1) theRecoil = G4Neutron::Neutron();
107 else if (Z == 1 && A == 1) theRecoil = G4Proton::Proton();
108 else if (Z == 1 && A == 2) theRecoil = G4Deuteron::Deuteron();
109 else if (Z == 1 && A == 3) theRecoil = G4Triton::Triton();
110 else if (Z == 2 && A == 3) theRecoil = G4He3::He3();
111 else if (Z == 2 && A == 4) theRecoil = G4Alpha::Alpha();
112 else {
114 ->GetIonTable()->GetIon(Z, A, 0.0);
115 }
116 if (nullptr == theRecoil) { return &theParticleChange; }
117
118 G4double mass2 = theSecondary->GetPDGMass();
119 G4double mass3 = theRecoil->GetPDGMass();
120
121 if (etot <= mass2 + mass3) {
122 // not possible kinematically
123 return &theParticleChange;
124 }
125
126 // sample kinematics
127 G4LorentzVector lv1(0.0, 0.0, 0.0, mass1);
128 G4LorentzVector lv = lv0 + lv1;
129 G4ThreeVector bst = lv.boostVector();
130 G4double ss = lv.mag2();
131
132 // tmax = 4*momCMS^2
133 G4double e2 = ss + mass2*mass2 - mass3*mass3;
134 G4double tmax = e2*e2/ss - 4*mass2*mass2;
135
136 G4double t = SampleT(theSecondary, A, tmax);
137
138 G4double phi = G4UniformRand()*CLHEP::twopi;
139 G4double cost = 1. - 2.0*t/tmax;
140
141 if (cost > 1.0) { cost = 1.0; }
142 else if(cost < -1.0) { cost = -1.0; }
143
144 G4double sint = std::sqrt((1.0-cost)*(1.0+cost));
145
146 if (verboseLevel>1) {
147 G4cout << " t= " << t << " tmax(GeV^2)= " << tmax/(GeV*GeV)
148 << " cos(t)=" << cost << " sin(t)=" << sint << G4endl;
149 }
150 G4double momentumCMS = 0.5*std::sqrt(tmax);
151 G4LorentzVector lv2(momentumCMS*sint*std::cos(phi),
152 momentumCMS*sint*std::sin(phi),
153 momentumCMS*cost,
154 std::sqrt(momentumCMS*momentumCMS + mass2*mass2));
155
156 // kinematics in the final stae, may be a warning should be added if
157 lv2.boost(bst);
158 if (lv2.e() < mass2) {
159 lv2.setE(mass2);
160 }
161 lv -= lv2;
162 if (lv.e() < mass3) {
163 lv.setE(mass3);
164 }
165
166 // prepare secondary particles
169
170 G4DynamicParticle * aSec = new G4DynamicParticle(theSecondary, lv2);
171 theParticleChange.AddSecondary(aSec, secID);
172
173 G4DynamicParticle * aRec = new G4DynamicParticle(theRecoil, lv);
174 theParticleChange.AddSecondary(aRec, secID);
175 return &theParticleChange;
176}
@ stopAndKill
double G4double
Definition G4Types.hh:83
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
const G4ParticleDefinition * SampleSecondaryType(const G4ParticleDefinition *, const G4int Z, const G4int A)
G4double SampleT(const G4ParticleDefinition *theSec, const G4int A, const G4double tmax) const
static G4Deuteron * Deuteron()
Definition G4Deuteron.cc:90
void SetStatusChange(G4HadFinalStateStatus aS)
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
void SetEnergyChange(G4double anEnergy)
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
G4IonTable * GetIonTable() const
static G4ParticleTable * GetParticleTable()
static G4Proton * Proton()
Definition G4Proton.cc:90
static G4Triton * Triton()
Definition G4Triton.cc:90

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

180{
181 G4double aa, bb, cc, dd;
182 G4Pow* g4pow = G4Pow::GetInstance();
183 if (A <= 62.) {
184 aa = g4pow->powZ(A, 1.63);
185 bb = 14.5*g4pow->powZ(A, 0.66);
186 cc = 1.4*g4pow->powZ(A, 0.33);
187 dd = 10.;
188 } else {
189 aa = g4pow->powZ(A, 1.33);
190 bb = 60.*g4pow->powZ(A, 0.33);
191 cc = 0.4*g4pow->powZ(A, 0.40);
192 dd = 10.;
193 }
194 G4double x1 = (1.0 - G4Exp(-tmax*bb))*aa/bb;
195 G4double x2 = (1.0 - G4Exp(-tmax*dd))*cc/dd;
196
197 G4double t;
198 G4double y = bb;
199 if(G4UniformRand()*(x1 + x2) < x2) y = dd;
200
201 const G4int maxN = 10000;
202 G4int count = 0;
203 do {
204 t = -G4Log(G4UniformRand())/y;
205 } while ( (t > tmax) && ++count < maxN );
206 /* Loop checking, 10.08.2015, A.Ribon */
207 if ( count >= maxN ) {
208 t = 0.0;
209 }
210 return t;
211}
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().


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