LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ dtrtri()

subroutine dtrtri ( character uplo,
character diag,
integer n,
double precision, dimension( lda, * ) a,
integer lda,
integer info )

DTRTRI

Download DTRTRI + dependencies [TGZ] [ZIP] [TXT]

Purpose:
!>
!> DTRTRI computes the inverse of a real upper or lower triangular
!> matrix A.
!>
!> This is the Level 3 BLAS version of the algorithm.
!> 
Parameters
[in]UPLO
!>          UPLO is CHARACTER*1
!>          = 'U':  A is upper triangular;
!>          = 'L':  A is lower triangular.
!> 
[in]DIAG
!>          DIAG is CHARACTER*1
!>          = 'N':  A is non-unit triangular;
!>          = 'U':  A is unit triangular.
!> 
[in]N
!>          N is INTEGER
!>          The order of the matrix A.  N >= 0.
!> 
[in,out]A
!>          A is DOUBLE PRECISION array, dimension (LDA,N)
!>          On entry, the triangular matrix A.  If UPLO = 'U', the
!>          leading N-by-N upper triangular part of the array A contains
!>          the upper triangular matrix, and the strictly lower
!>          triangular part of A is not referenced.  If UPLO = 'L', the
!>          leading N-by-N lower triangular part of the array A contains
!>          the lower triangular matrix, and the strictly upper
!>          triangular part of A is not referenced.  If DIAG = 'U', the
!>          diagonal elements of A are also not referenced and are
!>          assumed to be 1.
!>          On exit, the (triangular) inverse of the original matrix, in
!>          the same storage format.
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A.  LDA >= max(1,N).
!> 
[out]INFO
!>          INFO is INTEGER
!>          = 0: successful exit
!>          < 0: if INFO = -i, the i-th argument had an illegal value
!>          > 0: if INFO = i, A(i,i) is exactly zero.  The triangular
!>               matrix is singular and its inverse can not be computed.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 106 of file dtrtri.f.

107*
108* -- LAPACK computational routine --
109* -- LAPACK is a software package provided by Univ. of Tennessee, --
110* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
111*
112* .. Scalar Arguments ..
113 CHARACTER DIAG, UPLO
114 INTEGER INFO, LDA, N
115* ..
116* .. Array Arguments ..
117 DOUBLE PRECISION A( LDA, * )
118* ..
119*
120* =====================================================================
121*
122* .. Parameters ..
123 DOUBLE PRECISION ONE, ZERO
124 parameter( one = 1.0d+0, zero = 0.0d+0 )
125* ..
126* .. Local Scalars ..
127 LOGICAL NOUNIT, UPPER
128 INTEGER J, JB, NB, NN
129* ..
130* .. External Functions ..
131 LOGICAL LSAME
132 INTEGER ILAENV
133 EXTERNAL lsame, ilaenv
134* ..
135* .. External Subroutines ..
136 EXTERNAL dtrmm, dtrsm, dtrti2, xerbla
137* ..
138* .. Intrinsic Functions ..
139 INTRINSIC max, min
140* ..
141* .. Executable Statements ..
142*
143* Test the input parameters.
144*
145 info = 0
146 upper = lsame( uplo, 'U' )
147 nounit = lsame( diag, 'N' )
148 IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
149 info = -1
150 ELSE IF( .NOT.nounit .AND. .NOT.lsame( diag, 'U' ) ) THEN
151 info = -2
152 ELSE IF( n.LT.0 ) THEN
153 info = -3
154 ELSE IF( lda.LT.max( 1, n ) ) THEN
155 info = -5
156 END IF
157 IF( info.NE.0 ) THEN
158 CALL xerbla( 'DTRTRI', -info )
159 RETURN
160 END IF
161*
162* Quick return if possible
163*
164 IF( n.EQ.0 )
165 $ RETURN
166*
167* Check for singularity if non-unit.
168*
169 IF( nounit ) THEN
170 DO 10 info = 1, n
171 IF( a( info, info ).EQ.zero )
172 $ RETURN
173 10 CONTINUE
174 info = 0
175 END IF
176*
177* Determine the block size for this environment.
178*
179 nb = ilaenv( 1, 'DTRTRI', uplo // diag, n, -1, -1, -1 )
180 IF( nb.LE.1 .OR. nb.GE.n ) THEN
181*
182* Use unblocked code
183*
184 CALL dtrti2( uplo, diag, n, a, lda, info )
185 ELSE
186*
187* Use blocked code
188*
189 IF( upper ) THEN
190*
191* Compute inverse of upper triangular matrix
192*
193 DO 20 j = 1, n, nb
194 jb = min( nb, n-j+1 )
195*
196* Compute rows 1:j-1 of current block column
197*
198 CALL dtrmm( 'Left', 'Upper', 'No transpose', diag,
199 $ j-1,
200 $ jb, one, a, lda, a( 1, j ), lda )
201 CALL dtrsm( 'Right', 'Upper', 'No transpose', diag,
202 $ j-1,
203 $ jb, -one, a( j, j ), lda, a( 1, j ), lda )
204*
205* Compute inverse of current diagonal block
206*
207 CALL dtrti2( 'Upper', diag, jb, a( j, j ), lda, info )
208 20 CONTINUE
209 ELSE
210*
211* Compute inverse of lower triangular matrix
212*
213 nn = ( ( n-1 ) / nb )*nb + 1
214 DO 30 j = nn, 1, -nb
215 jb = min( nb, n-j+1 )
216 IF( j+jb.LE.n ) THEN
217*
218* Compute rows j+jb:n of current block column
219*
220 CALL dtrmm( 'Left', 'Lower', 'No transpose', diag,
221 $ n-j-jb+1, jb, one, a( j+jb, j+jb ), lda,
222 $ a( j+jb, j ), lda )
223 CALL dtrsm( 'Right', 'Lower', 'No transpose', diag,
224 $ n-j-jb+1, jb, -one, a( j, j ), lda,
225 $ a( j+jb, j ), lda )
226 END IF
227*
228* Compute inverse of current diagonal block
229*
230 CALL dtrti2( 'Lower', diag, jb, a( j, j ), lda, info )
231 30 CONTINUE
232 END IF
233 END IF
234*
235 RETURN
236*
237* End of DTRTRI
238*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:160
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine dtrmm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
DTRMM
Definition dtrmm.f:177
subroutine dtrsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
DTRSM
Definition dtrsm.f:181
subroutine dtrti2(uplo, diag, n, a, lda, info)
DTRTI2 computes the inverse of a triangular matrix (unblocked algorithm).
Definition dtrti2.f:108
Here is the call graph for this function:
Here is the caller graph for this function: