source: trunk/source/processes/electromagnetic/standard/src/G4CoulombScatteringModel.cc@ 1201

Last change on this file since 1201 was 1196, checked in by garnier, 16 years ago

update CVS release candidate geant4.9.3.01

File size: 7.1 KB
RevLine 
[819]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//
[1196]26// $Id: G4CoulombScatteringModel.cc,v 1.43 2009/10/28 10:14:13 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-03-cand-01 $
[819]28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class file
32//
33//
34// File name: G4CoulombScatteringModel
35//
36// Author: Vladimir Ivanchenko
37//
38// Creation date: 22.08.2005
39//
40// Modifications:
[1196]41//
[819]42// 01.08.06 V.Ivanchenko extend upper limit of table to TeV and review the
43// logic of building - only elements from G4ElementTable
44// 08.08.06 V.Ivanchenko build internal table in ekin scale, introduce faclim
45// 19.10.06 V.Ivanchenko use inheritance from G4eCoulombScatteringModel
46// 09.10.07 V.Ivanchenko reorganized methods, add cut dependence in scattering off e-
[961]47// 09.06.08 V.Ivanchenko SelectIsotope is moved to the base class
[1196]48// 16.06.09 Consolandi rows 109, 111-112, 183, 185-186
[819]49//
[1196]50//
[819]51// Class Description:
52//
53// -------------------------------------------------------------------
54//
55//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
56//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
57
58#include "G4CoulombScatteringModel.hh"
59#include "Randomize.hh"
60#include "G4ParticleChangeForGamma.hh"
61#include "G4ParticleTable.hh"
62#include "G4IonTable.hh"
63#include "G4Proton.hh"
64
65//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
66
67using namespace std;
68
[961]69G4CoulombScatteringModel::G4CoulombScatteringModel(const G4String& nam)
70 : G4eCoulombScatteringModel(nam)
71{}
[819]72
73//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
74
75G4CoulombScatteringModel::~G4CoulombScatteringModel()
76{}
77
78//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
79
80G4double G4CoulombScatteringModel::ComputeCrossSectionPerAtom(
81 const G4ParticleDefinition* p,
82 G4double kinEnergy,
83 G4double Z,
[961]84 G4double,
[819]85 G4double cutEnergy,
86 G4double)
87{
[961]88 SetupParticle(p);
[1196]89 if(kinEnergy < lowEnergyLimit) return 0.0;
90 SetupKinematic(kinEnergy, cutEnergy);
[819]91
[961]92 // save lab system kinematics
93 G4double xtkin = tkin;
94 G4double xmom2 = mom2;
95 G4double xinvb = invbeta2;
[819]96
97 // CM system
[1055]98 iz = G4int(Z);
[961]99 G4double m2 = fNistManager->GetAtomicMassAmu(iz)*amu_c2;
[819]100 G4double etot = tkin + mass;
101 G4double ptot = sqrt(mom2);
102
[961]103 G4double m12 = mass*mass;
104 G4double momCM= ptot*m2/sqrt(m12 + m2*m2 + 2.0*etot*m2);
[819]105
[961]106 mom2 = momCM*momCM;
107 tkin = sqrt(mom2 + m12) - mass;
[819]108
[1196]109 //invbeta2 = 1.0 + m12/mom2;
110
111 G4double fm = m2/(mass + m2);
112 invbeta2 = 1.0 + m12*fm*fm/mom2;
113
[961]114 SetupTarget(Z, tkin);
[819]115
[961]116 G4double xsec = CrossSectionPerAtom();
[819]117
[961]118 // restore Lab system kinematics
119 tkin = xtkin;
120 mom2 = xmom2;
121 invbeta2 = xinvb;
122
123 return xsec;
[819]124}
125
126//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
127
128void G4CoulombScatteringModel::SampleSecondaries(
129 std::vector<G4DynamicParticle*>* fvect,
130 const G4MaterialCutsCouple* couple,
131 const G4DynamicParticle* dp,
132 G4double cutEnergy,
[961]133 G4double)
[819]134{
135 G4double kinEnergy = dp->GetKineticEnergy();
[1196]136 if(kinEnergy < lowEnergyLimit) return;
[961]137 DefineMaterial(couple);
138 SetupParticle(dp->GetDefinition());
[1196]139 SetupKinematic(kinEnergy, cutEnergy);
[819]140
[961]141 // Choose nucleus
[1196]142 currentElement = SelectRandomAtom(couple,particle,
143 kinEnergy,ecut,kinEnergy);
[961]144
145 G4double Z = currentElement->GetZ();
[1055]146 iz = G4int(Z);
[961]147 G4int ia = SelectIsotopeNumber(currentElement);
148 G4double m2 = theParticleTable->GetIonTable()->GetNucleusMass(iz, ia);
[819]149
[961]150 // CM system
151 G4double etot = tkin + mass;
152 G4double ptot = sqrt(mom2);
[819]153
[961]154 G4double momCM= ptot*m2/sqrt(mass*mass + m2*m2 + 2.0*etot*m2);
155 mom2 = momCM*momCM;
156 G4double m12 = mass*mass;
157 G4double eCM = sqrt(mom2 + m12);
[819]158
[961]159 // a correction for heavy projectile
160 G4double fm = m2/(mass + m2);
161 invbeta2 = 1.0 + m12*fm*fm/mom2;
[819]162
[961]163 // sample scattering angle in CM system
164 SetupTarget(Z, eCM - mass);
[819]165
[961]166 G4double cost = SampleCosineTheta();
167 G4double z1 = 1.0 - cost;
168 if(z1 < 0.0) return;
169 G4double sint = sqrt(z1*(1.0 + cost));
170 G4double phi = twopi * G4UniformRand();
[819]171
[961]172 // kinematics in the Lab system
173 G4double bet = ptot/(etot + m2);
174 G4double gam = 1.0/sqrt((1.0 - bet)*(1.0 + bet));
175 G4double pzCM = momCM*cost;
[819]176
[961]177 G4ThreeVector v1(momCM*cos(phi)*sint,momCM*sin(phi)*sint,gam*(pzCM + bet*eCM));
[819]178 G4ThreeVector dir = dp->GetMomentumDirection();
179 G4ThreeVector newDirection = v1.unit();
180 newDirection.rotateUz(dir);
181 fParticleChange->ProposeMomentumDirection(newDirection);
[961]182
[1196]183 // G4double elab = gam*(eCM + bet*pzCM);
[819]184
[1196]185 G4double Ecm = sqrt(mass*mass + m2*m2 + 2.0*etot*m2);
186 G4double elab = etot - m2*(ptot/Ecm)*(ptot/Ecm)*(1.-cost) ;
187
188
189 G4double finalT = elab - mass;
190 if(finalT < 0.0) finalT = 0.0;
191 fParticleChange->SetProposedKineticEnergy(finalT);
192
[819]193 // recoil
[1196]194 G4double erec = kinEnergy - finalT;
[961]195
[1196]196 G4double tcut = recoilThreshold;
197 if(pCuts) { tcut= std::max(tcut,(*pCuts)[currentMaterialIndex]); }
198
199 if(erec > tcut) {
[819]200 G4ParticleDefinition* ion = theParticleTable->FindIon(iz, ia, 0, iz);
[1196]201 G4double plab = sqrt(finalT*(finalT + 2.0*mass));
[819]202 G4ThreeVector p2 = (ptot*dir - plab*newDirection).unit();
203 G4DynamicParticle* newdp = new G4DynamicParticle(ion, p2, erec);
204 fvect->push_back(newdp);
205 } else if(erec > 0.0) {
206 fParticleChange->ProposeLocalEnergyDeposit(erec);
207 fParticleChange->ProposeNonIonizingEnergyDeposit(erec);
208 }
209}
210
211//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
212
213
Note: See TracBrowser for help on using the repository browser.