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 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 954 955 956 957 958 959 960 961 962 963 964 965 966 967 968 969 970 971 972 973 974 975 976 977 978 979 980 981 982 983 984 985 986 987 988 989 990 991 992 993 994 995 996 997 998 999 1000 1001 1002 1003 1004 1005 1006 1007 1008 1009 1010 1011 1012 1013 1014 1015 1016 1017 1018 1019 1020 1021 1022 1023 1024 1025 1026 1027 1028 1029 1030 1031 1032 1033 1034 1035 1036 1037 1038 1039 1040 1041 1042 1043 1044 1045 1046 1047 1048 1049 1050 1051 1052 1053 1054 1055 1056 1057 1058 1059 1060 1061 1062 1063 1064 1065 1066 1067 1068 1069 1070 1071 1072 1073 1074 1075 1076 1077 1078 1079 1080 1081 1082 1083 1084 1085 1086 1087 1088 1089 1090 1091 1092 1093 1094 1095 1096 1097 1098 1099 1100 1101 1102 1103 1104 1105 1106 1107 1108 1109 1110 1111 1112 1113 1114 1115 1116 1117 1118 1119 1120 1121 1122 1123 1124 1125 1126 1127 1128 1129 1130 1131 1132 1133 1134 1135 1136 1137 1138 1139 1140 1141 1142 1143 1144 1145 1146 1147 1148 1149 1150 1151 1152 1153 1154 1155 1156 1157 1158 1159 1160 1161 1162 1163 1164 1165 1166 1167 1168 1169 1170 1171 1172 1173 1174 1175 1176 1177 1178 1179 1180 1181 1182 1183 1184 1185 1186 1187 1188 1189 1190 1191 1192 1193 1194 1195 1196 1197 1198 1199 1200 1201 1202 1203 1204 1205 1206 1207 1208 1209 1210 1211 1212 1213 1214 1215 1216 1217 1218 1219 1220 1221 1222 1223 1224 1225 1226 1227 1228 1229 1230 1231 1232 1233 1234 1235 1236 1237 1238 1239 1240 1241 1242 1243 1244 1245 1246 1247 1248 1249 1250 1251
|
SUBROUTINE PDDBTRF( N, BWL, BWU, A, JA, DESCA, AF, LAF, WORK,
$ LWORK, INFO )
*
* -- ScaLAPACK routine (version 2.0.2) --
* Univ. of Tennessee, Univ. of California Berkeley, Univ. of Colorado Denver
* May 1 2012
*
* .. Scalar Arguments ..
INTEGER BWL, BWU, INFO, JA, LAF, LWORK, N
* ..
* .. Array Arguments ..
INTEGER DESCA( * )
DOUBLE PRECISION A( * ), AF( * ), WORK( * )
* ..
*
* Purpose
* =======
*
* PDDBTRF computes a LU factorization
* of an N-by-N real banded
* diagonally dominant-like distributed matrix
* with bandwidth BWL, BWU: A(1:N, JA:JA+N-1).
* Reordering is used to increase parallelism in the factorization.
* This reordering results in factors that are DIFFERENT from those
* produced by equivalent sequential codes. These factors cannot
* be used directly by users; however, they can be used in
* subsequent calls to PDDBTRS to solve linear systems.
*
* The factorization has the form
*
* P A(1:N, JA:JA+N-1) P^T = L U
*
* where U is a banded upper triangular matrix and L is banded
* lower triangular, and P is a permutation matrix.
*
* =====================================================================
*
* Arguments
* =========
*
*
* N (global input) INTEGER
* The number of rows and columns to be operated on, i.e. the
* order of the distributed submatrix A(1:N, JA:JA+N-1). N >= 0.
*
* BWL (global input) INTEGER
* Number of subdiagonals. 0 <= BWL <= N-1
*
* BWU (global input) INTEGER
* Number of superdiagonals. 0 <= BWU <= N-1
*
* A (local input/local output) DOUBLE PRECISION pointer into
* local memory to an array with first dimension
* LLD_A >=(bwl+bwu+1) (stored in DESCA).
* On entry, this array contains the local pieces of the
* N-by-N unsymmetric banded distributed matrix
* A(1:N, JA:JA+N-1) to be factored.
* This local portion is stored in the packed banded format
* used in LAPACK. Please see the Notes below and the
* ScaLAPACK manual for more detail on the format of
* distributed matrices.
* On exit, this array contains information containing details
* of the factorization.
* Note that permutations are performed on the matrix, so that
* the factors returned are different from those returned
* by LAPACK.
*
* JA (global input) INTEGER
* The index in the global array A that points to the start of
* the matrix to be operated on (which may be either all of A
* or a submatrix of A).
*
* DESCA (global and local input) INTEGER array of dimension DLEN.
* if 1D type (DTYPE_A=501), DLEN >= 7;
* if 2D type (DTYPE_A=1), DLEN >= 9 .
* The array descriptor for the distributed matrix A.
* Contains information of mapping of A to memory. Please
* see NOTES below for full description and options.
*
* AF (local output) DOUBLE PRECISION array, dimension LAF.
* Auxiliary Fillin Space.
* Fillin is created during the factorization routine
* PDDBTRF and this is stored in AF. If a linear system
* is to be solved using PDDBTRS after the factorization
* routine, AF *must not be altered* after the factorization.
*
* LAF (local input) INTEGER
* Size of user-input Auxiliary Fillin space AF. Must be >=
* NB*(bwl+bwu)+6*max(bwl,bwu)*max(bwl,bwu)
* If LAF is not large enough, an error code will be returned
* and the minimum acceptable size will be returned in AF( 1 )
*
* WORK (local workspace/local output)
* DOUBLE PRECISION temporary workspace. This space may
* be overwritten in between calls to routines. WORK must be
* the size given in LWORK.
* On exit, WORK( 1 ) contains the minimal LWORK.
*
* LWORK (local input or global input) INTEGER
* Size of user-input workspace WORK.
* If LWORK is too small, the minimal acceptable size will be
* returned in WORK(1) and an error code is returned. LWORK>=
* max(bwl,bwu)*max(bwl,bwu)
*
* INFO (global output) INTEGER
* = 0: successful exit
* < 0: If the i-th argument is an array and the j-entry had
* an illegal value, then INFO = -(i*100+j), if the i-th
* argument is a scalar and had an illegal value, then
* INFO = -i.
* > 0: If INFO = K<=NPROCS, the submatrix stored on processor
* INFO and factored locally was not
* diagonally dominant-like, and
* the factorization was not completed.
* If INFO = K>NPROCS, the submatrix stored on processor
* INFO-NPROCS representing interactions with other
* processors was not
* stably factorable wo/interchanges,
* and the factorization was not completed.
*
* =====================================================================
*
* Restrictions
* ============
*
* The following are restrictions on the input parameters. Some of these
* are temporary and will be removed in future releases, while others
* may reflect fundamental technical limitations.
*
* Non-cyclic restriction: VERY IMPORTANT!
* P*NB>= mod(JA-1,NB)+N.
* The mapping for matrices must be blocked, reflecting the nature
* of the divide and conquer algorithm as a task-parallel algorithm.
* This formula in words is: no processor may have more than one
* chunk of the matrix.
*
* Blocksize cannot be too small:
* If the matrix spans more than one processor, the following
* restriction on NB, the size of each block on each processor,
* must hold:
* NB >= 2*MAX(BWL,BWU)
* The bulk of parallel computation is done on the matrix of size
* O(NB) on each processor. If this is too small, divide and conquer
* is a poor choice of algorithm.
*
* Submatrix reference:
* JA = IB
* Alignment restriction that prevents unnecessary communication.
*
*
* =====================================================================
*
* Notes
* =====
*
* If the factorization routine and the solve routine are to be called
* separately (to solve various sets of righthand sides using the same
* coefficient matrix), the auxiliary space AF *must not be altered*
* between calls to the factorization routine and the solve routine.
*
* The best algorithm for solving banded and tridiagonal linear systems
* depends on a variety of parameters, especially the bandwidth.
* Currently, only algorithms designed for the case N/P >> bw are
* implemented. These go by many names, including Divide and Conquer,
* Partitioning, domain decomposition-type, etc.
*
* Algorithm description: Divide and Conquer
*
* The Divide and Conqer algorithm assumes the matrix is narrowly
* banded compared with the number of equations. In this situation,
* it is best to distribute the input matrix A one-dimensionally,
* with columns atomic and rows divided amongst the processes.
* The basic algorithm divides the banded matrix up into
* P pieces with one stored on each processor,
* and then proceeds in 2 phases for the factorization or 3 for the
* solution of a linear system.
* 1) Local Phase:
* The individual pieces are factored independently and in
* parallel. These factors are applied to the matrix creating
* fillin, which is stored in a non-inspectable way in auxiliary
* space AF. Mathematically, this is equivalent to reordering
* the matrix A as P A P^T and then factoring the principal
* leading submatrix of size equal to the sum of the sizes of
* the matrices factored on each processor. The factors of
* these submatrices overwrite the corresponding parts of A
* in memory.
* 2) Reduced System Phase:
* A small (max(bwl,bwu)* (P-1)) system is formed representing
* interaction of the larger blocks, and is stored (as are its
* factors) in the space AF. A parallel Block Cyclic Reduction
* algorithm is used. For a linear system, a parallel front solve
* followed by an analagous backsolve, both using the structure
* of the factored matrix, are performed.
* 3) Backsubsitution Phase:
* For a linear system, a local backsubstitution is performed on
* each processor in parallel.
*
* Descriptors
* ===========
*
* Descriptors now have *types* and differ from ScaLAPACK 1.0.
*
* Note: banded codes can use either the old two dimensional
* or new one-dimensional descriptors, though the processor grid in
* both cases *must be one-dimensional*. We describe both types below.
*
* Each global data object is described by an associated description
* vector. This vector stores the information required to establish
* the mapping between an object element and its corresponding process
* and memory location.
*
* Let A be a generic term for any 2D block cyclicly distributed array.
* Such a global array has an associated description vector DESCA.
* In the following comments, the character _ should be read as
* "of the global array".
*
* NOTATION STORED IN EXPLANATION
* --------------- -------------- --------------------------------------
* DTYPE_A(global) DESCA( DTYPE_ )The descriptor type. In this case,
* DTYPE_A = 1.
* CTXT_A (global) DESCA( CTXT_ ) The BLACS context handle, indicating
* the BLACS process grid A is distribu-
* ted over. The context itself is glo-
* bal, but the handle (the integer
* value) may vary.
* M_A (global) DESCA( M_ ) The number of rows in the global
* array A.
* N_A (global) DESCA( N_ ) The number of columns in the global
* array A.
* MB_A (global) DESCA( MB_ ) The blocking factor used to distribute
* the rows of the array.
* NB_A (global) DESCA( NB_ ) The blocking factor used to distribute
* the columns of the array.
* RSRC_A (global) DESCA( RSRC_ ) The process row over which the first
* row of the array A is distributed.
* CSRC_A (global) DESCA( CSRC_ ) The process column over which the
* first column of the array A is
* distributed.
* LLD_A (local) DESCA( LLD_ ) The leading dimension of the local
* array. LLD_A >= MAX(1,LOCr(M_A)).
*
* Let K be the number of rows or columns of a distributed matrix,
* and assume that its process grid has dimension p x q.
* LOCr( K ) denotes the number of elements of K that a process
* would receive if K were distributed over the p processes of its
* process column.
* Similarly, LOCc( K ) denotes the number of elements of K that a
* process would receive if K were distributed over the q processes of
* its process row.
* The values of LOCr() and LOCc() may be determined via a call to the
* ScaLAPACK tool function, NUMROC:
* LOCr( M ) = NUMROC( M, MB_A, MYROW, RSRC_A, NPROW ),
* LOCc( N ) = NUMROC( N, NB_A, MYCOL, CSRC_A, NPCOL ).
* An upper bound for these quantities may be computed by:
* LOCr( M ) <= ceil( ceil(M/MB_A)/NPROW )*MB_A
* LOCc( N ) <= ceil( ceil(N/NB_A)/NPCOL )*NB_A
*
*
* One-dimensional descriptors:
*
* One-dimensional descriptors are a new addition to ScaLAPACK since
* version 1.0. They simplify and shorten the descriptor for 1D
* arrays.
*
* Since ScaLAPACK supports two-dimensional arrays as the fundamental
* object, we allow 1D arrays to be distributed either over the
* first dimension of the array (as if the grid were P-by-1) or the
* 2nd dimension (as if the grid were 1-by-P). This choice is
* indicated by the descriptor type (501 or 502)
* as described below.
*
* IMPORTANT NOTE: the actual BLACS grid represented by the
* CTXT entry in the descriptor may be *either* P-by-1 or 1-by-P
* irrespective of which one-dimensional descriptor type
* (501 or 502) is input.
* This routine will interpret the grid properly either way.
* ScaLAPACK routines *do not support intercontext operations* so that
* the grid passed to a single ScaLAPACK routine *must be the same*
* for all array descriptors passed to that routine.
*
* NOTE: In all cases where 1D descriptors are used, 2D descriptors
* may also be used, since a one-dimensional array is a special case
* of a two-dimensional array with one dimension of size unity.
* The two-dimensional array used in this case *must* be of the
* proper orientation:
* If the appropriate one-dimensional descriptor is DTYPEA=501
* (1 by P type), then the two dimensional descriptor must
* have a CTXT value that refers to a 1 by P BLACS grid;
* If the appropriate one-dimensional descriptor is DTYPEA=502
* (P by 1 type), then the two dimensional descriptor must
* have a CTXT value that refers to a P by 1 BLACS grid.
*
*
* Summary of allowed descriptors, types, and BLACS grids:
* DTYPE 501 502 1 1
* BLACS grid 1xP or Px1 1xP or Px1 1xP Px1
* -----------------------------------------------------
* A OK NO OK NO
* B NO OK NO OK
*
* Let A be a generic term for any 1D block cyclicly distributed array.
* Such a global array has an associated description vector DESCA.
* In the following comments, the character _ should be read as
* "of the global array".
*
* NOTATION STORED IN EXPLANATION
* --------------- ---------- ------------------------------------------
* DTYPE_A(global) DESCA( 1 ) The descriptor type. For 1D grids,
* TYPE_A = 501: 1-by-P grid.
* TYPE_A = 502: P-by-1 grid.
* CTXT_A (global) DESCA( 2 ) The BLACS context handle, indicating
* the BLACS process grid A is distribu-
* ted over. The context itself is glo-
* bal, but the handle (the integer
* value) may vary.
* N_A (global) DESCA( 3 ) The size of the array dimension being
* distributed.
* NB_A (global) DESCA( 4 ) The blocking factor used to distribute
* the distributed dimension of the array.
* SRC_A (global) DESCA( 5 ) The process row or column over which the
* first row or column of the array
* is distributed.
* LLD_A (local) DESCA( 6 ) The leading dimension of the local array
* storing the local blocks of the distri-
* buted array A. Minimum value of LLD_A
* depends on TYPE_A.
* TYPE_A = 501: LLD_A >=
* size of undistributed dimension, 1.
* TYPE_A = 502: LLD_A >=NB_A, 1.
* Reserved DESCA( 7 ) Reserved for future use.
*
* =====================================================================
*
* Code Developer: Andrew J. Cleary, University of Tennessee.
* Current address: Lawrence Livermore National Labs.
* Last modified by: Peter Arbenz, Institute of Scientific Computing,
* ETH, Zurich.
*
* =====================================================================
*
* .. Parameters ..
DOUBLE PRECISION ONE
PARAMETER ( ONE = 1.0D+0 )
DOUBLE PRECISION ZERO
PARAMETER ( ZERO = 0.0D+0 )
INTEGER INT_ONE
PARAMETER ( INT_ONE = 1 )
INTEGER DESCMULT, BIGNUM
PARAMETER ( DESCMULT = 100, BIGNUM = DESCMULT*DESCMULT )
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 )
* ..
* .. Local Scalars ..
INTEGER COMM_PROC, CSRC, FIRST_PROC, I, I1, I2, ICTXT,
$ ICTXT_NEW, ICTXT_SAVE, IDUM3, JA_NEW, LAF_MIN,
$ LEVEL_DIST, LLDA, MAX_BW, MBW2, MYCOL, MYROW,
$ MY_NUM_COLS, NB, NEXT_TRI_SIZE_M,
$ NEXT_TRI_SIZE_N, NP, NPCOL, NPROW, NP_SAVE,
$ ODD_SIZE, OFST, PART_OFFSET, PART_SIZE,
$ PREV_TRI_SIZE_M, PREV_TRI_SIZE_N, RETURN_CODE,
$ STORE_N_A, UP_PREV_TRI_SIZE_M,
$ UP_PREV_TRI_SIZE_N, WORK_SIZE_MIN, WORK_U
* ..
* .. Local Arrays ..
INTEGER DESCA_1XP( 7 ), PARAM_CHECK( 9, 3 )
* ..
* .. External Subroutines ..
EXTERNAL BLACS_GRIDEXIT, BLACS_GRIDINFO, DAXPY, DDBTRF,
$ DESC_CONVERT, DGEMM, DGEMV, DGERV2D, DGESD2D,
$ DLAMOV, DLATCPY, DTBTRS, DTRMM, DTRRV2D,
$ DTRSD2D, GLOBCHK, IGAMX2D, IGEBR2D, IGEBS2D,
$ PXERBLA, RESHAPE
* ..
* .. External Functions ..
INTEGER NUMROC
EXTERNAL NUMROC
* ..
* .. Intrinsic Functions ..
INTRINSIC MAX, MIN, MOD
* ..
* .. Executable Statements ..
*
* Test the input parameters
*
INFO = 0
*
* Convert descriptor into standard form for easy access to
* parameters, check that grid is of right shape.
*
DESCA_1XP( 1 ) = 501
*
CALL DESC_CONVERT( DESCA, DESCA_1XP, RETURN_CODE )
*
IF( RETURN_CODE.NE.0 ) THEN
INFO = -( 6*100+2 )
END IF
*
* Get values out of descriptor for use in code.
*
ICTXT = DESCA_1XP( 2 )
CSRC = DESCA_1XP( 5 )
NB = DESCA_1XP( 4 )
LLDA = DESCA_1XP( 6 )
STORE_N_A = DESCA_1XP( 3 )
*
* Get grid parameters
*
*
* Size of separator blocks is maximum of bandwidths
*
MAX_BW = MAX( BWL, BWU )
MBW2 = MAX_BW*MAX_BW
*
CALL BLACS_GRIDINFO( ICTXT, NPROW, NPCOL, MYROW, MYCOL )
NP = NPROW*NPCOL
*
*
*
IF( LWORK.LT.-1 ) THEN
INFO = -10
ELSE IF( LWORK.EQ.-1 ) THEN
IDUM3 = -1
ELSE
IDUM3 = 1
END IF
*
IF( N.LT.0 ) THEN
INFO = -1
END IF
*
IF( N+JA-1.GT.STORE_N_A ) THEN
INFO = -( 6*100+6 )
END IF
*
IF( ( BWL.GT.N-1 ) .OR. ( BWL.LT.0 ) ) THEN
INFO = -2
END IF
*
IF( ( BWU.GT.N-1 ) .OR. ( BWU.LT.0 ) ) THEN
INFO = -3
END IF
*
IF( LLDA.LT.( BWL+BWU+1 ) ) THEN
INFO = -( 6*100+6 )
END IF
*
IF( NB.LE.0 ) THEN
INFO = -( 6*100+4 )
END IF
*
* Argument checking that is specific to Divide & Conquer routine
*
IF( NPROW.NE.1 ) THEN
INFO = -( 6*100+2 )
END IF
*
IF( N.GT.NP*NB-MOD( JA-1, NB ) ) THEN
INFO = -( 1 )
CALL PXERBLA( ICTXT, 'PDDBTRF, D&C alg.: only 1 block per proc'
$ , -INFO )
RETURN
END IF
*
IF( ( JA+N-1.GT.NB ) .AND. ( NB.LT.2*MAX( BWL, BWU ) ) ) THEN
INFO = -( 6*100+4 )
CALL PXERBLA( ICTXT, 'PDDBTRF, D&C alg.: NB too small', -INFO )
RETURN
END IF
*
*
* Check auxiliary storage size
*
LAF_MIN = NB*( BWL+BWU ) + 6*MAX( BWL, BWU )*MAX( BWL, BWU )
*
IF( LAF.LT.LAF_MIN ) THEN
INFO = -8
* put minimum value of laf into AF( 1 )
AF( 1 ) = LAF_MIN
CALL PXERBLA( ICTXT, 'PDDBTRF: auxiliary storage error ',
$ -INFO )
RETURN
END IF
*
* Check worksize
*
WORK_SIZE_MIN = MAX( BWL, BWU )*MAX( BWL, BWU )
*
WORK( 1 ) = WORK_SIZE_MIN
*
IF( LWORK.LT.WORK_SIZE_MIN ) THEN
IF( LWORK.NE.-1 ) THEN
INFO = -10
CALL PXERBLA( ICTXT, 'PDDBTRF: worksize error ', -INFO )
END IF
RETURN
END IF
*
* Pack params and positions into arrays for global consistency check
*
PARAM_CHECK( 9, 1 ) = DESCA( 5 )
PARAM_CHECK( 8, 1 ) = DESCA( 4 )
PARAM_CHECK( 7, 1 ) = DESCA( 3 )
PARAM_CHECK( 6, 1 ) = DESCA( 1 )
PARAM_CHECK( 5, 1 ) = JA
PARAM_CHECK( 4, 1 ) = BWU
PARAM_CHECK( 3, 1 ) = BWL
PARAM_CHECK( 2, 1 ) = N
PARAM_CHECK( 1, 1 ) = IDUM3
*
PARAM_CHECK( 9, 2 ) = 605
PARAM_CHECK( 8, 2 ) = 604
PARAM_CHECK( 7, 2 ) = 603
PARAM_CHECK( 6, 2 ) = 601
PARAM_CHECK( 5, 2 ) = 5
PARAM_CHECK( 4, 2 ) = 3
PARAM_CHECK( 3, 2 ) = 2
PARAM_CHECK( 2, 2 ) = 1
PARAM_CHECK( 1, 2 ) = 10
*
* Want to find errors with MIN( ), so if no error, set it to a big
* number. If there already is an error, multiply by the the
* descriptor multiplier.
*
IF( INFO.GE.0 ) THEN
INFO = BIGNUM
ELSE IF( INFO.LT.-DESCMULT ) THEN
INFO = -INFO
ELSE
INFO = -INFO*DESCMULT
END IF
*
* Check consistency across processors
*
CALL GLOBCHK( ICTXT, 9, PARAM_CHECK, 9, PARAM_CHECK( 1, 3 ),
$ INFO )
*
* Prepare output: set info = 0 if no error, and divide by DESCMULT
* if error is not in a descriptor entry.
*
IF( INFO.EQ.BIGNUM ) THEN
INFO = 0
ELSE IF( MOD( INFO, DESCMULT ).EQ.0 ) THEN
INFO = -INFO / DESCMULT
ELSE
INFO = -INFO
END IF
*
IF( INFO.LT.0 ) THEN
CALL PXERBLA( ICTXT, 'PDDBTRF', -INFO )
RETURN
END IF
*
* Quick return if possible
*
IF( N.EQ.0 )
$ RETURN
*
*
* Adjust addressing into matrix space to properly get into
* the beginning part of the relevant data
*
PART_OFFSET = NB*( ( JA-1 ) / ( NPCOL*NB ) )
*
IF( ( MYCOL-CSRC ).LT.( JA-PART_OFFSET-1 ) / NB ) THEN
PART_OFFSET = PART_OFFSET + NB
END IF
*
IF( MYCOL.LT.CSRC ) THEN
PART_OFFSET = PART_OFFSET - NB
END IF
*
* Form a new BLACS grid (the "standard form" grid) with only procs
* holding part of the matrix, of size 1xNP where NP is adjusted,
* starting at csrc=0, with JA modified to reflect dropped procs.
*
* First processor to hold part of the matrix:
*
FIRST_PROC = MOD( ( JA-1 ) / NB+CSRC, NPCOL )
*
* Calculate new JA one while dropping off unused processors.
*
JA_NEW = MOD( JA-1, NB ) + 1
*
* Save and compute new value of NP
*
NP_SAVE = NP
NP = ( JA_NEW+N-2 ) / NB + 1
*
* Call utility routine that forms "standard-form" grid
*
CALL RESHAPE( ICTXT, INT_ONE, ICTXT_NEW, INT_ONE, FIRST_PROC,
$ INT_ONE, NP )
*
* Use new context from standard grid as context.
*
ICTXT_SAVE = ICTXT
ICTXT = ICTXT_NEW
DESCA_1XP( 2 ) = ICTXT_NEW
*
* Get information about new grid.
*
CALL BLACS_GRIDINFO( ICTXT, NPROW, NPCOL, MYROW, MYCOL )
*
* Drop out processors that do not have part of the matrix.
*
IF( MYROW.LT.0 ) THEN
GO TO 140
END IF
*
* ********************************
* Values reused throughout routine
*
* User-input value of partition size
*
PART_SIZE = NB
*
* Number of columns in each processor
*
MY_NUM_COLS = NUMROC( N, PART_SIZE, MYCOL, 0, NPCOL )
*
* Offset in columns to beginning of main partition in each proc
*
IF( MYCOL.EQ.0 ) THEN
PART_OFFSET = PART_OFFSET + MOD( JA_NEW-1, PART_SIZE )
MY_NUM_COLS = MY_NUM_COLS - MOD( JA_NEW-1, PART_SIZE )
END IF
*
* Offset in elements
*
OFST = PART_OFFSET*LLDA
*
* Size of main (or odd) partition in each processor
*
ODD_SIZE = MY_NUM_COLS
IF( MYCOL.LT.NP-1 ) THEN
ODD_SIZE = ODD_SIZE - MAX_BW
END IF
*
* Offset to workspace for Upper triangular factor
*
WORK_U = BWU*ODD_SIZE + 3*MBW2
*
*
* Zero out space for fillin
*
DO 10 I = 1, LAF_MIN
AF( I ) = ZERO
10 CONTINUE
*
* Zero out space for work
*
DO 20 I = 1, WORK_SIZE_MIN
WORK( I ) = ZERO
20 CONTINUE
*
* Begin main code
*
*
********************************************************************
* PHASE 1: Local computation phase.
********************************************************************
*
*
* Sizes of the extra triangles communicated bewtween processors
*
IF( MYCOL.GT.0 ) THEN
PREV_TRI_SIZE_M = MIN( BWL, NUMROC( N, PART_SIZE, MYCOL, 0,
$ NPCOL ) )
PREV_TRI_SIZE_N = MIN( BWL, NUMROC( N, PART_SIZE, MYCOL-1, 0,
$ NPCOL ) )
END IF
*
IF( MYCOL.GT.0 ) THEN
UP_PREV_TRI_SIZE_M = MIN( BWU,
$ NUMROC( N, PART_SIZE, MYCOL, 0, NPCOL ) )
UP_PREV_TRI_SIZE_N = MIN( BWU,
$ NUMROC( N, PART_SIZE, MYCOL-1, 0,
$ NPCOL ) )
END IF
*
IF( MYCOL.LT.NPCOL-1 ) THEN
NEXT_TRI_SIZE_M = MIN( BWL, NUMROC( N, PART_SIZE, MYCOL+1, 0,
$ NPCOL ) )
NEXT_TRI_SIZE_N = MIN( BWL, NUMROC( N, PART_SIZE, MYCOL, 0,
$ NPCOL ) )
END IF
*
IF( MYCOL.LT.NP-1 ) THEN
* Transfer last triangle D_i of local matrix to next processor
* which needs it to calculate fillin due to factorization of
* its main (odd) block A_i.
* Overlap the send with the factorization of A_i.
*
CALL DTRSD2D( ICTXT, 'U', 'N', NEXT_TRI_SIZE_M,
$ NEXT_TRI_SIZE_N, A( OFST+( MY_NUM_COLS-BWL )*
$ LLDA+( BWL+BWU+1 ) ), LLDA-1, 0, MYCOL+1 )
*
END IF
*
*
* Factor main partition A_i = L_i {U_i} in each processor
*
CALL DDBTRF( ODD_SIZE, ODD_SIZE, BWL, BWU, A( OFST+1 ), LLDA,
$ INFO )
*
IF( INFO.NE.0 ) THEN
INFO = MYCOL + 1
GO TO 30
END IF
*
*
IF( MYCOL.LT.NP-1 ) THEN
*
* Apply factorization to lower connection block BL_i
* transpose the connection block in preparation.
* Apply factorization to upper connection block BU_i
* Move the connection block in preparation.
*
CALL DLATCPY( 'U', BWL, BWL, A( ( OFST+( BWL+BWU+1 )+
$ ( ODD_SIZE-BWL )*LLDA ) ), LLDA-1,
$ AF( ODD_SIZE*BWU+2*MBW2+1+MAX_BW-BWL ), MAX_BW )
CALL DLAMOV( 'L', BWU, BWU, A( ( OFST+1+ODD_SIZE*LLDA ) ),
$ LLDA-1, AF( WORK_U+ODD_SIZE*BWL+2*MBW2+1+MAX_BW-
$ BWU ), MAX_BW )
*
* Perform the triangular system solve {L_i}{{BU'}_i} = {B_i}
*
CALL DTBTRS( 'L', 'N', 'U', BWU, BWL, BWU,
$ A( OFST+BWU+1+( ODD_SIZE-BWU )*LLDA ), LLDA,
$ AF( WORK_U+ODD_SIZE*BWL+2*MBW2+1+MAX_BW-BWU ),
$ MAX_BW, INFO )
*
* Perform the triangular solve {U_i}^T{BL'}_i^T = {BL_i}^T
*
CALL DTBTRS( 'U', 'T', 'N', BWL, BWU, BWL,
$ A( OFST+1+( ODD_SIZE-BWL )*LLDA ), LLDA,
$ AF( ODD_SIZE*BWU+2*MBW2+1+MAX_BW-BWL ), MAX_BW,
$ INFO )
*
* transpose resulting block to its location
* in main storage.
*
CALL DLATCPY( 'L', BWL, BWL, AF( ODD_SIZE*BWU+2*MBW2+1+MAX_BW-
$ BWL ), MAX_BW, A( ( OFST+( BWL+BWU+1 )+
$ ( ODD_SIZE-BWL )*LLDA ) ), LLDA-1 )
*
* Move the resulting block back to its location in main storage.
*
CALL DLAMOV( 'L', BWU, BWU, AF( WORK_U+ODD_SIZE*BWL+2*MBW2+1+
$ MAX_BW-BWU ), MAX_BW, A( ( OFST+1+ODD_SIZE*
$ LLDA ) ), LLDA-1 )
*
*
* Compute contribution to diagonal block(s) of reduced system.
* {C'}_i = {C_i}-{{BL'}_i}{{BU'}_i}
*
* The following method uses more flops than necessary but
* does not necessitate the writing of a new BLAS routine.
*
*
CALL DGEMM( 'T', 'N', MAX_BW, MAX_BW, MAX_BW, -ONE,
$ AF( ODD_SIZE*BWU+2*MBW2+1 ), MAX_BW,
$ AF( WORK_U+ODD_SIZE*BWL+2*MBW2+1 ), MAX_BW, ONE,
$ A( OFST+ODD_SIZE*LLDA+1+BWU ), LLDA-1 )
*
END IF
* End of "if ( MYCOL .lt. NP-1 )..." loop
*
*
30 CONTINUE
* If the processor could not locally factor, it jumps here.
*
IF( MYCOL.NE.0 ) THEN
* Discard temporary matrix stored beginning in
* AF( (odd_size+2*bwl, bwu)*bwl, bwu+1 ) and use for
* off_diagonal block of reduced system.
*
* Receive previously transmitted matrix section, which forms
* the right-hand-side for the triangular solve that calculates
* the "spike" fillin.
*
*
CALL DTRRV2D( ICTXT, 'U', 'N', PREV_TRI_SIZE_M,
$ PREV_TRI_SIZE_N, AF( WORK_U+1 ), BWL, 0,
$ MYCOL-1 )
*
IF( INFO.EQ.0 ) THEN
*
* Calculate the "spike" fillin, ${L_i} {{GU}_i} = {DL_i}$ .
*
* Transpose transmitted triangular matrix $DL_i$
*
DO 50 I1 = 1, BWL
DO 40 I2 = I1 + 1, BWL
AF( WORK_U+I2+( I1-1 )*BWL ) = AF( WORK_U+I1+( I2-1 )*
$ BWL )
AF( WORK_U+I1+( I2-1 )*BWL ) = ZERO
40 CONTINUE
50 CONTINUE
*
DO 60 I1 = 2, ODD_SIZE
I2 = MIN( I1-1, BWL )
CALL DGEMV( 'N', BWL, I2, -ONE,
$ AF( WORK_U+1+( I1-1-I2 )*BWL ), BWL,
$ A( OFST+BWU+1+I2+( I1-1-I2 )*LLDA ), LLDA-1,
$ ONE, AF( WORK_U+1+( I1-1 )*BWL ), 1 )
60 CONTINUE
*
*
* Calculate the "spike" fillin, ${U_i}^T {{GL}_i}^T = {DU_i}^T$
*
*
* Copy D block into AF storage for solve.
*
CALL DLAMOV( 'L', UP_PREV_TRI_SIZE_N, UP_PREV_TRI_SIZE_M,
$ A( OFST+1 ), LLDA-1, AF( 1 ), BWU )
*
DO 80 I1 = 1, ODD_SIZE
I2 = MIN( BWU, I1-1 )
CALL DGEMV( 'N', BWU, I2, -ONE, AF( ( I1-1-I2 )*BWU+1 ),
$ BWU, A( OFST+BWU+1-I2+( I1-1 )*LLDA ), 1,
$ ONE, AF( ( I1-1 )*BWU+1 ), 1 )
*
DO 70 I = 1, BWU
AF( ( I1-1 )*BWU+I ) = AF( ( I1-1 )*BWU+I ) /
$ A( ( I1-1 )*LLDA+BWU+1 )
70 CONTINUE
80 CONTINUE
*
* Calculate the update block for previous proc, E_i = GL_i{GU_i}
*
*
* Zero out space in case result is smaller than storage block
*
DO 90 I = 1, MBW2
AF( ODD_SIZE*BWU+2*MBW2+I ) = ZERO
90 CONTINUE
*
CALL DGEMM( 'N', 'T', BWU, BWL, ODD_SIZE, -ONE, AF( 1 ),
$ BWU, AF( WORK_U+1 ), BWL, ZERO,
$ AF( 1+MAX( 0, BWL-BWU )+ODD_SIZE*BWU+( 2*MAX_BW+
$ MAX( 0, BWU-BWL ) )*MAX_BW ), MAX_BW )
*
*
* Initiate send of E_i to previous processor to overlap
* with next computation.
*
CALL DGESD2D( ICTXT, MAX_BW, MAX_BW,
$ AF( ODD_SIZE*BWU+2*MBW2+1 ), MAX_BW, 0,
$ MYCOL-1 )
*
*
IF( MYCOL.LT.NP-1 ) THEN
*
* Calculate off-diagonal block(s) of reduced system.
* Note: for ease of use in solution of reduced system, store
* L's off-diagonal block in transpose form.
*
* Copy matrix HU_i (the last bwl rows of GU_i) to AFL storage
* as per requirements of BLAS routine DTRMM.
* Since we have GU_i stored,
* transpose HU_i to HU_i^T.
*
CALL DLAMOV( 'N', BWL, BWL,
$ AF( WORK_U+( ODD_SIZE-BWL )*BWL+1 ), BWL,
$ AF( ( ODD_SIZE )*BWU+1+( MAX_BW-BWL ) ),
$ MAX_BW )
*
CALL DTRMM( 'R', 'U', 'T', 'N', BWL, BWL, -ONE,
$ A( ( OFST+( BWL+BWU+1 )+( ODD_SIZE-BWL )*
$ LLDA ) ), LLDA-1, AF( ( ODD_SIZE )*BWU+1+
$ ( MAX_BW-BWL ) ), MAX_BW )
*
*
* Copy matrix HL_i (the last bwu rows of GL_i^T) to AFU store
* as per requirements of BLAS routine DTRMM.
* Since we have GL_i^T stored,
* transpose HL_i^T to HL_i.
*
CALL DLAMOV( 'N', BWU, BWU, AF( ( ODD_SIZE-BWU )*BWU+1 ),
$ BWU, AF( WORK_U+( ODD_SIZE )*BWL+1+MAX_BW-
$ BWU ), MAX_BW )
*
CALL DTRMM( 'R', 'L', 'N', 'N', BWU, BWU, -ONE,
$ A( ( OFST+1+ODD_SIZE*LLDA ) ), LLDA-1,
$ AF( WORK_U+( ODD_SIZE )*BWL+1+MAX_BW-BWU ),
$ MAX_BW )
*
END IF
*
END IF
* End of "if ( MYCOL .ne. 0 )..."
*
END IF
* End of "if (info.eq.0) then"
*
*
* Check to make sure no processors have found errors
*
CALL IGAMX2D( ICTXT, 'A', ' ', 1, 1, INFO, 1, INFO, INFO, -1, 0,
$ 0 )
*
IF( MYCOL.EQ.0 ) THEN
CALL IGEBS2D( ICTXT, 'A', ' ', 1, 1, INFO, 1 )
ELSE
CALL IGEBR2D( ICTXT, 'A', ' ', 1, 1, INFO, 1, 0, 0 )
END IF
*
IF( INFO.NE.0 ) THEN
GO TO 130
END IF
* No errors found, continue
*
*
********************************************************************
* PHASE 2: Formation and factorization of Reduced System.
********************************************************************
*
* Gather up local sections of reduced system
*
*
* The last processor does not participate in the factorization of
* the reduced system, having sent its E_i already.
IF( MYCOL.EQ.NPCOL-1 ) THEN
GO TO 120
END IF
*
* Initiate send of off-diag block(s) to overlap with next part.
* Off-diagonal block needed on neighboring processor to start
* algorithm.
*
IF( ( MOD( MYCOL+1, 2 ).EQ.0 ) .AND. ( MYCOL.GT.0 ) ) THEN
*
CALL DGESD2D( ICTXT, MAX_BW, MAX_BW, AF( ODD_SIZE*BWU+1 ),
$ MAX_BW, 0, MYCOL-1 )
*
CALL DGESD2D( ICTXT, MAX_BW, MAX_BW,
$ AF( WORK_U+ODD_SIZE*BWL+1 ), MAX_BW, 0, MYCOL-1 )
*
END IF
*
* Copy last diagonal block into AF storage for subsequent
* operations.
*
CALL DLAMOV( 'N', MAX_BW, MAX_BW, A( OFST+ODD_SIZE*LLDA+BWU+1 ),
$ LLDA-1, AF( ODD_SIZE*BWU+MBW2+1 ), MAX_BW )
*
* Receive cont. to diagonal block that is stored on this proc.
*
IF( MYCOL.LT.NPCOL-1 ) THEN
*
CALL DGERV2D( ICTXT, MAX_BW, MAX_BW,
$ AF( ODD_SIZE*BWU+2*MBW2+1 ), MAX_BW, 0, MYCOL+1 )
*
* Add contribution to diagonal block
*
CALL DAXPY( MBW2, ONE, AF( ODD_SIZE*BWU+2*MBW2+1 ), 1,
$ AF( ODD_SIZE*BWU+MBW2+1 ), 1 )
*
END IF
*
*
* *************************************
* Modification Loop
*
* The distance for sending and receiving for each level starts
* at 1 for the first level.
LEVEL_DIST = 1
*
* Do until this proc is needed to modify other procs' equations
*
100 CONTINUE
IF( MOD( ( MYCOL+1 ) / LEVEL_DIST, 2 ).NE.0 )
$ GO TO 110
*
* Receive and add contribution to diagonal block from the left
*
IF( MYCOL-LEVEL_DIST.GE.0 ) THEN
CALL DGERV2D( ICTXT, MAX_BW, MAX_BW, WORK( 1 ), MAX_BW, 0,
$ MYCOL-LEVEL_DIST )
*
CALL DAXPY( MBW2, ONE, WORK( 1 ), 1, AF( ODD_SIZE*BWU+MBW2+1 ),
$ 1 )
*
END IF
*
* Receive and add contribution to diagonal block from the right
*
IF( MYCOL+LEVEL_DIST.LT.NPCOL-1 ) THEN
CALL DGERV2D( ICTXT, MAX_BW, MAX_BW, WORK( 1 ), MAX_BW, 0,
$ MYCOL+LEVEL_DIST )
*
CALL DAXPY( MBW2, ONE, WORK( 1 ), 1, AF( ODD_SIZE*BWU+MBW2+1 ),
$ 1 )
*
END IF
*
LEVEL_DIST = LEVEL_DIST*2
*
GO TO 100
110 CONTINUE
* [End of GOTO Loop]
*
*
* *********************************
* Calculate and use this proc's blocks to modify other procs'...
*
* Factor diagonal block
*
CALL DDBTRF( MAX_BW, MAX_BW, MIN( MAX_BW-1, BWL ),
$ MIN( MAX_BW-1, BWU ), AF( ODD_SIZE*BWU+MBW2+1-
$ ( MIN( MAX_BW-1, BWU ) ) ), MAX_BW+1, INFO )
*
IF( INFO.NE.0 ) THEN
INFO = NPCOL + MYCOL
END IF
*
* ****************************************************************
* Receive offdiagonal block from processor to right.
* If this is the first group of processors, the receive comes
* from a different processor than otherwise.
*
IF( LEVEL_DIST.EQ.1 ) THEN
COMM_PROC = MYCOL + 1
*
* Move block into place that it will be expected to be for
* calcs.
*
CALL DLAMOV( 'N', MAX_BW, MAX_BW, AF( ODD_SIZE*BWU+1 ), MAX_BW,
$ AF( WORK_U+ODD_SIZE*BWL+2*MBW2+1 ), MAX_BW )
*
CALL DLAMOV( 'N', MAX_BW, MAX_BW, AF( WORK_U+ODD_SIZE*BWL+1 ),
$ MAX_BW, AF( ODD_SIZE*BWU+2*MBW2+1 ), MAX_BW )
*
ELSE
COMM_PROC = MYCOL + LEVEL_DIST / 2
END IF
*
IF( MYCOL / LEVEL_DIST.LE.( NPCOL-1 ) / LEVEL_DIST-2 ) THEN
*
CALL DGERV2D( ICTXT, MAX_BW, MAX_BW, AF( ODD_SIZE*BWU+1 ),
$ MAX_BW, 0, COMM_PROC )
*
CALL DGERV2D( ICTXT, MAX_BW, MAX_BW,
$ AF( WORK_U+ODD_SIZE*BWL+1 ), MAX_BW, 0,
$ COMM_PROC )
*
IF( INFO.EQ.0 ) THEN
*
*
* Modify upper off_diagonal block with diagonal block
*
*
CALL DTBTRS( 'L', 'N', 'U', BWU, MIN( BWL, BWU-1 ), BWU,
$ AF( ODD_SIZE*BWU+MBW2+1+( MAX_BW+1 )*( MAX_BW-
$ BWU ) ), MAX_BW+1, AF( WORK_U+ODD_SIZE*BWL+1+
$ MAX_BW-BWU ), MAX_BW, INFO )
*
* Modify lower off_diagonal block with diagonal block
*
*
CALL DTBTRS( 'U', 'T', 'N', BWL, MIN( BWU, BWL-1 ), BWL,
$ AF( ODD_SIZE*BWU+MBW2+1-MIN( BWU,
$ BWL-1 )+( MAX_BW+1 )*( MAX_BW-BWL ) ),
$ MAX_BW+1, AF( ODD_SIZE*BWU+1+MAX_BW-BWL ),
$ MAX_BW, INFO )
*
END IF
* End of "if ( info.eq.0 ) then"
*
* Calculate contribution from this block to next diagonal block
*
CALL DGEMM( 'T', 'N', MAX_BW, MAX_BW, MAX_BW, -ONE,
$ AF( ( ODD_SIZE )*BWU+1 ), MAX_BW,
$ AF( WORK_U+( ODD_SIZE )*BWL+1 ), MAX_BW, ZERO,
$ WORK( 1 ), MAX_BW )
*
* Send contribution to diagonal block's owning processor.
*
CALL DGESD2D( ICTXT, MAX_BW, MAX_BW, WORK( 1 ), MAX_BW, 0,
$ MYCOL+LEVEL_DIST )
*
END IF
* End of "if( mycol/level_dist .le. (npcol-1)/level_dist-2 )..."
*
*
* ****************************************************************
* Receive off_diagonal block from left and use to finish with this
* processor.
*
IF( ( MYCOL / LEVEL_DIST.GT.0 ) .AND.
$ ( MYCOL / LEVEL_DIST.LE.( NPCOL-1 ) / LEVEL_DIST-1 ) ) THEN
*
IF( LEVEL_DIST.GT.1 ) THEN
*
* Receive offdiagonal block(s) from proc level_dist/2 to the
* left
*
CALL DGERV2D( ICTXT, MAX_BW, MAX_BW,
$ AF( WORK_U+ODD_SIZE*BWL+2*MBW2+1 ), MAX_BW, 0,
$ MYCOL-LEVEL_DIST / 2 )
*
* Receive offdiagonal block(s) from proc level_dist/2 to the
* left
*
CALL DGERV2D( ICTXT, MAX_BW, MAX_BW,
$ AF( ODD_SIZE*BWU+2*MBW2+1 ), MAX_BW, 0,
$ MYCOL-LEVEL_DIST / 2 )
*
END IF
*
*
IF( INFO.EQ.0 ) THEN
*
* Use diagonal block(s) to modify this offdiagonal block
*
*
* Since DTBTRS has no "left-right" option, we must transpose
*
CALL DLATCPY( 'N', MAX_BW, MAX_BW,
$ AF( WORK_U+ODD_SIZE*BWL+2*MBW2+1 ), MAX_BW,
$ WORK( 1 ), MAX_BW )
*
CALL DTBTRS( 'L', 'N', 'U', MAX_BW, MIN( BWL, MAX_BW-1 ),
$ BWL, AF( ODD_SIZE*BWU+MBW2+1 ), MAX_BW+1,
$ WORK( 1+MAX_BW*( MAX_BW-BWL ) ), MAX_BW, INFO )
*
* Transpose back
*
CALL DLATCPY( 'N', MAX_BW, MAX_BW, WORK( 1 ), MAX_BW,
$ AF( WORK_U+ODD_SIZE*BWL+2*MBW2+1 ), MAX_BW )
*
*
*
* Since DTBTRS has no "left-right" option, we must transpose
*
CALL DLATCPY( 'N', MAX_BW, MAX_BW,
$ AF( ODD_SIZE*BWU+2*MBW2+1 ), MAX_BW,
$ WORK( 1 ), MAX_BW )
*
CALL DTBTRS( 'U', 'T', 'N', MAX_BW, MIN( BWU, MAX_BW-1 ),
$ BWU, AF( ODD_SIZE*BWU+MBW2+1-MIN( BWU,
$ MAX_BW-1 ) ), MAX_BW+1,
$ WORK( 1+MAX_BW*( MAX_BW-BWU ) ), MAX_BW, INFO )
*
* Transpose back
*
CALL DLATCPY( 'N', MAX_BW, MAX_BW, WORK( 1 ), MAX_BW,
$ AF( ODD_SIZE*BWU+2*MBW2+1 ), MAX_BW )
*
*
END IF
* End of "if( info.eq.0 ) then"
*
* Use offdiag block(s) to calculate modification to diag block
* of processor to the left
*
CALL DGEMM( 'N', 'T', MAX_BW, MAX_BW, MAX_BW, -ONE,
$ AF( ( ODD_SIZE )*BWU+2*MBW2+1 ), MAX_BW,
$ AF( WORK_U+( ODD_SIZE )*BWL+2*MBW2+1 ), MAX_BW,
$ ZERO, WORK( 1 ), MAX_BW )
*
* Send contribution to diagonal block's owning processor.
*
CALL DGESD2D( ICTXT, MAX_BW, MAX_BW, WORK( 1 ), MAX_BW, 0,
$ MYCOL-LEVEL_DIST )
*
* *******************************************************
*
IF( MYCOL / LEVEL_DIST.LE.( NPCOL-1 ) / LEVEL_DIST-2 ) THEN
*
* Decide which processor offdiagonal block(s) goes to
*
IF( ( MOD( MYCOL / ( 2*LEVEL_DIST ), 2 ) ).EQ.0 ) THEN
COMM_PROC = MYCOL + LEVEL_DIST
ELSE
COMM_PROC = MYCOL - LEVEL_DIST
END IF
*
* Use offdiagonal blocks to calculate offdiag
* block to send to neighboring processor. Depending
* on circumstances, may need to transpose the matrix.
*
CALL DGEMM( 'N', 'N', MAX_BW, MAX_BW, MAX_BW, -ONE,
$ AF( WORK_U+ODD_SIZE*BWL+2*MBW2+1 ), MAX_BW,
$ AF( ODD_SIZE*BWU+1 ), MAX_BW, ZERO, WORK( 1 ),
$ MAX_BW )
*
* Send contribution to offdiagonal block's owning processor.
*
CALL DGESD2D( ICTXT, MAX_BW, MAX_BW, WORK( 1 ), MAX_BW, 0,
$ COMM_PROC )
*
CALL DGEMM( 'N', 'N', MAX_BW, MAX_BW, MAX_BW, -ONE,
$ AF( ODD_SIZE*BWU+2*MBW2+1 ), MAX_BW,
$ AF( WORK_U+ODD_SIZE*BWL+1 ), MAX_BW, ZERO,
$ WORK( 1 ), MAX_BW )
*
* Send contribution to offdiagonal block's owning processor.
*
CALL DGESD2D( ICTXT, MAX_BW, MAX_BW, WORK( 1 ), MAX_BW, 0,
$ COMM_PROC )
*
END IF
*
END IF
* End of "if( mycol/level_dist.le. (npcol-1)/level_dist -1 )..."
*
120 CONTINUE
*
*
130 CONTINUE
*
*
* Free BLACS space used to hold standard-form grid.
*
IF( ICTXT_SAVE.NE.ICTXT_NEW ) THEN
CALL BLACS_GRIDEXIT( ICTXT_NEW )
END IF
*
140 CONTINUE
*
* Restore saved input parameters
*
ICTXT = ICTXT_SAVE
NP = NP_SAVE
*
* Output minimum worksize
*
WORK( 1 ) = WORK_SIZE_MIN
*
* Make INFO consistent across processors
*
CALL IGAMX2D( ICTXT, 'A', ' ', 1, 1, INFO, 1, INFO, INFO, -1, 0,
$ 0 )
*
IF( MYCOL.EQ.0 ) THEN
CALL IGEBS2D( ICTXT, 'A', ' ', 1, 1, INFO, 1 )
ELSE
CALL IGEBR2D( ICTXT, 'A', ' ', 1, 1, INFO, 1, 0, 0 )
END IF
*
*
RETURN
*
* End of PDDBTRF
*
END
|