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