1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911
|
PROGRAM PDLLTDRIVER
*
* -- ScaLAPACK testing driver (version 1.7) --
* University of Tennessee, Knoxville, Oak Ridge National Laboratory,
* and University of California, Berkeley.
* May 1, 1997
*
* Purpose
* =======
*
* PDLLTDRIVER is the main test program for the DOUBLE PRECISION
* ScaLAPACK Cholesky routines. This test driver performs an
* A = L*L**T or A = U**T*U factorization and solve, and optionally
* performs condition estimation and iterative refinement.
*
* The program must be driven by a short data file. An annotated
* example of a data file can be obtained by deleting the first 3
* characters from the following 18 lines:
* 'ScaLAPACK LLt factorization input file'
* 'Intel iPSC/860 hypercube, gamma model.'
* 'LLT.out' output file name (if any)
* 6 device out
* 'U' define Lower or Upper
* 1 number of problems sizes
* 31 100 200 values of N
* 1 number of NB's
* 2 10 24 values of NB
* 1 number of NRHS's
* 1 values of NRHS
* 1 Number of NBRHS's
* 1 values of NBRHS
* 1 number of process grids (ordered pairs of P & Q)
* 2 values of P
* 2 values of Q
* 1.0 threshold
* T (T or F) Test Cond. Est. and Iter. Ref. Routines
*
*
* Internal Parameters
* ===================
*
* TOTMEM INTEGER, default = 2000000
* TOTMEM is a machine-specific parameter indicating the
* maximum amount of available memory in bytes.
* The user should customize TOTMEM to his platform. Remember
* to leave room in memory for the operating system, the BLACS
* buffer, etc. For example, on a system with 8 MB of memory
* per process (e.g., one processor on an Intel iPSC/860), the
* parameters we use are TOTMEM=6200000 (leaving 1.8 MB for OS,
* code, BLACS buffer, etc). However, for PVM, we usually set
* TOTMEM = 2000000. Some experimenting with the maximum value
* of TOTMEM may be required.
*
* INTGSZ INTEGER, default = 4 bytes.
* DBLESZ INTEGER, default = 8 bytes.
* INTGSZ and DBLESZ indicate the length in bytes on the
* given platform for an integer and a double precision real.
* MEM DOUBLE PRECISION array, dimension ( TOTMEM / DBLESZ )
*
* All arrays used by SCALAPACK routines are allocated from
* this array and referenced by pointers. The integer IPA,
* for example, is a pointer to the starting element of MEM for
* the matrix A.
*
* =====================================================================
*
* .. Parameters ..
INTEGER BLOCK_CYCLIC_2D, CSRC_, CTXT_, DLEN_, DTYPE_,
$ LLD_, MB_, M_, NB_, N_, RSRC_
PARAMETER ( BLOCK_CYCLIC_2D = 1, DLEN_ = 9, DTYPE_ = 1,
$ CTXT_ = 2, M_ = 3, N_ = 4, MB_ = 5, NB_ = 6,
$ RSRC_ = 7, CSRC_ = 8, LLD_ = 9 )
INTEGER DBLESZ, INTGSZ, MEMSIZ, NTESTS, TOTMEM
DOUBLE PRECISION PADVAL, ZERO
PARAMETER ( DBLESZ = 8, TOTMEM = 2000000,
$ MEMSIZ = TOTMEM / DBLESZ, NTESTS = 20,
$ PADVAL = -9923.0D+0, ZERO = 0.0D+0 )
#ifdef TEST_INT64
PARAMETER ( INTGSZ = 8 )
#else
PARAMETER ( INTGSZ = 4 )
#endif
* ..
* .. Local Scalars ..
LOGICAL CHECK, EST
CHARACTER UPLO
CHARACTER*6 PASSED
CHARACTER*80 OUTFILE
INTEGER HH, I, IAM, IASEED, IBSEED, ICTXT, IMIDPAD,
$ INFO, IPA, IPA0, IPB, IPB0, IPBERR, IPFERR,
$ IPREPAD, IPOSTPAD, IPW, IPW2, ITEMP, J, K,
$ KFAIL, KK, KPASS, KSKIP, KTESTS, LCM, LCMQ,
$ LIWORK, LWORK, LW2, MYCOL, MYRHS, MYROW, N, NB,
$ NBRHS, NGRIDS, NMAT, NNB, NNBR, NNR, NOUT, NP,
$ NPCOL, NPROCS, NPROW, NQ, NRHS, WORKSIZ
REAL THRESH
DOUBLE PRECISION ANORM, ANORM1, FRESID, NOPS, RCOND,
$ SRESID, SRESID2, TMFLOPS
* ..
* .. Local Arrays ..
INTEGER DESCA( DLEN_ ), DESCB( DLEN_ ), IERR( 1 ),
$ NBRVAL( NTESTS ), NBVAL( NTESTS ),
$ NRVAL( NTESTS ), NVAL( NTESTS ),
$ PVAL( NTESTS ), QVAL( NTESTS )
DOUBLE PRECISION CTIME( 2 ), MEM( MEMSIZ ), WTIME( 2 )
* ..
* .. External Subroutines ..
EXTERNAL BLACS_BARRIER, BLACS_EXIT, BLACS_GRIDEXIT,
$ BLACS_GRIDINFO, BLACS_GRIDINIT, DESCINIT,
$ IGSUM2D, BLACS_PINFO, PDCHEKPAD, PDFILLPAD,
$ PDLAFCHK, PDLASCHK, PDLLTINFO,
$ PDMATGEN, PDPOCON, PDPORFS,
$ PDPOTRF, PDPOTRRV, PDPOTRS, SLBOOT,
$ SLCOMBINE, SLTIMER
* ..
* .. External Functions ..
LOGICAL LSAME
INTEGER ICEIL, ILCM, NUMROC
DOUBLE PRECISION PDLANSY
EXTERNAL ICEIL, ILCM, LSAME, NUMROC, PDLANSY
* ..
* .. Intrinsic Functions ..
INTRINSIC DBLE, MAX, MIN
* ..
* .. Data Statements ..
DATA KFAIL, KPASS, KSKIP, KTESTS / 4*0 /
* ..
* .. Executable Statements ..
*
* Get starting information
*
CALL BLACS_PINFO( IAM, NPROCS )
IASEED = 100
IBSEED = 200
CALL PDLLTINFO( OUTFILE, NOUT, UPLO, NMAT, NVAL, NTESTS, NNB,
$ NBVAL, NTESTS, NNR, NRVAL, NTESTS, NNBR, NBRVAL,
$ NTESTS, NGRIDS, PVAL, NTESTS, QVAL, NTESTS,
$ THRESH, EST, MEM, IAM, NPROCS )
CHECK = ( THRESH.GE.0.0E+0 )
*
* Print headings
*
IF( IAM.EQ.0 ) THEN
WRITE( NOUT, FMT = * )
WRITE( NOUT, FMT = 9995 )
WRITE( NOUT, FMT = 9994 )
WRITE( NOUT, FMT = * )
END IF
*
* Loop over different process grids
*
DO 50 I = 1, NGRIDS
*
NPROW = PVAL( I )
NPCOL = QVAL( I )
*
* Make sure grid information is correct
*
IERR( 1 ) = 0
IF( NPROW.LT.1 ) THEN
IF( IAM.EQ.0 )
$ WRITE( NOUT, FMT = 9999 ) 'GRID', 'nprow', NPROW
IERR( 1 ) = 1
ELSE IF( NPCOL.LT.1 ) THEN
IF( IAM.EQ.0 )
$ WRITE( NOUT, FMT = 9999 ) 'GRID', 'npcol', NPCOL
IERR( 1 ) = 1
ELSE IF( NPROW*NPCOL.GT.NPROCS ) THEN
IF( IAM.EQ.0 )
$ WRITE( NOUT, FMT = 9998 ) NPROW*NPCOL, NPROCS
IERR( 1 ) = 1
END IF
*
IF( IERR( 1 ).GT.0 ) THEN
IF( IAM.EQ.0 )
$ WRITE( NOUT, FMT = 9997 ) 'grid'
KSKIP = KSKIP + 1
GO TO 50
END IF
*
* Define process grid
*
CALL BLACS_GET( -1, 0, ICTXT )
CALL BLACS_GRIDINIT( ICTXT, 'Row-major', NPROW, NPCOL )
CALL BLACS_GRIDINFO( ICTXT, NPROW, NPCOL, MYROW, MYCOL )
*
* Go to bottom of process grid loop if this case doesn't use my
* process
*
IF( MYROW.GE.NPROW .OR. MYCOL.GE.NPCOL )
$ GO TO 50
*
DO 40 J = 1, NMAT
*
N = NVAL( J )
*
* Make sure matrix information is correct
*
IERR( 1 ) = 0
IF( N.LT.1 ) THEN
IF( IAM.EQ.0 )
$ WRITE( NOUT, FMT = 9999 ) 'MATRIX', 'N', N
IERR( 1 ) = 1
ELSE IF( N.LT.1 ) THEN
IF( IAM.EQ.0 )
$ WRITE( NOUT, FMT = 9999 ) 'MATRIX', 'N', N
IERR( 1 ) = 1
END IF
*
* Check all processes for an error
*
CALL IGSUM2D( ICTXT, 'All', ' ', 1, 1, IERR, 1, -1, 0 )
*
IF( IERR( 1 ).GT.0 ) THEN
IF( IAM.EQ.0 )
$ WRITE( NOUT, FMT = 9997 ) 'matrix'
KSKIP = KSKIP + 1
GO TO 40
END IF
*
DO 30 K = 1, NNB
*
NB = NBVAL( K )
*
* Make sure nb is legal
*
IERR( 1 ) = 0
IF( NB.LT.1 ) THEN
IERR( 1 ) = 1
IF( IAM.EQ.0 )
$ WRITE( NOUT, FMT = 9999 ) 'NB', 'NB', NB
END IF
*
* Check all processes for an error
*
CALL IGSUM2D( ICTXT, 'All', ' ', 1, 1, IERR, 1, -1, 0 )
*
IF( IERR( 1 ).GT.0 ) THEN
IF( IAM.EQ.0 )
$ WRITE( NOUT, FMT = 9997 ) 'NB'
KSKIP = KSKIP + 1
GO TO 30
END IF
*
* Padding constants
*
NP = NUMROC( N, NB, MYROW, 0, NPROW )
NQ = NUMROC( N, NB, MYCOL, 0, NPCOL )
IF( CHECK ) THEN
IPREPAD = MAX( NB, NP )
IMIDPAD = NB
IPOSTPAD = MAX( NB, NQ )
ELSE
IPREPAD = 0
IMIDPAD = 0
IPOSTPAD = 0
END IF
*
* Initialize the array descriptor for the matrix A
*
CALL DESCINIT( DESCA, N, N, NB, NB, 0, 0, ICTXT,
$ MAX( 1, NP )+IMIDPAD, IERR( 1 ) )
*
* Check all processes for an error
*
CALL IGSUM2D( ICTXT, 'All', ' ', 1, 1, IERR, 1, -1, 0 )
*
IF( IERR( 1 ).LT.0 ) THEN
IF( IAM.EQ.0 )
$ WRITE( NOUT, FMT = 9997 ) 'descriptor'
KSKIP = KSKIP + 1
GO TO 30
END IF
*
* Assign pointers into MEM for SCALAPACK arrays, A is
* allocated starting at position MEM( IPREPAD+1 )
*
IPA = IPREPAD+1
IF( EST ) THEN
IPA0 = IPA + DESCA( LLD_ )*NQ + IPOSTPAD + IPREPAD
IPW = IPA0 + DESCA( LLD_ )*NQ + IPOSTPAD + IPREPAD
ELSE
IPW = IPA + DESCA( LLD_ )*NQ + IPOSTPAD + IPREPAD
END IF
*
*
IF( CHECK ) THEN
*
* Calculate the amount of workspace required by
* the checking routines PDLAFCHK, PDPOTRRV, and
* PDLANSY
*
WORKSIZ = NP * DESCA( NB_ )
*
WORKSIZ = MAX( WORKSIZ, DESCA( MB_ ) * DESCA( NB_ ) )
*
LCM = ILCM( NPROW, NPCOL )
ITEMP = MAX( 2, 2 * NQ ) + NP
IF( NPROW.NE.NPCOL ) THEN
ITEMP = ITEMP +
$ NB * ICEIL( ICEIL( NP, NB ), LCM / NPROW )
END IF
WORKSIZ = MAX( WORKSIZ, ITEMP )
WORKSIZ = WORKSIZ + IPOSTPAD
*
ELSE
*
WORKSIZ = IPOSTPAD
*
END IF
*
* Check for adequate memory for problem size
*
IERR( 1 ) = 0
IF( IPW+WORKSIZ.GT.MEMSIZ ) THEN
IF( IAM.EQ.0 )
$ WRITE( NOUT, FMT = 9996 ) 'factorization',
$ ( IPW+WORKSIZ )*DBLESZ
IERR( 1 ) = 1
END IF
*
* Check all processes for an error
*
CALL IGSUM2D( ICTXT, 'All', ' ', 1, 1, IERR, 1, -1, 0 )
*
IF( IERR( 1 ).GT.0 ) THEN
IF( IAM.EQ.0 )
$ WRITE( NOUT, FMT = 9997 ) 'MEMORY'
KSKIP = KSKIP + 1
GO TO 30
END IF
*
* Generate a symmetric positive definite matrix A
*
CALL PDMATGEN( ICTXT, 'Symm', 'Diag', DESCA( M_ ),
$ DESCA( N_ ), DESCA( MB_ ), DESCA( NB_ ),
$ MEM( IPA ), DESCA( LLD_ ), DESCA( RSRC_ ),
$ DESCA( CSRC_ ), IASEED, 0, NP, 0, NQ,
$ MYROW, MYCOL, NPROW, NPCOL )
*
* Calculate inf-norm of A for residual error-checking
*
IF( CHECK ) THEN
CALL PDFILLPAD( ICTXT, NP, NQ, MEM( IPA-IPREPAD ),
$ DESCA( LLD_ ), IPREPAD, IPOSTPAD,
$ PADVAL )
CALL PDFILLPAD( ICTXT, WORKSIZ-IPOSTPAD, 1,
$ MEM( IPW-IPREPAD ), WORKSIZ-IPOSTPAD,
$ IPREPAD, IPOSTPAD, PADVAL )
ANORM = PDLANSY( 'I', UPLO, N, MEM( IPA ), 1, 1,
$ DESCA, MEM( IPW ) )
ANORM1 = PDLANSY( '1', UPLO, N, MEM( IPA ), 1, 1,
$ DESCA, MEM( IPW ) )
CALL PDCHEKPAD( ICTXT, 'PDLANSY', NP, NQ,
$ MEM( IPA-IPREPAD ), DESCA( LLD_ ),
$ IPREPAD, IPOSTPAD, PADVAL )
CALL PDCHEKPAD( ICTXT, 'PDLANSY',
$ WORKSIZ-IPOSTPAD, 1,
$ MEM( IPW-IPREPAD ), WORKSIZ-IPOSTPAD,
$ IPREPAD, IPOSTPAD, PADVAL )
END IF
*
IF( EST ) THEN
CALL PDMATGEN( ICTXT, 'Symm', 'Diag', DESCA( M_ ),
$ DESCA( N_ ), DESCA( MB_ ),
$ DESCA( NB_ ), MEM( IPA0 ),
$ DESCA( LLD_ ), DESCA( RSRC_ ),
$ DESCA( CSRC_ ), IASEED, 0, NP, 0, NQ,
$ MYROW, MYCOL, NPROW, NPCOL )
IF( CHECK )
$ CALL PDFILLPAD( ICTXT, NP, NQ,
$ MEM( IPA0-IPREPAD ), DESCA( LLD_ ),
$ IPREPAD, IPOSTPAD, PADVAL )
END IF
*
CALL SLBOOT()
CALL BLACS_BARRIER( ICTXT, 'All' )
*
* Perform LLt factorization
*
CALL SLTIMER( 1 )
*
CALL PDPOTRF( UPLO, N, MEM( IPA ), 1, 1, DESCA, INFO )
*
CALL SLTIMER( 1 )
*
IF( INFO.NE.0 ) THEN
IF( IAM.EQ.0 )
$ WRITE( NOUT, FMT = * ) 'PDPOTRF INFO=', INFO
KFAIL = KFAIL + 1
RCOND = ZERO
GO TO 60
END IF
*
IF( CHECK ) THEN
*
* Check for memory overwrite in LLt factorization
*
CALL PDCHEKPAD( ICTXT, 'PDPOTRF', NP, NQ,
$ MEM( IPA-IPREPAD ), DESCA( LLD_ ),
$ IPREPAD, IPOSTPAD, PADVAL )
END IF
*
IF( EST ) THEN
*
* Calculate workspace required for PDPOCON
*
LWORK = MAX( 1, 2*NP ) + MAX( 1, 2*NQ ) +
$ MAX( 2, DESCA( NB_ )*
$ MAX( 1, ICEIL( NPROW-1, NPCOL ) ),
$ NQ + DESCA( NB_ )*
$ MAX( 1, ICEIL( NPCOL-1, NPROW ) ) )
IPW2 = IPW + LWORK + IPOSTPAD + IPREPAD
LIWORK = MAX( 1, NP )
LW2 = ICEIL( LIWORK*INTGSZ, DBLESZ ) + IPOSTPAD
*
IERR( 1 ) = 0
IF( IPW2+LW2.GT.MEMSIZ ) THEN
IF( IAM.EQ.0 )
$ WRITE( NOUT, FMT = 9996 )'cond est',
$ ( IPW2+LW2 )*DBLESZ
IERR( 1 ) = 1
END IF
*
* Check all processes for an error
*
CALL IGSUM2D( ICTXT, 'All', ' ', 1, 1, IERR, 1,
$ -1, 0 )
*
IF( IERR( 1 ).GT.0 ) THEN
IF( IAM.EQ.0 )
$ WRITE( NOUT, FMT = 9997 ) 'MEMORY'
KSKIP = KSKIP + 1
GO TO 60
END IF
*
IF( CHECK ) THEN
CALL PDFILLPAD( ICTXT, LWORK, 1,
$ MEM( IPW-IPREPAD ), LWORK,
$ IPREPAD, IPOSTPAD, PADVAL )
CALL PDFILLPAD( ICTXT, LW2-IPOSTPAD, 1,
$ MEM( IPW2-IPREPAD ),
$ LW2-IPOSTPAD, IPREPAD,
$ IPOSTPAD, PADVAL )
END IF
*
* Compute condition number of the matrix
*
CALL PDPOCON( UPLO, N, MEM( IPA ), 1, 1, DESCA,
$ ANORM1, RCOND, MEM( IPW ), LWORK,
$ MEM( IPW2 ), LIWORK, INFO )
*
IF( CHECK ) THEN
CALL PDCHEKPAD( ICTXT, 'PDPOCON', NP, NQ,
$ MEM( IPA-IPREPAD ), DESCA( LLD_ ),
$ IPREPAD, IPOSTPAD, PADVAL )
CALL PDCHEKPAD( ICTXT, 'PDPOCON',
$ LWORK, 1, MEM( IPW-IPREPAD ),
$ LWORK, IPREPAD, IPOSTPAD,
$ PADVAL )
CALL PDCHEKPAD( ICTXT, 'PDPOCON',
$ LW2-IPOSTPAD, 1,
$ MEM( IPW2-IPREPAD ), LW2-IPOSTPAD,
$ IPREPAD, IPOSTPAD, PADVAL )
END IF
END IF
*
* Loop over the different values for NRHS
*
DO 20 HH = 1, NNR
*
NRHS = NRVAL( HH )
*
DO 10 KK = 1, NNBR
*
NBRHS = NBRVAL( KK )
*
* Initialize Array Descriptor for rhs
*
CALL DESCINIT( DESCB, N, NRHS, NB, NBRHS, 0, 0,
$ ICTXT, MAX( 1, NP )+IMIDPAD,
$ IERR( 1 ) )
*
* move IPW to allow room for RHS
*
MYRHS = NUMROC( DESCB( N_ ), DESCB( NB_ ), MYCOL,
$ DESCB( CSRC_ ), NPCOL )
IPB = IPW
*
IF( EST ) THEN
IPB0 = IPB + DESCB( LLD_ )*MYRHS + IPOSTPAD +
$ IPREPAD
IPFERR = IPB0 + DESCB( LLD_ )*MYRHS + IPOSTPAD
$ + IPREPAD
IPBERR = MYRHS + IPFERR + IPOSTPAD + IPREPAD
IPW = MYRHS + IPBERR + IPOSTPAD + IPREPAD
ELSE
IPW = IPB + DESCB( LLD_ )*MYRHS + IPOSTPAD +
$ IPREPAD
END IF
*
IF( CHECK ) THEN
*
* Calculate the amount of workspace required by
* the checking routines PDLASCHK
*
LCMQ = LCM / NPCOL
WORKSIZ = MAX( WORKSIZ-IPOSTPAD,
$ NQ * NBRHS + NP * NBRHS +
$ MAX( MAX( NQ*NB, 2*NBRHS ),
$ NBRHS * NUMROC( NUMROC(N,NB,0,0,NPCOL),NB,
$ 0,0,LCMQ ) ) )
WORKSIZ = IPOSTPAD + WORKSIZ
ELSE
WORKSIZ = IPOSTPAD
END IF
*
IERR( 1 ) = 0
IF( IPW+WORKSIZ.GT.MEMSIZ ) THEN
IF( IAM.EQ.0 )
$ WRITE( NOUT, FMT = 9996 )'solve',
$ ( IPW+WORKSIZ )*DBLESZ
IERR( 1 ) = 1
END IF
*
* Check all processes for an error
*
CALL IGSUM2D( ICTXT, 'All', ' ', 1, 1, IERR, 1,
$ -1, 0 )
*
IF( IERR( 1 ).GT.0 ) THEN
IF( IAM.EQ.0 )
$ WRITE( NOUT, FMT = 9997 ) 'MEMORY'
KSKIP = KSKIP + 1
GO TO 10
END IF
*
* Generate RHS
*
CALL PDMATGEN( ICTXT, 'No', 'No', DESCB( M_ ),
$ DESCB( N_ ), DESCB( MB_ ),
$ DESCB( NB_ ), MEM( IPB ),
$ DESCB( LLD_ ), DESCB( RSRC_ ),
$ DESCB( CSRC_ ), IBSEED, 0, NP, 0,
$ MYRHS, MYROW, MYCOL, NPROW, NPCOL )
*
IF( CHECK )
$ CALL PDFILLPAD( ICTXT, NP, MYRHS,
$ MEM( IPB-IPREPAD ),
$ DESCB( LLD_ ),
$ IPREPAD, IPOSTPAD, PADVAL )
*
IF( EST ) THEN
CALL PDMATGEN( ICTXT, 'No', 'No', DESCB( M_ ),
$ DESCB( N_ ), DESCB( MB_ ),
$ DESCB( NB_ ), MEM( IPB0 ),
$ DESCB( LLD_ ), DESCB( RSRC_ ),
$ DESCB( CSRC_ ), IBSEED, 0, NP, 0,
$ MYRHS, MYROW, MYCOL, NPROW,
$ NPCOL )
*
IF( CHECK ) THEN
CALL PDFILLPAD( ICTXT, NP, MYRHS,
$ MEM( IPB0-IPREPAD ),
$ DESCB( LLD_ ), IPREPAD,
$ IPOSTPAD, PADVAL )
CALL PDFILLPAD( ICTXT, 1, MYRHS,
$ MEM( IPFERR-IPREPAD ), 1,
$ IPREPAD, IPOSTPAD,
$ PADVAL )
CALL PDFILLPAD( ICTXT, 1, MYRHS,
$ MEM( IPBERR-IPREPAD ), 1,
$ IPREPAD, IPOSTPAD,
$ PADVAL )
END IF
END IF
*
CALL BLACS_BARRIER( ICTXT, 'All' )
CALL SLTIMER( 2 )
*
* Solve linear system via Cholesky factorization
*
CALL PDPOTRS( UPLO, N, NRHS, MEM( IPA ), 1, 1,
$ DESCA, MEM( IPB ), 1, 1, DESCB,
$ INFO )
*
CALL SLTIMER( 2 )
*
IF( CHECK ) THEN
*
* check for memory overwrite
*
CALL PDCHEKPAD( ICTXT, 'PDPOTRS', NP, NQ,
$ MEM( IPA-IPREPAD ),
$ DESCA( LLD_ ),
$ IPREPAD, IPOSTPAD, PADVAL )
CALL PDCHEKPAD( ICTXT, 'PDPOTRS', NP,
$ MYRHS, MEM( IPB-IPREPAD ),
$ DESCB( LLD_ ), IPREPAD,
$ IPOSTPAD, PADVAL )
*
CALL PDFILLPAD( ICTXT, WORKSIZ-IPOSTPAD, 1,
$ MEM( IPW-IPREPAD ),
$ WORKSIZ-IPOSTPAD, IPREPAD,
$ IPOSTPAD, PADVAL )
*
* check the solution to rhs
*
CALL PDLASCHK( 'Symm', 'Diag', N, NRHS,
$ MEM( IPB ), 1, 1, DESCB,
$ IASEED, 1, 1, DESCA, IBSEED,
$ ANORM, SRESID, MEM( IPW ) )
*
IF( IAM.EQ.0 .AND. SRESID.GT.THRESH )
$ WRITE( NOUT, FMT = 9985 ) SRESID
*
* check for memory overwrite
*
CALL PDCHEKPAD( ICTXT, 'PDLASCHK', NP,
$ MYRHS, MEM( IPB-IPREPAD ),
$ DESCB( LLD_ ), IPREPAD,
$ IPOSTPAD, PADVAL )
CALL PDCHEKPAD( ICTXT, 'PDLASCHK',
$ WORKSIZ-IPOSTPAD, 1,
$ MEM( IPW-IPREPAD ),
$ WORKSIZ-IPOSTPAD, IPREPAD,
$ IPOSTPAD, PADVAL )
*
* The second test is a NaN trap
*
IF( ( SRESID.LE.THRESH ).AND.
$ ( (SRESID-SRESID).EQ.0.0D+0 ) ) THEN
KPASS = KPASS + 1
PASSED = 'PASSED'
ELSE
KFAIL = KFAIL + 1
PASSED = 'FAILED'
END IF
ELSE
KPASS = KPASS + 1
SRESID = SRESID - SRESID
PASSED = 'BYPASS'
END IF
*
IF( EST ) THEN
*
* Calculate workspace required for PDPORFS
*
LWORK = MAX( 1, 3*NP )
IPW2 = IPW + LWORK + IPOSTPAD + IPREPAD
LIWORK = MAX( 1, NP )
LW2 = ICEIL( LIWORK*INTGSZ, DBLESZ ) +
$ IPOSTPAD
*
IERR( 1 ) = 0
IF( IPW2+LW2.GT.MEMSIZ ) THEN
IF( IAM.EQ.0 )
$ WRITE( NOUT, FMT = 9996 )
$ 'iter ref', ( IPW2+LW2 )*DBLESZ
IERR( 1 ) = 1
END IF
*
* Check all processes for an error
*
CALL IGSUM2D( ICTXT, 'All', ' ', 1, 1, IERR,
$ 1, -1, 0 )
*
IF( IERR( 1 ).GT.0 ) THEN
IF( IAM.EQ.0 )
$ WRITE( NOUT, FMT = 9997 )
$ 'MEMORY'
KSKIP = KSKIP + 1
GO TO 10
END IF
*
IF( CHECK ) THEN
CALL PDFILLPAD( ICTXT, LWORK, 1,
$ MEM( IPW-IPREPAD ),
$ LWORK, IPREPAD, IPOSTPAD,
$ PADVAL )
CALL PDFILLPAD( ICTXT, LW2-IPOSTPAD,
$ 1, MEM( IPW2-IPREPAD ),
$ LW2-IPOSTPAD,
$ IPREPAD, IPOSTPAD,
$ PADVAL )
END IF
*
* Use iterative refinement to improve the
* computed solution
*
CALL PDPORFS( UPLO, N, NRHS, MEM( IPA0 ),
$ 1, 1, DESCA, MEM( IPA ), 1, 1,
$ DESCA, MEM( IPB0 ), 1, 1,
$ DESCB, MEM( IPB ), 1, 1, DESCB,
$ MEM( IPFERR ), MEM( IPBERR ),
$ MEM( IPW ), LWORK, MEM( IPW2 ),
$ LIWORK, INFO )
*
* check for memory overwrite
*
IF( CHECK ) THEN
CALL PDCHEKPAD( ICTXT, 'PDPORFS', NP,
$ NQ, MEM( IPA0-IPREPAD ),
$ DESCA( LLD_ ), IPREPAD,
$ IPOSTPAD, PADVAL )
CALL PDCHEKPAD( ICTXT, 'PDPORFS', NP,
$ NQ, MEM( IPA-IPREPAD ),
$ DESCA( LLD_ ), IPREPAD,
$ IPOSTPAD, PADVAL )
CALL PDCHEKPAD( ICTXT, 'PDPORFS', NP,
$ MYRHS, MEM( IPB-IPREPAD ),
$ DESCB( LLD_ ), IPREPAD,
$ IPOSTPAD, PADVAL )
CALL PDCHEKPAD( ICTXT, 'PDPORFS', NP,
$ MYRHS,
$ MEM( IPB0-IPREPAD ),
$ DESCB( LLD_ ), IPREPAD,
$ IPOSTPAD, PADVAL )
CALL PDCHEKPAD( ICTXT, 'PDPORFS', 1,
$ MYRHS,
$ MEM( IPFERR-IPREPAD ), 1,
$ IPREPAD, IPOSTPAD,
$ PADVAL )
CALL PDCHEKPAD( ICTXT, 'PDPORFS', 1,
$ MYRHS,
$ MEM( IPBERR-IPREPAD ), 1,
$ IPREPAD, IPOSTPAD,
$ PADVAL )
CALL PDCHEKPAD( ICTXT, 'PDPORFS', LWORK,
$ 1, MEM( IPW-IPREPAD ),
$ LWORK, IPREPAD, IPOSTPAD,
$ PADVAL )
CALL PDCHEKPAD( ICTXT, 'PDPORFS',
$ LW2-IPOSTPAD, 1,
$ MEM( IPW2-IPREPAD ),
$ LW2-IPOSTPAD,
$ IPREPAD, IPOSTPAD,
$ PADVAL )
*
CALL PDFILLPAD( ICTXT, WORKSIZ-IPOSTPAD,
$ 1, MEM( IPW-IPREPAD ),
$ WORKSIZ-IPOSTPAD, IPREPAD,
$ IPOSTPAD, PADVAL )
*
* check the solution to rhs
*
CALL PDLASCHK( 'Symm', 'Diag', N, NRHS,
$ MEM( IPB ), 1, 1, DESCB,
$ IASEED, 1, 1, DESCA,
$ IBSEED, ANORM, SRESID2,
$ MEM( IPW ) )
*
IF( IAM.EQ.0 .AND. SRESID2.GT.THRESH )
$ WRITE( NOUT, FMT = 9985 ) SRESID2
*
* check for memory overwrite
*
CALL PDCHEKPAD( ICTXT, 'PDLASCHK', NP,
$ MYRHS, MEM( IPB-IPREPAD ),
$ DESCB( LLD_ ), IPREPAD,
$ IPOSTPAD, PADVAL )
CALL PDCHEKPAD( ICTXT, 'PDLASCHK',
$ WORKSIZ-IPOSTPAD, 1,
$ MEM( IPW-IPREPAD ),
$ WORKSIZ-IPOSTPAD,
$ IPREPAD, IPOSTPAD,
$ PADVAL )
END IF
END IF
*
* Gather maximum of all CPU and WALL clock timings
*
CALL SLCOMBINE( ICTXT, 'All', '>', 'W', 2, 1,
$ WTIME )
CALL SLCOMBINE( ICTXT, 'All', '>', 'C', 2, 1,
$ CTIME )
*
* Print results
*
IF( MYROW.EQ.0 .AND. MYCOL.EQ.0 ) THEN
*
* 1/3 N^3 + 1/2 N^2 flops for LLt factorization
*
NOPS = (DBLE(N)**3)/3.0D+0 +
$ (DBLE(N)**2)/2.0D+0
*
* nrhs * 2 N^2 flops for LLt solve.
*
NOPS = NOPS + 2.0D+0*(DBLE(N)**2)*DBLE(NRHS)
*
* Calculate total megaflops -- factorization and
* solve -- for WALL and CPU time, and print output
*
* Print WALL time if machine supports it
*
IF( WTIME( 1 ) + WTIME( 2 ) .GT. 0.0D+0 ) THEN
TMFLOPS = NOPS /
$ ( ( WTIME( 1 )+WTIME( 2 ) ) * 1.0D+6 )
ELSE
TMFLOPS = 0.0D+0
END IF
*
IF( WTIME( 2 ).GE.0.0D+0 )
$ WRITE( NOUT, FMT = 9993 ) 'WALL', UPLO, N,
$ NB, NRHS, NBRHS, NPROW, NPCOL,
$ WTIME( 1 ), WTIME( 2 ), TMFLOPS,
$ PASSED
*
* Print CPU time if machine supports it
*
IF( CTIME( 1 )+CTIME( 2 ).GT.0.0D+0 ) THEN
TMFLOPS = NOPS /
$ ( ( CTIME( 1 )+CTIME( 2 ) ) * 1.0D+6 )
ELSE
TMFLOPS = 0.0D+0
END IF
*
IF( CTIME( 2 ).GE.0.0D+0 )
$ WRITE( NOUT, FMT = 9993 ) 'CPU ', UPLO, N,
$ NB, NRHS, NBRHS, NPROW, NPCOL,
$ CTIME( 1 ), CTIME( 2 ), TMFLOPS,
$ PASSED
*
END IF
10 CONTINUE
20 CONTINUE
*
IF( CHECK .AND. SRESID.GT.THRESH ) THEN
*
* Compute FRESID = ||A - LL'|| / (||A|| * N * eps)
*
CALL PDPOTRRV( UPLO, N, MEM( IPA ), 1, 1, DESCA,
$ MEM( IPW ) )
CALL PDLAFCHK( 'Symm', 'Diag', N, N, MEM( IPA ), 1, 1,
$ DESCA, IASEED, ANORM, FRESID,
$ MEM( IPW ) )
*
* Check for memory overwrite
*
CALL PDCHEKPAD( ICTXT, 'PDPOTRRV', NP, NQ,
$ MEM( IPA-IPREPAD ), DESCA( LLD_ ),
$ IPREPAD, IPOSTPAD, PADVAL )
CALL PDCHEKPAD( ICTXT, 'PDGETRRV',
$ WORKSIZ-IPOSTPAD, 1,
$ MEM( IPW-IPREPAD ), WORKSIZ-IPOSTPAD,
$ IPREPAD, IPOSTPAD, PADVAL )
*
IF( IAM.EQ.0 ) THEN
IF( LSAME( UPLO, 'L' ) ) THEN
WRITE( NOUT, FMT = 9986 ) 'L*L''', FRESID
ELSE
WRITE( NOUT, FMT = 9986 ) 'U''*U', FRESID
END IF
END IF
END IF
*
30 CONTINUE
40 CONTINUE
CALL BLACS_GRIDEXIT( ICTXT )
50 CONTINUE
*
* Print ending messages and close output file
*
60 CONTINUE
IF( IAM.EQ.0 ) THEN
KTESTS = KPASS + KFAIL + KSKIP
WRITE( NOUT, FMT = * )
WRITE( NOUT, FMT = 9992 ) KTESTS
IF( CHECK ) THEN
WRITE( NOUT, FMT = 9991 ) KPASS
WRITE( NOUT, FMT = 9989 ) KFAIL
ELSE
WRITE( NOUT, FMT = 9990 ) KPASS
END IF
WRITE( NOUT, FMT = 9988 ) KSKIP
WRITE( NOUT, FMT = * )
WRITE( NOUT, FMT = * )
WRITE( NOUT, FMT = 9987 )
IF( NOUT.NE.6 .AND. NOUT.NE.0 )
$ CLOSE ( NOUT )
END IF
*
CALL BLACS_EXIT( 0 )
*
9999 FORMAT( 'ILLEGAL ', A6, ': ', A5, ' = ', I3,
$ '; It should be at least 1' )
9998 FORMAT( 'ILLEGAL GRID: nprow*npcol = ', I4, '. It can be at most',
$ I4 )
9997 FORMAT( 'Bad ', A6, ' parameters: going on to next test case.' )
9996 FORMAT( 'Unable to perform ', A, ': need TOTMEM of at least',
$ I11 )
9995 FORMAT( 'TIME UPLO N NB NRHS NBRHS P Q LLt Time ',
$ 'Slv Time MFLOPS CHECK' )
9994 FORMAT( '---- ---- ----- --- ---- ----- ---- ---- -------- ',
$ '-------- -------- ------' )
9993 FORMAT( A4, 4X, A1, 1X, I5, 1X, I3, 1X, I4, 1X, I5, 1X, I4, 1X,
$ I4, 1X, F8.2, 1X, F8.2, 1X, F8.2, 1X, A6 )
9992 FORMAT( 'Finished ', I6, ' tests, with the following results:' )
9991 FORMAT( I5, ' tests completed and passed residual checks.' )
9990 FORMAT( I5, ' tests completed without checking.' )
9989 FORMAT( I5, ' tests completed and failed residual checks.' )
9988 FORMAT( I5, ' tests skipped because of illegal input values.' )
9987 FORMAT( 'END OF TESTS.' )
9986 FORMAT( '||A - ', A4, '|| / (||A|| * N * eps) = ', G25.7 )
9985 FORMAT( '||Ax-b||/(||x||*||A||*eps*N) ', F25.7 )
*
STOP
*
* End of PDLLTDRIVER
*
END
|