Actual source code: ex1f90.F

petsc-3.3-p7 2013-05-11
  1: !
  2: !  Description: Creates an index set based on a set of integers. Views that index set
  3: !  and then destroys it.
  4: !
  5: !/*T
  6: !    Concepts: index sets^manipulating a general index set;
  7: !    Concepts: Fortran90^accessing indices of index set;
  8: !T*/
  9: !
 10: !  The following include statements are required for Fortran programs
 11: !  that use PETSc index sets:
 12: !     petscsys.h  - base PETSc routines
 13: !     petscis.h     - index sets (IS objects)
 14: !     petscis.h90   - to allow access to Fortran90 features of index sets
 15: !
 16:       program main
 17:       implicit none

 19: #include <finclude/petscsys.h>
 20: #include <finclude/petscis.h>
 21: #include <finclude/petscis.h90>

 23:       PetscErrorCode ierr
 24:       PetscInt indices(5),n
 25:       PetscInt five
 26:       PetscMPIInt rank
 27:       PetscInt, pointer :: idx(:)
 28:       IS      is

 30:       five = 5
 31:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
 32:       call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)

 34: !  Create an index set with 5 entries. Each processor creates
 35: !  its own index set with its own list of integers.
 36: 
 37:       indices(1) = rank + 1
 38:       indices(2) = rank + 2
 39:       indices(3) = rank + 3
 40:       indices(4) = rank + 4
 41:       indices(5) = rank + 5
 42:       call ISCreateGeneral(PETSC_COMM_SELF,five,indices,                   &
 43:      &                     PETSC_COPY_VALUES,is,ierr)

 45: !  Print the index set to stdout

 47:       call ISView(is,PETSC_VIEWER_STDOUT_SELF,ierr)

 49: !  Get the number of indices in the set

 51:       call ISGetLocalSize(is,n,ierr)

 53: !   Get the indices in the index set

 55:       call ISGetIndicesF90(is,idx,ierr)

 57:       if (associated(idx)) then
 58:          write (*,*) 'Association check passed'
 59:       else
 60:          write (*,*) 'Association check failed'
 61:       endif

 63: !   Now any code that needs access to the list of integers
 64: !   has access to it here

 66:       write(6,50) idx
 67:  50   format(5I3)

 69:       write(6,100) rank,idx(1),idx(5)
 70:  100  format('[',i5,'] First index = ',i5,' fifth index = ',i5)
 71: 
 72: !   Once we no longer need access to the indices they should
 73: !   returned to the system

 75:       call ISRestoreIndicesF90(is,idx,ierr)
 76: 
 77: !   All PETSc objects should be destroyed once they are
 78: !   no longer needed

 80:       call ISDestroy(is,ierr)
 81:       call PetscFinalize(ierr)
 82:       end

 84: