LAPACK  3.6.0
LAPACK: Linear Algebra PACKage
Collaboration diagram for complex16:

Functions

subroutine zlagge (M, N, KL, KU, D, A, LDA, ISEED, WORK, INFO)
 ZLAGGE More...
 
subroutine zlaghe (N, K, D, A, LDA, ISEED, WORK, INFO)
 ZLAGHE More...
 
subroutine zlagsy (N, K, D, A, LDA, ISEED, WORK, INFO)
 ZLAGSY More...
 
subroutine zlahilb (N, NRHS, A, LDA, X, LDX, B, LDB, WORK, INFO, PATH)
 ZLAHILB More...
 
subroutine zlakf2 (M, N, A, LDA, B, D, E, Z, LDZ)
 ZLAKF2 More...
 
subroutine zlarge (N, A, LDA, ISEED, WORK, INFO)
 ZLARGE More...
 
complex *16 function zlarnd (IDIST, ISEED)
 ZLARND More...
 
subroutine zlaror (SIDE, INIT, M, N, A, LDA, ISEED, X, INFO)
 ZLAROR More...
 
subroutine zlarot (LROWS, LLEFT, LRIGHT, NL, C, S, A, LDA, XLEFT, XRIGHT)
 ZLAROT More...
 
subroutine zlatm1 (MODE, COND, IRSIGN, IDIST, ISEED, D, N, INFO)
 ZLATM1 More...
 
complex *16 function zlatm2 (M, N, I, J, KL, KU, IDIST, ISEED, D, IGRADE, DL, DR, IPVTNG, IWORK, SPARSE)
 ZLATM2 More...
 
complex *16 function zlatm3 (M, N, I, J, ISUB, JSUB, KL, KU, IDIST, ISEED, D, IGRADE, DL, DR, IPVTNG, IWORK, SPARSE)
 ZLATM3 More...
 
subroutine zlatm5 (PRTYPE, M, N, A, LDA, B, LDB, C, LDC, D, LDD, E, LDE, F, LDF, R, LDR, L, LDL, ALPHA, QBLCKA, QBLCKB)
 ZLATM5 More...
 
subroutine zlatm6 (TYPE, N, A, LDA, B, X, LDX, Y, LDY, ALPHA, BETA, WX, WY, S, DIF)
 ZLATM6 More...
 
subroutine zlatme (N, DIST, ISEED, D, MODE, COND, DMAX, RSIGN, UPPER, SIM, DS, MODES, CONDS, KL, KU, ANORM, A, LDA, WORK, INFO)
 ZLATME More...
 
subroutine zlatmr (M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, RSIGN, GRADE, DL, MODEL, CONDL, DR, MODER, CONDR, PIVTNG, IPIVOT, KL, KU, SPARSE, ANORM, PACK, A, LDA, IWORK, INFO)
 ZLATMR More...
 
subroutine zlatms (M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
 ZLATMS More...
 
subroutine zlatmt (M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, RANK, KL, KU, PACK, A, LDA, WORK, INFO)
 ZLATMT More...
 

Detailed Description

This is the group of complex16 LAPACK TESTING MATGEN routines.

Function Documentation

subroutine zlagge ( integer  M,
integer  N,
integer  KL,
integer  KU,
double precision, dimension( * )  D,
complex*16, dimension( lda, * )  A,
integer  LDA,
integer, dimension( 4 )  ISEED,
complex*16, dimension( * )  WORK,
integer  INFO 
)

ZLAGGE

Purpose:
 ZLAGGE generates a complex general m by n matrix A, by pre- and post-
 multiplying a real diagonal matrix D with random unitary matrices:
 A = U*D*V. The lower and upper bandwidths may then be reduced to
 kl and ku by additional unitary transformations.
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.  N >= 0.
[in]KL
          KL is INTEGER
          The number of nonzero subdiagonals within the band of A.
          0 <= KL <= M-1.
[in]KU
          KU is INTEGER
          The number of nonzero superdiagonals within the band of A.
          0 <= KU <= N-1.
[in]D
          D is DOUBLE PRECISION array, dimension (min(M,N))
          The diagonal elements of the diagonal matrix D.
[out]A
          A is COMPLEX*16 array, dimension (LDA,N)
          The generated m by n matrix A.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= M.
[in,out]ISEED
          ISEED is INTEGER array, dimension (4)
          On entry, the seed of the random number generator; the array
          elements must be between 0 and 4095, and ISEED(4) must be
          odd.
          On exit, the seed is updated.
[out]WORK
          WORK is COMPLEX*16 array, dimension (M+N)
[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.
Date
November 2015

Definition at line 116 of file zlagge.f.

116 *
117 * -- LAPACK auxiliary routine (version 3.6.0) --
118 * -- LAPACK is a software package provided by Univ. of Tennessee, --
119 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
120 * November 2015
121 *
122 * .. Scalar Arguments ..
123  INTEGER info, kl, ku, lda, m, n
124 * ..
125 * .. Array Arguments ..
126  INTEGER iseed( 4 )
127  DOUBLE PRECISION d( * )
128  COMPLEX*16 a( lda, * ), work( * )
129 * ..
130 *
131 * =====================================================================
132 *
133 * .. Parameters ..
134  COMPLEX*16 zero, one
135  parameter( zero = ( 0.0d+0, 0.0d+0 ),
136  $ one = ( 1.0d+0, 0.0d+0 ) )
137 * ..
138 * .. Local Scalars ..
139  INTEGER i, j
140  DOUBLE PRECISION wn
141  COMPLEX*16 tau, wa, wb
142 * ..
143 * .. External Subroutines ..
144  EXTERNAL xerbla, zgemv, zgerc, zlacgv, zlarnv, zscal
145 * ..
146 * .. Intrinsic Functions ..
147  INTRINSIC abs, dble, max, min
148 * ..
149 * .. External Functions ..
150  DOUBLE PRECISION dznrm2
151  EXTERNAL dznrm2
152 * ..
153 * .. Executable Statements ..
154 *
155 * Test the input arguments
156 *
157  info = 0
158  IF( m.LT.0 ) THEN
159  info = -1
160  ELSE IF( n.LT.0 ) THEN
161  info = -2
162  ELSE IF( kl.LT.0 .OR. kl.GT.m-1 ) THEN
163  info = -3
164  ELSE IF( ku.LT.0 .OR. ku.GT.n-1 ) THEN
165  info = -4
166  ELSE IF( lda.LT.max( 1, m ) ) THEN
167  info = -7
168  END IF
169  IF( info.LT.0 ) THEN
170  CALL xerbla( 'ZLAGGE', -info )
171  RETURN
172  END IF
173 *
174 * initialize A to diagonal matrix
175 *
176  DO 20 j = 1, n
177  DO 10 i = 1, m
178  a( i, j ) = zero
179  10 CONTINUE
180  20 CONTINUE
181  DO 30 i = 1, min( m, n )
182  a( i, i ) = d( i )
183  30 CONTINUE
184 *
185 * Quick exit if the user wants a diagonal matrix
186 *
187  IF(( kl .EQ. 0 ).AND.( ku .EQ. 0)) RETURN
188 *
189 * pre- and post-multiply A by random unitary matrices
190 *
191  DO 40 i = min( m, n ), 1, -1
192  IF( i.LT.m ) THEN
193 *
194 * generate random reflection
195 *
196  CALL zlarnv( 3, iseed, m-i+1, work )
197  wn = dznrm2( m-i+1, work, 1 )
198  wa = ( wn / abs( work( 1 ) ) )*work( 1 )
199  IF( wn.EQ.zero ) THEN
200  tau = zero
201  ELSE
202  wb = work( 1 ) + wa
203  CALL zscal( m-i, one / wb, work( 2 ), 1 )
204  work( 1 ) = one
205  tau = dble( wb / wa )
206  END IF
207 *
208 * multiply A(i:m,i:n) by random reflection from the left
209 *
210  CALL zgemv( 'Conjugate transpose', m-i+1, n-i+1, one,
211  $ a( i, i ), lda, work, 1, zero, work( m+1 ), 1 )
212  CALL zgerc( m-i+1, n-i+1, -tau, work, 1, work( m+1 ), 1,
213  $ a( i, i ), lda )
214  END IF
215  IF( i.LT.n ) THEN
216 *
217 * generate random reflection
218 *
219  CALL zlarnv( 3, iseed, n-i+1, work )
220  wn = dznrm2( n-i+1, work, 1 )
221  wa = ( wn / abs( work( 1 ) ) )*work( 1 )
222  IF( wn.EQ.zero ) THEN
223  tau = zero
224  ELSE
225  wb = work( 1 ) + wa
226  CALL zscal( n-i, one / wb, work( 2 ), 1 )
227  work( 1 ) = one
228  tau = dble( wb / wa )
229  END IF
230 *
231 * multiply A(i:m,i:n) by random reflection from the right
232 *
233  CALL zgemv( 'No transpose', m-i+1, n-i+1, one, a( i, i ),
234  $ lda, work, 1, zero, work( n+1 ), 1 )
235  CALL zgerc( m-i+1, n-i+1, -tau, work( n+1 ), 1, work, 1,
236  $ a( i, i ), lda )
237  END IF
238  40 CONTINUE
239 *
240 * Reduce number of subdiagonals to KL and number of superdiagonals
241 * to KU
242 *
243  DO 70 i = 1, max( m-1-kl, n-1-ku )
244  IF( kl.LE.ku ) THEN
245 *
246 * annihilate subdiagonal elements first (necessary if KL = 0)
247 *
248  IF( i.LE.min( m-1-kl, n ) ) THEN
249 *
250 * generate reflection to annihilate A(kl+i+1:m,i)
251 *
252  wn = dznrm2( m-kl-i+1, a( kl+i, i ), 1 )
253  wa = ( wn / abs( a( kl+i, i ) ) )*a( kl+i, i )
254  IF( wn.EQ.zero ) THEN
255  tau = zero
256  ELSE
257  wb = a( kl+i, i ) + wa
258  CALL zscal( m-kl-i, one / wb, a( kl+i+1, i ), 1 )
259  a( kl+i, i ) = one
260  tau = dble( wb / wa )
261  END IF
262 *
263 * apply reflection to A(kl+i:m,i+1:n) from the left
264 *
265  CALL zgemv( 'Conjugate transpose', m-kl-i+1, n-i, one,
266  $ a( kl+i, i+1 ), lda, a( kl+i, i ), 1, zero,
267  $ work, 1 )
268  CALL zgerc( m-kl-i+1, n-i, -tau, a( kl+i, i ), 1, work,
269  $ 1, a( kl+i, i+1 ), lda )
270  a( kl+i, i ) = -wa
271  END IF
272 *
273  IF( i.LE.min( n-1-ku, m ) ) THEN
274 *
275 * generate reflection to annihilate A(i,ku+i+1:n)
276 *
277  wn = dznrm2( n-ku-i+1, a( i, ku+i ), lda )
278  wa = ( wn / abs( a( i, ku+i ) ) )*a( i, ku+i )
279  IF( wn.EQ.zero ) THEN
280  tau = zero
281  ELSE
282  wb = a( i, ku+i ) + wa
283  CALL zscal( n-ku-i, one / wb, a( i, ku+i+1 ), lda )
284  a( i, ku+i ) = one
285  tau = dble( wb / wa )
286  END IF
287 *
288 * apply reflection to A(i+1:m,ku+i:n) from the right
289 *
290  CALL zlacgv( n-ku-i+1, a( i, ku+i ), lda )
291  CALL zgemv( 'No transpose', m-i, n-ku-i+1, one,
292  $ a( i+1, ku+i ), lda, a( i, ku+i ), lda, zero,
293  $ work, 1 )
294  CALL zgerc( m-i, n-ku-i+1, -tau, work, 1, a( i, ku+i ),
295  $ lda, a( i+1, ku+i ), lda )
296  a( i, ku+i ) = -wa
297  END IF
298  ELSE
299 *
300 * annihilate superdiagonal elements first (necessary if
301 * KU = 0)
302 *
303  IF( i.LE.min( n-1-ku, m ) ) THEN
304 *
305 * generate reflection to annihilate A(i,ku+i+1:n)
306 *
307  wn = dznrm2( n-ku-i+1, a( i, ku+i ), lda )
308  wa = ( wn / abs( a( i, ku+i ) ) )*a( i, ku+i )
309  IF( wn.EQ.zero ) THEN
310  tau = zero
311  ELSE
312  wb = a( i, ku+i ) + wa
313  CALL zscal( n-ku-i, one / wb, a( i, ku+i+1 ), lda )
314  a( i, ku+i ) = one
315  tau = dble( wb / wa )
316  END IF
317 *
318 * apply reflection to A(i+1:m,ku+i:n) from the right
319 *
320  CALL zlacgv( n-ku-i+1, a( i, ku+i ), lda )
321  CALL zgemv( 'No transpose', m-i, n-ku-i+1, one,
322  $ a( i+1, ku+i ), lda, a( i, ku+i ), lda, zero,
323  $ work, 1 )
324  CALL zgerc( m-i, n-ku-i+1, -tau, work, 1, a( i, ku+i ),
325  $ lda, a( i+1, ku+i ), lda )
326  a( i, ku+i ) = -wa
327  END IF
328 *
329  IF( i.LE.min( m-1-kl, n ) ) THEN
330 *
331 * generate reflection to annihilate A(kl+i+1:m,i)
332 *
333  wn = dznrm2( m-kl-i+1, a( kl+i, i ), 1 )
334  wa = ( wn / abs( a( kl+i, i ) ) )*a( kl+i, i )
335  IF( wn.EQ.zero ) THEN
336  tau = zero
337  ELSE
338  wb = a( kl+i, i ) + wa
339  CALL zscal( m-kl-i, one / wb, a( kl+i+1, i ), 1 )
340  a( kl+i, i ) = one
341  tau = dble( wb / wa )
342  END IF
343 *
344 * apply reflection to A(kl+i:m,i+1:n) from the left
345 *
346  CALL zgemv( 'Conjugate transpose', m-kl-i+1, n-i, one,
347  $ a( kl+i, i+1 ), lda, a( kl+i, i ), 1, zero,
348  $ work, 1 )
349  CALL zgerc( m-kl-i+1, n-i, -tau, a( kl+i, i ), 1, work,
350  $ 1, a( kl+i, i+1 ), lda )
351  a( kl+i, i ) = -wa
352  END IF
353  END IF
354 *
355  IF (i .LE. n) THEN
356  DO 50 j = kl + i + 1, m
357  a( j, i ) = zero
358  50 CONTINUE
359  END IF
360 *
361  IF (i .LE. m) THEN
362  DO 60 j = ku + i + 1, n
363  a( i, j ) = zero
364  60 CONTINUE
365  END IF
366  70 CONTINUE
367  RETURN
368 *
369 * End of ZLAGGE
370 *
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine zlacgv(N, X, INCX)
ZLACGV conjugates a complex vector.
Definition: zlacgv.f:76
subroutine zscal(N, ZA, ZX, INCX)
ZSCAL
Definition: zscal.f:54
subroutine zgemv(TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
ZGEMV
Definition: zgemv.f:160
double precision function dznrm2(N, X, INCX)
DZNRM2
Definition: dznrm2.f:56
subroutine zlarnv(IDIST, ISEED, N, X)
ZLARNV returns a vector of random numbers from a uniform or normal distribution.
Definition: zlarnv.f:101
subroutine zgerc(M, N, ALPHA, X, INCX, Y, INCY, A, LDA)
ZGERC
Definition: zgerc.f:132

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine zlaghe ( integer  N,
integer  K,
double precision, dimension( * )  D,
complex*16, dimension( lda, * )  A,
integer  LDA,
integer, dimension( 4 )  ISEED,
complex*16, dimension( * )  WORK,
integer  INFO 
)

ZLAGHE

Purpose:
 ZLAGHE generates a complex hermitian matrix A, by pre- and post-
 multiplying a real diagonal matrix D with a random unitary matrix:
 A = U*D*U'. The semi-bandwidth may then be reduced to k by additional
 unitary transformations.
Parameters
[in]N
          N is INTEGER
          The order of the matrix A.  N >= 0.
[in]K
          K is INTEGER
          The number of nonzero subdiagonals within the band of A.
          0 <= K <= N-1.
[in]D
          D is DOUBLE PRECISION array, dimension (N)
          The diagonal elements of the diagonal matrix D.
[out]A
          A is COMPLEX*16 array, dimension (LDA,N)
          The generated n by n hermitian matrix A (the full matrix is
          stored).
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= N.
[in,out]ISEED
          ISEED is INTEGER array, dimension (4)
          On entry, the seed of the random number generator; the array
          elements must be between 0 and 4095, and ISEED(4) must be
          odd.
          On exit, the seed is updated.
[out]WORK
          WORK is COMPLEX*16 array, dimension (2*N)
[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.
Date
November 2011

Definition at line 104 of file zlaghe.f.

104 *
105 * -- LAPACK auxiliary routine (version 3.4.0) --
106 * -- LAPACK is a software package provided by Univ. of Tennessee, --
107 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
108 * November 2011
109 *
110 * .. Scalar Arguments ..
111  INTEGER info, k, lda, n
112 * ..
113 * .. Array Arguments ..
114  INTEGER iseed( 4 )
115  DOUBLE PRECISION d( * )
116  COMPLEX*16 a( lda, * ), work( * )
117 * ..
118 *
119 * =====================================================================
120 *
121 * .. Parameters ..
122  COMPLEX*16 zero, one, half
123  parameter( zero = ( 0.0d+0, 0.0d+0 ),
124  $ one = ( 1.0d+0, 0.0d+0 ),
125  $ half = ( 0.5d+0, 0.0d+0 ) )
126 * ..
127 * .. Local Scalars ..
128  INTEGER i, j
129  DOUBLE PRECISION wn
130  COMPLEX*16 alpha, tau, wa, wb
131 * ..
132 * .. External Subroutines ..
133  EXTERNAL xerbla, zaxpy, zgemv, zgerc, zhemv, zher2,
134  $ zlarnv, zscal
135 * ..
136 * .. External Functions ..
137  DOUBLE PRECISION dznrm2
138  COMPLEX*16 zdotc
139  EXTERNAL dznrm2, zdotc
140 * ..
141 * .. Intrinsic Functions ..
142  INTRINSIC abs, dble, dconjg, max
143 * ..
144 * .. Executable Statements ..
145 *
146 * Test the input arguments
147 *
148  info = 0
149  IF( n.LT.0 ) THEN
150  info = -1
151  ELSE IF( k.LT.0 .OR. k.GT.n-1 ) THEN
152  info = -2
153  ELSE IF( lda.LT.max( 1, n ) ) THEN
154  info = -5
155  END IF
156  IF( info.LT.0 ) THEN
157  CALL xerbla( 'ZLAGHE', -info )
158  RETURN
159  END IF
160 *
161 * initialize lower triangle of A to diagonal matrix
162 *
163  DO 20 j = 1, n
164  DO 10 i = j + 1, n
165  a( i, j ) = zero
166  10 CONTINUE
167  20 CONTINUE
168  DO 30 i = 1, n
169  a( i, i ) = d( i )
170  30 CONTINUE
171 *
172 * Generate lower triangle of hermitian matrix
173 *
174  DO 40 i = n - 1, 1, -1
175 *
176 * generate random reflection
177 *
178  CALL zlarnv( 3, iseed, n-i+1, work )
179  wn = dznrm2( n-i+1, work, 1 )
180  wa = ( wn / abs( work( 1 ) ) )*work( 1 )
181  IF( wn.EQ.zero ) THEN
182  tau = zero
183  ELSE
184  wb = work( 1 ) + wa
185  CALL zscal( n-i, one / wb, work( 2 ), 1 )
186  work( 1 ) = one
187  tau = dble( wb / wa )
188  END IF
189 *
190 * apply random reflection to A(i:n,i:n) from the left
191 * and the right
192 *
193 * compute y := tau * A * u
194 *
195  CALL zhemv( 'Lower', n-i+1, tau, a( i, i ), lda, work, 1, zero,
196  $ work( n+1 ), 1 )
197 *
198 * compute v := y - 1/2 * tau * ( y, u ) * u
199 *
200  alpha = -half*tau*zdotc( n-i+1, work( n+1 ), 1, work, 1 )
201  CALL zaxpy( n-i+1, alpha, work, 1, work( n+1 ), 1 )
202 *
203 * apply the transformation as a rank-2 update to A(i:n,i:n)
204 *
205  CALL zher2( 'Lower', n-i+1, -one, work, 1, work( n+1 ), 1,
206  $ a( i, i ), lda )
207  40 CONTINUE
208 *
209 * Reduce number of subdiagonals to K
210 *
211  DO 60 i = 1, n - 1 - k
212 *
213 * generate reflection to annihilate A(k+i+1:n,i)
214 *
215  wn = dznrm2( n-k-i+1, a( k+i, i ), 1 )
216  wa = ( wn / abs( a( k+i, i ) ) )*a( k+i, i )
217  IF( wn.EQ.zero ) THEN
218  tau = zero
219  ELSE
220  wb = a( k+i, i ) + wa
221  CALL zscal( n-k-i, one / wb, a( k+i+1, i ), 1 )
222  a( k+i, i ) = one
223  tau = dble( wb / wa )
224  END IF
225 *
226 * apply reflection to A(k+i:n,i+1:k+i-1) from the left
227 *
228  CALL zgemv( 'Conjugate transpose', n-k-i+1, k-1, one,
229  $ a( k+i, i+1 ), lda, a( k+i, i ), 1, zero, work, 1 )
230  CALL zgerc( n-k-i+1, k-1, -tau, a( k+i, i ), 1, work, 1,
231  $ a( k+i, i+1 ), lda )
232 *
233 * apply reflection to A(k+i:n,k+i:n) from the left and the right
234 *
235 * compute y := tau * A * u
236 *
237  CALL zhemv( 'Lower', n-k-i+1, tau, a( k+i, k+i ), lda,
238  $ a( k+i, i ), 1, zero, work, 1 )
239 *
240 * compute v := y - 1/2 * tau * ( y, u ) * u
241 *
242  alpha = -half*tau*zdotc( n-k-i+1, work, 1, a( k+i, i ), 1 )
243  CALL zaxpy( n-k-i+1, alpha, a( k+i, i ), 1, work, 1 )
244 *
245 * apply hermitian rank-2 update to A(k+i:n,k+i:n)
246 *
247  CALL zher2( 'Lower', n-k-i+1, -one, a( k+i, i ), 1, work, 1,
248  $ a( k+i, k+i ), lda )
249 *
250  a( k+i, i ) = -wa
251  DO 50 j = k + i + 1, n
252  a( j, i ) = zero
253  50 CONTINUE
254  60 CONTINUE
255 *
256 * Store full hermitian matrix
257 *
258  DO 80 j = 1, n
259  DO 70 i = j + 1, n
260  a( j, i ) = dconjg( a( i, j ) )
261  70 CONTINUE
262  80 CONTINUE
263  RETURN
264 *
265 * End of ZLAGHE
266 *
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine zher2(UPLO, N, ALPHA, X, INCX, Y, INCY, A, LDA)
ZHER2
Definition: zher2.f:152
subroutine zscal(N, ZA, ZX, INCX)
ZSCAL
Definition: zscal.f:54
subroutine zgemv(TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
ZGEMV
Definition: zgemv.f:160
double precision function dznrm2(N, X, INCX)
DZNRM2
Definition: dznrm2.f:56
subroutine zaxpy(N, ZA, ZX, INCX, ZY, INCY)
ZAXPY
Definition: zaxpy.f:53
subroutine zlarnv(IDIST, ISEED, N, X)
ZLARNV returns a vector of random numbers from a uniform or normal distribution.
Definition: zlarnv.f:101
subroutine zgerc(M, N, ALPHA, X, INCX, Y, INCY, A, LDA)
ZGERC
Definition: zgerc.f:132
subroutine zhemv(UPLO, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
ZHEMV
Definition: zhemv.f:156
complex *16 function zdotc(N, ZX, INCX, ZY, INCY)
ZDOTC
Definition: zdotc.f:54

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine zlagsy ( integer  N,
integer  K,
double precision, dimension( * )  D,
complex*16, dimension( lda, * )  A,
integer  LDA,
integer, dimension( 4 )  ISEED,
complex*16, dimension( * )  WORK,
integer  INFO 
)

ZLAGSY

Purpose:
 ZLAGSY generates a complex symmetric matrix A, by pre- and post-
 multiplying a real diagonal matrix D with a random unitary matrix:
 A = U*D*U**T. The semi-bandwidth may then be reduced to k by
 additional unitary transformations.
Parameters
[in]N
          N is INTEGER
          The order of the matrix A.  N >= 0.
[in]K
          K is INTEGER
          The number of nonzero subdiagonals within the band of A.
          0 <= K <= N-1.
[in]D
          D is DOUBLE PRECISION array, dimension (N)
          The diagonal elements of the diagonal matrix D.
[out]A
          A is COMPLEX*16 array, dimension (LDA,N)
          The generated n by n symmetric matrix A (the full matrix is
          stored).
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= N.
[in,out]ISEED
          ISEED is INTEGER array, dimension (4)
          On entry, the seed of the random number generator; the array
          elements must be between 0 and 4095, and ISEED(4) must be
          odd.
          On exit, the seed is updated.
[out]WORK
          WORK is COMPLEX*16 array, dimension (2*N)
[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.
Date
November 2011

Definition at line 104 of file zlagsy.f.

104 *
105 * -- LAPACK auxiliary routine (version 3.4.0) --
106 * -- LAPACK is a software package provided by Univ. of Tennessee, --
107 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
108 * November 2011
109 *
110 * .. Scalar Arguments ..
111  INTEGER info, k, lda, n
112 * ..
113 * .. Array Arguments ..
114  INTEGER iseed( 4 )
115  DOUBLE PRECISION d( * )
116  COMPLEX*16 a( lda, * ), work( * )
117 * ..
118 *
119 * =====================================================================
120 *
121 * .. Parameters ..
122  COMPLEX*16 zero, one, half
123  parameter( zero = ( 0.0d+0, 0.0d+0 ),
124  $ one = ( 1.0d+0, 0.0d+0 ),
125  $ half = ( 0.5d+0, 0.0d+0 ) )
126 * ..
127 * .. Local Scalars ..
128  INTEGER i, ii, j, jj
129  DOUBLE PRECISION wn
130  COMPLEX*16 alpha, tau, wa, wb
131 * ..
132 * .. External Subroutines ..
133  EXTERNAL xerbla, zaxpy, zgemv, zgerc, zlacgv, zlarnv,
134  $ zscal, zsymv
135 * ..
136 * .. External Functions ..
137  DOUBLE PRECISION dznrm2
138  COMPLEX*16 zdotc
139  EXTERNAL dznrm2, zdotc
140 * ..
141 * .. Intrinsic Functions ..
142  INTRINSIC abs, dble, max
143 * ..
144 * .. Executable Statements ..
145 *
146 * Test the input arguments
147 *
148  info = 0
149  IF( n.LT.0 ) THEN
150  info = -1
151  ELSE IF( k.LT.0 .OR. k.GT.n-1 ) THEN
152  info = -2
153  ELSE IF( lda.LT.max( 1, n ) ) THEN
154  info = -5
155  END IF
156  IF( info.LT.0 ) THEN
157  CALL xerbla( 'ZLAGSY', -info )
158  RETURN
159  END IF
160 *
161 * initialize lower triangle of A to diagonal matrix
162 *
163  DO 20 j = 1, n
164  DO 10 i = j + 1, n
165  a( i, j ) = zero
166  10 CONTINUE
167  20 CONTINUE
168  DO 30 i = 1, n
169  a( i, i ) = d( i )
170  30 CONTINUE
171 *
172 * Generate lower triangle of symmetric matrix
173 *
174  DO 60 i = n - 1, 1, -1
175 *
176 * generate random reflection
177 *
178  CALL zlarnv( 3, iseed, n-i+1, work )
179  wn = dznrm2( n-i+1, work, 1 )
180  wa = ( wn / abs( work( 1 ) ) )*work( 1 )
181  IF( wn.EQ.zero ) THEN
182  tau = zero
183  ELSE
184  wb = work( 1 ) + wa
185  CALL zscal( n-i, one / wb, work( 2 ), 1 )
186  work( 1 ) = one
187  tau = dble( wb / wa )
188  END IF
189 *
190 * apply random reflection to A(i:n,i:n) from the left
191 * and the right
192 *
193 * compute y := tau * A * conjg(u)
194 *
195  CALL zlacgv( n-i+1, work, 1 )
196  CALL zsymv( 'Lower', n-i+1, tau, a( i, i ), lda, work, 1, zero,
197  $ work( n+1 ), 1 )
198  CALL zlacgv( n-i+1, work, 1 )
199 *
200 * compute v := y - 1/2 * tau * ( u, y ) * u
201 *
202  alpha = -half*tau*zdotc( n-i+1, work, 1, work( n+1 ), 1 )
203  CALL zaxpy( n-i+1, alpha, work, 1, work( n+1 ), 1 )
204 *
205 * apply the transformation as a rank-2 update to A(i:n,i:n)
206 *
207 * CALL ZSYR2( 'Lower', N-I+1, -ONE, WORK, 1, WORK( N+1 ), 1,
208 * $ A( I, I ), LDA )
209 *
210  DO 50 jj = i, n
211  DO 40 ii = jj, n
212  a( ii, jj ) = a( ii, jj ) -
213  $ work( ii-i+1 )*work( n+jj-i+1 ) -
214  $ work( n+ii-i+1 )*work( jj-i+1 )
215  40 CONTINUE
216  50 CONTINUE
217  60 CONTINUE
218 *
219 * Reduce number of subdiagonals to K
220 *
221  DO 100 i = 1, n - 1 - k
222 *
223 * generate reflection to annihilate A(k+i+1:n,i)
224 *
225  wn = dznrm2( n-k-i+1, a( k+i, i ), 1 )
226  wa = ( wn / abs( a( k+i, i ) ) )*a( k+i, i )
227  IF( wn.EQ.zero ) THEN
228  tau = zero
229  ELSE
230  wb = a( k+i, i ) + wa
231  CALL zscal( n-k-i, one / wb, a( k+i+1, i ), 1 )
232  a( k+i, i ) = one
233  tau = dble( wb / wa )
234  END IF
235 *
236 * apply reflection to A(k+i:n,i+1:k+i-1) from the left
237 *
238  CALL zgemv( 'Conjugate transpose', n-k-i+1, k-1, one,
239  $ a( k+i, i+1 ), lda, a( k+i, i ), 1, zero, work, 1 )
240  CALL zgerc( n-k-i+1, k-1, -tau, a( k+i, i ), 1, work, 1,
241  $ a( k+i, i+1 ), lda )
242 *
243 * apply reflection to A(k+i:n,k+i:n) from the left and the right
244 *
245 * compute y := tau * A * conjg(u)
246 *
247  CALL zlacgv( n-k-i+1, a( k+i, i ), 1 )
248  CALL zsymv( 'Lower', n-k-i+1, tau, a( k+i, k+i ), lda,
249  $ a( k+i, i ), 1, zero, work, 1 )
250  CALL zlacgv( n-k-i+1, a( k+i, i ), 1 )
251 *
252 * compute v := y - 1/2 * tau * ( u, y ) * u
253 *
254  alpha = -half*tau*zdotc( n-k-i+1, a( k+i, i ), 1, work, 1 )
255  CALL zaxpy( n-k-i+1, alpha, a( k+i, i ), 1, work, 1 )
256 *
257 * apply symmetric rank-2 update to A(k+i:n,k+i:n)
258 *
259 * CALL ZSYR2( 'Lower', N-K-I+1, -ONE, A( K+I, I ), 1, WORK, 1,
260 * $ A( K+I, K+I ), LDA )
261 *
262  DO 80 jj = k + i, n
263  DO 70 ii = jj, n
264  a( ii, jj ) = a( ii, jj ) - a( ii, i )*work( jj-k-i+1 ) -
265  $ work( ii-k-i+1 )*a( jj, i )
266  70 CONTINUE
267  80 CONTINUE
268 *
269  a( k+i, i ) = -wa
270  DO 90 j = k + i + 1, n
271  a( j, i ) = zero
272  90 CONTINUE
273  100 CONTINUE
274 *
275 * Store full symmetric matrix
276 *
277  DO 120 j = 1, n
278  DO 110 i = j + 1, n
279  a( j, i ) = a( i, j )
280  110 CONTINUE
281  120 CONTINUE
282  RETURN
283 *
284 * End of ZLAGSY
285 *
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine zsymv(UPLO, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
ZSYMV computes a matrix-vector product for a complex symmetric matrix.
Definition: zsymv.f:159
subroutine zlacgv(N, X, INCX)
ZLACGV conjugates a complex vector.
Definition: zlacgv.f:76
subroutine zscal(N, ZA, ZX, INCX)
ZSCAL
Definition: zscal.f:54
subroutine zgemv(TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
ZGEMV
Definition: zgemv.f:160
double precision function dznrm2(N, X, INCX)
DZNRM2
Definition: dznrm2.f:56
subroutine zaxpy(N, ZA, ZX, INCX, ZY, INCY)
ZAXPY
Definition: zaxpy.f:53
subroutine zlarnv(IDIST, ISEED, N, X)
ZLARNV returns a vector of random numbers from a uniform or normal distribution.
Definition: zlarnv.f:101
subroutine zgerc(M, N, ALPHA, X, INCX, Y, INCY, A, LDA)
ZGERC
Definition: zgerc.f:132
complex *16 function zdotc(N, ZX, INCX, ZY, INCY)
ZDOTC
Definition: zdotc.f:54

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine zlahilb ( integer  N,
integer  NRHS,
complex*16, dimension(lda,n)  A,
integer  LDA,
complex*16, dimension(ldx, nrhs)  X,
integer  LDX,
complex*16, dimension(ldb, nrhs)  B,
integer  LDB,
double precision, dimension(n)  WORK,
integer  INFO,
character*3  PATH 
)

ZLAHILB

Purpose:
 ZLAHILB generates an N by N scaled Hilbert matrix in A along with
 NRHS right-hand sides in B and solutions in X such that A*X=B.

 The Hilbert matrix is scaled by M = LCM(1, 2, ..., 2*N-1) so that all
 entries are integers.  The right-hand sides are the first NRHS
 columns of M * the identity matrix, and the solutions are the
 first NRHS columns of the inverse Hilbert matrix.

 The condition number of the Hilbert matrix grows exponentially with
 its size, roughly as O(e ** (3.5*N)).  Additionally, the inverse
 Hilbert matrices beyond a relatively small dimension cannot be
 generated exactly without extra precision.  Precision is exhausted
 when the largest entry in the inverse Hilbert matrix is greater than
 2 to the power of the number of bits in the fraction of the data type
 used plus one, which is 24 for single precision.

 In single, the generated solution is exact for N <= 6 and has
 small componentwise error for 7 <= N <= 11.
Parameters
[in]N
          N is INTEGER
          The dimension of the matrix A.
[in]NRHS
          NRHS is INTEGER
          The requested number of right-hand sides.
[out]A
          A is COMPLEX array, dimension (LDA, N)
          The generated scaled Hilbert matrix.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= N.
[out]X
          X is COMPLEX array, dimension (LDX, NRHS)
          The generated exact solutions.  Currently, the first NRHS
          columns of the inverse Hilbert matrix.
[in]LDX
          LDX is INTEGER
          The leading dimension of the array X.  LDX >= N.
[out]B
          B is REAL array, dimension (LDB, NRHS)
          The generated right-hand sides.  Currently, the first NRHS
          columns of LCM(1, 2, ..., 2*N-1) * the identity matrix.
[in]LDB
          LDB is INTEGER
          The leading dimension of the array B.  LDB >= N.
[out]WORK
          WORK is REAL array, dimension (N)
[out]INFO
          INFO is INTEGER
          = 0: successful exit
          = 1: N is too large; the data is still generated but may not
               be not exact.
          < 0: if INFO = -i, the i-th argument had an illegal value
[in]PATH
          PATH is CHARACTER*3
          The LAPACK path name.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2015

Definition at line 136 of file zlahilb.f.

136 *
137 * -- LAPACK test routine (version 3.6.0) --
138 * -- LAPACK is a software package provided by Univ. of Tennessee, --
139 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
140 * November 2015
141 *
142 * .. Scalar Arguments ..
143  INTEGER n, nrhs, lda, ldx, ldb, info
144 * .. Array Arguments ..
145  DOUBLE PRECISION work(n)
146  COMPLEX*16 a(lda,n), x(ldx, nrhs), b(ldb, nrhs)
147  CHARACTER*3 path
148 * ..
149 *
150 * =====================================================================
151 * .. Local Scalars ..
152  INTEGER tm, ti, r
153  INTEGER m
154  INTEGER i, j
155  COMPLEX*16 tmp
156  CHARACTER*2 c2
157 
158 * .. Parameters ..
159 * NMAX_EXACT the largest dimension where the generated data is
160 * exact.
161 * NMAX_APPROX the largest dimension where the generated data has
162 * a small componentwise relative error.
163 * ??? complex uses how many bits ???
164  INTEGER nmax_exact, nmax_approx, size_d
165  parameter(nmax_exact = 6, nmax_approx = 11, size_d = 8)
166 
167 * d's are generated from random permuation of those eight elements.
168  COMPLEX*16 d1(8), d2(8), invd1(8), invd2(8)
169  DATA d1 /(-1,0),(0,1),(-1,-1),(0,-1),(1,0),(-1,1),(1,1),(1,-1)/
170  DATA d2 /(-1,0),(0,-1),(-1,1),(0,1),(1,0),(-1,-1),(1,-1),(1,1)/
171 
172  DATA invd1 /(-1,0),(0,-1),(-.5,.5),(0,1),(1,0),
173  $ (-.5,-.5),(.5,-.5),(.5,.5)/
174  DATA invd2 /(-1,0),(0,1),(-.5,-.5),(0,-1),(1,0),
175  $ (-.5,.5),(.5,.5),(.5,-.5)/
176 * ..
177 * .. External Functions
178  EXTERNAL zlaset, lsamen
179  INTRINSIC dble
180  LOGICAL lsamen
181 * ..
182 * .. Executable Statements ..
183  c2 = path( 2: 3 )
184 *
185 * Test the input arguments
186 *
187  info = 0
188  IF (n .LT. 0 .OR. n .GT. nmax_approx) THEN
189  info = -1
190  ELSE IF (nrhs .LT. 0) THEN
191  info = -2
192  ELSE IF (lda .LT. n) THEN
193  info = -4
194  ELSE IF (ldx .LT. n) THEN
195  info = -6
196  ELSE IF (ldb .LT. n) THEN
197  info = -8
198  END IF
199  IF (info .LT. 0) THEN
200  CALL xerbla('ZLAHILB', -info)
201  RETURN
202  END IF
203  IF (n .GT. nmax_exact) THEN
204  info = 1
205  END IF
206 
207 * Compute M = the LCM of the integers [1, 2*N-1]. The largest
208 * reasonable N is small enough that integers suffice (up to N = 11).
209  m = 1
210  DO i = 2, (2*n-1)
211  tm = m
212  ti = i
213  r = mod(tm, ti)
214  DO WHILE (r .NE. 0)
215  tm = ti
216  ti = r
217  r = mod(tm, ti)
218  END DO
219  m = (m / ti) * i
220  END DO
221 
222 * Generate the scaled Hilbert matrix in A
223 * If we are testing SY routines,
224 * take D1_i = D2_i, else, D1_i = D2_i*
225  IF ( lsamen( 2, c2, 'SY' ) ) THEN
226  DO j = 1, n
227  DO i = 1, n
228  a(i, j) = d1(mod(j,size_d)+1) * (dble(m) / (i + j - 1))
229  $ * d1(mod(i,size_d)+1)
230  END DO
231  END DO
232  ELSE
233  DO j = 1, n
234  DO i = 1, n
235  a(i, j) = d1(mod(j,size_d)+1) * (dble(m) / (i + j - 1))
236  $ * d2(mod(i,size_d)+1)
237  END DO
238  END DO
239  END IF
240 
241 * Generate matrix B as simply the first NRHS columns of M * the
242 * identity.
243  tmp = dble(m)
244  CALL zlaset('Full', n, nrhs, (0.0d+0,0.0d+0), tmp, b, ldb)
245 
246 * Generate the true solutions in X. Because B = the first NRHS
247 * columns of M*I, the true solutions are just the first NRHS columns
248 * of the inverse Hilbert matrix.
249  work(1) = n
250  DO j = 2, n
251  work(j) = ( ( (work(j-1)/(j-1)) * (j-1 - n) ) /(j-1) )
252  $ * (n +j -1)
253  END DO
254 
255 * If we are testing SY routines,
256 * take D1_i = D2_i, else, D1_i = D2_i*
257  IF ( lsamen( 2, c2, 'SY' ) ) THEN
258  DO j = 1, nrhs
259  DO i = 1, n
260  x(i, j) = invd1(mod(j,size_d)+1) *
261  $ ((work(i)*work(j)) / (i + j - 1))
262  $ * invd1(mod(i,size_d)+1)
263  END DO
264  END DO
265  ELSE
266  DO j = 1, nrhs
267  DO i = 1, n
268  x(i, j) = invd2(mod(j,size_d)+1) *
269  $ ((work(i)*work(j)) / (i + j - 1))
270  $ * invd1(mod(i,size_d)+1)
271  END DO
272  END DO
273  END IF
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine zlaset(UPLO, M, N, ALPHA, BETA, A, LDA)
ZLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values...
Definition: zlaset.f:108
logical function lsamen(N, CA, CB)
LSAMEN
Definition: lsamen.f:76

Here is the call graph for this function:

subroutine zlakf2 ( integer  M,
integer  N,
complex*16, dimension( lda, * )  A,
integer  LDA,
complex*16, dimension( lda, * )  B,
complex*16, dimension( lda, * )  D,
complex*16, dimension( lda, * )  E,
complex*16, dimension( ldz, * )  Z,
integer  LDZ 
)

ZLAKF2

Purpose:
 Form the 2*M*N by 2*M*N matrix

        Z = [ kron(In, A)  -kron(B', Im) ]
            [ kron(In, D)  -kron(E', Im) ],

 where In is the identity matrix of size n and X' is the transpose
 of X. kron(X, Y) is the Kronecker product between the matrices X
 and Y.
Parameters
[in]M
          M is INTEGER
          Size of matrix, must be >= 1.
[in]N
          N is INTEGER
          Size of matrix, must be >= 1.
[in]A
          A is COMPLEX*16, dimension ( LDA, M )
          The matrix A in the output matrix Z.
[in]LDA
          LDA is INTEGER
          The leading dimension of A, B, D, and E. ( LDA >= M+N )
[in]B
          B is COMPLEX*16, dimension ( LDA, N )
[in]D
          D is COMPLEX*16, dimension ( LDA, M )
[in]E
          E is COMPLEX*16, dimension ( LDA, N )

          The matrices used in forming the output matrix Z.
[out]Z
          Z is COMPLEX*16, dimension ( LDZ, 2*M*N )
          The resultant Kronecker M*N*2 by M*N*2 matrix (see above.)
[in]LDZ
          LDZ is INTEGER
          The leading dimension of Z. ( LDZ >= 2*M*N )
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 107 of file zlakf2.f.

107 *
108 * -- LAPACK computational routine (version 3.4.0) --
109 * -- LAPACK is a software package provided by Univ. of Tennessee, --
110 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
111 * November 2011
112 *
113 * .. Scalar Arguments ..
114  INTEGER lda, ldz, m, n
115 * ..
116 * .. Array Arguments ..
117  COMPLEX*16 a( lda, * ), b( lda, * ), d( lda, * ),
118  $ e( lda, * ), z( ldz, * )
119 * ..
120 *
121 * ====================================================================
122 *
123 * .. Parameters ..
124  COMPLEX*16 zero
125  parameter( zero = ( 0.0d+0, 0.0d+0 ) )
126 * ..
127 * .. Local Scalars ..
128  INTEGER i, ik, j, jk, l, mn, mn2
129 * ..
130 * .. External Subroutines ..
131  EXTERNAL zlaset
132 * ..
133 * .. Executable Statements ..
134 *
135 * Initialize Z
136 *
137  mn = m*n
138  mn2 = 2*mn
139  CALL zlaset( 'Full', mn2, mn2, zero, zero, z, ldz )
140 *
141  ik = 1
142  DO 50 l = 1, n
143 *
144 * form kron(In, A)
145 *
146  DO 20 i = 1, m
147  DO 10 j = 1, m
148  z( ik+i-1, ik+j-1 ) = a( i, j )
149  10 CONTINUE
150  20 CONTINUE
151 *
152 * form kron(In, D)
153 *
154  DO 40 i = 1, m
155  DO 30 j = 1, m
156  z( ik+mn+i-1, ik+j-1 ) = d( i, j )
157  30 CONTINUE
158  40 CONTINUE
159 *
160  ik = ik + m
161  50 CONTINUE
162 *
163  ik = 1
164  DO 90 l = 1, n
165  jk = mn + 1
166 *
167  DO 80 j = 1, n
168 *
169 * form -kron(B', Im)
170 *
171  DO 60 i = 1, m
172  z( ik+i-1, jk+i-1 ) = -b( j, l )
173  60 CONTINUE
174 *
175 * form -kron(E', Im)
176 *
177  DO 70 i = 1, m
178  z( ik+mn+i-1, jk+i-1 ) = -e( j, l )
179  70 CONTINUE
180 *
181  jk = jk + m
182  80 CONTINUE
183 *
184  ik = ik + m
185  90 CONTINUE
186 *
187  RETURN
188 *
189 * End of ZLAKF2
190 *
subroutine zlaset(UPLO, M, N, ALPHA, BETA, A, LDA)
ZLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values...
Definition: zlaset.f:108

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine zlarge ( integer  N,
complex*16, dimension( lda, * )  A,
integer  LDA,
integer, dimension( 4 )  ISEED,
complex*16, dimension( * )  WORK,
integer  INFO 
)

ZLARGE

Purpose:
 ZLARGE pre- and post-multiplies a complex general n by n matrix A
 with a random unitary matrix: A = U*D*U'.
Parameters
[in]N
          N is INTEGER
          The order of the matrix A.  N >= 0.
[in,out]A
          A is COMPLEX*16 array, dimension (LDA,N)
          On entry, the original n by n matrix A.
          On exit, A is overwritten by U*A*U' for some random
          unitary matrix U.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= N.
[in,out]ISEED
          ISEED is INTEGER array, dimension (4)
          On entry, the seed of the random number generator; the array
          elements must be between 0 and 4095, and ISEED(4) must be
          odd.
          On exit, the seed is updated.
[out]WORK
          WORK is COMPLEX*16 array, dimension (2*N)
[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.
Date
November 2011

Definition at line 89 of file zlarge.f.

89 *
90 * -- LAPACK auxiliary routine (version 3.4.0) --
91 * -- LAPACK is a software package provided by Univ. of Tennessee, --
92 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
93 * November 2011
94 *
95 * .. Scalar Arguments ..
96  INTEGER info, lda, n
97 * ..
98 * .. Array Arguments ..
99  INTEGER iseed( 4 )
100  COMPLEX*16 a( lda, * ), work( * )
101 * ..
102 *
103 * =====================================================================
104 *
105 * .. Parameters ..
106  COMPLEX*16 zero, one
107  parameter( zero = ( 0.0d+0, 0.0d+0 ),
108  $ one = ( 1.0d+0, 0.0d+0 ) )
109 * ..
110 * .. Local Scalars ..
111  INTEGER i
112  DOUBLE PRECISION wn
113  COMPLEX*16 tau, wa, wb
114 * ..
115 * .. External Subroutines ..
116  EXTERNAL xerbla, zgemv, zgerc, zlarnv, zscal
117 * ..
118 * .. Intrinsic Functions ..
119  INTRINSIC abs, dble, max
120 * ..
121 * .. External Functions ..
122  DOUBLE PRECISION dznrm2
123  EXTERNAL dznrm2
124 * ..
125 * .. Executable Statements ..
126 *
127 * Test the input arguments
128 *
129  info = 0
130  IF( n.LT.0 ) THEN
131  info = -1
132  ELSE IF( lda.LT.max( 1, n ) ) THEN
133  info = -3
134  END IF
135  IF( info.LT.0 ) THEN
136  CALL xerbla( 'ZLARGE', -info )
137  RETURN
138  END IF
139 *
140 * pre- and post-multiply A by random unitary matrix
141 *
142  DO 10 i = n, 1, -1
143 *
144 * generate random reflection
145 *
146  CALL zlarnv( 3, iseed, n-i+1, work )
147  wn = dznrm2( n-i+1, work, 1 )
148  wa = ( wn / abs( work( 1 ) ) )*work( 1 )
149  IF( wn.EQ.zero ) THEN
150  tau = zero
151  ELSE
152  wb = work( 1 ) + wa
153  CALL zscal( n-i, one / wb, work( 2 ), 1 )
154  work( 1 ) = one
155  tau = dble( wb / wa )
156  END IF
157 *
158 * multiply A(i:n,1:n) by random reflection from the left
159 *
160  CALL zgemv( 'Conjugate transpose', n-i+1, n, one, a( i, 1 ),
161  $ lda, work, 1, zero, work( n+1 ), 1 )
162  CALL zgerc( n-i+1, n, -tau, work, 1, work( n+1 ), 1, a( i, 1 ),
163  $ lda )
164 *
165 * multiply A(1:n,i:n) by random reflection from the right
166 *
167  CALL zgemv( 'No transpose', n, n-i+1, one, a( 1, i ), lda,
168  $ work, 1, zero, work( n+1 ), 1 )
169  CALL zgerc( n, n-i+1, -tau, work( n+1 ), 1, work, 1, a( 1, i ),
170  $ lda )
171  10 CONTINUE
172  RETURN
173 *
174 * End of ZLARGE
175 *
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine zscal(N, ZA, ZX, INCX)
ZSCAL
Definition: zscal.f:54
subroutine zgemv(TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
ZGEMV
Definition: zgemv.f:160
double precision function dznrm2(N, X, INCX)
DZNRM2
Definition: dznrm2.f:56
subroutine zlarnv(IDIST, ISEED, N, X)
ZLARNV returns a vector of random numbers from a uniform or normal distribution.
Definition: zlarnv.f:101
subroutine zgerc(M, N, ALPHA, X, INCX, Y, INCY, A, LDA)
ZGERC
Definition: zgerc.f:132

Here is the call graph for this function:

Here is the caller graph for this function:

complex*16 function zlarnd ( integer  IDIST,
integer, dimension( 4 )  ISEED 
)

ZLARND

Purpose:
 ZLARND returns a random complex number from a uniform or normal
 distribution.
Parameters
[in]IDIST
          IDIST is INTEGER
          Specifies the distribution of the random numbers:
          = 1:  real and imaginary parts each uniform (0,1)
          = 2:  real and imaginary parts each uniform (-1,1)
          = 3:  real and imaginary parts each normal (0,1)
          = 4:  uniformly distributed on the disc abs(z) <= 1
          = 5:  uniformly distributed on the circle abs(z) = 1
[in,out]ISEED
          ISEED is INTEGER array, dimension (4)
          On entry, the seed of the random number generator; the array
          elements must be between 0 and 4095, and ISEED(4) must be
          odd.
          On exit, the seed is updated.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011
Further Details:
  This routine calls the auxiliary routine DLARAN to generate a random
  real number from a uniform (0,1) distribution. The Box-Muller method
  is used to transform numbers from a uniform to a normal distribution.

Definition at line 77 of file zlarnd.f.

77 *
78 * -- LAPACK auxiliary routine (version 3.4.0) --
79 * -- LAPACK is a software package provided by Univ. of Tennessee, --
80 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
81 * November 2011
82 *
83 * .. Scalar Arguments ..
84  INTEGER idist
85 * ..
86 * .. Array Arguments ..
87  INTEGER iseed( 4 )
88 * ..
89 *
90 * =====================================================================
91 *
92 * .. Parameters ..
93  DOUBLE PRECISION zero, one, two
94  parameter( zero = 0.0d+0, one = 1.0d+0, two = 2.0d+0 )
95  DOUBLE PRECISION twopi
96  parameter( twopi = 6.2831853071795864769252867663d+0 )
97 * ..
98 * .. Local Scalars ..
99  DOUBLE PRECISION t1, t2
100 * ..
101 * .. External Functions ..
102  DOUBLE PRECISION dlaran
103  EXTERNAL dlaran
104 * ..
105 * .. Intrinsic Functions ..
106  INTRINSIC dcmplx, exp, log, sqrt
107 * ..
108 * .. Executable Statements ..
109 *
110 * Generate a pair of real random numbers from a uniform (0,1)
111 * distribution
112 *
113  t1 = dlaran( iseed )
114  t2 = dlaran( iseed )
115 *
116  IF( idist.EQ.1 ) THEN
117 *
118 * real and imaginary parts each uniform (0,1)
119 *
120  zlarnd = dcmplx( t1, t2 )
121  ELSE IF( idist.EQ.2 ) THEN
122 *
123 * real and imaginary parts each uniform (-1,1)
124 *
125  zlarnd = dcmplx( two*t1-one, two*t2-one )
126  ELSE IF( idist.EQ.3 ) THEN
127 *
128 * real and imaginary parts each normal (0,1)
129 *
130  zlarnd = sqrt( -two*log( t1 ) )*exp( dcmplx( zero, twopi*t2 ) )
131  ELSE IF( idist.EQ.4 ) THEN
132 *
133 * uniform distribution on the unit disc abs(z) <= 1
134 *
135  zlarnd = sqrt( t1 )*exp( dcmplx( zero, twopi*t2 ) )
136  ELSE IF( idist.EQ.5 ) THEN
137 *
138 * uniform distribution on the unit circle abs(z) = 1
139 *
140  zlarnd = exp( dcmplx( zero, twopi*t2 ) )
141  END IF
142  RETURN
143 *
144 * End of ZLARND
145 *
complex *16 function zlarnd(IDIST, ISEED)
ZLARND
Definition: zlarnd.f:77
double precision function dlaran(ISEED)
DLARAN
Definition: dlaran.f:69
subroutine zlaror ( character  SIDE,
character  INIT,
integer  M,
integer  N,
complex*16, dimension( lda, * )  A,
integer  LDA,
integer, dimension( 4 )  ISEED,
complex*16, dimension( * )  X,
integer  INFO 
)

ZLAROR

Purpose:
    ZLAROR pre- or post-multiplies an M by N matrix A by a random
    unitary matrix U, overwriting A. A may optionally be
    initialized to the identity matrix before multiplying by U.
    U is generated using the method of G.W. Stewart
    ( SIAM J. Numer. Anal. 17, 1980, pp. 403-409 ).
    (BLAS-2 version)
Parameters
[in]SIDE
          SIDE is CHARACTER*1
           SIDE specifies whether A is multiplied on the left or right
           by U.
       SIDE = 'L'   Multiply A on the left (premultiply) by U
       SIDE = 'R'   Multiply A on the right (postmultiply) by UC>       SIDE = 'C'   Multiply A on the left by U and the right by UC>       SIDE = 'T'   Multiply A on the left by U and the right by U'
           Not modified.
[in]INIT
          INIT is CHARACTER*1
           INIT specifies whether or not A should be initialized to
           the identity matrix.
              INIT = 'I'   Initialize A to (a section of) the
                           identity matrix before applying U.
              INIT = 'N'   No initialization.  Apply U to the
                           input matrix A.

           INIT = 'I' may be used to generate square (i.e., unitary)
           or rectangular orthogonal matrices (orthogonality being
           in the sense of ZDOTC):

           For square matrices, M=N, and SIDE many be either 'L' or
           'R'; the rows will be orthogonal to each other, as will the
           columns.
           For rectangular matrices where M < N, SIDE = 'R' will
           produce a dense matrix whose rows will be orthogonal and
           whose columns will not, while SIDE = 'L' will produce a
           matrix whose rows will be orthogonal, and whose first M
           columns will be orthogonal, the remaining columns being
           zero.
           For matrices where M > N, just use the previous
           explanation, interchanging 'L' and 'R' and "rows" and
           "columns".

           Not modified.
[in]M
          M is INTEGER
           Number of rows of A. Not modified.
[in]N
          N is INTEGER
           Number of columns of A. Not modified.
[in,out]A
           A is COMPLEX*16 array, dimension ( LDA, N )
           Input and output array. Overwritten by U A ( if SIDE = 'L' )
           or by A U ( if SIDE = 'R' )
           or by U A U* ( if SIDE = 'C')
           or by U A U' ( if SIDE = 'T') on exit.
[in]LDA
          LDA is INTEGER
           Leading dimension of A. Must be at least MAX ( 1, M ).
           Not modified.
[in,out]ISEED
          ISEED is INTEGER array, dimension ( 4 )
           On entry ISEED specifies the seed of the random number
           generator. The array elements should be between 0 and 4095;
           if not they will be reduced mod 4096.  Also, ISEED(4) must
           be odd.  The random number generator uses a linear
           congruential sequence limited to small integers, and so
           should produce machine independent random numbers. The
           values of ISEED are changed on exit, and can be used in the
           next call to ZLAROR to continue the same random number
           sequence.
           Modified.
[out]X
          X is COMPLEX*16 array, dimension ( 3*MAX( M, N ) )
           Workspace. Of length:
               2*M + N if SIDE = 'L',
               2*N + M if SIDE = 'R',
               3*N     if SIDE = 'C' or 'T'.
           Modified.
[out]INFO
          INFO is INTEGER
           An error flag.  It is set to:
            0  if no error.
            1  if ZLARND returned a bad random number (installation
               problem)
           -1  if SIDE is not L, R, C, or T.
           -3  if M is negative.
           -4  if N is negative or if SIDE is C or T and N is not equal
               to M.
           -6  if LDA is less than M.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 160 of file zlaror.f.

160 *
161 * -- LAPACK auxiliary routine (version 3.4.0) --
162 * -- LAPACK is a software package provided by Univ. of Tennessee, --
163 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
164 * November 2011
165 *
166 * .. Scalar Arguments ..
167  CHARACTER init, side
168  INTEGER info, lda, m, n
169 * ..
170 * .. Array Arguments ..
171  INTEGER iseed( 4 )
172  COMPLEX*16 a( lda, * ), x( * )
173 * ..
174 *
175 * =====================================================================
176 *
177 * .. Parameters ..
178  DOUBLE PRECISION zero, one, toosml
179  parameter( zero = 0.0d+0, one = 1.0d+0,
180  $ toosml = 1.0d-20 )
181  COMPLEX*16 czero, cone
182  parameter( czero = ( 0.0d+0, 0.0d+0 ),
183  $ cone = ( 1.0d+0, 0.0d+0 ) )
184 * ..
185 * .. Local Scalars ..
186  INTEGER irow, itype, ixfrm, j, jcol, kbeg, nxfrm
187  DOUBLE PRECISION factor, xabs, xnorm
188  COMPLEX*16 csign, xnorms
189 * ..
190 * .. External Functions ..
191  LOGICAL lsame
192  DOUBLE PRECISION dznrm2
193  COMPLEX*16 zlarnd
194  EXTERNAL lsame, dznrm2, zlarnd
195 * ..
196 * .. External Subroutines ..
197  EXTERNAL xerbla, zgemv, zgerc, zlacgv, zlaset, zscal
198 * ..
199 * .. Intrinsic Functions ..
200  INTRINSIC abs, dcmplx, dconjg
201 * ..
202 * .. Executable Statements ..
203 *
204  info = 0
205  IF( n.EQ.0 .OR. m.EQ.0 )
206  $ RETURN
207 *
208  itype = 0
209  IF( lsame( side, 'L' ) ) THEN
210  itype = 1
211  ELSE IF( lsame( side, 'R' ) ) THEN
212  itype = 2
213  ELSE IF( lsame( side, 'C' ) ) THEN
214  itype = 3
215  ELSE IF( lsame( side, 'T' ) ) THEN
216  itype = 4
217  END IF
218 *
219 * Check for argument errors.
220 *
221  IF( itype.EQ.0 ) THEN
222  info = -1
223  ELSE IF( m.LT.0 ) THEN
224  info = -3
225  ELSE IF( n.LT.0 .OR. ( itype.EQ.3 .AND. n.NE.m ) ) THEN
226  info = -4
227  ELSE IF( lda.LT.m ) THEN
228  info = -6
229  END IF
230  IF( info.NE.0 ) THEN
231  CALL xerbla( 'ZLAROR', -info )
232  RETURN
233  END IF
234 *
235  IF( itype.EQ.1 ) THEN
236  nxfrm = m
237  ELSE
238  nxfrm = n
239  END IF
240 *
241 * Initialize A to the identity matrix if desired
242 *
243  IF( lsame( init, 'I' ) )
244  $ CALL zlaset( 'Full', m, n, czero, cone, a, lda )
245 *
246 * If no rotation possible, still multiply by
247 * a random complex number from the circle |x| = 1
248 *
249 * 2) Compute Rotation by computing Householder
250 * Transformations H(2), H(3), ..., H(n). Note that the
251 * order in which they are computed is irrelevant.
252 *
253  DO 10 j = 1, nxfrm
254  x( j ) = czero
255  10 CONTINUE
256 *
257  DO 30 ixfrm = 2, nxfrm
258  kbeg = nxfrm - ixfrm + 1
259 *
260 * Generate independent normal( 0, 1 ) random numbers
261 *
262  DO 20 j = kbeg, nxfrm
263  x( j ) = zlarnd( 3, iseed )
264  20 CONTINUE
265 *
266 * Generate a Householder transformation from the random vector X
267 *
268  xnorm = dznrm2( ixfrm, x( kbeg ), 1 )
269  xabs = abs( x( kbeg ) )
270  IF( xabs.NE.czero ) THEN
271  csign = x( kbeg ) / xabs
272  ELSE
273  csign = cone
274  END IF
275  xnorms = csign*xnorm
276  x( nxfrm+kbeg ) = -csign
277  factor = xnorm*( xnorm+xabs )
278  IF( abs( factor ).LT.toosml ) THEN
279  info = 1
280  CALL xerbla( 'ZLAROR', -info )
281  RETURN
282  ELSE
283  factor = one / factor
284  END IF
285  x( kbeg ) = x( kbeg ) + xnorms
286 *
287 * Apply Householder transformation to A
288 *
289  IF( itype.EQ.1 .OR. itype.EQ.3 .OR. itype.EQ.4 ) THEN
290 *
291 * Apply H(k) on the left of A
292 *
293  CALL zgemv( 'C', ixfrm, n, cone, a( kbeg, 1 ), lda,
294  $ x( kbeg ), 1, czero, x( 2*nxfrm+1 ), 1 )
295  CALL zgerc( ixfrm, n, -dcmplx( factor ), x( kbeg ), 1,
296  $ x( 2*nxfrm+1 ), 1, a( kbeg, 1 ), lda )
297 *
298  END IF
299 *
300  IF( itype.GE.2 .AND. itype.LE.4 ) THEN
301 *
302 * Apply H(k)* (or H(k)') on the right of A
303 *
304  IF( itype.EQ.4 ) THEN
305  CALL zlacgv( ixfrm, x( kbeg ), 1 )
306  END IF
307 *
308  CALL zgemv( 'N', m, ixfrm, cone, a( 1, kbeg ), lda,
309  $ x( kbeg ), 1, czero, x( 2*nxfrm+1 ), 1 )
310  CALL zgerc( m, ixfrm, -dcmplx( factor ), x( 2*nxfrm+1 ), 1,
311  $ x( kbeg ), 1, a( 1, kbeg ), lda )
312 *
313  END IF
314  30 CONTINUE
315 *
316  x( 1 ) = zlarnd( 3, iseed )
317  xabs = abs( x( 1 ) )
318  IF( xabs.NE.zero ) THEN
319  csign = x( 1 ) / xabs
320  ELSE
321  csign = cone
322  END IF
323  x( 2*nxfrm ) = csign
324 *
325 * Scale the matrix A by D.
326 *
327  IF( itype.EQ.1 .OR. itype.EQ.3 .OR. itype.EQ.4 ) THEN
328  DO 40 irow = 1, m
329  CALL zscal( n, dconjg( x( nxfrm+irow ) ), a( irow, 1 ),
330  $ lda )
331  40 CONTINUE
332  END IF
333 *
334  IF( itype.EQ.2 .OR. itype.EQ.3 ) THEN
335  DO 50 jcol = 1, n
336  CALL zscal( m, x( nxfrm+jcol ), a( 1, jcol ), 1 )
337  50 CONTINUE
338  END IF
339 *
340  IF( itype.EQ.4 ) THEN
341  DO 60 jcol = 1, n
342  CALL zscal( m, dconjg( x( nxfrm+jcol ) ), a( 1, jcol ), 1 )
343  60 CONTINUE
344  END IF
345  RETURN
346 *
347 * End of ZLAROR
348 *
complex *16 function zlarnd(IDIST, ISEED)
ZLARND
Definition: zlarnd.f:77
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine zlacgv(N, X, INCX)
ZLACGV conjugates a complex vector.
Definition: zlacgv.f:76
subroutine zscal(N, ZA, ZX, INCX)
ZSCAL
Definition: zscal.f:54
subroutine zlaset(UPLO, M, N, ALPHA, BETA, A, LDA)
ZLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values...
Definition: zlaset.f:108
subroutine zgemv(TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
ZGEMV
Definition: zgemv.f:160
double precision function dznrm2(N, X, INCX)
DZNRM2
Definition: dznrm2.f:56
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine zgerc(M, N, ALPHA, X, INCX, Y, INCY, A, LDA)
ZGERC
Definition: zgerc.f:132

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine zlarot ( logical  LROWS,
logical  LLEFT,
logical  LRIGHT,
integer  NL,
complex*16  C,
complex*16  S,
complex*16, dimension( * )  A,
integer  LDA,
complex*16  XLEFT,
complex*16  XRIGHT 
)

ZLAROT

Purpose:
    ZLAROT applies a (Givens) rotation to two adjacent rows or
    columns, where one element of the first and/or last column/row
    for use on matrices stored in some format other than GE, so
    that elements of the matrix may be used or modified for which
    no array element is provided.

    One example is a symmetric matrix in SB format (bandwidth=4), for
    which UPLO='L':  Two adjacent rows will have the format:

    row j:     C> C> C> C> C> .  .  .  .
    row j+1:      C> C> C> C> C> .  .  .  .

    '*' indicates elements for which storage is provided,
    '.' indicates elements for which no storage is provided, but
    are not necessarily zero; their values are determined by
    symmetry.  ' ' indicates elements which are necessarily zero,
     and have no storage provided.

    Those columns which have two '*'s can be handled by DROT.
    Those columns which have no '*'s can be ignored, since as long
    as the Givens rotations are carefully applied to preserve
    symmetry, their values are determined.
    Those columns which have one '*' have to be handled separately,
    by using separate variables "p" and "q":

    row j:     C> C> C> C> C> p  .  .  .
    row j+1:   q  C> C> C> C> C> .  .  .  .

    The element p would have to be set correctly, then that column
    is rotated, setting p to its new value.  The next call to
    ZLAROT would rotate columns j and j+1, using p, and restore
    symmetry.  The element q would start out being zero, and be
    made non-zero by the rotation.  Later, rotations would presumably
    be chosen to zero q out.

    Typical Calling Sequences: rotating the i-th and (i+1)-st rows.
    ------- ------- ---------

      General dense matrix:

              CALL ZLAROT(.TRUE.,.FALSE.,.FALSE., N, C,S,
                      A(i,1),LDA, DUMMY, DUMMY)

      General banded matrix in GB format:

              j = MAX(1, i-KL )
              NL = MIN( N, i+KU+1 ) + 1-j
              CALL ZLAROT( .TRUE., i-KL.GE.1, i+KU.LT.N, NL, C,S,
                      A(KU+i+1-j,j),LDA-1, XLEFT, XRIGHT )

              [ note that i+1-j is just MIN(i,KL+1) ]

      Symmetric banded matrix in SY format, bandwidth K,
      lower triangle only:

              j = MAX(1, i-K )
              NL = MIN( K+1, i ) + 1
              CALL ZLAROT( .TRUE., i-K.GE.1, .TRUE., NL, C,S,
                      A(i,j), LDA, XLEFT, XRIGHT )

      Same, but upper triangle only:

              NL = MIN( K+1, N-i ) + 1
              CALL ZLAROT( .TRUE., .TRUE., i+K.LT.N, NL, C,S,
                      A(i,i), LDA, XLEFT, XRIGHT )

      Symmetric banded matrix in SB format, bandwidth K,
      lower triangle only:

              [ same as for SY, except:]
                  . . . .
                      A(i+1-j,j), LDA-1, XLEFT, XRIGHT )

              [ note that i+1-j is just MIN(i,K+1) ]

      Same, but upper triangle only:
                  . . .
                      A(K+1,i), LDA-1, XLEFT, XRIGHT )

      Rotating columns is just the transpose of rotating rows, except
      for GB and SB: (rotating columns i and i+1)

      GB:
              j = MAX(1, i-KU )
              NL = MIN( N, i+KL+1 ) + 1-j
              CALL ZLAROT( .TRUE., i-KU.GE.1, i+KL.LT.N, NL, C,S,
                      A(KU+j+1-i,i),LDA-1, XTOP, XBOTTM )

              [note that KU+j+1-i is just MAX(1,KU+2-i)]

      SB: (upper triangle)

                   . . . . . .
                      A(K+j+1-i,i),LDA-1, XTOP, XBOTTM )

      SB: (lower triangle)

                   . . . . . .
                      A(1,i),LDA-1, XTOP, XBOTTM )
  LROWS  - LOGICAL
           If .TRUE., then ZLAROT will rotate two rows.  If .FALSE.,
           then it will rotate two columns.
           Not modified.

  LLEFT  - LOGICAL
           If .TRUE., then XLEFT will be used instead of the
           corresponding element of A for the first element in the
           second row (if LROWS=.FALSE.) or column (if LROWS=.TRUE.)
           If .FALSE., then the corresponding element of A will be
           used.
           Not modified.

  LRIGHT - LOGICAL
           If .TRUE., then XRIGHT will be used instead of the
           corresponding element of A for the last element in the
           first row (if LROWS=.FALSE.) or column (if LROWS=.TRUE.) If
           .FALSE., then the corresponding element of A will be used.
           Not modified.

  NL     - INTEGER
           The length of the rows (if LROWS=.TRUE.) or columns (if
           LROWS=.FALSE.) to be rotated.  If XLEFT and/or XRIGHT are
           used, the columns/rows they are in should be included in
           NL, e.g., if LLEFT = LRIGHT = .TRUE., then NL must be at
           least 2.  The number of rows/columns to be rotated
           exclusive of those involving XLEFT and/or XRIGHT may
           not be negative, i.e., NL minus how many of LLEFT and
           LRIGHT are .TRUE. must be at least zero; if not, XERBLA
           will be called.
           Not modified.

  C, S   - COMPLEX*16
           Specify the Givens rotation to be applied.  If LROWS is
           true, then the matrix ( c  s )
                                 ( _  _ )
                                 (-s  c )  is applied from the left;
           if false, then the transpose (not conjugated) thereof is
           applied from the right.  Note that in contrast to the
           output of ZROTG or to most versions of ZROT, both C and S
           are complex.  For a Givens rotation, |C|**2 + |S|**2 should
           be 1, but this is not checked.
           Not modified.

  A      - COMPLEX*16 array.
           The array containing the rows/columns to be rotated.  The
           first element of A should be the upper left element to
           be rotated.
           Read and modified.

  LDA    - INTEGER
           The "effective" leading dimension of A.  If A contains
           a matrix stored in GE, HE, or SY format, then this is just
           the leading dimension of A as dimensioned in the calling
           routine.  If A contains a matrix stored in band (GB, HB, or
           SB) format, then this should be *one less* than the leading
           dimension used in the calling routine.  Thus, if A were
           dimensioned A(LDA,*) in ZLAROT, then A(1,j) would be the
           j-th element in the first of the two rows to be rotated,
           and A(2,j) would be the j-th in the second, regardless of
           how the array may be stored in the calling routine.  [A
           cannot, however, actually be dimensioned thus, since for
           band format, the row number may exceed LDA, which is not
           legal FORTRAN.]
           If LROWS=.TRUE., then LDA must be at least 1, otherwise
           it must be at least NL minus the number of .TRUE. values
           in XLEFT and XRIGHT.
           Not modified.

  XLEFT  - COMPLEX*16
           If LLEFT is .TRUE., then XLEFT will be used and modified
           instead of A(2,1) (if LROWS=.TRUE.) or A(1,2)
           (if LROWS=.FALSE.).
           Read and modified.

  XRIGHT - COMPLEX*16
           If LRIGHT is .TRUE., then XRIGHT will be used and modified
           instead of A(1,NL) (if LROWS=.TRUE.) or A(NL,1)
           (if LROWS=.FALSE.).
           Read and modified.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 231 of file zlarot.f.

231 *
232 * -- LAPACK auxiliary routine (version 3.4.0) --
233 * -- LAPACK is a software package provided by Univ. of Tennessee, --
234 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
235 * November 2011
236 *
237 * .. Scalar Arguments ..
238  LOGICAL lleft, lright, lrows
239  INTEGER lda, nl
240  COMPLEX*16 c, s, xleft, xright
241 * ..
242 * .. Array Arguments ..
243  COMPLEX*16 a( * )
244 * ..
245 *
246 * =====================================================================
247 *
248 * .. Local Scalars ..
249  INTEGER iinc, inext, ix, iy, iyt, j, nt
250  COMPLEX*16 tempx
251 * ..
252 * .. Local Arrays ..
253  COMPLEX*16 xt( 2 ), yt( 2 )
254 * ..
255 * .. External Subroutines ..
256  EXTERNAL xerbla
257 * ..
258 * .. Intrinsic Functions ..
259  INTRINSIC dconjg
260 * ..
261 * .. Executable Statements ..
262 *
263 * Set up indices, arrays for ends
264 *
265  IF( lrows ) THEN
266  iinc = lda
267  inext = 1
268  ELSE
269  iinc = 1
270  inext = lda
271  END IF
272 *
273  IF( lleft ) THEN
274  nt = 1
275  ix = 1 + iinc
276  iy = 2 + lda
277  xt( 1 ) = a( 1 )
278  yt( 1 ) = xleft
279  ELSE
280  nt = 0
281  ix = 1
282  iy = 1 + inext
283  END IF
284 *
285  IF( lright ) THEN
286  iyt = 1 + inext + ( nl-1 )*iinc
287  nt = nt + 1
288  xt( nt ) = xright
289  yt( nt ) = a( iyt )
290  END IF
291 *
292 * Check for errors
293 *
294  IF( nl.LT.nt ) THEN
295  CALL xerbla( 'ZLAROT', 4 )
296  RETURN
297  END IF
298  IF( lda.LE.0 .OR. ( .NOT.lrows .AND. lda.LT.nl-nt ) ) THEN
299  CALL xerbla( 'ZLAROT', 8 )
300  RETURN
301  END IF
302 *
303 * Rotate
304 *
305 * ZROT( NL-NT, A(IX),IINC, A(IY),IINC, C, S ) with complex C, S
306 *
307  DO 10 j = 0, nl - nt - 1
308  tempx = c*a( ix+j*iinc ) + s*a( iy+j*iinc )
309  a( iy+j*iinc ) = -dconjg( s )*a( ix+j*iinc ) +
310  $ dconjg( c )*a( iy+j*iinc )
311  a( ix+j*iinc ) = tempx
312  10 CONTINUE
313 *
314 * ZROT( NT, XT,1, YT,1, C, S ) with complex C, S
315 *
316  DO 20 j = 1, nt
317  tempx = c*xt( j ) + s*yt( j )
318  yt( j ) = -dconjg( s )*xt( j ) + dconjg( c )*yt( j )
319  xt( j ) = tempx
320  20 CONTINUE
321 *
322 * Stuff values back into XLEFT, XRIGHT, etc.
323 *
324  IF( lleft ) THEN
325  a( 1 ) = xt( 1 )
326  xleft = yt( 1 )
327  END IF
328 *
329  IF( lright ) THEN
330  xright = xt( nt )
331  a( iyt ) = yt( nt )
332  END IF
333 *
334  RETURN
335 *
336 * End of ZLAROT
337 *
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine zlatm1 ( integer  MODE,
double precision  COND,
integer  IRSIGN,
integer  IDIST,
integer, dimension( 4 )  ISEED,
complex*16, dimension( * )  D,
integer  N,
integer  INFO 
)

ZLATM1

Purpose:
    ZLATM1 computes the entries of D(1..N) as specified by
    MODE, COND and IRSIGN. IDIST and ISEED determine the generation
    of random numbers. ZLATM1 is called by ZLATMR to generate
    random test matrices for LAPACK programs.
Parameters
[in]MODE
          MODE is INTEGER
           On entry describes how D is to be computed:
           MODE = 0 means do not change D.
           MODE = 1 sets D(1)=1 and D(2:N)=1.0/COND
           MODE = 2 sets D(1:N-1)=1 and D(N)=1.0/COND
           MODE = 3 sets D(I)=COND**(-(I-1)/(N-1))
           MODE = 4 sets D(i)=1 - (i-1)/(N-1)*(1 - 1/COND)
           MODE = 5 sets D to random numbers in the range
                    ( 1/COND , 1 ) such that their logarithms
                    are uniformly distributed.
           MODE = 6 set D to random numbers from same distribution
                    as the rest of the matrix.
           MODE < 0 has the same meaning as ABS(MODE), except that
              the order of the elements of D is reversed.
           Thus if MODE is positive, D has entries ranging from
              1 to 1/COND, if negative, from 1/COND to 1,
           Not modified.
[in]COND
          COND is DOUBLE PRECISION
           On entry, used as described under MODE above.
           If used, it must be >= 1. Not modified.
[in]IRSIGN
          IRSIGN is INTEGER
           On entry, if MODE neither -6, 0 nor 6, determines sign of
           entries of D
           0 => leave entries of D unchanged
           1 => multiply each entry of D by random complex number
                uniformly distributed with absolute value 1
[in]IDIST
          IDIST is INTEGER
           On entry, IDIST specifies the type of distribution to be
           used to generate a random matrix .
           1 => real and imaginary parts each UNIFORM( 0, 1 )
           2 => real and imaginary parts each UNIFORM( -1, 1 )
           3 => real and imaginary parts each NORMAL( 0, 1 )
           4 => complex number uniform in DISK( 0, 1 )
           Not modified.
[in,out]ISEED
          ISEED is INTEGER array, dimension ( 4 )
           On entry ISEED specifies the seed of the random number
           generator. The random number generator uses a
           linear congruential sequence limited to small
           integers, and so should produce machine independent
           random numbers. The values of ISEED are changed on
           exit, and can be used in the next call to ZLATM1
           to continue the same random number sequence.
           Changed on exit.
[in,out]D
          D is COMPLEX*16 array, dimension ( N )
           Array to be computed according to MODE, COND and IRSIGN.
           May be changed on exit if MODE is nonzero.
[in]N
          N is INTEGER
           Number of entries of D. Not modified.
[out]INFO
          INFO is INTEGER
            0  => normal termination
           -1  => if MODE not in range -6 to 6
           -2  => if MODE neither -6, 0 nor 6, and
                  IRSIGN neither 0 nor 1
           -3  => if MODE neither -6, 0 nor 6 and COND less than 1
           -4  => if MODE equals 6 or -6 and IDIST not in range 1 to 4
           -7  => if N negative
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2015

Definition at line 139 of file zlatm1.f.

139 *
140 * -- LAPACK auxiliary routine (version 3.6.0) --
141 * -- LAPACK is a software package provided by Univ. of Tennessee, --
142 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
143 * November 2015
144 *
145 * .. Scalar Arguments ..
146  INTEGER idist, info, irsign, mode, n
147  DOUBLE PRECISION cond
148 * ..
149 * .. Array Arguments ..
150  INTEGER iseed( 4 )
151  COMPLEX*16 d( * )
152 * ..
153 *
154 * =====================================================================
155 *
156 * .. Parameters ..
157  DOUBLE PRECISION one
158  parameter( one = 1.0d0 )
159 * ..
160 * .. Local Scalars ..
161  INTEGER i
162  DOUBLE PRECISION alpha, temp
163  COMPLEX*16 ctemp
164 * ..
165 * .. External Functions ..
166  DOUBLE PRECISION dlaran
167  COMPLEX*16 zlarnd
168  EXTERNAL dlaran, zlarnd
169 * ..
170 * .. External Subroutines ..
171  EXTERNAL xerbla, zlarnv
172 * ..
173 * .. Intrinsic Functions ..
174  INTRINSIC abs, dble, exp, log
175 * ..
176 * .. Executable Statements ..
177 *
178 * Decode and Test the input parameters. Initialize flags & seed.
179 *
180  info = 0
181 *
182 * Quick return if possible
183 *
184  IF( n.EQ.0 )
185  $ RETURN
186 *
187 * Set INFO if an error
188 *
189  IF( mode.LT.-6 .OR. mode.GT.6 ) THEN
190  info = -1
191  ELSE IF( ( mode.NE.-6 .AND. mode.NE.0 .AND. mode.NE.6 ) .AND.
192  $ ( irsign.NE.0 .AND. irsign.NE.1 ) ) THEN
193  info = -2
194  ELSE IF( ( mode.NE.-6 .AND. mode.NE.0 .AND. mode.NE.6 ) .AND.
195  $ cond.LT.one ) THEN
196  info = -3
197  ELSE IF( ( mode.EQ.6 .OR. mode.EQ.-6 ) .AND.
198  $ ( idist.LT.1 .OR. idist.GT.4 ) ) THEN
199  info = -4
200  ELSE IF( n.LT.0 ) THEN
201  info = -7
202  END IF
203 *
204  IF( info.NE.0 ) THEN
205  CALL xerbla( 'ZLATM1', -info )
206  RETURN
207  END IF
208 *
209 * Compute D according to COND and MODE
210 *
211  IF( mode.NE.0 ) THEN
212  GO TO ( 10, 30, 50, 70, 90, 110 )abs( mode )
213 *
214 * One large D value:
215 *
216  10 CONTINUE
217  DO 20 i = 1, n
218  d( i ) = one / cond
219  20 CONTINUE
220  d( 1 ) = one
221  GO TO 120
222 *
223 * One small D value:
224 *
225  30 CONTINUE
226  DO 40 i = 1, n
227  d( i ) = one
228  40 CONTINUE
229  d( n ) = one / cond
230  GO TO 120
231 *
232 * Exponentially distributed D values:
233 *
234  50 CONTINUE
235  d( 1 ) = one
236  IF( n.GT.1 ) THEN
237  alpha = cond**( -one / dble( n-1 ) )
238  DO 60 i = 2, n
239  d( i ) = alpha**( i-1 )
240  60 CONTINUE
241  END IF
242  GO TO 120
243 *
244 * Arithmetically distributed D values:
245 *
246  70 CONTINUE
247  d( 1 ) = one
248  IF( n.GT.1 ) THEN
249  temp = one / cond
250  alpha = ( one-temp ) / dble( n-1 )
251  DO 80 i = 2, n
252  d( i ) = dble( n-i )*alpha + temp
253  80 CONTINUE
254  END IF
255  GO TO 120
256 *
257 * Randomly distributed D values on ( 1/COND , 1):
258 *
259  90 CONTINUE
260  alpha = log( one / cond )
261  DO 100 i = 1, n
262  d( i ) = exp( alpha*dlaran( iseed ) )
263  100 CONTINUE
264  GO TO 120
265 *
266 * Randomly distributed D values from IDIST
267 *
268  110 CONTINUE
269  CALL zlarnv( idist, iseed, n, d )
270 *
271  120 CONTINUE
272 *
273 * If MODE neither -6 nor 0 nor 6, and IRSIGN = 1, assign
274 * random signs to D
275 *
276  IF( ( mode.NE.-6 .AND. mode.NE.0 .AND. mode.NE.6 ) .AND.
277  $ irsign.EQ.1 ) THEN
278  DO 130 i = 1, n
279  ctemp = zlarnd( 3, iseed )
280  d( i ) = d( i )*( ctemp / abs( ctemp ) )
281  130 CONTINUE
282  END IF
283 *
284 * Reverse if MODE < 0
285 *
286  IF( mode.LT.0 ) THEN
287  DO 140 i = 1, n / 2
288  ctemp = d( i )
289  d( i ) = d( n+1-i )
290  d( n+1-i ) = ctemp
291  140 CONTINUE
292  END IF
293 *
294  END IF
295 *
296  RETURN
297 *
298 * End of ZLATM1
299 *
complex *16 function zlarnd(IDIST, ISEED)
ZLARND
Definition: zlarnd.f:77
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine zlarnv(IDIST, ISEED, N, X)
ZLARNV returns a vector of random numbers from a uniform or normal distribution.
Definition: zlarnv.f:101
double precision function dlaran(ISEED)
DLARAN
Definition: dlaran.f:69

Here is the call graph for this function:

Here is the caller graph for this function:

complex*16 function zlatm2 ( integer  M,
integer  N,
integer  I,
integer  J,
integer  KL,
integer  KU,
integer  IDIST,
integer, dimension( 4 )  ISEED,
complex*16, dimension( * )  D,
integer  IGRADE,
complex*16, dimension( * )  DL,
complex*16, dimension( * )  DR,
integer  IPVTNG,
integer, dimension( * )  IWORK,
double precision  SPARSE 
)

ZLATM2

Purpose:
    ZLATM2 returns the (I,J) entry of a random matrix of dimension
    (M, N) described by the other paramters. It is called by the
    ZLATMR routine in order to build random test matrices. No error
    checking on parameters is done, because this routine is called in
    a tight loop by ZLATMR which has already checked the parameters.

    Use of ZLATM2 differs from CLATM3 in the order in which the random
    number generator is called to fill in random matrix entries.
    With ZLATM2, the generator is called to fill in the pivoted matrix
    columnwise. With ZLATM3, the generator is called to fill in the
    matrix columnwise, after which it is pivoted. Thus, ZLATM3 can
    be used to construct random matrices which differ only in their
    order of rows and/or columns. ZLATM2 is used to construct band
    matrices while avoiding calling the random number generator for
    entries outside the band (and therefore generating random numbers

    The matrix whose (I,J) entry is returned is constructed as
    follows (this routine only computes one entry):

      If I is outside (1..M) or J is outside (1..N), return zero
         (this is convenient for generating matrices in band format).

      Generate a matrix A with random entries of distribution IDIST.

      Set the diagonal to D.

      Grade the matrix, if desired, from the left (by DL) and/or
         from the right (by DR or DL) as specified by IGRADE.

      Permute, if desired, the rows and/or columns as specified by
         IPVTNG and IWORK.

      Band the matrix to have lower bandwidth KL and upper
         bandwidth KU.

      Set random entries to zero as specified by SPARSE.
Parameters
[in]M
          M is INTEGER
           Number of rows of matrix. Not modified.
[in]N
          N is INTEGER
           Number of columns of matrix. Not modified.
[in]I
          I is INTEGER
           Row of entry to be returned. Not modified.
[in]J
          J is INTEGER
           Column of entry to be returned. Not modified.
[in]KL
          KL is INTEGER
           Lower bandwidth. Not modified.
[in]KU
          KU is INTEGER
           Upper bandwidth. Not modified.
[in]IDIST
          IDIST is INTEGER
           On entry, IDIST specifies the type of distribution to be
           used to generate a random matrix .
           1 => real and imaginary parts each UNIFORM( 0, 1 )
           2 => real and imaginary parts each UNIFORM( -1, 1 )
           3 => real and imaginary parts each NORMAL( 0, 1 )
           4 => complex number uniform in DISK( 0 , 1 )
           Not modified.
[in,out]ISEED
          ISEED is INTEGER array of dimension ( 4 )
           Seed for random number generator.
           Changed on exit.
[in]D
          D is COMPLEX*16 array of dimension ( MIN( I , J ) )
           Diagonal entries of matrix. Not modified.
[in]IGRADE
          IGRADE is INTEGER
           Specifies grading of matrix as follows:
           0  => no grading
           1  => matrix premultiplied by diag( DL )
           2  => matrix postmultiplied by diag( DR )
           3  => matrix premultiplied by diag( DL ) and
                         postmultiplied by diag( DR )
           4  => matrix premultiplied by diag( DL ) and
                         postmultiplied by inv( diag( DL ) )
           5  => matrix premultiplied by diag( DL ) and
                         postmultiplied by diag( CONJG(DL) )
           6  => matrix premultiplied by diag( DL ) and
                         postmultiplied by diag( DL )
           Not modified.
[in]DL
          DL is COMPLEX*16 array ( I or J, as appropriate )
           Left scale factors for grading matrix.  Not modified.
[in]DR
          DR is COMPLEX*16 array ( I or J, as appropriate )
           Right scale factors for grading matrix.  Not modified.
[in]IPVTNG
          IPVTNG is INTEGER
           On entry specifies pivoting permutations as follows:
           0 => none.
           1 => row pivoting.
           2 => column pivoting.
           3 => full pivoting, i.e., on both sides.
           Not modified.
[out]IWORK
          IWORK is INTEGER array ( I or J, as appropriate )
           This array specifies the permutation used. The
           row (or column) in position K was originally in
           position IWORK( K ).
           This differs from IWORK for ZLATM3. Not modified.
[in]SPARSE
          SPARSE is DOUBLE PRECISION between 0. and 1.
           On entry specifies the sparsity of the matrix
           if sparse matix is to be generated.
           SPARSE should lie between 0 and 1.
           A uniform ( 0, 1 ) random number x is generated and
           compared to SPARSE; if x is larger the matrix entry
           is unchanged and if x is smaller the entry is set
           to zero. Thus on the average a fraction SPARSE of the
           entries will be set to zero.
           Not modified.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 213 of file zlatm2.f.

213 *
214 * -- LAPACK auxiliary routine (version 3.4.0) --
215 * -- LAPACK is a software package provided by Univ. of Tennessee, --
216 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
217 * November 2011
218 *
219 * .. Scalar Arguments ..
220 *
221  INTEGER i, idist, igrade, ipvtng, j, kl, ku, m, n
222  DOUBLE PRECISION sparse
223 * ..
224 *
225 * .. Array Arguments ..
226 *
227  INTEGER iseed( 4 ), iwork( * )
228  COMPLEX*16 d( * ), dl( * ), dr( * )
229 * ..
230 *
231 * =====================================================================
232 *
233 * .. Parameters ..
234 *
235  COMPLEX*16 czero
236  parameter( czero = ( 0.0d0, 0.0d0 ) )
237  DOUBLE PRECISION zero
238  parameter( zero = 0.0d0 )
239 * ..
240 *
241 * .. Local Scalars ..
242 *
243  INTEGER isub, jsub
244  COMPLEX*16 ctemp
245 * ..
246 *
247 * .. External Functions ..
248 *
249  DOUBLE PRECISION dlaran
250  COMPLEX*16 zlarnd
251  EXTERNAL dlaran, zlarnd
252 * ..
253 *
254 * .. Intrinsic Functions ..
255 *
256  INTRINSIC dconjg
257 * ..
258 *
259 *-----------------------------------------------------------------------
260 *
261 * .. Executable Statements ..
262 *
263 *
264 * Check for I and J in range
265 *
266  IF( i.LT.1 .OR. i.GT.m .OR. j.LT.1 .OR. j.GT.n ) THEN
267  zlatm2 = czero
268  RETURN
269  END IF
270 *
271 * Check for banding
272 *
273  IF( j.GT.i+ku .OR. j.LT.i-kl ) THEN
274  zlatm2 = czero
275  RETURN
276  END IF
277 *
278 * Check for sparsity
279 *
280  IF( sparse.GT.zero ) THEN
281  IF( dlaran( iseed ).LT.sparse ) THEN
282  zlatm2 = czero
283  RETURN
284  END IF
285  END IF
286 *
287 * Compute subscripts depending on IPVTNG
288 *
289  IF( ipvtng.EQ.0 ) THEN
290  isub = i
291  jsub = j
292  ELSE IF( ipvtng.EQ.1 ) THEN
293  isub = iwork( i )
294  jsub = j
295  ELSE IF( ipvtng.EQ.2 ) THEN
296  isub = i
297  jsub = iwork( j )
298  ELSE IF( ipvtng.EQ.3 ) THEN
299  isub = iwork( i )
300  jsub = iwork( j )
301  END IF
302 *
303 * Compute entry and grade it according to IGRADE
304 *
305  IF( isub.EQ.jsub ) THEN
306  ctemp = d( isub )
307  ELSE
308  ctemp = zlarnd( idist, iseed )
309  END IF
310  IF( igrade.EQ.1 ) THEN
311  ctemp = ctemp*dl( isub )
312  ELSE IF( igrade.EQ.2 ) THEN
313  ctemp = ctemp*dr( jsub )
314  ELSE IF( igrade.EQ.3 ) THEN
315  ctemp = ctemp*dl( isub )*dr( jsub )
316  ELSE IF( igrade.EQ.4 .AND. isub.NE.jsub ) THEN
317  ctemp = ctemp*dl( isub ) / dl( jsub )
318  ELSE IF( igrade.EQ.5 ) THEN
319  ctemp = ctemp*dl( isub )*dconjg( dl( jsub ) )
320  ELSE IF( igrade.EQ.6 ) THEN
321  ctemp = ctemp*dl( isub )*dl( jsub )
322  END IF
323  zlatm2 = ctemp
324  RETURN
325 *
326 * End of ZLATM2
327 *
complex *16 function zlarnd(IDIST, ISEED)
ZLARND
Definition: zlarnd.f:77
complex *16 function zlatm2(M, N, I, J, KL, KU, IDIST, ISEED, D, IGRADE, DL, DR, IPVTNG, IWORK, SPARSE)
ZLATM2
Definition: zlatm2.f:213
double precision function dlaran(ISEED)
DLARAN
Definition: dlaran.f:69
complex*16 function zlatm3 ( integer  M,
integer  N,
integer  I,
integer  J,
integer  ISUB,
integer  JSUB,
integer  KL,
integer  KU,
integer  IDIST,
integer, dimension( 4 )  ISEED,
complex*16, dimension( * )  D,
integer  IGRADE,
complex*16, dimension( * )  DL,
complex*16, dimension( * )  DR,
integer  IPVTNG,
integer, dimension( * )  IWORK,
double precision  SPARSE 
)

ZLATM3

Purpose:
    ZLATM3 returns the (ISUB,JSUB) entry of a random matrix of
    dimension (M, N) described by the other paramters. (ISUB,JSUB)
    is the final position of the (I,J) entry after pivoting
    according to IPVTNG and IWORK. ZLATM3 is called by the
    ZLATMR routine in order to build random test matrices. No error
    checking on parameters is done, because this routine is called in
    a tight loop by ZLATMR which has already checked the parameters.

    Use of ZLATM3 differs from CLATM2 in the order in which the random
    number generator is called to fill in random matrix entries.
    With ZLATM2, the generator is called to fill in the pivoted matrix
    columnwise. With ZLATM3, the generator is called to fill in the
    matrix columnwise, after which it is pivoted. Thus, ZLATM3 can
    be used to construct random matrices which differ only in their
    order of rows and/or columns. ZLATM2 is used to construct band
    matrices while avoiding calling the random number generator for
    entries outside the band (and therefore generating random numbers
    in different orders for different pivot orders).

    The matrix whose (ISUB,JSUB) entry is returned is constructed as
    follows (this routine only computes one entry):

      If ISUB is outside (1..M) or JSUB is outside (1..N), return zero
         (this is convenient for generating matrices in band format).

      Generate a matrix A with random entries of distribution IDIST.

      Set the diagonal to D.

      Grade the matrix, if desired, from the left (by DL) and/or
         from the right (by DR or DL) as specified by IGRADE.

      Permute, if desired, the rows and/or columns as specified by
         IPVTNG and IWORK.

      Band the matrix to have lower bandwidth KL and upper
         bandwidth KU.

      Set random entries to zero as specified by SPARSE.
Parameters
[in]M
          M is INTEGER
           Number of rows of matrix. Not modified.
[in]N
          N is INTEGER
           Number of columns of matrix. Not modified.
[in]I
          I is INTEGER
           Row of unpivoted entry to be returned. Not modified.
[in]J
          J is INTEGER
           Column of unpivoted entry to be returned. Not modified.
[in,out]ISUB
          ISUB is INTEGER
           Row of pivoted entry to be returned. Changed on exit.
[in,out]JSUB
          JSUB is INTEGER
           Column of pivoted entry to be returned. Changed on exit.
[in]KL
          KL is INTEGER
           Lower bandwidth. Not modified.
[in]KU
          KU is INTEGER
           Upper bandwidth. Not modified.
[in]IDIST
          IDIST is INTEGER
           On entry, IDIST specifies the type of distribution to be
           used to generate a random matrix .
           1 => real and imaginary parts each UNIFORM( 0, 1 )
           2 => real and imaginary parts each UNIFORM( -1, 1 )
           3 => real and imaginary parts each NORMAL( 0, 1 )
           4 => complex number uniform in DISK( 0 , 1 )
           Not modified.
[in,out]ISEED
          ISEED is INTEGER array of dimension ( 4 )
           Seed for random number generator.
           Changed on exit.
[in]D
          D is COMPLEX*16 array of dimension ( MIN( I , J ) )
           Diagonal entries of matrix. Not modified.
[in]IGRADE
          IGRADE is INTEGER
           Specifies grading of matrix as follows:
           0  => no grading
           1  => matrix premultiplied by diag( DL )
           2  => matrix postmultiplied by diag( DR )
           3  => matrix premultiplied by diag( DL ) and
                         postmultiplied by diag( DR )
           4  => matrix premultiplied by diag( DL ) and
                         postmultiplied by inv( diag( DL ) )
           5  => matrix premultiplied by diag( DL ) and
                         postmultiplied by diag( CONJG(DL) )
           6  => matrix premultiplied by diag( DL ) and
                         postmultiplied by diag( DL )
           Not modified.
[in]DL
          DL is COMPLEX*16 array ( I or J, as appropriate )
           Left scale factors for grading matrix.  Not modified.
[in]DR
          DR is COMPLEX*16 array ( I or J, as appropriate )
           Right scale factors for grading matrix.  Not modified.
[in]IPVTNG
          IPVTNG is INTEGER
           On entry specifies pivoting permutations as follows:
           0 => none.
           1 => row pivoting.
           2 => column pivoting.
           3 => full pivoting, i.e., on both sides.
           Not modified.
[in]IWORK
          IWORK is INTEGER array ( I or J, as appropriate )
           This array specifies the permutation used. The
           row (or column) originally in position K is in
           position IWORK( K ) after pivoting.
           This differs from IWORK for ZLATM2. Not modified.
[in]SPARSE
          SPARSE is DOUBLE PRECISION between 0. and 1.
           On entry specifies the sparsity of the matrix
           if sparse matix is to be generated.
           SPARSE should lie between 0 and 1.
           A uniform ( 0, 1 ) random number x is generated and
           compared to SPARSE; if x is larger the matrix entry
           is unchanged and if x is smaller the entry is set
           to zero. Thus on the average a fraction SPARSE of the
           entries will be set to zero.
           Not modified.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 231 of file zlatm3.f.

231 *
232 * -- LAPACK auxiliary routine (version 3.4.0) --
233 * -- LAPACK is a software package provided by Univ. of Tennessee, --
234 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
235 * November 2011
236 *
237 * .. Scalar Arguments ..
238 *
239  INTEGER i, idist, igrade, ipvtng, isub, j, jsub, kl,
240  $ ku, m, n
241  DOUBLE PRECISION sparse
242 * ..
243 *
244 * .. Array Arguments ..
245 *
246  INTEGER iseed( 4 ), iwork( * )
247  COMPLEX*16 d( * ), dl( * ), dr( * )
248 * ..
249 *
250 * =====================================================================
251 *
252 * .. Parameters ..
253 *
254  DOUBLE PRECISION zero
255  parameter( zero = 0.0d0 )
256  COMPLEX*16 czero
257  parameter( czero = ( 0.0d0, 0.0d0 ) )
258 * ..
259 *
260 * .. Local Scalars ..
261 *
262  COMPLEX*16 ctemp
263 * ..
264 *
265 * .. External Functions ..
266 *
267  DOUBLE PRECISION dlaran
268  COMPLEX*16 zlarnd
269  EXTERNAL dlaran, zlarnd
270 * ..
271 *
272 * .. Intrinsic Functions ..
273 *
274  INTRINSIC dconjg
275 * ..
276 *
277 *-----------------------------------------------------------------------
278 *
279 * .. Executable Statements ..
280 *
281 *
282 * Check for I and J in range
283 *
284  IF( i.LT.1 .OR. i.GT.m .OR. j.LT.1 .OR. j.GT.n ) THEN
285  isub = i
286  jsub = j
287  zlatm3 = czero
288  RETURN
289  END IF
290 *
291 * Compute subscripts depending on IPVTNG
292 *
293  IF( ipvtng.EQ.0 ) THEN
294  isub = i
295  jsub = j
296  ELSE IF( ipvtng.EQ.1 ) THEN
297  isub = iwork( i )
298  jsub = j
299  ELSE IF( ipvtng.EQ.2 ) THEN
300  isub = i
301  jsub = iwork( j )
302  ELSE IF( ipvtng.EQ.3 ) THEN
303  isub = iwork( i )
304  jsub = iwork( j )
305  END IF
306 *
307 * Check for banding
308 *
309  IF( jsub.GT.isub+ku .OR. jsub.LT.isub-kl ) THEN
310  zlatm3 = czero
311  RETURN
312  END IF
313 *
314 * Check for sparsity
315 *
316  IF( sparse.GT.zero ) THEN
317  IF( dlaran( iseed ).LT.sparse ) THEN
318  zlatm3 = czero
319  RETURN
320  END IF
321  END IF
322 *
323 * Compute entry and grade it according to IGRADE
324 *
325  IF( i.EQ.j ) THEN
326  ctemp = d( i )
327  ELSE
328  ctemp = zlarnd( idist, iseed )
329  END IF
330  IF( igrade.EQ.1 ) THEN
331  ctemp = ctemp*dl( i )
332  ELSE IF( igrade.EQ.2 ) THEN
333  ctemp = ctemp*dr( j )
334  ELSE IF( igrade.EQ.3 ) THEN
335  ctemp = ctemp*dl( i )*dr( j )
336  ELSE IF( igrade.EQ.4 .AND. i.NE.j ) THEN
337  ctemp = ctemp*dl( i ) / dl( j )
338  ELSE IF( igrade.EQ.5 ) THEN
339  ctemp = ctemp*dl( i )*dconjg( dl( j ) )
340  ELSE IF( igrade.EQ.6 ) THEN
341  ctemp = ctemp*dl( i )*dl( j )
342  END IF
343  zlatm3 = ctemp
344  RETURN
345 *
346 * End of ZLATM3
347 *
complex *16 function zlarnd(IDIST, ISEED)
ZLARND
Definition: zlarnd.f:77
complex *16 function zlatm3(M, N, I, J, ISUB, JSUB, KL, KU, IDIST, ISEED, D, IGRADE, DL, DR, IPVTNG, IWORK, SPARSE)
ZLATM3
Definition: zlatm3.f:231
double precision function dlaran(ISEED)
DLARAN
Definition: dlaran.f:69
subroutine zlatm5 ( integer  PRTYPE,
integer  M,
integer  N,
complex*16, dimension( lda, * )  A,
integer  LDA,
complex*16, dimension( ldb, * )  B,
integer  LDB,
complex*16, dimension( ldc, * )  C,
integer  LDC,
complex*16, dimension( ldd, * )  D,
integer  LDD,
complex*16, dimension( lde, * )  E,
integer  LDE,
complex*16, dimension( ldf, * )  F,
integer  LDF,
complex*16, dimension( ldr, * )  R,
integer  LDR,
complex*16, dimension( ldl, * )  L,
integer  LDL,
double precision  ALPHA,
integer  QBLCKA,
integer  QBLCKB 
)

ZLATM5

Purpose:
 ZLATM5 generates matrices involved in the Generalized Sylvester
 equation:

     A * R - L * B = C
     D * R - L * E = F

 They also satisfy (the diagonalization condition)

  [ I -L ] ( [ A  -C ], [ D -F ] ) [ I  R ] = ( [ A    ], [ D    ] )
  [    I ] ( [     B ]  [    E ] ) [    I ]   ( [    B ]  [    E ] )
Parameters
[in]PRTYPE
          PRTYPE is INTEGER
          "Points" to a certian type of the matrices to generate
          (see futher details).
[in]M
          M is INTEGER
          Specifies the order of A and D and the number of rows in
          C, F,  R and L.
[in]N
          N is INTEGER
          Specifies the order of B and E and the number of columns in
          C, F, R and L.
[out]A
          A is COMPLEX*16 array, dimension (LDA, M).
          On exit A M-by-M is initialized according to PRTYPE.
[in]LDA
          LDA is INTEGER
          The leading dimension of A.
[out]B
          B is COMPLEX*16 array, dimension (LDB, N).
          On exit B N-by-N is initialized according to PRTYPE.
[in]LDB
          LDB is INTEGER
          The leading dimension of B.
[out]C
          C is COMPLEX*16 array, dimension (LDC, N).
          On exit C M-by-N is initialized according to PRTYPE.
[in]LDC
          LDC is INTEGER
          The leading dimension of C.
[out]D
          D is COMPLEX*16 array, dimension (LDD, M).
          On exit D M-by-M is initialized according to PRTYPE.
[in]LDD
          LDD is INTEGER
          The leading dimension of D.
[out]E
          E is COMPLEX*16 array, dimension (LDE, N).
          On exit E N-by-N is initialized according to PRTYPE.
[in]LDE
          LDE is INTEGER
          The leading dimension of E.
[out]F
          F is COMPLEX*16 array, dimension (LDF, N).
          On exit F M-by-N is initialized according to PRTYPE.
[in]LDF
          LDF is INTEGER
          The leading dimension of F.
[out]R
          R is COMPLEX*16 array, dimension (LDR, N).
          On exit R M-by-N is initialized according to PRTYPE.
[in]LDR
          LDR is INTEGER
          The leading dimension of R.
[out]L
          L is COMPLEX*16 array, dimension (LDL, N).
          On exit L M-by-N is initialized according to PRTYPE.
[in]LDL
          LDL is INTEGER
          The leading dimension of L.
[in]ALPHA
          ALPHA is DOUBLE PRECISION
          Parameter used in generating PRTYPE = 1 and 5 matrices.
[in]QBLCKA
          QBLCKA is INTEGER
          When PRTYPE = 3, specifies the distance between 2-by-2
          blocks on the diagonal in A. Otherwise, QBLCKA is not
          referenced. QBLCKA > 1.
[in]QBLCKB
          QBLCKB is INTEGER
          When PRTYPE = 3, specifies the distance between 2-by-2
          blocks on the diagonal in B. Otherwise, QBLCKB is not
          referenced. QBLCKB > 1.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011
Further Details:
  PRTYPE = 1: A and B are Jordan blocks, D and E are identity matrices

             A : if (i == j) then A(i, j) = 1.0
                 if (j == i + 1) then A(i, j) = -1.0
                 else A(i, j) = 0.0,            i, j = 1...M

             B : if (i == j) then B(i, j) = 1.0 - ALPHA
                 if (j == i + 1) then B(i, j) = 1.0
                 else B(i, j) = 0.0,            i, j = 1...N

             D : if (i == j) then D(i, j) = 1.0
                 else D(i, j) = 0.0,            i, j = 1...M

             E : if (i == j) then E(i, j) = 1.0
                 else E(i, j) = 0.0,            i, j = 1...N

             L =  R are chosen from [-10...10],
                  which specifies the right hand sides (C, F).

  PRTYPE = 2 or 3: Triangular and/or quasi- triangular.

             A : if (i <= j) then A(i, j) = [-1...1]
                 else A(i, j) = 0.0,             i, j = 1...M

                 if (PRTYPE = 3) then
                    A(k + 1, k + 1) = A(k, k)
                    A(k + 1, k) = [-1...1]
                    sign(A(k, k + 1) = -(sin(A(k + 1, k))
                        k = 1, M - 1, QBLCKA

             B : if (i <= j) then B(i, j) = [-1...1]
                 else B(i, j) = 0.0,            i, j = 1...N

                 if (PRTYPE = 3) then
                    B(k + 1, k + 1) = B(k, k)
                    B(k + 1, k) = [-1...1]
                    sign(B(k, k + 1) = -(sign(B(k + 1, k))
                        k = 1, N - 1, QBLCKB

             D : if (i <= j) then D(i, j) = [-1...1].
                 else D(i, j) = 0.0,            i, j = 1...M


             E : if (i <= j) then D(i, j) = [-1...1]
                 else E(i, j) = 0.0,            i, j = 1...N

                 L, R are chosen from [-10...10],
                 which specifies the right hand sides (C, F).

  PRTYPE = 4 Full
             A(i, j) = [-10...10]
             D(i, j) = [-1...1]    i,j = 1...M
             B(i, j) = [-10...10]
             E(i, j) = [-1...1]    i,j = 1...N
             R(i, j) = [-10...10]
             L(i, j) = [-1...1]    i = 1..M ,j = 1...N

             L, R specifies the right hand sides (C, F).

  PRTYPE = 5 special case common and/or close eigs.

Definition at line 270 of file zlatm5.f.

270 *
271 * -- LAPACK computational routine (version 3.4.0) --
272 * -- LAPACK is a software package provided by Univ. of Tennessee, --
273 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
274 * November 2011
275 *
276 * .. Scalar Arguments ..
277  INTEGER lda, ldb, ldc, ldd, lde, ldf, ldl, ldr, m, n,
278  $ prtype, qblcka, qblckb
279  DOUBLE PRECISION alpha
280 * ..
281 * .. Array Arguments ..
282  COMPLEX*16 a( lda, * ), b( ldb, * ), c( ldc, * ),
283  $ d( ldd, * ), e( lde, * ), f( ldf, * ),
284  $ l( ldl, * ), r( ldr, * )
285 * ..
286 *
287 * =====================================================================
288 *
289 * .. Parameters ..
290  COMPLEX*16 one, two, zero, half, twenty
291  parameter( one = ( 1.0d+0, 0.0d+0 ),
292  $ two = ( 2.0d+0, 0.0d+0 ),
293  $ zero = ( 0.0d+0, 0.0d+0 ),
294  $ half = ( 0.5d+0, 0.0d+0 ),
295  $ twenty = ( 2.0d+1, 0.0d+0 ) )
296 * ..
297 * .. Local Scalars ..
298  INTEGER i, j, k
299  COMPLEX*16 imeps, reeps
300 * ..
301 * .. Intrinsic Functions ..
302  INTRINSIC dcmplx, mod, sin
303 * ..
304 * .. External Subroutines ..
305  EXTERNAL zgemm
306 * ..
307 * .. Executable Statements ..
308 *
309  IF( prtype.EQ.1 ) THEN
310  DO 20 i = 1, m
311  DO 10 j = 1, m
312  IF( i.EQ.j ) THEN
313  a( i, j ) = one
314  d( i, j ) = one
315  ELSE IF( i.EQ.j-1 ) THEN
316  a( i, j ) = -one
317  d( i, j ) = zero
318  ELSE
319  a( i, j ) = zero
320  d( i, j ) = zero
321  END IF
322  10 CONTINUE
323  20 CONTINUE
324 *
325  DO 40 i = 1, n
326  DO 30 j = 1, n
327  IF( i.EQ.j ) THEN
328  b( i, j ) = one - alpha
329  e( i, j ) = one
330  ELSE IF( i.EQ.j-1 ) THEN
331  b( i, j ) = one
332  e( i, j ) = zero
333  ELSE
334  b( i, j ) = zero
335  e( i, j ) = zero
336  END IF
337  30 CONTINUE
338  40 CONTINUE
339 *
340  DO 60 i = 1, m
341  DO 50 j = 1, n
342  r( i, j ) = ( half-sin( dcmplx( i / j ) ) )*twenty
343  l( i, j ) = r( i, j )
344  50 CONTINUE
345  60 CONTINUE
346 *
347  ELSE IF( prtype.EQ.2 .OR. prtype.EQ.3 ) THEN
348  DO 80 i = 1, m
349  DO 70 j = 1, m
350  IF( i.LE.j ) THEN
351  a( i, j ) = ( half-sin( dcmplx( i ) ) )*two
352  d( i, j ) = ( half-sin( dcmplx( i*j ) ) )*two
353  ELSE
354  a( i, j ) = zero
355  d( i, j ) = zero
356  END IF
357  70 CONTINUE
358  80 CONTINUE
359 *
360  DO 100 i = 1, n
361  DO 90 j = 1, n
362  IF( i.LE.j ) THEN
363  b( i, j ) = ( half-sin( dcmplx( i+j ) ) )*two
364  e( i, j ) = ( half-sin( dcmplx( j ) ) )*two
365  ELSE
366  b( i, j ) = zero
367  e( i, j ) = zero
368  END IF
369  90 CONTINUE
370  100 CONTINUE
371 *
372  DO 120 i = 1, m
373  DO 110 j = 1, n
374  r( i, j ) = ( half-sin( dcmplx( i*j ) ) )*twenty
375  l( i, j ) = ( half-sin( dcmplx( i+j ) ) )*twenty
376  110 CONTINUE
377  120 CONTINUE
378 *
379  IF( prtype.EQ.3 ) THEN
380  IF( qblcka.LE.1 )
381  $ qblcka = 2
382  DO 130 k = 1, m - 1, qblcka
383  a( k+1, k+1 ) = a( k, k )
384  a( k+1, k ) = -sin( a( k, k+1 ) )
385  130 CONTINUE
386 *
387  IF( qblckb.LE.1 )
388  $ qblckb = 2
389  DO 140 k = 1, n - 1, qblckb
390  b( k+1, k+1 ) = b( k, k )
391  b( k+1, k ) = -sin( b( k, k+1 ) )
392  140 CONTINUE
393  END IF
394 *
395  ELSE IF( prtype.EQ.4 ) THEN
396  DO 160 i = 1, m
397  DO 150 j = 1, m
398  a( i, j ) = ( half-sin( dcmplx( i*j ) ) )*twenty
399  d( i, j ) = ( half-sin( dcmplx( i+j ) ) )*two
400  150 CONTINUE
401  160 CONTINUE
402 *
403  DO 180 i = 1, n
404  DO 170 j = 1, n
405  b( i, j ) = ( half-sin( dcmplx( i+j ) ) )*twenty
406  e( i, j ) = ( half-sin( dcmplx( i*j ) ) )*two
407  170 CONTINUE
408  180 CONTINUE
409 *
410  DO 200 i = 1, m
411  DO 190 j = 1, n
412  r( i, j ) = ( half-sin( dcmplx( j / i ) ) )*twenty
413  l( i, j ) = ( half-sin( dcmplx( i*j ) ) )*two
414  190 CONTINUE
415  200 CONTINUE
416 *
417  ELSE IF( prtype.GE.5 ) THEN
418  reeps = half*two*twenty / alpha
419  imeps = ( half-two ) / alpha
420  DO 220 i = 1, m
421  DO 210 j = 1, n
422  r( i, j ) = ( half-sin( dcmplx( i*j ) ) )*alpha / twenty
423  l( i, j ) = ( half-sin( dcmplx( i+j ) ) )*alpha / twenty
424  210 CONTINUE
425  220 CONTINUE
426 *
427  DO 230 i = 1, m
428  d( i, i ) = one
429  230 CONTINUE
430 *
431  DO 240 i = 1, m
432  IF( i.LE.4 ) THEN
433  a( i, i ) = one
434  IF( i.GT.2 )
435  $ a( i, i ) = one + reeps
436  IF( mod( i, 2 ).NE.0 .AND. i.LT.m ) THEN
437  a( i, i+1 ) = imeps
438  ELSE IF( i.GT.1 ) THEN
439  a( i, i-1 ) = -imeps
440  END IF
441  ELSE IF( i.LE.8 ) THEN
442  IF( i.LE.6 ) THEN
443  a( i, i ) = reeps
444  ELSE
445  a( i, i ) = -reeps
446  END IF
447  IF( mod( i, 2 ).NE.0 .AND. i.LT.m ) THEN
448  a( i, i+1 ) = one
449  ELSE IF( i.GT.1 ) THEN
450  a( i, i-1 ) = -one
451  END IF
452  ELSE
453  a( i, i ) = one
454  IF( mod( i, 2 ).NE.0 .AND. i.LT.m ) THEN
455  a( i, i+1 ) = imeps*2
456  ELSE IF( i.GT.1 ) THEN
457  a( i, i-1 ) = -imeps*2
458  END IF
459  END IF
460  240 CONTINUE
461 *
462  DO 250 i = 1, n
463  e( i, i ) = one
464  IF( i.LE.4 ) THEN
465  b( i, i ) = -one
466  IF( i.GT.2 )
467  $ b( i, i ) = one - reeps
468  IF( mod( i, 2 ).NE.0 .AND. i.LT.n ) THEN
469  b( i, i+1 ) = imeps
470  ELSE IF( i.GT.1 ) THEN
471  b( i, i-1 ) = -imeps
472  END IF
473  ELSE IF( i.LE.8 ) THEN
474  IF( i.LE.6 ) THEN
475  b( i, i ) = reeps
476  ELSE
477  b( i, i ) = -reeps
478  END IF
479  IF( mod( i, 2 ).NE.0 .AND. i.LT.n ) THEN
480  b( i, i+1 ) = one + imeps
481  ELSE IF( i.GT.1 ) THEN
482  b( i, i-1 ) = -one - imeps
483  END IF
484  ELSE
485  b( i, i ) = one - reeps
486  IF( mod( i, 2 ).NE.0 .AND. i.LT.n ) THEN
487  b( i, i+1 ) = imeps*2
488  ELSE IF( i.GT.1 ) THEN
489  b( i, i-1 ) = -imeps*2
490  END IF
491  END IF
492  250 CONTINUE
493  END IF
494 *
495 * Compute rhs (C, F)
496 *
497  CALL zgemm( 'N', 'N', m, n, m, one, a, lda, r, ldr, zero, c, ldc )
498  CALL zgemm( 'N', 'N', m, n, n, -one, l, ldl, b, ldb, one, c, ldc )
499  CALL zgemm( 'N', 'N', m, n, m, one, d, ldd, r, ldr, zero, f, ldf )
500  CALL zgemm( 'N', 'N', m, n, n, -one, l, ldl, e, lde, one, f, ldf )
501 *
502 * End of ZLATM5
503 *
subroutine zgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
ZGEMM
Definition: zgemm.f:189
logical function lde(RI, RJ, LR)
Definition: dblat2.f:2945

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine zlatm6 ( integer  TYPE,
integer  N,
complex*16, dimension( lda, * )  A,
integer  LDA,
complex*16, dimension( lda, * )  B,
complex*16, dimension( ldx, * )  X,
integer  LDX,
complex*16, dimension( ldy, * )  Y,
integer  LDY,
complex*16  ALPHA,
complex*16  BETA,
complex*16  WX,
complex*16  WY,
double precision, dimension( * )  S,
double precision, dimension( * )  DIF 
)

ZLATM6

Purpose:
 ZLATM6 generates test matrices for the generalized eigenvalue
 problem, their corresponding right and left eigenvector matrices,
 and also reciprocal condition numbers for all eigenvalues and
 the reciprocal condition numbers of eigenvectors corresponding to
 the 1th and 5th eigenvalues.

 Test Matrices
 =============

 Two kinds of test matrix pairs
          (A, B) = inverse(YH) * (Da, Db) * inverse(X)
 are used in the tests:

 Type 1:
    Da = 1+a   0    0    0    0    Db = 1   0   0   0   0
          0   2+a   0    0    0         0   1   0   0   0
          0    0   3+a   0    0         0   0   1   0   0
          0    0    0   4+a   0         0   0   0   1   0
          0    0    0    0   5+a ,      0   0   0   0   1
 and Type 2:
    Da = 1+i   0    0       0       0    Db = 1   0   0   0   0
          0   1-i   0       0       0         0   1   0   0   0
          0    0    1       0       0         0   0   1   0   0
          0    0    0 (1+a)+(1+b)i  0         0   0   0   1   0
          0    0    0       0 (1+a)-(1+b)i,   0   0   0   0   1 .

 In both cases the same inverse(YH) and inverse(X) are used to compute
 (A, B), giving the exact eigenvectors to (A,B) as (YH, X):

 YH:  =  1    0   -y    y   -y    X =  1   0  -x  -x   x
         0    1   -y    y   -y         0   1   x  -x  -x
         0    0    1    0    0         0   0   1   0   0
         0    0    0    1    0         0   0   0   1   0
         0    0    0    0    1,        0   0   0   0   1 , where

 a, b, x and y will have all values independently of each other.
Parameters
[in]TYPE
          TYPE is INTEGER
          Specifies the problem type (see futher details).
[in]N
          N is INTEGER
          Size of the matrices A and B.
[out]A
          A is COMPLEX*16 array, dimension (LDA, N).
          On exit A N-by-N is initialized according to TYPE.
[in]LDA
          LDA is INTEGER
          The leading dimension of A and of B.
[out]B
          B is COMPLEX*16 array, dimension (LDA, N).
          On exit B N-by-N is initialized according to TYPE.
[out]X
          X is COMPLEX*16 array, dimension (LDX, N).
          On exit X is the N-by-N matrix of right eigenvectors.
[in]LDX
          LDX is INTEGER
          The leading dimension of X.
[out]Y
          Y is COMPLEX*16 array, dimension (LDY, N).
          On exit Y is the N-by-N matrix of left eigenvectors.
[in]LDY
          LDY is INTEGER
          The leading dimension of Y.
[in]ALPHA
          ALPHA is COMPLEX*16
[in]BETA
          BETA is COMPLEX*16
 \verbatim
          Weighting constants for matrix A.
[in]WX
          WX is COMPLEX*16
          Constant for right eigenvector matrix.
[in]WY
          WY is COMPLEX*16
          Constant for left eigenvector matrix.
[out]S
          S is DOUBLE PRECISION array, dimension (N)
          S(i) is the reciprocal condition number for eigenvalue i.
[out]DIF
          DIF is DOUBLE PRECISION array, dimension (N)
          DIF(i) is the reciprocal condition number for eigenvector i.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 176 of file zlatm6.f.

176 *
177 * -- LAPACK computational routine (version 3.4.0) --
178 * -- LAPACK is a software package provided by Univ. of Tennessee, --
179 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
180 * November 2011
181 *
182 * .. Scalar Arguments ..
183  INTEGER lda, ldx, ldy, n, type
184  COMPLEX*16 alpha, beta, wx, wy
185 * ..
186 * .. Array Arguments ..
187  DOUBLE PRECISION dif( * ), s( * )
188  COMPLEX*16 a( lda, * ), b( lda, * ), x( ldx, * ),
189  $ y( ldy, * )
190 * ..
191 *
192 * =====================================================================
193 *
194 * .. Parameters ..
195  DOUBLE PRECISION rone, two, three
196  parameter( rone = 1.0d+0, two = 2.0d+0, three = 3.0d+0 )
197  COMPLEX*16 zero, one
198  parameter( zero = ( 0.0d+0, 0.0d+0 ),
199  $ one = ( 1.0d+0, 0.0d+0 ) )
200 * ..
201 * .. Local Scalars ..
202  INTEGER i, info, j
203 * ..
204 * .. Local Arrays ..
205  DOUBLE PRECISION rwork( 50 )
206  COMPLEX*16 work( 26 ), z( 8, 8 )
207 * ..
208 * .. Intrinsic Functions ..
209  INTRINSIC cdabs, dble, dcmplx, dconjg, sqrt
210 * ..
211 * .. External Subroutines ..
212  EXTERNAL zgesvd, zlacpy, zlakf2
213 * ..
214 * .. Executable Statements ..
215 *
216 * Generate test problem ...
217 * (Da, Db) ...
218 *
219  DO 20 i = 1, n
220  DO 10 j = 1, n
221 *
222  IF( i.EQ.j ) THEN
223  a( i, i ) = dcmplx( i ) + alpha
224  b( i, i ) = one
225  ELSE
226  a( i, j ) = zero
227  b( i, j ) = zero
228  END IF
229 *
230  10 CONTINUE
231  20 CONTINUE
232  IF( type.EQ.2 ) THEN
233  a( 1, 1 ) = dcmplx( rone, rone )
234  a( 2, 2 ) = dconjg( a( 1, 1 ) )
235  a( 3, 3 ) = one
236  a( 4, 4 ) = dcmplx( dble( one+alpha ), dble( one+beta ) )
237  a( 5, 5 ) = dconjg( a( 4, 4 ) )
238  END IF
239 *
240 * Form X and Y
241 *
242  CALL zlacpy( 'F', n, n, b, lda, y, ldy )
243  y( 3, 1 ) = -dconjg( wy )
244  y( 4, 1 ) = dconjg( wy )
245  y( 5, 1 ) = -dconjg( wy )
246  y( 3, 2 ) = -dconjg( wy )
247  y( 4, 2 ) = dconjg( wy )
248  y( 5, 2 ) = -dconjg( wy )
249 *
250  CALL zlacpy( 'F', n, n, b, lda, x, ldx )
251  x( 1, 3 ) = -wx
252  x( 1, 4 ) = -wx
253  x( 1, 5 ) = wx
254  x( 2, 3 ) = wx
255  x( 2, 4 ) = -wx
256  x( 2, 5 ) = -wx
257 *
258 * Form (A, B)
259 *
260  b( 1, 3 ) = wx + wy
261  b( 2, 3 ) = -wx + wy
262  b( 1, 4 ) = wx - wy
263  b( 2, 4 ) = wx - wy
264  b( 1, 5 ) = -wx + wy
265  b( 2, 5 ) = wx + wy
266  a( 1, 3 ) = wx*a( 1, 1 ) + wy*a( 3, 3 )
267  a( 2, 3 ) = -wx*a( 2, 2 ) + wy*a( 3, 3 )
268  a( 1, 4 ) = wx*a( 1, 1 ) - wy*a( 4, 4 )
269  a( 2, 4 ) = wx*a( 2, 2 ) - wy*a( 4, 4 )
270  a( 1, 5 ) = -wx*a( 1, 1 ) + wy*a( 5, 5 )
271  a( 2, 5 ) = wx*a( 2, 2 ) + wy*a( 5, 5 )
272 *
273 * Compute condition numbers
274 *
275  s( 1 ) = rone / sqrt( ( rone+three*cdabs( wy )*cdabs( wy ) ) /
276  $ ( rone+cdabs( a( 1, 1 ) )*cdabs( a( 1, 1 ) ) ) )
277  s( 2 ) = rone / sqrt( ( rone+three*cdabs( wy )*cdabs( wy ) ) /
278  $ ( rone+cdabs( a( 2, 2 ) )*cdabs( a( 2, 2 ) ) ) )
279  s( 3 ) = rone / sqrt( ( rone+two*cdabs( wx )*cdabs( wx ) ) /
280  $ ( rone+cdabs( a( 3, 3 ) )*cdabs( a( 3, 3 ) ) ) )
281  s( 4 ) = rone / sqrt( ( rone+two*cdabs( wx )*cdabs( wx ) ) /
282  $ ( rone+cdabs( a( 4, 4 ) )*cdabs( a( 4, 4 ) ) ) )
283  s( 5 ) = rone / sqrt( ( rone+two*cdabs( wx )*cdabs( wx ) ) /
284  $ ( rone+cdabs( a( 5, 5 ) )*cdabs( a( 5, 5 ) ) ) )
285 *
286  CALL zlakf2( 1, 4, a, lda, a( 2, 2 ), b, b( 2, 2 ), z, 8 )
287  CALL zgesvd( 'N', 'N', 8, 8, z, 8, rwork, work, 1, work( 2 ), 1,
288  $ work( 3 ), 24, rwork( 9 ), info )
289  dif( 1 ) = rwork( 8 )
290 *
291  CALL zlakf2( 4, 1, a, lda, a( 5, 5 ), b, b( 5, 5 ), z, 8 )
292  CALL zgesvd( 'N', 'N', 8, 8, z, 8, rwork, work, 1, work( 2 ), 1,
293  $ work( 3 ), 24, rwork( 9 ), info )
294  dif( 5 ) = rwork( 8 )
295 *
296  RETURN
297 *
298 * End of ZLATM6
299 *
subroutine zlacpy(UPLO, M, N, A, LDA, B, LDB)
ZLACPY copies all or part of one two-dimensional array to another.
Definition: zlacpy.f:105
subroutine zlakf2(M, N, A, LDA, B, D, E, Z, LDZ)
ZLAKF2
Definition: zlakf2.f:107
subroutine zgesvd(JOBU, JOBVT, M, N, A, LDA, S, U, LDU, VT, LDVT, WORK, LWORK, RWORK, INFO)
ZGESVD computes the singular value decomposition (SVD) for GE matrices
Definition: zgesvd.f:216

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine zlatme ( integer  N,
character  DIST,
integer, dimension( 4 )  ISEED,
complex*16, dimension( * )  D,
integer  MODE,
double precision  COND,
complex*16  DMAX,
character  RSIGN,
character  UPPER,
character  SIM,
double precision, dimension( * )  DS,
integer  MODES,
double precision  CONDS,
integer  KL,
integer  KU,
double precision  ANORM,
complex*16, dimension( lda, * )  A,
integer  LDA,
complex*16, dimension( * )  WORK,
integer  INFO 
)

ZLATME

Purpose:
    ZLATME generates random non-symmetric square matrices with
    specified eigenvalues for testing LAPACK programs.

    ZLATME operates by applying the following sequence of
    operations:

    1. Set the diagonal to D, where D may be input or
         computed according to MODE, COND, DMAX, and RSIGN
         as described below.

    2. If UPPER='T', the upper triangle of A is set to random values
         out of distribution DIST.

    3. If SIM='T', A is multiplied on the left by a random matrix
         X, whose singular values are specified by DS, MODES, and
         CONDS, and on the right by X inverse.

    4. If KL < N-1, the lower bandwidth is reduced to KL using
         Householder transformations.  If KU < N-1, the upper
         bandwidth is reduced to KU.

    5. If ANORM is not negative, the matrix is scaled to have
         maximum-element-norm ANORM.

    (Note: since the matrix cannot be reduced beyond Hessenberg form,
     no packing options are available.)
Parameters
[in]N
          N is INTEGER
           The number of columns (or rows) of A. Not modified.
[in]DIST
          DIST is CHARACTER*1
           On entry, DIST specifies the type of distribution to be used
           to generate the random eigen-/singular values, and on the
           upper triangle (see UPPER).
           'U' => UNIFORM( 0, 1 )  ( 'U' for uniform )
           'S' => UNIFORM( -1, 1 ) ( 'S' for symmetric )
           'N' => NORMAL( 0, 1 )   ( 'N' for normal )
           'D' => uniform on the complex disc |z| < 1.
           Not modified.
[in,out]ISEED
          ISEED is INTEGER array, dimension ( 4 )
           On entry ISEED specifies the seed of the random number
           generator. They should lie between 0 and 4095 inclusive,
           and ISEED(4) should be odd. The random number generator
           uses a linear congruential sequence limited to small
           integers, and so should produce machine independent
           random numbers. The values of ISEED are changed on
           exit, and can be used in the next call to ZLATME
           to continue the same random number sequence.
           Changed on exit.
[in,out]D
          D is COMPLEX*16 array, dimension ( N )
           This array is used to specify the eigenvalues of A.  If
           MODE=0, then D is assumed to contain the eigenvalues
           otherwise they will be computed according to MODE, COND,
           DMAX, and RSIGN and placed in D.
           Modified if MODE is nonzero.
[in]MODE
          MODE is INTEGER
           On entry this describes how the eigenvalues are to
           be specified:
           MODE = 0 means use D as input
           MODE = 1 sets D(1)=1 and D(2:N)=1.0/COND
           MODE = 2 sets D(1:N-1)=1 and D(N)=1.0/COND
           MODE = 3 sets D(I)=COND**(-(I-1)/(N-1))
           MODE = 4 sets D(i)=1 - (i-1)/(N-1)*(1 - 1/COND)
           MODE = 5 sets D to random numbers in the range
                    ( 1/COND , 1 ) such that their logarithms
                    are uniformly distributed.
           MODE = 6 set D to random numbers from same distribution
                    as the rest of the matrix.
           MODE < 0 has the same meaning as ABS(MODE), except that
              the order of the elements of D is reversed.
           Thus if MODE is between 1 and 4, D has entries ranging
              from 1 to 1/COND, if between -1 and -4, D has entries
              ranging from 1/COND to 1,
           Not modified.
[in]COND
          COND is DOUBLE PRECISION
           On entry, this is used as described under MODE above.
           If used, it must be >= 1. Not modified.
[in]DMAX
          DMAX is COMPLEX*16
           If MODE is neither -6, 0 nor 6, the contents of D, as
           computed according to MODE and COND, will be scaled by
           DMAX / max(abs(D(i))).  Note that DMAX need not be
           positive or real: if DMAX is negative or complex (or zero),
           D will be scaled by a negative or complex number (or zero).
           If RSIGN='F' then the largest (absolute) eigenvalue will be
           equal to DMAX.
           Not modified.
[in]RSIGN
          RSIGN is CHARACTER*1
           If MODE is not 0, 6, or -6, and RSIGN='T', then the
           elements of D, as computed according to MODE and COND, will
           be multiplied by a random complex number from the unit
           circle |z| = 1.  If RSIGN='F', they will not be.  RSIGN may
           only have the values 'T' or 'F'.
           Not modified.
[in]UPPER
          UPPER is CHARACTER*1
           If UPPER='T', then the elements of A above the diagonal
           will be set to random numbers out of DIST.  If UPPER='F',
           they will not.  UPPER may only have the values 'T' or 'F'.
           Not modified.
[in]SIM
          SIM is CHARACTER*1
           If SIM='T', then A will be operated on by a "similarity
           transform", i.e., multiplied on the left by a matrix X and
           on the right by X inverse.  X = U S V, where U and V are
           random unitary matrices and S is a (diagonal) matrix of
           singular values specified by DS, MODES, and CONDS.  If
           SIM='F', then A will not be transformed.
           Not modified.
[in,out]DS
          DS is DOUBLE PRECISION array, dimension ( N )
           This array is used to specify the singular values of X,
           in the same way that D specifies the eigenvalues of A.
           If MODE=0, the DS contains the singular values, which
           may not be zero.
           Modified if MODE is nonzero.
[in]MODES
          MODES is INTEGER
[in]CONDS
          CONDS is DOUBLE PRECISION
           Similar to MODE and COND, but for specifying the diagonal
           of S.  MODES=-6 and +6 are not allowed (since they would
           result in randomly ill-conditioned eigenvalues.)
[in]KL
          KL is INTEGER
           This specifies the lower bandwidth of the  matrix.  KL=1
           specifies upper Hessenberg form.  If KL is at least N-1,
           then A will have full lower bandwidth.
           Not modified.
[in]KU
          KU is INTEGER
           This specifies the upper bandwidth of the  matrix.  KU=1
           specifies lower Hessenberg form.  If KU is at least N-1,
           then A will have full upper bandwidth; if KU and KL
           are both at least N-1, then A will be dense.  Only one of
           KU and KL may be less than N-1.
           Not modified.
[in]ANORM
          ANORM is DOUBLE PRECISION
           If ANORM is not negative, then A will be scaled by a non-
           negative real number to make the maximum-element-norm of A
           to be ANORM.
           Not modified.
[out]A
          A is COMPLEX*16 array, dimension ( LDA, N )
           On exit A is the desired test matrix.
           Modified.
[in]LDA
          LDA is INTEGER
           LDA specifies the first dimension of A as declared in the
           calling program.  LDA must be at least M.
           Not modified.
[out]WORK
          WORK is COMPLEX*16 array, dimension ( 3*N )
           Workspace.
           Modified.
[out]INFO
          INFO is INTEGER
           Error code.  On exit, INFO will be set to one of the
           following values:
             0 => normal return
            -1 => N negative
            -2 => DIST illegal string
            -5 => MODE not in range -6 to 6
            -6 => COND less than 1.0, and MODE neither -6, 0 nor 6
            -9 => RSIGN is not 'T' or 'F'
           -10 => UPPER is not 'T' or 'F'
           -11 => SIM   is not 'T' or 'F'
           -12 => MODES=0 and DS has a zero singular value.
           -13 => MODES is not in the range -5 to 5.
           -14 => MODES is nonzero and CONDS is less than 1.
           -15 => KL is less than 1.
           -16 => KU is less than 1, or KL and KU are both less than
                  N-1.
           -19 => LDA is less than M.
            1  => Error return from ZLATM1 (computing D)
            2  => Cannot scale to DMAX (max. eigenvalue is 0)
            3  => Error return from DLATM1 (computing DS)
            4  => Error return from ZLARGE
            5  => Zero singular value from DLATM1.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 303 of file zlatme.f.

303 *
304 * -- LAPACK computational routine (version 3.4.0) --
305 * -- LAPACK is a software package provided by Univ. of Tennessee, --
306 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
307 * November 2011
308 *
309 * .. Scalar Arguments ..
310  CHARACTER dist, rsign, sim, upper
311  INTEGER info, kl, ku, lda, mode, modes, n
312  DOUBLE PRECISION anorm, cond, conds
313  COMPLEX*16 dmax
314 * ..
315 * .. Array Arguments ..
316  INTEGER iseed( 4 )
317  DOUBLE PRECISION ds( * )
318  COMPLEX*16 a( lda, * ), d( * ), work( * )
319 * ..
320 *
321 * =====================================================================
322 *
323 * .. Parameters ..
324  DOUBLE PRECISION zero
325  parameter( zero = 0.0d+0 )
326  DOUBLE PRECISION one
327  parameter( one = 1.0d+0 )
328  COMPLEX*16 czero
329  parameter( czero = ( 0.0d+0, 0.0d+0 ) )
330  COMPLEX*16 cone
331  parameter( cone = ( 1.0d+0, 0.0d+0 ) )
332 * ..
333 * .. Local Scalars ..
334  LOGICAL bads
335  INTEGER i, ic, icols, idist, iinfo, ir, irows, irsign,
336  $ isim, iupper, j, jc, jcr
337  DOUBLE PRECISION ralpha, temp
338  COMPLEX*16 alpha, tau, xnorms
339 * ..
340 * .. Local Arrays ..
341  DOUBLE PRECISION tempa( 1 )
342 * ..
343 * .. External Functions ..
344  LOGICAL lsame
345  DOUBLE PRECISION zlange
346  COMPLEX*16 zlarnd
347  EXTERNAL lsame, zlange, zlarnd
348 * ..
349 * .. External Subroutines ..
350  EXTERNAL dlatm1, xerbla, zcopy, zdscal, zgemv, zgerc,
352  $ zscal
353 * ..
354 * .. Intrinsic Functions ..
355  INTRINSIC abs, dconjg, max, mod
356 * ..
357 * .. Executable Statements ..
358 *
359 * 1) Decode and Test the input parameters.
360 * Initialize flags & seed.
361 *
362  info = 0
363 *
364 * Quick return if possible
365 *
366  IF( n.EQ.0 )
367  $ RETURN
368 *
369 * Decode DIST
370 *
371  IF( lsame( dist, 'U' ) ) THEN
372  idist = 1
373  ELSE IF( lsame( dist, 'S' ) ) THEN
374  idist = 2
375  ELSE IF( lsame( dist, 'N' ) ) THEN
376  idist = 3
377  ELSE IF( lsame( dist, 'D' ) ) THEN
378  idist = 4
379  ELSE
380  idist = -1
381  END IF
382 *
383 * Decode RSIGN
384 *
385  IF( lsame( rsign, 'T' ) ) THEN
386  irsign = 1
387  ELSE IF( lsame( rsign, 'F' ) ) THEN
388  irsign = 0
389  ELSE
390  irsign = -1
391  END IF
392 *
393 * Decode UPPER
394 *
395  IF( lsame( upper, 'T' ) ) THEN
396  iupper = 1
397  ELSE IF( lsame( upper, 'F' ) ) THEN
398  iupper = 0
399  ELSE
400  iupper = -1
401  END IF
402 *
403 * Decode SIM
404 *
405  IF( lsame( sim, 'T' ) ) THEN
406  isim = 1
407  ELSE IF( lsame( sim, 'F' ) ) THEN
408  isim = 0
409  ELSE
410  isim = -1
411  END IF
412 *
413 * Check DS, if MODES=0 and ISIM=1
414 *
415  bads = .false.
416  IF( modes.EQ.0 .AND. isim.EQ.1 ) THEN
417  DO 10 j = 1, n
418  IF( ds( j ).EQ.zero )
419  $ bads = .true.
420  10 CONTINUE
421  END IF
422 *
423 * Set INFO if an error
424 *
425  IF( n.LT.0 ) THEN
426  info = -1
427  ELSE IF( idist.EQ.-1 ) THEN
428  info = -2
429  ELSE IF( abs( mode ).GT.6 ) THEN
430  info = -5
431  ELSE IF( ( mode.NE.0 .AND. abs( mode ).NE.6 ) .AND. cond.LT.one )
432  $ THEN
433  info = -6
434  ELSE IF( irsign.EQ.-1 ) THEN
435  info = -9
436  ELSE IF( iupper.EQ.-1 ) THEN
437  info = -10
438  ELSE IF( isim.EQ.-1 ) THEN
439  info = -11
440  ELSE IF( bads ) THEN
441  info = -12
442  ELSE IF( isim.EQ.1 .AND. abs( modes ).GT.5 ) THEN
443  info = -13
444  ELSE IF( isim.EQ.1 .AND. modes.NE.0 .AND. conds.LT.one ) THEN
445  info = -14
446  ELSE IF( kl.LT.1 ) THEN
447  info = -15
448  ELSE IF( ku.LT.1 .OR. ( ku.LT.n-1 .AND. kl.LT.n-1 ) ) THEN
449  info = -16
450  ELSE IF( lda.LT.max( 1, n ) ) THEN
451  info = -19
452  END IF
453 *
454  IF( info.NE.0 ) THEN
455  CALL xerbla( 'ZLATME', -info )
456  RETURN
457  END IF
458 *
459 * Initialize random number generator
460 *
461  DO 20 i = 1, 4
462  iseed( i ) = mod( abs( iseed( i ) ), 4096 )
463  20 CONTINUE
464 *
465  IF( mod( iseed( 4 ), 2 ).NE.1 )
466  $ iseed( 4 ) = iseed( 4 ) + 1
467 *
468 * 2) Set up diagonal of A
469 *
470 * Compute D according to COND and MODE
471 *
472  CALL zlatm1( mode, cond, irsign, idist, iseed, d, n, iinfo )
473  IF( iinfo.NE.0 ) THEN
474  info = 1
475  RETURN
476  END IF
477  IF( mode.NE.0 .AND. abs( mode ).NE.6 ) THEN
478 *
479 * Scale by DMAX
480 *
481  temp = abs( d( 1 ) )
482  DO 30 i = 2, n
483  temp = max( temp, abs( d( i ) ) )
484  30 CONTINUE
485 *
486  IF( temp.GT.zero ) THEN
487  alpha = dmax / temp
488  ELSE
489  info = 2
490  RETURN
491  END IF
492 *
493  CALL zscal( n, alpha, d, 1 )
494 *
495  END IF
496 *
497  CALL zlaset( 'Full', n, n, czero, czero, a, lda )
498  CALL zcopy( n, d, 1, a, lda+1 )
499 *
500 * 3) If UPPER='T', set upper triangle of A to random numbers.
501 *
502  IF( iupper.NE.0 ) THEN
503  DO 40 jc = 2, n
504  CALL zlarnv( idist, iseed, jc-1, a( 1, jc ) )
505  40 CONTINUE
506  END IF
507 *
508 * 4) If SIM='T', apply similarity transformation.
509 *
510 * -1
511 * Transform is X A X , where X = U S V, thus
512 *
513 * it is U S V A V' (1/S) U'
514 *
515  IF( isim.NE.0 ) THEN
516 *
517 * Compute S (singular values of the eigenvector matrix)
518 * according to CONDS and MODES
519 *
520  CALL dlatm1( modes, conds, 0, 0, iseed, ds, n, iinfo )
521  IF( iinfo.NE.0 ) THEN
522  info = 3
523  RETURN
524  END IF
525 *
526 * Multiply by V and V'
527 *
528  CALL zlarge( n, a, lda, iseed, work, iinfo )
529  IF( iinfo.NE.0 ) THEN
530  info = 4
531  RETURN
532  END IF
533 *
534 * Multiply by S and (1/S)
535 *
536  DO 50 j = 1, n
537  CALL zdscal( n, ds( j ), a( j, 1 ), lda )
538  IF( ds( j ).NE.zero ) THEN
539  CALL zdscal( n, one / ds( j ), a( 1, j ), 1 )
540  ELSE
541  info = 5
542  RETURN
543  END IF
544  50 CONTINUE
545 *
546 * Multiply by U and U'
547 *
548  CALL zlarge( n, a, lda, iseed, work, iinfo )
549  IF( iinfo.NE.0 ) THEN
550  info = 4
551  RETURN
552  END IF
553  END IF
554 *
555 * 5) Reduce the bandwidth.
556 *
557  IF( kl.LT.n-1 ) THEN
558 *
559 * Reduce bandwidth -- kill column
560 *
561  DO 60 jcr = kl + 1, n - 1
562  ic = jcr - kl
563  irows = n + 1 - jcr
564  icols = n + kl - jcr
565 *
566  CALL zcopy( irows, a( jcr, ic ), 1, work, 1 )
567  xnorms = work( 1 )
568  CALL zlarfg( irows, xnorms, work( 2 ), 1, tau )
569  tau = dconjg( tau )
570  work( 1 ) = cone
571  alpha = zlarnd( 5, iseed )
572 *
573  CALL zgemv( 'C', irows, icols, cone, a( jcr, ic+1 ), lda,
574  $ work, 1, czero, work( irows+1 ), 1 )
575  CALL zgerc( irows, icols, -tau, work, 1, work( irows+1 ), 1,
576  $ a( jcr, ic+1 ), lda )
577 *
578  CALL zgemv( 'N', n, irows, cone, a( 1, jcr ), lda, work, 1,
579  $ czero, work( irows+1 ), 1 )
580  CALL zgerc( n, irows, -dconjg( tau ), work( irows+1 ), 1,
581  $ work, 1, a( 1, jcr ), lda )
582 *
583  a( jcr, ic ) = xnorms
584  CALL zlaset( 'Full', irows-1, 1, czero, czero,
585  $ a( jcr+1, ic ), lda )
586 *
587  CALL zscal( icols+1, alpha, a( jcr, ic ), lda )
588  CALL zscal( n, dconjg( alpha ), a( 1, jcr ), 1 )
589  60 CONTINUE
590  ELSE IF( ku.LT.n-1 ) THEN
591 *
592 * Reduce upper bandwidth -- kill a row at a time.
593 *
594  DO 70 jcr = ku + 1, n - 1
595  ir = jcr - ku
596  irows = n + ku - jcr
597  icols = n + 1 - jcr
598 *
599  CALL zcopy( icols, a( ir, jcr ), lda, work, 1 )
600  xnorms = work( 1 )
601  CALL zlarfg( icols, xnorms, work( 2 ), 1, tau )
602  tau = dconjg( tau )
603  work( 1 ) = cone
604  CALL zlacgv( icols-1, work( 2 ), 1 )
605  alpha = zlarnd( 5, iseed )
606 *
607  CALL zgemv( 'N', irows, icols, cone, a( ir+1, jcr ), lda,
608  $ work, 1, czero, work( icols+1 ), 1 )
609  CALL zgerc( irows, icols, -tau, work( icols+1 ), 1, work, 1,
610  $ a( ir+1, jcr ), lda )
611 *
612  CALL zgemv( 'C', icols, n, cone, a( jcr, 1 ), lda, work, 1,
613  $ czero, work( icols+1 ), 1 )
614  CALL zgerc( icols, n, -dconjg( tau ), work, 1,
615  $ work( icols+1 ), 1, a( jcr, 1 ), lda )
616 *
617  a( ir, jcr ) = xnorms
618  CALL zlaset( 'Full', 1, icols-1, czero, czero,
619  $ a( ir, jcr+1 ), lda )
620 *
621  CALL zscal( irows+1, alpha, a( ir, jcr ), 1 )
622  CALL zscal( n, dconjg( alpha ), a( jcr, 1 ), lda )
623  70 CONTINUE
624  END IF
625 *
626 * Scale the matrix to have norm ANORM
627 *
628  IF( anorm.GE.zero ) THEN
629  temp = zlange( 'M', n, n, a, lda, tempa )
630  IF( temp.GT.zero ) THEN
631  ralpha = anorm / temp
632  DO 80 j = 1, n
633  CALL zdscal( n, ralpha, a( 1, j ), 1 )
634  80 CONTINUE
635  END IF
636  END IF
637 *
638  RETURN
639 *
640 * End of ZLATME
641 *
subroutine dlatm1(MODE, COND, IRSIGN, IDIST, ISEED, D, N, INFO)
DLATM1
Definition: dlatm1.f:137
subroutine zcopy(N, ZX, INCX, ZY, INCY)
ZCOPY
Definition: zcopy.f:52
complex *16 function zlarnd(IDIST, ISEED)
ZLARND
Definition: zlarnd.f:77
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine zlatm1(MODE, COND, IRSIGN, IDIST, ISEED, D, N, INFO)
ZLATM1
Definition: zlatm1.f:139
subroutine zlacgv(N, X, INCX)
ZLACGV conjugates a complex vector.
Definition: zlacgv.f:76
subroutine zscal(N, ZA, ZX, INCX)
ZSCAL
Definition: zscal.f:54
subroutine zlaset(UPLO, M, N, ALPHA, BETA, A, LDA)
ZLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values...
Definition: zlaset.f:108
subroutine zgemv(TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
ZGEMV
Definition: zgemv.f:160
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine zlarnv(IDIST, ISEED, N, X)
ZLARNV returns a vector of random numbers from a uniform or normal distribution.
Definition: zlarnv.f:101
double precision function zlange(NORM, M, N, A, LDA, WORK)
ZLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition: zlange.f:117
subroutine zlarge(N, A, LDA, ISEED, WORK, INFO)
ZLARGE
Definition: zlarge.f:89
subroutine zlarfg(N, ALPHA, X, INCX, TAU)
ZLARFG generates an elementary reflector (Householder matrix).
Definition: zlarfg.f:108
subroutine zgerc(M, N, ALPHA, X, INCX, Y, INCY, A, LDA)
ZGERC
Definition: zgerc.f:132
subroutine zdscal(N, DA, ZX, INCX)
ZDSCAL
Definition: zdscal.f:54

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine zlatmr ( integer  M,
integer  N,
character  DIST,
integer, dimension( 4 )  ISEED,
character  SYM,
complex*16, dimension( * )  D,
integer  MODE,
double precision  COND,
complex*16  DMAX,
character  RSIGN,
character  GRADE,
complex*16, dimension( * )  DL,
integer  MODEL,
double precision  CONDL,
complex*16, dimension( * )  DR,
integer  MODER,
double precision  CONDR,
character  PIVTNG,
integer, dimension( * )  IPIVOT,
integer  KL,
integer  KU,
double precision  SPARSE,
double precision  ANORM,
character  PACK,
complex*16, dimension( lda, * )  A,
integer  LDA,
integer, dimension( * )  IWORK,
integer  INFO 
)

ZLATMR

Purpose:
    ZLATMR generates random matrices of various types for testing
    LAPACK programs.

    ZLATMR operates by applying the following sequence of
    operations:

      Generate a matrix A with random entries of distribution DIST
         which is symmetric if SYM='S', Hermitian if SYM='H', and
         nonsymmetric if SYM='N'.

      Set the diagonal to D, where D may be input or
         computed according to MODE, COND, DMAX and RSIGN
         as described below.

      Grade the matrix, if desired, from the left and/or right
         as specified by GRADE. The inputs DL, MODEL, CONDL, DR,
         MODER and CONDR also determine the grading as described
         below.

      Permute, if desired, the rows and/or columns as specified by
         PIVTNG and IPIVOT.

      Set random entries to zero, if desired, to get a random sparse
         matrix as specified by SPARSE.

      Make A a band matrix, if desired, by zeroing out the matrix
         outside a band of lower bandwidth KL and upper bandwidth KU.

      Scale A, if desired, to have maximum entry ANORM.

      Pack the matrix if desired. Options specified by PACK are:
         no packing
         zero out upper half (if symmetric or Hermitian)
         zero out lower half (if symmetric or Hermitian)
         store the upper half columnwise (if symmetric or Hermitian
             or square upper triangular)
         store the lower half columnwise (if symmetric or Hermitian
             or square lower triangular)
             same as upper half rowwise if symmetric
             same as conjugate upper half rowwise if Hermitian
         store the lower triangle in banded format
             (if symmetric or Hermitian)
         store the upper triangle in banded format
             (if symmetric or Hermitian)
         store the entire matrix in banded format

    Note: If two calls to ZLATMR differ only in the PACK parameter,
          they will generate mathematically equivalent matrices.

          If two calls to ZLATMR both have full bandwidth (KL = M-1
          and KU = N-1), and differ only in the PIVTNG and PACK
          parameters, then the matrices generated will differ only
          in the order of the rows and/or columns, and otherwise
          contain the same data. This consistency cannot be and
          is not maintained with less than full bandwidth.
Parameters
[in]M
          M is INTEGER
           Number of rows of A. Not modified.
[in]N
          N is INTEGER
           Number of columns of A. Not modified.
[in]DIST
          DIST is CHARACTER*1
           On entry, DIST specifies the type of distribution to be used
           to generate a random matrix .
           'U' => real and imaginary parts are independent
                  UNIFORM( 0, 1 )  ( 'U' for uniform )
           'S' => real and imaginary parts are independent
                  UNIFORM( -1, 1 ) ( 'S' for symmetric )
           'N' => real and imaginary parts are independent
                  NORMAL( 0, 1 )   ( 'N' for normal )
           'D' => uniform on interior of unit disk ( 'D' for disk )
           Not modified.
[in,out]ISEED
          ISEED is INTEGER array, dimension (4)
           On entry ISEED specifies the seed of the random number
           generator. They should lie between 0 and 4095 inclusive,
           and ISEED(4) should be odd. The random number generator
           uses a linear congruential sequence limited to small
           integers, and so should produce machine independent
           random numbers. The values of ISEED are changed on
           exit, and can be used in the next call to ZLATMR
           to continue the same random number sequence.
           Changed on exit.
[in]SYM
          SYM is CHARACTER*1
           If SYM='S', generated matrix is symmetric.
           If SYM='H', generated matrix is Hermitian.
           If SYM='N', generated matrix is nonsymmetric.
           Not modified.
[in,out]D
          D is COMPLEX*16 array, dimension (min(M,N))
           On entry this array specifies the diagonal entries
           of the diagonal of A.  D may either be specified
           on entry, or set according to MODE and COND as described
           below. If the matrix is Hermitian, the real part of D
           will be taken. May be changed on exit if MODE is nonzero.
[in]MODE
          MODE is INTEGER
           On entry describes how D is to be used:
           MODE = 0 means use D as input
           MODE = 1 sets D(1)=1 and D(2:N)=1.0/COND
           MODE = 2 sets D(1:N-1)=1 and D(N)=1.0/COND
           MODE = 3 sets D(I)=COND**(-(I-1)/(N-1))
           MODE = 4 sets D(i)=1 - (i-1)/(N-1)*(1 - 1/COND)
           MODE = 5 sets D to random numbers in the range
                    ( 1/COND , 1 ) such that their logarithms
                    are uniformly distributed.
           MODE = 6 set D to random numbers from same distribution
                    as the rest of the matrix.
           MODE < 0 has the same meaning as ABS(MODE), except that
              the order of the elements of D is reversed.
           Thus if MODE is positive, D has entries ranging from
              1 to 1/COND, if negative, from 1/COND to 1,
           Not modified.
[in]COND
          COND is DOUBLE PRECISION
           On entry, used as described under MODE above.
           If used, it must be >= 1. Not modified.
[in]DMAX
          DMAX is COMPLEX*16
           If MODE neither -6, 0 nor 6, the diagonal is scaled by
           DMAX / max(abs(D(i))), so that maximum absolute entry
           of diagonal is abs(DMAX). If DMAX is complex (or zero),
           diagonal will be scaled by a complex number (or zero).
[in]RSIGN
          RSIGN is CHARACTER*1
           If MODE neither -6, 0 nor 6, specifies sign of diagonal
           as follows:
           'T' => diagonal entries are multiplied by a random complex
                  number uniformly distributed with absolute value 1
           'F' => diagonal unchanged
           Not modified.
[in]GRADE
          GRADE is CHARACTER*1
           Specifies grading of matrix as follows:
           'N'  => no grading
           'L'  => matrix premultiplied by diag( DL )
                   (only if matrix nonsymmetric)
           'R'  => matrix postmultiplied by diag( DR )
                   (only if matrix nonsymmetric)
           'B'  => matrix premultiplied by diag( DL ) and
                         postmultiplied by diag( DR )
                   (only if matrix nonsymmetric)
           'H'  => matrix premultiplied by diag( DL ) and
                         postmultiplied by diag( CONJG(DL) )
                   (only if matrix Hermitian or nonsymmetric)
           'S'  => matrix premultiplied by diag( DL ) and
                         postmultiplied by diag( DL )
                   (only if matrix symmetric or nonsymmetric)
           'E'  => matrix premultiplied by diag( DL ) and
                         postmultiplied by inv( diag( DL ) )
                         ( 'S' for similarity )
                   (only if matrix nonsymmetric)
                   Note: if GRADE='S', then M must equal N.
           Not modified.
[in,out]DL
          DL is COMPLEX*16 array, dimension (M)
           If MODEL=0, then on entry this array specifies the diagonal
           entries of a diagonal matrix used as described under GRADE
           above. If MODEL is not zero, then DL will be set according
           to MODEL and CONDL, analogous to the way D is set according
           to MODE and COND (except there is no DMAX parameter for DL).
           If GRADE='E', then DL cannot have zero entries.
           Not referenced if GRADE = 'N' or 'R'. Changed on exit.
[in]MODEL
          MODEL is INTEGER
           This specifies how the diagonal array DL is to be computed,
           just as MODE specifies how D is to be computed.
           Not modified.
[in]CONDL
          CONDL is DOUBLE PRECISION
           When MODEL is not zero, this specifies the condition number
           of the computed DL.  Not modified.
[in,out]DR
          DR is COMPLEX*16 array, dimension (N)
           If MODER=0, then on entry this array specifies the diagonal
           entries of a diagonal matrix used as described under GRADE
           above. If MODER is not zero, then DR will be set according
           to MODER and CONDR, analogous to the way D is set according
           to MODE and COND (except there is no DMAX parameter for DR).
           Not referenced if GRADE = 'N', 'L', 'H' or 'S'.
           Changed on exit.
[in]MODER
          MODER is INTEGER
           This specifies how the diagonal array DR is to be computed,
           just as MODE specifies how D is to be computed.
           Not modified.
[in]CONDR
          CONDR is DOUBLE PRECISION
           When MODER is not zero, this specifies the condition number
           of the computed DR.  Not modified.
[in]PIVTNG
          PIVTNG is CHARACTER*1
           On entry specifies pivoting permutations as follows:
           'N' or ' ' => none.
           'L' => left or row pivoting (matrix must be nonsymmetric).
           'R' => right or column pivoting (matrix must be
                  nonsymmetric).
           'B' or 'F' => both or full pivoting, i.e., on both sides.
                         In this case, M must equal N

           If two calls to ZLATMR both have full bandwidth (KL = M-1
           and KU = N-1), and differ only in the PIVTNG and PACK
           parameters, then the matrices generated will differ only
           in the order of the rows and/or columns, and otherwise
           contain the same data. This consistency cannot be
           maintained with less than full bandwidth.
[in]IPIVOT
          IPIVOT is INTEGER array, dimension (N or M)
           This array specifies the permutation used.  After the
           basic matrix is generated, the rows, columns, or both
           are permuted.   If, say, row pivoting is selected, ZLATMR
           starts with the *last* row and interchanges the M-th and
           IPIVOT(M)-th rows, then moves to the next-to-last row,
           interchanging the (M-1)-th and the IPIVOT(M-1)-th rows,
           and so on.  In terms of "2-cycles", the permutation is
           (1 IPIVOT(1)) (2 IPIVOT(2)) ... (M IPIVOT(M))
           where the rightmost cycle is applied first.  This is the
           *inverse* of the effect of pivoting in LINPACK.  The idea
           is that factoring (with pivoting) an identity matrix
           which has been inverse-pivoted in this way should
           result in a pivot vector identical to IPIVOT.
           Not referenced if PIVTNG = 'N'. Not modified.
[in]SPARSE
          SPARSE is DOUBLE PRECISION
           On entry specifies the sparsity of the matrix if a sparse
           matrix is to be generated. SPARSE should lie between
           0 and 1. To generate a sparse matrix, for each matrix entry
           a uniform ( 0, 1 ) random number x is generated and
           compared to SPARSE; if x is larger the matrix entry
           is unchanged and if x is smaller the entry is set
           to zero. Thus on the average a fraction SPARSE of the
           entries will be set to zero.
           Not modified.
[in]KL
          KL is INTEGER
           On entry specifies the lower bandwidth of the  matrix. For
           example, KL=0 implies upper triangular, KL=1 implies upper
           Hessenberg, and KL at least M-1 implies the matrix is not
           banded. Must equal KU if matrix is symmetric or Hermitian.
           Not modified.
[in]KU
          KU is INTEGER
           On entry specifies the upper bandwidth of the  matrix. For
           example, KU=0 implies lower triangular, KU=1 implies lower
           Hessenberg, and KU at least N-1 implies the matrix is not
           banded. Must equal KL if matrix is symmetric or Hermitian.
           Not modified.
[in]ANORM
          ANORM is DOUBLE PRECISION
           On entry specifies maximum entry of output matrix
           (output matrix will by multiplied by a constant so that
           its largest absolute entry equal ANORM)
           if ANORM is nonnegative. If ANORM is negative no scaling
           is done. Not modified.
[in]PACK
          PACK is CHARACTER*1
           On entry specifies packing of matrix as follows:
           'N' => no packing
           'U' => zero out all subdiagonal entries
                  (if symmetric or Hermitian)
           'L' => zero out all superdiagonal entries
                  (if symmetric or Hermitian)
           'C' => store the upper triangle columnwise
                  (only if matrix symmetric or Hermitian or
                   square upper triangular)
           'R' => store the lower triangle columnwise
                  (only if matrix symmetric or Hermitian or
                   square lower triangular)
                  (same as upper half rowwise if symmetric)
                  (same as conjugate upper half rowwise if Hermitian)
           'B' => store the lower triangle in band storage scheme
                  (only if matrix symmetric or Hermitian)
           'Q' => store the upper triangle in band storage scheme
                  (only if matrix symmetric or Hermitian)
           'Z' => store the entire matrix in band storage scheme
                      (pivoting can be provided for by using this
                      option to store A in the trailing rows of
                      the allocated storage)

           Using these options, the various LAPACK packed and banded
           storage schemes can be obtained:
           GB               - use 'Z'
           PB, HB or TB     - use 'B' or 'Q'
           PP, HP or TP     - use 'C' or 'R'

           If two calls to ZLATMR differ only in the PACK parameter,
           they will generate mathematically equivalent matrices.
           Not modified.
[in,out]A
          A is COMPLEX*16 array, dimension (LDA,N)
           On exit A is the desired test matrix. Only those
           entries of A which are significant on output
           will be referenced (even if A is in packed or band
           storage format). The 'unoccupied corners' of A in
           band format will be zeroed out.
[in]LDA
          LDA is INTEGER
           on entry LDA specifies the first dimension of A as
           declared in the calling program.
           If PACK='N', 'U' or 'L', LDA must be at least max ( 1, M ).
           If PACK='C' or 'R', LDA must be at least 1.
           If PACK='B', or 'Q', LDA must be MIN ( KU+1, N )
           If PACK='Z', LDA must be at least KUU+KLL+1, where
           KUU = MIN ( KU, N-1 ) and KLL = MIN ( KL, N-1 )
           Not modified.
[out]IWORK
          IWORK is INTEGER array, dimension (N or M)
           Workspace. Not referenced if PIVTNG = 'N'. Changed on exit.
[out]INFO
          INFO is INTEGER
           Error parameter on exit:
             0 => normal return
            -1 => M negative or unequal to N and SYM='S' or 'H'
            -2 => N negative
            -3 => DIST illegal string
            -5 => SYM illegal string
            -7 => MODE not in range -6 to 6
            -8 => COND less than 1.0, and MODE neither -6, 0 nor 6
           -10 => MODE neither -6, 0 nor 6 and RSIGN illegal string
           -11 => GRADE illegal string, or GRADE='E' and
                  M not equal to N, or GRADE='L', 'R', 'B', 'S' or 'E'
                  and SYM = 'H', or GRADE='L', 'R', 'B', 'H' or 'E'
                  and SYM = 'S'
           -12 => GRADE = 'E' and DL contains zero
           -13 => MODEL not in range -6 to 6 and GRADE= 'L', 'B', 'H',
                  'S' or 'E'
           -14 => CONDL less than 1.0, GRADE='L', 'B', 'H', 'S' or 'E',
                  and MODEL neither -6, 0 nor 6
           -16 => MODER not in range -6 to 6 and GRADE= 'R' or 'B'
           -17 => CONDR less than 1.0, GRADE='R' or 'B', and
                  MODER neither -6, 0 nor 6
           -18 => PIVTNG illegal string, or PIVTNG='B' or 'F' and
                  M not equal to N, or PIVTNG='L' or 'R' and SYM='S'
                  or 'H'
           -19 => IPIVOT contains out of range number and
                  PIVTNG not equal to 'N'
           -20 => KL negative
           -21 => KU negative, or SYM='S' or 'H' and KU not equal to KL
           -22 => SPARSE not in range 0. to 1.
           -24 => PACK illegal string, or PACK='U', 'L', 'B' or 'Q'
                  and SYM='N', or PACK='C' and SYM='N' and either KL
                  not equal to 0 or N not equal to M, or PACK='R' and
                  SYM='N', and either KU not equal to 0 or N not equal
                  to M
           -26 => LDA too small
             1 => Error return from ZLATM1 (computing D)
             2 => Cannot scale diagonal to DMAX (max. entry is 0)
             3 => Error return from ZLATM1 (computing DL)
             4 => Error return from ZLATM1 (computing DR)
             5 => ANORM is positive, but matrix constructed prior to
                  attempting to scale it to have norm ANORM, is zero
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 492 of file zlatmr.f.

492 *
493 * -- LAPACK computational routine (version 3.4.0) --
494 * -- LAPACK is a software package provided by Univ. of Tennessee, --
495 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
496 * November 2011
497 *
498 * .. Scalar Arguments ..
499  CHARACTER dist, grade, pack, pivtng, rsign, sym
500  INTEGER info, kl, ku, lda, m, mode, model, moder, n
501  DOUBLE PRECISION anorm, cond, condl, condr, sparse
502  COMPLEX*16 dmax
503 * ..
504 * .. Array Arguments ..
505  INTEGER ipivot( * ), iseed( 4 ), iwork( * )
506  COMPLEX*16 a( lda, * ), d( * ), dl( * ), dr( * )
507 * ..
508 *
509 * =====================================================================
510 *
511 * .. Parameters ..
512  DOUBLE PRECISION zero
513  parameter( zero = 0.0d0 )
514  DOUBLE PRECISION one
515  parameter( one = 1.0d0 )
516  COMPLEX*16 cone
517  parameter( cone = ( 1.0d0, 0.0d0 ) )
518  COMPLEX*16 czero
519  parameter( czero = ( 0.0d0, 0.0d0 ) )
520 * ..
521 * .. Local Scalars ..
522  LOGICAL badpvt, dzero, fulbnd
523  INTEGER i, idist, igrade, iisub, ipack, ipvtng, irsign,
524  $ isub, isym, j, jjsub, jsub, k, kll, kuu, mnmin,
525  $ mnsub, mxsub, npvts
526  DOUBLE PRECISION onorm, temp
527  COMPLEX*16 calpha, ctemp
528 * ..
529 * .. Local Arrays ..
530  DOUBLE PRECISION tempa( 1 )
531 * ..
532 * .. External Functions ..
533  LOGICAL lsame
534  DOUBLE PRECISION zlangb, zlange, zlansb, zlansp, zlansy
535  COMPLEX*16 zlatm2, zlatm3
536  EXTERNAL lsame, zlangb, zlange, zlansb, zlansp, zlansy,
537  $ zlatm2, zlatm3
538 * ..
539 * .. External Subroutines ..
540  EXTERNAL xerbla, zdscal, zlatm1
541 * ..
542 * .. Intrinsic Functions ..
543  INTRINSIC abs, dble, dconjg, max, min, mod
544 * ..
545 * .. Executable Statements ..
546 *
547 * 1) Decode and Test the input parameters.
548 * Initialize flags & seed.
549 *
550  info = 0
551 *
552 * Quick return if possible
553 *
554  IF( m.EQ.0 .OR. n.EQ.0 )
555  $ RETURN
556 *
557 * Decode DIST
558 *
559  IF( lsame( dist, 'U' ) ) THEN
560  idist = 1
561  ELSE IF( lsame( dist, 'S' ) ) THEN
562  idist = 2
563  ELSE IF( lsame( dist, 'N' ) ) THEN
564  idist = 3
565  ELSE IF( lsame( dist, 'D' ) ) THEN
566  idist = 4
567  ELSE
568  idist = -1
569  END IF
570 *
571 * Decode SYM
572 *
573  IF( lsame( sym, 'H' ) ) THEN
574  isym = 0
575  ELSE IF( lsame( sym, 'N' ) ) THEN
576  isym = 1
577  ELSE IF( lsame( sym, 'S' ) ) THEN
578  isym = 2
579  ELSE
580  isym = -1
581  END IF
582 *
583 * Decode RSIGN
584 *
585  IF( lsame( rsign, 'F' ) ) THEN
586  irsign = 0
587  ELSE IF( lsame( rsign, 'T' ) ) THEN
588  irsign = 1
589  ELSE
590  irsign = -1
591  END IF
592 *
593 * Decode PIVTNG
594 *
595  IF( lsame( pivtng, 'N' ) ) THEN
596  ipvtng = 0
597  ELSE IF( lsame( pivtng, ' ' ) ) THEN
598  ipvtng = 0
599  ELSE IF( lsame( pivtng, 'L' ) ) THEN
600  ipvtng = 1
601  npvts = m
602  ELSE IF( lsame( pivtng, 'R' ) ) THEN
603  ipvtng = 2
604  npvts = n
605  ELSE IF( lsame( pivtng, 'B' ) ) THEN
606  ipvtng = 3
607  npvts = min( n, m )
608  ELSE IF( lsame( pivtng, 'F' ) ) THEN
609  ipvtng = 3
610  npvts = min( n, m )
611  ELSE
612  ipvtng = -1
613  END IF
614 *
615 * Decode GRADE
616 *
617  IF( lsame( grade, 'N' ) ) THEN
618  igrade = 0
619  ELSE IF( lsame( grade, 'L' ) ) THEN
620  igrade = 1
621  ELSE IF( lsame( grade, 'R' ) ) THEN
622  igrade = 2
623  ELSE IF( lsame( grade, 'B' ) ) THEN
624  igrade = 3
625  ELSE IF( lsame( grade, 'E' ) ) THEN
626  igrade = 4
627  ELSE IF( lsame( grade, 'H' ) ) THEN
628  igrade = 5
629  ELSE IF( lsame( grade, 'S' ) ) THEN
630  igrade = 6
631  ELSE
632  igrade = -1
633  END IF
634 *
635 * Decode PACK
636 *
637  IF( lsame( pack, 'N' ) ) THEN
638  ipack = 0
639  ELSE IF( lsame( pack, 'U' ) ) THEN
640  ipack = 1
641  ELSE IF( lsame( pack, 'L' ) ) THEN
642  ipack = 2
643  ELSE IF( lsame( pack, 'C' ) ) THEN
644  ipack = 3
645  ELSE IF( lsame( pack, 'R' ) ) THEN
646  ipack = 4
647  ELSE IF( lsame( pack, 'B' ) ) THEN
648  ipack = 5
649  ELSE IF( lsame( pack, 'Q' ) ) THEN
650  ipack = 6
651  ELSE IF( lsame( pack, 'Z' ) ) THEN
652  ipack = 7
653  ELSE
654  ipack = -1
655  END IF
656 *
657 * Set certain internal parameters
658 *
659  mnmin = min( m, n )
660  kll = min( kl, m-1 )
661  kuu = min( ku, n-1 )
662 *
663 * If inv(DL) is used, check to see if DL has a zero entry.
664 *
665  dzero = .false.
666  IF( igrade.EQ.4 .AND. model.EQ.0 ) THEN
667  DO 10 i = 1, m
668  IF( dl( i ).EQ.czero )
669  $ dzero = .true.
670  10 CONTINUE
671  END IF
672 *
673 * Check values in IPIVOT
674 *
675  badpvt = .false.
676  IF( ipvtng.GT.0 ) THEN
677  DO 20 j = 1, npvts
678  IF( ipivot( j ).LE.0 .OR. ipivot( j ).GT.npvts )
679  $ badpvt = .true.
680  20 CONTINUE
681  END IF
682 *
683 * Set INFO if an error
684 *
685  IF( m.LT.0 ) THEN
686  info = -1
687  ELSE IF( m.NE.n .AND. ( isym.EQ.0 .OR. isym.EQ.2 ) ) THEN
688  info = -1
689  ELSE IF( n.LT.0 ) THEN
690  info = -2
691  ELSE IF( idist.EQ.-1 ) THEN
692  info = -3
693  ELSE IF( isym.EQ.-1 ) THEN
694  info = -5
695  ELSE IF( mode.LT.-6 .OR. mode.GT.6 ) THEN
696  info = -7
697  ELSE IF( ( mode.NE.-6 .AND. mode.NE.0 .AND. mode.NE.6 ) .AND.
698  $ cond.LT.one ) THEN
699  info = -8
700  ELSE IF( ( mode.NE.-6 .AND. mode.NE.0 .AND. mode.NE.6 ) .AND.
701  $ irsign.EQ.-1 ) THEN
702  info = -10
703  ELSE IF( igrade.EQ.-1 .OR. ( igrade.EQ.4 .AND. m.NE.n ) .OR.
704  $ ( ( igrade.EQ.1 .OR. igrade.EQ.2 .OR. igrade.EQ.3 .OR.
705  $ igrade.EQ.4 .OR. igrade.EQ.6 ) .AND. isym.EQ.0 ) .OR.
706  $ ( ( igrade.EQ.1 .OR. igrade.EQ.2 .OR. igrade.EQ.3 .OR.
707  $ igrade.EQ.4 .OR. igrade.EQ.5 ) .AND. isym.EQ.2 ) ) THEN
708  info = -11
709  ELSE IF( igrade.EQ.4 .AND. dzero ) THEN
710  info = -12
711  ELSE IF( ( igrade.EQ.1 .OR. igrade.EQ.3 .OR. igrade.EQ.4 .OR.
712  $ igrade.EQ.5 .OR. igrade.EQ.6 ) .AND.
713  $ ( model.LT.-6 .OR. model.GT.6 ) ) THEN
714  info = -13
715  ELSE IF( ( igrade.EQ.1 .OR. igrade.EQ.3 .OR. igrade.EQ.4 .OR.
716  $ igrade.EQ.5 .OR. igrade.EQ.6 ) .AND.
717  $ ( model.NE.-6 .AND. model.NE.0 .AND. model.NE.6 ) .AND.
718  $ condl.LT.one ) THEN
719  info = -14
720  ELSE IF( ( igrade.EQ.2 .OR. igrade.EQ.3 ) .AND.
721  $ ( moder.LT.-6 .OR. moder.GT.6 ) ) THEN
722  info = -16
723  ELSE IF( ( igrade.EQ.2 .OR. igrade.EQ.3 ) .AND.
724  $ ( moder.NE.-6 .AND. moder.NE.0 .AND. moder.NE.6 ) .AND.
725  $ condr.LT.one ) THEN
726  info = -17
727  ELSE IF( ipvtng.EQ.-1 .OR. ( ipvtng.EQ.3 .AND. m.NE.n ) .OR.
728  $ ( ( ipvtng.EQ.1 .OR. ipvtng.EQ.2 ) .AND. ( isym.EQ.0 .OR.
729  $ isym.EQ.2 ) ) ) THEN
730  info = -18
731  ELSE IF( ipvtng.NE.0 .AND. badpvt ) THEN
732  info = -19
733  ELSE IF( kl.LT.0 ) THEN
734  info = -20
735  ELSE IF( ku.LT.0 .OR. ( ( isym.EQ.0 .OR. isym.EQ.2 ) .AND. kl.NE.
736  $ ku ) ) THEN
737  info = -21
738  ELSE IF( sparse.LT.zero .OR. sparse.GT.one ) THEN
739  info = -22
740  ELSE IF( ipack.EQ.-1 .OR. ( ( ipack.EQ.1 .OR. ipack.EQ.2 .OR.
741  $ ipack.EQ.5 .OR. ipack.EQ.6 ) .AND. isym.EQ.1 ) .OR.
742  $ ( ipack.EQ.3 .AND. isym.EQ.1 .AND. ( kl.NE.0 .OR. m.NE.
743  $ n ) ) .OR. ( ipack.EQ.4 .AND. isym.EQ.1 .AND. ( ku.NE.
744  $ 0 .OR. m.NE.n ) ) ) THEN
745  info = -24
746  ELSE IF( ( ( ipack.EQ.0 .OR. ipack.EQ.1 .OR. ipack.EQ.2 ) .AND.
747  $ lda.LT.max( 1, m ) ) .OR. ( ( ipack.EQ.3 .OR. ipack.EQ.
748  $ 4 ) .AND. lda.LT.1 ) .OR. ( ( ipack.EQ.5 .OR. ipack.EQ.
749  $ 6 ) .AND. lda.LT.kuu+1 ) .OR.
750  $ ( ipack.EQ.7 .AND. lda.LT.kll+kuu+1 ) ) THEN
751  info = -26
752  END IF
753 *
754  IF( info.NE.0 ) THEN
755  CALL xerbla( 'ZLATMR', -info )
756  RETURN
757  END IF
758 *
759 * Decide if we can pivot consistently
760 *
761  fulbnd = .false.
762  IF( kuu.EQ.n-1 .AND. kll.EQ.m-1 )
763  $ fulbnd = .true.
764 *
765 * Initialize random number generator
766 *
767  DO 30 i = 1, 4
768  iseed( i ) = mod( abs( iseed( i ) ), 4096 )
769  30 CONTINUE
770 *
771  iseed( 4 ) = 2*( iseed( 4 ) / 2 ) + 1
772 *
773 * 2) Set up D, DL, and DR, if indicated.
774 *
775 * Compute D according to COND and MODE
776 *
777  CALL zlatm1( mode, cond, irsign, idist, iseed, d, mnmin, info )
778  IF( info.NE.0 ) THEN
779  info = 1
780  RETURN
781  END IF
782  IF( mode.NE.0 .AND. mode.NE.-6 .AND. mode.NE.6 ) THEN
783 *
784 * Scale by DMAX
785 *
786  temp = abs( d( 1 ) )
787  DO 40 i = 2, mnmin
788  temp = max( temp, abs( d( i ) ) )
789  40 CONTINUE
790  IF( temp.EQ.zero .AND. dmax.NE.czero ) THEN
791  info = 2
792  RETURN
793  END IF
794  IF( temp.NE.zero ) THEN
795  calpha = dmax / temp
796  ELSE
797  calpha = cone
798  END IF
799  DO 50 i = 1, mnmin
800  d( i ) = calpha*d( i )
801  50 CONTINUE
802 *
803  END IF
804 *
805 * If matrix Hermitian, make D real
806 *
807  IF( isym.EQ.0 ) THEN
808  DO 60 i = 1, mnmin
809  d( i ) = dble( d( i ) )
810  60 CONTINUE
811  END IF
812 *
813 * Compute DL if grading set
814 *
815  IF( igrade.EQ.1 .OR. igrade.EQ.3 .OR. igrade.EQ.4 .OR. igrade.EQ.
816  $ 5 .OR. igrade.EQ.6 ) THEN
817  CALL zlatm1( model, condl, 0, idist, iseed, dl, m, info )
818  IF( info.NE.0 ) THEN
819  info = 3
820  RETURN
821  END IF
822  END IF
823 *
824 * Compute DR if grading set
825 *
826  IF( igrade.EQ.2 .OR. igrade.EQ.3 ) THEN
827  CALL zlatm1( moder, condr, 0, idist, iseed, dr, n, info )
828  IF( info.NE.0 ) THEN
829  info = 4
830  RETURN
831  END IF
832  END IF
833 *
834 * 3) Generate IWORK if pivoting
835 *
836  IF( ipvtng.GT.0 ) THEN
837  DO 70 i = 1, npvts
838  iwork( i ) = i
839  70 CONTINUE
840  IF( fulbnd ) THEN
841  DO 80 i = 1, npvts
842  k = ipivot( i )
843  j = iwork( i )
844  iwork( i ) = iwork( k )
845  iwork( k ) = j
846  80 CONTINUE
847  ELSE
848  DO 90 i = npvts, 1, -1
849  k = ipivot( i )
850  j = iwork( i )
851  iwork( i ) = iwork( k )
852  iwork( k ) = j
853  90 CONTINUE
854  END IF
855  END IF
856 *
857 * 4) Generate matrices for each kind of PACKing
858 * Always sweep matrix columnwise (if symmetric, upper
859 * half only) so that matrix generated does not depend
860 * on PACK
861 *
862  IF( fulbnd ) THEN
863 *
864 * Use ZLATM3 so matrices generated with differing PIVOTing only
865 * differ only in the order of their rows and/or columns.
866 *
867  IF( ipack.EQ.0 ) THEN
868  IF( isym.EQ.0 ) THEN
869  DO 110 j = 1, n
870  DO 100 i = 1, j
871  ctemp = zlatm3( m, n, i, j, isub, jsub, kl, ku,
872  $ idist, iseed, d, igrade, dl, dr, ipvtng,
873  $ iwork, sparse )
874  a( isub, jsub ) = ctemp
875  a( jsub, isub ) = dconjg( ctemp )
876  100 CONTINUE
877  110 CONTINUE
878  ELSE IF( isym.EQ.1 ) THEN
879  DO 130 j = 1, n
880  DO 120 i = 1, m
881  ctemp = zlatm3( m, n, i, j, isub, jsub, kl, ku,
882  $ idist, iseed, d, igrade, dl, dr, ipvtng,
883  $ iwork, sparse )
884  a( isub, jsub ) = ctemp
885  120 CONTINUE
886  130 CONTINUE
887  ELSE IF( isym.EQ.2 ) THEN
888  DO 150 j = 1, n
889  DO 140 i = 1, j
890  ctemp = zlatm3( m, n, i, j, isub, jsub, kl, ku,
891  $ idist, iseed, d, igrade, dl, dr, ipvtng,
892  $ iwork, sparse )
893  a( isub, jsub ) = ctemp
894  a( jsub, isub ) = ctemp
895  140 CONTINUE
896  150 CONTINUE
897  END IF
898 *
899  ELSE IF( ipack.EQ.1 ) THEN
900 *
901  DO 170 j = 1, n
902  DO 160 i = 1, j
903  ctemp = zlatm3( m, n, i, j, isub, jsub, kl, ku, idist,
904  $ iseed, d, igrade, dl, dr, ipvtng, iwork,
905  $ sparse )
906  mnsub = min( isub, jsub )
907  mxsub = max( isub, jsub )
908  IF( mxsub.EQ.isub .AND. isym.EQ.0 ) THEN
909  a( mnsub, mxsub ) = dconjg( ctemp )
910  ELSE
911  a( mnsub, mxsub ) = ctemp
912  END IF
913  IF( mnsub.NE.mxsub )
914  $ a( mxsub, mnsub ) = czero
915  160 CONTINUE
916  170 CONTINUE
917 *
918  ELSE IF( ipack.EQ.2 ) THEN
919 *
920  DO 190 j = 1, n
921  DO 180 i = 1, j
922  ctemp = zlatm3( m, n, i, j, isub, jsub, kl, ku, idist,
923  $ iseed, d, igrade, dl, dr, ipvtng, iwork,
924  $ sparse )
925  mnsub = min( isub, jsub )
926  mxsub = max( isub, jsub )
927  IF( mxsub.EQ.jsub .AND. isym.EQ.0 ) THEN
928  a( mxsub, mnsub ) = dconjg( ctemp )
929  ELSE
930  a( mxsub, mnsub ) = ctemp
931  END IF
932  IF( mnsub.NE.mxsub )
933  $ a( mnsub, mxsub ) = czero
934  180 CONTINUE
935  190 CONTINUE
936 *
937  ELSE IF( ipack.EQ.3 ) THEN
938 *
939  DO 210 j = 1, n
940  DO 200 i = 1, j
941  ctemp = zlatm3( m, n, i, j, isub, jsub, kl, ku, idist,
942  $ iseed, d, igrade, dl, dr, ipvtng, iwork,
943  $ sparse )
944 *
945 * Compute K = location of (ISUB,JSUB) entry in packed
946 * array
947 *
948  mnsub = min( isub, jsub )
949  mxsub = max( isub, jsub )
950  k = mxsub*( mxsub-1 ) / 2 + mnsub
951 *
952 * Convert K to (IISUB,JJSUB) location
953 *
954  jjsub = ( k-1 ) / lda + 1
955  iisub = k - lda*( jjsub-1 )
956 *
957  IF( mxsub.EQ.isub .AND. isym.EQ.0 ) THEN
958  a( iisub, jjsub ) = dconjg( ctemp )
959  ELSE
960  a( iisub, jjsub ) = ctemp
961  END IF
962  200 CONTINUE
963  210 CONTINUE
964 *
965  ELSE IF( ipack.EQ.4 ) THEN
966 *
967  DO 230 j = 1, n
968  DO 220 i = 1, j
969  ctemp = zlatm3( m, n, i, j, isub, jsub, kl, ku, idist,
970  $ iseed, d, igrade, dl, dr, ipvtng, iwork,
971  $ sparse )
972 *
973 * Compute K = location of (I,J) entry in packed array
974 *
975  mnsub = min( isub, jsub )
976  mxsub = max( isub, jsub )
977  IF( mnsub.EQ.1 ) THEN
978  k = mxsub
979  ELSE
980  k = n*( n+1 ) / 2 - ( n-mnsub+1 )*( n-mnsub+2 ) /
981  $ 2 + mxsub - mnsub + 1
982  END IF
983 *
984 * Convert K to (IISUB,JJSUB) location
985 *
986  jjsub = ( k-1 ) / lda + 1
987  iisub = k - lda*( jjsub-1 )
988 *
989  IF( mxsub.EQ.jsub .AND. isym.EQ.0 ) THEN
990  a( iisub, jjsub ) = dconjg( ctemp )
991  ELSE
992  a( iisub, jjsub ) = ctemp
993  END IF
994  220 CONTINUE
995  230 CONTINUE
996 *
997  ELSE IF( ipack.EQ.5 ) THEN
998 *
999  DO 250 j = 1, n
1000  DO 240 i = j - kuu, j
1001  IF( i.LT.1 ) THEN
1002  a( j-i+1, i+n ) = czero
1003  ELSE
1004  ctemp = zlatm3( m, n, i, j, isub, jsub, kl, ku,
1005  $ idist, iseed, d, igrade, dl, dr, ipvtng,
1006  $ iwork, sparse )
1007  mnsub = min( isub, jsub )
1008  mxsub = max( isub, jsub )
1009  IF( mxsub.EQ.jsub .AND. isym.EQ.0 ) THEN
1010  a( mxsub-mnsub+1, mnsub ) = dconjg( ctemp )
1011  ELSE
1012  a( mxsub-mnsub+1, mnsub ) = ctemp
1013  END IF
1014  END IF
1015  240 CONTINUE
1016  250 CONTINUE
1017 *
1018  ELSE IF( ipack.EQ.6 ) THEN
1019 *
1020  DO 270 j = 1, n
1021  DO 260 i = j - kuu, j
1022  ctemp = zlatm3( m, n, i, j, isub, jsub, kl, ku, idist,
1023  $ iseed, d, igrade, dl, dr, ipvtng, iwork,
1024  $ sparse )
1025  mnsub = min( isub, jsub )
1026  mxsub = max( isub, jsub )
1027  IF( mxsub.EQ.isub .AND. isym.EQ.0 ) THEN
1028  a( mnsub-mxsub+kuu+1, mxsub ) = dconjg( ctemp )
1029  ELSE
1030  a( mnsub-mxsub+kuu+1, mxsub ) = ctemp
1031  END IF
1032  260 CONTINUE
1033  270 CONTINUE
1034 *
1035  ELSE IF( ipack.EQ.7 ) THEN
1036 *
1037  IF( isym.NE.1 ) THEN
1038  DO 290 j = 1, n
1039  DO 280 i = j - kuu, j
1040  ctemp = zlatm3( m, n, i, j, isub, jsub, kl, ku,
1041  $ idist, iseed, d, igrade, dl, dr, ipvtng,
1042  $ iwork, sparse )
1043  mnsub = min( isub, jsub )
1044  mxsub = max( isub, jsub )
1045  IF( i.LT.1 )
1046  $ a( j-i+1+kuu, i+n ) = czero
1047  IF( mxsub.EQ.isub .AND. isym.EQ.0 ) THEN
1048  a( mnsub-mxsub+kuu+1, mxsub ) = dconjg( ctemp )
1049  ELSE
1050  a( mnsub-mxsub+kuu+1, mxsub ) = ctemp
1051  END IF
1052  IF( i.GE.1 .AND. mnsub.NE.mxsub ) THEN
1053  IF( mnsub.EQ.isub .AND. isym.EQ.0 ) THEN
1054  a( mxsub-mnsub+1+kuu,
1055  $ mnsub ) = dconjg( ctemp )
1056  ELSE
1057  a( mxsub-mnsub+1+kuu, mnsub ) = ctemp
1058  END IF
1059  END IF
1060  280 CONTINUE
1061  290 CONTINUE
1062  ELSE IF( isym.EQ.1 ) THEN
1063  DO 310 j = 1, n
1064  DO 300 i = j - kuu, j + kll
1065  ctemp = zlatm3( m, n, i, j, isub, jsub, kl, ku,
1066  $ idist, iseed, d, igrade, dl, dr, ipvtng,
1067  $ iwork, sparse )
1068  a( isub-jsub+kuu+1, jsub ) = ctemp
1069  300 CONTINUE
1070  310 CONTINUE
1071  END IF
1072 *
1073  END IF
1074 *
1075  ELSE
1076 *
1077 * Use ZLATM2
1078 *
1079  IF( ipack.EQ.0 ) THEN
1080  IF( isym.EQ.0 ) THEN
1081  DO 330 j = 1, n
1082  DO 320 i = 1, j
1083  a( i, j ) = zlatm2( m, n, i, j, kl, ku, idist,
1084  $ iseed, d, igrade, dl, dr, ipvtng,
1085  $ iwork, sparse )
1086  a( j, i ) = dconjg( a( i, j ) )
1087  320 CONTINUE
1088  330 CONTINUE
1089  ELSE IF( isym.EQ.1 ) THEN
1090  DO 350 j = 1, n
1091  DO 340 i = 1, m
1092  a( i, j ) = zlatm2( m, n, i, j, kl, ku, idist,
1093  $ iseed, d, igrade, dl, dr, ipvtng,
1094  $ iwork, sparse )
1095  340 CONTINUE
1096  350 CONTINUE
1097  ELSE IF( isym.EQ.2 ) THEN
1098  DO 370 j = 1, n
1099  DO 360 i = 1, j
1100  a( i, j ) = zlatm2( m, n, i, j, kl, ku, idist,
1101  $ iseed, d, igrade, dl, dr, ipvtng,
1102  $ iwork, sparse )
1103  a( j, i ) = a( i, j )
1104  360 CONTINUE
1105  370 CONTINUE
1106  END IF
1107 *
1108  ELSE IF( ipack.EQ.1 ) THEN
1109 *
1110  DO 390 j = 1, n
1111  DO 380 i = 1, j
1112  a( i, j ) = zlatm2( m, n, i, j, kl, ku, idist, iseed,
1113  $ d, igrade, dl, dr, ipvtng, iwork, sparse )
1114  IF( i.NE.j )
1115  $ a( j, i ) = czero
1116  380 CONTINUE
1117  390 CONTINUE
1118 *
1119  ELSE IF( ipack.EQ.2 ) THEN
1120 *
1121  DO 410 j = 1, n
1122  DO 400 i = 1, j
1123  IF( isym.EQ.0 ) THEN
1124  a( j, i ) = dconjg( zlatm2( m, n, i, j, kl, ku,
1125  $ idist, iseed, d, igrade, dl, dr,
1126  $ ipvtng, iwork, sparse ) )
1127  ELSE
1128  a( j, i ) = zlatm2( m, n, i, j, kl, ku, idist,
1129  $ iseed, d, igrade, dl, dr, ipvtng,
1130  $ iwork, sparse )
1131  END IF
1132  IF( i.NE.j )
1133  $ a( i, j ) = czero
1134  400 CONTINUE
1135  410 CONTINUE
1136 *
1137  ELSE IF( ipack.EQ.3 ) THEN
1138 *
1139  isub = 0
1140  jsub = 1
1141  DO 430 j = 1, n
1142  DO 420 i = 1, j
1143  isub = isub + 1
1144  IF( isub.GT.lda ) THEN
1145  isub = 1
1146  jsub = jsub + 1
1147  END IF
1148  a( isub, jsub ) = zlatm2( m, n, i, j, kl, ku, idist,
1149  $ iseed, d, igrade, dl, dr, ipvtng,
1150  $ iwork, sparse )
1151  420 CONTINUE
1152  430 CONTINUE
1153 *
1154  ELSE IF( ipack.EQ.4 ) THEN
1155 *
1156  IF( isym.EQ.0 .OR. isym.EQ.2 ) THEN
1157  DO 450 j = 1, n
1158  DO 440 i = 1, j
1159 *
1160 * Compute K = location of (I,J) entry in packed array
1161 *
1162  IF( i.EQ.1 ) THEN
1163  k = j
1164  ELSE
1165  k = n*( n+1 ) / 2 - ( n-i+1 )*( n-i+2 ) / 2 +
1166  $ j - i + 1
1167  END IF
1168 *
1169 * Convert K to (ISUB,JSUB) location
1170 *
1171  jsub = ( k-1 ) / lda + 1
1172  isub = k - lda*( jsub-1 )
1173 *
1174  a( isub, jsub ) = zlatm2( m, n, i, j, kl, ku,
1175  $ idist, iseed, d, igrade, dl, dr,
1176  $ ipvtng, iwork, sparse )
1177  IF( isym.EQ.0 )
1178  $ a( isub, jsub ) = dconjg( a( isub, jsub ) )
1179  440 CONTINUE
1180  450 CONTINUE
1181  ELSE
1182  isub = 0
1183  jsub = 1
1184  DO 470 j = 1, n
1185  DO 460 i = j, m
1186  isub = isub + 1
1187  IF( isub.GT.lda ) THEN
1188  isub = 1
1189  jsub = jsub + 1
1190  END IF
1191  a( isub, jsub ) = zlatm2( m, n, i, j, kl, ku,
1192  $ idist, iseed, d, igrade, dl, dr,
1193  $ ipvtng, iwork, sparse )
1194  460 CONTINUE
1195  470 CONTINUE
1196  END IF
1197 *
1198  ELSE IF( ipack.EQ.5 ) THEN
1199 *
1200  DO 490 j = 1, n
1201  DO 480 i = j - kuu, j
1202  IF( i.LT.1 ) THEN
1203  a( j-i+1, i+n ) = czero
1204  ELSE
1205  IF( isym.EQ.0 ) THEN
1206  a( j-i+1, i ) = dconjg( zlatm2( m, n, i, j, kl,
1207  $ ku, idist, iseed, d, igrade, dl,
1208  $ dr, ipvtng, iwork, sparse ) )
1209  ELSE
1210  a( j-i+1, i ) = zlatm2( m, n, i, j, kl, ku,
1211  $ idist, iseed, d, igrade, dl, dr,
1212  $ ipvtng, iwork, sparse )
1213  END IF
1214  END IF
1215  480 CONTINUE
1216  490 CONTINUE
1217 *
1218  ELSE IF( ipack.EQ.6 ) THEN
1219 *
1220  DO 510 j = 1, n
1221  DO 500 i = j - kuu, j
1222  a( i-j+kuu+1, j ) = zlatm2( m, n, i, j, kl, ku, idist,
1223  $ iseed, d, igrade, dl, dr, ipvtng,
1224  $ iwork, sparse )
1225  500 CONTINUE
1226  510 CONTINUE
1227 *
1228  ELSE IF( ipack.EQ.7 ) THEN
1229 *
1230  IF( isym.NE.1 ) THEN
1231  DO 530 j = 1, n
1232  DO 520 i = j - kuu, j
1233  a( i-j+kuu+1, j ) = zlatm2( m, n, i, j, kl, ku,
1234  $ idist, iseed, d, igrade, dl,
1235  $ dr, ipvtng, iwork, sparse )
1236  IF( i.LT.1 )
1237  $ a( j-i+1+kuu, i+n ) = czero
1238  IF( i.GE.1 .AND. i.NE.j ) THEN
1239  IF( isym.EQ.0 ) THEN
1240  a( j-i+1+kuu, i ) = dconjg( a( i-j+kuu+1,
1241  $ j ) )
1242  ELSE
1243  a( j-i+1+kuu, i ) = a( i-j+kuu+1, j )
1244  END IF
1245  END IF
1246  520 CONTINUE
1247  530 CONTINUE
1248  ELSE IF( isym.EQ.1 ) THEN
1249  DO 550 j = 1, n
1250  DO 540 i = j - kuu, j + kll
1251  a( i-j+kuu+1, j ) = zlatm2( m, n, i, j, kl, ku,
1252  $ idist, iseed, d, igrade, dl,
1253  $ dr, ipvtng, iwork, sparse )
1254  540 CONTINUE
1255  550 CONTINUE
1256  END IF
1257 *
1258  END IF
1259 *
1260  END IF
1261 *
1262 * 5) Scaling the norm
1263 *
1264  IF( ipack.EQ.0 ) THEN
1265  onorm = zlange( 'M', m, n, a, lda, tempa )
1266  ELSE IF( ipack.EQ.1 ) THEN
1267  onorm = zlansy( 'M', 'U', n, a, lda, tempa )
1268  ELSE IF( ipack.EQ.2 ) THEN
1269  onorm = zlansy( 'M', 'L', n, a, lda, tempa )
1270  ELSE IF( ipack.EQ.3 ) THEN
1271  onorm = zlansp( 'M', 'U', n, a, tempa )
1272  ELSE IF( ipack.EQ.4 ) THEN
1273  onorm = zlansp( 'M', 'L', n, a, tempa )
1274  ELSE IF( ipack.EQ.5 ) THEN
1275  onorm = zlansb( 'M', 'L', n, kll, a, lda, tempa )
1276  ELSE IF( ipack.EQ.6 ) THEN
1277  onorm = zlansb( 'M', 'U', n, kuu, a, lda, tempa )
1278  ELSE IF( ipack.EQ.7 ) THEN
1279  onorm = zlangb( 'M', n, kll, kuu, a, lda, tempa )
1280  END IF
1281 *
1282  IF( anorm.GE.zero ) THEN
1283 *
1284  IF( anorm.GT.zero .AND. onorm.EQ.zero ) THEN
1285 *
1286 * Desired scaling impossible
1287 *
1288  info = 5
1289  RETURN
1290 *
1291  ELSE IF( ( anorm.GT.one .AND. onorm.LT.one ) .OR.
1292  $ ( anorm.LT.one .AND. onorm.GT.one ) ) THEN
1293 *
1294 * Scale carefully to avoid over / underflow
1295 *
1296  IF( ipack.LE.2 ) THEN
1297  DO 560 j = 1, n
1298  CALL zdscal( m, one / onorm, a( 1, j ), 1 )
1299  CALL zdscal( m, anorm, a( 1, j ), 1 )
1300  560 CONTINUE
1301 *
1302  ELSE IF( ipack.EQ.3 .OR. ipack.EQ.4 ) THEN
1303 *
1304  CALL zdscal( n*( n+1 ) / 2, one / onorm, a, 1 )
1305  CALL zdscal( n*( n+1 ) / 2, anorm, a, 1 )
1306 *
1307  ELSE IF( ipack.GE.5 ) THEN
1308 *
1309  DO 570 j = 1, n
1310  CALL zdscal( kll+kuu+1, one / onorm, a( 1, j ), 1 )
1311  CALL zdscal( kll+kuu+1, anorm, a( 1, j ), 1 )
1312  570 CONTINUE
1313 *
1314  END IF
1315 *
1316  ELSE
1317 *
1318 * Scale straightforwardly
1319 *
1320  IF( ipack.LE.2 ) THEN
1321  DO 580 j = 1, n
1322  CALL zdscal( m, anorm / onorm, a( 1, j ), 1 )
1323  580 CONTINUE
1324 *
1325  ELSE IF( ipack.EQ.3 .OR. ipack.EQ.4 ) THEN
1326 *
1327  CALL zdscal( n*( n+1 ) / 2, anorm / onorm, a, 1 )
1328 *
1329  ELSE IF( ipack.GE.5 ) THEN
1330 *
1331  DO 590 j = 1, n
1332  CALL zdscal( kll+kuu+1, anorm / onorm, a( 1, j ), 1 )
1333  590 CONTINUE
1334  END IF
1335 *
1336  END IF
1337 *
1338  END IF
1339 *
1340 * End of ZLATMR
1341 *
double precision function zlansy(NORM, UPLO, N, A, LDA, WORK)
ZLANSY returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a complex symmetric matrix.
Definition: zlansy.f:125
double precision function zlansb(NORM, UPLO, N, K, AB, LDAB, WORK)
ZLANSB returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a symmetric band matrix.
Definition: zlansb.f:132
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
complex *16 function zlatm3(M, N, I, J, ISUB, JSUB, KL, KU, IDIST, ISEED, D, IGRADE, DL, DR, IPVTNG, IWORK, SPARSE)
ZLATM3
Definition: zlatm3.f:231
complex *16 function zlatm2(M, N, I, J, KL, KU, IDIST, ISEED, D, IGRADE, DL, DR, IPVTNG, IWORK, SPARSE)
ZLATM2
Definition: zlatm2.f:213
subroutine zlatm1(MODE, COND, IRSIGN, IDIST, ISEED, D, N, INFO)
ZLATM1
Definition: zlatm1.f:139
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
double precision function zlangb(NORM, N, KL, KU, AB, LDAB, WORK)
ZLANGB returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition: zlangb.f:127
double precision function zlange(NORM, M, N, A, LDA, WORK)
ZLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition: zlange.f:117
double precision function zlansp(NORM, UPLO, N, AP, WORK)
ZLANSP returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a symmetric matrix supplied in packed form.
Definition: zlansp.f:117
subroutine zdscal(N, DA, ZX, INCX)
ZDSCAL
Definition: zdscal.f:54

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine zlatms ( integer  M,
integer  N,
character  DIST,
integer, dimension( 4 )  ISEED,
character  SYM,
double precision, dimension( * )  D,
integer  MODE,
double precision  COND,
double precision  DMAX,
integer  KL,
integer  KU,
character  PACK,
complex*16, dimension( lda, * )  A,
integer  LDA,
complex*16, dimension( * )  WORK,
integer  INFO 
)

ZLATMS

Purpose:
    ZLATMS generates random matrices with specified singular values
    (or hermitian with specified eigenvalues)
    for testing LAPACK programs.

    ZLATMS operates by applying the following sequence of
    operations:

      Set the diagonal to D, where D may be input or
         computed according to MODE, COND, DMAX, and SYM
         as described below.

      Generate a matrix with the appropriate band structure, by one
         of two methods:

      Method A:
          Generate a dense M x N matrix by multiplying D on the left
              and the right by random unitary matrices, then:

          Reduce the bandwidth according to KL and KU, using
              Householder transformations.

      Method B:
          Convert the bandwidth-0 (i.e., diagonal) matrix to a
              bandwidth-1 matrix using Givens rotations, "chasing"
              out-of-band elements back, much as in QR; then convert
              the bandwidth-1 to a bandwidth-2 matrix, etc.  Note
              that for reasonably small bandwidths (relative to M and
              N) this requires less storage, as a dense matrix is not
              generated.  Also, for hermitian or symmetric matrices,
              only one triangle is generated.

      Method A is chosen if the bandwidth is a large fraction of the
          order of the matrix, and LDA is at least M (so a dense
          matrix can be stored.)  Method B is chosen if the bandwidth
          is small (< 1/2 N for hermitian or symmetric, < .3 N+M for
          non-symmetric), or LDA is less than M and not less than the
          bandwidth.

      Pack the matrix if desired. Options specified by PACK are:
         no packing
         zero out upper half (if hermitian)
         zero out lower half (if hermitian)
         store the upper half columnwise (if hermitian or upper
               triangular)
         store the lower half columnwise (if hermitian or lower
               triangular)
         store the lower triangle in banded format (if hermitian or
               lower triangular)
         store the upper triangle in banded format (if hermitian or
               upper triangular)
         store the entire matrix in banded format
      If Method B is chosen, and band format is specified, then the
         matrix will be generated in the band format, so no repacking
         will be necessary.
Parameters
[in]M
          M is INTEGER
           The number of rows of A. Not modified.
[in]N
          N is INTEGER
           The number of columns of A. N must equal M if the matrix
           is symmetric or hermitian (i.e., if SYM is not 'N')
           Not modified.
[in]DIST
          DIST is CHARACTER*1
           On entry, DIST specifies the type of distribution to be used
           to generate the random eigen-/singular values.
           'U' => UNIFORM( 0, 1 )  ( 'U' for uniform )
           'S' => UNIFORM( -1, 1 ) ( 'S' for symmetric )
           'N' => NORMAL( 0, 1 )   ( 'N' for normal )
           Not modified.
[in,out]ISEED
          ISEED is INTEGER array, dimension ( 4 )
           On entry ISEED specifies the seed of the random number
           generator. They should lie between 0 and 4095 inclusive,
           and ISEED(4) should be odd. The random number generator
           uses a linear congruential sequence limited to small
           integers, and so should produce machine independent
           random numbers. The values of ISEED are changed on
           exit, and can be used in the next call to ZLATMS
           to continue the same random number sequence.
           Changed on exit.
[in]SYM
          SYM is CHARACTER*1
           If SYM='H', the generated matrix is hermitian, with
             eigenvalues specified by D, COND, MODE, and DMAX; they
             may be positive, negative, or zero.
           If SYM='P', the generated matrix is hermitian, with
             eigenvalues (= singular values) specified by D, COND,
             MODE, and DMAX; they will not be negative.
           If SYM='N', the generated matrix is nonsymmetric, with
             singular values specified by D, COND, MODE, and DMAX;
             they will not be negative.
           If SYM='S', the generated matrix is (complex) symmetric,
             with singular values specified by D, COND, MODE, and
             DMAX; they will not be negative.
           Not modified.
[in,out]D
          D is DOUBLE PRECISION array, dimension ( MIN( M, N ) )
           This array is used to specify the singular values or
           eigenvalues of A (see SYM, above.)  If MODE=0, then D is
           assumed to contain the singular/eigenvalues, otherwise
           they will be computed according to MODE, COND, and DMAX,
           and placed in D.
           Modified if MODE is nonzero.
[in]MODE
          MODE is INTEGER
           On entry this describes how the singular/eigenvalues are to
           be specified:
           MODE = 0 means use D as input
           MODE = 1 sets D(1)=1 and D(2:N)=1.0/COND
           MODE = 2 sets D(1:N-1)=1 and D(N)=1.0/COND
           MODE = 3 sets D(I)=COND**(-(I-1)/(N-1))
           MODE = 4 sets D(i)=1 - (i-1)/(N-1)*(1 - 1/COND)
           MODE = 5 sets D to random numbers in the range
                    ( 1/COND , 1 ) such that their logarithms
                    are uniformly distributed.
           MODE = 6 set D to random numbers from same distribution
                    as the rest of the matrix.
           MODE < 0 has the same meaning as ABS(MODE), except that
              the order of the elements of D is reversed.
           Thus if MODE is positive, D has entries ranging from
              1 to 1/COND, if negative, from 1/COND to 1,
           If SYM='H', and MODE is neither 0, 6, nor -6, then
              the elements of D will also be multiplied by a random
              sign (i.e., +1 or -1.)
           Not modified.
[in]COND
          COND is DOUBLE PRECISION
           On entry, this is used as described under MODE above.
           If used, it must be >= 1. Not modified.
[in]DMAX
          DMAX is DOUBLE PRECISION
           If MODE is neither -6, 0 nor 6, the contents of D, as
           computed according to MODE and COND, will be scaled by
           DMAX / max(abs(D(i))); thus, the maximum absolute eigen- or
           singular value (which is to say the norm) will be abs(DMAX).
           Note that DMAX need not be positive: if DMAX is negative
           (or zero), D will be scaled by a negative number (or zero).
           Not modified.
[in]KL
          KL is INTEGER
           This specifies the lower bandwidth of the  matrix. For
           example, KL=0 implies upper triangular, KL=1 implies upper
           Hessenberg, and KL being at least M-1 means that the matrix
           has full lower bandwidth.  KL must equal KU if the matrix
           is symmetric or hermitian.
           Not modified.
[in]KU
          KU is INTEGER
           This specifies the upper bandwidth of the  matrix. For
           example, KU=0 implies lower triangular, KU=1 implies lower
           Hessenberg, and KU being at least N-1 means that the matrix
           has full upper bandwidth.  KL must equal KU if the matrix
           is symmetric or hermitian.
           Not modified.
[in]PACK
          PACK is CHARACTER*1
           This specifies packing of matrix as follows:
           'N' => no packing
           'U' => zero out all subdiagonal entries (if symmetric
                  or hermitian)
           'L' => zero out all superdiagonal entries (if symmetric
                  or hermitian)
           'C' => store the upper triangle columnwise (only if the
                  matrix is symmetric, hermitian, or upper triangular)
           'R' => store the lower triangle columnwise (only if the
                  matrix is symmetric, hermitian, or lower triangular)
           'B' => store the lower triangle in band storage scheme
                  (only if the matrix is symmetric, hermitian, or
                  lower triangular)
           'Q' => store the upper triangle in band storage scheme
                  (only if the matrix is symmetric, hermitian, or
                  upper triangular)
           'Z' => store the entire matrix in band storage scheme
                      (pivoting can be provided for by using this
                      option to store A in the trailing rows of
                      the allocated storage)

           Using these options, the various LAPACK packed and banded
           storage schemes can be obtained:
           GB                    - use 'Z'
           PB, SB, HB, or TB     - use 'B' or 'Q'
           PP, SP, HB, or TP     - use 'C' or 'R'

           If two calls to ZLATMS differ only in the PACK parameter,
           they will generate mathematically equivalent matrices.
           Not modified.
[in,out]A
          A is COMPLEX*16 array, dimension ( LDA, N )
           On exit A is the desired test matrix.  A is first generated
           in full (unpacked) form, and then packed, if so specified
           by PACK.  Thus, the first M elements of the first N
           columns will always be modified.  If PACK specifies a
           packed or banded storage scheme, all LDA elements of the
           first N columns will be modified; the elements of the
           array which do not correspond to elements of the generated
           matrix are set to zero.
           Modified.
[in]LDA
          LDA is INTEGER
           LDA specifies the first dimension of A as declared in the
           calling program.  If PACK='N', 'U', 'L', 'C', or 'R', then
           LDA must be at least M.  If PACK='B' or 'Q', then LDA must
           be at least MIN( KL, M-1) (which is equal to MIN(KU,N-1)).
           If PACK='Z', LDA must be large enough to hold the packed
           array: MIN( KU, N-1) + MIN( KL, M-1) + 1.
           Not modified.
[out]WORK
          WORK is COMPLEX*16 array, dimension ( 3*MAX( N, M ) )
           Workspace.
           Modified.
[out]INFO
          INFO is INTEGER
           Error code.  On exit, INFO will be set to one of the
           following values:
             0 => normal return
            -1 => M negative or unequal to N and SYM='S', 'H', or 'P'
            -2 => N negative
            -3 => DIST illegal string
            -5 => SYM illegal string
            -7 => MODE not in range -6 to 6
            -8 => COND less than 1.0, and MODE neither -6, 0 nor 6
           -10 => KL negative
           -11 => KU negative, or SYM is not 'N' and KU is not equal to
                  KL
           -12 => PACK illegal string, or PACK='U' or 'L', and SYM='N';
                  or PACK='C' or 'Q' and SYM='N' and KL is not zero;
                  or PACK='R' or 'B' and SYM='N' and KU is not zero;
                  or PACK='U', 'L', 'C', 'R', 'B', or 'Q', and M is not
                  N.
           -14 => LDA is less than M, or PACK='Z' and LDA is less than
                  MIN(KU,N-1) + MIN(KL,M-1) + 1.
            1  => Error return from DLATM1
            2  => Cannot scale to DMAX (max. sing. value is 0)
            3  => Error return from ZLAGGE, CLAGHE or CLAGSY
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 334 of file zlatms.f.

334 *
335 * -- LAPACK computational routine (version 3.4.0) --
336 * -- LAPACK is a software package provided by Univ. of Tennessee, --
337 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
338 * November 2011
339 *
340 * .. Scalar Arguments ..
341  CHARACTER dist, pack, sym
342  INTEGER info, kl, ku, lda, m, mode, n
343  DOUBLE PRECISION cond, dmax
344 * ..
345 * .. Array Arguments ..
346  INTEGER iseed( 4 )
347  DOUBLE PRECISION d( * )
348  COMPLEX*16 a( lda, * ), work( * )
349 * ..
350 *
351 * =====================================================================
352 *
353 * .. Parameters ..
354  DOUBLE PRECISION zero
355  parameter( zero = 0.0d+0 )
356  DOUBLE PRECISION one
357  parameter( one = 1.0d+0 )
358  COMPLEX*16 czero
359  parameter( czero = ( 0.0d+0, 0.0d+0 ) )
360  DOUBLE PRECISION twopi
361  parameter( twopi = 6.2831853071795864769252867663d+0 )
362 * ..
363 * .. Local Scalars ..
364  LOGICAL givens, ilextr, iltemp, topdwn, zsym
365  INTEGER i, ic, icol, idist, iendch, iinfo, il, ilda,
366  $ ioffg, ioffst, ipack, ipackg, ir, ir1, ir2,
367  $ irow, irsign, iskew, isym, isympk, j, jc, jch,
368  $ jkl, jku, jr, k, llb, minlda, mnmin, mr, nc,
369  $ uub
370  DOUBLE PRECISION alpha, angle, realc, temp
371  COMPLEX*16 c, ct, ctemp, dummy, extra, s, st
372 * ..
373 * .. External Functions ..
374  LOGICAL lsame
375  DOUBLE PRECISION dlarnd
376  COMPLEX*16 zlarnd
377  EXTERNAL lsame, dlarnd, zlarnd
378 * ..
379 * .. External Subroutines ..
380  EXTERNAL dlatm1, dscal, xerbla, zlagge, zlaghe, zlagsy,
381  $ zlarot, zlartg, zlaset
382 * ..
383 * .. Intrinsic Functions ..
384  INTRINSIC abs, cos, dble, dcmplx, dconjg, max, min, mod,
385  $ sin
386 * ..
387 * .. Executable Statements ..
388 *
389 * 1) Decode and Test the input parameters.
390 * Initialize flags & seed.
391 *
392  info = 0
393 *
394 * Quick return if possible
395 *
396  IF( m.EQ.0 .OR. n.EQ.0 )
397  $ RETURN
398 *
399 * Decode DIST
400 *
401  IF( lsame( dist, 'U' ) ) THEN
402  idist = 1
403  ELSE IF( lsame( dist, 'S' ) ) THEN
404  idist = 2
405  ELSE IF( lsame( dist, 'N' ) ) THEN
406  idist = 3
407  ELSE
408  idist = -1
409  END IF
410 *
411 * Decode SYM
412 *
413  IF( lsame( sym, 'N' ) ) THEN
414  isym = 1
415  irsign = 0
416  zsym = .false.
417  ELSE IF( lsame( sym, 'P' ) ) THEN
418  isym = 2
419  irsign = 0
420  zsym = .false.
421  ELSE IF( lsame( sym, 'S' ) ) THEN
422  isym = 2
423  irsign = 0
424  zsym = .true.
425  ELSE IF( lsame( sym, 'H' ) ) THEN
426  isym = 2
427  irsign = 1
428  zsym = .false.
429  ELSE
430  isym = -1
431  END IF
432 *
433 * Decode PACK
434 *
435  isympk = 0
436  IF( lsame( pack, 'N' ) ) THEN
437  ipack = 0
438  ELSE IF( lsame( pack, 'U' ) ) THEN
439  ipack = 1
440  isympk = 1
441  ELSE IF( lsame( pack, 'L' ) ) THEN
442  ipack = 2
443  isympk = 1
444  ELSE IF( lsame( pack, 'C' ) ) THEN
445  ipack = 3
446  isympk = 2
447  ELSE IF( lsame( pack, 'R' ) ) THEN
448  ipack = 4
449  isympk = 3
450  ELSE IF( lsame( pack, 'B' ) ) THEN
451  ipack = 5
452  isympk = 3
453  ELSE IF( lsame( pack, 'Q' ) ) THEN
454  ipack = 6
455  isympk = 2
456  ELSE IF( lsame( pack, 'Z' ) ) THEN
457  ipack = 7
458  ELSE
459  ipack = -1
460  END IF
461 *
462 * Set certain internal parameters
463 *
464  mnmin = min( m, n )
465  llb = min( kl, m-1 )
466  uub = min( ku, n-1 )
467  mr = min( m, n+llb )
468  nc = min( n, m+uub )
469 *
470  IF( ipack.EQ.5 .OR. ipack.EQ.6 ) THEN
471  minlda = uub + 1
472  ELSE IF( ipack.EQ.7 ) THEN
473  minlda = llb + uub + 1
474  ELSE
475  minlda = m
476  END IF
477 *
478 * Use Givens rotation method if bandwidth small enough,
479 * or if LDA is too small to store the matrix unpacked.
480 *
481  givens = .false.
482  IF( isym.EQ.1 ) THEN
483  IF( dble( llb+uub ).LT.0.3d0*dble( max( 1, mr+nc ) ) )
484  $ givens = .true.
485  ELSE
486  IF( 2*llb.LT.m )
487  $ givens = .true.
488  END IF
489  IF( lda.LT.m .AND. lda.GE.minlda )
490  $ givens = .true.
491 *
492 * Set INFO if an error
493 *
494  IF( m.LT.0 ) THEN
495  info = -1
496  ELSE IF( m.NE.n .AND. isym.NE.1 ) THEN
497  info = -1
498  ELSE IF( n.LT.0 ) THEN
499  info = -2
500  ELSE IF( idist.EQ.-1 ) THEN
501  info = -3
502  ELSE IF( isym.EQ.-1 ) THEN
503  info = -5
504  ELSE IF( abs( mode ).GT.6 ) THEN
505  info = -7
506  ELSE IF( ( mode.NE.0 .AND. abs( mode ).NE.6 ) .AND. cond.LT.one )
507  $ THEN
508  info = -8
509  ELSE IF( kl.LT.0 ) THEN
510  info = -10
511  ELSE IF( ku.LT.0 .OR. ( isym.NE.1 .AND. kl.NE.ku ) ) THEN
512  info = -11
513  ELSE IF( ipack.EQ.-1 .OR. ( isympk.EQ.1 .AND. isym.EQ.1 ) .OR.
514  $ ( isympk.EQ.2 .AND. isym.EQ.1 .AND. kl.GT.0 ) .OR.
515  $ ( isympk.EQ.3 .AND. isym.EQ.1 .AND. ku.GT.0 ) .OR.
516  $ ( isympk.NE.0 .AND. m.NE.n ) ) THEN
517  info = -12
518  ELSE IF( lda.LT.max( 1, minlda ) ) THEN
519  info = -14
520  END IF
521 *
522  IF( info.NE.0 ) THEN
523  CALL xerbla( 'ZLATMS', -info )
524  RETURN
525  END IF
526 *
527 * Initialize random number generator
528 *
529  DO 10 i = 1, 4
530  iseed( i ) = mod( abs( iseed( i ) ), 4096 )
531  10 CONTINUE
532 *
533  IF( mod( iseed( 4 ), 2 ).NE.1 )
534  $ iseed( 4 ) = iseed( 4 ) + 1
535 *
536 * 2) Set up D if indicated.
537 *
538 * Compute D according to COND and MODE
539 *
540  CALL dlatm1( mode, cond, irsign, idist, iseed, d, mnmin, iinfo )
541  IF( iinfo.NE.0 ) THEN
542  info = 1
543  RETURN
544  END IF
545 *
546 * Choose Top-Down if D is (apparently) increasing,
547 * Bottom-Up if D is (apparently) decreasing.
548 *
549  IF( abs( d( 1 ) ).LE.abs( d( mnmin ) ) ) THEN
550  topdwn = .true.
551  ELSE
552  topdwn = .false.
553  END IF
554 *
555  IF( mode.NE.0 .AND. abs( mode ).NE.6 ) THEN
556 *
557 * Scale by DMAX
558 *
559  temp = abs( d( 1 ) )
560  DO 20 i = 2, mnmin
561  temp = max( temp, abs( d( i ) ) )
562  20 CONTINUE
563 *
564  IF( temp.GT.zero ) THEN
565  alpha = dmax / temp
566  ELSE
567  info = 2
568  RETURN
569  END IF
570 *
571  CALL dscal( mnmin, alpha, d, 1 )
572 *
573  END IF
574 *
575  CALL zlaset( 'Full', lda, n, czero, czero, a, lda )
576 *
577 * 3) Generate Banded Matrix using Givens rotations.
578 * Also the special case of UUB=LLB=0
579 *
580 * Compute Addressing constants to cover all
581 * storage formats. Whether GE, HE, SY, GB, HB, or SB,
582 * upper or lower triangle or both,
583 * the (i,j)-th element is in
584 * A( i - ISKEW*j + IOFFST, j )
585 *
586  IF( ipack.GT.4 ) THEN
587  ilda = lda - 1
588  iskew = 1
589  IF( ipack.GT.5 ) THEN
590  ioffst = uub + 1
591  ELSE
592  ioffst = 1
593  END IF
594  ELSE
595  ilda = lda
596  iskew = 0
597  ioffst = 0
598  END IF
599 *
600 * IPACKG is the format that the matrix is generated in. If this is
601 * different from IPACK, then the matrix must be repacked at the
602 * end. It also signals how to compute the norm, for scaling.
603 *
604  ipackg = 0
605 *
606 * Diagonal Matrix -- We are done, unless it
607 * is to be stored HP/SP/PP/TP (PACK='R' or 'C')
608 *
609  IF( llb.EQ.0 .AND. uub.EQ.0 ) THEN
610  DO 30 j = 1, mnmin
611  a( ( 1-iskew )*j+ioffst, j ) = dcmplx( d( j ) )
612  30 CONTINUE
613 *
614  IF( ipack.LE.2 .OR. ipack.GE.5 )
615  $ ipackg = ipack
616 *
617  ELSE IF( givens ) THEN
618 *
619 * Check whether to use Givens rotations,
620 * Householder transformations, or nothing.
621 *
622  IF( isym.EQ.1 ) THEN
623 *
624 * Non-symmetric -- A = U D V
625 *
626  IF( ipack.GT.4 ) THEN
627  ipackg = ipack
628  ELSE
629  ipackg = 0
630  END IF
631 *
632  DO 40 j = 1, mnmin
633  a( ( 1-iskew )*j+ioffst, j ) = dcmplx( d( j ) )
634  40 CONTINUE
635 *
636  IF( topdwn ) THEN
637  jkl = 0
638  DO 70 jku = 1, uub
639 *
640 * Transform from bandwidth JKL, JKU-1 to JKL, JKU
641 *
642 * Last row actually rotated is M
643 * Last column actually rotated is MIN( M+JKU, N )
644 *
645  DO 60 jr = 1, min( m+jku, n ) + jkl - 1
646  extra = czero
647  angle = twopi*dlarnd( 1, iseed )
648  c = cos( angle )*zlarnd( 5, iseed )
649  s = sin( angle )*zlarnd( 5, iseed )
650  icol = max( 1, jr-jkl )
651  IF( jr.LT.m ) THEN
652  il = min( n, jr+jku ) + 1 - icol
653  CALL zlarot( .true., jr.GT.jkl, .false., il, c,
654  $ s, a( jr-iskew*icol+ioffst, icol ),
655  $ ilda, extra, dummy )
656  END IF
657 *
658 * Chase "EXTRA" back up
659 *
660  ir = jr
661  ic = icol
662  DO 50 jch = jr - jkl, 1, -jkl - jku
663  IF( ir.LT.m ) THEN
664  CALL zlartg( a( ir+1-iskew*( ic+1 )+ioffst,
665  $ ic+1 ), extra, realc, s, dummy )
666  dummy = zlarnd( 5, iseed )
667  c = dconjg( realc*dummy )
668  s = dconjg( -s*dummy )
669  END IF
670  irow = max( 1, jch-jku )
671  il = ir + 2 - irow
672  ctemp = czero
673  iltemp = jch.GT.jku
674  CALL zlarot( .false., iltemp, .true., il, c, s,
675  $ a( irow-iskew*ic+ioffst, ic ),
676  $ ilda, ctemp, extra )
677  IF( iltemp ) THEN
678  CALL zlartg( a( irow+1-iskew*( ic+1 )+ioffst,
679  $ ic+1 ), ctemp, realc, s, dummy )
680  dummy = zlarnd( 5, iseed )
681  c = dconjg( realc*dummy )
682  s = dconjg( -s*dummy )
683 *
684  icol = max( 1, jch-jku-jkl )
685  il = ic + 2 - icol
686  extra = czero
687  CALL zlarot( .true., jch.GT.jku+jkl, .true.,
688  $ il, c, s, a( irow-iskew*icol+
689  $ ioffst, icol ), ilda, extra,
690  $ ctemp )
691  ic = icol
692  ir = irow
693  END IF
694  50 CONTINUE
695  60 CONTINUE
696  70 CONTINUE
697 *
698  jku = uub
699  DO 100 jkl = 1, llb
700 *
701 * Transform from bandwidth JKL-1, JKU to JKL, JKU
702 *
703  DO 90 jc = 1, min( n+jkl, m ) + jku - 1
704  extra = czero
705  angle = twopi*dlarnd( 1, iseed )
706  c = cos( angle )*zlarnd( 5, iseed )
707  s = sin( angle )*zlarnd( 5, iseed )
708  irow = max( 1, jc-jku )
709  IF( jc.LT.n ) THEN
710  il = min( m, jc+jkl ) + 1 - irow
711  CALL zlarot( .false., jc.GT.jku, .false., il, c,
712  $ s, a( irow-iskew*jc+ioffst, jc ),
713  $ ilda, extra, dummy )
714  END IF
715 *
716 * Chase "EXTRA" back up
717 *
718  ic = jc
719  ir = irow
720  DO 80 jch = jc - jku, 1, -jkl - jku
721  IF( ic.LT.n ) THEN
722  CALL zlartg( a( ir+1-iskew*( ic+1 )+ioffst,
723  $ ic+1 ), extra, realc, s, dummy )
724  dummy = zlarnd( 5, iseed )
725  c = dconjg( realc*dummy )
726  s = dconjg( -s*dummy )
727  END IF
728  icol = max( 1, jch-jkl )
729  il = ic + 2 - icol
730  ctemp = czero
731  iltemp = jch.GT.jkl
732  CALL zlarot( .true., iltemp, .true., il, c, s,
733  $ a( ir-iskew*icol+ioffst, icol ),
734  $ ilda, ctemp, extra )
735  IF( iltemp ) THEN
736  CALL zlartg( a( ir+1-iskew*( icol+1 )+ioffst,
737  $ icol+1 ), ctemp, realc, s,
738  $ dummy )
739  dummy = zlarnd( 5, iseed )
740  c = dconjg( realc*dummy )
741  s = dconjg( -s*dummy )
742  irow = max( 1, jch-jkl-jku )
743  il = ir + 2 - irow
744  extra = czero
745  CALL zlarot( .false., jch.GT.jkl+jku, .true.,
746  $ il, c, s, a( irow-iskew*icol+
747  $ ioffst, icol ), ilda, extra,
748  $ ctemp )
749  ic = icol
750  ir = irow
751  END IF
752  80 CONTINUE
753  90 CONTINUE
754  100 CONTINUE
755 *
756  ELSE
757 *
758 * Bottom-Up -- Start at the bottom right.
759 *
760  jkl = 0
761  DO 130 jku = 1, uub
762 *
763 * Transform from bandwidth JKL, JKU-1 to JKL, JKU
764 *
765 * First row actually rotated is M
766 * First column actually rotated is MIN( M+JKU, N )
767 *
768  iendch = min( m, n+jkl ) - 1
769  DO 120 jc = min( m+jku, n ) - 1, 1 - jkl, -1
770  extra = czero
771  angle = twopi*dlarnd( 1, iseed )
772  c = cos( angle )*zlarnd( 5, iseed )
773  s = sin( angle )*zlarnd( 5, iseed )
774  irow = max( 1, jc-jku+1 )
775  IF( jc.GT.0 ) THEN
776  il = min( m, jc+jkl+1 ) + 1 - irow
777  CALL zlarot( .false., .false., jc+jkl.LT.m, il,
778  $ c, s, a( irow-iskew*jc+ioffst,
779  $ jc ), ilda, dummy, extra )
780  END IF
781 *
782 * Chase "EXTRA" back down
783 *
784  ic = jc
785  DO 110 jch = jc + jkl, iendch, jkl + jku
786  ilextr = ic.GT.0
787  IF( ilextr ) THEN
788  CALL zlartg( a( jch-iskew*ic+ioffst, ic ),
789  $ extra, realc, s, dummy )
790  dummy = zlarnd( 5, iseed )
791  c = realc*dummy
792  s = s*dummy
793  END IF
794  ic = max( 1, ic )
795  icol = min( n-1, jch+jku )
796  iltemp = jch + jku.LT.n
797  ctemp = czero
798  CALL zlarot( .true., ilextr, iltemp, icol+2-ic,
799  $ c, s, a( jch-iskew*ic+ioffst, ic ),
800  $ ilda, extra, ctemp )
801  IF( iltemp ) THEN
802  CALL zlartg( a( jch-iskew*icol+ioffst,
803  $ icol ), ctemp, realc, s, dummy )
804  dummy = zlarnd( 5, iseed )
805  c = realc*dummy
806  s = s*dummy
807  il = min( iendch, jch+jkl+jku ) + 2 - jch
808  extra = czero
809  CALL zlarot( .false., .true.,
810  $ jch+jkl+jku.LE.iendch, il, c, s,
811  $ a( jch-iskew*icol+ioffst,
812  $ icol ), ilda, ctemp, extra )
813  ic = icol
814  END IF
815  110 CONTINUE
816  120 CONTINUE
817  130 CONTINUE
818 *
819  jku = uub
820  DO 160 jkl = 1, llb
821 *
822 * Transform from bandwidth JKL-1, JKU to JKL, JKU
823 *
824 * First row actually rotated is MIN( N+JKL, M )
825 * First column actually rotated is N
826 *
827  iendch = min( n, m+jku ) - 1
828  DO 150 jr = min( n+jkl, m ) - 1, 1 - jku, -1
829  extra = czero
830  angle = twopi*dlarnd( 1, iseed )
831  c = cos( angle )*zlarnd( 5, iseed )
832  s = sin( angle )*zlarnd( 5, iseed )
833  icol = max( 1, jr-jkl+1 )
834  IF( jr.GT.0 ) THEN
835  il = min( n, jr+jku+1 ) + 1 - icol
836  CALL zlarot( .true., .false., jr+jku.LT.n, il,
837  $ c, s, a( jr-iskew*icol+ioffst,
838  $ icol ), ilda, dummy, extra )
839  END IF
840 *
841 * Chase "EXTRA" back down
842 *
843  ir = jr
844  DO 140 jch = jr + jku, iendch, jkl + jku
845  ilextr = ir.GT.0
846  IF( ilextr ) THEN
847  CALL zlartg( a( ir-iskew*jch+ioffst, jch ),
848  $ extra, realc, s, dummy )
849  dummy = zlarnd( 5, iseed )
850  c = realc*dummy
851  s = s*dummy
852  END IF
853  ir = max( 1, ir )
854  irow = min( m-1, jch+jkl )
855  iltemp = jch + jkl.LT.m
856  ctemp = czero
857  CALL zlarot( .false., ilextr, iltemp, irow+2-ir,
858  $ c, s, a( ir-iskew*jch+ioffst,
859  $ jch ), ilda, extra, ctemp )
860  IF( iltemp ) THEN
861  CALL zlartg( a( irow-iskew*jch+ioffst, jch ),
862  $ ctemp, realc, s, dummy )
863  dummy = zlarnd( 5, iseed )
864  c = realc*dummy
865  s = s*dummy
866  il = min( iendch, jch+jkl+jku ) + 2 - jch
867  extra = czero
868  CALL zlarot( .true., .true.,
869  $ jch+jkl+jku.LE.iendch, il, c, s,
870  $ a( irow-iskew*jch+ioffst, jch ),
871  $ ilda, ctemp, extra )
872  ir = irow
873  END IF
874  140 CONTINUE
875  150 CONTINUE
876  160 CONTINUE
877 *
878  END IF
879 *
880  ELSE
881 *
882 * Symmetric -- A = U D U'
883 * Hermitian -- A = U D U*
884 *
885  ipackg = ipack
886  ioffg = ioffst
887 *
888  IF( topdwn ) THEN
889 *
890 * Top-Down -- Generate Upper triangle only
891 *
892  IF( ipack.GE.5 ) THEN
893  ipackg = 6
894  ioffg = uub + 1
895  ELSE
896  ipackg = 1
897  END IF
898 *
899  DO 170 j = 1, mnmin
900  a( ( 1-iskew )*j+ioffg, j ) = dcmplx( d( j ) )
901  170 CONTINUE
902 *
903  DO 200 k = 1, uub
904  DO 190 jc = 1, n - 1
905  irow = max( 1, jc-k )
906  il = min( jc+1, k+2 )
907  extra = czero
908  ctemp = a( jc-iskew*( jc+1 )+ioffg, jc+1 )
909  angle = twopi*dlarnd( 1, iseed )
910  c = cos( angle )*zlarnd( 5, iseed )
911  s = sin( angle )*zlarnd( 5, iseed )
912  IF( zsym ) THEN
913  ct = c
914  st = s
915  ELSE
916  ctemp = dconjg( ctemp )
917  ct = dconjg( c )
918  st = dconjg( s )
919  END IF
920  CALL zlarot( .false., jc.GT.k, .true., il, c, s,
921  $ a( irow-iskew*jc+ioffg, jc ), ilda,
922  $ extra, ctemp )
923  CALL zlarot( .true., .true., .false.,
924  $ min( k, n-jc )+1, ct, st,
925  $ a( ( 1-iskew )*jc+ioffg, jc ), ilda,
926  $ ctemp, dummy )
927 *
928 * Chase EXTRA back up the matrix
929 *
930  icol = jc
931  DO 180 jch = jc - k, 1, -k
932  CALL zlartg( a( jch+1-iskew*( icol+1 )+ioffg,
933  $ icol+1 ), extra, realc, s, dummy )
934  dummy = zlarnd( 5, iseed )
935  c = dconjg( realc*dummy )
936  s = dconjg( -s*dummy )
937  ctemp = a( jch-iskew*( jch+1 )+ioffg, jch+1 )
938  IF( zsym ) THEN
939  ct = c
940  st = s
941  ELSE
942  ctemp = dconjg( ctemp )
943  ct = dconjg( c )
944  st = dconjg( s )
945  END IF
946  CALL zlarot( .true., .true., .true., k+2, c, s,
947  $ a( ( 1-iskew )*jch+ioffg, jch ),
948  $ ilda, ctemp, extra )
949  irow = max( 1, jch-k )
950  il = min( jch+1, k+2 )
951  extra = czero
952  CALL zlarot( .false., jch.GT.k, .true., il, ct,
953  $ st, a( irow-iskew*jch+ioffg, jch ),
954  $ ilda, extra, ctemp )
955  icol = jch
956  180 CONTINUE
957  190 CONTINUE
958  200 CONTINUE
959 *
960 * If we need lower triangle, copy from upper. Note that
961 * the order of copying is chosen to work for 'q' -> 'b'
962 *
963  IF( ipack.NE.ipackg .AND. ipack.NE.3 ) THEN
964  DO 230 jc = 1, n
965  irow = ioffst - iskew*jc
966  IF( zsym ) THEN
967  DO 210 jr = jc, min( n, jc+uub )
968  a( jr+irow, jc ) = a( jc-iskew*jr+ioffg, jr )
969  210 CONTINUE
970  ELSE
971  DO 220 jr = jc, min( n, jc+uub )
972  a( jr+irow, jc ) = dconjg( a( jc-iskew*jr+
973  $ ioffg, jr ) )
974  220 CONTINUE
975  END IF
976  230 CONTINUE
977  IF( ipack.EQ.5 ) THEN
978  DO 250 jc = n - uub + 1, n
979  DO 240 jr = n + 2 - jc, uub + 1
980  a( jr, jc ) = czero
981  240 CONTINUE
982  250 CONTINUE
983  END IF
984  IF( ipackg.EQ.6 ) THEN
985  ipackg = ipack
986  ELSE
987  ipackg = 0
988  END IF
989  END IF
990  ELSE
991 *
992 * Bottom-Up -- Generate Lower triangle only
993 *
994  IF( ipack.GE.5 ) THEN
995  ipackg = 5
996  IF( ipack.EQ.6 )
997  $ ioffg = 1
998  ELSE
999  ipackg = 2
1000  END IF
1001 *
1002  DO 260 j = 1, mnmin
1003  a( ( 1-iskew )*j+ioffg, j ) = dcmplx( d( j ) )
1004  260 CONTINUE
1005 *
1006  DO 290 k = 1, uub
1007  DO 280 jc = n - 1, 1, -1
1008  il = min( n+1-jc, k+2 )
1009  extra = czero
1010  ctemp = a( 1+( 1-iskew )*jc+ioffg, jc )
1011  angle = twopi*dlarnd( 1, iseed )
1012  c = cos( angle )*zlarnd( 5, iseed )
1013  s = sin( angle )*zlarnd( 5, iseed )
1014  IF( zsym ) THEN
1015  ct = c
1016  st = s
1017  ELSE
1018  ctemp = dconjg( ctemp )
1019  ct = dconjg( c )
1020  st = dconjg( s )
1021  END IF
1022  CALL zlarot( .false., .true., n-jc.GT.k, il, c, s,
1023  $ a( ( 1-iskew )*jc+ioffg, jc ), ilda,
1024  $ ctemp, extra )
1025  icol = max( 1, jc-k+1 )
1026  CALL zlarot( .true., .false., .true., jc+2-icol,
1027  $ ct, st, a( jc-iskew*icol+ioffg,
1028  $ icol ), ilda, dummy, ctemp )
1029 *
1030 * Chase EXTRA back down the matrix
1031 *
1032  icol = jc
1033  DO 270 jch = jc + k, n - 1, k
1034  CALL zlartg( a( jch-iskew*icol+ioffg, icol ),
1035  $ extra, realc, s, dummy )
1036  dummy = zlarnd( 5, iseed )
1037  c = realc*dummy
1038  s = s*dummy
1039  ctemp = a( 1+( 1-iskew )*jch+ioffg, jch )
1040  IF( zsym ) THEN
1041  ct = c
1042  st = s
1043  ELSE
1044  ctemp = dconjg( ctemp )
1045  ct = dconjg( c )
1046  st = dconjg( s )
1047  END IF
1048  CALL zlarot( .true., .true., .true., k+2, c, s,
1049  $ a( jch-iskew*icol+ioffg, icol ),
1050  $ ilda, extra, ctemp )
1051  il = min( n+1-jch, k+2 )
1052  extra = czero
1053  CALL zlarot( .false., .true., n-jch.GT.k, il,
1054  $ ct, st, a( ( 1-iskew )*jch+ioffg,
1055  $ jch ), ilda, ctemp, extra )
1056  icol = jch
1057  270 CONTINUE
1058  280 CONTINUE
1059  290 CONTINUE
1060 *
1061 * If we need upper triangle, copy from lower. Note that
1062 * the order of copying is chosen to work for 'b' -> 'q'
1063 *
1064  IF( ipack.NE.ipackg .AND. ipack.NE.4 ) THEN
1065  DO 320 jc = n, 1, -1
1066  irow = ioffst - iskew*jc
1067  IF( zsym ) THEN
1068  DO 300 jr = jc, max( 1, jc-uub ), -1
1069  a( jr+irow, jc ) = a( jc-iskew*jr+ioffg, jr )
1070  300 CONTINUE
1071  ELSE
1072  DO 310 jr = jc, max( 1, jc-uub ), -1
1073  a( jr+irow, jc ) = dconjg( a( jc-iskew*jr+
1074  $ ioffg, jr ) )
1075  310 CONTINUE
1076  END IF
1077  320 CONTINUE
1078  IF( ipack.EQ.6 ) THEN
1079  DO 340 jc = 1, uub
1080  DO 330 jr = 1, uub + 1 - jc
1081  a( jr, jc ) = czero
1082  330 CONTINUE
1083  340 CONTINUE
1084  END IF
1085  IF( ipackg.EQ.5 ) THEN
1086  ipackg = ipack
1087  ELSE
1088  ipackg = 0
1089  END IF
1090  END IF
1091  END IF
1092 *
1093 * Ensure that the diagonal is real if Hermitian
1094 *
1095  IF( .NOT.zsym ) THEN
1096  DO 350 jc = 1, n
1097  irow = ioffst + ( 1-iskew )*jc
1098  a( irow, jc ) = dcmplx( dble( a( irow, jc ) ) )
1099  350 CONTINUE
1100  END IF
1101 *
1102  END IF
1103 *
1104  ELSE
1105 *
1106 * 4) Generate Banded Matrix by first
1107 * Rotating by random Unitary matrices,
1108 * then reducing the bandwidth using Householder
1109 * transformations.
1110 *
1111 * Note: we should get here only if LDA .ge. N
1112 *
1113  IF( isym.EQ.1 ) THEN
1114 *
1115 * Non-symmetric -- A = U D V
1116 *
1117  CALL zlagge( mr, nc, llb, uub, d, a, lda, iseed, work,
1118  $ iinfo )
1119  ELSE
1120 *
1121 * Symmetric -- A = U D U' or
1122 * Hermitian -- A = U D U*
1123 *
1124  IF( zsym ) THEN
1125  CALL zlagsy( m, llb, d, a, lda, iseed, work, iinfo )
1126  ELSE
1127  CALL zlaghe( m, llb, d, a, lda, iseed, work, iinfo )
1128  END IF
1129  END IF
1130 *
1131  IF( iinfo.NE.0 ) THEN
1132  info = 3
1133  RETURN
1134  END IF
1135  END IF
1136 *
1137 * 5) Pack the matrix
1138 *
1139  IF( ipack.NE.ipackg ) THEN
1140  IF( ipack.EQ.1 ) THEN
1141 *
1142 * 'U' -- Upper triangular, not packed
1143 *
1144  DO 370 j = 1, m
1145  DO 360 i = j + 1, m
1146  a( i, j ) = czero
1147  360 CONTINUE
1148  370 CONTINUE
1149 *
1150  ELSE IF( ipack.EQ.2 ) THEN
1151 *
1152 * 'L' -- Lower triangular, not packed
1153 *
1154  DO 390 j = 2, m
1155  DO 380 i = 1, j - 1
1156  a( i, j ) = czero
1157  380 CONTINUE
1158  390 CONTINUE
1159 *
1160  ELSE IF( ipack.EQ.3 ) THEN
1161 *
1162 * 'C' -- Upper triangle packed Columnwise.
1163 *
1164  icol = 1
1165  irow = 0
1166  DO 410 j = 1, m
1167  DO 400 i = 1, j
1168  irow = irow + 1
1169  IF( irow.GT.lda ) THEN
1170  irow = 1
1171  icol = icol + 1
1172  END IF
1173  a( irow, icol ) = a( i, j )
1174  400 CONTINUE
1175  410 CONTINUE
1176 *
1177  ELSE IF( ipack.EQ.4 ) THEN
1178 *
1179 * 'R' -- Lower triangle packed Columnwise.
1180 *
1181  icol = 1
1182  irow = 0
1183  DO 430 j = 1, m
1184  DO 420 i = j, m
1185  irow = irow + 1
1186  IF( irow.GT.lda ) THEN
1187  irow = 1
1188  icol = icol + 1
1189  END IF
1190  a( irow, icol ) = a( i, j )
1191  420 CONTINUE
1192  430 CONTINUE
1193 *
1194  ELSE IF( ipack.GE.5 ) THEN
1195 *
1196 * 'B' -- The lower triangle is packed as a band matrix.
1197 * 'Q' -- The upper triangle is packed as a band matrix.
1198 * 'Z' -- The whole matrix is packed as a band matrix.
1199 *
1200  IF( ipack.EQ.5 )
1201  $ uub = 0
1202  IF( ipack.EQ.6 )
1203  $ llb = 0
1204 *
1205  DO 450 j = 1, uub
1206  DO 440 i = min( j+llb, m ), 1, -1
1207  a( i-j+uub+1, j ) = a( i, j )
1208  440 CONTINUE
1209  450 CONTINUE
1210 *
1211  DO 470 j = uub + 2, n
1212  DO 460 i = j - uub, min( j+llb, m )
1213  a( i-j+uub+1, j ) = a( i, j )
1214  460 CONTINUE
1215  470 CONTINUE
1216  END IF
1217 *
1218 * If packed, zero out extraneous elements.
1219 *
1220 * Symmetric/Triangular Packed --
1221 * zero out everything after A(IROW,ICOL)
1222 *
1223  IF( ipack.EQ.3 .OR. ipack.EQ.4 ) THEN
1224  DO 490 jc = icol, m
1225  DO 480 jr = irow + 1, lda
1226  a( jr, jc ) = czero
1227  480 CONTINUE
1228  irow = 0
1229  490 CONTINUE
1230 *
1231  ELSE IF( ipack.GE.5 ) THEN
1232 *
1233 * Packed Band --
1234 * 1st row is now in A( UUB+2-j, j), zero above it
1235 * m-th row is now in A( M+UUB-j,j), zero below it
1236 * last non-zero diagonal is now in A( UUB+LLB+1,j ),
1237 * zero below it, too.
1238 *
1239  ir1 = uub + llb + 2
1240  ir2 = uub + m + 2
1241  DO 520 jc = 1, n
1242  DO 500 jr = 1, uub + 1 - jc
1243  a( jr, jc ) = czero
1244  500 CONTINUE
1245  DO 510 jr = max( 1, min( ir1, ir2-jc ) ), lda
1246  a( jr, jc ) = czero
1247  510 CONTINUE
1248  520 CONTINUE
1249  END IF
1250  END IF
1251 *
1252  RETURN
1253 *
1254 * End of ZLATMS
1255 *
subroutine dlatm1(MODE, COND, IRSIGN, IDIST, ISEED, D, N, INFO)
DLATM1
Definition: dlatm1.f:137
double precision function dlarnd(IDIST, ISEED)
DLARND
Definition: dlarnd.f:75
complex *16 function zlarnd(IDIST, ISEED)
ZLARND
Definition: zlarnd.f:77
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine zlarot(LROWS, LLEFT, LRIGHT, NL, C, S, A, LDA, XLEFT, XRIGHT)
ZLAROT
Definition: zlarot.f:231
subroutine zlartg(F, G, CS, SN, R)
ZLARTG generates a plane rotation with real cosine and complex sine.
Definition: zlartg.f:105
subroutine zlagsy(N, K, D, A, LDA, ISEED, WORK, INFO)
ZLAGSY
Definition: zlagsy.f:104
subroutine zlaset(UPLO, M, N, ALPHA, BETA, A, LDA)
ZLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values...
Definition: zlaset.f:108
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine zlagge(M, N, KL, KU, D, A, LDA, ISEED, WORK, INFO)
ZLAGGE
Definition: zlagge.f:116
subroutine dscal(N, DA, DX, INCX)
DSCAL
Definition: dscal.f:55
subroutine zlaghe(N, K, D, A, LDA, ISEED, WORK, INFO)
ZLAGHE
Definition: zlaghe.f:104

Here is the call graph for this function:

subroutine zlatmt ( integer  M,
integer  N,
character  DIST,
integer, dimension( 4 )  ISEED,
character  SYM,
double precision, dimension( * )  D,
integer  MODE,
double precision  COND,
double precision  DMAX,
integer  RANK,
integer  KL,
integer  KU,
character  PACK,
complex*16, dimension( lda, * )  A,
integer  LDA,
complex*16, dimension( * )  WORK,
integer  INFO 
)

ZLATMT

Purpose:
    ZLATMT generates random matrices with specified singular values
    (or hermitian with specified eigenvalues)
    for testing LAPACK programs.

    ZLATMT operates by applying the following sequence of
    operations:

      Set the diagonal to D, where D may be input or
         computed according to MODE, COND, DMAX, and SYM
         as described below.

      Generate a matrix with the appropriate band structure, by one
         of two methods:

      Method A:
          Generate a dense M x N matrix by multiplying D on the left
              and the right by random unitary matrices, then:

          Reduce the bandwidth according to KL and KU, using
              Householder transformations.

      Method B:
          Convert the bandwidth-0 (i.e., diagonal) matrix to a
              bandwidth-1 matrix using Givens rotations, "chasing"
              out-of-band elements back, much as in QR; then convert
              the bandwidth-1 to a bandwidth-2 matrix, etc.  Note
              that for reasonably small bandwidths (relative to M and
              N) this requires less storage, as a dense matrix is not
              generated.  Also, for hermitian or symmetric matrices,
              only one triangle is generated.

      Method A is chosen if the bandwidth is a large fraction of the
          order of the matrix, and LDA is at least M (so a dense
          matrix can be stored.)  Method B is chosen if the bandwidth
          is small (< 1/2 N for hermitian or symmetric, < .3 N+M for
          non-symmetric), or LDA is less than M and not less than the
          bandwidth.

      Pack the matrix if desired. Options specified by PACK are:
         no packing
         zero out upper half (if hermitian)
         zero out lower half (if hermitian)
         store the upper half columnwise (if hermitian or upper
               triangular)
         store the lower half columnwise (if hermitian or lower
               triangular)
         store the lower triangle in banded format (if hermitian or
               lower triangular)
         store the upper triangle in banded format (if hermitian or
               upper triangular)
         store the entire matrix in banded format
      If Method B is chosen, and band format is specified, then the
         matrix will be generated in the band format, so no repacking
         will be necessary.
Parameters
[in]M
          M is INTEGER
           The number of rows of A. Not modified.
[in]N
          N is INTEGER
           The number of columns of A. N must equal M if the matrix
           is symmetric or hermitian (i.e., if SYM is not 'N')
           Not modified.
[in]DIST
          DIST is CHARACTER*1
           On entry, DIST specifies the type of distribution to be used
           to generate the random eigen-/singular values.
           'U' => UNIFORM( 0, 1 )  ( 'U' for uniform )
           'S' => UNIFORM( -1, 1 ) ( 'S' for symmetric )
           'N' => NORMAL( 0, 1 )   ( 'N' for normal )
           Not modified.
[in,out]ISEED
          ISEED is INTEGER array, dimension ( 4 )
           On entry ISEED specifies the seed of the random number
           generator. They should lie between 0 and 4095 inclusive,
           and ISEED(4) should be odd. The random number generator
           uses a linear congruential sequence limited to small
           integers, and so should produce machine independent
           random numbers. The values of ISEED are changed on
           exit, and can be used in the next call to ZLATMT
           to continue the same random number sequence.
           Changed on exit.
[in]SYM
          SYM is CHARACTER*1
           If SYM='H', the generated matrix is hermitian, with
             eigenvalues specified by D, COND, MODE, and DMAX; they
             may be positive, negative, or zero.
           If SYM='P', the generated matrix is hermitian, with
             eigenvalues (= singular values) specified by D, COND,
             MODE, and DMAX; they will not be negative.
           If SYM='N', the generated matrix is nonsymmetric, with
             singular values specified by D, COND, MODE, and DMAX;
             they will not be negative.
           If SYM='S', the generated matrix is (complex) symmetric,
             with singular values specified by D, COND, MODE, and
             DMAX; they will not be negative.
           Not modified.
[in,out]D
          D is DOUBLE PRECISION array, dimension ( MIN( M, N ) )
           This array is used to specify the singular values or
           eigenvalues of A (see SYM, above.)  If MODE=0, then D is
           assumed to contain the singular/eigenvalues, otherwise
           they will be computed according to MODE, COND, and DMAX,
           and placed in D.
           Modified if MODE is nonzero.
[in]MODE
          MODE is INTEGER
           On entry this describes how the singular/eigenvalues are to
           be specified:
           MODE = 0 means use D as input
           MODE = 1 sets D(1)=1 and D(2:RANK)=1.0/COND
           MODE = 2 sets D(1:RANK-1)=1 and D(RANK)=1.0/COND
           MODE = 3 sets D(I)=COND**(-(I-1)/(RANK-1))
           MODE = 4 sets D(i)=1 - (i-1)/(N-1)*(1 - 1/COND)
           MODE = 5 sets D to random numbers in the range
                    ( 1/COND , 1 ) such that their logarithms
                    are uniformly distributed.
           MODE = 6 set D to random numbers from same distribution
                    as the rest of the matrix.
           MODE < 0 has the same meaning as ABS(MODE), except that
              the order of the elements of D is reversed.
           Thus if MODE is positive, D has entries ranging from
              1 to 1/COND, if negative, from 1/COND to 1,
           If SYM='H', and MODE is neither 0, 6, nor -6, then
              the elements of D will also be multiplied by a random
              sign (i.e., +1 or -1.)
           Not modified.
[in]COND
          COND is DOUBLE PRECISION
           On entry, this is used as described under MODE above.
           If used, it must be >= 1. Not modified.
[in]DMAX
          DMAX is DOUBLE PRECISION
           If MODE is neither -6, 0 nor 6, the contents of D, as
           computed according to MODE and COND, will be scaled by
           DMAX / max(abs(D(i))); thus, the maximum absolute eigen- or
           singular value (which is to say the norm) will be abs(DMAX).
           Note that DMAX need not be positive: if DMAX is negative
           (or zero), D will be scaled by a negative number (or zero).
           Not modified.
[in]RANK
          RANK is INTEGER
           The rank of matrix to be generated for modes 1,2,3 only.
           D( RANK+1:N ) = 0.
           Not modified.
[in]KL
          KL is INTEGER
           This specifies the lower bandwidth of the  matrix. For
           example, KL=0 implies upper triangular, KL=1 implies upper
           Hessenberg, and KL being at least M-1 means that the matrix
           has full lower bandwidth.  KL must equal KU if the matrix
           is symmetric or hermitian.
           Not modified.
[in]KU
          KU is INTEGER
           This specifies the upper bandwidth of the  matrix. For
           example, KU=0 implies lower triangular, KU=1 implies lower
           Hessenberg, and KU being at least N-1 means that the matrix
           has full upper bandwidth.  KL must equal KU if the matrix
           is symmetric or hermitian.
           Not modified.
[in]PACK
          PACK is CHARACTER*1
           This specifies packing of matrix as follows:
           'N' => no packing
           'U' => zero out all subdiagonal entries (if symmetric
                  or hermitian)
           'L' => zero out all superdiagonal entries (if symmetric
                  or hermitian)
           'C' => store the upper triangle columnwise (only if the
                  matrix is symmetric, hermitian, or upper triangular)
           'R' => store the lower triangle columnwise (only if the
                  matrix is symmetric, hermitian, or lower triangular)
           'B' => store the lower triangle in band storage scheme
                  (only if the matrix is symmetric, hermitian, or
                  lower triangular)
           'Q' => store the upper triangle in band storage scheme
                  (only if the matrix is symmetric, hermitian, or
                  upper triangular)
           'Z' => store the entire matrix in band storage scheme
                      (pivoting can be provided for by using this
                      option to store A in the trailing rows of
                      the allocated storage)

           Using these options, the various LAPACK packed and banded
           storage schemes can be obtained:
           GB                    - use 'Z'
           PB, SB, HB, or TB     - use 'B' or 'Q'
           PP, SP, HB, or TP     - use 'C' or 'R'

           If two calls to ZLATMT differ only in the PACK parameter,
           they will generate mathematically equivalent matrices.
           Not modified.
[in,out]A
          A is COMPLEX*16 array, dimension ( LDA, N )
           On exit A is the desired test matrix.  A is first generated
           in full (unpacked) form, and then packed, if so specified
           by PACK.  Thus, the first M elements of the first N
           columns will always be modified.  If PACK specifies a
           packed or banded storage scheme, all LDA elements of the
           first N columns will be modified; the elements of the
           array which do not correspond to elements of the generated
           matrix are set to zero.
           Modified.
[in]LDA
          LDA is INTEGER
           LDA specifies the first dimension of A as declared in the
           calling program.  If PACK='N', 'U', 'L', 'C', or 'R', then
           LDA must be at least M.  If PACK='B' or 'Q', then LDA must
           be at least MIN( KL, M-1) (which is equal to MIN(KU,N-1)).
           If PACK='Z', LDA must be large enough to hold the packed
           array: MIN( KU, N-1) + MIN( KL, M-1) + 1.
           Not modified.
[out]WORK
          WORK is COMPLEX*16 array, dimension ( 3*MAX( N, M ) )
           Workspace.
           Modified.
[out]INFO
          INFO is INTEGER
           Error code.  On exit, INFO will be set to one of the
           following values:
             0 => normal return
            -1 => M negative or unequal to N and SYM='S', 'H', or 'P'
            -2 => N negative
            -3 => DIST illegal string
            -5 => SYM illegal string
            -7 => MODE not in range -6 to 6
            -8 => COND less than 1.0, and MODE neither -6, 0 nor 6
           -10 => KL negative
           -11 => KU negative, or SYM is not 'N' and KU is not equal to
                  KL
           -12 => PACK illegal string, or PACK='U' or 'L', and SYM='N';
                  or PACK='C' or 'Q' and SYM='N' and KL is not zero;
                  or PACK='R' or 'B' and SYM='N' and KU is not zero;
                  or PACK='U', 'L', 'C', 'R', 'B', or 'Q', and M is not
                  N.
           -14 => LDA is less than M, or PACK='Z' and LDA is less than
                  MIN(KU,N-1) + MIN(KL,M-1) + 1.
            1  => Error return from DLATM7
            2  => Cannot scale to DMAX (max. sing. value is 0)
            3  => Error return from ZLAGGE, ZLAGHE or ZLAGSY
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 342 of file zlatmt.f.

342 *
343 * -- LAPACK computational routine (version 3.4.0) --
344 * -- LAPACK is a software package provided by Univ. of Tennessee, --
345 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
346 * November 2011
347 *
348 * .. Scalar Arguments ..
349  DOUBLE PRECISION cond, dmax
350  INTEGER info, kl, ku, lda, m, mode, n, rank
351  CHARACTER dist, pack, sym
352 * ..
353 * .. Array Arguments ..
354  COMPLEX*16 a( lda, * ), work( * )
355  DOUBLE PRECISION d( * )
356  INTEGER iseed( 4 )
357 * ..
358 *
359 * =====================================================================
360 *
361 * .. Parameters ..
362  DOUBLE PRECISION zero
363  parameter( zero = 0.0d+0 )
364  DOUBLE PRECISION one
365  parameter( one = 1.0d+0 )
366  COMPLEX*16 czero
367  parameter( czero = ( 0.0d+0, 0.0d+0 ) )
368  DOUBLE PRECISION twopi
369  parameter( twopi = 6.2831853071795864769252867663d+0 )
370 * ..
371 * .. Local Scalars ..
372  COMPLEX*16 c, ct, dummy, extra, s, st, ztemp
373  DOUBLE PRECISION alpha, angle, realc, temp
374  INTEGER i, ic, icol, idist, iendch, iinfo, il, ilda,
375  $ ioffg, ioffst, ipack, ipackg, ir, ir1, ir2,
376  $ irow, irsign, iskew, isym, isympk, j, jc, jch,
377  $ jkl, jku, jr, k, llb, minlda, mnmin, mr, nc,
378  $ uub
379  LOGICAL csym, givens, ilextr, iltemp, topdwn
380 * ..
381 * .. External Functions ..
382  COMPLEX*16 zlarnd
383  DOUBLE PRECISION dlarnd
384  LOGICAL lsame
385  EXTERNAL zlarnd, dlarnd, lsame
386 * ..
387 * .. External Subroutines ..
388  EXTERNAL dlatm7, dscal, xerbla, zlagge, zlaghe,
390 * ..
391 * .. Intrinsic Functions ..
392  INTRINSIC abs, cos, dble, dcmplx, dconjg, max, min, mod,
393  $ sin
394 * ..
395 * .. Executable Statements ..
396 *
397 * 1) Decode and Test the input parameters.
398 * Initialize flags & seed.
399 *
400  info = 0
401 *
402 * Quick return if possible
403 *
404  IF( m.EQ.0 .OR. n.EQ.0 )
405  $ RETURN
406 *
407 * Decode DIST
408 *
409  IF( lsame( dist, 'U' ) ) THEN
410  idist = 1
411  ELSE IF( lsame( dist, 'S' ) ) THEN
412  idist = 2
413  ELSE IF( lsame( dist, 'N' ) ) THEN
414  idist = 3
415  ELSE
416  idist = -1
417  END IF
418 *
419 * Decode SYM
420 *
421  IF( lsame( sym, 'N' ) ) THEN
422  isym = 1
423  irsign = 0
424  csym = .false.
425  ELSE IF( lsame( sym, 'P' ) ) THEN
426  isym = 2
427  irsign = 0
428  csym = .false.
429  ELSE IF( lsame( sym, 'S' ) ) THEN
430  isym = 2
431  irsign = 0
432  csym = .true.
433  ELSE IF( lsame( sym, 'H' ) ) THEN
434  isym = 2
435  irsign = 1
436  csym = .false.
437  ELSE
438  isym = -1
439  END IF
440 *
441 * Decode PACK
442 *
443  isympk = 0
444  IF( lsame( pack, 'N' ) ) THEN
445  ipack = 0
446  ELSE IF( lsame( pack, 'U' ) ) THEN
447  ipack = 1
448  isympk = 1
449  ELSE IF( lsame( pack, 'L' ) ) THEN
450  ipack = 2
451  isympk = 1
452  ELSE IF( lsame( pack, 'C' ) ) THEN
453  ipack = 3
454  isympk = 2
455  ELSE IF( lsame( pack, 'R' ) ) THEN
456  ipack = 4
457  isympk = 3
458  ELSE IF( lsame( pack, 'B' ) ) THEN
459  ipack = 5
460  isympk = 3
461  ELSE IF( lsame( pack, 'Q' ) ) THEN
462  ipack = 6
463  isympk = 2
464  ELSE IF( lsame( pack, 'Z' ) ) THEN
465  ipack = 7
466  ELSE
467  ipack = -1
468  END IF
469 *
470 * Set certain internal parameters
471 *
472  mnmin = min( m, n )
473  llb = min( kl, m-1 )
474  uub = min( ku, n-1 )
475  mr = min( m, n+llb )
476  nc = min( n, m+uub )
477 *
478  IF( ipack.EQ.5 .OR. ipack.EQ.6 ) THEN
479  minlda = uub + 1
480  ELSE IF( ipack.EQ.7 ) THEN
481  minlda = llb + uub + 1
482  ELSE
483  minlda = m
484  END IF
485 *
486 * Use Givens rotation method if bandwidth small enough,
487 * or if LDA is too small to store the matrix unpacked.
488 *
489  givens = .false.
490  IF( isym.EQ.1 ) THEN
491  IF( dble( llb+uub ).LT.0.3d0*dble( max( 1, mr+nc ) ) )
492  $ givens = .true.
493  ELSE
494  IF( 2*llb.LT.m )
495  $ givens = .true.
496  END IF
497  IF( lda.LT.m .AND. lda.GE.minlda )
498  $ givens = .true.
499 *
500 * Set INFO if an error
501 *
502  IF( m.LT.0 ) THEN
503  info = -1
504  ELSE IF( m.NE.n .AND. isym.NE.1 ) THEN
505  info = -1
506  ELSE IF( n.LT.0 ) THEN
507  info = -2
508  ELSE IF( idist.EQ.-1 ) THEN
509  info = -3
510  ELSE IF( isym.EQ.-1 ) THEN
511  info = -5
512  ELSE IF( abs( mode ).GT.6 ) THEN
513  info = -7
514  ELSE IF( ( mode.NE.0 .AND. abs( mode ).NE.6 ) .AND. cond.LT.one )
515  $ THEN
516  info = -8
517  ELSE IF( kl.LT.0 ) THEN
518  info = -10
519  ELSE IF( ku.LT.0 .OR. ( isym.NE.1 .AND. kl.NE.ku ) ) THEN
520  info = -11
521  ELSE IF( ipack.EQ.-1 .OR. ( isympk.EQ.1 .AND. isym.EQ.1 ) .OR.
522  $ ( isympk.EQ.2 .AND. isym.EQ.1 .AND. kl.GT.0 ) .OR.
523  $ ( isympk.EQ.3 .AND. isym.EQ.1 .AND. ku.GT.0 ) .OR.
524  $ ( isympk.NE.0 .AND. m.NE.n ) ) THEN
525  info = -12
526  ELSE IF( lda.LT.max( 1, minlda ) ) THEN
527  info = -14
528  END IF
529 *
530  IF( info.NE.0 ) THEN
531  CALL xerbla( 'ZLATMT', -info )
532  RETURN
533  END IF
534 *
535 * Initialize random number generator
536 *
537  DO 100 i = 1, 4
538  iseed( i ) = mod( abs( iseed( i ) ), 4096 )
539  100 CONTINUE
540 *
541  IF( mod( iseed( 4 ), 2 ).NE.1 )
542  $ iseed( 4 ) = iseed( 4 ) + 1
543 *
544 * 2) Set up D if indicated.
545 *
546 * Compute D according to COND and MODE
547 *
548  CALL dlatm7( mode, cond, irsign, idist, iseed, d, mnmin, rank,
549  $ iinfo )
550  IF( iinfo.NE.0 ) THEN
551  info = 1
552  RETURN
553  END IF
554 *
555 * Choose Top-Down if D is (apparently) increasing,
556 * Bottom-Up if D is (apparently) decreasing.
557 *
558  IF( abs( d( 1 ) ).LE.abs( d( rank ) ) ) THEN
559  topdwn = .true.
560  ELSE
561  topdwn = .false.
562  END IF
563 *
564  IF( mode.NE.0 .AND. abs( mode ).NE.6 ) THEN
565 *
566 * Scale by DMAX
567 *
568  temp = abs( d( 1 ) )
569  DO 110 i = 2, rank
570  temp = max( temp, abs( d( i ) ) )
571  110 CONTINUE
572 *
573  IF( temp.GT.zero ) THEN
574  alpha = dmax / temp
575  ELSE
576  info = 2
577  RETURN
578  END IF
579 *
580  CALL dscal( rank, alpha, d, 1 )
581 *
582  END IF
583 *
584  CALL zlaset( 'Full', lda, n, czero, czero, a, lda )
585 *
586 * 3) Generate Banded Matrix using Givens rotations.
587 * Also the special case of UUB=LLB=0
588 *
589 * Compute Addressing constants to cover all
590 * storage formats. Whether GE, HE, SY, GB, HB, or SB,
591 * upper or lower triangle or both,
592 * the (i,j)-th element is in
593 * A( i - ISKEW*j + IOFFST, j )
594 *
595  IF( ipack.GT.4 ) THEN
596  ilda = lda - 1
597  iskew = 1
598  IF( ipack.GT.5 ) THEN
599  ioffst = uub + 1
600  ELSE
601  ioffst = 1
602  END IF
603  ELSE
604  ilda = lda
605  iskew = 0
606  ioffst = 0
607  END IF
608 *
609 * IPACKG is the format that the matrix is generated in. If this is
610 * different from IPACK, then the matrix must be repacked at the
611 * end. It also signals how to compute the norm, for scaling.
612 *
613  ipackg = 0
614 *
615 * Diagonal Matrix -- We are done, unless it
616 * is to be stored HP/SP/PP/TP (PACK='R' or 'C')
617 *
618  IF( llb.EQ.0 .AND. uub.EQ.0 ) THEN
619  DO 120 j = 1, mnmin
620  a( ( 1-iskew )*j+ioffst, j ) = dcmplx( d( j ) )
621  120 CONTINUE
622 *
623  IF( ipack.LE.2 .OR. ipack.GE.5 )
624  $ ipackg = ipack
625 *
626  ELSE IF( givens ) THEN
627 *
628 * Check whether to use Givens rotations,
629 * Householder transformations, or nothing.
630 *
631  IF( isym.EQ.1 ) THEN
632 *
633 * Non-symmetric -- A = U D V
634 *
635  IF( ipack.GT.4 ) THEN
636  ipackg = ipack
637  ELSE
638  ipackg = 0
639  END IF
640 *
641  DO 130 j = 1, mnmin
642  a( ( 1-iskew )*j+ioffst, j ) = dcmplx( d( j ) )
643  130 CONTINUE
644 *
645  IF( topdwn ) THEN
646  jkl = 0
647  DO 160 jku = 1, uub
648 *
649 * Transform from bandwidth JKL, JKU-1 to JKL, JKU
650 *
651 * Last row actually rotated is M
652 * Last column actually rotated is MIN( M+JKU, N )
653 *
654  DO 150 jr = 1, min( m+jku, n ) + jkl - 1
655  extra = czero
656  angle = twopi*dlarnd( 1, iseed )
657  c = cos( angle )*zlarnd( 5, iseed )
658  s = sin( angle )*zlarnd( 5, iseed )
659  icol = max( 1, jr-jkl )
660  IF( jr.LT.m ) THEN
661  il = min( n, jr+jku ) + 1 - icol
662  CALL zlarot( .true., jr.GT.jkl, .false., il, c,
663  $ s, a( jr-iskew*icol+ioffst, icol ),
664  $ ilda, extra, dummy )
665  END IF
666 *
667 * Chase "EXTRA" back up
668 *
669  ir = jr
670  ic = icol
671  DO 140 jch = jr - jkl, 1, -jkl - jku
672  IF( ir.LT.m ) THEN
673  CALL zlartg( a( ir+1-iskew*( ic+1 )+ioffst,
674  $ ic+1 ), extra, realc, s, dummy )
675  dummy = dlarnd( 5, iseed )
676  c = dconjg( realc*dummy )
677  s = dconjg( -s*dummy )
678  END IF
679  irow = max( 1, jch-jku )
680  il = ir + 2 - irow
681  ztemp = czero
682  iltemp = jch.GT.jku
683  CALL zlarot( .false., iltemp, .true., il, c, s,
684  $ a( irow-iskew*ic+ioffst, ic ),
685  $ ilda, ztemp, extra )
686  IF( iltemp ) THEN
687  CALL zlartg( a( irow+1-iskew*( ic+1 )+ioffst,
688  $ ic+1 ), ztemp, realc, s, dummy )
689  dummy = zlarnd( 5, iseed )
690  c = dconjg( realc*dummy )
691  s = dconjg( -s*dummy )
692 *
693  icol = max( 1, jch-jku-jkl )
694  il = ic + 2 - icol
695  extra = czero
696  CALL zlarot( .true., jch.GT.jku+jkl, .true.,
697  $ il, c, s, a( irow-iskew*icol+
698  $ ioffst, icol ), ilda, extra,
699  $ ztemp )
700  ic = icol
701  ir = irow
702  END IF
703  140 CONTINUE
704  150 CONTINUE
705  160 CONTINUE
706 *
707  jku = uub
708  DO 190 jkl = 1, llb
709 *
710 * Transform from bandwidth JKL-1, JKU to JKL, JKU
711 *
712  DO 180 jc = 1, min( n+jkl, m ) + jku - 1
713  extra = czero
714  angle = twopi*dlarnd( 1, iseed )
715  c = cos( angle )*zlarnd( 5, iseed )
716  s = sin( angle )*zlarnd( 5, iseed )
717  irow = max( 1, jc-jku )
718  IF( jc.LT.n ) THEN
719  il = min( m, jc+jkl ) + 1 - irow
720  CALL zlarot( .false., jc.GT.jku, .false., il, c,
721  $ s, a( irow-iskew*jc+ioffst, jc ),
722  $ ilda, extra, dummy )
723  END IF
724 *
725 * Chase "EXTRA" back up
726 *
727  ic = jc
728  ir = irow
729  DO 170 jch = jc - jku, 1, -jkl - jku
730  IF( ic.LT.n ) THEN
731  CALL zlartg( a( ir+1-iskew*( ic+1 )+ioffst,
732  $ ic+1 ), extra, realc, s, dummy )
733  dummy = zlarnd( 5, iseed )
734  c = dconjg( realc*dummy )
735  s = dconjg( -s*dummy )
736  END IF
737  icol = max( 1, jch-jkl )
738  il = ic + 2 - icol
739  ztemp = czero
740  iltemp = jch.GT.jkl
741  CALL zlarot( .true., iltemp, .true., il, c, s,
742  $ a( ir-iskew*icol+ioffst, icol ),
743  $ ilda, ztemp, extra )
744  IF( iltemp ) THEN
745  CALL zlartg( a( ir+1-iskew*( icol+1 )+ioffst,
746  $ icol+1 ), ztemp, realc, s,
747  $ dummy )
748  dummy = zlarnd( 5, iseed )
749  c = dconjg( realc*dummy )
750  s = dconjg( -s*dummy )
751  irow = max( 1, jch-jkl-jku )
752  il = ir + 2 - irow
753  extra = czero
754  CALL zlarot( .false., jch.GT.jkl+jku, .true.,
755  $ il, c, s, a( irow-iskew*icol+
756  $ ioffst, icol ), ilda, extra,
757  $ ztemp )
758  ic = icol
759  ir = irow
760  END IF
761  170 CONTINUE
762  180 CONTINUE
763  190 CONTINUE
764 *
765  ELSE
766 *
767 * Bottom-Up -- Start at the bottom right.
768 *
769  jkl = 0
770  DO 220 jku = 1, uub
771 *
772 * Transform from bandwidth JKL, JKU-1 to JKL, JKU
773 *
774 * First row actually rotated is M
775 * First column actually rotated is MIN( M+JKU, N )
776 *
777  iendch = min( m, n+jkl ) - 1
778  DO 210 jc = min( m+jku, n ) - 1, 1 - jkl, -1
779  extra = czero
780  angle = twopi*dlarnd( 1, iseed )
781  c = cos( angle )*zlarnd( 5, iseed )
782  s = sin( angle )*zlarnd( 5, iseed )
783  irow = max( 1, jc-jku+1 )
784  IF( jc.GT.0 ) THEN
785  il = min( m, jc+jkl+1 ) + 1 - irow
786  CALL zlarot( .false., .false., jc+jkl.LT.m, il,
787  $ c, s, a( irow-iskew*jc+ioffst,
788  $ jc ), ilda, dummy, extra )
789  END IF
790 *
791 * Chase "EXTRA" back down
792 *
793  ic = jc
794  DO 200 jch = jc + jkl, iendch, jkl + jku
795  ilextr = ic.GT.0
796  IF( ilextr ) THEN
797  CALL zlartg( a( jch-iskew*ic+ioffst, ic ),
798  $ extra, realc, s, dummy )
799  dummy = zlarnd( 5, iseed )
800  c = realc*dummy
801  s = s*dummy
802  END IF
803  ic = max( 1, ic )
804  icol = min( n-1, jch+jku )
805  iltemp = jch + jku.LT.n
806  ztemp = czero
807  CALL zlarot( .true., ilextr, iltemp, icol+2-ic,
808  $ c, s, a( jch-iskew*ic+ioffst, ic ),
809  $ ilda, extra, ztemp )
810  IF( iltemp ) THEN
811  CALL zlartg( a( jch-iskew*icol+ioffst,
812  $ icol ), ztemp, realc, s, dummy )
813  dummy = zlarnd( 5, iseed )
814  c = realc*dummy
815  s = s*dummy
816  il = min( iendch, jch+jkl+jku ) + 2 - jch
817  extra = czero
818  CALL zlarot( .false., .true.,
819  $ jch+jkl+jku.LE.iendch, il, c, s,
820  $ a( jch-iskew*icol+ioffst,
821  $ icol ), ilda, ztemp, extra )
822  ic = icol
823  END IF
824  200 CONTINUE
825  210 CONTINUE
826  220 CONTINUE
827 *
828  jku = uub
829  DO 250 jkl = 1, llb
830 *
831 * Transform from bandwidth JKL-1, JKU to JKL, JKU
832 *
833 * First row actually rotated is MIN( N+JKL, M )
834 * First column actually rotated is N
835 *
836  iendch = min( n, m+jku ) - 1
837  DO 240 jr = min( n+jkl, m ) - 1, 1 - jku, -1
838  extra = czero
839  angle = twopi*dlarnd( 1, iseed )
840  c = cos( angle )*zlarnd( 5, iseed )
841  s = sin( angle )*zlarnd( 5, iseed )
842  icol = max( 1, jr-jkl+1 )
843  IF( jr.GT.0 ) THEN
844  il = min( n, jr+jku+1 ) + 1 - icol
845  CALL zlarot( .true., .false., jr+jku.LT.n, il,
846  $ c, s, a( jr-iskew*icol+ioffst,
847  $ icol ), ilda, dummy, extra )
848  END IF
849 *
850 * Chase "EXTRA" back down
851 *
852  ir = jr
853  DO 230 jch = jr + jku, iendch, jkl + jku
854  ilextr = ir.GT.0
855  IF( ilextr ) THEN
856  CALL zlartg( a( ir-iskew*jch+ioffst, jch ),
857  $ extra, realc, s, dummy )
858  dummy = zlarnd( 5, iseed )
859  c = realc*dummy
860  s = s*dummy
861  END IF
862  ir = max( 1, ir )
863  irow = min( m-1, jch+jkl )
864  iltemp = jch + jkl.LT.m
865  ztemp = czero
866  CALL zlarot( .false., ilextr, iltemp, irow+2-ir,
867  $ c, s, a( ir-iskew*jch+ioffst,
868  $ jch ), ilda, extra, ztemp )
869  IF( iltemp ) THEN
870  CALL zlartg( a( irow-iskew*jch+ioffst, jch ),
871  $ ztemp, realc, s, dummy )
872  dummy = zlarnd( 5, iseed )
873  c = realc*dummy
874  s = s*dummy
875  il = min( iendch, jch+jkl+jku ) + 2 - jch
876  extra = czero
877  CALL zlarot( .true., .true.,
878  $ jch+jkl+jku.LE.iendch, il, c, s,
879  $ a( irow-iskew*jch+ioffst, jch ),
880  $ ilda, ztemp, extra )
881  ir = irow
882  END IF
883  230 CONTINUE
884  240 CONTINUE
885  250 CONTINUE
886 *
887  END IF
888 *
889  ELSE
890 *
891 * Symmetric -- A = U D U'
892 * Hermitian -- A = U D U*
893 *
894  ipackg = ipack
895  ioffg = ioffst
896 *
897  IF( topdwn ) THEN
898 *
899 * Top-Down -- Generate Upper triangle only
900 *
901  IF( ipack.GE.5 ) THEN
902  ipackg = 6
903  ioffg = uub + 1
904  ELSE
905  ipackg = 1
906  END IF
907 *
908  DO 260 j = 1, mnmin
909  a( ( 1-iskew )*j+ioffg, j ) = dcmplx( d( j ) )
910  260 CONTINUE
911 *
912  DO 290 k = 1, uub
913  DO 280 jc = 1, n - 1
914  irow = max( 1, jc-k )
915  il = min( jc+1, k+2 )
916  extra = czero
917  ztemp = a( jc-iskew*( jc+1 )+ioffg, jc+1 )
918  angle = twopi*dlarnd( 1, iseed )
919  c = cos( angle )*zlarnd( 5, iseed )
920  s = sin( angle )*zlarnd( 5, iseed )
921  IF( csym ) THEN
922  ct = c
923  st = s
924  ELSE
925  ztemp = dconjg( ztemp )
926  ct = dconjg( c )
927  st = dconjg( s )
928  END IF
929  CALL zlarot( .false., jc.GT.k, .true., il, c, s,
930  $ a( irow-iskew*jc+ioffg, jc ), ilda,
931  $ extra, ztemp )
932  CALL zlarot( .true., .true., .false.,
933  $ min( k, n-jc )+1, ct, st,
934  $ a( ( 1-iskew )*jc+ioffg, jc ), ilda,
935  $ ztemp, dummy )
936 *
937 * Chase EXTRA back up the matrix
938 *
939  icol = jc
940  DO 270 jch = jc - k, 1, -k
941  CALL zlartg( a( jch+1-iskew*( icol+1 )+ioffg,
942  $ icol+1 ), extra, realc, s, dummy )
943  dummy = zlarnd( 5, iseed )
944  c = dconjg( realc*dummy )
945  s = dconjg( -s*dummy )
946  ztemp = a( jch-iskew*( jch+1 )+ioffg, jch+1 )
947  IF( csym ) THEN
948  ct = c
949  st = s
950  ELSE
951  ztemp = dconjg( ztemp )
952  ct = dconjg( c )
953  st = dconjg( s )
954  END IF
955  CALL zlarot( .true., .true., .true., k+2, c, s,
956  $ a( ( 1-iskew )*jch+ioffg, jch ),
957  $ ilda, ztemp, extra )
958  irow = max( 1, jch-k )
959  il = min( jch+1, k+2 )
960  extra = czero
961  CALL zlarot( .false., jch.GT.k, .true., il, ct,
962  $ st, a( irow-iskew*jch+ioffg, jch ),
963  $ ilda, extra, ztemp )
964  icol = jch
965  270 CONTINUE
966  280 CONTINUE
967  290 CONTINUE
968 *
969 * If we need lower triangle, copy from upper. Note that
970 * the order of copying is chosen to work for 'q' -> 'b'
971 *
972  IF( ipack.NE.ipackg .AND. ipack.NE.3 ) THEN
973  DO 320 jc = 1, n
974  irow = ioffst - iskew*jc
975  IF( csym ) THEN
976  DO 300 jr = jc, min( n, jc+uub )
977  a( jr+irow, jc ) = a( jc-iskew*jr+ioffg, jr )
978  300 CONTINUE
979  ELSE
980  DO 310 jr = jc, min( n, jc+uub )
981  a( jr+irow, jc ) = dconjg( a( jc-iskew*jr+
982  $ ioffg, jr ) )
983  310 CONTINUE
984  END IF
985  320 CONTINUE
986  IF( ipack.EQ.5 ) THEN
987  DO 340 jc = n - uub + 1, n
988  DO 330 jr = n + 2 - jc, uub + 1
989  a( jr, jc ) = czero
990  330 CONTINUE
991  340 CONTINUE
992  END IF
993  IF( ipackg.EQ.6 ) THEN
994  ipackg = ipack
995  ELSE
996  ipackg = 0
997  END IF
998  END IF
999  ELSE
1000 *
1001 * Bottom-Up -- Generate Lower triangle only
1002 *
1003  IF( ipack.GE.5 ) THEN
1004  ipackg = 5
1005  IF( ipack.EQ.6 )
1006  $ ioffg = 1
1007  ELSE
1008  ipackg = 2
1009  END IF
1010 *
1011  DO 350 j = 1, mnmin
1012  a( ( 1-iskew )*j+ioffg, j ) = dcmplx( d( j ) )
1013  350 CONTINUE
1014 *
1015  DO 380 k = 1, uub
1016  DO 370 jc = n - 1, 1, -1
1017  il = min( n+1-jc, k+2 )
1018  extra = czero
1019  ztemp = a( 1+( 1-iskew )*jc+ioffg, jc )
1020  angle = twopi*dlarnd( 1, iseed )
1021  c = cos( angle )*zlarnd( 5, iseed )
1022  s = sin( angle )*zlarnd( 5, iseed )
1023  IF( csym ) THEN
1024  ct = c
1025  st = s
1026  ELSE
1027  ztemp = dconjg( ztemp )
1028  ct = dconjg( c )
1029  st = dconjg( s )
1030  END IF
1031  CALL zlarot( .false., .true., n-jc.GT.k, il, c, s,
1032  $ a( ( 1-iskew )*jc+ioffg, jc ), ilda,
1033  $ ztemp, extra )
1034  icol = max( 1, jc-k+1 )
1035  CALL zlarot( .true., .false., .true., jc+2-icol,
1036  $ ct, st, a( jc-iskew*icol+ioffg,
1037  $ icol ), ilda, dummy, ztemp )
1038 *
1039 * Chase EXTRA back down the matrix
1040 *
1041  icol = jc
1042  DO 360 jch = jc + k, n - 1, k
1043  CALL zlartg( a( jch-iskew*icol+ioffg, icol ),
1044  $ extra, realc, s, dummy )
1045  dummy = zlarnd( 5, iseed )
1046  c = realc*dummy
1047  s = s*dummy
1048  ztemp = a( 1+( 1-iskew )*jch+ioffg, jch )
1049  IF( csym ) THEN
1050  ct = c
1051  st = s
1052  ELSE
1053  ztemp = dconjg( ztemp )
1054  ct = dconjg( c )
1055  st = dconjg( s )
1056  END IF
1057  CALL zlarot( .true., .true., .true., k+2, c, s,
1058  $ a( jch-iskew*icol+ioffg, icol ),
1059  $ ilda, extra, ztemp )
1060  il = min( n+1-jch, k+2 )
1061  extra = czero
1062  CALL zlarot( .false., .true., n-jch.GT.k, il,
1063  $ ct, st, a( ( 1-iskew )*jch+ioffg,
1064  $ jch ), ilda, ztemp, extra )
1065  icol = jch
1066  360 CONTINUE
1067  370 CONTINUE
1068  380 CONTINUE
1069 *
1070 * If we need upper triangle, copy from lower. Note that
1071 * the order of copying is chosen to work for 'b' -> 'q'
1072 *
1073  IF( ipack.NE.ipackg .AND. ipack.NE.4 ) THEN
1074  DO 410 jc = n, 1, -1
1075  irow = ioffst - iskew*jc
1076  IF( csym ) THEN
1077  DO 390 jr = jc, max( 1, jc-uub ), -1
1078  a( jr+irow, jc ) = a( jc-iskew*jr+ioffg, jr )
1079  390 CONTINUE
1080  ELSE
1081  DO 400 jr = jc, max( 1, jc-uub ), -1
1082  a( jr+irow, jc ) = dconjg( a( jc-iskew*jr+
1083  $ ioffg, jr ) )
1084  400 CONTINUE
1085  END IF
1086  410 CONTINUE
1087  IF( ipack.EQ.6 ) THEN
1088  DO 430 jc = 1, uub
1089  DO 420 jr = 1, uub + 1 - jc
1090  a( jr, jc ) = czero
1091  420 CONTINUE
1092  430 CONTINUE
1093  END IF
1094  IF( ipackg.EQ.5 ) THEN
1095  ipackg = ipack
1096  ELSE
1097  ipackg = 0
1098  END IF
1099  END IF
1100  END IF
1101 *
1102 * Ensure that the diagonal is real if Hermitian
1103 *
1104  IF( .NOT.csym ) THEN
1105  DO 440 jc = 1, n
1106  irow = ioffst + ( 1-iskew )*jc
1107  a( irow, jc ) = dcmplx( dble( a( irow, jc ) ) )
1108  440 CONTINUE
1109  END IF
1110 *
1111  END IF
1112 *
1113  ELSE
1114 *
1115 * 4) Generate Banded Matrix by first
1116 * Rotating by random Unitary matrices,
1117 * then reducing the bandwidth using Householder
1118 * transformations.
1119 *
1120 * Note: we should get here only if LDA .ge. N
1121 *
1122  IF( isym.EQ.1 ) THEN
1123 *
1124 * Non-symmetric -- A = U D V
1125 *
1126  CALL zlagge( mr, nc, llb, uub, d, a, lda, iseed, work,
1127  $ iinfo )
1128  ELSE
1129 *
1130 * Symmetric -- A = U D U' or
1131 * Hermitian -- A = U D U*
1132 *
1133  IF( csym ) THEN
1134  CALL zlagsy( m, llb, d, a, lda, iseed, work, iinfo )
1135  ELSE
1136  CALL zlaghe( m, llb, d, a, lda, iseed, work, iinfo )
1137  END IF
1138  END IF
1139 *
1140  IF( iinfo.NE.0 ) THEN
1141  info = 3
1142  RETURN
1143  END IF
1144  END IF
1145 *
1146 * 5) Pack the matrix
1147 *
1148  IF( ipack.NE.ipackg ) THEN
1149  IF( ipack.EQ.1 ) THEN
1150 *
1151 * 'U' -- Upper triangular, not packed
1152 *
1153  DO 460 j = 1, m
1154  DO 450 i = j + 1, m
1155  a( i, j ) = czero
1156  450 CONTINUE
1157  460 CONTINUE
1158 *
1159  ELSE IF( ipack.EQ.2 ) THEN
1160 *
1161 * 'L' -- Lower triangular, not packed
1162 *
1163  DO 480 j = 2, m
1164  DO 470 i = 1, j - 1
1165  a( i, j ) = czero
1166  470 CONTINUE
1167  480 CONTINUE
1168 *
1169  ELSE IF( ipack.EQ.3 ) THEN
1170 *
1171 * 'C' -- Upper triangle packed Columnwise.
1172 *
1173  icol = 1
1174  irow = 0
1175  DO 500 j = 1, m
1176  DO 490 i = 1, j
1177  irow = irow + 1
1178  IF( irow.GT.lda ) THEN
1179  irow = 1
1180  icol = icol + 1
1181  END IF
1182  a( irow, icol ) = a( i, j )
1183  490 CONTINUE
1184  500 CONTINUE
1185 *
1186  ELSE IF( ipack.EQ.4 ) THEN
1187 *
1188 * 'R' -- Lower triangle packed Columnwise.
1189 *
1190  icol = 1
1191  irow = 0
1192  DO 520 j = 1, m
1193  DO 510 i = j, m
1194  irow = irow + 1
1195  IF( irow.GT.lda ) THEN
1196  irow = 1
1197  icol = icol + 1
1198  END IF
1199  a( irow, icol ) = a( i, j )
1200  510 CONTINUE
1201  520 CONTINUE
1202 *
1203  ELSE IF( ipack.GE.5 ) THEN
1204 *
1205 * 'B' -- The lower triangle is packed as a band matrix.
1206 * 'Q' -- The upper triangle is packed as a band matrix.
1207 * 'Z' -- The whole matrix is packed as a band matrix.
1208 *
1209  IF( ipack.EQ.5 )
1210  $ uub = 0
1211  IF( ipack.EQ.6 )
1212  $ llb = 0
1213 *
1214  DO 540 j = 1, uub
1215  DO 530 i = min( j+llb, m ), 1, -1
1216  a( i-j+uub+1, j ) = a( i, j )
1217  530 CONTINUE
1218  540 CONTINUE
1219 *
1220  DO 560 j = uub + 2, n
1221  DO 550 i = j - uub, min( j+llb, m )
1222  a( i-j+uub+1, j ) = a( i, j )
1223  550 CONTINUE
1224  560 CONTINUE
1225  END IF
1226 *
1227 * If packed, zero out extraneous elements.
1228 *
1229 * Symmetric/Triangular Packed --
1230 * zero out everything after A(IROW,ICOL)
1231 *
1232  IF( ipack.EQ.3 .OR. ipack.EQ.4 ) THEN
1233  DO 580 jc = icol, m
1234  DO 570 jr = irow + 1, lda
1235  a( jr, jc ) = czero
1236  570 CONTINUE
1237  irow = 0
1238  580 CONTINUE
1239 *
1240  ELSE IF( ipack.GE.5 ) THEN
1241 *
1242 * Packed Band --
1243 * 1st row is now in A( UUB+2-j, j), zero above it
1244 * m-th row is now in A( M+UUB-j,j), zero below it
1245 * last non-zero diagonal is now in A( UUB+LLB+1,j ),
1246 * zero below it, too.
1247 *
1248  ir1 = uub + llb + 2
1249  ir2 = uub + m + 2
1250  DO 610 jc = 1, n
1251  DO 590 jr = 1, uub + 1 - jc
1252  a( jr, jc ) = czero
1253  590 CONTINUE
1254  DO 600 jr = max( 1, min( ir1, ir2-jc ) ), lda
1255  a( jr, jc ) = czero
1256  600 CONTINUE
1257  610 CONTINUE
1258  END IF
1259  END IF
1260 *
1261  RETURN
1262 *
1263 * End of ZLATMT
1264 *
double precision function dlarnd(IDIST, ISEED)
DLARND
Definition: dlarnd.f:75
complex *16 function zlarnd(IDIST, ISEED)
ZLARND
Definition: zlarnd.f:77
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine zlarot(LROWS, LLEFT, LRIGHT, NL, C, S, A, LDA, XLEFT, XRIGHT)
ZLAROT
Definition: zlarot.f:231
subroutine zlartg(F, G, CS, SN, R)
ZLARTG generates a plane rotation with real cosine and complex sine.
Definition: zlartg.f:105
subroutine zlagsy(N, K, D, A, LDA, ISEED, WORK, INFO)
ZLAGSY
Definition: zlagsy.f:104
subroutine zlaset(UPLO, M, N, ALPHA, BETA, A, LDA)
ZLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values...
Definition: zlaset.f:108
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine zlagge(M, N, KL, KU, D, A, LDA, ISEED, WORK, INFO)
ZLAGGE
Definition: zlagge.f:116
subroutine dscal(N, DA, DX, INCX)
DSCAL
Definition: dscal.f:55
subroutine dlatm7(MODE, COND, IRSIGN, IDIST, ISEED, D, N, RANK, INFO)
DLATM7
Definition: dlatm7.f:124
subroutine zlaghe(N, K, D, A, LDA, ISEED, WORK, INFO)
ZLAGHE
Definition: zlaghe.f:104

Here is the call graph for this function:

Here is the caller graph for this function: