LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
clavsy.f
Go to the documentation of this file.
1*> \brief \b CLAVSY
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8* Definition:
9* ===========
10*
11* SUBROUTINE CLAVSY( UPLO, TRANS, DIAG, N, NRHS, A, LDA, IPIV, B,
12* LDB, INFO )
13*
14* .. Scalar Arguments ..
15* CHARACTER DIAG, TRANS, UPLO
16* INTEGER INFO, LDA, LDB, N, NRHS
17* ..
18* .. Array Arguments ..
19* INTEGER IPIV( * )
20* COMPLEX A( LDA, * ), B( LDB, * )
21* ..
22*
23*
24*> \par Purpose:
25* =============
26*>
27*> \verbatim
28*>
29*> CLAVSY performs one of the matrix-vector operations
30*> x := A*x or x := A'*x,
31*> where x is an N element vector and A is one of the factors
32*> from the block U*D*U' or L*D*L' factorization computed by CSYTRF.
33*>
34*> If TRANS = 'N', multiplies by U or U * D (or L or L * D)
35*> If TRANS = 'T', multiplies by U' or D * U' (or L' or D * L')
36*> \endverbatim
37*
38* Arguments:
39* ==========
40*
41*> \param[in] UPLO
42*> \verbatim
43*> UPLO is CHARACTER*1
44*> Specifies whether the factor stored in A is upper or lower
45*> triangular.
46*> = 'U': Upper triangular
47*> = 'L': Lower triangular
48*> \endverbatim
49*>
50*> \param[in] TRANS
51*> \verbatim
52*> TRANS is CHARACTER*1
53*> Specifies the operation to be performed:
54*> = 'N': x := A*x
55*> = 'T': x := A'*x
56*> \endverbatim
57*>
58*> \param[in] DIAG
59*> \verbatim
60*> DIAG is CHARACTER*1
61*> Specifies whether or not the diagonal blocks are unit
62*> matrices. If the diagonal blocks are assumed to be unit,
63*> then A = U or A = L, otherwise A = U*D or A = L*D.
64*> = 'U': Diagonal blocks are assumed to be unit matrices.
65*> = 'N': Diagonal blocks are assumed to be non-unit matrices.
66*> \endverbatim
67*>
68*> \param[in] N
69*> \verbatim
70*> N is INTEGER
71*> The number of rows and columns of the matrix A. N >= 0.
72*> \endverbatim
73*>
74*> \param[in] NRHS
75*> \verbatim
76*> NRHS is INTEGER
77*> The number of right hand sides, i.e., the number of vectors
78*> x to be multiplied by A. NRHS >= 0.
79*> \endverbatim
80*>
81*> \param[in] A
82*> \verbatim
83*> A is COMPLEX array, dimension (LDA,N)
84*> The block diagonal matrix D and the multipliers used to
85*> obtain the factor U or L as computed by CSYTRF.
86*> Stored as a 2-D triangular matrix.
87*> \endverbatim
88*>
89*> \param[in] LDA
90*> \verbatim
91*> LDA is INTEGER
92*> The leading dimension of the array A. LDA >= max(1,N).
93*> \endverbatim
94*>
95*> \param[in] IPIV
96*> \verbatim
97*> IPIV is INTEGER array, dimension (N)
98*> Details of the interchanges and the block structure of D,
99*> as determined by CSYTRF.
100*>
101*> If UPLO = 'U':
102*> If IPIV(k) > 0, then rows and columns k and IPIV(k)
103*> were interchanged and D(k,k) is a 1-by-1 diagonal block.
104*> (If IPIV( k ) = k, no interchange was done).
105*>
106*> If IPIV(k) = IPIV(k-1) < 0, then rows and
107*> columns k-1 and -IPIV(k) were interchanged,
108*> D(k-1:k,k-1:k) is a 2-by-2 diagonal block.
109*>
110*> If UPLO = 'L':
111*> If IPIV(k) > 0, then rows and columns k and IPIV(k)
112*> were interchanged and D(k,k) is a 1-by-1 diagonal block.
113*> (If IPIV( k ) = k, no interchange was done).
114*>
115*> If IPIV(k) = IPIV(k+1) < 0, then rows and
116*> columns k+1 and -IPIV(k) were interchanged,
117*> D(k:k+1,k:k+1) is a 2-by-2 diagonal block.
118*> \endverbatim
119*>
120*> \param[in,out] B
121*> \verbatim
122*> B is COMPLEX array, dimension (LDB,NRHS)
123*> On entry, B contains NRHS vectors of length N.
124*> On exit, B is overwritten with the product A * B.
125*> \endverbatim
126*>
127*> \param[in] LDB
128*> \verbatim
129*> LDB is INTEGER
130*> The leading dimension of the array B. LDB >= max(1,N).
131*> \endverbatim
132*>
133*> \param[out] INFO
134*> \verbatim
135*> INFO is INTEGER
136*> = 0: successful exit
137*> < 0: if INFO = -k, the k-th argument had an illegal value
138*> \endverbatim
139*
140* Authors:
141* ========
142*
143*> \author Univ. of Tennessee
144*> \author Univ. of California Berkeley
145*> \author Univ. of Colorado Denver
146*> \author NAG Ltd.
147*
148*> \ingroup complex_lin
149*
150* =====================================================================
151 SUBROUTINE clavsy( UPLO, TRANS, DIAG, N, NRHS, A, LDA, IPIV, B,
152 $ LDB, INFO )
153*
154* -- LAPACK test routine --
155* -- LAPACK is a software package provided by Univ. of Tennessee, --
156* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
157*
158* .. Scalar Arguments ..
159 CHARACTER DIAG, TRANS, UPLO
160 INTEGER INFO, LDA, LDB, N, NRHS
161* ..
162* .. Array Arguments ..
163 INTEGER IPIV( * )
164 COMPLEX A( LDA, * ), B( LDB, * )
165* ..
166*
167* =====================================================================
168*
169* .. Parameters ..
170 COMPLEX CONE
171 parameter( cone = ( 1.0e+0, 0.0e+0 ) )
172* ..
173* .. Local Scalars ..
174 LOGICAL NOUNIT
175 INTEGER J, K, KP
176 COMPLEX D11, D12, D21, D22, T1, T2
177* ..
178* .. External Functions ..
179 LOGICAL LSAME
180 EXTERNAL lsame
181* ..
182* .. External Subroutines ..
183 EXTERNAL cgemv, cgeru, cscal, cswap, xerbla
184* ..
185* .. Intrinsic Functions ..
186 INTRINSIC abs, max
187* ..
188* .. Executable Statements ..
189*
190* Test the input parameters.
191*
192 info = 0
193 IF( .NOT.lsame( uplo, 'U' ) .AND. .NOT.lsame( uplo, 'L' ) ) THEN
194 info = -1
195 ELSE IF( .NOT.lsame( trans, 'N' ) .AND. .NOT.lsame( trans, 'T' ) )
196 $ THEN
197 info = -2
198 ELSE IF( .NOT.lsame( diag, 'U' ) .AND. .NOT.lsame( diag, 'N' ) )
199 $ THEN
200 info = -3
201 ELSE IF( n.LT.0 ) THEN
202 info = -4
203 ELSE IF( lda.LT.max( 1, n ) ) THEN
204 info = -6
205 ELSE IF( ldb.LT.max( 1, n ) ) THEN
206 info = -9
207 END IF
208 IF( info.NE.0 ) THEN
209 CALL xerbla( 'CLAVSY ', -info )
210 RETURN
211 END IF
212*
213* Quick return if possible.
214*
215 IF( n.EQ.0 )
216 $ RETURN
217*
218 nounit = lsame( diag, 'N' )
219*------------------------------------------
220*
221* Compute B := A * B (No transpose)
222*
223*------------------------------------------
224 IF( lsame( trans, 'N' ) ) THEN
225*
226* Compute B := U*B
227* where U = P(m)*inv(U(m))* ... *P(1)*inv(U(1))
228*
229 IF( lsame( uplo, 'U' ) ) THEN
230*
231* Loop forward applying the transformations.
232*
233 k = 1
234 10 CONTINUE
235 IF( k.GT.n )
236 $ GO TO 30
237 IF( ipiv( k ).GT.0 ) THEN
238*
239* 1 x 1 pivot block
240*
241* Multiply by the diagonal element if forming U * D.
242*
243 IF( nounit )
244 $ CALL cscal( nrhs, a( k, k ), b( k, 1 ), ldb )
245*
246* Multiply by P(K) * inv(U(K)) if K > 1.
247*
248 IF( k.GT.1 ) THEN
249*
250* Apply the transformation.
251*
252 CALL cgeru( k-1, nrhs, cone, a( 1, k ), 1, b( k, 1 ),
253 $ ldb, b( 1, 1 ), ldb )
254*
255* Interchange if P(K) != I.
256*
257 kp = ipiv( k )
258 IF( kp.NE.k )
259 $ CALL cswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
260 END IF
261 k = k + 1
262 ELSE
263*
264* 2 x 2 pivot block
265*
266* Multiply by the diagonal block if forming U * D.
267*
268 IF( nounit ) THEN
269 d11 = a( k, k )
270 d22 = a( k+1, k+1 )
271 d12 = a( k, k+1 )
272 d21 = d12
273 DO 20 j = 1, nrhs
274 t1 = b( k, j )
275 t2 = b( k+1, j )
276 b( k, j ) = d11*t1 + d12*t2
277 b( k+1, j ) = d21*t1 + d22*t2
278 20 CONTINUE
279 END IF
280*
281* Multiply by P(K) * inv(U(K)) if K > 1.
282*
283 IF( k.GT.1 ) THEN
284*
285* Apply the transformations.
286*
287 CALL cgeru( k-1, nrhs, cone, a( 1, k ), 1, b( k, 1 ),
288 $ ldb, b( 1, 1 ), ldb )
289 CALL cgeru( k-1, nrhs, cone, a( 1, k+1 ), 1,
290 $ b( k+1, 1 ), ldb, b( 1, 1 ), ldb )
291*
292* Interchange if P(K) != I.
293*
294 kp = abs( ipiv( k ) )
295 IF( kp.NE.k )
296 $ CALL cswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
297 END IF
298 k = k + 2
299 END IF
300 GO TO 10
301 30 CONTINUE
302*
303* Compute B := L*B
304* where L = P(1)*inv(L(1))* ... *P(m)*inv(L(m)) .
305*
306 ELSE
307*
308* Loop backward applying the transformations to B.
309*
310 k = n
311 40 CONTINUE
312 IF( k.LT.1 )
313 $ GO TO 60
314*
315* Test the pivot index. If greater than zero, a 1 x 1
316* pivot was used, otherwise a 2 x 2 pivot was used.
317*
318 IF( ipiv( k ).GT.0 ) THEN
319*
320* 1 x 1 pivot block:
321*
322* Multiply by the diagonal element if forming L * D.
323*
324 IF( nounit )
325 $ CALL cscal( nrhs, a( k, k ), b( k, 1 ), ldb )
326*
327* Multiply by P(K) * inv(L(K)) if K < N.
328*
329 IF( k.NE.n ) THEN
330 kp = ipiv( k )
331*
332* Apply the transformation.
333*
334 CALL cgeru( n-k, nrhs, cone, a( k+1, k ), 1,
335 $ b( k, 1 ), ldb, b( k+1, 1 ), ldb )
336*
337* Interchange if a permutation was applied at the
338* K-th step of the factorization.
339*
340 IF( kp.NE.k )
341 $ CALL cswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
342 END IF
343 k = k - 1
344*
345 ELSE
346*
347* 2 x 2 pivot block:
348*
349* Multiply by the diagonal block if forming L * D.
350*
351 IF( nounit ) THEN
352 d11 = a( k-1, k-1 )
353 d22 = a( k, k )
354 d21 = a( k, k-1 )
355 d12 = d21
356 DO 50 j = 1, nrhs
357 t1 = b( k-1, j )
358 t2 = b( k, j )
359 b( k-1, j ) = d11*t1 + d12*t2
360 b( k, j ) = d21*t1 + d22*t2
361 50 CONTINUE
362 END IF
363*
364* Multiply by P(K) * inv(L(K)) if K < N.
365*
366 IF( k.NE.n ) THEN
367*
368* Apply the transformation.
369*
370 CALL cgeru( n-k, nrhs, cone, a( k+1, k ), 1,
371 $ b( k, 1 ), ldb, b( k+1, 1 ), ldb )
372 CALL cgeru( n-k, nrhs, cone, a( k+1, k-1 ), 1,
373 $ b( k-1, 1 ), ldb, b( k+1, 1 ), ldb )
374*
375* Interchange if a permutation was applied at the
376* K-th step of the factorization.
377*
378 kp = abs( ipiv( k ) )
379 IF( kp.NE.k )
380 $ CALL cswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
381 END IF
382 k = k - 2
383 END IF
384 GO TO 40
385 60 CONTINUE
386 END IF
387*----------------------------------------
388*
389* Compute B := A' * B (transpose)
390*
391*----------------------------------------
392 ELSE IF( lsame( trans, 'T' ) ) THEN
393*
394* Form B := U'*B
395* where U = P(m)*inv(U(m))* ... *P(1)*inv(U(1))
396* and U' = inv(U'(1))*P(1)* ... *inv(U'(m))*P(m)
397*
398 IF( lsame( uplo, 'U' ) ) THEN
399*
400* Loop backward applying the transformations.
401*
402 k = n
403 70 IF( k.LT.1 )
404 $ GO TO 90
405*
406* 1 x 1 pivot block.
407*
408 IF( ipiv( k ).GT.0 ) THEN
409 IF( k.GT.1 ) THEN
410*
411* Interchange if P(K) != I.
412*
413 kp = ipiv( k )
414 IF( kp.NE.k )
415 $ CALL cswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
416*
417* Apply the transformation
418*
419 CALL cgemv( 'Transpose', k-1, nrhs, cone, b, ldb,
420 $ a( 1, k ), 1, cone, b( k, 1 ), ldb )
421 END IF
422 IF( nounit )
423 $ CALL cscal( nrhs, a( k, k ), b( k, 1 ), ldb )
424 k = k - 1
425*
426* 2 x 2 pivot block.
427*
428 ELSE
429 IF( k.GT.2 ) THEN
430*
431* Interchange if P(K) != I.
432*
433 kp = abs( ipiv( k ) )
434 IF( kp.NE.k-1 )
435 $ CALL cswap( nrhs, b( k-1, 1 ), ldb, b( kp, 1 ),
436 $ ldb )
437*
438* Apply the transformations
439*
440 CALL cgemv( 'Transpose', k-2, nrhs, cone, b, ldb,
441 $ a( 1, k ), 1, cone, b( k, 1 ), ldb )
442 CALL cgemv( 'Transpose', k-2, nrhs, cone, b, ldb,
443 $ a( 1, k-1 ), 1, cone, b( k-1, 1 ), ldb )
444 END IF
445*
446* Multiply by the diagonal block if non-unit.
447*
448 IF( nounit ) THEN
449 d11 = a( k-1, k-1 )
450 d22 = a( k, k )
451 d12 = a( k-1, k )
452 d21 = d12
453 DO 80 j = 1, nrhs
454 t1 = b( k-1, j )
455 t2 = b( k, j )
456 b( k-1, j ) = d11*t1 + d12*t2
457 b( k, j ) = d21*t1 + d22*t2
458 80 CONTINUE
459 END IF
460 k = k - 2
461 END IF
462 GO TO 70
463 90 CONTINUE
464*
465* Form B := L'*B
466* where L = P(1)*inv(L(1))* ... *P(m)*inv(L(m))
467* and L' = 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 100 CONTINUE
475 IF( k.GT.n )
476 $ GO TO 120
477*
478* 1 x 1 pivot block
479*
480 IF( ipiv( k ).GT.0 ) THEN
481 IF( k.LT.n ) THEN
482*
483* Interchange if P(K) != I.
484*
485 kp = ipiv( k )
486 IF( kp.NE.k )
487 $ CALL cswap( nrhs, b( k, 1 ), ldb, b( kp, 1 ), ldb )
488*
489* Apply the transformation
490*
491 CALL cgemv( 'Transpose', n-k, nrhs, cone, b( k+1, 1 ),
492 $ ldb, a( k+1, k ), 1, cone, b( k, 1 ), ldb )
493 END IF
494 IF( nounit )
495 $ CALL cscal( nrhs, a( k, k ), b( k, 1 ), ldb )
496 k = k + 1
497*
498* 2 x 2 pivot block.
499*
500 ELSE
501 IF( k.LT.n-1 ) THEN
502*
503* Interchange if P(K) != I.
504*
505 kp = abs( ipiv( k ) )
506 IF( kp.NE.k+1 )
507 $ CALL cswap( nrhs, b( k+1, 1 ), ldb, b( kp, 1 ),
508 $ ldb )
509*
510* Apply the transformation
511*
512 CALL cgemv( 'Transpose', n-k-1, nrhs, cone,
513 $ b( k+2, 1 ), ldb, a( k+2, k+1 ), 1, cone,
514 $ b( k+1, 1 ), ldb )
515 CALL cgemv( 'Transpose', n-k-1, nrhs, cone,
516 $ b( k+2, 1 ), ldb, a( k+2, k ), 1, cone,
517 $ b( k, 1 ), ldb )
518 END IF
519*
520* Multiply by the diagonal block if non-unit.
521*
522 IF( nounit ) THEN
523 d11 = a( k, k )
524 d22 = a( k+1, k+1 )
525 d21 = a( k+1, k )
526 d12 = d21
527 DO 110 j = 1, nrhs
528 t1 = b( k, j )
529 t2 = b( k+1, j )
530 b( k, j ) = d11*t1 + d12*t2
531 b( k+1, j ) = d21*t1 + d22*t2
532 110 CONTINUE
533 END IF
534 k = k + 2
535 END IF
536 GO TO 100
537 120 CONTINUE
538 END IF
539 END IF
540 RETURN
541*
542* End of CLAVSY
543*
544 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine clavsy(uplo, trans, diag, n, nrhs, a, lda, ipiv, b, ldb, info)
CLAVSY
Definition clavsy.f:153
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 cscal(n, ca, cx, incx)
CSCAL
Definition cscal.f:78
subroutine cswap(n, cx, incx, cy, incy)
CSWAP
Definition cswap.f:81