LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine dgglse ( integer  M,
integer  N,
integer  P,
double precision, dimension( lda, * )  A,
integer  LDA,
double precision, dimension( ldb, * )  B,
integer  LDB,
double precision, dimension( * )  C,
double precision, dimension( * )  D,
double precision, dimension( * )  X,
double precision, dimension( * )  WORK,
integer  LWORK,
integer  INFO 
)

DGGLSE solves overdetermined or underdetermined systems for OTHER matrices

Download DGGLSE + dependencies [TGZ] [ZIP] [TXT]

Purpose:
 DGGLSE solves the linear equality-constrained least squares (LSE)
 problem:

         minimize || c - A*x ||_2   subject to   B*x = d

 where A is an M-by-N matrix, B is a P-by-N matrix, c is a given
 M-vector, and d is a given P-vector. It is assumed that
 P <= N <= M+P, and

          rank(B) = P and  rank( (A) ) = N.
                               ( (B) )

 These conditions ensure that the LSE problem has a unique solution,
 which is obtained using a generalized RQ factorization of the
 matrices (B, A) given by

    B = (0 R)*Q,   A = Z*T*Q.
Parameters
[in]M
          M is INTEGER
          The number of rows of the matrix A.  M >= 0.
[in]N
          N is INTEGER
          The number of columns of the matrices A and B. N >= 0.
[in]P
          P is INTEGER
          The number of rows of the matrix B. 0 <= P <= N <= M+P.
[in,out]A
          A is DOUBLE PRECISION array, dimension (LDA,N)
          On entry, the M-by-N matrix A.
          On exit, the elements on and above the diagonal of the array
          contain the min(M,N)-by-N upper trapezoidal matrix T.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A. LDA >= max(1,M).
[in,out]B
          B is DOUBLE PRECISION array, dimension (LDB,N)
          On entry, the P-by-N matrix B.
          On exit, the upper triangle of the subarray B(1:P,N-P+1:N)
          contains the P-by-P upper triangular matrix R.
[in]LDB
          LDB is INTEGER
          The leading dimension of the array B. LDB >= max(1,P).
[in,out]C
          C is DOUBLE PRECISION array, dimension (M)
          On entry, C contains the right hand side vector for the
          least squares part of the LSE problem.
          On exit, the residual sum of squares for the solution
          is given by the sum of squares of elements N-P+1 to M of
          vector C.
[in,out]D
          D is DOUBLE PRECISION array, dimension (P)
          On entry, D contains the right hand side vector for the
          constrained equation.
          On exit, D is destroyed.
[out]X
          X is DOUBLE PRECISION array, dimension (N)
          On exit, X is the solution of the LSE problem.
[out]WORK
          WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
[in]LWORK
          LWORK is INTEGER
          The dimension of the array WORK. LWORK >= max(1,M+N+P).
          For optimum performance LWORK >= P+min(M,N)+max(M,N)*NB,
          where NB is an upper bound for the optimal blocksizes for
          DGEQRF, SGERQF, DORMQR and SORMRQ.

          If LWORK = -1, then a workspace query is assumed; the routine
          only calculates the optimal size of the WORK array, returns
          this value as the first entry of the WORK array, and no error
          message related to LWORK is issued by XERBLA.
[out]INFO
          INFO is INTEGER
          = 0:  successful exit.
          < 0:  if INFO = -i, the i-th argument had an illegal value.
          = 1:  the upper triangular factor R associated with B in the
                generalized RQ factorization of the pair (B, A) is
                singular, so that rank(B) < P; the least squares
                solution could not be computed.
          = 2:  the (N-P) by (N-P) part of the upper trapezoidal factor
                T associated with A in the generalized RQ factorization
                of the pair (B, A) is singular, so that
                rank( (A) ) < N; the least squares solution could not
                    ( (B) )
                be computed.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 182 of file dgglse.f.

182 *
183 * -- LAPACK driver routine (version 3.4.0) --
184 * -- LAPACK is a software package provided by Univ. of Tennessee, --
185 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
186 * November 2011
187 *
188 * .. Scalar Arguments ..
189  INTEGER info, lda, ldb, lwork, m, n, p
190 * ..
191 * .. Array Arguments ..
192  DOUBLE PRECISION a( lda, * ), b( ldb, * ), c( * ), d( * ),
193  $ work( * ), x( * )
194 * ..
195 *
196 * =====================================================================
197 *
198 * .. Parameters ..
199  DOUBLE PRECISION one
200  parameter ( one = 1.0d+0 )
201 * ..
202 * .. Local Scalars ..
203  LOGICAL lquery
204  INTEGER lopt, lwkmin, lwkopt, mn, nb, nb1, nb2, nb3,
205  $ nb4, nr
206 * ..
207 * .. External Subroutines ..
208  EXTERNAL daxpy, dcopy, dgemv, dggrqf, dormqr, dormrq,
209  $ dtrmv, dtrtrs, xerbla
210 * ..
211 * .. External Functions ..
212  INTEGER ilaenv
213  EXTERNAL ilaenv
214 * ..
215 * .. Intrinsic Functions ..
216  INTRINSIC int, max, min
217 * ..
218 * .. Executable Statements ..
219 *
220 * Test the input parameters
221 *
222  info = 0
223  mn = min( m, n )
224  lquery = ( lwork.EQ.-1 )
225  IF( m.LT.0 ) THEN
226  info = -1
227  ELSE IF( n.LT.0 ) THEN
228  info = -2
229  ELSE IF( p.LT.0 .OR. p.GT.n .OR. p.LT.n-m ) THEN
230  info = -3
231  ELSE IF( lda.LT.max( 1, m ) ) THEN
232  info = -5
233  ELSE IF( ldb.LT.max( 1, p ) ) THEN
234  info = -7
235  END IF
236 *
237 * Calculate workspace
238 *
239  IF( info.EQ.0) THEN
240  IF( n.EQ.0 ) THEN
241  lwkmin = 1
242  lwkopt = 1
243  ELSE
244  nb1 = ilaenv( 1, 'DGEQRF', ' ', m, n, -1, -1 )
245  nb2 = ilaenv( 1, 'DGERQF', ' ', m, n, -1, -1 )
246  nb3 = ilaenv( 1, 'DORMQR', ' ', m, n, p, -1 )
247  nb4 = ilaenv( 1, 'DORMRQ', ' ', m, n, p, -1 )
248  nb = max( nb1, nb2, nb3, nb4 )
249  lwkmin = m + n + p
250  lwkopt = p + mn + max( m, n )*nb
251  END IF
252  work( 1 ) = lwkopt
253 *
254  IF( lwork.LT.lwkmin .AND. .NOT.lquery ) THEN
255  info = -12
256  END IF
257  END IF
258 *
259  IF( info.NE.0 ) THEN
260  CALL xerbla( 'DGGLSE', -info )
261  RETURN
262  ELSE IF( lquery ) THEN
263  RETURN
264  END IF
265 *
266 * Quick return if possible
267 *
268  IF( n.EQ.0 )
269  $ RETURN
270 *
271 * Compute the GRQ factorization of matrices B and A:
272 *
273 * B*Q**T = ( 0 T12 ) P Z**T*A*Q**T = ( R11 R12 ) N-P
274 * N-P P ( 0 R22 ) M+P-N
275 * N-P P
276 *
277 * where T12 and R11 are upper triangular, and Q and Z are
278 * orthogonal.
279 *
280  CALL dggrqf( p, m, n, b, ldb, work, a, lda, work( p+1 ),
281  $ work( p+mn+1 ), lwork-p-mn, info )
282  lopt = work( p+mn+1 )
283 *
284 * Update c = Z**T *c = ( c1 ) N-P
285 * ( c2 ) M+P-N
286 *
287  CALL dormqr( 'Left', 'Transpose', m, 1, mn, a, lda, work( p+1 ),
288  $ c, max( 1, m ), work( p+mn+1 ), lwork-p-mn, info )
289  lopt = max( lopt, int( work( p+mn+1 ) ) )
290 *
291 * Solve T12*x2 = d for x2
292 *
293  IF( p.GT.0 ) THEN
294  CALL dtrtrs( 'Upper', 'No transpose', 'Non-unit', p, 1,
295  $ b( 1, n-p+1 ), ldb, d, p, info )
296 *
297  IF( info.GT.0 ) THEN
298  info = 1
299  RETURN
300  END IF
301 *
302 * Put the solution in X
303 *
304  CALL dcopy( p, d, 1, x( n-p+1 ), 1 )
305 *
306 * Update c1
307 *
308  CALL dgemv( 'No transpose', n-p, p, -one, a( 1, n-p+1 ), lda,
309  $ d, 1, one, c, 1 )
310  END IF
311 *
312 * Solve R11*x1 = c1 for x1
313 *
314  IF( n.GT.p ) THEN
315  CALL dtrtrs( 'Upper', 'No transpose', 'Non-unit', n-p, 1,
316  $ a, lda, c, n-p, info )
317 *
318  IF( info.GT.0 ) THEN
319  info = 2
320  RETURN
321  END IF
322 *
323 * Put the solutions in X
324 *
325  CALL dcopy( n-p, c, 1, x, 1 )
326  END IF
327 *
328 * Compute the residual vector:
329 *
330  IF( m.LT.n ) THEN
331  nr = m + p - n
332  IF( nr.GT.0 )
333  $ CALL dgemv( 'No transpose', nr, n-m, -one, a( n-p+1, m+1 ),
334  $ lda, d( nr+1 ), 1, one, c( n-p+1 ), 1 )
335  ELSE
336  nr = p
337  END IF
338  IF( nr.GT.0 ) THEN
339  CALL dtrmv( 'Upper', 'No transpose', 'Non unit', nr,
340  $ a( n-p+1, n-p+1 ), lda, d, 1 )
341  CALL daxpy( nr, -one, d, 1, c( n-p+1 ), 1 )
342  END IF
343 *
344 * Backward transformation x = Q**T*x
345 *
346  CALL dormrq( 'Left', 'Transpose', n, 1, p, b, ldb, work( 1 ), x,
347  $ n, work( p+mn+1 ), lwork-p-mn, info )
348  work( 1 ) = p + mn + max( lopt, int( work( p+mn+1 ) ) )
349 *
350  RETURN
351 *
352 * End of DGGLSE
353 *
subroutine dcopy(N, DX, INCX, DY, INCY)
DCOPY
Definition: dcopy.f:53
subroutine dtrtrs(UPLO, TRANS, DIAG, N, NRHS, A, LDA, B, LDB, INFO)
DTRTRS
Definition: dtrtrs.f:142
subroutine dgemv(TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
DGEMV
Definition: dgemv.f:158
subroutine daxpy(N, DA, DX, INCX, DY, INCY)
DAXPY
Definition: daxpy.f:54
subroutine dormqr(SIDE, TRANS, M, N, K, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
DORMQR
Definition: dormqr.f:169
subroutine dormrq(SIDE, TRANS, M, N, K, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
DORMRQ
Definition: dormrq.f:169
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
Definition: tstiee.f:83
subroutine dggrqf(M, P, N, A, LDA, TAUA, B, LDB, TAUB, WORK, LWORK, INFO)
DGGRQF
Definition: dggrqf.f:216
subroutine dtrmv(UPLO, TRANS, DIAG, N, A, LDA, X, INCX)
DTRMV
Definition: dtrmv.f:149

Here is the call graph for this function:

Here is the caller graph for this function: