xref: /petsc/src/mat/tests/ex133.c (revision 40badf4fbc550ac1f60bd080eaff6de6d55b946d)
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   PetscErrorCode ierr;
11   PetscMPIInt    size;
12   PetscScalar    x[6][9];
13 
14   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
15   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
16   PetscCheckFalse(size > 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"Test is only for sequential");
17   CHKERRQ(MatCreateSeqSBAIJ(PETSC_COMM_SELF,bs,m*bs,m*bs,1,NULL,&A));
18   CHKERRQ(MatSetOption(A,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE));
19   rstart = 0;
20 
21   row[0] =rstart+0;  row[1] =rstart+2;
22   col[0] =rstart+0;  col[1] =rstart+1;  col[2] =rstart+3;
23   for (i=0; i<6; i++) {
24     for (j =0; j< 9; j++) x[i][j] = (PetscScalar)val++;
25   }
26   CHKERRQ(MatSetValuesBlocked(A,2,row,3,col,&x[0][0],INSERT_VALUES));
27   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
28   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
29   CHKERRQ(MatView(A,PETSC_VIEWER_BINARY_WORLD));
30 
31   CHKERRQ(MatDestroy(&A));
32   ierr = PetscFinalize();
33   return ierr;
34 }
35