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