LAPACK  3.9.0
LAPACK: Linear Algebra PACKage

◆ 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] \par Purpose: @verbatim 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. \endverbatim @param [in] M @verbatim M is INTEGER The number of rows of the matrix A. M >= 0. \endverbatim @param [in] N @verbatim N is INTEGER The number of columns of the matrix A. M >= N >= 0. \endverbatim @param [in] MB @verbatim 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). \endverbatim @param [in] NB @verbatim 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). \endverbatim @param [in,out] A @verbatim 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. \endverbatim @param [in] LDA @verbatim LDA is INTEGER The leading dimension of the array A. LDA >= max(1,M). \endverbatim @param [in] T @verbatim 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). \endverbatim @param [in] LDT @verbatim LDT is INTEGER The leading dimension of the array T. LDT >= max(1,min(NB1,N)). \endverbatim @param [out] WORK @verbatim (workspace) REAL array, dimension (MAX(2,LWORK)) On exit, if INFO = 0, WORK(1) returns the optimal LWORK. \endverbatim @param [in] LWORK @verbatim 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. \endverbatim @param [out] INFO @verbatim INFO is INTEGER = 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value \endverbatim \author Univ. of Tennessee \author Univ. of California Berkeley \author Univ. of Colorado Denver \author NAG Ltd. \date November 2019 \par Contributors: @verbatim November 2019, Igor Kozachenko, Computer Science Division, University of California, Berkeley \endverbatim

Definition at line 176 of file sorgtsqr.f.

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