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 PetscErrorCode ierr; 16c4762a1bSJed Brown PetscScalar v,diag=-4.0; 17c4762a1bSJed Brown IS is; 18c4762a1bSJed Brown 19c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 20*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 21*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 22c4762a1bSJed Brown n = 2*size; 23c4762a1bSJed Brown 24c4762a1bSJed Brown /* create A Square matrix for the five point stencil,YET AGAIN*/ 25*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 26*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n)); 27*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 28*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A)); 29c4762a1bSJed Brown for (i=0; i<m; i++) { 30c4762a1bSJed Brown for (j=2*rank; j<2*rank+2; j++) { 31c4762a1bSJed Brown v = -1.0; Ii = j + n*i; 32*5f80ce2aSJacob Faibussowitsch if (i>0) {J = Ii - n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 33*5f80ce2aSJacob Faibussowitsch if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 34*5f80ce2aSJacob Faibussowitsch if (j>0) {J = Ii - 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 35*5f80ce2aSJacob Faibussowitsch if (j<n-1) {J = Ii + 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 36*5f80ce2aSJacob Faibussowitsch v = 4.0; CHKERRQ(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 37c4762a1bSJed Brown } 38c4762a1bSJed Brown } 39*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 40*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 41c4762a1bSJed Brown 42c4762a1bSJed Brown /* Create AN IS required by MatZeroRows() */ 43c4762a1bSJed Brown Imax = n*rank; if (Imax>= n*m -m - 1) Imax = m*n - m - 1; 44*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(PETSC_COMM_SELF,m,Imax,1,&is)); 45c4762a1bSJed Brown 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(TestMatZeroRows_Basic(A,is,0.0)); 47*5f80ce2aSJacob Faibussowitsch CHKERRQ(TestMatZeroRows_Basic(A,is,diag)); 48c4762a1bSJed Brown 49*5f80ce2aSJacob Faibussowitsch CHKERRQ(TestMatZeroRows_with_no_allocation(A,is,0.0)); 50*5f80ce2aSJacob Faibussowitsch CHKERRQ(TestMatZeroRows_with_no_allocation(A,is,diag)); 51c4762a1bSJed Brown 52*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 53c4762a1bSJed Brown 54c4762a1bSJed Brown /* Now Create a rectangular matrix with five point stencil (app) 55c4762a1bSJed Brown n+size is used so that this dimension is always divisible by size. 56c4762a1bSJed Brown This way, we can always use bs = size for any number of procs */ 57*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 58*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*(n+size))); 59*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 60*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A)); 61c4762a1bSJed Brown for (i=0; i<m; i++) { 62c4762a1bSJed Brown for (j=2*rank; j<2*rank+2; j++) { 63c4762a1bSJed Brown v = -1.0; Ii = j + n*i; 64*5f80ce2aSJacob Faibussowitsch if (i>0) {J = Ii - n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 65*5f80ce2aSJacob Faibussowitsch if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 66*5f80ce2aSJacob Faibussowitsch if (j>0) {J = Ii - 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 67*5f80ce2aSJacob Faibussowitsch if (j<n+size-1) {J = Ii + 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 68*5f80ce2aSJacob Faibussowitsch v = 4.0; CHKERRQ(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 69c4762a1bSJed Brown } 70c4762a1bSJed Brown } 71*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 72*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 73c4762a1bSJed Brown 74*5f80ce2aSJacob Faibussowitsch CHKERRQ(TestMatZeroRows_Basic(A,is,0.0)); 75*5f80ce2aSJacob Faibussowitsch CHKERRQ(TestMatZeroRows_Basic(A,is,diag)); 76c4762a1bSJed Brown 77*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 78*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&is)); 79c4762a1bSJed Brown ierr = PetscFinalize(); 80c4762a1bSJed Brown return ierr; 81c4762a1bSJed Brown } 82c4762a1bSJed Brown 83c4762a1bSJed Brown PetscErrorCode TestMatZeroRows_Basic(Mat A,IS is,PetscScalar diag) 84c4762a1bSJed Brown { 85c4762a1bSJed Brown Mat B; 86c4762a1bSJed Brown PetscBool keepnonzeropattern; 87c4762a1bSJed Brown 88c4762a1bSJed Brown /* Now copy A into B, and test it with MatZeroRows() */ 89*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(A,MAT_COPY_VALUES,&B)); 90c4762a1bSJed Brown 91*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-keep_nonzero_pattern",&keepnonzeropattern)); 92c4762a1bSJed Brown if (keepnonzeropattern) { 93*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(B,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE)); 94c4762a1bSJed Brown } 95c4762a1bSJed Brown 96*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatZeroRowsIS(B,is,diag,0,0)); 97*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_WORLD)); 98*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 99c4762a1bSJed Brown return 0; 100c4762a1bSJed Brown } 101c4762a1bSJed Brown 102c4762a1bSJed Brown PetscErrorCode TestMatZeroRows_with_no_allocation(Mat A,IS is,PetscScalar diag) 103c4762a1bSJed Brown { 104c4762a1bSJed Brown Mat B; 105c4762a1bSJed Brown 106c4762a1bSJed Brown /* Now copy A into B, and test it with MatZeroRows() */ 107*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(A,MAT_COPY_VALUES,&B)); 108c4762a1bSJed Brown /* Set this flag after assembly. This way, it affects only MatZeroRows() */ 109*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE)); 110c4762a1bSJed Brown 111*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatZeroRowsIS(B,is,diag,0,0)); 112*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_WORLD)); 113*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 114c4762a1bSJed Brown return 0; 115c4762a1bSJed Brown } 116c4762a1bSJed Brown 117c4762a1bSJed Brown /*TEST 118c4762a1bSJed Brown 119c4762a1bSJed Brown test: 120c4762a1bSJed Brown nsize: 2 121c4762a1bSJed Brown filter: grep -v "MPI processes" 122c4762a1bSJed Brown 123c4762a1bSJed Brown test: 124c4762a1bSJed Brown suffix: 2 125c4762a1bSJed Brown nsize: 3 126c4762a1bSJed Brown args: -mat_type mpibaij -mat_block_size 3 127c4762a1bSJed Brown filter: grep -v "MPI processes" 128c4762a1bSJed Brown 129c4762a1bSJed Brown test: 130c4762a1bSJed Brown suffix: 3 131c4762a1bSJed Brown nsize: 3 132c4762a1bSJed Brown args: -mat_type mpiaij -keep_nonzero_pattern 133c4762a1bSJed Brown filter: grep -v "MPI processes" 134c4762a1bSJed Brown 135c4762a1bSJed Brown test: 136c4762a1bSJed Brown suffix: 4 137c4762a1bSJed Brown nsize: 3 138c4762a1bSJed Brown args: -keep_nonzero_pattern -mat_type mpibaij -mat_block_size 3 139c4762a1bSJed Brown filter: grep -v "MPI processes" 140c4762a1bSJed Brown 141c4762a1bSJed Brown TEST*/ 142