197929ea7SJunchao Zhang static char help[] = "Scatters from a sequential vector to a parallel 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 PetscViewer sviewer;
1697929ea7SJunchao Zhang
17327415f7SBarry Smith PetscFunctionBeginUser;
18*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
2197929ea7SJunchao Zhang
2208401ef6SPierre Jolivet PetscCheck(size >= 2, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Must run more than one processor");
2397929ea7SJunchao Zhang
249566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-bs", &bs, NULL));
2597929ea7SJunchao Zhang n = bs * n;
2697929ea7SJunchao Zhang
2797929ea7SJunchao Zhang /* Create vector x over shared memory */
289566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
299566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x, n, PETSC_DECIDE));
309566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x));
3197929ea7SJunchao Zhang
329566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(x, &low, NULL));
339566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &array));
34ad540459SPierre Jolivet for (i = 0; i < n; i++) array[i] = (PetscScalar)(i + low);
359566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &array));
3697929ea7SJunchao Zhang
3797929ea7SJunchao Zhang /* Create a sequential vector y */
389566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &y));
399566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &array));
40ad540459SPierre Jolivet for (i = 0; i < n; i++) array[i] = (PetscScalar)(i + 100 * rank);
419566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &array));
4297929ea7SJunchao Zhang
4397929ea7SJunchao Zhang /* Create two index sets */
44dd400576SPatrick Sanan if (rank == 0) {
459566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, 3, ix0, PETSC_COPY_VALUES, &isx));
469566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, 3, iy0, PETSC_COPY_VALUES, &isy));
4797929ea7SJunchao Zhang } else {
489566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, 3, ix1, PETSC_COPY_VALUES, &isx));
499566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, 3, iy1, PETSC_COPY_VALUES, &isy));
5097929ea7SJunchao Zhang }
5197929ea7SJunchao Zhang
5297929ea7SJunchao Zhang if (rank == 10) {
539566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n[%d] isx:\n", rank));
549566063dSJacob Faibussowitsch PetscCall(ISView(isx, PETSC_VIEWER_STDOUT_SELF));
5597929ea7SJunchao Zhang }
5697929ea7SJunchao Zhang
579566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(y, isy, x, isx, &ctx));
589566063dSJacob Faibussowitsch PetscCall(VecScatterSetFromOptions(ctx));
5997929ea7SJunchao Zhang
6097929ea7SJunchao Zhang /* Test forward vecscatter */
619566063dSJacob Faibussowitsch PetscCall(VecSet(x, 0.0));
629566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ctx, y, x, ADD_VALUES, SCATTER_FORWARD));
639566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ctx, y, x, ADD_VALUES, SCATTER_FORWARD));
649566063dSJacob Faibussowitsch PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
6597929ea7SJunchao Zhang
6697929ea7SJunchao Zhang /* Test reverse vecscatter */
679566063dSJacob Faibussowitsch PetscCall(VecScale(x, -1.0));
689566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ctx, x, y, ADD_VALUES, SCATTER_REVERSE));
699566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ctx, x, y, ADD_VALUES, SCATTER_REVERSE));
709566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &sviewer));
7148a46eb9SPierre Jolivet if (rank == 1) PetscCall(VecView(y, sviewer));
729566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &sviewer));
7397929ea7SJunchao Zhang
7497929ea7SJunchao Zhang /* Free spaces */
759566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ctx));
769566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isx));
779566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isy));
789566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y));
809566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
81b122ec5aSJacob Faibussowitsch return 0;
8297929ea7SJunchao Zhang }
8397929ea7SJunchao Zhang
8497929ea7SJunchao Zhang /*TEST
8597929ea7SJunchao Zhang
8697929ea7SJunchao Zhang test:
8797929ea7SJunchao Zhang nsize: 3
8897929ea7SJunchao Zhang
8997929ea7SJunchao Zhang TEST*/
90