xref: /petsc/src/mat/tests/ex133.c (revision 98d129c30f3ee9fdddc40fdbc5a989b7be64f888)
1 static char help[] = "Test saving SeqSBAIJ matrix that is missing diagonal entries.";
2 
3 #include <petscmat.h>
4 
5 int main(int argc, char **args)
6 {
7   Mat         A;
8   PetscInt    bs = 3, m = 4, i, j, val = 10, row[2], col[3], rstart;
9   PetscMPIInt size;
10   PetscScalar x[6][9];
11 
12   PetscFunctionBeginUser;
13   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
14   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
15   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Test is only for sequential");
16   PetscCall(MatCreateSeqSBAIJ(PETSC_COMM_SELF, bs, m * bs, m * bs, 1, NULL, &A));
17   PetscCall(MatSetOption(A, MAT_IGNORE_LOWER_TRIANGULAR, PETSC_TRUE));
18   rstart = 0;
19 
20   row[0] = rstart + 0;
21   row[1] = rstart + 2;
22   col[0] = rstart + 0;
23   col[1] = rstart + 1;
24   col[2] = rstart + 3;
25   for (i = 0; i < 6; i++) {
26     for (j = 0; j < 9; j++) x[i][j] = (PetscScalar)val++;
27   }
28   PetscCall(MatSetValuesBlocked(A, 2, row, 3, col, &x[0][0], INSERT_VALUES));
29   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
30   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
31   PetscCall(MatView(A, PETSC_VIEWER_BINARY_WORLD));
32 
33   PetscCall(MatDestroy(&A));
34   PetscCall(PetscFinalize());
35   return 0;
36 }
37