1c4762a1bSJed Brown static char help[] = "Tests the use of MatZeroRows() for parallel matrices.\n\ 2c4762a1bSJed Brown This example also tests the use of MatDuplicate() for both MPIAIJ and MPIBAIJ matrices"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown extern PetscErrorCode TestMatZeroRows_Basic(Mat, IS, PetscScalar); 7c4762a1bSJed Brown extern PetscErrorCode TestMatZeroRows_with_no_allocation(Mat, IS, PetscScalar); 8c4762a1bSJed Brown 9d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 10d71ae5a4SJacob Faibussowitsch { 11c4762a1bSJed Brown Mat A; 12c4762a1bSJed Brown PetscInt i, j, m = 3, n, Ii, J, Imax; 13c4762a1bSJed Brown PetscMPIInt rank, size; 14c4762a1bSJed Brown PetscScalar v, diag = -4.0; 15c4762a1bSJed Brown IS is; 16c4762a1bSJed Brown 17327415f7SBarry Smith PetscFunctionBeginUser; 18*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help)); 199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 209566063dSJacob Faibussowitsch PetscCallMPI(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*/ 249566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 259566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n)); 269566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 279566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 28c4762a1bSJed Brown for (i = 0; i < m; i++) { 29c4762a1bSJed Brown for (j = 2 * rank; j < 2 * rank + 2; j++) { 309371c9d4SSatish Balay v = -1.0; 319371c9d4SSatish Balay Ii = j + n * i; 329371c9d4SSatish Balay if (i > 0) { 339371c9d4SSatish Balay J = Ii - n; 349371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 359371c9d4SSatish Balay } 369371c9d4SSatish Balay if (i < m - 1) { 379371c9d4SSatish Balay J = Ii + n; 389371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 399371c9d4SSatish Balay } 409371c9d4SSatish Balay if (j > 0) { 419371c9d4SSatish Balay J = Ii - 1; 429371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 439371c9d4SSatish Balay } 449371c9d4SSatish Balay if (j < n - 1) { 459371c9d4SSatish Balay J = Ii + 1; 469371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 479371c9d4SSatish Balay } 489371c9d4SSatish Balay v = 4.0; 499371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &Ii, &v, INSERT_VALUES)); 50c4762a1bSJed Brown } 51c4762a1bSJed Brown } 529566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 539566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* Create AN IS required by MatZeroRows() */ 569371c9d4SSatish Balay Imax = n * rank; 579371c9d4SSatish Balay if (Imax >= n * m - m - 1) Imax = m * n - m - 1; 589566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, m, Imax, 1, &is)); 59c4762a1bSJed Brown 609566063dSJacob Faibussowitsch PetscCall(TestMatZeroRows_Basic(A, is, 0.0)); 619566063dSJacob Faibussowitsch PetscCall(TestMatZeroRows_Basic(A, is, diag)); 62c4762a1bSJed Brown 639566063dSJacob Faibussowitsch PetscCall(TestMatZeroRows_with_no_allocation(A, is, 0.0)); 649566063dSJacob Faibussowitsch PetscCall(TestMatZeroRows_with_no_allocation(A, is, diag)); 65c4762a1bSJed Brown 669566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 67c4762a1bSJed Brown 68c4762a1bSJed Brown /* Now Create a rectangular matrix with five point stencil (app) 69c4762a1bSJed Brown n+size is used so that this dimension is always divisible by size. 70c4762a1bSJed Brown This way, we can always use bs = size for any number of procs */ 719566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 729566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, m * n, m * (n + size))); 739566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 749566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 75c4762a1bSJed Brown for (i = 0; i < m; i++) { 76c4762a1bSJed Brown for (j = 2 * rank; j < 2 * rank + 2; j++) { 779371c9d4SSatish Balay v = -1.0; 789371c9d4SSatish Balay Ii = j + n * i; 799371c9d4SSatish Balay if (i > 0) { 809371c9d4SSatish Balay J = Ii - n; 819371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 829371c9d4SSatish Balay } 839371c9d4SSatish Balay if (i < m - 1) { 849371c9d4SSatish Balay J = Ii + n; 859371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 869371c9d4SSatish Balay } 879371c9d4SSatish Balay if (j > 0) { 889371c9d4SSatish Balay J = Ii - 1; 899371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 909371c9d4SSatish Balay } 919371c9d4SSatish Balay if (j < n + size - 1) { 929371c9d4SSatish Balay J = Ii + 1; 939371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 949371c9d4SSatish Balay } 959371c9d4SSatish Balay v = 4.0; 969371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &Ii, &v, INSERT_VALUES)); 97c4762a1bSJed Brown } 98c4762a1bSJed Brown } 999566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1009566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 101c4762a1bSJed Brown 1029566063dSJacob Faibussowitsch PetscCall(TestMatZeroRows_Basic(A, is, 0.0)); 1039566063dSJacob Faibussowitsch PetscCall(TestMatZeroRows_Basic(A, is, diag)); 104c4762a1bSJed Brown 1059566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1069566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 1079566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 108b122ec5aSJacob Faibussowitsch return 0; 109c4762a1bSJed Brown } 110c4762a1bSJed Brown 111d71ae5a4SJacob Faibussowitsch PetscErrorCode TestMatZeroRows_Basic(Mat A, IS is, PetscScalar diag) 112d71ae5a4SJacob Faibussowitsch { 113c4762a1bSJed Brown Mat B; 114c4762a1bSJed Brown PetscBool keepnonzeropattern; 115c4762a1bSJed Brown 116c4762a1bSJed Brown /* Now copy A into B, and test it with MatZeroRows() */ 1173ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 1189566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &B)); 119c4762a1bSJed Brown 1209566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-keep_nonzero_pattern", &keepnonzeropattern)); 1211baa6e33SBarry Smith if (keepnonzeropattern) PetscCall(MatSetOption(B, MAT_KEEP_NONZERO_PATTERN, PETSC_TRUE)); 122c4762a1bSJed Brown 1239566063dSJacob Faibussowitsch PetscCall(MatZeroRowsIS(B, is, diag, 0, 0)); 1249566063dSJacob Faibussowitsch PetscCall(MatView(B, PETSC_VIEWER_STDOUT_WORLD)); 1259566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 1263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 127c4762a1bSJed Brown } 128c4762a1bSJed Brown 129d71ae5a4SJacob Faibussowitsch PetscErrorCode TestMatZeroRows_with_no_allocation(Mat A, IS is, PetscScalar diag) 130d71ae5a4SJacob Faibussowitsch { 131c4762a1bSJed Brown Mat B; 132c4762a1bSJed Brown 133c4762a1bSJed Brown /* Now copy A into B, and test it with MatZeroRows() */ 1343ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 1359566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &B)); 136c4762a1bSJed Brown /* Set this flag after assembly. This way, it affects only MatZeroRows() */ 1379566063dSJacob Faibussowitsch PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE)); 138c4762a1bSJed Brown 1399566063dSJacob Faibussowitsch PetscCall(MatZeroRowsIS(B, is, diag, 0, 0)); 1409566063dSJacob Faibussowitsch PetscCall(MatView(B, PETSC_VIEWER_STDOUT_WORLD)); 1419566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 1423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 143c4762a1bSJed Brown } 144c4762a1bSJed Brown 145c4762a1bSJed Brown /*TEST 146c4762a1bSJed Brown 147c4762a1bSJed Brown test: 148c4762a1bSJed Brown nsize: 2 1498cc725e6SPierre Jolivet filter: grep -v " MPI process" 150c4762a1bSJed Brown 151c4762a1bSJed Brown test: 152c4762a1bSJed Brown suffix: 2 153c4762a1bSJed Brown nsize: 3 154c4762a1bSJed Brown args: -mat_type mpibaij -mat_block_size 3 1558cc725e6SPierre Jolivet filter: grep -v " MPI process" 156c4762a1bSJed Brown 157c4762a1bSJed Brown test: 158c4762a1bSJed Brown suffix: 3 159c4762a1bSJed Brown nsize: 3 160c4762a1bSJed Brown args: -mat_type mpiaij -keep_nonzero_pattern 1618cc725e6SPierre Jolivet filter: grep -v " MPI process" 162c4762a1bSJed Brown 163c4762a1bSJed Brown test: 164c4762a1bSJed Brown suffix: 4 165c4762a1bSJed Brown nsize: 3 166c4762a1bSJed Brown args: -keep_nonzero_pattern -mat_type mpibaij -mat_block_size 3 1678cc725e6SPierre Jolivet filter: grep -v " MPI process" 168c4762a1bSJed Brown 169e9c27b8aSBarry Smith test: 170e9c27b8aSBarry Smith requires: !defined(PETSC_HAVE_THREADSAFETY) 171e9c27b8aSBarry Smith suffix: 5 172e9c27b8aSBarry Smith nsize: 3 173e9c27b8aSBarry Smith args: -mat_type mpibaij -mat_block_size 3 -mat_view 174e9c27b8aSBarry Smith filter: grep -v " MPI process" 175e9c27b8aSBarry Smith 176c4762a1bSJed Brown TEST*/ 177