LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ sorgtsqr()

subroutine sorgtsqr ( integer  m,
integer  n,
integer  mb,
integer  nb,
real, dimension( lda, * )  a,
integer  lda,
real, dimension( ldt, * )  t,
integer  ldt,
real, dimension( * )  work,
integer  lwork,
integer  info 
)

SORGTSQR

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

Purpose:
 SORGTSQR generates an M-by-N real matrix Q_out with orthonormal columns,
 which are the first N columns of a product of real orthogonal
 matrices of order M which are returned by SLATSQR

      Q_out = first_N_columns_of( Q(1)_in * Q(2)_in * ... * Q(k)_in ).

 See the documentation for SLATSQR.
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 matrix A. M >= N >= 0.
[in]MB
          MB is INTEGER
          The row block size used by SLATSQR to return
          arrays A and T. MB > N.
          (Note that if MB > M, then M is used instead of MB
          as the row block size).
[in]NB
          NB is INTEGER
          The column block size used by SLATSQR to return
          arrays A and T. NB >= 1.
          (Note that if NB > N, then N is used instead of NB
          as the column block size).
[in,out]A
          A is REAL array, dimension (LDA,N)

          On entry:

             The elements on and above the diagonal are not accessed.
             The elements below the diagonal represent the unit
             lower-trapezoidal blocked matrix V computed by SLATSQR
             that defines the input matrices Q_in(k) (ones on the
             diagonal are not stored) (same format as the output A
             below the diagonal in SLATSQR).

          On exit:

             The array A contains an M-by-N orthonormal matrix Q_out,
             i.e the columns of A are orthogonal unit vectors.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,M).
[in]T
          T is REAL array,
          dimension (LDT, N * NIRB)
          where NIRB = Number_of_input_row_blocks
                     = MAX( 1, CEIL((M-N)/(MB-N)) )
          Let NICB = Number_of_input_col_blocks
                   = CEIL(N/NB)

          The upper-triangular block reflectors used to define the
          input matrices Q_in(k), k=(1:NIRB*NICB). The block
          reflectors are stored in compact form in NIRB block
          reflector sequences. Each of NIRB block reflector sequences
          is stored in a larger NB-by-N column block of T and consists
          of NICB smaller NB-by-NB upper-triangular column blocks.
          (same format as the output T in SLATSQR).
[in]LDT
          LDT is INTEGER
          The leading dimension of the array T.
          LDT >= max(1,min(NB1,N)).
[out]WORK
          (workspace) REAL array, dimension (MAX(2,LWORK))
          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
[in]LWORK
          LWORK is INTEGER
          The dimension of the array WORK.  LWORK >= (M+NB)*N.
          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
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Contributors:
 November 2019, Igor Kozachenko,
                Computer Science Division,
                University of California, Berkeley

Definition at line 174 of file sorgtsqr.f.

176 IMPLICIT NONE
177*
178* -- LAPACK computational routine --
179* -- LAPACK is a software package provided by Univ. of Tennessee, --
180* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
181*
182* .. Scalar Arguments ..
183 INTEGER INFO, LDA, LDT, LWORK, M, N, MB, NB
184* ..
185* .. Array Arguments ..
186 REAL A( LDA, * ), T( LDT, * ), WORK( * )
187* ..
188*
189* =====================================================================
190*
191* .. Parameters ..
192 REAL ONE, ZERO
193 parameter( one = 1.0e+0, zero = 0.0e+0 )
194* ..
195* .. Local Scalars ..
196 LOGICAL LQUERY
197 INTEGER IINFO, LDC, LWORKOPT, LC, LW, NBLOCAL, J
198* ..
199* .. External Functions ..
200 REAL SROUNDUP_LWORK
201 EXTERNAL sroundup_lwork
202* ..
203* .. External Subroutines ..
204 EXTERNAL scopy, slamtsqr, slaset, xerbla
205* ..
206* .. Intrinsic Functions ..
207 INTRINSIC max, min
208* ..
209* .. Executable Statements ..
210*
211* Test the input parameters
212*
213 lquery = lwork.EQ.-1
214 info = 0
215 IF( m.LT.0 ) THEN
216 info = -1
217 ELSE IF( n.LT.0 .OR. m.LT.n ) THEN
218 info = -2
219 ELSE IF( mb.LE.n ) THEN
220 info = -3
221 ELSE IF( nb.LT.1 ) THEN
222 info = -4
223 ELSE IF( lda.LT.max( 1, m ) ) THEN
224 info = -6
225 ELSE IF( ldt.LT.max( 1, min( nb, n ) ) ) THEN
226 info = -8
227 ELSE
228*
229* Test the input LWORK for the dimension of the array WORK.
230* This workspace is used to store array C(LDC, N) and WORK(LWORK)
231* in the call to SLAMTSQR. See the documentation for SLAMTSQR.
232*
233 IF( lwork.LT.2 .AND. (.NOT.lquery) ) THEN
234 info = -10
235 ELSE
236*
237* Set block size for column blocks
238*
239 nblocal = min( nb, n )
240*
241* LWORK = -1, then set the size for the array C(LDC,N)
242* in SLAMTSQR call and set the optimal size of the work array
243* WORK(LWORK) in SLAMTSQR call.
244*
245 ldc = m
246 lc = ldc*n
247 lw = n * nblocal
248*
249 lworkopt = lc+lw
250*
251 IF( ( lwork.LT.max( 1, lworkopt ) ).AND.(.NOT.lquery) ) THEN
252 info = -10
253 END IF
254 END IF
255*
256 END IF
257*
258* Handle error in the input parameters and return workspace query.
259*
260 IF( info.NE.0 ) THEN
261 CALL xerbla( 'SORGTSQR', -info )
262 RETURN
263 ELSE IF ( lquery ) THEN
264 work( 1 ) = sroundup_lwork( lworkopt )
265 RETURN
266 END IF
267*
268* Quick return if possible
269*
270 IF( min( m, n ).EQ.0 ) THEN
271 work( 1 ) = sroundup_lwork( lworkopt )
272 RETURN
273 END IF
274*
275* (1) Form explicitly the tall-skinny M-by-N left submatrix Q1_in
276* of M-by-M orthogonal matrix Q_in, which is implicitly stored in
277* the subdiagonal part of input array A and in the input array T.
278* Perform by the following operation using the routine SLAMTSQR.
279*
280* Q1_in = Q_in * ( I ), where I is a N-by-N identity matrix,
281* ( 0 ) 0 is a (M-N)-by-N zero matrix.
282*
283* (1a) Form M-by-N matrix in the array WORK(1:LDC*N) with ones
284* on the diagonal and zeros elsewhere.
285*
286 CALL slaset( 'F', m, n, zero, one, work, ldc )
287*
288* (1b) On input, WORK(1:LDC*N) stores ( I );
289* ( 0 )
290*
291* On output, WORK(1:LDC*N) stores Q1_in.
292*
293 CALL slamtsqr( 'L', 'N', m, n, n, mb, nblocal, a, lda, t, ldt,
294 $ work, ldc, work( lc+1 ), lw, iinfo )
295*
296* (2) Copy the result from the part of the work array (1:M,1:N)
297* with the leading dimension LDC that starts at WORK(1) into
298* the output array A(1:M,1:N) column-by-column.
299*
300 DO j = 1, n
301 CALL scopy( m, work( (j-1)*ldc + 1 ), 1, a( 1, j ), 1 )
302 END DO
303*
304 work( 1 ) = sroundup_lwork( lworkopt )
305 RETURN
306*
307* End of SORGTSQR
308*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine scopy(n, sx, incx, sy, incy)
SCOPY
Definition scopy.f:82
subroutine slamtsqr(side, trans, m, n, k, mb, nb, a, lda, t, ldt, c, ldc, work, lwork, info)
SLAMTSQR
Definition slamtsqr.f:199
subroutine slaset(uplo, m, n, alpha, beta, a, lda)
SLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition slaset.f:110
real function sroundup_lwork(lwork)
SROUNDUP_LWORK
Here is the call graph for this function:
Here is the caller graph for this function: