Actual source code: ex17f.F90

  1: #include <petsc/finclude/petscmat.h>
  2: program main
  3:   use petscmat

  5:   implicit none

  7:   Mat A
  8:   MatPartitioning part
  9:   IS is
 10:   PetscInt   ::     i, m, N
 11:   PetscInt   ::     rstart, rend
 12:   PetscInt, pointer, dimension(:) ::   emptyranks, bigranks, cols
 13:   PetscScalar, pointer, dimension(:) :: vals
 14:   PetscInt :: &
 15:     nbigranks = 10, &
 16:     nemptyranks = 10
 17:   PetscMPIInt   ::  rank, sizef
 18:   PetscErrorCode ierr
 19:   PetscBool set
 20:   PetscInt, parameter :: zero = 0, one = 1, two = 2, three = 3

 22:   PetscCallA(PetscInitialize(ierr))

 24:   PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD, rank, ierr))
 25:   PetscCallMPIA(MPI_Comm_size(PETSC_COMM_WORLD, sizef, ierr))

 27:   allocate (emptyranks(nemptyranks))
 28:   allocate (bigranks(nbigranks))

 30:   PetscCallA(PetscOptionsGetIntArray(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-emptyranks', emptyranks, nemptyranks, set, ierr))
 31:   PetscCallA(PetscOptionsGetIntArray(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-bigranks', bigranks, nbigranks, set, ierr))

 33:   m = 1
 34:   do i = 1, nemptyranks
 35:     if (rank == emptyranks(i)) m = 0
 36:   end do
 37:   do i = 1, nbigranks
 38:     if (rank == bigranks(i)) m = 5
 39:   end do

 41:   deallocate (emptyranks)
 42:   deallocate (bigranks)

 44:   PetscCallA(MatCreate(PETSC_COMM_WORLD, A, ierr))
 45:   PetscCallA(MatSetsizes(A, m, m, PETSC_DECIDE, PETSC_DECIDE, ierr))
 46:   PetscCallA(MatSetFromOptions(A, ierr))
 47:   PetscCallA(MatSeqAIJSetPreallocation(A, three, PETSC_NULL_INTEGER_ARRAY, ierr))
 48:   PetscCallA(MatMPIAIJSetPreallocation(A, three, PETSC_NULL_INTEGER_ARRAY, two, PETSC_NULL_INTEGER_ARRAY, ierr))
 49:   PetscCallA(MatSeqBAIJSetPreallocation(A, one, three, PETSC_NULL_INTEGER_ARRAY, ierr))
 50:   PetscCallA(MatMPIBAIJSetPreallocation(A, one, three, PETSC_NULL_INTEGER_ARRAY, two, PETSC_NULL_INTEGER_ARRAY, ierr))
 51:   PetscCallA(MatSeqSBAIJSetPreallocation(A, one, two, PETSC_NULL_INTEGER_ARRAY, ierr))
 52:   PetscCallA(MatMPISBAIJSetPreallocation(A, one, two, PETSC_NULL_INTEGER_ARRAY, one, PETSC_NULL_INTEGER_ARRAY, ierr))

 54:   PetscCallA(MatGetSize(A, PETSC_NULL_INTEGER, N, ierr))
 55:   PetscCallA(MatGetOwnershipRange(A, rstart, rend, ierr))

 57:   allocate (cols(1:3))
 58:   allocate (vals(1:3))
 59:   do i = rstart, rend - 1

 61:     cols = (/mod((i + N - 1), N), i, mod((i + 1), N)/)
 62:     vals = [1.0, 1.0, 1.0]
 63:     PetscCallA(MatSetValues(A, one, [i], three, cols, vals, INSERT_VALUES, ierr))
 64:   end do
 65:   deallocate (cols)
 66:   deallocate (vals)
 67:   PetscCallA(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY, ierr))
 68:   PetscCallA(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY, ierr))
 69:   PetscCallA(MatView(A, PETSC_VIEWER_STDOUT_WORLD, ierr))

 71:   PetscCallA(MatPartitioningCreate(PETSC_COMM_WORLD, part, ierr))
 72:   PetscCallA(MatPartitioningSetAdjacency(part, A, ierr))
 73:   PetscCallA(MatPartitioningSetFromOptions(part, ierr))
 74:   PetscCallA(MatPartitioningApply(part, is, ierr))
 75:   PetscCallA(ISView(is, PETSC_VIEWER_STDOUT_WORLD, ierr))
 76:   PetscCallA(ISDestroy(is, ierr))
 77:   PetscCallA(MatPartitioningDestroy(part, ierr))
 78:   PetscCallA(MatDestroy(A, ierr))
 79:   PetscCallA(PetscFinalize(ierr))

 81: end program

 83: !/*TEST
 84: !
 85: !   test:
 86: !      nsize: 8
 87: !      args: -emptyranks 0,2,4 -bigranks 1,3,7 -mat_partitioning_type average
 88: !      output_file: output/ex17_1.out
 89: !      # cannot test with external package partitioners since they produce different results on different systems
 90: !
 91: !TEST*/