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

SGGGLM

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

Purpose:
 SGGGLM solves a general Gauss-Markov linear model (GLM) problem:

         minimize || y ||_2   subject to   d = A*x + B*y
             x

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

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

 Under these assumptions, the constrained equation is always
 consistent, and there is a unique solution x and a minimal 2-norm
 solution y, which is obtained using a generalized QR factorization
 of the matrices (A, B) given by

    A = Q*(R),   B = Q*T*Z.
          (0)

 In particular, if matrix B is square nonsingular, then the problem
 GLM is equivalent to the following weighted linear least squares
 problem

              minimize || inv(B)*(d-A*x) ||_2
                  x

 where inv(B) denotes the inverse of B.
Parameters
[in]N
          N is INTEGER
          The number of rows of the matrices A and B.  N >= 0.
[in]M
          M is INTEGER
          The number of columns of the matrix A.  0 <= M <= N.
[in]P
          P is INTEGER
          The number of columns of the matrix B.  P >= N-M.
[in,out]A
          A is REAL array, dimension (LDA,M)
          On entry, the N-by-M matrix A.
          On exit, the upper triangular part of the array A contains
          the M-by-M upper triangular matrix R.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A. LDA >= max(1,N).
[in,out]B
          B is REAL array, dimension (LDB,P)
          On entry, the N-by-P matrix B.
          On exit, if N <= P, the upper triangle of the subarray
          B(1:N,P-N+1:P) contains the N-by-N upper triangular matrix T;
          if N > P, the elements on and above the (N-P)th subdiagonal
          contain the N-by-P upper trapezoidal matrix T.
[in]LDB
          LDB is INTEGER
          The leading dimension of the array B. LDB >= max(1,N).
[in,out]D
          D is REAL array, dimension (N)
          On entry, D is the left hand side of the GLM equation.
          On exit, D is destroyed.
[out]X
          X is REAL array, dimension (M)
[out]Y
          Y is REAL array, dimension (P)

          On exit, X and Y are the solutions of the GLM problem.
[out]WORK
          WORK is REAL 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,N+M+P).
          For optimum performance, LWORK >= M+min(N,P)+max(N,P)*NB,
          where NB is an upper bound for the optimal blocksizes for
          SGEQRF, SGERQF, SORMQR 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 A in the
                generalized QR factorization of the pair (A, B) is
                singular, so that rank(A) < M; the least squares
                solution could not be computed.
          = 2:  the bottom (N-M) by (N-M) part of the upper trapezoidal
                factor T associated with B in the generalized QR
                factorization of the pair (A, B) is singular, so that
                rank( A B ) < N; the least squares solution could not
                be computed.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2015

Definition at line 187 of file sggglm.f.

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  REAL a( lda, * ), b( ldb, * ), d( * ), work( * ),
198  $ x( * ), y( * )
199 * ..
200 *
201 * ===================================================================
202 *
203 * .. Parameters ..
204  REAL zero, one
205  parameter ( zero = 0.0e+0, one = 1.0e+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 scopy, sgemv, sggqrf, sormqr, sormrq, strtrs,
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, 'SGEQRF', ' ', n, m, -1, -1 )
250  nb2 = ilaenv( 1, 'SGERQF', ' ', n, m, -1, -1 )
251  nb3 = ilaenv( 1, 'SORMQR', ' ', n, m, p, -1 )
252  nb4 = ilaenv( 1, 'SORMRQ', ' ', 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( 'SGGGLM', -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 sggqrf( 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 sormqr( '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 strtrs( '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 scopy( 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 sgemv( '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 strtrs( '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 scopy( m, d, 1, x, 1 )
335  END IF
336 *
337 * Backward transformation y = Z**T *y
338 *
339  CALL sormrq( '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 SGGGLM
347 *
subroutine sormqr(SIDE, TRANS, M, N, K, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
SORMQR
Definition: sormqr.f:170
subroutine strtrs(UPLO, TRANS, DIAG, N, NRHS, A, LDA, B, LDB, INFO)
STRTRS
Definition: strtrs.f:142
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine sgemv(TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
SGEMV
Definition: sgemv.f:158
subroutine sormrq(SIDE, TRANS, M, N, K, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
SORMRQ
Definition: sormrq.f:170
subroutine sggqrf(N, M, P, A, LDA, TAUA, B, LDB, TAUB, WORK, LWORK, INFO)
SGGQRF
Definition: sggqrf.f:217
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
Definition: tstiee.f:83
subroutine scopy(N, SX, INCX, SY, INCY)
SCOPY
Definition: scopy.f:53

Here is the call graph for this function:

Here is the caller graph for this function: