LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ zpotrf2()

recursive subroutine zpotrf2 ( character  uplo,
integer  n,
complex*16, dimension( lda, * )  a,
integer  lda,
integer  info 
)

ZPOTRF2

Purpose:
 ZPOTRF2 computes the Cholesky factorization of a Hermitian
 positive definite matrix A using the recursive algorithm.

 The factorization has the form
    A = U**H * U,  if UPLO = 'U', or
    A = L  * L**H,  if UPLO = 'L',
 where U is an upper triangular matrix and L is lower triangular.

 This is the recursive version of the algorithm. It divides
 the matrix into four submatrices:

        [  A11 | A12  ]  where A11 is n1 by n1 and A22 is n2 by n2
    A = [ -----|----- ]  with n1 = n/2
        [  A21 | A22  ]       n2 = n-n1

 The subroutine calls itself to factor A11. Update and scale A21
 or A12, update A22 then call itself to factor A22.
Parameters
[in]UPLO
          UPLO is CHARACTER*1
          = 'U':  Upper triangle of A is stored;
          = 'L':  Lower triangle of A is stored.
[in]N
          N is INTEGER
          The order of the matrix A.  N >= 0.
[in,out]A
          A is COMPLEX*16 array, dimension (LDA,N)
          On entry, the Hermitian matrix A.  If UPLO = 'U', the leading
          N-by-N upper triangular part of A contains the upper
          triangular part of the matrix A, and the strictly lower
          triangular part of A is not referenced.  If UPLO = 'L', the
          leading N-by-N lower triangular part of A contains the lower
          triangular part of the matrix A, and the strictly upper
          triangular part of A is not referenced.

          On exit, if INFO = 0, the factor U or L from the Cholesky
          factorization A = U**H*U or A = L*L**H.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,N).
[out]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
          > 0:  if INFO = i, the leading principal minor of order i
                is not positive, and the factorization could not be
                completed.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 105 of file zpotrf2.f.

106*
107* -- LAPACK computational routine --
108* -- LAPACK is a software package provided by Univ. of Tennessee, --
109* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
110*
111* .. Scalar Arguments ..
112 CHARACTER UPLO
113 INTEGER INFO, LDA, N
114* ..
115* .. Array Arguments ..
116 COMPLEX*16 A( LDA, * )
117* ..
118*
119* =====================================================================
120*
121* .. Parameters ..
122 DOUBLE PRECISION ONE, ZERO
123 parameter( one = 1.0d+0, zero = 0.0d+0 )
124 COMPLEX*16 CONE
125 parameter( cone = (1.0d+0, 0.0d+0) )
126* ..
127* .. Local Scalars ..
128 LOGICAL UPPER
129 INTEGER N1, N2, IINFO
130 DOUBLE PRECISION AJJ
131* ..
132* .. External Functions ..
133 LOGICAL LSAME, DISNAN
134 EXTERNAL lsame, disnan
135* ..
136* .. External Subroutines ..
137 EXTERNAL zherk, ztrsm, xerbla
138* ..
139* .. Intrinsic Functions ..
140 INTRINSIC max, dble, sqrt
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( 'ZPOTRF2', -info )
157 RETURN
158 END IF
159*
160* Quick return if possible
161*
162 IF( n.EQ.0 )
163 $ RETURN
164*
165* N=1 case
166*
167 IF( n.EQ.1 ) THEN
168*
169* Test for non-positive-definiteness
170*
171 ajj = dble( a( 1, 1 ) )
172 IF( ajj.LE.zero.OR.disnan( ajj ) ) THEN
173 info = 1
174 RETURN
175 END IF
176*
177* Factor
178*
179 a( 1, 1 ) = sqrt( ajj )
180*
181* Use recursive code
182*
183 ELSE
184 n1 = n/2
185 n2 = n-n1
186*
187* Factor A11
188*
189 CALL zpotrf2( uplo, n1, a( 1, 1 ), lda, iinfo )
190 IF ( iinfo.NE.0 ) THEN
191 info = iinfo
192 RETURN
193 END IF
194*
195* Compute the Cholesky factorization A = U**H*U
196*
197 IF( upper ) THEN
198*
199* Update and scale A12
200*
201 CALL ztrsm( 'L', 'U', 'C', 'N', n1, n2, cone,
202 $ a( 1, 1 ), lda, a( 1, n1+1 ), lda )
203*
204* Update and factor A22
205*
206 CALL zherk( uplo, 'C', n2, n1, -one, a( 1, n1+1 ), lda,
207 $ one, a( n1+1, n1+1 ), lda )
208 CALL zpotrf2( uplo, n2, a( n1+1, n1+1 ), lda, iinfo )
209 IF ( iinfo.NE.0 ) THEN
210 info = iinfo + n1
211 RETURN
212 END IF
213*
214* Compute the Cholesky factorization A = L*L**H
215*
216 ELSE
217*
218* Update and scale A21
219*
220 CALL ztrsm( 'R', 'L', 'C', 'N', n2, n1, cone,
221 $ a( 1, 1 ), lda, a( n1+1, 1 ), lda )
222*
223* Update and factor A22
224*
225 CALL zherk( uplo, 'N', n2, n1, -one, a( n1+1, 1 ), lda,
226 $ one, a( n1+1, n1+1 ), lda )
227 CALL zpotrf2( uplo, n2, a( n1+1, n1+1 ), lda, iinfo )
228 IF ( iinfo.NE.0 ) THEN
229 info = iinfo + n1
230 RETURN
231 END IF
232 END IF
233 END IF
234 RETURN
235*
236* End of ZPOTRF2
237*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zherk(uplo, trans, n, k, alpha, a, lda, beta, c, ldc)
ZHERK
Definition zherk.f:173
logical function disnan(din)
DISNAN tests input for NaN.
Definition disnan.f:59
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
recursive subroutine zpotrf2(uplo, n, a, lda, info)
ZPOTRF2
Definition zpotrf2.f:106
subroutine ztrsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
ZTRSM
Definition ztrsm.f:180
Here is the call graph for this function:
Here is the caller graph for this function: