G4GaussLegendreQ.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 //
00027 // $Id$
00028 //
00029 
00030 #include "G4GaussLegendreQ.hh"
00031 #include "G4PhysicalConstants.hh"
00032 
00033 G4GaussLegendreQ::G4GaussLegendreQ( function pFunction )
00034    : G4VGaussianQuadrature(pFunction)
00035 {
00036 }
00037 
00038 // --------------------------------------------------------------------------
00039 //
00040 // Constructor for GaussLegendre quadrature method. The value nLegendre sets
00041 // the accuracy required, i.e the number of points where the function pFunction
00042 // will be evaluated during integration. The constructor creates the arrays for
00043 // abscissas and weights that are used in Gauss-Legendre quadrature method. 
00044 // The values a and b are the limits of integration of the pFunction.
00045 // nLegendre MUST BE EVEN !!!
00046 
00047 G4GaussLegendreQ::G4GaussLegendreQ( function pFunction,
00048                                     G4int nLegendre     )
00049    : G4VGaussianQuadrature(pFunction)
00050 {
00051    const G4double tolerance = 1.6e-10 ;
00052    G4int k = nLegendre ;
00053    fNumber = (nLegendre + 1)/2 ;
00054    if(2*fNumber != k)
00055    {
00056       G4Exception("G4GaussLegendreQ::G4GaussLegendreQ()", "InvalidCall",
00057                   FatalException, "Invalid nLegendre argument !") ;
00058    }
00059    G4double newton0=0.0, newton1=0.0,
00060             temp1=0.0, temp2=0.0, temp3=0.0, temp=0.0 ;
00061 
00062    fAbscissa = new G4double[fNumber] ;
00063    fWeight   = new G4double[fNumber] ;
00064       
00065    for(G4int i=1;i<=fNumber;i++)      // Loop over the desired roots
00066    {
00067       newton0 = std::cos(pi*(i - 0.25)/(k + 0.5)) ;  // Initial root
00068       do                                            // approximation
00069       {               // loop of Newton's method               
00070          temp1 = 1.0 ;
00071          temp2 = 0.0 ;
00072          for(G4int j=1;j<=k;j++)
00073          {
00074             temp3 = temp2 ;
00075             temp2 = temp1 ;
00076             temp1 = ((2.0*j - 1.0)*newton0*temp2 - (j - 1.0)*temp3)/j ;
00077          }
00078          temp = k*(newton0*temp1 - temp2)/(newton0*newton0 - 1.0) ;
00079          newton1 = newton0 ;
00080          newton0 = newton1 - temp1/temp ;       // Newton's method
00081       }
00082       while(std::fabs(newton0 - newton1) > tolerance) ;
00083 
00084       fAbscissa[fNumber-i] = newton0 ;
00085       fWeight[fNumber-i] = 2.0/((1.0 - newton0*newton0)*temp*temp) ;
00086    }
00087 }
00088 
00089 // --------------------------------------------------------------------------
00090 //
00091 // Returns the integral of the function to be pointed by fFunction between a
00092 // and b, by 2*fNumber point Gauss-Legendre integration: the function is
00093 // evaluated exactly 2*fNumber times at interior points in the range of
00094 // integration. Since the weights and abscissas are, in this case, symmetric
00095 // around the midpoint of the range of integration, there are actually only
00096 // fNumber distinct values of each.
00097 
00098 G4double 
00099 G4GaussLegendreQ::Integral(G4double a, G4double b) const 
00100 {
00101    G4double xMean = 0.5*(a + b),
00102             xDiff = 0.5*(b - a),
00103             integral = 0.0, dx = 0.0 ;
00104    for(G4int i=0;i<fNumber;i++)
00105    {
00106       dx = xDiff*fAbscissa[i] ;
00107       integral += fWeight[i]*(fFunction(xMean + dx) + fFunction(xMean - dx)) ;
00108    }
00109    return integral *= xDiff ;
00110 }
00111 
00112 // --------------------------------------------------------------------------
00113 //
00114 // Returns the integral of the function to be pointed by fFunction between a
00115 // and b, by ten point Gauss-Legendre integration: the function is evaluated
00116 // exactly ten times at interior points in the range of integration. Since the
00117 // weights and abscissas are, in this case, symmetric around the midpoint of
00118 // the range of integration, there are actually only five distinct values of
00119 // each.
00120 
00121 G4double 
00122    G4GaussLegendreQ::QuickIntegral(G4double a, G4double b) const 
00123 {
00124    // From Abramowitz M., Stegan I.A. 1964 , Handbook of Math... , p. 916
00125 
00126    static G4double abscissa[] = { 0.148874338981631, 0.433395394129247,
00127                                   0.679409568299024, 0.865063366688985,
00128                                   0.973906528517172                      } ;
00129    
00130    static G4double weight[] =   { 0.295524224714753, 0.269266719309996, 
00131                                   0.219086362515982, 0.149451349150581,
00132                                   0.066671344308688                      } ;
00133    G4double xMean = 0.5*(a + b),
00134             xDiff = 0.5*(b - a),
00135             integral = 0.0, dx = 0.0 ;
00136    for(G4int i=0;i<5;i++)
00137    {
00138       dx = xDiff*abscissa[i] ;
00139       integral += weight[i]*(fFunction(xMean + dx) + fFunction(xMean - dx)) ;
00140    }
00141    return integral *= xDiff ;
00142 }
00143 
00144 // -------------------------------------------------------------------------
00145 //
00146 // Returns the integral of the function to be pointed by fFunction between a
00147 // and b, by 96 point Gauss-Legendre integration: the function is evaluated
00148 // exactly ten times at interior points in the range of integration. Since the
00149 // weights and abscissas are, in this case, symmetric around the midpoint of
00150 // the range of integration, there are actually only five distinct values of
00151 // each.
00152 
00153 G4double 
00154    G4GaussLegendreQ::AccurateIntegral(G4double a, G4double b) const 
00155 {   
00156    // From Abramowitz M., Stegan I.A. 1964 , Handbook of Math... , p. 919
00157    
00158    static 
00159    G4double abscissa[] = { 
00160                            0.016276744849602969579, 0.048812985136049731112,
00161                            0.081297495464425558994, 0.113695850110665920911,
00162                            0.145973714654896941989, 0.178096882367618602759,  // 6
00163                            
00164                            0.210031310460567203603, 0.241743156163840012328,
00165                            0.273198812591049141487, 0.304364944354496353024,
00166                            0.335208522892625422616, 0.365696861472313635031,  // 12
00167    
00168                            0.395797649828908603285, 0.425478988407300545365,
00169                            0.454709422167743008636, 0.483457973920596359768,
00170                            0.511694177154667673586, 0.539388108324357436227,  // 18
00171    
00172                            0.566510418561397168404, 0.593032364777572080684,
00173                            0.618925840125468570386, 0.644163403784967106798,
00174                            0.668718310043916153953, 0.692564536642171561344,  // 24
00175    
00176                            0.715676812348967626225, 0.738030643744400132851,
00177                            0.759602341176647498703, 0.780369043867433217604,
00178                            0.800308744139140817229, 0.819400310737931675539,  // 30
00179    
00180                            0.837623511228187121494, 0.854959033434601455463,
00181                            0.871388505909296502874, 0.886894517402420416057,
00182                            0.901460635315852341319, 0.915071423120898074206,  // 36
00183    
00184                            0.927712456722308690965, 0.939370339752755216932,
00185                            0.950032717784437635756, 0.959688291448742539300,
00186                            0.968326828463264212174, 0.975939174585136466453,  // 42
00187    
00188                            0.982517263563014677447, 0.988054126329623799481,
00189                            0.992543900323762624572, 0.995981842987209290650,
00190                            0.998364375863181677724, 0.999689503883230766828   // 48
00191                                                                             } ;
00192    
00193    static 
00194    G4double weight[] =   {  
00195                            0.032550614492363166242, 0.032516118713868835987,
00196                            0.032447163714064269364, 0.032343822568575928429,
00197                            0.032206204794030250669, 0.032034456231992663218,  // 6
00198    
00199                            0.031828758894411006535, 0.031589330770727168558,
00200                            0.031316425596862355813, 0.031010332586313837423,
00201                            0.030671376123669149014, 0.030299915420827593794,  // 12
00202    
00203                            0.029896344136328385984, 0.029461089958167905970,
00204                            0.028994614150555236543, 0.028497411065085385646,
00205                            0.027970007616848334440, 0.027412962726029242823,  // 18
00206    
00207                            0.026826866725591762198, 0.026212340735672413913,
00208                            0.025570036005349361499, 0.024900633222483610288,
00209                            0.024204841792364691282, 0.023483399085926219842,  // 24
00210    
00211                            0.022737069658329374001, 0.021966644438744349195,
00212                            0.021172939892191298988, 0.020356797154333324595,
00213                            0.019519081140145022410, 0.018660679627411467385,  // 30
00214    
00215                            0.017782502316045260838, 0.016885479864245172450,
00216                            0.015970562902562291381, 0.015038721026994938006,
00217                            0.014090941772314860916, 0.013128229566961572637,  // 36
00218    
00219                            0.012151604671088319635, 0.011162102099838498591,
00220                            0.010160770535008415758, 0.009148671230783386633,
00221                            0.008126876925698759217, 0.007096470791153865269,  // 42
00222    
00223                            0.006058545504235961683, 0.005014202742927517693,
00224                            0.003964554338444686674, 0.002910731817934946408,
00225                            0.001853960788946921732, 0.000796792065552012429   // 48
00226                                                                             } ;
00227    G4double xMean = 0.5*(a + b),
00228             xDiff = 0.5*(b - a),
00229             integral = 0.0, dx = 0.0 ;
00230    for(G4int i=0;i<48;i++)
00231    {
00232       dx = xDiff*abscissa[i] ;
00233       integral += weight[i]*(fFunction(xMean + dx) + fFunction(xMean - dx)) ;
00234    }
00235    return integral *= xDiff ;
00236 }

Generated on Mon May 27 17:48:19 2013 for Geant4 by  doxygen 1.4.7