LAPACK 3.11.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ dlattr()

subroutine dlattr ( integer  IMAT,
character  UPLO,
character  TRANS,
character  DIAG,
integer, dimension( 4 )  ISEED,
integer  N,
double precision, dimension( lda, * )  A,
integer  LDA,
double precision, dimension( * )  B,
double precision, dimension( * )  WORK,
integer  INFO 
)

DLATTR

Purpose:
 DLATTR generates a triangular test matrix.
 IMAT and UPLO uniquely specify the properties of the test
 matrix, which is returned in the array A.
Parameters
[in]IMAT
          IMAT is INTEGER
          An integer key describing which matrix to generate for this
          path.
[in]UPLO
          UPLO is CHARACTER*1
          Specifies whether the matrix A will be upper or lower
          triangular.
          = 'U':  Upper triangular
          = 'L':  Lower triangular
[in]TRANS
          TRANS is CHARACTER*1
          Specifies whether the matrix or its transpose will be used.
          = 'N':  No transpose
          = 'T':  Transpose
          = 'C':  Conjugate transpose (= Transpose)
[out]DIAG
          DIAG is CHARACTER*1
          Specifies whether or not the matrix A is unit triangular.
          = 'N':  Non-unit triangular
          = 'U':  Unit triangular
[in,out]ISEED
          ISEED is INTEGER array, dimension (4)
          The seed vector for the random number generator (used in
          DLATMS).  Modified on exit.
[in]N
          N is INTEGER
          The order of the matrix to be generated.
[out]A
          A is DOUBLE PRECISION array, dimension (LDA,N)
          The triangular matrix A.  If UPLO = 'U', the leading n by n
          upper triangular part of the array A contains the upper
          triangular matrix, and the strictly lower triangular part of
          A is not referenced.  If UPLO = 'L', the leading n by n lower
          triangular part of the array A contains the lower triangular
          matrix, and the strictly upper triangular part of A is not
          referenced.  If DIAG = 'U', the diagonal elements of A are
          set so that A(k,k) = k for 1 <= k <= n.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,N).
[out]B
          B is DOUBLE PRECISION array, dimension (N)
          The right hand side vector, if IMAT > 10.
[out]WORK
          WORK is DOUBLE PRECISION array, dimension (3*N)
[out]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0: if INFO = -k, the k-th argument had an illegal value
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 131 of file dlattr.f.

133*
134* -- LAPACK test routine --
135* -- LAPACK is a software package provided by Univ. of Tennessee, --
136* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
137*
138* .. Scalar Arguments ..
139 CHARACTER DIAG, TRANS, UPLO
140 INTEGER IMAT, INFO, LDA, N
141* ..
142* .. Array Arguments ..
143 INTEGER ISEED( 4 )
144 DOUBLE PRECISION A( LDA, * ), B( * ), WORK( * )
145* ..
146*
147* =====================================================================
148*
149* .. Parameters ..
150 DOUBLE PRECISION ONE, TWO, ZERO
151 parameter( one = 1.0d+0, two = 2.0d+0, zero = 0.0d+0 )
152* ..
153* .. Local Scalars ..
154 LOGICAL UPPER
155 CHARACTER DIST, TYPE
156 CHARACTER*3 PATH
157 INTEGER I, IY, J, JCOUNT, KL, KU, MODE
158 DOUBLE PRECISION ANORM, BIGNUM, BNORM, BSCAL, C, CNDNUM, PLUS1,
159 $ PLUS2, RA, RB, REXP, S, SFAC, SMLNUM, STAR1,
160 $ TEXP, TLEFT, TSCAL, ULP, UNFL, X, Y, Z
161* ..
162* .. External Functions ..
163 LOGICAL LSAME
164 INTEGER IDAMAX
165 DOUBLE PRECISION DLAMCH, DLARND
166 EXTERNAL lsame, idamax, dlamch, dlarnd
167* ..
168* .. External Subroutines ..
169 EXTERNAL dcopy, dlabad, dlarnv, dlatb4, dlatms, drot,
170 $ drotg, dscal, dswap
171* ..
172* .. Intrinsic Functions ..
173 INTRINSIC abs, dble, max, sign, sqrt
174* ..
175* .. Executable Statements ..
176*
177 path( 1: 1 ) = 'Double precision'
178 path( 2: 3 ) = 'TR'
179 unfl = dlamch( 'Safe minimum' )
180 ulp = dlamch( 'Epsilon' )*dlamch( 'Base' )
181 smlnum = unfl
182 bignum = ( one-ulp ) / smlnum
183 CALL dlabad( smlnum, bignum )
184 IF( ( imat.GE.7 .AND. imat.LE.10 ) .OR. imat.EQ.18 ) THEN
185 diag = 'U'
186 ELSE
187 diag = 'N'
188 END IF
189 info = 0
190*
191* Quick return if N.LE.0.
192*
193 IF( n.LE.0 )
194 $ RETURN
195*
196* Call DLATB4 to set parameters for DLATMS.
197*
198 upper = lsame( uplo, 'U' )
199 IF( upper ) THEN
200 CALL dlatb4( path, imat, n, n, TYPE, KL, KU, ANORM, MODE,
201 $ CNDNUM, DIST )
202 ELSE
203 CALL dlatb4( path, -imat, n, n, TYPE, KL, KU, ANORM, MODE,
204 $ CNDNUM, DIST )
205 END IF
206*
207* IMAT <= 6: Non-unit triangular matrix
208*
209 IF( imat.LE.6 ) THEN
210 CALL dlatms( n, n, dist, iseed, TYPE, B, MODE, CNDNUM, ANORM,
211 $ KL, KU, 'No packing', A, LDA, WORK, INFO )
212*
213* IMAT > 6: Unit triangular matrix
214* The diagonal is deliberately set to something other than 1.
215*
216* IMAT = 7: Matrix is the identity
217*
218 ELSE IF( imat.EQ.7 ) THEN
219 IF( upper ) THEN
220 DO 20 j = 1, n
221 DO 10 i = 1, j - 1
222 a( i, j ) = zero
223 10 CONTINUE
224 a( j, j ) = j
225 20 CONTINUE
226 ELSE
227 DO 40 j = 1, n
228 a( j, j ) = j
229 DO 30 i = j + 1, n
230 a( i, j ) = zero
231 30 CONTINUE
232 40 CONTINUE
233 END IF
234*
235* IMAT > 7: Non-trivial unit triangular matrix
236*
237* Generate a unit triangular matrix T with condition CNDNUM by
238* forming a triangular matrix with known singular values and
239* filling in the zero entries with Givens rotations.
240*
241 ELSE IF( imat.LE.10 ) THEN
242 IF( upper ) THEN
243 DO 60 j = 1, n
244 DO 50 i = 1, j - 1
245 a( i, j ) = zero
246 50 CONTINUE
247 a( j, j ) = j
248 60 CONTINUE
249 ELSE
250 DO 80 j = 1, n
251 a( j, j ) = j
252 DO 70 i = j + 1, n
253 a( i, j ) = zero
254 70 CONTINUE
255 80 CONTINUE
256 END IF
257*
258* Since the trace of a unit triangular matrix is 1, the product
259* of its singular values must be 1. Let s = sqrt(CNDNUM),
260* x = sqrt(s) - 1/sqrt(s), y = sqrt(2/(n-2))*x, and z = x**2.
261* The following triangular matrix has singular values s, 1, 1,
262* ..., 1, 1/s:
263*
264* 1 y y y ... y y z
265* 1 0 0 ... 0 0 y
266* 1 0 ... 0 0 y
267* . ... . . .
268* . . . .
269* 1 0 y
270* 1 y
271* 1
272*
273* To fill in the zeros, we first multiply by a matrix with small
274* condition number of the form
275*
276* 1 0 0 0 0 ...
277* 1 + * 0 0 ...
278* 1 + 0 0 0
279* 1 + * 0 0
280* 1 + 0 0
281* ...
282* 1 + 0
283* 1 0
284* 1
285*
286* Each element marked with a '*' is formed by taking the product
287* of the adjacent elements marked with '+'. The '*'s can be
288* chosen freely, and the '+'s are chosen so that the inverse of
289* T will have elements of the same magnitude as T. If the *'s in
290* both T and inv(T) have small magnitude, T is well conditioned.
291* The two offdiagonals of T are stored in WORK.
292*
293* The product of these two matrices has the form
294*
295* 1 y y y y y . y y z
296* 1 + * 0 0 . 0 0 y
297* 1 + 0 0 . 0 0 y
298* 1 + * . . . .
299* 1 + . . . .
300* . . . . .
301* . . . .
302* 1 + y
303* 1 y
304* 1
305*
306* Now we multiply by Givens rotations, using the fact that
307*
308* [ c s ] [ 1 w ] [ -c -s ] = [ 1 -w ]
309* [ -s c ] [ 0 1 ] [ s -c ] [ 0 1 ]
310* and
311* [ -c -s ] [ 1 0 ] [ c s ] = [ 1 0 ]
312* [ s -c ] [ w 1 ] [ -s c ] [ -w 1 ]
313*
314* where c = w / sqrt(w**2+4) and s = 2 / sqrt(w**2+4).
315*
316 star1 = 0.25d0
317 sfac = 0.5d0
318 plus1 = sfac
319 DO 90 j = 1, n, 2
320 plus2 = star1 / plus1
321 work( j ) = plus1
322 work( n+j ) = star1
323 IF( j+1.LE.n ) THEN
324 work( j+1 ) = plus2
325 work( n+j+1 ) = zero
326 plus1 = star1 / plus2
327 rexp = dlarnd( 2, iseed )
328 star1 = star1*( sfac**rexp )
329 IF( rexp.LT.zero ) THEN
330 star1 = -sfac**( one-rexp )
331 ELSE
332 star1 = sfac**( one+rexp )
333 END IF
334 END IF
335 90 CONTINUE
336*
337 x = sqrt( cndnum ) - 1 / sqrt( cndnum )
338 IF( n.GT.2 ) THEN
339 y = sqrt( 2.d0 / ( n-2 ) )*x
340 ELSE
341 y = zero
342 END IF
343 z = x*x
344*
345 IF( upper ) THEN
346 IF( n.GT.3 ) THEN
347 CALL dcopy( n-3, work, 1, a( 2, 3 ), lda+1 )
348 IF( n.GT.4 )
349 $ CALL dcopy( n-4, work( n+1 ), 1, a( 2, 4 ), lda+1 )
350 END IF
351 DO 100 j = 2, n - 1
352 a( 1, j ) = y
353 a( j, n ) = y
354 100 CONTINUE
355 a( 1, n ) = z
356 ELSE
357 IF( n.GT.3 ) THEN
358 CALL dcopy( n-3, work, 1, a( 3, 2 ), lda+1 )
359 IF( n.GT.4 )
360 $ CALL dcopy( n-4, work( n+1 ), 1, a( 4, 2 ), lda+1 )
361 END IF
362 DO 110 j = 2, n - 1
363 a( j, 1 ) = y
364 a( n, j ) = y
365 110 CONTINUE
366 a( n, 1 ) = z
367 END IF
368*
369* Fill in the zeros using Givens rotations.
370*
371 IF( upper ) THEN
372 DO 120 j = 1, n - 1
373 ra = a( j, j+1 )
374 rb = 2.0d0
375 CALL drotg( ra, rb, c, s )
376*
377* Multiply by [ c s; -s c] on the left.
378*
379 IF( n.GT.j+1 )
380 $ CALL drot( n-j-1, a( j, j+2 ), lda, a( j+1, j+2 ),
381 $ lda, c, s )
382*
383* Multiply by [-c -s; s -c] on the right.
384*
385 IF( j.GT.1 )
386 $ CALL drot( j-1, a( 1, j+1 ), 1, a( 1, j ), 1, -c, -s )
387*
388* Negate A(J,J+1).
389*
390 a( j, j+1 ) = -a( j, j+1 )
391 120 CONTINUE
392 ELSE
393 DO 130 j = 1, n - 1
394 ra = a( j+1, j )
395 rb = 2.0d0
396 CALL drotg( ra, rb, c, s )
397*
398* Multiply by [ c -s; s c] on the right.
399*
400 IF( n.GT.j+1 )
401 $ CALL drot( n-j-1, a( j+2, j+1 ), 1, a( j+2, j ), 1, c,
402 $ -s )
403*
404* Multiply by [-c s; -s -c] on the left.
405*
406 IF( j.GT.1 )
407 $ CALL drot( j-1, a( j, 1 ), lda, a( j+1, 1 ), lda, -c,
408 $ s )
409*
410* Negate A(J+1,J).
411*
412 a( j+1, j ) = -a( j+1, j )
413 130 CONTINUE
414 END IF
415*
416* IMAT > 10: Pathological test cases. These triangular matrices
417* are badly scaled or badly conditioned, so when used in solving a
418* triangular system they may cause overflow in the solution vector.
419*
420 ELSE IF( imat.EQ.11 ) THEN
421*
422* Type 11: Generate a triangular matrix with elements between
423* -1 and 1. Give the diagonal norm 2 to make it well-conditioned.
424* Make the right hand side large so that it requires scaling.
425*
426 IF( upper ) THEN
427 DO 140 j = 1, n
428 CALL dlarnv( 2, iseed, j, a( 1, j ) )
429 a( j, j ) = sign( two, a( j, j ) )
430 140 CONTINUE
431 ELSE
432 DO 150 j = 1, n
433 CALL dlarnv( 2, iseed, n-j+1, a( j, j ) )
434 a( j, j ) = sign( two, a( j, j ) )
435 150 CONTINUE
436 END IF
437*
438* Set the right hand side so that the largest value is BIGNUM.
439*
440 CALL dlarnv( 2, iseed, n, b )
441 iy = idamax( n, b, 1 )
442 bnorm = abs( b( iy ) )
443 bscal = bignum / max( one, bnorm )
444 CALL dscal( n, bscal, b, 1 )
445*
446 ELSE IF( imat.EQ.12 ) THEN
447*
448* Type 12: Make the first diagonal element in the solve small to
449* cause immediate overflow when dividing by T(j,j).
450* In type 12, the offdiagonal elements are small (CNORM(j) < 1).
451*
452 CALL dlarnv( 2, iseed, n, b )
453 tscal = one / max( one, dble( n-1 ) )
454 IF( upper ) THEN
455 DO 160 j = 1, n
456 CALL dlarnv( 2, iseed, j, a( 1, j ) )
457 CALL dscal( j-1, tscal, a( 1, j ), 1 )
458 a( j, j ) = sign( one, a( j, j ) )
459 160 CONTINUE
460 a( n, n ) = smlnum*a( n, n )
461 ELSE
462 DO 170 j = 1, n
463 CALL dlarnv( 2, iseed, n-j+1, a( j, j ) )
464 IF( n.GT.j )
465 $ CALL dscal( n-j, tscal, a( j+1, j ), 1 )
466 a( j, j ) = sign( one, a( j, j ) )
467 170 CONTINUE
468 a( 1, 1 ) = smlnum*a( 1, 1 )
469 END IF
470*
471 ELSE IF( imat.EQ.13 ) THEN
472*
473* Type 13: Make the first diagonal element in the solve small to
474* cause immediate overflow when dividing by T(j,j).
475* In type 13, the offdiagonal elements are O(1) (CNORM(j) > 1).
476*
477 CALL dlarnv( 2, iseed, n, b )
478 IF( upper ) THEN
479 DO 180 j = 1, n
480 CALL dlarnv( 2, iseed, j, a( 1, j ) )
481 a( j, j ) = sign( one, a( j, j ) )
482 180 CONTINUE
483 a( n, n ) = smlnum*a( n, n )
484 ELSE
485 DO 190 j = 1, n
486 CALL dlarnv( 2, iseed, n-j+1, a( j, j ) )
487 a( j, j ) = sign( one, a( j, j ) )
488 190 CONTINUE
489 a( 1, 1 ) = smlnum*a( 1, 1 )
490 END IF
491*
492 ELSE IF( imat.EQ.14 ) THEN
493*
494* Type 14: T is diagonal with small numbers on the diagonal to
495* make the growth factor underflow, but a small right hand side
496* chosen so that the solution does not overflow.
497*
498 IF( upper ) THEN
499 jcount = 1
500 DO 210 j = n, 1, -1
501 DO 200 i = 1, j - 1
502 a( i, j ) = zero
503 200 CONTINUE
504 IF( jcount.LE.2 ) THEN
505 a( j, j ) = smlnum
506 ELSE
507 a( j, j ) = one
508 END IF
509 jcount = jcount + 1
510 IF( jcount.GT.4 )
511 $ jcount = 1
512 210 CONTINUE
513 ELSE
514 jcount = 1
515 DO 230 j = 1, n
516 DO 220 i = j + 1, n
517 a( i, j ) = zero
518 220 CONTINUE
519 IF( jcount.LE.2 ) THEN
520 a( j, j ) = smlnum
521 ELSE
522 a( j, j ) = one
523 END IF
524 jcount = jcount + 1
525 IF( jcount.GT.4 )
526 $ jcount = 1
527 230 CONTINUE
528 END IF
529*
530* Set the right hand side alternately zero and small.
531*
532 IF( upper ) THEN
533 b( 1 ) = zero
534 DO 240 i = n, 2, -2
535 b( i ) = zero
536 b( i-1 ) = smlnum
537 240 CONTINUE
538 ELSE
539 b( n ) = zero
540 DO 250 i = 1, n - 1, 2
541 b( i ) = zero
542 b( i+1 ) = smlnum
543 250 CONTINUE
544 END IF
545*
546 ELSE IF( imat.EQ.15 ) THEN
547*
548* Type 15: Make the diagonal elements small to cause gradual
549* overflow when dividing by T(j,j). To control the amount of
550* scaling needed, the matrix is bidiagonal.
551*
552 texp = one / max( one, dble( n-1 ) )
553 tscal = smlnum**texp
554 CALL dlarnv( 2, iseed, n, b )
555 IF( upper ) THEN
556 DO 270 j = 1, n
557 DO 260 i = 1, j - 2
558 a( i, j ) = 0.d0
559 260 CONTINUE
560 IF( j.GT.1 )
561 $ a( j-1, j ) = -one
562 a( j, j ) = tscal
563 270 CONTINUE
564 b( n ) = one
565 ELSE
566 DO 290 j = 1, n
567 DO 280 i = j + 2, n
568 a( i, j ) = 0.d0
569 280 CONTINUE
570 IF( j.LT.n )
571 $ a( j+1, j ) = -one
572 a( j, j ) = tscal
573 290 CONTINUE
574 b( 1 ) = one
575 END IF
576*
577 ELSE IF( imat.EQ.16 ) THEN
578*
579* Type 16: One zero diagonal element.
580*
581 iy = n / 2 + 1
582 IF( upper ) THEN
583 DO 300 j = 1, n
584 CALL dlarnv( 2, iseed, j, a( 1, j ) )
585 IF( j.NE.iy ) THEN
586 a( j, j ) = sign( two, a( j, j ) )
587 ELSE
588 a( j, j ) = zero
589 END IF
590 300 CONTINUE
591 ELSE
592 DO 310 j = 1, n
593 CALL dlarnv( 2, iseed, n-j+1, a( j, j ) )
594 IF( j.NE.iy ) THEN
595 a( j, j ) = sign( two, a( j, j ) )
596 ELSE
597 a( j, j ) = zero
598 END IF
599 310 CONTINUE
600 END IF
601 CALL dlarnv( 2, iseed, n, b )
602 CALL dscal( n, two, b, 1 )
603*
604 ELSE IF( imat.EQ.17 ) THEN
605*
606* Type 17: Make the offdiagonal elements large to cause overflow
607* when adding a column of T. In the non-transposed case, the
608* matrix is constructed to cause overflow when adding a column in
609* every other step.
610*
611 tscal = unfl / ulp
612 tscal = ( one-ulp ) / tscal
613 DO 330 j = 1, n
614 DO 320 i = 1, n
615 a( i, j ) = 0.d0
616 320 CONTINUE
617 330 CONTINUE
618 texp = one
619 IF( upper ) THEN
620 DO 340 j = n, 2, -2
621 a( 1, j ) = -tscal / dble( n+1 )
622 a( j, j ) = one
623 b( j ) = texp*( one-ulp )
624 a( 1, j-1 ) = -( tscal / dble( n+1 ) ) / dble( n+2 )
625 a( j-1, j-1 ) = one
626 b( j-1 ) = texp*dble( n*n+n-1 )
627 texp = texp*2.d0
628 340 CONTINUE
629 b( 1 ) = ( dble( n+1 ) / dble( n+2 ) )*tscal
630 ELSE
631 DO 350 j = 1, n - 1, 2
632 a( n, j ) = -tscal / dble( n+1 )
633 a( j, j ) = one
634 b( j ) = texp*( one-ulp )
635 a( n, j+1 ) = -( tscal / dble( n+1 ) ) / dble( n+2 )
636 a( j+1, j+1 ) = one
637 b( j+1 ) = texp*dble( n*n+n-1 )
638 texp = texp*2.d0
639 350 CONTINUE
640 b( n ) = ( dble( n+1 ) / dble( n+2 ) )*tscal
641 END IF
642*
643 ELSE IF( imat.EQ.18 ) THEN
644*
645* Type 18: Generate a unit triangular matrix with elements
646* between -1 and 1, and make the right hand side large so that it
647* requires scaling.
648*
649 IF( upper ) THEN
650 DO 360 j = 1, n
651 CALL dlarnv( 2, iseed, j-1, a( 1, j ) )
652 a( j, j ) = zero
653 360 CONTINUE
654 ELSE
655 DO 370 j = 1, n
656 IF( j.LT.n )
657 $ CALL dlarnv( 2, iseed, n-j, a( j+1, j ) )
658 a( j, j ) = zero
659 370 CONTINUE
660 END IF
661*
662* Set the right hand side so that the largest value is BIGNUM.
663*
664 CALL dlarnv( 2, iseed, n, b )
665 iy = idamax( n, b, 1 )
666 bnorm = abs( b( iy ) )
667 bscal = bignum / max( one, bnorm )
668 CALL dscal( n, bscal, b, 1 )
669*
670 ELSE IF( imat.EQ.19 ) THEN
671*
672* Type 19: Generate a triangular matrix with elements between
673* BIGNUM/(n-1) and BIGNUM so that at least one of the column
674* norms will exceed BIGNUM.
675* 1/3/91: DLATRS no longer can handle this case
676*
677 tleft = bignum / max( one, dble( n-1 ) )
678 tscal = bignum*( dble( n-1 ) / max( one, dble( n ) ) )
679 IF( upper ) THEN
680 DO 390 j = 1, n
681 CALL dlarnv( 2, iseed, j, a( 1, j ) )
682 DO 380 i = 1, j
683 a( i, j ) = sign( tleft, a( i, j ) ) + tscal*a( i, j )
684 380 CONTINUE
685 390 CONTINUE
686 ELSE
687 DO 410 j = 1, n
688 CALL dlarnv( 2, iseed, n-j+1, a( j, j ) )
689 DO 400 i = j, n
690 a( i, j ) = sign( tleft, a( i, j ) ) + tscal*a( i, j )
691 400 CONTINUE
692 410 CONTINUE
693 END IF
694 CALL dlarnv( 2, iseed, n, b )
695 CALL dscal( n, two, b, 1 )
696 END IF
697*
698* Flip the matrix if the transpose will be used.
699*
700 IF( .NOT.lsame( trans, 'N' ) ) THEN
701 IF( upper ) THEN
702 DO 420 j = 1, n / 2
703 CALL dswap( n-2*j+1, a( j, j ), lda, a( j+1, n-j+1 ),
704 $ -1 )
705 420 CONTINUE
706 ELSE
707 DO 430 j = 1, n / 2
708 CALL dswap( n-2*j+1, a( j, j ), 1, a( n-j+1, j+1 ),
709 $ -lda )
710 430 CONTINUE
711 END IF
712 END IF
713*
714 RETURN
715*
716* End of DLATTR
717*
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:69
subroutine dlabad(SMALL, LARGE)
DLABAD
Definition: dlabad.f:74
subroutine dlarnv(IDIST, ISEED, N, X)
DLARNV returns a vector of random numbers from a uniform or normal distribution.
Definition: dlarnv.f:97
integer function idamax(N, DX, INCX)
IDAMAX
Definition: idamax.f:71
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:53
subroutine dcopy(N, DX, INCX, DY, INCY)
DCOPY
Definition: dcopy.f:82
subroutine drot(N, DX, INCX, DY, INCY, C, S)
DROT
Definition: drot.f:92
subroutine dscal(N, DA, DX, INCX)
DSCAL
Definition: dscal.f:79
subroutine dswap(N, DX, INCX, DY, INCY)
DSWAP
Definition: dswap.f:82
subroutine dlatb4(PATH, IMAT, M, N, TYPE, KL, KU, ANORM, MODE, CNDNUM, DIST)
DLATB4
Definition: dlatb4.f:120
double precision function dlarnd(IDIST, ISEED)
DLARND
Definition: dlarnd.f:73
subroutine dlatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
DLATMS
Definition: dlatms.f:321
subroutine drotg(a, b, c, s)
DROTG
Definition: drotg.f90:93
Here is the call graph for this function:
Here is the caller graph for this function: