KTH framework for Nek5000 toolboxes; testing version  0.0.1
dpotrf.f
Go to the documentation of this file.
1  SUBROUTINE dpotrf( UPLO, N, A, LDA, INFO )
2 *
3 * -- LAPACK routine (version 3.0) --
4 * Univ. of Tennessee, Univ. of California Berkeley, NAG Ltd.,
5 * Courant Institute, Argonne National Lab, and Rice University
6 * March 31, 1993
7 *
8 * .. Scalar Arguments ..
9  CHARACTER UPLO
10  INTEGER INFO, LDA, N
11 * ..
12 * .. Array Arguments ..
13  DOUBLE PRECISION A( LDA, * )
14 * ..
15 *
16 * Purpose
17 * =======
18 *
19 * DPOTRF computes the Cholesky factorization of a real symmetric
20 * positive definite matrix A.
21 *
22 * The factorization has the form
23 * A = U**T * U, if UPLO = 'U', or
24 * A = L * L**T, if UPLO = 'L',
25 * where U is an upper triangular matrix and L is lower triangular.
26 *
27 * This is the block version of the algorithm, calling Level 3 BLAS.
28 *
29 * Arguments
30 * =========
31 *
32 * UPLO (input) CHARACTER*1
33 * = 'U': Upper triangle of A is stored;
34 * = 'L': Lower triangle of A is stored.
35 *
36 * N (input) INTEGER
37 * The order of the matrix A. N >= 0.
38 *
39 * A (input/output) DOUBLE PRECISION array, dimension (LDA,N)
40 * On entry, the symmetric matrix A. If UPLO = 'U', the leading
41 * N-by-N upper triangular part of A contains the upper
42 * triangular part of the matrix A, and the strictly lower
43 * triangular part of A is not referenced. If UPLO = 'L', the
44 * leading N-by-N lower triangular part of A contains the lower
45 * triangular part of the matrix A, and the strictly upper
46 * triangular part of A is not referenced.
47 *
48 * On exit, if INFO = 0, the factor U or L from the Cholesky
49 * factorization A = U**T*U or A = L*L**T.
50 *
51 * LDA (input) INTEGER
52 * The leading dimension of the array A. LDA >= max(1,N).
53 *
54 * INFO (output) INTEGER
55 * = 0: successful exit
56 * < 0: if INFO = -i, the i-th argument had an illegal value
57 * > 0: if INFO = i, the leading minor of order i is not
58 * positive definite, and the factorization could not be
59 * completed.
60 *
61 * =====================================================================
62 *
63 * .. Parameters ..
64  DOUBLE PRECISION ONE
65  parameter( one = 1.0d+0 )
66 * ..
67 * .. Local Scalars ..
68  LOGICAL UPPER
69  INTEGER J, JB, NB
70 * ..
71 * .. External Functions ..
72  LOGICAL LSAME
73  INTEGER ILAENV
74  EXTERNAL lsame, ilaenv
75 * ..
76 * .. External Subroutines ..
77  EXTERNAL dgemm, dpotf2, dsyrk, dtrsm, xerbla
78 * ..
79 * .. Intrinsic Functions ..
80  INTRINSIC max, min
81 * ..
82 * .. Executable Statements ..
83 *
84 * Test the input parameters.
85 *
86  info = 0
87  upper = lsame( uplo, 'U' )
88  IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
89  info = -1
90  ELSE IF( n.LT.0 ) THEN
91  info = -2
92  ELSE IF( lda.LT.max( 1, n ) ) THEN
93  info = -4
94  END IF
95  IF( info.NE.0 ) THEN
96  CALL xerbla( 'DPOTRF', -info )
97  RETURN
98  END IF
99 *
100 * Quick return if possible
101 *
102  IF( n.EQ.0 )
103  $ RETURN
104 *
105 * Determine the block size for this environment.
106 *
107  nb = ilaenv( 1, 'DPOTRF', uplo, n, -1, -1, -1 )
108  IF( nb.LE.1 .OR. nb.GE.n ) THEN
109 *
110 * Use unblocked code.
111 *
112  CALL dpotf2( uplo, n, a, lda, info )
113  ELSE
114 *
115 * Use blocked code.
116 *
117  IF( upper ) THEN
118 *
119 * Compute the Cholesky factorization A = U'*U.
120 *
121  DO 10 j = 1, n, nb
122 *
123 * Update and factorize the current diagonal block and test
124 * for non-positive-definiteness.
125 *
126  jb = min( nb, n-j+1 )
127  CALL dsyrk( 'Upper', 'Transpose', jb, j-1, -one,
128  $ a( 1, j ), lda, one, a( j, j ), lda )
129  CALL dpotf2( 'Upper', jb, a( j, j ), lda, info )
130  IF( info.NE.0 )
131  $ GO TO 30
132  IF( j+jb.LE.n ) THEN
133 *
134 * Compute the current block row.
135 *
136  CALL dgemm( 'Transpose', 'No transpose', jb, n-j-jb+1,
137  $ j-1, -one, a( 1, j ), lda, a( 1, j+jb ),
138  $ lda, one, a( j, j+jb ), lda )
139  CALL dtrsm( 'Left', 'Upper', 'Transpose', 'Non-unit',
140  $ jb, n-j-jb+1, one, a( j, j ), lda,
141  $ a( j, j+jb ), lda )
142  END IF
143  10 CONTINUE
144 *
145  ELSE
146 *
147 * Compute the Cholesky factorization A = L*L'.
148 *
149  DO 20 j = 1, n, nb
150 *
151 * Update and factorize the current diagonal block and test
152 * for non-positive-definiteness.
153 *
154  jb = min( nb, n-j+1 )
155  CALL dsyrk( 'Lower', 'No transpose', jb, j-1, -one,
156  $ a( j, 1 ), lda, one, a( j, j ), lda )
157  CALL dpotf2( 'Lower', jb, a( j, j ), lda, info )
158  IF( info.NE.0 )
159  $ GO TO 30
160  IF( j+jb.LE.n ) THEN
161 *
162 * Compute the current block column.
163 *
164  CALL dgemm( 'No transpose', 'Transpose', n-j-jb+1, jb,
165  $ j-1, -one, a( j+jb, 1 ), lda, a( j, 1 ),
166  $ lda, one, a( j+jb, j ), lda )
167  CALL dtrsm( 'Right', 'Lower', 'Transpose', 'Non-unit',
168  $ n-j-jb+1, jb, one, a( j, j ), lda,
169  $ a( j+jb, j ), lda )
170  END IF
171  20 CONTINUE
172  END IF
173  END IF
174  GO TO 40
175 *
176  30 CONTINUE
177  info = info + j - 1
178 *
179  40 CONTINUE
180  RETURN
181 *
182 * End of DPOTRF
183 *
184  END
subroutine dgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
Definition: dgemm.f:3
subroutine dpotf2(UPLO, N, A, LDA, INFO)
Definition: dpotf2.f:2
subroutine dpotrf(UPLO, N, A, LDA, INFO)
Definition: dpotrf.f:2
subroutine dsyrk(UPLO, TRANS, N, K, ALPHA, A, LDA, BETA, C, LDC)
Definition: dsyrk.f:3
subroutine dtrsm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
Definition: dtrsm.f:3
subroutine xerbla(SRNAME, INFO)
Definition: xerbla.f:2