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
G4LEAntiSigmaPlusInelastic.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: AntiSigmaPlus 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 G4LEAntiSigmaPlusInelastic::ModelDescription(std::ostream& outFile) const
38{
39 outFile << "G4LEAntiSigmaPlusInelastic is one of the Low Energy\n"
40 << "Parameterized (LEP) models used to implement inelastic\n"
41 << "antiSigma+ scattering from nuclei. It is a re-engineered\n"
42 << "version of the GHEISHA code of H. Fesefeldt. It divides the\n"
43 << "initial collision products into backward- and forward-going\n"
44 << "clusters which are then decayed into final state hadrons. The\n"
45 << "model does not conserve energy on an event-by-event basis. It\n"
46 << "may be applied to antiSigma+ with initial energies between 0\n"
47 << "and 25 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 << "G4LEAntiSigmaPlusInelastic::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 const G4double anni = std::min( 1.3*currentParticle.GetTotalMomentum()/GeV, 0.4 );
118 if ((currentParticle.GetKineticEnergy()/MeV > cutOff) ||
119 (G4UniformRand() > anni) )
120 Cascade(vec, vecLen, originalIncident, currentParticle, targetParticle,
121 incidentHasChanged, targetHasChanged, quasiElastic);
122
123 CalculateMomenta(vec, vecLen, originalIncident, originalTarget,
124 modifiedOriginal, targetNucleus, currentParticle,
125 targetParticle, incidentHasChanged, targetHasChanged,
126 quasiElastic);
127
128 SetUpChange(vec, vecLen, currentParticle, targetParticle, incidentHasChanged);
129
130 if (isotopeProduction) DoIsotopeCounting(originalIncident, targetNucleus);
131
132 delete originalTarget;
133 return &theParticleChange;
134}
135
136void G4LEAntiSigmaPlusInelastic::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 CASASP by H. Fesefeldt (13-Sep-1987)
147 //
148 // AntiSigmaPlus 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()/MeV;
157 const G4double etOriginal = originalIncident->GetTotalEnergy()/MeV;
158 const G4double pOriginal = originalIncident->GetTotalMomentum()/MeV;
159 const G4double targetMass = targetParticle.GetMass()/MeV;
160 G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
161 targetMass*targetMass +
162 2.0*targetMass*etOriginal );
163 G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
164
165 static G4bool first = true;
166 const G4int numMul = 1200;
167 const G4int numMulA = 400;
168 const G4int numSec = 60;
169 static G4double protmul[numMul], protnorm[numSec]; // proton constants
170 static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
171 static G4double protmulA[numMulA], protnormA[numSec]; // proton constants
172 static G4double neutmulA[numMulA], neutnormA[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.7, 0.7 };
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 && nt<=numSec )
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 //
230 // do the same for annihilation channels
231 //
232 for( i=0; i<numMulA; ++i )protmulA[i] = 0.0;
233 for( i=0; i<numSec; ++i )protnormA[i] = 0.0;
234 counter = -1;
235 for( npos=1; npos<(numSec/3); ++npos )
236 {
237 nneg = npos;
238 for( nzero=0; nzero<numSec/3; ++nzero )
239 {
240 if( ++counter < numMulA )
241 {
242 nt = npos+nneg+nzero;
243 if( nt>1 && nt<=numSec )
244 {
245 protmulA[counter] = Pmltpc(npos,nneg,nzero,nt,b[0],c);
246 protnormA[nt-1] += protmulA[counter];
247 }
248 }
249 }
250 }
251 for( i=0; i<numMulA; ++i )neutmulA[i] = 0.0;
252 for( i=0; i<numSec; ++i )neutnormA[i] = 0.0;
253 counter = -1;
254 for( npos=0; npos<numSec/3; ++npos )
255 {
256 nneg = npos+1;
257 for( nzero=0; nzero<numSec/3; ++nzero )
258 {
259 if( ++counter < numMulA )
260 {
261 nt = npos+nneg+nzero;
262 if( nt>1 && nt<=numSec )
263 {
264 neutmulA[counter] = Pmltpc(npos,nneg,nzero,nt,b[1],c);
265 neutnormA[nt-1] += neutmulA[counter];
266 }
267 }
268 }
269 }
270 for( i=0; i<numSec; ++i )
271 {
272 if( protnormA[i] > 0.0 )protnormA[i] = 1.0/protnormA[i];
273 if( neutnormA[i] > 0.0 )neutnormA[i] = 1.0/neutnormA[i];
274 }
275 } // end of initialization
276
277 const G4double expxu = 82.; // upper bound for arg. of exp
278 const G4double expxl = -expxu; // lower bound for arg. of exp
287 const G4double anhl[] = {1.00,1.00,1.00,1.00,1.00,1.00,1.00,1.00,0.97,0.88,
288 0.85,0.81,0.75,0.64,0.64,0.55,0.55,0.45,0.47,0.40,
289 0.39,0.36,0.33,0.10,0.01};
290 G4int iplab = G4int( pOriginal/GeV*10.0 );
291 if( iplab > 9 )iplab = G4int( (pOriginal/GeV- 1.0)*5.0 ) + 10;
292 if( iplab > 14 )iplab = G4int( pOriginal/GeV- 2.0 ) + 15;
293 if( iplab > 22 )iplab = G4int( (pOriginal/GeV-10.0)/10.0 ) + 23;
294 if( iplab > 24 )iplab = 24;
295 if( G4UniformRand() > anhl[iplab] )
296 {
297 if( availableEnergy <= aPiPlus->GetPDGMass()/MeV )
298 {
299 quasiElastic = true;
300 return;
301 }
302 G4double n, anpn;
303 GetNormalizationConstant( availableEnergy, n, anpn );
304 G4double ran = G4UniformRand();
305 G4double dum, excs = 0.0;
306 if( targetParticle.GetDefinition() == aProton )
307 {
308 counter = -1;
309 for( npos=0; npos<numSec/3 && ran>=excs; ++npos )
310 {
311 for( nneg=std::max(0,npos-1); nneg<=(npos+1) && ran>=excs; ++nneg )
312 {
313 for( nzero=0; nzero<numSec/3 && ran>=excs; ++nzero )
314 {
315 if( ++counter < numMul )
316 {
317 nt = npos+nneg+nzero;
318 if( (nt>0) && (nt<=numSec) )
319 {
320 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
321 dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
322 if( std::fabs(dum) < 1.0 )
323 {
324 if( test >= 1.0e-10 )excs += dum*test;
325 }
326 else
327 excs += dum*test;
328 }
329 }
330 }
331 }
332 }
333 if( ran >= excs ) // 3 previous loops continued to the end
334 {
335 quasiElastic = true;
336 return;
337 }
338 npos--; nneg--; nzero--;
339 G4int ncht = std::min( 3, std::max( 1, npos-nneg+2 ) );
340 switch( ncht )
341 {
342 case 1:
343 if( G4UniformRand() < 0.5 )
344 currentParticle.SetDefinitionAndUpdateE( anAntiLambda );
345 else
346 currentParticle.SetDefinitionAndUpdateE( anAntiSigmaZero );
347 incidentHasChanged = true;
348 break;
349 case 2:
350 if( G4UniformRand() >= 0.5 )
351 {
352 if( G4UniformRand() < 0.5 )
353 currentParticle.SetDefinitionAndUpdateE( anAntiLambda );
354 else
355 currentParticle.SetDefinitionAndUpdateE( anAntiSigmaZero );
356 incidentHasChanged = true;
357 }
358 targetParticle.SetDefinitionAndUpdateE( aNeutron );
359 targetHasChanged = true;
360 break;
361 case 3:
362 targetParticle.SetDefinitionAndUpdateE( aNeutron );
363 targetHasChanged = true;
364 break;
365 }
366 }
367 else // target must be a neutron
368 {
369 counter = -1;
370 for( npos=0; npos<numSec/3 && ran>=excs; ++npos )
371 {
372 for( nneg=npos; nneg<=(npos+2) && ran>=excs; ++nneg )
373 {
374 for( nzero=0; nzero<numSec/3 && ran>=excs; ++nzero )
375 {
376 if( ++counter < numMul )
377 {
378 nt = npos+nneg+nzero;
379 if( (nt>0) && (nt<=numSec) )
380 {
381 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
382 dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
383 if( std::fabs(dum) < 1.0 )
384 {
385 if( test >= 1.0e-10 )excs += dum*test;
386 }
387 else
388 excs += dum*test;
389 }
390 }
391 }
392 }
393 }
394 if( ran >= excs ) // 3 previous loops continued to the end
395 {
396 quasiElastic = true;
397 return;
398 }
399 npos--; nneg--; nzero--;
400 G4int ncht = std::min( 3, std::max( 1, npos-nneg+3 ) );
401 switch( ncht )
402 {
403 case 1:
404 if( G4UniformRand() < 0.5 )
405 currentParticle.SetDefinitionAndUpdateE( anAntiLambda );
406 else
407 currentParticle.SetDefinitionAndUpdateE( anAntiSigmaZero );
408 incidentHasChanged = true;
409 targetParticle.SetDefinitionAndUpdateE( aProton );
410 targetHasChanged = true;
411 break;
412 case 2:
413 if( G4UniformRand() < 0.5 )
414 {
415 if( G4UniformRand() < 0.5 )
416 {
417 currentParticle.SetDefinitionAndUpdateE( anAntiLambda );
418 incidentHasChanged = true;
419 }
420 else
421 {
422 targetParticle.SetDefinitionAndUpdateE( aProton );
423 targetHasChanged = true;
424 }
425 }
426 else
427 {
428 if( G4UniformRand() < 0.5 )
429 {
430 currentParticle.SetDefinitionAndUpdateE( anAntiSigmaZero );
431 incidentHasChanged = true;
432 }
433 else
434 {
435 targetParticle.SetDefinitionAndUpdateE( aProton );
436 targetHasChanged = true;
437 }
438 }
439 break;
440 case 3:
441 break;
442 }
443 }
444 }
445 else // random number <= anhl[iplab]
446 {
447 if( centerofmassEnergy <= aPiPlus->GetPDGMass()/MeV+aKaonPlus->GetPDGMass()/MeV )
448 {
449 quasiElastic = true;
450 return;
451 }
452 G4double n, anpn;
453 GetNormalizationConstant( -centerofmassEnergy, n, anpn );
454 G4double ran = G4UniformRand();
455 G4double dum, excs = 0.0;
456 if( targetParticle.GetDefinition() == aProton )
457 {
458 counter = -1;
459 for( npos=1; npos<numSec/3 && ran>=excs; ++npos )
460 {
461 nneg = npos;
462 for( nzero=0; nzero<numSec/3 && ran>=excs; ++nzero )
463 {
464 if( ++counter < numMulA )
465 {
466 nt = npos+nneg+nzero;
467 if( nt>1 && nt<=numSec )
468 {
469 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
470 dum = (pi/anpn)*nt*protmulA[counter]*protnormA[nt-1]/(2.0*n*n);
471 if( std::fabs(dum) < 1.0 )
472 {
473 if( test >= 1.0e-10 )excs += dum*test;
474 }
475 else
476 excs += dum*test;
477 }
478 }
479 }
480 }
481 if( ran >= excs ) // 3 previous loops continued to the end
482 {
483 quasiElastic = true;
484 return;
485 }
486 npos--; nzero--;
487 }
488 else // target must be a neutron
489 {
490 counter = -1;
491 for( npos=0; npos<numSec/3 && ran>=excs; ++npos )
492 {
493 nneg = npos+1;
494 for( nzero=0; nzero<numSec/3 && ran>=excs; ++nzero )
495 {
496 if( ++counter < numMulA )
497 {
498 nt = npos+nneg+nzero;
499 if( nt>1 && nt<=numSec )
500 {
501 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
502 dum = (pi/anpn)*nt*neutmulA[counter]*neutnormA[nt-1]/(2.0*n*n);
503 if( std::fabs(dum) < 1.0 )
504 {
505 if( test >= 1.0e-10 )excs += dum*test;
506 }
507 else
508 excs += dum*test;
509 }
510 }
511 }
512 }
513 if( ran >= excs ) // 3 previous loops continued to the end
514 {
515 quasiElastic = true;
516 return;
517 }
518 npos--; nzero--;
519 }
520 if( nzero > 0 )
521 {
522 if( nneg > 0 )
523 {
524 if( G4UniformRand() < 0.5 )
525 {
526 vec.Initialize( 1 );
528 p->SetDefinition( aKaonMinus );
529 (G4UniformRand() < 0.5) ? p->SetSide( -1 ) : p->SetSide( 1 );
530 vec.SetElement( vecLen++, p );
531 --nneg;
532 }
533 else
534 {
535 vec.Initialize( 1 );
537 p->SetDefinition( aKaonZL );
538 (G4UniformRand() < 0.5) ? p->SetSide( -1 ) : p->SetSide( 1 );
539 vec.SetElement( vecLen++, p );
540 --nzero;
541 }
542 }
543 else // nneg == 0
544 {
545 vec.Initialize( 1 );
547 p->SetDefinition( aKaonZL );
548 (G4UniformRand() < 0.5) ? p->SetSide( -1 ) : p->SetSide( 1 );
549 vec.SetElement( vecLen++, p );
550 --nzero;
551 }
552 }
553 else // nzero == 0
554 {
555 if( nneg > 0 )
556 {
557 vec.Initialize( 1 );
559 p->SetDefinition( aKaonMinus );
560 (G4UniformRand() < 0.5) ? p->SetSide( -1 ) : p->SetSide( 1 );
561 vec.SetElement( vecLen++, p );
562 --nneg;
563 }
564 }
565 currentParticle.SetMass( 0.0 );
566 targetParticle.SetMass( 0.0 );
567 }
568 SetUpPions( npos, nneg, nzero, vec, vecLen );
569 return;
570}
571
572 /* end of file */
573
@ 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
static G4AntiLambda * AntiLambda()
static G4AntiSigmaZero * AntiSigmaZero()
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
G4double GetTotalMomentum() 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
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 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
void SetDefinition(G4ParticleDefinition *aParticleDefinition)
G4double GetMass() const
void SetMass(const G4double mas)
const G4double pi