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:       program main
  6: #include <petsc/finclude/petscis.h>
  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(1) = 0
 20:       inputindices(2) = 1
 21:       inputindices(3) = 3
 22:       inputindices(4) = 4

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

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

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

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

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

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

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

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