:orphan: # PetscTimSort Sorts an array in place in increasing order using Tim Peters adaptive sorting algorithm. ## Synopsis ``` #include "petscsys.h" PetscErrorCode PetscTimSort(PetscInt n, void *arr, size_t size, int (*cmp)(const void *, const void *, void *), void *ctx) ``` Not Collective ## Input Parameters - ***n -*** number of values - ***arr -*** array to be sorted - ***size -*** size in bytes of the datatype held in arr - ***cmp -*** function pointer to comparison function - ***ctx -*** optional context to be passed to comparison function, NULL if not needed ## Output Parameter - ***arr -*** sorted array ## Notes Timsort makes the assumption that input data is already likely partially ordered, or that it contains contiguous sections (termed 'runs') where the data is locally ordered (but not necessarily globally ordered). It therefore aims to select slices of the array in such a way that resulting mergesorts operate on near perfectly length-balanced arrays. To do so it repeatedly triggers attempts throughout to merge adjacent runs. Should one run continuously "win" a comparison the algorithm begins the "gallop" phase. It will aggressively search the "winner" for the location of the "losers" next entry (and vice versa) to copy all preceding elements into place in bulk. However if the data is truly unordered (as is the case with random data) the immense gains possible from these searches are expected __not__ to repay their costs. While adjacent arrays are almost all nearly the same size, they likely all contain similar data. ## Sample usage The comparison function must follow the `qsort()` comparison function paradigm, returning the sign of the difference between its arguments. If left < right : return -1, if left == right : return 0, if left > right : return 1. The user may also change or reverse the order of the sort by flipping the above. Note that stability of the sort is only guaranteed if the comparison function forms a valid trigraph. For example when sorting an array of type "my_type" in increasing order ```none int my_increasing_comparison_function(const void *left, const void *right, void *ctx) { my_type l = *(my_type *) left, r = *(my_type *) right; return (l < r) ? -1 : (l > r); } ``` Note the context is unused here but you may use it to pass and subsequently access whatever information required inside the comparison function. The context pointer will unaltered except for any changes made inside the comparison function. Then pass the function ```none PetscTimSort(n, arr, sizeof(arr[0]), my_increasing_comparison_function, ctx) ``` ## Fortran Notes To use this from fortran you must write a comparison subroutine with 4 arguments which accepts left, right, ctx and returns result. For example ```none subroutine CompareIntegers(left,right,ctx,result) implicit none PetscInt,intent(in) :: left, right type(UserCtx) :: ctx integer,intent(out) :: result if (left < right) then result = -1 else if (left == right) then result = 0 else result = 1 end if return end subroutine CompareIntegers ``` ## References - **** -*** Tim Peters. https://bugs.python.org/file4451/timsort.txt ## See Also `PetscTimSortWithArray()`, `PetscIntSortSemiOrdered()`, `PetscRealSortSemiOrdered()`, `PetscMPIIntSortSemiOrdered()` ## Level developer ## Location src/sys/utils/sortso.c --- [Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/sys/utils/sortso.c) [Index of all Sys routines](index.md) [Table of Contents for all manual pages](/manualpages/index.md) [Index of all manual pages](/manualpages/singleindex.md)