129 SUBROUTINE dsytrs_aa( UPLO, N, NRHS, A, LDA, IPIV, B, LDB,
130 $ WORK, LWORK, INFO )
140 INTEGER N, NRHS, LDA, LDB, LWORK, INFO
144 DOUBLE PRECISION A( LDA, * ), B( LDB, * ), WORK( * )
150 parameter( one = 1.0d+0 )
153 LOGICAL LQUERY, UPPER
154 INTEGER K, KP, LWKOPT
169 upper = lsame( uplo,
'U' )
170 lquery = ( lwork.EQ.-1 )
171 IF( .NOT.upper .AND. .NOT.lsame( uplo,
'L' ) )
THEN
173 ELSE IF( n.LT.0 )
THEN
175 ELSE IF( nrhs.LT.0 )
THEN
177 ELSE IF( lda.LT.max( 1, n ) )
THEN
179 ELSE IF( ldb.LT.max( 1, n ) )
THEN
181 ELSE IF( lwork.LT.max( 1, 3*n-2 ) .AND. .NOT.lquery )
THEN
185 CALL xerbla(
'DSYTRS_AA', -info )
187 ELSE IF( lquery )
THEN
195 IF( n.EQ.0 .OR. nrhs.EQ.0 )
211 $
CALL dswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
216 CALL dtrsm(
'L',
'U',
'T',
'U', n-1, nrhs, one, a( 1, 2 ),
217 $ lda, b( 2, 1 ), ldb)
224 CALL dlacpy(
'F', 1, n, a( 1, 1 ), lda+1, work( n ), 1)
226 CALL dlacpy(
'F', 1, n-1, a( 1, 2 ), lda+1, work( 1 ), 1 )
227 CALL dlacpy(
'F', 1, n-1, a( 1, 2 ), lda+1, work( 2*n ), 1 )
229 CALL dgtsv( n, nrhs, work( 1 ), work( n ), work( 2*n ), b, ldb,
238 CALL dtrsm(
'L',
'U',
'N',
'U', n-1, nrhs, one, a( 1, 2 ),
239 $ lda, b( 2, 1 ), ldb)
246 $
CALL dswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
263 $
CALL dswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
268 CALL dtrsm(
'L',
'L',
'N',
'U', n-1, nrhs, one, a( 2, 1 ),
269 $ lda, b( 2, 1 ), ldb)
276 CALL dlacpy(
'F', 1, n, a(1, 1), lda+1, work(n), 1)
278 CALL dlacpy(
'F', 1, n-1, a( 2, 1 ), lda+1, work( 1 ), 1 )
279 CALL dlacpy(
'F', 1, n-1, a( 2, 1 ), lda+1, work( 2*n ), 1 )
281 CALL dgtsv( n, nrhs, work( 1 ), work(n), work( 2*n ), b, ldb,
290 CALL dtrsm(
'L',
'L',
'T',
'U', n-1, nrhs, one, a( 2, 1 ),
291 $ lda, b( 2, 1 ), ldb)
298 $
CALL dswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
subroutine dlacpy(UPLO, M, N, A, LDA, B, LDB)
DLACPY copies all or part of one two-dimensional array to another.
subroutine xerbla(SRNAME, INFO)
XERBLA
subroutine dswap(N, DX, INCX, DY, INCY)
DSWAP
subroutine dtrsm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
DTRSM
subroutine dgtsv(N, NRHS, DL, D, DU, B, LDB, INFO)
DGTSV computes the solution to system of linear equations A * X = B for GT matrices
subroutine dsytrs_aa(UPLO, N, NRHS, A, LDA, IPIV, B, LDB, WORK, LWORK, INFO)
DSYTRS_AA