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 10d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 11d71ae5a4SJacob Faibussowitsch { 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 18327415f7SBarry Smith PetscFunctionBeginUser; 199566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 219566063dSJacob Faibussowitsch PetscCallMPI(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*/ 259566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 269566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n)); 279566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 289566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 29c4762a1bSJed Brown for (i = 0; i < m; i++) { 30c4762a1bSJed Brown for (j = 2 * rank; j < 2 * rank + 2; j++) { 319371c9d4SSatish Balay v = -1.0; 329371c9d4SSatish Balay Ii = j + n * i; 339371c9d4SSatish Balay if (i > 0) { 349371c9d4SSatish Balay J = Ii - n; 359371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 369371c9d4SSatish Balay } 379371c9d4SSatish Balay if (i < m - 1) { 389371c9d4SSatish Balay J = Ii + n; 399371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 409371c9d4SSatish Balay } 419371c9d4SSatish Balay if (j > 0) { 429371c9d4SSatish Balay J = Ii - 1; 439371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 449371c9d4SSatish Balay } 459371c9d4SSatish Balay if (j < n - 1) { 469371c9d4SSatish Balay J = Ii + 1; 479371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 489371c9d4SSatish Balay } 499371c9d4SSatish Balay v = 4.0; 509371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &Ii, &v, INSERT_VALUES)); 51c4762a1bSJed Brown } 52c4762a1bSJed Brown } 539566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 549566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 55c4762a1bSJed Brown 56c4762a1bSJed Brown /* Create AN IS required by MatZeroRows() */ 579371c9d4SSatish Balay Imax = n * rank; 589371c9d4SSatish Balay if (Imax >= n * m - m - 1) Imax = m * n - m - 1; 599566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, m, Imax, 1, &is)); 60c4762a1bSJed Brown 619566063dSJacob Faibussowitsch PetscCall(TestMatZeroRows_Basic(A, is, 0.0)); 629566063dSJacob Faibussowitsch PetscCall(TestMatZeroRows_Basic(A, is, diag)); 63c4762a1bSJed Brown 649566063dSJacob Faibussowitsch PetscCall(TestMatZeroRows_with_no_allocation(A, is, 0.0)); 659566063dSJacob Faibussowitsch PetscCall(TestMatZeroRows_with_no_allocation(A, is, diag)); 66c4762a1bSJed Brown 679566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 68c4762a1bSJed Brown 69c4762a1bSJed Brown /* Now Create a rectangular matrix with five point stencil (app) 70c4762a1bSJed Brown n+size is used so that this dimension is always divisible by size. 71c4762a1bSJed Brown This way, we can always use bs = size for any number of procs */ 729566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 739566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, m * n, m * (n + size))); 749566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 759566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 76c4762a1bSJed Brown for (i = 0; i < m; i++) { 77c4762a1bSJed Brown for (j = 2 * rank; j < 2 * rank + 2; j++) { 789371c9d4SSatish Balay v = -1.0; 799371c9d4SSatish Balay Ii = j + n * i; 809371c9d4SSatish Balay if (i > 0) { 819371c9d4SSatish Balay J = Ii - n; 829371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 839371c9d4SSatish Balay } 849371c9d4SSatish Balay if (i < m - 1) { 859371c9d4SSatish Balay J = Ii + n; 869371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 879371c9d4SSatish Balay } 889371c9d4SSatish Balay if (j > 0) { 899371c9d4SSatish Balay J = Ii - 1; 909371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 919371c9d4SSatish Balay } 929371c9d4SSatish Balay if (j < n + size - 1) { 939371c9d4SSatish Balay J = Ii + 1; 949371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 959371c9d4SSatish Balay } 969371c9d4SSatish Balay v = 4.0; 979371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &Ii, &v, INSERT_VALUES)); 98c4762a1bSJed Brown } 99c4762a1bSJed Brown } 1009566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1019566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 102c4762a1bSJed Brown 1039566063dSJacob Faibussowitsch PetscCall(TestMatZeroRows_Basic(A, is, 0.0)); 1049566063dSJacob Faibussowitsch PetscCall(TestMatZeroRows_Basic(A, is, diag)); 105c4762a1bSJed Brown 1069566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1079566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 1089566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 109b122ec5aSJacob Faibussowitsch return 0; 110c4762a1bSJed Brown } 111c4762a1bSJed Brown 112d71ae5a4SJacob Faibussowitsch PetscErrorCode TestMatZeroRows_Basic(Mat A, IS is, PetscScalar diag) 113d71ae5a4SJacob Faibussowitsch { 114c4762a1bSJed Brown Mat B; 115c4762a1bSJed Brown PetscBool keepnonzeropattern; 116c4762a1bSJed Brown 117c4762a1bSJed Brown /* Now copy A into B, and test it with MatZeroRows() */ 118*3ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 1199566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &B)); 120c4762a1bSJed Brown 1219566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-keep_nonzero_pattern", &keepnonzeropattern)); 1221baa6e33SBarry Smith if (keepnonzeropattern) PetscCall(MatSetOption(B, MAT_KEEP_NONZERO_PATTERN, PETSC_TRUE)); 123c4762a1bSJed Brown 1249566063dSJacob Faibussowitsch PetscCall(MatZeroRowsIS(B, is, diag, 0, 0)); 1259566063dSJacob Faibussowitsch PetscCall(MatView(B, PETSC_VIEWER_STDOUT_WORLD)); 1269566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 127*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 128c4762a1bSJed Brown } 129c4762a1bSJed Brown 130d71ae5a4SJacob Faibussowitsch PetscErrorCode TestMatZeroRows_with_no_allocation(Mat A, IS is, PetscScalar diag) 131d71ae5a4SJacob Faibussowitsch { 132c4762a1bSJed Brown Mat B; 133c4762a1bSJed Brown 134c4762a1bSJed Brown /* Now copy A into B, and test it with MatZeroRows() */ 135*3ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 1369566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &B)); 137c4762a1bSJed Brown /* Set this flag after assembly. This way, it affects only MatZeroRows() */ 1389566063dSJacob Faibussowitsch PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE)); 139c4762a1bSJed Brown 1409566063dSJacob Faibussowitsch PetscCall(MatZeroRowsIS(B, is, diag, 0, 0)); 1419566063dSJacob Faibussowitsch PetscCall(MatView(B, PETSC_VIEWER_STDOUT_WORLD)); 1429566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 143*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 144c4762a1bSJed Brown } 145c4762a1bSJed Brown 146c4762a1bSJed Brown /*TEST 147c4762a1bSJed Brown 148c4762a1bSJed Brown test: 149c4762a1bSJed Brown nsize: 2 1508cc725e6SPierre Jolivet filter: grep -v " MPI process" 151c4762a1bSJed Brown 152c4762a1bSJed Brown test: 153c4762a1bSJed Brown suffix: 2 154c4762a1bSJed Brown nsize: 3 155c4762a1bSJed Brown args: -mat_type mpibaij -mat_block_size 3 1568cc725e6SPierre Jolivet filter: grep -v " MPI process" 157c4762a1bSJed Brown 158c4762a1bSJed Brown test: 159c4762a1bSJed Brown suffix: 3 160c4762a1bSJed Brown nsize: 3 161c4762a1bSJed Brown args: -mat_type mpiaij -keep_nonzero_pattern 1628cc725e6SPierre Jolivet filter: grep -v " MPI process" 163c4762a1bSJed Brown 164c4762a1bSJed Brown test: 165c4762a1bSJed Brown suffix: 4 166c4762a1bSJed Brown nsize: 3 167c4762a1bSJed Brown args: -keep_nonzero_pattern -mat_type mpibaij -mat_block_size 3 1688cc725e6SPierre Jolivet filter: grep -v " MPI process" 169c4762a1bSJed Brown 170c4762a1bSJed Brown TEST*/ 171