Geant4 11.4.0
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),
99 nSrcPair(rhs.nSrcPair),
104 theNKaon(rhs.theNKaon),
108#endif
109 theHelicity(rhs.theHelicity),
110 emissionTime(rhs.emissionTime),
111 outOfWell(rhs.outOfWell),
112 theSrcPartner(rhs.theSrcPartner),
113 theMass(rhs.theMass)
114 {
115 if(rhs.thePropagationEnergy == &(rhs.theFrozenEnergy))
117 else
121 else
123 // ID intentionally not copied
124 ID = nextID++;
125
126 theBiasCollisionVector = rhs.theBiasCollisionVector;
127 }
128
129 protected:
130 /// \brief Helper method for the assignment operator
131 void swap(Particle &rhs) {
132 std::swap(theZ, rhs.theZ);
133 std::swap(theA, rhs.theA);
134 std::swap(theS, rhs.theS);
136 std::swap(theType, rhs.theType);
137 if(rhs.thePropagationEnergy == &(rhs.theFrozenEnergy))
139 else
141 std::swap(theEnergy, rhs.theEnergy);
142 std::swap(theFrozenEnergy, rhs.theFrozenEnergy);
145 else
147 std::swap(theMomentum, rhs.theMomentum);
148 std::swap(theFrozenMomentum, rhs.theFrozenMomentum);
149 std::swap(thePosition, rhs.thePosition);
150 std::swap(nCollisions, rhs.nCollisions);
151 std::swap(nDecays, rhs.nDecays);
152 std::swap(nSrcPair, rhs.nSrcPair),
154 // ID intentionally not swapped
155
156#ifdef INCLXX_IN_GEANT4_MODE
159#endif
160
161 std::swap(theHelicity, rhs.theHelicity);
162 std::swap(emissionTime, rhs.emissionTime);
163 std::swap(outOfWell, rhs.outOfWell);
164 std::swap(theSrcPartner, rhs.theSrcPartner);
165
166 std::swap(theMass, rhs.theMass);
167 std::swap(rpCorrelated, rhs.rpCorrelated);
169
170 std::swap(theParticleBias, rhs.theParticleBias);
171 std::swap(theBiasCollisionVector, rhs.theBiasCollisionVector);
172
173 }
174
175 public:
176
177 /** \brief Assignment operator
178 *
179 * Does not copy the particle ID.
180 */
182 Particle temporaryParticle(rhs);
183 swap(temporaryParticle);
184 return *this;
185 }
186
187 /**
188 * Get the particle type.
189 * @see G4INCL::ParticleType
190 */
192 return theType;
193 };
194
195 /// \brief Get the particle species
197 return ParticleSpecies(theType);
198 };
199
201 theType = t;
202 switch(theType)
203 {
204 case DeltaPlusPlus:
205 theA = 1;
206 theZ = 2;
207 theS = 0;
208 break;
209 case Proton:
210 case DeltaPlus:
211 theA = 1;
212 theZ = 1;
213 theS = 0;
214 break;
215 case Neutron:
216 case DeltaZero:
217 theA = 1;
218 theZ = 0;
219 theS = 0;
220 break;
221 case DeltaMinus:
222 theA = 1;
223 theZ = -1;
224 theS = 0;
225 break;
226 case PiPlus:
227 theA = 0;
228 theZ = 1;
229 theS = 0;
230 break;
231 case PiZero:
232 case Eta:
233 case Omega:
234 case EtaPrime:
235 case Photon:
236 theA = 0;
237 theZ = 0;
238 theS = 0;
239 break;
240 case PiMinus:
241 theA = 0;
242 theZ = -1;
243 theS = 0;
244 break;
245 case Lambda:
246 theA = 1;
247 theZ = 0;
248 theS = -1;
249 break;
250 case SigmaPlus:
251 theA = 1;
252 theZ = 1;
253 theS = -1;
254 break;
255 case SigmaZero:
256 theA = 1;
257 theZ = 0;
258 theS = -1;
259 break;
260 case SigmaMinus:
261 theA = 1;
262 theZ = -1;
263 theS = -1;
264 break;
265 case antiProton:
266 theA = -1;
267 theZ = -1;
268 theS = 0;
269 break;
270 case XiMinus:
271 theA = 1;
272 theZ = -1;
273 theS = -2;
274 break;
275 case XiZero:
276 theA = 1;
277 theZ = 0;
278 theS = -2;
279 break;
280 case antiNeutron:
281 theA = -1;
282 theZ = 0;
283 theS = 0;
284 break;
285 case antiLambda:
286 theA = -1;
287 theZ = 0;
288 theS = 1;
289 break;
290 case antiSigmaMinus:
291 theA = -1;
292 theZ = 1;
293 theS = 1;
294 break;
295 case antiSigmaPlus:
296 theA = -1;
297 theZ = -1;
298 theS = 1;
299 break;
300 case antiSigmaZero:
301 theA = -1;
302 theZ = 0;
303 theS = 1;
304 break;
305 case antiXiMinus:
306 theA = -1;
307 theZ = 1;
308 theS = 2;
309 break;
310 case antiXiZero:
311 theA = -1;
312 theZ = 0;
313 theS = 2;
314 break;
315 case KPlus:
316 theA = 0;
317 theZ = 1;
318 theS = 1;
319 break;
320 case KZero:
321 theA = 0;
322 theZ = 0;
323 theS = 1;
324 break;
325 case KZeroBar:
326 theA = 0;
327 theZ = 0;
328 theS = -1;
329 break;
330 case KShort:
331 theA = 0;
332 theZ = 0;
333// theS should not be defined
334 break;
335 case KLong:
336 theA = 0;
337 theZ = 0;
338// theS should not be defined
339 break;
340 case KMinus:
341 theA = 0;
342 theZ = -1;
343 theS = -1;
344 break;
345 case Composite:
346 // INCL_ERROR("Trying to set particle type to Composite! Construct a Cluster object instead" << '\n');
347 theA = 0;
348 theZ = 0;
349 theS = 0;
350 break;
351 case antiComposite:
352 theA = 0;
353 theZ = 0;
354 theS = 0;
355 break;
356 case UnknownParticle:
357 theA = 0;
358 theZ = 0;
359 theS = 0;
360 INCL_ERROR("Trying to set particle type to Unknown!" << '\n');
361 break;
362 }
363
364 if( !isResonance() && t!=Composite && t!=antiComposite )
365 setINCLMass();
366 }
367
368 /**
369 * Is this a nucleon?
370 */
373 return true;
374 else
375 return false;
376 };
377
381
385
388 }
389
393
397
398 virtual void makeParticipant() {
400 }
401
405
409
410 /** \brief Is this a pion? */
411 G4bool isPion() const { return (theType == PiPlus || theType == PiZero || theType == PiMinus); }
412
413 /** \brief Is this an eta? */
414 G4bool isEta() const { return (theType == Eta); }
415
416 /** \brief Is this an omega? */
417 G4bool isOmega() const { return (theType == Omega); }
418
419 /** \brief Is this an etaprime? */
420 G4bool isEtaPrime() const { return (theType == EtaPrime); }
421
422 /** \brief Is this a photon? */
423 G4bool isPhoton() const { return (theType == Photon); }
424
425 /** \brief Is it a resonance? */
426 inline G4bool isResonance() const { return isDelta(); }
427
428 /** \brief Is it a Delta? */
429 inline G4bool isDelta() const {
430 return (theType==DeltaPlusPlus || theType==DeltaPlus ||
432
433 /** \brief Is this a Sigma? */
434 G4bool isSigma() const { return (theType == SigmaPlus || theType == SigmaZero || theType == SigmaMinus); }
435
436 /** \brief Is this a Kaon? */
437 G4bool isKaon() const { return (theType == KPlus || theType == KZero); }
438
439 /** \brief Is this an antiKaon? */
440 G4bool isAntiKaon() const { return (theType == KZeroBar || theType == KMinus); }
441
442 /** \brief Is this a Lambda? */
443 G4bool isLambda() const { return (theType == Lambda); }
444
445 /** \brief Is this a Nucleon or a Lambda? */
446 G4bool isNucleonorLambda() const { return (isNucleon() || isLambda()); }
447
448 /** \brief Is this an Hyperon? */
449 G4bool isHyperon() const { return (isLambda() || isSigma() ); } //|| isXi()
450
451 /** \brief Is this a Meson? */
452 G4bool isMeson() const { return (isPion() || isKaon() || isAntiKaon() || isEta() || isEtaPrime() || isOmega()); }
453
454 /** \brief Is this a Baryon? */
455 G4bool isBaryon() const { return (isNucleon() || isResonance() || isHyperon()); }
456
457 /** \brief Is this a Strange? */
458 G4bool isStrange() const { return (isKaon() || isAntiKaon() || isHyperon()); }
459
460 /** \brief Is this a Xi? */
461 G4bool isXi() const { return (theType == XiZero || theType == XiMinus); }
462
463 /** \brief Is this an antinucleon? */
465
466 /** \brief Is this an antiSigma? */
468
469 /** \brief Is this an antiXi? */
470 G4bool isAntiXi() const { return (theType == antiXiZero || theType == antiXiMinus); }
471
472 /** \brief Is this an antiLambda? */
473 G4bool isAntiLambda() const { return (theType == antiLambda); }
474
475 /** \brief Is this an antiHyperon? */
476 G4bool isAntiHyperon() const { return (isAntiLambda() || isAntiSigma() || isAntiXi()); }
477
478 /** \brief Is this an antiBaryon? */
479 G4bool isAntiBaryon() const { return (isAntiNucleon() || isAntiHyperon()); }
480
481 /** \brief Is this an antiNucleon or an antiLambda? */
483
484 /** \brief Returns the baryon number. */
485 G4int getA() const { return theA; }
486
487 /** \brief Returns the charge number. */
488 G4int getZ() const { return theZ; }
489
490 /** \brief Returns the strangeness number. */
491 G4int getS() const { return theS; }
492
493 /** \brief Returns the strangeness number. */
494 G4int getSrcPair() const { return nSrcPair; }
495
497 const G4double P = theMomentum.mag();
498 return P/theEnergy;
499 }
500
501 /**
502 * Returns a three vector we can give to the boost() -method.
503 *
504 * In order to go to the particle rest frame you need to multiply
505 * the boost vector by -1.0.
506 */
508 return theMomentum / theEnergy;
509 }
510
511 /**
512 * Boost the particle using a boost vector.
513 *
514 * Example (go to the particle rest frame):
515 * particle->boost(particle->boostVector());
516 */
517 void boost(const ThreeVector &aBoostVector) {
518 const G4double beta2 = aBoostVector.mag2();
519 const G4double gamma = 1.0 / std::sqrt(1.0 - beta2);
520 const G4double bp = theMomentum.dot(aBoostVector);
521 const G4double alpha = (gamma*gamma)/(1.0 + gamma);
522
523 theMomentum = theMomentum + aBoostVector * (alpha * bp - gamma * theEnergy);
524 theEnergy = gamma * (theEnergy - bp);
525 }
526
527 /** \brief Lorentz-contract the particle position around some center
528 *
529 * Apply Lorentz contraction to the position component along the
530 * direction of the boost vector.
531 *
532 * \param aBoostVector the boost vector (velocity) [c]
533 * \param refPos the reference position
534 */
535 void lorentzContract(const ThreeVector &aBoostVector, const ThreeVector &refPos) {
536 const G4double beta2 = aBoostVector.mag2();
537 const G4double gamma = 1.0 / std::sqrt(1.0 - beta2);
538 const ThreeVector theRelativePosition = thePosition - refPos;
539 const ThreeVector transversePosition = theRelativePosition - aBoostVector * (theRelativePosition.dot(aBoostVector) / aBoostVector.mag2());
540 const ThreeVector longitudinalPosition = theRelativePosition - transversePosition;
541
542 thePosition = refPos + transversePosition + longitudinalPosition / gamma;
543 }
544
545 /** \brief Get the cached particle mass. */
546 inline G4double getMass() const { return theMass; }
547
548 /** \brief Get the INCL particle mass. */
549 inline G4double getINCLMass() const {
550 switch(theType) {
551 case Proton:
552 case Neutron:
553 case PiPlus:
554 case PiMinus:
555 case PiZero:
556 case Lambda:
557 case SigmaPlus:
558 case SigmaZero:
559 case SigmaMinus:
560 case antiProton:
561 case XiZero:
562 case XiMinus:
563 case antiNeutron:
564 case antiLambda:
565 case antiSigmaPlus:
566 case antiSigmaZero:
567 case antiSigmaMinus:
568 case antiXiZero:
569 case antiXiMinus:
570 case KPlus:
571 case KZero:
572 case KZeroBar:
573 case KShort:
574 case KLong:
575 case KMinus:
576 case Eta:
577 case Omega:
578 case EtaPrime:
579 case Photon:
581 break;
582
583 case DeltaPlusPlus:
584 case DeltaPlus:
585 case DeltaZero:
586 case DeltaMinus:
587 return theMass;
588 break;
589
590 case Composite:
592 break;
593 case antiComposite:
595 break;
596
597 default:
598 INCL_ERROR("Particle::getINCLMass: Unknown particle type." << '\n');
599 return 0.0;
600 break;
601 }
602 }
603
604 /** \brief Get the tabulated particle mass. */
605 inline virtual G4double getTableMass() const {
606 switch(theType) {
607 case Proton:
608 case Neutron:
609 case PiPlus:
610 case PiMinus:
611 case PiZero:
612 case Lambda:
613 case SigmaPlus:
614 case SigmaZero:
615 case SigmaMinus:
616 case antiProton:
617 case XiZero:
618 case XiMinus:
619 case antiNeutron:
620 case antiLambda:
621 case antiSigmaPlus:
622 case antiSigmaZero:
623 case antiSigmaMinus:
624 case antiXiZero:
625 case antiXiMinus:
626 case KPlus:
627 case KZero:
628 case KZeroBar:
629 case KShort:
630 case KLong:
631 case KMinus:
632 case Eta:
633 case Omega:
634 case EtaPrime:
635 case Photon:
637 break;
638
639 case DeltaPlusPlus:
640 case DeltaPlus:
641 case DeltaZero:
642 case DeltaMinus:
643 return theMass;
644 break;
645
646 case Composite:
648 break;
649 case antiComposite:
651 break;
652
653 default:
654 INCL_ERROR("Particle::getTableMass: Unknown particle type." << '\n');
655 return 0.0;
656 break;
657 }
658 }
659
660 /** \brief Get the real particle mass. */
661 inline G4double getRealMass() const {
662 switch(theType) {
663 case Proton:
664 case Neutron:
665 case PiPlus:
666 case PiMinus:
667 case PiZero:
668 case Lambda:
669 case SigmaPlus:
670 case SigmaZero:
671 case SigmaMinus:
672 case antiProton:
673 case XiZero:
674 case XiMinus:
675 case antiNeutron:
676 case antiLambda:
677 case antiSigmaPlus:
678 case antiSigmaZero:
679 case antiSigmaMinus:
680 case antiXiZero:
681 case antiXiMinus:
682 case KPlus:
683 case KZero:
684 case KZeroBar:
685 case KShort:
686 case KLong:
687 case KMinus:
688 case Eta:
689 case Omega:
690 case EtaPrime:
691 case Photon:
693 break;
694
695 case DeltaPlusPlus:
696 case DeltaPlus:
697 case DeltaZero:
698 case DeltaMinus:
699 return theMass;
700 break;
701
702 case Composite:
704 break;
705 case antiComposite:
707 break;
708
709 default:
710 INCL_ERROR("Particle::getRealMass: Unknown particle type." << '\n');
711 return 0.0;
712 break;
713 }
714 }
715
716 /// \brief Set the mass of the Particle to its real mass
718
719 /// \brief Set the mass of the Particle to its table mass
721
722 /// \brief Set the mass of the Particle to its table mass
724
725 /**\brief Computes correction on the emission Q-value
726 *
727 * Computes the correction that must be applied to INCL particles in
728 * order to obtain the correct Q-value for particle emission from a given
729 * nucleus. For absorption, the correction is obviously equal to minus
730 * the value returned by this function.
731 *
732 * \param AParent the mass number of the emitting nucleus
733 * \param ZParent the charge number of the emitting nucleus
734 * \return the correction
735 */
736 G4double getEmissionQValueCorrection(const G4int AParent, const G4int ZParent) const {
737 const G4int SParent = 0;
738 const G4int ADaughter = AParent - theA;
739 const G4int ZDaughter = ZParent - theZ;
740 const G4int SDaughter = 0;
741
742 // Note the minus sign here
743 G4double theQValue;
744 if(isCluster())
745 theQValue = -ParticleTable::getTableQValue(theA, theZ, theS, ADaughter, ZDaughter, SDaughter);
746 else {
747 const G4double massTableParent = ParticleTable::getTableMass(AParent,ZParent,SParent);
748 const G4double massTableDaughter = ParticleTable::getTableMass(ADaughter,ZDaughter,SDaughter);
749 const G4double massTableParticle = getTableMass();
750 theQValue = massTableParent - massTableDaughter - massTableParticle;
751 }
752
753 const G4double massINCLParent = ParticleTable::getINCLMass(AParent,ZParent,SParent);
754 const G4double massINCLDaughter = ParticleTable::getINCLMass(ADaughter,ZDaughter,SDaughter);
755 const G4double massINCLParticle = getINCLMass();
756
757 // The rhs corresponds to the INCL Q-value
758 return theQValue - (massINCLParent-massINCLDaughter-massINCLParticle);
759 }
760
761 /**\brief Computes correction on the transfer Q-value
762 *
763 * Computes the correction that must be applied to INCL particles in
764 * order to obtain the correct Q-value for particle transfer from a given
765 * nucleus to another.
766 *
767 * Assumes that the receving nucleus is INCL's target nucleus, with the
768 * INCL separation energy.
769 *
770 * \param AFrom the mass number of the donating nucleus
771 * \param ZFrom the charge number of the donating nucleus
772 * \param ATo the mass number of the receiving nucleus
773 * \param ZTo the charge number of the receiving nucleus
774 * \return the correction
775 */
776 G4double getTransferQValueCorrection(const G4int AFrom, const G4int ZFrom, const G4int ATo, const G4int ZTo) const {
777 const G4int SFrom = 0;
778 const G4int STo = 0;
779 const G4int AFromDaughter = AFrom - theA;
780 const G4int ZFromDaughter = ZFrom - theZ;
781 const G4int SFromDaughter = 0;
782 const G4int AToDaughter = ATo + theA;
783 const G4int ZToDaughter = ZTo + theZ;
784 const G4int SToDaughter = 0;
785 const G4double theQValue = ParticleTable::getTableQValue(AToDaughter,ZToDaughter,SToDaughter,AFromDaughter,ZFromDaughter,SFromDaughter,AFrom,ZFrom,SFrom);
786
787 const G4double massINCLTo = ParticleTable::getINCLMass(ATo,ZTo,STo);
788 const G4double massINCLToDaughter = ParticleTable::getINCLMass(AToDaughter,ZToDaughter,SToDaughter);
789 /* Note that here we have to use the table mass in the INCL Q-value. We
790 * cannot use theMass, because at this stage the particle is probably
791 * still off-shell; and we cannot use getINCLMass(), because it leads to
792 * violations of global energy conservation.
793 */
794 const G4double massINCLParticle = getTableMass();
795
796 // The rhs corresponds to the INCL Q-value for particle absorption
797 return theQValue - (massINCLToDaughter-massINCLTo-massINCLParticle);
798 }
799
800 /**\brief Computes correction on the emission Q-value for hypernuclei
801 *
802 * Computes the correction that must be applied to INCL particles in
803 * order to obtain the correct Q-value for particle emission from a given
804 * nucleus. For absorption, the correction is obviously equal to minus
805 * the value returned by this function.
806 *
807 * \param AParent the mass number of the emitting nucleus
808 * \param ZParent the charge number of the emitting nucleus
809 * \param SParent the strangess number of the emitting nucleus
810 * \return the correction
811 */
812 G4double getEmissionQValueCorrection(const G4int AParent, const G4int ZParent, const G4int SParent) const {
813 const G4int ADaughter = AParent - theA;
814 const G4int ZDaughter = ZParent - theZ;
815 const G4int SDaughter = SParent - theS;
816
817 // Note the minus sign here
818 G4double theQValue;
819 if(isCluster())
820 theQValue = -ParticleTable::getTableQValue(theA, theZ, theS, ADaughter, ZDaughter, SDaughter);
821 else {
822 const G4double massTableParent = ParticleTable::getTableMass(AParent,ZParent,SParent);
823 const G4double massTableDaughter = ParticleTable::getTableMass(ADaughter,ZDaughter,SDaughter);
824 const G4double massTableParticle = getTableMass();
825 theQValue = massTableParent - massTableDaughter - massTableParticle;
826 }
827
828 const G4double massINCLParent = ParticleTable::getINCLMass(AParent,ZParent,SParent);
829 const G4double massINCLDaughter = ParticleTable::getINCLMass(ADaughter,ZDaughter,SDaughter);
830 const G4double massINCLParticle = getINCLMass();
831
832 // The rhs corresponds to the INCL Q-value
833 return theQValue - (massINCLParent-massINCLDaughter-massINCLParticle);
834 }
835
836 /**\brief Computes correction on the transfer Q-value for hypernuclei
837 *
838 * Computes the correction that must be applied to INCL particles in
839 * order to obtain the correct Q-value for particle transfer from a given
840 * nucleus to another.
841 *
842 * Assumes that the receving nucleus is INCL's target nucleus, with the
843 * INCL separation energy.
844 *
845 * \param AFrom the mass number of the donating nucleus
846 * \param ZFrom the charge number of the donating nucleus
847 * \param SFrom the strangess number of the donating nucleus
848 * \param ATo the mass number of the receiving nucleus
849 * \param ZTo the charge number of the receiving nucleus
850 * \param STo the strangess number of the receiving nucleus
851 * \return the correction
852 */
853 G4double getTransferQValueCorrection(const G4int AFrom, const G4int ZFrom, const G4int SFrom, const G4int ATo, const G4int ZTo , const G4int STo) const {
854 const G4int AFromDaughter = AFrom - theA;
855 const G4int ZFromDaughter = ZFrom - theZ;
856 const G4int SFromDaughter = SFrom - theS;
857 const G4int AToDaughter = ATo + theA;
858 const G4int ZToDaughter = ZTo + theZ;
859 const G4int SToDaughter = STo + theS;
860 const G4double theQValue = ParticleTable::getTableQValue(AToDaughter,ZToDaughter,SFromDaughter,AFromDaughter,ZFromDaughter,SToDaughter,AFrom,ZFrom,SFrom);
861
862 const G4double massINCLTo = ParticleTable::getINCLMass(ATo,ZTo,STo);
863 const G4double massINCLToDaughter = ParticleTable::getINCLMass(AToDaughter,ZToDaughter,SToDaughter);
864 /* Note that here we have to use the table mass in the INCL Q-value. We
865 * cannot use theMass, because at this stage the particle is probably
866 * still off-shell; and we cannot use getINCLMass(), because it leads to
867 * violations of global energy conservation.
868 */
869 const G4double massINCLParticle = getTableMass();
870
871 // The rhs corresponds to the INCL Q-value for particle absorption
872 return theQValue - (massINCLToDaughter-massINCLTo-massINCLParticle);
873 }
874
875
876
877 /** \brief Get the the particle invariant mass.
878 *
879 * Uses the relativistic invariant
880 * \f[ m = \sqrt{E^2 - {\vec p}^2}\f]
881 **/
883 const G4double mass = std::pow(theEnergy, 2) - theMomentum.dot(theMomentum);
884 if(mass < 0.0) {
885 INCL_ERROR("E*E - p*p is negative." << '\n');
886 return 0.0;
887 } else {
888 return std::sqrt(mass);
889 }
890 };
891
892 /// \brief Get the particle kinetic energy.
893 inline G4double getKineticEnergy() const { return theEnergy - theMass; }
894
895 /// \brief Get the particle potential energy.
897
898 /// \brief Set the particle potential energy.
900
901 /**
902 * Get the energy of the particle in MeV.
903 */
905 {
906 return theEnergy;
907 };
908
909 /**
910 * Set the mass of the particle in MeV/c^2.
911 */
912 void setMass(G4double mass)
913 {
914 this->theMass = mass;
915 }
916
917 /**
918 * Set the energy of the particle in MeV.
919 */
920 void setEnergy(G4double energy)
921 {
922 this->theEnergy = energy;
923 };
924
925 /**
926 * Get the momentum vector.
927 */
929 {
930 return theMomentum;
931 };
932
933 /** Get the angular momentum w.r.t. the origin */
935 {
936 return thePosition.vector(theMomentum);
937 };
938
939 /**
940 * Set the momentum vector.
941 */
942 virtual void setMomentum(const G4INCL::ThreeVector &momentum)
943 {
944 this->theMomentum = momentum;
945 };
946
947 /**
948 * Set the position vector.
949 */
951 {
952 return thePosition;
953 };
954
956 {
957 this->thePosition = position;
958 };
959
960 G4double getHelicity() { return theHelicity; };
961 void setHelicity(G4double h) { theHelicity = h; };
962
963 void propagate(G4double step) {
964 thePosition += ((*thePropagationMomentum)*(step/(*thePropagationEnergy)));
965 };
966
967 /** \brief Return the number of collisions undergone by the particle. **/
969
970 /** \brief Set the number of collisions undergone by the particle. **/
972
973 /** \brief Increment the number of collisions undergone by the particle. **/
975
976 /** \brief Return the number of decays undergone by the particle. **/
977 G4int getNumberOfDecays() const { return nDecays; }
978
979 /** \brief Set the number of decays undergone by the particle. **/
981
982 /** \brief Increment the number of decays undergone by the particle. **/
984
985 /** \brief Set the number of srcpairs. **/
986 void setNumberOfSrcPair(int n) { nSrcPair = n; }
987
988 /** \brief Mark the particle as out of its potential well
989 *
990 * This flag is used to control pions created outside their potential well
991 * in delta decay. The pion potential checks it and returns zero if it is
992 * true (necessary in order to correctly enforce energy conservation). The
993 * Nucleus::applyFinalState() method uses it to determine whether new
994 * avatars should be generated for the particle.
995 */
996 void setOutOfWell() { outOfWell = true; }
997
998 /// \brief Check if the particle is out of its potential well
999 G4bool isOutOfWell() const { return outOfWell; }
1000
1001 /// \brief Set and reset src partner
1002 void setSrcPartner() { theSrcPartner = true; }
1003 void resetSrcPartner() { theSrcPartner = false; nSrcPair=0; }
1004
1005 /// \brief Check if the particle is a src partner
1006 G4bool isSrcPartner() const { return theSrcPartner; }
1007
1008 void setEmissionTime(G4double t) { emissionTime = t; }
1009 G4double getEmissionTime() { return emissionTime; };
1010
1011 /** \brief Transverse component of the position w.r.t. the momentum. */
1015
1016 /** \brief Longitudinal component of the position w.r.t. the momentum. */
1020
1021 /** \brief Rescale the momentum to match the total energy. */
1023
1024 /** \brief Recompute the energy to match the momentum. */
1026
1028 return ((theType == Composite || theType == antiComposite));
1029 }
1030
1031 /// \brief Set the frozen particle momentum
1032 void setFrozenMomentum(const ThreeVector &momentum) { theFrozenMomentum = momentum; }
1033
1034 /// \brief Set the frozen particle momentum
1035 void setFrozenEnergy(const G4double energy) { theFrozenEnergy = energy; }
1036
1037 /// \brief Get the frozen particle momentum
1039
1040 /// \brief Get the frozen particle momentum
1042
1043 /// \brief Get the propagation velocity of the particle
1044 ThreeVector getPropagationVelocity() const { return (*thePropagationMomentum)/(*thePropagationEnergy); }
1045
1046 /** \brief Freeze particle propagation
1047 *
1048 * Make the particle use theFrozenMomentum and theFrozenEnergy for
1049 * propagation. The normal state can be restored by calling the
1050 * thawPropagation() method.
1051 */
1056
1057 /** \brief Unfreeze particle propagation
1058 *
1059 * Make the particle use theMomentum and theEnergy for propagation. Call
1060 * this method to restore the normal propagation if the
1061 * freezePropagation() method has been called.
1062 */
1067
1068 /** \brief Rotate the particle position and momentum
1069 *
1070 * \param angle the rotation angle
1071 * \param axis a unit vector representing the rotation axis
1072 */
1073 virtual void rotatePositionAndMomentum(const G4double angle, const ThreeVector &axis) {
1074 rotatePosition(angle, axis);
1075 rotateMomentum(angle, axis);
1076 }
1077
1078 /** \brief Rotate the particle position
1079 *
1080 * \param angle the rotation angle
1081 * \param axis a unit vector representing the rotation axis
1082 */
1083 virtual void rotatePosition(const G4double angle, const ThreeVector &axis) {
1084 thePosition.rotate(angle, axis);
1085 }
1086
1087 /** \brief Rotate the particle momentum
1088 *
1089 * \param angle the rotation angle
1090 * \param axis a unit vector representing the rotation axis
1091 */
1092 virtual void rotateMomentum(const G4double angle, const ThreeVector &axis) {
1093 theMomentum.rotate(angle, axis);
1094 theFrozenMomentum.rotate(angle, axis);
1095 }
1096
1097 std::string print() const {
1098 std::stringstream ss;
1099 ss << "Particle (ID = " << ID << ") type = ";
1101 ss << ", SRC pair = " << nSrcPair;
1102 ss << ", Potential energy = " << thePotentialEnergy;
1103 ss << '\n'
1104 << " energy = " << theEnergy << '\n'
1105 << " momentum = "
1106 << theMomentum.print()
1107 << '\n'
1108 << " position = "
1109 << thePosition.print()
1110 << '\n';
1111 return ss.str();
1112 };
1113
1114 std::string dump() const {
1115 std::stringstream ss;
1116 ss << "(particle " << ID << " ";
1118 ss << nSrcPair << " ";
1119 ss << '\n'
1120 << thePosition.dump()
1121 << '\n'
1122 << theMomentum.dump()
1123 << '\n'
1124 << theEnergy << ")" << '\n';
1125 return ss.str();
1126 };
1127
1128 long getID() const { return ID; };
1129
1130 /**
1131 * Return a NULL pointer
1132 */
1134 INCL_WARN("Particle::getParticles() method was called on a Particle object" << '\n');
1135 return 0;
1136 }
1137
1138 /** \brief Return the reflection momentum
1139 *
1140 * The reflection momentum is used by calls to getSurfaceRadius to compute
1141 * the radius of the sphere where the nucleon moves. It is necessary to
1142 * introduce fuzzy r-p correlations.
1143 */
1145 if(rpCorrelated)
1146 return theMomentum.mag();
1147 else
1148 return uncorrelatedMomentum;
1149 }
1150
1151 /// \brief Set the uncorrelated momentum
1153
1154 /// \brief Make the particle follow a strict r-p correlation
1155 void rpCorrelate() { rpCorrelated = true; }
1156
1157 /// \brief Make the particle not follow a strict r-p correlation
1158 void rpDecorrelate() { rpCorrelated = false; }
1159
1160 /// \brief Get the cosine of the angle between position and momentum
1162 const G4double norm = thePosition.mag2()*thePropagationMomentum->mag2();
1163 if(norm>0.)
1164 return thePosition.dot(*thePropagationMomentum) / std::sqrt(norm);
1165 else
1166 return 1.;
1167 }
1168
1169 /// \brief General bias vector function
1170 static G4double getTotalBias();
1171 static void setINCLBiasVector(std::vector<G4double> NewVector);
1172 static void FillINCLBiasVector(G4double newBias);
1173 static G4double getBiasFromVector(std::vector<G4int> VectorBias);
1174
1175 static std::vector<G4int> MergeVectorBias(Particle const * const p1, Particle const * const p2);
1176 static std::vector<G4int> MergeVectorBias(std::vector<G4int> p1, Particle const * const p2);
1177
1178 /// \brief Get the particle bias.
1180
1181 /// \brief Set the particle bias.
1182 void setParticleBias(G4double ParticleBias) { this->theParticleBias = ParticleBias; }
1183
1184 /// \brief Get the vector list of biased vertices on the particle path.
1185 std::vector<G4int> getBiasCollisionVector() const { return theBiasCollisionVector; }
1186
1187 /// \brief Set the vector list of biased vertices on the particle path.
1188 void setBiasCollisionVector(std::vector<G4int> BiasCollisionVector) {
1189 this->theBiasCollisionVector = BiasCollisionVector;
1190 this->setParticleBias(Particle::getBiasFromVector(std::move(BiasCollisionVector)));
1191 }
1192
1193 /** \brief Number of Kaon inside de nucleus
1194 *
1195 * Put in the Particle class in order to calculate the
1196 * "correct" mass of composit particle.
1197 *
1198 */
1199
1200 G4int getNumberOfKaon() const { return theNKaon; };
1201 void setNumberOfKaon(const G4int NK) { theNKaon = NK; }
1202
1203#ifdef INCLXX_IN_GEANT4_MODE
1205 void setParentResonancePDGCode(const G4int parentPDGCode) { theParentResonancePDGCode = parentPDGCode; };
1207 void setParentResonanceID(const G4int parentID) { theParentResonanceID = parentID; };
1208#endif
1209
1210 public:
1211 /** \brief Time ordered vector of all bias applied
1212 *
1213 * /!\ Caution /!\
1214 * methods Assotiated to G4VectorCache<T> are:
1215 * Push_back(…),
1216 * operator[],
1217 * Begin(),
1218 * End(),
1219 * Clear(),
1220 * Size() and
1221 * Pop_back()
1222 *
1223 */
1224#ifdef INCLXX_IN_GEANT4_MODE
1225 static std::vector<G4double> INCLBiasVector;
1226 //static G4VectorCache<G4double> INCLBiasVector;
1227#else
1228 static G4ThreadLocal std::vector<G4double> INCLBiasVector;
1229 //static G4VectorCache<G4double> INCLBiasVector;
1230#endif
1232
1233 protected:
1248 long ID;
1249
1252
1254 /// \brief The number of Kaons inside the nucleus (update during the cascade)
1256
1257#ifdef INCLXX_IN_GEANT4_MODE
1260#endif
1261
1262 private:
1263 G4double theHelicity;
1264 G4double emissionTime;
1265 G4bool outOfWell;
1266 G4bool theSrcPartner;
1267
1268 /// \brief Time ordered vector of all biased vertices on the particle path
1269 std::vector<G4int> theBiasCollisionVector;
1270
1271 G4double theMass;
1272 static G4ThreadLocal long nextID;
1273
1275 };
1276}
1277
1278#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?
G4bool isAntiBaryon() const
Is this an antiBaryon?
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 isAntiSigma() const
Is this an antiSigma?
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 setNumberOfSrcPair(int n)
Set the number of srcpairs.
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 a Strange?
G4int getSrcPair() const
Returns the strangeness number.
void incrementNumberOfDecays()
Increment the number of decays undergone by the particle.
G4bool isEtaPrime() const
Is this an etaprime?
void setSrcPartner()
Set and reset src partner.
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.
G4bool isAntiHyperon() const
Is this an antiHyperon?
ThreeVector getPropagationVelocity() const
Get the propagation velocity of the particle.
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.
G4bool isAntiXi() const
Is this an antiXi?
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.
G4bool isAntiNucleon() const
Is this an antinucleon?
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?
G4bool isSrcPartner() const
Check if the particle is a src partner.
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.
G4bool isXi() const
Is this a Xi?
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.
G4bool isAntiNucleonorAntiLambda() const
Is this an antiNucleon or an antiLambda?
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)
G4bool isAntiLambda() const
Is this an antiLambda?
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.
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.
G4double dot(const ThreeVector &v) const
G4double mag2() 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
const axis_t axis_to_type< N >::axis
Definition pugixml.cc:9668
#define G4ThreadLocal
Definition tls.hh:77