LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
sorgtsqr.f
Go to the documentation of this file.
1*> \brief \b SORGTSQR
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download SORGTSQR + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/sorgtsqr.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/sorgtsqr.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/sorgtsqr.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE SORGTSQR( M, N, MB, NB, A, LDA, T, LDT, WORK, LWORK,
20* $ INFO )
21*
22* .. Scalar Arguments ..
23* INTEGER INFO, LDA, LDT, LWORK, M, N, MB, NB
24* ..
25* .. Array Arguments ..
26* REAL A( LDA, * ), T( LDT, * ), WORK( * )
27* ..
28*
29*> \par Purpose:
30* =============
31*>
32*> \verbatim
33*>
34*> SORGTSQR generates an M-by-N real matrix Q_out with orthonormal columns,
35*> which are the first N columns of a product of real orthogonal
36*> matrices of order M which are returned by SLATSQR
37*>
38*> Q_out = first_N_columns_of( Q(1)_in * Q(2)_in * ... * Q(k)_in ).
39*>
40*> See the documentation for SLATSQR.
41*> \endverbatim
42*
43* Arguments:
44* ==========
45*
46*> \param[in] M
47*> \verbatim
48*> M is INTEGER
49*> The number of rows of the matrix A. M >= 0.
50*> \endverbatim
51*>
52*> \param[in] N
53*> \verbatim
54*> N is INTEGER
55*> The number of columns of the matrix A. M >= N >= 0.
56*> \endverbatim
57*>
58*> \param[in] MB
59*> \verbatim
60*> MB is INTEGER
61*> The row block size used by SLATSQR to return
62*> arrays A and T. MB > N.
63*> (Note that if MB > M, then M is used instead of MB
64*> as the row block size).
65*> \endverbatim
66*>
67*> \param[in] NB
68*> \verbatim
69*> NB is INTEGER
70*> The column block size used by SLATSQR to return
71*> arrays A and T. NB >= 1.
72*> (Note that if NB > N, then N is used instead of NB
73*> as the column block size).
74*> \endverbatim
75*>
76*> \param[in,out] A
77*> \verbatim
78*> A is REAL array, dimension (LDA,N)
79*>
80*> On entry:
81*>
82*> The elements on and above the diagonal are not accessed.
83*> The elements below the diagonal represent the unit
84*> lower-trapezoidal blocked matrix V computed by SLATSQR
85*> that defines the input matrices Q_in(k) (ones on the
86*> diagonal are not stored) (same format as the output A
87*> below the diagonal in SLATSQR).
88*>
89*> On exit:
90*>
91*> The array A contains an M-by-N orthonormal matrix Q_out,
92*> i.e the columns of A are orthogonal unit vectors.
93*> \endverbatim
94*>
95*> \param[in] LDA
96*> \verbatim
97*> LDA is INTEGER
98*> The leading dimension of the array A. LDA >= max(1,M).
99*> \endverbatim
100*>
101*> \param[in] T
102*> \verbatim
103*> T is REAL array,
104*> dimension (LDT, N * NIRB)
105*> where NIRB = Number_of_input_row_blocks
106*> = MAX( 1, CEIL((M-N)/(MB-N)) )
107*> Let NICB = Number_of_input_col_blocks
108*> = CEIL(N/NB)
109*>
110*> The upper-triangular block reflectors used to define the
111*> input matrices Q_in(k), k=(1:NIRB*NICB). The block
112*> reflectors are stored in compact form in NIRB block
113*> reflector sequences. Each of NIRB block reflector sequences
114*> is stored in a larger NB-by-N column block of T and consists
115*> of NICB smaller NB-by-NB upper-triangular column blocks.
116*> (same format as the output T in SLATSQR).
117*> \endverbatim
118*>
119*> \param[in] LDT
120*> \verbatim
121*> LDT is INTEGER
122*> The leading dimension of the array T.
123*> LDT >= max(1,min(NB1,N)).
124*> \endverbatim
125*>
126*> \param[out] WORK
127*> \verbatim
128*> (workspace) REAL array, dimension (MAX(2,LWORK))
129*> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
130*> \endverbatim
131*>
132*> \param[in] LWORK
133*> \verbatim
134*> LWORK is INTEGER
135*> The dimension of the array WORK. LWORK >= (M+NB)*N.
136*> If LWORK = -1, then a workspace query is assumed.
137*> The routine only calculates the optimal size of the WORK
138*> array, returns this value as the first entry of the WORK
139*> array, and no error message related to LWORK is issued
140*> by XERBLA.
141*> \endverbatim
142*>
143*> \param[out] INFO
144*> \verbatim
145*> INFO is INTEGER
146*> = 0: successful exit
147*> < 0: if INFO = -i, the i-th argument had an illegal value
148*> \endverbatim
149*>
150* Authors:
151* ========
152*
153*> \author Univ. of Tennessee
154*> \author Univ. of California Berkeley
155*> \author Univ. of Colorado Denver
156*> \author NAG Ltd.
157*
158*> \ingroup ungtsqr
159*
160*> \par Contributors:
161* ==================
162*>
163*> \verbatim
164*>
165*> November 2019, Igor Kozachenko,
166*> Computer Science Division,
167*> University of California, Berkeley
168*>
169*> \endverbatim
170*
171* =====================================================================
172 SUBROUTINE sorgtsqr( M, N, MB, NB, A, LDA, T, LDT, WORK, LWORK,
173 $ INFO )
174 IMPLICIT NONE
175*
176* -- LAPACK computational routine --
177* -- LAPACK is a software package provided by Univ. of Tennessee, --
178* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
179*
180* .. Scalar Arguments ..
181 INTEGER INFO, LDA, LDT, LWORK, M, N, MB, NB
182* ..
183* .. Array Arguments ..
184 REAL A( LDA, * ), T( LDT, * ), WORK( * )
185* ..
186*
187* =====================================================================
188*
189* .. Parameters ..
190 REAL ONE, ZERO
191 parameter( one = 1.0e+0, zero = 0.0e+0 )
192* ..
193* .. Local Scalars ..
194 LOGICAL LQUERY
195 INTEGER IINFO, LDC, LWORKOPT, LC, LW, NBLOCAL, J
196* ..
197* .. External Functions ..
198 REAL SROUNDUP_LWORK
199 EXTERNAL sroundup_lwork
200* ..
201* .. External Subroutines ..
202 EXTERNAL scopy, slamtsqr, slaset, xerbla
203* ..
204* .. Intrinsic Functions ..
205 INTRINSIC max, min
206* ..
207* .. Executable Statements ..
208*
209* Test the input parameters
210*
211 lquery = lwork.EQ.-1
212 info = 0
213 IF( m.LT.0 ) THEN
214 info = -1
215 ELSE IF( n.LT.0 .OR. m.LT.n ) THEN
216 info = -2
217 ELSE IF( mb.LE.n ) THEN
218 info = -3
219 ELSE IF( nb.LT.1 ) THEN
220 info = -4
221 ELSE IF( lda.LT.max( 1, m ) ) THEN
222 info = -6
223 ELSE IF( ldt.LT.max( 1, min( nb, n ) ) ) THEN
224 info = -8
225 ELSE
226*
227* Test the input LWORK for the dimension of the array WORK.
228* This workspace is used to store array C(LDC, N) and WORK(LWORK)
229* in the call to SLAMTSQR. See the documentation for SLAMTSQR.
230*
231 IF( lwork.LT.2 .AND. (.NOT.lquery) ) THEN
232 info = -10
233 ELSE
234*
235* Set block size for column blocks
236*
237 nblocal = min( nb, n )
238*
239* LWORK = -1, then set the size for the array C(LDC,N)
240* in SLAMTSQR call and set the optimal size of the work array
241* WORK(LWORK) in SLAMTSQR call.
242*
243 ldc = m
244 lc = ldc*n
245 lw = n * nblocal
246*
247 lworkopt = lc+lw
248*
249 IF( ( lwork.LT.max( 1, lworkopt ) ).AND.(.NOT.lquery) ) THEN
250 info = -10
251 END IF
252 END IF
253*
254 END IF
255*
256* Handle error in the input parameters and return workspace query.
257*
258 IF( info.NE.0 ) THEN
259 CALL xerbla( 'SORGTSQR', -info )
260 RETURN
261 ELSE IF ( lquery ) THEN
262 work( 1 ) = sroundup_lwork( lworkopt )
263 RETURN
264 END IF
265*
266* Quick return if possible
267*
268 IF( min( m, n ).EQ.0 ) THEN
269 work( 1 ) = sroundup_lwork( lworkopt )
270 RETURN
271 END IF
272*
273* (1) Form explicitly the tall-skinny M-by-N left submatrix Q1_in
274* of M-by-M orthogonal matrix Q_in, which is implicitly stored in
275* the subdiagonal part of input array A and in the input array T.
276* Perform by the following operation using the routine SLAMTSQR.
277*
278* Q1_in = Q_in * ( I ), where I is a N-by-N identity matrix,
279* ( 0 ) 0 is a (M-N)-by-N zero matrix.
280*
281* (1a) Form M-by-N matrix in the array WORK(1:LDC*N) with ones
282* on the diagonal and zeros elsewhere.
283*
284 CALL slaset( 'F', m, n, zero, one, work, ldc )
285*
286* (1b) On input, WORK(1:LDC*N) stores ( I );
287* ( 0 )
288*
289* On output, WORK(1:LDC*N) stores Q1_in.
290*
291 CALL slamtsqr( 'L', 'N', m, n, n, mb, nblocal, a, lda, t, ldt,
292 $ work, ldc, work( lc+1 ), lw, iinfo )
293*
294* (2) Copy the result from the part of the work array (1:M,1:N)
295* with the leading dimension LDC that starts at WORK(1) into
296* the output array A(1:M,1:N) column-by-column.
297*
298 DO j = 1, n
299 CALL scopy( m, work( (j-1)*ldc + 1 ), 1, a( 1, j ), 1 )
300 END DO
301*
302 work( 1 ) = sroundup_lwork( lworkopt )
303 RETURN
304*
305* End of SORGTSQR
306*
307 END
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:201
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:108
subroutine sorgtsqr(m, n, mb, nb, a, lda, t, ldt, work, lwork, info)
SORGTSQR
Definition sorgtsqr.f:174