Actual source code: ex9f.F90

petsc-3.8.4 2018-03-24
Report Typos and Errors
  1: !
  2: !
  3: ! Description: Illustrates the use of VecCreateGhost()
  4: !
  5: !/*T
  6: !   Concepts: vectors^assembling vectors;
  7: !   Concepts: vectors^ghost padding;
  8: !   Processors: n
  9: !
 10: !   Description: 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:  #include <petsc/finclude/petscvec.h>
 19:       use petscvec
 20:       implicit none

 22:       PetscMPIInt rank,size
 23:       PetscInt nlocal,nghost,ifrom(2)
 24:       PetscErrorCode ierr
 25:       PetscInt i,rstart,rend,ione
 26:       PetscBool   flag
 27:       PetscScalar  value,tarray(20)
 28:       Vec          lx,gx,gxs

 30:       nlocal = 6
 31:       nghost = 2

 33:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
 34:       if (ierr .ne. 0) then
 35:         print*,'PetscInitialize failed'
 36:         stop
 37:       endif
 38:       call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
 39:       call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)

 41:       if (size .ne. 2) then
 42:        SETERRA(PETSC_COMM_WORLD,1,'Requires 2 processors')
 43:       endif

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


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

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

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


 89: !      Test VecDuplicate

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

 94: !      Access the local Form

 96:        call VecGhostGetLocalForm(gx,lx,ierr)

 98: !     Set the values from 0 to 12 into the 'global' vector

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

102:        ione = 1
103:        do 10, i=rstart,rend-1
104:          value = i
105:          call VecSetValues(gx,ione,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:        if (rank .eq. 0) then
117:           call VecView(lx,PETSC_VIEWER_STDOUT_SELF,ierr)
118:        endif
119:        call MPI_Barrier(PETSC_COMM_WORLD,ierr)
120:        if (rank .eq. 1) then
121:           call VecView(lx,PETSC_VIEWER_STDOUT_SELF,ierr)
122:        endif

124:        call VecGhostRestoreLocalForm(gx,lx,ierr)
125:        call VecDestroy(gx,ierr)
126:        call PetscFinalize(ierr)
127:        end