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