Actual source code: ex3f90.F90

petsc-3.10.5 2019-03-28
Report Typos and Errors
  1: !
  2: !    Description:  Creates an index set based on blocks of integers. Views that index set
  3: !    and then destroys it.
  4: !
  5: !/*T
  6: !    Concepts: index sets^manipulating a block index set;
  7: !    Concepts: Fortran90^accessing indices in index set;
  8: !
  9: !T*/
 10: !
 11:       program main
 12:  #include <petsc/finclude/petscis.h>
 13:       use petscis
 14:       implicit none

 16:       PetscErrorCode ierr
 17:       PetscInt n,bs,issize
 18:       PetscInt inputindices(4)
 19:       PetscInt, pointer :: indices(:)
 20:       IS       set
 21:       PetscBool  isablock;

 23:       n               = 4
 24:       bs              = 3
 25:       inputindices(1) = 0
 26:       inputindices(2) = 1
 27:       inputindices(3) = 3
 28:       inputindices(4) = 4

 30:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
 31:       if (ierr .ne. 0) then
 32:          print*,'Unable to initialize PETSc'
 33:          stop
 34:        endif
 35: !
 36: !    Create a block index set. The index set has 4 blocks each of size 3.
 37: !    The indices are {0,1,2,3,4,5,9,10,11,12,13,14}
 38: !    Note each processor is generating its own index set
 39: !    (in this case they are all identical)
 40: !
 41:       call ISCreateBlock(PETSC_COMM_SELF,bs,n,inputindices,PETSC_COPY_VALUES,set,ierr);CHKERRA(ierr)
 42:       call ISView(set,PETSC_VIEWER_STDOUT_SELF,ierr);CHKERRA(ierr)

 44: !
 45: !    Extract indices from set.
 46: !
 47:       call ISGetLocalSize(set,issize,ierr);CHKERRA(ierr)
 48:       call ISGetIndicesF90(set,indices,ierr);CHKERRA(ierr)
 49:       write(6,100)indices
 50:  100  format(12I3)
 51:       call ISRestoreIndicesF90(set,indices,ierr);CHKERRA(ierr)

 53: !
 54: !    Extract the block indices. This returns one index per block.
 55: !
 56:       call ISBlockGetIndicesF90(set,indices,ierr);CHKERRA(ierr)
 57:       write(6,200)indices
 58:  200  format(4I3)
 59:       call ISBlockRestoreIndicesF90(set,indices,ierr);CHKERRA(ierr)

 61: !
 62: !    Check if this is really a block index set
 63: !
 64:       call PetscObjectTypeCompare(set,ISBLOCK,isablock,ierr);CHKERRA(ierr)
 65:       if (.not. isablock) then
 66:         write(6,*) 'Index set is not blocked!'
 67:       endif

 69: !
 70: !    Determine the block size of the index set
 71: !
 72:       call ISGetBlockSize(set,bs,ierr);CHKERRA(ierr)
 73:       if (bs .ne. 3) then
 74:         write(6,*) 'Blocksize != 3'
 75:       endif

 77: !
 78: !    Get the number of blocks
 79: !
 80:       call ISBlockGetLocalSize(set,n,ierr);CHKERRA(ierr)
 81:       if (n .ne. 4) then
 82:         write(6,*) 'Number of blocks != 4'
 83:       endif

 85:       call ISDestroy(set,ierr);CHKERRA(ierr)
 86:       call PetscFinalize(ierr)
 87:       end

 89: !/*TEST
 90: !
 91: !   test:
 92: !      filter: sort -b
 93: !      filter_output: sort -b
 94: !
 95: !TEST*/