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