1 2 static char help[]= "Scatters between parallel vectors. \n\ 3 uses block index sets\n\n"; 4 5 #include <petscvec.h> 6 7 int main(int argc,char **argv) 8 { 9 PetscErrorCode ierr; 10 PetscInt bs=1,n=5,N,i,low; 11 PetscInt ix0[3] = {5,7,9},iy0[3] = {1,2,4},ix1[3] = {2,3,1},iy1[3] = {0,3,9}; 12 PetscMPIInt size,rank; 13 PetscScalar *array; 14 Vec x,x1,y; 15 IS isx,isy; 16 VecScatter ctx; 17 VecScatterType type; 18 PetscBool flg; 19 20 ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 21 ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 22 ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 23 24 if (size <2) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Must run more than one processor"); 25 26 ierr = PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL);CHKERRQ(ierr); 27 n = bs*n; 28 29 /* Create vector x over shared memory */ 30 ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr); 31 ierr = VecSetSizes(x,n,PETSC_DECIDE);CHKERRQ(ierr); 32 ierr = VecSetFromOptions(x);CHKERRQ(ierr); 33 34 ierr = VecGetOwnershipRange(x,&low,NULL);CHKERRQ(ierr); 35 ierr = VecGetArray(x,&array);CHKERRQ(ierr); 36 for (i=0; i<n; i++) { 37 array[i] = (PetscScalar)(i + low); 38 } 39 ierr = VecRestoreArray(x,&array);CHKERRQ(ierr); 40 /* ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 41 42 /* Test some vector functions */ 43 ierr = VecAssemblyBegin(x);CHKERRQ(ierr); 44 ierr = VecAssemblyEnd(x);CHKERRQ(ierr); 45 46 ierr = VecGetSize(x,&N);CHKERRQ(ierr); 47 ierr = VecGetLocalSize(x,&n);CHKERRQ(ierr); 48 49 ierr = VecDuplicate(x,&x1);CHKERRQ(ierr); 50 ierr = VecCopy(x,x1);CHKERRQ(ierr); 51 ierr = VecEqual(x,x1,&flg);CHKERRQ(ierr); 52 if (!flg) SETERRQ(PetscObjectComm((PetscObject)x),PETSC_ERR_ARG_WRONG,"x1 != x"); 53 54 ierr = VecScale(x1,2.0);CHKERRQ(ierr); 55 ierr = VecSet(x1,10.0);CHKERRQ(ierr); 56 /* ierr = VecView(x1,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 57 58 /* Create vector y over shared memory */ 59 ierr = VecCreate(PETSC_COMM_WORLD,&y);CHKERRQ(ierr); 60 ierr = VecSetSizes(y,n,PETSC_DECIDE);CHKERRQ(ierr); 61 ierr = VecSetFromOptions(y);CHKERRQ(ierr); 62 ierr = VecGetArray(y,&array);CHKERRQ(ierr); 63 for (i=0; i<n; i++) { 64 array[i] = -(PetscScalar) (i + 100*rank); 65 } 66 ierr = VecRestoreArray(y,&array);CHKERRQ(ierr); 67 ierr = VecAssemblyBegin(y);CHKERRQ(ierr); 68 ierr = VecAssemblyEnd(y);CHKERRQ(ierr); 69 /* ierr = VecView(y,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 70 71 /* Create two index sets */ 72 if (rank == 0) { 73 ierr = ISCreateBlock(PETSC_COMM_SELF,bs,3,ix0,PETSC_COPY_VALUES,&isx);CHKERRQ(ierr); 74 ierr = ISCreateBlock(PETSC_COMM_SELF,bs,3,iy0,PETSC_COPY_VALUES,&isy);CHKERRQ(ierr); 75 } else { 76 ierr = ISCreateBlock(PETSC_COMM_SELF,bs,3,ix1,PETSC_COPY_VALUES,&isx);CHKERRQ(ierr); 77 ierr = ISCreateBlock(PETSC_COMM_SELF,bs,3,iy1,PETSC_COPY_VALUES,&isy);CHKERRQ(ierr); 78 } 79 80 if (rank == 10) { 81 ierr = PetscPrintf(PETSC_COMM_SELF,"\n[%d] isx:\n",rank);CHKERRQ(ierr); 82 ierr = ISView(isx,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 83 ierr = PetscPrintf(PETSC_COMM_SELF,"\n[%d] isy:\n",rank);CHKERRQ(ierr); 84 ierr = ISView(isy,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 85 } 86 87 /* Create Vector scatter */ 88 ierr = VecScatterCreate(x,isx,y,isy,&ctx);CHKERRQ(ierr); 89 ierr = VecScatterSetFromOptions(ctx);CHKERRQ(ierr); 90 ierr = VecScatterGetType(ctx,&type);CHKERRQ(ierr); 91 ierr = PetscPrintf(PETSC_COMM_WORLD,"scatter type %s\n",type);CHKERRQ(ierr); 92 93 /* Test forward vecscatter */ 94 ierr = VecSet(y,0.0);CHKERRQ(ierr); 95 ierr = VecScatterBegin(ctx,x,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 96 ierr = VecScatterEnd(ctx,x,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 97 ierr = PetscPrintf(PETSC_COMM_WORLD,"\nSCATTER_FORWARD y:\n");CHKERRQ(ierr); 98 ierr = VecView(y,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 99 100 /* Test reverse vecscatter */ 101 ierr = VecSet(x,0.0);CHKERRQ(ierr); 102 ierr = VecSet(y,0.0);CHKERRQ(ierr); 103 ierr = VecGetOwnershipRange(y,&low,NULL);CHKERRQ(ierr); 104 ierr = VecGetArray(y,&array);CHKERRQ(ierr); 105 for (i=0; i<n; i++) { 106 array[i] = (PetscScalar)(i+ low); 107 } 108 ierr = VecRestoreArray(y,&array);CHKERRQ(ierr); 109 ierr = VecScatterBegin(ctx,y,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 110 ierr = VecScatterEnd(ctx,y,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 111 ierr = PetscPrintf(PETSC_COMM_WORLD,"\nSCATTER_REVERSE x:\n");CHKERRQ(ierr); 112 ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 113 114 /* Free objects */ 115 ierr = VecScatterDestroy(&ctx);CHKERRQ(ierr); 116 ierr = ISDestroy(&isx);CHKERRQ(ierr); 117 ierr = ISDestroy(&isy);CHKERRQ(ierr); 118 ierr = VecDestroy(&x);CHKERRQ(ierr); 119 ierr = VecDestroy(&x1);CHKERRQ(ierr); 120 ierr = VecDestroy(&y);CHKERRQ(ierr); 121 ierr = PetscFinalize(); 122 return ierr; 123 } 124 125 /*TEST 126 127 test: 128 nsize: 2 129 130 TEST*/ 131