xref: /petsc/src/ksp/pc/impls/tfs/ivec.c (revision 57508ece14a6b1339c0bbf016ecd72f673a062b0)
1827bd09bSSatish Balay /**********************************ivec.c**************************************
2827bd09bSSatish Balay 
3827bd09bSSatish Balay Author: Henry M. Tufo III
4827bd09bSSatish Balay 
5827bd09bSSatish Balay e-mail: hmt@cs.brown.edu
6827bd09bSSatish Balay 
7827bd09bSSatish Balay snail-mail:
8827bd09bSSatish Balay Division of Applied Mathematics
9827bd09bSSatish Balay Brown University
10827bd09bSSatish Balay Providence, RI 02912
11827bd09bSSatish Balay 
12827bd09bSSatish Balay Last Modification:
13827bd09bSSatish Balay 6.21.97
14827bd09bSSatish Balay ***********************************ivec.c*************************************/
15827bd09bSSatish Balay 
16c6db04a5SJed Brown #include <../src/ksp/pc/impls/tfs/tfs.h>
17827bd09bSSatish Balay 
18827bd09bSSatish Balay /* sorting args ivec.c ivec.c ... */
19827bd09bSSatish Balay #define SORT_OPT   6
20827bd09bSSatish Balay #define SORT_STACK 50000
21827bd09bSSatish Balay 
22827bd09bSSatish Balay /* allocate an address and size stack for sorter(s) */
23827bd09bSSatish Balay static void    *offset_stack[2 * SORT_STACK];
2452f87cdaSBarry Smith static PetscInt size_stack[SORT_STACK];
25827bd09bSSatish Balay 
267b1ae94cSBarry Smith /***********************************ivec.c*************************************/
27d71ae5a4SJacob Faibussowitsch PetscInt *PCTFS_ivec_copy(PetscInt *arg1, PetscInt *arg2, PetscInt n)
28d71ae5a4SJacob Faibussowitsch {
292fa5cd67SKarl Rupp   while (n--) *arg1++ = *arg2++;
304ad8454bSPierre Jolivet   return arg1;
31827bd09bSSatish Balay }
32827bd09bSSatish Balay 
337b1ae94cSBarry Smith /***********************************ivec.c*************************************/
34d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_ivec_zero(PetscInt *arg1, PetscInt n)
35d71ae5a4SJacob Faibussowitsch {
363fdc5746SBarry Smith   PetscFunctionBegin;
372fa5cd67SKarl Rupp   while (n--) *arg1++ = 0;
383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
39827bd09bSSatish Balay }
40827bd09bSSatish Balay 
417b1ae94cSBarry Smith /***********************************ivec.c*************************************/
42d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_ivec_set(PetscInt *arg1, PetscInt arg2, PetscInt n)
43d71ae5a4SJacob Faibussowitsch {
443fdc5746SBarry Smith   PetscFunctionBegin;
452fa5cd67SKarl Rupp   while (n--) *arg1++ = arg2;
463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
47827bd09bSSatish Balay }
48827bd09bSSatish Balay 
497b1ae94cSBarry Smith /***********************************ivec.c*************************************/
50d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_ivec_max(PetscInt *arg1, PetscInt *arg2, PetscInt n)
51d71ae5a4SJacob Faibussowitsch {
523fdc5746SBarry Smith   PetscFunctionBegin;
539371c9d4SSatish Balay   while (n--) {
549371c9d4SSatish Balay     *arg1 = PetscMax(*arg1, *arg2);
559371c9d4SSatish Balay     arg1++;
569371c9d4SSatish Balay     arg2++;
579371c9d4SSatish Balay   }
583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
59827bd09bSSatish Balay }
60827bd09bSSatish Balay 
617b1ae94cSBarry Smith /***********************************ivec.c*************************************/
62d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_ivec_min(PetscInt *arg1, PetscInt *arg2, PetscInt n)
63d71ae5a4SJacob Faibussowitsch {
643fdc5746SBarry Smith   PetscFunctionBegin;
652fa5cd67SKarl Rupp   while (n--) {
662fa5cd67SKarl Rupp     *(arg1) = PetscMin(*arg1, *arg2);
672fa5cd67SKarl Rupp     arg1++;
682fa5cd67SKarl Rupp     arg2++;
692fa5cd67SKarl Rupp   }
703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
71827bd09bSSatish Balay }
72827bd09bSSatish Balay 
737b1ae94cSBarry Smith /***********************************ivec.c*************************************/
74d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_ivec_mult(PetscInt *arg1, PetscInt *arg2, PetscInt n)
75d71ae5a4SJacob Faibussowitsch {
763fdc5746SBarry Smith   PetscFunctionBegin;
772fa5cd67SKarl Rupp   while (n--) *arg1++ *= *arg2++;
783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
79827bd09bSSatish Balay }
80827bd09bSSatish Balay 
817b1ae94cSBarry Smith /***********************************ivec.c*************************************/
82d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_ivec_add(PetscInt *arg1, PetscInt *arg2, PetscInt n)
83d71ae5a4SJacob Faibussowitsch {
843fdc5746SBarry Smith   PetscFunctionBegin;
852fa5cd67SKarl Rupp   while (n--) *arg1++ += *arg2++;
863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
87827bd09bSSatish Balay }
88827bd09bSSatish Balay 
897b1ae94cSBarry Smith /***********************************ivec.c*************************************/
90d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_ivec_lxor(PetscInt *arg1, PetscInt *arg2, PetscInt n)
91d71ae5a4SJacob Faibussowitsch {
923fdc5746SBarry Smith   PetscFunctionBegin;
932fa5cd67SKarl Rupp   while (n--) {
94*57508eceSPierre Jolivet     *arg1 = (*arg1 || *arg2) && !(*arg1 && *arg2);
952fa5cd67SKarl Rupp     arg1++;
962fa5cd67SKarl Rupp     arg2++;
972fa5cd67SKarl Rupp   }
983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
99827bd09bSSatish Balay }
100827bd09bSSatish Balay 
1017b1ae94cSBarry Smith /***********************************ivec.c*************************************/
102d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_ivec_xor(PetscInt *arg1, PetscInt *arg2, PetscInt n)
103d71ae5a4SJacob Faibussowitsch {
1043fdc5746SBarry Smith   PetscFunctionBegin;
1052fa5cd67SKarl Rupp   while (n--) *arg1++ ^= *arg2++;
1063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
107827bd09bSSatish Balay }
108827bd09bSSatish Balay 
1097b1ae94cSBarry Smith /***********************************ivec.c*************************************/
110d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_ivec_or(PetscInt *arg1, PetscInt *arg2, PetscInt n)
111d71ae5a4SJacob Faibussowitsch {
1123fdc5746SBarry Smith   PetscFunctionBegin;
1132fa5cd67SKarl Rupp   while (n--) *arg1++ |= *arg2++;
1143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
115827bd09bSSatish Balay }
116827bd09bSSatish Balay 
1177b1ae94cSBarry Smith /***********************************ivec.c*************************************/
118d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_ivec_lor(PetscInt *arg1, PetscInt *arg2, PetscInt n)
119d71ae5a4SJacob Faibussowitsch {
1203fdc5746SBarry Smith   PetscFunctionBegin;
1212fa5cd67SKarl Rupp   while (n--) {
1222fa5cd67SKarl Rupp     *arg1 = (*arg1 || *arg2);
1232fa5cd67SKarl Rupp     arg1++;
1242fa5cd67SKarl Rupp     arg2++;
1252fa5cd67SKarl Rupp   }
1263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
127827bd09bSSatish Balay }
128827bd09bSSatish Balay 
1297b1ae94cSBarry Smith /***********************************ivec.c*************************************/
130d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_ivec_and(PetscInt *arg1, PetscInt *arg2, PetscInt n)
131d71ae5a4SJacob Faibussowitsch {
1323fdc5746SBarry Smith   PetscFunctionBegin;
1332fa5cd67SKarl Rupp   while (n--) *arg1++ &= *arg2++;
1343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
135827bd09bSSatish Balay }
136827bd09bSSatish Balay 
1377b1ae94cSBarry Smith /***********************************ivec.c*************************************/
138d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_ivec_land(PetscInt *arg1, PetscInt *arg2, PetscInt n)
139d71ae5a4SJacob Faibussowitsch {
1403fdc5746SBarry Smith   PetscFunctionBegin;
1412fa5cd67SKarl Rupp   while (n--) {
1422fa5cd67SKarl Rupp     *arg1 = (*arg1 && *arg2);
1432fa5cd67SKarl Rupp     arg1++;
1442fa5cd67SKarl Rupp     arg2++;
1452fa5cd67SKarl Rupp   }
1463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
147827bd09bSSatish Balay }
148827bd09bSSatish Balay 
1497b1ae94cSBarry Smith /***********************************ivec.c*************************************/
150d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_ivec_and3(PetscInt *arg1, PetscInt *arg2, PetscInt *arg3, PetscInt n)
151d71ae5a4SJacob Faibussowitsch {
1523fdc5746SBarry Smith   PetscFunctionBegin;
1532fa5cd67SKarl Rupp   while (n--) *arg1++ = (*arg2++ & *arg3++);
1543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
155827bd09bSSatish Balay }
156827bd09bSSatish Balay 
1577b1ae94cSBarry Smith /***********************************ivec.c*************************************/
158d71ae5a4SJacob Faibussowitsch PetscInt PCTFS_ivec_sum(PetscInt *arg1, PetscInt n)
159d71ae5a4SJacob Faibussowitsch {
16052f87cdaSBarry Smith   PetscInt tmp = 0;
1612fa5cd67SKarl Rupp   while (n--) tmp += *arg1++;
1624ad8454bSPierre Jolivet   return tmp;
163827bd09bSSatish Balay }
164827bd09bSSatish Balay 
1657b1ae94cSBarry Smith /***********************************ivec.c*************************************/
166d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_ivec_non_uniform(PetscInt *arg1, PetscInt *arg2, PetscInt n, ...)
167d71ae5a4SJacob Faibussowitsch {
16852f87cdaSBarry Smith   PetscInt  i, j, type;
1696ae68b59SBarry Smith   PetscInt *arg3;
1706ae68b59SBarry Smith   va_list   ap;
171827bd09bSSatish Balay 
1723fdc5746SBarry Smith   PetscFunctionBegin;
1736ae68b59SBarry Smith   va_start(ap, n);
1746ae68b59SBarry Smith   arg3 = va_arg(ap, PetscInt *);
1756ae68b59SBarry Smith   va_end(ap);
1766ae68b59SBarry Smith 
177827bd09bSSatish Balay   /* LATER: if we're really motivated we can sort and then unsort */
178db4deed7SKarl Rupp   for (i = 0; i < n;) {
179827bd09bSSatish Balay     /* clump 'em for now */
180827bd09bSSatish Balay     j    = i + 1;
181827bd09bSSatish Balay     type = arg3[i];
1822fa5cd67SKarl Rupp     while ((j < n) && (arg3[j] == type)) j++;
183827bd09bSSatish Balay 
184827bd09bSSatish Balay     /* how many together */
185827bd09bSSatish Balay     j -= i;
186827bd09bSSatish Balay 
187827bd09bSSatish Balay     /* call appropriate ivec function */
1883ba16761SJacob Faibussowitsch     if (type == GL_MAX) PetscCall(PCTFS_ivec_max(arg1, arg2, j));
1893ba16761SJacob Faibussowitsch     else if (type == GL_MIN) PetscCall(PCTFS_ivec_min(arg1, arg2, j));
1903ba16761SJacob Faibussowitsch     else if (type == GL_MULT) PetscCall(PCTFS_ivec_mult(arg1, arg2, j));
1913ba16761SJacob Faibussowitsch     else if (type == GL_ADD) PetscCall(PCTFS_ivec_add(arg1, arg2, j));
1923ba16761SJacob Faibussowitsch     else if (type == GL_B_XOR) PetscCall(PCTFS_ivec_xor(arg1, arg2, j));
1933ba16761SJacob Faibussowitsch     else if (type == GL_B_OR) PetscCall(PCTFS_ivec_or(arg1, arg2, j));
1943ba16761SJacob Faibussowitsch     else if (type == GL_B_AND) PetscCall(PCTFS_ivec_and(arg1, arg2, j));
1953ba16761SJacob Faibussowitsch     else if (type == GL_L_XOR) PetscCall(PCTFS_ivec_lxor(arg1, arg2, j));
1963ba16761SJacob Faibussowitsch     else if (type == GL_L_OR) PetscCall(PCTFS_ivec_lor(arg1, arg2, j));
1973ba16761SJacob Faibussowitsch     else if (type == GL_L_AND) PetscCall(PCTFS_ivec_land(arg1, arg2, j));
198db4deed7SKarl Rupp     else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "unrecognized type passed to PCTFS_ivec_non_uniform()!");
199827bd09bSSatish Balay 
2009371c9d4SSatish Balay     arg1 += j;
2019371c9d4SSatish Balay     arg2 += j;
2029371c9d4SSatish Balay     i += j;
203827bd09bSSatish Balay   }
2043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
205827bd09bSSatish Balay }
206827bd09bSSatish Balay 
2077b1ae94cSBarry Smith /***********************************ivec.c*************************************/
208d71ae5a4SJacob Faibussowitsch vfp PCTFS_ivec_fct_addr(PetscInt type)
209d71ae5a4SJacob Faibussowitsch {
2104ad8454bSPierre Jolivet   if (type == NON_UNIFORM) return (PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_non_uniform;
2114ad8454bSPierre Jolivet   else if (type == GL_MAX) return (PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_max;
2124ad8454bSPierre Jolivet   else if (type == GL_MIN) return (PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_min;
2134ad8454bSPierre Jolivet   else if (type == GL_MULT) return (PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_mult;
2144ad8454bSPierre Jolivet   else if (type == GL_ADD) return (PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_add;
2154ad8454bSPierre Jolivet   else if (type == GL_B_XOR) return (PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_xor;
2164ad8454bSPierre Jolivet   else if (type == GL_B_OR) return (PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_or;
2174ad8454bSPierre Jolivet   else if (type == GL_B_AND) return (PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_and;
2184ad8454bSPierre Jolivet   else if (type == GL_L_XOR) return (PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_lxor;
2194ad8454bSPierre Jolivet   else if (type == GL_L_OR) return (PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_lor;
2204ad8454bSPierre Jolivet   else if (type == GL_L_AND) return (PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_land;
221827bd09bSSatish Balay 
222827bd09bSSatish Balay   /* catch all ... not good if we get here */
2234ad8454bSPierre Jolivet   return NULL;
224827bd09bSSatish Balay }
225827bd09bSSatish Balay 
2267b1ae94cSBarry Smith /******************************************************************************/
227d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_ivec_sort(PetscInt *ar, PetscInt size)
228d71ae5a4SJacob Faibussowitsch {
22952f87cdaSBarry Smith   PetscInt  *pi, *pj, temp;
23052f87cdaSBarry Smith   PetscInt **top_a = (PetscInt **)offset_stack;
23152f87cdaSBarry Smith   PetscInt  *top_s = size_stack, *bottom_s = size_stack;
232827bd09bSSatish Balay 
233b458e8f1SJose E. Roman   PetscFunctionBegin;
234827bd09bSSatish Balay   /* we're really interested in the offset of the last element */
235827bd09bSSatish Balay   /* ==> length of the list is now size + 1                    */
236827bd09bSSatish Balay   size--;
237827bd09bSSatish Balay 
238827bd09bSSatish Balay   /* do until we're done ... return when stack is exhausted */
239db4deed7SKarl Rupp   for (;;) {
240827bd09bSSatish Balay     /* if list is large enough use quicksort partition exchange code */
241db4deed7SKarl Rupp     if (size > SORT_OPT) {
242827bd09bSSatish Balay       /* start up pointer at element 1 and down at size     */
243827bd09bSSatish Balay       pi = ar + 1;
244827bd09bSSatish Balay       pj = ar + size;
245827bd09bSSatish Balay 
246827bd09bSSatish Balay       /* find middle element in list and swap w/ element 1 */
247827bd09bSSatish Balay       SWAP(*(ar + (size >> 1)), *pi)
248827bd09bSSatish Balay 
249827bd09bSSatish Balay       /* order element 0,1,size-1 st {M,L,...,U} w/L<=M<=U */
250827bd09bSSatish Balay       /* note ==> pivot_value in index 0                   */
251db4deed7SKarl Rupp       if (*pi > *pj) { SWAP(*pi, *pj) }
2529371c9d4SSatish Balay       if (*ar > *pj) {
2539371c9d4SSatish Balay         SWAP(*ar, *pj)
2549371c9d4SSatish Balay       } else if (*pi > *ar) {
2559371c9d4SSatish Balay         SWAP(*(ar), *(ar + 1))
2569371c9d4SSatish Balay       }
257827bd09bSSatish Balay 
258827bd09bSSatish Balay       /* partition about pivot_value ...                              */
259827bd09bSSatish Balay       /* note lists of length 2 are not guaranteed to be sorted */
260db4deed7SKarl Rupp       for (;;) {
261827bd09bSSatish Balay         /* walk up ... and down ... swap if equal to pivot! */
2629371c9d4SSatish Balay         do pi++;
2639371c9d4SSatish Balay         while (*pi < *ar);
2649371c9d4SSatish Balay         do pj--;
2659371c9d4SSatish Balay         while (*pj > *ar);
266827bd09bSSatish Balay 
267827bd09bSSatish Balay         /* if we've crossed we're done */
268827bd09bSSatish Balay         if (pj < pi) break;
269827bd09bSSatish Balay 
270827bd09bSSatish Balay         /* else swap */
271827bd09bSSatish Balay         SWAP(*pi, *pj)
272827bd09bSSatish Balay       }
273827bd09bSSatish Balay 
274827bd09bSSatish Balay       /* place pivot_value in it's correct location */
275827bd09bSSatish Balay       SWAP(*ar, *pj)
276827bd09bSSatish Balay 
277827bd09bSSatish Balay       /* test stack_size to see if we've exhausted our stack */
27808401ef6SPierre Jolivet       PetscCheck(top_s - bottom_s < SORT_STACK, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCTFS_ivec_sort() :: STACK EXHAUSTED!!!");
279827bd09bSSatish Balay 
280827bd09bSSatish Balay       /* push right hand child iff length > 1 */
281db4deed7SKarl Rupp       if ((*top_s = size - ((PetscInt)(pi - ar)))) {
282827bd09bSSatish Balay         *(top_a++) = pi;
283827bd09bSSatish Balay         size -= *top_s + 2;
284827bd09bSSatish Balay         top_s++;
2859371c9d4SSatish Balay       } else if (size -= *top_s + 2)
2869371c9d4SSatish Balay         ; /* set up for next loop iff there is something to do */
2879371c9d4SSatish Balay       else { /* might as well pop - note NR_OPT >=2 ==> we're ok! */ ar = *(--top_a);
288827bd09bSSatish Balay         size                                                            = *(--top_s);
289827bd09bSSatish Balay       }
290db4deed7SKarl Rupp     } else { /* else sort small list directly then pop another off stack */
291827bd09bSSatish Balay 
292827bd09bSSatish Balay       /* insertion sort for bottom */
293db4deed7SKarl Rupp       for (pj = ar + 1; pj <= ar + size; pj++) {
294827bd09bSSatish Balay         temp = *pj;
295db4deed7SKarl Rupp         for (pi = pj - 1; pi >= ar; pi--) {
296827bd09bSSatish Balay           if (*pi <= temp) break;
297827bd09bSSatish Balay           *(pi + 1) = *pi;
298827bd09bSSatish Balay         }
299827bd09bSSatish Balay         *(pi + 1) = temp;
300827bd09bSSatish Balay       }
301827bd09bSSatish Balay 
302827bd09bSSatish Balay       /* check to see if stack is exhausted ==> DONE */
3033ba16761SJacob Faibussowitsch       if (top_s == bottom_s) PetscFunctionReturn(PETSC_SUCCESS);
304827bd09bSSatish Balay 
305827bd09bSSatish Balay       /* else pop another list from the stack */
306827bd09bSSatish Balay       ar   = *(--top_a);
307827bd09bSSatish Balay       size = *(--top_s);
308827bd09bSSatish Balay     }
309827bd09bSSatish Balay   }
310827bd09bSSatish Balay }
311827bd09bSSatish Balay 
3127b1ae94cSBarry Smith /******************************************************************************/
313d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_ivec_sort_companion(PetscInt *ar, PetscInt *ar2, PetscInt size)
314d71ae5a4SJacob Faibussowitsch {
31552f87cdaSBarry Smith   PetscInt  *pi, *pj, temp, temp2;
31652f87cdaSBarry Smith   PetscInt **top_a = (PetscInt **)offset_stack;
31752f87cdaSBarry Smith   PetscInt  *top_s = size_stack, *bottom_s = size_stack;
31852f87cdaSBarry Smith   PetscInt  *pi2, *pj2;
31952f87cdaSBarry Smith   PetscInt   mid;
320827bd09bSSatish Balay 
3213fdc5746SBarry Smith   PetscFunctionBegin;
322827bd09bSSatish Balay   /* we're really interested in the offset of the last element */
323827bd09bSSatish Balay   /* ==> length of the list is now size + 1                    */
324827bd09bSSatish Balay   size--;
325827bd09bSSatish Balay 
326827bd09bSSatish Balay   /* do until we're done ... return when stack is exhausted */
327db4deed7SKarl Rupp   for (;;) {
328827bd09bSSatish Balay     /* if list is large enough use quicksort partition exchange code */
329db4deed7SKarl Rupp     if (size > SORT_OPT) {
330827bd09bSSatish Balay       /* start up pointer at element 1 and down at size     */
331827bd09bSSatish Balay       mid = size >> 1;
332827bd09bSSatish Balay       pi  = ar + 1;
333827bd09bSSatish Balay       pj  = ar + mid;
334827bd09bSSatish Balay       pi2 = ar2 + 1;
335827bd09bSSatish Balay       pj2 = ar2 + mid;
336827bd09bSSatish Balay 
337827bd09bSSatish Balay       /* find middle element in list and swap w/ element 1 */
338827bd09bSSatish Balay       SWAP(*pi, *pj)
339827bd09bSSatish Balay       SWAP(*pi2, *pj2)
340827bd09bSSatish Balay 
341827bd09bSSatish Balay       /* order element 0,1,size-1 st {M,L,...,U} w/L<=M<=U */
342827bd09bSSatish Balay       /* note ==> pivot_value in index 0                   */
343827bd09bSSatish Balay       pj  = ar + size;
344827bd09bSSatish Balay       pj2 = ar2 + size;
345db4deed7SKarl Rupp       if (*pi > *pj) { SWAP(*pi, *pj) SWAP(*pi2, *pj2) }
3469371c9d4SSatish Balay       if (*ar > *pj) {
3479371c9d4SSatish Balay         SWAP(*ar, *pj) SWAP(*ar2, *pj2)
3489371c9d4SSatish Balay       } else if (*pi > *ar) {
3499371c9d4SSatish Balay         SWAP(*(ar), *(ar + 1)) SWAP(*(ar2), *(ar2 + 1))
3509371c9d4SSatish Balay       }
351827bd09bSSatish Balay 
352827bd09bSSatish Balay       /* partition about pivot_value ...                              */
353827bd09bSSatish Balay       /* note lists of length 2 are not guaranteed to be sorted */
354db4deed7SKarl Rupp       for (;;) {
355827bd09bSSatish Balay         /* walk up ... and down ... swap if equal to pivot! */
3569371c9d4SSatish Balay         do {
3579371c9d4SSatish Balay           pi++;
3589371c9d4SSatish Balay           pi2++;
3599371c9d4SSatish Balay         } while (*pi < *ar);
3609371c9d4SSatish Balay         do {
3619371c9d4SSatish Balay           pj--;
3629371c9d4SSatish Balay           pj2--;
3639371c9d4SSatish Balay         } while (*pj > *ar);
364827bd09bSSatish Balay 
365827bd09bSSatish Balay         /* if we've crossed we're done */
366827bd09bSSatish Balay         if (pj < pi) break;
367827bd09bSSatish Balay 
368827bd09bSSatish Balay         /* else swap */
369827bd09bSSatish Balay         SWAP(*pi, *pj)
370827bd09bSSatish Balay         SWAP(*pi2, *pj2)
371827bd09bSSatish Balay       }
372827bd09bSSatish Balay 
373827bd09bSSatish Balay       /* place pivot_value in it's correct location */
374827bd09bSSatish Balay       SWAP(*ar, *pj)
375827bd09bSSatish Balay       SWAP(*ar2, *pj2)
376827bd09bSSatish Balay 
377827bd09bSSatish Balay       /* test stack_size to see if we've exhausted our stack */
37808401ef6SPierre Jolivet       PetscCheck(top_s - bottom_s < SORT_STACK, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCTFS_ivec_sort_companion() :: STACK EXHAUSTED!!!");
379827bd09bSSatish Balay 
380827bd09bSSatish Balay       /* push right hand child iff length > 1 */
381db4deed7SKarl Rupp       if ((*top_s = size - ((PetscInt)(pi - ar)))) {
382827bd09bSSatish Balay         *(top_a++) = pi;
383827bd09bSSatish Balay         *(top_a++) = pi2;
384827bd09bSSatish Balay         size -= *top_s + 2;
385827bd09bSSatish Balay         top_s++;
3869371c9d4SSatish Balay       } else if (size -= *top_s + 2)
3879371c9d4SSatish Balay         ; /* set up for next loop iff there is something to do */
3889371c9d4SSatish Balay       else { /* might as well pop - note NR_OPT >=2 ==> we're ok! */ ar2 = *(--top_a);
389827bd09bSSatish Balay         ar                                                               = *(--top_a);
390827bd09bSSatish Balay         size                                                             = *(--top_s);
391827bd09bSSatish Balay       }
392db4deed7SKarl Rupp     } else { /* else sort small list directly then pop another off stack */
393827bd09bSSatish Balay 
394827bd09bSSatish Balay       /* insertion sort for bottom */
395db4deed7SKarl Rupp       for (pj = ar + 1, pj2 = ar2 + 1; pj <= ar + size; pj++, pj2++) {
396827bd09bSSatish Balay         temp  = *pj;
397827bd09bSSatish Balay         temp2 = *pj2;
398db4deed7SKarl Rupp         for (pi = pj - 1, pi2 = pj2 - 1; pi >= ar; pi--, pi2--) {
399827bd09bSSatish Balay           if (*pi <= temp) break;
400827bd09bSSatish Balay           *(pi + 1)  = *pi;
401827bd09bSSatish Balay           *(pi2 + 1) = *pi2;
402827bd09bSSatish Balay         }
403827bd09bSSatish Balay         *(pi + 1)  = temp;
404827bd09bSSatish Balay         *(pi2 + 1) = temp2;
405827bd09bSSatish Balay       }
406827bd09bSSatish Balay 
407827bd09bSSatish Balay       /* check to see if stack is exhausted ==> DONE */
4083ba16761SJacob Faibussowitsch       if (top_s == bottom_s) PetscFunctionReturn(PETSC_SUCCESS);
409827bd09bSSatish Balay 
410827bd09bSSatish Balay       /* else pop another list from the stack */
411827bd09bSSatish Balay       ar2  = *(--top_a);
412827bd09bSSatish Balay       ar   = *(--top_a);
413827bd09bSSatish Balay       size = *(--top_s);
414827bd09bSSatish Balay     }
415827bd09bSSatish Balay   }
416827bd09bSSatish Balay }
417827bd09bSSatish Balay 
4187b1ae94cSBarry Smith /******************************************************************************/
419d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_ivec_sort_companion_hack(PetscInt *ar, PetscInt **ar2, PetscInt size)
420d71ae5a4SJacob Faibussowitsch {
42152f87cdaSBarry Smith   PetscInt  *pi, *pj, temp, *ptr;
42252f87cdaSBarry Smith   PetscInt **top_a = (PetscInt **)offset_stack;
42352f87cdaSBarry Smith   PetscInt  *top_s = size_stack, *bottom_s = size_stack;
42452f87cdaSBarry Smith   PetscInt **pi2, **pj2;
42552f87cdaSBarry Smith   PetscInt   mid;
426827bd09bSSatish Balay 
4273fdc5746SBarry Smith   PetscFunctionBegin;
428827bd09bSSatish Balay   /* we're really interested in the offset of the last element */
429827bd09bSSatish Balay   /* ==> length of the list is now size + 1                    */
430827bd09bSSatish Balay   size--;
431827bd09bSSatish Balay 
432827bd09bSSatish Balay   /* do until we're done ... return when stack is exhausted */
433db4deed7SKarl Rupp   for (;;) {
434827bd09bSSatish Balay     /* if list is large enough use quicksort partition exchange code */
435db4deed7SKarl Rupp     if (size > SORT_OPT) {
436827bd09bSSatish Balay       /* start up pointer at element 1 and down at size     */
437827bd09bSSatish Balay       mid = size >> 1;
438827bd09bSSatish Balay       pi  = ar + 1;
439827bd09bSSatish Balay       pj  = ar + mid;
440827bd09bSSatish Balay       pi2 = ar2 + 1;
441827bd09bSSatish Balay       pj2 = ar2 + mid;
442827bd09bSSatish Balay 
443827bd09bSSatish Balay       /* find middle element in list and swap w/ element 1 */
444827bd09bSSatish Balay       SWAP(*pi, *pj)
445827bd09bSSatish Balay       P_SWAP(*pi2, *pj2)
446827bd09bSSatish Balay 
447827bd09bSSatish Balay       /* order element 0,1,size-1 st {M,L,...,U} w/L<=M<=U */
448827bd09bSSatish Balay       /* note ==> pivot_value in index 0                   */
449827bd09bSSatish Balay       pj  = ar + size;
450827bd09bSSatish Balay       pj2 = ar2 + size;
451db4deed7SKarl Rupp       if (*pi > *pj) { SWAP(*pi, *pj) P_SWAP(*pi2, *pj2) }
4529371c9d4SSatish Balay       if (*ar > *pj) {
4539371c9d4SSatish Balay         SWAP(*ar, *pj) P_SWAP(*ar2, *pj2)
4549371c9d4SSatish Balay       } else if (*pi > *ar) {
4559371c9d4SSatish Balay         SWAP(*(ar), *(ar + 1)) P_SWAP(*(ar2), *(ar2 + 1))
4569371c9d4SSatish Balay       }
457827bd09bSSatish Balay 
458827bd09bSSatish Balay       /* partition about pivot_value ...                              */
459827bd09bSSatish Balay       /* note lists of length 2 are not guaranteed to be sorted */
460db4deed7SKarl Rupp       for (;;) {
461827bd09bSSatish Balay         /* walk up ... and down ... swap if equal to pivot! */
4629371c9d4SSatish Balay         do {
4639371c9d4SSatish Balay           pi++;
4649371c9d4SSatish Balay           pi2++;
4659371c9d4SSatish Balay         } while (*pi < *ar);
4669371c9d4SSatish Balay         do {
4679371c9d4SSatish Balay           pj--;
4689371c9d4SSatish Balay           pj2--;
4699371c9d4SSatish Balay         } while (*pj > *ar);
470827bd09bSSatish Balay 
471827bd09bSSatish Balay         /* if we've crossed we're done */
472827bd09bSSatish Balay         if (pj < pi) break;
473827bd09bSSatish Balay 
474827bd09bSSatish Balay         /* else swap */
475827bd09bSSatish Balay         SWAP(*pi, *pj)
476827bd09bSSatish Balay         P_SWAP(*pi2, *pj2)
477827bd09bSSatish Balay       }
478827bd09bSSatish Balay 
479827bd09bSSatish Balay       /* place pivot_value in it's correct location */
480827bd09bSSatish Balay       SWAP(*ar, *pj)
481827bd09bSSatish Balay       P_SWAP(*ar2, *pj2)
482827bd09bSSatish Balay 
483827bd09bSSatish Balay       /* test stack_size to see if we've exhausted our stack */
48408401ef6SPierre Jolivet       PetscCheck(top_s - bottom_s < SORT_STACK, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCTFS_ivec_sort_companion_hack() :: STACK EXHAUSTED!!!");
485827bd09bSSatish Balay 
486827bd09bSSatish Balay       /* push right hand child iff length > 1 */
487db4deed7SKarl Rupp       if ((*top_s = size - ((PetscInt)(pi - ar)))) {
488827bd09bSSatish Balay         *(top_a++) = pi;
48952f87cdaSBarry Smith         *(top_a++) = (PetscInt *)pi2;
490827bd09bSSatish Balay         size -= *top_s + 2;
491827bd09bSSatish Balay         top_s++;
4929371c9d4SSatish Balay       } else if (size -= *top_s + 2)
4939371c9d4SSatish Balay         ; /* set up for next loop iff there is something to do */
4949371c9d4SSatish Balay       else { /* might as well pop - note NR_OPT >=2 ==> we're ok! */ ar2 = (PetscInt **)*(--top_a);
495827bd09bSSatish Balay         ar                                                               = *(--top_a);
496827bd09bSSatish Balay         size                                                             = *(--top_s);
497827bd09bSSatish Balay       }
4982fa5cd67SKarl Rupp     } else { /* else sort small list directly then pop another off stack */
499827bd09bSSatish Balay       /* insertion sort for bottom */
500db4deed7SKarl Rupp       for (pj = ar + 1, pj2 = ar2 + 1; pj <= ar + size; pj++, pj2++) {
501827bd09bSSatish Balay         temp = *pj;
502827bd09bSSatish Balay         ptr  = *pj2;
503db4deed7SKarl Rupp         for (pi = pj - 1, pi2 = pj2 - 1; pi >= ar; pi--, pi2--) {
504827bd09bSSatish Balay           if (*pi <= temp) break;
505827bd09bSSatish Balay           *(pi + 1)  = *pi;
506827bd09bSSatish Balay           *(pi2 + 1) = *pi2;
507827bd09bSSatish Balay         }
508827bd09bSSatish Balay         *(pi + 1)  = temp;
509827bd09bSSatish Balay         *(pi2 + 1) = ptr;
510827bd09bSSatish Balay       }
511827bd09bSSatish Balay 
512827bd09bSSatish Balay       /* check to see if stack is exhausted ==> DONE */
5133ba16761SJacob Faibussowitsch       if (top_s == bottom_s) PetscFunctionReturn(PETSC_SUCCESS);
514827bd09bSSatish Balay 
515827bd09bSSatish Balay       /* else pop another list from the stack */
51652f87cdaSBarry Smith       ar2  = (PetscInt **)*(--top_a);
517827bd09bSSatish Balay       ar   = *(--top_a);
518827bd09bSSatish Balay       size = *(--top_s);
519827bd09bSSatish Balay     }
520827bd09bSSatish Balay   }
521827bd09bSSatish Balay }
522827bd09bSSatish Balay 
5237b1ae94cSBarry Smith /******************************************************************************/
524d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_SMI_sort(void *ar1, void *ar2, PetscInt size, PetscInt type)
525d71ae5a4SJacob Faibussowitsch {
5263fdc5746SBarry Smith   PetscFunctionBegin;
527e7e72b3dSBarry Smith   if (type == SORT_INTEGER) {
5283ba16761SJacob Faibussowitsch     if (ar2) PetscCall(PCTFS_ivec_sort_companion((PetscInt *)ar1, (PetscInt *)ar2, size));
5293ba16761SJacob Faibussowitsch     else PetscCall(PCTFS_ivec_sort((PetscInt *)ar1, size));
530e7e72b3dSBarry Smith   } else if (type == SORT_INT_PTR) {
5313ba16761SJacob Faibussowitsch     if (ar2) PetscCall(PCTFS_ivec_sort_companion_hack((PetscInt *)ar1, (PetscInt **)ar2, size));
5323ba16761SJacob Faibussowitsch     else PetscCall(PCTFS_ivec_sort((PetscInt *)ar1, size));
533ca8e9878SJed Brown   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCTFS_SMI_sort only does SORT_INTEGER!");
5343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
535827bd09bSSatish Balay }
536827bd09bSSatish Balay 
5377b1ae94cSBarry Smith /***********************************ivec.c*************************************/
538d71ae5a4SJacob Faibussowitsch PetscInt PCTFS_ivec_linear_search(PetscInt item, PetscInt *list, PetscInt n)
539d71ae5a4SJacob Faibussowitsch {
54052f87cdaSBarry Smith   PetscInt tmp = n - 1;
5415fd66863SKarl Rupp 
5422fa5cd67SKarl Rupp   while (n--) {
5434ad8454bSPierre Jolivet     if (*list++ == item) return tmp - n;
5442fa5cd67SKarl Rupp   }
5454ad8454bSPierre Jolivet   return -1;
546827bd09bSSatish Balay }
547827bd09bSSatish Balay 
5487b1ae94cSBarry Smith /***********************************ivec.c*************************************/
549d71ae5a4SJacob Faibussowitsch PetscInt PCTFS_ivec_binary_search(PetscInt item, PetscInt *list, PetscInt rh)
550d71ae5a4SJacob Faibussowitsch {
55152f87cdaSBarry Smith   PetscInt mid, lh = 0;
552827bd09bSSatish Balay 
553827bd09bSSatish Balay   rh--;
554db4deed7SKarl Rupp   while (lh <= rh) {
555827bd09bSSatish Balay     mid = (lh + rh) >> 1;
5564ad8454bSPierre Jolivet     if (*(list + mid) == item) return mid;
5572fa5cd67SKarl Rupp     if (*(list + mid) > item) rh = mid - 1;
5582fa5cd67SKarl Rupp     else lh = mid + 1;
559827bd09bSSatish Balay   }
5604ad8454bSPierre Jolivet   return -1;
561827bd09bSSatish Balay }
562827bd09bSSatish Balay 
5637b1ae94cSBarry Smith /*********************************ivec.c*************************************/
564d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_rvec_copy(PetscScalar *arg1, PetscScalar *arg2, PetscInt n)
565d71ae5a4SJacob Faibussowitsch {
5663fdc5746SBarry Smith   PetscFunctionBegin;
5672fa5cd67SKarl Rupp   while (n--) *arg1++ = *arg2++;
5683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
569827bd09bSSatish Balay }
570827bd09bSSatish Balay 
5717b1ae94cSBarry Smith /*********************************ivec.c*************************************/
572d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_rvec_zero(PetscScalar *arg1, PetscInt n)
573d71ae5a4SJacob Faibussowitsch {
5743fdc5746SBarry Smith   PetscFunctionBegin;
5752fa5cd67SKarl Rupp   while (n--) *arg1++ = 0.0;
5763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
577827bd09bSSatish Balay }
578827bd09bSSatish Balay 
5797b1ae94cSBarry Smith /***********************************ivec.c*************************************/
580d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_rvec_one(PetscScalar *arg1, PetscInt n)
581d71ae5a4SJacob Faibussowitsch {
5823fdc5746SBarry Smith   PetscFunctionBegin;
5832fa5cd67SKarl Rupp   while (n--) *arg1++ = 1.0;
5843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
585827bd09bSSatish Balay }
586827bd09bSSatish Balay 
5877b1ae94cSBarry Smith /***********************************ivec.c*************************************/
588d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_rvec_set(PetscScalar *arg1, PetscScalar arg2, PetscInt n)
589d71ae5a4SJacob Faibussowitsch {
5903fdc5746SBarry Smith   PetscFunctionBegin;
5912fa5cd67SKarl Rupp   while (n--) *arg1++ = arg2;
5923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
593827bd09bSSatish Balay }
594827bd09bSSatish Balay 
5957b1ae94cSBarry Smith /***********************************ivec.c*************************************/
596d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_rvec_scale(PetscScalar *arg1, PetscScalar arg2, PetscInt n)
597d71ae5a4SJacob Faibussowitsch {
5983fdc5746SBarry Smith   PetscFunctionBegin;
5992fa5cd67SKarl Rupp   while (n--) *arg1++ *= arg2;
6003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
601827bd09bSSatish Balay }
602827bd09bSSatish Balay 
6037b1ae94cSBarry Smith /*********************************ivec.c*************************************/
604d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_rvec_add(PetscScalar *arg1, PetscScalar *arg2, PetscInt n)
605d71ae5a4SJacob Faibussowitsch {
6063fdc5746SBarry Smith   PetscFunctionBegin;
6072fa5cd67SKarl Rupp   while (n--) *arg1++ += *arg2++;
6083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
609827bd09bSSatish Balay }
610827bd09bSSatish Balay 
6117b1ae94cSBarry Smith /*********************************ivec.c*************************************/
612d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_rvec_mult(PetscScalar *arg1, PetscScalar *arg2, PetscInt n)
613d71ae5a4SJacob Faibussowitsch {
6143fdc5746SBarry Smith   PetscFunctionBegin;
6152fa5cd67SKarl Rupp   while (n--) *arg1++ *= *arg2++;
6163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
617827bd09bSSatish Balay }
618827bd09bSSatish Balay 
6197b1ae94cSBarry Smith /*********************************ivec.c*************************************/
620d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_rvec_max(PetscScalar *arg1, PetscScalar *arg2, PetscInt n)
621d71ae5a4SJacob Faibussowitsch {
6223fdc5746SBarry Smith   PetscFunctionBegin;
6232fa5cd67SKarl Rupp   while (n--) {
6242fa5cd67SKarl Rupp     *arg1 = PetscMax(*arg1, *arg2);
6252fa5cd67SKarl Rupp     arg1++;
6262fa5cd67SKarl Rupp     arg2++;
6272fa5cd67SKarl Rupp   }
6283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
629827bd09bSSatish Balay }
630827bd09bSSatish Balay 
6317b1ae94cSBarry Smith /*********************************ivec.c*************************************/
632d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_rvec_max_abs(PetscScalar *arg1, PetscScalar *arg2, PetscInt n)
633d71ae5a4SJacob Faibussowitsch {
6343fdc5746SBarry Smith   PetscFunctionBegin;
6352fa5cd67SKarl Rupp   while (n--) {
6362fa5cd67SKarl Rupp     *arg1 = MAX_FABS(*arg1, *arg2);
6372fa5cd67SKarl Rupp     arg1++;
6382fa5cd67SKarl Rupp     arg2++;
6392fa5cd67SKarl Rupp   }
6403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
641827bd09bSSatish Balay }
642827bd09bSSatish Balay 
6437b1ae94cSBarry Smith /*********************************ivec.c*************************************/
644d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_rvec_min(PetscScalar *arg1, PetscScalar *arg2, PetscInt n)
645d71ae5a4SJacob Faibussowitsch {
6463fdc5746SBarry Smith   PetscFunctionBegin;
6472fa5cd67SKarl Rupp   while (n--) {
6482fa5cd67SKarl Rupp     *arg1 = PetscMin(*arg1, *arg2);
6492fa5cd67SKarl Rupp     arg1++;
6502fa5cd67SKarl Rupp     arg2++;
6512fa5cd67SKarl Rupp   }
6523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
653827bd09bSSatish Balay }
654827bd09bSSatish Balay 
6557b1ae94cSBarry Smith /*********************************ivec.c*************************************/
656d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_rvec_min_abs(PetscScalar *arg1, PetscScalar *arg2, PetscInt n)
657d71ae5a4SJacob Faibussowitsch {
6583fdc5746SBarry Smith   PetscFunctionBegin;
6592fa5cd67SKarl Rupp   while (n--) {
6602fa5cd67SKarl Rupp     *arg1 = MIN_FABS(*arg1, *arg2);
6612fa5cd67SKarl Rupp     arg1++;
6622fa5cd67SKarl Rupp     arg2++;
6632fa5cd67SKarl Rupp   }
6643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
665827bd09bSSatish Balay }
666827bd09bSSatish Balay 
6677b1ae94cSBarry Smith /*********************************ivec.c*************************************/
66866976f2fSJacob Faibussowitsch static PetscErrorCode PCTFS_rvec_exists(PetscScalar *arg1, PetscScalar *arg2, PetscInt n)
669d71ae5a4SJacob Faibussowitsch {
6703fdc5746SBarry Smith   PetscFunctionBegin;
6712fa5cd67SKarl Rupp   while (n--) {
6722fa5cd67SKarl Rupp     *arg1 = EXISTS(*arg1, *arg2);
6732fa5cd67SKarl Rupp     arg1++;
6742fa5cd67SKarl Rupp     arg2++;
6752fa5cd67SKarl Rupp   }
6763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
677827bd09bSSatish Balay }
678827bd09bSSatish Balay 
6797b1ae94cSBarry Smith /***********************************ivec.c*************************************/
68066976f2fSJacob Faibussowitsch static PetscErrorCode PCTFS_rvec_non_uniform(PetscScalar *arg1, PetscScalar *arg2, PetscInt n, PetscInt *arg3)
681d71ae5a4SJacob Faibussowitsch {
68252f87cdaSBarry Smith   PetscInt i, j, type;
683827bd09bSSatish Balay 
6843fdc5746SBarry Smith   PetscFunctionBegin;
685827bd09bSSatish Balay   /* LATER: if we're really motivated we can sort and then unsort */
686db4deed7SKarl Rupp   for (i = 0; i < n;) {
687827bd09bSSatish Balay     /* clump 'em for now */
688827bd09bSSatish Balay     j    = i + 1;
689827bd09bSSatish Balay     type = arg3[i];
6902fa5cd67SKarl Rupp     while ((j < n) && (arg3[j] == type)) j++;
691827bd09bSSatish Balay 
692827bd09bSSatish Balay     /* how many together */
693827bd09bSSatish Balay     j -= i;
694827bd09bSSatish Balay 
695827bd09bSSatish Balay     /* call appropriate ivec function */
6963ba16761SJacob Faibussowitsch     if (type == GL_MAX) PetscCall(PCTFS_rvec_max(arg1, arg2, j));
6973ba16761SJacob Faibussowitsch     else if (type == GL_MIN) PetscCall(PCTFS_rvec_min(arg1, arg2, j));
6983ba16761SJacob Faibussowitsch     else if (type == GL_MULT) PetscCall(PCTFS_rvec_mult(arg1, arg2, j));
6993ba16761SJacob Faibussowitsch     else if (type == GL_ADD) PetscCall(PCTFS_rvec_add(arg1, arg2, j));
7003ba16761SJacob Faibussowitsch     else if (type == GL_MAX_ABS) PetscCall(PCTFS_rvec_max_abs(arg1, arg2, j));
7013ba16761SJacob Faibussowitsch     else if (type == GL_MIN_ABS) PetscCall(PCTFS_rvec_min_abs(arg1, arg2, j));
7023ba16761SJacob Faibussowitsch     else if (type == GL_EXISTS) PetscCall(PCTFS_rvec_exists(arg1, arg2, j));
703ca8e9878SJed Brown     else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "unrecognized type passed to PCTFS_rvec_non_uniform()!");
704827bd09bSSatish Balay 
7059371c9d4SSatish Balay     arg1 += j;
7069371c9d4SSatish Balay     arg2 += j;
7079371c9d4SSatish Balay     i += j;
708827bd09bSSatish Balay   }
7093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
710827bd09bSSatish Balay }
711827bd09bSSatish Balay 
7127b1ae94cSBarry Smith /***********************************ivec.c*************************************/
713d71ae5a4SJacob Faibussowitsch vfp PCTFS_rvec_fct_addr(PetscInt type)
714d71ae5a4SJacob Faibussowitsch {
7154ad8454bSPierre Jolivet   if (type == NON_UNIFORM) return (PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_rvec_non_uniform;
7164ad8454bSPierre Jolivet   else if (type == GL_MAX) return (PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_rvec_max;
7174ad8454bSPierre Jolivet   else if (type == GL_MIN) return (PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_rvec_min;
7184ad8454bSPierre Jolivet   else if (type == GL_MULT) return (PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_rvec_mult;
7194ad8454bSPierre Jolivet   else if (type == GL_ADD) return (PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_rvec_add;
7204ad8454bSPierre Jolivet   else if (type == GL_MAX_ABS) return (PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_rvec_max_abs;
7214ad8454bSPierre Jolivet   else if (type == GL_MIN_ABS) return (PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_rvec_min_abs;
7224ad8454bSPierre Jolivet   else if (type == GL_EXISTS) return (PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_rvec_exists;
723827bd09bSSatish Balay 
724827bd09bSSatish Balay   /* catch all ... not good if we get here */
7254ad8454bSPierre Jolivet   return NULL;
726827bd09bSSatish Balay }
727