Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4ChargeExchange Class Reference

#include <G4ChargeExchange.hh>

Inheritance diagram for G4ChargeExchange:
G4HadronicInteraction

Public Member Functions

 G4ChargeExchange ()
 
virtual ~G4ChargeExchange ()
 
virtual G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
void SetLowestEnergyLimit (G4double value)
 
void SetRecoilKinEnergyLimit (G4double value)
 
G4double SampleT (G4double p, G4double A)
 
- 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 &, G4Nucleus &)
 
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)
 
const G4HadronicInteractionGetMyPointer () const
 
virtual G4int GetVerboseLevel () const
 
virtual 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
 
G4bool operator== (const G4HadronicInteraction &right) const
 
G4bool operator!= (const G4HadronicInteraction &right) const
 
virtual const std::pair
< G4double, G4double
GetFatalEnergyCheckLevels () const
 
virtual std::pair< G4double,
G4double
GetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
virtual void ModelDescription (std::ostream &outFile) const
 

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 52 of file G4ChargeExchange.hh.

Constructor & Destructor Documentation

G4ChargeExchange::G4ChargeExchange ( )

Definition at line 48 of file G4ChargeExchange.cc.

References G4Alpha::Alpha(), G4AntiLambda::AntiLambda(), G4AntiNeutron::AntiNeutron(), G4AntiOmegaMinus::AntiOmegaMinus(), G4AntiProton::AntiProton(), G4AntiSigmaMinus::AntiSigmaMinus(), G4AntiSigmaPlus::AntiSigmaPlus(), G4AntiSigmaZero::AntiSigmaZero(), G4AntiXiMinus::AntiXiMinus(), G4AntiXiZero::AntiXiZero(), G4Deuteron::Deuteron(), python.hepunit::GeV, G4He3::He3(), G4KaonMinus::KaonMinus(), G4KaonPlus::KaonPlus(), G4KaonZeroLong::KaonZeroLong(), G4KaonZeroShort::KaonZeroShort(), python.hepunit::keV, G4Lambda::Lambda(), python.hepunit::MeV, G4Neutron::Neutron(), G4OmegaMinus::OmegaMinus(), G4PionMinus::PionMinus(), G4PionPlus::PionPlus(), G4PionZero::PionZero(), G4Proton::Proton(), G4HadronicInteraction::SetMaxEnergy(), G4HadronicInteraction::SetMinEnergy(), G4SigmaMinus::SigmaMinus(), G4SigmaPlus::SigmaPlus(), G4SigmaZero::SigmaZero(), python.hepunit::TeV, G4Triton::Triton(), G4XiMinus::XiMinus(), and G4XiZero::XiZero().

48  : G4HadronicInteraction("Charge Exchange")
49 {
50  SetMinEnergy( 0.0*GeV );
51  SetMaxEnergy( 100.*TeV );
52 
53  lowEnergyRecoilLimit = 100.*keV;
54  lowestEnergyLimit = 1.*MeV;
55 
56  theProton = G4Proton::Proton();
57  theNeutron = G4Neutron::Neutron();
58  theAProton = G4AntiProton::AntiProton();
59  theANeutron = G4AntiNeutron::AntiNeutron();
60  thePiPlus = G4PionPlus::PionPlus();
61  thePiMinus = G4PionMinus::PionMinus();
62  thePiZero = G4PionZero::PionZero();
63  theKPlus = G4KaonPlus::KaonPlus();
64  theKMinus = G4KaonMinus::KaonMinus();
67  theL = G4Lambda::Lambda();
68  theAntiL = G4AntiLambda::AntiLambda();
69  theSPlus = G4SigmaPlus::SigmaPlus();
70  theASPlus = G4AntiSigmaPlus::AntiSigmaPlus();
71  theSMinus = G4SigmaMinus::SigmaMinus();
72  theASMinus = G4AntiSigmaMinus::AntiSigmaMinus();
73  theS0 = G4SigmaZero::SigmaZero();
75  theXiMinus = G4XiMinus::XiMinus();
76  theXi0 = G4XiZero::XiZero();
77  theAXiMinus = G4AntiXiMinus::AntiXiMinus();
78  theAXi0 = G4AntiXiZero::AntiXiZero();
79  theOmega = G4OmegaMinus::OmegaMinus();
81  theD = G4Deuteron::Deuteron();
82  theT = G4Triton::Triton();
83  theA = G4Alpha::Alpha();
84  theHe3 = G4He3::He3();
85 }
static G4AntiOmegaMinus * AntiOmegaMinus()
static G4OmegaMinus * OmegaMinus()
static G4KaonZeroLong * KaonZeroLong()
static G4AntiSigmaPlus * AntiSigmaPlus()
static G4SigmaZero * SigmaZero()
Definition: G4SigmaZero.cc:99
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:113
static G4AntiSigmaMinus * AntiSigmaMinus()
static G4XiZero * XiZero()
Definition: G4XiZero.cc:106
void SetMinEnergy(G4double anEnergy)
static G4KaonZeroShort * KaonZeroShort()
static G4AntiProton * AntiProton()
Definition: G4AntiProton.cc:93
static G4XiMinus * XiMinus()
Definition: G4XiMinus.cc:106
static G4AntiXiMinus * AntiXiMinus()
static G4Triton * Triton()
Definition: G4Triton.cc:95
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
G4HadronicInteraction(const G4String &modelName="HadronicModel")
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4PionZero * PionZero()
Definition: G4PionZero.cc:104
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
static G4SigmaMinus * SigmaMinus()
static G4AntiLambda * AntiLambda()
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
static G4AntiSigmaZero * AntiSigmaZero()
void SetMaxEnergy(const G4double anEnergy)
static G4AntiXiZero * AntiXiZero()
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
static G4SigmaPlus * SigmaPlus()
Definition: G4SigmaPlus.cc:108
static G4Lambda * Lambda()
Definition: G4Lambda.cc:108
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:113
static G4He3 * He3()
Definition: G4He3.cc:94
static G4AntiNeutron * AntiNeutron()
G4ChargeExchange::~G4ChargeExchange ( )
virtual

Definition at line 87 of file G4ChargeExchange.cc.

88 {}

Member Function Documentation

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

Implements G4HadronicInteraction.

Definition at line 90 of file G4ChargeExchange.cc.

References G4HadFinalState::AddSecondary(), CLHEP::HepLorentzVector::boost(), CLHEP::HepLorentzVector::boostVector(), G4HadFinalState::Clear(), CLHEP::HepLorentzVector::e(), G4cout, G4endl, G4UniformRand, G4HadProjectile::Get4Momentum(), G4Nucleus::GetA_asInt(), G4HadProjectile::GetDefinition(), G4IonTable::GetIon(), G4ParticleTable::GetIonTable(), G4HadProjectile::GetKineticEnergy(), G4NucleiProperties::GetNuclearMass(), G4ParticleDefinition::GetParticleName(), G4ParticleTable::GetParticleTable(), G4ParticleDefinition::GetPDGEncoding(), G4ParticleDefinition::GetPDGMass(), G4HadProjectile::GetTotalMomentum(), G4Nucleus::GetZ_asInt(), python.hepunit::GeV, python.hepunit::MeV, N, SampleT(), G4DynamicParticle::SetDefinition(), G4HadFinalState::SetEnergyChange(), G4DynamicParticle::SetKineticEnergy(), G4HadFinalState::SetLocalEnergyDeposit(), G4HadFinalState::SetMomentumChange(), G4HadFinalState::SetStatusChange(), stopAndKill, G4HadronicInteraction::theParticleChange, python.hepunit::twopi, CLHEP::Hep3Vector::unit(), CLHEP::HepLorentzVector::vect(), G4HadronicInteraction::verboseLevel, test::x, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

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

Definition at line 303 of file G4ChargeExchange.cc.

References G4UniformRand.

Referenced by ApplyYourself().

304 {
305  G4double aa, bb, cc, dd;
306  if (A <= 62.) {
307  aa = std::pow(A, 1.63);
308  bb = 14.5*std::pow(A, 0.66);
309  cc = 1.4*std::pow(A, 0.33);
310  dd = 10.;
311  } else {
312  aa = std::pow(A, 1.33);
313  bb = 60.*std::pow(A, 0.33);
314  cc = 0.4*std::pow(A, 0.40);
315  dd = 10.;
316  }
317  G4double x1 = (1.0 - std::exp(-tmax*bb))*aa/bb;
318  G4double x2 = (1.0 - std::exp(-tmax*dd))*cc/dd;
319 
320  G4double t;
321  G4double y = bb;
322  if(G4UniformRand()*(x1 + x2) < x2) y = dd;
323 
324  do {t = -std::log(G4UniformRand())/y;} while (t > tmax);
325 
326  return t;
327 }
#define G4UniformRand()
Definition: Randomize.hh:87
double G4double
Definition: G4Types.hh:76
void G4ChargeExchange::SetLowestEnergyLimit ( G4double  value)
inline

Definition at line 112 of file G4ChargeExchange.hh.

113 {
114  lowestEnergyLimit = value;
115 }
const XML_Char int const XML_Char * value
void G4ChargeExchange::SetRecoilKinEnergyLimit ( G4double  value)
inline

Definition at line 107 of file G4ChargeExchange.hh.

108 {
109  lowEnergyRecoilLimit = value;
110 }
const XML_Char int const XML_Char * value

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