16dd63270SBarry Smith #include <petsc/private/ftnimpl.h> 24d3610e3SJacob Faibussowitsch #include <petscsys.h> 34d3610e3SJacob Faibussowitsch 44d3610e3SJacob Faibussowitsch #if defined(PETSC_HAVE_FORTRAN_CAPS) 54d3610e3SJacob Faibussowitsch #define petsctimsort_ PETSCTIMSORT 64d3610e3SJacob Faibussowitsch #define petsctimsortwitharray_ PETSCTIMSORTWITHARRAY 74d3610e3SJacob Faibussowitsch #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) 84d3610e3SJacob Faibussowitsch #define petsctimsort_ petsctimsort 94d3610e3SJacob Faibussowitsch #define petsctimsortwitharray_ petsctimsortwitharray 104d3610e3SJacob Faibussowitsch #endif 114d3610e3SJacob Faibussowitsch 124d3610e3SJacob Faibussowitsch struct fc_c { 134d3610e3SJacob Faibussowitsch void (*fcmp)(const void *a, const void *b, void *c, int *res); 144d3610e3SJacob Faibussowitsch void *fctx; 154d3610e3SJacob Faibussowitsch } fc_c; 164d3610e3SJacob Faibussowitsch 17*2a8381b2SBarry Smith int cmp_via_fortran(const void *a, const void *b, PetscCtx ctx) 184d3610e3SJacob Faibussowitsch { 194d3610e3SJacob Faibussowitsch int result; 204d3610e3SJacob Faibussowitsch struct fc_c *fc = (struct fc_c *)ctx; 214d3610e3SJacob Faibussowitsch fc->fcmp(a, b, fc->fctx, &result); 224d3610e3SJacob Faibussowitsch return result; 234d3610e3SJacob Faibussowitsch } 244d3610e3SJacob Faibussowitsch 25*2a8381b2SBarry Smith PETSC_EXTERN void petsctimsort_(PetscInt *n, void *arr, size_t *size, void (*cmp)(const void *, const void *, void *, int *), PetscCtx ctx, PetscErrorCode *ierr) 264d3610e3SJacob Faibussowitsch { 274d3610e3SJacob Faibussowitsch struct fc_c fc = {cmp, ctx}; 284d3610e3SJacob Faibussowitsch *ierr = PetscTimSort(*n, arr, *size, cmp_via_fortran, &fc); 294d3610e3SJacob Faibussowitsch } 304d3610e3SJacob Faibussowitsch 31*2a8381b2SBarry Smith PETSC_EXTERN void petsctimsortwitharray_(PetscInt *n, void *arr, size_t *asize, void *barr, size_t *bsize, void (*cmp)(const void *, const void *, void *, int *), PetscCtx ctx, PetscErrorCode *ierr) 324d3610e3SJacob Faibussowitsch { 334d3610e3SJacob Faibussowitsch struct fc_c fc = {cmp, ctx}; 344d3610e3SJacob Faibussowitsch *ierr = PetscTimSortWithArray(*n, arr, *asize, barr, *bsize, cmp_via_fortran, &fc); 354d3610e3SJacob Faibussowitsch } 36