LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
dpotrf.f
Go to the documentation of this file.
1*> \brief \b DPOTRF
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download DPOTRF + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dpotrf.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dpotrf.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dpotrf.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE DPOTRF( UPLO, N, A, LDA, INFO )
20*
21* .. Scalar Arguments ..
22* CHARACTER UPLO
23* INTEGER INFO, LDA, N
24* ..
25* .. Array Arguments ..
26* DOUBLE PRECISION A( LDA, * )
27* ..
28*
29*
30*> \par Purpose:
31* =============
32*>
33*> \verbatim
34*>
35*> DPOTRF computes the Cholesky factorization of a real symmetric
36*> positive definite matrix A.
37*>
38*> The factorization has the form
39*> A = U**T * U, if UPLO = 'U', or
40*> A = L * L**T, if UPLO = 'L',
41*> where U is an upper triangular matrix and L is lower triangular.
42*>
43*> This is the block version of the algorithm, calling Level 3 BLAS.
44*> \endverbatim
45*
46* Arguments:
47* ==========
48*
49*> \param[in] UPLO
50*> \verbatim
51*> UPLO is CHARACTER*1
52*> = 'U': Upper triangle of A is stored;
53*> = 'L': Lower triangle of A is stored.
54*> \endverbatim
55*>
56*> \param[in] N
57*> \verbatim
58*> N is INTEGER
59*> The order of the matrix A. N >= 0.
60*> \endverbatim
61*>
62*> \param[in,out] A
63*> \verbatim
64*> A is DOUBLE PRECISION array, dimension (LDA,N)
65*> On entry, the symmetric matrix A. If UPLO = 'U', the leading
66*> N-by-N upper triangular part of A contains the upper
67*> triangular part of the matrix A, and the strictly lower
68*> triangular part of A is not referenced. If UPLO = 'L', the
69*> leading N-by-N lower triangular part of A contains the lower
70*> triangular part of the matrix A, and the strictly upper
71*> triangular part of A is not referenced.
72*>
73*> On exit, if INFO = 0, the factor U or L from the Cholesky
74*> factorization A = U**T*U or A = L*L**T.
75*> \endverbatim
76*>
77*> \param[in] LDA
78*> \verbatim
79*> LDA is INTEGER
80*> The leading dimension of the array A. LDA >= max(1,N).
81*> \endverbatim
82*>
83*> \param[out] INFO
84*> \verbatim
85*> INFO is INTEGER
86*> = 0: successful exit
87*> < 0: if INFO = -i, the i-th argument had an illegal value
88*> > 0: if INFO = i, the leading principal minor of order i
89*> is not positive, and the factorization could not be
90*> completed.
91*> \endverbatim
92*
93* Authors:
94* ========
95*
96*> \author Univ. of Tennessee
97*> \author Univ. of California Berkeley
98*> \author Univ. of Colorado Denver
99*> \author NAG Ltd.
100*
101*> \ingroup potrf
102*
103* =====================================================================
104 SUBROUTINE dpotrf( UPLO, N, A, LDA, INFO )
105*
106* -- LAPACK computational routine --
107* -- LAPACK is a software package provided by Univ. of Tennessee, --
108* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
109*
110* .. Scalar Arguments ..
111 CHARACTER UPLO
112 INTEGER INFO, LDA, N
113* ..
114* .. Array Arguments ..
115 DOUBLE PRECISION A( LDA, * )
116* ..
117*
118* =====================================================================
119*
120* .. Parameters ..
121 DOUBLE PRECISION ONE
122 parameter( one = 1.0d+0 )
123* ..
124* .. Local Scalars ..
125 LOGICAL UPPER
126 INTEGER J, JB, NB
127* ..
128* .. External Functions ..
129 LOGICAL LSAME
130 INTEGER ILAENV
131 EXTERNAL lsame, ilaenv
132* ..
133* .. External Subroutines ..
134 EXTERNAL dgemm, dpotrf2, dsyrk, dtrsm,
135 $ xerbla
136* ..
137* .. Intrinsic Functions ..
138 INTRINSIC max, min
139* ..
140* .. Executable Statements ..
141*
142* Test the input parameters.
143*
144 info = 0
145 upper = lsame( uplo, 'U' )
146 IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
147 info = -1
148 ELSE IF( n.LT.0 ) THEN
149 info = -2
150 ELSE IF( lda.LT.max( 1, n ) ) THEN
151 info = -4
152 END IF
153 IF( info.NE.0 ) THEN
154 CALL xerbla( 'DPOTRF', -info )
155 RETURN
156 END IF
157*
158* Quick return if possible
159*
160 IF( n.EQ.0 )
161 $ RETURN
162*
163* Determine the block size for this environment.
164*
165 nb = ilaenv( 1, 'DPOTRF', uplo, n, -1, -1, -1 )
166 IF( nb.LE.1 .OR. nb.GE.n ) THEN
167*
168* Use unblocked code.
169*
170 CALL dpotrf2( uplo, n, a, lda, info )
171 ELSE
172*
173* Use blocked code.
174*
175 IF( upper ) THEN
176*
177* Compute the Cholesky factorization A = U**T*U.
178*
179 DO 10 j = 1, n, nb
180*
181* Update and factorize the current diagonal block and test
182* for non-positive-definiteness.
183*
184 jb = min( nb, n-j+1 )
185 CALL dsyrk( 'Upper', 'Transpose', jb, j-1, -one,
186 $ a( 1, j ), lda, one, a( j, j ), lda )
187 CALL dpotrf2( 'Upper', jb, a( j, j ), lda, info )
188 IF( info.NE.0 )
189 $ GO TO 30
190 IF( j+jb.LE.n ) THEN
191*
192* Compute the current block row.
193*
194 CALL dgemm( 'Transpose', 'No transpose', jb,
195 $ n-j-jb+1,
196 $ j-1, -one, a( 1, j ), lda, a( 1, j+jb ),
197 $ lda, one, a( j, j+jb ), lda )
198 CALL dtrsm( 'Left', 'Upper', 'Transpose',
199 $ 'Non-unit',
200 $ jb, n-j-jb+1, one, a( j, j ), lda,
201 $ a( j, j+jb ), lda )
202 END IF
203 10 CONTINUE
204*
205 ELSE
206*
207* Compute the Cholesky factorization A = L*L**T.
208*
209 DO 20 j = 1, n, nb
210*
211* Update and factorize the current diagonal block and test
212* for non-positive-definiteness.
213*
214 jb = min( nb, n-j+1 )
215 CALL dsyrk( 'Lower', 'No transpose', jb, j-1, -one,
216 $ a( j, 1 ), lda, one, a( j, j ), lda )
217 CALL dpotrf2( 'Lower', jb, a( j, j ), lda, info )
218 IF( info.NE.0 )
219 $ GO TO 30
220 IF( j+jb.LE.n ) THEN
221*
222* Compute the current block column.
223*
224 CALL dgemm( 'No transpose', 'Transpose', n-j-jb+1,
225 $ jb,
226 $ j-1, -one, a( j+jb, 1 ), lda, a( j, 1 ),
227 $ lda, one, a( j+jb, j ), lda )
228 CALL dtrsm( 'Right', 'Lower', 'Transpose',
229 $ 'Non-unit',
230 $ n-j-jb+1, jb, one, a( j, j ), lda,
231 $ a( j+jb, j ), lda )
232 END IF
233 20 CONTINUE
234 END IF
235 END IF
236 GO TO 40
237*
238 30 CONTINUE
239 info = info + j - 1
240*
241 40 CONTINUE
242 RETURN
243*
244* End of DPOTRF
245*
246 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
DGEMM
Definition dgemm.f:188
subroutine dsyrk(uplo, trans, n, k, alpha, a, lda, beta, c, ldc)
DSYRK
Definition dsyrk.f:169
recursive subroutine dpotrf2(uplo, n, a, lda, info)
DPOTRF2
Definition dpotrf2.f:106
subroutine dpotrf(uplo, n, a, lda, info)
DPOTRF
Definition dpotrf.f:105
subroutine dtrsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
DTRSM
Definition dtrsm.f:181