Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LELambdaInelastic.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: Lambda Inelastic Process
29// J.L. Chuma, TRIUMF, 18-Feb-1997
30// Modified by J.L.Chuma 30-Apr-97: added originalTarget for CalculateMomenta
31
32#include <iostream>
33
36#include "G4SystemOfUnits.hh"
37#include "Randomize.hh"
38
41{
42 SetMinEnergy(0.0);
43 SetMaxEnergy(25.*GeV);
44 G4cout << "WARNING: model G4LELambdaInelastic is being deprecated and will\n"
45 << "disappear in Geant4 version 10.0" << G4endl;
46}
47
48
49void G4LELambdaInelastic::ModelDescription(std::ostream& outFile) const
50{
51 outFile << "G4LELambdaInelastic is one of the Low Energy Parameterized\n"
52 << "(LEP) models used to implement inelastic Lambda scattering\n"
53 << "from nuclei. It is a re-engineered version of the GHEISHA\n"
54 << "code of H. Fesefeldt. It divides the initial collision\n"
55 << "products into backward- and forward-going clusters which are\n"
56 << "then decayed into final state hadrons. The model does not\n"
57 << "conserve energy on an event-by-event basis. It may be\n"
58 << "applied to lambdas with initial energies between 0 and 25\n"
59 << "GeV.\n";
60}
61
62
65 G4Nucleus& targetNucleus)
66{
67 const G4HadProjectile *originalIncident = &aTrack;
68
69 // create the target particle
70
71 G4DynamicParticle* originalTarget = targetNucleus.ReturnTargetParticle();
72
73 if (verboseLevel > 1) {
74 const G4Material *targetMaterial = aTrack.GetMaterial();
75 G4cout << "G4LELambdaInelastic::ApplyYourself called" << G4endl;
76 G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy()/MeV << "MeV, ";
77 G4cout << "target material = " << targetMaterial->GetName() << ", ";
78 G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName()
79 << G4endl;
80 }
81
82 // Fermi motion and evaporation
83 // As of Geant3, the Fermi energy calculation had not been done
84 G4double ek = originalIncident->GetKineticEnergy()/MeV;
85 G4double amas = originalIncident->GetDefinition()->GetPDGMass()/MeV;
86 G4ReactionProduct modifiedOriginal;
87 modifiedOriginal = *originalIncident;
88
89 G4double tkin = targetNucleus.Cinema( ek );
90 ek += tkin;
91 modifiedOriginal.SetKineticEnergy( ek*MeV );
92 G4double et = ek + amas;
93 G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
94 G4double pp = modifiedOriginal.GetMomentum().mag()/MeV;
95 if (pp > 0.0) {
96 G4ThreeVector momentum = modifiedOriginal.GetMomentum();
97 modifiedOriginal.SetMomentum( momentum * (p/pp) );
98 }
99
100 // calculate black track energies
101 tkin = targetNucleus.EvaporationEffects(ek);
102 ek -= tkin;
103 modifiedOriginal.SetKineticEnergy(ek*MeV);
104 et = ek + amas;
105 p = std::sqrt( std::abs((et-amas)*(et+amas)) );
106 pp = modifiedOriginal.GetMomentum().mag()/MeV;
107 if (pp > 0.0) {
108 G4ThreeVector momentum = modifiedOriginal.GetMomentum();
109 modifiedOriginal.SetMomentum( momentum * (p/pp) );
110 }
111
112 G4ReactionProduct currentParticle = modifiedOriginal;
113 G4ReactionProduct targetParticle;
114 targetParticle = *originalTarget;
115 currentParticle.SetSide(1); // incident always goes in forward hemisphere
116 targetParticle.SetSide(-1); // target always goes in backward hemisphere
117 G4bool incidentHasChanged = false;
118 G4bool targetHasChanged = false;
119 G4bool quasiElastic = false;
120 G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec; // vec will contain the secondary particles
121 G4int vecLen = 0;
122 vec.Initialize(0);
123
124 const G4double cutOff = 0.1;
125 if (currentParticle.GetKineticEnergy()/MeV > cutOff)
126 Cascade(vec, vecLen, originalIncident, currentParticle, targetParticle,
127 incidentHasChanged, targetHasChanged, quasiElastic);
128
129 CalculateMomenta(vec, vecLen, originalIncident, originalTarget,
130 modifiedOriginal, targetNucleus, currentParticle,
131 targetParticle, incidentHasChanged, targetHasChanged,
132 quasiElastic);
133
134 SetUpChange(vec, vecLen, currentParticle, targetParticle, incidentHasChanged);
135
136 if (isotopeProduction) DoIsotopeCounting(originalIncident, targetNucleus);
137
138 delete originalTarget;
139 return &theParticleChange;
140}
141
142
143void G4LELambdaInelastic::Cascade(
145 G4int& vecLen,
146 const G4HadProjectile* originalIncident,
147 G4ReactionProduct& currentParticle,
148 G4ReactionProduct& targetParticle,
149 G4bool& incidentHasChanged,
150 G4bool& targetHasChanged,
151 G4bool& quasiElastic)
152{
153 // derived from original FORTRAN code CASL0 by H. Fesefeldt (13-Sep-1987)
154 //
155 // Lambda undergoes interaction with nucleon within a nucleus. Check if it
156 // is energetically possible to produce pions/kaons. In not, assume
157 // nuclear excitation occurs and input particle is degraded in energy. No
158 // other particles are produced. If reaction is possible, find the correct
159 // number of pions/protons/neutrons produced using an interpolation to
160 // multiplicity data. Replace some pions or protons/neutrons by kaons or
161 // strange baryons according to the average multiplicity per inelastic
162 // reaction.
163
164 const G4double mOriginal = originalIncident->GetDefinition()->GetPDGMass()/MeV;
165 const G4double etOriginal = originalIncident->GetTotalEnergy()/MeV;
166 const G4double targetMass = targetParticle.GetMass()/MeV;
167 G4double centerofmassEnergy = std::sqrt(mOriginal*mOriginal +
168 targetMass*targetMass +
169 2.0*targetMass*etOriginal );
170 G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
171 if (availableEnergy <= G4PionPlus::PionPlus()->GetPDGMass()/MeV) {
172 quasiElastic = true;
173 return;
174 }
175 static G4bool first = true;
176 const G4int numMul = 1200;
177 const G4int numSec = 60;
178 static G4double protmul[numMul], protnorm[numSec]; // proton constants
179 static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
180
181 // npos = number of pi+, nneg = number of pi-, nzero = number of pi0
182 G4int counter, nt=0, npos=0, nneg=0, nzero=0;
183 G4double test;
184 const G4double c = 1.25;
185 const G4double b[] = { 0.70, 0.35 };
186 if( first ) { // compute normalization constants, this will only be Done once
187 first = false;
188 G4int i;
189 for( i=0; i<numMul; ++i )protmul[i] = 0.0;
190 for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
191 counter = -1;
192 for( npos=0; npos<(numSec/3); ++npos ) {
193 for( nneg=std::max(0,npos-2); nneg<=(npos+1); ++nneg ) {
194 for( nzero=0; nzero<numSec/3; ++nzero ) {
195 if( ++counter < numMul ) {
196 nt = npos+nneg+nzero;
197 if( nt>0 && nt<=numSec ) {
198 protmul[counter] = Pmltpc(npos,nneg,nzero,nt,b[0],c);
199 protnorm[nt-1] += protmul[counter];
200 }
201 }
202 }
203 }
204 }
205 for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
206 for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
207 counter = -1;
208 for( npos=0; npos<numSec/3; ++npos ) {
209 for( nneg=std::max(0,npos-1); nneg<=(npos+2); ++nneg ) {
210 for( nzero=0; nzero<numSec/3; ++nzero ) {
211 if( ++counter < numMul ) {
212 nt = npos+nneg+nzero;
213 if( nt>0 && nt<=numSec ) {
214 neutmul[counter] = Pmltpc(npos,nneg,nzero,nt,b[1],c);
215 neutnorm[nt-1] += neutmul[counter];
216 }
217 }
218 }
219 }
220 }
221 for( i=0; i<numSec; ++i ) {
222 if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
223 if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
224 }
225 } // end of initialization
226
227 const G4double expxu = 82.; // upper bound for arg. of exp
228 const G4double expxl = -expxu; // lower bound for arg. of exp
234
235 // energetically possible to produce pion(s) --> inelastic scattering
236 // otherwise quasi-elastic scattering
237
238 G4double n, anpn;
239 GetNormalizationConstant( availableEnergy, n, anpn );
240 G4double ran = G4UniformRand();
241 G4double dum, excs = 0.0;
242 if( targetParticle.GetDefinition() == aProton ) {
243 counter = -1;
244 for( npos=0; npos<numSec/3 && ran>=excs; ++npos ) {
245 for( nneg=std::max(0,npos-2); nneg<=(npos+1) && ran>=excs; ++nneg ) {
246 for( nzero=0; nzero<numSec/3 && ran>=excs; ++nzero ) {
247 if( ++counter < numMul ) {
248 nt = npos+nneg+nzero;
249 if( nt>0 && nt<=numSec ) {
250 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
251 dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
252 if( std::fabs(dum) < 1.0 ) {
253 if( test >= 1.0e-10 )excs += dum*test;
254 } else {
255 excs += dum*test;
256 }
257 }
258 }
259 }
260 }
261 }
262 if( ran >= excs ) // 3 previous loops continued to the end
263 {
264 quasiElastic = true;
265 return;
266 }
267 npos--; nneg--; nzero--;
268 G4int ncht = std::max( 1, npos-nneg );
269 switch( ncht ) {
270 case 1:
271 currentParticle.SetDefinitionAndUpdateE( aSigmaPlus );
272 incidentHasChanged = true;
273 break;
274 case 2:
275 if( G4UniformRand() < 0.5 ) {
276 if( G4UniformRand() < 0.5 ) {
277 currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
278 incidentHasChanged = true;
279 }
280 } else {
281 currentParticle.SetDefinitionAndUpdateE( aSigmaPlus );
282 incidentHasChanged = true;
283 targetParticle.SetDefinitionAndUpdateE( aNeutron );
284 targetHasChanged = true;
285 }
286 break;
287 case 3:
288 if( G4UniformRand() < 0.5 ) {
289 if( G4UniformRand() < 0.5 ) {
290 currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
291 incidentHasChanged = true;
292 }
293 targetParticle.SetDefinitionAndUpdateE( aNeutron );
294 targetHasChanged = true;
295 } else {
296 currentParticle.SetDefinitionAndUpdateE( aSigmaMinus );
297 incidentHasChanged = true;
298 }
299 break;
300 default:
301 currentParticle.SetDefinitionAndUpdateE( aSigmaMinus );
302 incidentHasChanged = true;
303 targetParticle.SetDefinitionAndUpdateE( aNeutron );
304 targetHasChanged = true;
305 break;
306 }
307 }
308 else // target must be a neutron
309 {
310 counter = -1;
311 for( npos=0; npos<numSec/3 && ran>=excs; ++npos ) {
312 for( nneg=std::max(0,npos-1); nneg<=(npos+2) && ran>=excs; ++nneg ) {
313 for( nzero=0; nzero<numSec/3 && ran>=excs; ++nzero ) {
314 if( ++counter < numMul ) {
315 nt = npos+nneg+nzero;
316 if( nt>0 && nt<=numSec ) {
317 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
318 dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
319 if( std::fabs(dum) < 1.0 ) {
320 if( test >= 1.0e-10 )excs += dum*test;
321 } else {
322 excs += dum*test;
323 }
324 }
325 }
326 }
327 }
328 }
329 if( ran >= excs ) // 3 previous loops continued to the end
330 {
331 quasiElastic = true;
332 return;
333 }
334 npos--; nneg--; nzero--;
335 G4int ncht = std::max( 1, npos-nneg+3 );
336 switch( ncht ) {
337 case 1:
338 currentParticle.SetDefinitionAndUpdateE( aSigmaPlus );
339 incidentHasChanged = true;
340 targetParticle.SetDefinitionAndUpdateE( aProton );
341 targetHasChanged = true;
342 break;
343 case 2:
344 if( G4UniformRand() < 0.5 ) {
345 if( G4UniformRand() < 0.5 ) {
346 currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
347 incidentHasChanged = true;
348 }
349 targetParticle.SetDefinitionAndUpdateE( aProton );
350 targetHasChanged = true;
351 } else {
352 currentParticle.SetDefinitionAndUpdateE( aSigmaPlus );
353 incidentHasChanged = true;
354 }
355 break;
356 case 3:
357 if( G4UniformRand() < 0.5 ) {
358 if( G4UniformRand() < 0.5 ) {
359 currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
360 incidentHasChanged = true;
361 }
362 } else {
363 currentParticle.SetDefinitionAndUpdateE( aSigmaMinus );
364 incidentHasChanged = true;
365 targetParticle.SetDefinitionAndUpdateE( aProton );
366 targetHasChanged = true;
367 }
368 break;
369 default:
370 currentParticle.SetDefinitionAndUpdateE( aSigmaMinus );
371 incidentHasChanged = true;
372 break;
373 }
374 }
375
376 SetUpPions( npos, nneg, nzero, vec, vecLen );
377 return;
378}
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
double mag() const
G4ParticleDefinition * GetDefinition() const
void Initialize(G4int items)
Definition: G4FastVector.hh:63
const G4Material * GetMaterial() const
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetTotalEnergy() const
void SetMinEnergy(G4double anEnergy)
void SetMaxEnergy(const G4double anEnergy)
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)
G4LELambdaInelastic(const G4String &name="G4LELambdaInelastic")
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
G4double GetMass() const
static G4SigmaMinus * SigmaMinus()
static G4SigmaPlus * SigmaPlus()
Definition: G4SigmaPlus.cc:108
static G4SigmaZero * SigmaZero()
Definition: G4SigmaZero.cc:99
const G4double pi