LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
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*> \htmlonly
9*> Download DPOTRF + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dpotrf.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dpotrf.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dpotrf.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE DPOTRF( UPLO, N, A, LDA, INFO )
22*
23* .. Scalar Arguments ..
24* CHARACTER 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*> DPOTRF computes the Cholesky factorization of a real symmetric
38*> positive definite matrix A.
39*>
40*> The factorization has the form
41*> A = U**T * U, if UPLO = 'U', or
42*> A = L * L**T, if UPLO = 'L',
43*> where U is an upper triangular matrix and L is lower triangular.
44*>
45*> This is the block version of the algorithm, calling Level 3 BLAS.
46*> \endverbatim
47*
48* Arguments:
49* ==========
50*
51*> \param[in] UPLO
52*> \verbatim
53*> UPLO is CHARACTER*1
54*> = 'U': Upper triangle of A is stored;
55*> = 'L': Lower triangle of A is stored.
56*> \endverbatim
57*>
58*> \param[in] N
59*> \verbatim
60*> N is INTEGER
61*> The order of the matrix A. N >= 0.
62*> \endverbatim
63*>
64*> \param[in,out] A
65*> \verbatim
66*> A is DOUBLE PRECISION array, dimension (LDA,N)
67*> On entry, the symmetric matrix A. If UPLO = 'U', the leading
68*> N-by-N upper triangular part of A contains the upper
69*> triangular part of the matrix A, and the strictly lower
70*> triangular part of A is not referenced. If UPLO = 'L', the
71*> leading N-by-N lower triangular part of A contains the lower
72*> triangular part of the matrix A, and the strictly upper
73*> triangular part of A is not referenced.
74*>
75*> On exit, if INFO = 0, the factor U or L from the Cholesky
76*> factorization A = U**T*U or A = L*L**T.
77*> \endverbatim
78*>
79*> \param[in] LDA
80*> \verbatim
81*> LDA is INTEGER
82*> The leading dimension of the array A. LDA >= max(1,N).
83*> \endverbatim
84*>
85*> \param[out] INFO
86*> \verbatim
87*> INFO is INTEGER
88*> = 0: successful exit
89*> < 0: if INFO = -i, the i-th argument had an illegal value
90*> > 0: if INFO = i, the leading principal minor of order i
91*> is not positive, and the factorization could not be
92*> completed.
93*> \endverbatim
94*
95* Authors:
96* ========
97*
98*> \author Univ. of Tennessee
99*> \author Univ. of California Berkeley
100*> \author Univ. of Colorado Denver
101*> \author NAG Ltd.
102*
103*> \ingroup potrf
104*
105* =====================================================================
106 SUBROUTINE dpotrf( UPLO, N, A, LDA, INFO )
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 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
124 parameter( one = 1.0d+0 )
125* ..
126* .. Local Scalars ..
127 LOGICAL UPPER
128 INTEGER J, JB, NB
129* ..
130* .. External Functions ..
131 LOGICAL LSAME
132 INTEGER ILAENV
133 EXTERNAL lsame, ilaenv
134* ..
135* .. External Subroutines ..
136 EXTERNAL dgemm, dpotrf2, dsyrk, dtrsm, 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 IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
148 info = -1
149 ELSE IF( n.LT.0 ) THEN
150 info = -2
151 ELSE IF( lda.LT.max( 1, n ) ) THEN
152 info = -4
153 END IF
154 IF( info.NE.0 ) THEN
155 CALL xerbla( 'DPOTRF', -info )
156 RETURN
157 END IF
158*
159* Quick return if possible
160*
161 IF( n.EQ.0 )
162 $ RETURN
163*
164* Determine the block size for this environment.
165*
166 nb = ilaenv( 1, 'DPOTRF', uplo, n, -1, -1, -1 )
167 IF( nb.LE.1 .OR. nb.GE.n ) THEN
168*
169* Use unblocked code.
170*
171 CALL dpotrf2( uplo, n, a, lda, info )
172 ELSE
173*
174* Use blocked code.
175*
176 IF( upper ) THEN
177*
178* Compute the Cholesky factorization A = U**T*U.
179*
180 DO 10 j = 1, n, nb
181*
182* Update and factorize the current diagonal block and test
183* for non-positive-definiteness.
184*
185 jb = min( nb, n-j+1 )
186 CALL dsyrk( 'Upper', 'Transpose', jb, j-1, -one,
187 $ a( 1, j ), lda, one, a( j, j ), lda )
188 CALL dpotrf2( 'Upper', jb, a( j, j ), lda, info )
189 IF( info.NE.0 )
190 $ GO TO 30
191 IF( j+jb.LE.n ) THEN
192*
193* Compute the current block row.
194*
195 CALL dgemm( 'Transpose', 'No transpose', jb, 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', 'Non-unit',
199 $ jb, n-j-jb+1, one, a( j, j ), lda,
200 $ a( j, j+jb ), lda )
201 END IF
202 10 CONTINUE
203*
204 ELSE
205*
206* Compute the Cholesky factorization A = L*L**T.
207*
208 DO 20 j = 1, n, nb
209*
210* Update and factorize the current diagonal block and test
211* for non-positive-definiteness.
212*
213 jb = min( nb, n-j+1 )
214 CALL dsyrk( 'Lower', 'No transpose', jb, j-1, -one,
215 $ a( j, 1 ), lda, one, a( j, j ), lda )
216 CALL dpotrf2( 'Lower', jb, a( j, j ), lda, info )
217 IF( info.NE.0 )
218 $ GO TO 30
219 IF( j+jb.LE.n ) THEN
220*
221* Compute the current block column.
222*
223 CALL dgemm( 'No transpose', 'Transpose', n-j-jb+1, jb,
224 $ j-1, -one, a( j+jb, 1 ), lda, a( j, 1 ),
225 $ lda, one, a( j+jb, j ), lda )
226 CALL dtrsm( 'Right', 'Lower', 'Transpose', 'Non-unit',
227 $ n-j-jb+1, jb, one, a( j, j ), lda,
228 $ a( j+jb, j ), lda )
229 END IF
230 20 CONTINUE
231 END IF
232 END IF
233 GO TO 40
234*
235 30 CONTINUE
236 info = info + j - 1
237*
238 40 CONTINUE
239 RETURN
240*
241* End of DPOTRF
242*
243 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:107
subroutine dtrsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
DTRSM
Definition dtrsm.f:181