Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4INCLStrangeAbsorbtionChannel.cc
이 파일의 문서화 페이지로 가기
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 
39 #include "G4INCLKinematicsUtils.hh"
41 #include "G4INCLRandom.hh"
42 #include "G4INCLGlobals.hh"
43 #include "G4INCLLogger.hh"
44 
45 namespace G4INCL {
46 
48  : particle1(p1), particle2(p2)
49  {
50 
51  }
52 
54 
55  void StrangeAbsorbtionChannel::sampleAngles(G4double *ctet_par, G4double *stet_par, G4double *phi_par) {
56  (*ctet_par) = -1.0 + 2.0*Random::shoot();
57  if(std::abs(*ctet_par) > 1.0) (*ctet_par) = Math::sign(*ctet_par); // needed?
58  (*stet_par) = std::sqrt(1.-(*ctet_par)*(*ctet_par));
59  (*phi_par) = Math::twoPi * Random::shoot();
60  }
61 
63  Particle * nucleon;
64  Particle * strange;
65 
66  ThreeVector const incidentDirection = particle1->getMomentum() + particle2->getMomentum();
67 
68  if(particle1->isNucleon()) {
69  nucleon = particle1;
70  strange = particle2;
71  } else {
72  nucleon = particle2;
73  strange = particle1;
74  }
75 // assert(strange->isSigma() || strange->isAntiKaon());
76 
77  ParticleType finalType = Neutron;
78  if(ParticleConfig::isPair(nucleon, strange, Neutron, KZeroBar)) {
79  finalType = PiZero;
80  } else if(ParticleConfig::isPair(nucleon, strange, Proton, KZeroBar)) {
81  finalType = PiPlus;
82  } else if(ParticleConfig::isPair(nucleon, strange, Neutron, KMinus)) {
83  finalType = PiMinus;
84  } else if(ParticleConfig::isPair(nucleon, strange, Proton, KMinus)) {
85  finalType = PiZero;
86  } else if(ParticleConfig::isPair(nucleon, strange, Proton, SigmaMinus)) {
87  finalType = Neutron;
88  } else if(ParticleConfig::isPair(nucleon, strange, Neutron, SigmaZero)) {
89  finalType = Neutron;
90  } else if(ParticleConfig::isPair(nucleon, strange, Proton, SigmaZero)) {
91  finalType = Proton;
92  } else if(ParticleConfig::isPair(nucleon, strange, Neutron, SigmaPlus)) {
93  finalType = Proton;
94  } else {
95  INCL_ERROR("Unknown particle pair in Strange-N absorbtion: " << nucleon << '\t' << strange << '\n');
96  return;
97  }
98 
99  G4double energycm = KinematicsUtils::totalEnergyInCM(nucleon, strange);
100 
101  G4double finalTypemass = ParticleTable::getINCLMass(finalType);
102  nucleon->setType(Lambda); // nucleon becomes the lambda
103  G4double lambdamass = nucleon->getMass();
104 
105  G4double fi, ctet, stet;
106  sampleAngles(&ctet, &stet, &fi);
107 
108  G4double cfi = std::cos(fi);
109  G4double sfi = std::sin(fi);
110  G4double beta = incidentDirection.mag();
111 
112  G4double q1, q2, q3;
113  G4double sal=0.0;
114  if (beta >= 1.0e-10)
115  sal = incidentDirection.perp()/beta;
116  if (sal >= 1.0e-6) {
117  G4double b1 = incidentDirection.getX();
118  G4double b2 = incidentDirection.getY();
119  G4double b3 = incidentDirection.getZ();
120  G4double cal = b3/beta;
121  G4double t1 = ctet+cal*stet*sfi/sal;
122  G4double t2 = stet/sal;
123  q1=(b1*t1+b2*t2*cfi)/beta;
124  q2=(b2*t1-b1*t2*cfi)/beta;
125  q3=(b3*t1/beta-t2*sfi);
126  } else {
127  q1 = stet*cfi;
128  q2 = stet*sfi;
129  q3 = ctet;
130  }
131 
133  lambdamass,
134  finalTypemass);
135 
136  q1 *= xq;
137  q2 *= xq;
138  q3 *= xq;
139 
140  ThreeVector finalMomentum(q1, q2, q3);
141 
142  strange->setType(finalType);
143  strange->setMomentum(finalMomentum);
144  strange->adjustEnergyFromMomentum();
145  nucleon->setMomentum(-finalMomentum);
146  nucleon->adjustEnergyFromMomentum();
147 
148  fs->addModifiedParticle(nucleon); // nucleon became a lambda
149  fs->addModifiedParticle(strange); // the strange particle became an unstange particle
150  }
151 
152 }
void addModifiedParticle(Particle *p)
TTree * t1
Definition: plottest35.C:26
G4double totalEnergyInCM(Particle const *const p1, Particle const *const p2)
const G4INCL::ThreeVector & getMomentum() const
G4bool nucleon(G4int ityp)
G4double getMass() const
Get the cached particle mass.
Double_t beta
G4double getZ() const
G4double momentumInCM(Particle const *const p1, Particle const *const p2)
gives the momentum in the CM frame of two particles.
double G4double
Definition: G4Types.hh:76
G4double shoot()
Definition: G4INCLRandom.cc:93
#define INCL_ERROR(x)
G4double getINCLMass(const G4int A, const G4int Z)
Get INCL nuclear mass (in MeV/c^2)
void sampleAngles(G4double *, G4double *, G4double *)
G4double getX() const
G4double perp() const
void setType(ParticleType t)
G4int sign(const T t)
const G4double twoPi
TTree * t2
Definition: plottest35.C:36
virtual void setMomentum(const G4INCL::ThreeVector &momentum)
G4bool isPair(Particle const *const p1, Particle const *const p2, ParticleType t1, ParticleType t2)
G4double mag() const
G4double adjustEnergyFromMomentum()
Recompute the energy to match the momentum.
G4bool isNucleon() const
G4double getY() const