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