xref: /petsc/src/ksp/pc/impls/tfs/ivec.c (revision 7b1ae94c5cf0701b5dfcb9679471e0ea993d816f)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2827bd09bSSatish Balay 
3827bd09bSSatish Balay /**********************************ivec.c**************************************
4827bd09bSSatish Balay 
5827bd09bSSatish Balay Author: Henry M. Tufo III
6827bd09bSSatish Balay 
7827bd09bSSatish Balay e-mail: hmt@cs.brown.edu
8827bd09bSSatish Balay 
9827bd09bSSatish Balay snail-mail:
10827bd09bSSatish Balay Division of Applied Mathematics
11827bd09bSSatish Balay Brown University
12827bd09bSSatish Balay Providence, RI 02912
13827bd09bSSatish Balay 
14827bd09bSSatish Balay Last Modification:
15827bd09bSSatish Balay 6.21.97
16827bd09bSSatish Balay ***********************************ivec.c*************************************/
17827bd09bSSatish Balay 
187758a8cdSBarry Smith #include "src/ksp/pc/impls/tfs/tfs.h"
19827bd09bSSatish Balay 
20827bd09bSSatish Balay /* sorting args ivec.c ivec.c ... */
21827bd09bSSatish Balay #define   SORT_OPT	6
22827bd09bSSatish Balay #define   SORT_STACK	50000
23827bd09bSSatish Balay 
24827bd09bSSatish Balay 
25827bd09bSSatish Balay /* allocate an address and size stack for sorter(s) */
26827bd09bSSatish Balay static void *offset_stack[2*SORT_STACK];
27827bd09bSSatish Balay static int   size_stack[SORT_STACK];
28827bd09bSSatish Balay 
29*7b1ae94cSBarry Smith /***********************************ivec.c*************************************/
30a501084fSBarry Smith int *ivec_copy( int *arg1,  int *arg2,  int n)
31827bd09bSSatish Balay {
32827bd09bSSatish Balay   while (n--)  {*arg1++ = *arg2++;}
33827bd09bSSatish Balay   return(arg1);
34827bd09bSSatish Balay }
35827bd09bSSatish Balay 
36*7b1ae94cSBarry Smith /***********************************ivec.c*************************************/
370924e98cSBarry Smith PetscErrorCode ivec_zero( int *arg1,  int n)
38827bd09bSSatish Balay {
393fdc5746SBarry Smith   PetscFunctionBegin;
40827bd09bSSatish Balay   while (n--)  {*arg1++ = 0;}
413fdc5746SBarry Smith   PetscFunctionReturn(0);
42827bd09bSSatish Balay }
43827bd09bSSatish Balay 
44*7b1ae94cSBarry Smith /***********************************ivec.c*************************************/
450924e98cSBarry Smith PetscErrorCode ivec_set( int *arg1,  int arg2,  int n)
46827bd09bSSatish Balay {
473fdc5746SBarry Smith   PetscFunctionBegin;
48827bd09bSSatish Balay   while (n--)  {*arg1++ = arg2;}
493fdc5746SBarry Smith   PetscFunctionReturn(0);
50827bd09bSSatish Balay }
51827bd09bSSatish Balay 
52*7b1ae94cSBarry Smith /***********************************ivec.c*************************************/
530924e98cSBarry Smith PetscErrorCode ivec_max( int *arg1,  int *arg2,  int n)
54827bd09bSSatish Balay {
553fdc5746SBarry Smith   PetscFunctionBegin;
5639945688SSatish Balay   while (n--)  {*arg1 = PetscMax(*arg1,*arg2); arg1++; arg2++;}
573fdc5746SBarry Smith   PetscFunctionReturn(0);
58827bd09bSSatish Balay }
59827bd09bSSatish Balay 
60*7b1ae94cSBarry Smith /***********************************ivec.c*************************************/
610924e98cSBarry Smith PetscErrorCode ivec_min( int *arg1,  int *arg2,  int n)
62827bd09bSSatish Balay {
633fdc5746SBarry Smith   PetscFunctionBegin;
6439945688SSatish Balay   while (n--)  {*(arg1) = PetscMin(*arg1,*arg2); arg1++; arg2++;}
653fdc5746SBarry Smith   PetscFunctionReturn(0);
66827bd09bSSatish Balay }
67827bd09bSSatish Balay 
68*7b1ae94cSBarry Smith /***********************************ivec.c*************************************/
690924e98cSBarry Smith PetscErrorCode ivec_mult( int *arg1,  int *arg2,  int n)
70827bd09bSSatish Balay {
713fdc5746SBarry Smith   PetscFunctionBegin;
72827bd09bSSatish Balay   while (n--)  {*arg1++ *= *arg2++;}
733fdc5746SBarry Smith   PetscFunctionReturn(0);
74827bd09bSSatish Balay }
75827bd09bSSatish Balay 
76*7b1ae94cSBarry Smith /***********************************ivec.c*************************************/
770924e98cSBarry Smith PetscErrorCode ivec_add( int *arg1,  int *arg2,  int n)
78827bd09bSSatish Balay {
793fdc5746SBarry Smith   PetscFunctionBegin;
80827bd09bSSatish Balay   while (n--)  {*arg1++ += *arg2++;}
813fdc5746SBarry Smith   PetscFunctionReturn(0);
82827bd09bSSatish Balay }
83827bd09bSSatish Balay 
84*7b1ae94cSBarry Smith /***********************************ivec.c*************************************/
850924e98cSBarry Smith PetscErrorCode ivec_lxor( int *arg1,  int *arg2,  int n)
86827bd09bSSatish Balay {
873fdc5746SBarry Smith   PetscFunctionBegin;
88827bd09bSSatish Balay   while (n--) {*arg1=((*arg1 || *arg2) && !(*arg1 && *arg2)); arg1++; arg2++;}
893fdc5746SBarry Smith   PetscFunctionReturn(0);
90827bd09bSSatish Balay }
91827bd09bSSatish Balay 
92*7b1ae94cSBarry Smith /***********************************ivec.c*************************************/
930924e98cSBarry Smith PetscErrorCode ivec_xor( int *arg1,  int *arg2,  int n)
94827bd09bSSatish Balay {
953fdc5746SBarry Smith   PetscFunctionBegin;
96827bd09bSSatish Balay   while (n--)  {*arg1++ ^= *arg2++;}
973fdc5746SBarry Smith   PetscFunctionReturn(0);
98827bd09bSSatish Balay }
99827bd09bSSatish Balay 
100*7b1ae94cSBarry Smith /***********************************ivec.c*************************************/
1010924e98cSBarry Smith PetscErrorCode ivec_or( int *arg1,  int *arg2,  int n)
102827bd09bSSatish Balay {
1033fdc5746SBarry Smith   PetscFunctionBegin;
104827bd09bSSatish Balay   while (n--)  {*arg1++ |= *arg2++;}
1053fdc5746SBarry Smith   PetscFunctionReturn(0);
106827bd09bSSatish Balay }
107827bd09bSSatish Balay 
108*7b1ae94cSBarry Smith /***********************************ivec.c*************************************/
1090924e98cSBarry Smith PetscErrorCode ivec_lor( int *arg1,  int *arg2,  int n)
110827bd09bSSatish Balay {
1113fdc5746SBarry Smith   PetscFunctionBegin;
112827bd09bSSatish Balay   while (n--)  {*arg1 = (*arg1 || *arg2); arg1++; arg2++;}
1133fdc5746SBarry Smith   PetscFunctionReturn(0);
114827bd09bSSatish Balay }
115827bd09bSSatish Balay 
116*7b1ae94cSBarry Smith /***********************************ivec.c*************************************/
1170924e98cSBarry Smith PetscErrorCode ivec_and( int *arg1,  int *arg2,  int n)
118827bd09bSSatish Balay {
1193fdc5746SBarry Smith   PetscFunctionBegin;
120827bd09bSSatish Balay   while (n--)  {*arg1++ &= *arg2++;}
1213fdc5746SBarry Smith   PetscFunctionReturn(0);
122827bd09bSSatish Balay }
123827bd09bSSatish Balay 
124*7b1ae94cSBarry Smith /***********************************ivec.c*************************************/
1250924e98cSBarry Smith PetscErrorCode ivec_land( int *arg1,  int *arg2,  int n)
126827bd09bSSatish Balay {
1273fdc5746SBarry Smith   PetscFunctionBegin;
128827bd09bSSatish Balay   while (n--) {*arg1 = (*arg1 && *arg2); arg1++; arg2++;}
1293fdc5746SBarry Smith   PetscFunctionReturn(0);
130827bd09bSSatish Balay }
131827bd09bSSatish Balay 
132*7b1ae94cSBarry Smith /***********************************ivec.c*************************************/
1330924e98cSBarry Smith PetscErrorCode ivec_and3( int *arg1,  int *arg2,  int *arg3, int n)
134827bd09bSSatish Balay {
1353fdc5746SBarry Smith   PetscFunctionBegin;
136827bd09bSSatish Balay   while (n--)  {*arg1++ = (*arg2++ & *arg3++);}
1373fdc5746SBarry Smith   PetscFunctionReturn(0);
138827bd09bSSatish Balay }
139827bd09bSSatish Balay 
140*7b1ae94cSBarry Smith /***********************************ivec.c*************************************/
141a501084fSBarry Smith int ivec_sum( int *arg1,  int n)
142827bd09bSSatish Balay {
143a501084fSBarry Smith    int tmp = 0;
144827bd09bSSatish Balay 
145827bd09bSSatish Balay 
146827bd09bSSatish Balay   while (n--) {tmp += *arg1++;}
147827bd09bSSatish Balay   return(tmp);
148827bd09bSSatish Balay }
149827bd09bSSatish Balay 
150*7b1ae94cSBarry Smith /***********************************ivec.c*************************************/
1510924e98cSBarry Smith PetscErrorCode ivec_non_uniform(int *arg1, int *arg2,  int n,  int *arg3)
152827bd09bSSatish Balay {
153a501084fSBarry Smith    int i, j, type;
154827bd09bSSatish Balay 
155827bd09bSSatish Balay 
1563fdc5746SBarry Smith   PetscFunctionBegin;
157827bd09bSSatish Balay   /* LATER: if we're really motivated we can sort and then unsort */
158827bd09bSSatish Balay   for (i=0;i<n;)
159827bd09bSSatish Balay     {
160827bd09bSSatish Balay       /* clump 'em for now */
161827bd09bSSatish Balay       j=i+1;
162827bd09bSSatish Balay       type = arg3[i];
163827bd09bSSatish Balay       while ((j<n)&&(arg3[j]==type))
164827bd09bSSatish Balay 	{j++;}
165827bd09bSSatish Balay 
166827bd09bSSatish Balay       /* how many together */
167827bd09bSSatish Balay       j -= i;
168827bd09bSSatish Balay 
169827bd09bSSatish Balay       /* call appropriate ivec function */
170827bd09bSSatish Balay       if (type == GL_MAX)
171827bd09bSSatish Balay 	{ivec_max(arg1,arg2,j);}
172827bd09bSSatish Balay       else if (type == GL_MIN)
173827bd09bSSatish Balay 	{ivec_min(arg1,arg2,j);}
174827bd09bSSatish Balay       else if (type == GL_MULT)
175827bd09bSSatish Balay 	{ivec_mult(arg1,arg2,j);}
176827bd09bSSatish Balay       else if (type == GL_ADD)
177827bd09bSSatish Balay 	{ivec_add(arg1,arg2,j);}
178827bd09bSSatish Balay       else if (type == GL_B_XOR)
179827bd09bSSatish Balay 	{ivec_xor(arg1,arg2,j);}
180827bd09bSSatish Balay       else if (type == GL_B_OR)
181827bd09bSSatish Balay 	{ivec_or(arg1,arg2,j);}
182827bd09bSSatish Balay       else if (type == GL_B_AND)
183827bd09bSSatish Balay 	{ivec_and(arg1,arg2,j);}
184827bd09bSSatish Balay       else if (type == GL_L_XOR)
185827bd09bSSatish Balay 	{ivec_lxor(arg1,arg2,j);}
186827bd09bSSatish Balay       else if (type == GL_L_OR)
187827bd09bSSatish Balay 	{ivec_lor(arg1,arg2,j);}
188827bd09bSSatish Balay       else if (type == GL_L_AND)
189827bd09bSSatish Balay 	{ivec_land(arg1,arg2,j);}
190827bd09bSSatish Balay       else
191827bd09bSSatish Balay 	{error_msg_fatal("unrecognized type passed to ivec_non_uniform()!");}
192827bd09bSSatish Balay 
193827bd09bSSatish Balay       arg1+=j; arg2+=j; i+=j;
194827bd09bSSatish Balay     }
1953fdc5746SBarry Smith   PetscFunctionReturn(0);
196827bd09bSSatish Balay }
197827bd09bSSatish Balay 
198*7b1ae94cSBarry Smith /***********************************ivec.c*************************************/
199a501084fSBarry Smith vfp ivec_fct_addr( int type)
200827bd09bSSatish Balay {
2013fdc5746SBarry Smith   PetscFunctionBegin;
202827bd09bSSatish Balay   if (type == NON_UNIFORM)
2033fdc5746SBarry Smith     {return((PetscErrorCode (*)(void*, void *, int, ...))&ivec_non_uniform);}
204827bd09bSSatish Balay   else if (type == GL_MAX)
2053fdc5746SBarry Smith     {return((PetscErrorCode (*)(void*, void *, int, ...))&ivec_max);}
206827bd09bSSatish Balay   else if (type == GL_MIN)
2073fdc5746SBarry Smith     {return((PetscErrorCode (*)(void*, void *, int, ...))&ivec_min);}
208827bd09bSSatish Balay   else if (type == GL_MULT)
2093fdc5746SBarry Smith     {return((PetscErrorCode (*)(void*, void *, int, ...))&ivec_mult);}
210827bd09bSSatish Balay   else if (type == GL_ADD)
2113fdc5746SBarry Smith     {return((PetscErrorCode (*)(void*, void *, int, ...))&ivec_add);}
212827bd09bSSatish Balay   else if (type == GL_B_XOR)
2133fdc5746SBarry Smith     {return((PetscErrorCode (*)(void*, void *, int, ...))&ivec_xor);}
214827bd09bSSatish Balay   else if (type == GL_B_OR)
2153fdc5746SBarry Smith     {return((PetscErrorCode (*)(void*, void *, int, ...))&ivec_or);}
216827bd09bSSatish Balay   else if (type == GL_B_AND)
2173fdc5746SBarry Smith     {return((PetscErrorCode (*)(void*, void *, int, ...))&ivec_and);}
218827bd09bSSatish Balay   else if (type == GL_L_XOR)
2193fdc5746SBarry Smith     {return((PetscErrorCode (*)(void*, void *, int, ...))&ivec_lxor);}
220827bd09bSSatish Balay   else if (type == GL_L_OR)
2213fdc5746SBarry Smith     {return((PetscErrorCode (*)(void*, void *, int, ...))&ivec_lor);}
222827bd09bSSatish Balay   else if (type == GL_L_AND)
2233fdc5746SBarry Smith     {return((PetscErrorCode (*)(void*, void *, int, ...))&ivec_land);}
224827bd09bSSatish Balay 
225827bd09bSSatish Balay   /* catch all ... not good if we get here */
226827bd09bSSatish Balay   return(NULL);
227827bd09bSSatish Balay }
228827bd09bSSatish Balay 
229*7b1ae94cSBarry Smith /******************************************************************************/
2300924e98cSBarry Smith PetscErrorCode ivec_sort( int *ar,  int size)
231827bd09bSSatish Balay {
232a501084fSBarry Smith    int *pi, *pj, temp;
233a501084fSBarry Smith    int **top_a = (int **) offset_stack;
234a501084fSBarry Smith    int *top_s = size_stack, *bottom_s = size_stack;
235827bd09bSSatish Balay 
236827bd09bSSatish Balay 
237827bd09bSSatish Balay   /* we're really interested in the offset of the last element */
238827bd09bSSatish Balay   /* ==> length of the list is now size + 1                    */
239827bd09bSSatish Balay   size--;
240827bd09bSSatish Balay 
241827bd09bSSatish Balay   /* do until we're done ... return when stack is exhausted */
242827bd09bSSatish Balay   for (;;)
243827bd09bSSatish Balay     {
244827bd09bSSatish Balay       /* if list is large enough use quicksort partition exchange code */
245827bd09bSSatish Balay       if (size > SORT_OPT)
246827bd09bSSatish Balay 	{
247827bd09bSSatish Balay 	  /* start up pointer at element 1 and down at size     */
248827bd09bSSatish Balay 	  pi = ar+1;
249827bd09bSSatish Balay 	  pj = ar+size;
250827bd09bSSatish Balay 
251827bd09bSSatish Balay 	  /* find middle element in list and swap w/ element 1 */
252827bd09bSSatish Balay 	  SWAP(*(ar+(size>>1)),*pi)
253827bd09bSSatish Balay 
254827bd09bSSatish Balay 	  /* order element 0,1,size-1 st {M,L,...,U} w/L<=M<=U */
255827bd09bSSatish Balay 	  /* note ==> pivot_value in index 0                   */
256827bd09bSSatish Balay 	  if (*pi > *pj)
257827bd09bSSatish Balay 	    {SWAP(*pi,*pj)}
258827bd09bSSatish Balay 	  if (*ar > *pj)
259827bd09bSSatish Balay 	    {SWAP(*ar,*pj)}
260827bd09bSSatish Balay 	  else if (*pi > *ar)
261827bd09bSSatish Balay 	    {SWAP(*(ar),*(ar+1))}
262827bd09bSSatish Balay 
263827bd09bSSatish Balay 	  /* partition about pivot_value ...  	                    */
264827bd09bSSatish Balay 	  /* note lists of length 2 are not guaranteed to be sorted */
265827bd09bSSatish Balay 	  for(;;)
266827bd09bSSatish Balay 	    {
267827bd09bSSatish Balay 	      /* walk up ... and down ... swap if equal to pivot! */
268827bd09bSSatish Balay 	      do pi++; while (*pi<*ar);
269827bd09bSSatish Balay 	      do pj--; while (*pj>*ar);
270827bd09bSSatish Balay 
271827bd09bSSatish Balay 	      /* if we've crossed we're done */
272827bd09bSSatish Balay 	      if (pj<pi) break;
273827bd09bSSatish Balay 
274827bd09bSSatish Balay 	      /* else swap */
275827bd09bSSatish Balay 	      SWAP(*pi,*pj)
276827bd09bSSatish Balay 	    }
277827bd09bSSatish Balay 
278827bd09bSSatish Balay 	  /* place pivot_value in it's correct location */
279827bd09bSSatish Balay 	  SWAP(*ar,*pj)
280827bd09bSSatish Balay 
281827bd09bSSatish Balay 	  /* test stack_size to see if we've exhausted our stack */
282827bd09bSSatish Balay 	  if (top_s-bottom_s >= SORT_STACK)
283827bd09bSSatish Balay 	    {error_msg_fatal("ivec_sort() :: STACK EXHAUSTED!!!");}
284827bd09bSSatish Balay 
285827bd09bSSatish Balay 	  /* push right hand child iff length > 1 */
286827bd09bSSatish Balay 	  if ((*top_s = size-((int) (pi-ar))))
287827bd09bSSatish Balay 	    {
288827bd09bSSatish Balay 	      *(top_a++) = pi;
289827bd09bSSatish Balay 	      size -= *top_s+2;
290827bd09bSSatish Balay 	      top_s++;
291827bd09bSSatish Balay 	    }
292827bd09bSSatish Balay 	  /* set up for next loop iff there is something to do */
293827bd09bSSatish Balay 	  else if (size -= *top_s+2)
294827bd09bSSatish Balay 	    {;}
295827bd09bSSatish Balay 	  /* might as well pop - note NR_OPT >=2 ==> we're ok! */
296827bd09bSSatish Balay 	  else
297827bd09bSSatish Balay 	    {
298827bd09bSSatish Balay 	      ar = *(--top_a);
299827bd09bSSatish Balay 	      size = *(--top_s);
300827bd09bSSatish Balay 	    }
301827bd09bSSatish Balay 	}
302827bd09bSSatish Balay 
303827bd09bSSatish Balay       /* else sort small list directly then pop another off stack */
304827bd09bSSatish Balay       else
305827bd09bSSatish Balay 	{
306827bd09bSSatish Balay 	  /* insertion sort for bottom */
307827bd09bSSatish Balay           for (pj=ar+1;pj<=ar+size;pj++)
308827bd09bSSatish Balay             {
309827bd09bSSatish Balay               temp = *pj;
310827bd09bSSatish Balay               for (pi=pj-1;pi>=ar;pi--)
311827bd09bSSatish Balay                 {
312827bd09bSSatish Balay                   if (*pi <= temp) break;
313827bd09bSSatish Balay                   *(pi+1)=*pi;
314827bd09bSSatish Balay                 }
315827bd09bSSatish Balay               *(pi+1)=temp;
316827bd09bSSatish Balay 	    }
317827bd09bSSatish Balay 
318827bd09bSSatish Balay 	  /* check to see if stack is exhausted ==> DONE */
3193fdc5746SBarry Smith 	  if (top_s==bottom_s)   PetscFunctionReturn(0);
320827bd09bSSatish Balay 
321827bd09bSSatish Balay 	  /* else pop another list from the stack */
322827bd09bSSatish Balay 	  ar = *(--top_a);
323827bd09bSSatish Balay 	  size = *(--top_s);
324827bd09bSSatish Balay 	}
325827bd09bSSatish Balay     }
3263fdc5746SBarry Smith   PetscFunctionReturn(0);
327827bd09bSSatish Balay }
328827bd09bSSatish Balay 
329*7b1ae94cSBarry Smith /******************************************************************************/
3300924e98cSBarry Smith PetscErrorCode ivec_sort_companion( int *ar,  int *ar2,  int size)
331827bd09bSSatish Balay {
332a501084fSBarry Smith    int *pi, *pj, temp, temp2;
333a501084fSBarry Smith    int **top_a = (int **)offset_stack;
334a501084fSBarry Smith    int *top_s = size_stack, *bottom_s = size_stack;
335a501084fSBarry Smith    int *pi2, *pj2;
336a501084fSBarry Smith    int mid;
337827bd09bSSatish Balay 
3383fdc5746SBarry Smith   PetscFunctionBegin;
339827bd09bSSatish Balay   /* we're really interested in the offset of the last element */
340827bd09bSSatish Balay   /* ==> length of the list is now size + 1                    */
341827bd09bSSatish Balay   size--;
342827bd09bSSatish Balay 
343827bd09bSSatish Balay   /* do until we're done ... return when stack is exhausted */
344827bd09bSSatish Balay   for (;;)
345827bd09bSSatish Balay     {
346827bd09bSSatish Balay       /* if list is large enough use quicksort partition exchange code */
347827bd09bSSatish Balay       if (size > SORT_OPT)
348827bd09bSSatish Balay 	{
349827bd09bSSatish Balay 	  /* start up pointer at element 1 and down at size     */
350827bd09bSSatish Balay 	  mid = size>>1;
351827bd09bSSatish Balay 	  pi = ar+1;
352827bd09bSSatish Balay 	  pj = ar+mid;
353827bd09bSSatish Balay 	  pi2 = ar2+1;
354827bd09bSSatish Balay 	  pj2 = ar2+mid;
355827bd09bSSatish Balay 
356827bd09bSSatish Balay 	  /* find middle element in list and swap w/ element 1 */
357827bd09bSSatish Balay 	  SWAP(*pi,*pj)
358827bd09bSSatish Balay 	  SWAP(*pi2,*pj2)
359827bd09bSSatish Balay 
360827bd09bSSatish Balay 	  /* order element 0,1,size-1 st {M,L,...,U} w/L<=M<=U */
361827bd09bSSatish Balay 	  /* note ==> pivot_value in index 0                   */
362827bd09bSSatish Balay 	  pj = ar+size;
363827bd09bSSatish Balay 	  pj2 = ar2+size;
364827bd09bSSatish Balay 	  if (*pi > *pj)
365827bd09bSSatish Balay 	    {SWAP(*pi,*pj) SWAP(*pi2,*pj2)}
366827bd09bSSatish Balay 	  if (*ar > *pj)
367827bd09bSSatish Balay 	    {SWAP(*ar,*pj) SWAP(*ar2,*pj2)}
368827bd09bSSatish Balay 	  else if (*pi > *ar)
369827bd09bSSatish Balay 	    {SWAP(*(ar),*(ar+1)) SWAP(*(ar2),*(ar2+1))}
370827bd09bSSatish Balay 
371827bd09bSSatish Balay 	  /* partition about pivot_value ...  	                    */
372827bd09bSSatish Balay 	  /* note lists of length 2 are not guaranteed to be sorted */
373827bd09bSSatish Balay 	  for(;;)
374827bd09bSSatish Balay 	    {
375827bd09bSSatish Balay 	      /* walk up ... and down ... swap if equal to pivot! */
376827bd09bSSatish Balay 	      do {pi++; pi2++;} while (*pi<*ar);
377827bd09bSSatish Balay 	      do {pj--; pj2--;} while (*pj>*ar);
378827bd09bSSatish Balay 
379827bd09bSSatish Balay 	      /* if we've crossed we're done */
380827bd09bSSatish Balay 	      if (pj<pi) break;
381827bd09bSSatish Balay 
382827bd09bSSatish Balay 	      /* else swap */
383827bd09bSSatish Balay 	      SWAP(*pi,*pj)
384827bd09bSSatish Balay 	      SWAP(*pi2,*pj2)
385827bd09bSSatish Balay 	    }
386827bd09bSSatish Balay 
387827bd09bSSatish Balay 	  /* place pivot_value in it's correct location */
388827bd09bSSatish Balay 	  SWAP(*ar,*pj)
389827bd09bSSatish Balay 	  SWAP(*ar2,*pj2)
390827bd09bSSatish Balay 
391827bd09bSSatish Balay 	  /* test stack_size to see if we've exhausted our stack */
392827bd09bSSatish Balay 	  if (top_s-bottom_s >= SORT_STACK)
393827bd09bSSatish Balay 	    {error_msg_fatal("ivec_sort_companion() :: STACK EXHAUSTED!!!");}
394827bd09bSSatish Balay 
395827bd09bSSatish Balay 	  /* push right hand child iff length > 1 */
396827bd09bSSatish Balay 	  if ((*top_s = size-((int) (pi-ar))))
397827bd09bSSatish Balay 	    {
398827bd09bSSatish Balay 	      *(top_a++) = pi;
399827bd09bSSatish Balay 	      *(top_a++) = pi2;
400827bd09bSSatish Balay 	      size -= *top_s+2;
401827bd09bSSatish Balay 	      top_s++;
402827bd09bSSatish Balay 	    }
403827bd09bSSatish Balay 	  /* set up for next loop iff there is something to do */
404827bd09bSSatish Balay 	  else if (size -= *top_s+2)
405827bd09bSSatish Balay 	    {;}
406827bd09bSSatish Balay 	  /* might as well pop - note NR_OPT >=2 ==> we're ok! */
407827bd09bSSatish Balay 	  else
408827bd09bSSatish Balay 	    {
409827bd09bSSatish Balay 	      ar2 = *(--top_a);
410827bd09bSSatish Balay 	      ar  = *(--top_a);
411827bd09bSSatish Balay 	      size = *(--top_s);
412827bd09bSSatish Balay 	    }
413827bd09bSSatish Balay 	}
414827bd09bSSatish Balay 
415827bd09bSSatish Balay       /* else sort small list directly then pop another off stack */
416827bd09bSSatish Balay       else
417827bd09bSSatish Balay 	{
418827bd09bSSatish Balay 	  /* insertion sort for bottom */
419827bd09bSSatish Balay           for (pj=ar+1, pj2=ar2+1;pj<=ar+size;pj++,pj2++)
420827bd09bSSatish Balay             {
421827bd09bSSatish Balay               temp = *pj;
422827bd09bSSatish Balay               temp2 = *pj2;
423827bd09bSSatish Balay               for (pi=pj-1,pi2=pj2-1;pi>=ar;pi--,pi2--)
424827bd09bSSatish Balay                 {
425827bd09bSSatish Balay                   if (*pi <= temp) break;
426827bd09bSSatish Balay                   *(pi+1)=*pi;
427827bd09bSSatish Balay                   *(pi2+1)=*pi2;
428827bd09bSSatish Balay                 }
429827bd09bSSatish Balay               *(pi+1)=temp;
430827bd09bSSatish Balay               *(pi2+1)=temp2;
431827bd09bSSatish Balay 	    }
432827bd09bSSatish Balay 
433827bd09bSSatish Balay 	  /* check to see if stack is exhausted ==> DONE */
4343fdc5746SBarry Smith 	  if (top_s==bottom_s)   PetscFunctionReturn(0);
435827bd09bSSatish Balay 
436827bd09bSSatish Balay 	  /* else pop another list from the stack */
437827bd09bSSatish Balay 	  ar2 = *(--top_a);
438827bd09bSSatish Balay 	  ar  = *(--top_a);
439827bd09bSSatish Balay 	  size = *(--top_s);
440827bd09bSSatish Balay 	}
441827bd09bSSatish Balay     }
4423fdc5746SBarry Smith   PetscFunctionReturn(0);
443827bd09bSSatish Balay }
444827bd09bSSatish Balay 
445*7b1ae94cSBarry Smith /******************************************************************************/
4461d7d0905SBarry Smith PetscErrorCode ivec_sort_companion_hack( int *ar,  int **ar2, int size)
447827bd09bSSatish Balay {
448a501084fSBarry Smith    int *pi, *pj, temp, *ptr;
449a501084fSBarry Smith    int **top_a = (int **)offset_stack;
450a501084fSBarry Smith    int *top_s = size_stack, *bottom_s = size_stack;
451a501084fSBarry Smith    int **pi2, **pj2;
452a501084fSBarry Smith    int mid;
453827bd09bSSatish Balay 
4543fdc5746SBarry Smith   PetscFunctionBegin;
455827bd09bSSatish Balay   /* we're really interested in the offset of the last element */
456827bd09bSSatish Balay   /* ==> length of the list is now size + 1                    */
457827bd09bSSatish Balay   size--;
458827bd09bSSatish Balay 
459827bd09bSSatish Balay   /* do until we're done ... return when stack is exhausted */
460827bd09bSSatish Balay   for (;;)
461827bd09bSSatish Balay     {
462827bd09bSSatish Balay       /* if list is large enough use quicksort partition exchange code */
463827bd09bSSatish Balay       if (size > SORT_OPT)
464827bd09bSSatish Balay 	{
465827bd09bSSatish Balay 	  /* start up pointer at element 1 and down at size     */
466827bd09bSSatish Balay 	  mid = size>>1;
467827bd09bSSatish Balay 	  pi = ar+1;
468827bd09bSSatish Balay 	  pj = ar+mid;
469827bd09bSSatish Balay 	  pi2 = ar2+1;
470827bd09bSSatish Balay 	  pj2 = ar2+mid;
471827bd09bSSatish Balay 
472827bd09bSSatish Balay 	  /* find middle element in list and swap w/ element 1 */
473827bd09bSSatish Balay 	  SWAP(*pi,*pj)
474827bd09bSSatish Balay 	  P_SWAP(*pi2,*pj2)
475827bd09bSSatish Balay 
476827bd09bSSatish Balay 	  /* order element 0,1,size-1 st {M,L,...,U} w/L<=M<=U */
477827bd09bSSatish Balay 	  /* note ==> pivot_value in index 0                   */
478827bd09bSSatish Balay 	  pj = ar+size;
479827bd09bSSatish Balay 	  pj2 = ar2+size;
480827bd09bSSatish Balay 	  if (*pi > *pj)
481827bd09bSSatish Balay 	    {SWAP(*pi,*pj) P_SWAP(*pi2,*pj2)}
482827bd09bSSatish Balay 	  if (*ar > *pj)
483827bd09bSSatish Balay 	    {SWAP(*ar,*pj) P_SWAP(*ar2,*pj2)}
484827bd09bSSatish Balay 	  else if (*pi > *ar)
485827bd09bSSatish Balay 	    {SWAP(*(ar),*(ar+1)) P_SWAP(*(ar2),*(ar2+1))}
486827bd09bSSatish Balay 
487827bd09bSSatish Balay 	  /* partition about pivot_value ...  	                    */
488827bd09bSSatish Balay 	  /* note lists of length 2 are not guaranteed to be sorted */
489827bd09bSSatish Balay 	  for(;;)
490827bd09bSSatish Balay 	    {
491827bd09bSSatish Balay 	      /* walk up ... and down ... swap if equal to pivot! */
492827bd09bSSatish Balay 	      do {pi++; pi2++;} while (*pi<*ar);
493827bd09bSSatish Balay 	      do {pj--; pj2--;} while (*pj>*ar);
494827bd09bSSatish Balay 
495827bd09bSSatish Balay 	      /* if we've crossed we're done */
496827bd09bSSatish Balay 	      if (pj<pi) break;
497827bd09bSSatish Balay 
498827bd09bSSatish Balay 	      /* else swap */
499827bd09bSSatish Balay 	      SWAP(*pi,*pj)
500827bd09bSSatish Balay 	      P_SWAP(*pi2,*pj2)
501827bd09bSSatish Balay 	    }
502827bd09bSSatish Balay 
503827bd09bSSatish Balay 	  /* place pivot_value in it's correct location */
504827bd09bSSatish Balay 	  SWAP(*ar,*pj)
505827bd09bSSatish Balay 	  P_SWAP(*ar2,*pj2)
506827bd09bSSatish Balay 
507827bd09bSSatish Balay 	  /* test stack_size to see if we've exhausted our stack */
508827bd09bSSatish Balay 	  if (top_s-bottom_s >= SORT_STACK)
509827bd09bSSatish Balay          {error_msg_fatal("ivec_sort_companion_hack() :: STACK EXHAUSTED!!!");}
510827bd09bSSatish Balay 
511827bd09bSSatish Balay 	  /* push right hand child iff length > 1 */
512827bd09bSSatish Balay 	  if ((*top_s = size-((int) (pi-ar))))
513827bd09bSSatish Balay 	    {
514827bd09bSSatish Balay 	      *(top_a++) = pi;
515827bd09bSSatish Balay 	      *(top_a++) = (int*) pi2;
516827bd09bSSatish Balay 	      size -= *top_s+2;
517827bd09bSSatish Balay 	      top_s++;
518827bd09bSSatish Balay 	    }
519827bd09bSSatish Balay 	  /* set up for next loop iff there is something to do */
520827bd09bSSatish Balay 	  else if (size -= *top_s+2)
521827bd09bSSatish Balay 	    {;}
522827bd09bSSatish Balay 	  /* might as well pop - note NR_OPT >=2 ==> we're ok! */
523827bd09bSSatish Balay 	  else
524827bd09bSSatish Balay 	    {
525827bd09bSSatish Balay 	      ar2 = (int **) *(--top_a);
526827bd09bSSatish Balay 	      ar  = *(--top_a);
527827bd09bSSatish Balay 	      size = *(--top_s);
528827bd09bSSatish Balay 	    }
529827bd09bSSatish Balay 	}
530827bd09bSSatish Balay 
531827bd09bSSatish Balay       /* else sort small list directly then pop another off stack */
532827bd09bSSatish Balay       else
533827bd09bSSatish Balay 	{
534827bd09bSSatish Balay 	  /* insertion sort for bottom */
535827bd09bSSatish Balay           for (pj=ar+1, pj2=ar2+1;pj<=ar+size;pj++,pj2++)
536827bd09bSSatish Balay             {
537827bd09bSSatish Balay               temp = *pj;
538827bd09bSSatish Balay               ptr = *pj2;
539827bd09bSSatish Balay               for (pi=pj-1,pi2=pj2-1;pi>=ar;pi--,pi2--)
540827bd09bSSatish Balay                 {
541827bd09bSSatish Balay                   if (*pi <= temp) break;
542827bd09bSSatish Balay                   *(pi+1)=*pi;
543827bd09bSSatish Balay                   *(pi2+1)=*pi2;
544827bd09bSSatish Balay                 }
545827bd09bSSatish Balay               *(pi+1)=temp;
546827bd09bSSatish Balay               *(pi2+1)=ptr;
547827bd09bSSatish Balay 	    }
548827bd09bSSatish Balay 
549827bd09bSSatish Balay 	  /* check to see if stack is exhausted ==> DONE */
5503fdc5746SBarry Smith 	  if (top_s==bottom_s)   PetscFunctionReturn(0);
551827bd09bSSatish Balay 
552827bd09bSSatish Balay 	  /* else pop another list from the stack */
553827bd09bSSatish Balay 	  ar2 = (int **)*(--top_a);
554827bd09bSSatish Balay 	  ar  = *(--top_a);
555827bd09bSSatish Balay 	  size = *(--top_s);
556827bd09bSSatish Balay 	}
557827bd09bSSatish Balay     }
5583fdc5746SBarry Smith   PetscFunctionReturn(0);
559827bd09bSSatish Balay }
560827bd09bSSatish Balay 
561*7b1ae94cSBarry Smith /******************************************************************************/
5620924e98cSBarry Smith PetscErrorCode SMI_sort(void *ar1, void *ar2, int size, int type)
563827bd09bSSatish Balay {
5643fdc5746SBarry Smith   PetscFunctionBegin;
565827bd09bSSatish Balay   if (type == SORT_INTEGER)
566827bd09bSSatish Balay     {
567827bd09bSSatish Balay       if (ar2)
568827bd09bSSatish Balay 	{ivec_sort_companion((int*)ar1,(int*)ar2,size);}
569827bd09bSSatish Balay       else
570827bd09bSSatish Balay 	{ivec_sort((int*)ar1,size);}
571827bd09bSSatish Balay     }
572827bd09bSSatish Balay   else if (type == SORT_INT_PTR)
573827bd09bSSatish Balay     {
574827bd09bSSatish Balay       if (ar2)
575827bd09bSSatish Balay 	{ivec_sort_companion_hack((int*)ar1,(int **)ar2,size);}
576827bd09bSSatish Balay       else
577827bd09bSSatish Balay 	{ivec_sort((int*)ar1,size);}
578827bd09bSSatish Balay     }
579827bd09bSSatish Balay 
580827bd09bSSatish Balay   else
581827bd09bSSatish Balay     {
582827bd09bSSatish Balay       error_msg_fatal("SMI_sort only does SORT_INTEGER!");
583827bd09bSSatish Balay     }
5843fdc5746SBarry Smith   PetscFunctionReturn(0);
585827bd09bSSatish Balay }
586827bd09bSSatish Balay 
587*7b1ae94cSBarry Smith /***********************************ivec.c*************************************/
588*7b1ae94cSBarry Smith int ivec_linear_search( int item,  int *list,  int n)
589827bd09bSSatish Balay {
590a501084fSBarry Smith    int tmp = n-1;
5913fdc5746SBarry Smith   PetscFunctionBegin;
592827bd09bSSatish Balay   while (n--)  {if (*list++ == item) {return(tmp-n);}}
593827bd09bSSatish Balay   return(-1);
594827bd09bSSatish Balay }
595827bd09bSSatish Balay 
596*7b1ae94cSBarry Smith /***********************************ivec.c*************************************/
597*7b1ae94cSBarry Smith int ivec_binary_search( int item,  int *list,  int rh)
598827bd09bSSatish Balay {
599a501084fSBarry Smith    int mid, lh=0;
600827bd09bSSatish Balay 
601827bd09bSSatish Balay   rh--;
602827bd09bSSatish Balay   while (lh<=rh)
603827bd09bSSatish Balay     {
604827bd09bSSatish Balay       mid = (lh+rh)>>1;
605827bd09bSSatish Balay       if (*(list+mid) == item)
606827bd09bSSatish Balay 	{return(mid);}
607827bd09bSSatish Balay       if (*(list+mid) > item)
608827bd09bSSatish Balay 	{rh = mid-1;}
609827bd09bSSatish Balay       else
610827bd09bSSatish Balay 	{lh = mid+1;}
611827bd09bSSatish Balay     }
612827bd09bSSatish Balay   return(-1);
613827bd09bSSatish Balay }
614827bd09bSSatish Balay 
615*7b1ae94cSBarry Smith /*********************************ivec.c*************************************/
6160924e98cSBarry Smith PetscErrorCode rvec_copy( PetscScalar *arg1,  PetscScalar *arg2,  int n)
617827bd09bSSatish Balay {
6183fdc5746SBarry Smith   PetscFunctionBegin;
619827bd09bSSatish Balay   while (n--)  {*arg1++ = *arg2++;}
6203fdc5746SBarry Smith   PetscFunctionReturn(0);
621827bd09bSSatish Balay }
622827bd09bSSatish Balay 
623*7b1ae94cSBarry Smith /*********************************ivec.c*************************************/
6240924e98cSBarry Smith PetscErrorCode rvec_zero( PetscScalar *arg1,  int n)
625827bd09bSSatish Balay {
6263fdc5746SBarry Smith   PetscFunctionBegin;
627827bd09bSSatish Balay   while (n--)  {*arg1++ = 0.0;}
6283fdc5746SBarry Smith   PetscFunctionReturn(0);
629827bd09bSSatish Balay }
630827bd09bSSatish Balay 
631*7b1ae94cSBarry Smith /***********************************ivec.c*************************************/
6320924e98cSBarry Smith PetscErrorCode rvec_one( PetscScalar *arg1,  int n)
633827bd09bSSatish Balay {
6343fdc5746SBarry Smith   PetscFunctionBegin;
635827bd09bSSatish Balay   while (n--)  {*arg1++ = 1.0;}
6363fdc5746SBarry Smith   PetscFunctionReturn(0);
637827bd09bSSatish Balay }
638827bd09bSSatish Balay 
639*7b1ae94cSBarry Smith /***********************************ivec.c*************************************/
6400924e98cSBarry Smith PetscErrorCode rvec_set( PetscScalar *arg1,  PetscScalar arg2,  int n)
641827bd09bSSatish Balay {
6423fdc5746SBarry Smith   PetscFunctionBegin;
643827bd09bSSatish Balay   while (n--)  {*arg1++ = arg2;}
6443fdc5746SBarry Smith   PetscFunctionReturn(0);
645827bd09bSSatish Balay }
646827bd09bSSatish Balay 
647*7b1ae94cSBarry Smith /***********************************ivec.c*************************************/
6480924e98cSBarry Smith PetscErrorCode rvec_scale( PetscScalar *arg1,  PetscScalar arg2,  int n)
649827bd09bSSatish Balay {
6503fdc5746SBarry Smith   PetscFunctionBegin;
651827bd09bSSatish Balay   while (n--)  {*arg1++ *= arg2;}
6523fdc5746SBarry Smith   PetscFunctionReturn(0);
653827bd09bSSatish Balay }
654827bd09bSSatish Balay 
655*7b1ae94cSBarry Smith /*********************************ivec.c*************************************/
6560924e98cSBarry Smith PetscErrorCode rvec_add( PetscScalar *arg1,  PetscScalar *arg2,  int n)
657827bd09bSSatish Balay {
6583fdc5746SBarry Smith   PetscFunctionBegin;
659827bd09bSSatish Balay   while (n--)  {*arg1++ += *arg2++;}
6603fdc5746SBarry Smith   PetscFunctionReturn(0);
661827bd09bSSatish Balay }
662827bd09bSSatish Balay 
663*7b1ae94cSBarry Smith /*********************************ivec.c*************************************/
6640924e98cSBarry Smith PetscErrorCode rvec_mult( PetscScalar *arg1,  PetscScalar *arg2,  int n)
665827bd09bSSatish Balay {
6663fdc5746SBarry Smith   PetscFunctionBegin;
667827bd09bSSatish Balay   while (n--)  {*arg1++ *= *arg2++;}
6683fdc5746SBarry Smith   PetscFunctionReturn(0);
669827bd09bSSatish Balay }
670827bd09bSSatish Balay 
671*7b1ae94cSBarry Smith /*********************************ivec.c*************************************/
6720924e98cSBarry Smith PetscErrorCode rvec_max( PetscScalar *arg1,  PetscScalar *arg2,  int n)
673827bd09bSSatish Balay {
6743fdc5746SBarry Smith   PetscFunctionBegin;
67539945688SSatish Balay   while (n--)  {*arg1 = PetscMax(*arg1,*arg2); arg1++; arg2++;}
6763fdc5746SBarry Smith   PetscFunctionReturn(0);
677827bd09bSSatish Balay }
678827bd09bSSatish Balay 
679*7b1ae94cSBarry Smith /*********************************ivec.c*************************************/
6800924e98cSBarry Smith PetscErrorCode rvec_max_abs( PetscScalar *arg1,  PetscScalar *arg2,  int n)
681827bd09bSSatish Balay {
6823fdc5746SBarry Smith   PetscFunctionBegin;
683827bd09bSSatish Balay   while (n--)  {*arg1 = MAX_FABS(*arg1,*arg2); arg1++; arg2++;}
6843fdc5746SBarry Smith   PetscFunctionReturn(0);
685827bd09bSSatish Balay }
686827bd09bSSatish Balay 
687*7b1ae94cSBarry Smith /*********************************ivec.c*************************************/
6880924e98cSBarry Smith PetscErrorCode rvec_min( PetscScalar *arg1,  PetscScalar *arg2,  int n)
689827bd09bSSatish Balay {
6903fdc5746SBarry Smith   PetscFunctionBegin;
69139945688SSatish Balay   while (n--)  {*arg1 = PetscMin(*arg1,*arg2); arg1++; arg2++;}
6923fdc5746SBarry Smith   PetscFunctionReturn(0);
693827bd09bSSatish Balay }
694827bd09bSSatish Balay 
695*7b1ae94cSBarry Smith /*********************************ivec.c*************************************/
6960924e98cSBarry Smith PetscErrorCode rvec_min_abs( PetscScalar *arg1,  PetscScalar *arg2,  int n)
697827bd09bSSatish Balay {
6983fdc5746SBarry Smith   PetscFunctionBegin;
699827bd09bSSatish Balay   while (n--)  {*arg1 = MIN_FABS(*arg1,*arg2); arg1++; arg2++;}
7003fdc5746SBarry Smith   PetscFunctionReturn(0);
701827bd09bSSatish Balay }
702827bd09bSSatish Balay 
703*7b1ae94cSBarry Smith /*********************************ivec.c*************************************/
7040924e98cSBarry Smith PetscErrorCode rvec_exists( PetscScalar *arg1,  PetscScalar *arg2,  int n)
705827bd09bSSatish Balay {
7063fdc5746SBarry Smith   PetscFunctionBegin;
707827bd09bSSatish Balay   while (n--)  {*arg1 = EXISTS(*arg1,*arg2); arg1++; arg2++;}
7083fdc5746SBarry Smith   PetscFunctionReturn(0);
709827bd09bSSatish Balay }
710827bd09bSSatish Balay 
711*7b1ae94cSBarry Smith /***********************************ivec.c*************************************/
7120924e98cSBarry Smith PetscErrorCode rvec_non_uniform(PetscScalar *arg1, PetscScalar *arg2,  int n,  int *arg3)
713827bd09bSSatish Balay {
714a501084fSBarry Smith    int i, j, type;
715827bd09bSSatish Balay 
7163fdc5746SBarry Smith   PetscFunctionBegin;
717827bd09bSSatish Balay   /* LATER: if we're really motivated we can sort and then unsort */
718827bd09bSSatish Balay   for (i=0;i<n;)
719827bd09bSSatish Balay     {
720827bd09bSSatish Balay       /* clump 'em for now */
721827bd09bSSatish Balay       j=i+1;
722827bd09bSSatish Balay       type = arg3[i];
723827bd09bSSatish Balay       while ((j<n)&&(arg3[j]==type))
724827bd09bSSatish Balay 	{j++;}
725827bd09bSSatish Balay 
726827bd09bSSatish Balay       /* how many together */
727827bd09bSSatish Balay       j -= i;
728827bd09bSSatish Balay 
729827bd09bSSatish Balay       /* call appropriate ivec function */
730827bd09bSSatish Balay       if (type == GL_MAX)
731827bd09bSSatish Balay 	{rvec_max(arg1,arg2,j);}
732827bd09bSSatish Balay       else if (type == GL_MIN)
733827bd09bSSatish Balay 	{rvec_min(arg1,arg2,j);}
734827bd09bSSatish Balay       else if (type == GL_MULT)
735827bd09bSSatish Balay 	{rvec_mult(arg1,arg2,j);}
736827bd09bSSatish Balay       else if (type == GL_ADD)
737827bd09bSSatish Balay 	{rvec_add(arg1,arg2,j);}
738827bd09bSSatish Balay       else if (type == GL_MAX_ABS)
739827bd09bSSatish Balay 	{rvec_max_abs(arg1,arg2,j);}
740827bd09bSSatish Balay       else if (type == GL_MIN_ABS)
741827bd09bSSatish Balay 	{rvec_min_abs(arg1,arg2,j);}
742827bd09bSSatish Balay       else if (type == GL_EXISTS)
743827bd09bSSatish Balay 	{rvec_exists(arg1,arg2,j);}
744827bd09bSSatish Balay       else
745827bd09bSSatish Balay 	{error_msg_fatal("unrecognized type passed to rvec_non_uniform()!");}
746827bd09bSSatish Balay 
747827bd09bSSatish Balay       arg1+=j; arg2+=j; i+=j;
748827bd09bSSatish Balay     }
7493fdc5746SBarry Smith   PetscFunctionReturn(0);
750827bd09bSSatish Balay }
751827bd09bSSatish Balay 
752*7b1ae94cSBarry Smith /***********************************ivec.c*************************************/
753a501084fSBarry Smith vfp rvec_fct_addr( int type)
754827bd09bSSatish Balay {
755827bd09bSSatish Balay   if (type == NON_UNIFORM)
7563fdc5746SBarry Smith     {return((PetscErrorCode (*)(void*, void *, int, ...))&rvec_non_uniform);}
757827bd09bSSatish Balay   else if (type == GL_MAX)
7583fdc5746SBarry Smith     {return((PetscErrorCode (*)(void*, void *, int, ...))&rvec_max);}
759827bd09bSSatish Balay   else if (type == GL_MIN)
7603fdc5746SBarry Smith     {return((PetscErrorCode (*)(void*, void *, int, ...))&rvec_min);}
761827bd09bSSatish Balay   else if (type == GL_MULT)
7623fdc5746SBarry Smith     {return((PetscErrorCode (*)(void*, void *, int, ...))&rvec_mult);}
763827bd09bSSatish Balay   else if (type == GL_ADD)
7643fdc5746SBarry Smith     {return((PetscErrorCode (*)(void*, void *, int, ...))&rvec_add);}
765827bd09bSSatish Balay   else if (type == GL_MAX_ABS)
7663fdc5746SBarry Smith     {return((PetscErrorCode (*)(void*, void *, int, ...))&rvec_max_abs);}
767827bd09bSSatish Balay   else if (type == GL_MIN_ABS)
7683fdc5746SBarry Smith     {return((PetscErrorCode (*)(void*, void *, int, ...))&rvec_min_abs);}
769827bd09bSSatish Balay   else if (type == GL_EXISTS)
7703fdc5746SBarry Smith     {return((PetscErrorCode (*)(void*, void *, int, ...))&rvec_exists);}
771827bd09bSSatish Balay 
772827bd09bSSatish Balay   /* catch all ... not good if we get here */
773827bd09bSSatish Balay   return(NULL);
774827bd09bSSatish Balay }
775827bd09bSSatish Balay 
776827bd09bSSatish Balay 
777827bd09bSSatish Balay 
778827bd09bSSatish Balay 
779827bd09bSSatish Balay 
780827bd09bSSatish Balay 
781