Lots of trailing whitespaces in the files of Syd. Cleaning this. No big deal.
[platform/upstream/lapack.git] / SRC / dsytrf_aasen.f
1 *> \brief \b DSYTRF_AASEN
2 *
3 *  =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 *            http://www.netlib.org/lapack/explore-html/
7 *
8 *> \htmlonly
9 *> Download DSYTRF_AASEN + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dsytrf_aasen.f">
11 *> [TGZ]</a>
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dsytrf_aasen.f">
13 *> [ZIP]</a>
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dsytrf_aasen.f">
15 *> [TXT]</a>
16 *> \endhtmlonly
17 *
18 *  Definition:
19 *  ===========
20 *
21 *       SUBROUTINE DSYTRF_AASEN( UPLO, N, A, LDA, IPIV, WORK, LWORK, INFO )
22 *
23 *       .. Scalar Arguments ..
24 *       CHARACTER          UPLO
25 *       INTEGER            N, LDA, LWORK, INFO
26 *       ..
27 *       .. Array Arguments ..
28 *       INTEGER            IPIV( * )
29 *       DOUBLE PRECISION   A( LDA, * ), WORK( * )
30 *       ..
31 *
32 *> \par Purpose:
33 *  =============
34 *>
35 *> \verbatim
36 *>
37 *> DSYTRF_AASEN computes the factorization of a real symmetric matrix A
38 *> using the Aasen's algorithm.  The form of the factorization is
39 *>
40 *>    A = U*T*U**T  or  A = L*T*L**T
41 *>
42 *> where U (or L) is a product of permutation and unit upper (lower)
43 *> triangular matrices, and T is a symmetric tridiagonal matrix.
44 *>
45 *> This is the blocked version of the algorithm, calling Level 3 BLAS.
46 *> \endverbatim
47 *
48 *  Arguments:
49 *  ==========
50 *
51 *> \param[in] UPLO
52 *> \verbatim
53 *>          UPLO is CHARACTER*1
54 *>          = 'U':  Upper triangle of A is stored;
55 *>          = 'L':  Lower triangle of A is stored.
56 *> \endverbatim
57 *>
58 *> \param[in] N
59 *> \verbatim
60 *>          N is INTEGER
61 *>          The order of the matrix A.  N >= 0.
62 *> \endverbatim
63 *>
64 *> \param[in,out] A
65 *> \verbatim
66 *>          A is DOUBLE PRECISION array, dimension (LDA,N)
67 *>          On entry, the symmetric matrix A.  If UPLO = 'U', the leading
68 *>          N-by-N upper triangular part of A contains the upper
69 *>          triangular part of the matrix A, and the strictly lower
70 *>          triangular part of A is not referenced.  If UPLO = 'L', the
71 *>          leading N-by-N lower triangular part of A contains the lower
72 *>          triangular part of the matrix A, and the strictly upper
73 *>          triangular part of A is not referenced.
74 *>
75 *>          On exit, the tridiagonal matrix is stored in the diagonals
76 *>          and the subdiagonals of A just below (or above) the diagonals,
77 *>          and L is stored below (or above) the subdiaonals, when UPLO
78 *>          is 'L' (or 'U').
79 *> \endverbatim
80 *>
81 *> \param[in] LDA
82 *> \verbatim
83 *>          LDA is INTEGER
84 *>          The leading dimension of the array A.  LDA >= max(1,N).
85 *> \endverbatim
86 *>
87 *> \param[out] IPIV
88 *> \verbatim
89 *>          IPIV is INTEGER array, dimension (N)
90 *>          On exit, it contains the details of the interchanges, i.e.,
91 *>          the row and column k of A were interchanged with the
92 *>          row and column IPIV(k).
93 *> \endverbatim
94 *>
95 *> \param[out] WORK
96 *> \verbatim
97 *>          WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
98 *>          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
99 *> \endverbatim
100 *>
101 *> \param[in] LWORK
102 *> \verbatim
103 *>          LWORK is INTEGER
104 *>          The length of WORK.  LWORK >=2*N. For optimum performance
105 *>          LWORK >= N*(1+NB), where NB is the optimal blocksize.
106 *>
107 *>          If LWORK = -1, then a workspace query is assumed; the routine
108 *>          only calculates the optimal size of the WORK array, returns
109 *>          this value as the first entry of the WORK array, and no error
110 *>          message related to LWORK is issued by XERBLA.
111 *> \endverbatim
112 *>
113 *> \param[out] INFO
114 *> \verbatim
115 *>          INFO is INTEGER
116 *>          = 0:  successful exit
117 *>          < 0:  if INFO = -i, the i-th argument had an illegal value
118 *>          > 0:  if INFO = i, D(i,i) is exactly zero.  The factorization
119 *>                has been completed, but the block diagonal matrix D is
120 *>                exactly singular, and division by zero will occur if it
121 *>                is used to solve a system of equations.
122 *> \endverbatim
123 *
124 *  Authors:
125 *  ========
126 *
127 *> \author Univ. of Tennessee
128 *> \author Univ. of California Berkeley
129 *> \author Univ. of Colorado Denver
130 *> \author NAG Ltd.
131 *
132 *> \date November 2016
133 *
134 *> \ingroup doubleSYcomputational
135 *
136 *  @precisions fortran d -> s
137 *
138 *  =====================================================================
139       SUBROUTINE DSYTRF_AASEN( UPLO, N, A, LDA, IPIV, WORK, LWORK, INFO)
140 *
141 *  -- LAPACK computational routine (version 3.4.0) --
142 *  -- LAPACK is a software package provided by Univ. of Tennessee,    --
143 *  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
144 *     November 2016
145 *
146       IMPLICIT NONE
147 *
148 *     .. Scalar Arguments ..
149       CHARACTER          UPLO
150       INTEGER            N, LDA, LWORK, INFO
151 *     ..
152 *     .. Array Arguments ..
153       INTEGER            IPIV( * )
154       DOUBLE PRECISION   A( LDA, * ), WORK( * )
155 *     ..
156 *
157 *  =====================================================================
158 *     .. Parameters ..
159       DOUBLE PRECISION   ZERO, ONE
160       PARAMETER          ( ZERO = 0.0D+0, ONE = 1.0D+0 )
161 *
162 *     .. Local Scalars ..
163       LOGICAL            LQUERY, UPPER
164       INTEGER            J, LWKOPT, IINFO
165       INTEGER            NB, MJ, NJ, K1, K2, J1, J2, J3, JB
166       DOUBLE PRECISION   ALPHA
167 *     ..
168 *     .. External Functions ..
169       LOGICAL            LSAME
170       INTEGER            ILAENV
171       EXTERNAL           LSAME, ILAENV
172 *     ..
173 *     .. External Subroutines ..
174       EXTERNAL           XERBLA
175 *     ..
176 *     .. Intrinsic Functions ..
177       INTRINSIC          MAX
178 *     ..
179 *     .. Executable Statements ..
180 *
181 *     Determine the block size
182 *
183       NB = ILAENV( 1, 'DSYTRF', UPLO, N, -1, -1, -1 )
184 *
185 *     Test the input parameters.
186 *
187       INFO = 0
188       UPPER = LSAME( UPLO, 'U' )
189       LQUERY = ( LWORK.EQ.-1 )
190       IF( .NOT.UPPER .AND. .NOT.LSAME( UPLO, 'L' ) ) THEN
191          INFO = -1
192       ELSE IF( N.LT.0 ) THEN
193          INFO = -2
194       ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
195          INFO = -4
196       ELSE IF( LWORK.LT.( 2*N ) .AND. .NOT.LQUERY ) THEN
197          INFO = -7
198       END IF
199 *
200       IF( INFO.EQ.0 ) THEN
201          LWKOPT = (NB+1)*N
202          WORK( 1 ) = LWKOPT
203       END IF
204 *
205       IF( INFO.NE.0 ) THEN
206          CALL XERBLA( 'DSYTRF_AASEN', -INFO )
207          RETURN
208       ELSE IF( LQUERY ) THEN
209          RETURN
210       END IF
211 *
212 *     Quick return
213 *
214       IF ( N.EQ.0 ) THEN
215           RETURN
216       ENDIF
217       IPIV( 1 ) = 1
218       IF ( N.EQ.1 ) THEN
219          IF ( A( 1, 1 ).EQ.ZERO ) THEN
220             INFO = 1
221          END IF
222          RETURN
223       END IF
224 *
225 *     Adjubst block size based on the workspace size
226 *
227       IF( LWORK.LT.((1+NB)*N) ) THEN
228          NB = ( LWORK-N ) / N
229       END IF
230 *
231       IF( UPPER ) THEN
232 *
233 *        .....................................................
234 *        Factorize A as L*D*L**T using the upper triangle of A
235 *        .....................................................
236 *
237 *        Copy first row A(1, 1:N) into H(1:n) (stored in WORK(1:N))
238 *
239          CALL DCOPY( N, A( 1, 1 ), LDA, WORK( 1 ), 1 )
240 *
241 *        J is the main loop index, increasing from 1 to N in steps of
242 *        JB, where JB is the number of columns factorized by DLASYF;
243 *        JB is either NB, or N-J+1 for the last block
244 *
245          J = 0
246  10      CONTINUE
247          IF( J.GE.N )
248      $      GO TO 20
249 *
250 *        each step of the main loop
251 *         J is the last column of the previous panel
252 *         J1 is the first column of the current panel
253 *         K1 identifies if the previous column of the panel has been
254 *          explicitly stored, e.g., K1=1 for the first panel, and
255 *          K1=0 for the rest
256 *
257          J1 = J + 1
258          JB = MIN( N-J1+1, NB )
259          K1 = MAX(1, J)-J
260 *
261 *        Panel factorization
262 *
263          CALL DLASYF_AASEN( UPLO, 2-K1, N-J, JB,
264      $                      A( MAX(1, J), J+1 ), LDA,
265      $                      IPIV( J+1 ), WORK, N, WORK( N*NB+1 ),
266      $                      IINFO )
267          IF( (IINFO.GT.0) .AND. (INFO.EQ.0) ) THEN
268              INFO = IINFO+J
269          ENDIF
270 *
271 *        Ajust IPIV and apply it back (J-th step picks (J+1)-th pivot)
272 *
273          DO J2 = J+2, MIN(N, J+JB+1)
274             IPIV( J2 ) = IPIV( J2 ) + J
275             IF( (J2.NE.IPIV(J2)) .AND. ((J1-K1).GT.2) ) THEN
276                CALL DSWAP( J1-K1-2, A( 1, J2 ), 1,
277      $                              A( 1, IPIV(J2) ), 1 )
278             END IF
279          END DO
280          J = J + JB
281 *
282 *        Trailing submatrix update, where
283 *         the row A(J1-1, J2-1:N) stores U(J1, J2+1:N) and
284 *         WORK stores the current block of the auxiriarly matrix H
285 *
286          IF( J.LT.N ) THEN
287 *
288 *           If first panel and JB=1 (NB=1), then nothing to do
289 *
290             IF( J1.GT.1 .OR. JB.GT.1 ) THEN
291 *
292 *              Merge rank-1 update with BLAS-3 update
293 *
294                ALPHA = A( J, J+1 )
295                A( J, J+1 ) = ONE
296                CALL DCOPY( N-J, A( J-1, J+1 ), LDA,
297      $                          WORK( (J+1-J1+1)+JB*N ), 1 )
298                CALL DSCAL( N-J, ALPHA, WORK( (J+1-J1+1)+JB*N ), 1 )
299 *
300 *              K1 identifies if the previous column of the panel has been
301 *               explicitly stored, e.g., K1=1 and K2= 0 for the first panel,
302 *               while K1=0 and K2=1 for the rest
303 *
304                IF( J1.GT.1 ) THEN
305 *
306 *                 Not first panel
307 *
308                   K2 = 1
309                ELSE
310 *
311 *                 First panel
312 *
313                   K2 = 0
314 *
315 *                 First update skips the first column
316 *
317                   JB = JB - 1
318                END IF
319 *
320                DO J2 = J+1, N, NB
321                   NJ = MIN( NB, N-J2+1 )
322 *
323 *                 Update (J2, J2) diagonal block with DGEMV
324 *
325                   J3 = J2
326                   DO MJ = NJ-1, 1, -1
327                      CALL DGEMV( 'No transpose', MJ, JB+1,
328      $                          -ONE, WORK( J3-J1+1+K1*N ), N,
329      $                                A( J1-K2, J3 ), 1,
330      $                           ONE, A( J3, J3 ), LDA )
331                      J3 = J3 + 1
332                   END DO
333 *
334 *                 Update off-diagonal block of J2-th block row with DGEMM
335 *
336                   CALL DGEMM( 'Transpose', 'Transpose',
337      $                        NJ, N-J3+1, JB+1,
338      $                       -ONE, A( J1-K2, J2 ), LDA,
339      $                             WORK( J3-J1+1+K1*N ), N,
340      $                        ONE, A( J2, J3 ), LDA )
341                END DO
342 *
343 *              Recover T( J, J+1 )
344 *
345                A( J, J+1 ) = ALPHA
346             END IF
347 *
348 *           WORK(J+1, 1) stores H(J+1, 1)
349 *
350             CALL DCOPY( N-J, A( J+1, J+1 ), LDA, WORK( 1 ), 1 )
351          END IF
352          GO TO 10
353       ELSE
354 *
355 *        .....................................................
356 *        Factorize A as L*D*L**T using the lower triangle of A
357 *        .....................................................
358 *
359 *        copy first column A(1:N, 1) into H(1:N, 1)
360 *         (stored in WORK(1:N))
361 *
362          CALL DCOPY( N, A( 1, 1 ), 1, WORK( 1 ), 1 )
363 *
364 *        J is the main loop index, increasing from 1 to N in steps of
365 *        JB, where JB is the number of columns factorized by DLASYF;
366 *        JB is either NB, or N-J+1 for the last block
367 *
368          J = 0
369  11      CONTINUE
370          IF( J.GE.N )
371      $      GO TO 20
372 *
373 *        each step of the main loop
374 *         J is the last column of the previous panel
375 *         J1 is the first column of the current panel
376 *         K1 identifies if the previous column of the panel has been
377 *          explicitly stored, e.g., K1=1 for the first panel, and
378 *          K1=0 for the rest
379 *
380          J1 = J+1
381          JB = MIN( N-J1+1, NB )
382          K1 = MAX(1, J)-J
383 *
384 *        Panel factorization
385 *
386          CALL DLASYF_AASEN( UPLO, 2-K1, N-J, JB,
387      $                      A( J+1, MAX(1, J) ), LDA,
388      $                      IPIV( J+1 ), WORK, N, WORK( N*NB+1 ), IINFO)
389          IF( (IINFO.GT.0) .AND. (INFO.EQ.0) ) THEN
390             INFO = IINFO+J
391          ENDIF
392 *
393 *        Ajust IPIV and apply it back (J-th step picks (J+1)-th pivot)
394 *
395          DO J2 = J+2, MIN(N, J+JB+1)
396             IPIV( J2 ) = IPIV( J2 ) + J
397             IF( (J2.NE.IPIV(J2)) .AND. ((J1-K1).GT.2) ) THEN
398                CALL DSWAP( J1-K1-2, A( J2, 1 ), LDA,
399      $                              A( IPIV(J2), 1 ), LDA )
400             END IF
401          END DO
402          J = J + JB
403 *
404 *        Trailing submatrix update, where
405 *          A(J2+1, J1-1) stores L(J2+1, J1) and
406 *          WORK(J2+1, 1) stores H(J2+1, 1)
407 *
408          IF( J.LT.N ) THEN
409 *
410 *           if first panel and JB=1 (NB=1), then nothing to do
411 *
412             IF( J1.GT.1 .OR. JB.GT.1 ) THEN
413 *
414 *              Merge rank-1 update with BLAS-3 update
415 *
416                ALPHA = A( J+1, J )
417                A( J+1, J ) = ONE
418                CALL DCOPY( N-J, A( J+1, J-1 ), 1,
419      $                          WORK( (J+1-J1+1)+JB*N ), 1 )
420                CALL DSCAL( N-J, ALPHA, WORK( (J+1-J1+1)+JB*N ), 1 )
421 *
422 *              K1 identifies if the previous column of the panel has been
423 *               explicitly stored, e.g., K1=1 and K2= 0 for the first panel,
424 *               while K1=0 and K2=1 for the rest
425 *
426                IF( J1.GT.1 ) THEN
427 *
428 *                 Not first panel
429 *
430                   K2 = 1
431                ELSE
432 *
433 *                 First panel
434 *
435                   K2 = 0
436 *
437 *                 First update skips the first column
438 *
439                   JB = JB - 1
440                END IF
441 *
442                DO J2 = J+1, N, NB
443                   NJ = MIN( NB, N-J2+1 )
444 *
445 *                 Update (J2, J2) diagonal block with DGEMV
446 *
447                   J3 = J2
448                   DO MJ = NJ-1, 1, -1
449                      CALL DGEMV( 'No transpose', MJ, JB+1,
450      $                          -ONE, WORK( J3-J1+1+K1*N ), N,
451      $                                A( J3, J1-K2 ), LDA,
452      $                           ONE, A( J3, J3 ), 1 )
453                      J3 = J3 + 1
454                   END DO
455 *
456 *                 Update off-diagonal block in J2-th block column with DGEMM
457 *
458                   CALL DGEMM( 'No transpose', 'Transpose',
459      $                        N-J3+1, NJ, JB+1,
460      $                       -ONE, WORK( J3-J1+1+K1*N ), N,
461      $                             A( J2, J1-K2 ), LDA,
462      $                        ONE, A( J3, J2 ), LDA )
463                END DO
464 *
465 *              Recover T( J+1, J )
466 *
467                A( J+1, J ) = ALPHA
468             END IF
469 *
470 *           WORK(J+1, 1) stores H(J+1, 1)
471 *
472             CALL DCOPY( N-J, A( J+1, J+1 ), 1, WORK( 1 ), 1 )
473          END IF
474          GO TO 11
475       END IF
476 *
477    20 CONTINUE
478       RETURN
479 *
480 *     End of DSYTRF_AASEN
481 *
482       END