STRSM(3S)STRSM(3S)NAME
STRSM, DTRSM, CTRSM, ZTRSM - Solves a real or complex triangular system
of equations with multiple right-hand sides
SYNOPSIS
Single precision
Fortran:
CALL STRSM (side, uplo, transa, diag, m, n, alpha, a, lda, b,
ldb)
C/C++:
#include <scsl_blas.h>
void strsm (char *side, char *uplo, char *transA, char *diag,
int m, int n, float alpha, float *a, int lda, float *b, int
ldb);
Double precision
Fortran:
CALL DTRSM (side, uplo, transa, diag, m, n, alpha, a, lda, b,
ldb)
C/C++:
#include <scsl_blas.h>
void dtrsm (char *side, char *uplo, char *transA, char *diag,
int m, int n, double alpha, double *a, int lda, double *b, int
ldb);
Single precision complex
Fortran:
CALL CTRSM (side, uplo, transa, diag, m, n, alpha, a, lda, b,
ldb)
C/C++:
#include <scsl_blas.h>
void ctrsm (char *side, char *uplo, char *transA, char *diag,
int m, int n, scsl_complex *alpha, scsl_complex *a, int lda,
scsl_complex *b, int ldb);
C++ STL:
#include <complex.h>
#include <scsl_blas.h>
void ctrsm (char *side, char *uplo, char *transA, char *diag,
int m, int n, complex<float> *alpha, complex<float> *a, int
lda, complex<float> *b, int ldb);
Double precision complex
Page 1
STRSM(3S)STRSM(3S)
Fortran:
CALL ZTRSM (side, uplo, transa, diag, m, n, alpha, a, lda, b,
ldb)
C/C++:
#include <scsl_blas.h>
void ztrsm (char *side, char *uplo, char *transA, char *diag,
int m, int n, scsl_zomplex *alpha, scsl_zomplex *a, int lda,
scsl_zomplex *b, int ldb);
C++ STL:
#include <complex.h>
#include <scsl_blas.h>
void ztrsm (char *side, char *uplo, char *transA, char *diag,
int m, int n, complex<double> *alpha, complex<double> *a, int
lda, complex<double> *b, int ldb);
IMPLEMENTATION
These routines are part of the SCSL Scientific Library and can be loaded
using either the -lscs or the -lscs_mp option. The -lscs_mp option
directs the linker to use the multi-processor version of the library.
When linking to SCSL with -lscs or -lscs_mp, the default integer size is
4 bytes (32 bits). Another version of SCSL is available in which integers
are 8 bytes (64 bits). This version allows the user access to larger
memory sizes and helps when porting legacy Cray codes. It can be loaded
by using the -lscs_i8 option or the -lscs_i8_mp option. A program may use
only one of the two versions; 4-byte integer and 8-byte integer library
calls cannot be mixed.
The C and C++ prototypes shown above are appropriate for the 4-byte
integer version of SCSL. When using the 8-byte integer version, the
variables of type int become long long and the <scsl_blas_i8.h> header
file should be included.
DESCRIPTION
STRSM and DTRSM solve a real triangular syste m of equations with
multiple right-hand sides.
CTRSM and ZTRSM solve a complex triangular system of equations with
multiple right-hand sides.
These routines solve one of the following matrix equations, using the
operation associated with each:
op(A) X = alpha B
B <- alpha op(A-1)B
or
Page 2
STRSM(3S)STRSM(3S)
X op(A) = alpha B
B <- alpha B op(A-1)
where
* alpha is a scalar
* X and B are m-by-n matrices
* A is either a unit or nonunit upper or lower triangular matrix
* A-1 is the inverse of A
* op(A) is one of the following:
op(A) = A
op(A) = AT
op(A) = AH (CTRSM and ZTRSM only)
where AT is the transpose of A, and AH is the conjugate transpose of A.
See the NOTES section of this man page for information about the
interpretation of the data types described in the following arguments.
These routines have the following arguments:
side Character. (input)
Specifies whether op(A) appears on the left or right of X, as
follows:
side = 'L' or 'l': op(A)*X = alpha*B
side = 'R' or 'r': X*op(A) = alpha*B
For C/C++, a pointer to this character is passed.
uplo Character. (input)
Specifies whether matrix A is an upper or lower triangular
matrix, as follows:
uplo = 'U' or 'u': A is an upper triangular matrix.
uplo = 'L' or 'l': A is a lower triangular matrix.
For C/C++, a pointer to this character is passed.
transa Character*1. (input)
Specifies the form of op(A) to be used in the matrix
multiplication, as follows:
Page 3
STRSM(3S)STRSM(3S)
transa = 'N' or 'n': op(A) = A
transa = 'T' or 't': op(A) = AT
transa = 'C' or 'c': op(A) = AT (STRSM, DTRSM),
or op(A) = AH (CTRSM, ZTRSM)
For C/C++, a pointer to this character is passed.
diag Character. (input)
Specifies whether A is unit triangular, as follows:
diag = 'U' or 'u': A is assumed to be unit triangular.
diag = 'N' or 'n': A is not assumed to be unit triangular.
For C/C++, a pointer to this character is passed.
m Integer. (input)
Specifies the number of rows in B. m must be >= 0.
n Integer. (input)
Specifies the number of columns in B. n must be >= 0.
alpha Scalar factor. (input)
STRSM: Single precision.
DTRSM: Double precision.
CTRSM: Single precision complex.
ZTRSM: Double precision complex.
When alpha is 0, a is not referenced, and b need not be set
before entry.
For C/C++, a pointer to this scalar is passed when alpha is
complex; otherwise, alpha is passed by value.
a Array of dimension (lda,k). (input)
STRSM: Single precision array.
DTRSM: Double precision array.
CTRSM: Single precision complex array.
ZTRSM: Double precision complex array.
When side = 'L' or 'l', k is m; when side = 'R' or 'r', it is
n.
Contains the matrix A.
Before entry with uplo = 'U' or 'u', the leading k-by-k upper
triangular part of array a must contain the upper triangular
matrix. The strictly lower triangular part of a is not
referenced.
Before entry with uplo = 'L' or 'l', the leading k-by-k lower
triangular part of array a must contain the lower triangular
matrix. The strictly upper triangular part of a is not
Page 4
STRSM(3S)STRSM(3S)
referenced.
When diag = 'U' or 'u', the diagonal elements of a are not
referenced, but they are assumed to be unity.
lda Integer. (input)
Specifies the first dimension of a as declared in the calling
program.
When side = 'L' or 'l', lda >= MAX(1,m).
When side = 'R' or 'r', lda >= MAX(1,n).
b Array of dimension (ldb,n). (input)
STRSM: Real array.
DTRSM: Double precision array.
CTRSM: Complex array.
ZTRSM: Double complex array.
Contains the matrix B.
Before entry, the leading m-by-n part of array b must contain
the right-hand side matrix B. On exit, the solution matrix X
overwrites array b.
ldb Integer. (input)
Specifies the first dimension of b as declared in the calling
program. ldb >= MAX(1,m).
NOTES
These routines are Level 3 Basic Linear Algebra Subprograms (Level 3
BLAS).
Data Types
The following data types are described in this documentation:
Term Used Data type
Fortran:
Array dimensioned n x(n)
Array of dimensions (m,n) x(m,n)
Character CHARACTER
Integer INTEGER (INTEGER*8 for -lscs_i8[_mp])
Single precision REAL
Double precision DOUBLE PRECISION
Page 5
STRSM(3S)STRSM(3S)
Single precision complex COMPLEX
Double precision complex DOUBLE COMPLEX
C/C++:
Array dimensioned n x[n]
Array of dimensions (m,n) x[m*n]
Character char
Integer int (long long for -lscs_i8[_mp])
Single precision float
Double precision double
Single precision complex scsl_complex
Double precision complex scsl_zomplex
C++ STL:
Array dimensioned n x[n]
Array of dimensions (m,n) x[m*n]
Character char
Integer int (long long for -lscs_i8[_mp])
Single precision float
Double precision double
Single precision complex complex<float>
Double precision complex complex<double>
Note that you can explicitly declare multidimensional C/C++ arrays
provided that the array dimensions are swapped with respect to the
Fortran declaration (e.g., x[n][m] in C/C++ versus x(m,n) in Fortran).
To avoid a compiler type mismatch error in C++ (or a compiler warning
message in C), however, the array should be cast to a pointer of the
appropriate type when passed as an argument to a SCSL routine.
SEE ALSOINTRO_SCSL(3S), INTRO_BLAS3(3S)
Page 6
STRSM(3S)STRSM(3S)INTRO_CBLAS(3S) for information about using the C interface to Fortran 77
Basic Linear Algebra Subprograms (legacy BLAS) set forth by the Basic
Linear Algebra Subprograms Technical Forum.
Page 7