source: trunk/source/processes/electromagnetic/lowenergy/src/G4PenelopeBremsstrahlungAngular.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.9 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//
[1192]26// $Id: G4PenelopeBremsstrahlungAngular.cc,v 1.8 2009/06/10 13:32:36 mantero Exp $
[1196]27// GEANT4 tag $Name: geant4-09-03-cand-01 $
[819]28//
29// --------------------------------------------------------------
30//
31// File name: G4PenelopeBremsstrahlungAngular
32//
33// Author: Luciano Pandola
34//
35// Creation date: February 2003
36//
37// History:
38// -----------
39// 04 Feb 2003 L. Pandola 1st implementation
40// 19 Mar 2003 L. Pandola Bugs fixed
41// 07 Nov 2003 L. Pandola Added GetAtomicNumber method for testing
42// purposes
43//----------------------------------------------------------------
44
45#include "G4PenelopeBremsstrahlungAngular.hh"
46#include "G4PenelopeInterpolator.hh"
47#include "Randomize.hh"
48#include "globals.hh"
49
50G4PenelopeBremsstrahlungAngular::G4PenelopeBremsstrahlungAngular (G4int Zed)
51 : Zmat(Zed)
52{
53 InterpolationTableForZ();
54 InterpolationForK();
55}
56
57
58G4PenelopeBremsstrahlungAngular::~G4PenelopeBremsstrahlungAngular()
59{
60}
61
62G4int G4PenelopeBremsstrahlungAngular::GetAtomicNumber()
63{
64 return Zmat;
65}
66
67void G4PenelopeBremsstrahlungAngular::InterpolationTableForZ()
68{
69 G4double pZ[NumberofZPoints] = {2.0,8.0,13.0,47.0,79.0,92.0};
70 G4double pX[NumberofZPoints],pY[NumberofZPoints];
71 G4double QQ1[NumberofZPoints][NumberofEPoints][NumberofKPoints];
72 G4double QQ2[NumberofZPoints][NumberofEPoints][NumberofKPoints];
73
74 //Read information from DataBase file
75 char* path = getenv("G4LEDATA");
76 if (!path)
77 {
78 G4String excep = "G4PenelopeBremsstrahlungAngular - G4LEDATA environment variable not set!";
79 G4Exception(excep);
80 }
81 G4String pathString(path);
82 G4String pathFile = pathString + "/penelope/br-ang-pen.dat";
83 std::ifstream file(pathFile);
84 std::filebuf* lsdp = file.rdbuf();
85
86 if (!(lsdp->is_open()))
87 {
88 G4String excep = "G4PenelopeBremsstrahlungAngular - data file " + pathFile + " not found!";
89 G4Exception(excep);
90 }
91 G4int i=0,j=0,k=0; // i=index for Z, j=index for E, k=index for K
92 G4double a1,a2;
93 while(i != -1) {
94 file >> i >> j >> k >> a1 >> a2;
95 if (i > -1){
96 QQ1[i][j][k]=a1;
97 QQ2[i][j][k]=a2;
98 }
99 }
100 file.close();
101
102
103 //Interpolation in Z
104 for (i=0;i<NumberofEPoints;i++){
105 for (j=0;j<NumberofKPoints;j++){
106 for (k=0;k<NumberofZPoints;k++){
107 pX[k]=std::log(QQ1[k][i][j]);
108 pY[k]=QQ2[k][i][j];
109 }
110 G4PenelopeInterpolator* interpolator1 = new G4PenelopeInterpolator(pZ,pX,NumberofZPoints);
111 Q1[i][j]=std::exp(interpolator1->CubicSplineInterpolation((G4double) Zmat));
112 delete interpolator1;
113 G4PenelopeInterpolator* interpolator2 = new G4PenelopeInterpolator(pZ,pY,NumberofZPoints);
114 Q2[i][j]=interpolator2->CubicSplineInterpolation((G4double) Zmat);
115 delete interpolator2;
116 }
117 }
118
119
120 //std::ofstream fil("matrice.dat",std::ios::app);
121 //fil << "Numero atomico: " << Zmat << G4endl;
122 //for (i=0;i<NumberofEPoints;i++)
123 //{
124 // fil << Q1[i][0] << " " << Q1[i][1] << " " << Q1[i][2] << " " << Q1[i][3] << G4endl;
125 //}
126 //fil.close();
127
128}
129
130void G4PenelopeBremsstrahlungAngular::InterpolationForK()
131{
132 G4double pE[NumberofEPoints] = {1.0e-03,5.0e-03,1.0e-02,5.0e-02,1.0e-01,5.0e-01};
133 G4double pK[NumberofKPoints] = {0.0,0.6,0.8,0.95};
134 G4double ppK[reducedEnergyGrid];
135 G4double pX[NumberofKPoints];
136 G4int i,j;
137
138 for(i=0;i<reducedEnergyGrid;i++){
139 ppK[i]=((G4double) i) * 0.05;
140 }
141
142 for(i=0;i<NumberofEPoints;i++){
143 betas[i]=std::sqrt(pE[i]*(pE[i]+2*electron_mass_c2))/(pE[i]+electron_mass_c2);
144 }
145
146 for (i=0;i<NumberofEPoints;i++){
147 for (j=0;j<NumberofKPoints;j++){
148 Q1[i][j]=Q1[i][j]/((G4double) Zmat);
149 }
150 }
151
152 //Expanded table of distribution parameters
153 for (i=0;i<NumberofEPoints;i++){
154 for (j=0;j<NumberofKPoints;j++){
155 pX[j]=std::log(Q1[i][j]); //logarithmic
156 }
157 G4PenelopeInterpolator* interpolator = new G4PenelopeInterpolator(pK,pX,NumberofKPoints);
158 for (j=0;j<reducedEnergyGrid;j++){
159 Q1E[i][j]=interpolator->CubicSplineInterpolation(ppK[j]);
160 }
161 delete interpolator;
162 for (j=0;j<NumberofKPoints;j++){
163 pX[j]=Q2[i][j];
164 }
165 G4PenelopeInterpolator* interpolator2 = new G4PenelopeInterpolator(pK,pX,NumberofKPoints);
166 for (j=0;j<reducedEnergyGrid;j++){
167 Q2E[i][j]=interpolator2->CubicSplineInterpolation(ppK[j]);
168 }
169 delete interpolator2;
170 }
171}
172
173G4double G4PenelopeBremsstrahlungAngular::ExtractCosTheta(G4double e1,G4double e2)
174{
175 //e1 = kinetic energy of the electron
176 //e2 = energy of the bremsstrahlung photon
177
178 G4double beta = std::sqrt(e1*(e1+2*electron_mass_c2))/(e1+electron_mass_c2);
179
180
181
182 G4double RK=20.0*e2/e1;
183 G4int ik=std::min((G4int) RK,19);
184
185 G4double P10=0,P11=0,P1=0;
186 G4double P20=0,P21=0,P2=0;
187 G4double pX[NumberofEPoints];
188 //First coefficient
189 G4int i;
190 G4int j = ik;
191 for (i=0;i<NumberofEPoints;i++){
192 pX[i]=Q1E[i][j];
193 }
194 G4PenelopeInterpolator* interpolator = new G4PenelopeInterpolator(betas,pX,NumberofEPoints);
195 P10=interpolator->CubicSplineInterpolation(beta);
196 delete interpolator;
197 j++; //(j=ik+1)
198 for (i=0;i<NumberofEPoints;i++){
199 pX[i]=Q1E[i][j];
200 }
201 G4PenelopeInterpolator* interpolator2 = new G4PenelopeInterpolator(betas,pX,NumberofEPoints);
202 P11=interpolator2->CubicSplineInterpolation(beta);
203 delete interpolator2;
204 P1=P10+(RK-(G4double) ik)*(P11-P10);
205
206 //Second coefficient
207 j = ik;
208 for (i=0;i<NumberofEPoints;i++){
209 pX[i]=Q2E[i][j];
210 }
211 G4PenelopeInterpolator* interpolator3 = new G4PenelopeInterpolator(betas,pX,NumberofEPoints);
212 P20=interpolator3->CubicSplineInterpolation(beta);
213 delete interpolator3;
214 j++; //(j=ik+1)
215 for (i=0;i<NumberofEPoints;i++){
216 pX[i]=Q2E[i][j];
217 }
218 G4PenelopeInterpolator* interpolator4 = new G4PenelopeInterpolator(betas,pX,NumberofEPoints);
219 P21=interpolator4->CubicSplineInterpolation(beta);
220 delete interpolator4;
221 P2=P20+(RK-(G4double) ik)*(P21-P20);
222
223 //Sampling from the Lorenz-trasformed dipole distributions
224 P1=std::min(std::exp(P1)/beta,1.0);
225 G4double betap = std::min(std::max(beta*(1.0+P2/beta),0.0),0.9999);
226
227 G4double cdt=0,testf=0;
228
229 if (G4UniformRand() < P1){
230 do{
231 cdt = 2.0*G4UniformRand()-1.0;
232 testf=2.0*G4UniformRand()-(1.0+cdt*cdt);
233 }while(testf>0);
234 }
235 else{
236 do{
237 cdt = 2.0*G4UniformRand()-1.0;
238 testf=G4UniformRand()-(1.0-cdt*cdt);
239 }while(testf>0);
240 }
241 cdt = (cdt+betap)/(1.0+betap*cdt);
242 return cdt;
243}
Note: See TracBrowser for help on using the repository browser.