LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine derrqr ( character*3  PATH,
integer  NUNIT 
)

DERRQR

Purpose:
 DERRQR tests the error exits for the DOUBLE PRECISION routines
 that use the QR decomposition of a general matrix.
Parameters
[in]PATH
          PATH is CHARACTER*3
          The LAPACK path name for the routines to be tested.
[in]NUNIT
          NUNIT is INTEGER
          The unit number for output.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 57 of file derrqr.f.

57 *
58 * -- LAPACK test routine (version 3.4.0) --
59 * -- LAPACK is a software package provided by Univ. of Tennessee, --
60 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
61 * November 2011
62 *
63 * .. Scalar Arguments ..
64  CHARACTER*3 path
65  INTEGER nunit
66 * ..
67 *
68 * =====================================================================
69 *
70 * .. Parameters ..
71  INTEGER nmax
72  parameter ( nmax = 2 )
73 * ..
74 * .. Local Scalars ..
75  INTEGER i, info, j
76 * ..
77 * .. Local Arrays ..
78  DOUBLE PRECISION a( nmax, nmax ), af( nmax, nmax ), b( nmax ),
79  $ w( nmax ), x( nmax )
80 * ..
81 * .. External Subroutines ..
82  EXTERNAL alaesm, chkxer, dgeqr2, dgeqr2p, dgeqrf,
84  $ dormqr
85 * ..
86 * .. Scalars in Common ..
87  LOGICAL lerr, ok
88  CHARACTER*32 srnamt
89  INTEGER infot, nout
90 * ..
91 * .. Common blocks ..
92  COMMON / infoc / infot, nout, ok, lerr
93  COMMON / srnamc / srnamt
94 * ..
95 * .. Intrinsic Functions ..
96  INTRINSIC dble
97 * ..
98 * .. Executable Statements ..
99 *
100  nout = nunit
101  WRITE( nout, fmt = * )
102 *
103 * Set the variables to innocuous values.
104 *
105  DO 20 j = 1, nmax
106  DO 10 i = 1, nmax
107  a( i, j ) = 1.d0 / dble( i+j )
108  af( i, j ) = 1.d0 / dble( i+j )
109  10 CONTINUE
110  b( j ) = 0.d0
111  w( j ) = 0.d0
112  x( j ) = 0.d0
113  20 CONTINUE
114  ok = .true.
115 *
116 * Error exits for QR factorization
117 *
118 * DGEQRF
119 *
120  srnamt = 'DGEQRF'
121  infot = 1
122  CALL dgeqrf( -1, 0, a, 1, b, w, 1, info )
123  CALL chkxer( 'DGEQRF', infot, nout, lerr, ok )
124  infot = 2
125  CALL dgeqrf( 0, -1, a, 1, b, w, 1, info )
126  CALL chkxer( 'DGEQRF', infot, nout, lerr, ok )
127  infot = 4
128  CALL dgeqrf( 2, 1, a, 1, b, w, 1, info )
129  CALL chkxer( 'DGEQRF', infot, nout, lerr, ok )
130  infot = 7
131  CALL dgeqrf( 1, 2, a, 1, b, w, 1, info )
132  CALL chkxer( 'DGEQRF', infot, nout, lerr, ok )
133 *
134 * DGEQRFP
135 *
136  srnamt = 'DGEQRFP'
137  infot = 1
138  CALL dgeqrfp( -1, 0, a, 1, b, w, 1, info )
139  CALL chkxer( 'DGEQRFP', infot, nout, lerr, ok )
140  infot = 2
141  CALL dgeqrfp( 0, -1, a, 1, b, w, 1, info )
142  CALL chkxer( 'DGEQRFP', infot, nout, lerr, ok )
143  infot = 4
144  CALL dgeqrfp( 2, 1, a, 1, b, w, 1, info )
145  CALL chkxer( 'DGEQRFP', infot, nout, lerr, ok )
146  infot = 7
147  CALL dgeqrfp( 1, 2, a, 1, b, w, 1, info )
148  CALL chkxer( 'DGEQRFP', infot, nout, lerr, ok )
149 *
150 * DGEQR2
151 *
152  srnamt = 'DGEQR2'
153  infot = 1
154  CALL dgeqr2( -1, 0, a, 1, b, w, info )
155  CALL chkxer( 'DGEQR2', infot, nout, lerr, ok )
156  infot = 2
157  CALL dgeqr2( 0, -1, a, 1, b, w, info )
158  CALL chkxer( 'DGEQR2', infot, nout, lerr, ok )
159  infot = 4
160  CALL dgeqr2( 2, 1, a, 1, b, w, info )
161  CALL chkxer( 'DGEQR2', infot, nout, lerr, ok )
162 *
163 * DGEQR2P
164 *
165  srnamt = 'DGEQR2P'
166  infot = 1
167  CALL dgeqr2p( -1, 0, a, 1, b, w, info )
168  CALL chkxer( 'DGEQR2P', infot, nout, lerr, ok )
169  infot = 2
170  CALL dgeqr2p( 0, -1, a, 1, b, w, info )
171  CALL chkxer( 'DGEQR2P', infot, nout, lerr, ok )
172  infot = 4
173  CALL dgeqr2p( 2, 1, a, 1, b, w, info )
174  CALL chkxer( 'DGEQR2P', infot, nout, lerr, ok )
175 *
176 * DGEQRS
177 *
178  srnamt = 'DGEQRS'
179  infot = 1
180  CALL dgeqrs( -1, 0, 0, a, 1, x, b, 1, w, 1, info )
181  CALL chkxer( 'DGEQRS', infot, nout, lerr, ok )
182  infot = 2
183  CALL dgeqrs( 0, -1, 0, a, 1, x, b, 1, w, 1, info )
184  CALL chkxer( 'DGEQRS', infot, nout, lerr, ok )
185  infot = 2
186  CALL dgeqrs( 1, 2, 0, a, 2, x, b, 2, w, 1, info )
187  CALL chkxer( 'DGEQRS', infot, nout, lerr, ok )
188  infot = 3
189  CALL dgeqrs( 0, 0, -1, a, 1, x, b, 1, w, 1, info )
190  CALL chkxer( 'DGEQRS', infot, nout, lerr, ok )
191  infot = 5
192  CALL dgeqrs( 2, 1, 0, a, 1, x, b, 2, w, 1, info )
193  CALL chkxer( 'DGEQRS', infot, nout, lerr, ok )
194  infot = 8
195  CALL dgeqrs( 2, 1, 0, a, 2, x, b, 1, w, 1, info )
196  CALL chkxer( 'DGEQRS', infot, nout, lerr, ok )
197  infot = 10
198  CALL dgeqrs( 1, 1, 2, a, 1, x, b, 1, w, 1, info )
199  CALL chkxer( 'DGEQRS', infot, nout, lerr, ok )
200 *
201 * DORGQR
202 *
203  srnamt = 'DORGQR'
204  infot = 1
205  CALL dorgqr( -1, 0, 0, a, 1, x, w, 1, info )
206  CALL chkxer( 'DORGQR', infot, nout, lerr, ok )
207  infot = 2
208  CALL dorgqr( 0, -1, 0, a, 1, x, w, 1, info )
209  CALL chkxer( 'DORGQR', infot, nout, lerr, ok )
210  infot = 2
211  CALL dorgqr( 1, 2, 0, a, 1, x, w, 2, info )
212  CALL chkxer( 'DORGQR', infot, nout, lerr, ok )
213  infot = 3
214  CALL dorgqr( 0, 0, -1, a, 1, x, w, 1, info )
215  CALL chkxer( 'DORGQR', infot, nout, lerr, ok )
216  infot = 3
217  CALL dorgqr( 1, 1, 2, a, 1, x, w, 1, info )
218  CALL chkxer( 'DORGQR', infot, nout, lerr, ok )
219  infot = 5
220  CALL dorgqr( 2, 2, 0, a, 1, x, w, 2, info )
221  CALL chkxer( 'DORGQR', infot, nout, lerr, ok )
222  infot = 8
223  CALL dorgqr( 2, 2, 0, a, 2, x, w, 1, info )
224  CALL chkxer( 'DORGQR', infot, nout, lerr, ok )
225 *
226 * DORG2R
227 *
228  srnamt = 'DORG2R'
229  infot = 1
230  CALL dorg2r( -1, 0, 0, a, 1, x, w, info )
231  CALL chkxer( 'DORG2R', infot, nout, lerr, ok )
232  infot = 2
233  CALL dorg2r( 0, -1, 0, a, 1, x, w, info )
234  CALL chkxer( 'DORG2R', infot, nout, lerr, ok )
235  infot = 2
236  CALL dorg2r( 1, 2, 0, a, 1, x, w, info )
237  CALL chkxer( 'DORG2R', infot, nout, lerr, ok )
238  infot = 3
239  CALL dorg2r( 0, 0, -1, a, 1, x, w, info )
240  CALL chkxer( 'DORG2R', infot, nout, lerr, ok )
241  infot = 3
242  CALL dorg2r( 2, 1, 2, a, 2, x, w, info )
243  CALL chkxer( 'DORG2R', infot, nout, lerr, ok )
244  infot = 5
245  CALL dorg2r( 2, 1, 0, a, 1, x, w, info )
246  CALL chkxer( 'DORG2R', infot, nout, lerr, ok )
247 *
248 * DORMQR
249 *
250  srnamt = 'DORMQR'
251  infot = 1
252  CALL dormqr( '/', 'N', 0, 0, 0, a, 1, x, af, 1, w, 1, info )
253  CALL chkxer( 'DORMQR', infot, nout, lerr, ok )
254  infot = 2
255  CALL dormqr( 'L', '/', 0, 0, 0, a, 1, x, af, 1, w, 1, info )
256  CALL chkxer( 'DORMQR', infot, nout, lerr, ok )
257  infot = 3
258  CALL dormqr( 'L', 'N', -1, 0, 0, a, 1, x, af, 1, w, 1, info )
259  CALL chkxer( 'DORMQR', infot, nout, lerr, ok )
260  infot = 4
261  CALL dormqr( 'L', 'N', 0, -1, 0, a, 1, x, af, 1, w, 1, info )
262  CALL chkxer( 'DORMQR', infot, nout, lerr, ok )
263  infot = 5
264  CALL dormqr( 'L', 'N', 0, 0, -1, a, 1, x, af, 1, w, 1, info )
265  CALL chkxer( 'DORMQR', infot, nout, lerr, ok )
266  infot = 5
267  CALL dormqr( 'L', 'N', 0, 1, 1, a, 1, x, af, 1, w, 1, info )
268  CALL chkxer( 'DORMQR', infot, nout, lerr, ok )
269  infot = 5
270  CALL dormqr( 'R', 'N', 1, 0, 1, a, 1, x, af, 1, w, 1, info )
271  CALL chkxer( 'DORMQR', infot, nout, lerr, ok )
272  infot = 7
273  CALL dormqr( 'L', 'N', 2, 1, 0, a, 1, x, af, 2, w, 1, info )
274  CALL chkxer( 'DORMQR', infot, nout, lerr, ok )
275  infot = 7
276  CALL dormqr( 'R', 'N', 1, 2, 0, a, 1, x, af, 1, w, 1, info )
277  CALL chkxer( 'DORMQR', infot, nout, lerr, ok )
278  infot = 10
279  CALL dormqr( 'L', 'N', 2, 1, 0, a, 2, x, af, 1, w, 1, info )
280  CALL chkxer( 'DORMQR', infot, nout, lerr, ok )
281  infot = 12
282  CALL dormqr( 'L', 'N', 1, 2, 0, a, 1, x, af, 1, w, 1, info )
283  CALL chkxer( 'DORMQR', infot, nout, lerr, ok )
284  infot = 12
285  CALL dormqr( 'R', 'N', 2, 1, 0, a, 1, x, af, 2, w, 1, info )
286  CALL chkxer( 'DORMQR', infot, nout, lerr, ok )
287 *
288 * DORM2R
289 *
290  srnamt = 'DORM2R'
291  infot = 1
292  CALL dorm2r( '/', 'N', 0, 0, 0, a, 1, x, af, 1, w, info )
293  CALL chkxer( 'DORM2R', infot, nout, lerr, ok )
294  infot = 2
295  CALL dorm2r( 'L', '/', 0, 0, 0, a, 1, x, af, 1, w, info )
296  CALL chkxer( 'DORM2R', infot, nout, lerr, ok )
297  infot = 3
298  CALL dorm2r( 'L', 'N', -1, 0, 0, a, 1, x, af, 1, w, info )
299  CALL chkxer( 'DORM2R', infot, nout, lerr, ok )
300  infot = 4
301  CALL dorm2r( 'L', 'N', 0, -1, 0, a, 1, x, af, 1, w, info )
302  CALL chkxer( 'DORM2R', infot, nout, lerr, ok )
303  infot = 5
304  CALL dorm2r( 'L', 'N', 0, 0, -1, a, 1, x, af, 1, w, info )
305  CALL chkxer( 'DORM2R', infot, nout, lerr, ok )
306  infot = 5
307  CALL dorm2r( 'L', 'N', 0, 1, 1, a, 1, x, af, 1, w, info )
308  CALL chkxer( 'DORM2R', infot, nout, lerr, ok )
309  infot = 5
310  CALL dorm2r( 'R', 'N', 1, 0, 1, a, 1, x, af, 1, w, info )
311  CALL chkxer( 'DORM2R', infot, nout, lerr, ok )
312  infot = 7
313  CALL dorm2r( 'L', 'N', 2, 1, 0, a, 1, x, af, 2, w, info )
314  CALL chkxer( 'DORM2R', infot, nout, lerr, ok )
315  infot = 7
316  CALL dorm2r( 'R', 'N', 1, 2, 0, a, 1, x, af, 1, w, info )
317  CALL chkxer( 'DORM2R', infot, nout, lerr, ok )
318  infot = 10
319  CALL dorm2r( 'L', 'N', 2, 1, 0, a, 2, x, af, 1, w, info )
320  CALL chkxer( 'DORM2R', infot, nout, lerr, ok )
321 *
322 * Print a summary line.
323 *
324  CALL alaesm( path, ok, nout )
325 *
326  RETURN
327 *
328 * End of DERRQR
329 *
subroutine alaesm(PATH, OK, NOUT)
ALAESM
Definition: alaesm.f:65
subroutine dgeqr2(M, N, A, LDA, TAU, WORK, INFO)
DGEQR2 computes the QR factorization of a general rectangular matrix using an unblocked algorithm...
Definition: dgeqr2.f:123
subroutine dormqr(SIDE, TRANS, M, N, K, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
DORMQR
Definition: dormqr.f:169
subroutine chkxer(SRNAMT, INFOT, NOUT, LERR, OK)
Definition: cblat2.f:3199
subroutine dgeqr2p(M, N, A, LDA, TAU, WORK, INFO)
DGEQR2P computes the QR factorization of a general rectangular matrix with non-negative diagonal elem...
Definition: dgeqr2p.f:126
subroutine dgeqrs(M, N, NRHS, A, LDA, TAU, B, LDB, WORK, LWORK, INFO)
DGEQRS
Definition: dgeqrs.f:123
subroutine dorm2r(SIDE, TRANS, M, N, K, A, LDA, TAU, C, LDC, WORK, INFO)
DORM2R multiplies a general matrix by the orthogonal matrix from a QR factorization determined by sge...
Definition: dorm2r.f:161
subroutine dgeqrf(M, N, A, LDA, TAU, WORK, LWORK, INFO)
DGEQRF
Definition: dgeqrf.f:138
subroutine dorg2r(M, N, K, A, LDA, TAU, WORK, INFO)
DORG2R generates all or part of the orthogonal matrix Q from a QR factorization determined by sgeqrf ...
Definition: dorg2r.f:116
subroutine dgeqrfp(M, N, A, LDA, TAU, WORK, LWORK, INFO)
DGEQRFP
Definition: dgeqrfp.f:141
subroutine dorgqr(M, N, K, A, LDA, TAU, WORK, LWORK, INFO)
DORGQR
Definition: dorgqr.f:130

Here is the call graph for this function:

Here is the caller graph for this function: