xref: /petsc/src/vec/is/sf/tests/ex1f.F90 (revision 030f984af8d8bb4c203755d35bded3c05b3d83ce)
1!
2!  Tests VecScatterCreateToAll Fortran stub
3      program main
4#include <petsc/finclude/petscvec.h>
5      use petscvec
6      implicit none
7
8      PetscErrorCode ierr
9      PetscInt  nlocal, row
10      PetscScalar num
11      PetscMPIInt rank
12      Vec v1, v2
13      VecScatter toall
14
15      call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
16      if (ierr .ne. 0) then
17        print*,'Unable to initialize PETSc'
18        stop
19      endif
20      call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
21
22      nlocal = 1
23      call VecCreateMPI(PETSC_COMM_WORLD,nlocal,PETSC_DECIDE,v1,ierr)
24
25      row = rank
26      num = rank
27      call VecSetValue(v1,row,num,INSERT_VALUES,ierr)
28      call VecAssemblyBegin(v1,ierr)
29      call VecAssemblyEnd(v1,ierr)
30
31      call VecScatterCreateToAll(v1,toall,v2,ierr)
32
33      call VecScatterBegin(toall,v1,v2,INSERT_VALUES,SCATTER_FORWARD,ierr)
34      call VecScatterEnd(toall,v1,v2,INSERT_VALUES,SCATTER_FORWARD,ierr)
35
36      call VecScatterDestroy(toall,ierr)
37! Destroy v2 and then re-create it in VecScatterCreateToAll() to test if petsc can differentiate NULL projects with destroyed objects
38      call VecDestroy(v2,ierr)
39
40      call VecScatterCreateToAll(v1,toall,v2,ierr)
41      call VecScatterBegin(toall,v1,v2,INSERT_VALUES,SCATTER_FORWARD,ierr)
42      call VecScatterEnd(toall,v1,v2,INSERT_VALUES,SCATTER_FORWARD,ierr)
43
44      if (rank.eq.2) then
45         call PetscObjectSetName(v2, 'v2',ierr)
46         call VecView(v2,PETSC_VIEWER_STDOUT_SELF,ierr)
47      end if
48
49      call VecScatterDestroy(toall,ierr)
50      call VecDestroy(v1,ierr)
51      call VecDestroy(v2,ierr)
52! It is OK to destroy again
53      call VecDestroy(v2,ierr)
54
55      call PetscFinalize(ierr)
56      end
57
58!/*TEST
59!
60!     test:
61!       nsize: 4
62!
63!TEST*/
64