LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zpotrf.f
Go to the documentation of this file.
1*> \brief \b ZPOTRF
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> \htmlonly
9*> Download ZPOTRF + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zpotrf.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zpotrf.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zpotrf.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE ZPOTRF( UPLO, N, A, LDA, INFO )
22*
23* .. Scalar Arguments ..
24* CHARACTER UPLO
25* INTEGER INFO, LDA, N
26* ..
27* .. Array Arguments ..
28* COMPLEX*16 A( LDA, * )
29* ..
30*
31*
32*> \par Purpose:
33* =============
34*>
35*> \verbatim
36*>
37*> ZPOTRF computes the Cholesky factorization of a complex Hermitian
38*> positive definite matrix A.
39*>
40*> The factorization has the form
41*> A = U**H * U, if UPLO = 'U', or
42*> A = L * L**H, 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 COMPLEX*16 array, dimension (LDA,N)
67*> On entry, the Hermitian 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**H *U or A = L*L**H.
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 zpotrf( 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 COMPLEX*16 A( LDA, * )
118* ..
119*
120* =====================================================================
121*
122* .. Parameters ..
123 DOUBLE PRECISION ONE
124 COMPLEX*16 CONE
125 parameter( one = 1.0d+0, cone = ( 1.0d+0, 0.0d+0 ) )
126* ..
127* .. Local Scalars ..
128 LOGICAL UPPER
129 INTEGER J, JB, NB
130* ..
131* .. External Functions ..
132 LOGICAL LSAME
133 INTEGER ILAENV
134 EXTERNAL lsame, ilaenv
135* ..
136* .. External Subroutines ..
137 EXTERNAL xerbla, zgemm, zherk, zpotrf2, ztrsm
138* ..
139* .. Intrinsic Functions ..
140 INTRINSIC max, min
141* ..
142* .. Executable Statements ..
143*
144* Test the input parameters.
145*
146 info = 0
147 upper = lsame( uplo, 'U' )
148 IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
149 info = -1
150 ELSE IF( n.LT.0 ) THEN
151 info = -2
152 ELSE IF( lda.LT.max( 1, n ) ) THEN
153 info = -4
154 END IF
155 IF( info.NE.0 ) THEN
156 CALL xerbla( 'ZPOTRF', -info )
157 RETURN
158 END IF
159*
160* Quick return if possible
161*
162 IF( n.EQ.0 )
163 $ RETURN
164*
165* Determine the block size for this environment.
166*
167 nb = ilaenv( 1, 'ZPOTRF', uplo, n, -1, -1, -1 )
168 IF( nb.LE.1 .OR. nb.GE.n ) THEN
169*
170* Use unblocked code.
171*
172 CALL zpotrf2( uplo, n, a, lda, info )
173 ELSE
174*
175* Use blocked code.
176*
177 IF( upper ) THEN
178*
179* Compute the Cholesky factorization A = U**H *U.
180*
181 DO 10 j = 1, n, nb
182*
183* Update and factorize the current diagonal block and test
184* for non-positive-definiteness.
185*
186 jb = min( nb, n-j+1 )
187 CALL zherk( 'Upper', 'Conjugate transpose', jb, j-1,
188 $ -one, a( 1, j ), lda, one, a( j, j ), lda )
189 CALL zpotrf2( 'Upper', jb, a( j, j ), lda, info )
190 IF( info.NE.0 )
191 $ GO TO 30
192 IF( j+jb.LE.n ) THEN
193*
194* Compute the current block row.
195*
196 CALL zgemm( 'Conjugate transpose', 'No transpose', jb,
197 $ n-j-jb+1, j-1, -cone, a( 1, j ), lda,
198 $ a( 1, j+jb ), lda, cone, a( j, j+jb ),
199 $ lda )
200 CALL ztrsm( 'Left', 'Upper', 'Conjugate transpose',
201 $ 'Non-unit', jb, n-j-jb+1, cone, a( j, j ),
202 $ lda, a( j, j+jb ), lda )
203 END IF
204 10 CONTINUE
205*
206 ELSE
207*
208* Compute the Cholesky factorization A = L*L**H.
209*
210 DO 20 j = 1, n, nb
211*
212* Update and factorize the current diagonal block and test
213* for non-positive-definiteness.
214*
215 jb = min( nb, n-j+1 )
216 CALL zherk( 'Lower', 'No transpose', jb, j-1, -one,
217 $ a( j, 1 ), lda, one, a( j, j ), lda )
218 CALL zpotrf2( 'Lower', jb, a( j, j ), lda, info )
219 IF( info.NE.0 )
220 $ GO TO 30
221 IF( j+jb.LE.n ) THEN
222*
223* Compute the current block column.
224*
225 CALL zgemm( 'No transpose', 'Conjugate transpose',
226 $ n-j-jb+1, jb, j-1, -cone, a( j+jb, 1 ),
227 $ lda, a( j, 1 ), lda, cone, a( j+jb, j ),
228 $ lda )
229 CALL ztrsm( 'Right', 'Lower', 'Conjugate transpose',
230 $ 'Non-unit', n-j-jb+1, jb, cone, a( j, j ),
231 $ lda, 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 ZPOTRF
245*
246 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
ZGEMM
Definition zgemm.f:188
subroutine zherk(uplo, trans, n, k, alpha, a, lda, beta, c, ldc)
ZHERK
Definition zherk.f:173
recursive subroutine zpotrf2(uplo, n, a, lda, info)
ZPOTRF2
Definition zpotrf2.f:106
subroutine zpotrf(uplo, n, a, lda, info)
ZPOTRF
Definition zpotrf.f:107
subroutine ztrsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
ZTRSM
Definition ztrsm.f:180