static char help[] = "Tests the use of MatZeroRows() for parallel matrices.\n\ This example also tests the use of MatDuplicate() for both MPIAIJ and MPIBAIJ matrices"; #include extern PetscErrorCode TestMatZeroRows_Basic(Mat,IS,PetscScalar); extern PetscErrorCode TestMatZeroRows_with_no_allocation(Mat,IS,PetscScalar); int main(int argc,char **args) { Mat A; PetscInt i,j,m = 3,n,Ii,J,Imax; PetscMPIInt rank,size; PetscScalar v,diag=-4.0; IS is; CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); n = 2*size; /* create A Square matrix for the five point stencil,YET AGAIN*/ CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n)); CHKERRQ(MatSetFromOptions(A)); CHKERRQ(MatSetUp(A)); for (i=0; i0) {J = Ii - n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} if (i0) {J = Ii - 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} if (j= n*m -m - 1) Imax = m*n - m - 1; CHKERRQ(ISCreateStride(PETSC_COMM_SELF,m,Imax,1,&is)); CHKERRQ(TestMatZeroRows_Basic(A,is,0.0)); CHKERRQ(TestMatZeroRows_Basic(A,is,diag)); CHKERRQ(TestMatZeroRows_with_no_allocation(A,is,0.0)); CHKERRQ(TestMatZeroRows_with_no_allocation(A,is,diag)); CHKERRQ(MatDestroy(&A)); /* Now Create a rectangular matrix with five point stencil (app) n+size is used so that this dimension is always divisible by size. This way, we can always use bs = size for any number of procs */ CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*(n+size))); CHKERRQ(MatSetFromOptions(A)); CHKERRQ(MatSetUp(A)); for (i=0; i0) {J = Ii - n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} if (i0) {J = Ii - 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} if (j