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