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