S3L_gen_band_factor performs the LU factorization of an n x n general banded array with lower bandwidth bl and upper bandwidth bl. The non-zero diagonals of the array should be stored in an S3L array a of size [2*bl+2*bu+1,n].
In the more general case, a can be a multidimensional array, where axis_r and axis_d denote the array axes whose extents are 2*bl+2*bu+1 and n respectively. The format of the array a is described in the following example:
Consider a 7 x 7 (n=7) banded array with bl = 1, bu = 2. c is the main diagonal, b is the first superdiagonal and a the second. d is the first subdiagonal. The contents of the composite array a used as input to S3L_gen_band_factor should have the following organization:
* * * * * * * * * * * * * * * * * * * * * * * a0 a1 a2 a3 a4 * b0 b1 b2 b3 b4 b5 c0 c1 c2 c3 c4 c5 c6 d0 d1 d2 d3 d4 d5 * |
Note that, items denoted by '*' are not referenced.
If a is two-dimensional, S3L_gen_band_factor is more efficient when axis_r is the first axis, axis_d is the second axis, and array a is block-distributed along the second axis. For C programs, the indices of the first and second axes are 0 and 1, respectively. For Fortran programs, the corresponding indices are 1 and 2.
If a has more than two dimensions, S3L_gen_band_factor is most efficient when axes axis_r and axis_d of a are local (that is, are not distributed).
The C and Fortran syntax for S3L_gen_band_factor are shown below.
#include <s3l/s3l-c.h> #include <s3l/s3l_errno-c.h> int S3L_gen_band_factor(a, bl, bu, factors, axis_r, axis_d) S3L_array_t a int bl int bu int *factors int axis_r int axis_d |
include `s3l/s3l-f.h' include `s3l/s3l_errno-f.h' subroutine S3L_gen_band_factor(a, bl, bu, factors, axis_r, axis_d, ier) integer*4 a integer*4 bl integer*4 bu integer*4 factors integer*4 axis_r integer*4 axis_d integer*4 ier |
a - S3L array handle for a real or complex parallel array of size [1+2*bl+2*bl,N].
bl - Lower bandwidth of a.
bu - Upper bandwidth of a.
axis_r - Specifies the row axis along which factorization will occur.
axis_d - Specifies the column axis along which factorization will occur.
This function uses the following arguments for output:
a - Upon successful completion, S3L_gen_band_factor stores the factorization results in a.
factors - Pointer to an internal structure that holds the factorization.
ier (Fortran only) - When called from a Fortran program, this function returns error status in ier.
On success, S3L_gen_band_factor returns S3L_SUCCESS.
S3L_gen_band_factor performs generic checking of the arrays it accepts as arguments. If an array argument contains an invalid or corrupted value, the function terminates and an error code indicating which value of the array handle was invalid is returned. See Appendix A of this manual for a detailed list of these error codes.
In addition, the following conditions will cause the function to terminate and return the associated error code:
S3L_ERR_ARG_DTYPE - The type of a is not one of: real, double, complex or double complex.
S3L_ERR_INDX_INVALID - bl or bu value is invalid for either of the following reasons:
Less than 0 (C/C++) or less than 1 (F77/F90).
Greater than the extent of a along axis_d.
S3L_ERR_ARG_EXTENTS - The extent of a along axis axis_r is not equal to 2*bl+2*bu+1.
S3L_ERR_ARRTOOSMALL - The extents of a along axis axis_d are such that the block size in a block distribution is less than bu + bl + 1.
S3L_ERR_ARG_AXISNUM - An axis argument is invalid; that is, it is either:
It is less than 0 (C/C++) or less than 1 (F77/F90).
It is greater than the rank of the referenced array.
axis_d is equal to axis_r.
S3L_ERR_BAND_FFAIL - The factorization could not be completed.
../examples/s3l/band/ex_band.c ../examples/s3l/band-f/ex_band.f
S3L_gen_band_solve(3) S3L_gen_band_free_factors(3)
S3L_gen_band_free_factors frees internal memory associated with a banded matrix factorization.
The C and Fortran syntax for S3L_gen_band_free_factors are shown below.
#include <s3l/s3l-c.h> #include <s3l/s3l_errno-c.h> int S3L_gen_band_free_factors(factors) int *factors |
include `s3l/s3l-f.h' include `s3l/s3l_errno-f.h' subroutine S3L_gen_band_free_factors(factors, ier) integer*4 factors integer*4 ier |
factors - Pointer to the internal structure that will be freed.
This function uses the following argument for output:
ier (Fortran only) - When called from a Fortran program, this function returns error status in ier.
On success, S3L_gen_band_free_factors returns S3L_SUCCESS.
The following condition will cause S3L_gen_band_free_factors to terminate and return the associated error code:
S3L_ERR_ARG_SETUP - The factors value is invalid.
../examples/s3l/band/ex_band.c ../examples/s3l/band-f/ex_band.f
S3L_gen_band_solve(3) S3L_gen_band_factor(3)
S3L_gen_band_solve solves a banded system whose factorization has been computed by a prior call to S3L_gen_band_factor.
The factored banded matrix is stored in array a, whose dimensions are 2*bu + 2*bl + 1 x n. The right-hand-side is stored in array b, whose dimensions are n x nrhs.
If a and b have more than two dimensions, axis_r and axis_d refer to those axes of a whose extents are 2*bu + 2*bl + 1 and n, respectively. Likewise, axis_row and axis_col refer to the axes of b with extents n and nrhs.
Two-Dimensional Arrays: If a and b are two-dimensional, S3L_gen_band_solve is more efficient when axis_r = 0, axis_d = 1, array a is block distributed along axis 1, axis_row = 0, axis_col = 1 and array b is block distributed along axis 0.
Note that the values cited in the previous paragraph apply to programs using the C/C++ interface--that is, they assume zero-based array indexing. When S3L_gen_band_solve is called from F77 or F90 applications, these values must be increased by one. Therefore, when a and b are two-dimensional and S3L_gen_band_solve is called by a Fortran program, the solver is more efficient when axis_r = 1, axis_d = 2, array a is block distributed along axis 2, axis_row = 1, axis_col = 2 and array b is block distributed along axis 1.
When a and b are two-dimensional and nrhs is greater than 1, the size of a must be such that n is divisible by the number of processors.
Arrays With More Than Two Dimensions: If a and b have more than two dimensions, S3L_gen_band_solve is more efficient when axes axis_r and axis_d of a and axes axis_row and axis_col are local (not distributed).
The C and Fortran syntax for S3L_gen_band_solve are shown below.
#include <s3l/s3l-c.h> #include <s3l/s3l_errno-c.h> int S3L_gen_band_solve(a, bl, bu, factors, axis_r, axis_d, b, axis_row, axix_col) S3L_array_t a int bl int bu int *factors int axis_r int axis_d S3L_array_t b int axis_row int axis_col |
include `s3l/s3l-f.h' include `s3l/s3l_errno-f.h' subroutine S3L_gen_band_solve(a, bl, bu, factors, axis_r, axis_d, b, axis_row, axis_col, ier) integer*4 a integer*4 bl integer*4 bu integer*4 factors integer*4 axis_r integer*4 axis_d integer*8 b integer*4 axis_row integer*4 axis_col integer*4 ier |
a - S3L array handle for a real or complex parallel array of size [1+2*bl+2*bu,n].
bl - Lower bandwidth of a.
bu - Upper bandwidth of a.
factors - Pointer to an internal structure that holds the factorization results.
axis_r - Specifies the axis of array a whose extent is 1+2*bl+2*bu+1
axis_d - Specifies the axis of array a whose extent is n.
b - S3L array handle containing the right-hand side of the matrix equation ax=b.
This function uses the following argument for output:
b - On output, b is overwritten by the solution to the matrix equation ax=b.
ier (Fortran only) - When called from a Fortran program, this function returns error status in ier.
On success, S3L_gen_band_solve returns S3L_SUCCESS.
S3L_gen_band_solve performs generic checking of the arrays it accepts as arguments. If an array argument contains an invalid or corrupted value, the function terminates and an error code indicating which value of the array handle was invalid is returned. See Appendix A of this manual for a detailed list of these error codes.
In addition, the following conditions will cause the function to terminate and return the associated error code:
S3L_ERR_ARG_DTYPE - The type of a is not one of: real, double, complex or double complex.
S3L_ERR_INDX_INVALID - bl or bu value is invalid for either of the following reasons:
It is less than 0 (C/C++) or less than 1 (F77/F90).
It is greater than the extent of a along axis_d.
S3L_ERR_ARG_EXTENTS - The extent of a along axis axis_r is not equal to 2*bl+2*bu+1.
S3L_ERR_ARRTOOSMALL - The extents of a along axis axis_d are such that the block size in a block distribution is less than bu + bl + 1.
S3L_ERR_ARG_AXISNUM - An axis argument is invalid; that is, it is either:
Less than 0 (C/C++) or less than 1 (F77/F90).
Greater than the rank of the referenced array
axis_d is equal to axis_r.
S3L_ERR_MATCH_RANK - The rank of a is not the same as that of b.
S3L_ERR_ARG_SETUP - The factors value does not correspond to a valid setup.
S3L_ERR_MATCH_EXTENTS - The extents of a along axis_d do not equal the extents of b along axis_row or some of the other extents of a and b do not match.
../examples/s3l/band/ex_band.c ../examples/s3l/band-f/ex_band.f
S3L_gen_band_factor(3) S3L_gen_band_free_factors(3)
S3L_gen_trid_factor factors a tridiagonal matrix, whose diagonal is stored in vector D. The first upper subdiagonal is stored in U, and the first lower subdiagonal in L.
On return, the integer factors contains a pointer to an internal setup structure that holds the factorization. Subsequent calls to S3L_gen_trid_solve use the value in factors to access the factorization results.
The contents of the vectors D, U, and L may be altered. These altered vectors, along with the factors parameter, have to be passed to a subsequent call to S3L_gen_trid_solve to produce the solution to a tridiagonal system.
D, U, and L must have the same extents and type. If they are one-dimensional, all three must be of length n. The first n-1 entries of U contain the elements of the superdiagonal. The last n-1 entries of L contain the elements of the first subdiagonal. The last element of U and the first element of L are not referenced and can be initialized arbitrarily.
If D, U and L have more than one dimension, axis_d is the axis along which the multidimensional arrays are factored. If they are one-dimensional, axis_d must be 0 in C/C++ programs and 1 in F77/F90 programs.
S3L_gen_trid_factor is based on the ScaLAPACK routines pxdttrf, where x is single, double, complex, or double complex. It does no pivoting; consequently, the matrix has to be positive definite for the factorization to be stable.
For one-dimensional arrays, the routine is more efficient when D, U, and L are block distributed. For multiple dimensions, the routine is more efficient when axis_d is a local axis.
The C and Fortran syntax for S3L_gen_trid_factor are shown below.
#include <s3l/s3l-c.h> #include <s3l/s3l_errno-c.h> int S3L_gen_trid_factor(D, U, L, factors, axis_d) S3L_array_t D S3L_array_t U S3L_array_t L int *factors int axis_d |
include `s3l/s3l-f.h' include `s3l/s3l_errno-f.h' subroutine S3L_gen_trid_factor(D, U, L, factors, axis_d, ier) integer*8 D integer*8 U integer*8 L integer*4 factors integer*4 axis_d integer*4 ier |
D - Vector containing the diagonal for the matrix being factored.
U - Vector containing the first upper diagonal for the matrix being factored.
L - Vector containing the first lower diagonal for the matrix being factored.
axis_d - When D, U, and L are one-dimensional, axis_d must be 0 (C/C++ programs) or 1 (F77/F90 programs). For multidimensional arrays, axis_d specifies the axis along which the arrays are factored.
This function uses the following arguments for output:
D - On output, D is overwritten with the partial result of the factorization.
U - On output, U is overwritten with the partial result of the factorization.
L - On output, L is overwritten with the partial result of the factorization.
factors - Upon completion, factors points to the internal data structure containing the factorization results.
ier (Fortran only) - When called from a Fortran program, this function returns error status in ier.
On success, S3L_gen_trid_factor returns S3L_SUCCESS.
S3L_gen_trid_factor performs generic checking of the arrays it accepts as arguments. If an array argument contains an invalid or corrupted value, the function terminates and an error code indicating which value of the array handle was invalid is returned. See Appendix A of this manual for a detailed list of these error codes.
In addition, the following conditions will cause the function to terminate and return the associated error code:
S3L_ERR_MATCH_DTYPE - The arrays are not the same data type.
S3L_ERR_MATCH_RANK - The arrays do not have the same rank.
S3L_ERR_MATCH_EXTENTS - The arrays do not have the same extents.
S3L_ERR_ARG_DTYPE - The array type cannot be operated on by the routine (that is, it is integer or long long).
S3L_ERR_ARRTOOSMALL - The array extent is too small, making the length of the main diagonal less than two times the number of processes.
S3L_ERR_ARG_AXISNUM - An axis argument is invalid; that is, it is either:
Less than 0 (C/C++) or less than 1 (F77/F90).
Greater than the rank of the referenced array.
S3L_ERR_FACTOR_FAIL - The tridiagonal matrix could not be factored for some reason. For example, it might not be diagonally dominant.
../examples/s3l/trid/ex_trid.c ../examples/s3l/trid-f/ex_trid.f
S3L_gen_trid_solve(3) S3L_gen_trid_free_factors(3)
S3L_gen_trid_free_factors frees the internal memory setup that was reserved by a prior call to S3L_gen_trid_factor. The factors argument contains the value returned by the earlier S3L_gen_trid_factor call.
The C and Fortran syntax for S3L_gen_trid_free_factors are shown below.
#include <s3l/s3l-c.h> #include <s3l/s3l_errno-c.h> int S3L_gen_band_free_factors(factors) int *factors |
include `s3l/s3l-f.h' include `s3l/s3l_errno-f.h' subroutine S3L_gen_band_free_factors(factors, ier) integer*4 factors integer*4 ier |
factors - Pointer to the internal structure that will be freed.
This function uses the following argument for output:
ier (Fortran only) - When called from a Fortran program, this function returns error status in ier.
On success, S3L_gen_trid_free_factors returns S3L_SUCCESS.
The following condition will cause S3L_gen_trid_free_factors to terminate and return the associated error code:
S3L_ERR_ARG_SETUP - The factors value is invalid.
../examples/s3l/trid/ex_trid.c ../examples/s3l/trid-f/ex_trid.f
S3L_gen_trid_solve(3) S3L_gen_trid_factor(3)
S3L_gen_trid_solve solves a tridiagonal system that has been previously factored via a call to S3L_gen_trid_factor.
If D, U, and L are of length n, B (the right-hand side of the tridiagonal system) must be of size n x nrhs. If D, U, and L are multidimensional, axis_d is the axis along which the system is solved. The rank of B must be one greater than the rank of D, U, and L.
If the rank of B is greater than 2, row_b and col_b specify the axes whose dimensions are n and nrhs, respectively. The extents of all other axes must be the same as the corresponding axes of D, U, and L.
When computing multiple tridiagonal systems in which only the right-hand-side matrix changes, the factorization routine S3L_gen_trid_factor need only be called once, before the first call to S3l_gen_trid_solve. Then, S3L_gen_trid_solve can be called repeatedly without calling S3L_gen_trid_factor again.
The C and Fortran syntax for S3L_gen_trid_solve are shown below.
#include <s3l/s3l-c.h> #include <s3l/s3l_errno-c.h> int S3L_gen_trid_solve(D, U, L, factors, B, row_b, col_b) S3L_array_t D S3L_array_t U S3L_array_t L int *factors S3L_array_t B int axis_d int axis_d int row_b |
include `s3l/s3l-f.h' include `s3l/s3l_errno-f.h' subroutine S3L_gen_trid_solve(D, U, L, factors, B, axis_d, row_b, col_b, ier) integer*8 D integer*8 U integer*8 L integer*4 factors integer*8 B integer*4 axis_d integer*4 row_b integer*4 col_b integer*4 ier |
D - Vector containing the diagonal for the matrix being factored.
U - Vector containing the first upper subdiagonal for the matrix being factored.
L - Vector containing the first lower subdiagonal for the matrix being factored.
factors - Pointer to an internal structure that holds the factorization results.
B - The right-hand side of the tridiagonal system to be solved.
axis_d - When D, U, and L are one-dimensional, axis_d must be 0 (C/C++ programs) or 1 (F77/F90 programs). For multidimensional arrays, axis_d specifies the axis along which factorization was carried out.
row_b - Indicates the row axis of the right-hand side array, B. The value of row_b depends on the following:
When B is two-dimensional and its sides are n x nrhs, row_b is 0 (C/C++) or 1 (F77/F90).
When B is two-dimensional and its sides are nrhs x n, row_b is 1 (C/C++) or 2 (F77/F90).
When B has more than two dimensions, row_b identifies the side of B with an extent of n. For C/C++ programs, the row_b value is zero-based and for F77/F90 programs, it is one-based.
col_b - Indicates the column axis of the right-hand side array, B that has an extent of nrhs. The value of col_b is determined as follows:
When B is two-dimensional and its sides are n x nrhs, col_b is 1 (C/C++) or 2 (F77/F90).
When B is two-dimensional and its sides are nrhs x n, col_b is 0 (C/C++) or 1 (F77/F90).
When B has more than two dimensions, col_b identifies the side of B with an extent of nhrs. For C/C++ programs, the col_b value is zero-based and for F77/F90 programs, it is one-based.
This function uses the following argument for output:
B - On output, B is overwritten with the solution to the tridiagonal system.
ier (Fortran only) - When called from a Fortran program, this function returns error status in ier.
On success, S3L_gen_trid_solve returns S3L_SUCCESS.
S3L_gen_trid_solve performs generic checking of the arrays it accepts as arguments. If an array argument contains an invalid or corrupted value, the function terminates and an error code indicating which value of the array handle was invalid is returned. See Appendix A of this manual for a detailed list of these error codes.
In addition, the following conditions will cause the function to terminate and return the associated error code.
S3L_ERR_MATCH_DTYPE - The arrays are not the same data type.
S3L_ERR_MATCH_RANK - The arrays do not have compatible rank.
S3L_ERR_MATCH_EXTENTS - The arrays do not have compatible extents.
S3L_ERR_ARG_DTYPE - The array type cannot be operated on by the routine (that is, it is integer or long long).
S3L_ERR_ARRTOOSMALL - The array extent is too small, making the length of the main diagonal less than two times the number of processes.
S3L_ERR_ARG_AXISNUM - An axis argument is invalid; that is, it is either:
Less than 0 (C/C++) or less than 1 (F77/F90).
Greater than the rank of the referenced array.
row_b is equal to col_b.
S3L_ERR_ARG_SETUP - The factors value does not correspond to a valid setup.
../examples/s3l/trid/ex_trid.c ../examples/s3l/trid-f/ex_trid.f
S3L_gen_trid_factor(3) S3L_gen_trid_free_factors(3)