10dd791a8SStefano Zampini static char help[] = "Test VecScatter of different block sizes across processes\n\n"; 297929ea7SJunchao Zhang 397929ea7SJunchao Zhang #include <petscvec.h> 4d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 5d71ae5a4SJacob Faibussowitsch { 697929ea7SJunchao Zhang PetscInt i, bs, n, low, high; 797929ea7SJunchao Zhang PetscMPIInt nproc, rank; 897929ea7SJunchao Zhang Vec x, y, z; 997929ea7SJunchao Zhang IS ix, iy; 1097929ea7SJunchao Zhang VecScatter vscat; 1197929ea7SJunchao Zhang const PetscScalar *yv; 1297929ea7SJunchao Zhang 13327415f7SBarry Smith PetscFunctionBeginUser; 14*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 159566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &nproc)); 169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 1708401ef6SPierre Jolivet PetscCheck(nproc == 2, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This test can only run on two MPI ranks"); 1897929ea7SJunchao Zhang 1997929ea7SJunchao Zhang /* Create an MPI vector x of size 12 on two processes, and set x = {0, 1, 2, .., 11} */ 209566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(PETSC_COMM_WORLD, 6, PETSC_DECIDE, &x)); 219566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(x, &low, &high)); 229566063dSJacob Faibussowitsch for (i = low; i < high; i++) PetscCall(VecSetValue(x, i, (PetscScalar)i, INSERT_VALUES)); 239566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(x)); 249566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(x)); 2597929ea7SJunchao Zhang 2697929ea7SJunchao Zhang /* Create a seq vector y, and a parallel to sequential (PtoS) vecscatter to scatter x to y */ 27dd400576SPatrick Sanan if (rank == 0) { 2897929ea7SJunchao Zhang /* On rank 0, seq y is of size 6. We will scatter x[0,1,2,6,7,8] to y[0,1,2,3,4,5] using IS with bs=3 */ 2997929ea7SJunchao Zhang PetscInt idx[2] = {0, 2}; 3097929ea7SJunchao Zhang PetscInt idy[2] = {0, 1}; 3197929ea7SJunchao Zhang n = 6; 3297929ea7SJunchao Zhang bs = 3; 339566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &y)); 349566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, 2, idx, PETSC_COPY_VALUES, &ix)); 359566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, 2, idy, PETSC_COPY_VALUES, &iy)); 3697929ea7SJunchao Zhang } else { 3797929ea7SJunchao Zhang /* On rank 1, seq y is of size 4. We will scatter x[4,5,10,11] to y[0,1,2,3] using IS with bs=2 */ 3897929ea7SJunchao Zhang PetscInt idx[2] = {2, 5}; 3997929ea7SJunchao Zhang PetscInt idy[2] = {0, 1}; 4097929ea7SJunchao Zhang n = 4; 4197929ea7SJunchao Zhang bs = 2; 429566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &y)); 439566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, 2, idx, PETSC_COPY_VALUES, &ix)); 449566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, 2, idy, PETSC_COPY_VALUES, &iy)); 4597929ea7SJunchao Zhang } 469566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(x, ix, y, iy, &vscat)); 4797929ea7SJunchao Zhang 4897929ea7SJunchao Zhang /* Do the vecscatter */ 499566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(vscat, x, y, INSERT_VALUES, SCATTER_FORWARD)); 509566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(vscat, x, y, INSERT_VALUES, SCATTER_FORWARD)); 5197929ea7SJunchao Zhang 5297929ea7SJunchao Zhang /* Print y. Since y is sequential, we put y in a parallel z to print its value on both ranks */ 539566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(y, &yv)); 549566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, n, PETSC_DECIDE, yv, &z)); 559566063dSJacob Faibussowitsch PetscCall(VecView(z, PETSC_VIEWER_STDOUT_WORLD)); 569566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(y, &yv)); 5797929ea7SJunchao Zhang 589566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ix)); 599566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iy)); 609566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&z)); 639566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&vscat)); 6497929ea7SJunchao Zhang 659566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 66b122ec5aSJacob Faibussowitsch return 0; 6797929ea7SJunchao Zhang } 6897929ea7SJunchao Zhang 6997929ea7SJunchao Zhang /*TEST 7097929ea7SJunchao Zhang 7197929ea7SJunchao Zhang test: 7297929ea7SJunchao Zhang nsize: 2 7397929ea7SJunchao Zhang args: 7497929ea7SJunchao Zhang requires: 7597929ea7SJunchao Zhang TEST*/ 76