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