Geant4-11
G4LEpp.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
27 // G4 Low energy model: n-n or p-p scattering
28 // F.W. Jones, L.G. Greeniaus, H.P. Wellisch
29
30// FWJ 27-AUG-2010: extended Coulomb-suppressed data to 5 GeV
31
32#include "G4LEpp.hh"
34#include "G4SystemOfUnits.hh"
35#include "Randomize.hh"
36#include "G4ios.hh"
37
38// Initialization of static data arrays:
39#include "G4LEppData.hh"
40
42
43
45{
47 SetMinEnergy(0.);
48 SetMaxEnergy(5.*GeV);
49}
50
52{}
53
55G4LEpp::ApplyYourself(const G4HadProjectile& aTrack, G4Nucleus& targetNucleus)
56{
58 const G4HadProjectile* aParticle = &aTrack;
59
60 G4double P = aParticle->GetTotalMomentum();
61 G4double Px = aParticle->Get4Momentum().x();
62 G4double Py = aParticle->Get4Momentum().y();
63 G4double Pz = aParticle->Get4Momentum().z();
64 G4double E = aParticle->GetTotalEnergy();
65 G4ThreeVector theInitial = aParticle->Get4Momentum().vect().unit();
66
67 if (verboseLevel > 1) {
68 G4double ek = aParticle->GetKineticEnergy();
69 G4double E0 = aParticle->GetDefinition()->GetPDGMass();
70 G4double Q = aParticle->GetDefinition()->GetPDGCharge();
71 G4int A = targetNucleus.GetA_asInt();
72 G4int Z = targetNucleus.GetZ_asInt();
73 G4cout << "G4LEpp:ApplyYourself: incident particle: "
74 << aParticle->GetDefinition()->GetParticleName() << G4endl;
75 G4cout << "P = " << P/GeV << " GeV/c"
76 << ", Px = " << Px/GeV << " GeV/c"
77 << ", Py = " << Py/GeV << " GeV/c"
78 << ", Pz = " << Pz/GeV << " GeV/c" << G4endl;
79 G4cout << "E = " << E/GeV << " GeV"
80 << ", kinetic energy = " << ek/GeV << " GeV"
81 << ", mass = " << E0/GeV << " GeV"
82 << ", charge = " << Q << G4endl;
83 G4cout << "G4LEpp:ApplyYourself: material:" << G4endl;
84 G4cout << "A = " << A
85 << ", Z = " << Z
86 << ", atomic mass "
87 << G4Proton::Proton()->GetPDGMass()/GeV << "GeV"
88 << G4endl;
89 //
90 // GHEISHA ADD operation to get total energy, mass, charge
91 //
92 E += proton_mass_c2;
93 G4double E02 = E*E - P*P;
94 E0 = std::sqrt(std::fabs(E02));
95 if (E02 < 0)E0 *= -1;
96 Q += Z;
97 G4cout << "G4LEpp:ApplyYourself: total:" << G4endl;
98 G4cout << "E = " << E/GeV << " GeV"
99 << ", mass = " << E0/GeV << " GeV"
100 << ", charge = " << Q << G4endl;
101 }
102 G4double t = SampleInvariantT(aParticle->GetDefinition(), P, 0, 0);
103 G4double cost = 1.0 - 2*t/(P*P);
104 if(cost > 1.0) { cost = 1.0; }
105 if(cost <-1.0) { cost =-1.0; }
106 G4double sint = std::sqrt((1.0 - cost)*(1.0 + cost));
108 // Get the target particle
109 G4DynamicParticle* targetParticle = targetNucleus.ReturnTargetParticle();
110
111 G4double E1 = aParticle->GetTotalEnergy();
112 G4double M1 = aParticle->GetDefinition()->GetPDGMass();
113 G4double E2 = targetParticle->GetTotalEnergy();
114 G4double M2 = targetParticle->GetDefinition()->GetPDGMass();
115 G4double totalEnergy = E1 + E2;
116 G4double pseudoMass = std::sqrt(totalEnergy*totalEnergy - P*P);
117
118 // Transform into centre of mass system
119
120 G4double px = (M2/pseudoMass)*Px;
121 G4double py = (M2/pseudoMass)*Py;
122 G4double pz = (M2/pseudoMass)*Pz;
123 G4double p = std::sqrt(px*px + py*py + pz*pz);
124
125 if (verboseLevel > 1) {
126 G4cout << " E1, M1 (GeV) " << E1/GeV << " " << M1/GeV << G4endl;
127 G4cout << " E2, M2 (GeV) " << E2/GeV << " " << M2/GeV << G4endl;
128 G4cout << " particle 1 momentum in CM " << px/GeV
129 << " " << py/GeV << " "
130 << pz/GeV << " " << p/GeV << G4endl;
131 }
132
133 // First scatter w.r.t. Z axis
134 G4double pxnew = p*sint*std::cos(phi);
135 G4double pynew = p*sint*std::sin(phi);
136 G4double pznew = p*cost;
137
138 // Rotate according to the direction of the incident particle
139 if (px*px + py*py > 0) {
140 G4double ph, cosp, sinp;
141 cost = pz/p;
142 sint = (std::sqrt((1-cost)*(1+cost)) + std::sqrt(px*px+py*py)/p)/2;
143 py < 0 ? ph = 3*halfpi : ph = halfpi;
144 if (std::fabs(px) > 0.000001*GeV) ph = std::atan2(py,px);
145 cosp = std::cos(ph);
146 sinp = std::sin(ph);
147 px = (cost*cosp*pxnew - sinp*pynew + sint*cosp*pznew);
148 py = (cost*sinp*pxnew + cosp*pynew + sint*sinp*pznew);
149 pz = (-sint*pxnew + cost*pznew);
150 }
151 else {
152 px = pxnew;
153 py = pynew;
154 pz = pznew;
155 }
156
157 if (verboseLevel > 1) {
158 G4cout << " AFTER SCATTER..." << G4endl;
159 G4cout << " particle 1 momentum in CM " << px/GeV << " " << py/GeV << " "
160 << pz/GeV << " " << p/GeV << G4endl;
161 }
162
163 // Transform to lab system
164
165 G4double E1pM2 = E1 + M2;
166 G4double betaCM = P/E1pM2;
167 G4double betaCMx = Px/E1pM2;
168 G4double betaCMy = Py/E1pM2;
169 G4double betaCMz = Pz/E1pM2;
170 G4double gammaCM = E1pM2/std::sqrt(E1pM2*E1pM2 - P*P);
171
172 if (verboseLevel > 1) {
173 G4cout << " betaCM " << betaCMx << " " << betaCMy << " "
174 << betaCMz << " " << betaCM << G4endl;
175 G4cout << " gammaCM " << gammaCM << G4endl;
176 }
177
178 // Now following GLOREN...
179
180 G4double BETA[5], PA[5], PB[5];
181 BETA[1] = -betaCMx;
182 BETA[2] = -betaCMy;
183 BETA[3] = -betaCMz;
184 BETA[4] = gammaCM;
185
186 //The incident particle...
187
188 PA[1] = px;
189 PA[2] = py;
190 PA[3] = pz;
191 PA[4] = std::sqrt(M1*M1 + p*p);
192
193 G4double BETPA = BETA[1]*PA[1] + BETA[2]*PA[2] + BETA[3]*PA[3];
194 G4double BPGAM = (BETPA * BETA[4]/(BETA[4] + 1.) - PA[4]) * BETA[4];
195
196 PB[1] = PA[1] + BPGAM * BETA[1];
197 PB[2] = PA[2] + BPGAM * BETA[2];
198 PB[3] = PA[3] + BPGAM * BETA[3];
199 PB[4] = (PA[4] - BETPA) * BETA[4];
200
202 newP->SetDefinition(aParticle->GetDefinition());
203 newP->SetMomentum(G4ThreeVector(PB[1], PB[2], PB[3]));
204
205 //The target particle...
206
207 PA[1] = -px;
208 PA[2] = -py;
209 PA[3] = -pz;
210 PA[4] = std::sqrt(M2*M2 + p*p);
211
212 BETPA = BETA[1]*PA[1] + BETA[2]*PA[2] + BETA[3]*PA[3];
213 BPGAM = (BETPA * BETA[4]/(BETA[4] + 1.) - PA[4]) * BETA[4];
214
215 PB[1] = PA[1] + BPGAM * BETA[1];
216 PB[2] = PA[2] + BPGAM * BETA[2];
217 PB[3] = PA[3] + BPGAM * BETA[3];
218 PB[4] = (PA[4] - BETPA) * BETA[4];
219
220 targetParticle->SetMomentum(G4ThreeVector(PB[1], PB[2], PB[3]));
221
222 if (verboseLevel > 1) {
223 G4cout << " particle 1 momentum in LAB "
224 << newP->GetMomentum()/GeV
225 << " " << newP->GetTotalMomentum()/GeV << G4endl;
226 G4cout << " particle 2 momentum in LAB "
227 << targetParticle->GetMomentum()/GeV
228 << " " << targetParticle->GetTotalMomentum()/GeV << G4endl;
229 G4cout << " TOTAL momentum in LAB "
230 << (newP->GetMomentum()+targetParticle->GetMomentum())/GeV
231 << " "
232 << (newP->GetMomentum()+targetParticle->GetMomentum()).mag()/GeV
233 << G4endl;
234 }
235
238 delete newP;
239
240 // Recoil particle
241 theParticleChange.AddSecondary(targetParticle, secID);
242 return &theParticleChange;
243}
244
246//
247// sample momentum transfer using Lab. momentum
248
250 G4double plab, G4int , G4int )
251{
252 G4double nMass = p->GetPDGMass(); // 939.565346*MeV;
253 G4double ek = std::sqrt(plab*plab+nMass*nMass) - nMass;
254
255 // Find energy bin
256
257 G4int je1 = 0;
258 G4int je2 = NENERGY - 1;
259 ek /= GeV;
260
261 do
262 {
263 G4int midBin = (je1 + je2)/2;
264
265 if (ek < elab[midBin]) je2 = midBin;
266 else je1 = midBin;
267 }
268 while (je2 - je1 > 1); /* Loop checking, 10.08.2015, A.Ribon */
269
270 G4double delab = elab[je2] - elab[je1];
271
272 // Sample the angle
273
274 G4double sample = G4UniformRand();
275 G4int ke1 = 0;
276 G4int ke2 = NANGLE - 1;
277 G4double dsig, b, rc;
278
279 dsig = Sig[je2][0] - Sig[je1][0];
280 rc = dsig/delab;
281 b = Sig[je1][0] - rc*elab[je1];
282
283 G4double sigint1 = rc*ek + b;
284 G4double sigint2 = 0.;
285
286 do
287 {
288 G4int midBin = (ke1 + ke2)/2;
289 dsig = Sig[je2][midBin] - Sig[je1][midBin];
290 rc = dsig/delab;
291 b = Sig[je1][midBin] - rc*elab[je1];
292 G4double sigint = rc*ek + b;
293
294 if (sample < sigint)
295 {
296 ke2 = midBin;
297 sigint2 = sigint;
298 }
299 else
300 {
301 ke1 = midBin;
302 sigint1 = sigint;
303 }
304 }
305 while (ke2 - ke1 > 1); /* Loop checking, 10.08.2015, A.Ribon */
306
307 dsig = sigint2 - sigint1;
308 rc = 1./dsig;
309 b = ke1 - rc*sigint1;
310
311 G4double kint = rc*sample + b;
312 G4double theta = (0.5 + kint)*pi/180.;
313 G4double t = 0.5*plab*plab*(1 - std::cos(theta));
314
315 return t;
316}
317// end of file
static constexpr double twopi
Definition: G4SIunits.hh:56
static constexpr double GeV
Definition: G4SIunits.hh:203
static constexpr double pi
Definition: G4SIunits.hh:55
static constexpr double halfpi
Definition: G4SIunits.hh:57
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector unit() const
Hep3Vector vect() const
const G4ThreeVector & GetMomentumDirection() const
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetTotalEnergy() const
void SetMomentum(const G4ThreeVector &momentum)
G4ThreeVector GetMomentum() const
G4double GetTotalMomentum() const
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
void SetEnergyChange(G4double anEnergy)
void SetMomentumChange(const G4ThreeVector &aV)
G4double GetTotalMomentum() const
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
G4double GetTotalEnergy() const
void SetMinEnergy(G4double anEnergy)
const G4String & GetModelName() const
void SetMaxEnergy(const G4double anEnergy)
G4double SampleInvariantT(const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A) override
Definition: G4LEpp.cc:249
~G4LEpp() override
Definition: G4LEpp.cc:51
@ NENERGY
Definition: G4LEpp.hh:45
@ NANGLE
Definition: G4LEpp.hh:45
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus) override
Definition: G4LEpp.cc:55
static const G4float Sig[NENERGY][NANGLE]
Definition: G4LEpp.hh:65
G4LEpp()
Definition: G4LEpp.cc:44
static const G4float elab[NENERGY]
Definition: G4LEpp.hh:66
G4int GetA_asInt() const
Definition: G4Nucleus.hh:99
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:105
G4DynamicParticle * ReturnTargetParticle() const
Definition: G4Nucleus.cc:340
G4double GetPDGCharge() const
const G4String & GetParticleName() const
static G4int GetModelID(const G4int modelIndex)
static G4Proton * Proton()
Definition: G4Proton.cc:92
float proton_mass_c2
Definition: hepunit.py:274
static double Q[]
static double P[]