1 *> \brief \b DTRTI2 computes the inverse of a triangular matrix (unblocked algorithm).
3 * =========== DOCUMENTATION ===========
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
9 *> Download DTRTI2 + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dtrti2.f">
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dtrti2.f">
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dtrti2.f">
21 * SUBROUTINE DTRTI2( UPLO, DIAG, N, A, LDA, INFO )
23 * .. Scalar Arguments ..
24 * CHARACTER DIAG, UPLO
25 * INTEGER INFO, LDA, N
27 * .. Array Arguments ..
28 * DOUBLE PRECISION A( LDA, * )
37 *> DTRTI2 computes the inverse of a real upper or lower triangular
40 *> This is the Level 2 BLAS version of the algorithm.
48 *> UPLO is CHARACTER*1
49 *> Specifies whether the matrix A is upper or lower triangular.
50 *> = 'U': Upper triangular
51 *> = 'L': Lower triangular
56 *> DIAG is CHARACTER*1
57 *> Specifies whether or not the matrix A is unit triangular.
58 *> = 'N': Non-unit triangular
59 *> = 'U': Unit triangular
65 *> The order of the matrix A. N >= 0.
70 *> A is DOUBLE PRECISION array, dimension (LDA,N)
71 *> On entry, the triangular matrix A. If UPLO = 'U', the
72 *> leading n by n upper triangular part of the array A contains
73 *> the upper triangular matrix, and the strictly lower
74 *> triangular part of A is not referenced. If UPLO = 'L', the
75 *> leading n by n lower triangular part of the array A contains
76 *> the lower triangular matrix, and the strictly upper
77 *> triangular part of A is not referenced. If DIAG = 'U', the
78 *> diagonal elements of A are also not referenced and are
81 *> On exit, the (triangular) inverse of the original matrix, in
82 *> the same storage format.
88 *> The leading dimension of the array A. LDA >= max(1,N).
94 *> = 0: successful exit
95 *> < 0: if INFO = -k, the k-th argument had an illegal value
101 *> \author Univ. of Tennessee
102 *> \author Univ. of California Berkeley
103 *> \author Univ. of Colorado Denver
106 *> \date September 2012
108 *> \ingroup doubleOTHERcomputational
110 * =====================================================================
111 SUBROUTINE DTRTI2( UPLO, DIAG, N, A, LDA, INFO )
113 * -- LAPACK computational routine (version 3.4.2) --
114 * -- LAPACK is a software package provided by Univ. of Tennessee, --
115 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
118 * .. Scalar Arguments ..
122 * .. Array Arguments ..
123 DOUBLE PRECISION A( LDA, * )
126 * =====================================================================
130 PARAMETER ( ONE = 1.0D+0 )
132 * .. Local Scalars ..
133 LOGICAL NOUNIT, UPPER
137 * .. External Functions ..
141 * .. External Subroutines ..
142 EXTERNAL DSCAL, DTRMV, XERBLA
144 * .. Intrinsic Functions ..
147 * .. Executable Statements ..
149 * Test the input parameters.
152 UPPER = LSAME( UPLO, 'U' )
153 NOUNIT = LSAME( DIAG, 'N' )
154 IF( .NOT.UPPER .AND. .NOT.LSAME( UPLO, 'L' ) ) THEN
156 ELSE IF( .NOT.NOUNIT .AND. .NOT.LSAME( DIAG, 'U' ) ) THEN
158 ELSE IF( N.LT.0 ) THEN
160 ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
164 CALL XERBLA( 'DTRTI2', -INFO )
170 * Compute inverse of upper triangular matrix.
174 A( J, J ) = ONE / A( J, J )
180 * Compute elements 1:j-1 of j-th column.
182 CALL DTRMV( 'Upper', 'No transpose', DIAG, J-1, A, LDA,
184 CALL DSCAL( J-1, AJJ, A( 1, J ), 1 )
188 * Compute inverse of lower triangular matrix.
192 A( J, J ) = ONE / A( J, J )
199 * Compute elements j+1:n of j-th column.
201 CALL DTRMV( 'Lower', 'No transpose', DIAG, N-J,
202 $ A( J+1, J+1 ), LDA, A( J+1, J ), 1 )
203 CALL DSCAL( N-J, AJJ, A( J+1, J ), 1 )