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
G4LEXiMinusInelastic.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: XiMinus 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 G4LEXiMinusInelastic::ModelDescription(std::ostream& outFile) const
38{
39 outFile << "G4LEXiMinusInelastic is one of the Low Energy Parameterized\n"
40 << "(LEP) models used to implement inelastic Xi- 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 Xi- 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 << "G4LEXiMinusInelastic::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 G4LEXiMinusInelastic::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 CASXM by H. Fesefeldt (17-Jan-1989)
146 //
147 // XiMinus 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-1); 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=npos; 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 //
240 G4double n, anpn;
241 GetNormalizationConstant( availableEnergy, n, anpn );
242 G4double ran = G4UniformRand();
243 G4double dum, excs = 0.0;
244 if( targetParticle.GetDefinition() == aProton )
245 {
246 counter = -1;
247 for( npos=0; npos<numSec/3 && ran>=excs; ++npos )
248 {
249 for( nneg=std::max(0,npos-1); nneg<=(npos+1) && ran>=excs; ++nneg )
250 {
251 for( nzero=0; nzero<numSec/3 && ran>=excs; ++nzero )
252 {
253 if( ++counter < numMul )
254 {
255 nt = npos+nneg+nzero;
256 if( nt>0 && nt<=numSec )
257 {
258 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
259 dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
260 if( std::fabs(dum) < 1.0 )
261 {
262 if( test >= 1.0e-10 )excs += dum*test;
263 }
264 else
265 excs += dum*test;
266 }
267 }
268 }
269 }
270 }
271 if( ran >= excs ) // 3 previous loops continued to the end
272 {
273 quasiElastic = true;
274 return;
275 }
276 npos--; nneg--; nzero--;
277 //
278 // number of secondary mesons determined by kno distribution
279 // check for total charge of final state mesons to determine
280 // the kind of baryons to be produced, taking into account
281 // charge and strangeness conservation
282 //
283 if( npos < nneg )
284 {
285 if( npos+1 == nneg )
286 {
287 currentParticle.SetDefinitionAndUpdateE( aXiZero );
288 incidentHasChanged = true;
289 }
290 else // charge mismatch
291 {
292 currentParticle.SetDefinitionAndUpdateE( aSigmaPlus );
293 incidentHasChanged = true;
294 //
295 // correct the strangeness by replacing a pi- by a kaon-
296 //
297 vec.Initialize( 1 );
299 p->SetDefinition( aKaonMinus );
300 (G4UniformRand() < 0.5) ? p->SetSide( -1 ) : p->SetSide( 1 );
301 vec.SetElement( vecLen++, p );
302 --nneg;
303 }
304 }
305 else if( npos == nneg )
306 {
307 if( G4UniformRand() >= 0.5 )
308 {
309 currentParticle.SetDefinitionAndUpdateE( aXiZero );
310 incidentHasChanged = true;
311 targetParticle.SetDefinitionAndUpdateE( aNeutron );
312 targetHasChanged = true;
313 }
314 }
315 else
316 {
317 targetParticle.SetDefinitionAndUpdateE( aNeutron );
318 targetHasChanged = true;
319 }
320 }
321 else // target must be a neutron
322 {
323 counter = -1;
324 for( npos=0; npos<numSec/3 && ran>=excs; ++npos )
325 {
326 for( nneg=npos; nneg<=(npos+2) && ran>=excs; ++nneg )
327 {
328 for( nzero=0; nzero<numSec/3 && ran>=excs; ++nzero )
329 {
330 if( ++counter < numMul )
331 {
332 nt = npos+nneg+nzero;
333 if( nt>0 && nt<=numSec )
334 {
335 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
336 dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
337 if( std::fabs(dum) < 1.0 )
338 {
339 if( test >= 1.0e-10 )excs += dum*test;
340 }
341 else
342 excs += dum*test;
343 }
344 }
345 }
346 }
347 }
348 if( ran >= excs ) // 3 previous loops continued to the end
349 {
350 quasiElastic = true;
351 return;
352 }
353 npos--; nneg--; nzero--;
354 if( npos+1 < nneg )
355 {
356 if( npos+2 == nneg )
357 {
358 currentParticle.SetDefinitionAndUpdateE( aXiZero );
359 incidentHasChanged = true;
360 targetParticle.SetDefinitionAndUpdateE( aProton );
361 targetHasChanged = true;
362 }
363 else // charge mismatch
364 {
365 currentParticle.SetDefinitionAndUpdateE( aSigmaPlus );
366 incidentHasChanged = true;
367 targetParticle.SetDefinitionAndUpdateE( aProton );
368 targetHasChanged = true;
369 //
370 // correct the strangeness by replacing a pi- by a kaon-
371 //
372 vec.Initialize( 1 );
374 p->SetDefinition( aKaonMinus );
375 (G4UniformRand() < 0.5) ? p->SetSide( -1 ) : p->SetSide( 1 );
376 vec.SetElement( vecLen++, p );
377 --nneg;
378 }
379 }
380 else if( npos+1 == nneg )
381 {
382 if( G4UniformRand() < 0.5 )
383 {
384 currentParticle.SetDefinitionAndUpdateE( aXiZero );
385 incidentHasChanged = true;
386 }
387 else
388 {
389 targetParticle.SetDefinitionAndUpdateE( aProton );
390 targetHasChanged = true;
391 }
392 }
393 }
394 SetUpPions(npos, nneg, nzero, vec, vecLen);
395 return;
396}
397
398 /* end of file */
399
@ 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 G4XiZero * XiZero()
Definition: G4XiZero.cc:106
const G4double pi