Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LEAntiLambdaInelastic.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: AntiLambda Inelastic Process
29// J.L. Chuma, TRIUMF, 19-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 G4LEAntiLambdaInelastic::ModelDescription(std::ostream& outFile) const
38{
39 outFile << "G4LEAntiLambdaInelastic is one of the Low Energy Parameterized\n"
40 << "(LEP) models used to implement inelastic anti-lambda\n"
41 << "scattering from nuclei. It is a re-engineered version of the\n"
42 << "GHEISHA code of H. Fesefeldt. It divides the initial\n"
43 << "collision products into backward- and forward-going clusters\n"
44 << "which are then decayed into final state hadrons. The model\n"
45 << "does not conserve energy on an event-by-event basis. It may\n"
46 << "be applied to anti-lambdas with initial energies between 0 and\n"
47 << "25 GeV.\n";
48}
49
50
53 G4Nucleus& targetNucleus)
54{
55 const G4HadProjectile *originalIncident = &aTrack;
56
57 // create the target particle
58 G4DynamicParticle* originalTarget = targetNucleus.ReturnTargetParticle();
59
60 if (verboseLevel > 1) {
61 const G4Material *targetMaterial = aTrack.GetMaterial();
62 G4cout << "G4LEAntiLambdaInelastic::ApplyYourself called" << G4endl;
63 G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy()/MeV << "MeV, ";
64 G4cout << "target material = " << targetMaterial->GetName() << ", ";
65 G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName()
66 << G4endl;
67 }
68
69 // Fermi motion and evaporation
70 // As of Geant3, the Fermi energy calculation had not been Done
71 G4double ek = originalIncident->GetKineticEnergy()/MeV;
72 G4double amas = originalIncident->GetDefinition()->GetPDGMass()/MeV;
73 G4ReactionProduct modifiedOriginal;
74 modifiedOriginal = *originalIncident;
75
76 G4double tkin = targetNucleus.Cinema( ek );
77 ek += tkin;
78 modifiedOriginal.SetKineticEnergy( ek*MeV );
79 G4double et = ek + amas;
80 G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
81 G4double pp = modifiedOriginal.GetMomentum().mag()/MeV;
82 if (pp > 0.0) {
83 G4ThreeVector momentum = modifiedOriginal.GetMomentum();
84 modifiedOriginal.SetMomentum( momentum * (p/pp) );
85 }
86
87 // calculate black track energies
88 tkin = targetNucleus.EvaporationEffects( ek );
89 ek -= tkin;
90 modifiedOriginal.SetKineticEnergy( ek*MeV );
91 et = ek + amas;
92 p = std::sqrt( std::abs((et-amas)*(et+amas)) );
93 pp = modifiedOriginal.GetMomentum().mag()/MeV;
94 if (pp > 0.0) {
95 G4ThreeVector momentum = modifiedOriginal.GetMomentum();
96 modifiedOriginal.SetMomentum( momentum * (p/pp) );
97 }
98
99 G4ReactionProduct currentParticle = modifiedOriginal;
100 G4ReactionProduct targetParticle;
101 targetParticle = *originalTarget;
102 currentParticle.SetSide(1); // incident always goes in forward hemisphere
103 targetParticle.SetSide(-1); // target always goes in backward hemisphere
104 G4bool incidentHasChanged = false;
105 G4bool targetHasChanged = false;
106 G4bool quasiElastic = false;
107 G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec; // vec will contain the secondary particles
108 G4int vecLen = 0;
109 vec.Initialize(0);
110
111 const G4double cutOff = 0.1;
112 const G4double anni = std::min( 1.3*currentParticle.GetTotalMomentum()/GeV, 0.4 );
113 if ( (originalIncident->GetKineticEnergy()/MeV > cutOff) || (G4UniformRand() > anni) )
114 Cascade(vec, vecLen, originalIncident, currentParticle, targetParticle,
115 incidentHasChanged, targetHasChanged, quasiElastic);
116
117 CalculateMomenta(vec, vecLen, originalIncident, originalTarget,
118 modifiedOriginal, targetNucleus, currentParticle,
119 targetParticle, incidentHasChanged, targetHasChanged,
120 quasiElastic);
121
122 SetUpChange(vec, vecLen, currentParticle, targetParticle, incidentHasChanged);
123
124 if (isotopeProduction) DoIsotopeCounting(originalIncident, targetNucleus);
125
126 delete originalTarget;
127 return &theParticleChange;
128}
129
130
131void G4LEAntiLambdaInelastic::Cascade(
133 G4int &vecLen,
134 const G4HadProjectile *originalIncident,
135 G4ReactionProduct &currentParticle,
136 G4ReactionProduct &targetParticle,
137 G4bool &incidentHasChanged,
138 G4bool &targetHasChanged,
139 G4bool &quasiElastic )
140{
141 // derived from original FORTRAN code CASAL0 by H. Fesefeldt (13-Sep-1987)
142 //
143 // AntiLambda undergoes interaction with nucleon within a nucleus. Check if it is
144 // energetically possible to produce pions/kaons. In not, assume nuclear excitation
145 // occurs and input particle is degraded in energy. No other particles are produced.
146 // If reaction is possible, find the correct number of pions/protons/neutrons
147 // produced using an interpolation to multiplicity data. Replace some pions or
148 // protons/neutrons by kaons or strange baryons according to the average
149 // multiplicity per Inelastic reaction.
150
151 const G4double mOriginal = originalIncident->GetDefinition()->GetPDGMass()/MeV;
152 const G4double etOriginal = originalIncident->GetTotalEnergy()/MeV;
153 const G4double targetMass = targetParticle.GetMass()/MeV;
154 const G4double pOriginal = originalIncident->GetTotalMomentum()/GeV;
155 G4double centerofmassEnergy = std::sqrt(mOriginal*mOriginal +
156 targetMass*targetMass +
157 2.0*targetMass*etOriginal);
158 G4double availableEnergy = centerofmassEnergy - (targetMass+mOriginal);
159
160 static G4bool first = true;
161 const G4int numMul = 1200;
162 const G4int numMulA = 400;
163 const G4int numSec = 60;
164 static G4double protmul[numMul], protnorm[numSec]; // proton constants
165 static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
166 static G4double protmulA[numMulA], protnormA[numSec]; // proton constants
167 static G4double neutmulA[numMulA], neutnormA[numSec]; // neutron constants
168
169 // npos = number of pi+, nneg = number of pi-, nzero = number of pi0
170 G4int nt=0;
171 G4int npos = 0, nneg = 0, nzero = 0;
172 G4double test;
173 const G4double c = 1.25;
174 const G4double b[] = { 0.7, 0.7 };
175 if (first) { // Computation of normalization constants will only be done once
176 first = false;
177 G4int i;
178 for (i = 0; i < numMul; ++i) protmul[i] = 0.0;
179 for (i = 0; i < numSec; ++i) protnorm[i] = 0.0;
180 G4int counter = -1;
181 for (npos = 0; npos < (numSec/3); ++npos) {
182 for (nneg = std::max(0, npos-2); nneg <= (npos+1); ++nneg) {
183 for (nzero = 0; nzero < numSec/3; ++nzero) {
184 if (++counter < numMul) {
185 nt = npos + nneg + nzero;
186 if (nt > 0 && nt <= numSec) {
187 protmul[counter] = Pmltpc(npos,nneg,nzero,nt,b[0],c);
188 protnorm[nt-1] += protmul[counter];
189 }
190 }
191 }
192 }
193 }
194
195 for (i = 0; i < numMul; ++i) neutmul[i] = 0.0;
196 for (i = 0; i < numSec; ++i) neutnorm[i] = 0.0;
197 counter = -1;
198 for (npos = 0; npos < numSec/3; ++npos) {
199 for (nneg = std::max(0,npos-1); nneg <= (npos+2); ++nneg) {
200 for (nzero = 0; nzero < numSec/3; ++nzero) {
201 if (++counter < numMul) {
202 nt = npos + nneg + nzero;
203 if (nt > 0 && nt <= numSec) {
204 neutmul[counter] = Pmltpc(npos,nneg,nzero,nt,b[1],c);
205 neutnorm[nt-1] += neutmul[counter];
206 }
207 }
208 }
209 }
210 }
211
212 for (i = 0; i < numSec; ++i) {
213 if (protnorm[i] > 0.0) protnorm[i] = 1.0/protnorm[i];
214 if (neutnorm[i] > 0.0) neutnorm[i] = 1.0/neutnorm[i];
215 }
216
217 // do the same for annihilation channels
218 for (i = 0; i < numMulA; ++i) protmulA[i] = 0.0;
219 for (i = 0; i < numSec; ++i) protnormA[i] = 0.0;
220 counter = -1;
221 for (npos = 1; npos < (numSec/3); ++npos) {
222 nneg = npos - 1;
223 for (nzero = 0; nzero < numSec/3; ++nzero) {
224 if (++counter < numMulA) {
225 nt = npos + nneg + nzero;
226 if (nt > 1 && nt <= numSec) {
227 protmulA[counter] = Pmltpc(npos,nneg,nzero,nt,b[0],c);
228 protnormA[nt-1] += protmulA[counter];
229 }
230 }
231 }
232 }
233
234 for (i = 0; i < numMulA; ++i) neutmulA[i] = 0.0;
235 for (i = 0; i < numSec; ++i) neutnormA[i] = 0.0;
236 counter = -1;
237 for (npos = 0; npos < numSec/3; ++npos) {
238 nneg = npos;
239 for (nzero = 0; nzero < numSec/3; ++nzero) {
240 if (++counter < numMulA) {
241 nt = npos + nneg + nzero;
242 if (nt > 1 && nt <= numSec) {
243 neutmulA[counter] = Pmltpc(npos,nneg,nzero,nt,b[1],c);
244 neutnormA[nt-1] += neutmulA[counter];
245 }
246 }
247 }
248 }
249 for (i = 0; i < numSec; ++i) {
250 if (protnormA[i] > 0.0) protnormA[i] = 1.0/protnormA[i];
251 if (neutnormA[i] > 0.0) neutnormA[i] = 1.0/neutnormA[i];
252 }
253 } // end of initialization
254
255 const G4double expxu = 82.; // upper bound for arg. of exp
256 const G4double expxl = -expxu; // lower bound for arg. of exp
257
269 const G4double anhl[] = {1.00,1.00,1.00,1.00,1.00,1.00,1.00,1.00,0.97,0.88,
270 0.85,0.81,0.75,0.64,0.64,0.55,0.55,0.45,0.47,0.40,
271 0.39,0.36,0.33,0.10,0.01};
272 G4int iplab = G4int( pOriginal*10.0 );
273 if (iplab > 9) iplab = G4int( (pOriginal- 1.0)*5.0 ) + 10;
274 if (iplab > 14) iplab = G4int( pOriginal- 2.0 ) + 15;
275 if (iplab > 22) iplab = G4int( (pOriginal-10.0)/10.0 ) + 23;
276 if (iplab > 24) iplab = 24;
277
278 if (G4UniformRand() > anhl[iplab]) {
279 if (availableEnergy <= aPiPlus->GetPDGMass()/MeV) {
280 // not energetically possible to produce pion(s)
281 quasiElastic = true;
282 return;
283 }
284 G4double n, anpn;
285 GetNormalizationConstant (availableEnergy, n, anpn);
286 G4double ran = G4UniformRand();
287 G4double dum, excs = 0.0;
288 if (targetParticle.GetDefinition() == aProton) {
289 G4int counter = -1;
290 for (npos = 0; npos < numSec/3 && ran>=excs; ++npos) {
291 for (nneg = std::max(0,npos-2); nneg <= (npos+1) && ran >= excs; ++nneg) {
292 for (nzero = 0; nzero < numSec/3 && ran >= excs; ++nzero) {
293 if (++counter < numMul) {
294 nt = npos + nneg + nzero;
295 if (nt > 0 && nt <= numSec) {
296 test = std::exp(std::min(expxu,
297 std::max(expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
298 dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
299 if (std::fabs(dum) < 1.0) {
300 if (test >= 1.0e-10 )excs += dum*test;
301 } else {
302 excs += dum*test;
303 }
304 }
305 }
306 }
307 }
308 }
309
310 if (ran >= excs) {
311 // 3 previous loops continued to the end
312 quasiElastic = true;
313 return;
314 }
315 npos--; nneg--; nzero--;
316 G4int ncht = std::min(4, std::max(1, npos-nneg+2 ) );
317 switch (ncht)
318 {
319 case 1:
320 currentParticle.SetDefinitionAndUpdateE(anAntiSigmaMinus);
321 incidentHasChanged = true;
322 break;
323 case 2:
324 if (G4UniformRand() < 0.5) {
325 if (G4UniformRand() < 0.5) {
326 currentParticle.SetDefinitionAndUpdateE(anAntiSigmaZero);
327 incidentHasChanged = true;
328 } else {
329 currentParticle.SetDefinitionAndUpdateE(anAntiSigmaMinus);
330 incidentHasChanged = true;
331 targetParticle.SetDefinitionAndUpdateE(aNeutron);
332 targetHasChanged = true;
333 }
334 } else {
335 if (G4UniformRand() >= 0.5) {
336 currentParticle.SetDefinitionAndUpdateE(anAntiSigmaMinus);
337 incidentHasChanged = true;
338 targetParticle.SetDefinitionAndUpdateE(aNeutron);
339 targetHasChanged = true;
340 }
341 }
342 break;
343 case 3:
344 if (G4UniformRand() < 0.5) {
345 if (G4UniformRand() < 0.5) {
346 currentParticle.SetDefinitionAndUpdateE(anAntiSigmaZero);
347 incidentHasChanged = true;
348 targetParticle.SetDefinitionAndUpdateE(aNeutron);
349 targetHasChanged = true;
350 } else {
351 currentParticle.SetDefinitionAndUpdateE(anAntiSigmaPlus);
352 incidentHasChanged = true;
353 }
354 } else {
355 if (G4UniformRand() < 0.5) {
356 targetParticle.SetDefinitionAndUpdateE(aNeutron);
357 targetHasChanged = true;
358 } else {
359 currentParticle.SetDefinitionAndUpdateE(anAntiSigmaPlus);
360 incidentHasChanged = true;
361 }
362 }
363 break;
364 case 4:
365 currentParticle.SetDefinitionAndUpdateE(anAntiSigmaPlus);
366 incidentHasChanged = true;
367 targetParticle.SetDefinitionAndUpdateE(aNeutron);
368 targetHasChanged = true;
369 break;
370 } // end switch
371
372 } else {
373 // target must be a neutron
374 G4int counter = -1;
375 for (npos = 0; npos < numSec/3 && ran>=excs; ++npos) {
376 for (nneg = std::max(0,npos-1); nneg <= (npos+2) && ran>=excs; ++nneg) {
377 for (nzero = 0; nzero < numSec/3 && ran >= excs; ++nzero) {
378 if (++counter < numMul) {
379 nt = npos + nneg + nzero;
380 if (nt > 0 && nt <= numSec) {
381 test = std::exp(std::min(expxu,
382 std::max(expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
383 dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
384 if (std::fabs(dum) < 1.0) {
385 if (test >= 1.0e-10) excs += dum*test;
386 } else {
387 excs += dum*test;
388 }
389 }
390 }
391 }
392 }
393 }
394 if (ran >= excs) {
395 // 3 previous loops continued to the end
396 quasiElastic = true;
397 return;
398 }
399
400 npos--; nneg--; nzero--;
401 G4int ncht = std::min(4, std::max(1, npos-nneg+3) );
402 switch (ncht)
403 {
404 case 1:
405 currentParticle.SetDefinitionAndUpdateE(anAntiSigmaMinus);
406 incidentHasChanged = true;
407 targetParticle.SetDefinitionAndUpdateE(aProton);
408 targetHasChanged = true;
409 break;
410 case 2:
411 if (G4UniformRand() < 0.5) {
412 if (G4UniformRand() < 0.5) {
413 currentParticle.SetDefinitionAndUpdateE(anAntiSigmaZero);
414 incidentHasChanged = true;
415 targetParticle.SetDefinitionAndUpdateE(aProton);
416 targetHasChanged = true;
417 } else {
418 currentParticle.SetDefinitionAndUpdateE(anAntiSigmaMinus);
419 incidentHasChanged = true;
420 }
421 } else {
422 if (G4UniformRand() < 0.5) {
423 targetParticle.SetDefinitionAndUpdateE(aProton);
424 targetHasChanged = true;
425 } else {
426 currentParticle.SetDefinitionAndUpdateE(anAntiSigmaMinus);
427 incidentHasChanged = true;
428 }
429 }
430 break;
431 case 3:
432 if (G4UniformRand() < 0.5) {
433 if (G4UniformRand() < 0.5) {
434 currentParticle.SetDefinitionAndUpdateE(anAntiSigmaZero);
435 incidentHasChanged = true;
436 } else {
437 currentParticle.SetDefinitionAndUpdateE(anAntiSigmaPlus);
438 incidentHasChanged = true;
439 targetParticle.SetDefinitionAndUpdateE( aProton );
440 targetHasChanged = true;
441 }
442 } else {
443 if (G4UniformRand() >= 0.5) {
444 currentParticle.SetDefinitionAndUpdateE(anAntiSigmaPlus);
445 incidentHasChanged = true;
446 targetParticle.SetDefinitionAndUpdateE(aProton);
447 targetHasChanged = true;
448 }
449 }
450 break;
451 default:
452 currentParticle.SetDefinitionAndUpdateE(anAntiSigmaPlus);
453 incidentHasChanged = true;
454 break;
455 } // end of switch
456 }
457 } else { // random number <= anhl[iplab]
458 if (centerofmassEnergy <= aPiPlus->GetPDGMass()/MeV+aKaonPlus->GetPDGMass()/MeV) {
459 quasiElastic = true;
460 return;
461 }
462
463 // annihilation channels
464 G4double n, anpn;
465 GetNormalizationConstant(-centerofmassEnergy, n, anpn);
466 G4double ran = G4UniformRand();
467 G4double dum, excs = 0.0;
468 if (targetParticle.GetDefinition() == aProton) {
469 G4int counter = -1;
470 for (npos = 1; npos < numSec/3 && ran>=excs; ++npos) {
471 nneg = npos - 1;
472 for (nzero = 0; nzero < numSec/3 && ran >= excs; ++nzero) {
473 if (++counter < numMulA) {
474 nt = npos + nneg + nzero;
475 if (nt > 1 && nt <= numSec) {
476 test = std::exp(std::min(expxu,
477 std::max(expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
478 dum = (pi/anpn)*nt*protmulA[counter]*protnormA[nt-1]/(2.0*n*n);
479 if (std::fabs(dum) < 1.0) {
480 if (test >= 1.0e-10) excs += dum*test;
481 } else {
482 excs += dum*test;
483 }
484 }
485 }
486 }
487 }
488 } else {
489 // target must be a neutron
490 G4int counter = -1;
491 for (npos = 0; npos < numSec/3 && ran >= excs; ++npos) {
492 nneg = npos;
493 for (nzero = 0; nzero < numSec/3 && ran >= excs; ++nzero) {
494 if (++counter < numMulA) {
495 nt = npos + nneg + nzero;
496 if (nt > 1 && nt <= numSec) {
497 test = std::exp(std::min(expxu,
498 std::max(expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
499 dum = (pi/anpn)*nt*neutmulA[counter]*neutnormA[nt-1]/(2.0*n*n);
500 if (std::fabs(dum) < 1.0) {
501 if (test >= 1.0e-10) excs += dum*test;
502 } else {
503 excs += dum*test;
504 }
505 }
506 }
507 }
508 }
509 }
510 if (ran >= excs) {
511 // 3 previous loops continued to the end
512 quasiElastic = true;
513 return;
514 }
515 npos--; nzero--;
516 currentParticle.SetMass( 0.0 );
517 targetParticle.SetMass( 0.0 );
518 }
519
520 SetUpPions(npos, nneg, nzero, vec, vecLen);
521 if (currentParticle.GetMass() == 0.0) {
522 if (nzero == 0) {
523 if (nneg > 0) {
524 for (G4int i = 0; i < vecLen; ++i) {
525 if (vec[i]->GetDefinition() == aPiMinus) {
526 vec[i]->SetDefinitionAndUpdateE( aKaonMinus );
527 break;
528 }
529 }
530 }
531 } else {
532 // nzero > 0
533 if (nneg == 0) {
534 for (G4int i = 0; i < vecLen; ++i) {
535 if (vec[i]->GetDefinition() == aPiZero) {
536 vec[i]->SetDefinitionAndUpdateE(aKaonZL);
537 break;
538 }
539 }
540 } else {
541 // nneg > 0
542 if (G4UniformRand() < 0.5) {
543 if (nneg > 0) {
544 for (G4int i = 0; i < vecLen; ++i) {
545 if (vec[i]->GetDefinition() == aPiMinus) {
546 vec[i]->SetDefinitionAndUpdateE(aKaonMinus);
547 break;
548 }
549 }
550 }
551 } else {
552 // random number >= 0.5
553 for (G4int i = 0; i < vecLen; ++i) {
554 if (vec[i]->GetDefinition() == aPiZero) {
555 vec[i]->SetDefinitionAndUpdateE( aKaonZL );
556 break;
557 }
558 }
559 }
560 }
561 }
562 }
563 return;
564}
565
566 /* end of file */
567
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
static G4AntiSigmaMinus * AntiSigmaMinus()
static G4AntiSigmaPlus * AntiSigmaPlus()
static G4AntiSigmaZero * AntiSigmaZero()
G4ParticleDefinition * GetDefinition() const
void Initialize(G4int items)
Definition: G4FastVector.hh:63
const G4Material * GetMaterial() const
G4double GetTotalMomentum() const
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() 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
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:113
static G4KaonZeroLong * KaonZeroLong()
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 G4PionZero * PionZero()
Definition: G4PionZero.cc:104
static G4Proton * Proton()
Definition: G4Proton.cc:93
void SetMomentum(const G4double x, const G4double y, const G4double z)
G4double GetTotalMomentum() const
G4ThreeVector GetMomentum() const
void SetSide(const G4int sid)
void SetDefinitionAndUpdateE(G4ParticleDefinition *aParticleDefinition)
void SetKineticEnergy(const G4double en)
G4ParticleDefinition * GetDefinition() const
G4double GetMass() const
void SetMass(const G4double mas)
const G4double pi