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