1 static char help[] = "Tests the use of MatZeroRows() for uniprocessor matrices.\n\n"; 2 3 #include <petscmat.h> 4 5 int main(int argc, char **args) 6 { 7 Mat C; 8 PetscInt i, j, m = 5, n = 5, Ii, J; 9 PetscScalar v, five = 5.0; 10 IS isrow; 11 PetscBool keepnonzeropattern; 12 13 PetscFunctionBeginUser; 14 PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 15 /* create the matrix for the five point stencil, YET AGAIN*/ 16 PetscCall(MatCreate(PETSC_COMM_SELF, &C)); 17 PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n)); 18 PetscCall(MatSetFromOptions(C)); 19 PetscCall(MatSetUp(C)); 20 for (i = 0; i < m; i++) { 21 for (j = 0; j < n; j++) { 22 v = -1.0; 23 Ii = j + n * i; 24 if (i > 0) { 25 J = Ii - n; 26 PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 27 } 28 if (i < m - 1) { 29 J = Ii + n; 30 PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 31 } 32 if (j > 0) { 33 J = Ii - 1; 34 PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 35 } 36 if (j < n - 1) { 37 J = Ii + 1; 38 PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 39 } 40 v = 4.0; 41 PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES)); 42 } 43 } 44 PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 45 PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 46 47 PetscCall(ISCreateStride(PETSC_COMM_SELF, (m * n) / 2, 0, 2, &isrow)); 48 49 PetscCall(PetscOptionsHasName(NULL, NULL, "-keep_nonzero_pattern", &keepnonzeropattern)); 50 if (keepnonzeropattern) PetscCall(MatSetOption(C, MAT_KEEP_NONZERO_PATTERN, PETSC_TRUE)); 51 52 PetscCall(MatZeroRowsIS(C, isrow, five, 0, 0)); 53 54 PetscCall(MatView(C, PETSC_VIEWER_STDOUT_SELF)); 55 56 PetscCall(ISDestroy(&isrow)); 57 PetscCall(MatDestroy(&C)); 58 PetscCall(PetscFinalize()); 59 return 0; 60 } 61 62 /*TEST 63 64 test: 65 66 test: 67 suffix: 2 68 args: -mat_type seqbaij -mat_block_size 5 69 70 test: 71 suffix: 3 72 args: -keep_nonzero_pattern 73 74 test: 75 suffix: 4 76 args: -keep_nonzero_pattern -mat_type seqbaij -mat_block_size 5 77 78 TEST*/ 79