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

Last change on this file since 1315 was 1315, checked in by garnier, 14 years ago

update geant4-09-04-beta-cand-01 interfaces-V09-03-09 vis-V09-03-08

File size: 8.1 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// $Id: G4CoulombScatteringModel.cc,v 1.49 2010/05/27 14:22:05 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-04-beta-cand-01 $
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:
41//
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-
47// 09.06.08 V.Ivanchenko SelectIsotope is moved to the base class
48// 16.06.09 Consolandi rows 109, 111-112, 183, 185-186
49// 27.05.10 V.Ivanchenko added G4WentzelOKandVIxSection class to
50//              compute cross sections and sample scattering angle
51//
52//
53// Class Description:
54//
55// -------------------------------------------------------------------
56//
57//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
58//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
59
60#include "G4CoulombScatteringModel.hh"
61#include "Randomize.hh"
62#include "G4ParticleChangeForGamma.hh"
63#include "G4ParticleTable.hh"
64#include "G4IonTable.hh"
65#include "G4Proton.hh"
66#include "G4NucleiProperties.hh"
67
68//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
69
70using namespace std;
71
72G4CoulombScatteringModel::G4CoulombScatteringModel(const G4String& nam)
73  : G4eCoulombScatteringModel(nam)
74{}
75
76//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
77
78G4CoulombScatteringModel::~G4CoulombScatteringModel()
79{}
80
81//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
82
83G4double G4CoulombScatteringModel::ComputeCrossSectionPerAtom(
84                                const G4ParticleDefinition* p,
85                                G4double kinEnergy, 
86                                G4double Z, 
87                                G4double, 
88                                G4double cutEnergy,
89                                G4double)
90{
91  //G4cout << "### G4CoulombScatteringModel::ComputeCrossSectionPerAtom  for "
92  //     << p->GetParticleName()<<" Z= "<<Z<<" e(MeV)= "<< kinEnergy/MeV
93  //     <<" cut(MeV)= " << cutEnergy<< G4endl;
94  G4double xsec = 0.0;
95  if(p != particle) { SetupParticle(p); }
96  if(kinEnergy < lowEnergyLimit) { return 0.0; }
97  DefineMaterial(CurrentCouple());
98
99  // Lab system
100  G4int iz = G4int(Z);
101  G4double etot = kinEnergy + mass;
102  G4double m2 = fNistManager->GetAtomicMassAmu(iz)*amu_c2;
103
104  // 03.09.2009 C.Consaldi suggested to use relativistic reduced mass
105  //            from publucation
106  // A.P. Martynenko, R.N. Faustov, Teoret. mat. Fiz. 64 (1985) 179
107  G4double Ecm  = sqrt(mass*mass + m2*m2 + 2.0*etot*m2);
108  G4double mu_rel = mass*m2/Ecm;
109  G4double tkin = Ecm - mu_rel;
110  wokvi->SetRelativisticMass(mu_rel);
111 
112  cosTetMinNuc = wokvi->SetupKinematic(tkin, currentMaterial);
113  if(cosThetaMax < cosTetMinNuc) {
114    cosTetMinNuc = wokvi->SetupTarget(iz, cutEnergy);
115    cosTetMaxNuc = cosThetaMax; 
116    if(iz == 1 && cosTetMaxNuc < 0.0 && particle == theProton) { 
117      cosTetMaxNuc = 0.0; 
118    }
119    xsec =  wokvi->ComputeNuclearCrossSection(cosTetMinNuc, cosTetMaxNuc);
120    elecRatio = wokvi->ComputeElectronCrossSection(cosTetMinNuc, cosThetaMax);
121    xsec += elecRatio;
122    if(xsec > 0.0) { elecRatio /= xsec; } 
123  }
124  /*
125  G4cout << "e(MeV)= " << kinEnergy/MeV << " xsec(b)= " << xsec/barn 
126         << "cosTetMinNuc= " << cosTetMinNuc
127         << " cosTetMaxNuc= " << cosTetMaxNuc
128         << " cosTetMaxElec= " << cosTetMaxElec
129         << " screenZ= " << screenZ
130         << " formfactA= " << formfactA << G4endl;
131  */
132  return xsec;
133}
134
135//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
136
137void G4CoulombScatteringModel::SampleSecondaries(
138                               std::vector<G4DynamicParticle*>* fvect,
139                               const G4MaterialCutsCouple* couple,
140                               const G4DynamicParticle* dp,
141                               G4double cutEnergy, 
142                               G4double)
143{
144  G4double kinEnergy = dp->GetKineticEnergy();
145  if(kinEnergy < lowEnergyLimit) { return; }
146  DefineMaterial(couple);
147  SetupParticle(dp->GetDefinition());
148
149  // Choose nucleus
150  currentElement = SelectRandomAtom(couple,particle,
151                                    kinEnergy,cutEnergy,kinEnergy);
152
153  G4double Z = currentElement->GetZ();
154  G4int iz = G4int(Z);
155  G4int ia = SelectIsotopeNumber(currentElement);
156  G4double targetMass = G4NucleiProperties::GetNuclearMass(ia, iz);
157 
158  if(ComputeCrossSectionPerAtom(particle,kinEnergy, Z,
159                                kinEnergy, cutEnergy, kinEnergy) == 0.0) 
160    { return; }
161
162  G4ThreeVector newDirection = 
163    wokvi->SampleSingleScattering(cosTetMinNuc, cosTetMaxNuc, elecRatio);
164
165  // kinematics in the Lab system
166  G4double etot = mass + kinEnergy;
167  G4double ptot = sqrt(kinEnergy*(etot + mass));
168  G4double bet  = ptot/(etot + targetMass);
169  G4double gam  = 1.0/sqrt((1.0 - bet)*(1.0 + bet));
170  G4double eCM  = sqrt(mass*mass + targetMass*targetMass + 2*targetMass*etot);
171  G4double pCM  = ptot*targetMass/eCM;
172  G4double e1   = sqrt(mass*mass + pCM*pCM);
173
174  newDirection *= pCM;
175
176  G4ThreeVector v1(newDirection.x(),newDirection.y(),gam*(newDirection.z() + bet*e1));
177  G4double finalT = gam*(e1 + bet*newDirection.z()) - mass;
178  newDirection = v1.unit();
179
180  G4ThreeVector dir = dp->GetMomentumDirection(); 
181  newDirection.rotateUz(dir);   
182  fParticleChange->ProposeMomentumDirection(newDirection);   
183 
184  // recoil
185  G4double trec = kinEnergy - finalT;
186  if(finalT <= lowEnergyLimit) { 
187    trec = kinEnergy; 
188    finalT = 0.0;
189  } 
190   
191  fParticleChange->SetProposedKineticEnergy(finalT);
192
193  //  G4cout << "sint= " << sint << " Erec(eV)= " << erec/eV << G4endl;
194
195  G4double tcut = recoilThreshold;
196  if(pCuts) { tcut= std::max(tcut,(*pCuts)[currentMaterialIndex]); } 
197  /* 
198  G4cout << "sint= " << sint << " Erec(eV)= " << erec/eV
199         << " tcut(eV)= " << tcut/eV << " th(eV)= " << recoilThreshold/eV
200         << " cut(eV)= " << (*pCuts)[currentMaterialIndex]/eV
201         << "  "  << fvect->size()
202         << G4endl;
203  */
204  if(trec > tcut) {
205    G4ParticleDefinition* ion = theParticleTable->FindIon(iz, ia, 0, iz);
206    G4double plab = sqrt(finalT*(finalT + 2.0*mass));
207    G4ThreeVector p2 = (ptot*dir - plab*newDirection).unit();
208    G4DynamicParticle* newdp  = new G4DynamicParticle(ion, p2, trec);
209    fvect->push_back(newdp);
210  } else if(trec > 0.0) {
211    fParticleChange->ProposeLocalEnergyDeposit(trec);
212    fParticleChange->ProposeNonIonizingEnergyDeposit(trec);
213  }
214}
215
216//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
217
Note: See TracBrowser for help on using the repository browser.