! !======================================================================= ! MODULE INTEGRATION5 ! ! This module contains Gauss quadrature routines in order to integrate ! a function F over the interval [A,B]. ! ! These routines are: ! ! * SUBROUTINE GAUSSQ(KIND,N,ALPHA,BETA,KPTS,ENDPTS,B,T,W) ! USE ACCURACY_REAL ! CONTAINS ! !======================================================================= ! SUBROUTINE GAUSSQ(KIND,N,ALPHA,BETA,KPTS,ENDPTS,B,T,W) ! ! This set of routines computes the nodes T(J) and weights ! W(J) for Gaussian-type quadrature rules with pre-assigned ! nodes. these are used when one wishes to approximate ! ! integral (from a to b) f(x) w(x) dx ! ! n ! by sum w f(t ) ! j=1 j j ! ! (note w(x) and W(J) have no connection with each other.) ! here w(x) is one of six possible non-negative weight ! functions (listed below), and f(x) is the ! function to be integrated. Gaussian quadrature is particularly ! useful on infinite intervals (with appropriate weight ! functions), since then other techniques often fail. ! ! Associated with each weight function w(x) is a set of ! orthogonal polynomials. The nodes T(J) are just the zeroes ! of the proper n-th degree polynomial. ! ! input parameters (all real numbers are in double precision) ! ! KIND an integer between 1 and 6 giving the type of ! quadrature rule: ! ! KIND = 1: Legendre quadrature, w(x) = 1 on (-1, 1) ! ! KIND = 2: Chebyshev quadrature of the first kind ! w(x) = 1/sqrt(1 - x*x) on (-1, +1) ! ! KIND = 3: Chebyshev quadrature of the second kind ! w(x) = sqrt(1 - x*x) on (-1, 1) ! ! KIND = 4: Hermite quadrature, w(x) = exp(-x*x) on ! (-infinity, +infinity) ! ! KIND = 5: Jacobi quadrature, w(x) = (1-x)**ALPHA * (1+x)** ! BETA on (-1, 1), ALPHA, BETA .GT. -1. ! note: KIND=2 and 3 are a special case of this. ! ! KIND = 6: generalized Laguerre quadrature, w(x) = exp(-x)* ! x**ALPHA on (0, +infinity), ALPHA .GT. -1 ! ! N the number of points used for the quadrature rule ! ALPHA real parameter used only for Gauss-Jacobi and Gauss- ! Laguerre quadrature (otherwise use 0.d0). ! BETA real parameter used only for Gauss-Jacobi quadrature-- ! (otherwise use 0.d0) ! KPTS (integer) normally 0, unless the left or right end- ! point (or both) of the interval is required to be a ! node (this is called gauss-radau or Gauss-Lobatto ! quadrature). Then KPTS is the number of fixed ! endpoints (1 or 2). ! ENDPTS real array of length 2. Contains the values of ! any fixed endpoints, if KPTS = 1 or 2. ! B real scratch array of length N ! ! Output parameters (both double precision arrays of length N) ! ! T will contain the desired nodes. ! W will contain the desired weights W(J). ! ! Underflow may sometimes occur, but is harmless. ! ! References ! 1. Golub, G. H., and Welsch, J. H., "Calculation of Gaussian ! quadrature rules," Mathematics of Computation 23 (April, ! 1969), pp. 221-230. ! 2. Golub, G. H., "Some modified matrix eigenvalue problems," ! SIAM Review 15 (April, 1973), pp. 318-334 (section 7). ! 3. Stroud and Secrest, Gaussian Quadrature Formulas, Prentice- ! Hall, Englewood Cliffs, N.J., 1966. ! ! Original version 20 Jan 1975 from Stanford ! Modified 21 Dec 1983 by Eric Grosse ! IMTQL2 => GAUSQ2 ! hex constant => D1MACH (from core library) ! compute pi using datan ! removed accuracy claims, description of method ! added single precision version ! ! Last modified (DS) : 7 Aug 2020 ! ! USE REAL_NUMBERS, ONLY : ZERO,ONE USE SPECFUNC_SLATEC, ONLY : DGAMMA ! IMPLICIT NONE ! REAL (WP) :: B(N),T(N),W(N) REAL (WP) :: ENDPTS(2),MUZERO,T1,GAM,ALPHA,BETA ! REAL (WP) :: DSQRT ! INTEGER :: KIND,KPTS INTEGER :: N,I,IERR ! CALL CLASS(KIND,N,ALPHA,BETA,B,T,MUZERO) ! ! ! The matrix of coefficients is assumed to be symmetric. ! The array T contains the diagonal elements, the array ! B the off-diagonal elements. ! Make appropriate changes in the lower right 2 by 2 ! submatrix. ! IF (KPTS == 0) GO TO 100 ! IF (KPTS == 2) GO TO 50 ! ! ! if KPTS=1, only T(N) must be changed ! T(N) = SOLVE(ENDPTS(1),N,T,B)*B(N-1)**2 + ENDPTS(1) ! GO TO 100 ! ! ! if KPTS=2, T(N) and B(N-1) must be recomputed ! 50 GAM = SOLVE(ENDPTS(1),N,T,B) ! T1 = ((ENDPTS(1) - ENDPTS(2))/(SOLVE(ENDPTS(2),N,T,B) - GAM)) ! B(N-1) = DSQRT(T1) ! T(N) = ENDPTS(1) + GAM*T1 ! ! ! Note that the indices of the elements of B run from 1 to n-1 ! and thus the value of B(N) is arbitrary. ! Now compute the eigenvalues of the symmetric tridiagonal ! matrix, which has been modified as necessary. ! the method used is a ql-type method with origin shifting ! 100 W(1) = ONE ! ! DO I=2,N ! W(I) = ZERO ! END DO ! ! CALL GAUSQ2(N,T,B,W,IERR) ! ! DO I=1,N ! W(I) = MUZERO * W(I) * W(I) ! END DO ! ! END SUBROUTINE GAUSSQ ! !======================================================================= ! FUNCTION SOLVE(SHIFT,N,A,B) ! ! This procedure performs elimination to solve for the ! n-th component of the solution delta to the equation ! ! (jn - shift*identity) * delta = en, ! ! where en is the vector of all zeroes except for 1 in ! the n-th position. ! ! The matrix jn is symmetric tridiagonal, with diagonal ! elements A(I), off-diagonal elements B(I). This equation ! must be solved to obtain the appropriate changes in the lower ! 2 by 2 submatrix of coefficients for orthogonal polynomials. ! ! Last modified (DS) : 7 Aug 2020 ! USE REAL_NUMBERS, ONLY : ONE ! IMPLICIT NONE ! REAL (WP) :: SOLVE REAL (WP) :: SHIFT,A(N),B(N),ALPHA ! INTEGER :: N,NM1,I ! ALPHA = A(1) - SHIFT ! NM1 = N - 1 ! ! DO I=2,NM1 ! ALPHA = A(I) - SHIFT - B(I-1)**2/ALPHA ! END DO ! SOLVE = ONE/ALPHA ! ! END FUNCTION SOLVE ! !======================================================================= ! SUBROUTINE CLASS(KIND,N,ALPHA,BETA,B,A,MUZERO) ! ! This procedure supplies the coefficients A(J), B(J) of the ! recurrence relation ! ! b p (x) = (x - a ) p (x) - b p (x) ! j j j j-1 j-1 j-2 ! ! for the various classical (normalized) orthogonal polynomials, ! and the zero-th moment ! ! muzero = integral w(x) dx ! ! of the given polynomial's weight function w(x). Since the ! polynomials are orthonormalized, the tridiagonal matrix is ! guaranteed to be symmetric. ! ! The input parameter ALPHA is used only for Laguerre and ! Jacobi polynomials, and the parameter BETA is used only for ! Jacobi polynomials. The Laguerre and Jacobi polynomials ! require the Gamma function. ! ! Last modified (DS) : 4 Jun 2020 ! USE REAL_NUMBERS, ONLY : ZERO,ONE,TWO,FOUR, & HALF ! IMPLICIT NONE ! REAL (WP) :: A(N),B(N) REAL (WP) :: MUZERO,ALPHA,BETA REAL (WP) :: ABI,A2B2,DGAMMA,PI,AB ! REAL (WP) :: DATAN,DFLOAT,DSQRT ! INTEGER :: KIND,N INTEGER :: NM1,I ! PI = FOUR * DATAN(ONE) ! NM1 = N - 1 ! ! IF(KIND == 1) THEN ! ! ! KIND = 1: Legendre polynomials p(x) ! on (-1, +1), w(x) = 1. ! MUZERO = TWO ! DO I=1,NM1 ! A(I) = ZERO ! ABI = DFLOAT(I) ! B(I) = ABI/DSQRT(FOUR*ABI*ABI - ONE) END DO ! A(N) = ZERO ! RETURN ! ! ELSE IF(KIND == 2) THEN ! ! ! KIND = 2: Chebyshev polynomials of the first kind t(x) ! on (-1, +1), w(x) = 1 / sqrt(1 - x*x) ! MUZERO = PI ! ! DO I=1,NM1 ! A(I) = ZERO ! B(I) = HALF ! END DO ! ! B(1) = DSQRT(HALF) ! A(N) = ZERO ! ! RETURN ! ! ELSE IF(KIND == 3) THEN ! ! ! KIND = 3: Chebyshev polynomials of the second kind u(x) ! on (-1, +1), w(x) = sqrt(1 - x*x) ! MUZERO = PI/TWO ! ! DO I=1,NM1 ! A(I) = ZERO ! B(I) = HALF ! END DO ! ! A(N) = ZERO ! ! RETURN ! ! ELSE IF(KIND == 4) THEN ! ! ! KIND = 4: Hermite polynomials h(x) on (-infinity, ! +infinity), w(x) = exp(-x**2) ! MUZERO = DSQRT(PI) ! ! DO I=1,NM1 ! A(I) = ZERO ! B(I) = DSQRT( DFLOAT(I)/TWO ) ! END DO ! ! A(N) = ZERO ! ! RETURN ! ! ELSE IF(KIND == 5) THEN ! ! ! KIND = 5: Jacobi polynomials p(alpha, beta)(x) on ! (-1, +1), w(x) = (1-x)**alpha + (1+x)**beta, alpha and ! beta greater than -1 ! AB = ALPHA + BETA ! ABI = TWO + AB ! MUZERO = TWO ** (AB + ONE) * DGAMMA(ALPHA + ONE) * & ! DGAMMA(BETA + ONE) / DGAMMA(ABI) ! A(1) = (BETA - ALPHA) / ABI ! B(1) = DSQRT( FOUR*(ONE + ALPHA)*(ONE + BETA) / & ! ((ABI + ONE)* ABI*ABI) & ! ) ! A2B2 = BETA*BETA - ALPHA*ALPHA ! ! DO I=2,NM1 ! ABI = TWO*I + AB ! A(I) = A2B2/((ABI - TWO)*ABI) ! B(I) = DSQRT( FOUR*DFLOAT(I)*(DFLOAT(I) + ALPHA) * & ! (DFLOAT(I) + BETA)*(DFLOAT(I) + AB) / & ! ((ABI*ABI - ONE)*ABI*ABI) & ! ) ! END DO ! ! ABI = TWO*N + AB ! A(N) = A2B2/((ABI - TWO)*ABI) ! ! RETURN ! ! ELSE IF(KIND == 6) THEN ! ! ! KIND = 6: Laguerre polynomials l(alpha)(x) on ! (0, +infinity), w(x) = exp(-x) * x**alpha, alpha greater ! than -1. ! MUZERO = DGAMMA(ALPHA + ONE) ! ! DO I = 1, NM1 ! A(I) = TWO*DFLOAT(I) - ONE + ALPHA ! B(I) = DSQRT( DFLOAT(I)*(DFLOAT(I) + ALPHA) ) ! END DO ! ! A(N) = TWO*DFLOAT(N) - ONE + ALPHA ! ! RETURN ! ! END IF ! ! END SUBROUTINE CLASS ! !======================================================================= ! SUBROUTINE GAUSQ2(N,D,E,Z,IERR) ! ! This subroutine is a translation of an Algol procedure, ! Num. Math. 12, 377-383(1968) by Martin and Wilkinson, ! as modified in Num. Math. 15, 450(1970) by Dubrulle. ! Handbook for Auto. Comp., vol.ii-Linear Algebra, 241-248(1971). ! This is a modified version of the 'EISPACK' routine IMTQL2. ! ! This subroutine finds the eigenvalues and first components of the ! eigenvectors of a symmetric tridiagonal matrix by the implicit ql ! method. ! ! On input: ! ! N is the order of the matrix; ! ! D contains the diagonal elements of the input matrix; ! ! E contains the subdiagonal elements of the input matrix ! in its first N-1 positions. E(N) is arbitrary; ! ! Z contains the first row of the identity matrix. ! ! On output: ! ! D contains the eigenvalues in ascending order. If an ! error exit is made, the eigenvalues are correct but ! unordered for indices 1, 2, ..., ierr-1; ! ! E has been destroyed; ! ! Z contains the first components of the orthonormal eigenvectors ! of the symmetric tridiagonal matrix. If an error exit is ! made, Z contains the eigenvectors associated with the stored ! eigenvalues; ! ! IERR is set to ! ZERO for normal return, ! J if the j-th eigenvalue has not been ! determined after 30 iterations. ! ! ------------------------------------------------------------------ ! ! Last modified (DS) : 11 Aug 2020 ! ! USE REAL_NUMBERS, ONLY : ZERO,ONE,TWO USE MACHINE_ACCURACY, ONLY : D1MACH ! IMPLICIT NONE ! REAL (WP) :: D(N),E(N),Z(N) REAL (WP) :: B,C,F,G,P,R,S,MACHEP ! REAL (WP) :: DSQRT,DABS,DSIGN ! INTEGER :: I,J,K,L,M,N,II,MML,IERR ! MACHEP=D1MACH(4) ! ! IERR = 0 ! IF(N == 1) GO TO 1001 ! ! E(N) = ZERO ! DO L=1,N ! ! J = 0 ! ! ! :::::::::: look for small sub-diagonal element :::::::::: ! 105 DO M=L,N ! IF(M == N) GO TO 120 ! IF( DABS(E(M)) <= MACHEP * (DABS(D(M)) + DABS(D(M+1)))& ! ) GO TO 120 ! END DO ! ! 120 P = D(L) ! IF(M == L) GO TO 240 ! IF(J == 30) GO TO 1000 ! J = J + 1 ! ! :::::::::: form shift :::::::::: ! G = (D(L+1) - P) / (TWO * E(L)) ! R = DSQRT(G*G+ONE) ! G = D(M) - P + E(L) / (G + DSIGN(R,G)) ! S = ONE ! C = ONE ! P = ZERO ! MML = M - L ! ! ! :::::::::: for i=m-1 step -1 until l do -- :::::::::: ! DO II=1,MML I = M - II ! F = S * E(I) ! B = C * E(I) ! ! IF(DABS(F) < DABS(G)) GO TO 150 ! ! C = G / F ! R = DSQRT(C*C+ONE) ! E(I+1) = F * R ! S = ONE / R ! C = C * S ! GO TO 160 ! ! 150 S = F / G ! R = DSQRT(S*S+ONE) ! E(I+1) = G * R ! C = ONE / R ! S = S * C ! ! 160 G = D(I+1) - P ! R = (D(I) - G) * S + TWO * C * B ! P = S * R ! D(I+1) = G + P ! G = C * R - B ! ! ! :::::::::: form first component of vector :::::::::: ! F = Z(I+1) ! Z(I+1) = S * Z(I) + C * F ! Z(I) = C * Z(I) - S * F ! ! END DO ! ! D(L) = D(L) - P ! ! E(L) = G ! E(M) = ZERO ! GO TO 105 ! ! 240 CONTINUE ! ! END DO ! ! ! :::::::::: order eigenvalues and eigenvectors :::::::::: ! DO II=2,N ! ! I = II - 1 ! K = I ! P = D(I) ! ! DO J=II,N ! IF (D(J) >= P) GO TO 260 ! K = J ! P = D(J) ! 260 CONTINUE ! END DO ! ! IF(K == I) GO TO 300 ! D(K) = D(I) ! D(I) = P ! P = Z(I) ! Z(I) = Z(K) ! Z(K) = P ! ! 300 CONTINUE ! ! END DO ! ! GO TO 1001 ! ! ! :::::::::: set error -- no convergence to an ! eigenvalue after 30 iterations :::::::::: ! 1000 IERR = L ! 1001 RETURN ! ! ! :::::::::: last card of GAUSQ2 :::::::::: ! END SUBROUTINE GAUSQ2 ! END MODULE INTEGRATION5