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

CGGLSE solves overdetermined or underdetermined systems for OTHER matrices

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

Purpose:
 CGGLSE 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 COMPLEX 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 COMPLEX 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 COMPLEX 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 COMPLEX 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 COMPLEX array, dimension (N)
          On exit, X is the solution of the LSE problem.
[out]WORK
          WORK is COMPLEX 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
          CGEQRF, CGERQF, CUNMQR and CUNMRQ.

          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 cgglse.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  COMPLEX a( lda, * ), b( ldb, * ), c( * ), d( * ),
193  $ work( * ), x( * )
194 * ..
195 *
196 * =====================================================================
197 *
198 * .. Parameters ..
199  COMPLEX cone
200  parameter ( cone = ( 1.0e+0, 0.0e+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 caxpy, ccopy, cgemv, cggrqf, ctrmv, ctrtrs,
209  $ cunmqr, cunmrq, 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, 'CGEQRF', ' ', m, n, -1, -1 )
245  nb2 = ilaenv( 1, 'CGERQF', ' ', m, n, -1, -1 )
246  nb3 = ilaenv( 1, 'CUNMQR', ' ', m, n, p, -1 )
247  nb4 = ilaenv( 1, 'CUNMRQ', ' ', 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( 'CGGLSE', -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**H = ( 0 T12 ) P Z**H*A*Q**H = ( 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 * unitary.
279 *
280  CALL cggrqf( 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**H *c = ( c1 ) N-P
285 * ( c2 ) M+P-N
286 *
287  CALL cunmqr( 'Left', 'Conjugate Transpose', m, 1, mn, a, lda,
288  $ work( p+1 ), c, max( 1, m ), work( p+mn+1 ),
289  $ lwork-p-mn, info )
290  lopt = max( lopt, int( work( p+mn+1 ) ) )
291 *
292 * Solve T12*x2 = d for x2
293 *
294  IF( p.GT.0 ) THEN
295  CALL ctrtrs( 'Upper', 'No transpose', 'Non-unit', p, 1,
296  $ b( 1, n-p+1 ), ldb, d, p, info )
297 *
298  IF( info.GT.0 ) THEN
299  info = 1
300  RETURN
301  END IF
302 *
303 * Put the solution in X
304 *
305  CALL ccopy( p, d, 1, x( n-p+1 ), 1 )
306 *
307 * Update c1
308 *
309  CALL cgemv( 'No transpose', n-p, p, -cone, a( 1, n-p+1 ), lda,
310  $ d, 1, cone, c, 1 )
311  END IF
312 *
313 * Solve R11*x1 = c1 for x1
314 *
315  IF( n.GT.p ) THEN
316  CALL ctrtrs( 'Upper', 'No transpose', 'Non-unit', n-p, 1,
317  $ a, lda, c, n-p, info )
318 *
319  IF( info.GT.0 ) THEN
320  info = 2
321  RETURN
322  END IF
323 *
324 * Put the solutions in X
325 *
326  CALL ccopy( n-p, c, 1, x, 1 )
327  END IF
328 *
329 * Compute the residual vector:
330 *
331  IF( m.LT.n ) THEN
332  nr = m + p - n
333  IF( nr.GT.0 )
334  $ CALL cgemv( 'No transpose', nr, n-m, -cone, a( n-p+1, m+1 ),
335  $ lda, d( nr+1 ), 1, cone, c( n-p+1 ), 1 )
336  ELSE
337  nr = p
338  END IF
339  IF( nr.GT.0 ) THEN
340  CALL ctrmv( 'Upper', 'No transpose', 'Non unit', nr,
341  $ a( n-p+1, n-p+1 ), lda, d, 1 )
342  CALL caxpy( nr, -cone, d, 1, c( n-p+1 ), 1 )
343  END IF
344 *
345 * Backward transformation x = Q**H*x
346 *
347  CALL cunmrq( 'Left', 'Conjugate Transpose', n, 1, p, b, ldb,
348  $ work( 1 ), x, n, work( p+mn+1 ), lwork-p-mn, info )
349  work( 1 ) = p + mn + max( lopt, int( work( p+mn+1 ) ) )
350 *
351  RETURN
352 *
353 * End of CGGLSE
354 *
subroutine cggrqf(M, P, N, A, LDA, TAUA, B, LDB, TAUB, WORK, LWORK, INFO)
CGGRQF
Definition: cggrqf.f:216
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine cunmrq(SIDE, TRANS, M, N, K, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
CUNMRQ
Definition: cunmrq.f:170
subroutine cunmqr(SIDE, TRANS, M, N, K, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
CUNMQR
Definition: cunmqr.f:170
subroutine cgemv(TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
CGEMV
Definition: cgemv.f:160
subroutine ctrtrs(UPLO, TRANS, DIAG, N, NRHS, A, LDA, B, LDB, INFO)
CTRTRS
Definition: ctrtrs.f:142
subroutine ctrmv(UPLO, TRANS, DIAG, N, A, LDA, X, INCX)
CTRMV
Definition: ctrmv.f:149
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
Definition: tstiee.f:83
subroutine ccopy(N, CX, INCX, CY, INCY)
CCOPY
Definition: ccopy.f:52
subroutine caxpy(N, CA, CX, INCX, CY, INCY)
CAXPY
Definition: caxpy.f:53

Here is the call graph for this function:

Here is the caller graph for this function: