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