1*> \brief \b SORHR_COL01
2*
3*  =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6*            http://www.netlib.org/lapack/explore-html/
7*
8*  Definition:
9*  ===========
10*
11*       SUBROUTINE SORHR_COL01( M, N, MB1, NB1, NB2, RESULT )
12*
13*       .. Scalar Arguments ..
14*       INTEGER           M, N, MB1, NB1, NB2
15*       .. Return values ..
16*       REAL              RESULT(6)
17*
18*
19*> \par Purpose:
20*  =============
21*>
22*> \verbatim
23*>
24*> SORHR_COL01 tests SORGTSQR and SORHR_COL using SLATSQR, SGEMQRT.
25*> Therefore, SLATSQR (part of SGEQR), SGEMQRT (part of SGEMQR)
26*> have to be tested before this test.
27*>
28*> \endverbatim
29*
30*  Arguments:
31*  ==========
32*
33*> \param[in] M
34*> \verbatim
35*>          M is INTEGER
36*>          Number of rows in test matrix.
37*> \endverbatim
38*> \param[in] N
39*> \verbatim
40*>          N is INTEGER
41*>          Number of columns in test matrix.
42*> \endverbatim
43*> \param[in] MB1
44*> \verbatim
45*>          MB1 is INTEGER
46*>          Number of row in row block in an input test matrix.
47*> \endverbatim
48*>
49*> \param[in] NB1
50*> \verbatim
51*>          NB1 is INTEGER
52*>          Number of columns in column block an input test matrix.
53*> \endverbatim
54*>
55*> \param[in] NB2
56*> \verbatim
57*>          NB2 is INTEGER
58*>          Number of columns in column block in an output test matrix.
59*> \endverbatim
60*>
61*> \param[out] RESULT
62*> \verbatim
63*>          RESULT is REAL array, dimension (6)
64*>          Results of each of the six tests below.
65*>
66*>            A is a m-by-n test input matrix to be factored.
67*>            so that A = Q_gr * ( R )
68*>                               ( 0 ),
69*>
70*>            Q_qr is an implicit m-by-m orthogonal Q matrix, the result
71*>            of factorization in blocked WY-representation,
72*>            stored in SGEQRT output format.
73*>
74*>            R is a n-by-n upper-triangular matrix,
75*>
76*>            0 is a (m-n)-by-n zero matrix,
77*>
78*>            Q is an explicit m-by-m orthogonal matrix Q = Q_gr * I
79*>
80*>            C is an m-by-n random matrix,
81*>
82*>            D is an n-by-m random matrix.
83*>
84*>          The six tests are:
85*>
86*>          RESULT(1) = |R - (Q**H) * A| / ( eps * m * |A| )
87*>            is equivalent to test for | A - Q * R | / (eps * m * |A|),
88*>
89*>          RESULT(2) = |I - (Q**H) * Q| / ( eps * m ),
90*>
91*>          RESULT(3) = | Q_qr * C - Q * C | / (eps * m * |C|),
92*>
93*>          RESULT(4) = | (Q_gr**H) * C - (Q**H) * C | / (eps * m * |C|)
94*>
95*>          RESULT(5) = | D * Q_qr - D * Q | / (eps * m * |D|)
96*>
97*>          RESULT(6) = | D * (Q_qr**H) - D * (Q**H) | / (eps * m * |D|),
98*>
99*>          where:
100*>            Q_qr * C, (Q_gr**H) * C, D * Q_qr, D * (Q_qr**H) are
101*>            computed using SGEMQRT,
102*>
103*>            Q * C, (Q**H) * C, D * Q, D * (Q**H)  are
104*>            computed using SGEMM.
105*> \endverbatim
106*
107*  Authors:
108*  ========
109*
110*> \author Univ. of Tennessee
111*> \author Univ. of California Berkeley
112*> \author Univ. of Colorado Denver
113*> \author NAG Ltd.
114*
115*> \ingroup single_lin
116*
117*  =====================================================================
118      SUBROUTINE SORHR_COL01( M, N, MB1, NB1, NB2, RESULT )
119      IMPLICIT NONE
120*
121*  -- LAPACK test routine --
122*  -- LAPACK is a software package provided by Univ. of Tennessee,    --
123*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
124*
125*     .. Scalar Arguments ..
126      INTEGER           M, N, MB1, NB1, NB2
127*     .. Return values ..
128      REAL              RESULT(6)
129*
130*  =====================================================================
131*
132*     ..
133*     .. Local allocatable arrays
134      REAL            , ALLOCATABLE ::  A(:,:), AF(:,:), Q(:,:), R(:,:),
135     $                   RWORK(:), WORK( : ), T1(:,:), T2(:,:), DIAG(:),
136     $                   C(:,:), CF(:,:), D(:,:), DF(:,:)
137*
138*     .. Parameters ..
139      REAL               ONE, ZERO
140      PARAMETER          ( ZERO = 0.0E+0, ONE = 1.0E+0 )
141*     ..
142*     .. Local Scalars ..
143      LOGICAL            TESTZEROS
144      INTEGER            INFO, I, J, K, L, LWORK, NB1_UB, NB2_UB, NRB
145      REAL               ANORM, EPS, RESID, CNORM, DNORM
146*     ..
147*     .. Local Arrays ..
148      INTEGER            ISEED( 4 )
149      REAL               WORKQUERY( 1 )
150*     ..
151*     .. External Functions ..
152      REAL               SLAMCH, SLANGE, SLANSY
153      EXTERNAL           SLAMCH, SLANGE, SLANSY
154*     ..
155*     .. External Subroutines ..
156      EXTERNAL           SLACPY, SLARNV, SLASET, SLATSQR, SORHR_COL,
157     $                   SORGTSQR, SSCAL, SGEMM, SGEMQRT, SSYRK
158*     ..
159*     .. Intrinsic Functions ..
160      INTRINSIC          CEILING, REAL, MAX, MIN
161*     ..
162*     .. Scalars in Common ..
163      CHARACTER(LEN=32)  SRNAMT
164*     ..
165*     .. Common blocks ..
166      COMMON             / SRMNAMC / SRNAMT
167*     ..
168*     .. Data statements ..
169      DATA ISEED / 1988, 1989, 1990, 1991 /
170*
171*     TEST MATRICES WITH HALF OF MATRIX BEING ZEROS
172*
173      TESTZEROS = .FALSE.
174*
175      EPS = SLAMCH( 'Epsilon' )
176      K = MIN( M, N )
177      L = MAX( M, N, 1)
178*
179*     Dynamically allocate local arrays
180*
181      ALLOCATE ( A(M,N), AF(M,N), Q(L,L), R(M,L), RWORK(L),
182     $           C(M,N), CF(M,N),
183     $           D(N,M), DF(N,M) )
184*
185*     Put random numbers into A and copy to AF
186*
187      DO J = 1, N
188         CALL SLARNV( 2, ISEED, M, A( 1, J ) )
189      END DO
190      IF( TESTZEROS ) THEN
191         IF( M.GE.4 ) THEN
192            DO J = 1, N
193               CALL SLARNV( 2, ISEED, M/2, A( M/4, J ) )
194            END DO
195         END IF
196      END IF
197      CALL SLACPY( 'Full', M, N, A, M, AF, M )
198*
199*     Number of row blocks in SLATSQR
200*
201      NRB = MAX( 1, CEILING( REAL( M - N ) / REAL( MB1 - N ) ) )
202*
203      ALLOCATE ( T1( NB1, N * NRB ) )
204      ALLOCATE ( T2( NB2, N ) )
205      ALLOCATE ( DIAG( N ) )
206*
207*     Begin determine LWORK for the array WORK and allocate memory.
208*
209*     SLATSQR requires NB1 to be bounded by N.
210*
211      NB1_UB = MIN( NB1, N)
212*
213*     SGEMQRT requires NB2 to be bounded by N.
214*
215      NB2_UB = MIN( NB2, N)
216*
217      CALL SLATSQR( M, N, MB1, NB1_UB, AF, M, T1, NB1,
218     $              WORKQUERY, -1, INFO )
219      LWORK = INT( WORKQUERY( 1 ) )
220      CALL SORGTSQR( M, N, MB1, NB1, AF, M, T1, NB1, WORKQUERY, -1,
221     $               INFO )
222
223      LWORK = MAX( LWORK, INT( WORKQUERY( 1 ) ) )
224*
225*     In SGEMQRT, WORK is N*NB2_UB if SIDE = 'L',
226*                or  M*NB2_UB if SIDE = 'R'.
227*
228      LWORK = MAX( LWORK, NB2_UB * N, NB2_UB * M )
229*
230      ALLOCATE ( WORK( LWORK ) )
231*
232*     End allocate memory for WORK.
233*
234*
235*     Begin Householder reconstruction routines
236*
237*     Factor the matrix A in the array AF.
238*
239      SRNAMT = 'SLATSQR'
240      CALL SLATSQR( M, N, MB1, NB1_UB, AF, M, T1, NB1, WORK, LWORK,
241     $              INFO )
242*
243*     Copy the factor R into the array R.
244*
245      SRNAMT = 'SLACPY'
246      CALL SLACPY( 'U', N, N, AF, M, R, M )
247*
248*     Reconstruct the orthogonal matrix Q.
249*
250      SRNAMT = 'SORGTSQR'
251      CALL SORGTSQR( M, N, MB1, NB1, AF, M, T1, NB1, WORK, LWORK,
252     $               INFO )
253*
254*     Perform the Householder reconstruction, the result is stored
255*     the arrays AF and T2.
256*
257      SRNAMT = 'SORHR_COL'
258      CALL SORHR_COL( M, N, NB2, AF, M, T2, NB2, DIAG, INFO )
259*
260*     Compute the factor R_hr corresponding to the Householder
261*     reconstructed Q_hr and place it in the upper triangle of AF to
262*     match the Q storage format in SGEQRT. R_hr = R_tsqr * S,
263*     this means changing the sign of I-th row of the matrix R_tsqr
264*     according to sign of of I-th diagonal element DIAG(I) of the
265*     matrix S.
266*
267      SRNAMT = 'SLACPY'
268      CALL SLACPY( 'U', N, N, R, M, AF, M )
269*
270      DO I = 1, N
271         IF( DIAG( I ).EQ.-ONE ) THEN
272            CALL SSCAL( N+1-I, -ONE, AF( I, I ), M )
273         END IF
274      END DO
275*
276*     End Householder reconstruction routines.
277*
278*
279*     Generate the m-by-m matrix Q
280*
281      CALL SLASET( 'Full', M, M, ZERO, ONE, Q, M )
282*
283      SRNAMT = 'SGEMQRT'
284      CALL SGEMQRT( 'L', 'N', M, M, K, NB2_UB, AF, M, T2, NB2, Q, M,
285     $              WORK, INFO )
286*
287*     Copy R
288*
289      CALL SLASET( 'Full', M, N, ZERO, ZERO, R, M )
290*
291      CALL SLACPY( 'Upper', M, N, AF, M, R, M )
292*
293*     TEST 1
294*     Compute |R - (Q**T)*A| / ( eps * m * |A| ) and store in RESULT(1)
295*
296      CALL SGEMM( 'T', 'N', M, N, M, -ONE, Q, M, A, M, ONE, R, M )
297*
298      ANORM = SLANGE( '1', M, N, A, M, RWORK )
299      RESID = SLANGE( '1', M, N, R, M, RWORK )
300      IF( ANORM.GT.ZERO ) THEN
301         RESULT( 1 ) = RESID / ( EPS * MAX( 1, M ) * ANORM )
302      ELSE
303         RESULT( 1 ) = ZERO
304      END IF
305*
306*     TEST 2
307*     Compute |I - (Q**T)*Q| / ( eps * m ) and store in RESULT(2)
308*
309      CALL SLASET( 'Full', M, M, ZERO, ONE, R, M )
310      CALL SSYRK( 'U', 'T', M, M, -ONE, Q, M, ONE, R, M )
311      RESID = SLANSY( '1', 'Upper', M, R, M, RWORK )
312      RESULT( 2 ) = RESID / ( EPS * MAX( 1, M ) )
313*
314*     Generate random m-by-n matrix C
315*
316      DO J = 1, N
317         CALL SLARNV( 2, ISEED, M, C( 1, J ) )
318      END DO
319      CNORM = SLANGE( '1', M, N, C, M, RWORK )
320      CALL SLACPY( 'Full', M, N, C, M, CF, M )
321*
322*     Apply Q to C as Q*C = CF
323*
324      SRNAMT = 'SGEMQRT'
325      CALL SGEMQRT( 'L', 'N', M, N, K, NB2_UB, AF, M, T2, NB2, CF, M,
326     $               WORK, INFO )
327*
328*     TEST 3
329*     Compute |CF - Q*C| / ( eps *  m * |C| )
330*
331      CALL SGEMM( 'N', 'N', M, N, M, -ONE, Q, M, C, M, ONE, CF, M )
332      RESID = SLANGE( '1', M, N, CF, M, RWORK )
333      IF( CNORM.GT.ZERO ) THEN
334         RESULT( 3 ) = RESID / ( EPS * MAX( 1, M ) * CNORM )
335      ELSE
336         RESULT( 3 ) = ZERO
337      END IF
338*
339*     Copy C into CF again
340*
341      CALL SLACPY( 'Full', M, N, C, M, CF, M )
342*
343*     Apply Q to C as (Q**T)*C = CF
344*
345      SRNAMT = 'SGEMQRT'
346      CALL SGEMQRT( 'L', 'T', M, N, K, NB2_UB, AF, M, T2, NB2, CF, M,
347     $               WORK, INFO )
348*
349*     TEST 4
350*     Compute |CF - (Q**T)*C| / ( eps * m * |C|)
351*
352      CALL SGEMM( 'T', 'N', M, N, M, -ONE, Q, M, C, M, ONE, CF, M )
353      RESID = SLANGE( '1', M, N, CF, M, RWORK )
354      IF( CNORM.GT.ZERO ) THEN
355         RESULT( 4 ) = RESID / ( EPS * MAX( 1, M ) * CNORM )
356      ELSE
357         RESULT( 4 ) = ZERO
358      END IF
359*
360*     Generate random n-by-m matrix D and a copy DF
361*
362      DO J = 1, M
363         CALL SLARNV( 2, ISEED, N, D( 1, J ) )
364      END DO
365      DNORM = SLANGE( '1', N, M, D, N, RWORK )
366      CALL SLACPY( 'Full', N, M, D, N, DF, N )
367*
368*     Apply Q to D as D*Q = DF
369*
370      SRNAMT = 'SGEMQRT'
371      CALL SGEMQRT( 'R', 'N', N, M, K, NB2_UB, AF, M, T2, NB2, DF, N,
372     $               WORK, INFO )
373*
374*     TEST 5
375*     Compute |DF - D*Q| / ( eps * m * |D| )
376*
377      CALL SGEMM( 'N', 'N', N, M, M, -ONE, D, N, Q, M, ONE, DF, N )
378      RESID = SLANGE( '1', N, M, DF, N, RWORK )
379      IF( DNORM.GT.ZERO ) THEN
380         RESULT( 5 ) = RESID / ( EPS * MAX( 1, M ) * DNORM )
381      ELSE
382         RESULT( 5 ) = ZERO
383      END IF
384*
385*     Copy D into DF again
386*
387      CALL SLACPY( 'Full', N, M, D, N, DF, N )
388*
389*     Apply Q to D as D*QT = DF
390*
391      SRNAMT = 'SGEMQRT'
392      CALL SGEMQRT( 'R', 'T', N, M, K, NB2_UB, AF, M, T2, NB2, DF, N,
393     $               WORK, INFO )
394*
395*     TEST 6
396*     Compute |DF - D*(Q**T)| / ( eps * m * |D| )
397*
398      CALL SGEMM( 'N', 'T', N, M, M, -ONE, D, N, Q, M, ONE, DF, N )
399      RESID = SLANGE( '1', N, M, DF, N, RWORK )
400      IF( DNORM.GT.ZERO ) THEN
401         RESULT( 6 ) = RESID / ( EPS * MAX( 1, M ) * DNORM )
402      ELSE
403         RESULT( 6 ) = ZERO
404      END IF
405*
406*     Deallocate all arrays
407*
408      DEALLOCATE ( A, AF, Q, R, RWORK, WORK, T1, T2, DIAG,
409     $             C, D, CF, DF )
410*
411      RETURN
412*
413*     End of SORHR_COL01
414*
415      END
416