xref: /petsc/src/ksp/pc/impls/tfs/ivec.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
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*************************************/
28d71ae5a4SJacob Faibussowitsch PetscInt *PCTFS_ivec_copy(PetscInt *arg1, PetscInt *arg2, PetscInt n)
29d71ae5a4SJacob Faibussowitsch {
302fa5cd67SKarl Rupp   while (n--) *arg1++ = *arg2++;
31827bd09bSSatish Balay   return (arg1);
32827bd09bSSatish Balay }
33827bd09bSSatish Balay 
347b1ae94cSBarry Smith /***********************************ivec.c*************************************/
35d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_ivec_zero(PetscInt *arg1, PetscInt n)
36d71ae5a4SJacob Faibussowitsch {
373fdc5746SBarry Smith   PetscFunctionBegin;
382fa5cd67SKarl Rupp   while (n--) *arg1++ = 0;
39*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
40827bd09bSSatish Balay }
41827bd09bSSatish Balay 
427b1ae94cSBarry Smith /***********************************ivec.c*************************************/
43d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_ivec_set(PetscInt *arg1, PetscInt arg2, PetscInt n)
44d71ae5a4SJacob Faibussowitsch {
453fdc5746SBarry Smith   PetscFunctionBegin;
462fa5cd67SKarl Rupp   while (n--) *arg1++ = arg2;
47*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
48827bd09bSSatish Balay }
49827bd09bSSatish Balay 
507b1ae94cSBarry Smith /***********************************ivec.c*************************************/
51d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_ivec_max(PetscInt *arg1, PetscInt *arg2, PetscInt n)
52d71ae5a4SJacob Faibussowitsch {
533fdc5746SBarry Smith   PetscFunctionBegin;
549371c9d4SSatish Balay   while (n--) {
559371c9d4SSatish Balay     *arg1 = PetscMax(*arg1, *arg2);
569371c9d4SSatish Balay     arg1++;
579371c9d4SSatish Balay     arg2++;
589371c9d4SSatish Balay   }
59*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
60827bd09bSSatish Balay }
61827bd09bSSatish Balay 
627b1ae94cSBarry Smith /***********************************ivec.c*************************************/
63d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_ivec_min(PetscInt *arg1, PetscInt *arg2, PetscInt n)
64d71ae5a4SJacob Faibussowitsch {
653fdc5746SBarry Smith   PetscFunctionBegin;
662fa5cd67SKarl Rupp   while (n--) {
672fa5cd67SKarl Rupp     *(arg1) = PetscMin(*arg1, *arg2);
682fa5cd67SKarl Rupp     arg1++;
692fa5cd67SKarl Rupp     arg2++;
702fa5cd67SKarl Rupp   }
71*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
72827bd09bSSatish Balay }
73827bd09bSSatish Balay 
747b1ae94cSBarry Smith /***********************************ivec.c*************************************/
75d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_ivec_mult(PetscInt *arg1, PetscInt *arg2, PetscInt n)
76d71ae5a4SJacob Faibussowitsch {
773fdc5746SBarry Smith   PetscFunctionBegin;
782fa5cd67SKarl Rupp   while (n--) *arg1++ *= *arg2++;
79*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
80827bd09bSSatish Balay }
81827bd09bSSatish Balay 
827b1ae94cSBarry Smith /***********************************ivec.c*************************************/
83d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_ivec_add(PetscInt *arg1, PetscInt *arg2, PetscInt n)
84d71ae5a4SJacob Faibussowitsch {
853fdc5746SBarry Smith   PetscFunctionBegin;
862fa5cd67SKarl Rupp   while (n--) *arg1++ += *arg2++;
87*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
88827bd09bSSatish Balay }
89827bd09bSSatish Balay 
907b1ae94cSBarry Smith /***********************************ivec.c*************************************/
91d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_ivec_lxor(PetscInt *arg1, PetscInt *arg2, PetscInt n)
92d71ae5a4SJacob Faibussowitsch {
933fdc5746SBarry Smith   PetscFunctionBegin;
942fa5cd67SKarl Rupp   while (n--) {
952fa5cd67SKarl Rupp     *arg1 = ((*arg1 || *arg2) && !(*arg1 && *arg2));
962fa5cd67SKarl Rupp     arg1++;
972fa5cd67SKarl Rupp     arg2++;
982fa5cd67SKarl Rupp   }
99*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
100827bd09bSSatish Balay }
101827bd09bSSatish Balay 
1027b1ae94cSBarry Smith /***********************************ivec.c*************************************/
103d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_ivec_xor(PetscInt *arg1, PetscInt *arg2, PetscInt n)
104d71ae5a4SJacob Faibussowitsch {
1053fdc5746SBarry Smith   PetscFunctionBegin;
1062fa5cd67SKarl Rupp   while (n--) *arg1++ ^= *arg2++;
107*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
108827bd09bSSatish Balay }
109827bd09bSSatish Balay 
1107b1ae94cSBarry Smith /***********************************ivec.c*************************************/
111d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_ivec_or(PetscInt *arg1, PetscInt *arg2, PetscInt n)
112d71ae5a4SJacob Faibussowitsch {
1133fdc5746SBarry Smith   PetscFunctionBegin;
1142fa5cd67SKarl Rupp   while (n--) *arg1++ |= *arg2++;
115*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
116827bd09bSSatish Balay }
117827bd09bSSatish Balay 
1187b1ae94cSBarry Smith /***********************************ivec.c*************************************/
119d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_ivec_lor(PetscInt *arg1, PetscInt *arg2, PetscInt n)
120d71ae5a4SJacob Faibussowitsch {
1213fdc5746SBarry Smith   PetscFunctionBegin;
1222fa5cd67SKarl Rupp   while (n--) {
1232fa5cd67SKarl Rupp     *arg1 = (*arg1 || *arg2);
1242fa5cd67SKarl Rupp     arg1++;
1252fa5cd67SKarl Rupp     arg2++;
1262fa5cd67SKarl Rupp   }
127*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
128827bd09bSSatish Balay }
129827bd09bSSatish Balay 
1307b1ae94cSBarry Smith /***********************************ivec.c*************************************/
131d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_ivec_and(PetscInt *arg1, PetscInt *arg2, PetscInt n)
132d71ae5a4SJacob Faibussowitsch {
1333fdc5746SBarry Smith   PetscFunctionBegin;
1342fa5cd67SKarl Rupp   while (n--) *arg1++ &= *arg2++;
135*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
136827bd09bSSatish Balay }
137827bd09bSSatish Balay 
1387b1ae94cSBarry Smith /***********************************ivec.c*************************************/
139d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_ivec_land(PetscInt *arg1, PetscInt *arg2, PetscInt n)
140d71ae5a4SJacob Faibussowitsch {
1413fdc5746SBarry Smith   PetscFunctionBegin;
1422fa5cd67SKarl Rupp   while (n--) {
1432fa5cd67SKarl Rupp     *arg1 = (*arg1 && *arg2);
1442fa5cd67SKarl Rupp     arg1++;
1452fa5cd67SKarl Rupp     arg2++;
1462fa5cd67SKarl Rupp   }
147*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
148827bd09bSSatish Balay }
149827bd09bSSatish Balay 
1507b1ae94cSBarry Smith /***********************************ivec.c*************************************/
151d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_ivec_and3(PetscInt *arg1, PetscInt *arg2, PetscInt *arg3, PetscInt n)
152d71ae5a4SJacob Faibussowitsch {
1533fdc5746SBarry Smith   PetscFunctionBegin;
1542fa5cd67SKarl Rupp   while (n--) *arg1++ = (*arg2++ & *arg3++);
155*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
156827bd09bSSatish Balay }
157827bd09bSSatish Balay 
1587b1ae94cSBarry Smith /***********************************ivec.c*************************************/
159d71ae5a4SJacob Faibussowitsch PetscInt PCTFS_ivec_sum(PetscInt *arg1, PetscInt n)
160d71ae5a4SJacob Faibussowitsch {
16152f87cdaSBarry Smith   PetscInt tmp = 0;
1622fa5cd67SKarl Rupp   while (n--) tmp += *arg1++;
163827bd09bSSatish Balay   return (tmp);
164827bd09bSSatish Balay }
165827bd09bSSatish Balay 
1667b1ae94cSBarry Smith /***********************************ivec.c*************************************/
167d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_ivec_non_uniform(PetscInt *arg1, PetscInt *arg2, PetscInt n, ...)
168d71ae5a4SJacob Faibussowitsch {
16952f87cdaSBarry Smith   PetscInt  i, j, type;
1706ae68b59SBarry Smith   PetscInt *arg3;
1716ae68b59SBarry Smith   va_list   ap;
172827bd09bSSatish Balay 
1733fdc5746SBarry Smith   PetscFunctionBegin;
1746ae68b59SBarry Smith   va_start(ap, n);
1756ae68b59SBarry Smith   arg3 = va_arg(ap, PetscInt *);
1766ae68b59SBarry Smith   va_end(ap);
1776ae68b59SBarry Smith 
178827bd09bSSatish Balay   /* LATER: if we're really motivated we can sort and then unsort */
179db4deed7SKarl Rupp   for (i = 0; i < n;) {
180827bd09bSSatish Balay     /* clump 'em for now */
181827bd09bSSatish Balay     j    = i + 1;
182827bd09bSSatish Balay     type = arg3[i];
1832fa5cd67SKarl Rupp     while ((j < n) && (arg3[j] == type)) j++;
184827bd09bSSatish Balay 
185827bd09bSSatish Balay     /* how many together */
186827bd09bSSatish Balay     j -= i;
187827bd09bSSatish Balay 
188827bd09bSSatish Balay     /* call appropriate ivec function */
189*3ba16761SJacob Faibussowitsch     if (type == GL_MAX) PetscCall(PCTFS_ivec_max(arg1, arg2, j));
190*3ba16761SJacob Faibussowitsch     else if (type == GL_MIN) PetscCall(PCTFS_ivec_min(arg1, arg2, j));
191*3ba16761SJacob Faibussowitsch     else if (type == GL_MULT) PetscCall(PCTFS_ivec_mult(arg1, arg2, j));
192*3ba16761SJacob Faibussowitsch     else if (type == GL_ADD) PetscCall(PCTFS_ivec_add(arg1, arg2, j));
193*3ba16761SJacob Faibussowitsch     else if (type == GL_B_XOR) PetscCall(PCTFS_ivec_xor(arg1, arg2, j));
194*3ba16761SJacob Faibussowitsch     else if (type == GL_B_OR) PetscCall(PCTFS_ivec_or(arg1, arg2, j));
195*3ba16761SJacob Faibussowitsch     else if (type == GL_B_AND) PetscCall(PCTFS_ivec_and(arg1, arg2, j));
196*3ba16761SJacob Faibussowitsch     else if (type == GL_L_XOR) PetscCall(PCTFS_ivec_lxor(arg1, arg2, j));
197*3ba16761SJacob Faibussowitsch     else if (type == GL_L_OR) PetscCall(PCTFS_ivec_lor(arg1, arg2, j));
198*3ba16761SJacob Faibussowitsch     else if (type == GL_L_AND) PetscCall(PCTFS_ivec_land(arg1, arg2, j));
199db4deed7SKarl Rupp     else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "unrecognized type passed to PCTFS_ivec_non_uniform()!");
200827bd09bSSatish Balay 
2019371c9d4SSatish Balay     arg1 += j;
2029371c9d4SSatish Balay     arg2 += j;
2039371c9d4SSatish Balay     i += j;
204827bd09bSSatish Balay   }
205*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
206827bd09bSSatish Balay }
207827bd09bSSatish Balay 
2087b1ae94cSBarry Smith /***********************************ivec.c*************************************/
209d71ae5a4SJacob Faibussowitsch vfp PCTFS_ivec_fct_addr(PetscInt type)
210d71ae5a4SJacob Faibussowitsch {
2112fa5cd67SKarl Rupp   if (type == NON_UNIFORM) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_non_uniform);
2122fa5cd67SKarl Rupp   else if (type == GL_MAX) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_max);
2132fa5cd67SKarl Rupp   else if (type == GL_MIN) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_min);
2142fa5cd67SKarl Rupp   else if (type == GL_MULT) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_mult);
2152fa5cd67SKarl Rupp   else if (type == GL_ADD) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_add);
2162fa5cd67SKarl Rupp   else if (type == GL_B_XOR) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_xor);
2172fa5cd67SKarl Rupp   else if (type == GL_B_OR) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_or);
2182fa5cd67SKarl Rupp   else if (type == GL_B_AND) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_and);
2192fa5cd67SKarl Rupp   else if (type == GL_L_XOR) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_lxor);
2202fa5cd67SKarl Rupp   else if (type == GL_L_OR) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_lor);
2212fa5cd67SKarl Rupp   else if (type == GL_L_AND) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_land);
222827bd09bSSatish Balay 
223827bd09bSSatish Balay   /* catch all ... not good if we get here */
224827bd09bSSatish Balay   return (NULL);
225827bd09bSSatish Balay }
226827bd09bSSatish Balay 
2277b1ae94cSBarry Smith /******************************************************************************/
228d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_ivec_sort(PetscInt *ar, PetscInt size)
229d71ae5a4SJacob Faibussowitsch {
23052f87cdaSBarry Smith   PetscInt  *pi, *pj, temp;
23152f87cdaSBarry Smith   PetscInt **top_a = (PetscInt **)offset_stack;
23252f87cdaSBarry Smith   PetscInt  *top_s = size_stack, *bottom_s = size_stack;
233827bd09bSSatish Balay 
234b458e8f1SJose E. Roman   PetscFunctionBegin;
235827bd09bSSatish Balay   /* we're really interested in the offset of the last element */
236827bd09bSSatish Balay   /* ==> length of the list is now size + 1                    */
237827bd09bSSatish Balay   size--;
238827bd09bSSatish Balay 
239827bd09bSSatish Balay   /* do until we're done ... return when stack is exhausted */
240db4deed7SKarl Rupp   for (;;) {
241827bd09bSSatish Balay     /* if list is large enough use quicksort partition exchange code */
242db4deed7SKarl Rupp     if (size > SORT_OPT) {
243827bd09bSSatish Balay       /* start up pointer at element 1 and down at size     */
244827bd09bSSatish Balay       pi = ar + 1;
245827bd09bSSatish Balay       pj = ar + size;
246827bd09bSSatish Balay 
247827bd09bSSatish Balay       /* find middle element in list and swap w/ element 1 */
248827bd09bSSatish Balay       SWAP(*(ar + (size >> 1)), *pi)
249827bd09bSSatish Balay 
250827bd09bSSatish Balay       /* order element 0,1,size-1 st {M,L,...,U} w/L<=M<=U */
251827bd09bSSatish Balay       /* note ==> pivot_value in index 0                   */
252db4deed7SKarl Rupp       if (*pi > *pj) { SWAP(*pi, *pj) }
2539371c9d4SSatish Balay       if (*ar > *pj) {
2549371c9d4SSatish Balay         SWAP(*ar, *pj)
2559371c9d4SSatish Balay       } else if (*pi > *ar) {
2569371c9d4SSatish Balay         SWAP(*(ar), *(ar + 1))
2579371c9d4SSatish Balay       }
258827bd09bSSatish Balay 
259827bd09bSSatish Balay       /* partition about pivot_value ...                              */
260827bd09bSSatish Balay       /* note lists of length 2 are not guaranteed to be sorted */
261db4deed7SKarl Rupp       for (;;) {
262827bd09bSSatish Balay         /* walk up ... and down ... swap if equal to pivot! */
2639371c9d4SSatish Balay         do pi++;
2649371c9d4SSatish Balay         while (*pi < *ar);
2659371c9d4SSatish Balay         do pj--;
2669371c9d4SSatish Balay         while (*pj > *ar);
267827bd09bSSatish Balay 
268827bd09bSSatish Balay         /* if we've crossed we're done */
269827bd09bSSatish Balay         if (pj < pi) break;
270827bd09bSSatish Balay 
271827bd09bSSatish Balay         /* else swap */
272827bd09bSSatish Balay         SWAP(*pi, *pj)
273827bd09bSSatish Balay       }
274827bd09bSSatish Balay 
275827bd09bSSatish Balay       /* place pivot_value in it's correct location */
276827bd09bSSatish Balay       SWAP(*ar, *pj)
277827bd09bSSatish Balay 
278827bd09bSSatish Balay       /* test stack_size to see if we've exhausted our stack */
27908401ef6SPierre Jolivet       PetscCheck(top_s - bottom_s < SORT_STACK, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCTFS_ivec_sort() :: STACK EXHAUSTED!!!");
280827bd09bSSatish Balay 
281827bd09bSSatish Balay       /* push right hand child iff length > 1 */
282db4deed7SKarl Rupp       if ((*top_s = size - ((PetscInt)(pi - ar)))) {
283827bd09bSSatish Balay         *(top_a++) = pi;
284827bd09bSSatish Balay         size -= *top_s + 2;
285827bd09bSSatish Balay         top_s++;
2869371c9d4SSatish Balay       } else if (size -= *top_s + 2)
2879371c9d4SSatish Balay         ; /* set up for next loop iff there is something to do */
2889371c9d4SSatish Balay       else { /* might as well pop - note NR_OPT >=2 ==> we're ok! */ ar = *(--top_a);
289827bd09bSSatish Balay         size                                                            = *(--top_s);
290827bd09bSSatish Balay       }
291db4deed7SKarl Rupp     } else { /* else sort small list directly then pop another off stack */
292827bd09bSSatish Balay 
293827bd09bSSatish Balay       /* insertion sort for bottom */
294db4deed7SKarl Rupp       for (pj = ar + 1; pj <= ar + size; pj++) {
295827bd09bSSatish Balay         temp = *pj;
296db4deed7SKarl Rupp         for (pi = pj - 1; pi >= ar; pi--) {
297827bd09bSSatish Balay           if (*pi <= temp) break;
298827bd09bSSatish Balay           *(pi + 1) = *pi;
299827bd09bSSatish Balay         }
300827bd09bSSatish Balay         *(pi + 1) = temp;
301827bd09bSSatish Balay       }
302827bd09bSSatish Balay 
303827bd09bSSatish Balay       /* check to see if stack is exhausted ==> DONE */
304*3ba16761SJacob Faibussowitsch       if (top_s == bottom_s) PetscFunctionReturn(PETSC_SUCCESS);
305827bd09bSSatish Balay 
306827bd09bSSatish Balay       /* else pop another list from the stack */
307827bd09bSSatish Balay       ar   = *(--top_a);
308827bd09bSSatish Balay       size = *(--top_s);
309827bd09bSSatish Balay     }
310827bd09bSSatish Balay   }
311827bd09bSSatish Balay }
312827bd09bSSatish Balay 
3137b1ae94cSBarry Smith /******************************************************************************/
314d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_ivec_sort_companion(PetscInt *ar, PetscInt *ar2, PetscInt size)
315d71ae5a4SJacob Faibussowitsch {
31652f87cdaSBarry Smith   PetscInt  *pi, *pj, temp, temp2;
31752f87cdaSBarry Smith   PetscInt **top_a = (PetscInt **)offset_stack;
31852f87cdaSBarry Smith   PetscInt  *top_s = size_stack, *bottom_s = size_stack;
31952f87cdaSBarry Smith   PetscInt  *pi2, *pj2;
32052f87cdaSBarry Smith   PetscInt   mid;
321827bd09bSSatish Balay 
3223fdc5746SBarry Smith   PetscFunctionBegin;
323827bd09bSSatish Balay   /* we're really interested in the offset of the last element */
324827bd09bSSatish Balay   /* ==> length of the list is now size + 1                    */
325827bd09bSSatish Balay   size--;
326827bd09bSSatish Balay 
327827bd09bSSatish Balay   /* do until we're done ... return when stack is exhausted */
328db4deed7SKarl Rupp   for (;;) {
329827bd09bSSatish Balay     /* if list is large enough use quicksort partition exchange code */
330db4deed7SKarl Rupp     if (size > SORT_OPT) {
331827bd09bSSatish Balay       /* start up pointer at element 1 and down at size     */
332827bd09bSSatish Balay       mid = size >> 1;
333827bd09bSSatish Balay       pi  = ar + 1;
334827bd09bSSatish Balay       pj  = ar + mid;
335827bd09bSSatish Balay       pi2 = ar2 + 1;
336827bd09bSSatish Balay       pj2 = ar2 + mid;
337827bd09bSSatish Balay 
338827bd09bSSatish Balay       /* find middle element in list and swap w/ element 1 */
339827bd09bSSatish Balay       SWAP(*pi, *pj)
340827bd09bSSatish Balay       SWAP(*pi2, *pj2)
341827bd09bSSatish Balay 
342827bd09bSSatish Balay       /* order element 0,1,size-1 st {M,L,...,U} w/L<=M<=U */
343827bd09bSSatish Balay       /* note ==> pivot_value in index 0                   */
344827bd09bSSatish Balay       pj  = ar + size;
345827bd09bSSatish Balay       pj2 = ar2 + size;
346db4deed7SKarl Rupp       if (*pi > *pj) { SWAP(*pi, *pj) SWAP(*pi2, *pj2) }
3479371c9d4SSatish Balay       if (*ar > *pj) {
3489371c9d4SSatish Balay         SWAP(*ar, *pj) SWAP(*ar2, *pj2)
3499371c9d4SSatish Balay       } else if (*pi > *ar) {
3509371c9d4SSatish Balay         SWAP(*(ar), *(ar + 1)) SWAP(*(ar2), *(ar2 + 1))
3519371c9d4SSatish Balay       }
352827bd09bSSatish Balay 
353827bd09bSSatish Balay       /* partition about pivot_value ...                              */
354827bd09bSSatish Balay       /* note lists of length 2 are not guaranteed to be sorted */
355db4deed7SKarl Rupp       for (;;) {
356827bd09bSSatish Balay         /* walk up ... and down ... swap if equal to pivot! */
3579371c9d4SSatish Balay         do {
3589371c9d4SSatish Balay           pi++;
3599371c9d4SSatish Balay           pi2++;
3609371c9d4SSatish Balay         } while (*pi < *ar);
3619371c9d4SSatish Balay         do {
3629371c9d4SSatish Balay           pj--;
3639371c9d4SSatish Balay           pj2--;
3649371c9d4SSatish Balay         } while (*pj > *ar);
365827bd09bSSatish Balay 
366827bd09bSSatish Balay         /* if we've crossed we're done */
367827bd09bSSatish Balay         if (pj < pi) break;
368827bd09bSSatish Balay 
369827bd09bSSatish Balay         /* else swap */
370827bd09bSSatish Balay         SWAP(*pi, *pj)
371827bd09bSSatish Balay         SWAP(*pi2, *pj2)
372827bd09bSSatish Balay       }
373827bd09bSSatish Balay 
374827bd09bSSatish Balay       /* place pivot_value in it's correct location */
375827bd09bSSatish Balay       SWAP(*ar, *pj)
376827bd09bSSatish Balay       SWAP(*ar2, *pj2)
377827bd09bSSatish Balay 
378827bd09bSSatish Balay       /* test stack_size to see if we've exhausted our stack */
37908401ef6SPierre Jolivet       PetscCheck(top_s - bottom_s < SORT_STACK, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCTFS_ivec_sort_companion() :: STACK EXHAUSTED!!!");
380827bd09bSSatish Balay 
381827bd09bSSatish Balay       /* push right hand child iff length > 1 */
382db4deed7SKarl Rupp       if ((*top_s = size - ((PetscInt)(pi - ar)))) {
383827bd09bSSatish Balay         *(top_a++) = pi;
384827bd09bSSatish Balay         *(top_a++) = pi2;
385827bd09bSSatish Balay         size -= *top_s + 2;
386827bd09bSSatish Balay         top_s++;
3879371c9d4SSatish Balay       } else if (size -= *top_s + 2)
3889371c9d4SSatish Balay         ; /* set up for next loop iff there is something to do */
3899371c9d4SSatish Balay       else { /* might as well pop - note NR_OPT >=2 ==> we're ok! */ ar2 = *(--top_a);
390827bd09bSSatish Balay         ar                                                               = *(--top_a);
391827bd09bSSatish Balay         size                                                             = *(--top_s);
392827bd09bSSatish Balay       }
393db4deed7SKarl Rupp     } else { /* else sort small list directly then pop another off stack */
394827bd09bSSatish Balay 
395827bd09bSSatish Balay       /* insertion sort for bottom */
396db4deed7SKarl Rupp       for (pj = ar + 1, pj2 = ar2 + 1; pj <= ar + size; pj++, pj2++) {
397827bd09bSSatish Balay         temp  = *pj;
398827bd09bSSatish Balay         temp2 = *pj2;
399db4deed7SKarl Rupp         for (pi = pj - 1, pi2 = pj2 - 1; pi >= ar; pi--, pi2--) {
400827bd09bSSatish Balay           if (*pi <= temp) break;
401827bd09bSSatish Balay           *(pi + 1)  = *pi;
402827bd09bSSatish Balay           *(pi2 + 1) = *pi2;
403827bd09bSSatish Balay         }
404827bd09bSSatish Balay         *(pi + 1)  = temp;
405827bd09bSSatish Balay         *(pi2 + 1) = temp2;
406827bd09bSSatish Balay       }
407827bd09bSSatish Balay 
408827bd09bSSatish Balay       /* check to see if stack is exhausted ==> DONE */
409*3ba16761SJacob Faibussowitsch       if (top_s == bottom_s) PetscFunctionReturn(PETSC_SUCCESS);
410827bd09bSSatish Balay 
411827bd09bSSatish Balay       /* else pop another list from the stack */
412827bd09bSSatish Balay       ar2  = *(--top_a);
413827bd09bSSatish Balay       ar   = *(--top_a);
414827bd09bSSatish Balay       size = *(--top_s);
415827bd09bSSatish Balay     }
416827bd09bSSatish Balay   }
417827bd09bSSatish Balay }
418827bd09bSSatish Balay 
4197b1ae94cSBarry Smith /******************************************************************************/
420d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_ivec_sort_companion_hack(PetscInt *ar, PetscInt **ar2, PetscInt size)
421d71ae5a4SJacob Faibussowitsch {
42252f87cdaSBarry Smith   PetscInt  *pi, *pj, temp, *ptr;
42352f87cdaSBarry Smith   PetscInt **top_a = (PetscInt **)offset_stack;
42452f87cdaSBarry Smith   PetscInt  *top_s = size_stack, *bottom_s = size_stack;
42552f87cdaSBarry Smith   PetscInt **pi2, **pj2;
42652f87cdaSBarry Smith   PetscInt   mid;
427827bd09bSSatish Balay 
4283fdc5746SBarry Smith   PetscFunctionBegin;
429827bd09bSSatish Balay   /* we're really interested in the offset of the last element */
430827bd09bSSatish Balay   /* ==> length of the list is now size + 1                    */
431827bd09bSSatish Balay   size--;
432827bd09bSSatish Balay 
433827bd09bSSatish Balay   /* do until we're done ... return when stack is exhausted */
434db4deed7SKarl Rupp   for (;;) {
435827bd09bSSatish Balay     /* if list is large enough use quicksort partition exchange code */
436db4deed7SKarl Rupp     if (size > SORT_OPT) {
437827bd09bSSatish Balay       /* start up pointer at element 1 and down at size     */
438827bd09bSSatish Balay       mid = size >> 1;
439827bd09bSSatish Balay       pi  = ar + 1;
440827bd09bSSatish Balay       pj  = ar + mid;
441827bd09bSSatish Balay       pi2 = ar2 + 1;
442827bd09bSSatish Balay       pj2 = ar2 + mid;
443827bd09bSSatish Balay 
444827bd09bSSatish Balay       /* find middle element in list and swap w/ element 1 */
445827bd09bSSatish Balay       SWAP(*pi, *pj)
446827bd09bSSatish Balay       P_SWAP(*pi2, *pj2)
447827bd09bSSatish Balay 
448827bd09bSSatish Balay       /* order element 0,1,size-1 st {M,L,...,U} w/L<=M<=U */
449827bd09bSSatish Balay       /* note ==> pivot_value in index 0                   */
450827bd09bSSatish Balay       pj  = ar + size;
451827bd09bSSatish Balay       pj2 = ar2 + size;
452db4deed7SKarl Rupp       if (*pi > *pj) { SWAP(*pi, *pj) P_SWAP(*pi2, *pj2) }
4539371c9d4SSatish Balay       if (*ar > *pj) {
4549371c9d4SSatish Balay         SWAP(*ar, *pj) P_SWAP(*ar2, *pj2)
4559371c9d4SSatish Balay       } else if (*pi > *ar) {
4569371c9d4SSatish Balay         SWAP(*(ar), *(ar + 1)) P_SWAP(*(ar2), *(ar2 + 1))
4579371c9d4SSatish Balay       }
458827bd09bSSatish Balay 
459827bd09bSSatish Balay       /* partition about pivot_value ...                              */
460827bd09bSSatish Balay       /* note lists of length 2 are not guaranteed to be sorted */
461db4deed7SKarl Rupp       for (;;) {
462827bd09bSSatish Balay         /* walk up ... and down ... swap if equal to pivot! */
4639371c9d4SSatish Balay         do {
4649371c9d4SSatish Balay           pi++;
4659371c9d4SSatish Balay           pi2++;
4669371c9d4SSatish Balay         } while (*pi < *ar);
4679371c9d4SSatish Balay         do {
4689371c9d4SSatish Balay           pj--;
4699371c9d4SSatish Balay           pj2--;
4709371c9d4SSatish Balay         } while (*pj > *ar);
471827bd09bSSatish Balay 
472827bd09bSSatish Balay         /* if we've crossed we're done */
473827bd09bSSatish Balay         if (pj < pi) break;
474827bd09bSSatish Balay 
475827bd09bSSatish Balay         /* else swap */
476827bd09bSSatish Balay         SWAP(*pi, *pj)
477827bd09bSSatish Balay         P_SWAP(*pi2, *pj2)
478827bd09bSSatish Balay       }
479827bd09bSSatish Balay 
480827bd09bSSatish Balay       /* place pivot_value in it's correct location */
481827bd09bSSatish Balay       SWAP(*ar, *pj)
482827bd09bSSatish Balay       P_SWAP(*ar2, *pj2)
483827bd09bSSatish Balay 
484827bd09bSSatish Balay       /* test stack_size to see if we've exhausted our stack */
48508401ef6SPierre Jolivet       PetscCheck(top_s - bottom_s < SORT_STACK, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCTFS_ivec_sort_companion_hack() :: STACK EXHAUSTED!!!");
486827bd09bSSatish Balay 
487827bd09bSSatish Balay       /* push right hand child iff length > 1 */
488db4deed7SKarl Rupp       if ((*top_s = size - ((PetscInt)(pi - ar)))) {
489827bd09bSSatish Balay         *(top_a++) = pi;
49052f87cdaSBarry Smith         *(top_a++) = (PetscInt *)pi2;
491827bd09bSSatish Balay         size -= *top_s + 2;
492827bd09bSSatish Balay         top_s++;
4939371c9d4SSatish Balay       } else if (size -= *top_s + 2)
4949371c9d4SSatish Balay         ; /* set up for next loop iff there is something to do */
4959371c9d4SSatish Balay       else { /* might as well pop - note NR_OPT >=2 ==> we're ok! */ ar2 = (PetscInt **)*(--top_a);
496827bd09bSSatish Balay         ar                                                               = *(--top_a);
497827bd09bSSatish Balay         size                                                             = *(--top_s);
498827bd09bSSatish Balay       }
4992fa5cd67SKarl Rupp     } else { /* else sort small list directly then pop another off stack */
500827bd09bSSatish Balay       /* insertion sort for bottom */
501db4deed7SKarl Rupp       for (pj = ar + 1, pj2 = ar2 + 1; pj <= ar + size; pj++, pj2++) {
502827bd09bSSatish Balay         temp = *pj;
503827bd09bSSatish Balay         ptr  = *pj2;
504db4deed7SKarl Rupp         for (pi = pj - 1, pi2 = pj2 - 1; pi >= ar; pi--, pi2--) {
505827bd09bSSatish Balay           if (*pi <= temp) break;
506827bd09bSSatish Balay           *(pi + 1)  = *pi;
507827bd09bSSatish Balay           *(pi2 + 1) = *pi2;
508827bd09bSSatish Balay         }
509827bd09bSSatish Balay         *(pi + 1)  = temp;
510827bd09bSSatish Balay         *(pi2 + 1) = ptr;
511827bd09bSSatish Balay       }
512827bd09bSSatish Balay 
513827bd09bSSatish Balay       /* check to see if stack is exhausted ==> DONE */
514*3ba16761SJacob Faibussowitsch       if (top_s == bottom_s) PetscFunctionReturn(PETSC_SUCCESS);
515827bd09bSSatish Balay 
516827bd09bSSatish Balay       /* else pop another list from the stack */
51752f87cdaSBarry Smith       ar2  = (PetscInt **)*(--top_a);
518827bd09bSSatish Balay       ar   = *(--top_a);
519827bd09bSSatish Balay       size = *(--top_s);
520827bd09bSSatish Balay     }
521827bd09bSSatish Balay   }
522827bd09bSSatish Balay }
523827bd09bSSatish Balay 
5247b1ae94cSBarry Smith /******************************************************************************/
525d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_SMI_sort(void *ar1, void *ar2, PetscInt size, PetscInt type)
526d71ae5a4SJacob Faibussowitsch {
5273fdc5746SBarry Smith   PetscFunctionBegin;
528e7e72b3dSBarry Smith   if (type == SORT_INTEGER) {
529*3ba16761SJacob Faibussowitsch     if (ar2) PetscCall(PCTFS_ivec_sort_companion((PetscInt *)ar1, (PetscInt *)ar2, size));
530*3ba16761SJacob Faibussowitsch     else PetscCall(PCTFS_ivec_sort((PetscInt *)ar1, size));
531e7e72b3dSBarry Smith   } else if (type == SORT_INT_PTR) {
532*3ba16761SJacob Faibussowitsch     if (ar2) PetscCall(PCTFS_ivec_sort_companion_hack((PetscInt *)ar1, (PetscInt **)ar2, size));
533*3ba16761SJacob Faibussowitsch     else PetscCall(PCTFS_ivec_sort((PetscInt *)ar1, size));
534ca8e9878SJed Brown   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCTFS_SMI_sort only does SORT_INTEGER!");
535*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
536827bd09bSSatish Balay }
537827bd09bSSatish Balay 
5387b1ae94cSBarry Smith /***********************************ivec.c*************************************/
539d71ae5a4SJacob Faibussowitsch PetscInt PCTFS_ivec_linear_search(PetscInt item, PetscInt *list, PetscInt n)
540d71ae5a4SJacob Faibussowitsch {
54152f87cdaSBarry Smith   PetscInt tmp = n - 1;
5425fd66863SKarl Rupp 
5432fa5cd67SKarl Rupp   while (n--) {
5442fa5cd67SKarl Rupp     if (*list++ == item) return (tmp - n);
5452fa5cd67SKarl Rupp   }
546827bd09bSSatish Balay   return (-1);
547827bd09bSSatish Balay }
548827bd09bSSatish Balay 
5497b1ae94cSBarry Smith /***********************************ivec.c*************************************/
550d71ae5a4SJacob Faibussowitsch PetscInt PCTFS_ivec_binary_search(PetscInt item, PetscInt *list, PetscInt rh)
551d71ae5a4SJacob Faibussowitsch {
55252f87cdaSBarry Smith   PetscInt mid, lh = 0;
553827bd09bSSatish Balay 
554827bd09bSSatish Balay   rh--;
555db4deed7SKarl Rupp   while (lh <= rh) {
556827bd09bSSatish Balay     mid = (lh + rh) >> 1;
5572fa5cd67SKarl Rupp     if (*(list + mid) == item) return (mid);
5582fa5cd67SKarl Rupp     if (*(list + mid) > item) rh = mid - 1;
5592fa5cd67SKarl Rupp     else lh = mid + 1;
560827bd09bSSatish Balay   }
561827bd09bSSatish Balay   return (-1);
562827bd09bSSatish Balay }
563827bd09bSSatish Balay 
5647b1ae94cSBarry Smith /*********************************ivec.c*************************************/
565d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_rvec_copy(PetscScalar *arg1, PetscScalar *arg2, PetscInt n)
566d71ae5a4SJacob Faibussowitsch {
5673fdc5746SBarry Smith   PetscFunctionBegin;
5682fa5cd67SKarl Rupp   while (n--) *arg1++ = *arg2++;
569*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
570827bd09bSSatish Balay }
571827bd09bSSatish Balay 
5727b1ae94cSBarry Smith /*********************************ivec.c*************************************/
573d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_rvec_zero(PetscScalar *arg1, PetscInt n)
574d71ae5a4SJacob Faibussowitsch {
5753fdc5746SBarry Smith   PetscFunctionBegin;
5762fa5cd67SKarl Rupp   while (n--) *arg1++ = 0.0;
577*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
578827bd09bSSatish Balay }
579827bd09bSSatish Balay 
5807b1ae94cSBarry Smith /***********************************ivec.c*************************************/
581d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_rvec_one(PetscScalar *arg1, PetscInt n)
582d71ae5a4SJacob Faibussowitsch {
5833fdc5746SBarry Smith   PetscFunctionBegin;
5842fa5cd67SKarl Rupp   while (n--) *arg1++ = 1.0;
585*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
586827bd09bSSatish Balay }
587827bd09bSSatish Balay 
5887b1ae94cSBarry Smith /***********************************ivec.c*************************************/
589d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_rvec_set(PetscScalar *arg1, PetscScalar arg2, PetscInt n)
590d71ae5a4SJacob Faibussowitsch {
5913fdc5746SBarry Smith   PetscFunctionBegin;
5922fa5cd67SKarl Rupp   while (n--) *arg1++ = arg2;
593*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
594827bd09bSSatish Balay }
595827bd09bSSatish Balay 
5967b1ae94cSBarry Smith /***********************************ivec.c*************************************/
597d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_rvec_scale(PetscScalar *arg1, PetscScalar arg2, PetscInt n)
598d71ae5a4SJacob Faibussowitsch {
5993fdc5746SBarry Smith   PetscFunctionBegin;
6002fa5cd67SKarl Rupp   while (n--) *arg1++ *= arg2;
601*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
602827bd09bSSatish Balay }
603827bd09bSSatish Balay 
6047b1ae94cSBarry Smith /*********************************ivec.c*************************************/
605d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_rvec_add(PetscScalar *arg1, PetscScalar *arg2, PetscInt n)
606d71ae5a4SJacob Faibussowitsch {
6073fdc5746SBarry Smith   PetscFunctionBegin;
6082fa5cd67SKarl Rupp   while (n--) *arg1++ += *arg2++;
609*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
610827bd09bSSatish Balay }
611827bd09bSSatish Balay 
6127b1ae94cSBarry Smith /*********************************ivec.c*************************************/
613d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_rvec_mult(PetscScalar *arg1, PetscScalar *arg2, PetscInt n)
614d71ae5a4SJacob Faibussowitsch {
6153fdc5746SBarry Smith   PetscFunctionBegin;
6162fa5cd67SKarl Rupp   while (n--) *arg1++ *= *arg2++;
617*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
618827bd09bSSatish Balay }
619827bd09bSSatish Balay 
6207b1ae94cSBarry Smith /*********************************ivec.c*************************************/
621d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_rvec_max(PetscScalar *arg1, PetscScalar *arg2, PetscInt n)
622d71ae5a4SJacob Faibussowitsch {
6233fdc5746SBarry Smith   PetscFunctionBegin;
6242fa5cd67SKarl Rupp   while (n--) {
6252fa5cd67SKarl Rupp     *arg1 = PetscMax(*arg1, *arg2);
6262fa5cd67SKarl Rupp     arg1++;
6272fa5cd67SKarl Rupp     arg2++;
6282fa5cd67SKarl Rupp   }
629*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
630827bd09bSSatish Balay }
631827bd09bSSatish Balay 
6327b1ae94cSBarry Smith /*********************************ivec.c*************************************/
633d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_rvec_max_abs(PetscScalar *arg1, PetscScalar *arg2, PetscInt n)
634d71ae5a4SJacob Faibussowitsch {
6353fdc5746SBarry Smith   PetscFunctionBegin;
6362fa5cd67SKarl Rupp   while (n--) {
6372fa5cd67SKarl Rupp     *arg1 = MAX_FABS(*arg1, *arg2);
6382fa5cd67SKarl Rupp     arg1++;
6392fa5cd67SKarl Rupp     arg2++;
6402fa5cd67SKarl Rupp   }
641*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
642827bd09bSSatish Balay }
643827bd09bSSatish Balay 
6447b1ae94cSBarry Smith /*********************************ivec.c*************************************/
645d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_rvec_min(PetscScalar *arg1, PetscScalar *arg2, PetscInt n)
646d71ae5a4SJacob Faibussowitsch {
6473fdc5746SBarry Smith   PetscFunctionBegin;
6482fa5cd67SKarl Rupp   while (n--) {
6492fa5cd67SKarl Rupp     *arg1 = PetscMin(*arg1, *arg2);
6502fa5cd67SKarl Rupp     arg1++;
6512fa5cd67SKarl Rupp     arg2++;
6522fa5cd67SKarl Rupp   }
653*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
654827bd09bSSatish Balay }
655827bd09bSSatish Balay 
6567b1ae94cSBarry Smith /*********************************ivec.c*************************************/
657d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_rvec_min_abs(PetscScalar *arg1, PetscScalar *arg2, PetscInt n)
658d71ae5a4SJacob Faibussowitsch {
6593fdc5746SBarry Smith   PetscFunctionBegin;
6602fa5cd67SKarl Rupp   while (n--) {
6612fa5cd67SKarl Rupp     *arg1 = MIN_FABS(*arg1, *arg2);
6622fa5cd67SKarl Rupp     arg1++;
6632fa5cd67SKarl Rupp     arg2++;
6642fa5cd67SKarl Rupp   }
665*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
666827bd09bSSatish Balay }
667827bd09bSSatish Balay 
6687b1ae94cSBarry Smith /*********************************ivec.c*************************************/
669d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_rvec_exists(PetscScalar *arg1, PetscScalar *arg2, PetscInt n)
670d71ae5a4SJacob Faibussowitsch {
6713fdc5746SBarry Smith   PetscFunctionBegin;
6722fa5cd67SKarl Rupp   while (n--) {
6732fa5cd67SKarl Rupp     *arg1 = EXISTS(*arg1, *arg2);
6742fa5cd67SKarl Rupp     arg1++;
6752fa5cd67SKarl Rupp     arg2++;
6762fa5cd67SKarl Rupp   }
677*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
678827bd09bSSatish Balay }
679827bd09bSSatish Balay 
6807b1ae94cSBarry Smith /***********************************ivec.c*************************************/
681d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_rvec_non_uniform(PetscScalar *arg1, PetscScalar *arg2, PetscInt n, PetscInt *arg3)
682d71ae5a4SJacob Faibussowitsch {
68352f87cdaSBarry Smith   PetscInt i, j, type;
684827bd09bSSatish Balay 
6853fdc5746SBarry Smith   PetscFunctionBegin;
686827bd09bSSatish Balay   /* LATER: if we're really motivated we can sort and then unsort */
687db4deed7SKarl Rupp   for (i = 0; i < n;) {
688827bd09bSSatish Balay     /* clump 'em for now */
689827bd09bSSatish Balay     j    = i + 1;
690827bd09bSSatish Balay     type = arg3[i];
6912fa5cd67SKarl Rupp     while ((j < n) && (arg3[j] == type)) j++;
692827bd09bSSatish Balay 
693827bd09bSSatish Balay     /* how many together */
694827bd09bSSatish Balay     j -= i;
695827bd09bSSatish Balay 
696827bd09bSSatish Balay     /* call appropriate ivec function */
697*3ba16761SJacob Faibussowitsch     if (type == GL_MAX) PetscCall(PCTFS_rvec_max(arg1, arg2, j));
698*3ba16761SJacob Faibussowitsch     else if (type == GL_MIN) PetscCall(PCTFS_rvec_min(arg1, arg2, j));
699*3ba16761SJacob Faibussowitsch     else if (type == GL_MULT) PetscCall(PCTFS_rvec_mult(arg1, arg2, j));
700*3ba16761SJacob Faibussowitsch     else if (type == GL_ADD) PetscCall(PCTFS_rvec_add(arg1, arg2, j));
701*3ba16761SJacob Faibussowitsch     else if (type == GL_MAX_ABS) PetscCall(PCTFS_rvec_max_abs(arg1, arg2, j));
702*3ba16761SJacob Faibussowitsch     else if (type == GL_MIN_ABS) PetscCall(PCTFS_rvec_min_abs(arg1, arg2, j));
703*3ba16761SJacob Faibussowitsch     else if (type == GL_EXISTS) PetscCall(PCTFS_rvec_exists(arg1, arg2, j));
704ca8e9878SJed Brown     else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "unrecognized type passed to PCTFS_rvec_non_uniform()!");
705827bd09bSSatish Balay 
7069371c9d4SSatish Balay     arg1 += j;
7079371c9d4SSatish Balay     arg2 += j;
7089371c9d4SSatish Balay     i += j;
709827bd09bSSatish Balay   }
710*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
711827bd09bSSatish Balay }
712827bd09bSSatish Balay 
7137b1ae94cSBarry Smith /***********************************ivec.c*************************************/
714d71ae5a4SJacob Faibussowitsch vfp PCTFS_rvec_fct_addr(PetscInt type)
715d71ae5a4SJacob Faibussowitsch {
7162fa5cd67SKarl Rupp   if (type == NON_UNIFORM) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_rvec_non_uniform);
7172fa5cd67SKarl Rupp   else if (type == GL_MAX) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_rvec_max);
7182fa5cd67SKarl Rupp   else if (type == GL_MIN) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_rvec_min);
7192fa5cd67SKarl Rupp   else if (type == GL_MULT) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_rvec_mult);
7202fa5cd67SKarl Rupp   else if (type == GL_ADD) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_rvec_add);
7212fa5cd67SKarl Rupp   else if (type == GL_MAX_ABS) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_rvec_max_abs);
7222fa5cd67SKarl Rupp   else if (type == GL_MIN_ABS) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_rvec_min_abs);
7232fa5cd67SKarl Rupp   else if (type == GL_EXISTS) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_rvec_exists);
724827bd09bSSatish Balay 
725827bd09bSSatish Balay   /* catch all ... not good if we get here */
726827bd09bSSatish Balay   return (NULL);
727827bd09bSSatish Balay }
728