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