Geant4  10.02.p03
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),
87  theParticipantType(rhs.theParticipantType),
88  theType(rhs.theType),
89  theEnergy(rhs.theEnergy),
90  theFrozenEnergy(rhs.theFrozenEnergy),
91  theMomentum(rhs.theMomentum),
92  theFrozenMomentum(rhs.theFrozenMomentum),
93  thePosition(rhs.thePosition),
94  nCollisions(rhs.nCollisions),
95  nDecays(rhs.nDecays),
96  thePotentialEnergy(rhs.thePotentialEnergy),
97  rpCorrelated(rhs.rpCorrelated),
98  uncorrelatedMomentum(rhs.uncorrelatedMomentum),
99  theHelicity(rhs.theHelicity),
100  emissionTime(rhs.emissionTime),
101  outOfWell(rhs.outOfWell),
102  theMass(rhs.theMass)
103  {
104  if(rhs.thePropagationEnergy == &(rhs.theFrozenEnergy))
105  thePropagationEnergy = &theFrozenEnergy;
106  else
107  thePropagationEnergy = &theEnergy;
108  if(rhs.thePropagationMomentum == &(rhs.theFrozenMomentum))
109  thePropagationMomentum = &theFrozenMomentum;
110  else
111  thePropagationMomentum = &theMomentum;
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))
124  thePropagationEnergy = &theFrozenEnergy;
125  else
126  thePropagationEnergy = &theEnergy;
127  std::swap(theEnergy, rhs.theEnergy);
128  std::swap(theFrozenEnergy, rhs.theFrozenEnergy);
129  if(rhs.thePropagationMomentum == &(rhs.theFrozenMomentum))
130  thePropagationMomentum = &theFrozenMomentum;
131  else
132  thePropagationMomentum = &theMomentum;
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);
147  std::swap(uncorrelatedMomentum, rhs.uncorrelatedMomentum);
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  theA = 0;
203  theZ = 0;
204  break;
205  case PiMinus:
206  theA = 0;
207  theZ = -1;
208  break;
209  case Composite:
210  // INCL_ERROR("Trying to set particle type to Composite! Construct a Cluster object instead" << '\n');
211  theA = 0;
212  theZ = 0;
213  break;
214  case UnknownParticle:
215  theA = 0;
216  theZ = 0;
217  INCL_ERROR("Trying to set particle type to Unknown!" << '\n');
218  break;
219  }
220 
221  if( !isResonance() && t!=Composite )
222  setINCLMass();
223  }
224 
228  G4bool isNucleon() const {
229  if(theType == G4INCL::Proton || theType == G4INCL::Neutron)
230  return true;
231  else
232  return false;
233  };
234 
236  return theParticipantType;
237  }
238 
240  theParticipantType = p;
241  }
242 
244  return (theParticipantType==Participant);
245  }
246 
248  return (theParticipantType==TargetSpectator);
249  }
250 
252  return (theParticipantType==ProjectileSpectator);
253  }
254 
255  virtual void makeParticipant() {
256  theParticipantType = Participant;
257  }
258 
259  virtual void makeTargetSpectator() {
260  theParticipantType = TargetSpectator;
261  }
262 
263  virtual void makeProjectileSpectator() {
264  theParticipantType = ProjectileSpectator;
265  }
266 
268  G4bool isPion() const { return (theType == PiPlus || theType == PiZero || theType == PiMinus); }
269 
271  inline G4bool isResonance() const { return isDelta(); }
272 
274  inline G4bool isDelta() const {
275  return (theType==DeltaPlusPlus || theType==DeltaPlus ||
276  theType==DeltaZero || theType==DeltaMinus);
277  }
278 
280  G4int getA() const { return theA; }
281 
283  G4int getZ() const { return theZ; }
284 
285  G4double getBeta() const {
286  const G4double P = theMomentum.mag();
287  return P/theEnergy;
288  }
289 
297  return theMomentum / theEnergy;
298  }
299 
306  void boost(const ThreeVector &aBoostVector) {
307  const G4double beta2 = aBoostVector.mag2();
308  const G4double gamma = 1.0 / std::sqrt(1.0 - beta2);
309  const G4double bp = theMomentum.dot(aBoostVector);
310  const G4double alpha = (gamma*gamma)/(1.0 + gamma);
311 
312  theMomentum = theMomentum + aBoostVector * (alpha * bp - gamma * theEnergy);
313  theEnergy = gamma * (theEnergy - bp);
314  }
315 
324  void lorentzContract(const ThreeVector &aBoostVector, const ThreeVector &refPos) {
325  const G4double beta2 = aBoostVector.mag2();
326  const G4double gamma = 1.0 / std::sqrt(1.0 - beta2);
327  const ThreeVector theRelativePosition = thePosition - refPos;
328  const ThreeVector transversePosition = theRelativePosition - aBoostVector * (theRelativePosition.dot(aBoostVector) / aBoostVector.mag2());
329  const ThreeVector longitudinalPosition = theRelativePosition - transversePosition;
330 
331  thePosition = refPos + transversePosition + longitudinalPosition / gamma;
332  }
333 
335  inline G4double getMass() const { return theMass; }
336 
338  inline G4double getINCLMass() const {
339  switch(theType) {
340  case Proton:
341  case Neutron:
342  case PiPlus:
343  case PiMinus:
344  case PiZero:
345  return ParticleTable::getINCLMass(theType);
346  break;
347 
348  case DeltaPlusPlus:
349  case DeltaPlus:
350  case DeltaZero:
351  case DeltaMinus:
352  return theMass;
353  break;
354 
355  case Composite:
356  return ParticleTable::getINCLMass(theA,theZ);
357  break;
358 
359  default:
360  INCL_ERROR("Particle::getINCLMass: Unknown particle type." << '\n');
361  return 0.0;
362  break;
363  }
364  }
365 
367  inline virtual G4double getTableMass() const {
368  switch(theType) {
369  case Proton:
370  case Neutron:
371  case PiPlus:
372  case PiMinus:
373  case PiZero:
374  return ParticleTable::getTableParticleMass(theType);
375  break;
376 
377  case DeltaPlusPlus:
378  case DeltaPlus:
379  case DeltaZero:
380  case DeltaMinus:
381  return theMass;
382  break;
383 
384  case Composite:
385  return ParticleTable::getTableMass(theA,theZ);
386  break;
387 
388  default:
389  INCL_ERROR("Particle::getTableMass: Unknown particle type." << '\n');
390  return 0.0;
391  break;
392  }
393  }
394 
396  inline G4double getRealMass() const {
397  switch(theType) {
398  case Proton:
399  case Neutron:
400  case PiPlus:
401  case PiMinus:
402  case PiZero:
403  return ParticleTable::getRealMass(theType);
404  break;
405 
406  case DeltaPlusPlus:
407  case DeltaPlus:
408  case DeltaZero:
409  case DeltaMinus:
410  return theMass;
411  break;
412 
413  case Composite:
414  return ParticleTable::getRealMass(theA,theZ);
415  break;
416 
417  default:
418  INCL_ERROR("Particle::getRealMass: Unknown particle type." << '\n');
419  return 0.0;
420  break;
421  }
422  }
423 
425  void setRealMass() { setMass(getRealMass()); }
426 
428  void setTableMass() { setMass(getTableMass()); }
429 
431  void setINCLMass() { setMass(getINCLMass()); }
432 
444  G4double getEmissionQValueCorrection(const G4int AParent, const G4int ZParent) const {
445  const G4int ADaughter = AParent - theA;
446  const G4int ZDaughter = ZParent - theZ;
447 
448  // Note the minus sign here
449  G4double theQValue;
450  if(isCluster())
451  theQValue = -ParticleTable::getTableQValue(theA, theZ, ADaughter, ZDaughter);
452  else {
453  const G4double massTableParent = ParticleTable::getTableMass(AParent,ZParent);
454  const G4double massTableDaughter = ParticleTable::getTableMass(ADaughter,ZDaughter);
455  const G4double massTableParticle = getTableMass();
456  theQValue = massTableParent - massTableDaughter - massTableParticle;
457  }
458 
459  const G4double massINCLParent = ParticleTable::getINCLMass(AParent,ZParent);
460  const G4double massINCLDaughter = ParticleTable::getINCLMass(ADaughter,ZDaughter);
461  const G4double massINCLParticle = getINCLMass();
462 
463  // The rhs corresponds to the INCL Q-value
464  return theQValue - (massINCLParent-massINCLDaughter-massINCLParticle);
465  }
466 
482  G4double getTransferQValueCorrection(const G4int AFrom, const G4int ZFrom, const G4int ATo, const G4int ZTo) const {
483  const G4int AFromDaughter = AFrom - theA;
484  const G4int ZFromDaughter = ZFrom - theZ;
485  const G4int AToDaughter = ATo + theA;
486  const G4int ZToDaughter = ZTo + theZ;
487  const G4double theQValue = ParticleTable::getTableQValue(AToDaughter,ZToDaughter,AFromDaughter,ZFromDaughter,AFrom,ZFrom);
488 
489  const G4double massINCLTo = ParticleTable::getINCLMass(ATo,ZTo);
490  const G4double massINCLToDaughter = ParticleTable::getINCLMass(AToDaughter,ZToDaughter);
491  /* Note that here we have to use the table mass in the INCL Q-value. We
492  * cannot use theMass, because at this stage the particle is probably
493  * still off-shell; and we cannot use getINCLMass(), because it leads to
494  * violations of global energy conservation.
495  */
496  const G4double massINCLParticle = getTableMass();
497 
498  // The rhs corresponds to the INCL Q-value for particle absorption
499  return theQValue - (massINCLToDaughter-massINCLTo-massINCLParticle);
500  }
501 
508  const G4double mass = std::pow(theEnergy, 2) - theMomentum.dot(theMomentum);
509  if(mass < 0.0) {
510  INCL_ERROR("E*E - p*p is negative." << '\n');
511  return 0.0;
512  } else {
513  return std::sqrt(mass);
514  }
515  };
516 
518  inline G4double getKineticEnergy() const { return theEnergy - theMass; }
519 
521  inline G4double getPotentialEnergy() const { return thePotentialEnergy; }
522 
524  inline void setPotentialEnergy(G4double v) { thePotentialEnergy = v; }
525 
530  {
531  return theEnergy;
532  };
533 
537  void setMass(G4double mass)
538  {
539  this->theMass = mass;
540  }
541 
545  void setEnergy(G4double energy)
546  {
547  this->theEnergy = energy;
548  };
549 
554  {
555  return theMomentum;
556  };
557 
560  {
561  return thePosition.vector(theMomentum);
562  };
563 
567  virtual void setMomentum(const G4INCL::ThreeVector &momentum)
568  {
569  this->theMomentum = momentum;
570  };
571 
576  {
577  return thePosition;
578  };
579 
580  virtual void setPosition(const G4INCL::ThreeVector &position)
581  {
582  this->thePosition = position;
583  };
584 
585  G4double getHelicity() { return theHelicity; };
586  void setHelicity(G4double h) { theHelicity = h; };
587 
588  void propagate(G4double step) {
589  thePosition += ((*thePropagationMomentum)*(step/(*thePropagationEnergy)));
590  };
591 
593  G4int getNumberOfCollisions() const { return nCollisions; }
594 
596  void setNumberOfCollisions(G4int n) { nCollisions = n; }
597 
599  void incrementNumberOfCollisions() { nCollisions++; }
600 
602  G4int getNumberOfDecays() const { return nDecays; }
603 
605  void setNumberOfDecays(G4int n) { nDecays = n; }
606 
608  void incrementNumberOfDecays() { nDecays++; }
609 
618  void setOutOfWell() { outOfWell = true; }
619 
621  G4bool isOutOfWell() const { return outOfWell; }
622 
623  void setEmissionTime(G4double t) { emissionTime = t; }
624  G4double getEmissionTime() { return emissionTime; };
625 
628  return thePosition - getLongitudinalPosition();
629  }
630 
633  return *thePropagationMomentum * (thePosition.dot(*thePropagationMomentum)/thePropagationMomentum->mag2());
634  }
635 
637  const ThreeVector &adjustMomentumFromEnergy();
638 
640  G4double adjustEnergyFromMomentum();
641 
642  G4bool isCluster() const {
643  return (theType == Composite);
644  }
645 
647  void setFrozenMomentum(const ThreeVector &momentum) { theFrozenMomentum = momentum; }
648 
650  void setFrozenEnergy(const G4double energy) { theFrozenEnergy = energy; }
651 
653  ThreeVector getFrozenMomentum() const { return theFrozenMomentum; }
654 
656  G4double getFrozenEnergy() const { return theFrozenEnergy; }
657 
659  ThreeVector getPropagationVelocity() const { return (*thePropagationMomentum)/(*thePropagationEnergy); }
660 
668  thePropagationMomentum = &theFrozenMomentum;
669  thePropagationEnergy = &theFrozenEnergy;
670  }
671 
679  thePropagationMomentum = &theMomentum;
680  thePropagationEnergy = &theEnergy;
681  }
682 
688  virtual void rotatePositionAndMomentum(const G4double angle, const ThreeVector &axis) {
689  rotatePosition(angle, axis);
690  rotateMomentum(angle, axis);
691  }
692 
698  virtual void rotatePosition(const G4double angle, const ThreeVector &axis) {
699  thePosition.rotate(angle, axis);
700  }
701 
707  virtual void rotateMomentum(const G4double angle, const ThreeVector &axis) {
708  theMomentum.rotate(angle, axis);
709  theFrozenMomentum.rotate(angle, axis);
710  }
711 
712  std::string print() const {
713  std::stringstream ss;
714  ss << "Particle (ID = " << ID << ") type = ";
715  ss << ParticleTable::getName(theType);
716  ss << '\n'
717  << " energy = " << theEnergy << '\n'
718  << " momentum = "
719  << theMomentum.print()
720  << '\n'
721  << " position = "
722  << thePosition.print()
723  << '\n';
724  return ss.str();
725  };
726 
727  std::string dump() const {
728  std::stringstream ss;
729  ss << "(particle " << ID << " ";
730  ss << ParticleTable::getName(theType);
731  ss << '\n'
732  << thePosition.dump()
733  << '\n'
734  << theMomentum.dump()
735  << '\n'
736  << theEnergy << ")" << '\n';
737  return ss.str();
738  };
739 
740  long getID() const { return ID; };
741 
745  ParticleList const *getParticles() const {
746  INCL_WARN("Particle::getParticles() method was called on a Particle object" << '\n');
747  return 0;
748  }
749 
757  if(rpCorrelated)
758  return theMomentum.mag();
759  else
760  return uncorrelatedMomentum;
761  }
762 
764  void setUncorrelatedMomentum(const G4double p) { uncorrelatedMomentum = p; }
765 
767  void rpCorrelate() { rpCorrelated = true; }
768 
770  void rpDecorrelate() { rpCorrelated = false; }
771 
774  const G4double norm = thePosition.mag2()*thePropagationMomentum->mag2();
775  if(norm>0.)
776  return thePosition.dot(*thePropagationMomentum) / std::sqrt(norm);
777  else
778  return 1.;
779  }
780 
781  protected:
782  G4int theZ, theA;
795  long ID;
796 
799 
800  private:
804 
806  static G4ThreadLocal long nextID;
807 
809  };
810 }
811 
812 #endif /* PARTICLE_HH_ */
ParticipantType theParticipantType
G4bool isNucleon() const
#define INCL_DECLARE_ALLOCATION_POOL(T)
void setParticipantType(ParticipantType const p)
G4double getBeta() const
void setMass(G4double mass)
ThreeVector boostVector() const
G4bool isOutOfWell() const
Check if the particle is out of its potential well.
virtual void makeTargetSpectator()
G4double * thePropagationEnergy
G4ThreadLocal ParticleMassFn getTableParticleMass
Static pointer to the mass function for particles.
G4double getEnergy() const
G4bool isPion() const
Is this a pion?
ThreeVector getFrozenMomentum() const
Get the frozen particle momentum.
G4double mag2() const
G4double getINCLMass() const
Get the INCL particle mass.
const G4INCL::ThreeVector & getPosition() const
G4double getEmissionTime()
G4double getMass() const
Get the cached particle mass.
void setEmissionTime(G4double t)
G4INCL::ThreeVector theFrozenMomentum
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)
G4double getReflectionMomentum() const
Return the reflection momentum.
virtual void rotatePosition(const G4double angle, const ThreeVector &axis)
Rotate the particle position.
void boost(const ThreeVector &aBoostVector)
Float_t norm
static G4double angle[DIM]
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.
void rotateMomentum(const G4double angle, const ThreeVector &axis) const
void incrementNumberOfCollisions()
Increment the number of collisions undergone by the particle.
virtual void rotateMomentum(const G4double angle, const ThreeVector &axis)
Rotate the particle momentum.
void setUncorrelatedMomentum(const G4double p)
Set the uncorrelated momentum.
G4int getA() const
Returns the baryon number.
G4double uncorrelatedMomentum
G4int getZ() const
Returns the charge number.
G4double getKineticEnergy() const
Get the particle kinetic energy.
#define G4ThreadLocal
Definition: tls.hh:89
int G4int
Definition: G4Types.hh:78
std::string print() const
static double P[]
void incrementNumberOfDecays()
Increment the number of decays undergone by the particle.
void rpCorrelate()
Make the particle follow a strict r-p correlation.
void setEnergy(G4double energy)
G4double getFrozenEnergy() const
Get the frozen particle momentum.
ThreeVector vector(const ThreeVector &v) const
G4double getRealMass(const G4INCL::ParticleType t)
Get particle mass (in MeV/c^2)
G4int getNumberOfCollisions() const
Return the number of collisions undergone by the particle.
void boost(const ThreeVector &b) const
Char_t n[5]
G4bool isTargetSpectator() const
#define position
Definition: xmlparse.cc:622
double energy
Definition: plottest35.C:25
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.
G4INCL::ParticleType getType() const
ParticleList const * getParticles() const
G4bool isCluster() const
bool G4bool
Definition: G4Types.hh:79
void lorentzContract(const ThreeVector &aBoostVector, const ThreeVector &refPos)
Lorentz-contract the particle position around some center.
virtual void setPosition(const G4INCL::ThreeVector &position)
G4bool isDelta() const
Is it a Delta?
void setPotentialEnergy(G4double v)
Set the particle potential energy.
static G4ThreadLocal long nextID
virtual G4double getTableMass() const
Get the tabulated particle mass.
G4double dot(const ThreeVector &v) const
G4int getNumberOfDecays() const
Return the number of decays undergone by the particle.
G4INCL::ThreeVector thePosition
std::string dump() const
G4double getTransferQValueCorrection(const G4int AFrom, const G4int ZFrom, const G4int ATo, const G4int ZTo) const
Computes correction on the transfer Q-value.
Particle(const Particle &rhs)
Copy constructor.
void rotatePosition(const G4double angle, const ThreeVector &axis) const
G4INCL::ParticleType theType
void setFrozenMomentum(const ThreeVector &momentum)
Set the frozen particle momentum.
G4double getPotentialEnergy() const
Get the particle potential energy.
G4double getINCLMass(const G4int A, const G4int Z)
Get INCL nuclear mass (in MeV/c^2)
G4ThreadLocal NuclearMassFn getTableMass
Static pointer to the mass function for nuclei.
void setTableMass()
Set the mass of the Particle to its table mass.
void setRealMass()
Set the mass of the Particle to its real mass.
void swap(Particle &rhs)
Helper method for the assignment operator.
ThreeVector getTransversePosition() const
Transverse component of the position w.r.t. the momentum.
ThreeVector getLongitudinalPosition() const
Longitudinal component of the position w.r.t. the momentum.
void rotatePositionAndMomentum(const G4double angle, const ThreeVector &axis) const
void setType(ParticleType t)
G4bool isProjectileSpectator() const
static const G4double bp
void setNumberOfCollisions(G4int n)
Set the number of collisions undergone by the particle.
G4double getRealMass() const
Get the real particle mass.
G4double getHelicity()
G4double getCosRPAngle() const
Get the cosine of the angle between position and momentum.
Particle & operator=(const Particle &rhs)
Assignment operator.
G4INCL::ThreeVector theMomentum
virtual G4INCL::ParticleSpecies getSpecies() const
Get the particle species.
virtual void makeParticipant()
double G4double
Definition: G4Types.hh:76
G4double theFrozenEnergy
const G4INCL::ThreeVector & getMomentum() const
G4INCL::ThreeVector * thePropagationMomentum
static const G4double alpha
G4double getInvariantMass() const
Get the the particle invariant mass.
ParticleList::iterator ParticleMutableIter
virtual void rotatePositionAndMomentum(const G4double angle, const ThreeVector &axis)
Rotate the particle position and momentum.
virtual G4INCL::ThreeVector getAngularMomentum() const
G4bool isResonance() const
Is it a resonance?
ThreeVector getPropagationVelocity() const
Get the propagation velocity of the particle.
G4double thePotentialEnergy
ParticipantType getParticipantType() const
ParticleList::const_iterator ParticleIter
G4bool isParticipant() const
G4double getEmissionQValueCorrection(const G4int AParent, const G4int ZParent) const
Computes correction on the emission Q-value.
long getID() const
std::string getName(const ParticleType t)
Get the native INCL name of the particle.
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.