Actual source code: ex33.c
2: static char help[] = "Tests the routines VecScatterCreateToAll(), VecScatterCreateToZero()\n\n";
4: #include <petscvec.h>
6: int main(int argc, char **argv)
7: {
8: PetscInt n = 3, i, len, start, end;
9: PetscMPIInt size, rank;
10: PetscScalar value, *yy;
11: Vec x, y, z, y_t;
12: VecScatter toall, tozero;
14: PetscFunctionBeginUser;
15: PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
16: PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
17: PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
19: /* create two vectors */
20: PetscCall(VecCreateMPI(PETSC_COMM_WORLD, PETSC_DECIDE, size * n, &x));
22: /* each processor inserts its values */
24: PetscCall(VecGetOwnershipRange(x, &start, &end));
25: for (i = start; i < end; i++) {
26: value = (PetscScalar)i;
27: PetscCall(VecSetValues(x, 1, &i, &value, INSERT_VALUES));
28: }
29: PetscCall(VecAssemblyBegin(x));
30: PetscCall(VecAssemblyEnd(x));
31: PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
33: PetscCall(VecScatterCreateToAll(x, &toall, &y));
34: PetscCall(VecScatterBegin(toall, x, y, INSERT_VALUES, SCATTER_FORWARD));
35: PetscCall(VecScatterEnd(toall, x, y, INSERT_VALUES, SCATTER_FORWARD));
36: PetscCall(VecScatterDestroy(&toall));
38: /* Cannot view the above vector with VecView(), so place it in an MPI Vec
39: and do a VecView() */
40: PetscCall(VecGetArray(y, &yy));
41: PetscCall(VecGetLocalSize(y, &len));
42: PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, len, PETSC_DECIDE, yy, &y_t));
43: PetscCall(VecView(y_t, PETSC_VIEWER_STDOUT_WORLD));
44: PetscCall(VecDestroy(&y_t));
45: PetscCall(VecRestoreArray(y, &yy));
47: PetscCall(VecScatterCreateToAll(x, &tozero, &z));
48: PetscCall(VecScatterBegin(tozero, x, z, INSERT_VALUES, SCATTER_FORWARD));
49: PetscCall(VecScatterEnd(tozero, x, z, INSERT_VALUES, SCATTER_FORWARD));
50: PetscCall(VecScatterDestroy(&tozero));
51: if (rank == 0) PetscCall(VecView(z, PETSC_VIEWER_STDOUT_SELF));
52: PetscCall(VecDestroy(&z));
54: PetscCall(VecScatterCreateToZero(x, &tozero, &z));
55: PetscCall(VecScatterBegin(tozero, x, z, INSERT_VALUES, SCATTER_FORWARD));
56: PetscCall(VecScatterEnd(tozero, x, z, INSERT_VALUES, SCATTER_FORWARD));
57: PetscCall(VecScatterDestroy(&tozero));
58: PetscCall(VecDestroy(&z));
60: PetscCall(VecDestroy(&x));
61: PetscCall(VecDestroy(&y));
63: PetscCall(PetscFinalize());
64: return 0;
65: }
67: /*TEST
69: test:
70: nsize: 4
72: TEST*/