source: trunk/source/processes/hadronic/models/coherent_elastic/src/G4LEpp.cc@ 1228

Last change on this file since 1228 was 819, checked in by garnier, 17 years ago

import all except CVS

File size: 12.3 KB
Line 
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
31#include "G4LEpp.hh"
32#include "Randomize.hh"
33#include "G4ios.hh"
34
35// Initialization of static data arrays:
36#include "G4LEppData.hh"
37
38G4LEpp::G4LEpp():G4HadronicInteraction("G4LEpp")
39{
40 // theParticleChange.SetNumberOfSecondaries(1);
41 // SetMinEnergy(10.*MeV);
42 // SetMaxEnergy(1200.*MeV);
43
44 SetCoulombEffects(0);
45
46 SetMinEnergy(0.);
47 SetMaxEnergy(1200.*GeV);
48}
49
50G4LEpp::~G4LEpp()
51{
52 // theParticleChange.Clear();
53}
54
55
56void
57G4LEpp::SetCoulombEffects(G4int State)
58{
59 if (State) {
60 for(G4int i=0; i<NANGLE; i++)
61 {
62 sig[i] = SigCoul[i];
63 }
64 elab = ElabCoul;
65 }
66 else {
67 for(G4int i=0; i<NANGLE; i++)
68 {
69 sig[i] = Sig[i];
70 }
71 elab = Elab;
72 }
73}
74
75
76G4HadFinalState*
77G4LEpp::ApplyYourself(const G4HadProjectile& aTrack, G4Nucleus& targetNucleus)
78{
79 theParticleChange.Clear();
80 const G4HadProjectile* aParticle = &aTrack;
81
82 G4double P = aParticle->GetTotalMomentum();
83 G4double Px = aParticle->Get4Momentum().x();
84 G4double Py = aParticle->Get4Momentum().y();
85 G4double Pz = aParticle->Get4Momentum().z();
86 G4double ek = aParticle->GetKineticEnergy();
87 G4ThreeVector theInitial = aParticle->Get4Momentum().vect();
88
89// if (verboseLevel > 1)
90 {
91 G4double E = aParticle->GetTotalEnergy();
92 G4double E0 = aParticle->GetDefinition()->GetPDGMass();
93 G4double Q = aParticle->GetDefinition()->GetPDGCharge();
94 G4double N = targetNucleus.GetN();
95 G4double Z = targetNucleus.GetZ();
96 G4cout << "G4LEpp:ApplyYourself: incident particle: "
97 << aParticle->GetDefinition()->GetParticleName() << G4endl;
98 G4cout << "P = " << P/GeV << " GeV/c"
99 << ", Px = " << Px/GeV << " GeV/c"
100 << ", Py = " << Py/GeV << " GeV/c"
101 << ", Pz = " << Pz/GeV << " GeV/c" << G4endl;
102 G4cout << "E = " << E/GeV << " GeV"
103 << ", kinetic energy = " << ek/GeV << " GeV"
104 << ", mass = " << E0/GeV << " GeV"
105 << ", charge = " << Q << G4endl;
106 G4cout << "G4LEpp:ApplyYourself: material:" << G4endl;
107 G4cout << "A = " << N
108 << ", Z = " << Z
109 << ", atomic mass "
110 << G4Proton::Proton()->GetPDGMass()/GeV << "GeV"
111 << G4endl;
112 //
113 // GHEISHA ADD operation to get total energy, mass, charge
114 //
115 E += G4Proton::Proton()->GetPDGMass();
116 G4double E02 = E*E - P*P;
117 E0 = std::sqrt(std::abs(E02));
118 if (E02 < 0)E0 *= -1;
119 Q += Z;
120 G4cout << "G4LEpp:ApplyYourself: total:" << G4endl;
121 G4cout << "E = " << E/GeV << " GeV"
122 << ", mass = " << E0/GeV << " GeV"
123 << ", charge = " << Q << G4endl;
124 }
125
126 // Find energy bin
127
128 G4int je1 = 0;
129 G4int je2 = NENERGY - 1;
130 ek = ek/GeV;
131 do {
132 G4int midBin = (je1 + je2)/2;
133 if (ek < elab[midBin])
134 je2 = midBin;
135 else
136 je1 = midBin;
137 } while (je2 - je1 > 1);
138 // G4int j;
139 //std::abs(ek-elab[je1]) < std::abs(ek-elab[je2]) ? j = je1 : j = je2;
140 G4double delab = elab[je2] - elab[je1];
141
142 // Sample the angle
143
144 G4float sample = G4UniformRand();
145 G4int ke1 = 0;
146 G4int ke2 = NANGLE - 1;
147 G4double dsig = sig[je2][0] - sig[je1][0];
148 G4double rc = dsig/delab;
149 G4double b = sig[je1][0] - rc*elab[je1];
150 G4double sigint1 = rc*ek + b;
151 G4double sigint2 = 0.;
152
153 if (verboseLevel > 1) G4cout << "sample=" << sample << G4endl
154 << ke1 << " " << ke2 << " "
155 << sigint1 << " " << sigint2 << G4endl;
156
157 do {
158 G4int midBin = (ke1 + ke2)/2;
159 dsig = sig[je2][midBin] - sig[je1][midBin];
160 rc = dsig/delab;
161 b = sig[je1][midBin] - rc*elab[je1];
162 G4double sigint = rc*ek + b;
163 if (sample < sigint) {
164 ke2 = midBin;
165 sigint2 = sigint;
166 }
167 else {
168 ke1 = midBin;
169 sigint1 = sigint;
170 }
171 if (verboseLevel > 1)G4cout << ke1 << " " << ke2 << " "
172 << sigint1 << " " << sigint2 << G4endl;
173 } while (ke2 - ke1 > 1);
174
175 // sigint1 and sigint2 should be recoverable from above loop
176
177 // G4double dsig = sig[je2][ke1] - sig[je1][ke1];
178 // G4double rc = dsig/delab;
179 // G4double b = sig[je1][ke1] - rc*elab[je1];
180 // G4double sigint1 = rc*ek + b;
181
182 // G4double dsig = sig[je2][ke2] - sig[je1][ke2];
183 // G4double rc = dsig/delab;
184 // G4double b = sig[je1][ke2] - rc*elab[je1];
185 // G4double sigint2 = rc*ek + b;
186
187 dsig = sigint2 - sigint1;
188 rc = 1./dsig;
189 b = ke1 - rc*sigint1;
190 G4double kint = rc*sample + b;
191 G4double theta = (0.5 + kint)*pi/180.;
192 if (theta < 0.) theta = 0.;
193
194 // G4int k;
195 //std::abs(sample-sig[j][ke1]) < std::abs(sample-sig[j][ke2]) ? k = ke1 : k = ke2;
196 // G4double theta = (0.5 + k)*pi/180.;
197
198 if (verboseLevel > 1) {
199 G4cout << " energy bin " << je1 << " energy=" << elab[je1] << G4endl;
200 G4cout << " angle bin " << kint << " angle=" << theta/degree << G4endl;
201 }
202
203
204 // Get the target particle
205
206 G4DynamicParticle* targetParticle = targetNucleus.ReturnTargetParticle();
207
208 G4double E1 = aParticle->GetTotalEnergy();
209 G4double M1 = aParticle->GetDefinition()->GetPDGMass();
210 G4double E2 = targetParticle->GetTotalEnergy();
211 G4double M2 = targetParticle->GetDefinition()->GetPDGMass();
212 G4double totalEnergy = E1 + E2;
213 G4double pseudoMass = std::sqrt(totalEnergy*totalEnergy - P*P);
214 // pseudoMass also = std::sqrt(M1*M1 + M2*M2 + 2*M2*E1)
215
216 // Transform into centre of mass system
217
218 G4double px = (M2/pseudoMass)*Px;
219 G4double py = (M2/pseudoMass)*Py;
220 G4double pz = (M2/pseudoMass)*Pz;
221 G4double p = std::sqrt(px*px + py*py + pz*pz);
222
223 if (verboseLevel > 1) {
224 G4cout << " E1, M1 (GeV) " << E1/GeV << " " << M1/GeV << G4endl;
225 G4cout << " E2, M2 (GeV) " << E2/GeV << " " << M2/GeV << G4endl;
226 G4cout << " particle 1 momentum in CM " << px/GeV << " " << py/GeV << " "
227 << pz/GeV << " " << p/GeV << G4endl;
228 }
229
230 // First scatter w.r.t. Z axis
231 G4double phi = G4UniformRand()*twopi;
232 G4double pxnew = p*std::sin(theta)*std::cos(phi);
233 G4double pynew = p*std::sin(theta)*std::sin(phi);
234 G4double pznew = p*std::cos(theta);
235
236 // Rotate according to the direction of the incident particle
237 if (px*px + py*py > 0) {
238 G4double cost, sint, ph, cosp, sinp;
239 cost = pz/p;
240 sint = (std::sqrt(std::abs((1-cost)*(1+cost))) + std::sqrt(px*px+py*py)/p)/2;
241 py < 0 ? ph = 3*halfpi : ph = halfpi;
242 if (std::abs(px) > 0.000001*GeV) ph = std::atan2(py,px);
243 cosp = std::cos(ph);
244 sinp = std::sin(ph);
245 px = (cost*cosp*pxnew - sinp*pynew + sint*cosp*pznew);
246 py = (cost*sinp*pxnew + cosp*pynew + sint*sinp*pznew);
247 pz = (-sint*pxnew + cost*pznew);
248 // G4ThreeVector it(a,b,c);
249 // p0->SetMomentum(it);
250 // G4ThreeVector aTargetMom = theInitial - it;
251 // targetParticle->SetMomentum(aTargetMom);
252 }
253 else {
254 px = pxnew;
255 py = pynew;
256 pz = pznew;
257 }
258
259 if (verboseLevel > 1) {
260 G4cout << " AFTER SCATTER..." << G4endl;
261 G4cout << " particle 1 momentum in CM " << px/GeV << " " << py/GeV << " "
262 << pz/GeV << " " << p/GeV << G4endl;
263 }
264
265 // Transform to lab system
266
267 G4double E1pM2 = E1 + M2;
268 G4double betaCM = P/E1pM2;
269 G4double betaCMx = Px/E1pM2;
270 G4double betaCMy = Py/E1pM2;
271 G4double betaCMz = Pz/E1pM2;
272 G4double gammaCM = E1pM2/std::sqrt(E1pM2*E1pM2 - P*P);
273
274 if (verboseLevel > 1) {
275 G4cout << " betaCM " << betaCMx << " " << betaCMy << " "
276 << betaCMz << " " << betaCM << G4endl;
277 G4cout << " gammaCM " << gammaCM << G4endl;
278 }
279
280 // Now following GLOREN...
281
282 G4double BETA[5], PA[5], PB[5];
283 BETA[1] = -betaCMx;
284 BETA[2] = -betaCMy;
285 BETA[3] = -betaCMz;
286 BETA[4] = gammaCM;
287
288 //The incident particle...
289
290 PA[1] = px;
291 PA[2] = py;
292 PA[3] = pz;
293 PA[4] = std::sqrt(M1*M1 + p*p);
294
295 G4double BETPA = BETA[1]*PA[1] + BETA[2]*PA[2] + BETA[3]*PA[3];
296 G4double BPGAM = (BETPA * BETA[4]/(BETA[4] + 1.) - PA[4]) * BETA[4];
297
298 PB[1] = PA[1] + BPGAM * BETA[1];
299 PB[2] = PA[2] + BPGAM * BETA[2];
300 PB[3] = PA[3] + BPGAM * BETA[3];
301 PB[4] = (PA[4] - BETPA) * BETA[4];
302
303 G4DynamicParticle* newP = new G4DynamicParticle;
304 newP->SetDefinition(const_cast<G4ParticleDefinition *>(aParticle->GetDefinition()) );
305 newP->SetMomentum(G4ThreeVector(PB[1], PB[2], PB[3]));
306
307
308 //The target particle...
309
310 PA[1] = -px;
311 PA[2] = -py;
312 PA[3] = -pz;
313 PA[4] = std::sqrt(M2*M2 + p*p);
314
315 BETPA = BETA[1]*PA[1] + BETA[2]*PA[2] + BETA[3]*PA[3];
316 BPGAM = (BETPA * BETA[4]/(BETA[4] + 1.) - PA[4]) * BETA[4];
317
318 PB[1] = PA[1] + BPGAM * BETA[1];
319 PB[2] = PA[2] + BPGAM * BETA[2];
320 PB[3] = PA[3] + BPGAM * BETA[3];
321 PB[4] = (PA[4] - BETPA) * BETA[4];
322
323 targetParticle->SetMomentum(G4ThreeVector(PB[1], PB[2], PB[3]));
324
325 // G4double ektotal = newP->GetKineticEnergy() +
326 // targetParticle->GetKineticEnergy();
327
328 if (verboseLevel > 1) {
329 G4cout << " particle 1 momentum in LAB "
330 << newP->GetMomentum()*(1./GeV)
331 << " " << newP->GetTotalMomentum()/GeV << G4endl;
332 G4cout << " particle 2 momentum in LAB "
333 << targetParticle->GetMomentum()*(1./GeV)
334 << " " << targetParticle->GetTotalMomentum()/GeV << G4endl;
335 G4cout << " TOTAL momentum in LAB "
336 << (newP->GetMomentum()+targetParticle->GetMomentum())*(1./GeV)
337 << " "
338 << (newP->GetMomentum()+targetParticle->GetMomentum()).mag()/GeV
339 << G4endl;
340 }
341
342 // if (theta < pi/2.) {
343 // G4double p = newP->GetMomentum().mag();
344 // G4ThreeVector m = newP->GetMomentum();
345 // if (p > DBL_MIN)
346 // theParticleChange.SetMomentumChange(m.x()/p, m.y()/p, m.z()/p);
347 // else
348 // theParticleChange.SetMomentumChange(0., 0., 0.);
349
350 theParticleChange.SetMomentumChange( newP->GetMomentumDirection());
351 theParticleChange.SetEnergyChange(newP->GetKineticEnergy());
352 delete newP;
353
354 // }
355 // else {
356 // // charge exchange
357 // theParticleChange.SetNumberOfSecondaries(2);
358 // theParticleChange.AddSecondary(newP);
359 // theParticleChange.SetStatusChange(fStopAndKill);
360 // // theParticleChange.SetEnergyChange(0.0);
361 // }
362
363 // Recoil particle
364 G4DynamicParticle* p1 = new G4DynamicParticle;
365 p1->SetDefinition(targetParticle->GetDefinition());
366 p1->SetMomentum(targetParticle->GetMomentum());
367 theParticleChange.AddSecondary(p1);
368
369
370 return &theParticleChange;
371}
372
373 // end of file
Note: See TracBrowser for help on using the repository browser.