LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ clavhp()

subroutine clavhp ( character  uplo,
character  trans,
character  diag,
integer  n,
integer  nrhs,
complex, dimension( * )  a,
integer, dimension( * )  ipiv,
complex, dimension( ldb, * )  b,
integer  ldb,
integer  info 
)

CLAVHP

Purpose:
    CLAVHP  performs one of the matrix-vector operations
       x := A*x  or  x := A^H*x,
    where x is an N element vector and  A is one of the factors
    from the symmetric factorization computed by CHPTRF.
    CHPTRF produces a factorization of the form
         U * D * U^H     or     L * D * L^H,
    where U (or L) is a product of permutation and unit upper (lower)
    triangular matrices, U^H (or L^H) is the conjugate transpose of
    U (or L), and D is Hermitian and block diagonal with 1 x 1 and
    2 x 2 diagonal blocks.  The multipliers for the transformations
    and the upper or lower triangular parts of the diagonal blocks
    are stored columnwise in packed format in the linear array A.

    If TRANS = 'N' or 'n', CLAVHP multiplies either by U or U * D
    (or L or L * D).
    If TRANS = 'C' or 'c', CLAVHP multiplies either by U^H or D * U^H
    (or L^H or D * L^H ).
  UPLO   - CHARACTER*1
           On entry, UPLO specifies whether the triangular matrix
           stored in A is upper or lower triangular.
              UPLO = 'U' or 'u'   The matrix is upper triangular.
              UPLO = 'L' or 'l'   The matrix is lower triangular.
           Unchanged on exit.

  TRANS  - CHARACTER*1
           On entry, TRANS specifies the operation to be performed as
           follows:
              TRANS = 'N' or 'n'   x := A*x.
              TRANS = 'C' or 'c'   x := A^H*x.
           Unchanged on exit.

  DIAG   - CHARACTER*1
           On entry, DIAG specifies whether the diagonal blocks are
           assumed to be unit matrices, as follows:
              DIAG = 'U' or 'u'   Diagonal blocks are unit matrices.
              DIAG = 'N' or 'n'   Diagonal blocks are non-unit.
           Unchanged on exit.

  N      - INTEGER
           On entry, N specifies the order of the matrix A.
           N must be at least zero.
           Unchanged on exit.

  NRHS   - INTEGER
           On entry, NRHS specifies the number of right hand sides,
           i.e., the number of vectors x to be multiplied by A.
           NRHS must be at least zero.
           Unchanged on exit.

  A      - COMPLEX array, dimension( N*(N+1)/2 )
           On entry, A contains a block diagonal matrix and the
           multipliers of the transformations used to obtain it,
           stored as a packed triangular matrix.
           Unchanged on exit.

  IPIV   - INTEGER array, dimension( N )
           On entry, IPIV contains the vector of pivot indices as
           determined by CSPTRF or CHPTRF.
           If IPIV( K ) = K, no interchange was done.
           If IPIV( K ) <> K but IPIV( K ) > 0, then row K was inter-
           changed with row IPIV( K ) and a 1 x 1 pivot block was used.
           If IPIV( K ) < 0 and UPLO = 'U', then row K-1 was exchanged
           with row | IPIV( K ) | and a 2 x 2 pivot block was used.
           If IPIV( K ) < 0 and UPLO = 'L', then row K+1 was exchanged
           with row | IPIV( K ) | and a 2 x 2 pivot block was used.

  B      - COMPLEX array, dimension( LDB, NRHS )
           On entry, B contains NRHS vectors of length N.
           On exit, B is overwritten with the product A * B.

  LDB    - INTEGER
           On entry, LDB contains the leading dimension of B as
           declared in the calling program.  LDB must be at least
           max( 1, N ).
           Unchanged on exit.

  INFO   - INTEGER
           INFO is the error flag.
           On exit, a value of 0 indicates a successful exit.
           A negative value, say -K, indicates that the K-th argument
           has an illegal value.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 129 of file clavhp.f.

131*
132* -- LAPACK test routine --
133* -- LAPACK is a software package provided by Univ. of Tennessee, --
134* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
135*
136* .. Scalar Arguments ..
137 CHARACTER DIAG, TRANS, UPLO
138 INTEGER INFO, LDB, N, NRHS
139* ..
140* .. Array Arguments ..
141 INTEGER IPIV( * )
142 COMPLEX A( * ), B( LDB, * )
143* ..
144*
145* =====================================================================
146*
147* .. Parameters ..
148 COMPLEX ONE
149 parameter( one = ( 1.0e+0, 0.0e+0 ) )
150* ..
151* .. Local Scalars ..
152 LOGICAL NOUNIT
153 INTEGER J, K, KC, KCNEXT, KP
154 COMPLEX D11, D12, D21, D22, T1, T2
155* ..
156* .. External Functions ..
157 LOGICAL LSAME
158 EXTERNAL lsame
159* ..
160* .. External Subroutines ..
161 EXTERNAL cgemv, cgeru, clacgv, cscal, cswap, xerbla
162* ..
163* .. Intrinsic Functions ..
164 INTRINSIC abs, conjg, max
165* ..
166* .. Executable Statements ..
167*
168* Test the input parameters.
169*
170 info = 0
171 IF( .NOT.lsame( uplo, 'U' ) .AND. .NOT.lsame( uplo, 'L' ) ) THEN
172 info = -1
173 ELSE IF( .NOT.lsame( trans, 'N' ) .AND. .NOT.lsame( trans, 'C' ) )
174 $ THEN
175 info = -2
176 ELSE IF( .NOT.lsame( diag, 'U' ) .AND. .NOT.lsame( diag, 'N' ) )
177 $ THEN
178 info = -3
179 ELSE IF( n.LT.0 ) THEN
180 info = -4
181 ELSE IF( ldb.LT.max( 1, n ) ) THEN
182 info = -8
183 END IF
184 IF( info.NE.0 ) THEN
185 CALL xerbla( 'CLAVHP ', -info )
186 RETURN
187 END IF
188*
189* Quick return if possible.
190*
191 IF( n.EQ.0 )
192 $ RETURN
193*
194 nounit = lsame( diag, 'N' )
195*------------------------------------------
196*
197* Compute B := A * B (No transpose)
198*
199*------------------------------------------
200 IF( lsame( trans, 'N' ) ) THEN
201*
202* Compute B := U*B
203* where U = P(m)*inv(U(m))* ... *P(1)*inv(U(1))
204*
205 IF( lsame( uplo, 'U' ) ) THEN
206*
207* Loop forward applying the transformations.
208*
209 k = 1
210 kc = 1
211 10 CONTINUE
212 IF( k.GT.n )
213 $ GO TO 30
214*
215* 1 x 1 pivot block
216*
217 IF( ipiv( k ).GT.0 ) THEN
218*
219* Multiply by the diagonal element if forming U * D.
220*
221 IF( nounit )
222 $ CALL cscal( nrhs, a( kc+k-1 ), b( k, 1 ), ldb )
223*
224* Multiply by P(K) * inv(U(K)) if K > 1.
225*
226 IF( k.GT.1 ) THEN
227*
228* Apply the transformation.
229*
230 CALL cgeru( k-1, nrhs, one, a( kc ), 1, b( k, 1 ),
231 $ ldb, b( 1, 1 ), ldb )
232*
233* Interchange if P(K) != I.
234*
235 kp = ipiv( k )
236 IF( kp.NE.k )
237 $ CALL cswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
238 END IF
239 kc = kc + k
240 k = k + 1
241 ELSE
242*
243* 2 x 2 pivot block
244*
245 kcnext = kc + k
246*
247* Multiply by the diagonal block if forming U * D.
248*
249 IF( nounit ) THEN
250 d11 = a( kcnext-1 )
251 d22 = a( kcnext+k )
252 d12 = a( kcnext+k-1 )
253 d21 = conjg( d12 )
254 DO 20 j = 1, nrhs
255 t1 = b( k, j )
256 t2 = b( k+1, j )
257 b( k, j ) = d11*t1 + d12*t2
258 b( k+1, j ) = d21*t1 + d22*t2
259 20 CONTINUE
260 END IF
261*
262* Multiply by P(K) * inv(U(K)) if K > 1.
263*
264 IF( k.GT.1 ) THEN
265*
266* Apply the transformations.
267*
268 CALL cgeru( k-1, nrhs, one, a( kc ), 1, b( k, 1 ),
269 $ ldb, b( 1, 1 ), ldb )
270 CALL cgeru( k-1, nrhs, one, a( kcnext ), 1,
271 $ b( k+1, 1 ), ldb, b( 1, 1 ), ldb )
272*
273* Interchange if P(K) != I.
274*
275 kp = abs( ipiv( k ) )
276 IF( kp.NE.k )
277 $ CALL cswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
278 END IF
279 kc = kcnext + k + 1
280 k = k + 2
281 END IF
282 GO TO 10
283 30 CONTINUE
284*
285* Compute B := L*B
286* where L = P(1)*inv(L(1))* ... *P(m)*inv(L(m)) .
287*
288 ELSE
289*
290* Loop backward applying the transformations to B.
291*
292 k = n
293 kc = n*( n+1 ) / 2 + 1
294 40 CONTINUE
295 IF( k.LT.1 )
296 $ GO TO 60
297 kc = kc - ( n-k+1 )
298*
299* Test the pivot index. If greater than zero, a 1 x 1
300* pivot was used, otherwise a 2 x 2 pivot was used.
301*
302 IF( ipiv( k ).GT.0 ) THEN
303*
304* 1 x 1 pivot block:
305*
306* Multiply by the diagonal element if forming L * D.
307*
308 IF( nounit )
309 $ CALL cscal( nrhs, a( kc ), b( k, 1 ), ldb )
310*
311* Multiply by P(K) * inv(L(K)) if K < N.
312*
313 IF( k.NE.n ) THEN
314 kp = ipiv( k )
315*
316* Apply the transformation.
317*
318 CALL cgeru( n-k, nrhs, one, a( kc+1 ), 1, b( k, 1 ),
319 $ ldb, b( k+1, 1 ), ldb )
320*
321* Interchange if a permutation was applied at the
322* K-th step of the factorization.
323*
324 IF( kp.NE.k )
325 $ CALL cswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
326 END IF
327 k = k - 1
328*
329 ELSE
330*
331* 2 x 2 pivot block:
332*
333 kcnext = kc - ( n-k+2 )
334*
335* Multiply by the diagonal block if forming L * D.
336*
337 IF( nounit ) THEN
338 d11 = a( kcnext )
339 d22 = a( kc )
340 d21 = a( kcnext+1 )
341 d12 = conjg( d21 )
342 DO 50 j = 1, nrhs
343 t1 = b( k-1, j )
344 t2 = b( k, j )
345 b( k-1, j ) = d11*t1 + d12*t2
346 b( k, j ) = d21*t1 + d22*t2
347 50 CONTINUE
348 END IF
349*
350* Multiply by P(K) * inv(L(K)) if K < N.
351*
352 IF( k.NE.n ) THEN
353*
354* Apply the transformation.
355*
356 CALL cgeru( n-k, nrhs, one, a( kc+1 ), 1, b( k, 1 ),
357 $ ldb, b( k+1, 1 ), ldb )
358 CALL cgeru( n-k, nrhs, one, a( kcnext+2 ), 1,
359 $ b( k-1, 1 ), ldb, b( k+1, 1 ), ldb )
360*
361* Interchange if a permutation was applied at the
362* K-th step of the factorization.
363*
364 kp = abs( ipiv( k ) )
365 IF( kp.NE.k )
366 $ CALL cswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
367 END IF
368 kc = kcnext
369 k = k - 2
370 END IF
371 GO TO 40
372 60 CONTINUE
373 END IF
374*-------------------------------------------------
375*
376* Compute B := A^H * B (conjugate transpose)
377*
378*-------------------------------------------------
379 ELSE
380*
381* Form B := U^H*B
382* where U = P(m)*inv(U(m))* ... *P(1)*inv(U(1))
383* and U^H = inv(U^H(1))*P(1)* ... *inv(U^H(m))*P(m)
384*
385 IF( lsame( uplo, 'U' ) ) THEN
386*
387* Loop backward applying the transformations.
388*
389 k = n
390 kc = n*( n+1 ) / 2 + 1
391 70 IF( k.LT.1 )
392 $ GO TO 90
393 kc = kc - k
394*
395* 1 x 1 pivot block.
396*
397 IF( ipiv( k ).GT.0 ) THEN
398 IF( k.GT.1 ) THEN
399*
400* Interchange if P(K) != I.
401*
402 kp = ipiv( k )
403 IF( kp.NE.k )
404 $ CALL cswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
405*
406* Apply the transformation:
407* y := y - B' * conjg(x)
408* where x is a column of A and y is a row of B.
409*
410 CALL clacgv( nrhs, b( k, 1 ), ldb )
411 CALL cgemv( 'Conjugate', k-1, nrhs, one, b, ldb,
412 $ a( kc ), 1, one, b( k, 1 ), ldb )
413 CALL clacgv( nrhs, b( k, 1 ), ldb )
414 END IF
415 IF( nounit )
416 $ CALL cscal( nrhs, a( kc+k-1 ), b( k, 1 ), ldb )
417 k = k - 1
418*
419* 2 x 2 pivot block.
420*
421 ELSE
422 kcnext = kc - ( k-1 )
423 IF( k.GT.2 ) THEN
424*
425* Interchange if P(K) != I.
426*
427 kp = abs( ipiv( k ) )
428 IF( kp.NE.k-1 )
429 $ CALL cswap( nrhs, b( k-1, 1 ), ldb, b( kp, 1 ),
430 $ ldb )
431*
432* Apply the transformations.
433*
434 CALL clacgv( nrhs, b( k, 1 ), ldb )
435 CALL cgemv( 'Conjugate', k-2, nrhs, one, b, ldb,
436 $ a( kc ), 1, one, b( k, 1 ), ldb )
437 CALL clacgv( nrhs, b( k, 1 ), ldb )
438*
439 CALL clacgv( nrhs, b( k-1, 1 ), ldb )
440 CALL cgemv( 'Conjugate', k-2, nrhs, one, b, ldb,
441 $ a( kcnext ), 1, one, b( k-1, 1 ), ldb )
442 CALL clacgv( nrhs, b( k-1, 1 ), ldb )
443 END IF
444*
445* Multiply by the diagonal block if non-unit.
446*
447 IF( nounit ) THEN
448 d11 = a( kc-1 )
449 d22 = a( kc+k-1 )
450 d12 = a( kc+k-2 )
451 d21 = conjg( d12 )
452 DO 80 j = 1, nrhs
453 t1 = b( k-1, j )
454 t2 = b( k, j )
455 b( k-1, j ) = d11*t1 + d12*t2
456 b( k, j ) = d21*t1 + d22*t2
457 80 CONTINUE
458 END IF
459 kc = kcnext
460 k = k - 2
461 END IF
462 GO TO 70
463 90 CONTINUE
464*
465* Form B := L^H*B
466* where L = P(1)*inv(L(1))* ... *P(m)*inv(L(m))
467* and L^H = inv(L(m))*P(m)* ... *inv(L(1))*P(1)
468*
469 ELSE
470*
471* Loop forward applying the L-transformations.
472*
473 k = 1
474 kc = 1
475 100 CONTINUE
476 IF( k.GT.n )
477 $ GO TO 120
478*
479* 1 x 1 pivot block
480*
481 IF( ipiv( k ).GT.0 ) THEN
482 IF( k.LT.n ) THEN
483*
484* Interchange if P(K) != I.
485*
486 kp = ipiv( k )
487 IF( kp.NE.k )
488 $ CALL cswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
489*
490* Apply the transformation
491*
492 CALL clacgv( nrhs, b( k, 1 ), ldb )
493 CALL cgemv( 'Conjugate', n-k, nrhs, one, b( k+1, 1 ),
494 $ ldb, a( kc+1 ), 1, one, b( k, 1 ), ldb )
495 CALL clacgv( nrhs, b( k, 1 ), ldb )
496 END IF
497 IF( nounit )
498 $ CALL cscal( nrhs, a( kc ), b( k, 1 ), ldb )
499 kc = kc + n - k + 1
500 k = k + 1
501*
502* 2 x 2 pivot block.
503*
504 ELSE
505 kcnext = kc + n - k + 1
506 IF( k.LT.n-1 ) THEN
507*
508* Interchange if P(K) != I.
509*
510 kp = abs( ipiv( k ) )
511 IF( kp.NE.k+1 )
512 $ CALL cswap( nrhs, b( k+1, 1 ), ldb, b( kp, 1 ),
513 $ ldb )
514*
515* Apply the transformation
516*
517 CALL clacgv( nrhs, b( k+1, 1 ), ldb )
518 CALL cgemv( 'Conjugate', n-k-1, nrhs, one,
519 $ b( k+2, 1 ), ldb, a( kcnext+1 ), 1, one,
520 $ b( k+1, 1 ), ldb )
521 CALL clacgv( nrhs, b( k+1, 1 ), ldb )
522*
523 CALL clacgv( nrhs, b( k, 1 ), ldb )
524 CALL cgemv( 'Conjugate', n-k-1, nrhs, one,
525 $ b( k+2, 1 ), ldb, a( kc+2 ), 1, one,
526 $ b( k, 1 ), ldb )
527 CALL clacgv( nrhs, b( k, 1 ), ldb )
528 END IF
529*
530* Multiply by the diagonal block if non-unit.
531*
532 IF( nounit ) THEN
533 d11 = a( kc )
534 d22 = a( kcnext )
535 d21 = a( kc+1 )
536 d12 = conjg( d21 )
537 DO 110 j = 1, nrhs
538 t1 = b( k, j )
539 t2 = b( k+1, j )
540 b( k, j ) = d11*t1 + d12*t2
541 b( k+1, j ) = d21*t1 + d22*t2
542 110 CONTINUE
543 END IF
544 kc = kcnext + ( n-k )
545 k = k + 2
546 END IF
547 GO TO 100
548 120 CONTINUE
549 END IF
550*
551 END IF
552 RETURN
553*
554* End of CLAVHP
555*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine cgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
CGEMV
Definition cgemv.f:160
subroutine cgeru(m, n, alpha, x, incx, y, incy, a, lda)
CGERU
Definition cgeru.f:130
subroutine clacgv(n, x, incx)
CLACGV conjugates a complex vector.
Definition clacgv.f:74
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine cscal(n, ca, cx, incx)
CSCAL
Definition cscal.f:78
subroutine cswap(n, cx, incx, cy, incy)
CSWAP
Definition cswap.f:81
Here is the call graph for this function:
Here is the caller graph for this function: