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