Actual source code: ex9f.F

  1: !
  2: !     "$Id: ex9f.F,v 1.10 2001/01/22 23:27:17 balay Exp $";
  3: !
  4: ! Description: Illustrates the use of VecCreateGhost()
  5: !
  6: !/*T
  7: !   Concepts: vectors^assembling vectors; ghost padding
  8: !   Processors: n
  9: !
 10: !   Comment: Ghost padding is one way to handle local calculations that
 11: !      involve values from other processors. VecCreateGhost() provides
 12: !      a way to create vectors with extra room at the end of the vector 
 13: !      array to contain the needed ghost values from other processors, 
 14: !      vector computations are otherwise unaffected.
 15: !T*/

 17:       program main
 18:       implicit none

 20: !
 21: !  The following include statements are required for Fortran programs
 22: !  that use PETSc vectors:
 23: !     petsc.h       - base PETSc routines
 24: !     petscvec.h    - vectors
 25: !

 27: #include "include/finclude/petsc.h"
 28: #include "include/finclude/petscvec.h"

 30:       integer    rank,nlocal,nghost,ifrom(2),size,ierr,i,rstart,rend
 31:       integer    flag
 32:       Scalar     value,tarray(20)
 33:       Vec        lx,gx,gxs

 35:       nlocal = 6
 36:       nghost = 2

 38:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
 39:       call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
 40:       call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)

 42:       if (size .ne. 2) then
 43:        SETERRQ(1,'Must run with two processors',ierr)
 44:       endif

 46: !
 47: !     Construct a two dimensional graph connecting nlocal degrees of
 48: !     freedom per processor. From this we will generate the global
 49: !     indices of needed ghost values
 50: !
 51: !     For simplicity we generate the entire graph on each processor:
 52: !     in real application the graph would stored in parallel, but this
 53: !     example is only to demonstrate the management of ghost padding
 54: !     with VecCreateGhost().
 55: !
 56: !     In this example we consider the vector as representing
 57: !     degrees of freedom in a one dimensional grid with periodic
 58: !     boundary conditions.
 59: !
 60: !        ----Processor  1---------  ----Processor 2 --------
 61: !         0    1   2   3   4    5    6    7   8   9   10   11
 62: !                               |----|
 63: !         |-------------------------------------------------|
 64: !


 67:       if (rank .eq. 0) then
 68:         ifrom(1) = 11
 69:         ifrom(2) = 6
 70:       else
 71:         ifrom(1) = 0
 72:         ifrom(2) = 5
 73:       endif

 75: !     Create the vector with two slots for ghost points. Note that both
 76: !     the local vector (lx) and the global vector (gx) share the same
 77: !     array for storing vector values.

 79:       call PetscOptionsHasName(PETSC_NULL_CHARACTER,'-allocate',         &
 80:      &     flag,ierr)
 81:       if (flag .ne. 0) then
 82:         call VecCreateGhostWithArray(PETSC_COMM_WORLD,nlocal,            &
 83:      &        PETSC_DECIDE,nghost,ifrom,tarray,gxs,ierr)
 84:       else
 85:         call VecCreateGhost(PETSC_COMM_WORLD,nlocal,PETSC_DECIDE,        &
 86:      &       nghost,ifrom,gxs,ierr)
 87:       endif


 90: !      Test VecDuplicate

 92:        call VecDuplicate(gxs,gx,ierr)
 93:        call VecDestroy(gxs,ierr)

 95: !      Access the local Form

 97:        call VecGhostGetLocalForm(gx,lx,ierr)

 99: !     Set the values from 0 to 12 into the "global" vector

101:        call VecGetOwnershipRange(gx,rstart,rend,ierr)

103:        do 10, i=rstart,rend-1
104:          value = i
105:          call VecSetValues(gx,1,i,value,INSERT_VALUES,ierr)
106:  10    continue

108:        call VecAssemblyBegin(gx,ierr)
109:        call VecAssemblyEnd(gx,ierr)

111:        call VecGhostUpdateBegin(gx,INSERT_VALUES,SCATTER_FORWARD,ierr)
112:        call VecGhostUpdateEnd(gx,INSERT_VALUES,SCATTER_FORWARD,ierr)

114: !     Print out each vector, including the ghost padding region.

116:        call VecView(lx,PETSC_VIEWER_STDOUT_SELF,ierr)

118:        call VecGhostRestoreLocalForm(gx,lx,ierr)
119:        call VecDestroy(gx,ierr)
120:        call PetscFinalize(ierr)
121:        end
122: