source: trunk/source/global/HEPNumerics/src/G4GaussLaguerreQ.cc@ 1036

Last change on this file since 1036 was 850, checked in by garnier, 17 years ago

geant4.8.2 beta

File size: 4.6 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//
27// $Id: G4GaussLaguerreQ.cc,v 1.8 2007/11/13 17:35:06 gcosmo Exp $
28// GEANT4 tag $Name: HEAD $
29//
30#include "G4GaussLaguerreQ.hh"
31
32
33
34// ------------------------------------------------------------
35//
36// Constructor for Gauss-Laguerre quadrature method: integral from zero to
37// infinity of std::pow(x,alpha)*std::exp(-x)*f(x).
38// The value of nLaguerre sets the accuracy.
39// The constructor creates arrays fAbscissa[0,..,nLaguerre-1] and
40// fWeight[0,..,nLaguerre-1] .
41//
42
43G4GaussLaguerreQ::G4GaussLaguerreQ( function pFunction,
44 G4double alpha,
45 G4int nLaguerre )
46 : G4VGaussianQuadrature(pFunction)
47{
48 const G4double tolerance = 1.0e-10 ;
49 const G4int maxNumber = 12 ;
50 G4int i=1, k=1 ;
51 G4double newton0=0.0, newton1=0.0,
52 temp1=0.0, temp2=0.0, temp3=0.0, temp=0.0, cofi=0.0 ;
53
54 fNumber = nLaguerre ;
55 fAbscissa = new G4double[fNumber] ;
56 fWeight = new G4double[fNumber] ;
57
58 for(i=1;i<=fNumber;i++) // Loop over the desired roots
59 {
60 if(i == 1)
61 {
62 newton0 = (1.0 + alpha)*(3.0 + 0.92*alpha)
63 / (1.0 + 2.4*fNumber + 1.8*alpha) ;
64 }
65 else if(i == 2)
66 {
67 newton0 += (15.0 + 6.25*alpha)/(1.0 + 0.9*alpha + 2.5*fNumber) ;
68 }
69 else
70 {
71 cofi = i - 2 ;
72 newton0 += ((1.0+2.55*cofi)/(1.9*cofi)
73 + 1.26*cofi*alpha/(1.0+3.5*cofi))
74 * (newton0 - fAbscissa[i-3])/(1.0 + 0.3*alpha) ;
75 }
76 for(k=1;k<=maxNumber;k++)
77 {
78 temp1 = 1.0 ;
79 temp2 = 0.0 ;
80 for(G4int j=1;j<=fNumber;j++)
81 {
82 temp3 = temp2 ;
83 temp2 = temp1 ;
84 temp1 = ((2*j - 1 + alpha - newton0)*temp2
85 - (j - 1 + alpha)*temp3)/j ;
86 }
87 temp = (fNumber*temp1 - (fNumber +alpha)*temp2)/newton0 ;
88 newton1 = newton0 ;
89 newton0 = newton1 - temp1/temp ;
90 if(std::fabs(newton0 - newton1) <= tolerance)
91 {
92 break ;
93 }
94 }
95 if(k > maxNumber)
96 {
97 G4Exception("G4GaussLaguerreQ::G4GaussLaguerreQ()",
98 "OutOfRange", FatalException,
99 "Too many iterations in Gauss-Laguerre constructor") ;
100 }
101
102 fAbscissa[i-1] = newton0 ;
103 fWeight[i-1] = -std::exp(GammaLogarithm(alpha + fNumber)
104 - GammaLogarithm((G4double)fNumber))/(temp*fNumber*temp2) ;
105 }
106}
107
108// -----------------------------------------------------------------
109//
110// Gauss-Laguerre method for integration of
111// std::pow(x,alpha)*std::exp(-x)*pFunction(x)
112// from zero up to infinity. pFunction is evaluated in fNumber points
113// for which fAbscissa[i] and fWeight[i] arrays were created in
114// G4VGaussianQuadrature(double,int) constructor
115
116G4double
117G4GaussLaguerreQ::Integral() const
118{
119 G4double integral = 0.0 ;
120 for(G4int i=0;i<fNumber;i++)
121 {
122 integral += fWeight[i]*fFunction(fAbscissa[i]) ;
123 }
124 return integral ;
125}
Note: See TracBrowser for help on using the repository browser.