1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests the use of MatZeroRows() for parallel matrices.\n\ 3c4762a1bSJed Brown This example also tests the use of MatDuplicate() for both MPIAIJ and MPIBAIJ matrices"; 4c4762a1bSJed Brown 5c4762a1bSJed Brown #include <petscmat.h> 6c4762a1bSJed Brown 7c4762a1bSJed Brown extern PetscErrorCode TestMatZeroRows_Basic(Mat,IS,PetscScalar); 8c4762a1bSJed Brown extern PetscErrorCode TestMatZeroRows_with_no_allocation(Mat,IS,PetscScalar); 9c4762a1bSJed Brown 10c4762a1bSJed Brown int main(int argc,char **args) 11c4762a1bSJed Brown { 12c4762a1bSJed Brown Mat A; 13c4762a1bSJed Brown PetscInt i,j,m = 3,n,Ii,J,Imax; 14c4762a1bSJed Brown PetscMPIInt rank,size; 15c4762a1bSJed Brown PetscScalar v,diag=-4.0; 16c4762a1bSJed Brown IS is; 17c4762a1bSJed Brown 18*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 195f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 205f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 21c4762a1bSJed Brown n = 2*size; 22c4762a1bSJed Brown 23c4762a1bSJed Brown /* create A Square matrix for the five point stencil,YET AGAIN*/ 245f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 255f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n)); 265f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 275f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A)); 28c4762a1bSJed Brown for (i=0; i<m; i++) { 29c4762a1bSJed Brown for (j=2*rank; j<2*rank+2; j++) { 30c4762a1bSJed Brown v = -1.0; Ii = j + n*i; 315f80ce2aSJacob Faibussowitsch if (i>0) {J = Ii - n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 325f80ce2aSJacob Faibussowitsch if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 335f80ce2aSJacob Faibussowitsch if (j>0) {J = Ii - 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 345f80ce2aSJacob Faibussowitsch if (j<n-1) {J = Ii + 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 355f80ce2aSJacob Faibussowitsch v = 4.0; CHKERRQ(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 36c4762a1bSJed Brown } 37c4762a1bSJed Brown } 385f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 395f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 40c4762a1bSJed Brown 41c4762a1bSJed Brown /* Create AN IS required by MatZeroRows() */ 42c4762a1bSJed Brown Imax = n*rank; if (Imax>= n*m -m - 1) Imax = m*n - m - 1; 435f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(PETSC_COMM_SELF,m,Imax,1,&is)); 44c4762a1bSJed Brown 455f80ce2aSJacob Faibussowitsch CHKERRQ(TestMatZeroRows_Basic(A,is,0.0)); 465f80ce2aSJacob Faibussowitsch CHKERRQ(TestMatZeroRows_Basic(A,is,diag)); 47c4762a1bSJed Brown 485f80ce2aSJacob Faibussowitsch CHKERRQ(TestMatZeroRows_with_no_allocation(A,is,0.0)); 495f80ce2aSJacob Faibussowitsch CHKERRQ(TestMatZeroRows_with_no_allocation(A,is,diag)); 50c4762a1bSJed Brown 515f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 52c4762a1bSJed Brown 53c4762a1bSJed Brown /* Now Create a rectangular matrix with five point stencil (app) 54c4762a1bSJed Brown n+size is used so that this dimension is always divisible by size. 55c4762a1bSJed Brown This way, we can always use bs = size for any number of procs */ 565f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 575f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*(n+size))); 585f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 595f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A)); 60c4762a1bSJed Brown for (i=0; i<m; i++) { 61c4762a1bSJed Brown for (j=2*rank; j<2*rank+2; j++) { 62c4762a1bSJed Brown v = -1.0; Ii = j + n*i; 635f80ce2aSJacob Faibussowitsch if (i>0) {J = Ii - n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 645f80ce2aSJacob Faibussowitsch if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 655f80ce2aSJacob Faibussowitsch if (j>0) {J = Ii - 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 665f80ce2aSJacob Faibussowitsch if (j<n+size-1) {J = Ii + 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 675f80ce2aSJacob Faibussowitsch v = 4.0; CHKERRQ(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 68c4762a1bSJed Brown } 69c4762a1bSJed Brown } 705f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 715f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 72c4762a1bSJed Brown 735f80ce2aSJacob Faibussowitsch CHKERRQ(TestMatZeroRows_Basic(A,is,0.0)); 745f80ce2aSJacob Faibussowitsch CHKERRQ(TestMatZeroRows_Basic(A,is,diag)); 75c4762a1bSJed Brown 765f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 775f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&is)); 78*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 79*b122ec5aSJacob Faibussowitsch return 0; 80c4762a1bSJed Brown } 81c4762a1bSJed Brown 82c4762a1bSJed Brown PetscErrorCode TestMatZeroRows_Basic(Mat A,IS is,PetscScalar diag) 83c4762a1bSJed Brown { 84c4762a1bSJed Brown Mat B; 85c4762a1bSJed Brown PetscBool keepnonzeropattern; 86c4762a1bSJed Brown 87c4762a1bSJed Brown /* Now copy A into B, and test it with MatZeroRows() */ 885f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(A,MAT_COPY_VALUES,&B)); 89c4762a1bSJed Brown 905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-keep_nonzero_pattern",&keepnonzeropattern)); 91c4762a1bSJed Brown if (keepnonzeropattern) { 925f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(B,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE)); 93c4762a1bSJed Brown } 94c4762a1bSJed Brown 955f80ce2aSJacob Faibussowitsch CHKERRQ(MatZeroRowsIS(B,is,diag,0,0)); 965f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_WORLD)); 975f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 98c4762a1bSJed Brown return 0; 99c4762a1bSJed Brown } 100c4762a1bSJed Brown 101c4762a1bSJed Brown PetscErrorCode TestMatZeroRows_with_no_allocation(Mat A,IS is,PetscScalar diag) 102c4762a1bSJed Brown { 103c4762a1bSJed Brown Mat B; 104c4762a1bSJed Brown 105c4762a1bSJed Brown /* Now copy A into B, and test it with MatZeroRows() */ 1065f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(A,MAT_COPY_VALUES,&B)); 107c4762a1bSJed Brown /* Set this flag after assembly. This way, it affects only MatZeroRows() */ 1085f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE)); 109c4762a1bSJed Brown 1105f80ce2aSJacob Faibussowitsch CHKERRQ(MatZeroRowsIS(B,is,diag,0,0)); 1115f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_WORLD)); 1125f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 113c4762a1bSJed Brown return 0; 114c4762a1bSJed Brown } 115c4762a1bSJed Brown 116c4762a1bSJed Brown /*TEST 117c4762a1bSJed Brown 118c4762a1bSJed Brown test: 119c4762a1bSJed Brown nsize: 2 120c4762a1bSJed Brown filter: grep -v "MPI processes" 121c4762a1bSJed Brown 122c4762a1bSJed Brown test: 123c4762a1bSJed Brown suffix: 2 124c4762a1bSJed Brown nsize: 3 125c4762a1bSJed Brown args: -mat_type mpibaij -mat_block_size 3 126c4762a1bSJed Brown filter: grep -v "MPI processes" 127c4762a1bSJed Brown 128c4762a1bSJed Brown test: 129c4762a1bSJed Brown suffix: 3 130c4762a1bSJed Brown nsize: 3 131c4762a1bSJed Brown args: -mat_type mpiaij -keep_nonzero_pattern 132c4762a1bSJed Brown filter: grep -v "MPI processes" 133c4762a1bSJed Brown 134c4762a1bSJed Brown test: 135c4762a1bSJed Brown suffix: 4 136c4762a1bSJed Brown nsize: 3 137c4762a1bSJed Brown args: -keep_nonzero_pattern -mat_type mpibaij -mat_block_size 3 138c4762a1bSJed Brown filter: grep -v "MPI processes" 139c4762a1bSJed Brown 140c4762a1bSJed Brown TEST*/ 141