xref: /petsc/src/vec/is/sf/tests/ex12.c (revision 732aec7a18f2199fb53bb9a2f3aef439a834ce31)
197929ea7SJunchao Zhang static char help[] = "Scatters from a parallel vector to a sequential vector. \n\
297929ea7SJunchao Zhang uses block index sets\n\n";
397929ea7SJunchao Zhang 
497929ea7SJunchao Zhang #include <petscvec.h>
597929ea7SJunchao Zhang 
main(int argc,char ** argv)6d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
7d71ae5a4SJacob Faibussowitsch {
897929ea7SJunchao Zhang   PetscInt     bs = 1, n = 5, i, low;
997929ea7SJunchao Zhang   PetscInt     ix0[3] = {5, 7, 9}, iy0[3] = {1, 2, 4}, ix1[3] = {2, 3, 4}, iy1[3] = {0, 1, 3};
1097929ea7SJunchao Zhang   PetscMPIInt  size, rank;
1197929ea7SJunchao Zhang   PetscScalar *array;
1297929ea7SJunchao Zhang   Vec          x, y;
1397929ea7SJunchao Zhang   IS           isx, isy;
1497929ea7SJunchao Zhang   VecScatter   ctx;
1597929ea7SJunchao Zhang 
16327415f7SBarry Smith   PetscFunctionBeginUser;
17*c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
199566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
2097929ea7SJunchao Zhang 
2108401ef6SPierre Jolivet   PetscCheck(size >= 2, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Must run more than one processor");
2297929ea7SJunchao Zhang 
239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-bs", &bs, NULL));
2497929ea7SJunchao Zhang   n = bs * n;
2597929ea7SJunchao Zhang 
2697929ea7SJunchao Zhang   /* Create vector x over shared memory */
279566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
289566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x, n, PETSC_DECIDE));
299566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
3097929ea7SJunchao Zhang 
319566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(x, &low, NULL));
329566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x, &array));
33ad540459SPierre Jolivet   for (i = 0; i < n; i++) array[i] = (PetscScalar)(i + low);
349566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(x, &array));
3597929ea7SJunchao Zhang 
3697929ea7SJunchao Zhang   /* Create a sequential vector y */
379566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &y));
389566063dSJacob Faibussowitsch   PetscCall(VecSet(y, 0.0));
3997929ea7SJunchao Zhang 
4097929ea7SJunchao Zhang   /* Create two index sets */
41dd400576SPatrick Sanan   if (rank == 0) {
429566063dSJacob Faibussowitsch     PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, 3, ix0, PETSC_COPY_VALUES, &isx));
439566063dSJacob Faibussowitsch     PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, 3, iy0, PETSC_COPY_VALUES, &isy));
4497929ea7SJunchao Zhang   } else {
459566063dSJacob Faibussowitsch     PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, 3, ix1, PETSC_COPY_VALUES, &isx));
469566063dSJacob Faibussowitsch     PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, 3, iy1, PETSC_COPY_VALUES, &isy));
4797929ea7SJunchao Zhang   }
4897929ea7SJunchao Zhang 
4997929ea7SJunchao Zhang   if (rank == 10) {
509566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n[%d] isx:\n", rank));
519566063dSJacob Faibussowitsch     PetscCall(ISView(isx, PETSC_VIEWER_STDOUT_SELF));
5297929ea7SJunchao Zhang   }
5397929ea7SJunchao Zhang 
549566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(x, isx, y, isy, &ctx));
559566063dSJacob Faibussowitsch   PetscCall(VecScatterSetFromOptions(ctx));
5697929ea7SJunchao Zhang 
5797929ea7SJunchao Zhang   /* Test forward vecscatter */
589566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(ctx, x, y, ADD_VALUES, SCATTER_FORWARD));
599566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(ctx, x, y, ADD_VALUES, SCATTER_FORWARD));
6097929ea7SJunchao Zhang   if (rank == 0) {
619566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] y:\n", rank));
629566063dSJacob Faibussowitsch     PetscCall(VecView(y, PETSC_VIEWER_STDOUT_SELF));
6397929ea7SJunchao Zhang   }
6497929ea7SJunchao Zhang 
6597929ea7SJunchao Zhang   /* Test reverse vecscatter */
669566063dSJacob Faibussowitsch   PetscCall(VecScale(y, -1.0));
6748a46eb9SPierre Jolivet   if (rank) PetscCall(VecScale(y, 1.0 / (size - 1)));
6897929ea7SJunchao Zhang 
699566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(ctx, y, x, ADD_VALUES, SCATTER_REVERSE));
709566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(ctx, y, x, ADD_VALUES, SCATTER_REVERSE));
719566063dSJacob Faibussowitsch   PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
7297929ea7SJunchao Zhang 
7397929ea7SJunchao Zhang   /* Free spaces */
749566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&ctx));
759566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isx));
769566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isy));
779566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
789566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
799566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
80b122ec5aSJacob Faibussowitsch   return 0;
8197929ea7SJunchao Zhang }
8297929ea7SJunchao Zhang 
8397929ea7SJunchao Zhang /*TEST
8497929ea7SJunchao Zhang 
8597929ea7SJunchao Zhang    test:
8697929ea7SJunchao Zhang       nsize: 3
8797929ea7SJunchao Zhang 
8897929ea7SJunchao Zhang TEST*/
89