Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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 // theParticleChange.SetNumberOfSecondaries(1);
44 // SetMinEnergy(10.*MeV);
45 // SetMaxEnergy(1200.*MeV);
46
48
49 SetMinEnergy(0.);
50 SetMaxEnergy(5.*GeV);
51}
52
54{
55 // theParticleChange.Clear();
56}
57
58
59void
61{
62 if (State) {
63 for(G4int i=0; i<NANGLE; i++)
64 {
65 sig[i] = SigCoul[i];
66 }
67 elab = ElabCoul;
68 SetMaxEnergy(1.2*GeV);
69 }
70 else {
71 for(G4int i=0; i<NANGLE; i++)
72 {
73 sig[i] = Sig[i];
74 }
75 elab = Elab;
76 SetMaxEnergy(5.*GeV);
77 }
78}
79
80
82G4LEpp::ApplyYourself(const G4HadProjectile& aTrack, G4Nucleus& targetNucleus)
83{
85 const G4HadProjectile* aParticle = &aTrack;
86
87 G4double P = aParticle->GetTotalMomentum();
88 G4double Px = aParticle->Get4Momentum().x();
89 G4double Py = aParticle->Get4Momentum().y();
90 G4double Pz = aParticle->Get4Momentum().z();
91 G4double ek = aParticle->GetKineticEnergy();
92 G4ThreeVector theInitial = aParticle->Get4Momentum().vect();
93
94 if (verboseLevel > 1) {
95 G4double E = aParticle->GetTotalEnergy();
96 G4double E0 = aParticle->GetDefinition()->GetPDGMass();
97 G4double Q = aParticle->GetDefinition()->GetPDGCharge();
98 G4int A = targetNucleus.GetA_asInt();
99 G4int Z = targetNucleus.GetZ_asInt();
100 G4cout << "G4LEpp:ApplyYourself: incident particle: "
101 << aParticle->GetDefinition()->GetParticleName() << G4endl;
102 G4cout << "P = " << P/GeV << " GeV/c"
103 << ", Px = " << Px/GeV << " GeV/c"
104 << ", Py = " << Py/GeV << " GeV/c"
105 << ", Pz = " << Pz/GeV << " GeV/c" << G4endl;
106 G4cout << "E = " << E/GeV << " GeV"
107 << ", kinetic energy = " << ek/GeV << " GeV"
108 << ", mass = " << E0/GeV << " GeV"
109 << ", charge = " << Q << G4endl;
110 G4cout << "G4LEpp:ApplyYourself: material:" << G4endl;
111 G4cout << "A = " << A
112 << ", Z = " << Z
113 << ", atomic mass "
114 << G4Proton::Proton()->GetPDGMass()/GeV << "GeV"
115 << G4endl;
116 //
117 // GHEISHA ADD operation to get total energy, mass, charge
118 //
119 E += proton_mass_c2;
120 G4double E02 = E*E - P*P;
121 E0 = std::sqrt(std::fabs(E02));
122 if (E02 < 0)E0 *= -1;
123 Q += Z;
124 G4cout << "G4LEpp:ApplyYourself: total:" << G4endl;
125 G4cout << "E = " << E/GeV << " GeV"
126 << ", mass = " << E0/GeV << " GeV"
127 << ", charge = " << Q << G4endl;
128 }
129
130 // Find energy bin
131
132 G4int je1 = 0;
133 G4int je2 = NENERGY - 1;
134 ek = ek/GeV;
135 do {
136 G4int midBin = (je1 + je2)/2;
137 if (ek < elab[midBin])
138 je2 = midBin;
139 else
140 je1 = midBin;
141 } while (je2 - je1 > 1);
142 G4double delab = elab[je2] - elab[je1];
143
144 // Sample the angle
145
146 G4float sample = G4UniformRand();
147 G4int ke1 = 0;
148 G4int ke2 = NANGLE - 1;
149 G4double dsig = sig[je2][0] - sig[je1][0];
150 G4double rc = dsig/delab;
151 G4double b = sig[je1][0] - rc*elab[je1];
152 G4double sigint1 = rc*ek + b;
153 G4double sigint2 = 0.;
154
155 if (verboseLevel > 1) G4cout << "sample=" << sample << G4endl
156 << ke1 << " " << ke2 << " "
157 << sigint1 << " " << sigint2 << G4endl;
158
159 do {
160 G4int midBin = (ke1 + ke2)/2;
161 dsig = sig[je2][midBin] - sig[je1][midBin];
162 rc = dsig/delab;
163 b = sig[je1][midBin] - rc*elab[je1];
164 G4double sigint = rc*ek + b;
165 if (sample < sigint) {
166 ke2 = midBin;
167 sigint2 = sigint;
168 }
169 else {
170 ke1 = midBin;
171 sigint1 = sigint;
172 }
173 if (verboseLevel > 1)G4cout << ke1 << " " << ke2 << " "
174 << sigint1 << " " << sigint2 << G4endl;
175 } while (ke2 - ke1 > 1);
176
177 dsig = sigint2 - sigint1;
178 rc = 1./dsig;
179 b = ke1 - rc*sigint1;
180 G4double kint = rc*sample + b;
181 G4double theta = (0.5 + kint)*pi/180.;
182 if (theta < 0.) { theta = 0.; }
183
184 if (verboseLevel > 1) {
185 G4cout << " energy bin " << je1 << " energy=" << elab[je1] << G4endl;
186 G4cout << " angle bin " << kint << " angle=" << theta/degree << G4endl;
187 }
188
189 // Get the target particle
190 G4DynamicParticle* targetParticle = targetNucleus.ReturnTargetParticle();
191
192 G4double E1 = aParticle->GetTotalEnergy();
193 G4double M1 = aParticle->GetDefinition()->GetPDGMass();
194 G4double E2 = targetParticle->GetTotalEnergy();
195 G4double M2 = targetParticle->GetDefinition()->GetPDGMass();
196 G4double totalEnergy = E1 + E2;
197 G4double pseudoMass = std::sqrt(totalEnergy*totalEnergy - P*P);
198
199 // Transform into centre of mass system
200
201 G4double px = (M2/pseudoMass)*Px;
202 G4double py = (M2/pseudoMass)*Py;
203 G4double pz = (M2/pseudoMass)*Pz;
204 G4double p = std::sqrt(px*px + py*py + pz*pz);
205
206 if (verboseLevel > 1) {
207 G4cout << " E1, M1 (GeV) " << E1/GeV << " " << M1/GeV << G4endl;
208 G4cout << " E2, M2 (GeV) " << E2/GeV << " " << M2/GeV << G4endl;
209 G4cout << " particle 1 momentum in CM " << px/GeV << " " << py/GeV << " "
210 << pz/GeV << " " << p/GeV << G4endl;
211 }
212
213 // First scatter w.r.t. Z axis
214 G4double phi = G4UniformRand()*twopi;
215 G4double pxnew = p*std::sin(theta)*std::cos(phi);
216 G4double pynew = p*std::sin(theta)*std::sin(phi);
217 G4double pznew = p*std::cos(theta);
218
219 // Rotate according to the direction of the incident particle
220 if (px*px + py*py > 0) {
221 G4double cost, sint, ph, cosp, sinp;
222 cost = pz/p;
223 sint = (std::sqrt(std::fabs((1-cost)*(1+cost))) + std::sqrt(px*px+py*py)/p)/2;
224 py < 0 ? ph = 3*halfpi : ph = halfpi;
225 if (std::fabs(px) > 0.000001*GeV) ph = std::atan2(py,px);
226 cosp = std::cos(ph);
227 sinp = std::sin(ph);
228 px = (cost*cosp*pxnew - sinp*pynew + sint*cosp*pznew);
229 py = (cost*sinp*pxnew + cosp*pynew + sint*sinp*pznew);
230 pz = (-sint*pxnew + cost*pznew);
231 }
232 else {
233 px = pxnew;
234 py = pynew;
235 pz = pznew;
236 }
237
238 if (verboseLevel > 1) {
239 G4cout << " AFTER SCATTER..." << G4endl;
240 G4cout << " particle 1 momentum in CM " << px/GeV << " " << py/GeV << " "
241 << pz/GeV << " " << p/GeV << G4endl;
242 }
243
244 // Transform to lab system
245
246 G4double E1pM2 = E1 + M2;
247 G4double betaCM = P/E1pM2;
248 G4double betaCMx = Px/E1pM2;
249 G4double betaCMy = Py/E1pM2;
250 G4double betaCMz = Pz/E1pM2;
251 G4double gammaCM = E1pM2/std::sqrt(E1pM2*E1pM2 - P*P);
252
253 if (verboseLevel > 1) {
254 G4cout << " betaCM " << betaCMx << " " << betaCMy << " "
255 << betaCMz << " " << betaCM << G4endl;
256 G4cout << " gammaCM " << gammaCM << G4endl;
257 }
258
259 // Now following GLOREN...
260
261 G4double BETA[5], PA[5], PB[5];
262 BETA[1] = -betaCMx;
263 BETA[2] = -betaCMy;
264 BETA[3] = -betaCMz;
265 BETA[4] = gammaCM;
266
267 //The incident particle...
268
269 PA[1] = px;
270 PA[2] = py;
271 PA[3] = pz;
272 PA[4] = std::sqrt(M1*M1 + p*p);
273
274 G4double BETPA = BETA[1]*PA[1] + BETA[2]*PA[2] + BETA[3]*PA[3];
275 G4double BPGAM = (BETPA * BETA[4]/(BETA[4] + 1.) - PA[4]) * BETA[4];
276
277 PB[1] = PA[1] + BPGAM * BETA[1];
278 PB[2] = PA[2] + BPGAM * BETA[2];
279 PB[3] = PA[3] + BPGAM * BETA[3];
280 PB[4] = (PA[4] - BETPA) * BETA[4];
281
283 newP->SetDefinition(const_cast<G4ParticleDefinition *>(aParticle->GetDefinition()) );
284 newP->SetMomentum(G4ThreeVector(PB[1], PB[2], PB[3]));
285
286 //The target particle...
287
288 PA[1] = -px;
289 PA[2] = -py;
290 PA[3] = -pz;
291 PA[4] = std::sqrt(M2*M2 + p*p);
292
293 BETPA = BETA[1]*PA[1] + BETA[2]*PA[2] + BETA[3]*PA[3];
294 BPGAM = (BETPA * BETA[4]/(BETA[4] + 1.) - PA[4]) * BETA[4];
295
296 PB[1] = PA[1] + BPGAM * BETA[1];
297 PB[2] = PA[2] + BPGAM * BETA[2];
298 PB[3] = PA[3] + BPGAM * BETA[3];
299 PB[4] = (PA[4] - BETPA) * BETA[4];
300
301 targetParticle->SetMomentum(G4ThreeVector(PB[1], PB[2], PB[3]));
302
303 if (verboseLevel > 1) {
304 G4cout << " particle 1 momentum in LAB "
305 << newP->GetMomentum()/GeV
306 << " " << newP->GetTotalMomentum()/GeV << G4endl;
307 G4cout << " particle 2 momentum in LAB "
308 << targetParticle->GetMomentum()/GeV
309 << " " << targetParticle->GetTotalMomentum()/GeV << G4endl;
310 G4cout << " TOTAL momentum in LAB "
311 << (newP->GetMomentum()+targetParticle->GetMomentum())/GeV
312 << " "
313 << (newP->GetMomentum()+targetParticle->GetMomentum()).mag()/GeV
314 << G4endl;
315 }
316
319 delete newP;
320
321 // Recoil particle
322 theParticleChange.AddSecondary(targetParticle);
323 return &theParticleChange;
324}
325
326 // end of file
#define State(theXInfo)
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:64
float G4float
Definition: G4Types.hh:65
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
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)
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)
void SetMaxEnergy(const G4double anEnergy)
~G4LEpp()
Definition: G4LEpp.cc:53
G4LEpp()
Definition: G4LEpp.cc:41
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
Definition: G4LEpp.cc:82
void SetCoulombEffects(G4int State)
Definition: G4LEpp.cc:60
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
G4DynamicParticle * ReturnTargetParticle() const
Definition: G4Nucleus.cc:227
G4double GetPDGCharge() const
const G4String & GetParticleName() const
static G4Proton * Proton()
Definition: G4Proton.cc:93