xref: /petsc/src/mat/tests/ex133.c (revision daa037dfd3c3bec8dc8659548d2b20b07c1dc6de)
1 
2 static char help[] = "Test saving SeqSBAIJ matrix that is missing diagonal entries.";
3 
4 #include <petscmat.h>
5 
6 int main(int argc,char **args)
7 {
8   Mat            A;
9   PetscInt       bs=3,m=4,i,j,val = 10,row[2],col[3],rstart;
10   PetscMPIInt    size;
11   PetscScalar    x[6][9];
12 
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;  row[1] =rstart+2;
21   col[0] =rstart+0;  col[1] =rstart+1;  col[2] =rstart+3;
22   for (i=0; i<6; i++) {
23     for (j =0; j< 9; j++) x[i][j] = (PetscScalar)val++;
24   }
25   PetscCall(MatSetValuesBlocked(A,2,row,3,col,&x[0][0],INSERT_VALUES));
26   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
27   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
28   PetscCall(MatView(A,PETSC_VIEWER_BINARY_WORLD));
29 
30   PetscCall(MatDestroy(&A));
31   PetscCall(PetscFinalize());
32   return 0;
33 }
34