Actual source code: ex1f.F90

  1: !
  2: ! Test the workaround for a bug in Open MPI 2.1.1 on Ubuntu 18.04.2
  3: ! See https://lists.mcs.anl.gov/pipermail/petsc-dev/2019-July/024803.html
  4: !
  5: ! Contributed-by:       Fabian Jakub  <Fabian.Jakub@physik.uni-muenchen.de>
  6: program main
  7: #include "petsc/finclude/petsc.h"

  9:   use petsc
 10:   implicit none

 12:   PetscInt, parameter :: Ndof=1, stencil_size=1
 13:   PetscInt, parameter :: Nx=3, Ny=3
 14:   PetscErrorCode :: myid, commsize, ierr
 15:   PetscScalar, pointer :: xv1d(:)

 17:   type(tDM) :: da
 18:   type(tVec) :: gVec!, naturalVec

 20:   PetscCallA(PetscInitialize(PETSC_NULL_CHARACTER, ierr))
 21:   PetscCallA(mpi_comm_rank(PETSC_COMM_WORLD, myid, ierr))
 22:   PetscCallA(mpi_comm_size(PETSC_COMM_WORLD, commsize, ierr))

 24:   PetscCallA(DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC, DM_BOUNDARY_PERIODIC,DMDA_STENCIL_STAR,Nx, Ny, PETSC_DECIDE, PETSC_DECIDE, Ndof, stencil_size,PETSC_NULL_INTEGER, PETSC_NULL_INTEGER, da, ierr))
 25:   PetscCallA(DMSetup(da, ierr))
 26:   PetscCallA(DMSetFromOptions(da, ierr))

 28:   PetscCallA(DMCreateGlobalVector(da, gVec, ierr))
 29:   PetscCallA(VecGetArrayF90(gVec, xv1d, ierr))
 30:   xv1d(:) = real(myid, kind(xv1d))
 31:   !print *,myid, 'xv1d', xv1d, ':', xv1d
 32:   PetscCallA(VecRestoreArrayF90(gVec, xv1d, ierr))

 34:   PetscCallA(PetscObjectViewFromOptions(gVec, PETSC_NULL_VEC, '-show_gVec', ierr))

 36:   PetscCallA(VecDestroy(gVec, ierr))
 37:   PetscCallA(DMDestroy(da, ierr))
 38:   PetscCallA(PetscFinalize(ierr))
 39: end program

 41: !/*TEST
 42: !
 43: !   test:
 44: !      nsize: 9
 45: !      args: -show_gVec
 46: !TEST*/