1: *> \brief \b DGGGLM
2: *
3: * =========== DOCUMENTATION ===========
4: *
5: * Online html documentation available at
6: * http://www.netlib.org/lapack/explore-html/
7: *
8: *> \htmlonly
9: *> Download DGGGLM + dependencies
10: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dggglm.f">
11: *> [TGZ]</a>
12: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dggglm.f">
13: *> [ZIP]</a>
14: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dggglm.f">
15: *> [TXT]</a>
16: *> \endhtmlonly
17: *
18: * Definition:
19: * ===========
20: *
21: * SUBROUTINE DGGGLM( N, M, P, A, LDA, B, LDB, D, X, Y, WORK, LWORK,
22: * INFO )
23: *
24: * .. Scalar Arguments ..
25: * INTEGER INFO, LDA, LDB, LWORK, M, N, P
26: * ..
27: * .. Array Arguments ..
28: * DOUBLE PRECISION A( LDA, * ), B( LDB, * ), D( * ), WORK( * ),
29: * $ X( * ), Y( * )
30: * ..
31: *
32: *
33: *> \par Purpose:
34: * =============
35: *>
36: *> \verbatim
37: *>
38: *> DGGGLM solves a general Gauss-Markov linear model (GLM) problem:
39: *>
40: *> minimize || y ||_2 subject to d = A*x + B*y
41: *> x
42: *>
43: *> where A is an N-by-M matrix, B is an N-by-P matrix, and d is a
44: *> given N-vector. It is assumed that M <= N <= M+P, and
45: *>
46: *> rank(A) = M and rank( A B ) = N.
47: *>
48: *> Under these assumptions, the constrained equation is always
49: *> consistent, and there is a unique solution x and a minimal 2-norm
50: *> solution y, which is obtained using a generalized QR factorization
51: *> of the matrices (A, B) given by
52: *>
53: *> A = Q*(R), B = Q*T*Z.
54: *> (0)
55: *>
56: *> In particular, if matrix B is square nonsingular, then the problem
57: *> GLM is equivalent to the following weighted linear least squares
58: *> problem
59: *>
60: *> minimize || inv(B)*(d-A*x) ||_2
61: *> x
62: *>
63: *> where inv(B) denotes the inverse of B.
64: *> \endverbatim
65: *
66: * Arguments:
67: * ==========
68: *
69: *> \param[in] N
70: *> \verbatim
71: *> N is INTEGER
72: *> The number of rows of the matrices A and B. N >= 0.
73: *> \endverbatim
74: *>
75: *> \param[in] M
76: *> \verbatim
77: *> M is INTEGER
78: *> The number of columns of the matrix A. 0 <= M <= N.
79: *> \endverbatim
80: *>
81: *> \param[in] P
82: *> \verbatim
83: *> P is INTEGER
84: *> The number of columns of the matrix B. P >= N-M.
85: *> \endverbatim
86: *>
87: *> \param[in,out] A
88: *> \verbatim
89: *> A is DOUBLE PRECISION array, dimension (LDA,M)
90: *> On entry, the N-by-M matrix A.
91: *> On exit, the upper triangular part of the array A contains
92: *> the M-by-M upper triangular matrix R.
93: *> \endverbatim
94: *>
95: *> \param[in] LDA
96: *> \verbatim
97: *> LDA is INTEGER
98: *> The leading dimension of the array A. LDA >= max(1,N).
99: *> \endverbatim
100: *>
101: *> \param[in,out] B
102: *> \verbatim
103: *> B is DOUBLE PRECISION array, dimension (LDB,P)
104: *> On entry, the N-by-P matrix B.
105: *> On exit, if N <= P, the upper triangle of the subarray
106: *> B(1:N,P-N+1:P) contains the N-by-N upper triangular matrix T;
107: *> if N > P, the elements on and above the (N-P)th subdiagonal
108: *> contain the N-by-P upper trapezoidal matrix T.
109: *> \endverbatim
110: *>
111: *> \param[in] LDB
112: *> \verbatim
113: *> LDB is INTEGER
114: *> The leading dimension of the array B. LDB >= max(1,N).
115: *> \endverbatim
116: *>
117: *> \param[in,out] D
118: *> \verbatim
119: *> D is DOUBLE PRECISION array, dimension (N)
120: *> On entry, D is the left hand side of the GLM equation.
121: *> On exit, D is destroyed.
122: *> \endverbatim
123: *>
124: *> \param[out] X
125: *> \verbatim
126: *> X is DOUBLE PRECISION array, dimension (M)
127: *> \endverbatim
128: *>
129: *> \param[out] Y
130: *> \verbatim
131: *> Y is DOUBLE PRECISION array, dimension (P)
132: *>
133: *> On exit, X and Y are the solutions of the GLM problem.
134: *> \endverbatim
135: *>
136: *> \param[out] WORK
137: *> \verbatim
138: *> WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
139: *> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
140: *> \endverbatim
141: *>
142: *> \param[in] LWORK
143: *> \verbatim
144: *> LWORK is INTEGER
145: *> The dimension of the array WORK. LWORK >= max(1,N+M+P).
146: *> For optimum performance, LWORK >= M+min(N,P)+max(N,P)*NB,
147: *> where NB is an upper bound for the optimal blocksizes for
148: *> DGEQRF, SGERQF, DORMQR and SORMRQ.
149: *>
150: *> If LWORK = -1, then a workspace query is assumed; the routine
151: *> only calculates the optimal size of the WORK array, returns
152: *> this value as the first entry of the WORK array, and no error
153: *> message related to LWORK is issued by XERBLA.
154: *> \endverbatim
155: *>
156: *> \param[out] INFO
157: *> \verbatim
158: *> INFO is INTEGER
159: *> = 0: successful exit.
160: *> < 0: if INFO = -i, the i-th argument had an illegal value.
161: *> = 1: the upper triangular factor R associated with A in the
162: *> generalized QR factorization of the pair (A, B) is
163: *> singular, so that rank(A) < M; the least squares
164: *> solution could not be computed.
165: *> = 2: the bottom (N-M) by (N-M) part of the upper trapezoidal
166: *> factor T associated with B in the generalized QR
167: *> factorization of the pair (A, B) is singular, so that
168: *> rank( A B ) < N; the least squares solution could not
169: *> be computed.
170: *> \endverbatim
171: *
172: * Authors:
173: * ========
174: *
175: *> \author Univ. of Tennessee
176: *> \author Univ. of California Berkeley
177: *> \author Univ. of Colorado Denver
178: *> \author NAG Ltd.
179: *
180: *> \date November 2015
181: *
182: *> \ingroup doubleOTHEReigen
183: *
184: * =====================================================================
185: SUBROUTINE DGGGLM( N, M, P, A, LDA, B, LDB, D, X, Y, WORK, LWORK,
186: $ INFO )
187: *
188: * -- LAPACK driver routine (version 3.6.0) --
189: * -- LAPACK is a software package provided by Univ. of Tennessee, --
190: * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
191: * November 2015
192: *
193: * .. Scalar Arguments ..
194: INTEGER INFO, LDA, LDB, LWORK, M, N, P
195: * ..
196: * .. Array Arguments ..
197: DOUBLE PRECISION A( LDA, * ), B( LDB, * ), D( * ), WORK( * ),
198: $ X( * ), Y( * )
199: * ..
200: *
201: * ===================================================================
202: *
203: * .. Parameters ..
204: DOUBLE PRECISION ZERO, ONE
205: PARAMETER ( ZERO = 0.0D+0, ONE = 1.0D+0 )
206: * ..
207: * .. Local Scalars ..
208: LOGICAL LQUERY
209: INTEGER I, LOPT, LWKMIN, LWKOPT, NB, NB1, NB2, NB3,
210: $ NB4, NP
211: * ..
212: * .. External Subroutines ..
213: EXTERNAL DCOPY, DGEMV, DGGQRF, DORMQR, DORMRQ, DTRTRS,
214: $ XERBLA
215: * ..
216: * .. External Functions ..
217: INTEGER ILAENV
218: EXTERNAL ILAENV
219: * ..
220: * .. Intrinsic Functions ..
221: INTRINSIC INT, MAX, MIN
222: * ..
223: * .. Executable Statements ..
224: *
225: * Test the input parameters
226: *
227: INFO = 0
228: NP = MIN( N, P )
229: LQUERY = ( LWORK.EQ.-1 )
230: IF( N.LT.0 ) THEN
231: INFO = -1
232: ELSE IF( M.LT.0 .OR. M.GT.N ) THEN
233: INFO = -2
234: ELSE IF( P.LT.0 .OR. P.LT.N-M ) THEN
235: INFO = -3
236: ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
237: INFO = -5
238: ELSE IF( LDB.LT.MAX( 1, N ) ) THEN
239: INFO = -7
240: END IF
241: *
242: * Calculate workspace
243: *
244: IF( INFO.EQ.0) THEN
245: IF( N.EQ.0 ) THEN
246: LWKMIN = 1
247: LWKOPT = 1
248: ELSE
249: NB1 = ILAENV( 1, 'DGEQRF', ' ', N, M, -1, -1 )
250: NB2 = ILAENV( 1, 'DGERQF', ' ', N, M, -1, -1 )
251: NB3 = ILAENV( 1, 'DORMQR', ' ', N, M, P, -1 )
252: NB4 = ILAENV( 1, 'DORMRQ', ' ', N, M, P, -1 )
253: NB = MAX( NB1, NB2, NB3, NB4 )
254: LWKMIN = M + N + P
255: LWKOPT = M + NP + MAX( N, P )*NB
256: END IF
257: WORK( 1 ) = LWKOPT
258: *
259: IF( LWORK.LT.LWKMIN .AND. .NOT.LQUERY ) THEN
260: INFO = -12
261: END IF
262: END IF
263: *
264: IF( INFO.NE.0 ) THEN
265: CALL XERBLA( 'DGGGLM', -INFO )
266: RETURN
267: ELSE IF( LQUERY ) THEN
268: RETURN
269: END IF
270: *
271: * Quick return if possible
272: *
273: IF( N.EQ.0 )
274: $ RETURN
275: *
276: * Compute the GQR factorization of matrices A and B:
277: *
278: * Q**T*A = ( R11 ) M, Q**T*B*Z**T = ( T11 T12 ) M
279: * ( 0 ) N-M ( 0 T22 ) N-M
280: * M M+P-N N-M
281: *
282: * where R11 and T22 are upper triangular, and Q and Z are
283: * orthogonal.
284: *
285: CALL DGGQRF( N, M, P, A, LDA, WORK, B, LDB, WORK( M+1 ),
286: $ WORK( M+NP+1 ), LWORK-M-NP, INFO )
287: LOPT = WORK( M+NP+1 )
288: *
289: * Update left-hand-side vector d = Q**T*d = ( d1 ) M
290: * ( d2 ) N-M
291: *
292: CALL DORMQR( 'Left', 'Transpose', N, 1, M, A, LDA, WORK, D,
293: $ MAX( 1, N ), WORK( M+NP+1 ), LWORK-M-NP, INFO )
294: LOPT = MAX( LOPT, INT( WORK( M+NP+1 ) ) )
295: *
296: * Solve T22*y2 = d2 for y2
297: *
298: IF( N.GT.M ) THEN
299: CALL DTRTRS( 'Upper', 'No transpose', 'Non unit', N-M, 1,
300: $ B( M+1, M+P-N+1 ), LDB, D( M+1 ), N-M, INFO )
301: *
302: IF( INFO.GT.0 ) THEN
303: INFO = 1
304: RETURN
305: END IF
306: *
307: CALL DCOPY( N-M, D( M+1 ), 1, Y( M+P-N+1 ), 1 )
308: END IF
309: *
310: * Set y1 = 0
311: *
312: DO 10 I = 1, M + P - N
313: Y( I ) = ZERO
314: 10 CONTINUE
315: *
316: * Update d1 = d1 - T12*y2
317: *
318: CALL DGEMV( 'No transpose', M, N-M, -ONE, B( 1, M+P-N+1 ), LDB,
319: $ Y( M+P-N+1 ), 1, ONE, D, 1 )
320: *
321: * Solve triangular system: R11*x = d1
322: *
323: IF( M.GT.0 ) THEN
324: CALL DTRTRS( 'Upper', 'No Transpose', 'Non unit', M, 1, A, LDA,
325: $ D, M, INFO )
326: *
327: IF( INFO.GT.0 ) THEN
328: INFO = 2
329: RETURN
330: END IF
331: *
332: * Copy D to X
333: *
334: CALL DCOPY( M, D, 1, X, 1 )
335: END IF
336: *
337: * Backward transformation y = Z**T *y
338: *
339: CALL DORMRQ( 'Left', 'Transpose', P, 1, NP,
340: $ B( MAX( 1, N-P+1 ), 1 ), LDB, WORK( M+1 ), Y,
341: $ MAX( 1, P ), WORK( M+NP+1 ), LWORK-M-NP, INFO )
342: WORK( 1 ) = M + NP + MAX( LOPT, INT( WORK( M+NP+1 ) ) )
343: *
344: RETURN
345: *
346: * End of DGGGLM
347: *
348: END
CVSweb interface <joel.bertrand@systella.fr>