LAPACK  3.10.0
LAPACK: Linear Algebra PACKage
cpotrf.f
Go to the documentation of this file.
1 C> \brief \b CPOTRF VARIANT: top-looking block version of the algorithm, calling Level 3 BLAS.
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 * Definition:
9 * ===========
10 *
11 * SUBROUTINE CPOTRF ( UPLO, N, A, LDA, INFO )
12 *
13 * .. Scalar Arguments ..
14 * CHARACTER UPLO
15 * INTEGER INFO, LDA, N
16 * ..
17 * .. Array Arguments ..
18 * COMPLEX A( LDA, * )
19 * ..
20 *
21 * Purpose
22 * =======
23 *
24 C>\details \b Purpose:
25 C>\verbatim
26 C>
27 C> CPOTRF computes the Cholesky factorization of a real symmetric
28 C> positive definite matrix A.
29 C>
30 C> The factorization has the form
31 C> A = U**H * U, if UPLO = 'U', or
32 C> A = L * L**H, if UPLO = 'L',
33 C> where U is an upper triangular matrix and L is lower triangular.
34 C>
35 C> This is the top-looking block version of the algorithm, calling Level 3 BLAS.
36 C>
37 C>\endverbatim
38 *
39 * Arguments:
40 * ==========
41 *
42 C> \param[in] UPLO
43 C> \verbatim
44 C> UPLO is CHARACTER*1
45 C> = 'U': Upper triangle of A is stored;
46 C> = 'L': Lower triangle of A is stored.
47 C> \endverbatim
48 C>
49 C> \param[in] N
50 C> \verbatim
51 C> N is INTEGER
52 C> The order of the matrix A. N >= 0.
53 C> \endverbatim
54 C>
55 C> \param[in,out] A
56 C> \verbatim
57 C> A is COMPLEX array, dimension (LDA,N)
58 C> On entry, the symmetric matrix A. If UPLO = 'U', the leading
59 C> N-by-N upper triangular part of A contains the upper
60 C> triangular part of the matrix A, and the strictly lower
61 C> triangular part of A is not referenced. If UPLO = 'L', the
62 C> leading N-by-N lower triangular part of A contains the lower
63 C> triangular part of the matrix A, and the strictly upper
64 C> triangular part of A is not referenced.
65 C> \endverbatim
66 C> \verbatim
67 C> On exit, if INFO = 0, the factor U or L from the Cholesky
68 C> factorization A = U**H*U or A = L*L**H.
69 C> \endverbatim
70 C>
71 C> \param[in] LDA
72 C> \verbatim
73 C> LDA is INTEGER
74 C> The leading dimension of the array A. LDA >= max(1,N).
75 C> \endverbatim
76 C>
77 C> \param[out] INFO
78 C> \verbatim
79 C> INFO is INTEGER
80 C> = 0: successful exit
81 C> < 0: if INFO = -i, the i-th argument had an illegal value
82 C> > 0: if INFO = i, the leading minor of order i is not
83 C> positive definite, and the factorization could not be
84 C> completed.
85 C> \endverbatim
86 C>
87 *
88 * Authors:
89 * ========
90 *
91 C> \author Univ. of Tennessee
92 C> \author Univ. of California Berkeley
93 C> \author Univ. of Colorado Denver
94 C> \author NAG Ltd.
95 *
96 C> \date December 2016
97 *
98 C> \ingroup variantsPOcomputational
99 *
100 * =====================================================================
101  SUBROUTINE cpotrf ( UPLO, N, A, LDA, INFO )
102 *
103 * -- LAPACK computational routine --
104 * -- LAPACK is a software package provided by Univ. of Tennessee, --
105 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
106 *
107 * .. Scalar Arguments ..
108  CHARACTER UPLO
109  INTEGER INFO, LDA, N
110 * ..
111 * .. Array Arguments ..
112  COMPLEX A( LDA, * )
113 * ..
114 *
115 * =====================================================================
116 *
117 * .. Parameters ..
118  REAL ONE
119  COMPLEX CONE
120  parameter( one = 1.0e+0, cone = ( 1.0e+0, 0.0e+0 ) )
121 * ..
122 * .. Local Scalars ..
123  LOGICAL UPPER
124  INTEGER J, JB, NB
125 * ..
126 * .. External Functions ..
127  LOGICAL LSAME
128  INTEGER ILAENV
129  EXTERNAL lsame, ilaenv
130 * ..
131 * .. External Subroutines ..
132  EXTERNAL cgemm, cpotf2, cherk, ctrsm, xerbla
133 * ..
134 * .. Intrinsic Functions ..
135  INTRINSIC max, min
136 * ..
137 * .. Executable Statements ..
138 *
139 * Test the input parameters.
140 *
141  info = 0
142  upper = lsame( uplo, 'U' )
143  IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
144  info = -1
145  ELSE IF( n.LT.0 ) THEN
146  info = -2
147  ELSE IF( lda.LT.max( 1, n ) ) THEN
148  info = -4
149  END IF
150  IF( info.NE.0 ) THEN
151  CALL xerbla( 'CPOTRF', -info )
152  RETURN
153  END IF
154 *
155 * Quick return if possible
156 *
157  IF( n.EQ.0 )
158  $ RETURN
159 *
160 * Determine the block size for this environment.
161 *
162  nb = ilaenv( 1, 'CPOTRF', uplo, n, -1, -1, -1 )
163  IF( nb.LE.1 .OR. nb.GE.n ) THEN
164 *
165 * Use unblocked code.
166 *
167  CALL cpotf2( uplo, n, a, lda, info )
168  ELSE
169 *
170 * Use blocked code.
171 *
172  IF( upper ) THEN
173 *
174 * Compute the Cholesky factorization A = U'*U.
175 *
176  DO 10 j = 1, n, nb
177 
178  jb = min( nb, n-j+1 )
179 *
180 * Compute the current block.
181 *
182  CALL ctrsm( 'Left', 'Upper', 'Conjugate Transpose',
183  $ 'Non-unit', j-1, jb, cone, a( 1, 1 ), lda,
184  $ a( 1, j ), lda )
185 
186  CALL cherk( 'Upper', 'Conjugate Transpose', jb, j-1,
187  $ -one, a( 1, j ), lda, one, a( j, j ), lda )
188 *
189 * Update and factorize the current diagonal block and test
190 * for non-positive-definiteness.
191 *
192  CALL cpotf2( 'Upper', jb, a( j, j ), lda, info )
193  IF( info.NE.0 )
194  $ GO TO 30
195 
196  10 CONTINUE
197 *
198  ELSE
199 *
200 * Compute the Cholesky factorization A = L*L'.
201 *
202  DO 20 j = 1, n, nb
203 
204  jb = min( nb, n-j+1 )
205 *
206 * Compute the current block.
207 *
208  CALL ctrsm( 'Right', 'Lower', 'Conjugate Transpose',
209  $ 'Non-unit', jb, j-1, cone, a( 1, 1 ), lda,
210  $ a( j, 1 ), lda )
211 
212  CALL cherk( 'Lower', 'No Transpose', jb, j-1,
213  $ -one, a( j, 1 ), lda,
214  $ one, a( j, j ), lda )
215 *
216 * Update and factorize the current diagonal block and test
217 * for non-positive-definiteness.
218 *
219  CALL cpotf2( 'Lower', jb, a( j, j ), lda, info )
220  IF( info.NE.0 )
221  $ GO TO 30
222 
223  20 CONTINUE
224  END IF
225  END IF
226  GO TO 40
227 *
228  30 CONTINUE
229  info = info + j - 1
230 *
231  40 CONTINUE
232  RETURN
233 *
234 * End of CPOTRF
235 *
236  END
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:60
subroutine cgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
CGEMM
Definition: cgemm.f:187
subroutine cherk(UPLO, TRANS, N, K, ALPHA, A, LDA, BETA, C, LDC)
CHERK
Definition: cherk.f:173
subroutine ctrsm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
CTRSM
Definition: ctrsm.f:180
subroutine cpotrf(UPLO, N, A, LDA, INFO)
CPOTRF
Definition: cpotrf.f:107
subroutine cpotf2(UPLO, N, A, LDA, INFO)
CPOTF2 computes the Cholesky factorization of a symmetric/Hermitian positive definite matrix (unblock...
Definition: cpotf2.f:109