Actual source code: ex3f90.F90

  1: !
  2: !    Description:  Creates an index set based on blocks of integers. Views that index set
  3: !    and then destroys it.
  4: !
  5: #include <petsc/finclude/petscis.h>
  6: program main
  7:   use petscis
  8:   implicit none

 10:   PetscErrorCode ierr
 11:   PetscInt n, bs, issize
 12:   PetscInt inputindices(4)
 13:   PetscInt, pointer :: indices(:)
 14:   IS set
 15:   PetscBool isablock

 17:   n = 4
 18:   bs = 3
 19:   inputindices = [0, 1, 3, 4]

 21:   PetscCallA(PetscInitialize(ierr))
 22: !
 23: ! Create a block index set. The index set has 4 blocks each of size 3.
 24: ! The indices are {0,1,2,3,4,5,9,10,11,12,13,14}
 25: ! Note each processor is generating its own index set
 26: ! (in this case they are all identical)
 27: !
 28:   PetscCallA(ISCreateBlock(PETSC_COMM_SELF, bs, n, inputindices, PETSC_COPY_VALUES, set, ierr))
 29:   PetscCallA(ISView(set, PETSC_VIEWER_STDOUT_SELF, ierr))

 31: !
 32: ! Extract indices from set.
 33: !
 34:   PetscCallA(ISGetLocalSize(set, issize, ierr))
 35:   PetscCallA(ISGetIndices(set, indices, ierr))
 36:   write (6, 100) indices
 37: 100 format(12I3)
 38:   PetscCallA(ISRestoreIndices(set, indices, ierr))

 40: !
 41: ! Extract the block indices. This returns one index per block.
 42: !
 43:   PetscCallA(ISBlockGetIndices(set, indices, ierr))
 44:   write (6, 200) indices
 45: 200 format(4I3)
 46:   PetscCallA(ISBlockRestoreIndices(set, indices, ierr))

 48: !
 49: ! Check if this is really a block index set
 50: !
 51:   PetscCallA(PetscObjectTypeCompare(PetscObjectCast(set), ISBLOCK, isablock, ierr))
 52:   if (.not. isablock) then
 53:     write (6, *) 'Index set is not blocked!'
 54:   end if

 56: !
 57: ! Determine the block size of the index set
 58: !
 59:   PetscCallA(ISGetBlockSize(set, bs, ierr))
 60:   if (bs /= 3) then
 61:     write (6, *) 'Blocksize != 3'
 62:   end if

 64: !
 65: ! Get the number of blocks
 66: !
 67:   PetscCallA(ISBlockGetLocalSize(set, n, ierr))
 68:   if (n /= 4) then
 69:     write (6, *) 'Number of blocks != 4'
 70:   end if

 72:   PetscCallA(ISDestroy(set, ierr))
 73:   PetscCallA(PetscFinalize(ierr))
 74: end

 76: !/*TEST
 77: !
 78: !   test:
 79: !      filter: sort -b
 80: !      filter_output: sort -b
 81: !
 82: !TEST*/