Geant4  10.00.p02
G4GaussLegendreQ.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 // $Id: G4GaussLegendreQ.cc 69546 2013-05-08 09:50:34Z gcosmo $
28 //
29 
30 #include "G4GaussLegendreQ.hh"
31 #include "G4PhysicalConstants.hh"
32 
34  : G4VGaussianQuadrature(pFunction)
35 {
36 }
37 
38 // --------------------------------------------------------------------------
39 //
40 // Constructor for GaussLegendre quadrature method. The value nLegendre sets
41 // the accuracy required, i.e the number of points where the function pFunction
42 // will be evaluated during integration. The constructor creates the arrays for
43 // abscissas and weights that are used in Gauss-Legendre quadrature method.
44 // The values a and b are the limits of integration of the pFunction.
45 // nLegendre MUST BE EVEN !!!
46 
48  G4int nLegendre )
49  : G4VGaussianQuadrature(pFunction)
50 {
51  const G4double tolerance = 1.6e-10 ;
52  G4int k = nLegendre ;
53  fNumber = (nLegendre + 1)/2 ;
54  if(2*fNumber != k)
55  {
56  G4Exception("G4GaussLegendreQ::G4GaussLegendreQ()", "InvalidCall",
57  FatalException, "Invalid nLegendre argument !") ;
58  }
59  G4double newton0=0.0, newton1=0.0,
60  temp1=0.0, temp2=0.0, temp3=0.0, temp=0.0 ;
61 
62  fAbscissa = new G4double[fNumber] ;
63  fWeight = new G4double[fNumber] ;
64 
65  for(G4int i=1;i<=fNumber;i++) // Loop over the desired roots
66  {
67  newton0 = std::cos(pi*(i - 0.25)/(k + 0.5)) ; // Initial root
68  do // approximation
69  { // loop of Newton's method
70  temp1 = 1.0 ;
71  temp2 = 0.0 ;
72  for(G4int j=1;j<=k;j++)
73  {
74  temp3 = temp2 ;
75  temp2 = temp1 ;
76  temp1 = ((2.0*j - 1.0)*newton0*temp2 - (j - 1.0)*temp3)/j ;
77  }
78  temp = k*(newton0*temp1 - temp2)/(newton0*newton0 - 1.0) ;
79  newton1 = newton0 ;
80  newton0 = newton1 - temp1/temp ; // Newton's method
81  }
82  while(std::fabs(newton0 - newton1) > tolerance) ;
83 
84  fAbscissa[fNumber-i] = newton0 ;
85  fWeight[fNumber-i] = 2.0/((1.0 - newton0*newton0)*temp*temp) ;
86  }
87 }
88 
89 // --------------------------------------------------------------------------
90 //
91 // Returns the integral of the function to be pointed by fFunction between a
92 // and b, by 2*fNumber point Gauss-Legendre integration: the function is
93 // evaluated exactly 2*fNumber times at interior points in the range of
94 // integration. Since the weights and abscissas are, in this case, symmetric
95 // around the midpoint of the range of integration, there are actually only
96 // fNumber distinct values of each.
97 
98 G4double
100 {
101  G4double xMean = 0.5*(a + b),
102  xDiff = 0.5*(b - a),
103  integral = 0.0, dx = 0.0 ;
104  for(G4int i=0;i<fNumber;i++)
105  {
106  dx = xDiff*fAbscissa[i] ;
107  integral += fWeight[i]*(fFunction(xMean + dx) + fFunction(xMean - dx)) ;
108  }
109  return integral *= xDiff ;
110 }
111 
112 // --------------------------------------------------------------------------
113 //
114 // Returns the integral of the function to be pointed by fFunction between a
115 // and b, by ten point Gauss-Legendre integration: the function is evaluated
116 // exactly ten times at interior points in the range of integration. Since the
117 // weights and abscissas are, in this case, symmetric around the midpoint of
118 // the range of integration, there are actually only five distinct values of
119 // each.
120 
121 G4double
123 {
124  // From Abramowitz M., Stegan I.A. 1964 , Handbook of Math... , p. 916
125 
126  static const G4double abscissa[] = { 0.148874338981631, 0.433395394129247,
127  0.679409568299024, 0.865063366688985,
128  0.973906528517172 } ;
129 
130  static const G4double weight[] = { 0.295524224714753, 0.269266719309996,
131  0.219086362515982, 0.149451349150581,
132  0.066671344308688 } ;
133  G4double xMean = 0.5*(a + b),
134  xDiff = 0.5*(b - a),
135  integral = 0.0, dx = 0.0 ;
136  for(G4int i=0;i<5;i++)
137  {
138  dx = xDiff*abscissa[i] ;
139  integral += weight[i]*(fFunction(xMean + dx) + fFunction(xMean - dx)) ;
140  }
141  return integral *= xDiff ;
142 }
143 
144 // -------------------------------------------------------------------------
145 //
146 // Returns the integral of the function to be pointed by fFunction between a
147 // and b, by 96 point Gauss-Legendre integration: the function is evaluated
148 // exactly ten times at interior points in the range of integration. Since the
149 // weights and abscissas are, in this case, symmetric around the midpoint of
150 // the range of integration, there are actually only five distinct values of
151 // each.
152 
153 G4double
155 {
156  // From Abramowitz M., Stegan I.A. 1964 , Handbook of Math... , p. 919
157 
158  static const
159  G4double abscissa[] = {
160  0.016276744849602969579, 0.048812985136049731112,
161  0.081297495464425558994, 0.113695850110665920911,
162  0.145973714654896941989, 0.178096882367618602759, // 6
163 
164  0.210031310460567203603, 0.241743156163840012328,
165  0.273198812591049141487, 0.304364944354496353024,
166  0.335208522892625422616, 0.365696861472313635031, // 12
167 
168  0.395797649828908603285, 0.425478988407300545365,
169  0.454709422167743008636, 0.483457973920596359768,
170  0.511694177154667673586, 0.539388108324357436227, // 18
171 
172  0.566510418561397168404, 0.593032364777572080684,
173  0.618925840125468570386, 0.644163403784967106798,
174  0.668718310043916153953, 0.692564536642171561344, // 24
175 
176  0.715676812348967626225, 0.738030643744400132851,
177  0.759602341176647498703, 0.780369043867433217604,
178  0.800308744139140817229, 0.819400310737931675539, // 30
179 
180  0.837623511228187121494, 0.854959033434601455463,
181  0.871388505909296502874, 0.886894517402420416057,
182  0.901460635315852341319, 0.915071423120898074206, // 36
183 
184  0.927712456722308690965, 0.939370339752755216932,
185  0.950032717784437635756, 0.959688291448742539300,
186  0.968326828463264212174, 0.975939174585136466453, // 42
187 
188  0.982517263563014677447, 0.988054126329623799481,
189  0.992543900323762624572, 0.995981842987209290650,
190  0.998364375863181677724, 0.999689503883230766828 // 48
191  } ;
192 
193  static const
194  G4double weight[] = {
195  0.032550614492363166242, 0.032516118713868835987,
196  0.032447163714064269364, 0.032343822568575928429,
197  0.032206204794030250669, 0.032034456231992663218, // 6
198 
199  0.031828758894411006535, 0.031589330770727168558,
200  0.031316425596862355813, 0.031010332586313837423,
201  0.030671376123669149014, 0.030299915420827593794, // 12
202 
203  0.029896344136328385984, 0.029461089958167905970,
204  0.028994614150555236543, 0.028497411065085385646,
205  0.027970007616848334440, 0.027412962726029242823, // 18
206 
207  0.026826866725591762198, 0.026212340735672413913,
208  0.025570036005349361499, 0.024900633222483610288,
209  0.024204841792364691282, 0.023483399085926219842, // 24
210 
211  0.022737069658329374001, 0.021966644438744349195,
212  0.021172939892191298988, 0.020356797154333324595,
213  0.019519081140145022410, 0.018660679627411467385, // 30
214 
215  0.017782502316045260838, 0.016885479864245172450,
216  0.015970562902562291381, 0.015038721026994938006,
217  0.014090941772314860916, 0.013128229566961572637, // 36
218 
219  0.012151604671088319635, 0.011162102099838498591,
220  0.010160770535008415758, 0.009148671230783386633,
221  0.008126876925698759217, 0.007096470791153865269, // 42
222 
223  0.006058545504235961683, 0.005014202742927517693,
224  0.003964554338444686674, 0.002910731817934946408,
225  0.001853960788946921732, 0.000796792065552012429 // 48
226  } ;
227  G4double xMean = 0.5*(a + b),
228  xDiff = 0.5*(b - a),
229  integral = 0.0, dx = 0.0 ;
230  for(G4int i=0;i<48;i++)
231  {
232  dx = xDiff*abscissa[i] ;
233  integral += weight[i]*(fFunction(xMean + dx) + fFunction(xMean - dx)) ;
234  }
235  return integral *= xDiff ;
236 }
G4double QuickIntegral(G4double a, G4double b) const
G4double Integral(G4double a, G4double b) const
const G4double pi
G4double a
Definition: TRTMaterials.hh:39
int G4int
Definition: G4Types.hh:78
G4GaussLegendreQ(function pFunction)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
double G4double
Definition: G4Types.hh:76
G4double AccurateIntegral(G4double a, G4double b) const