xref: /petsc/src/ksp/pc/impls/tfs/ivec.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1827bd09bSSatish Balay 
2827bd09bSSatish Balay /**********************************ivec.c**************************************
3827bd09bSSatish Balay 
4827bd09bSSatish Balay Author: Henry M. Tufo III
5827bd09bSSatish Balay 
6827bd09bSSatish Balay e-mail: hmt@cs.brown.edu
7827bd09bSSatish Balay 
8827bd09bSSatish Balay snail-mail:
9827bd09bSSatish Balay Division of Applied Mathematics
10827bd09bSSatish Balay Brown University
11827bd09bSSatish Balay Providence, RI 02912
12827bd09bSSatish Balay 
13827bd09bSSatish Balay Last Modification:
14827bd09bSSatish Balay 6.21.97
15827bd09bSSatish Balay ***********************************ivec.c*************************************/
16827bd09bSSatish Balay 
17c6db04a5SJed Brown #include <../src/ksp/pc/impls/tfs/tfs.h>
18827bd09bSSatish Balay 
19827bd09bSSatish Balay /* sorting args ivec.c ivec.c ... */
20827bd09bSSatish Balay #define SORT_OPT   6
21827bd09bSSatish Balay #define SORT_STACK 50000
22827bd09bSSatish Balay 
23827bd09bSSatish Balay /* allocate an address and size stack for sorter(s) */
24827bd09bSSatish Balay static void    *offset_stack[2 * SORT_STACK];
2552f87cdaSBarry Smith static PetscInt size_stack[SORT_STACK];
26827bd09bSSatish Balay 
277b1ae94cSBarry Smith /***********************************ivec.c*************************************/
28*9371c9d4SSatish Balay PetscInt *PCTFS_ivec_copy(PetscInt *arg1, PetscInt *arg2, PetscInt n) {
292fa5cd67SKarl Rupp   while (n--) *arg1++ = *arg2++;
30827bd09bSSatish Balay   return (arg1);
31827bd09bSSatish Balay }
32827bd09bSSatish Balay 
337b1ae94cSBarry Smith /***********************************ivec.c*************************************/
34*9371c9d4SSatish Balay PetscErrorCode PCTFS_ivec_zero(PetscInt *arg1, PetscInt n) {
353fdc5746SBarry Smith   PetscFunctionBegin;
362fa5cd67SKarl Rupp   while (n--) *arg1++ = 0;
373fdc5746SBarry Smith   PetscFunctionReturn(0);
38827bd09bSSatish Balay }
39827bd09bSSatish Balay 
407b1ae94cSBarry Smith /***********************************ivec.c*************************************/
41*9371c9d4SSatish Balay PetscErrorCode PCTFS_ivec_set(PetscInt *arg1, PetscInt arg2, PetscInt n) {
423fdc5746SBarry Smith   PetscFunctionBegin;
432fa5cd67SKarl Rupp   while (n--) *arg1++ = arg2;
443fdc5746SBarry Smith   PetscFunctionReturn(0);
45827bd09bSSatish Balay }
46827bd09bSSatish Balay 
477b1ae94cSBarry Smith /***********************************ivec.c*************************************/
48*9371c9d4SSatish Balay PetscErrorCode PCTFS_ivec_max(PetscInt *arg1, PetscInt *arg2, PetscInt n) {
493fdc5746SBarry Smith   PetscFunctionBegin;
50*9371c9d4SSatish Balay   while (n--) {
51*9371c9d4SSatish Balay     *arg1 = PetscMax(*arg1, *arg2);
52*9371c9d4SSatish Balay     arg1++;
53*9371c9d4SSatish Balay     arg2++;
54*9371c9d4SSatish Balay   }
553fdc5746SBarry Smith   PetscFunctionReturn(0);
56827bd09bSSatish Balay }
57827bd09bSSatish Balay 
587b1ae94cSBarry Smith /***********************************ivec.c*************************************/
59*9371c9d4SSatish Balay PetscErrorCode PCTFS_ivec_min(PetscInt *arg1, PetscInt *arg2, PetscInt n) {
603fdc5746SBarry Smith   PetscFunctionBegin;
612fa5cd67SKarl Rupp   while (n--) {
622fa5cd67SKarl Rupp     *(arg1) = PetscMin(*arg1, *arg2);
632fa5cd67SKarl Rupp     arg1++;
642fa5cd67SKarl Rupp     arg2++;
652fa5cd67SKarl Rupp   }
663fdc5746SBarry Smith   PetscFunctionReturn(0);
67827bd09bSSatish Balay }
68827bd09bSSatish Balay 
697b1ae94cSBarry Smith /***********************************ivec.c*************************************/
70*9371c9d4SSatish Balay PetscErrorCode PCTFS_ivec_mult(PetscInt *arg1, PetscInt *arg2, PetscInt n) {
713fdc5746SBarry Smith   PetscFunctionBegin;
722fa5cd67SKarl Rupp   while (n--) *arg1++ *= *arg2++;
733fdc5746SBarry Smith   PetscFunctionReturn(0);
74827bd09bSSatish Balay }
75827bd09bSSatish Balay 
767b1ae94cSBarry Smith /***********************************ivec.c*************************************/
77*9371c9d4SSatish Balay PetscErrorCode PCTFS_ivec_add(PetscInt *arg1, PetscInt *arg2, PetscInt n) {
783fdc5746SBarry Smith   PetscFunctionBegin;
792fa5cd67SKarl Rupp   while (n--) *arg1++ += *arg2++;
803fdc5746SBarry Smith   PetscFunctionReturn(0);
81827bd09bSSatish Balay }
82827bd09bSSatish Balay 
837b1ae94cSBarry Smith /***********************************ivec.c*************************************/
84*9371c9d4SSatish Balay PetscErrorCode PCTFS_ivec_lxor(PetscInt *arg1, PetscInt *arg2, PetscInt n) {
853fdc5746SBarry Smith   PetscFunctionBegin;
862fa5cd67SKarl Rupp   while (n--) {
872fa5cd67SKarl Rupp     *arg1 = ((*arg1 || *arg2) && !(*arg1 && *arg2));
882fa5cd67SKarl Rupp     arg1++;
892fa5cd67SKarl Rupp     arg2++;
902fa5cd67SKarl Rupp   }
913fdc5746SBarry Smith   PetscFunctionReturn(0);
92827bd09bSSatish Balay }
93827bd09bSSatish Balay 
947b1ae94cSBarry Smith /***********************************ivec.c*************************************/
95*9371c9d4SSatish Balay PetscErrorCode PCTFS_ivec_xor(PetscInt *arg1, PetscInt *arg2, PetscInt n) {
963fdc5746SBarry Smith   PetscFunctionBegin;
972fa5cd67SKarl Rupp   while (n--) *arg1++ ^= *arg2++;
983fdc5746SBarry Smith   PetscFunctionReturn(0);
99827bd09bSSatish Balay }
100827bd09bSSatish Balay 
1017b1ae94cSBarry Smith /***********************************ivec.c*************************************/
102*9371c9d4SSatish Balay PetscErrorCode PCTFS_ivec_or(PetscInt *arg1, PetscInt *arg2, PetscInt n) {
1033fdc5746SBarry Smith   PetscFunctionBegin;
1042fa5cd67SKarl Rupp   while (n--) *arg1++ |= *arg2++;
1053fdc5746SBarry Smith   PetscFunctionReturn(0);
106827bd09bSSatish Balay }
107827bd09bSSatish Balay 
1087b1ae94cSBarry Smith /***********************************ivec.c*************************************/
109*9371c9d4SSatish Balay PetscErrorCode PCTFS_ivec_lor(PetscInt *arg1, PetscInt *arg2, PetscInt n) {
1103fdc5746SBarry Smith   PetscFunctionBegin;
1112fa5cd67SKarl Rupp   while (n--) {
1122fa5cd67SKarl Rupp     *arg1 = (*arg1 || *arg2);
1132fa5cd67SKarl Rupp     arg1++;
1142fa5cd67SKarl Rupp     arg2++;
1152fa5cd67SKarl Rupp   }
1163fdc5746SBarry Smith   PetscFunctionReturn(0);
117827bd09bSSatish Balay }
118827bd09bSSatish Balay 
1197b1ae94cSBarry Smith /***********************************ivec.c*************************************/
120*9371c9d4SSatish Balay PetscErrorCode PCTFS_ivec_and(PetscInt *arg1, PetscInt *arg2, PetscInt n) {
1213fdc5746SBarry Smith   PetscFunctionBegin;
1222fa5cd67SKarl Rupp   while (n--) *arg1++ &= *arg2++;
1233fdc5746SBarry Smith   PetscFunctionReturn(0);
124827bd09bSSatish Balay }
125827bd09bSSatish Balay 
1267b1ae94cSBarry Smith /***********************************ivec.c*************************************/
127*9371c9d4SSatish Balay PetscErrorCode PCTFS_ivec_land(PetscInt *arg1, PetscInt *arg2, PetscInt n) {
1283fdc5746SBarry Smith   PetscFunctionBegin;
1292fa5cd67SKarl Rupp   while (n--) {
1302fa5cd67SKarl Rupp     *arg1 = (*arg1 && *arg2);
1312fa5cd67SKarl Rupp     arg1++;
1322fa5cd67SKarl Rupp     arg2++;
1332fa5cd67SKarl Rupp   }
1343fdc5746SBarry Smith   PetscFunctionReturn(0);
135827bd09bSSatish Balay }
136827bd09bSSatish Balay 
1377b1ae94cSBarry Smith /***********************************ivec.c*************************************/
138*9371c9d4SSatish Balay PetscErrorCode PCTFS_ivec_and3(PetscInt *arg1, PetscInt *arg2, PetscInt *arg3, PetscInt n) {
1393fdc5746SBarry Smith   PetscFunctionBegin;
1402fa5cd67SKarl Rupp   while (n--) *arg1++ = (*arg2++ & *arg3++);
1413fdc5746SBarry Smith   PetscFunctionReturn(0);
142827bd09bSSatish Balay }
143827bd09bSSatish Balay 
1447b1ae94cSBarry Smith /***********************************ivec.c*************************************/
145*9371c9d4SSatish Balay PetscInt PCTFS_ivec_sum(PetscInt *arg1, PetscInt n) {
14652f87cdaSBarry Smith   PetscInt tmp = 0;
1472fa5cd67SKarl Rupp   while (n--) tmp += *arg1++;
148827bd09bSSatish Balay   return (tmp);
149827bd09bSSatish Balay }
150827bd09bSSatish Balay 
1517b1ae94cSBarry Smith /***********************************ivec.c*************************************/
152*9371c9d4SSatish Balay PetscErrorCode PCTFS_ivec_non_uniform(PetscInt *arg1, PetscInt *arg2, PetscInt n, ...) {
15352f87cdaSBarry Smith   PetscInt  i, j, type;
1546ae68b59SBarry Smith   PetscInt *arg3;
1556ae68b59SBarry Smith   va_list   ap;
156827bd09bSSatish Balay 
1573fdc5746SBarry Smith   PetscFunctionBegin;
1586ae68b59SBarry Smith   va_start(ap, n);
1596ae68b59SBarry Smith   arg3 = va_arg(ap, PetscInt *);
1606ae68b59SBarry Smith   va_end(ap);
1616ae68b59SBarry Smith 
162827bd09bSSatish Balay   /* LATER: if we're really motivated we can sort and then unsort */
163db4deed7SKarl Rupp   for (i = 0; i < n;) {
164827bd09bSSatish Balay     /* clump 'em for now */
165827bd09bSSatish Balay     j    = i + 1;
166827bd09bSSatish Balay     type = arg3[i];
1672fa5cd67SKarl Rupp     while ((j < n) && (arg3[j] == type)) j++;
168827bd09bSSatish Balay 
169827bd09bSSatish Balay     /* how many together */
170827bd09bSSatish Balay     j -= i;
171827bd09bSSatish Balay 
172827bd09bSSatish Balay     /* call appropriate ivec function */
1732fa5cd67SKarl Rupp     if (type == GL_MAX) PCTFS_ivec_max(arg1, arg2, j);
1742fa5cd67SKarl Rupp     else if (type == GL_MIN) PCTFS_ivec_min(arg1, arg2, j);
1752fa5cd67SKarl Rupp     else if (type == GL_MULT) PCTFS_ivec_mult(arg1, arg2, j);
1762fa5cd67SKarl Rupp     else if (type == GL_ADD) PCTFS_ivec_add(arg1, arg2, j);
1772fa5cd67SKarl Rupp     else if (type == GL_B_XOR) PCTFS_ivec_xor(arg1, arg2, j);
1782fa5cd67SKarl Rupp     else if (type == GL_B_OR) PCTFS_ivec_or(arg1, arg2, j);
1792fa5cd67SKarl Rupp     else if (type == GL_B_AND) PCTFS_ivec_and(arg1, arg2, j);
1802fa5cd67SKarl Rupp     else if (type == GL_L_XOR) PCTFS_ivec_lxor(arg1, arg2, j);
1812fa5cd67SKarl Rupp     else if (type == GL_L_OR) PCTFS_ivec_lor(arg1, arg2, j);
1822fa5cd67SKarl Rupp     else if (type == GL_L_AND) PCTFS_ivec_land(arg1, arg2, j);
183db4deed7SKarl Rupp     else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "unrecognized type passed to PCTFS_ivec_non_uniform()!");
184827bd09bSSatish Balay 
185*9371c9d4SSatish Balay     arg1 += j;
186*9371c9d4SSatish Balay     arg2 += j;
187*9371c9d4SSatish Balay     i += j;
188827bd09bSSatish Balay   }
1893fdc5746SBarry Smith   PetscFunctionReturn(0);
190827bd09bSSatish Balay }
191827bd09bSSatish Balay 
1927b1ae94cSBarry Smith /***********************************ivec.c*************************************/
193*9371c9d4SSatish Balay vfp PCTFS_ivec_fct_addr(PetscInt type) {
1942fa5cd67SKarl Rupp   if (type == NON_UNIFORM) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_non_uniform);
1952fa5cd67SKarl Rupp   else if (type == GL_MAX) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_max);
1962fa5cd67SKarl Rupp   else if (type == GL_MIN) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_min);
1972fa5cd67SKarl Rupp   else if (type == GL_MULT) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_mult);
1982fa5cd67SKarl Rupp   else if (type == GL_ADD) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_add);
1992fa5cd67SKarl Rupp   else if (type == GL_B_XOR) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_xor);
2002fa5cd67SKarl Rupp   else if (type == GL_B_OR) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_or);
2012fa5cd67SKarl Rupp   else if (type == GL_B_AND) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_and);
2022fa5cd67SKarl Rupp   else if (type == GL_L_XOR) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_lxor);
2032fa5cd67SKarl Rupp   else if (type == GL_L_OR) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_lor);
2042fa5cd67SKarl Rupp   else if (type == GL_L_AND) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_land);
205827bd09bSSatish Balay 
206827bd09bSSatish Balay   /* catch all ... not good if we get here */
207827bd09bSSatish Balay   return (NULL);
208827bd09bSSatish Balay }
209827bd09bSSatish Balay 
2107b1ae94cSBarry Smith /******************************************************************************/
211*9371c9d4SSatish Balay PetscErrorCode PCTFS_ivec_sort(PetscInt *ar, PetscInt size) {
21252f87cdaSBarry Smith   PetscInt  *pi, *pj, temp;
21352f87cdaSBarry Smith   PetscInt **top_a = (PetscInt **)offset_stack;
21452f87cdaSBarry Smith   PetscInt  *top_s = size_stack, *bottom_s = size_stack;
215827bd09bSSatish Balay 
216b458e8f1SJose E. Roman   PetscFunctionBegin;
217827bd09bSSatish Balay   /* we're really interested in the offset of the last element */
218827bd09bSSatish Balay   /* ==> length of the list is now size + 1                    */
219827bd09bSSatish Balay   size--;
220827bd09bSSatish Balay 
221827bd09bSSatish Balay   /* do until we're done ... return when stack is exhausted */
222db4deed7SKarl Rupp   for (;;) {
223827bd09bSSatish Balay     /* if list is large enough use quicksort partition exchange code */
224db4deed7SKarl Rupp     if (size > SORT_OPT) {
225827bd09bSSatish Balay       /* start up pointer at element 1 and down at size     */
226827bd09bSSatish Balay       pi = ar + 1;
227827bd09bSSatish Balay       pj = ar + size;
228827bd09bSSatish Balay 
229827bd09bSSatish Balay       /* find middle element in list and swap w/ element 1 */
230827bd09bSSatish Balay       SWAP(*(ar + (size >> 1)), *pi)
231827bd09bSSatish Balay 
232827bd09bSSatish Balay       /* order element 0,1,size-1 st {M,L,...,U} w/L<=M<=U */
233827bd09bSSatish Balay       /* note ==> pivot_value in index 0                   */
234db4deed7SKarl Rupp       if (*pi > *pj) { SWAP(*pi, *pj) }
235*9371c9d4SSatish Balay       if (*ar > *pj) {
236*9371c9d4SSatish Balay         SWAP(*ar, *pj)
237*9371c9d4SSatish Balay       } else if (*pi > *ar) {
238*9371c9d4SSatish Balay         SWAP(*(ar), *(ar + 1))
239*9371c9d4SSatish Balay       }
240827bd09bSSatish Balay 
241827bd09bSSatish Balay       /* partition about pivot_value ...                              */
242827bd09bSSatish Balay       /* note lists of length 2 are not guaranteed to be sorted */
243db4deed7SKarl Rupp       for (;;) {
244827bd09bSSatish Balay         /* walk up ... and down ... swap if equal to pivot! */
245*9371c9d4SSatish Balay         do pi++;
246*9371c9d4SSatish Balay         while (*pi < *ar);
247*9371c9d4SSatish Balay         do pj--;
248*9371c9d4SSatish Balay         while (*pj > *ar);
249827bd09bSSatish Balay 
250827bd09bSSatish Balay         /* if we've crossed we're done */
251827bd09bSSatish Balay         if (pj < pi) break;
252827bd09bSSatish Balay 
253827bd09bSSatish Balay         /* else swap */
254827bd09bSSatish Balay         SWAP(*pi, *pj)
255827bd09bSSatish Balay       }
256827bd09bSSatish Balay 
257827bd09bSSatish Balay       /* place pivot_value in it's correct location */
258827bd09bSSatish Balay       SWAP(*ar, *pj)
259827bd09bSSatish Balay 
260827bd09bSSatish Balay       /* test stack_size to see if we've exhausted our stack */
26108401ef6SPierre Jolivet       PetscCheck(top_s - bottom_s < SORT_STACK, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCTFS_ivec_sort() :: STACK EXHAUSTED!!!");
262827bd09bSSatish Balay 
263827bd09bSSatish Balay       /* push right hand child iff length > 1 */
264db4deed7SKarl Rupp       if ((*top_s = size - ((PetscInt)(pi - ar)))) {
265827bd09bSSatish Balay         *(top_a++) = pi;
266827bd09bSSatish Balay         size -= *top_s + 2;
267827bd09bSSatish Balay         top_s++;
268*9371c9d4SSatish Balay       } else if (size -= *top_s + 2)
269*9371c9d4SSatish Balay         ; /* set up for next loop iff there is something to do */
270*9371c9d4SSatish Balay       else { /* might as well pop - note NR_OPT >=2 ==> we're ok! */ ar = *(--top_a);
271827bd09bSSatish Balay         size                                                            = *(--top_s);
272827bd09bSSatish Balay       }
273db4deed7SKarl Rupp     } else { /* else sort small list directly then pop another off stack */
274827bd09bSSatish Balay 
275827bd09bSSatish Balay       /* insertion sort for bottom */
276db4deed7SKarl Rupp       for (pj = ar + 1; pj <= ar + size; pj++) {
277827bd09bSSatish Balay         temp = *pj;
278db4deed7SKarl Rupp         for (pi = pj - 1; pi >= ar; pi--) {
279827bd09bSSatish Balay           if (*pi <= temp) break;
280827bd09bSSatish Balay           *(pi + 1) = *pi;
281827bd09bSSatish Balay         }
282827bd09bSSatish Balay         *(pi + 1) = temp;
283827bd09bSSatish Balay       }
284827bd09bSSatish Balay 
285827bd09bSSatish Balay       /* check to see if stack is exhausted ==> DONE */
2863fdc5746SBarry Smith       if (top_s == bottom_s) PetscFunctionReturn(0);
287827bd09bSSatish Balay 
288827bd09bSSatish Balay       /* else pop another list from the stack */
289827bd09bSSatish Balay       ar   = *(--top_a);
290827bd09bSSatish Balay       size = *(--top_s);
291827bd09bSSatish Balay     }
292827bd09bSSatish Balay   }
293827bd09bSSatish Balay }
294827bd09bSSatish Balay 
2957b1ae94cSBarry Smith /******************************************************************************/
296*9371c9d4SSatish Balay PetscErrorCode PCTFS_ivec_sort_companion(PetscInt *ar, PetscInt *ar2, PetscInt size) {
29752f87cdaSBarry Smith   PetscInt  *pi, *pj, temp, temp2;
29852f87cdaSBarry Smith   PetscInt **top_a = (PetscInt **)offset_stack;
29952f87cdaSBarry Smith   PetscInt  *top_s = size_stack, *bottom_s = size_stack;
30052f87cdaSBarry Smith   PetscInt  *pi2, *pj2;
30152f87cdaSBarry Smith   PetscInt   mid;
302827bd09bSSatish Balay 
3033fdc5746SBarry Smith   PetscFunctionBegin;
304827bd09bSSatish Balay   /* we're really interested in the offset of the last element */
305827bd09bSSatish Balay   /* ==> length of the list is now size + 1                    */
306827bd09bSSatish Balay   size--;
307827bd09bSSatish Balay 
308827bd09bSSatish Balay   /* do until we're done ... return when stack is exhausted */
309db4deed7SKarl Rupp   for (;;) {
310827bd09bSSatish Balay     /* if list is large enough use quicksort partition exchange code */
311db4deed7SKarl Rupp     if (size > SORT_OPT) {
312827bd09bSSatish Balay       /* start up pointer at element 1 and down at size     */
313827bd09bSSatish Balay       mid = size >> 1;
314827bd09bSSatish Balay       pi  = ar + 1;
315827bd09bSSatish Balay       pj  = ar + mid;
316827bd09bSSatish Balay       pi2 = ar2 + 1;
317827bd09bSSatish Balay       pj2 = ar2 + mid;
318827bd09bSSatish Balay 
319827bd09bSSatish Balay       /* find middle element in list and swap w/ element 1 */
320827bd09bSSatish Balay       SWAP(*pi, *pj)
321827bd09bSSatish Balay       SWAP(*pi2, *pj2)
322827bd09bSSatish Balay 
323827bd09bSSatish Balay       /* order element 0,1,size-1 st {M,L,...,U} w/L<=M<=U */
324827bd09bSSatish Balay       /* note ==> pivot_value in index 0                   */
325827bd09bSSatish Balay       pj  = ar + size;
326827bd09bSSatish Balay       pj2 = ar2 + size;
327db4deed7SKarl Rupp       if (*pi > *pj) { SWAP(*pi, *pj) SWAP(*pi2, *pj2) }
328*9371c9d4SSatish Balay       if (*ar > *pj) {
329*9371c9d4SSatish Balay         SWAP(*ar, *pj) SWAP(*ar2, *pj2)
330*9371c9d4SSatish Balay       } else if (*pi > *ar) {
331*9371c9d4SSatish Balay         SWAP(*(ar), *(ar + 1)) SWAP(*(ar2), *(ar2 + 1))
332*9371c9d4SSatish Balay       }
333827bd09bSSatish Balay 
334827bd09bSSatish Balay       /* partition about pivot_value ...                              */
335827bd09bSSatish Balay       /* note lists of length 2 are not guaranteed to be sorted */
336db4deed7SKarl Rupp       for (;;) {
337827bd09bSSatish Balay         /* walk up ... and down ... swap if equal to pivot! */
338*9371c9d4SSatish Balay         do {
339*9371c9d4SSatish Balay           pi++;
340*9371c9d4SSatish Balay           pi2++;
341*9371c9d4SSatish Balay         } while (*pi < *ar);
342*9371c9d4SSatish Balay         do {
343*9371c9d4SSatish Balay           pj--;
344*9371c9d4SSatish Balay           pj2--;
345*9371c9d4SSatish Balay         } while (*pj > *ar);
346827bd09bSSatish Balay 
347827bd09bSSatish Balay         /* if we've crossed we're done */
348827bd09bSSatish Balay         if (pj < pi) break;
349827bd09bSSatish Balay 
350827bd09bSSatish Balay         /* else swap */
351827bd09bSSatish Balay         SWAP(*pi, *pj)
352827bd09bSSatish Balay         SWAP(*pi2, *pj2)
353827bd09bSSatish Balay       }
354827bd09bSSatish Balay 
355827bd09bSSatish Balay       /* place pivot_value in it's correct location */
356827bd09bSSatish Balay       SWAP(*ar, *pj)
357827bd09bSSatish Balay       SWAP(*ar2, *pj2)
358827bd09bSSatish Balay 
359827bd09bSSatish Balay       /* test stack_size to see if we've exhausted our stack */
36008401ef6SPierre Jolivet       PetscCheck(top_s - bottom_s < SORT_STACK, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCTFS_ivec_sort_companion() :: STACK EXHAUSTED!!!");
361827bd09bSSatish Balay 
362827bd09bSSatish Balay       /* push right hand child iff length > 1 */
363db4deed7SKarl Rupp       if ((*top_s = size - ((PetscInt)(pi - ar)))) {
364827bd09bSSatish Balay         *(top_a++) = pi;
365827bd09bSSatish Balay         *(top_a++) = pi2;
366827bd09bSSatish Balay         size -= *top_s + 2;
367827bd09bSSatish Balay         top_s++;
368*9371c9d4SSatish Balay       } else if (size -= *top_s + 2)
369*9371c9d4SSatish Balay         ; /* set up for next loop iff there is something to do */
370*9371c9d4SSatish Balay       else { /* might as well pop - note NR_OPT >=2 ==> we're ok! */ ar2 = *(--top_a);
371827bd09bSSatish Balay         ar                                                               = *(--top_a);
372827bd09bSSatish Balay         size                                                             = *(--top_s);
373827bd09bSSatish Balay       }
374db4deed7SKarl Rupp     } else { /* else sort small list directly then pop another off stack */
375827bd09bSSatish Balay 
376827bd09bSSatish Balay       /* insertion sort for bottom */
377db4deed7SKarl Rupp       for (pj = ar + 1, pj2 = ar2 + 1; pj <= ar + size; pj++, pj2++) {
378827bd09bSSatish Balay         temp  = *pj;
379827bd09bSSatish Balay         temp2 = *pj2;
380db4deed7SKarl Rupp         for (pi = pj - 1, pi2 = pj2 - 1; pi >= ar; pi--, pi2--) {
381827bd09bSSatish Balay           if (*pi <= temp) break;
382827bd09bSSatish Balay           *(pi + 1)  = *pi;
383827bd09bSSatish Balay           *(pi2 + 1) = *pi2;
384827bd09bSSatish Balay         }
385827bd09bSSatish Balay         *(pi + 1)  = temp;
386827bd09bSSatish Balay         *(pi2 + 1) = temp2;
387827bd09bSSatish Balay       }
388827bd09bSSatish Balay 
389827bd09bSSatish Balay       /* check to see if stack is exhausted ==> DONE */
3903fdc5746SBarry Smith       if (top_s == bottom_s) PetscFunctionReturn(0);
391827bd09bSSatish Balay 
392827bd09bSSatish Balay       /* else pop another list from the stack */
393827bd09bSSatish Balay       ar2  = *(--top_a);
394827bd09bSSatish Balay       ar   = *(--top_a);
395827bd09bSSatish Balay       size = *(--top_s);
396827bd09bSSatish Balay     }
397827bd09bSSatish Balay   }
398827bd09bSSatish Balay }
399827bd09bSSatish Balay 
4007b1ae94cSBarry Smith /******************************************************************************/
401*9371c9d4SSatish Balay PetscErrorCode PCTFS_ivec_sort_companion_hack(PetscInt *ar, PetscInt **ar2, PetscInt size) {
40252f87cdaSBarry Smith   PetscInt  *pi, *pj, temp, *ptr;
40352f87cdaSBarry Smith   PetscInt **top_a = (PetscInt **)offset_stack;
40452f87cdaSBarry Smith   PetscInt  *top_s = size_stack, *bottom_s = size_stack;
40552f87cdaSBarry Smith   PetscInt **pi2, **pj2;
40652f87cdaSBarry Smith   PetscInt   mid;
407827bd09bSSatish Balay 
4083fdc5746SBarry Smith   PetscFunctionBegin;
409827bd09bSSatish Balay   /* we're really interested in the offset of the last element */
410827bd09bSSatish Balay   /* ==> length of the list is now size + 1                    */
411827bd09bSSatish Balay   size--;
412827bd09bSSatish Balay 
413827bd09bSSatish Balay   /* do until we're done ... return when stack is exhausted */
414db4deed7SKarl Rupp   for (;;) {
415827bd09bSSatish Balay     /* if list is large enough use quicksort partition exchange code */
416db4deed7SKarl Rupp     if (size > SORT_OPT) {
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) }
433*9371c9d4SSatish Balay       if (*ar > *pj) {
434*9371c9d4SSatish Balay         SWAP(*ar, *pj) P_SWAP(*ar2, *pj2)
435*9371c9d4SSatish Balay       } else if (*pi > *ar) {
436*9371c9d4SSatish Balay         SWAP(*(ar), *(ar + 1)) P_SWAP(*(ar2), *(ar2 + 1))
437*9371c9d4SSatish Balay       }
438827bd09bSSatish Balay 
439827bd09bSSatish Balay       /* partition about pivot_value ...                              */
440827bd09bSSatish Balay       /* note lists of length 2 are not guaranteed to be sorted */
441db4deed7SKarl Rupp       for (;;) {
442827bd09bSSatish Balay         /* walk up ... and down ... swap if equal to pivot! */
443*9371c9d4SSatish Balay         do {
444*9371c9d4SSatish Balay           pi++;
445*9371c9d4SSatish Balay           pi2++;
446*9371c9d4SSatish Balay         } while (*pi < *ar);
447*9371c9d4SSatish Balay         do {
448*9371c9d4SSatish Balay           pj--;
449*9371c9d4SSatish Balay           pj2--;
450*9371c9d4SSatish Balay         } while (*pj > *ar);
451827bd09bSSatish Balay 
452827bd09bSSatish Balay         /* if we've crossed we're done */
453827bd09bSSatish Balay         if (pj < pi) break;
454827bd09bSSatish Balay 
455827bd09bSSatish Balay         /* else swap */
456827bd09bSSatish Balay         SWAP(*pi, *pj)
457827bd09bSSatish Balay         P_SWAP(*pi2, *pj2)
458827bd09bSSatish Balay       }
459827bd09bSSatish Balay 
460827bd09bSSatish Balay       /* place pivot_value in it's correct location */
461827bd09bSSatish Balay       SWAP(*ar, *pj)
462827bd09bSSatish Balay       P_SWAP(*ar2, *pj2)
463827bd09bSSatish Balay 
464827bd09bSSatish Balay       /* test stack_size to see if we've exhausted our stack */
46508401ef6SPierre Jolivet       PetscCheck(top_s - bottom_s < SORT_STACK, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCTFS_ivec_sort_companion_hack() :: STACK EXHAUSTED!!!");
466827bd09bSSatish Balay 
467827bd09bSSatish Balay       /* push right hand child iff length > 1 */
468db4deed7SKarl Rupp       if ((*top_s = size - ((PetscInt)(pi - ar)))) {
469827bd09bSSatish Balay         *(top_a++) = pi;
47052f87cdaSBarry Smith         *(top_a++) = (PetscInt *)pi2;
471827bd09bSSatish Balay         size -= *top_s + 2;
472827bd09bSSatish Balay         top_s++;
473*9371c9d4SSatish Balay       } else if (size -= *top_s + 2)
474*9371c9d4SSatish Balay         ; /* set up for next loop iff there is something to do */
475*9371c9d4SSatish Balay       else { /* might as well pop - note NR_OPT >=2 ==> we're ok! */ ar2 = (PetscInt **)*(--top_a);
476827bd09bSSatish Balay         ar                                                               = *(--top_a);
477827bd09bSSatish Balay         size                                                             = *(--top_s);
478827bd09bSSatish Balay       }
4792fa5cd67SKarl Rupp     } else { /* else sort small list directly then pop another off stack */
480827bd09bSSatish Balay       /* insertion sort for bottom */
481db4deed7SKarl Rupp       for (pj = ar + 1, pj2 = ar2 + 1; pj <= ar + size; pj++, pj2++) {
482827bd09bSSatish Balay         temp = *pj;
483827bd09bSSatish Balay         ptr  = *pj2;
484db4deed7SKarl Rupp         for (pi = pj - 1, pi2 = pj2 - 1; pi >= ar; pi--, pi2--) {
485827bd09bSSatish Balay           if (*pi <= temp) break;
486827bd09bSSatish Balay           *(pi + 1)  = *pi;
487827bd09bSSatish Balay           *(pi2 + 1) = *pi2;
488827bd09bSSatish Balay         }
489827bd09bSSatish Balay         *(pi + 1)  = temp;
490827bd09bSSatish Balay         *(pi2 + 1) = ptr;
491827bd09bSSatish Balay       }
492827bd09bSSatish Balay 
493827bd09bSSatish Balay       /* check to see if stack is exhausted ==> DONE */
4943fdc5746SBarry Smith       if (top_s == bottom_s) PetscFunctionReturn(0);
495827bd09bSSatish Balay 
496827bd09bSSatish Balay       /* else pop another list from the stack */
49752f87cdaSBarry Smith       ar2  = (PetscInt **)*(--top_a);
498827bd09bSSatish Balay       ar   = *(--top_a);
499827bd09bSSatish Balay       size = *(--top_s);
500827bd09bSSatish Balay     }
501827bd09bSSatish Balay   }
502827bd09bSSatish Balay }
503827bd09bSSatish Balay 
5047b1ae94cSBarry Smith /******************************************************************************/
505*9371c9d4SSatish Balay PetscErrorCode PCTFS_SMI_sort(void *ar1, void *ar2, PetscInt size, PetscInt type) {
5063fdc5746SBarry Smith   PetscFunctionBegin;
507e7e72b3dSBarry Smith   if (type == SORT_INTEGER) {
508ca8e9878SJed Brown     if (ar2) PCTFS_ivec_sort_companion((PetscInt *)ar1, (PetscInt *)ar2, size);
509ca8e9878SJed Brown     else PCTFS_ivec_sort((PetscInt *)ar1, size);
510e7e72b3dSBarry Smith   } else if (type == SORT_INT_PTR) {
511ca8e9878SJed Brown     if (ar2) PCTFS_ivec_sort_companion_hack((PetscInt *)ar1, (PetscInt **)ar2, size);
512ca8e9878SJed Brown     else PCTFS_ivec_sort((PetscInt *)ar1, size);
513ca8e9878SJed Brown   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCTFS_SMI_sort only does SORT_INTEGER!");
5143fdc5746SBarry Smith   PetscFunctionReturn(0);
515827bd09bSSatish Balay }
516827bd09bSSatish Balay 
5177b1ae94cSBarry Smith /***********************************ivec.c*************************************/
518*9371c9d4SSatish Balay PetscInt PCTFS_ivec_linear_search(PetscInt item, PetscInt *list, PetscInt n) {
51952f87cdaSBarry Smith   PetscInt tmp = n - 1;
5205fd66863SKarl Rupp 
5212fa5cd67SKarl Rupp   while (n--) {
5222fa5cd67SKarl Rupp     if (*list++ == item) return (tmp - n);
5232fa5cd67SKarl Rupp   }
524827bd09bSSatish Balay   return (-1);
525827bd09bSSatish Balay }
526827bd09bSSatish Balay 
5277b1ae94cSBarry Smith /***********************************ivec.c*************************************/
528*9371c9d4SSatish Balay PetscInt PCTFS_ivec_binary_search(PetscInt item, PetscInt *list, PetscInt rh) {
52952f87cdaSBarry Smith   PetscInt mid, lh = 0;
530827bd09bSSatish Balay 
531827bd09bSSatish Balay   rh--;
532db4deed7SKarl Rupp   while (lh <= rh) {
533827bd09bSSatish Balay     mid = (lh + rh) >> 1;
5342fa5cd67SKarl Rupp     if (*(list + mid) == item) return (mid);
5352fa5cd67SKarl Rupp     if (*(list + mid) > item) rh = mid - 1;
5362fa5cd67SKarl Rupp     else lh = mid + 1;
537827bd09bSSatish Balay   }
538827bd09bSSatish Balay   return (-1);
539827bd09bSSatish Balay }
540827bd09bSSatish Balay 
5417b1ae94cSBarry Smith /*********************************ivec.c*************************************/
542*9371c9d4SSatish Balay PetscErrorCode PCTFS_rvec_copy(PetscScalar *arg1, PetscScalar *arg2, PetscInt n) {
5433fdc5746SBarry Smith   PetscFunctionBegin;
5442fa5cd67SKarl Rupp   while (n--) *arg1++ = *arg2++;
5453fdc5746SBarry Smith   PetscFunctionReturn(0);
546827bd09bSSatish Balay }
547827bd09bSSatish Balay 
5487b1ae94cSBarry Smith /*********************************ivec.c*************************************/
549*9371c9d4SSatish Balay PetscErrorCode PCTFS_rvec_zero(PetscScalar *arg1, PetscInt n) {
5503fdc5746SBarry Smith   PetscFunctionBegin;
5512fa5cd67SKarl Rupp   while (n--) *arg1++ = 0.0;
5523fdc5746SBarry Smith   PetscFunctionReturn(0);
553827bd09bSSatish Balay }
554827bd09bSSatish Balay 
5557b1ae94cSBarry Smith /***********************************ivec.c*************************************/
556*9371c9d4SSatish Balay PetscErrorCode PCTFS_rvec_one(PetscScalar *arg1, PetscInt n) {
5573fdc5746SBarry Smith   PetscFunctionBegin;
5582fa5cd67SKarl Rupp   while (n--) *arg1++ = 1.0;
5593fdc5746SBarry Smith   PetscFunctionReturn(0);
560827bd09bSSatish Balay }
561827bd09bSSatish Balay 
5627b1ae94cSBarry Smith /***********************************ivec.c*************************************/
563*9371c9d4SSatish Balay PetscErrorCode PCTFS_rvec_set(PetscScalar *arg1, PetscScalar arg2, PetscInt n) {
5643fdc5746SBarry Smith   PetscFunctionBegin;
5652fa5cd67SKarl Rupp   while (n--) *arg1++ = arg2;
5663fdc5746SBarry Smith   PetscFunctionReturn(0);
567827bd09bSSatish Balay }
568827bd09bSSatish Balay 
5697b1ae94cSBarry Smith /***********************************ivec.c*************************************/
570*9371c9d4SSatish Balay PetscErrorCode PCTFS_rvec_scale(PetscScalar *arg1, PetscScalar arg2, PetscInt n) {
5713fdc5746SBarry Smith   PetscFunctionBegin;
5722fa5cd67SKarl Rupp   while (n--) *arg1++ *= arg2;
5733fdc5746SBarry Smith   PetscFunctionReturn(0);
574827bd09bSSatish Balay }
575827bd09bSSatish Balay 
5767b1ae94cSBarry Smith /*********************************ivec.c*************************************/
577*9371c9d4SSatish Balay PetscErrorCode PCTFS_rvec_add(PetscScalar *arg1, PetscScalar *arg2, PetscInt n) {
5783fdc5746SBarry Smith   PetscFunctionBegin;
5792fa5cd67SKarl Rupp   while (n--) *arg1++ += *arg2++;
5803fdc5746SBarry Smith   PetscFunctionReturn(0);
581827bd09bSSatish Balay }
582827bd09bSSatish Balay 
5837b1ae94cSBarry Smith /*********************************ivec.c*************************************/
584*9371c9d4SSatish Balay PetscErrorCode PCTFS_rvec_mult(PetscScalar *arg1, PetscScalar *arg2, PetscInt n) {
5853fdc5746SBarry Smith   PetscFunctionBegin;
5862fa5cd67SKarl Rupp   while (n--) *arg1++ *= *arg2++;
5873fdc5746SBarry Smith   PetscFunctionReturn(0);
588827bd09bSSatish Balay }
589827bd09bSSatish Balay 
5907b1ae94cSBarry Smith /*********************************ivec.c*************************************/
591*9371c9d4SSatish Balay PetscErrorCode PCTFS_rvec_max(PetscScalar *arg1, PetscScalar *arg2, PetscInt n) {
5923fdc5746SBarry Smith   PetscFunctionBegin;
5932fa5cd67SKarl Rupp   while (n--) {
5942fa5cd67SKarl Rupp     *arg1 = PetscMax(*arg1, *arg2);
5952fa5cd67SKarl Rupp     arg1++;
5962fa5cd67SKarl Rupp     arg2++;
5972fa5cd67SKarl Rupp   }
5983fdc5746SBarry Smith   PetscFunctionReturn(0);
599827bd09bSSatish Balay }
600827bd09bSSatish Balay 
6017b1ae94cSBarry Smith /*********************************ivec.c*************************************/
602*9371c9d4SSatish Balay PetscErrorCode PCTFS_rvec_max_abs(PetscScalar *arg1, PetscScalar *arg2, PetscInt n) {
6033fdc5746SBarry Smith   PetscFunctionBegin;
6042fa5cd67SKarl Rupp   while (n--) {
6052fa5cd67SKarl Rupp     *arg1 = MAX_FABS(*arg1, *arg2);
6062fa5cd67SKarl Rupp     arg1++;
6072fa5cd67SKarl Rupp     arg2++;
6082fa5cd67SKarl Rupp   }
6093fdc5746SBarry Smith   PetscFunctionReturn(0);
610827bd09bSSatish Balay }
611827bd09bSSatish Balay 
6127b1ae94cSBarry Smith /*********************************ivec.c*************************************/
613*9371c9d4SSatish Balay PetscErrorCode PCTFS_rvec_min(PetscScalar *arg1, PetscScalar *arg2, PetscInt n) {
6143fdc5746SBarry Smith   PetscFunctionBegin;
6152fa5cd67SKarl Rupp   while (n--) {
6162fa5cd67SKarl Rupp     *arg1 = PetscMin(*arg1, *arg2);
6172fa5cd67SKarl Rupp     arg1++;
6182fa5cd67SKarl Rupp     arg2++;
6192fa5cd67SKarl Rupp   }
6203fdc5746SBarry Smith   PetscFunctionReturn(0);
621827bd09bSSatish Balay }
622827bd09bSSatish Balay 
6237b1ae94cSBarry Smith /*********************************ivec.c*************************************/
624*9371c9d4SSatish Balay PetscErrorCode PCTFS_rvec_min_abs(PetscScalar *arg1, PetscScalar *arg2, PetscInt n) {
6253fdc5746SBarry Smith   PetscFunctionBegin;
6262fa5cd67SKarl Rupp   while (n--) {
6272fa5cd67SKarl Rupp     *arg1 = MIN_FABS(*arg1, *arg2);
6282fa5cd67SKarl Rupp     arg1++;
6292fa5cd67SKarl Rupp     arg2++;
6302fa5cd67SKarl Rupp   }
6313fdc5746SBarry Smith   PetscFunctionReturn(0);
632827bd09bSSatish Balay }
633827bd09bSSatish Balay 
6347b1ae94cSBarry Smith /*********************************ivec.c*************************************/
635*9371c9d4SSatish Balay PetscErrorCode PCTFS_rvec_exists(PetscScalar *arg1, PetscScalar *arg2, PetscInt n) {
6363fdc5746SBarry Smith   PetscFunctionBegin;
6372fa5cd67SKarl Rupp   while (n--) {
6382fa5cd67SKarl Rupp     *arg1 = EXISTS(*arg1, *arg2);
6392fa5cd67SKarl Rupp     arg1++;
6402fa5cd67SKarl Rupp     arg2++;
6412fa5cd67SKarl Rupp   }
6423fdc5746SBarry Smith   PetscFunctionReturn(0);
643827bd09bSSatish Balay }
644827bd09bSSatish Balay 
6457b1ae94cSBarry Smith /***********************************ivec.c*************************************/
646*9371c9d4SSatish Balay PetscErrorCode PCTFS_rvec_non_uniform(PetscScalar *arg1, PetscScalar *arg2, PetscInt n, PetscInt *arg3) {
64752f87cdaSBarry Smith   PetscInt i, j, type;
648827bd09bSSatish Balay 
6493fdc5746SBarry Smith   PetscFunctionBegin;
650827bd09bSSatish Balay   /* LATER: if we're really motivated we can sort and then unsort */
651db4deed7SKarl Rupp   for (i = 0; i < n;) {
652827bd09bSSatish Balay     /* clump 'em for now */
653827bd09bSSatish Balay     j    = i + 1;
654827bd09bSSatish Balay     type = arg3[i];
6552fa5cd67SKarl Rupp     while ((j < n) && (arg3[j] == type)) j++;
656827bd09bSSatish Balay 
657827bd09bSSatish Balay     /* how many together */
658827bd09bSSatish Balay     j -= i;
659827bd09bSSatish Balay 
660827bd09bSSatish Balay     /* call appropriate ivec function */
6612fa5cd67SKarl Rupp     if (type == GL_MAX) PCTFS_rvec_max(arg1, arg2, j);
6622fa5cd67SKarl Rupp     else if (type == GL_MIN) PCTFS_rvec_min(arg1, arg2, j);
6632fa5cd67SKarl Rupp     else if (type == GL_MULT) PCTFS_rvec_mult(arg1, arg2, j);
6642fa5cd67SKarl Rupp     else if (type == GL_ADD) PCTFS_rvec_add(arg1, arg2, j);
6652fa5cd67SKarl Rupp     else if (type == GL_MAX_ABS) PCTFS_rvec_max_abs(arg1, arg2, j);
6662fa5cd67SKarl Rupp     else if (type == GL_MIN_ABS) PCTFS_rvec_min_abs(arg1, arg2, j);
6672fa5cd67SKarl Rupp     else if (type == GL_EXISTS) PCTFS_rvec_exists(arg1, arg2, j);
668ca8e9878SJed Brown     else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "unrecognized type passed to PCTFS_rvec_non_uniform()!");
669827bd09bSSatish Balay 
670*9371c9d4SSatish Balay     arg1 += j;
671*9371c9d4SSatish Balay     arg2 += j;
672*9371c9d4SSatish Balay     i += j;
673827bd09bSSatish Balay   }
6743fdc5746SBarry Smith   PetscFunctionReturn(0);
675827bd09bSSatish Balay }
676827bd09bSSatish Balay 
6777b1ae94cSBarry Smith /***********************************ivec.c*************************************/
678*9371c9d4SSatish Balay vfp PCTFS_rvec_fct_addr(PetscInt type) {
6792fa5cd67SKarl Rupp   if (type == NON_UNIFORM) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_rvec_non_uniform);
6802fa5cd67SKarl Rupp   else if (type == GL_MAX) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_rvec_max);
6812fa5cd67SKarl Rupp   else if (type == GL_MIN) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_rvec_min);
6822fa5cd67SKarl Rupp   else if (type == GL_MULT) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_rvec_mult);
6832fa5cd67SKarl Rupp   else if (type == GL_ADD) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_rvec_add);
6842fa5cd67SKarl Rupp   else if (type == GL_MAX_ABS) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_rvec_max_abs);
6852fa5cd67SKarl Rupp   else if (type == GL_MIN_ABS) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_rvec_min_abs);
6862fa5cd67SKarl Rupp   else if (type == GL_EXISTS) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_rvec_exists);
687827bd09bSSatish Balay 
688827bd09bSSatish Balay   /* catch all ... not good if we get here */
689827bd09bSSatish Balay   return (NULL);
690827bd09bSSatish Balay }
691