// // ******************************************************************** // * License and Disclaimer * // * * // * The Geant4 software is copyright of the Copyright Holders of * // * the Geant4 Collaboration. It is provided under the terms and * // * conditions of the Geant4 Software License, included in the file * // * LICENSE and available at http://cern.ch/geant4/license . These * // * include a list of copyright holders. * // * * // * Neither the authors of this software system, nor their employing * // * institutes,nor the agencies providing financial support for this * // * work make any representation or warranty, express or implied, * // * regarding this software system or assume any liability for its * // * use. Please see the license in the file LICENSE and URL above * // * for the full disclaimer and the limitation of liability. * // * * // * This code implementation is the result of the scientific and * // * technical work of the GEANT4 collaboration. * // * By using, copying, modifying or distributing the software (or * // * any work based on the software) you agree to acknowledge its * // * use in resulting scientific publications, and indicate your * // * acceptance of all terms of the Geant4 Software license. * // ******************************************************************** // // // $Id: G4GaussLegendreQ.cc,v 1.8 2007/11/13 17:35:06 gcosmo Exp $ // GEANT4 tag $Name: geant4-09-04-beta-01 $ // #include "G4GaussLegendreQ.hh" G4GaussLegendreQ::G4GaussLegendreQ( function pFunction ) : G4VGaussianQuadrature(pFunction) { } // -------------------------------------------------------------------------- // // Constructor for GaussLegendre quadrature method. The value nLegendre sets // the accuracy required, i.e the number of points where the function pFunction // will be evaluated during integration. The constructor creates the arrays for // abscissas and weights that are used in Gauss-Legendre quadrature method. // The values a and b are the limits of integration of the pFunction. // nLegendre MUST BE EVEN !!! G4GaussLegendreQ::G4GaussLegendreQ( function pFunction, G4int nLegendre ) : G4VGaussianQuadrature(pFunction) { const G4double tolerance = 1.6e-10 ; G4int k = nLegendre ; fNumber = (nLegendre + 1)/2 ; if(2*fNumber != k) { G4Exception("G4GaussLegendreQ::G4GaussLegendreQ()", "InvalidCall", FatalException, "Invalid nLegendre argument !") ; } G4double newton0=0.0, newton1=0.0, temp1=0.0, temp2=0.0, temp3=0.0, temp=0.0 ; fAbscissa = new G4double[fNumber] ; fWeight = new G4double[fNumber] ; for(G4int i=1;i<=fNumber;i++) // Loop over the desired roots { newton0 = std::cos(pi*(i - 0.25)/(k + 0.5)) ; // Initial root do // approximation { // loop of Newton's method temp1 = 1.0 ; temp2 = 0.0 ; for(G4int j=1;j<=k;j++) { temp3 = temp2 ; temp2 = temp1 ; temp1 = ((2.0*j - 1.0)*newton0*temp2 - (j - 1.0)*temp3)/j ; } temp = k*(newton0*temp1 - temp2)/(newton0*newton0 - 1.0) ; newton1 = newton0 ; newton0 = newton1 - temp1/temp ; // Newton's method } while(std::fabs(newton0 - newton1) > tolerance) ; fAbscissa[fNumber-i] = newton0 ; fWeight[fNumber-i] = 2.0/((1.0 - newton0*newton0)*temp*temp) ; } } // -------------------------------------------------------------------------- // // Returns the integral of the function to be pointed by fFunction between a // and b, by 2*fNumber point Gauss-Legendre integration: the function is // evaluated exactly 2*fNumber times at interior points in the range of // integration. Since the weights and abscissas are, in this case, symmetric // around the midpoint of the range of integration, there are actually only // fNumber distinct values of each. G4double G4GaussLegendreQ::Integral(G4double a, G4double b) const { G4double xMean = 0.5*(a + b), xDiff = 0.5*(b - a), integral = 0.0, dx = 0.0 ; for(G4int i=0;i