LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
ssytrf_aa_2stage.f
Go to the documentation of this file.
1*> \brief \b SSYTRF_AA_2STAGE
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> \htmlonly
9*> Download SSYTRF_AA_2STAGE + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/ssytrf_aa_2stage.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/ssytrf_aa_2stage.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/ssytrf_aa_2stage.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE SSYTRF_AA_2STAGE( UPLO, N, A, LDA, TB, LTB, IPIV,
22* IPIV2, WORK, LWORK, INFO )
23*
24* .. Scalar Arguments ..
25* CHARACTER UPLO
26* INTEGER N, LDA, LTB, LWORK, INFO
27* ..
28* .. Array Arguments ..
29* INTEGER IPIV( * ), IPIV2( * )
30* REAL A( LDA, * ), TB( * ), WORK( * )
31* ..
32*
33*> \par Purpose:
34* =============
35*>
36*> \verbatim
37*>
38*> SSYTRF_AA_2STAGE computes the factorization of a real symmetric matrix A
39*> using the Aasen's algorithm. The form of the factorization is
40*>
41*> A = U**T*T*U or A = L*T*L**T
42*>
43*> where U (or L) is a product of permutation and unit upper (lower)
44*> triangular matrices, and T is a symmetric band matrix with the
45*> bandwidth of NB (NB is internally selected and stored in TB( 1 ), and T is
46*> LU factorized with partial pivoting).
47*>
48*> This is the blocked version of the algorithm, calling Level 3 BLAS.
49*> \endverbatim
50*
51* Arguments:
52* ==========
53*
54*> \param[in] UPLO
55*> \verbatim
56*> UPLO is CHARACTER*1
57*> = 'U': Upper triangle of A is stored;
58*> = 'L': Lower triangle of A is stored.
59*> \endverbatim
60*>
61*> \param[in] N
62*> \verbatim
63*> N is INTEGER
64*> The order of the matrix A. N >= 0.
65*> \endverbatim
66*>
67*> \param[in,out] A
68*> \verbatim
69*> A is REAL array, dimension (LDA,N)
70*> On entry, the symmetric matrix A. If UPLO = 'U', the leading
71*> N-by-N upper triangular part of A contains the upper
72*> triangular part of the matrix A, and the strictly lower
73*> triangular part of A is not referenced. If UPLO = 'L', the
74*> leading N-by-N lower triangular part of A contains the lower
75*> triangular part of the matrix A, and the strictly upper
76*> triangular part of A is not referenced.
77*>
78*> On exit, L is stored below (or above) the subdiagonal blocks,
79*> when UPLO is 'L' (or 'U').
80*> \endverbatim
81*>
82*> \param[in] LDA
83*> \verbatim
84*> LDA is INTEGER
85*> The leading dimension of the array A. LDA >= max(1,N).
86*> \endverbatim
87*>
88*> \param[out] TB
89*> \verbatim
90*> TB is REAL array, dimension (LTB)
91*> On exit, details of the LU factorization of the band matrix.
92*> \endverbatim
93*>
94*> \param[in] LTB
95*> \verbatim
96*> LTB is INTEGER
97*> The size of the array TB. LTB >= 4*N, internally
98*> used to select NB such that LTB >= (3*NB+1)*N.
99*>
100*> If LTB = -1, then a workspace query is assumed; the
101*> routine only calculates the optimal size of LTB,
102*> returns this value as the first entry of TB, and
103*> no error message related to LTB is issued by XERBLA.
104*> \endverbatim
105*>
106*> \param[out] IPIV
107*> \verbatim
108*> IPIV is INTEGER array, dimension (N)
109*> On exit, it contains the details of the interchanges, i.e.,
110*> the row and column k of A were interchanged with the
111*> row and column IPIV(k).
112*> \endverbatim
113*>
114*> \param[out] IPIV2
115*> \verbatim
116*> IPIV2 is INTEGER array, dimension (N)
117*> On exit, it contains the details of the interchanges, i.e.,
118*> the row and column k of T were interchanged with the
119*> row and column IPIV(k).
120*> \endverbatim
121*>
122*> \param[out] WORK
123*> \verbatim
124*> WORK is REAL workspace of size LWORK
125*> \endverbatim
126*>
127*> \param[in] LWORK
128*> \verbatim
129*> LWORK is INTEGER
130*> The size of WORK. LWORK >= N, internally used to select NB
131*> such that LWORK >= N*NB.
132*>
133*> If LWORK = -1, then a workspace query is assumed; the
134*> routine only calculates the optimal size of the WORK array,
135*> returns this value as the first entry of the WORK array, and
136*> no error message related to LWORK is issued by XERBLA.
137*> \endverbatim
138*>
139*> \param[out] INFO
140*> \verbatim
141*> INFO is INTEGER
142*> = 0: successful exit
143*> < 0: if INFO = -i, the i-th argument had an illegal value.
144*> > 0: if INFO = i, band LU factorization failed on i-th column
145*> \endverbatim
146*
147* Authors:
148* ========
149*
150*> \author Univ. of Tennessee
151*> \author Univ. of California Berkeley
152*> \author Univ. of Colorado Denver
153*> \author NAG Ltd.
154*
155*> \ingroup hetrf_aa_2stage
156*
157* =====================================================================
158 SUBROUTINE ssytrf_aa_2stage( UPLO, N, A, LDA, TB, LTB, IPIV,
159 $ IPIV2, WORK, LWORK, INFO )
160*
161* -- LAPACK computational routine --
162* -- LAPACK is a software package provided by Univ. of Tennessee, --
163* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
164*
165 IMPLICIT NONE
166*
167* .. Scalar Arguments ..
168 CHARACTER UPLO
169 INTEGER N, LDA, LTB, LWORK, INFO
170* ..
171* .. Array Arguments ..
172 INTEGER IPIV( * ), IPIV2( * )
173 REAL A( LDA, * ), TB( * ), WORK( * )
174* ..
175*
176* =====================================================================
177* .. Parameters ..
178 REAL ZERO, ONE
179 parameter( zero = 0.0e+0, one = 1.0e+0 )
180*
181* .. Local Scalars ..
182 LOGICAL UPPER, TQUERY, WQUERY
183 INTEGER I, J, K, I1, I2, TD
184 INTEGER LDTB, NB, KB, JB, NT, IINFO
185 REAL PIV
186* ..
187* .. External Functions ..
188 LOGICAL LSAME
189 INTEGER ILAENV
190 REAL SROUNDUP_LWORK
191 EXTERNAL lsame, ilaenv, sroundup_lwork
192* ..
193* .. External Subroutines ..
194 EXTERNAL xerbla, scopy, slacpy,
196 $ ssygst, sswap, strsm
197* ..
198* .. Intrinsic Functions ..
199 INTRINSIC min, max
200* ..
201* .. Executable Statements ..
202*
203* Test the input parameters.
204*
205 info = 0
206 upper = lsame( uplo, 'U' )
207 wquery = ( lwork.EQ.-1 )
208 tquery = ( ltb.EQ.-1 )
209 IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
210 info = -1
211 ELSE IF( n.LT.0 ) THEN
212 info = -2
213 ELSE IF( lda.LT.max( 1, n ) ) THEN
214 info = -4
215 ELSE IF ( ltb .LT. 4*n .AND. .NOT.tquery ) THEN
216 info = -6
217 ELSE IF ( lwork .LT. n .AND. .NOT.wquery ) THEN
218 info = -10
219 END IF
220*
221 IF( info.NE.0 ) THEN
222 CALL xerbla( 'SSYTRF_AA_2STAGE', -info )
223 RETURN
224 END IF
225*
226* Answer the query
227*
228 nb = ilaenv( 1, 'SSYTRF_AA_2STAGE', uplo, n, -1, -1, -1 )
229 IF( info.EQ.0 ) THEN
230 IF( tquery ) THEN
231 tb( 1 ) = (3*nb+1)*n
232 END IF
233 IF( wquery ) THEN
234 work( 1 ) = sroundup_lwork(n*nb)
235 END IF
236 END IF
237 IF( tquery .OR. wquery ) THEN
238 RETURN
239 END IF
240*
241* Quick return
242*
243 IF ( n.EQ.0 ) THEN
244 RETURN
245 ENDIF
246*
247* Determine the number of the block size
248*
249 ldtb = ltb/n
250 IF( ldtb .LT. 3*nb+1 ) THEN
251 nb = (ldtb-1)/3
252 END IF
253 IF( lwork .LT. nb*n ) THEN
254 nb = lwork/n
255 END IF
256*
257* Determine the number of the block columns
258*
259 nt = (n+nb-1)/nb
260 td = 2*nb
261 kb = min(nb, n)
262*
263* Initialize vectors/matrices
264*
265 DO j = 1, kb
266 ipiv( j ) = j
267 END DO
268*
269* Save NB
270*
271 tb( 1 ) = nb
272*
273 IF( upper ) THEN
274*
275* .....................................................
276* Factorize A as U**T*D*U using the upper triangle of A
277* .....................................................
278*
279 DO j = 0, nt-1
280*
281* Generate Jth column of W and H
282*
283 kb = min(nb, n-j*nb)
284 DO i = 1, j-1
285 IF( i.EQ.1 ) THEN
286* H(I,J) = T(I,I)*U(I,J) + T(I+1,I)*U(I+1,J)
287 IF( i .EQ. (j-1) ) THEN
288 jb = nb+kb
289 ELSE
290 jb = 2*nb
291 END IF
292 CALL sgemm( 'NoTranspose', 'NoTranspose',
293 $ nb, kb, jb,
294 $ one, tb( td+1 + (i*nb)*ldtb ), ldtb-1,
295 $ a( (i-1)*nb+1, j*nb+1 ), lda,
296 $ zero, work( i*nb+1 ), n )
297 ELSE
298* H(I,J) = T(I,I-1)*U(I-1,J) + T(I,I)*U(I,J) + T(I,I+1)*U(I+1,J)
299 IF( i .EQ. j-1) THEN
300 jb = 2*nb+kb
301 ELSE
302 jb = 3*nb
303 END IF
304 CALL sgemm( 'NoTranspose', 'NoTranspose',
305 $ nb, kb, jb,
306 $ one, tb( td+nb+1 + ((i-1)*nb)*ldtb ),
307 $ ldtb-1,
308 $ a( (i-2)*nb+1, j*nb+1 ), lda,
309 $ zero, work( i*nb+1 ), n )
310 END IF
311 END DO
312*
313* Compute T(J,J)
314*
315 CALL slacpy( 'Upper', kb, kb, a( j*nb+1, j*nb+1 ), lda,
316 $ tb( td+1 + (j*nb)*ldtb ), ldtb-1 )
317 IF( j.GT.1 ) THEN
318* T(J,J) = U(1:J,J)'*H(1:J)
319 CALL sgemm( 'Transpose', 'NoTranspose',
320 $ kb, kb, (j-1)*nb,
321 $ -one, a( 1, j*nb+1 ), lda,
322 $ work( nb+1 ), n,
323 $ one, tb( td+1 + (j*nb)*ldtb ), ldtb-1 )
324* T(J,J) += U(J,J)'*T(J,J-1)*U(J-1,J)
325 CALL sgemm( 'Transpose', 'NoTranspose',
326 $ kb, nb, kb,
327 $ one, a( (j-1)*nb+1, j*nb+1 ), lda,
328 $ tb( td+nb+1 + ((j-1)*nb)*ldtb ), ldtb-1,
329 $ zero, work( 1 ), n )
330 CALL sgemm( 'NoTranspose', 'NoTranspose',
331 $ kb, kb, nb,
332 $ -one, work( 1 ), n,
333 $ a( (j-2)*nb+1, j*nb+1 ), lda,
334 $ one, tb( td+1 + (j*nb)*ldtb ), ldtb-1 )
335 END IF
336 IF( j.GT.0 ) THEN
337 CALL ssygst( 1, 'Upper', kb,
338 $ tb( td+1 + (j*nb)*ldtb ), ldtb-1,
339 $ a( (j-1)*nb+1, j*nb+1 ), lda, iinfo )
340 END IF
341*
342* Expand T(J,J) into full format
343*
344 DO i = 1, kb
345 DO k = i+1, kb
346 tb( td+(k-i)+1 + (j*nb+i-1)*ldtb )
347 $ = tb( td-(k-(i+1)) + (j*nb+k-1)*ldtb )
348 END DO
349 END DO
350*
351 IF( j.LT.nt-1 ) THEN
352 IF( j.GT.0 ) THEN
353*
354* Compute H(J,J)
355*
356 IF( j.EQ.1 ) THEN
357 CALL sgemm( 'NoTranspose', 'NoTranspose',
358 $ kb, kb, kb,
359 $ one, tb( td+1 + (j*nb)*ldtb ), ldtb-1,
360 $ a( (j-1)*nb+1, j*nb+1 ), lda,
361 $ zero, work( j*nb+1 ), n )
362 ELSE
363 CALL sgemm( 'NoTranspose', 'NoTranspose',
364 $ kb, kb, nb+kb,
365 $ one, tb( td+nb+1 + ((j-1)*nb)*ldtb ),
366 $ ldtb-1,
367 $ a( (j-2)*nb+1, j*nb+1 ), lda,
368 $ zero, work( j*nb+1 ), n )
369 END IF
370*
371* Update with the previous column
372*
373 CALL sgemm( 'Transpose', 'NoTranspose',
374 $ nb, n-(j+1)*nb, j*nb,
375 $ -one, work( nb+1 ), n,
376 $ a( 1, (j+1)*nb+1 ), lda,
377 $ one, a( j*nb+1, (j+1)*nb+1 ), lda )
378 END IF
379*
380* Copy panel to workspace to call SGETRF
381*
382 DO k = 1, nb
383 CALL scopy( n-(j+1)*nb,
384 $ a( j*nb+k, (j+1)*nb+1 ), lda,
385 $ work( 1+(k-1)*n ), 1 )
386 END DO
387*
388* Factorize panel
389*
390 CALL sgetrf( n-(j+1)*nb, nb,
391 $ work, n,
392 $ ipiv( (j+1)*nb+1 ), iinfo )
393c IF (IINFO.NE.0 .AND. INFO.EQ.0) THEN
394c INFO = IINFO+(J+1)*NB
395c END IF
396*
397* Copy panel back
398*
399 DO k = 1, nb
400 CALL scopy( n-(j+1)*nb,
401 $ work( 1+(k-1)*n ), 1,
402 $ a( j*nb+k, (j+1)*nb+1 ), lda )
403 END DO
404*
405* Compute T(J+1, J), zero out for GEMM update
406*
407 kb = min(nb, n-(j+1)*nb)
408 CALL slaset( 'Full', kb, nb, zero, zero,
409 $ tb( td+nb+1 + (j*nb)*ldtb), ldtb-1 )
410 CALL slacpy( 'Upper', kb, nb,
411 $ work, n,
412 $ tb( td+nb+1 + (j*nb)*ldtb ), ldtb-1 )
413 IF( j.GT.0 ) THEN
414 CALL strsm( 'R', 'U', 'N', 'U', kb, nb, one,
415 $ a( (j-1)*nb+1, j*nb+1 ), lda,
416 $ tb( td+nb+1 + (j*nb)*ldtb ), ldtb-1 )
417 END IF
418*
419* Copy T(J,J+1) into T(J+1, J), both upper/lower for GEMM
420* updates
421*
422 DO k = 1, nb
423 DO i = 1, kb
424 tb( td-nb+k-i+1 + (j*nb+nb+i-1)*ldtb )
425 $ = tb( td+nb+i-k+1 + (j*nb+k-1)*ldtb )
426 END DO
427 END DO
428 CALL slaset( 'Lower', kb, nb, zero, one,
429 $ a( j*nb+1, (j+1)*nb+1), lda )
430*
431* Apply pivots to trailing submatrix of A
432*
433 DO k = 1, kb
434* > Adjust ipiv
435 ipiv( (j+1)*nb+k ) = ipiv( (j+1)*nb+k ) + (j+1)*nb
436*
437 i1 = (j+1)*nb+k
438 i2 = ipiv( (j+1)*nb+k )
439 IF( i1.NE.i2 ) THEN
440* > Apply pivots to previous columns of L
441 CALL sswap( k-1, a( (j+1)*nb+1, i1 ), 1,
442 $ a( (j+1)*nb+1, i2 ), 1 )
443* > Swap A(I1+1:M, I1) with A(I2, I1+1:M)
444 IF( i2.GT.(i1+1) )
445 $ CALL sswap( i2-i1-1, a( i1, i1+1 ), lda,
446 $ a( i1+1, i2 ), 1 )
447* > Swap A(I2+1:M, I1) with A(I2+1:M, I2)
448 IF( i2.LT.n )
449 $ CALL sswap( n-i2, a( i1, i2+1 ), lda,
450 $ a( i2, i2+1 ), lda )
451* > Swap A(I1, I1) with A(I2, I2)
452 piv = a( i1, i1 )
453 a( i1, i1 ) = a( i2, i2 )
454 a( i2, i2 ) = piv
455* > Apply pivots to previous columns of L
456 IF( j.GT.0 ) THEN
457 CALL sswap( j*nb, a( 1, i1 ), 1,
458 $ a( 1, i2 ), 1 )
459 END IF
460 ENDIF
461 END DO
462 END IF
463 END DO
464 ELSE
465*
466* .....................................................
467* Factorize A as L*D*L**T using the lower triangle of A
468* .....................................................
469*
470 DO j = 0, nt-1
471*
472* Generate Jth column of W and H
473*
474 kb = min(nb, n-j*nb)
475 DO i = 1, j-1
476 IF( i.EQ.1 ) THEN
477* H(I,J) = T(I,I)*L(J,I)' + T(I+1,I)'*L(J,I+1)'
478 IF( i .EQ. (j-1) ) THEN
479 jb = nb+kb
480 ELSE
481 jb = 2*nb
482 END IF
483 CALL sgemm( 'NoTranspose', 'Transpose',
484 $ nb, kb, jb,
485 $ one, tb( td+1 + (i*nb)*ldtb ), ldtb-1,
486 $ a( j*nb+1, (i-1)*nb+1 ), lda,
487 $ zero, work( i*nb+1 ), n )
488 ELSE
489* H(I,J) = T(I,I-1)*L(J,I-1)' + T(I,I)*L(J,I)' + T(I,I+1)*L(J,I+1)'
490 IF( i .EQ. j-1) THEN
491 jb = 2*nb+kb
492 ELSE
493 jb = 3*nb
494 END IF
495 CALL sgemm( 'NoTranspose', 'Transpose',
496 $ nb, kb, jb,
497 $ one, tb( td+nb+1 + ((i-1)*nb)*ldtb ),
498 $ ldtb-1,
499 $ a( j*nb+1, (i-2)*nb+1 ), lda,
500 $ zero, work( i*nb+1 ), n )
501 END IF
502 END DO
503*
504* Compute T(J,J)
505*
506 CALL slacpy( 'Lower', kb, kb, a( j*nb+1, j*nb+1 ), lda,
507 $ tb( td+1 + (j*nb)*ldtb ), ldtb-1 )
508 IF( j.GT.1 ) THEN
509* T(J,J) = L(J,1:J)*H(1:J)
510 CALL sgemm( 'NoTranspose', 'NoTranspose',
511 $ kb, kb, (j-1)*nb,
512 $ -one, a( j*nb+1, 1 ), lda,
513 $ work( nb+1 ), n,
514 $ one, tb( td+1 + (j*nb)*ldtb ), ldtb-1 )
515* T(J,J) += L(J,J)*T(J,J-1)*L(J,J-1)'
516 CALL sgemm( 'NoTranspose', 'NoTranspose',
517 $ kb, nb, kb,
518 $ one, a( j*nb+1, (j-1)*nb+1 ), lda,
519 $ tb( td+nb+1 + ((j-1)*nb)*ldtb ), ldtb-1,
520 $ zero, work( 1 ), n )
521 CALL sgemm( 'NoTranspose', 'Transpose',
522 $ kb, kb, nb,
523 $ -one, work( 1 ), n,
524 $ a( j*nb+1, (j-2)*nb+1 ), lda,
525 $ one, tb( td+1 + (j*nb)*ldtb ), ldtb-1 )
526 END IF
527 IF( j.GT.0 ) THEN
528 CALL ssygst( 1, 'Lower', kb,
529 $ tb( td+1 + (j*nb)*ldtb ), ldtb-1,
530 $ a( j*nb+1, (j-1)*nb+1 ), lda, iinfo )
531 END IF
532*
533* Expand T(J,J) into full format
534*
535 DO i = 1, kb
536 DO k = i+1, kb
537 tb( td-(k-(i+1)) + (j*nb+k-1)*ldtb )
538 $ = tb( td+(k-i)+1 + (j*nb+i-1)*ldtb )
539 END DO
540 END DO
541*
542 IF( j.LT.nt-1 ) THEN
543 IF( j.GT.0 ) THEN
544*
545* Compute H(J,J)
546*
547 IF( j.EQ.1 ) THEN
548 CALL sgemm( 'NoTranspose', 'Transpose',
549 $ kb, kb, kb,
550 $ one, tb( td+1 + (j*nb)*ldtb ), ldtb-1,
551 $ a( j*nb+1, (j-1)*nb+1 ), lda,
552 $ zero, work( j*nb+1 ), n )
553 ELSE
554 CALL sgemm( 'NoTranspose', 'Transpose',
555 $ kb, kb, nb+kb,
556 $ one, tb( td+nb+1 + ((j-1)*nb)*ldtb ),
557 $ ldtb-1,
558 $ a( j*nb+1, (j-2)*nb+1 ), lda,
559 $ zero, work( j*nb+1 ), n )
560 END IF
561*
562* Update with the previous column
563*
564 CALL sgemm( 'NoTranspose', 'NoTranspose',
565 $ n-(j+1)*nb, nb, j*nb,
566 $ -one, a( (j+1)*nb+1, 1 ), lda,
567 $ work( nb+1 ), n,
568 $ one, a( (j+1)*nb+1, j*nb+1 ), lda )
569 END IF
570*
571* Factorize panel
572*
573 CALL sgetrf( n-(j+1)*nb, nb,
574 $ a( (j+1)*nb+1, j*nb+1 ), lda,
575 $ ipiv( (j+1)*nb+1 ), iinfo )
576c IF (IINFO.NE.0 .AND. INFO.EQ.0) THEN
577c INFO = IINFO+(J+1)*NB
578c END IF
579*
580* Compute T(J+1, J), zero out for GEMM update
581*
582 kb = min(nb, n-(j+1)*nb)
583 CALL slaset( 'Full', kb, nb, zero, zero,
584 $ tb( td+nb+1 + (j*nb)*ldtb), ldtb-1 )
585 CALL slacpy( 'Upper', kb, nb,
586 $ a( (j+1)*nb+1, j*nb+1 ), lda,
587 $ tb( td+nb+1 + (j*nb)*ldtb ), ldtb-1 )
588 IF( j.GT.0 ) THEN
589 CALL strsm( 'R', 'L', 'T', 'U', kb, nb, one,
590 $ a( j*nb+1, (j-1)*nb+1 ), lda,
591 $ tb( td+nb+1 + (j*nb)*ldtb ), ldtb-1 )
592 END IF
593*
594* Copy T(J+1,J) into T(J, J+1), both upper/lower for GEMM
595* updates
596*
597 DO k = 1, nb
598 DO i = 1, kb
599 tb( td-nb+k-i+1 + (j*nb+nb+i-1)*ldtb ) =
600 $ tb( td+nb+i-k+1 + (j*nb+k-1)*ldtb )
601 END DO
602 END DO
603 CALL slaset( 'Upper', kb, nb, zero, one,
604 $ a( (j+1)*nb+1, j*nb+1), lda )
605*
606* Apply pivots to trailing submatrix of A
607*
608 DO k = 1, kb
609* > Adjust ipiv
610 ipiv( (j+1)*nb+k ) = ipiv( (j+1)*nb+k ) + (j+1)*nb
611*
612 i1 = (j+1)*nb+k
613 i2 = ipiv( (j+1)*nb+k )
614 IF( i1.NE.i2 ) THEN
615* > Apply pivots to previous columns of L
616 CALL sswap( k-1, a( i1, (j+1)*nb+1 ), lda,
617 $ a( i2, (j+1)*nb+1 ), lda )
618* > Swap A(I1+1:M, I1) with A(I2, I1+1:M)
619 IF( i2.GT.(i1+1) )
620 $ CALL sswap( i2-i1-1, a( i1+1, i1 ), 1,
621 $ a( i2, i1+1 ), lda )
622* > Swap A(I2+1:M, I1) with A(I2+1:M, I2)
623 IF( i2.LT.n )
624 $ CALL sswap( n-i2, a( i2+1, i1 ), 1,
625 $ a( i2+1, i2 ), 1 )
626* > Swap A(I1, I1) with A(I2, I2)
627 piv = a( i1, i1 )
628 a( i1, i1 ) = a( i2, i2 )
629 a( i2, i2 ) = piv
630* > Apply pivots to previous columns of L
631 IF( j.GT.0 ) THEN
632 CALL sswap( j*nb, a( i1, 1 ), lda,
633 $ a( i2, 1 ), lda )
634 END IF
635 ENDIF
636 END DO
637*
638* Apply pivots to previous columns of L
639*
640c CALL SLASWP( J*NB, A( 1, 1 ), LDA,
641c $ (J+1)*NB+1, (J+1)*NB+KB, IPIV, 1 )
642 END IF
643 END DO
644 END IF
645*
646* Factor the band matrix
647 CALL sgbtrf( n, n, nb, nb, tb, ldtb, ipiv2, info )
648*
649 RETURN
650*
651* End of SSYTRF_AA_2STAGE
652*
653 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine scopy(n, sx, incx, sy, incy)
SCOPY
Definition scopy.f:82
subroutine sgbtrf(m, n, kl, ku, ab, ldab, ipiv, info)
SGBTRF
Definition sgbtrf.f:144
subroutine sgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
SGEMM
Definition sgemm.f:188
subroutine sgetrf(m, n, a, lda, ipiv, info)
SGETRF
Definition sgetrf.f:108
subroutine ssygst(itype, uplo, n, a, lda, b, ldb, info)
SSYGST
Definition ssygst.f:127
subroutine ssytrf_aa_2stage(uplo, n, a, lda, tb, ltb, ipiv, ipiv2, work, lwork, info)
SSYTRF_AA_2STAGE
subroutine slacpy(uplo, m, n, a, lda, b, ldb)
SLACPY copies all or part of one two-dimensional array to another.
Definition slacpy.f:103
subroutine slaset(uplo, m, n, alpha, beta, a, lda)
SLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition slaset.f:110
subroutine sswap(n, sx, incx, sy, incy)
SSWAP
Definition sswap.f:82
subroutine strsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
STRSM
Definition strsm.f:181