xref: /petsc/src/mat/tests/ex141.c (revision ebead697dbf761eb322f829370bbe90b3bd93fa3)
1 
2 static char help[] = "Tests converting a SBAIJ matrix to BAIJ format with MatConvert. Modified from ex55.c\n\n";
3 
4 #include <petscmat.h>
5 /* Usage: ./ex141 -mat_view */
6 
7 int main(int argc,char **args)
8 {
9   Mat            C,B;
10   PetscInt       i,bs=2,mbs,m,block,d_nz=6,col[3];
11   PetscMPIInt    size;
12   char           file[PETSC_MAX_PATH_LEN];
13   PetscViewer    fd;
14   PetscBool      equal,loadmat;
15   PetscScalar    value[4];
16   PetscInt       ridx[2],cidx[2];
17 
18   PetscFunctionBeginUser;
19   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
20   PetscCall(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&loadmat));
21   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
22   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
23 
24   /* input matrix C */
25   if (loadmat) {
26     /* Open binary file. Load a sbaij matrix, then destroy the viewer. */
27     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd));
28     PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
29     PetscCall(MatSetType(C,MATSEQSBAIJ));
30     PetscCall(MatLoad(C,fd));
31     PetscCall(PetscViewerDestroy(&fd));
32   } else { /* Create a sbaij mat with bs>1  */
33     mbs  =8;
34     PetscCall(PetscOptionsGetInt(NULL,NULL,"-mbs",&mbs,NULL));
35     m    = mbs*bs;
36     PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
37     PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m,m));
38     PetscCall(MatSetType(C,MATSBAIJ));
39     PetscCall(MatSetFromOptions(C));
40     PetscCall(MatSeqSBAIJSetPreallocation(C,bs,d_nz,NULL));
41     PetscCall(MatSetUp(C));
42     PetscCall(MatSetOption(C,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE));
43 
44     for (block=0; block<mbs; block++) {
45       /* diagonal blocks */
46       value[0] = -1.0; value[1] = 4.0; value[2] = -1.0;
47       for (i=1+block*bs; i<bs-1+block*bs; i++) {
48         col[0] = i-1; col[1] = i; col[2] = i+1;
49         PetscCall(MatSetValues(C,1,&i,3,col,value,INSERT_VALUES));
50       }
51       i = bs - 1+block*bs; col[0] = bs - 2+block*bs; col[1] = bs - 1+block*bs;
52 
53       value[0]=-1.0; value[1]=4.0;
54       PetscCall(MatSetValues(C,1,&i,2,col,value,INSERT_VALUES));
55 
56       i = 0+block*bs; col[0] = 0+block*bs; col[1] = 1+block*bs;
57 
58       value[0]=4.0; value[1] = -1.0;
59       PetscCall(MatSetValues(C,1,&i,2,col,value,INSERT_VALUES));
60     }
61     /* off-diagonal blocks */
62     value[0]=-1.0; value[1] = -0.1; value[2] = 0.0; value[3] = -1.0; /* row-oriented */
63     for (block=0; block<mbs-1; block++) {
64       for (i=0; i<bs; i++) {
65         ridx[i] = block*bs+i; cidx[i] = (block+1)*bs+i;
66       }
67       PetscCall(MatSetValues(C,bs,ridx,bs,cidx,value,INSERT_VALUES));
68     }
69     PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
70     PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
71   }
72 
73   /* convert C to BAIJ format */
74   PetscCall(MatConvert(C,MATSEQBAIJ,MAT_INITIAL_MATRIX,&B));
75   PetscCall(MatMultEqual(B,C,10,&equal));
76   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatConvert fails!");
77 
78   PetscCall(MatDestroy(&B));
79   PetscCall(MatDestroy(&C));
80   PetscCall(PetscFinalize());
81   return 0;
82 }
83 
84 /*TEST
85 
86    test:
87      output_file: output/ex141.out
88 
89 TEST*/
90