xref: /petsc/src/dm/tests/ex1f.F90 (revision bfe80ac4a46d58cb7760074b25f5e81b2f541d8a)
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>
6program main
7#include "petsc/finclude/petscdmda.h"
8  use petscvec
9  use petscdm
10  use petscdmda
11  implicit none
12
13  PetscInt, parameter :: Ndof=1, stencil_size=1
14  PetscInt, parameter :: Nx=3, Ny=3
15  PetscErrorCode :: myid, commsize, ierr
16  PetscScalar, pointer :: xv1d(:)
17  PetscInt, pointer :: lx(:), ly(:)
18  PetscMPIInt, pointer :: nb(:)
19
20  type(tDM) :: da
21  type(tVec) :: gVec!, naturalVec
22
23  PetscCallA(PetscInitialize(PETSC_NULL_CHARACTER, ierr))
24  PetscCallA(mpi_comm_rank(PETSC_COMM_WORLD, myid, ierr))
25  PetscCallA(mpi_comm_size(PETSC_COMM_WORLD, commsize, ierr))
26
27  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_ARRAY, PETSC_NULL_INTEGER_ARRAY, da, ierr))
28  PetscCallA(DMSetFromOptions(da, ierr))
29  PetscCallA(DMSetup(da, ierr))
30
31  PetscCallA(DMCreateGlobalVector(da, gVec, ierr))
32  PetscCallA(VecGetArray(gVec, xv1d, ierr))
33  xv1d(:) = real(myid, kind(xv1d))
34  !print *,myid, 'xv1d', xv1d, ':', xv1d
35  PetscCallA(VecRestoreArray(gVec, xv1d, ierr))
36
37  PetscCallA(PetscObjectViewFromOptions(PetscObjectCast(gVec), PETSC_NULL_OBJECT, '-show_gVec', ierr))
38
39  PetscCallA(DMDAGetOwnershipRanges(da, lx, ly, PETSC_NULL_INTEGER_POINTER, ierr))
40  PetscCallA(DMDARestoreOwnershipRanges(da, lx, ly, PETSC_NULL_INTEGER_POINTER, ierr))
41  PetscCallA(DMDAGetNeighbors(da, nb, ierr))
42  PetscCallA(DMDARestoreNeighbors(da, nb, ierr))
43
44  PetscCallA(VecDestroy(gVec, ierr))
45  PetscCallA(DMDestroy(da, ierr))
46  PetscCallA(PetscFinalize(ierr))
47end program
48
49!/*TEST
50!
51!   test:
52!      nsize: 9
53!      args: -show_gVec
54!TEST*/
55