Actual source code: ex49f.F90

  1: !
  2: !  Test Fortran binding of sort routines
  3: !
  4: #include "petsc/finclude/petsc.h"
  5: module ex49fmodule
  6:   use petsc
  7:   implicit none
  8:   type uctx
  9:     PetscInt myint
 10:   end type uctx
 11: contains
 12:   subroutine CompareIntegers(a, b, ctx, res)
 13:     PetscInt, intent(in) :: a, b
 14:     type(uctx) :: ctx
 15:     integer, intent(out)  :: res

 17:     if (a < b) then
 18:       res = -1
 19:     else if (a == b) then
 20:       res = 0
 21:     else
 22:       res = 1
 23:     end if
 24:   end subroutine CompareIntegers
 25: end module ex49fmodule

 27: program main
 28:   use ex49fmodule
 29:   implicit none

 31:   PetscErrorCode ierr
 32:   PetscCount, parameter::  iN = 3
 33:   PetscInt, parameter ::  N = 3
 34:   PetscInt x(N), x1(N), y(N), z(N)
 35:   PetscMPIInt mx(N), my(N)
 36:   PetscScalar s(N)
 37:   PetscReal r(N)
 38:   PetscMPIInt, parameter:: two = 2, five = 5, seven = 7
 39:   type(uctx)::            ctx
 40:   PetscInt i
 41:   PetscSizeT sizeofentry

 43:   PetscCallA(PetscInitialize(ierr))

 45:   x = [3, 2, 1]
 46:   x1 = [3, 2, 1]
 47:   y = [6, 5, 4]
 48:   z = [3, 5, 2]
 49:   mx = [five, seven, two]
 50:   my = [five, seven, two]
 51:   s = [1.0, 2.0, 3.0]
 52:   r = [1.0, 2.0, 3.0]
 53: #if defined(PETSC_USE_64BIT_INDICES)
 54:   sizeofentry = 8
 55: #else
 56:   sizeofentry = 4
 57: #endif
 58:   ctx%myint = 1
 59:   PetscCallA(PetscSortInt(iN, x, ierr))
 60:   PetscCallA(PetscTimSort(N, x1, sizeofentry, CompareIntegers, ctx, ierr))
 61:   do i = 1, N
 62:     PetscCheckA(x1(i) == x(i), PETSC_COMM_SELF, PETSC_ERR_PLIB, 'PetscTimSort and PetscSortInt arrays did not match')
 63:   end do
 64:   PetscCallA(PetscSortIntWithArray(iN, y, x, ierr))
 65:   PetscCallA(PetscSortIntWithArrayPair(iN, x, y, z, ierr))

 67:   PetscCallA(PetscSortMPIInt(iN, mx, ierr))
 68:   PetscCallA(PetscSortMPIIntWithArray(iN, mx, my, ierr))
 69:   PetscCallA(PetscSortMPIIntWithIntArray(iN, mx, y, ierr))

 71:   PetscCallA(PetscSortIntWithScalarArray(iN, x, s, ierr))

 73:   PetscCallA(PetscSortReal(iN, r, ierr))
 74:   PetscCallA(PetscSortRealWithArrayInt(iN, r, x, ierr))

 76:   PetscCallA(PetscFinalize(ierr))
 77: end program main

 79: !/*TEST
 80: !
 81: !   test:
 82: !      output_file: output/empty.out
 83: !
 84: !TEST*/