LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zggsvp.f
Go to the documentation of this file.
1*> \brief \b ZGGSVP
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> \htmlonly
9*> Download ZGGSVP + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zggsvp.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zggsvp.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zggsvp.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE ZGGSVP( JOBU, JOBV, JOBQ, M, P, N, A, LDA, B, LDB,
22* TOLA, TOLB, K, L, U, LDU, V, LDV, Q, LDQ,
23* IWORK, RWORK, TAU, WORK, INFO )
24*
25* .. Scalar Arguments ..
26* CHARACTER JOBQ, JOBU, JOBV
27* INTEGER INFO, K, L, LDA, LDB, LDQ, LDU, LDV, M, N, P
28* DOUBLE PRECISION TOLA, TOLB
29* ..
30* .. Array Arguments ..
31* INTEGER IWORK( * )
32* DOUBLE PRECISION RWORK( * )
33* COMPLEX*16 A( LDA, * ), B( LDB, * ), Q( LDQ, * ),
34* $ TAU( * ), U( LDU, * ), V( LDV, * ), WORK( * )
35* ..
36*
37*
38*> \par Purpose:
39* =============
40*>
41*> \verbatim
42*>
43*> This routine is deprecated and has been replaced by routine ZGGSVP3.
44*>
45*> ZGGSVP computes unitary matrices U, V and Q such that
46*>
47*> N-K-L K L
48*> U**H*A*Q = K ( 0 A12 A13 ) if M-K-L >= 0;
49*> L ( 0 0 A23 )
50*> M-K-L ( 0 0 0 )
51*>
52*> N-K-L K L
53*> = K ( 0 A12 A13 ) if M-K-L < 0;
54*> M-K ( 0 0 A23 )
55*>
56*> N-K-L K L
57*> V**H*B*Q = L ( 0 0 B13 )
58*> P-L ( 0 0 0 )
59*>
60*> where the K-by-K matrix A12 and L-by-L matrix B13 are nonsingular
61*> upper triangular; A23 is L-by-L upper triangular if M-K-L >= 0,
62*> otherwise A23 is (M-K)-by-L upper trapezoidal. K+L = the effective
63*> numerical rank of the (M+P)-by-N matrix (A**H,B**H)**H.
64*>
65*> This decomposition is the preprocessing step for computing the
66*> Generalized Singular Value Decomposition (GSVD), see subroutine
67*> ZGGSVD.
68*> \endverbatim
69*
70* Arguments:
71* ==========
72*
73*> \param[in] JOBU
74*> \verbatim
75*> JOBU is CHARACTER*1
76*> = 'U': Unitary matrix U is computed;
77*> = 'N': U is not computed.
78*> \endverbatim
79*>
80*> \param[in] JOBV
81*> \verbatim
82*> JOBV is CHARACTER*1
83*> = 'V': Unitary matrix V is computed;
84*> = 'N': V is not computed.
85*> \endverbatim
86*>
87*> \param[in] JOBQ
88*> \verbatim
89*> JOBQ is CHARACTER*1
90*> = 'Q': Unitary matrix Q is computed;
91*> = 'N': Q is not computed.
92*> \endverbatim
93*>
94*> \param[in] M
95*> \verbatim
96*> M is INTEGER
97*> The number of rows of the matrix A. M >= 0.
98*> \endverbatim
99*>
100*> \param[in] P
101*> \verbatim
102*> P is INTEGER
103*> The number of rows of the matrix B. P >= 0.
104*> \endverbatim
105*>
106*> \param[in] N
107*> \verbatim
108*> N is INTEGER
109*> The number of columns of the matrices A and B. N >= 0.
110*> \endverbatim
111*>
112*> \param[in,out] A
113*> \verbatim
114*> A is COMPLEX*16 array, dimension (LDA,N)
115*> On entry, the M-by-N matrix A.
116*> On exit, A contains the triangular (or trapezoidal) matrix
117*> described in the Purpose section.
118*> \endverbatim
119*>
120*> \param[in] LDA
121*> \verbatim
122*> LDA is INTEGER
123*> The leading dimension of the array A. LDA >= max(1,M).
124*> \endverbatim
125*>
126*> \param[in,out] B
127*> \verbatim
128*> B is COMPLEX*16 array, dimension (LDB,N)
129*> On entry, the P-by-N matrix B.
130*> On exit, B contains the triangular matrix described in
131*> the Purpose section.
132*> \endverbatim
133*>
134*> \param[in] LDB
135*> \verbatim
136*> LDB is INTEGER
137*> The leading dimension of the array B. LDB >= max(1,P).
138*> \endverbatim
139*>
140*> \param[in] TOLA
141*> \verbatim
142*> TOLA is DOUBLE PRECISION
143*> \endverbatim
144*>
145*> \param[in] TOLB
146*> \verbatim
147*> TOLB is DOUBLE PRECISION
148*>
149*> TOLA and TOLB are the thresholds to determine the effective
150*> numerical rank of matrix B and a subblock of A. Generally,
151*> they are set to
152*> TOLA = MAX(M,N)*norm(A)*MAZHEPS,
153*> TOLB = MAX(P,N)*norm(B)*MAZHEPS.
154*> The size of TOLA and TOLB may affect the size of backward
155*> errors of the decomposition.
156*> \endverbatim
157*>
158*> \param[out] K
159*> \verbatim
160*> K is INTEGER
161*> \endverbatim
162*>
163*> \param[out] L
164*> \verbatim
165*> L is INTEGER
166*>
167*> On exit, K and L specify the dimension of the subblocks
168*> described in Purpose section.
169*> K + L = effective numerical rank of (A**H,B**H)**H.
170*> \endverbatim
171*>
172*> \param[out] U
173*> \verbatim
174*> U is COMPLEX*16 array, dimension (LDU,M)
175*> If JOBU = 'U', U contains the unitary matrix U.
176*> If JOBU = 'N', U is not referenced.
177*> \endverbatim
178*>
179*> \param[in] LDU
180*> \verbatim
181*> LDU is INTEGER
182*> The leading dimension of the array U. LDU >= max(1,M) if
183*> JOBU = 'U'; LDU >= 1 otherwise.
184*> \endverbatim
185*>
186*> \param[out] V
187*> \verbatim
188*> V is COMPLEX*16 array, dimension (LDV,P)
189*> If JOBV = 'V', V contains the unitary matrix V.
190*> If JOBV = 'N', V is not referenced.
191*> \endverbatim
192*>
193*> \param[in] LDV
194*> \verbatim
195*> LDV is INTEGER
196*> The leading dimension of the array V. LDV >= max(1,P) if
197*> JOBV = 'V'; LDV >= 1 otherwise.
198*> \endverbatim
199*>
200*> \param[out] Q
201*> \verbatim
202*> Q is COMPLEX*16 array, dimension (LDQ,N)
203*> If JOBQ = 'Q', Q contains the unitary matrix Q.
204*> If JOBQ = 'N', Q is not referenced.
205*> \endverbatim
206*>
207*> \param[in] LDQ
208*> \verbatim
209*> LDQ is INTEGER
210*> The leading dimension of the array Q. LDQ >= max(1,N) if
211*> JOBQ = 'Q'; LDQ >= 1 otherwise.
212*> \endverbatim
213*>
214*> \param[out] IWORK
215*> \verbatim
216*> IWORK is INTEGER array, dimension (N)
217*> \endverbatim
218*>
219*> \param[out] RWORK
220*> \verbatim
221*> RWORK is DOUBLE PRECISION array, dimension (2*N)
222*> \endverbatim
223*>
224*> \param[out] TAU
225*> \verbatim
226*> TAU is COMPLEX*16 array, dimension (N)
227*> \endverbatim
228*>
229*> \param[out] WORK
230*> \verbatim
231*> WORK is COMPLEX*16 array, dimension (max(3*N,M,P))
232*> \endverbatim
233*>
234*> \param[out] INFO
235*> \verbatim
236*> INFO is INTEGER
237*> = 0: successful exit
238*> < 0: if INFO = -i, the i-th argument had an illegal value.
239*> \endverbatim
240*
241* Authors:
242* ========
243*
244*> \author Univ. of Tennessee
245*> \author Univ. of California Berkeley
246*> \author Univ. of Colorado Denver
247*> \author NAG Ltd.
248*
249*> \ingroup ggsvp
250*
251*> \par Further Details:
252* =====================
253*>
254*> \verbatim
255*>
256*> The subroutine uses LAPACK subroutine ZGEQPF for the QR factorization
257*> with column pivoting to detect the effective numerical rank of the
258*> a matrix. It may be replaced by a better rank determination strategy.
259*> \endverbatim
260*>
261* =====================================================================
262 SUBROUTINE zggsvp( JOBU, JOBV, JOBQ, M, P, N, A, LDA, B, LDB,
263 $ TOLA, TOLB, K, L, U, LDU, V, LDV, Q, LDQ,
264 $ IWORK, RWORK, TAU, WORK, INFO )
265*
266* -- LAPACK computational routine --
267* -- LAPACK is a software package provided by Univ. of Tennessee, --
268* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
269*
270* .. Scalar Arguments ..
271 CHARACTER JOBQ, JOBU, JOBV
272 INTEGER INFO, K, L, LDA, LDB, LDQ, LDU, LDV, M, N, P
273 DOUBLE PRECISION TOLA, TOLB
274* ..
275* .. Array Arguments ..
276 INTEGER IWORK( * )
277 DOUBLE PRECISION RWORK( * )
278 COMPLEX*16 A( LDA, * ), B( LDB, * ), Q( LDQ, * ),
279 $ tau( * ), u( ldu, * ), v( ldv, * ), work( * )
280* ..
281*
282* =====================================================================
283*
284* .. Parameters ..
285 COMPLEX*16 CZERO, CONE
286 PARAMETER ( CZERO = ( 0.0d+0, 0.0d+0 ),
287 $ cone = ( 1.0d+0, 0.0d+0 ) )
288* ..
289* .. Local Scalars ..
290 LOGICAL FORWRD, WANTQ, WANTU, WANTV
291 INTEGER I, J
292 COMPLEX*16 T
293* ..
294* .. External Functions ..
295 LOGICAL LSAME
296 EXTERNAL LSAME
297* ..
298* .. External Subroutines ..
299 EXTERNAL xerbla, zgeqpf, zgeqr2, zgerq2, zlacpy, zlapmt,
301* ..
302* .. Intrinsic Functions ..
303 INTRINSIC abs, dble, dimag, max, min
304* ..
305* .. Statement Functions ..
306 DOUBLE PRECISION CABS1
307* ..
308* .. Statement Function definitions ..
309 cabs1( t ) = abs( dble( t ) ) + abs( dimag( t ) )
310* ..
311* .. Executable Statements ..
312*
313* Test the input parameters
314*
315 wantu = lsame( jobu, 'U' )
316 wantv = lsame( jobv, 'V' )
317 wantq = lsame( jobq, 'Q' )
318 forwrd = .true.
319*
320 info = 0
321 IF( .NOT.( wantu .OR. lsame( jobu, 'N' ) ) ) THEN
322 info = -1
323 ELSE IF( .NOT.( wantv .OR. lsame( jobv, 'N' ) ) ) THEN
324 info = -2
325 ELSE IF( .NOT.( wantq .OR. lsame( jobq, 'N' ) ) ) THEN
326 info = -3
327 ELSE IF( m.LT.0 ) THEN
328 info = -4
329 ELSE IF( p.LT.0 ) THEN
330 info = -5
331 ELSE IF( n.LT.0 ) THEN
332 info = -6
333 ELSE IF( lda.LT.max( 1, m ) ) THEN
334 info = -8
335 ELSE IF( ldb.LT.max( 1, p ) ) THEN
336 info = -10
337 ELSE IF( ldu.LT.1 .OR. ( wantu .AND. ldu.LT.m ) ) THEN
338 info = -16
339 ELSE IF( ldv.LT.1 .OR. ( wantv .AND. ldv.LT.p ) ) THEN
340 info = -18
341 ELSE IF( ldq.LT.1 .OR. ( wantq .AND. ldq.LT.n ) ) THEN
342 info = -20
343 END IF
344 IF( info.NE.0 ) THEN
345 CALL xerbla( 'ZGGSVP', -info )
346 RETURN
347 END IF
348*
349* QR with column pivoting of B: B*P = V*( S11 S12 )
350* ( 0 0 )
351*
352 DO 10 i = 1, n
353 iwork( i ) = 0
354 10 CONTINUE
355 CALL zgeqpf( p, n, b, ldb, iwork, tau, work, rwork, info )
356*
357* Update A := A*P
358*
359 CALL zlapmt( forwrd, m, n, a, lda, iwork )
360*
361* Determine the effective rank of matrix B.
362*
363 l = 0
364 DO 20 i = 1, min( p, n )
365 IF( cabs1( b( i, i ) ).GT.tolb )
366 $ l = l + 1
367 20 CONTINUE
368*
369 IF( wantv ) THEN
370*
371* Copy the details of V, and form V.
372*
373 CALL zlaset( 'Full', p, p, czero, czero, v, ldv )
374 IF( p.GT.1 )
375 $ CALL zlacpy( 'Lower', p-1, n, b( 2, 1 ), ldb, v( 2, 1 ),
376 $ ldv )
377 CALL zung2r( p, p, min( p, n ), v, ldv, tau, work, info )
378 END IF
379*
380* Clean up B
381*
382 DO 40 j = 1, l - 1
383 DO 30 i = j + 1, l
384 b( i, j ) = czero
385 30 CONTINUE
386 40 CONTINUE
387 IF( p.GT.l )
388 $ CALL zlaset( 'Full', p-l, n, czero, czero, b( l+1, 1 ), ldb )
389*
390 IF( wantq ) THEN
391*
392* Set Q = I and Update Q := Q*P
393*
394 CALL zlaset( 'Full', n, n, czero, cone, q, ldq )
395 CALL zlapmt( forwrd, n, n, q, ldq, iwork )
396 END IF
397*
398 IF( p.GE.l .AND. n.NE.l ) THEN
399*
400* RQ factorization of ( S11 S12 ) = ( 0 S12 )*Z
401*
402 CALL zgerq2( l, n, b, ldb, tau, work, info )
403*
404* Update A := A*Z**H
405*
406 CALL zunmr2( 'Right', 'Conjugate transpose', m, n, l, b, ldb,
407 $ tau, a, lda, work, info )
408 IF( wantq ) THEN
409*
410* Update Q := Q*Z**H
411*
412 CALL zunmr2( 'Right', 'Conjugate transpose', n, n, l, b,
413 $ ldb, tau, q, ldq, work, info )
414 END IF
415*
416* Clean up B
417*
418 CALL zlaset( 'Full', l, n-l, czero, czero, b, ldb )
419 DO 60 j = n - l + 1, n
420 DO 50 i = j - n + l + 1, l
421 b( i, j ) = czero
422 50 CONTINUE
423 60 CONTINUE
424*
425 END IF
426*
427* Let N-L L
428* A = ( A11 A12 ) M,
429*
430* then the following does the complete QR decomposition of A11:
431*
432* A11 = U*( 0 T12 )*P1**H
433* ( 0 0 )
434*
435 DO 70 i = 1, n - l
436 iwork( i ) = 0
437 70 CONTINUE
438 CALL zgeqpf( m, n-l, a, lda, iwork, tau, work, rwork, info )
439*
440* Determine the effective rank of A11
441*
442 k = 0
443 DO 80 i = 1, min( m, n-l )
444 IF( cabs1( a( i, i ) ).GT.tola )
445 $ k = k + 1
446 80 CONTINUE
447*
448* Update A12 := U**H*A12, where A12 = A( 1:M, N-L+1:N )
449*
450 CALL zunm2r( 'Left', 'Conjugate transpose', m, l, min( m, n-l ),
451 $ a, lda, tau, a( 1, n-l+1 ), lda, work, info )
452*
453 IF( wantu ) THEN
454*
455* Copy the details of U, and form U
456*
457 CALL zlaset( 'Full', m, m, czero, czero, u, ldu )
458 IF( m.GT.1 )
459 $ CALL zlacpy( 'Lower', m-1, n-l, a( 2, 1 ), lda, u( 2, 1 ),
460 $ ldu )
461 CALL zung2r( m, m, min( m, n-l ), u, ldu, tau, work, info )
462 END IF
463*
464 IF( wantq ) THEN
465*
466* Update Q( 1:N, 1:N-L ) = Q( 1:N, 1:N-L )*P1
467*
468 CALL zlapmt( forwrd, n, n-l, q, ldq, iwork )
469 END IF
470*
471* Clean up A: set the strictly lower triangular part of
472* A(1:K, 1:K) = 0, and A( K+1:M, 1:N-L ) = 0.
473*
474 DO 100 j = 1, k - 1
475 DO 90 i = j + 1, k
476 a( i, j ) = czero
477 90 CONTINUE
478 100 CONTINUE
479 IF( m.GT.k )
480 $ CALL zlaset( 'Full', m-k, n-l, czero, czero, a( k+1, 1 ), lda )
481*
482 IF( n-l.GT.k ) THEN
483*
484* RQ factorization of ( T11 T12 ) = ( 0 T12 )*Z1
485*
486 CALL zgerq2( k, n-l, a, lda, tau, work, info )
487*
488 IF( wantq ) THEN
489*
490* Update Q( 1:N,1:N-L ) = Q( 1:N,1:N-L )*Z1**H
491*
492 CALL zunmr2( 'Right', 'Conjugate transpose', n, n-l, k, a,
493 $ lda, tau, q, ldq, work, info )
494 END IF
495*
496* Clean up A
497*
498 CALL zlaset( 'Full', k, n-l-k, czero, czero, a, lda )
499 DO 120 j = n - l - k + 1, n - l
500 DO 110 i = j - n + l + k + 1, k
501 a( i, j ) = czero
502 110 CONTINUE
503 120 CONTINUE
504*
505 END IF
506*
507 IF( m.GT.k ) THEN
508*
509* QR factorization of A( K+1:M,N-L+1:N )
510*
511 CALL zgeqr2( m-k, l, a( k+1, n-l+1 ), lda, tau, work, info )
512*
513 IF( wantu ) THEN
514*
515* Update U(:,K+1:M) := U(:,K+1:M)*U1
516*
517 CALL zunm2r( 'Right', 'No transpose', m, m-k, min( m-k, l ),
518 $ a( k+1, n-l+1 ), lda, tau, u( 1, k+1 ), ldu,
519 $ work, info )
520 END IF
521*
522* Clean up
523*
524 DO 140 j = n - l + 1, n
525 DO 130 i = j - n + k + l + 1, m
526 a( i, j ) = czero
527 130 CONTINUE
528 140 CONTINUE
529*
530 END IF
531*
532 RETURN
533*
534* End of ZGGSVP
535*
536 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zgeqr2(m, n, a, lda, tau, work, info)
ZGEQR2 computes the QR factorization of a general rectangular matrix using an unblocked algorithm.
Definition zgeqr2.f:130
subroutine zgerq2(m, n, a, lda, tau, work, info)
ZGERQ2 computes the RQ factorization of a general rectangular matrix using an unblocked algorithm.
Definition zgerq2.f:123
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
Definition zlacpy.f:103
subroutine zlapmt(forwrd, m, n, x, ldx, k)
ZLAPMT performs a forward or backward permutation of the columns of a matrix.
Definition zlapmt.f:104
subroutine zlaset(uplo, m, n, alpha, beta, a, lda)
ZLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition zlaset.f:106
subroutine zung2r(m, n, k, a, lda, tau, work, info)
ZUNG2R
Definition zung2r.f:114
subroutine zunm2r(side, trans, m, n, k, a, lda, tau, c, ldc, work, info)
ZUNM2R multiplies a general matrix by the unitary matrix from a QR factorization determined by cgeqrf...
Definition zunm2r.f:159
subroutine zunmr2(side, trans, m, n, k, a, lda, tau, c, ldc, work, info)
ZUNMR2 multiplies a general matrix by the unitary matrix from a RQ factorization determined by cgerqf...
Definition zunmr2.f:159
subroutine zgeqpf(m, n, a, lda, jpvt, tau, work, rwork, info)
ZGEQPF
Definition zgeqpf.f:148
subroutine zggsvp(jobu, jobv, jobq, m, p, n, a, lda, b, ldb, tola, tolb, k, l, u, ldu, v, ldv, q, ldq, iwork, rwork, tau, work, info)
ZGGSVP
Definition zggsvp.f:265