These functions approximate the integral of a real valued function of two
variables over a rectangular region, using the Gauss-Legendre Quadrature
method of specified order in each variable.
a, b, c, d, Z, Result | Subroutine |
Short-precision real | SGLNQ2 |
Long-precision real | DGLNQ2 |
Fortran | SGLNQ2 | DGLNQ2 (subf, a, b, n1, c, d, n2, z, ldz) |
C and C++ | sglnq2 | dglnq2 (subf, a, b, n1, c, d, n2, z, ldz); |
PL/I | SGLNQ2 | DGLNQ2 (subf, a, b, n1, c, d, n2, z, ldz); |
Specified as: subf must be declared as an external subroutine in your application program. It can be whatever name you choose.
The integral:
is approximated for a real valued function of two variables s and t, over a rectangular region, using the Gauss-Legendre Quadrature method of specified order in each variable. The region of integration is:
The method gives a good approximation when your integrand is closely approximated by a function of the form f(s, t), where f is a polynomial of degree less than 2(n1) for s and 2(n2) for t. See the function description for SGLNQ and DGLNQ--Numerical Quadrature Performed on a Function Using Gauss-Legendre Quadrature and references [26] and [92]. The result is returned as the function value.
To achieve optimal performance in this subroutine and in the functional evaluation, specify the first variable integrated in this subroutine as the variable having more points. The first variable integrated is the variable in the inner integral. For example, in the following integration, x is the first variable integrated:
This is the suggested order of integration if the x variable has more points than the y variable. On the other hand, if the y variable has more points, you make y the first variable integrated.
Because the order of integration does not matter to the resulting approximation, you may be able to reverse the order that x and y are integrated and get better performance. This can be expressed as:
Results are mathematically equivalent. However, because the algorithm is computed in a different way, results may not be bitwise identical.
Table 163 shows how to assign your variables to the _GLNQ2 and
subf arguments for the x-y integration shown on the
left and for the y-x integration shown on the right.
For examples of how to do each of these, see Example 1 and Example 2.
Table 163. How to Assign Your Variables for x-y Integration Versus y-x Integration
_GLNQ2 and SUBF Arguments |
Variables for x-y Integration |
Variables for y-x Integration |
---|---|---|
For _GLNQ2: a b n1 c d n2 For subf: s t n1 n2 |
|
|
None
This example shows how to compute the integral of the function f given by:
over the intervals (0.0, 2.0) for the first variable x and (-2.0, -1.0) for the second variable y, using the Gauss-Legendre method with 10 points in the x variable and 5 points in the y variable:
Because the variable x has more points, it is the first variable integrated. This allows the SGLNQ2 subroutine and the FUN1 evaluation to achieve optimal performance. Therefore, the x and y variables correspond to S and T in the FUN1 subroutine. Also, the x and y variables correspond to the A, B, N1 and C, D, N2 sets of arguments, respectively, for SGLNQ2.
The user-supplied subroutine FUN1, which evaluates the integrand function, is coded in Fortran as follows:
SUBROUTINE FUN1 (S,N1,T,N2,Z,LDZ) INTEGER*4 N1,N2,LDZ REAL*4 S(*),T(*),Z(LDZ,*) DO 1 J=1,N2 DO 2 I=1,N1 2 Z(I,J)=EXP(S(I))*SIN(T(J)) 1 CONTINUE RETURN END
. . . DO 1 I=1,N1 1 T1(I)=EXP(S(I)) DO 2 J=1,N2 2 T2(J)=SIN(T(J)) DO 3 J=1,N2 DO 4 I=1,N1 4 Z(I,J)=T1(I)*T2(J) 3 CONTINUE . . .
When coding your application, this is the preferred technique. It reduces the number of evaluations performed and, therefore, provides better performance.
The user-supplied subroutine FUN1, which evaluates the integrand function, is coded in C as follows:
void fun1(s, n1, t, n2, z, ldz) float *s, *t, *z; int *n1, *n2, *ldz; { int i, j; for(j = 0; j < *n2; ++j, z += *ldz) { for(i = 0; i < *n1; ++i) z[i] = exp(s[i]) * sin(t[j]); } }
The user-supplied subroutine FUN1, which evaluates the integrand function, is coded in C++ as follows:
void fun1(float *s, int *n1, float *t, int *n2, float *z, int *ldz) { int i, j; for(j = 0; j < *n2; ++j, z += *ldz) { for(i = 0; i < *n1; ++i) z[i] = exp(s[i]) * sin(t[j]); } }
The user-supplied subroutine FUN1, which evaluates the integrand function, is coded in PL/I as follows:
FUN1: PROCEDURE(S,N1,T,N2,Z,LDZ) OPTIONS(FORTRAN,NOMAP); DCL (N1,N2,LDZ,I,J) REAL FIXED BINARY(31,0); DCL (S(10),T(10),Z(5,10)) REAL FLOAT DEC(16) ALIGNED CONNECTED; DO J=1 TO N1; DO I=1 TO N2; Z(I,J)=EXP(S(J))*SIN(T(I)); END; END; RETURN; END FUN1;
EXTERNAL FUN1 . . . SUBF A B N1 C D N2 Z LDZ | | | | | | | | | XYINT = SGLNQ2( FUN1 , 0.0 , 2.0 , 10 , -2.0 , -1.0 , 5 , Z , 10 ) . . .
XYINT = -6.1108
This example shows how to reverse the order of integration of the variables x and y. It computes the integral of the function f given by:
over the intervals (0.0, 1.0) for the variable x and (0.0, 20.0) for the variable y, using the Gauss-Legendre method with 5 points in the x variable and 48 points in the y variable. Because the order of integration does not matter to the approximation:
the variable y, having more points, is the first variable integrated (performing the integration shown on the right.) This allows the SGLNQ2 subroutine and the FUN1 evaluation to achieve optimal performance. Therefore, the x and y variables correspond to T and S in the FUN2 subroutine. Also, the x and y variables correspond to the C, D, N2 and A, B, N1 sets of arguments, respectively, for SGLNQ2.
The user-supplied subroutine FUN2, which evaluates the integrand function, is coded in Fortran as follows:
SUBROUTINE FUN2 (S,N1,T,N2,Z,LDZ) INTEGER*4 N1,N2,LDZ REAL*4 S(*),T(*),Z(LDZ,*) DO 1 J=1,N2 DO 2 I=1,N1 2 Z(I,J)=COS(T(J))*SIN(S(I)) 1 CONTINUE RETURN END
EXTERNAL FUN2. . . . SUBF A B N1 C D N2 Z LDZ | | | | | | | | | YXINT = SGLNQ2( FUN2 , 0.0 , 20.0 , 48 , 0.0 , 1.0 , 5 , Z , 48 ) . . .
YXINT = 0.4981