xref: /petsc/src/ksp/pc/impls/tfs/ivec.c (revision 58d68138c660dfb4e9f5b03334792cd4f2ffd7cc)
1 
2 /**********************************ivec.c**************************************
3 
4 Author: Henry M. Tufo III
5 
6 e-mail: hmt@cs.brown.edu
7 
8 snail-mail:
9 Division of Applied Mathematics
10 Brown University
11 Providence, RI 02912
12 
13 Last Modification:
14 6.21.97
15 ***********************************ivec.c*************************************/
16 
17 #include <../src/ksp/pc/impls/tfs/tfs.h>
18 
19 /* sorting args ivec.c ivec.c ... */
20 #define SORT_OPT   6
21 #define SORT_STACK 50000
22 
23 /* allocate an address and size stack for sorter(s) */
24 static void    *offset_stack[2 * SORT_STACK];
25 static PetscInt size_stack[SORT_STACK];
26 
27 /***********************************ivec.c*************************************/
28 PetscInt *PCTFS_ivec_copy(PetscInt *arg1, PetscInt *arg2, PetscInt n) {
29   while (n--) *arg1++ = *arg2++;
30   return (arg1);
31 }
32 
33 /***********************************ivec.c*************************************/
34 PetscErrorCode PCTFS_ivec_zero(PetscInt *arg1, PetscInt n) {
35   PetscFunctionBegin;
36   while (n--) *arg1++ = 0;
37   PetscFunctionReturn(0);
38 }
39 
40 /***********************************ivec.c*************************************/
41 PetscErrorCode PCTFS_ivec_set(PetscInt *arg1, PetscInt arg2, PetscInt n) {
42   PetscFunctionBegin;
43   while (n--) *arg1++ = arg2;
44   PetscFunctionReturn(0);
45 }
46 
47 /***********************************ivec.c*************************************/
48 PetscErrorCode PCTFS_ivec_max(PetscInt *arg1, PetscInt *arg2, PetscInt n) {
49   PetscFunctionBegin;
50   while (n--) {
51     *arg1 = PetscMax(*arg1, *arg2);
52     arg1++;
53     arg2++;
54   }
55   PetscFunctionReturn(0);
56 }
57 
58 /***********************************ivec.c*************************************/
59 PetscErrorCode PCTFS_ivec_min(PetscInt *arg1, PetscInt *arg2, PetscInt n) {
60   PetscFunctionBegin;
61   while (n--) {
62     *(arg1) = PetscMin(*arg1, *arg2);
63     arg1++;
64     arg2++;
65   }
66   PetscFunctionReturn(0);
67 }
68 
69 /***********************************ivec.c*************************************/
70 PetscErrorCode PCTFS_ivec_mult(PetscInt *arg1, PetscInt *arg2, PetscInt n) {
71   PetscFunctionBegin;
72   while (n--) *arg1++ *= *arg2++;
73   PetscFunctionReturn(0);
74 }
75 
76 /***********************************ivec.c*************************************/
77 PetscErrorCode PCTFS_ivec_add(PetscInt *arg1, PetscInt *arg2, PetscInt n) {
78   PetscFunctionBegin;
79   while (n--) *arg1++ += *arg2++;
80   PetscFunctionReturn(0);
81 }
82 
83 /***********************************ivec.c*************************************/
84 PetscErrorCode PCTFS_ivec_lxor(PetscInt *arg1, PetscInt *arg2, PetscInt n) {
85   PetscFunctionBegin;
86   while (n--) {
87     *arg1 = ((*arg1 || *arg2) && !(*arg1 && *arg2));
88     arg1++;
89     arg2++;
90   }
91   PetscFunctionReturn(0);
92 }
93 
94 /***********************************ivec.c*************************************/
95 PetscErrorCode PCTFS_ivec_xor(PetscInt *arg1, PetscInt *arg2, PetscInt n) {
96   PetscFunctionBegin;
97   while (n--) *arg1++ ^= *arg2++;
98   PetscFunctionReturn(0);
99 }
100 
101 /***********************************ivec.c*************************************/
102 PetscErrorCode PCTFS_ivec_or(PetscInt *arg1, PetscInt *arg2, PetscInt n) {
103   PetscFunctionBegin;
104   while (n--) *arg1++ |= *arg2++;
105   PetscFunctionReturn(0);
106 }
107 
108 /***********************************ivec.c*************************************/
109 PetscErrorCode PCTFS_ivec_lor(PetscInt *arg1, PetscInt *arg2, PetscInt n) {
110   PetscFunctionBegin;
111   while (n--) {
112     *arg1 = (*arg1 || *arg2);
113     arg1++;
114     arg2++;
115   }
116   PetscFunctionReturn(0);
117 }
118 
119 /***********************************ivec.c*************************************/
120 PetscErrorCode PCTFS_ivec_and(PetscInt *arg1, PetscInt *arg2, PetscInt n) {
121   PetscFunctionBegin;
122   while (n--) *arg1++ &= *arg2++;
123   PetscFunctionReturn(0);
124 }
125 
126 /***********************************ivec.c*************************************/
127 PetscErrorCode PCTFS_ivec_land(PetscInt *arg1, PetscInt *arg2, PetscInt n) {
128   PetscFunctionBegin;
129   while (n--) {
130     *arg1 = (*arg1 && *arg2);
131     arg1++;
132     arg2++;
133   }
134   PetscFunctionReturn(0);
135 }
136 
137 /***********************************ivec.c*************************************/
138 PetscErrorCode PCTFS_ivec_and3(PetscInt *arg1, PetscInt *arg2, PetscInt *arg3, PetscInt n) {
139   PetscFunctionBegin;
140   while (n--) *arg1++ = (*arg2++ & *arg3++);
141   PetscFunctionReturn(0);
142 }
143 
144 /***********************************ivec.c*************************************/
145 PetscInt PCTFS_ivec_sum(PetscInt *arg1, PetscInt n) {
146   PetscInt tmp = 0;
147   while (n--) tmp += *arg1++;
148   return (tmp);
149 }
150 
151 /***********************************ivec.c*************************************/
152 PetscErrorCode PCTFS_ivec_non_uniform(PetscInt *arg1, PetscInt *arg2, PetscInt n, ...) {
153   PetscInt  i, j, type;
154   PetscInt *arg3;
155   va_list   ap;
156 
157   PetscFunctionBegin;
158   va_start(ap, n);
159   arg3 = va_arg(ap, PetscInt *);
160   va_end(ap);
161 
162   /* LATER: if we're really motivated we can sort and then unsort */
163   for (i = 0; i < n;) {
164     /* clump 'em for now */
165     j    = i + 1;
166     type = arg3[i];
167     while ((j < n) && (arg3[j] == type)) j++;
168 
169     /* how many together */
170     j -= i;
171 
172     /* call appropriate ivec function */
173     if (type == GL_MAX) PCTFS_ivec_max(arg1, arg2, j);
174     else if (type == GL_MIN) PCTFS_ivec_min(arg1, arg2, j);
175     else if (type == GL_MULT) PCTFS_ivec_mult(arg1, arg2, j);
176     else if (type == GL_ADD) PCTFS_ivec_add(arg1, arg2, j);
177     else if (type == GL_B_XOR) PCTFS_ivec_xor(arg1, arg2, j);
178     else if (type == GL_B_OR) PCTFS_ivec_or(arg1, arg2, j);
179     else if (type == GL_B_AND) PCTFS_ivec_and(arg1, arg2, j);
180     else if (type == GL_L_XOR) PCTFS_ivec_lxor(arg1, arg2, j);
181     else if (type == GL_L_OR) PCTFS_ivec_lor(arg1, arg2, j);
182     else if (type == GL_L_AND) PCTFS_ivec_land(arg1, arg2, j);
183     else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "unrecognized type passed to PCTFS_ivec_non_uniform()!");
184 
185     arg1 += j;
186     arg2 += j;
187     i += j;
188   }
189   PetscFunctionReturn(0);
190 }
191 
192 /***********************************ivec.c*************************************/
193 vfp PCTFS_ivec_fct_addr(PetscInt type) {
194   if (type == NON_UNIFORM) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_non_uniform);
195   else if (type == GL_MAX) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_max);
196   else if (type == GL_MIN) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_min);
197   else if (type == GL_MULT) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_mult);
198   else if (type == GL_ADD) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_add);
199   else if (type == GL_B_XOR) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_xor);
200   else if (type == GL_B_OR) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_or);
201   else if (type == GL_B_AND) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_and);
202   else if (type == GL_L_XOR) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_lxor);
203   else if (type == GL_L_OR) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_lor);
204   else if (type == GL_L_AND) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_ivec_land);
205 
206   /* catch all ... not good if we get here */
207   return (NULL);
208 }
209 
210 /******************************************************************************/
211 PetscErrorCode PCTFS_ivec_sort(PetscInt *ar, PetscInt size) {
212   PetscInt  *pi, *pj, temp;
213   PetscInt **top_a = (PetscInt **)offset_stack;
214   PetscInt  *top_s = size_stack, *bottom_s = size_stack;
215 
216   PetscFunctionBegin;
217   /* we're really interested in the offset of the last element */
218   /* ==> length of the list is now size + 1                    */
219   size--;
220 
221   /* do until we're done ... return when stack is exhausted */
222   for (;;) {
223     /* if list is large enough use quicksort partition exchange code */
224     if (size > SORT_OPT) {
225       /* start up pointer at element 1 and down at size     */
226       pi = ar + 1;
227       pj = ar + size;
228 
229       /* find middle element in list and swap w/ element 1 */
230       SWAP(*(ar + (size >> 1)), *pi)
231 
232       /* order element 0,1,size-1 st {M,L,...,U} w/L<=M<=U */
233       /* note ==> pivot_value in index 0                   */
234       if (*pi > *pj) { SWAP(*pi, *pj) }
235       if (*ar > *pj) {
236         SWAP(*ar, *pj)
237       } else if (*pi > *ar) {
238         SWAP(*(ar), *(ar + 1))
239       }
240 
241       /* partition about pivot_value ...                              */
242       /* note lists of length 2 are not guaranteed to be sorted */
243       for (;;) {
244         /* walk up ... and down ... swap if equal to pivot! */
245         do pi++;
246         while (*pi < *ar);
247         do pj--;
248         while (*pj > *ar);
249 
250         /* if we've crossed we're done */
251         if (pj < pi) break;
252 
253         /* else swap */
254         SWAP(*pi, *pj)
255       }
256 
257       /* place pivot_value in it's correct location */
258       SWAP(*ar, *pj)
259 
260       /* test stack_size to see if we've exhausted our stack */
261       PetscCheck(top_s - bottom_s < SORT_STACK, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCTFS_ivec_sort() :: STACK EXHAUSTED!!!");
262 
263       /* push right hand child iff length > 1 */
264       if ((*top_s = size - ((PetscInt)(pi - ar)))) {
265         *(top_a++) = pi;
266         size -= *top_s + 2;
267         top_s++;
268       } else if (size -= *top_s + 2)
269         ; /* set up for next loop iff there is something to do */
270       else { /* might as well pop - note NR_OPT >=2 ==> we're ok! */ ar = *(--top_a);
271         size                                                            = *(--top_s);
272       }
273     } else { /* else sort small list directly then pop another off stack */
274 
275       /* insertion sort for bottom */
276       for (pj = ar + 1; pj <= ar + size; pj++) {
277         temp = *pj;
278         for (pi = pj - 1; pi >= ar; pi--) {
279           if (*pi <= temp) break;
280           *(pi + 1) = *pi;
281         }
282         *(pi + 1) = temp;
283       }
284 
285       /* check to see if stack is exhausted ==> DONE */
286       if (top_s == bottom_s) PetscFunctionReturn(0);
287 
288       /* else pop another list from the stack */
289       ar   = *(--top_a);
290       size = *(--top_s);
291     }
292   }
293 }
294 
295 /******************************************************************************/
296 PetscErrorCode PCTFS_ivec_sort_companion(PetscInt *ar, PetscInt *ar2, PetscInt size) {
297   PetscInt  *pi, *pj, temp, temp2;
298   PetscInt **top_a = (PetscInt **)offset_stack;
299   PetscInt  *top_s = size_stack, *bottom_s = size_stack;
300   PetscInt  *pi2, *pj2;
301   PetscInt   mid;
302 
303   PetscFunctionBegin;
304   /* we're really interested in the offset of the last element */
305   /* ==> length of the list is now size + 1                    */
306   size--;
307 
308   /* do until we're done ... return when stack is exhausted */
309   for (;;) {
310     /* if list is large enough use quicksort partition exchange code */
311     if (size > SORT_OPT) {
312       /* start up pointer at element 1 and down at size     */
313       mid = size >> 1;
314       pi  = ar + 1;
315       pj  = ar + mid;
316       pi2 = ar2 + 1;
317       pj2 = ar2 + mid;
318 
319       /* find middle element in list and swap w/ element 1 */
320       SWAP(*pi, *pj)
321       SWAP(*pi2, *pj2)
322 
323       /* order element 0,1,size-1 st {M,L,...,U} w/L<=M<=U */
324       /* note ==> pivot_value in index 0                   */
325       pj  = ar + size;
326       pj2 = ar2 + size;
327       if (*pi > *pj) { SWAP(*pi, *pj) SWAP(*pi2, *pj2) }
328       if (*ar > *pj) {
329         SWAP(*ar, *pj) SWAP(*ar2, *pj2)
330       } else if (*pi > *ar) {
331         SWAP(*(ar), *(ar + 1)) SWAP(*(ar2), *(ar2 + 1))
332       }
333 
334       /* partition about pivot_value ...                              */
335       /* note lists of length 2 are not guaranteed to be sorted */
336       for (;;) {
337         /* walk up ... and down ... swap if equal to pivot! */
338         do {
339           pi++;
340           pi2++;
341         } while (*pi < *ar);
342         do {
343           pj--;
344           pj2--;
345         } while (*pj > *ar);
346 
347         /* if we've crossed we're done */
348         if (pj < pi) break;
349 
350         /* else swap */
351         SWAP(*pi, *pj)
352         SWAP(*pi2, *pj2)
353       }
354 
355       /* place pivot_value in it's correct location */
356       SWAP(*ar, *pj)
357       SWAP(*ar2, *pj2)
358 
359       /* test stack_size to see if we've exhausted our stack */
360       PetscCheck(top_s - bottom_s < SORT_STACK, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCTFS_ivec_sort_companion() :: STACK EXHAUSTED!!!");
361 
362       /* push right hand child iff length > 1 */
363       if ((*top_s = size - ((PetscInt)(pi - ar)))) {
364         *(top_a++) = pi;
365         *(top_a++) = pi2;
366         size -= *top_s + 2;
367         top_s++;
368       } else if (size -= *top_s + 2)
369         ; /* set up for next loop iff there is something to do */
370       else { /* might as well pop - note NR_OPT >=2 ==> we're ok! */ ar2 = *(--top_a);
371         ar                                                               = *(--top_a);
372         size                                                             = *(--top_s);
373       }
374     } else { /* else sort small list directly then pop another off stack */
375 
376       /* insertion sort for bottom */
377       for (pj = ar + 1, pj2 = ar2 + 1; pj <= ar + size; pj++, pj2++) {
378         temp  = *pj;
379         temp2 = *pj2;
380         for (pi = pj - 1, pi2 = pj2 - 1; pi >= ar; pi--, pi2--) {
381           if (*pi <= temp) break;
382           *(pi + 1)  = *pi;
383           *(pi2 + 1) = *pi2;
384         }
385         *(pi + 1)  = temp;
386         *(pi2 + 1) = temp2;
387       }
388 
389       /* check to see if stack is exhausted ==> DONE */
390       if (top_s == bottom_s) PetscFunctionReturn(0);
391 
392       /* else pop another list from the stack */
393       ar2  = *(--top_a);
394       ar   = *(--top_a);
395       size = *(--top_s);
396     }
397   }
398 }
399 
400 /******************************************************************************/
401 PetscErrorCode PCTFS_ivec_sort_companion_hack(PetscInt *ar, PetscInt **ar2, PetscInt size) {
402   PetscInt  *pi, *pj, temp, *ptr;
403   PetscInt **top_a = (PetscInt **)offset_stack;
404   PetscInt  *top_s = size_stack, *bottom_s = size_stack;
405   PetscInt **pi2, **pj2;
406   PetscInt   mid;
407 
408   PetscFunctionBegin;
409   /* we're really interested in the offset of the last element */
410   /* ==> length of the list is now size + 1                    */
411   size--;
412 
413   /* do until we're done ... return when stack is exhausted */
414   for (;;) {
415     /* if list is large enough use quicksort partition exchange code */
416     if (size > SORT_OPT) {
417       /* start up pointer at element 1 and down at size     */
418       mid = size >> 1;
419       pi  = ar + 1;
420       pj  = ar + mid;
421       pi2 = ar2 + 1;
422       pj2 = ar2 + mid;
423 
424       /* find middle element in list and swap w/ element 1 */
425       SWAP(*pi, *pj)
426       P_SWAP(*pi2, *pj2)
427 
428       /* order element 0,1,size-1 st {M,L,...,U} w/L<=M<=U */
429       /* note ==> pivot_value in index 0                   */
430       pj  = ar + size;
431       pj2 = ar2 + size;
432       if (*pi > *pj) { SWAP(*pi, *pj) P_SWAP(*pi2, *pj2) }
433       if (*ar > *pj) {
434         SWAP(*ar, *pj) P_SWAP(*ar2, *pj2)
435       } else if (*pi > *ar) {
436         SWAP(*(ar), *(ar + 1)) P_SWAP(*(ar2), *(ar2 + 1))
437       }
438 
439       /* partition about pivot_value ...                              */
440       /* note lists of length 2 are not guaranteed to be sorted */
441       for (;;) {
442         /* walk up ... and down ... swap if equal to pivot! */
443         do {
444           pi++;
445           pi2++;
446         } while (*pi < *ar);
447         do {
448           pj--;
449           pj2--;
450         } while (*pj > *ar);
451 
452         /* if we've crossed we're done */
453         if (pj < pi) break;
454 
455         /* else swap */
456         SWAP(*pi, *pj)
457         P_SWAP(*pi2, *pj2)
458       }
459 
460       /* place pivot_value in it's correct location */
461       SWAP(*ar, *pj)
462       P_SWAP(*ar2, *pj2)
463 
464       /* test stack_size to see if we've exhausted our stack */
465       PetscCheck(top_s - bottom_s < SORT_STACK, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCTFS_ivec_sort_companion_hack() :: STACK EXHAUSTED!!!");
466 
467       /* push right hand child iff length > 1 */
468       if ((*top_s = size - ((PetscInt)(pi - ar)))) {
469         *(top_a++) = pi;
470         *(top_a++) = (PetscInt *)pi2;
471         size -= *top_s + 2;
472         top_s++;
473       } else if (size -= *top_s + 2)
474         ; /* set up for next loop iff there is something to do */
475       else { /* might as well pop - note NR_OPT >=2 ==> we're ok! */ ar2 = (PetscInt **)*(--top_a);
476         ar                                                               = *(--top_a);
477         size                                                             = *(--top_s);
478       }
479     } else { /* else sort small list directly then pop another off stack */
480       /* insertion sort for bottom */
481       for (pj = ar + 1, pj2 = ar2 + 1; pj <= ar + size; pj++, pj2++) {
482         temp = *pj;
483         ptr  = *pj2;
484         for (pi = pj - 1, pi2 = pj2 - 1; pi >= ar; pi--, pi2--) {
485           if (*pi <= temp) break;
486           *(pi + 1)  = *pi;
487           *(pi2 + 1) = *pi2;
488         }
489         *(pi + 1)  = temp;
490         *(pi2 + 1) = ptr;
491       }
492 
493       /* check to see if stack is exhausted ==> DONE */
494       if (top_s == bottom_s) PetscFunctionReturn(0);
495 
496       /* else pop another list from the stack */
497       ar2  = (PetscInt **)*(--top_a);
498       ar   = *(--top_a);
499       size = *(--top_s);
500     }
501   }
502 }
503 
504 /******************************************************************************/
505 PetscErrorCode PCTFS_SMI_sort(void *ar1, void *ar2, PetscInt size, PetscInt type) {
506   PetscFunctionBegin;
507   if (type == SORT_INTEGER) {
508     if (ar2) PCTFS_ivec_sort_companion((PetscInt *)ar1, (PetscInt *)ar2, size);
509     else PCTFS_ivec_sort((PetscInt *)ar1, size);
510   } else if (type == SORT_INT_PTR) {
511     if (ar2) PCTFS_ivec_sort_companion_hack((PetscInt *)ar1, (PetscInt **)ar2, size);
512     else PCTFS_ivec_sort((PetscInt *)ar1, size);
513   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCTFS_SMI_sort only does SORT_INTEGER!");
514   PetscFunctionReturn(0);
515 }
516 
517 /***********************************ivec.c*************************************/
518 PetscInt PCTFS_ivec_linear_search(PetscInt item, PetscInt *list, PetscInt n) {
519   PetscInt tmp = n - 1;
520 
521   while (n--) {
522     if (*list++ == item) return (tmp - n);
523   }
524   return (-1);
525 }
526 
527 /***********************************ivec.c*************************************/
528 PetscInt PCTFS_ivec_binary_search(PetscInt item, PetscInt *list, PetscInt rh) {
529   PetscInt mid, lh = 0;
530 
531   rh--;
532   while (lh <= rh) {
533     mid = (lh + rh) >> 1;
534     if (*(list + mid) == item) return (mid);
535     if (*(list + mid) > item) rh = mid - 1;
536     else lh = mid + 1;
537   }
538   return (-1);
539 }
540 
541 /*********************************ivec.c*************************************/
542 PetscErrorCode PCTFS_rvec_copy(PetscScalar *arg1, PetscScalar *arg2, PetscInt n) {
543   PetscFunctionBegin;
544   while (n--) *arg1++ = *arg2++;
545   PetscFunctionReturn(0);
546 }
547 
548 /*********************************ivec.c*************************************/
549 PetscErrorCode PCTFS_rvec_zero(PetscScalar *arg1, PetscInt n) {
550   PetscFunctionBegin;
551   while (n--) *arg1++ = 0.0;
552   PetscFunctionReturn(0);
553 }
554 
555 /***********************************ivec.c*************************************/
556 PetscErrorCode PCTFS_rvec_one(PetscScalar *arg1, PetscInt n) {
557   PetscFunctionBegin;
558   while (n--) *arg1++ = 1.0;
559   PetscFunctionReturn(0);
560 }
561 
562 /***********************************ivec.c*************************************/
563 PetscErrorCode PCTFS_rvec_set(PetscScalar *arg1, PetscScalar arg2, PetscInt n) {
564   PetscFunctionBegin;
565   while (n--) *arg1++ = arg2;
566   PetscFunctionReturn(0);
567 }
568 
569 /***********************************ivec.c*************************************/
570 PetscErrorCode PCTFS_rvec_scale(PetscScalar *arg1, PetscScalar arg2, PetscInt n) {
571   PetscFunctionBegin;
572   while (n--) *arg1++ *= arg2;
573   PetscFunctionReturn(0);
574 }
575 
576 /*********************************ivec.c*************************************/
577 PetscErrorCode PCTFS_rvec_add(PetscScalar *arg1, PetscScalar *arg2, PetscInt n) {
578   PetscFunctionBegin;
579   while (n--) *arg1++ += *arg2++;
580   PetscFunctionReturn(0);
581 }
582 
583 /*********************************ivec.c*************************************/
584 PetscErrorCode PCTFS_rvec_mult(PetscScalar *arg1, PetscScalar *arg2, PetscInt n) {
585   PetscFunctionBegin;
586   while (n--) *arg1++ *= *arg2++;
587   PetscFunctionReturn(0);
588 }
589 
590 /*********************************ivec.c*************************************/
591 PetscErrorCode PCTFS_rvec_max(PetscScalar *arg1, PetscScalar *arg2, PetscInt n) {
592   PetscFunctionBegin;
593   while (n--) {
594     *arg1 = PetscMax(*arg1, *arg2);
595     arg1++;
596     arg2++;
597   }
598   PetscFunctionReturn(0);
599 }
600 
601 /*********************************ivec.c*************************************/
602 PetscErrorCode PCTFS_rvec_max_abs(PetscScalar *arg1, PetscScalar *arg2, PetscInt n) {
603   PetscFunctionBegin;
604   while (n--) {
605     *arg1 = MAX_FABS(*arg1, *arg2);
606     arg1++;
607     arg2++;
608   }
609   PetscFunctionReturn(0);
610 }
611 
612 /*********************************ivec.c*************************************/
613 PetscErrorCode PCTFS_rvec_min(PetscScalar *arg1, PetscScalar *arg2, PetscInt n) {
614   PetscFunctionBegin;
615   while (n--) {
616     *arg1 = PetscMin(*arg1, *arg2);
617     arg1++;
618     arg2++;
619   }
620   PetscFunctionReturn(0);
621 }
622 
623 /*********************************ivec.c*************************************/
624 PetscErrorCode PCTFS_rvec_min_abs(PetscScalar *arg1, PetscScalar *arg2, PetscInt n) {
625   PetscFunctionBegin;
626   while (n--) {
627     *arg1 = MIN_FABS(*arg1, *arg2);
628     arg1++;
629     arg2++;
630   }
631   PetscFunctionReturn(0);
632 }
633 
634 /*********************************ivec.c*************************************/
635 PetscErrorCode PCTFS_rvec_exists(PetscScalar *arg1, PetscScalar *arg2, PetscInt n) {
636   PetscFunctionBegin;
637   while (n--) {
638     *arg1 = EXISTS(*arg1, *arg2);
639     arg1++;
640     arg2++;
641   }
642   PetscFunctionReturn(0);
643 }
644 
645 /***********************************ivec.c*************************************/
646 PetscErrorCode PCTFS_rvec_non_uniform(PetscScalar *arg1, PetscScalar *arg2, PetscInt n, PetscInt *arg3) {
647   PetscInt i, j, type;
648 
649   PetscFunctionBegin;
650   /* LATER: if we're really motivated we can sort and then unsort */
651   for (i = 0; i < n;) {
652     /* clump 'em for now */
653     j    = i + 1;
654     type = arg3[i];
655     while ((j < n) && (arg3[j] == type)) j++;
656 
657     /* how many together */
658     j -= i;
659 
660     /* call appropriate ivec function */
661     if (type == GL_MAX) PCTFS_rvec_max(arg1, arg2, j);
662     else if (type == GL_MIN) PCTFS_rvec_min(arg1, arg2, j);
663     else if (type == GL_MULT) PCTFS_rvec_mult(arg1, arg2, j);
664     else if (type == GL_ADD) PCTFS_rvec_add(arg1, arg2, j);
665     else if (type == GL_MAX_ABS) PCTFS_rvec_max_abs(arg1, arg2, j);
666     else if (type == GL_MIN_ABS) PCTFS_rvec_min_abs(arg1, arg2, j);
667     else if (type == GL_EXISTS) PCTFS_rvec_exists(arg1, arg2, j);
668     else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "unrecognized type passed to PCTFS_rvec_non_uniform()!");
669 
670     arg1 += j;
671     arg2 += j;
672     i += j;
673   }
674   PetscFunctionReturn(0);
675 }
676 
677 /***********************************ivec.c*************************************/
678 vfp PCTFS_rvec_fct_addr(PetscInt type) {
679   if (type == NON_UNIFORM) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_rvec_non_uniform);
680   else if (type == GL_MAX) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_rvec_max);
681   else if (type == GL_MIN) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_rvec_min);
682   else if (type == GL_MULT) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_rvec_mult);
683   else if (type == GL_ADD) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_rvec_add);
684   else if (type == GL_MAX_ABS) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_rvec_max_abs);
685   else if (type == GL_MIN_ABS) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_rvec_min_abs);
686   else if (type == GL_EXISTS) return ((PetscErrorCode(*)(void *, void *, PetscInt, ...)) & PCTFS_rvec_exists);
687 
688   /* catch all ... not good if we get here */
689   return (NULL);
690 }
691