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