Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4INCLParticle.hh
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 // INCL++ intra-nuclear cascade model
27 // Alain Boudard, CEA-Saclay, France
28 // Joseph Cugnon, University of Liege, Belgium
29 // Jean-Christophe David, CEA-Saclay, France
30 // Pekka Kaitaniemi, CEA-Saclay, France, and Helsinki Institute of Physics, Finland
31 // Sylvie Leray, CEA-Saclay, France
32 // Davide Mancusi, CEA-Saclay, France
33 //
34 #define INCLXX_IN_GEANT4_MODE 1
35 
36 #include "globals.hh"
37 
38 /*
39  * G4INCLParticle.hh
40  *
41  * \date Jun 5, 2009
42  * \author Pekka Kaitaniemi
43  */
44 
45 #ifndef PARTICLE_HH_
46 #define PARTICLE_HH_
47 
48 #include "G4INCLThreeVector.hh"
49 #include "G4INCLParticleTable.hh"
50 #include "G4INCLParticleType.hh"
51 #include "G4INCLParticleSpecies.hh"
52 #include "G4INCLLogger.hh"
53 #include "G4INCLUnorderedVector.hh"
54 #include "G4INCLAllocationPool.hh"
55 #include <sstream>
56 #include <string>
57 
58 namespace G4INCL {
59 
60  class Particle;
61 
62  class ParticleList : public UnorderedVector<Particle*> {
63  public:
64  void rotatePositionAndMomentum(const G4double angle, const ThreeVector &axis) const;
65  void rotatePosition(const G4double angle, const ThreeVector &axis) const;
66  void rotateMomentum(const G4double angle, const ThreeVector &axis) const;
67  void boost(const ThreeVector &b) const;
68  };
69 
70  typedef ParticleList::const_iterator ParticleIter;
71  typedef ParticleList::iterator ParticleMutableIter;
72 
73  class Particle {
74  public:
75  Particle();
77  Particle(ParticleType t, ThreeVector const &momentum, ThreeVector const &position);
78  virtual ~Particle() {}
79 
84  Particle(const Particle &rhs) :
85  theZ(rhs.theZ),
86  theA(rhs.theA),
88  theType(rhs.theType),
89  theEnergy(rhs.theEnergy),
95  nDecays(rhs.nDecays),
99  theHelicity(rhs.theHelicity),
100  emissionTime(rhs.emissionTime),
101  outOfWell(rhs.outOfWell),
102  theMass(rhs.theMass)
103  {
104  if(rhs.thePropagationEnergy == &(rhs.theFrozenEnergy))
106  else
108  if(rhs.thePropagationMomentum == &(rhs.theFrozenMomentum))
110  else
112  // ID intentionally not copied
113  ID = nextID++;
114  }
115 
116  protected:
118  void swap(Particle &rhs) {
119  std::swap(theZ, rhs.theZ);
120  std::swap(theA, rhs.theA);
121  std::swap(theParticipantType, rhs.theParticipantType);
122  std::swap(theType, rhs.theType);
123  if(rhs.thePropagationEnergy == &(rhs.theFrozenEnergy))
125  else
127  std::swap(theEnergy, rhs.theEnergy);
128  std::swap(theFrozenEnergy, rhs.theFrozenEnergy);
129  if(rhs.thePropagationMomentum == &(rhs.theFrozenMomentum))
131  else
133  std::swap(theMomentum, rhs.theMomentum);
134  std::swap(theFrozenMomentum, rhs.theFrozenMomentum);
135  std::swap(thePosition, rhs.thePosition);
136  std::swap(nCollisions, rhs.nCollisions);
137  std::swap(nDecays, rhs.nDecays);
138  std::swap(thePotentialEnergy, rhs.thePotentialEnergy);
139  // ID intentionally not swapped
140 
141  std::swap(theHelicity, rhs.theHelicity);
142  std::swap(emissionTime, rhs.emissionTime);
143  std::swap(outOfWell, rhs.outOfWell);
144 
145  std::swap(theMass, rhs.theMass);
146  std::swap(rpCorrelated, rhs.rpCorrelated);
148  }
149 
150  public:
151 
156  Particle &operator=(const Particle &rhs) {
157  Particle temporaryParticle(rhs);
158  swap(temporaryParticle);
159  return *this;
160  }
161 
167  return theType;
168  };
169 
172  return ParticleSpecies(theType);
173  };
174 
176  theType = t;
177  switch(theType)
178  {
179  case DeltaPlusPlus:
180  theA = 1;
181  theZ = 2;
182  break;
183  case Proton:
184  case DeltaPlus:
185  theA = 1;
186  theZ = 1;
187  break;
188  case Neutron:
189  case DeltaZero:
190  theA = 1;
191  theZ = 0;
192  break;
193  case DeltaMinus:
194  theA = 1;
195  theZ = -1;
196  break;
197  case PiPlus:
198  theA = 0;
199  theZ = 1;
200  break;
201  case PiZero:
202  case Eta:
203  case Omega:
204  case EtaPrime:
205  case Photon:
206  theA = 0;
207  theZ = 0;
208  break;
209  case PiMinus:
210  theA = 0;
211  theZ = -1;
212  break;
213  case Composite:
214  // INCL_ERROR("Trying to set particle type to Composite! Construct a Cluster object instead" << '\n');
215  theA = 0;
216  theZ = 0;
217  break;
218  case UnknownParticle:
219  theA = 0;
220  theZ = 0;
221  INCL_ERROR("Trying to set particle type to Unknown!" << '\n');
222  break;
223  }
224 
225  if( !isResonance() && t!=Composite )
226  setINCLMass();
227  }
228 
232  G4bool isNucleon() const {
234  return true;
235  else
236  return false;
237  };
238 
240  return theParticipantType;
241  }
242 
245  }
246 
249  }
250 
253  }
254 
257  }
258 
259  virtual void makeParticipant() {
261  }
262 
263  virtual void makeTargetSpectator() {
265  }
266 
267  virtual void makeProjectileSpectator() {
269  }
270 
272  G4bool isPion() const { return (theType == PiPlus || theType == PiZero || theType == PiMinus); }
273 
275  G4bool isEta() const { return (theType == Eta); }
276 
278  G4bool isOmega() const { return (theType == Omega); }
279 
281  G4bool isEtaPrime() const { return (theType == EtaPrime); }
282 
284  inline G4bool isResonance() const { return isDelta(); }
285 
287  inline G4bool isDelta() const {
288  return (theType==DeltaPlusPlus || theType==DeltaPlus ||
290  }
291 
293  G4int getA() const { return theA; }
294 
296  G4int getZ() const { return theZ; }
297 
298  G4double getBeta() const {
299  const G4double P = theMomentum.mag();
300  return P/theEnergy;
301  }
302 
310  return theMomentum / theEnergy;
311  }
312 
319  void boost(const ThreeVector &aBoostVector) {
320  const G4double beta2 = aBoostVector.mag2();
321  const G4double gamma = 1.0 / std::sqrt(1.0 - beta2);
322  const G4double bp = theMomentum.dot(aBoostVector);
323  const G4double alpha = (gamma*gamma)/(1.0 + gamma);
324 
325  theMomentum = theMomentum + aBoostVector * (alpha * bp - gamma * theEnergy);
326  theEnergy = gamma * (theEnergy - bp);
327  }
328 
337  void lorentzContract(const ThreeVector &aBoostVector, const ThreeVector &refPos) {
338  const G4double beta2 = aBoostVector.mag2();
339  const G4double gamma = 1.0 / std::sqrt(1.0 - beta2);
340  const ThreeVector theRelativePosition = thePosition - refPos;
341  const ThreeVector transversePosition = theRelativePosition - aBoostVector * (theRelativePosition.dot(aBoostVector) / aBoostVector.mag2());
342  const ThreeVector longitudinalPosition = theRelativePosition - transversePosition;
343 
344  thePosition = refPos + transversePosition + longitudinalPosition / gamma;
345  }
346 
348  inline G4double getMass() const { return theMass; }
349 
351  inline G4double getINCLMass() const {
352  switch(theType) {
353  case Proton:
354  case Neutron:
355  case PiPlus:
356  case PiMinus:
357  case PiZero:
358  case Eta:
359  case Omega:
360  case EtaPrime:
361  case Photon:
363  break;
364 
365  case DeltaPlusPlus:
366  case DeltaPlus:
367  case DeltaZero:
368  case DeltaMinus:
369  return theMass;
370  break;
371 
372  case Composite:
374  break;
375 
376  default:
377  INCL_ERROR("Particle::getINCLMass: Unknown particle type." << '\n');
378  return 0.0;
379  break;
380  }
381  }
382 
384  inline virtual G4double getTableMass() const {
385  switch(theType) {
386  case Proton:
387  case Neutron:
388  case PiPlus:
389  case PiMinus:
390  case PiZero:
391  case Eta:
392  case Omega:
393  case EtaPrime:
394  case Photon:
396  break;
397 
398  case DeltaPlusPlus:
399  case DeltaPlus:
400  case DeltaZero:
401  case DeltaMinus:
402  return theMass;
403  break;
404 
405  case Composite:
407  break;
408 
409  default:
410  INCL_ERROR("Particle::getTableMass: Unknown particle type." << '\n');
411  return 0.0;
412  break;
413  }
414  }
415 
417  inline G4double getRealMass() const {
418  switch(theType) {
419  case Proton:
420  case Neutron:
421  case PiPlus:
422  case PiMinus:
423  case PiZero:
424  case Eta:
425  case Omega:
426  case EtaPrime:
427  case Photon:
429  break;
430 
431  case DeltaPlusPlus:
432  case DeltaPlus:
433  case DeltaZero:
434  case DeltaMinus:
435  return theMass;
436  break;
437 
438  case Composite:
440  break;
441 
442  default:
443  INCL_ERROR("Particle::getRealMass: Unknown particle type." << '\n');
444  return 0.0;
445  break;
446  }
447  }
448 
451 
454 
457 
469  G4double getEmissionQValueCorrection(const G4int AParent, const G4int ZParent) const {
470  const G4int ADaughter = AParent - theA;
471  const G4int ZDaughter = ZParent - theZ;
472 
473  // Note the minus sign here
474  G4double theQValue;
475  if(isCluster())
476  theQValue = -ParticleTable::getTableQValue(theA, theZ, ADaughter, ZDaughter);
477  else {
478  const G4double massTableParent = ParticleTable::getTableMass(AParent,ZParent);
479  const G4double massTableDaughter = ParticleTable::getTableMass(ADaughter,ZDaughter);
480  const G4double massTableParticle = getTableMass();
481  theQValue = massTableParent - massTableDaughter - massTableParticle;
482  }
483 
484  const G4double massINCLParent = ParticleTable::getINCLMass(AParent,ZParent);
485  const G4double massINCLDaughter = ParticleTable::getINCLMass(ADaughter,ZDaughter);
486  const G4double massINCLParticle = getINCLMass();
487 
488  // The rhs corresponds to the INCL Q-value
489  return theQValue - (massINCLParent-massINCLDaughter-massINCLParticle);
490  }
491 
507  G4double getTransferQValueCorrection(const G4int AFrom, const G4int ZFrom, const G4int ATo, const G4int ZTo) const {
508  const G4int AFromDaughter = AFrom - theA;
509  const G4int ZFromDaughter = ZFrom - theZ;
510  const G4int AToDaughter = ATo + theA;
511  const G4int ZToDaughter = ZTo + theZ;
512  const G4double theQValue = ParticleTable::getTableQValue(AToDaughter,ZToDaughter,AFromDaughter,ZFromDaughter,AFrom,ZFrom);
513 
514  const G4double massINCLTo = ParticleTable::getINCLMass(ATo,ZTo);
515  const G4double massINCLToDaughter = ParticleTable::getINCLMass(AToDaughter,ZToDaughter);
516  /* Note that here we have to use the table mass in the INCL Q-value. We
517  * cannot use theMass, because at this stage the particle is probably
518  * still off-shell; and we cannot use getINCLMass(), because it leads to
519  * violations of global energy conservation.
520  */
521  const G4double massINCLParticle = getTableMass();
522 
523  // The rhs corresponds to the INCL Q-value for particle absorption
524  return theQValue - (massINCLToDaughter-massINCLTo-massINCLParticle);
525  }
526 
533  const G4double mass = std::pow(theEnergy, 2) - theMomentum.dot(theMomentum);
534  if(mass < 0.0) {
535  INCL_ERROR("E*E - p*p is negative." << '\n');
536  return 0.0;
537  } else {
538  return std::sqrt(mass);
539  }
540  };
541 
543  inline G4double getKineticEnergy() const { return theEnergy - theMass; }
544 
546  inline G4double getPotentialEnergy() const { return thePotentialEnergy; }
547 
550 
555  {
556  return theEnergy;
557  };
558 
562  void setMass(G4double mass)
563  {
564  this->theMass = mass;
565  }
566 
570  void setEnergy(G4double energy)
571  {
572  this->theEnergy = energy;
573  };
574 
579  {
580  return theMomentum;
581  };
582 
585  {
587  };
588 
592  virtual void setMomentum(const G4INCL::ThreeVector &momentum)
593  {
594  this->theMomentum = momentum;
595  };
596 
601  {
602  return thePosition;
603  };
604 
605  virtual void setPosition(const G4INCL::ThreeVector &position)
606  {
607  this->thePosition = position;
608  };
609 
610  G4double getHelicity() { return theHelicity; };
611  void setHelicity(G4double h) { theHelicity = h; };
612 
613  void propagate(G4double step) {
614  thePosition += ((*thePropagationMomentum)*(step/(*thePropagationEnergy)));
615  };
616 
619 
622 
625 
627  G4int getNumberOfDecays() const { return nDecays; }
628 
631 
634 
643  void setOutOfWell() { outOfWell = true; }
644 
646  G4bool isOutOfWell() const { return outOfWell; }
647 
648  void setEmissionTime(G4double t) { emissionTime = t; }
649  G4double getEmissionTime() { return emissionTime; };
650 
654  }
655 
659  }
660 
663 
666 
667  G4bool isCluster() const {
668  return (theType == Composite);
669  }
670 
672  void setFrozenMomentum(const ThreeVector &momentum) { theFrozenMomentum = momentum; }
673 
675  void setFrozenEnergy(const G4double energy) { theFrozenEnergy = energy; }
676 
679 
682 
684  ThreeVector getPropagationVelocity() const { return (*thePropagationMomentum)/(*thePropagationEnergy); }
685 
695  }
696 
706  }
707 
713  virtual void rotatePositionAndMomentum(const G4double angle, const ThreeVector &axis) {
714  rotatePosition(angle, axis);
715  rotateMomentum(angle, axis);
716  }
717 
723  virtual void rotatePosition(const G4double angle, const ThreeVector &axis) {
724  thePosition.rotate(angle, axis);
725  }
726 
732  virtual void rotateMomentum(const G4double angle, const ThreeVector &axis) {
733  theMomentum.rotate(angle, axis);
734  theFrozenMomentum.rotate(angle, axis);
735  }
736 
737  std::string print() const {
738  std::stringstream ss;
739  ss << "Particle (ID = " << ID << ") type = ";
741  ss << '\n'
742  << " energy = " << theEnergy << '\n'
743  << " momentum = "
744  << theMomentum.print()
745  << '\n'
746  << " position = "
747  << thePosition.print()
748  << '\n';
749  return ss.str();
750  };
751 
752  std::string dump() const {
753  std::stringstream ss;
754  ss << "(particle " << ID << " ";
756  ss << '\n'
757  << thePosition.dump()
758  << '\n'
759  << theMomentum.dump()
760  << '\n'
761  << theEnergy << ")" << '\n';
762  return ss.str();
763  };
764 
765  long getID() const { return ID; };
766 
770  ParticleList const *getParticles() const {
771  INCL_WARN("Particle::getParticles() method was called on a Particle object" << '\n');
772  return 0;
773  }
774 
782  if(rpCorrelated)
783  return theMomentum.mag();
784  else
785  return uncorrelatedMomentum;
786  }
787 
790 
792  void rpCorrelate() { rpCorrelated = true; }
793 
795  void rpDecorrelate() { rpCorrelated = false; }
796 
800  if(norm>0.)
801  return thePosition.dot(*thePropagationMomentum) / std::sqrt(norm);
802  else
803  return 1.;
804  }
805 
806  protected:
820  long ID;
821 
824 
825  private:
826  G4double theHelicity;
827  G4double emissionTime;
828  G4bool outOfWell;
829 
830  G4double theMass;
831  static G4ThreadLocal long nextID;
832 
834  };
835 }
836 
837 #endif /* PARTICLE_HH_ */
G4bool isEta() const
Is this a eta?
G4int getA() const
Returns the baryon number.
ParticipantType theParticipantType
void rotatePosition(const G4double angle, const ThreeVector &axis) const
#define INCL_DECLARE_ALLOCATION_POOL(T)
void setParticipantType(ParticipantType const p)
G4double getCosRPAngle() const
Get the cosine of the angle between position and momentum.
G4int getNumberOfCollisions() const
Return the number of collisions undergone by the particle.
void setMass(G4double mass)
G4bool isResonance() const
Is it a resonance?
virtual void makeTargetSpectator()
G4double getReflectionMomentum() const
Return the reflection momentum.
G4double getMass() const
Get the cached particle mass.
G4double * thePropagationEnergy
G4double getBeta() const
G4ThreadLocal ParticleMassFn getTableParticleMass
Static pointer to the mass function for particles.
G4double dot(const ThreeVector &v) const
void rotate(const G4double angle, const ThreeVector &axis)
Rotate the vector by a given angle around a given axis.
G4double getEmissionTime()
ThreeVector getPropagationVelocity() const
Get the propagation velocity of the particle.
void setEmissionTime(G4double t)
G4INCL::ThreeVector theFrozenMomentum
const char * p
Definition: xmltok.h:285
G4double getTableQValue(const G4int A1, const G4int Z1, const G4int A2, const G4int Z2)
Get Q-value (in MeV/c^2)
#define INCL_ERROR(x)
G4bool isTargetSpectator() const
virtual void rotatePosition(const G4double angle, const ThreeVector &axis)
Rotate the particle position.
void boost(const ThreeVector &aBoostVector)
const G4INCL::ThreeVector & getMomentum() const
G4double adjustEnergyFromMomentum()
Recompute the energy to match the momentum.
G4bool isDelta() const
Is it a Delta?
std::string print() const
static G4double angle[DIM]
G4double getEnergy() const
G4bool isEtaPrime() const
Is this a etaprime?
void setINCLMass()
Set the mass of the Particle to its table mass.
void thawPropagation()
Unfreeze particle propagation.
#define INCL_WARN(x)
virtual void makeProjectileSpectator()
void propagate(G4double step)
Singleton for recycling allocation of instances of a given class.
void setOutOfWell()
Mark the particle as out of its potential well.
G4double getINCLMass() const
Get the INCL particle mass.
void incrementNumberOfCollisions()
Increment the number of collisions undergone by the particle.
virtual void rotateMomentum(const G4double angle, const ThreeVector &axis)
Rotate the particle momentum.
G4bool isParticipant() const
void setUncorrelatedMomentum(const G4double p)
Set the uncorrelated momentum.
G4double uncorrelatedMomentum
#define G4ThreadLocal
Definition: tls.hh:89
ThreeVector boostVector() const
int G4int
Definition: G4Types.hh:78
ThreeVector vector(const ThreeVector &v) const
G4double mag2() const
G4double getPotentialEnergy() const
Get the particle potential energy.
ParticleList const * getParticles() const
static double P[]
void incrementNumberOfDecays()
Increment the number of decays undergone by the particle.
G4double getInvariantMass() const
Get the the particle invariant mass.
void rpCorrelate()
Make the particle follow a strict r-p correlation.
void setEnergy(G4double energy)
void boost(const ThreeVector &b) const
G4double getRealMass(const G4INCL::ParticleType t)
Get particle mass (in MeV/c^2)
G4bool isOmega() const
Is this a omega?
#define position
Definition: xmlparse.cc:622
void rpDecorrelate()
Make the particle not follow a strict r-p correlation.
void setNumberOfDecays(G4int n)
Set the number of decays undergone by the particle.
G4double getTransferQValueCorrection(const G4int AFrom, const G4int ZFrom, const G4int ATo, const G4int ZTo) const
Computes correction on the transfer Q-value.
ParticipantType getParticipantType() const
bool G4bool
Definition: G4Types.hh:79
void lorentzContract(const ThreeVector &aBoostVector, const ThreeVector &refPos)
Lorentz-contract the particle position around some center.
virtual G4INCL::ThreeVector getAngularMomentum() const
virtual void setPosition(const G4INCL::ThreeVector &position)
G4int getZ() const
Returns the charge number.
void setPotentialEnergy(G4double v)
Set the particle potential energy.
std::string dump() const
G4bool isOutOfWell() const
Check if the particle is out of its potential well.
G4INCL::ThreeVector thePosition
ThreeVector getFrozenMomentum() const
Get the frozen particle momentum.
ThreeVector getTransversePosition() const
Transverse component of the position w.r.t. the momentum.
Particle(const Particle &rhs)
Copy constructor.
G4INCL::ParticleType theType
void setFrozenMomentum(const ThreeVector &momentum)
Set the frozen particle momentum.
G4double getINCLMass(const G4int A, const G4int Z)
Get INCL nuclear mass (in MeV/c^2)
G4int getNumberOfDecays() const
Return the number of decays undergone by the particle.
G4INCL::ParticleType getType() const
G4double energy(const ThreeVector &p, const G4double m)
G4ThreadLocal NuclearMassFn getTableMass
Static pointer to the mass function for nuclei.
void setTableMass()
Set the mass of the Particle to its table mass.
const G4INCL::ThreeVector & getPosition() const
void setRealMass()
Set the mass of the Particle to its real mass.
void swap(Particle &rhs)
Helper method for the assignment operator.
G4bool isProjectileSpectator() const
void setType(ParticleType t)
static const G4double bp
std::string dump() const
void setNumberOfCollisions(G4int n)
Set the number of collisions undergone by the particle.
G4double getHelicity()
Particle & operator=(const Particle &rhs)
Assignment operator.
void rotatePositionAndMomentum(const G4double angle, const ThreeVector &axis) const
G4bool isNucleon() const
G4bool isCluster() const
G4INCL::ThreeVector theMomentum
virtual G4double getTableMass() const
Get the tabulated particle mass.
std::string print() const
G4double getKineticEnergy() const
Get the particle kinetic energy.
virtual G4INCL::ParticleSpecies getSpecies() const
Get the particle species.
G4double getRealMass() const
Get the real particle mass.
virtual void makeParticipant()
G4double mag() const
double G4double
Definition: G4Types.hh:76
G4double theFrozenEnergy
G4double getEmissionQValueCorrection(const G4int AParent, const G4int ZParent) const
Computes correction on the emission Q-value.
void rotateMomentum(const G4double angle, const ThreeVector &axis) const
G4INCL::ThreeVector * thePropagationMomentum
static const G4double alpha
G4bool isPion() const
Is this a pion?
const ThreeVector & adjustMomentumFromEnergy()
Rescale the momentum to match the total energy.
ParticleList::iterator ParticleMutableIter
virtual void rotatePositionAndMomentum(const G4double angle, const ThreeVector &axis)
Rotate the particle position and momentum.
G4double thePotentialEnergy
ThreeVector getLongitudinalPosition() const
Longitudinal component of the position w.r.t. the momentum.
ParticleList::const_iterator ParticleIter
std::string getName(const ParticleType t)
Get the native INCL name of the particle.
long getID() const
void setHelicity(G4double h)
virtual void setMomentum(const G4INCL::ThreeVector &momentum)
void setFrozenEnergy(const G4double energy)
Set the frozen particle momentum.
void freezePropagation()
Freeze particle propagation.
G4double getFrozenEnergy() const
Get the frozen particle momentum.