LAPACK 3.11.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
spotrf.f
Go to the documentation of this file.
1C> \brief \b SPOTRF 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 SPOTRF ( UPLO, N, A, LDA, INFO )
12*
13* .. Scalar Arguments ..
14* CHARACTER UPLO
15* INTEGER INFO, LDA, N
16* ..
17* .. Array Arguments ..
18* REAL A( LDA, * )
19* ..
20*
21* Purpose
22* =======
23*
24C>\details \b Purpose:
25C>\verbatim
26C>
27C> SPOTRF computes the Cholesky factorization of a real symmetric
28C> positive definite matrix A.
29C>
30C> The factorization has the form
31C> A = U**T * U, if UPLO = 'U', or
32C> A = L * L**T, if UPLO = 'L',
33C> where U is an upper triangular matrix and L is lower triangular.
34C>
35C> This is the top-looking block version of the algorithm, calling Level 3 BLAS.
36C>
37C>\endverbatim
38*
39* Arguments:
40* ==========
41*
42C> \param[in] UPLO
43C> \verbatim
44C> UPLO is CHARACTER*1
45C> = 'U': Upper triangle of A is stored;
46C> = 'L': Lower triangle of A is stored.
47C> \endverbatim
48C>
49C> \param[in] N
50C> \verbatim
51C> N is INTEGER
52C> The order of the matrix A. N >= 0.
53C> \endverbatim
54C>
55C> \param[in,out] A
56C> \verbatim
57C> A is REAL array, dimension (LDA,N)
58C> On entry, the symmetric matrix A. If UPLO = 'U', the leading
59C> N-by-N upper triangular part of A contains the upper
60C> triangular part of the matrix A, and the strictly lower
61C> triangular part of A is not referenced. If UPLO = 'L', the
62C> leading N-by-N lower triangular part of A contains the lower
63C> triangular part of the matrix A, and the strictly upper
64C> triangular part of A is not referenced.
65C> \endverbatim
66C> \verbatim
67C> On exit, if INFO = 0, the factor U or L from the Cholesky
68C> factorization A = U**T*U or A = L*L**T.
69C> \endverbatim
70C>
71C> \param[in] LDA
72C> \verbatim
73C> LDA is INTEGER
74C> The leading dimension of the array A. LDA >= max(1,N).
75C> \endverbatim
76C>
77C> \param[out] INFO
78C> \verbatim
79C> INFO is INTEGER
80C> = 0: successful exit
81C> < 0: if INFO = -i, the i-th argument had an illegal value
82C> > 0: if INFO = i, the leading minor of order i is not
83C> positive definite, and the factorization could not be
84C> completed.
85C> \endverbatim
86C>
87*
88* Authors:
89* ========
90*
91C> \author Univ. of Tennessee
92C> \author Univ. of California Berkeley
93C> \author Univ. of Colorado Denver
94C> \author NAG Ltd.
95*
96C> \date December 2016
97*
98C> \ingroup variantsPOcomputational
99*
100* =====================================================================
101 SUBROUTINE spotrf ( 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 REAL A( LDA, * )
113* ..
114*
115* =====================================================================
116*
117* .. Parameters ..
118 REAL ONE
119 parameter( one = 1.0e+0 )
120* ..
121* .. Local Scalars ..
122 LOGICAL UPPER
123 INTEGER J, JB, NB
124* ..
125* .. External Functions ..
126 LOGICAL LSAME
127 INTEGER ILAENV
128 EXTERNAL lsame, ilaenv
129* ..
130* .. External Subroutines ..
131 EXTERNAL sgemm, spotf2, ssyrk, strsm, xerbla
132* ..
133* .. Intrinsic Functions ..
134 INTRINSIC max, min
135* ..
136* .. Executable Statements ..
137*
138* Test the input parameters.
139*
140 info = 0
141 upper = lsame( uplo, 'U' )
142 IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
143 info = -1
144 ELSE IF( n.LT.0 ) THEN
145 info = -2
146 ELSE IF( lda.LT.max( 1, n ) ) THEN
147 info = -4
148 END IF
149 IF( info.NE.0 ) THEN
150 CALL xerbla( 'SPOTRF', -info )
151 RETURN
152 END IF
153*
154* Quick return if possible
155*
156 IF( n.EQ.0 )
157 $ RETURN
158*
159* Determine the block size for this environment.
160*
161 nb = ilaenv( 1, 'SPOTRF', uplo, n, -1, -1, -1 )
162 IF( nb.LE.1 .OR. nb.GE.n ) THEN
163*
164* Use unblocked code.
165*
166 CALL spotf2( uplo, n, a, lda, info )
167 ELSE
168*
169* Use blocked code.
170*
171 IF( upper ) THEN
172*
173* Compute the Cholesky factorization A = U'*U.
174*
175 DO 10 j = 1, n, nb
176
177 jb = min( nb, n-j+1 )
178*
179* Compute the current block.
180*
181 CALL strsm( 'Left', 'Upper', 'Transpose', 'Non-unit',
182 $ j-1, jb, one, a( 1, 1 ), lda,
183 $ a( 1, j ), lda )
184
185 CALL ssyrk( 'Upper', 'Transpose', jb, j-1, -one,
186 $ a( 1, j ), lda,
187 $ one, a( j, j ), lda )
188*
189* Update and factorize the current diagonal block and test
190* for non-positive-definiteness.
191*
192 CALL spotf2( '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 strsm( 'Right', 'Lower', 'Transpose', 'Non-unit',
209 $ jb, j-1, one, a( 1, 1 ), lda,
210 $ a( j, 1 ), lda )
211
212 CALL ssyrk( '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 spotf2( '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 SPOTRF
235*
236 END
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:60
subroutine spotf2(UPLO, N, A, LDA, INFO)
SPOTF2 computes the Cholesky factorization of a symmetric/Hermitian positive definite matrix (unblock...
Definition: spotf2.f:109
subroutine spotrf(UPLO, N, A, LDA, INFO)
SPOTRF
Definition: spotrf.f:107
subroutine strsm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
STRSM
Definition: strsm.f:181
subroutine ssyrk(UPLO, TRANS, N, K, ALPHA, A, LDA, BETA, C, LDC)
SSYRK
Definition: ssyrk.f:169
subroutine sgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
SGEMM
Definition: sgemm.f:187