1 #define PETSCKSP_DLL 2 3 /**********************************ivec.c************************************** 4 5 Author: Henry M. Tufo III 6 7 e-mail: hmt@cs.brown.edu 8 9 snail-mail: 10 Division of Applied Mathematics 11 Brown University 12 Providence, RI 02912 13 14 Last Modification: 15 6.21.97 16 ***********************************ivec.c*************************************/ 17 18 /**********************************ivec.c************************************** 19 File Description: 20 ----------------- 21 22 ***********************************ivec.c*************************************/ 23 #include "src/ksp/pc/impls/tfs/tfs.h" 24 25 /* sorting args ivec.c ivec.c ... */ 26 #define SORT_OPT 6 27 #define SORT_STACK 50000 28 29 30 /* allocate an address and size stack for sorter(s) */ 31 static void *offset_stack[2*SORT_STACK]; 32 static int size_stack[SORT_STACK]; 33 static long psize_stack[SORT_STACK]; 34 35 36 37 /**********************************ivec.c************************************** 38 Function ivec_dump() 39 40 Input : 41 Output: 42 Return: 43 Description: 44 ***********************************ivec.c*************************************/ 45 PetscErrorCode ivec_dump(int *v, int n, int tag, int tag2, char * s) 46 { 47 int i; 48 PetscFunctionBegin; 49 printf("%2d %2d %s %2d :: ",tag,tag2,s,my_id); 50 for (i=0;i<n;i++) 51 {printf("%2d ",v[i]);} 52 printf("\n"); 53 fflush(stdout); 54 PetscFunctionReturn(0); 55 } 56 57 58 59 /**********************************ivec.c************************************** 60 Function ivec_lb_ub() 61 62 Input : 63 Output: 64 Return: 65 Description: 66 ***********************************ivec.c*************************************/ 67 PetscErrorCode ivec_lb_ub( int *arg1, int n, int *lb, int *ub) 68 { 69 int min = INT_MAX; 70 int max = INT_MIN; 71 72 PetscFunctionBegin; 73 while (n--) 74 { 75 min = PetscMin(min,*arg1); 76 max = PetscMax(max,*arg1); 77 arg1++; 78 } 79 80 *lb=min; 81 *ub=max; 82 PetscFunctionReturn(0); 83 } 84 85 86 87 /**********************************ivec.c************************************** 88 Function ivec_copy() 89 90 Input : 91 Output: 92 Return: 93 Description: 94 ***********************************ivec.c*************************************/ 95 int *ivec_copy( int *arg1, int *arg2, int n) 96 { 97 while (n--) {*arg1++ = *arg2++;} 98 return(arg1); 99 } 100 101 102 103 /**********************************ivec.c************************************** 104 Function ivec_zero() 105 106 Input : 107 Output: 108 Return: 109 Description: 110 ***********************************ivec.c*************************************/ 111 PetscErrorCode ivec_zero( int *arg1, int n) 112 { 113 PetscFunctionBegin; 114 while (n--) {*arg1++ = 0;} 115 PetscFunctionReturn(0); 116 } 117 118 119 120 /**********************************ivec.c************************************** 121 Function ivec_comp() 122 123 Input : 124 Output: 125 Return: 126 Description: 127 ***********************************ivec.c*************************************/ 128 PetscErrorCode ivec_comp( int *arg1, int n) 129 { 130 PetscFunctionBegin; 131 while (n--) {*arg1 = ~*arg1; arg1++;} 132 PetscFunctionReturn(0); 133 } 134 135 136 137 /**********************************ivec.c************************************** 138 Function ivec_neg_one() 139 140 Input : 141 Output: 142 Return: 143 Description: 144 ***********************************ivec.c*************************************/ 145 PetscErrorCode ivec_neg_one( int *arg1, int n) 146 { 147 PetscFunctionBegin; 148 while (n--) {*arg1++ = -1;} 149 PetscFunctionReturn(0); 150 } 151 152 153 154 /**********************************ivec.c************************************** 155 Function ivec_pos_one() 156 157 Input : 158 Output: 159 Return: 160 Description: 161 ***********************************ivec.c*************************************/ 162 PetscErrorCode ivec_pos_one( int *arg1, int n) 163 { 164 PetscFunctionBegin; 165 while (n--) {*arg1++ = 1;} 166 PetscFunctionReturn(0); 167 } 168 169 170 171 /**********************************ivec.c************************************** 172 Function ivec_c_index() 173 174 Input : 175 Output: 176 Return: 177 Description: 178 ***********************************ivec.c*************************************/ 179 PetscErrorCode ivec_c_index( int *arg1, int n) 180 { 181 int i=0; 182 183 PetscFunctionBegin; 184 while (n--) {*arg1++ = i++;} 185 PetscFunctionReturn(0); 186 } 187 188 189 /**********************************ivec.c************************************** 190 Function ivec_set() 191 192 Input : 193 Output: 194 Return: 195 Description: 196 ***********************************ivec.c*************************************/ 197 PetscErrorCode ivec_set( int *arg1, int arg2, int n) 198 { 199 PetscFunctionBegin; 200 while (n--) {*arg1++ = arg2;} 201 PetscFunctionReturn(0); 202 } 203 204 205 206 /**********************************ivec.c************************************** 207 Function ivec_cmp() 208 209 Input : 210 Output: 211 Return: 212 Description: 213 ***********************************ivec.c*************************************/ 214 int 215 ivec_cmp( int *arg1, int *arg2, int n) 216 { 217 PetscFunctionBegin; 218 while (n--) {if (*arg1++ != *arg2++) {return(FALSE);}} 219 return(TRUE); 220 } 221 222 223 224 /**********************************ivec.c************************************** 225 Function ivec_max() 226 227 Input : 228 Output: 229 Return: 230 Description: 231 ***********************************ivec.c*************************************/ 232 PetscErrorCode ivec_max( int *arg1, int *arg2, int n) 233 { 234 PetscFunctionBegin; 235 while (n--) {*arg1 = PetscMax(*arg1,*arg2); arg1++; arg2++;} 236 PetscFunctionReturn(0); 237 } 238 239 240 241 /**********************************ivec.c************************************** 242 Function ivec_min() 243 244 Input : 245 Output: 246 Return: 247 Description: 248 ***********************************ivec.c*************************************/ 249 PetscErrorCode ivec_min( int *arg1, int *arg2, int n) 250 { 251 PetscFunctionBegin; 252 while (n--) {*(arg1) = PetscMin(*arg1,*arg2); arg1++; arg2++;} 253 PetscFunctionReturn(0); 254 } 255 256 257 258 /**********************************ivec.c************************************** 259 Function ivec_mult() 260 261 Input : 262 Output: 263 Return: 264 Description: 265 ***********************************ivec.c*************************************/ 266 PetscErrorCode ivec_mult( int *arg1, int *arg2, int n) 267 { 268 PetscFunctionBegin; 269 while (n--) {*arg1++ *= *arg2++;} 270 PetscFunctionReturn(0); 271 } 272 273 274 275 /**********************************ivec.c************************************** 276 Function ivec_add() 277 278 Input : 279 Output: 280 Return: 281 Description: 282 ***********************************ivec.c*************************************/ 283 PetscErrorCode ivec_add( int *arg1, int *arg2, int n) 284 { 285 PetscFunctionBegin; 286 while (n--) {*arg1++ += *arg2++;} 287 PetscFunctionReturn(0); 288 } 289 290 291 292 /**********************************ivec.c************************************** 293 Function ivec_lxor() 294 295 Input : 296 Output: 297 Return: 298 Description: 299 ***********************************ivec.c*************************************/ 300 PetscErrorCode ivec_lxor( int *arg1, int *arg2, int n) 301 { 302 PetscFunctionBegin; 303 while (n--) {*arg1=((*arg1 || *arg2) && !(*arg1 && *arg2)); arg1++; arg2++;} 304 PetscFunctionReturn(0); 305 } 306 307 308 309 /**********************************ivec.c************************************** 310 Function ivec_xor() 311 312 Input : 313 Output: 314 Return: 315 Description: 316 ***********************************ivec.c*************************************/ 317 PetscErrorCode ivec_xor( int *arg1, int *arg2, int n) 318 { 319 PetscFunctionBegin; 320 while (n--) {*arg1++ ^= *arg2++;} 321 PetscFunctionReturn(0); 322 } 323 324 325 326 /**********************************ivec.c************************************** 327 Function ivec_or() 328 329 Input : 330 Output: 331 Return: 332 Description: 333 ***********************************ivec.c*************************************/ 334 PetscErrorCode ivec_or( int *arg1, int *arg2, int n) 335 { 336 PetscFunctionBegin; 337 while (n--) {*arg1++ |= *arg2++;} 338 PetscFunctionReturn(0); 339 } 340 341 342 343 /**********************************ivec.c************************************** 344 Function ivec_lor() 345 346 Input : 347 Output: 348 Return: 349 Description: 350 ***********************************ivec.c*************************************/ 351 PetscErrorCode ivec_lor( int *arg1, int *arg2, int n) 352 { 353 PetscFunctionBegin; 354 while (n--) {*arg1 = (*arg1 || *arg2); arg1++; arg2++;} 355 PetscFunctionReturn(0); 356 } 357 358 359 360 /**********************************ivec.c************************************** 361 Function ivec_or3() 362 363 Input : 364 Output: 365 Return: 366 Description: 367 ***********************************ivec.c*************************************/ 368 PetscErrorCode ivec_or3( int *arg1, int *arg2, int *arg3, 369 int n) 370 { 371 PetscFunctionBegin; 372 while (n--) {*arg1++ = (*arg2++ | *arg3++);} 373 PetscFunctionReturn(0); 374 } 375 376 377 378 /**********************************ivec.c************************************** 379 Function ivec_and() 380 381 Input : 382 Output: 383 Return: 384 Description: 385 ***********************************ivec.c*************************************/ 386 PetscErrorCode ivec_and( int *arg1, int *arg2, int n) 387 { 388 PetscFunctionBegin; 389 while (n--) {*arg1++ &= *arg2++;} 390 PetscFunctionReturn(0); 391 } 392 393 394 395 /**********************************ivec.c************************************** 396 Function ivec_land() 397 398 Input : 399 Output: 400 Return: 401 Description: 402 ***********************************ivec.c*************************************/ 403 PetscErrorCode ivec_land( int *arg1, int *arg2, int n) 404 { 405 PetscFunctionBegin; 406 while (n--) {*arg1 = (*arg1 && *arg2); arg1++; arg2++;} 407 PetscFunctionReturn(0); 408 } 409 410 411 412 /**********************************ivec.c************************************** 413 Function ivec_and3() 414 415 Input : 416 Output: 417 Return: 418 Description: 419 ***********************************ivec.c*************************************/ 420 PetscErrorCode ivec_and3( int *arg1, int *arg2, int *arg3, int n) 421 { 422 PetscFunctionBegin; 423 while (n--) {*arg1++ = (*arg2++ & *arg3++);} 424 PetscFunctionReturn(0); 425 } 426 427 428 429 /**********************************ivec.c************************************** 430 Function ivec_sum 431 432 Input : 433 Output: 434 Return: 435 Description: 436 ***********************************ivec.c*************************************/ 437 int ivec_sum( int *arg1, int n) 438 { 439 int tmp = 0; 440 441 442 while (n--) {tmp += *arg1++;} 443 return(tmp); 444 } 445 446 447 448 /**********************************ivec.c************************************** 449 Function ivec_reduce_and 450 451 Input : 452 Output: 453 Return: 454 Description: 455 ***********************************ivec.c*************************************/ 456 int ivec_reduce_and( int *arg1, int n) 457 { 458 int tmp = ALL_ONES; 459 460 461 while (n--) {tmp &= *arg1++;} 462 return(tmp); 463 } 464 465 466 467 /**********************************ivec.c************************************** 468 Function ivec_reduce_or 469 470 Input : 471 Output: 472 Return: 473 Description: 474 ***********************************ivec.c*************************************/ 475 int ivec_reduce_or( int *arg1, int n) 476 { 477 int tmp = 0; 478 479 480 while (n--) {tmp |= *arg1++;} 481 return(tmp); 482 } 483 484 485 486 /**********************************ivec.c************************************** 487 Function ivec_prod 488 489 Input : 490 Output: 491 Return: 492 Description: 493 ***********************************ivec.c*************************************/ 494 int ivec_prod( int *arg1, int n) 495 { 496 int tmp = 1; 497 498 499 while (n--) {tmp *= *arg1++;} 500 return(tmp); 501 } 502 503 504 505 /**********************************ivec.c************************************** 506 Function ivec_u_sum 507 508 Input : 509 Output: 510 Return: 511 Description: 512 ***********************************ivec.c*************************************/ 513 int ivec_u_sum( unsigned *arg1, int n) 514 { 515 unsigned tmp = 0; 516 517 518 while (n--) {tmp += *arg1++;} 519 return(tmp); 520 } 521 522 523 524 /**********************************ivec.c************************************** 525 Function ivec_lb() 526 527 Input : 528 Output: 529 Return: 530 Description: 531 ***********************************ivec.c*************************************/ 532 int 533 ivec_lb( int *arg1, int n) 534 { 535 int min = INT_MAX; 536 537 538 while (n--) {min = PetscMin(min,*arg1); arg1++;} 539 return(min); 540 } 541 542 543 544 /**********************************ivec.c************************************** 545 Function ivec_ub() 546 547 Input : 548 Output: 549 Return: 550 Description: 551 ***********************************ivec.c*************************************/ 552 int 553 ivec_ub( int *arg1, int n) 554 { 555 int max = INT_MIN; 556 557 558 while (n--) {max = PetscMax(max,*arg1); arg1++;} 559 return(max); 560 } 561 562 563 564 /**********************************ivec.c************************************** 565 Function split_buf() 566 567 Input : 568 Output: 569 Return: 570 Description: 571 572 assumes that sizeof(int) == 4bytes!!! 573 ***********************************ivec.c*************************************/ 574 int 575 ivec_split_buf(int *buf1, int **buf2, int size) 576 { 577 *buf2 = (buf1 + (size>>3)); 578 return(size); 579 } 580 581 582 583 /**********************************ivec.c************************************** 584 Function ivec_non_uniform() 585 586 Input : 587 Output: 588 Return: 589 Description: 590 ***********************************ivec.c*************************************/ 591 PetscErrorCode ivec_non_uniform(int *arg1, int *arg2, int n, int *arg3) 592 { 593 int i, j, type; 594 595 596 PetscFunctionBegin; 597 /* LATER: if we're really motivated we can sort and then unsort */ 598 for (i=0;i<n;) 599 { 600 /* clump 'em for now */ 601 j=i+1; 602 type = arg3[i]; 603 while ((j<n)&&(arg3[j]==type)) 604 {j++;} 605 606 /* how many together */ 607 j -= i; 608 609 /* call appropriate ivec function */ 610 if (type == GL_MAX) 611 {ivec_max(arg1,arg2,j);} 612 else if (type == GL_MIN) 613 {ivec_min(arg1,arg2,j);} 614 else if (type == GL_MULT) 615 {ivec_mult(arg1,arg2,j);} 616 else if (type == GL_ADD) 617 {ivec_add(arg1,arg2,j);} 618 else if (type == GL_B_XOR) 619 {ivec_xor(arg1,arg2,j);} 620 else if (type == GL_B_OR) 621 {ivec_or(arg1,arg2,j);} 622 else if (type == GL_B_AND) 623 {ivec_and(arg1,arg2,j);} 624 else if (type == GL_L_XOR) 625 {ivec_lxor(arg1,arg2,j);} 626 else if (type == GL_L_OR) 627 {ivec_lor(arg1,arg2,j);} 628 else if (type == GL_L_AND) 629 {ivec_land(arg1,arg2,j);} 630 else 631 {error_msg_fatal("unrecognized type passed to ivec_non_uniform()!");} 632 633 arg1+=j; arg2+=j; i+=j; 634 } 635 PetscFunctionReturn(0); 636 } 637 638 639 640 /**********************************ivec.c************************************** 641 Function ivec_addr() 642 643 Input : 644 Output: 645 Return: 646 Description: 647 ***********************************ivec.c*************************************/ 648 vfp ivec_fct_addr( int type) 649 { 650 PetscFunctionBegin; 651 if (type == NON_UNIFORM) 652 {return((PetscErrorCode (*)(void*, void *, int, ...))&ivec_non_uniform);} 653 else if (type == GL_MAX) 654 {return((PetscErrorCode (*)(void*, void *, int, ...))&ivec_max);} 655 else if (type == GL_MIN) 656 {return((PetscErrorCode (*)(void*, void *, int, ...))&ivec_min);} 657 else if (type == GL_MULT) 658 {return((PetscErrorCode (*)(void*, void *, int, ...))&ivec_mult);} 659 else if (type == GL_ADD) 660 {return((PetscErrorCode (*)(void*, void *, int, ...))&ivec_add);} 661 else if (type == GL_B_XOR) 662 {return((PetscErrorCode (*)(void*, void *, int, ...))&ivec_xor);} 663 else if (type == GL_B_OR) 664 {return((PetscErrorCode (*)(void*, void *, int, ...))&ivec_or);} 665 else if (type == GL_B_AND) 666 {return((PetscErrorCode (*)(void*, void *, int, ...))&ivec_and);} 667 else if (type == GL_L_XOR) 668 {return((PetscErrorCode (*)(void*, void *, int, ...))&ivec_lxor);} 669 else if (type == GL_L_OR) 670 {return((PetscErrorCode (*)(void*, void *, int, ...))&ivec_lor);} 671 else if (type == GL_L_AND) 672 {return((PetscErrorCode (*)(void*, void *, int, ...))&ivec_land);} 673 674 /* catch all ... not good if we get here */ 675 return(NULL); 676 } 677 678 679 /**********************************ivec.c************************************** 680 Function ct_bits() 681 682 Input : 683 Output: 684 Return: 685 Description: MUST FIX THIS!!! 686 ***********************************ivec.c*************************************/ 687 #if defined(notusing) 688 static 689 int 690 ivec_ct_bits( int *ptr, int n) 691 { 692 int tmp=0; 693 694 695 /* should expand to full 32 bit */ 696 while (n--) 697 { 698 if (*ptr&128) {tmp++;} 699 if (*ptr&64) {tmp++;} 700 if (*ptr&32) {tmp++;} 701 if (*ptr&16) {tmp++;} 702 if (*ptr&8) {tmp++;} 703 if (*ptr&4) {tmp++;} 704 if (*ptr&2) {tmp++;} 705 if (*ptr&1) {tmp++;} 706 ptr++; 707 } 708 709 return(tmp); 710 } 711 #endif 712 713 714 /****************************************************************************** 715 Function: ivec_sort(). 716 717 Input : offset of list to be sorted, number of elements to be sorted. 718 Output: sorted list (in ascending order). 719 Return: none. 720 Description: stack based (nonrecursive) quicksort w/brute-shell bottom. 721 ******************************************************************************/ 722 PetscErrorCode ivec_sort( int *ar, int size) 723 { 724 int *pi, *pj, temp; 725 int **top_a = (int **) offset_stack; 726 int *top_s = size_stack, *bottom_s = size_stack; 727 728 729 /* we're really interested in the offset of the last element */ 730 /* ==> length of the list is now size + 1 */ 731 size--; 732 733 /* do until we're done ... return when stack is exhausted */ 734 for (;;) 735 { 736 /* if list is large enough use quicksort partition exchange code */ 737 if (size > SORT_OPT) 738 { 739 /* start up pointer at element 1 and down at size */ 740 pi = ar+1; 741 pj = ar+size; 742 743 /* find middle element in list and swap w/ element 1 */ 744 SWAP(*(ar+(size>>1)),*pi) 745 746 /* order element 0,1,size-1 st {M,L,...,U} w/L<=M<=U */ 747 /* note ==> pivot_value in index 0 */ 748 if (*pi > *pj) 749 {SWAP(*pi,*pj)} 750 if (*ar > *pj) 751 {SWAP(*ar,*pj)} 752 else if (*pi > *ar) 753 {SWAP(*(ar),*(ar+1))} 754 755 /* partition about pivot_value ... */ 756 /* note lists of length 2 are not guaranteed to be sorted */ 757 for(;;) 758 { 759 /* walk up ... and down ... swap if equal to pivot! */ 760 do pi++; while (*pi<*ar); 761 do pj--; while (*pj>*ar); 762 763 /* if we've crossed we're done */ 764 if (pj<pi) break; 765 766 /* else swap */ 767 SWAP(*pi,*pj) 768 } 769 770 /* place pivot_value in it's correct location */ 771 SWAP(*ar,*pj) 772 773 /* test stack_size to see if we've exhausted our stack */ 774 if (top_s-bottom_s >= SORT_STACK) 775 {error_msg_fatal("ivec_sort() :: STACK EXHAUSTED!!!");} 776 777 /* push right hand child iff length > 1 */ 778 if ((*top_s = size-((int) (pi-ar)))) 779 { 780 *(top_a++) = pi; 781 size -= *top_s+2; 782 top_s++; 783 } 784 /* set up for next loop iff there is something to do */ 785 else if (size -= *top_s+2) 786 {;} 787 /* might as well pop - note NR_OPT >=2 ==> we're ok! */ 788 else 789 { 790 ar = *(--top_a); 791 size = *(--top_s); 792 } 793 } 794 795 /* else sort small list directly then pop another off stack */ 796 else 797 { 798 /* insertion sort for bottom */ 799 for (pj=ar+1;pj<=ar+size;pj++) 800 { 801 temp = *pj; 802 for (pi=pj-1;pi>=ar;pi--) 803 { 804 if (*pi <= temp) break; 805 *(pi+1)=*pi; 806 } 807 *(pi+1)=temp; 808 } 809 810 /* check to see if stack is exhausted ==> DONE */ 811 if (top_s==bottom_s) PetscFunctionReturn(0); 812 813 /* else pop another list from the stack */ 814 ar = *(--top_a); 815 size = *(--top_s); 816 } 817 } 818 PetscFunctionReturn(0); 819 } 820 821 822 823 /****************************************************************************** 824 Function: ivec_sort_companion(). 825 826 Input : offset of list to be sorted, number of elements to be sorted. 827 Output: sorted list (in ascending order). 828 Return: none. 829 Description: stack based (nonrecursive) quicksort w/brute-shell bottom. 830 ******************************************************************************/ 831 PetscErrorCode ivec_sort_companion( int *ar, int *ar2, int size) 832 { 833 int *pi, *pj, temp, temp2; 834 int **top_a = (int **)offset_stack; 835 int *top_s = size_stack, *bottom_s = size_stack; 836 int *pi2, *pj2; 837 int mid; 838 839 PetscFunctionBegin; 840 /* we're really interested in the offset of the last element */ 841 /* ==> length of the list is now size + 1 */ 842 size--; 843 844 /* do until we're done ... return when stack is exhausted */ 845 for (;;) 846 { 847 /* if list is large enough use quicksort partition exchange code */ 848 if (size > SORT_OPT) 849 { 850 /* start up pointer at element 1 and down at size */ 851 mid = size>>1; 852 pi = ar+1; 853 pj = ar+mid; 854 pi2 = ar2+1; 855 pj2 = ar2+mid; 856 857 /* find middle element in list and swap w/ element 1 */ 858 SWAP(*pi,*pj) 859 SWAP(*pi2,*pj2) 860 861 /* order element 0,1,size-1 st {M,L,...,U} w/L<=M<=U */ 862 /* note ==> pivot_value in index 0 */ 863 pj = ar+size; 864 pj2 = ar2+size; 865 if (*pi > *pj) 866 {SWAP(*pi,*pj) SWAP(*pi2,*pj2)} 867 if (*ar > *pj) 868 {SWAP(*ar,*pj) SWAP(*ar2,*pj2)} 869 else if (*pi > *ar) 870 {SWAP(*(ar),*(ar+1)) SWAP(*(ar2),*(ar2+1))} 871 872 /* partition about pivot_value ... */ 873 /* note lists of length 2 are not guaranteed to be sorted */ 874 for(;;) 875 { 876 /* walk up ... and down ... swap if equal to pivot! */ 877 do {pi++; pi2++;} while (*pi<*ar); 878 do {pj--; pj2--;} while (*pj>*ar); 879 880 /* if we've crossed we're done */ 881 if (pj<pi) break; 882 883 /* else swap */ 884 SWAP(*pi,*pj) 885 SWAP(*pi2,*pj2) 886 } 887 888 /* place pivot_value in it's correct location */ 889 SWAP(*ar,*pj) 890 SWAP(*ar2,*pj2) 891 892 /* test stack_size to see if we've exhausted our stack */ 893 if (top_s-bottom_s >= SORT_STACK) 894 {error_msg_fatal("ivec_sort_companion() :: STACK EXHAUSTED!!!");} 895 896 /* push right hand child iff length > 1 */ 897 if ((*top_s = size-((int) (pi-ar)))) 898 { 899 *(top_a++) = pi; 900 *(top_a++) = pi2; 901 size -= *top_s+2; 902 top_s++; 903 } 904 /* set up for next loop iff there is something to do */ 905 else if (size -= *top_s+2) 906 {;} 907 /* might as well pop - note NR_OPT >=2 ==> we're ok! */ 908 else 909 { 910 ar2 = *(--top_a); 911 ar = *(--top_a); 912 size = *(--top_s); 913 } 914 } 915 916 /* else sort small list directly then pop another off stack */ 917 else 918 { 919 /* insertion sort for bottom */ 920 for (pj=ar+1, pj2=ar2+1;pj<=ar+size;pj++,pj2++) 921 { 922 temp = *pj; 923 temp2 = *pj2; 924 for (pi=pj-1,pi2=pj2-1;pi>=ar;pi--,pi2--) 925 { 926 if (*pi <= temp) break; 927 *(pi+1)=*pi; 928 *(pi2+1)=*pi2; 929 } 930 *(pi+1)=temp; 931 *(pi2+1)=temp2; 932 } 933 934 /* check to see if stack is exhausted ==> DONE */ 935 if (top_s==bottom_s) PetscFunctionReturn(0); 936 937 /* else pop another list from the stack */ 938 ar2 = *(--top_a); 939 ar = *(--top_a); 940 size = *(--top_s); 941 } 942 } 943 PetscFunctionReturn(0); 944 } 945 946 947 948 /****************************************************************************** 949 Function: ivec_sort_companion_hack(). 950 951 Input : offset of list to be sorted, number of elements to be sorted. 952 Output: sorted list (in ascending order). 953 Return: none. 954 Description: stack based (nonrecursive) quicksort w/brute-shell bottom. 955 ******************************************************************************/ 956 PetscErrorCode ivec_sort_companion_hack( int *ar, int **ar2, 957 int size) 958 { 959 int *pi, *pj, temp, *ptr; 960 int **top_a = (int **)offset_stack; 961 int *top_s = size_stack, *bottom_s = size_stack; 962 int **pi2, **pj2; 963 int mid; 964 965 PetscFunctionBegin; 966 /* we're really interested in the offset of the last element */ 967 /* ==> length of the list is now size + 1 */ 968 size--; 969 970 /* do until we're done ... return when stack is exhausted */ 971 for (;;) 972 { 973 /* if list is large enough use quicksort partition exchange code */ 974 if (size > SORT_OPT) 975 { 976 /* start up pointer at element 1 and down at size */ 977 mid = size>>1; 978 pi = ar+1; 979 pj = ar+mid; 980 pi2 = ar2+1; 981 pj2 = ar2+mid; 982 983 /* find middle element in list and swap w/ element 1 */ 984 SWAP(*pi,*pj) 985 P_SWAP(*pi2,*pj2) 986 987 /* order element 0,1,size-1 st {M,L,...,U} w/L<=M<=U */ 988 /* note ==> pivot_value in index 0 */ 989 pj = ar+size; 990 pj2 = ar2+size; 991 if (*pi > *pj) 992 {SWAP(*pi,*pj) P_SWAP(*pi2,*pj2)} 993 if (*ar > *pj) 994 {SWAP(*ar,*pj) P_SWAP(*ar2,*pj2)} 995 else if (*pi > *ar) 996 {SWAP(*(ar),*(ar+1)) P_SWAP(*(ar2),*(ar2+1))} 997 998 /* partition about pivot_value ... */ 999 /* note lists of length 2 are not guaranteed to be sorted */ 1000 for(;;) 1001 { 1002 /* walk up ... and down ... swap if equal to pivot! */ 1003 do {pi++; pi2++;} while (*pi<*ar); 1004 do {pj--; pj2--;} while (*pj>*ar); 1005 1006 /* if we've crossed we're done */ 1007 if (pj<pi) break; 1008 1009 /* else swap */ 1010 SWAP(*pi,*pj) 1011 P_SWAP(*pi2,*pj2) 1012 } 1013 1014 /* place pivot_value in it's correct location */ 1015 SWAP(*ar,*pj) 1016 P_SWAP(*ar2,*pj2) 1017 1018 /* test stack_size to see if we've exhausted our stack */ 1019 if (top_s-bottom_s >= SORT_STACK) 1020 {error_msg_fatal("ivec_sort_companion_hack() :: STACK EXHAUSTED!!!");} 1021 1022 /* push right hand child iff length > 1 */ 1023 if ((*top_s = size-((int) (pi-ar)))) 1024 { 1025 *(top_a++) = pi; 1026 *(top_a++) = (int*) pi2; 1027 size -= *top_s+2; 1028 top_s++; 1029 } 1030 /* set up for next loop iff there is something to do */ 1031 else if (size -= *top_s+2) 1032 {;} 1033 /* might as well pop - note NR_OPT >=2 ==> we're ok! */ 1034 else 1035 { 1036 ar2 = (int **) *(--top_a); 1037 ar = *(--top_a); 1038 size = *(--top_s); 1039 } 1040 } 1041 1042 /* else sort small list directly then pop another off stack */ 1043 else 1044 { 1045 /* insertion sort for bottom */ 1046 for (pj=ar+1, pj2=ar2+1;pj<=ar+size;pj++,pj2++) 1047 { 1048 temp = *pj; 1049 ptr = *pj2; 1050 for (pi=pj-1,pi2=pj2-1;pi>=ar;pi--,pi2--) 1051 { 1052 if (*pi <= temp) break; 1053 *(pi+1)=*pi; 1054 *(pi2+1)=*pi2; 1055 } 1056 *(pi+1)=temp; 1057 *(pi2+1)=ptr; 1058 } 1059 1060 /* check to see if stack is exhausted ==> DONE */ 1061 if (top_s==bottom_s) PetscFunctionReturn(0); 1062 1063 /* else pop another list from the stack */ 1064 ar2 = (int **)*(--top_a); 1065 ar = *(--top_a); 1066 size = *(--top_s); 1067 } 1068 } 1069 PetscFunctionReturn(0); 1070 } 1071 1072 1073 1074 /****************************************************************************** 1075 Function: SMI_sort(). 1076 Input : offset of list to be sorted, number of elements to be sorted. 1077 Output: sorted list (in ascending order). 1078 Return: none. 1079 Description: stack based (nonrecursive) quicksort w/brute-shell bottom. 1080 ******************************************************************************/ 1081 PetscErrorCode SMI_sort(void *ar1, void *ar2, int size, int type) 1082 { 1083 PetscFunctionBegin; 1084 if (type == SORT_INTEGER) 1085 { 1086 if (ar2) 1087 {ivec_sort_companion((int*)ar1,(int*)ar2,size);} 1088 else 1089 {ivec_sort((int*)ar1,size);} 1090 } 1091 else if (type == SORT_INT_PTR) 1092 { 1093 if (ar2) 1094 {ivec_sort_companion_hack((int*)ar1,(int **)ar2,size);} 1095 else 1096 {ivec_sort((int*)ar1,size);} 1097 } 1098 1099 else 1100 { 1101 error_msg_fatal("SMI_sort only does SORT_INTEGER!"); 1102 } 1103 PetscFunctionReturn(0); 1104 } 1105 1106 1107 1108 /**********************************ivec.c************************************** 1109 Function ivec_linear_search() 1110 1111 Input : 1112 Output: 1113 Return: 1114 Description: 1115 ***********************************ivec.c*************************************/ 1116 int 1117 ivec_linear_search( int item, int *list, int n) 1118 { 1119 int tmp = n-1; 1120 PetscFunctionBegin; 1121 while (n--) {if (*list++ == item) {return(tmp-n);}} 1122 return(-1); 1123 } 1124 1125 1126 1127 /**********************************ivec.c************************************** 1128 Function ivec_binary_search() 1129 1130 Input : 1131 Output: 1132 Return: 1133 Description: 1134 ***********************************ivec.c*************************************/ 1135 int 1136 ivec_binary_search( int item, int *list, int rh) 1137 { 1138 int mid, lh=0; 1139 1140 rh--; 1141 while (lh<=rh) 1142 { 1143 mid = (lh+rh)>>1; 1144 if (*(list+mid) == item) 1145 {return(mid);} 1146 if (*(list+mid) > item) 1147 {rh = mid-1;} 1148 else 1149 {lh = mid+1;} 1150 } 1151 return(-1); 1152 } 1153 1154 1155 1156 /**********************************ivec.c************************************** 1157 Function rvec_dump 1158 1159 Input : 1160 Output: 1161 Return: 1162 Description: 1163 ***********************************ivec.c*************************************/ 1164 PetscErrorCode rvec_dump(PetscScalar *v, int n, int tag, int tag2, char * s) 1165 { 1166 int i; 1167 PetscFunctionBegin; 1168 printf("%2d %2d %s %2d :: ",tag,tag2,s,my_id); 1169 for (i=0;i<n;i++) 1170 {printf("%f ",v[i]);} 1171 printf("\n"); 1172 fflush(stdout); 1173 PetscFunctionReturn(0); 1174 } 1175 1176 1177 1178 /**********************************ivec.c************************************** 1179 Function rvec_lb_ub() 1180 1181 Input : 1182 Output: 1183 Return: 1184 Description: 1185 ***********************************ivec.c*************************************/ 1186 PetscErrorCode rvec_lb_ub( PetscScalar *arg1, int n, PetscScalar *lb, PetscScalar *ub) 1187 { 1188 PetscScalar min = REAL_MAX; 1189 PetscScalar max = -REAL_MAX; 1190 1191 PetscFunctionBegin; 1192 while (n--) 1193 { 1194 min = PetscMin(min,*arg1); 1195 max = PetscMax(max,*arg1); 1196 arg1++; 1197 } 1198 1199 *lb=min; 1200 *ub=max; 1201 PetscFunctionReturn(0); 1202 } 1203 1204 1205 1206 /********************************ivec.c************************************** 1207 Function rvec_copy() 1208 1209 Input : 1210 Output: 1211 Return: 1212 Description: 1213 *********************************ivec.c*************************************/ 1214 PetscErrorCode rvec_copy( PetscScalar *arg1, PetscScalar *arg2, int n) 1215 { 1216 PetscFunctionBegin; 1217 while (n--) {*arg1++ = *arg2++;} 1218 PetscFunctionReturn(0); 1219 } 1220 1221 1222 1223 /********************************ivec.c************************************** 1224 Function rvec_zero() 1225 1226 Input : 1227 Output: 1228 Return: 1229 Description: 1230 *********************************ivec.c*************************************/ 1231 PetscErrorCode rvec_zero( PetscScalar *arg1, int n) 1232 { 1233 PetscFunctionBegin; 1234 while (n--) {*arg1++ = 0.0;} 1235 PetscFunctionReturn(0); 1236 } 1237 1238 1239 1240 /**********************************ivec.c************************************** 1241 Function rvec_one() 1242 1243 Input : 1244 Output: 1245 Return: 1246 Description: 1247 ***********************************ivec.c*************************************/ 1248 PetscErrorCode rvec_one( PetscScalar *arg1, int n) 1249 { 1250 PetscFunctionBegin; 1251 while (n--) {*arg1++ = 1.0;} 1252 PetscFunctionReturn(0); 1253 } 1254 1255 1256 1257 /**********************************ivec.c************************************** 1258 Function rvec_neg_one() 1259 1260 Input : 1261 Output: 1262 Return: 1263 Description: 1264 ***********************************ivec.c*************************************/ 1265 PetscErrorCode rvec_neg_one( PetscScalar *arg1, int n) 1266 { 1267 PetscFunctionBegin; 1268 while (n--) {*arg1++ = -1.0;} 1269 PetscFunctionReturn(0); 1270 } 1271 1272 1273 1274 /**********************************ivec.c************************************** 1275 Function rvec_set() 1276 1277 Input : 1278 Output: 1279 Return: 1280 Description: 1281 ***********************************ivec.c*************************************/ 1282 PetscErrorCode rvec_set( PetscScalar *arg1, PetscScalar arg2, int n) 1283 { 1284 PetscFunctionBegin; 1285 while (n--) {*arg1++ = arg2;} 1286 PetscFunctionReturn(0); 1287 } 1288 1289 1290 1291 /**********************************ivec.c************************************** 1292 Function rvec_scale() 1293 1294 Input : 1295 Output: 1296 Return: 1297 Description: 1298 ***********************************ivec.c*************************************/ 1299 PetscErrorCode rvec_scale( PetscScalar *arg1, PetscScalar arg2, int n) 1300 { 1301 PetscFunctionBegin; 1302 while (n--) {*arg1++ *= arg2;} 1303 PetscFunctionReturn(0); 1304 } 1305 1306 1307 1308 /********************************ivec.c************************************** 1309 Function rvec_add() 1310 1311 Input : 1312 Output: 1313 Return: 1314 Description: 1315 *********************************ivec.c*************************************/ 1316 PetscErrorCode rvec_add( PetscScalar *arg1, PetscScalar *arg2, int n) 1317 { 1318 PetscFunctionBegin; 1319 while (n--) {*arg1++ += *arg2++;} 1320 PetscFunctionReturn(0); 1321 } 1322 1323 1324 1325 /********************************ivec.c************************************** 1326 Function rvec_dot() 1327 1328 Input : 1329 Output: 1330 Return: 1331 Description: 1332 *********************************ivec.c*************************************/ 1333 PetscScalar 1334 rvec_dot( PetscScalar *arg1, PetscScalar *arg2, int n) 1335 { 1336 PetscScalar dot=0.0; 1337 1338 while (n--) {dot+= *arg1++ * *arg2++;} 1339 1340 return(dot); 1341 } 1342 1343 1344 1345 /********************************ivec.c************************************** 1346 Function rvec_axpy() 1347 1348 Input : 1349 Output: 1350 Return: 1351 Description: 1352 *********************************ivec.c*************************************/ 1353 PetscErrorCode rvec_axpy( PetscScalar *arg1, PetscScalar *arg2, PetscScalar scale, 1354 int n) 1355 { 1356 PetscFunctionBegin; 1357 while (n--) {*arg1++ += scale * *arg2++;} 1358 PetscFunctionReturn(0); 1359 } 1360 1361 1362 /********************************ivec.c************************************** 1363 Function rvec_mult() 1364 1365 Input : 1366 Output: 1367 Return: 1368 Description: 1369 *********************************ivec.c*************************************/ 1370 PetscErrorCode rvec_mult( PetscScalar *arg1, PetscScalar *arg2, int n) 1371 { 1372 PetscFunctionBegin; 1373 while (n--) {*arg1++ *= *arg2++;} 1374 PetscFunctionReturn(0); 1375 } 1376 1377 1378 1379 /********************************ivec.c************************************** 1380 Function rvec_max() 1381 1382 Input : 1383 Output: 1384 Return: 1385 Description: 1386 *********************************ivec.c*************************************/ 1387 PetscErrorCode rvec_max( PetscScalar *arg1, PetscScalar *arg2, int n) 1388 { 1389 PetscFunctionBegin; 1390 while (n--) {*arg1 = PetscMax(*arg1,*arg2); arg1++; arg2++;} 1391 PetscFunctionReturn(0); 1392 } 1393 1394 1395 1396 /********************************ivec.c************************************** 1397 Function rvec_max_abs() 1398 1399 Input : 1400 Output: 1401 Return: 1402 Description: 1403 *********************************ivec.c*************************************/ 1404 PetscErrorCode rvec_max_abs( PetscScalar *arg1, PetscScalar *arg2, int n) 1405 { 1406 PetscFunctionBegin; 1407 while (n--) {*arg1 = MAX_FABS(*arg1,*arg2); arg1++; arg2++;} 1408 PetscFunctionReturn(0); 1409 } 1410 1411 1412 1413 /********************************ivec.c************************************** 1414 Function rvec_min() 1415 1416 Input : 1417 Output: 1418 Return: 1419 Description: 1420 *********************************ivec.c*************************************/ 1421 PetscErrorCode rvec_min( PetscScalar *arg1, PetscScalar *arg2, int n) 1422 { 1423 PetscFunctionBegin; 1424 while (n--) {*arg1 = PetscMin(*arg1,*arg2); arg1++; arg2++;} 1425 PetscFunctionReturn(0); 1426 } 1427 1428 1429 1430 /********************************ivec.c************************************** 1431 Function rvec_min_abs() 1432 1433 Input : 1434 Output: 1435 Return: 1436 Description: 1437 *********************************ivec.c*************************************/ 1438 PetscErrorCode rvec_min_abs( PetscScalar *arg1, PetscScalar *arg2, int n) 1439 { 1440 PetscFunctionBegin; 1441 while (n--) {*arg1 = MIN_FABS(*arg1,*arg2); arg1++; arg2++;} 1442 PetscFunctionReturn(0); 1443 } 1444 1445 1446 1447 /********************************ivec.c************************************** 1448 Function rvec_exists() 1449 1450 Input : 1451 Output: 1452 Return: 1453 Description: 1454 *********************************ivec.c*************************************/ 1455 PetscErrorCode rvec_exists( PetscScalar *arg1, PetscScalar *arg2, int n) 1456 { 1457 PetscFunctionBegin; 1458 while (n--) {*arg1 = EXISTS(*arg1,*arg2); arg1++; arg2++;} 1459 PetscFunctionReturn(0); 1460 } 1461 1462 1463 1464 /**********************************ivec.c************************************** 1465 Function rvec_non_uniform() 1466 1467 Input : 1468 Output: 1469 Return: 1470 Description: 1471 ***********************************ivec.c*************************************/ 1472 PetscErrorCode rvec_non_uniform(PetscScalar *arg1, PetscScalar *arg2, int n, int *arg3) 1473 { 1474 int i, j, type; 1475 1476 PetscFunctionBegin; 1477 /* LATER: if we're really motivated we can sort and then unsort */ 1478 for (i=0;i<n;) 1479 { 1480 /* clump 'em for now */ 1481 j=i+1; 1482 type = arg3[i]; 1483 while ((j<n)&&(arg3[j]==type)) 1484 {j++;} 1485 1486 /* how many together */ 1487 j -= i; 1488 1489 /* call appropriate ivec function */ 1490 if (type == GL_MAX) 1491 {rvec_max(arg1,arg2,j);} 1492 else if (type == GL_MIN) 1493 {rvec_min(arg1,arg2,j);} 1494 else if (type == GL_MULT) 1495 {rvec_mult(arg1,arg2,j);} 1496 else if (type == GL_ADD) 1497 {rvec_add(arg1,arg2,j);} 1498 else if (type == GL_MAX_ABS) 1499 {rvec_max_abs(arg1,arg2,j);} 1500 else if (type == GL_MIN_ABS) 1501 {rvec_min_abs(arg1,arg2,j);} 1502 else if (type == GL_EXISTS) 1503 {rvec_exists(arg1,arg2,j);} 1504 else 1505 {error_msg_fatal("unrecognized type passed to rvec_non_uniform()!");} 1506 1507 arg1+=j; arg2+=j; i+=j; 1508 } 1509 PetscFunctionReturn(0); 1510 } 1511 1512 1513 1514 /**********************************ivec.c************************************** 1515 Function rvec_fct_addr() 1516 1517 Input : 1518 Output: 1519 Return: 1520 Description: 1521 ***********************************ivec.c*************************************/ 1522 vfp rvec_fct_addr( int type) 1523 { 1524 if (type == NON_UNIFORM) 1525 {return((PetscErrorCode (*)(void*, void *, int, ...))&rvec_non_uniform);} 1526 else if (type == GL_MAX) 1527 {return((PetscErrorCode (*)(void*, void *, int, ...))&rvec_max);} 1528 else if (type == GL_MIN) 1529 {return((PetscErrorCode (*)(void*, void *, int, ...))&rvec_min);} 1530 else if (type == GL_MULT) 1531 {return((PetscErrorCode (*)(void*, void *, int, ...))&rvec_mult);} 1532 else if (type == GL_ADD) 1533 {return((PetscErrorCode (*)(void*, void *, int, ...))&rvec_add);} 1534 else if (type == GL_MAX_ABS) 1535 {return((PetscErrorCode (*)(void*, void *, int, ...))&rvec_max_abs);} 1536 else if (type == GL_MIN_ABS) 1537 {return((PetscErrorCode (*)(void*, void *, int, ...))&rvec_min_abs);} 1538 else if (type == GL_EXISTS) 1539 {return((PetscErrorCode (*)(void*, void *, int, ...))&rvec_exists);} 1540 1541 /* catch all ... not good if we get here */ 1542 return(NULL); 1543 } 1544 1545 1546 /****************************************************************************** 1547 Function: my_sort(). 1548 Input : offset of list to be sorted, number of elements to be sorted. 1549 Output: sorted list (in ascending order). 1550 Return: none. 1551 Description: stack based (nonrecursive) quicksort w/brute-shell bottom. 1552 ******************************************************************************/ 1553 PetscErrorCode rvec_sort( PetscScalar *ar, int Size) 1554 { 1555 PetscScalar *pi, *pj, temp; 1556 PetscScalar **top_a = (PetscScalar **)offset_stack; 1557 long *top_s = psize_stack, *bottom_s = psize_stack; 1558 long size = (long) Size; 1559 1560 PetscFunctionBegin; 1561 /* we're really interested in the offset of the last element */ 1562 /* ==> length of the list is now size + 1 */ 1563 size--; 1564 1565 /* do until we're done ... return when stack is exhausted */ 1566 for (;;) 1567 { 1568 /* if list is large enough use quicksort partition exchange code */ 1569 if (size > SORT_OPT) 1570 { 1571 /* start up pointer at element 1 and down at size */ 1572 pi = ar+1; 1573 pj = ar+size; 1574 1575 /* find middle element in list and swap w/ element 1 */ 1576 SWAP(*(ar+(size>>1)),*pi) 1577 1578 pj = ar+size; 1579 1580 /* order element 0,1,size-1 st {M,L,...,U} w/L<=M<=U */ 1581 /* note ==> pivot_value in index 0 */ 1582 if (*pi > *pj) 1583 {SWAP(*pi,*pj)} 1584 if (*ar > *pj) 1585 {SWAP(*ar,*pj)} 1586 else if (*pi > *ar) 1587 {SWAP(*(ar),*(ar+1))} 1588 1589 /* partition about pivot_value ... */ 1590 /* note lists of length 2 are not guaranteed to be sorted */ 1591 for(;;) 1592 { 1593 /* walk up ... and down ... swap if equal to pivot! */ 1594 do pi++; while (*pi<*ar); 1595 do pj--; while (*pj>*ar); 1596 1597 /* if we've crossed we're done */ 1598 if (pj<pi) break; 1599 1600 /* else swap */ 1601 SWAP(*pi,*pj) 1602 } 1603 1604 /* place pivot_value in it's correct location */ 1605 SWAP(*ar,*pj) 1606 1607 /* test stack_size to see if we've exhausted our stack */ 1608 if (top_s-bottom_s >= SORT_STACK) 1609 {error_msg_fatal("\nSTACK EXHAUSTED!!!\n");} 1610 1611 /* push right hand child iff length > 1 */ 1612 if ((*top_s = size-(pi-ar))) 1613 { 1614 *(top_a++) = pi; 1615 size -= *top_s+2; 1616 top_s++; 1617 } 1618 /* set up for next loop iff there is something to do */ 1619 else if (size -= *top_s+2) 1620 {;} 1621 /* might as well pop - note NR_OPT >=2 ==> we're ok! */ 1622 else 1623 { 1624 ar = *(--top_a); 1625 size = *(--top_s); 1626 } 1627 } 1628 1629 /* else sort small list directly then pop another off stack */ 1630 else 1631 { 1632 /* insertion sort for bottom */ 1633 for (pj=ar+1;pj<=ar+size;pj++) 1634 { 1635 temp = *pj; 1636 for (pi=pj-1;pi>=ar;pi--) 1637 { 1638 if (*pi <= temp) break; 1639 *(pi+1)=*pi; 1640 } 1641 *(pi+1)=temp; 1642 } 1643 1644 /* check to see if stack is exhausted ==> DONE */ 1645 if (top_s==bottom_s) PetscFunctionReturn(0); 1646 1647 /* else pop another list from the stack */ 1648 ar = *(--top_a); 1649 size = *(--top_s); 1650 } 1651 } 1652 PetscFunctionReturn(0); 1653 } 1654 1655 1656 1657 /****************************************************************************** 1658 Function: my_sort(). 1659 Input : offset of list to be sorted, number of elements to be sorted. 1660 Output: sorted list (in ascending order). 1661 Return: none. 1662 Description: stack based (nonrecursive) quicksort w/brute-shell bottom. 1663 ******************************************************************************/ 1664 PetscErrorCode rvec_sort_companion( PetscScalar *ar, int *ar2, int Size) 1665 { 1666 PetscScalar *pi, *pj, temp; 1667 PetscScalar **top_a = (PetscScalar **)offset_stack; 1668 long *top_s = psize_stack, *bottom_s = psize_stack; 1669 long size = (long) Size; 1670 1671 int *pi2, *pj2; 1672 int ptr; 1673 long mid; 1674 1675 PetscFunctionBegin; 1676 /* we're really interested in the offset of the last element */ 1677 /* ==> length of the list is now size + 1 */ 1678 size--; 1679 1680 /* do until we're done ... return when stack is exhausted */ 1681 for (;;) 1682 { 1683 /* if list is large enough use quicksort partition exchange code */ 1684 if (size > SORT_OPT) 1685 { 1686 /* start up pointer at element 1 and down at size */ 1687 mid = size>>1; 1688 pi = ar+1; 1689 pj = ar+mid; 1690 pi2 = ar2+1; 1691 pj2 = ar2+mid; 1692 1693 /* find middle element in list and swap w/ element 1 */ 1694 SWAP(*pi,*pj) 1695 P_SWAP(*pi2,*pj2) 1696 1697 /* order element 0,1,size-1 st {M,L,...,U} w/L<=M<=U */ 1698 /* note ==> pivot_value in index 0 */ 1699 pj = ar+size; 1700 pj2 = ar2+size; 1701 if (*pi > *pj) 1702 {SWAP(*pi,*pj) P_SWAP(*pi2,*pj2)} 1703 if (*ar > *pj) 1704 {SWAP(*ar,*pj) P_SWAP(*ar2,*pj2)} 1705 else if (*pi > *ar) 1706 {SWAP(*(ar),*(ar+1)) P_SWAP(*(ar2),*(ar2+1))} 1707 1708 /* partition about pivot_value ... */ 1709 /* note lists of length 2 are not guaranteed to be sorted */ 1710 for(;;) 1711 { 1712 /* walk up ... and down ... swap if equal to pivot! */ 1713 do {pi++; pi2++;} while (*pi<*ar); 1714 do {pj--; pj2--;} while (*pj>*ar); 1715 1716 /* if we've crossed we're done */ 1717 if (pj<pi) break; 1718 1719 /* else swap */ 1720 SWAP(*pi,*pj) 1721 P_SWAP(*pi2,*pj2) 1722 } 1723 1724 /* place pivot_value in it's correct location */ 1725 SWAP(*ar,*pj) 1726 P_SWAP(*ar2,*pj2) 1727 1728 /* test stack_size to see if we've exhausted our stack */ 1729 if (top_s-bottom_s >= SORT_STACK) 1730 {error_msg_fatal("\nSTACK EXHAUSTED!!!\n");} 1731 1732 /* push right hand child iff length > 1 */ 1733 if ((*top_s = size-(pi-ar))) 1734 { 1735 *(top_a++) = pi; 1736 *(top_a++) = (PetscScalar *) pi2; 1737 size -= *top_s+2; 1738 top_s++; 1739 } 1740 /* set up for next loop iff there is something to do */ 1741 else if (size -= *top_s+2) 1742 {;} 1743 /* might as well pop - note NR_OPT >=2 ==> we're ok! */ 1744 else 1745 { 1746 ar2 = (int*) *(--top_a); 1747 ar = *(--top_a); 1748 size = *(--top_s); 1749 } 1750 } 1751 1752 /* else sort small list directly then pop another off stack */ 1753 else 1754 { 1755 /* insertion sort for bottom */ 1756 for (pj=ar+1, pj2=ar2+1;pj<=ar+size;pj++,pj2++) 1757 { 1758 temp = *pj; 1759 ptr = *pj2; 1760 for (pi=pj-1,pi2=pj2-1;pi>=ar;pi--,pi2--) 1761 { 1762 if (*pi <= temp) break; 1763 *(pi+1)=*pi; 1764 *(pi2+1)=*pi2; 1765 } 1766 *(pi+1)=temp; 1767 *(pi2+1)=ptr; 1768 } 1769 1770 /* check to see if stack is exhausted ==> DONE */ 1771 if (top_s==bottom_s) PetscFunctionReturn(0); 1772 1773 /* else pop another list from the stack */ 1774 ar2 = (int*) *(--top_a); 1775 ar = *(--top_a); 1776 size = *(--top_s); 1777 } 1778 } 1779 PetscFunctionReturn(0); 1780 } 1781 1782 1783 1784 1785 1786 /**********************************ivec.c************************************** 1787 Function ivec_binary_search() 1788 1789 Input : 1790 Output: 1791 Return: 1792 Description: 1793 ***********************************ivec.c*************************************/ 1794 int 1795 rvec_binary_search( PetscScalar item, PetscScalar *list, int rh) 1796 { 1797 int mid, lh=0; 1798 PetscFunctionBegin; 1799 rh--; 1800 while (lh<=rh) 1801 { 1802 mid = (lh+rh)>>1; 1803 if (*(list+mid) == item) 1804 {return(mid);} 1805 if (*(list+mid) > item) 1806 {rh = mid-1;} 1807 else 1808 {lh = mid+1;} 1809 } 1810 return(-1); 1811 } 1812 1813 1814 1815 1816 1817