xref: /petsc/src/ksp/pc/impls/tfs/ivec.c (revision b458e8f169278db94fa1d489e1d3db410fc8a4c8)
1827bd09bSSatish Balay 
29895aa37SBarry Smith 
3827bd09bSSatish Balay /**********************************ivec.c**************************************
4827bd09bSSatish Balay 
5827bd09bSSatish Balay Author: Henry M. Tufo III
6827bd09bSSatish Balay 
7827bd09bSSatish Balay e-mail: hmt@cs.brown.edu
8827bd09bSSatish Balay 
9827bd09bSSatish Balay snail-mail:
10827bd09bSSatish Balay Division of Applied Mathematics
11827bd09bSSatish Balay Brown University
12827bd09bSSatish Balay Providence, RI 02912
13827bd09bSSatish Balay 
14827bd09bSSatish Balay Last Modification:
15827bd09bSSatish Balay 6.21.97
16827bd09bSSatish Balay ***********************************ivec.c*************************************/
17827bd09bSSatish Balay 
18c6db04a5SJed Brown #include <../src/ksp/pc/impls/tfs/tfs.h>
19827bd09bSSatish Balay 
20827bd09bSSatish Balay /* sorting args ivec.c ivec.c ... */
21827bd09bSSatish Balay #define   SORT_OPT     6
22827bd09bSSatish Balay #define   SORT_STACK   50000
23827bd09bSSatish Balay 
24827bd09bSSatish Balay 
25827bd09bSSatish Balay /* allocate an address and size stack for sorter(s) */
26827bd09bSSatish Balay static void     *offset_stack[2*SORT_STACK];
2752f87cdaSBarry Smith static PetscInt size_stack[SORT_STACK];
28827bd09bSSatish Balay 
297b1ae94cSBarry Smith /***********************************ivec.c*************************************/
30ca8e9878SJed Brown PetscInt *PCTFS_ivec_copy(PetscInt *arg1, PetscInt *arg2, PetscInt n)
31827bd09bSSatish Balay {
322fa5cd67SKarl Rupp   while (n--) *arg1++ = *arg2++;
33827bd09bSSatish Balay   return(arg1);
34827bd09bSSatish Balay }
35827bd09bSSatish Balay 
367b1ae94cSBarry Smith /***********************************ivec.c*************************************/
37ca8e9878SJed Brown PetscErrorCode PCTFS_ivec_zero(PetscInt *arg1, PetscInt n)
38827bd09bSSatish Balay {
393fdc5746SBarry Smith   PetscFunctionBegin;
402fa5cd67SKarl Rupp   while (n--) *arg1++ = 0;
413fdc5746SBarry Smith   PetscFunctionReturn(0);
42827bd09bSSatish Balay }
43827bd09bSSatish Balay 
447b1ae94cSBarry Smith /***********************************ivec.c*************************************/
45ca8e9878SJed Brown PetscErrorCode PCTFS_ivec_set(PetscInt *arg1, PetscInt arg2, PetscInt n)
46827bd09bSSatish Balay {
473fdc5746SBarry Smith   PetscFunctionBegin;
482fa5cd67SKarl Rupp   while (n--) *arg1++ = arg2;
493fdc5746SBarry Smith   PetscFunctionReturn(0);
50827bd09bSSatish Balay }
51827bd09bSSatish Balay 
527b1ae94cSBarry Smith /***********************************ivec.c*************************************/
53ca8e9878SJed Brown PetscErrorCode PCTFS_ivec_max(PetscInt *arg1, PetscInt *arg2, PetscInt n)
54827bd09bSSatish Balay {
553fdc5746SBarry Smith   PetscFunctionBegin;
5639945688SSatish Balay   while (n--) { *arg1 = PetscMax(*arg1,*arg2); arg1++; arg2++; }
573fdc5746SBarry Smith   PetscFunctionReturn(0);
58827bd09bSSatish Balay }
59827bd09bSSatish Balay 
607b1ae94cSBarry Smith /***********************************ivec.c*************************************/
61ca8e9878SJed Brown PetscErrorCode PCTFS_ivec_min(PetscInt *arg1,  PetscInt *arg2,  PetscInt n)
62827bd09bSSatish Balay {
633fdc5746SBarry Smith   PetscFunctionBegin;
642fa5cd67SKarl Rupp   while (n--) {
652fa5cd67SKarl Rupp     *(arg1) = PetscMin(*arg1,*arg2);
662fa5cd67SKarl Rupp     arg1++;
672fa5cd67SKarl Rupp     arg2++;
682fa5cd67SKarl Rupp   }
693fdc5746SBarry Smith   PetscFunctionReturn(0);
70827bd09bSSatish Balay }
71827bd09bSSatish Balay 
727b1ae94cSBarry Smith /***********************************ivec.c*************************************/
73ca8e9878SJed Brown PetscErrorCode PCTFS_ivec_mult(PetscInt *arg1,  PetscInt *arg2,  PetscInt n)
74827bd09bSSatish Balay {
753fdc5746SBarry Smith   PetscFunctionBegin;
762fa5cd67SKarl Rupp   while (n--) *arg1++ *= *arg2++;
773fdc5746SBarry Smith   PetscFunctionReturn(0);
78827bd09bSSatish Balay }
79827bd09bSSatish Balay 
807b1ae94cSBarry Smith /***********************************ivec.c*************************************/
81ca8e9878SJed Brown PetscErrorCode PCTFS_ivec_add(PetscInt *arg1,  PetscInt *arg2,  PetscInt n)
82827bd09bSSatish Balay {
833fdc5746SBarry Smith   PetscFunctionBegin;
842fa5cd67SKarl Rupp   while (n--) *arg1++ += *arg2++;
853fdc5746SBarry Smith   PetscFunctionReturn(0);
86827bd09bSSatish Balay }
87827bd09bSSatish Balay 
887b1ae94cSBarry Smith /***********************************ivec.c*************************************/
89ca8e9878SJed Brown PetscErrorCode PCTFS_ivec_lxor(PetscInt *arg1,  PetscInt *arg2,  PetscInt n)
90827bd09bSSatish Balay {
913fdc5746SBarry Smith   PetscFunctionBegin;
922fa5cd67SKarl Rupp   while (n--) {
932fa5cd67SKarl Rupp     *arg1=((*arg1 || *arg2) && !(*arg1 && *arg2));
942fa5cd67SKarl Rupp     arg1++;
952fa5cd67SKarl Rupp     arg2++;
962fa5cd67SKarl Rupp   }
973fdc5746SBarry Smith   PetscFunctionReturn(0);
98827bd09bSSatish Balay }
99827bd09bSSatish Balay 
1007b1ae94cSBarry Smith /***********************************ivec.c*************************************/
101ca8e9878SJed Brown PetscErrorCode PCTFS_ivec_xor(PetscInt *arg1,  PetscInt *arg2,  PetscInt n)
102827bd09bSSatish Balay {
1033fdc5746SBarry Smith   PetscFunctionBegin;
1042fa5cd67SKarl Rupp   while (n--) *arg1++ ^= *arg2++;
1053fdc5746SBarry Smith   PetscFunctionReturn(0);
106827bd09bSSatish Balay }
107827bd09bSSatish Balay 
1087b1ae94cSBarry Smith /***********************************ivec.c*************************************/
109ca8e9878SJed Brown PetscErrorCode PCTFS_ivec_or(PetscInt *arg1,  PetscInt *arg2,  PetscInt n)
110827bd09bSSatish Balay {
1113fdc5746SBarry Smith   PetscFunctionBegin;
1122fa5cd67SKarl Rupp   while (n--) *arg1++ |= *arg2++;
1133fdc5746SBarry Smith   PetscFunctionReturn(0);
114827bd09bSSatish Balay }
115827bd09bSSatish Balay 
1167b1ae94cSBarry Smith /***********************************ivec.c*************************************/
117ca8e9878SJed Brown PetscErrorCode PCTFS_ivec_lor(PetscInt *arg1,  PetscInt *arg2,  PetscInt n)
118827bd09bSSatish Balay {
1193fdc5746SBarry Smith   PetscFunctionBegin;
1202fa5cd67SKarl Rupp   while (n--) {
1212fa5cd67SKarl Rupp     *arg1 = (*arg1 || *arg2);
1222fa5cd67SKarl Rupp     arg1++;
1232fa5cd67SKarl Rupp     arg2++;
1242fa5cd67SKarl Rupp   }
1253fdc5746SBarry Smith   PetscFunctionReturn(0);
126827bd09bSSatish Balay }
127827bd09bSSatish Balay 
1287b1ae94cSBarry Smith /***********************************ivec.c*************************************/
129ca8e9878SJed Brown PetscErrorCode PCTFS_ivec_and(PetscInt *arg1,  PetscInt *arg2,  PetscInt n)
130827bd09bSSatish Balay {
1313fdc5746SBarry Smith   PetscFunctionBegin;
1322fa5cd67SKarl Rupp   while (n--) *arg1++ &= *arg2++;
1333fdc5746SBarry Smith   PetscFunctionReturn(0);
134827bd09bSSatish Balay }
135827bd09bSSatish Balay 
1367b1ae94cSBarry Smith /***********************************ivec.c*************************************/
137ca8e9878SJed Brown PetscErrorCode PCTFS_ivec_land(PetscInt *arg1,  PetscInt *arg2,  PetscInt n)
138827bd09bSSatish Balay {
1393fdc5746SBarry Smith   PetscFunctionBegin;
1402fa5cd67SKarl Rupp   while (n--) {
1412fa5cd67SKarl Rupp     *arg1 = (*arg1 && *arg2);
1422fa5cd67SKarl Rupp     arg1++;
1432fa5cd67SKarl Rupp     arg2++;
1442fa5cd67SKarl Rupp   }
1453fdc5746SBarry Smith   PetscFunctionReturn(0);
146827bd09bSSatish Balay }
147827bd09bSSatish Balay 
1487b1ae94cSBarry Smith /***********************************ivec.c*************************************/
149ca8e9878SJed Brown PetscErrorCode PCTFS_ivec_and3(PetscInt *arg1,  PetscInt *arg2,  PetscInt *arg3, PetscInt n)
150827bd09bSSatish Balay {
1513fdc5746SBarry Smith   PetscFunctionBegin;
1522fa5cd67SKarl Rupp   while (n--) *arg1++ = (*arg2++ & *arg3++);
1533fdc5746SBarry Smith   PetscFunctionReturn(0);
154827bd09bSSatish Balay }
155827bd09bSSatish Balay 
1567b1ae94cSBarry Smith /***********************************ivec.c*************************************/
157ca8e9878SJed Brown PetscInt PCTFS_ivec_sum(PetscInt *arg1,  PetscInt n)
158827bd09bSSatish Balay {
15952f87cdaSBarry Smith   PetscInt tmp = 0;
1602fa5cd67SKarl Rupp   while (n--) tmp += *arg1++;
161827bd09bSSatish Balay   return(tmp);
162827bd09bSSatish Balay }
163827bd09bSSatish Balay 
1647b1ae94cSBarry Smith /***********************************ivec.c*************************************/
165ca8e9878SJed Brown PetscErrorCode PCTFS_ivec_non_uniform(PetscInt *arg1, PetscInt *arg2,  PetscInt n,  PetscInt *arg3)
166827bd09bSSatish Balay {
16752f87cdaSBarry Smith   PetscInt i, j, type;
168827bd09bSSatish Balay 
1693fdc5746SBarry Smith   PetscFunctionBegin;
170827bd09bSSatish Balay   /* LATER: if we're really motivated we can sort and then unsort */
171db4deed7SKarl Rupp   for (i=0; i<n;) {
172827bd09bSSatish Balay     /* clump 'em for now */
173827bd09bSSatish Balay     j    =i+1;
174827bd09bSSatish Balay     type = arg3[i];
1752fa5cd67SKarl Rupp     while ((j<n)&&(arg3[j]==type)) j++;
176827bd09bSSatish Balay 
177827bd09bSSatish Balay     /* how many together */
178827bd09bSSatish Balay     j -= i;
179827bd09bSSatish Balay 
180827bd09bSSatish Balay     /* call appropriate ivec function */
1812fa5cd67SKarl Rupp     if (type == GL_MAX)        PCTFS_ivec_max(arg1,arg2,j);
1822fa5cd67SKarl Rupp     else if (type == GL_MIN)   PCTFS_ivec_min(arg1,arg2,j);
1832fa5cd67SKarl Rupp     else if (type == GL_MULT)  PCTFS_ivec_mult(arg1,arg2,j);
1842fa5cd67SKarl Rupp     else if (type == GL_ADD)   PCTFS_ivec_add(arg1,arg2,j);
1852fa5cd67SKarl Rupp     else if (type == GL_B_XOR) PCTFS_ivec_xor(arg1,arg2,j);
1862fa5cd67SKarl Rupp     else if (type == GL_B_OR)  PCTFS_ivec_or(arg1,arg2,j);
1872fa5cd67SKarl Rupp     else if (type == GL_B_AND) PCTFS_ivec_and(arg1,arg2,j);
1882fa5cd67SKarl Rupp     else if (type == GL_L_XOR) PCTFS_ivec_lxor(arg1,arg2,j);
1892fa5cd67SKarl Rupp     else if (type == GL_L_OR)  PCTFS_ivec_lor(arg1,arg2,j);
1902fa5cd67SKarl Rupp     else if (type == GL_L_AND) PCTFS_ivec_land(arg1,arg2,j);
191db4deed7SKarl Rupp     else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"unrecognized type passed to PCTFS_ivec_non_uniform()!");
192827bd09bSSatish Balay 
193827bd09bSSatish Balay     arg1+=j; arg2+=j; i+=j;
194827bd09bSSatish Balay   }
1953fdc5746SBarry Smith   PetscFunctionReturn(0);
196827bd09bSSatish Balay }
197827bd09bSSatish Balay 
1987b1ae94cSBarry Smith /***********************************ivec.c*************************************/
199ca8e9878SJed Brown vfp PCTFS_ivec_fct_addr(PetscInt type)
200827bd09bSSatish Balay {
2012fa5cd67SKarl Rupp   if (type == NON_UNIFORM)   return((PetscErrorCode (*)(void*, void*, PetscInt, ...))&PCTFS_ivec_non_uniform);
2022fa5cd67SKarl Rupp   else if (type == GL_MAX)   return((PetscErrorCode (*)(void*, void*, PetscInt, ...))&PCTFS_ivec_max);
2032fa5cd67SKarl Rupp   else if (type == GL_MIN)   return((PetscErrorCode (*)(void*, void*, PetscInt, ...))&PCTFS_ivec_min);
2042fa5cd67SKarl Rupp   else if (type == GL_MULT)  return((PetscErrorCode (*)(void*, void*, PetscInt, ...))&PCTFS_ivec_mult);
2052fa5cd67SKarl Rupp   else if (type == GL_ADD)   return((PetscErrorCode (*)(void*, void*, PetscInt, ...))&PCTFS_ivec_add);
2062fa5cd67SKarl Rupp   else if (type == GL_B_XOR) return((PetscErrorCode (*)(void*, void*, PetscInt, ...))&PCTFS_ivec_xor);
2072fa5cd67SKarl Rupp   else if (type == GL_B_OR)  return((PetscErrorCode (*)(void*, void*, PetscInt, ...))&PCTFS_ivec_or);
2082fa5cd67SKarl Rupp   else if (type == GL_B_AND) return((PetscErrorCode (*)(void*, void*, PetscInt, ...))&PCTFS_ivec_and);
2092fa5cd67SKarl Rupp   else if (type == GL_L_XOR) return((PetscErrorCode (*)(void*, void*, PetscInt, ...))&PCTFS_ivec_lxor);
2102fa5cd67SKarl Rupp   else if (type == GL_L_OR)  return((PetscErrorCode (*)(void*, void*, PetscInt, ...))&PCTFS_ivec_lor);
2112fa5cd67SKarl Rupp   else if (type == GL_L_AND) return((PetscErrorCode (*)(void*, void*, PetscInt, ...))&PCTFS_ivec_land);
212827bd09bSSatish Balay 
213827bd09bSSatish Balay   /* catch all ... not good if we get here */
214827bd09bSSatish Balay   return(NULL);
215827bd09bSSatish Balay }
216827bd09bSSatish Balay 
2177b1ae94cSBarry Smith /******************************************************************************/
218ca8e9878SJed Brown PetscErrorCode PCTFS_ivec_sort(PetscInt *ar,  PetscInt size)
219827bd09bSSatish Balay {
22052f87cdaSBarry Smith   PetscInt *pi, *pj, temp;
22152f87cdaSBarry Smith   PetscInt **top_a = (PetscInt**) offset_stack;
22252f87cdaSBarry Smith   PetscInt *top_s  = size_stack, *bottom_s = size_stack;
223827bd09bSSatish Balay 
224*b458e8f1SJose E. Roman   PetscFunctionBegin;
225827bd09bSSatish Balay   /* we're really interested in the offset of the last element */
226827bd09bSSatish Balay   /* ==> length of the list is now size + 1                    */
227827bd09bSSatish Balay   size--;
228827bd09bSSatish Balay 
229827bd09bSSatish Balay   /* do until we're done ... return when stack is exhausted */
230db4deed7SKarl Rupp   for (;;) {
231827bd09bSSatish Balay     /* if list is large enough use quicksort partition exchange code */
232db4deed7SKarl Rupp     if (size > SORT_OPT) {
233827bd09bSSatish Balay       /* start up pointer at element 1 and down at size     */
234827bd09bSSatish Balay       pi = ar+1;
235827bd09bSSatish Balay       pj = ar+size;
236827bd09bSSatish Balay 
237827bd09bSSatish Balay       /* find middle element in list and swap w/ element 1 */
238827bd09bSSatish Balay       SWAP(*(ar+(size>>1)),*pi)
239827bd09bSSatish Balay 
240827bd09bSSatish Balay       /* order element 0,1,size-1 st {M,L,...,U} w/L<=M<=U */
241827bd09bSSatish Balay       /* note ==> pivot_value in index 0                   */
242db4deed7SKarl Rupp       if (*pi > *pj) { SWAP(*pi,*pj) }
243db4deed7SKarl Rupp       if (*ar > *pj) { SWAP(*ar,*pj) }
244db4deed7SKarl Rupp       else if (*pi > *ar) { SWAP(*(ar),*(ar+1)) }
245827bd09bSSatish Balay 
246827bd09bSSatish Balay       /* partition about pivot_value ...                              */
247827bd09bSSatish Balay       /* note lists of length 2 are not guaranteed to be sorted */
248db4deed7SKarl Rupp       for (;;) {
249827bd09bSSatish Balay         /* walk up ... and down ... swap if equal to pivot! */
250827bd09bSSatish Balay         do pi++; while (*pi<*ar);
251827bd09bSSatish Balay         do pj--; while (*pj>*ar);
252827bd09bSSatish Balay 
253827bd09bSSatish Balay         /* if we've crossed we're done */
254827bd09bSSatish Balay         if (pj<pi) break;
255827bd09bSSatish Balay 
256827bd09bSSatish Balay         /* else swap */
257827bd09bSSatish Balay         SWAP(*pi,*pj)
258827bd09bSSatish Balay       }
259827bd09bSSatish Balay 
260827bd09bSSatish Balay       /* place pivot_value in it's correct location */
261827bd09bSSatish Balay       SWAP(*ar,*pj)
262827bd09bSSatish Balay 
263827bd09bSSatish Balay       /* test stack_size to see if we've exhausted our stack */
264ca8e9878SJed Brown       if (top_s-bottom_s >= SORT_STACK) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_ivec_sort() :: STACK EXHAUSTED!!!");
265827bd09bSSatish Balay 
266827bd09bSSatish Balay       /* push right hand child iff length > 1 */
267db4deed7SKarl Rupp       if ((*top_s = size-((PetscInt) (pi-ar)))) {
268827bd09bSSatish Balay         *(top_a++) = pi;
269827bd09bSSatish Balay         size      -= *top_s+2;
270827bd09bSSatish Balay         top_s++;
2712fa5cd67SKarl Rupp       } else if (size -= *top_s+2) ;   /* set up for next loop iff there is something to do */
272db4deed7SKarl Rupp       else { /* might as well pop - note NR_OPT >=2 ==> we're ok! */
273827bd09bSSatish Balay         ar   = *(--top_a);
274827bd09bSSatish Balay         size = *(--top_s);
275827bd09bSSatish Balay       }
276db4deed7SKarl Rupp     } else { /* else sort small list directly then pop another off stack */
277827bd09bSSatish Balay 
278827bd09bSSatish Balay       /* insertion sort for bottom */
279db4deed7SKarl Rupp       for (pj=ar+1; pj<=ar+size; pj++) {
280827bd09bSSatish Balay         temp = *pj;
281db4deed7SKarl Rupp         for (pi=pj-1; pi>=ar; pi--) {
282827bd09bSSatish Balay           if (*pi <= temp) break;
283827bd09bSSatish Balay           *(pi+1)=*pi;
284827bd09bSSatish Balay         }
285827bd09bSSatish Balay         *(pi+1)=temp;
286827bd09bSSatish Balay       }
287827bd09bSSatish Balay 
288827bd09bSSatish Balay       /* check to see if stack is exhausted ==> DONE */
2893fdc5746SBarry Smith       if (top_s==bottom_s) PetscFunctionReturn(0);
290827bd09bSSatish Balay 
291827bd09bSSatish Balay       /* else pop another list from the stack */
292827bd09bSSatish Balay       ar   = *(--top_a);
293827bd09bSSatish Balay       size = *(--top_s);
294827bd09bSSatish Balay     }
295827bd09bSSatish Balay   }
296827bd09bSSatish Balay }
297827bd09bSSatish Balay 
2987b1ae94cSBarry Smith /******************************************************************************/
299ca8e9878SJed Brown PetscErrorCode PCTFS_ivec_sort_companion(PetscInt *ar,  PetscInt *ar2,  PetscInt size)
300827bd09bSSatish Balay {
30152f87cdaSBarry Smith   PetscInt *pi, *pj, temp, temp2;
30252f87cdaSBarry Smith   PetscInt **top_a = (PetscInt**)offset_stack;
30352f87cdaSBarry Smith   PetscInt *top_s  = size_stack, *bottom_s = size_stack;
30452f87cdaSBarry Smith   PetscInt *pi2, *pj2;
30552f87cdaSBarry Smith   PetscInt mid;
306827bd09bSSatish Balay 
3073fdc5746SBarry Smith   PetscFunctionBegin;
308827bd09bSSatish Balay   /* we're really interested in the offset of the last element */
309827bd09bSSatish Balay   /* ==> length of the list is now size + 1                    */
310827bd09bSSatish Balay   size--;
311827bd09bSSatish Balay 
312827bd09bSSatish Balay   /* do until we're done ... return when stack is exhausted */
313db4deed7SKarl Rupp   for (;;) {
314db4deed7SKarl Rupp 
315827bd09bSSatish Balay     /* if list is large enough use quicksort partition exchange code */
316db4deed7SKarl Rupp     if (size > SORT_OPT) {
317db4deed7SKarl Rupp 
318827bd09bSSatish Balay       /* start up pointer at element 1 and down at size     */
319827bd09bSSatish Balay       mid = size>>1;
320827bd09bSSatish Balay       pi  = ar+1;
321827bd09bSSatish Balay       pj  = ar+mid;
322827bd09bSSatish Balay       pi2 = ar2+1;
323827bd09bSSatish Balay       pj2 = ar2+mid;
324827bd09bSSatish Balay 
325827bd09bSSatish Balay       /* find middle element in list and swap w/ element 1 */
326827bd09bSSatish Balay       SWAP(*pi,*pj)
327827bd09bSSatish Balay       SWAP(*pi2,*pj2)
328827bd09bSSatish Balay 
329827bd09bSSatish Balay       /* order element 0,1,size-1 st {M,L,...,U} w/L<=M<=U */
330827bd09bSSatish Balay       /* note ==> pivot_value in index 0                   */
331827bd09bSSatish Balay       pj  = ar+size;
332827bd09bSSatish Balay       pj2 = ar2+size;
333db4deed7SKarl Rupp       if (*pi > *pj) { SWAP(*pi,*pj) SWAP(*pi2,*pj2) }
334db4deed7SKarl Rupp       if (*ar > *pj) { SWAP(*ar,*pj) SWAP(*ar2,*pj2) }
335db4deed7SKarl Rupp       else if (*pi > *ar) { SWAP(*(ar),*(ar+1)) SWAP(*(ar2),*(ar2+1)) }
336827bd09bSSatish Balay 
337827bd09bSSatish Balay       /* partition about pivot_value ...                              */
338827bd09bSSatish Balay       /* note lists of length 2 are not guaranteed to be sorted */
339db4deed7SKarl Rupp       for (;;) {
340827bd09bSSatish Balay         /* walk up ... and down ... swap if equal to pivot! */
341827bd09bSSatish Balay         do { pi++; pi2++; } while (*pi<*ar);
342827bd09bSSatish Balay         do { pj--; pj2--; } while (*pj>*ar);
343827bd09bSSatish Balay 
344827bd09bSSatish Balay         /* if we've crossed we're done */
345827bd09bSSatish Balay         if (pj<pi) break;
346827bd09bSSatish Balay 
347827bd09bSSatish Balay         /* else swap */
348827bd09bSSatish Balay         SWAP(*pi,*pj)
349827bd09bSSatish Balay         SWAP(*pi2,*pj2)
350827bd09bSSatish Balay       }
351827bd09bSSatish Balay 
352827bd09bSSatish Balay       /* place pivot_value in it's correct location */
353827bd09bSSatish Balay       SWAP(*ar,*pj)
354827bd09bSSatish Balay       SWAP(*ar2,*pj2)
355827bd09bSSatish Balay 
356827bd09bSSatish Balay       /* test stack_size to see if we've exhausted our stack */
357ca8e9878SJed Brown       if (top_s-bottom_s >= SORT_STACK) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_ivec_sort_companion() :: STACK EXHAUSTED!!!");
358827bd09bSSatish Balay 
359827bd09bSSatish Balay       /* push right hand child iff length > 1 */
360db4deed7SKarl Rupp       if ((*top_s = size-((PetscInt) (pi-ar)))) {
361827bd09bSSatish Balay         *(top_a++) = pi;
362827bd09bSSatish Balay         *(top_a++) = pi2;
363827bd09bSSatish Balay         size      -= *top_s+2;
364827bd09bSSatish Balay         top_s++;
3652fa5cd67SKarl Rupp       } else if (size -= *top_s+2) ;   /* set up for next loop iff there is something to do */
366db4deed7SKarl Rupp       else {  /* might as well pop - note NR_OPT >=2 ==> we're ok! */
367827bd09bSSatish Balay         ar2  = *(--top_a);
368827bd09bSSatish Balay         ar   = *(--top_a);
369827bd09bSSatish Balay         size = *(--top_s);
370827bd09bSSatish Balay       }
371db4deed7SKarl Rupp     } else { /* else sort small list directly then pop another off stack */
372827bd09bSSatish Balay 
373827bd09bSSatish Balay       /* insertion sort for bottom */
374db4deed7SKarl Rupp       for (pj=ar+1, pj2=ar2+1; pj<=ar+size; pj++,pj2++) {
375827bd09bSSatish Balay         temp  = *pj;
376827bd09bSSatish Balay         temp2 = *pj2;
377db4deed7SKarl Rupp         for (pi=pj-1,pi2=pj2-1; pi>=ar; pi--,pi2--) {
378827bd09bSSatish Balay           if (*pi <= temp) break;
379827bd09bSSatish Balay           *(pi+1) =*pi;
380827bd09bSSatish Balay           *(pi2+1)=*pi2;
381827bd09bSSatish Balay         }
382827bd09bSSatish Balay         *(pi+1) =temp;
383827bd09bSSatish Balay         *(pi2+1)=temp2;
384827bd09bSSatish Balay       }
385827bd09bSSatish Balay 
386827bd09bSSatish Balay       /* check to see if stack is exhausted ==> DONE */
3873fdc5746SBarry Smith       if (top_s==bottom_s) PetscFunctionReturn(0);
388827bd09bSSatish Balay 
389827bd09bSSatish Balay       /* else pop another list from the stack */
390827bd09bSSatish Balay       ar2  = *(--top_a);
391827bd09bSSatish Balay       ar   = *(--top_a);
392827bd09bSSatish Balay       size = *(--top_s);
393827bd09bSSatish Balay     }
394827bd09bSSatish Balay   }
395827bd09bSSatish Balay }
396827bd09bSSatish Balay 
3977b1ae94cSBarry Smith /******************************************************************************/
398ca8e9878SJed Brown PetscErrorCode PCTFS_ivec_sort_companion_hack(PetscInt *ar,  PetscInt **ar2, PetscInt size)
399827bd09bSSatish Balay {
40052f87cdaSBarry Smith   PetscInt *pi, *pj, temp, *ptr;
40152f87cdaSBarry Smith   PetscInt **top_a = (PetscInt**)offset_stack;
40252f87cdaSBarry Smith   PetscInt *top_s  = size_stack, *bottom_s = size_stack;
40352f87cdaSBarry Smith   PetscInt **pi2, **pj2;
40452f87cdaSBarry Smith   PetscInt mid;
405827bd09bSSatish Balay 
4063fdc5746SBarry Smith   PetscFunctionBegin;
407827bd09bSSatish Balay   /* we're really interested in the offset of the last element */
408827bd09bSSatish Balay   /* ==> length of the list is now size + 1                    */
409827bd09bSSatish Balay   size--;
410827bd09bSSatish Balay 
411827bd09bSSatish Balay   /* do until we're done ... return when stack is exhausted */
412db4deed7SKarl Rupp   for (;;) {
413db4deed7SKarl Rupp 
414827bd09bSSatish Balay     /* if list is large enough use quicksort partition exchange code */
415db4deed7SKarl Rupp     if (size > SORT_OPT) {
416db4deed7SKarl Rupp 
417827bd09bSSatish Balay       /* start up pointer at element 1 and down at size     */
418827bd09bSSatish Balay       mid = size>>1;
419827bd09bSSatish Balay       pi  = ar+1;
420827bd09bSSatish Balay       pj  = ar+mid;
421827bd09bSSatish Balay       pi2 = ar2+1;
422827bd09bSSatish Balay       pj2 = ar2+mid;
423827bd09bSSatish Balay 
424827bd09bSSatish Balay       /* find middle element in list and swap w/ element 1 */
425827bd09bSSatish Balay       SWAP(*pi,*pj)
426827bd09bSSatish Balay       P_SWAP(*pi2,*pj2)
427827bd09bSSatish Balay 
428827bd09bSSatish Balay       /* order element 0,1,size-1 st {M,L,...,U} w/L<=M<=U */
429827bd09bSSatish Balay       /* note ==> pivot_value in index 0                   */
430827bd09bSSatish Balay       pj  = ar+size;
431827bd09bSSatish Balay       pj2 = ar2+size;
432db4deed7SKarl Rupp       if (*pi > *pj) { SWAP(*pi,*pj) P_SWAP(*pi2,*pj2) }
433db4deed7SKarl Rupp       if (*ar > *pj) { SWAP(*ar,*pj) P_SWAP(*ar2,*pj2) }
434db4deed7SKarl Rupp       else if (*pi > *ar) { SWAP(*(ar),*(ar+1)) P_SWAP(*(ar2),*(ar2+1)) }
435827bd09bSSatish Balay 
436827bd09bSSatish Balay       /* partition about pivot_value ...                              */
437827bd09bSSatish Balay       /* note lists of length 2 are not guaranteed to be sorted */
438db4deed7SKarl Rupp       for (;;) {
439db4deed7SKarl Rupp 
440827bd09bSSatish Balay         /* walk up ... and down ... swap if equal to pivot! */
441827bd09bSSatish Balay         do {pi++; pi2++;} while (*pi<*ar);
442827bd09bSSatish Balay         do {pj--; pj2--;} while (*pj>*ar);
443827bd09bSSatish Balay 
444827bd09bSSatish Balay         /* if we've crossed we're done */
445827bd09bSSatish Balay         if (pj<pi) break;
446827bd09bSSatish Balay 
447827bd09bSSatish Balay         /* else swap */
448827bd09bSSatish Balay         SWAP(*pi,*pj)
449827bd09bSSatish Balay         P_SWAP(*pi2,*pj2)
450827bd09bSSatish Balay       }
451827bd09bSSatish Balay 
452827bd09bSSatish Balay       /* place pivot_value in it's correct location */
453827bd09bSSatish Balay       SWAP(*ar,*pj)
454827bd09bSSatish Balay       P_SWAP(*ar2,*pj2)
455827bd09bSSatish Balay 
456827bd09bSSatish Balay       /* test stack_size to see if we've exhausted our stack */
457ca8e9878SJed Brown       if (top_s-bottom_s >= SORT_STACK) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_ivec_sort_companion_hack() :: STACK EXHAUSTED!!!");
458827bd09bSSatish Balay 
459827bd09bSSatish Balay       /* push right hand child iff length > 1 */
460db4deed7SKarl Rupp       if ((*top_s = size-((PetscInt) (pi-ar)))) {
461827bd09bSSatish Balay         *(top_a++) = pi;
46252f87cdaSBarry Smith         *(top_a++) = (PetscInt*) pi2;
463827bd09bSSatish Balay         size      -= *top_s+2;
464827bd09bSSatish Balay         top_s++;
4652fa5cd67SKarl Rupp       } else if (size -= *top_s+2) ;   /* set up for next loop iff there is something to do */
466db4deed7SKarl Rupp       else { /* might as well pop - note NR_OPT >=2 ==> we're ok! */
46752f87cdaSBarry Smith         ar2  = (PetscInt**) *(--top_a);
468827bd09bSSatish Balay         ar   = *(--top_a);
469827bd09bSSatish Balay         size = *(--top_s);
470827bd09bSSatish Balay       }
4712fa5cd67SKarl Rupp     } else  { /* else sort small list directly then pop another off stack */
472827bd09bSSatish Balay       /* insertion sort for bottom */
473db4deed7SKarl Rupp       for (pj=ar+1, pj2=ar2+1; pj<=ar+size; pj++,pj2++) {
474827bd09bSSatish Balay         temp = *pj;
475827bd09bSSatish Balay         ptr  = *pj2;
476db4deed7SKarl Rupp         for (pi=pj-1,pi2=pj2-1; pi>=ar; pi--,pi2--) {
477827bd09bSSatish Balay           if (*pi <= temp) break;
478827bd09bSSatish Balay           *(pi+1) =*pi;
479827bd09bSSatish Balay           *(pi2+1)=*pi2;
480827bd09bSSatish Balay         }
481827bd09bSSatish Balay         *(pi+1) =temp;
482827bd09bSSatish Balay         *(pi2+1)=ptr;
483827bd09bSSatish Balay       }
484827bd09bSSatish Balay 
485827bd09bSSatish Balay       /* check to see if stack is exhausted ==> DONE */
4863fdc5746SBarry Smith       if (top_s==bottom_s) PetscFunctionReturn(0);
487827bd09bSSatish Balay 
488827bd09bSSatish Balay       /* else pop another list from the stack */
48952f87cdaSBarry Smith       ar2  = (PetscInt**)*(--top_a);
490827bd09bSSatish Balay       ar   = *(--top_a);
491827bd09bSSatish Balay       size = *(--top_s);
492827bd09bSSatish Balay     }
493827bd09bSSatish Balay   }
494827bd09bSSatish Balay }
495827bd09bSSatish Balay 
4967b1ae94cSBarry Smith /******************************************************************************/
497ca8e9878SJed Brown PetscErrorCode PCTFS_SMI_sort(void *ar1, void *ar2, PetscInt size, PetscInt type)
498827bd09bSSatish Balay {
4993fdc5746SBarry Smith   PetscFunctionBegin;
500e7e72b3dSBarry Smith   if (type == SORT_INTEGER) {
501ca8e9878SJed Brown     if (ar2) PCTFS_ivec_sort_companion((PetscInt*)ar1,(PetscInt*)ar2,size);
502ca8e9878SJed Brown     else PCTFS_ivec_sort((PetscInt*)ar1,size);
503e7e72b3dSBarry Smith   } else if (type == SORT_INT_PTR) {
504ca8e9878SJed Brown     if (ar2) PCTFS_ivec_sort_companion_hack((PetscInt*)ar1,(PetscInt**)ar2,size);
505ca8e9878SJed Brown     else PCTFS_ivec_sort((PetscInt*)ar1,size);
506ca8e9878SJed Brown   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_SMI_sort only does SORT_INTEGER!");
5073fdc5746SBarry Smith   PetscFunctionReturn(0);
508827bd09bSSatish Balay }
509827bd09bSSatish Balay 
5107b1ae94cSBarry Smith /***********************************ivec.c*************************************/
511ca8e9878SJed Brown PetscInt PCTFS_ivec_linear_search(PetscInt item,  PetscInt *list,  PetscInt n)
512827bd09bSSatish Balay {
51352f87cdaSBarry Smith   PetscInt tmp = n-1;
5145fd66863SKarl Rupp 
5153fdc5746SBarry Smith   PetscFunctionBegin;
5162fa5cd67SKarl Rupp   while (n--) {
5172fa5cd67SKarl Rupp     if (*list++ == item) return(tmp-n);
5182fa5cd67SKarl Rupp   }
519827bd09bSSatish Balay   return(-1);
520827bd09bSSatish Balay }
521827bd09bSSatish Balay 
5227b1ae94cSBarry Smith /***********************************ivec.c*************************************/
523ca8e9878SJed Brown PetscInt PCTFS_ivec_binary_search(PetscInt item,  PetscInt *list,  PetscInt rh)
524827bd09bSSatish Balay {
52552f87cdaSBarry Smith   PetscInt mid, lh=0;
526827bd09bSSatish Balay 
527827bd09bSSatish Balay   rh--;
528db4deed7SKarl Rupp   while (lh<=rh) {
529827bd09bSSatish Balay     mid = (lh+rh)>>1;
5302fa5cd67SKarl Rupp     if (*(list+mid) == item) return(mid);
5312fa5cd67SKarl Rupp     if (*(list+mid) > item) rh = mid-1;
5322fa5cd67SKarl Rupp     else lh = mid+1;
533827bd09bSSatish Balay   }
534827bd09bSSatish Balay   return(-1);
535827bd09bSSatish Balay }
536827bd09bSSatish Balay 
5377b1ae94cSBarry Smith /*********************************ivec.c*************************************/
538ca8e9878SJed Brown PetscErrorCode PCTFS_rvec_copy(PetscScalar *arg1,  PetscScalar *arg2,  PetscInt n)
539827bd09bSSatish Balay {
5403fdc5746SBarry Smith   PetscFunctionBegin;
5412fa5cd67SKarl Rupp   while (n--) *arg1++ = *arg2++;
5423fdc5746SBarry Smith   PetscFunctionReturn(0);
543827bd09bSSatish Balay }
544827bd09bSSatish Balay 
5457b1ae94cSBarry Smith /*********************************ivec.c*************************************/
546ca8e9878SJed Brown PetscErrorCode PCTFS_rvec_zero(PetscScalar *arg1,  PetscInt n)
547827bd09bSSatish Balay {
5483fdc5746SBarry Smith   PetscFunctionBegin;
5492fa5cd67SKarl Rupp   while (n--) *arg1++ = 0.0;
5503fdc5746SBarry Smith   PetscFunctionReturn(0);
551827bd09bSSatish Balay }
552827bd09bSSatish Balay 
5537b1ae94cSBarry Smith /***********************************ivec.c*************************************/
554ca8e9878SJed Brown PetscErrorCode PCTFS_rvec_one(PetscScalar *arg1,  PetscInt n)
555827bd09bSSatish Balay {
5563fdc5746SBarry Smith   PetscFunctionBegin;
5572fa5cd67SKarl Rupp   while (n--) *arg1++ = 1.0;
5583fdc5746SBarry Smith   PetscFunctionReturn(0);
559827bd09bSSatish Balay }
560827bd09bSSatish Balay 
5617b1ae94cSBarry Smith /***********************************ivec.c*************************************/
562ca8e9878SJed Brown PetscErrorCode PCTFS_rvec_set(PetscScalar *arg1,  PetscScalar arg2,  PetscInt n)
563827bd09bSSatish Balay {
5643fdc5746SBarry Smith   PetscFunctionBegin;
5652fa5cd67SKarl Rupp   while (n--) *arg1++ = arg2;
5663fdc5746SBarry Smith   PetscFunctionReturn(0);
567827bd09bSSatish Balay }
568827bd09bSSatish Balay 
5697b1ae94cSBarry Smith /***********************************ivec.c*************************************/
570ca8e9878SJed Brown PetscErrorCode PCTFS_rvec_scale(PetscScalar *arg1,  PetscScalar arg2,  PetscInt n)
571827bd09bSSatish Balay {
5723fdc5746SBarry Smith   PetscFunctionBegin;
5732fa5cd67SKarl Rupp   while (n--) *arg1++ *= arg2;
5743fdc5746SBarry Smith   PetscFunctionReturn(0);
575827bd09bSSatish Balay }
576827bd09bSSatish Balay 
5777b1ae94cSBarry Smith /*********************************ivec.c*************************************/
578ca8e9878SJed Brown PetscErrorCode PCTFS_rvec_add(PetscScalar *arg1,  PetscScalar *arg2,  PetscInt n)
579827bd09bSSatish Balay {
5803fdc5746SBarry Smith   PetscFunctionBegin;
5812fa5cd67SKarl Rupp   while (n--) *arg1++ += *arg2++;
5823fdc5746SBarry Smith   PetscFunctionReturn(0);
583827bd09bSSatish Balay }
584827bd09bSSatish Balay 
5857b1ae94cSBarry Smith /*********************************ivec.c*************************************/
586ca8e9878SJed Brown PetscErrorCode PCTFS_rvec_mult(PetscScalar *arg1,  PetscScalar *arg2,  PetscInt n)
587827bd09bSSatish Balay {
5883fdc5746SBarry Smith   PetscFunctionBegin;
5892fa5cd67SKarl Rupp   while (n--) *arg1++ *= *arg2++;
5903fdc5746SBarry Smith   PetscFunctionReturn(0);
591827bd09bSSatish Balay }
592827bd09bSSatish Balay 
5937b1ae94cSBarry Smith /*********************************ivec.c*************************************/
594ca8e9878SJed Brown PetscErrorCode PCTFS_rvec_max(PetscScalar *arg1,  PetscScalar *arg2,  PetscInt n)
595827bd09bSSatish Balay {
5963fdc5746SBarry Smith   PetscFunctionBegin;
5972fa5cd67SKarl Rupp   while (n--) {
5982fa5cd67SKarl Rupp     *arg1 = PetscMax(*arg1,*arg2);
5992fa5cd67SKarl Rupp     arg1++;
6002fa5cd67SKarl Rupp     arg2++;
6012fa5cd67SKarl Rupp   }
6023fdc5746SBarry Smith   PetscFunctionReturn(0);
603827bd09bSSatish Balay }
604827bd09bSSatish Balay 
6057b1ae94cSBarry Smith /*********************************ivec.c*************************************/
606ca8e9878SJed Brown PetscErrorCode PCTFS_rvec_max_abs(PetscScalar *arg1,  PetscScalar *arg2,  PetscInt n)
607827bd09bSSatish Balay {
6083fdc5746SBarry Smith   PetscFunctionBegin;
6092fa5cd67SKarl Rupp   while (n--) {
6102fa5cd67SKarl Rupp     *arg1 = MAX_FABS(*arg1,*arg2);
6112fa5cd67SKarl Rupp     arg1++;
6122fa5cd67SKarl Rupp     arg2++;
6132fa5cd67SKarl Rupp   }
6143fdc5746SBarry Smith   PetscFunctionReturn(0);
615827bd09bSSatish Balay }
616827bd09bSSatish Balay 
6177b1ae94cSBarry Smith /*********************************ivec.c*************************************/
618ca8e9878SJed Brown PetscErrorCode PCTFS_rvec_min(PetscScalar *arg1,  PetscScalar *arg2,  PetscInt n)
619827bd09bSSatish Balay {
6203fdc5746SBarry Smith   PetscFunctionBegin;
6212fa5cd67SKarl Rupp   while (n--) {
6222fa5cd67SKarl Rupp     *arg1 = PetscMin(*arg1,*arg2);
6232fa5cd67SKarl Rupp     arg1++;
6242fa5cd67SKarl Rupp     arg2++;
6252fa5cd67SKarl Rupp   }
6263fdc5746SBarry Smith   PetscFunctionReturn(0);
627827bd09bSSatish Balay }
628827bd09bSSatish Balay 
6297b1ae94cSBarry Smith /*********************************ivec.c*************************************/
630ca8e9878SJed Brown PetscErrorCode PCTFS_rvec_min_abs(PetscScalar *arg1,  PetscScalar *arg2,  PetscInt n)
631827bd09bSSatish Balay {
6323fdc5746SBarry Smith   PetscFunctionBegin;
6332fa5cd67SKarl Rupp   while (n--) {
6342fa5cd67SKarl Rupp     *arg1 = MIN_FABS(*arg1,*arg2);
6352fa5cd67SKarl Rupp     arg1++;
6362fa5cd67SKarl Rupp     arg2++;
6372fa5cd67SKarl Rupp   }
6383fdc5746SBarry Smith   PetscFunctionReturn(0);
639827bd09bSSatish Balay }
640827bd09bSSatish Balay 
6417b1ae94cSBarry Smith /*********************************ivec.c*************************************/
642ca8e9878SJed Brown PetscErrorCode PCTFS_rvec_exists(PetscScalar *arg1,  PetscScalar *arg2,  PetscInt n)
643827bd09bSSatish Balay {
6443fdc5746SBarry Smith   PetscFunctionBegin;
6452fa5cd67SKarl Rupp   while (n--) {
6462fa5cd67SKarl Rupp     *arg1 = EXISTS(*arg1,*arg2);
6472fa5cd67SKarl Rupp     arg1++;
6482fa5cd67SKarl Rupp     arg2++;
6492fa5cd67SKarl Rupp   }
6503fdc5746SBarry Smith   PetscFunctionReturn(0);
651827bd09bSSatish Balay }
652827bd09bSSatish Balay 
6537b1ae94cSBarry Smith /***********************************ivec.c*************************************/
654ca8e9878SJed Brown PetscErrorCode PCTFS_rvec_non_uniform(PetscScalar *arg1, PetscScalar *arg2,  PetscInt n,  PetscInt *arg3)
655827bd09bSSatish Balay {
65652f87cdaSBarry Smith   PetscInt i, j, type;
657827bd09bSSatish Balay 
6583fdc5746SBarry Smith   PetscFunctionBegin;
659827bd09bSSatish Balay   /* LATER: if we're really motivated we can sort and then unsort */
660db4deed7SKarl Rupp   for (i=0; i<n;) {
661db4deed7SKarl Rupp 
662827bd09bSSatish Balay     /* clump 'em for now */
663827bd09bSSatish Balay     j    =i+1;
664827bd09bSSatish Balay     type = arg3[i];
6652fa5cd67SKarl Rupp     while ((j<n)&&(arg3[j]==type)) j++;
666827bd09bSSatish Balay 
667827bd09bSSatish Balay     /* how many together */
668827bd09bSSatish Balay     j -= i;
669827bd09bSSatish Balay 
670827bd09bSSatish Balay     /* call appropriate ivec function */
6712fa5cd67SKarl Rupp     if (type == GL_MAX)          PCTFS_rvec_max(arg1,arg2,j);
6722fa5cd67SKarl Rupp     else if (type == GL_MIN)     PCTFS_rvec_min(arg1,arg2,j);
6732fa5cd67SKarl Rupp     else if (type == GL_MULT)    PCTFS_rvec_mult(arg1,arg2,j);
6742fa5cd67SKarl Rupp     else if (type == GL_ADD)     PCTFS_rvec_add(arg1,arg2,j);
6752fa5cd67SKarl Rupp     else if (type == GL_MAX_ABS) PCTFS_rvec_max_abs(arg1,arg2,j);
6762fa5cd67SKarl Rupp     else if (type == GL_MIN_ABS) PCTFS_rvec_min_abs(arg1,arg2,j);
6772fa5cd67SKarl Rupp     else if (type == GL_EXISTS)  PCTFS_rvec_exists(arg1,arg2,j);
678ca8e9878SJed Brown     else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"unrecognized type passed to PCTFS_rvec_non_uniform()!");
679827bd09bSSatish Balay 
680827bd09bSSatish Balay     arg1+=j; arg2+=j; i+=j;
681827bd09bSSatish Balay   }
6823fdc5746SBarry Smith   PetscFunctionReturn(0);
683827bd09bSSatish Balay }
684827bd09bSSatish Balay 
6857b1ae94cSBarry Smith /***********************************ivec.c*************************************/
686ca8e9878SJed Brown vfp PCTFS_rvec_fct_addr(PetscInt type)
687827bd09bSSatish Balay {
6882fa5cd67SKarl Rupp   if (type == NON_UNIFORM)     return((PetscErrorCode (*)(void*, void*, PetscInt, ...))&PCTFS_rvec_non_uniform);
6892fa5cd67SKarl Rupp   else if (type == GL_MAX)     return((PetscErrorCode (*)(void*, void*, PetscInt, ...))&PCTFS_rvec_max);
6902fa5cd67SKarl Rupp   else if (type == GL_MIN)     return((PetscErrorCode (*)(void*, void*, PetscInt, ...))&PCTFS_rvec_min);
6912fa5cd67SKarl Rupp   else if (type == GL_MULT)    return((PetscErrorCode (*)(void*, void*, PetscInt, ...))&PCTFS_rvec_mult);
6922fa5cd67SKarl Rupp   else if (type == GL_ADD)     return((PetscErrorCode (*)(void*, void*, PetscInt, ...))&PCTFS_rvec_add);
6932fa5cd67SKarl Rupp   else if (type == GL_MAX_ABS) return((PetscErrorCode (*)(void*, void*, PetscInt, ...))&PCTFS_rvec_max_abs);
6942fa5cd67SKarl Rupp   else if (type == GL_MIN_ABS) return((PetscErrorCode (*)(void*, void*, PetscInt, ...))&PCTFS_rvec_min_abs);
6952fa5cd67SKarl Rupp   else if (type == GL_EXISTS)  return((PetscErrorCode (*)(void*, void*, PetscInt, ...))&PCTFS_rvec_exists);
696827bd09bSSatish Balay 
697827bd09bSSatish Balay   /* catch all ... not good if we get here */
698827bd09bSSatish Balay   return(NULL);
699827bd09bSSatish Balay }
700827bd09bSSatish Balay 
701827bd09bSSatish Balay 
702827bd09bSSatish Balay 
703827bd09bSSatish Balay 
704827bd09bSSatish Balay 
705827bd09bSSatish Balay 
706