1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests the use of MatZeroRows() for uniprocessor matrices.\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6*d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 7*d71ae5a4SJacob Faibussowitsch { 8c4762a1bSJed Brown Mat C; 9c4762a1bSJed Brown PetscInt i, j, m = 5, n = 5, Ii, J; 10c4762a1bSJed Brown PetscScalar v, five = 5.0; 11c4762a1bSJed Brown IS isrow; 12c4762a1bSJed Brown PetscBool keepnonzeropattern; 13c4762a1bSJed Brown 14327415f7SBarry Smith PetscFunctionBeginUser; 159566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 16c4762a1bSJed Brown /* create the matrix for the five point stencil, YET AGAIN*/ 179566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &C)); 189566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n)); 199566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 209566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 21c4762a1bSJed Brown for (i = 0; i < m; i++) { 22c4762a1bSJed Brown for (j = 0; j < n; j++) { 239371c9d4SSatish Balay v = -1.0; 249371c9d4SSatish Balay Ii = j + n * i; 259371c9d4SSatish Balay if (i > 0) { 269371c9d4SSatish Balay J = Ii - n; 279371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 289371c9d4SSatish Balay } 299371c9d4SSatish Balay if (i < m - 1) { 309371c9d4SSatish Balay J = Ii + n; 319371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 329371c9d4SSatish Balay } 339371c9d4SSatish Balay if (j > 0) { 349371c9d4SSatish Balay J = Ii - 1; 359371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 369371c9d4SSatish Balay } 379371c9d4SSatish Balay if (j < n - 1) { 389371c9d4SSatish Balay J = Ii + 1; 399371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 409371c9d4SSatish Balay } 419371c9d4SSatish Balay v = 4.0; 429371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES)); 43c4762a1bSJed Brown } 44c4762a1bSJed Brown } 459566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 469566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 47c4762a1bSJed Brown 489566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, (m * n) / 2, 0, 2, &isrow)); 49c4762a1bSJed Brown 509566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-keep_nonzero_pattern", &keepnonzeropattern)); 511baa6e33SBarry Smith if (keepnonzeropattern) PetscCall(MatSetOption(C, MAT_KEEP_NONZERO_PATTERN, PETSC_TRUE)); 52c4762a1bSJed Brown 539566063dSJacob Faibussowitsch PetscCall(MatZeroRowsIS(C, isrow, five, 0, 0)); 54c4762a1bSJed Brown 559566063dSJacob Faibussowitsch PetscCall(MatView(C, PETSC_VIEWER_STDOUT_SELF)); 56c4762a1bSJed Brown 579566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrow)); 589566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 599566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 60b122ec5aSJacob Faibussowitsch return 0; 61c4762a1bSJed Brown } 62c4762a1bSJed Brown 63c4762a1bSJed Brown /*TEST 64c4762a1bSJed Brown 65c4762a1bSJed Brown test: 66c4762a1bSJed Brown 67c4762a1bSJed Brown test: 68c4762a1bSJed Brown suffix: 2 69c4762a1bSJed Brown args: -mat_type seqbaij -mat_block_size 5 70c4762a1bSJed Brown 71c4762a1bSJed Brown test: 72c4762a1bSJed Brown suffix: 3 73c4762a1bSJed Brown args: -keep_nonzero_pattern 74c4762a1bSJed Brown 75c4762a1bSJed Brown test: 76c4762a1bSJed Brown suffix: 4 77c4762a1bSJed Brown args: -keep_nonzero_pattern -mat_type seqbaij -mat_block_size 5 78c4762a1bSJed Brown 79c4762a1bSJed Brown TEST*/ 80