Geant4  10.01.p01
G4LEHadronProtonElastic.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 
35 #include "G4PhysicalConstants.hh"
36 #include "G4SystemOfUnits.hh"
37 #include "Randomize.hh"
38 #include "G4ios.hh"
39 
40 #include "Randomize.hh"
41 
42 
44  G4HadronElastic("G4LEHadronProtonElastic")
45 {
46  SetMinEnergy(0.);
47  SetMaxEnergy(20.*MeV);
48 }
49 
51 {
53 }
54 
57  G4Nucleus& targetNucleus)
58 {
60  const G4HadProjectile* aParticle = &aTrack;
61 
62  G4double P = aParticle->GetTotalMomentum();
63  G4double Px = aParticle->Get4Momentum().x();
64  G4double Py = aParticle->Get4Momentum().y();
65  G4double Pz = aParticle->Get4Momentum().z();
66  G4double ek = aParticle->GetKineticEnergy();
67  G4ThreeVector theInitial = aParticle->Get4Momentum().vect();
68 
69  if (verboseLevel > 1)
70  {
71  G4double E = aParticle->GetTotalEnergy();
72  G4double E0 = aParticle->GetDefinition()->GetPDGMass();
73  G4double Q = aParticle->GetDefinition()->GetPDGCharge();
74  G4int A = targetNucleus.GetA_asInt();
75  G4int Z = targetNucleus.GetZ_asInt();
76  G4cout << "G4LEHadronProtonElastic:ApplyYourself: incident particle: "
77  << aParticle->GetDefinition()->GetParticleName() << G4endl;
78  G4cout << "P = " << P/GeV << " GeV/c"
79  << ", Px = " << Px/GeV << " GeV/c"
80  << ", Py = " << Py/GeV << " GeV/c"
81  << ", Pz = " << Pz/GeV << " GeV/c" << G4endl;
82  G4cout << "E = " << E/GeV << " GeV"
83  << ", kinetic energy = " << ek/GeV << " GeV"
84  << ", mass = " << E0/GeV << " GeV"
85  << ", charge = " << Q << G4endl;
86  G4cout << "G4LEHadronProtonElastic:ApplyYourself: material:" << G4endl;
87  G4cout << "A = " << A
88  << ", Z = " << Z
89  << ", atomic mass "
90  << G4Proton::Proton()->GetPDGMass()/GeV << "GeV"
91  << G4endl;
92  //
93  // GHEISHA ADD operation to get total energy, mass, charge
94  //
95  E += proton_mass_c2;
96  G4double E02 = E*E - P*P;
97  E0 = std::sqrt(std::abs(E02));
98  if (E02 < 0)E0 *= -1;
99  Q += Z;
100  G4cout << "G4LEHadronProtonElastic:ApplyYourself: total:" << G4endl;
101  G4cout << "E = " << E/GeV << " GeV"
102  << ", mass = " << E0/GeV << " GeV"
103  << ", charge = " << Q << G4endl;
104  }
105 
106  G4double theta = (0.5)*pi/180.;
107 
108  // Get the target particle
109 
110  G4DynamicParticle* targetParticle = targetNucleus.ReturnTargetParticle();
111 
112  G4double E1 = aParticle->GetTotalEnergy();
113  G4double M1 = aParticle->GetDefinition()->GetPDGMass();
114  G4double E2 = targetParticle->GetTotalEnergy();
115  G4double M2 = targetParticle->GetDefinition()->GetPDGMass();
116  G4double totalEnergy = E1 + E2;
117  G4double pseudoMass = std::sqrt(totalEnergy*totalEnergy - P*P);
118 
119  // Transform into centre of mass system
120 
121  G4double px = (M2/pseudoMass)*Px;
122  G4double py = (M2/pseudoMass)*Py;
123  G4double pz = (M2/pseudoMass)*Pz;
124  G4double p = std::sqrt(px*px + py*py + pz*pz);
125 
126  if (verboseLevel > 1) {
127  G4cout << " E1, M1 (GeV) " << E1/GeV << " " << M1/GeV << G4endl;
128  G4cout << " E2, M2 (GeV) " << E2/GeV << " " << M2/GeV << G4endl;
129  G4cout << " particle 1 momentum in CM " << px/GeV << " " << py/GeV << " "
130  << pz/GeV << " " << p/GeV << G4endl;
131  }
132 
133  // First scatter w.r.t. Z axis
134  G4double phi = G4UniformRand()*twopi;
135  G4double pxnew = p*std::sin(theta)*std::cos(phi);
136  G4double pynew = p*std::sin(theta)*std::sin(phi);
137  G4double pznew = p*std::cos(theta);
138 
139  // Rotate according to the direction of the incident particle
140  if (px*px + py*py > 0)
141  {
142  G4double cost, sint, ph, cosp, sinp;
143  cost = pz/p;
144  sint = (std::sqrt(std::fabs((1-cost)*(1+cost)))
145  + std::sqrt(px*px+py*py)/p)/2;
146  py < 0 ? ph = 3*halfpi : ph = halfpi;
147  if (std::abs(px) > 0.000001*GeV) ph = std::atan2(py,px);
148  cosp = std::cos(ph);
149  sinp = std::sin(ph);
150  px = (cost*cosp*pxnew - sinp*pynew + sint*cosp*pznew);
151  py = (cost*sinp*pxnew + cosp*pynew + sint*sinp*pznew);
152  pz = (-sint*pxnew + cost*pznew);
153  }
154  else {
155  px = pxnew;
156  py = pynew;
157  pz = pznew;
158  }
159 
160  if (verboseLevel > 1) {
161  G4cout << " AFTER SCATTER..." << G4endl;
162  G4cout << " particle 1 momentum in CM " << px/GeV
163  << " " << py/GeV << " " << pz/GeV << " " << p/GeV
164  << G4endl;
165  }
166 
167  // Transform to lab system
168 
169  G4double E1pM2 = E1 + M2;
170  G4double betaCM = P/E1pM2;
171  G4double betaCMx = Px/E1pM2;
172  G4double betaCMy = Py/E1pM2;
173  G4double betaCMz = Pz/E1pM2;
174  G4double gammaCM = E1pM2/std::sqrt(E1pM2*E1pM2 - P*P);
175 
176  if (verboseLevel > 1) {
177  G4cout << " betaCM " << betaCMx << " " << betaCMy << " "
178  << betaCMz << " " << betaCM << G4endl;
179  G4cout << " gammaCM " << gammaCM << G4endl;
180  }
181 
182  // Now following GLOREN...
183 
184  G4double BETA[5], PA[5], PB[5];
185  BETA[1] = -betaCMx;
186  BETA[2] = -betaCMy;
187  BETA[3] = -betaCMz;
188  BETA[4] = gammaCM;
189 
190  //The incident particle...
191 
192  PA[1] = px;
193  PA[2] = py;
194  PA[3] = pz;
195  PA[4] = std::sqrt(M1*M1 + p*p);
196 
197  G4double BETPA = BETA[1]*PA[1] + BETA[2]*PA[2] + BETA[3]*PA[3];
198  G4double BPGAM = (BETPA * BETA[4]/(BETA[4] + 1.) - PA[4]) * BETA[4];
199 
200  PB[1] = PA[1] + BPGAM * BETA[1];
201  PB[2] = PA[2] + BPGAM * BETA[2];
202  PB[3] = PA[3] + BPGAM * BETA[3];
203  PB[4] = (PA[4] - BETPA) * BETA[4];
204 
206  newP->SetDefinition(aParticle->GetDefinition());
207  newP->SetMomentum(G4ThreeVector(PB[1], PB[2], PB[3]));
208 
209  //The target particle...
210 
211  PA[1] = -px;
212  PA[2] = -py;
213  PA[3] = -pz;
214  PA[4] = std::sqrt(M2*M2 + p*p);
215 
216  BETPA = BETA[1]*PA[1] + BETA[2]*PA[2] + BETA[3]*PA[3];
217  BPGAM = (BETPA * BETA[4]/(BETA[4] + 1.) - PA[4]) * BETA[4];
218 
219  PB[1] = PA[1] + BPGAM * BETA[1];
220  PB[2] = PA[2] + BPGAM * BETA[2];
221  PB[3] = PA[3] + BPGAM * BETA[3];
222  PB[4] = (PA[4] - BETPA) * BETA[4];
223 
224  targetParticle->SetMomentum(G4ThreeVector(PB[1], PB[2], PB[3]));
225 
226  if (verboseLevel > 1) {
227  G4cout << " particle 1 momentum in LAB "
228  << newP->GetMomentum()*(1./GeV)
229  << " " << newP->GetTotalMomentum()/GeV << G4endl;
230  G4cout << " particle 2 momentum in LAB "
231  << targetParticle->GetMomentum()*(1./GeV)
232  << " " << targetParticle->GetTotalMomentum()/GeV << G4endl;
233  G4cout << " TOTAL momentum in LAB "
234  << (newP->GetMomentum()+targetParticle->GetMomentum())*(1./GeV)
235  << " "
236  << (newP->GetMomentum()+targetParticle->GetMomentum()).mag()/GeV
237  << G4endl;
238  }
239 
242  delete newP;
243  theParticleChange.AddSecondary(targetParticle);
244 
245  return &theParticleChange;
246 }
247 
249 //
250 // sample momentum transfer using Lab. momentum
251 
252 G4double
254  G4double plab, G4int , G4int )
255 {
256  G4double hMass = p->GetPDGMass();
257  G4double pCMS = 0.5*plab;
258  // pCMS *= 50;
259  G4double hEcms = std::sqrt(pCMS*pCMS+hMass*hMass);
260  // G4double gamma = hEcms/hMass;
261  // gamma *= 15;
262  G4double beta = pCMS/hEcms; // std::sqrt(1-1./gamma/gamma); //
263  // beta /= 0.8; // 0.95; // 1.0; // 1.1 // 0.5*pi; // pi; twopi;
264  G4double cosDipole = RandCosThetaDipPen();
265  G4double cosTheta = cosDipole + beta;
266  cosTheta /= 1. + cosDipole*beta;
267  G4double t = 2.*pCMS*pCMS*(1.-cosTheta);
268  return t;
269 
270 }
271 
273 //
274 // 1 + cos^2(theta) random distribution in the projectile rest frame, Penelope algorithm
275 
277 {
278  G4double x, cosTheta, signX, modX, power = 1./3.;
279 
280  if( G4UniformRand() > 0.25)
281  {
282  cosTheta = 2.*G4UniformRand()-1.;
283  }
284  else
285  {
286  x = 2.*G4UniformRand()-1.;
287 
288  if ( x < 0. )
289  {
290  modX = -x;
291  signX = -1.;
292  }
293  else
294  {
295  modX = x;
296  signX = 1.;
297  }
298  cosTheta = signX*std::pow(modX,power);
299  }
300  return cosTheta;
301 }
302 
303  // end of file
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
static const double MeV
Definition: G4SIunits.hh:193
G4double SampleInvariantT(const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
void SetMomentum(const G4ThreeVector &momentum)
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
G4double GetTotalEnergy() const
const G4double pi
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:95
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 G4double A[nN]
const G4LorentzVector & Get4Momentum() const
void SetEnergyChange(G4double anEnergy)
G4double GetPDGMass() const
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
void SetMaxEnergy(const G4double anEnergy)
#define G4endl
Definition: G4ios.hh:61
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
double G4double
Definition: G4Types.hh:76
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
G4double GetPDGCharge() const
void SetMomentumChange(const G4ThreeVector &aV)
G4ThreeVector GetMomentum() const
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
G4double GetTotalMomentum() const
G4double GetTotalEnergy() const