197929ea7SJunchao Zhang 297929ea7SJunchao Zhang static char help[] = "Scatters from a sequential vector to a parallel vector. \n\ 397929ea7SJunchao Zhang uses block index sets\n\n"; 497929ea7SJunchao Zhang 597929ea7SJunchao Zhang #include <petscvec.h> 697929ea7SJunchao Zhang 7*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 8*d71ae5a4SJacob Faibussowitsch { 997929ea7SJunchao Zhang PetscInt bs = 1, n = 5, i, low; 1097929ea7SJunchao Zhang PetscInt ix0[3] = {5, 7, 9}, iy0[3] = {1, 2, 4}, ix1[3] = {2, 3, 4}, iy1[3] = {0, 1, 3}; 1197929ea7SJunchao Zhang PetscMPIInt size, rank; 1297929ea7SJunchao Zhang PetscScalar *array; 1397929ea7SJunchao Zhang Vec x, y; 1497929ea7SJunchao Zhang IS isx, isy; 1597929ea7SJunchao Zhang VecScatter ctx; 1697929ea7SJunchao Zhang PetscViewer sviewer; 1797929ea7SJunchao Zhang 18327415f7SBarry Smith PetscFunctionBeginUser; 199566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 2297929ea7SJunchao Zhang 2308401ef6SPierre Jolivet PetscCheck(size >= 2, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Must run more than one processor"); 2497929ea7SJunchao Zhang 259566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-bs", &bs, NULL)); 2697929ea7SJunchao Zhang n = bs * n; 2797929ea7SJunchao Zhang 2897929ea7SJunchao Zhang /* Create vector x over shared memory */ 299566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &x)); 309566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x, n, PETSC_DECIDE)); 319566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x)); 3297929ea7SJunchao Zhang 339566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(x, &low, NULL)); 349566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &array)); 35ad540459SPierre Jolivet for (i = 0; i < n; i++) array[i] = (PetscScalar)(i + low); 369566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &array)); 3797929ea7SJunchao Zhang 3897929ea7SJunchao Zhang /* Create a sequential vector y */ 399566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &y)); 409566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &array)); 41ad540459SPierre Jolivet for (i = 0; i < n; i++) array[i] = (PetscScalar)(i + 100 * rank); 429566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &array)); 4397929ea7SJunchao Zhang 4497929ea7SJunchao Zhang /* Create two index sets */ 45dd400576SPatrick Sanan if (rank == 0) { 469566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, 3, ix0, PETSC_COPY_VALUES, &isx)); 479566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, 3, iy0, PETSC_COPY_VALUES, &isy)); 4897929ea7SJunchao Zhang } else { 499566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, 3, ix1, PETSC_COPY_VALUES, &isx)); 509566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, 3, iy1, PETSC_COPY_VALUES, &isy)); 5197929ea7SJunchao Zhang } 5297929ea7SJunchao Zhang 5397929ea7SJunchao Zhang if (rank == 10) { 549566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n[%d] isx:\n", rank)); 559566063dSJacob Faibussowitsch PetscCall(ISView(isx, PETSC_VIEWER_STDOUT_SELF)); 5697929ea7SJunchao Zhang } 5797929ea7SJunchao Zhang 589566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(y, isy, x, isx, &ctx)); 599566063dSJacob Faibussowitsch PetscCall(VecScatterSetFromOptions(ctx)); 6097929ea7SJunchao Zhang 6197929ea7SJunchao Zhang /* Test forward vecscatter */ 629566063dSJacob Faibussowitsch PetscCall(VecSet(x, 0.0)); 639566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ctx, y, x, ADD_VALUES, SCATTER_FORWARD)); 649566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ctx, y, x, ADD_VALUES, SCATTER_FORWARD)); 659566063dSJacob Faibussowitsch PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD)); 6697929ea7SJunchao Zhang 6797929ea7SJunchao Zhang /* Test reverse vecscatter */ 689566063dSJacob Faibussowitsch PetscCall(VecScale(x, -1.0)); 699566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ctx, x, y, ADD_VALUES, SCATTER_REVERSE)); 709566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ctx, x, y, ADD_VALUES, SCATTER_REVERSE)); 719566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &sviewer)); 7248a46eb9SPierre Jolivet if (rank == 1) PetscCall(VecView(y, sviewer)); 739566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &sviewer)); 7497929ea7SJunchao Zhang 7597929ea7SJunchao Zhang /* Free spaces */ 769566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ctx)); 779566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isx)); 789566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isy)); 799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 819566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 82b122ec5aSJacob Faibussowitsch return 0; 8397929ea7SJunchao Zhang } 8497929ea7SJunchao Zhang 8597929ea7SJunchao Zhang /*TEST 8697929ea7SJunchao Zhang 8797929ea7SJunchao Zhang test: 8897929ea7SJunchao Zhang nsize: 3 8997929ea7SJunchao Zhang 9097929ea7SJunchao Zhang TEST*/ 91