Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/global/HEPNumerics/include/G4Integrator.hh

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 // G4Integrator
 27 //
 28 // Class description:
 29 //
 30 // Template class collecting integrator methods for generic funtions.
 31 
 32 // Author: V.Grichine, 04.09.1999 - First implementation based on
 33 //         G4SimpleIntegration class with H.P.Wellisch, G.Cosmo, and
 34 //         E.TCherniaev advises
 35 // --------------------------------------------------------------------
 36 #ifndef G4INTEGRATOR_HH
 37 #define G4INTEGRATOR_HH 1
 38 
 39 #include "G4Types.hh"
 40 #include <CLHEP/Units/PhysicalConstants.h>
 41 #include <cmath>
 42 
 43 template <class T, class F>
 44 class G4Integrator
 45 {
 46  public:
 47   G4Integrator() { ; }
 48   ~G4Integrator() { ; }
 49 
 50   G4double Simpson(T& typeT, F f, G4double a, G4double b, G4int n);
 51   G4double Simpson(T* ptrT, F f, G4double a, G4double b, G4int n);
 52   G4double Simpson(G4double (*f)(G4double), G4double a, G4double b, G4int n);
 53   // Simpson integration method
 54 
 55   G4double AdaptiveGauss(T& typeT, F f, G4double a, G4double b, G4double e);
 56   G4double AdaptiveGauss(T* ptrT, F f, G4double a, G4double b, G4double e);
 57   G4double AdaptiveGauss(G4double (*f)(G4double), G4double a, G4double b,
 58                          G4double e);
 59   // Adaptive Gauss method
 60 
 61   // Integration methods involving orthogohol polynomials
 62 
 63   G4double Legendre(T& typeT, F f, G4double a, G4double b, G4int n);
 64   G4double Legendre(T* ptrT, F f, G4double a, G4double b, G4int n);
 65   G4double Legendre(G4double (*f)(G4double), G4double a, G4double b, G4int n);
 66   //
 67   // Methods involving Legendre polynomials
 68 
 69   G4double Legendre10(T& typeT, F f, G4double a, G4double b);
 70   G4double Legendre10(T* ptrT, F f, G4double a, G4double b);
 71   G4double Legendre10(G4double (*f)(G4double), G4double a, G4double b);
 72   //
 73   // Legendre10 is very fast and accurate enough
 74 
 75   G4double Legendre96(T& typeT, F f, G4double a, G4double b);
 76   G4double Legendre96(T* ptrT, F f, G4double a, G4double b);
 77   G4double Legendre96(G4double (*f)(G4double), G4double a, G4double b);
 78   //
 79   // Legendre96 is very accurate and fast enough
 80 
 81   G4double Chebyshev(T& typeT, F f, G4double a, G4double b, G4int n);
 82   G4double Chebyshev(T* ptrT, F f, G4double a, G4double b, G4int n);
 83   G4double Chebyshev(G4double (*f)(G4double), G4double a, G4double b, G4int n);
 84   //
 85   // Methods involving Chebyshev  polynomials
 86 
 87   G4double Laguerre(T& typeT, F f, G4double alpha, G4int n);
 88   G4double Laguerre(T* ptrT, F f, G4double alpha, G4int n);
 89   G4double Laguerre(G4double (*f)(G4double), G4double alpha, G4int n);
 90   //
 91   // Method involving Laguerre polynomials
 92 
 93   G4double Hermite(T& typeT, F f, G4int n);
 94   G4double Hermite(T* ptrT, F f, G4int n);
 95   G4double Hermite(G4double (*f)(G4double), G4int n);
 96   //
 97   // Method involving Hermite polynomials
 98 
 99   G4double Jacobi(T& typeT, F f, G4double alpha, G4double beta, G4int n);
100   G4double Jacobi(T* ptrT, F f, G4double alpha, G4double beta, G4int n);
101   G4double Jacobi(G4double (*f)(G4double), G4double alpha, G4double beta,
102                   G4int n);
103   // Method involving Jacobi polynomials
104 
105  protected:
106   // Auxiliary functions for adaptive Gauss method
107 
108   G4double Gauss(T& typeT, F f, G4double a, G4double b);
109   G4double Gauss(T* ptrT, F f, G4double a, G4double b);
110   G4double Gauss(G4double (*f)(G4double), G4double a, G4double b);
111 
112   void AdaptGauss(T& typeT, F f, G4double a, G4double b, G4double e,
113                   G4double& sum, G4int& n);
114   void AdaptGauss(T* typeT, F f, G4double a, G4double b, G4double e,
115                   G4double& sum, G4int& n);
116   void AdaptGauss(G4double (*f)(G4double), G4double a, G4double b, G4double e,
117                   G4double& sum, G4int& n);
118 
119   G4double GammaLogarithm(G4double xx);
120 };
121 
122 #include "G4Integrator.icc"
123 
124 #endif
125