Actual source code: ex17f.F
1: !
2: !
3: ! "Scatters from a parallel vector to a sequential vector. In
4: ! this case each local vector is as long as the entire parallel vector.
5: !
6: implicit none
7: #include "include/finclude/petsc.h"
8: #include "include/finclude/petscis.h"
9: #include "include/finclude/petscvec.h"
10: #include "include/finclude/petscsys.h"
12: PetscErrorCode ierr
13: PetscMPIInt size,rank
14: PetscInt n,NN,low,high,iglobal,i,ione,first,stride
15: PetscScalar value,zero
16: Vec x,y
17: IS is1,is2
18: VecScatter ctx
20: n = 5
21: zero = 0.d0
22: call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
24: call MPI_Comm_size(MPI_COMM_WORLD,size,ierr)
25: call MPI_Comm_rank(MPI_COMM_WORLD,rank,ierr)
27: ! create two vectors
28: ! one parallel and one sequential. The sequential one on each processor
29: ! is as long as the entire parallel one.
31: NN = size*n
33: call VecCreateMPI(MPI_COMM_WORLD,PETSC_DECIDE,NN,y,ierr)
34: call VecCreateSeq(MPI_COMM_SELF,NN,x,ierr)
36: call VecSet(x,zero,ierr)
37: call VecGetOwnershipRange(y,low,high,ierr)
38: ione = 1
39: do 10, i=0,n-1
40: iglobal = i + low
41: value = i + 10*rank
42: call VecSetValues(y,ione,iglobal,value,INSERT_VALUES,ierr)
43: 10 continue
45: call VecAssemblyBegin(y,ierr)
46: call VecAssemblyEnd(y,ierr)
47: !
48: ! View the parallel vector
49: !
50: call VecView(y,PETSC_VIEWER_STDOUT_WORLD,ierr)
52: ! create two index sets and the scatter context to move the contents of
53: ! of the parallel vector to each sequential vector. If you want the
54: ! parallel vector delivered to only one processor then create a is2
55: ! of length zero on all processors except the one to receive the parallel vector
57: first = 0
58: stride = 1
59: call ISCreateStride(PETSC_COMM_SELF,NN,first,stride,is1,ierr)
60: call ISCreateStride(PETSC_COMM_SELF,NN,first,stride,is2,ierr)
61: call VecScatterCreate(y,is2,x,is1,ctx,ierr)
62: call VecScatterBegin(y,x,ADD_VALUES,SCATTER_FORWARD,ctx,ierr)
63: call VecScatterEnd(y,x,ADD_VALUES,SCATTER_FORWARD,ctx,ierr)
64: call VecScatterDestroy(ctx,ierr)
65: !
66: ! View the sequential vector on the 0th processor
67: !
68: if (rank .eq. 0) then
69: call VecView(x,PETSC_VIEWER_STDOUT_SELF,ierr)
70: endif
72: call VecDestroy(x,ierr)
73: call VecDestroy(y,ierr)
74: call ISDestroy(is1,ierr)
75: call ISDestroy(is2,ierr)
77: call PetscFinalize(ierr)
78: end
79: