Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4INCLSigmaZeroDecayChannel.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// 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
41#include "G4INCLRandom.hh"
42#include "G4INCLGlobals.hh"
43
44namespace G4INCL {
45
47 :theParticle(p), incidentDirection(dir)
48 {}
49
51
52
54// assert(p->getType() == SigmaZero);
55 const G4double gamma = std::sqrt(1+std::pow(p->getMomentum().mag()/p->getMass(),2));
56 const G4double tau = ParticleTable::getWidth(SigmaZero)*3E8*1E15*gamma; // fm
57 const G4double t = -tau * std::log(Random::shoot());
58 return t;
59 }
60
61 void SigmaZeroDecayChannel::sampleAngles(G4double *ctet_par, G4double *stet_par, G4double *phi_par) {
62
63 (*ctet_par) = -1.0 + 2.0*Random::shoot();
64 if(std::abs(*ctet_par) > 1.0) (*ctet_par) = Math::sign(*ctet_par); // needed?
65 (*stet_par) = std::sqrt(1.-(*ctet_par)*(*ctet_par));
66 (*phi_par) = Math::twoPi * Random::shoot();
67 }
68
70
71// assert( theParticle->getType() == SigmaZero);
72 ParticleType createdType = Photon;
73
74 const G4double sqrtS = theParticle->getMass();
75
76 theParticle->setType(Lambda);
77 G4double phi, c_tet, s_tet;
78 sampleAngles(&c_tet, &s_tet, &phi);
79
80 G4double c_phi = std::cos(phi);
81 G4double s_phi = std::sin(phi);
82 G4double beta = incidentDirection.mag();
83
84 G4double q1, q2, q3;
85 G4double sal=0.0;
86 if (beta >= 1.0e-10)
87 sal = incidentDirection.perp()/beta;
88 if (sal >= 1.0e-6) {
89 G4double b1 = incidentDirection.getX();
90 G4double b2 = incidentDirection.getY();
91 G4double b3 = incidentDirection.getZ();
92 G4double cal = b3/beta;
93 G4double t1 = c_tet+cal*s_tet*s_phi/sal;
94 G4double t2 = s_tet/sal;
95 q1=(b1*t1+b2*t2*c_phi)/beta;
96 q2=(b2*t1-b1*t2*c_phi)/beta;
97 q3=(b3*t1/beta-t2*s_phi);
98 } else {
99 q1 = s_tet*c_phi;
100 q2 = s_tet*s_phi;
101 q3 = c_tet;
102 }
103
105 theParticle->getMass(),
106 ParticleTable::getINCLMass(createdType));
107 q1 *= xq;
108 q2 *= xq;
109 q3 *= xq;
110
111 ThreeVector createdMomentum(q1, q2, q3);
112 ThreeVector createdPosition(theParticle->getPosition());
113 Particle *createdParticle = new Particle(createdType, createdMomentum, createdPosition);
114 theParticle->setMomentum(-createdMomentum);
115 theParticle->adjustEnergyFromMomentum();
116
117 fs->addModifiedParticle(theParticle);
118 fs->addCreatedParticle(createdParticle);
119
120 }
121}
double G4double
Definition: G4Types.hh:83
void addModifiedParticle(Particle *p)
void addCreatedParticle(Particle *p)
const G4INCL::ThreeVector & getPosition() const
G4double adjustEnergyFromMomentum()
Recompute the energy to match the momentum.
const G4INCL::ThreeVector & getMomentum() const
virtual void setMomentum(const G4INCL::ThreeVector &momentum)
void setType(ParticleType t)
G4double getMass() const
Get the cached particle mass.
static G4double computeDecayTime(Particle *p)
SigmaZeroDecayChannel(Particle *, ThreeVector const &)
G4double getY() const
G4double getZ() const
G4double mag() const
G4double perp() const
G4double getX() const
G4double momentumInCM(Particle const *const p1, Particle const *const p2)
gives the momentum in the CM frame of two particles.
const G4double twoPi
G4int sign(const T t)
G4double getINCLMass(const G4int A, const G4int Z, const G4int S)
Get INCL nuclear mass (in MeV/c^2)
G4double getWidth(const ParticleType t)
Get particle width (in s)
G4double shoot()
Definition: G4INCLRandom.cc:93