Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4LEXiZeroInelastic.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// $Id$
27//
28// Hadronic Process: XiZero Inelastic Process
29// J.L. Chuma, TRIUMF, 20-Feb-1997
30// Modified by J.L.Chuma 30-Apr-97: added originalTarget for CalculateMomenta
31
34#include "G4SystemOfUnits.hh"
35#include "Randomize.hh"
36
37void G4LEXiZeroInelastic::ModelDescription(std::ostream& outFile) const
38{
39 outFile << "G4LEXiZeroInelastic is one of the Low Energy Parameterized\n"
40 << "(LEP) models used to implement inelastic X0 scattering\n"
41 << "from nuclei. It is a re-engineered version of the GHEISHA\n"
42 << "code of H. Fesefeldt. It divides the initial collision\n"
43 << "products into backward- and forward-going clusters which are\n"
44 << "then decayed into final state hadrons. The model does not\n"
45 << "conserve energy on an event-by-event basis. It may be\n"
46 << "applied to X0 with initial energies between 0 and 25\n"
47 << "GeV.\n";
48}
49
50
53 G4Nucleus& targetNucleus)
54{
55 const G4HadProjectile *originalIncident = &aTrack;
56 if (originalIncident->GetKineticEnergy()<= 0.1*MeV) {
60 return &theParticleChange;
61 }
62
63 // create the target particle
64 G4DynamicParticle* originalTarget = targetNucleus.ReturnTargetParticle();
65
66 if (verboseLevel > 1) {
67 const G4Material *targetMaterial = aTrack.GetMaterial();
68 G4cout << "G4LEXiZeroInelastic::ApplyYourself called" << G4endl;
69 G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy()/MeV << "MeV, ";
70 G4cout << "target material = " << targetMaterial->GetName() << ", ";
71 G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName()
72 << G4endl;
73 }
74
75 // Fermi motion and evaporation
76 // As of Geant3, the Fermi energy calculation had not been done
77 G4double ek = originalIncident->GetKineticEnergy()/MeV;
78 G4double amas = originalIncident->GetDefinition()->GetPDGMass()/MeV;
79 G4ReactionProduct modifiedOriginal;
80 modifiedOriginal = *originalIncident;
81
82 G4double tkin = targetNucleus.Cinema(ek);
83 ek += tkin;
84 modifiedOriginal.SetKineticEnergy(ek*MeV);
85 G4double et = ek + amas;
86 G4double p = std::sqrt(std::abs((et-amas)*(et+amas)) );
87 G4double pp = modifiedOriginal.GetMomentum().mag()/MeV;
88 if (pp > 0.0) {
89 G4ThreeVector momentum = modifiedOriginal.GetMomentum();
90 modifiedOriginal.SetMomentum( momentum * (p/pp) );
91 }
92
93 // calculate black track energies
94 tkin = targetNucleus.EvaporationEffects(ek);
95 ek -= tkin;
96 modifiedOriginal.SetKineticEnergy(ek*MeV);
97 et = ek + amas;
98 p = std::sqrt( std::abs((et-amas)*(et+amas)) );
99 pp = modifiedOriginal.GetMomentum().mag()/MeV;
100 if (pp > 0.0) {
101 G4ThreeVector momentum = modifiedOriginal.GetMomentum();
102 modifiedOriginal.SetMomentum( momentum * (p/pp) );
103 }
104 G4ReactionProduct currentParticle = modifiedOriginal;
105 G4ReactionProduct targetParticle;
106 targetParticle = *originalTarget;
107 currentParticle.SetSide(1); // incident always goes in forward hemisphere
108 targetParticle.SetSide(-1); // target always goes in backward hemisphere
109 G4bool incidentHasChanged = false;
110 G4bool targetHasChanged = false;
111 G4bool quasiElastic = false;
112 G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec; // vec will contain the secondary particles
113 G4int vecLen = 0;
114 vec.Initialize(0);
115
116 const G4double cutOff = 0.1;
117 if (currentParticle.GetKineticEnergy()/MeV > cutOff)
118 Cascade(vec, vecLen, originalIncident, currentParticle, targetParticle,
119 incidentHasChanged, targetHasChanged, quasiElastic);
120
121 CalculateMomenta(vec, vecLen, originalIncident, originalTarget,
122 modifiedOriginal, targetNucleus, currentParticle,
123 targetParticle, incidentHasChanged, targetHasChanged,
124 quasiElastic);
125
126 SetUpChange(vec, vecLen, currentParticle, targetParticle, incidentHasChanged);
127
128 if (isotopeProduction) DoIsotopeCounting(originalIncident, targetNucleus);
129
130 delete originalTarget;
131 return &theParticleChange;
132}
133
134
135void G4LEXiZeroInelastic::Cascade(
137 G4int& vecLen,
138 const G4HadProjectile *originalIncident,
139 G4ReactionProduct &currentParticle,
140 G4ReactionProduct &targetParticle,
141 G4bool &incidentHasChanged,
142 G4bool &targetHasChanged,
143 G4bool &quasiElastic)
144{
145 // derived from original FORTRAN code CASX0 by H. Fesefeldt (20-Jan-1989)
146 //
147 // XiZero undergoes interaction with nucleon within a nucleus. Check if it is
148 // energetically possible to produce pions/kaons. In not, assume nuclear excitation
149 // occurs and input particle is degraded in energy. No other particles are produced.
150 // If reaction is possible, find the correct number of pions/protons/neutrons
151 // produced using an interpolation to multiplicity data. Replace some pions or
152 // protons/neutrons by kaons or strange baryons according to the average
153 // multiplicity per inelastic reaction.
154 //
155 const G4double mOriginal = originalIncident->GetDefinition()->GetPDGMass()/MeV;
156 const G4double etOriginal = originalIncident->GetTotalEnergy()/MeV;
157 const G4double targetMass = targetParticle.GetMass()/MeV;
158 G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
159 targetMass*targetMass +
160 2.0*targetMass*etOriginal );
161 G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
162 if( availableEnergy <= G4PionPlus::PionPlus()->GetPDGMass()/MeV )
163 {
164 quasiElastic = true;
165 return;
166 }
167 static G4bool first = true;
168 const G4int numMul = 1200;
169 const G4int numSec = 60;
170 static G4double protmul[numMul], protnorm[numSec]; // proton constants
171 static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
172 // npos = number of pi+, nneg = number of pi-, nzero = number of pi0
173 G4int counter, nt=0, npos=0, nneg=0, nzero=0;
174 G4double test;
175 const G4double c = 1.25;
176 const G4double b[] = { 0.7, 0.7 };
177 if( first ) // compute normalization constants, this will only be Done once
178 {
179 first = false;
180 G4int i;
181 for( i=0; i<numMul; ++i )protmul[i] = 0.0;
182 for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
183 counter = -1;
184 for( npos=0; npos<(numSec/3); ++npos )
185 {
186 for( nneg=std::max(0,npos-2); nneg<=(npos+1); ++nneg )
187 {
188 for( nzero=0; nzero<numSec/3; ++nzero )
189 {
190 if( ++counter < numMul )
191 {
192 nt = npos+nneg+nzero;
193 if( nt>0 && nt<=numSec )
194 {
195 protmul[counter] = Pmltpc(npos,nneg,nzero,nt,b[0],c);
196 protnorm[nt-1] += protmul[counter];
197 }
198 }
199 }
200 }
201 }
202 for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
203 for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
204 counter = -1;
205 for( npos=0; npos<numSec/3; ++npos )
206 {
207 for( nneg=std::max(0,npos-1); nneg<=(npos+2); ++nneg )
208 {
209 for( nzero=0; nzero<numSec/3; ++nzero )
210 {
211 if( ++counter < numMul )
212 {
213 nt = npos+nneg+nzero;
214 if( nt>0 && nt<=numSec )
215 {
216 neutmul[counter] = Pmltpc(npos,nneg,nzero,nt,b[1],c);
217 neutnorm[nt-1] += neutmul[counter];
218 }
219 }
220 }
221 }
222 }
223 for( i=0; i<numSec; ++i )
224 {
225 if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
226 if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
227 }
228 } // end of initialization
229
230 const G4double expxu = 82.; // upper bound for arg. of exp
231 const G4double expxl = -expxu; // lower bound for arg. of exp
237
238 // energetically possible to produce pion(s) --> inelastic scattering
239 G4double n, anpn;
240 GetNormalizationConstant(availableEnergy, n, anpn);
241 G4double ran = G4UniformRand();
242 G4double dum, excs = 0.0;
243 if( targetParticle.GetDefinition() == aProton )
244 {
245 counter = -1;
246 for( npos=0; npos<numSec/3 && ran>=excs; ++npos )
247 {
248 for( nneg=std::max(0,npos-2); nneg<=(npos+1) && ran>=excs; ++nneg )
249 {
250 for( nzero=0; nzero<numSec/3 && ran>=excs; ++nzero )
251 {
252 if( ++counter < numMul )
253 {
254 nt = npos+nneg+nzero;
255 if( nt>0 && nt<=numSec )
256 {
257 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
258 dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
259 if( std::fabs(dum) < 1.0 )
260 {
261 if( test >= 1.0e-10 )excs += dum*test;
262 }
263 else
264 excs += dum*test;
265 }
266 }
267 }
268 }
269 }
270 if( ran >= excs ) // 3 previous loops continued to the end
271 {
272 quasiElastic = true;
273 return;
274 }
275 npos--; nneg--; nzero--;
276 //
277 // number of secondary mesons determined by kno distribution
278 // check for total charge of final state mesons to determine
279 // the kind of baryons to be produced, taking into account
280 // charge and strangeness conservation
281 //
282 if( npos < nneg+1 )
283 {
284 if( npos != nneg ) // charge mismatch
285 {
286 currentParticle.SetDefinitionAndUpdateE( aSigmaPlus );
287 incidentHasChanged = true;
288 //
289 // correct the strangeness by replacing a pi- by a kaon-
290 //
291 vec.Initialize( 1 );
293 p->SetDefinition( aKaonMinus );
294 (G4UniformRand() < 0.5) ? p->SetSide( -1 ) : p->SetSide( 1 );
295 vec.SetElement( vecLen++, p );
296 --nneg;
297 }
298 }
299 else if( npos == nneg+1 )
300 {
301 if( G4UniformRand() < 0.5 )
302 {
303 targetParticle.SetDefinitionAndUpdateE( aNeutron );
304 targetHasChanged = true;
305 }
306 else
307 {
308 currentParticle.SetDefinitionAndUpdateE( aXiMinus );
309 incidentHasChanged = true;
310 }
311 }
312 else
313 {
314 currentParticle.SetDefinitionAndUpdateE( aXiMinus );
315 incidentHasChanged = true;
316 targetParticle.SetDefinitionAndUpdateE( aNeutron );
317 targetHasChanged = true;
318 }
319 }
320 else // target must be a neutron
321 {
322 counter = -1;
323 for( npos=0; npos<numSec/3 && ran>=excs; ++npos )
324 {
325 for( nneg=std::max(0,npos-1); nneg<=(npos+2) && ran>=excs; ++nneg )
326 {
327 for( nzero=0; nzero<numSec/3 && ran>=excs; ++nzero )
328 {
329 if( ++counter < numMul )
330 {
331 nt = npos+nneg+nzero;
332 if( nt>0 && nt<=numSec )
333 {
334 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
335 dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
336 if( std::fabs(dum) < 1.0 )
337 {
338 if( test >= 1.0e-10 )excs += dum*test;
339 }
340 else
341 excs += dum*test;
342 }
343 }
344 }
345 }
346 }
347 if( ran >= excs ) // 3 previous loops continued to the end
348 {
349 quasiElastic = true;
350 return;
351 }
352 npos--; nneg--; nzero--;
353 if( npos < nneg )
354 {
355 if( npos+1 == nneg )
356 {
357 targetParticle.SetDefinitionAndUpdateE( aProton );
358 targetHasChanged = true;
359 }
360 else // charge mismatch
361 {
362 currentParticle.SetDefinitionAndUpdateE( aSigmaPlus );
363 incidentHasChanged = true;
364 targetParticle.SetDefinitionAndUpdateE( aProton );
365 targetHasChanged = true;
366 //
367 // correct the strangeness by replacing a pi- by a kaon-
368 //
369 vec.Initialize( 1 );
371 p->SetDefinition( aKaonMinus );
372 (G4UniformRand() < 0.5) ? p->SetSide( -1 ) : p->SetSide( 1 );
373 vec.SetElement( vecLen++, p );
374 --nneg;
375 }
376 }
377 else if( npos == nneg )
378 {
379 if( G4UniformRand() >= 0.5 )
380 {
381 currentParticle.SetDefinitionAndUpdateE( aXiMinus );
382 incidentHasChanged = true;
383 targetParticle.SetDefinitionAndUpdateE( aProton );
384 targetHasChanged = true;
385 }
386 }
387 else
388 {
389 currentParticle.SetDefinitionAndUpdateE( aXiMinus );
390 incidentHasChanged = true;
391 }
392 }
393 SetUpPions(npos, nneg, nzero, vec, vecLen);
394 return;
395}
396
397 /* end of file */
398
@ isAlive
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
Hep3Vector unit() const
double mag() const
Hep3Vector vect() const
G4ParticleDefinition * GetDefinition() const
void SetElement(G4int anIndex, Type *anElement)
Definition: G4FastVector.hh:76
void Initialize(G4int items)
Definition: G4FastVector.hh:63
void SetStatusChange(G4HadFinalStateStatus aS)
void SetEnergyChange(G4double anEnergy)
void SetMomentumChange(const G4ThreeVector &aV)
const G4Material * GetMaterial() const
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
G4double GetTotalEnergy() const
G4double Pmltpc(G4int np, G4int nm, G4int nz, G4int n, G4double b, G4double c)
void GetNormalizationConstant(const G4double availableEnergy, G4double &n, G4double &anpn)
void SetUpPions(const G4int np, const G4int nm, const G4int nz, G4FastVector< G4ReactionProduct, GHADLISTSIZE > &vec, G4int &vecLen)
void CalculateMomenta(G4FastVector< G4ReactionProduct, GHADLISTSIZE > &vec, G4int &vecLen, const G4HadProjectile *originalIncident, const G4DynamicParticle *originalTarget, G4ReactionProduct &modifiedOriginal, G4Nucleus &targetNucleus, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged, G4bool &targetHasChanged, G4bool quasiElastic)
void DoIsotopeCounting(const G4HadProjectile *theProjectile, const G4Nucleus &aNucleus)
void SetUpChange(G4FastVector< G4ReactionProduct, GHADLISTSIZE > &vec, G4int &vecLen, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged)
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:113
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
virtual void ModelDescription(std::ostream &outFile) const
const G4String & GetName() const
Definition: G4Material.hh:177
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4double EvaporationEffects(G4double kineticEnergy)
Definition: G4Nucleus.cc:264
G4double Cinema(G4double kineticEnergy)
Definition: G4Nucleus.cc:368
G4DynamicParticle * ReturnTargetParticle() const
Definition: G4Nucleus.cc:227
const G4String & GetParticleName() const
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
static G4Proton * Proton()
Definition: G4Proton.cc:93
void SetMomentum(const G4double x, const G4double y, const G4double z)
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
void SetSide(const G4int sid)
void SetDefinitionAndUpdateE(G4ParticleDefinition *aParticleDefinition)
void SetKineticEnergy(const G4double en)
G4ParticleDefinition * GetDefinition() const
void SetDefinition(G4ParticleDefinition *aParticleDefinition)
G4double GetMass() const
static G4SigmaPlus * SigmaPlus()
Definition: G4SigmaPlus.cc:108
static G4XiMinus * XiMinus()
Definition: G4XiMinus.cc:106
const G4double pi