#define PETSC_DLL /* This file contains routines for sorting doubles. Values are sorted in place. These are provided because the general sort routines incur a great deal of overhead in calling the comparision routines. */ #include "petsc.h" /*I "petsc.h" I*/ #include "petscsys.h" /*I "petscsys.h" I*/ #define SWAP(a,b,t) {t=a;a=b;b=t;} #undef __FUNCT__ #define __FUNCT__ "PetscSortReal_Private" /* A simple version of quicksort; taken from Kernighan and Ritchie, page 87 */ static PetscErrorCode PetscSortReal_Private(PetscReal *v,PetscInt right) { PetscInt i,last; PetscReal vl,tmp; PetscFunctionBegin; if (right <= 1) { if (right == 1) { if (v[0] > v[1]) SWAP(v[0],v[1],tmp); } PetscFunctionReturn(0); } SWAP(v[0],v[right/2],tmp); vl = v[0]; last = 0; for (i=1; i<=right; i++) { if (v[i] < vl) {last++; SWAP(v[last],v[i],tmp);} } SWAP(v[0],v[last],tmp); PetscSortReal_Private(v,last-1); PetscSortReal_Private(v+last+1,right-(last+1)); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "PetscSortReal" /*@ PetscSortReal - Sorts an array of doubles in place in increasing order. Not Collective Input Parameters: + n - number of values - v - array of doubles Level: intermediate Concepts: sorting^doubles .seealso: PetscSortInt(), PetscSortRealWithPermutation() @*/ PetscErrorCode PETSC_DLLEXPORT PetscSortReal(PetscInt n,PetscReal v[]) { PetscInt j,k; PetscReal tmp,vk; PetscFunctionBegin; if (n<8) { for (k=0; k v[j]) { SWAP(v[k],v[j],tmp); vk = v[k]; } } } } else { PetscSortReal_Private(v,n-1); } PetscFunctionReturn(0); }