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: !
  6:       program main
  7: #include <petsc/finclude/petscis.h>
  8:       use petscis
  9:       implicit none

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

 18:       n               = 4
 19:       bs              = 3
 20:       inputindices(1) = 0
 21:       inputindices(2) = 1
 22:       inputindices(3) = 3
 23:       inputindices(4) = 4

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

 35: !
 36: !    Extract indices from set.
 37: !
 38:       PetscCallA(ISGetLocalSize(set,issize,ierr))
 39:       PetscCallA(ISGetIndicesF90(set,indices,ierr))
 40:       write(6,100)indices
 41:  100  format(12I3)
 42:       PetscCallA(ISRestoreIndicesF90(set,indices,ierr))

 44: !
 45: !    Extract the block indices. This returns one index per block.
 46: !
 47:       PetscCallA(ISBlockGetIndicesF90(set,indices,ierr))
 48:       write(6,200)indices
 49:  200  format(4I3)
 50:       PetscCallA(ISBlockRestoreIndicesF90(set,indices,ierr))

 52: !
 53: !    Check if this is really a block index set
 54: !
 55:       PetscCallA(PetscObjectTypeCompare(set,ISBLOCK,isablock,ierr))
 56:       if (.not. isablock) then
 57:         write(6,*) 'Index set is not blocked!'
 58:       endif

 60: !
 61: !    Determine the block size of the index set
 62: !
 63:       PetscCallA(ISGetBlockSize(set,bs,ierr))
 64:       if (bs .ne. 3) then
 65:         write(6,*) 'Blocksize != 3'
 66:       endif

 68: !
 69: !    Get the number of blocks
 70: !
 71:       PetscCallA(ISBlockGetLocalSize(set,n,ierr))
 72:       if (n .ne. 4) then
 73:         write(6,*) 'Number of blocks != 4'
 74:       endif

 76:       PetscCallA(ISDestroy(set,ierr))
 77:       PetscCallA(PetscFinalize(ierr))
 78:       end

 80: !/*TEST
 81: !
 82: !   test:
 83: !      filter: sort -b
 84: !      filter_output: sort -b
 85: !
 86: !TEST*/