Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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"
50#include "G4INCLParticleType.hh"
52#include "G4INCLLogger.hh"
55#include <sstream>
56#include <string>
57
58namespace 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;
69 std::vector<G4int> getParticleListBiasVector() const;
70 };
71
72 typedef ParticleList::const_iterator ParticleIter;
73 typedef ParticleList::iterator ParticleMutableIter;
74
75 class Particle {
76 public:
77 Particle();
78 Particle(ParticleType t, G4double energy, ThreeVector const &momentum, ThreeVector const &position);
79 Particle(ParticleType t, ThreeVector const &momentum, ThreeVector const &position);
80 virtual ~Particle() {}
81
82 /** \brief Copy constructor
83 *
84 * Does not copy the particle ID.
85 */
86 Particle(const Particle &rhs) :
87 theZ(rhs.theZ),
88 theA(rhs.theA),
89 theS(rhs.theS),
91 theType(rhs.theType),
98 nDecays(rhs.nDecays),
103 theNKaon(rhs.theNKaon),
106 theHelicity(rhs.theHelicity),
107 emissionTime(rhs.emissionTime),
108 outOfWell(rhs.outOfWell),
109 theMass(rhs.theMass)
110 {
111 if(rhs.thePropagationEnergy == &(rhs.theFrozenEnergy))
113 else
117 else
119 // ID intentionally not copied
120 ID = nextID++;
121
122 theBiasCollisionVector = rhs.theBiasCollisionVector;
123 }
124
125 protected:
126 /// \brief Helper method for the assignment operator
127 void swap(Particle &rhs) {
128 std::swap(theZ, rhs.theZ);
129 std::swap(theA, rhs.theA);
130 std::swap(theS, rhs.theS);
132 std::swap(theType, rhs.theType);
133 if(rhs.thePropagationEnergy == &(rhs.theFrozenEnergy))
135 else
137 std::swap(theEnergy, rhs.theEnergy);
138 std::swap(theFrozenEnergy, rhs.theFrozenEnergy);
141 else
143 std::swap(theMomentum, rhs.theMomentum);
144 std::swap(theFrozenMomentum, rhs.theFrozenMomentum);
145 std::swap(thePosition, rhs.thePosition);
146 std::swap(nCollisions, rhs.nCollisions);
147 std::swap(nDecays, rhs.nDecays);
149 // ID intentionally not swapped
150
153
154 std::swap(theHelicity, rhs.theHelicity);
155 std::swap(emissionTime, rhs.emissionTime);
156 std::swap(outOfWell, rhs.outOfWell);
157
158 std::swap(theMass, rhs.theMass);
159 std::swap(rpCorrelated, rhs.rpCorrelated);
161
162 std::swap(theParticleBias, rhs.theParticleBias);
163 std::swap(theBiasCollisionVector, rhs.theBiasCollisionVector);
164 }
165
166 public:
167
168 /** \brief Assignment operator
169 *
170 * Does not copy the particle ID.
171 */
173 Particle temporaryParticle(rhs);
174 swap(temporaryParticle);
175 return *this;
176 }
177
178 /**
179 * Get the particle type.
180 * @see G4INCL::ParticleType
181 */
183 return theType;
184 };
185
186 /// \brief Get the particle species
188 return ParticleSpecies(theType);
189 };
190
192 theType = t;
193 switch(theType)
194 {
195 case DeltaPlusPlus:
196 theA = 1;
197 theZ = 2;
198 theS = 0;
199 break;
200 case Proton:
201 case DeltaPlus:
202 theA = 1;
203 theZ = 1;
204 theS = 0;
205 break;
206 case Neutron:
207 case DeltaZero:
208 theA = 1;
209 theZ = 0;
210 theS = 0;
211 break;
212 case DeltaMinus:
213 theA = 1;
214 theZ = -1;
215 theS = 0;
216 break;
217 case PiPlus:
218 theA = 0;
219 theZ = 1;
220 theS = 0;
221 break;
222 case PiZero:
223 case Eta:
224 case Omega:
225 case EtaPrime:
226 case Photon:
227 theA = 0;
228 theZ = 0;
229 theS = 0;
230 break;
231 case PiMinus:
232 theA = 0;
233 theZ = -1;
234 theS = 0;
235 break;
236 case Lambda:
237 theA = 1;
238 theZ = 0;
239 theS = -1;
240 break;
241 case SigmaPlus:
242 theA = 1;
243 theZ = 1;
244 theS = -1;
245 break;
246 case SigmaZero:
247 theA = 1;
248 theZ = 0;
249 theS = -1;
250 break;
251 case SigmaMinus:
252 theA = 1;
253 theZ = -1;
254 theS = -1;
255 break;
256 case KPlus:
257 theA = 0;
258 theZ = 1;
259 theS = 1;
260 break;
261 case KZero:
262 theA = 0;
263 theZ = 0;
264 theS = 1;
265 break;
266 case KZeroBar:
267 theA = 0;
268 theZ = 0;
269 theS = -1;
270 break;
271 case KShort:
272 theA = 0;
273 theZ = 0;
274// theS should not be defined
275 break;
276 case KLong:
277 theA = 0;
278 theZ = 0;
279// theS should not be defined
280 break;
281 case KMinus:
282 theA = 0;
283 theZ = -1;
284 theS = -1;
285 break;
286 case Composite:
287 // INCL_ERROR("Trying to set particle type to Composite! Construct a Cluster object instead" << '\n');
288 theA = 0;
289 theZ = 0;
290 theS = 0;
291 break;
292 case UnknownParticle:
293 theA = 0;
294 theZ = 0;
295 theS = 0;
296 INCL_ERROR("Trying to set particle type to Unknown!" << '\n');
297 break;
298 }
299
300 if( !isResonance() && t!=Composite )
301 setINCLMass();
302 }
303
304 /**
305 * Is this a nucleon?
306 */
309 return true;
310 else
311 return false;
312 };
313
315 return theParticipantType;
316 }
317
320 }
321
324 }
325
328 }
329
332 }
333
334 virtual void makeParticipant() {
336 }
337
338 virtual void makeTargetSpectator() {
340 }
341
342 virtual void makeProjectileSpectator() {
344 }
345
346 /** \brief Is this a pion? */
347 G4bool isPion() const { return (theType == PiPlus || theType == PiZero || theType == PiMinus); }
348
349 /** \brief Is this an eta? */
350 G4bool isEta() const { return (theType == Eta); }
351
352 /** \brief Is this an omega? */
353 G4bool isOmega() const { return (theType == Omega); }
354
355 /** \brief Is this an etaprime? */
356 G4bool isEtaPrime() const { return (theType == EtaPrime); }
357
358 /** \brief Is this a photon? */
359 G4bool isPhoton() const { return (theType == Photon); }
360
361 /** \brief Is it a resonance? */
362 inline G4bool isResonance() const { return isDelta(); }
363
364 /** \brief Is it a Delta? */
365 inline G4bool isDelta() const {
366 return (theType==DeltaPlusPlus || theType==DeltaPlus ||
368
369 /** \brief Is this a Sigma? */
370 G4bool isSigma() const { return (theType == SigmaPlus || theType == SigmaZero || theType == SigmaMinus); }
371
372 /** \brief Is this a Kaon? */
373 G4bool isKaon() const { return (theType == KPlus || theType == KZero); }
374
375 /** \brief Is this an antiKaon? */
376 G4bool isAntiKaon() const { return (theType == KZeroBar || theType == KMinus); }
377
378 /** \brief Is this a Lambda? */
379 G4bool isLambda() const { return (theType == Lambda); }
380
381 /** \brief Is this a Nucleon or a Lambda? */
382 G4bool isNucleonorLambda() const { return (isNucleon() || isLambda()); }
383
384 /** \brief Is this an Hyperon? */
385 G4bool isHyperon() const { return (isLambda() || isSigma()); }
386
387 /** \brief Is this a Meson? */
388 G4bool isMeson() const { return (isPion() || isKaon() || isAntiKaon() || isEta() || isEtaPrime() || isOmega()); }
389
390 /** \brief Is this a Baryon? */
391 G4bool isBaryon() const { return (isNucleon() || isResonance() || isHyperon()); }
392
393 /** \brief Is this an Strange? */
394 G4bool isStrange() const { return (isKaon() || isAntiKaon() || isHyperon()); }
395
396 /** \brief Returns the baryon number. */
397 G4int getA() const { return theA; }
398
399 /** \brief Returns the charge number. */
400 G4int getZ() const { return theZ; }
401
402 /** \brief Returns the strangeness number. */
403 G4int getS() const { return theS; }
404
406 const G4double P = theMomentum.mag();
407 return P/theEnergy;
408 }
409
410 /**
411 * Returns a three vector we can give to the boost() -method.
412 *
413 * In order to go to the particle rest frame you need to multiply
414 * the boost vector by -1.0.
415 */
417 return theMomentum / theEnergy;
418 }
419
420 /**
421 * Boost the particle using a boost vector.
422 *
423 * Example (go to the particle rest frame):
424 * particle->boost(particle->boostVector());
425 */
426 void boost(const ThreeVector &aBoostVector) {
427 const G4double beta2 = aBoostVector.mag2();
428 const G4double gamma = 1.0 / std::sqrt(1.0 - beta2);
429 const G4double bp = theMomentum.dot(aBoostVector);
430 const G4double alpha = (gamma*gamma)/(1.0 + gamma);
431
432 theMomentum = theMomentum + aBoostVector * (alpha * bp - gamma * theEnergy);
433 theEnergy = gamma * (theEnergy - bp);
434 }
435
436 /** \brief Lorentz-contract the particle position around some center
437 *
438 * Apply Lorentz contraction to the position component along the
439 * direction of the boost vector.
440 *
441 * \param aBoostVector the boost vector (velocity) [c]
442 * \param refPos the reference position
443 */
444 void lorentzContract(const ThreeVector &aBoostVector, const ThreeVector &refPos) {
445 const G4double beta2 = aBoostVector.mag2();
446 const G4double gamma = 1.0 / std::sqrt(1.0 - beta2);
447 const ThreeVector theRelativePosition = thePosition - refPos;
448 const ThreeVector transversePosition = theRelativePosition - aBoostVector * (theRelativePosition.dot(aBoostVector) / aBoostVector.mag2());
449 const ThreeVector longitudinalPosition = theRelativePosition - transversePosition;
450
451 thePosition = refPos + transversePosition + longitudinalPosition / gamma;
452 }
453
454 /** \brief Get the cached particle mass. */
455 inline G4double getMass() const { return theMass; }
456
457 /** \brief Get the INCL particle mass. */
458 inline G4double getINCLMass() const {
459 switch(theType) {
460 case Proton:
461 case Neutron:
462 case PiPlus:
463 case PiMinus:
464 case PiZero:
465 case Lambda:
466 case SigmaPlus:
467 case SigmaZero:
468 case SigmaMinus:
469 case KPlus:
470 case KZero:
471 case KZeroBar:
472 case KShort:
473 case KLong:
474 case KMinus:
475 case Eta:
476 case Omega:
477 case EtaPrime:
478 case Photon:
480 break;
481
482 case DeltaPlusPlus:
483 case DeltaPlus:
484 case DeltaZero:
485 case DeltaMinus:
486 return theMass;
487 break;
488
489 case Composite:
491 break;
492
493 default:
494 INCL_ERROR("Particle::getINCLMass: Unknown particle type." << '\n');
495 return 0.0;
496 break;
497 }
498 }
499
500 /** \brief Get the tabulated particle mass. */
501 inline virtual G4double getTableMass() const {
502 switch(theType) {
503 case Proton:
504 case Neutron:
505 case PiPlus:
506 case PiMinus:
507 case PiZero:
508 case Lambda:
509 case SigmaPlus:
510 case SigmaZero:
511 case SigmaMinus:
512 case KPlus:
513 case KZero:
514 case KZeroBar:
515 case KShort:
516 case KLong:
517 case KMinus:
518 case Eta:
519 case Omega:
520 case EtaPrime:
521 case Photon:
523 break;
524
525 case DeltaPlusPlus:
526 case DeltaPlus:
527 case DeltaZero:
528 case DeltaMinus:
529 return theMass;
530 break;
531
532 case Composite:
534 break;
535
536 default:
537 INCL_ERROR("Particle::getTableMass: Unknown particle type." << '\n');
538 return 0.0;
539 break;
540 }
541 }
542
543 /** \brief Get the real particle mass. */
544 inline G4double getRealMass() const {
545 switch(theType) {
546 case Proton:
547 case Neutron:
548 case PiPlus:
549 case PiMinus:
550 case PiZero:
551 case Lambda:
552 case SigmaPlus:
553 case SigmaZero:
554 case SigmaMinus:
555 case KPlus:
556 case KZero:
557 case KZeroBar:
558 case KShort:
559 case KLong:
560 case KMinus:
561 case Eta:
562 case Omega:
563 case EtaPrime:
564 case Photon:
566 break;
567
568 case DeltaPlusPlus:
569 case DeltaPlus:
570 case DeltaZero:
571 case DeltaMinus:
572 return theMass;
573 break;
574
575 case Composite:
577 break;
578
579 default:
580 INCL_ERROR("Particle::getRealMass: Unknown particle type." << '\n');
581 return 0.0;
582 break;
583 }
584 }
585
586 /// \brief Set the mass of the Particle to its real mass
588
589 /// \brief Set the mass of the Particle to its table mass
591
592 /// \brief Set the mass of the Particle to its table mass
594
595 /**\brief Computes correction on the emission Q-value
596 *
597 * Computes the correction that must be applied to INCL particles in
598 * order to obtain the correct Q-value for particle emission from a given
599 * nucleus. For absorption, the correction is obviously equal to minus
600 * the value returned by this function.
601 *
602 * \param AParent the mass number of the emitting nucleus
603 * \param ZParent the charge number of the emitting nucleus
604 * \return the correction
605 */
606 G4double getEmissionQValueCorrection(const G4int AParent, const G4int ZParent) const {
607 const G4int SParent = 0;
608 const G4int ADaughter = AParent - theA;
609 const G4int ZDaughter = ZParent - theZ;
610 const G4int SDaughter = 0;
611
612 // Note the minus sign here
613 G4double theQValue;
614 if(isCluster())
615 theQValue = -ParticleTable::getTableQValue(theA, theZ, theS, ADaughter, ZDaughter, SDaughter);
616 else {
617 const G4double massTableParent = ParticleTable::getTableMass(AParent,ZParent,SParent);
618 const G4double massTableDaughter = ParticleTable::getTableMass(ADaughter,ZDaughter,SDaughter);
619 const G4double massTableParticle = getTableMass();
620 theQValue = massTableParent - massTableDaughter - massTableParticle;
621 }
622
623 const G4double massINCLParent = ParticleTable::getINCLMass(AParent,ZParent,SParent);
624 const G4double massINCLDaughter = ParticleTable::getINCLMass(ADaughter,ZDaughter,SDaughter);
625 const G4double massINCLParticle = getINCLMass();
626
627 // The rhs corresponds to the INCL Q-value
628 return theQValue - (massINCLParent-massINCLDaughter-massINCLParticle);
629 }
630
631 /**\brief Computes correction on the transfer Q-value
632 *
633 * Computes the correction that must be applied to INCL particles in
634 * order to obtain the correct Q-value for particle transfer from a given
635 * nucleus to another.
636 *
637 * Assumes that the receving nucleus is INCL's target nucleus, with the
638 * INCL separation energy.
639 *
640 * \param AFrom the mass number of the donating nucleus
641 * \param ZFrom the charge number of the donating nucleus
642 * \param ATo the mass number of the receiving nucleus
643 * \param ZTo the charge number of the receiving nucleus
644 * \return the correction
645 */
646 G4double getTransferQValueCorrection(const G4int AFrom, const G4int ZFrom, const G4int ATo, const G4int ZTo) const {
647 const G4int SFrom = 0;
648 const G4int STo = 0;
649 const G4int AFromDaughter = AFrom - theA;
650 const G4int ZFromDaughter = ZFrom - theZ;
651 const G4int SFromDaughter = 0;
652 const G4int AToDaughter = ATo + theA;
653 const G4int ZToDaughter = ZTo + theZ;
654 const G4int SToDaughter = 0;
655 const G4double theQValue = ParticleTable::getTableQValue(AToDaughter,ZToDaughter,SToDaughter,AFromDaughter,ZFromDaughter,SFromDaughter,AFrom,ZFrom,SFrom);
656
657 const G4double massINCLTo = ParticleTable::getINCLMass(ATo,ZTo,STo);
658 const G4double massINCLToDaughter = ParticleTable::getINCLMass(AToDaughter,ZToDaughter,SToDaughter);
659 /* Note that here we have to use the table mass in the INCL Q-value. We
660 * cannot use theMass, because at this stage the particle is probably
661 * still off-shell; and we cannot use getINCLMass(), because it leads to
662 * violations of global energy conservation.
663 */
664 const G4double massINCLParticle = getTableMass();
665
666 // The rhs corresponds to the INCL Q-value for particle absorption
667 return theQValue - (massINCLToDaughter-massINCLTo-massINCLParticle);
668 }
669
670 /**\brief Computes correction on the emission Q-value for hypernuclei
671 *
672 * Computes the correction that must be applied to INCL particles in
673 * order to obtain the correct Q-value for particle emission from a given
674 * nucleus. For absorption, the correction is obviously equal to minus
675 * the value returned by this function.
676 *
677 * \param AParent the mass number of the emitting nucleus
678 * \param ZParent the charge number of the emitting nucleus
679 * \param SParent the strangess number of the emitting nucleus
680 * \return the correction
681 */
682 G4double getEmissionQValueCorrection(const G4int AParent, const G4int ZParent, const G4int SParent) const {
683 const G4int ADaughter = AParent - theA;
684 const G4int ZDaughter = ZParent - theZ;
685 const G4int SDaughter = SParent - theS;
686
687 // Note the minus sign here
688 G4double theQValue;
689 if(isCluster())
690 theQValue = -ParticleTable::getTableQValue(theA, theZ, theS, ADaughter, ZDaughter, SDaughter);
691 else {
692 const G4double massTableParent = ParticleTable::getTableMass(AParent,ZParent,SParent);
693 const G4double massTableDaughter = ParticleTable::getTableMass(ADaughter,ZDaughter,SDaughter);
694 const G4double massTableParticle = getTableMass();
695 theQValue = massTableParent - massTableDaughter - massTableParticle;
696 }
697
698 const G4double massINCLParent = ParticleTable::getINCLMass(AParent,ZParent,SParent);
699 const G4double massINCLDaughter = ParticleTable::getINCLMass(ADaughter,ZDaughter,SDaughter);
700 const G4double massINCLParticle = getINCLMass();
701
702 // The rhs corresponds to the INCL Q-value
703 return theQValue - (massINCLParent-massINCLDaughter-massINCLParticle);
704 }
705
706 /**\brief Computes correction on the transfer Q-value for hypernuclei
707 *
708 * Computes the correction that must be applied to INCL particles in
709 * order to obtain the correct Q-value for particle transfer from a given
710 * nucleus to another.
711 *
712 * Assumes that the receving nucleus is INCL's target nucleus, with the
713 * INCL separation energy.
714 *
715 * \param AFrom the mass number of the donating nucleus
716 * \param ZFrom the charge number of the donating nucleus
717 * \param SFrom the strangess number of the donating nucleus
718 * \param ATo the mass number of the receiving nucleus
719 * \param ZTo the charge number of the receiving nucleus
720 * \param STo the strangess number of the receiving nucleus
721 * \return the correction
722 */
723 G4double getTransferQValueCorrection(const G4int AFrom, const G4int ZFrom, const G4int SFrom, const G4int ATo, const G4int ZTo , const G4int STo) const {
724 const G4int AFromDaughter = AFrom - theA;
725 const G4int ZFromDaughter = ZFrom - theZ;
726 const G4int SFromDaughter = SFrom - theS;
727 const G4int AToDaughter = ATo + theA;
728 const G4int ZToDaughter = ZTo + theZ;
729 const G4int SToDaughter = STo + theS;
730 const G4double theQValue = ParticleTable::getTableQValue(AToDaughter,ZToDaughter,SFromDaughter,AFromDaughter,ZFromDaughter,SToDaughter,AFrom,ZFrom,SFrom);
731
732 const G4double massINCLTo = ParticleTable::getINCLMass(ATo,ZTo,STo);
733 const G4double massINCLToDaughter = ParticleTable::getINCLMass(AToDaughter,ZToDaughter,SToDaughter);
734 /* Note that here we have to use the table mass in the INCL Q-value. We
735 * cannot use theMass, because at this stage the particle is probably
736 * still off-shell; and we cannot use getINCLMass(), because it leads to
737 * violations of global energy conservation.
738 */
739 const G4double massINCLParticle = getTableMass();
740
741 // The rhs corresponds to the INCL Q-value for particle absorption
742 return theQValue - (massINCLToDaughter-massINCLTo-massINCLParticle);
743 }
744
745
746
747 /** \brief Get the the particle invariant mass.
748 *
749 * Uses the relativistic invariant
750 * \f[ m = \sqrt{E^2 - {\vec p}^2}\f]
751 **/
753 const G4double mass = std::pow(theEnergy, 2) - theMomentum.dot(theMomentum);
754 if(mass < 0.0) {
755 INCL_ERROR("E*E - p*p is negative." << '\n');
756 return 0.0;
757 } else {
758 return std::sqrt(mass);
759 }
760 };
761
762 /// \brief Get the particle kinetic energy.
763 inline G4double getKineticEnergy() const { return theEnergy - theMass; }
764
765 /// \brief Get the particle potential energy.
767
768 /// \brief Set the particle potential energy.
770
771 /**
772 * Get the energy of the particle in MeV.
773 */
775 {
776 return theEnergy;
777 };
778
779 /**
780 * Set the mass of the particle in MeV/c^2.
781 */
782 void setMass(G4double mass)
783 {
784 this->theMass = mass;
785 }
786
787 /**
788 * Set the energy of the particle in MeV.
789 */
790 void setEnergy(G4double energy)
791 {
792 this->theEnergy = energy;
793 };
794
795 /**
796 * Get the momentum vector.
797 */
799 {
800 return theMomentum;
801 };
802
803 /** Get the angular momentum w.r.t. the origin */
805 {
807 };
808
809 /**
810 * Set the momentum vector.
811 */
812 virtual void setMomentum(const G4INCL::ThreeVector &momentum)
813 {
814 this->theMomentum = momentum;
815 };
816
817 /**
818 * Set the position vector.
819 */
821 {
822 return thePosition;
823 };
824
826 {
827 this->thePosition = position;
828 };
829
830 G4double getHelicity() { return theHelicity; };
831 void setHelicity(G4double h) { theHelicity = h; };
832
833 void propagate(G4double step) {
834 thePosition += ((*thePropagationMomentum)*(step/(*thePropagationEnergy)));
835 };
836
837 /** \brief Return the number of collisions undergone by the particle. **/
839
840 /** \brief Set the number of collisions undergone by the particle. **/
842
843 /** \brief Increment the number of collisions undergone by the particle. **/
845
846 /** \brief Return the number of decays undergone by the particle. **/
847 G4int getNumberOfDecays() const { return nDecays; }
848
849 /** \brief Set the number of decays undergone by the particle. **/
851
852 /** \brief Increment the number of decays undergone by the particle. **/
854
855 /** \brief Mark the particle as out of its potential well
856 *
857 * This flag is used to control pions created outside their potential well
858 * in delta decay. The pion potential checks it and returns zero if it is
859 * true (necessary in order to correctly enforce energy conservation). The
860 * Nucleus::applyFinalState() method uses it to determine whether new
861 * avatars should be generated for the particle.
862 */
863 void setOutOfWell() { outOfWell = true; }
864
865 /// \brief Check if the particle is out of its potential well
866 G4bool isOutOfWell() const { return outOfWell; }
867
868 void setEmissionTime(G4double t) { emissionTime = t; }
869 G4double getEmissionTime() { return emissionTime; };
870
871 /** \brief Transverse component of the position w.r.t. the momentum. */
874 }
875
876 /** \brief Longitudinal component of the position w.r.t. the momentum. */
879 }
880
881 /** \brief Rescale the momentum to match the total energy. */
883
884 /** \brief Recompute the energy to match the momentum. */
886
888 return (theType == Composite);
889 }
890
891 /// \brief Set the frozen particle momentum
892 void setFrozenMomentum(const ThreeVector &momentum) { theFrozenMomentum = momentum; }
893
894 /// \brief Set the frozen particle momentum
895 void setFrozenEnergy(const G4double energy) { theFrozenEnergy = energy; }
896
897 /// \brief Get the frozen particle momentum
899
900 /// \brief Get the frozen particle momentum
902
903 /// \brief Get the propagation velocity of the particle
904 ThreeVector getPropagationVelocity() const { return (*thePropagationMomentum)/(*thePropagationEnergy); }
905
906 /** \brief Freeze particle propagation
907 *
908 * Make the particle use theFrozenMomentum and theFrozenEnergy for
909 * propagation. The normal state can be restored by calling the
910 * thawPropagation() method.
911 */
915 }
916
917 /** \brief Unfreeze particle propagation
918 *
919 * Make the particle use theMomentum and theEnergy for propagation. Call
920 * this method to restore the normal propagation if the
921 * freezePropagation() method has been called.
922 */
926 }
927
928 /** \brief Rotate the particle position and momentum
929 *
930 * \param angle the rotation angle
931 * \param axis a unit vector representing the rotation axis
932 */
933 virtual void rotatePositionAndMomentum(const G4double angle, const ThreeVector &axis) {
934 rotatePosition(angle, axis);
935 rotateMomentum(angle, axis);
936 }
937
938 /** \brief Rotate the particle position
939 *
940 * \param angle the rotation angle
941 * \param axis a unit vector representing the rotation axis
942 */
943 virtual void rotatePosition(const G4double angle, const ThreeVector &axis) {
944 thePosition.rotate(angle, axis);
945 }
946
947 /** \brief Rotate the particle momentum
948 *
949 * \param angle the rotation angle
950 * \param axis a unit vector representing the rotation axis
951 */
952 virtual void rotateMomentum(const G4double angle, const ThreeVector &axis) {
953 theMomentum.rotate(angle, axis);
954 theFrozenMomentum.rotate(angle, axis);
955 }
956
957 std::string print() const {
958 std::stringstream ss;
959 ss << "Particle (ID = " << ID << ") type = ";
961 ss << '\n'
962 << " energy = " << theEnergy << '\n'
963 << " momentum = "
964 << theMomentum.print()
965 << '\n'
966 << " position = "
967 << thePosition.print()
968 << '\n';
969 return ss.str();
970 };
971
972 std::string dump() const {
973 std::stringstream ss;
974 ss << "(particle " << ID << " ";
976 ss << '\n'
977 << thePosition.dump()
978 << '\n'
979 << theMomentum.dump()
980 << '\n'
981 << theEnergy << ")" << '\n';
982 return ss.str();
983 };
984
985 long getID() const { return ID; };
986
987 /**
988 * Return a NULL pointer
989 */
990 ParticleList const *getParticles() const {
991 INCL_WARN("Particle::getParticles() method was called on a Particle object" << '\n');
992 return 0;
993 }
994
995 /** \brief Return the reflection momentum
996 *
997 * The reflection momentum is used by calls to getSurfaceRadius to compute
998 * the radius of the sphere where the nucleon moves. It is necessary to
999 * introduce fuzzy r-p correlations.
1000 */
1002 if(rpCorrelated)
1003 return theMomentum.mag();
1004 else
1005 return uncorrelatedMomentum;
1006 }
1007
1008 /// \brief Set the uncorrelated momentum
1010
1011 /// \brief Make the particle follow a strict r-p correlation
1012 void rpCorrelate() { rpCorrelated = true; }
1013
1014 /// \brief Make the particle not follow a strict r-p correlation
1015 void rpDecorrelate() { rpCorrelated = false; }
1016
1017 /// \brief Get the cosine of the angle between position and momentum
1020 if(norm>0.)
1021 return thePosition.dot(*thePropagationMomentum) / std::sqrt(norm);
1022 else
1023 return 1.;
1024 }
1025
1026 /// \brief General bias vector function
1027 static G4double getTotalBias();
1028 static void setINCLBiasVector(std::vector<G4double> NewVector);
1029 static void FillINCLBiasVector(G4double newBias);
1030 static G4double getBiasFromVector(std::vector<G4int> VectorBias);
1031
1032 static std::vector<G4int> MergeVectorBias(Particle const * const p1, Particle const * const p2);
1033 static std::vector<G4int> MergeVectorBias(std::vector<G4int> p1, Particle const * const p2);
1034
1035 /// \brief Get the particle bias.
1037
1038 /// \brief Set the particle bias.
1039 void setParticleBias(G4double ParticleBias) { this->theParticleBias = ParticleBias; }
1040
1041 /// \brief Get the vector list of biased vertices on the particle path.
1042 std::vector<G4int> getBiasCollisionVector() const { return theBiasCollisionVector; }
1043
1044 /// \brief Set the vector list of biased vertices on the particle path.
1045 void setBiasCollisionVector(std::vector<G4int> BiasCollisionVector) {
1046 this->theBiasCollisionVector = BiasCollisionVector;
1047 this->setParticleBias(Particle::getBiasFromVector(BiasCollisionVector));
1048 }
1049
1050 /** \brief Number of Kaon inside de nucleus
1051 *
1052 * Put in the Particle class in order to calculate the
1053 * "correct" mass of composit particle.
1054 *
1055 */
1056
1057 G4int getNumberOfKaon() const { return theNKaon; };
1058 void setNumberOfKaon(const G4int NK) { theNKaon = NK; }
1059
1061 void setParentResonancePDGCode(const G4int parentPDGCode) { theParentResonancePDGCode = parentPDGCode; };
1063 void setParentResonanceID(const G4int parentID) { theParentResonanceID = parentID; };
1064
1065 public:
1066 /** \brief Time ordered vector of all bias applied
1067 *
1068 * /!\ Caution /!\
1069 * methods Assotiated to G4VectorCache<T> are:
1070 * Push_back(…),
1071 * operator[],
1072 * Begin(),
1073 * End(),
1074 * Clear(),
1075 * Size() and
1076 * Pop_back()
1077 *
1078 */
1079#ifdef INCLXX_IN_GEANT4_MODE
1080 static std::vector<G4double> INCLBiasVector;
1081 //static G4VectorCache<G4double> INCLBiasVector;
1082#else
1083 static G4ThreadLocal std::vector<G4double> INCLBiasVector;
1084 //static G4VectorCache<G4double> INCLBiasVector;
1085#endif
1087
1088 protected:
1102 long ID;
1103
1106
1108 /// \brief The number of Kaons inside the nucleus (update during the cascade)
1110
1113
1114 private:
1115 G4double theHelicity;
1116 G4double emissionTime;
1117 G4bool outOfWell;
1118
1119 /// \brief Time ordered vector of all biased vertices on the particle path
1120 std::vector<G4int> theBiasCollisionVector;
1121
1122 G4double theMass;
1123 static G4ThreadLocal long nextID;
1124
1126 };
1127}
1128
1129#endif /* PARTICLE_HH_ */
Singleton for recycling allocation of instances of a given class.
#define INCL_DECLARE_ALLOCATION_POOL(T)
#define INCL_ERROR(x)
#define INCL_WARN(x)
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4double getParticleListBias() const
std::vector< G4int > getParticleListBiasVector() const
void rotateMomentum(const G4double angle, const ThreeVector &axis) const
void boost(const ThreeVector &b) const
void rotatePosition(const G4double angle, const ThreeVector &axis) const
void rotatePositionAndMomentum(const G4double angle, const ThreeVector &axis) const
G4bool isCluster() const
ThreeVector boostVector() const
G4double getINCLMass() const
Get the INCL particle mass.
G4bool isBaryon() const
Is this a Baryon?
virtual G4INCL::ParticleSpecies getSpecies() const
Get the particle species.
void setPotentialEnergy(G4double v)
Set the particle potential energy.
ParticleList const * getParticles() const
G4bool isPhoton() const
Is this a photon?
G4INCL::ThreeVector * thePropagationMomentum
G4bool isLambda() const
Is this a Lambda?
G4int getS() const
Returns the strangeness number.
G4bool isOutOfWell() const
Check if the particle is out of its potential well.
G4double getEmissionQValueCorrection(const G4int AParent, const G4int ZParent) const
Computes correction on the emission Q-value.
G4INCL::ThreeVector theMomentum
void setBiasCollisionVector(std::vector< G4int > BiasCollisionVector)
Set the vector list of biased vertices on the particle path.
void setNumberOfKaon(const G4int NK)
void setParentResonancePDGCode(const G4int parentPDGCode)
G4bool isMeson() const
Is this a Meson?
virtual G4INCL::ThreeVector getAngularMomentum() const
void setFrozenMomentum(const ThreeVector &momentum)
Set the frozen particle momentum.
G4double getFrozenEnergy() const
Get the frozen particle momentum.
void setUncorrelatedMomentum(const G4double p)
Set the uncorrelated momentum.
G4double getEnergy() const
G4double getPotentialEnergy() const
Get the particle potential energy.
G4bool isStrange() const
Is this an Strange?
void incrementNumberOfDecays()
Increment the number of decays undergone by the particle.
G4bool isEtaPrime() const
Is this an etaprime?
G4double getHelicity()
static std::vector< G4double > INCLBiasVector
Time ordered vector of all bias applied.
void rpCorrelate()
Make the particle follow a strict r-p correlation.
void setMass(G4double mass)
G4bool isNucleonorLambda() const
Is this a Nucleon or a Lambda?
G4bool isOmega() const
Is this an omega?
void lorentzContract(const ThreeVector &aBoostVector, const ThreeVector &refPos)
Lorentz-contract the particle position around some center.
virtual void makeTargetSpectator()
ParticipantType getParticipantType() const
void setHelicity(G4double h)
void setFrozenEnergy(const G4double energy)
Set the frozen particle momentum.
ThreeVector getPropagationVelocity() const
Get the propagation velocity of the particle.
G4double thePotentialEnergy
virtual void rotatePositionAndMomentum(const G4double angle, const ThreeVector &axis)
Rotate the particle position and momentum.
void setParticipantType(ParticipantType const p)
ParticipantType theParticipantType
std::string dump() const
void propagate(G4double step)
G4int getZ() const
Returns the charge number.
static void FillINCLBiasVector(G4double newBias)
G4bool isSigma() const
Is this a Sigma?
const G4INCL::ThreeVector & getPosition() const
void setParticleBias(G4double ParticleBias)
Set the particle bias.
G4double getBeta() const
G4bool isParticipant() const
void setParentResonanceID(const G4int parentID)
G4double getKineticEnergy() const
Get the particle kinetic energy.
G4INCL::ParticleType theType
Particle(const Particle &rhs)
Copy constructor.
G4double getReflectionMomentum() const
Return the reflection momentum.
G4bool isTargetSpectator() const
static std::vector< G4int > MergeVectorBias(Particle const *const p1, Particle const *const p2)
virtual void rotateMomentum(const G4double angle, const ThreeVector &axis)
Rotate the particle momentum.
G4int theNKaon
The number of Kaons inside the nucleus (update during the cascade)
G4double getCosRPAngle() const
Get the cosine of the angle between position and momentum.
std::vector< G4int > getBiasCollisionVector() const
Get the vector list of biased vertices on the particle path.
void setEmissionTime(G4double t)
ThreeVector getLongitudinalPosition() const
Longitudinal component of the position w.r.t. the momentum.
G4int getNumberOfCollisions() const
Return the number of collisions undergone by the particle.
G4bool isEta() const
Is this an eta?
virtual void makeProjectileSpectator()
G4double adjustEnergyFromMomentum()
Recompute the energy to match the momentum.
static void setINCLBiasVector(std::vector< G4double > NewVector)
static G4double getTotalBias()
General bias vector function.
virtual void rotatePosition(const G4double angle, const ThreeVector &axis)
Rotate the particle position.
G4double getParticleBias() const
Get the particle bias.
void setNumberOfDecays(G4int n)
Set the number of decays undergone by the particle.
G4int getNumberOfDecays() const
Return the number of decays undergone by the particle.
const ThreeVector & adjustMomentumFromEnergy()
Rescale the momentum to match the total energy.
G4bool isPion() const
Is this a pion?
void setINCLMass()
Set the mass of the Particle to its table mass.
void incrementNumberOfCollisions()
Increment the number of collisions 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.
ThreeVector getFrozenMomentum() const
Get the frozen particle momentum.
G4double getRealMass() const
Get the real particle mass.
G4double uncorrelatedMomentum
G4bool isProjectileSpectator() const
G4bool isAntiKaon() const
Is this an antiKaon?
const G4INCL::ThreeVector & getMomentum() const
void rpDecorrelate()
Make the particle not follow a strict r-p correlation.
Particle & operator=(const Particle &rhs)
Assignment operator.
virtual void setMomentum(const G4INCL::ThreeVector &momentum)
G4INCL::ParticleType getType() const
G4bool isKaon() const
Is this a Kaon?
G4int getNumberOfKaon() const
Number of Kaon inside de nucleus.
void thawPropagation()
Unfreeze particle propagation.
virtual G4double getTableMass() const
Get the tabulated particle mass.
G4double getEmissionQValueCorrection(const G4int AParent, const G4int ZParent, const G4int SParent) const
Computes correction on the emission Q-value for hypernuclei.
G4double getInvariantMass() const
Get the the particle invariant mass.
G4double getEmissionTime()
virtual void makeParticipant()
G4bool isResonance() const
Is it a resonance?
void setEnergy(G4double energy)
void setType(ParticleType t)
std::string print() const
ThreeVector getTransversePosition() const
Transverse component of the position w.r.t. the momentum.
G4double getMass() const
Get the cached particle mass.
void boost(const ThreeVector &aBoostVector)
static G4ThreadLocal G4int nextBiasedCollisionID
virtual void setPosition(const G4INCL::ThreeVector &position)
void setTableMass()
Set the mass of the Particle to its table mass.
void setOutOfWell()
Mark the particle as out of its potential well.
void swap(Particle &rhs)
Helper method for the assignment operator.
long getID() const
G4double * thePropagationEnergy
G4bool isDelta() const
Is it a Delta?
void freezePropagation()
Freeze particle propagation.
G4int getParentResonancePDGCode() const
void setRealMass()
Set the mass of the Particle to its real mass.
static G4double getBiasFromVector(std::vector< G4int > VectorBias)
G4int getA() const
Returns the baryon number.
G4bool isHyperon() const
Is this an Hyperon?
G4int getParentResonanceID() const
G4INCL::ThreeVector thePosition
G4INCL::ThreeVector theFrozenMomentum
void setNumberOfCollisions(G4int n)
Set the number of collisions undergone by the particle.
G4bool isNucleon() const
G4double getTransferQValueCorrection(const G4int AFrom, const G4int ZFrom, const G4int SFrom, const G4int ATo, const G4int ZTo, const G4int STo) const
Computes correction on the transfer Q-value for hypernuclei.
std::string dump() const
G4double mag() const
std::string print() const
void rotate(const G4double angle, const ThreeVector &axis)
Rotate the vector by a given angle around a given axis.
G4double dot(const ThreeVector &v) const
G4double mag2() const
ThreeVector vector(const ThreeVector &v) const
G4ThreadLocal NuclearMassFn getTableMass
Static pointer to the mass function for nuclei.
G4double getTableQValue(const G4int A1, const G4int Z1, const G4int S1, const G4int A2, const G4int Z2, const G4int S2)
Get Q-value (in MeV/c^2)
G4ThreadLocal ParticleMassFn getTableParticleMass
Static pointer to the mass function for particles.
G4double getRealMass(const G4INCL::ParticleType t)
Get particle mass (in MeV/c^2)
G4double getINCLMass(const G4int A, const G4int Z, const G4int S)
Get INCL nuclear mass (in MeV/c^2)
std::string getName(const ParticleType t)
Get the native INCL name of the particle.
ParticleList::iterator ParticleMutableIter
ParticleList::const_iterator ParticleIter
#define G4ThreadLocal
Definition: tls.hh:77