LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dtrtri.f
Go to the documentation of this file.
1*> \brief \b DTRTRI
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> \htmlonly
9*> Download DTRTRI + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dtrtri.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dtrtri.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dtrtri.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE DTRTRI( UPLO, DIAG, N, A, LDA, INFO )
22*
23* .. Scalar Arguments ..
24* CHARACTER DIAG, UPLO
25* INTEGER INFO, LDA, N
26* ..
27* .. Array Arguments ..
28* DOUBLE PRECISION A( LDA, * )
29* ..
30*
31*
32*> \par Purpose:
33* =============
34*>
35*> \verbatim
36*>
37*> DTRTRI computes the inverse of a real upper or lower triangular
38*> matrix A.
39*>
40*> This is the Level 3 BLAS version of the algorithm.
41*> \endverbatim
42*
43* Arguments:
44* ==========
45*
46*> \param[in] UPLO
47*> \verbatim
48*> UPLO is CHARACTER*1
49*> = 'U': A is upper triangular;
50*> = 'L': A is lower triangular.
51*> \endverbatim
52*>
53*> \param[in] DIAG
54*> \verbatim
55*> DIAG is CHARACTER*1
56*> = 'N': A is non-unit triangular;
57*> = 'U': A is unit triangular.
58*> \endverbatim
59*>
60*> \param[in] N
61*> \verbatim
62*> N is INTEGER
63*> The order of the matrix A. N >= 0.
64*> \endverbatim
65*>
66*> \param[in,out] A
67*> \verbatim
68*> A is DOUBLE PRECISION array, dimension (LDA,N)
69*> On entry, the triangular matrix A. If UPLO = 'U', the
70*> leading N-by-N upper triangular part of the array A contains
71*> the upper triangular matrix, and the strictly lower
72*> triangular part of A is not referenced. If UPLO = 'L', the
73*> leading N-by-N lower triangular part of the array A contains
74*> the lower triangular matrix, and the strictly upper
75*> triangular part of A is not referenced. If DIAG = 'U', the
76*> diagonal elements of A are also not referenced and are
77*> assumed to be 1.
78*> On exit, the (triangular) inverse of the original matrix, in
79*> the same storage format.
80*> \endverbatim
81*>
82*> \param[in] LDA
83*> \verbatim
84*> LDA is INTEGER
85*> The leading dimension of the array A. LDA >= max(1,N).
86*> \endverbatim
87*>
88*> \param[out] INFO
89*> \verbatim
90*> INFO is INTEGER
91*> = 0: successful exit
92*> < 0: if INFO = -i, the i-th argument had an illegal value
93*> > 0: if INFO = i, A(i,i) is exactly zero. The triangular
94*> matrix is singular and its inverse can not be computed.
95*> \endverbatim
96*
97* Authors:
98* ========
99*
100*> \author Univ. of Tennessee
101*> \author Univ. of California Berkeley
102*> \author Univ. of Colorado Denver
103*> \author NAG Ltd.
104*
105*> \ingroup trtri
106*
107* =====================================================================
108 SUBROUTINE dtrtri( UPLO, DIAG, N, A, LDA, INFO )
109*
110* -- LAPACK computational routine --
111* -- LAPACK is a software package provided by Univ. of Tennessee, --
112* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
113*
114* .. Scalar Arguments ..
115 CHARACTER DIAG, UPLO
116 INTEGER INFO, LDA, N
117* ..
118* .. Array Arguments ..
119 DOUBLE PRECISION A( LDA, * )
120* ..
121*
122* =====================================================================
123*
124* .. Parameters ..
125 DOUBLE PRECISION ONE, ZERO
126 parameter( one = 1.0d+0, zero = 0.0d+0 )
127* ..
128* .. Local Scalars ..
129 LOGICAL NOUNIT, UPPER
130 INTEGER J, JB, NB, NN
131* ..
132* .. External Functions ..
133 LOGICAL LSAME
134 INTEGER ILAENV
135 EXTERNAL lsame, ilaenv
136* ..
137* .. External Subroutines ..
138 EXTERNAL dtrmm, dtrsm, dtrti2, xerbla
139* ..
140* .. Intrinsic Functions ..
141 INTRINSIC max, min
142* ..
143* .. Executable Statements ..
144*
145* Test the input parameters.
146*
147 info = 0
148 upper = lsame( uplo, 'U' )
149 nounit = lsame( diag, 'N' )
150 IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
151 info = -1
152 ELSE IF( .NOT.nounit .AND. .NOT.lsame( diag, 'U' ) ) THEN
153 info = -2
154 ELSE IF( n.LT.0 ) THEN
155 info = -3
156 ELSE IF( lda.LT.max( 1, n ) ) THEN
157 info = -5
158 END IF
159 IF( info.NE.0 ) THEN
160 CALL xerbla( 'DTRTRI', -info )
161 RETURN
162 END IF
163*
164* Quick return if possible
165*
166 IF( n.EQ.0 )
167 $ RETURN
168*
169* Check for singularity if non-unit.
170*
171 IF( nounit ) THEN
172 DO 10 info = 1, n
173 IF( a( info, info ).EQ.zero )
174 $ RETURN
175 10 CONTINUE
176 info = 0
177 END IF
178*
179* Determine the block size for this environment.
180*
181 nb = ilaenv( 1, 'DTRTRI', uplo // diag, n, -1, -1, -1 )
182 IF( nb.LE.1 .OR. nb.GE.n ) THEN
183*
184* Use unblocked code
185*
186 CALL dtrti2( uplo, diag, n, a, lda, info )
187 ELSE
188*
189* Use blocked code
190*
191 IF( upper ) THEN
192*
193* Compute inverse of upper triangular matrix
194*
195 DO 20 j = 1, n, nb
196 jb = min( nb, n-j+1 )
197*
198* Compute rows 1:j-1 of current block column
199*
200 CALL dtrmm( 'Left', 'Upper', 'No transpose', diag, j-1,
201 $ jb, one, a, lda, a( 1, j ), lda )
202 CALL dtrsm( 'Right', 'Upper', 'No transpose', diag, 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*
239 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
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:110
subroutine dtrtri(uplo, diag, n, a, lda, info)
DTRTRI
Definition dtrtri.f:109