xref: /petsc/src/mat/tests/ex141.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests converting a SBAIJ matrix to BAIJ format with MatConvert. Modified from ex55.c\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown /* Usage: ./ex141 -mat_view */
6c4762a1bSJed Brown 
7c4762a1bSJed Brown int main(int argc,char **args)
8c4762a1bSJed Brown {
9c4762a1bSJed Brown   Mat            C,B;
10c4762a1bSJed Brown   PetscInt       i,bs=2,mbs,m,block,d_nz=6,col[3];
11c4762a1bSJed Brown   PetscMPIInt    size;
12c4762a1bSJed Brown   char           file[PETSC_MAX_PATH_LEN];
13c4762a1bSJed Brown   PetscViewer    fd;
14c4762a1bSJed Brown   PetscBool      equal,loadmat;
15c4762a1bSJed Brown   PetscScalar    value[4];
16c4762a1bSJed Brown   PetscInt       ridx[2],cidx[2];
17c4762a1bSJed Brown 
18*327415f7SBarry Smith   PetscFunctionBeginUser;
199566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&loadmat));
219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
22be096a46SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
23c4762a1bSJed Brown 
24c4762a1bSJed Brown   /* input matrix C */
25c4762a1bSJed Brown   if (loadmat) {
26c4762a1bSJed Brown     /* Open binary file. Load a sbaij matrix, then destroy the viewer. */
279566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd));
289566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
299566063dSJacob Faibussowitsch     PetscCall(MatSetType(C,MATSEQSBAIJ));
309566063dSJacob Faibussowitsch     PetscCall(MatLoad(C,fd));
319566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&fd));
32c4762a1bSJed Brown   } else { /* Create a sbaij mat with bs>1  */
33c4762a1bSJed Brown     mbs  =8;
349566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetInt(NULL,NULL,"-mbs",&mbs,NULL));
35c4762a1bSJed Brown     m    = mbs*bs;
369566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
379566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m,m));
389566063dSJacob Faibussowitsch     PetscCall(MatSetType(C,MATSBAIJ));
399566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(C));
409566063dSJacob Faibussowitsch     PetscCall(MatSeqSBAIJSetPreallocation(C,bs,d_nz,NULL));
419566063dSJacob Faibussowitsch     PetscCall(MatSetUp(C));
429566063dSJacob Faibussowitsch     PetscCall(MatSetOption(C,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE));
43c4762a1bSJed Brown 
44c4762a1bSJed Brown     for (block=0; block<mbs; block++) {
45c4762a1bSJed Brown       /* diagonal blocks */
46c4762a1bSJed Brown       value[0] = -1.0; value[1] = 4.0; value[2] = -1.0;
47c4762a1bSJed Brown       for (i=1+block*bs; i<bs-1+block*bs; i++) {
48c4762a1bSJed Brown         col[0] = i-1; col[1] = i; col[2] = i+1;
499566063dSJacob Faibussowitsch         PetscCall(MatSetValues(C,1,&i,3,col,value,INSERT_VALUES));
50c4762a1bSJed Brown       }
51c4762a1bSJed Brown       i = bs - 1+block*bs; col[0] = bs - 2+block*bs; col[1] = bs - 1+block*bs;
52c4762a1bSJed Brown 
53c4762a1bSJed Brown       value[0]=-1.0; value[1]=4.0;
549566063dSJacob Faibussowitsch       PetscCall(MatSetValues(C,1,&i,2,col,value,INSERT_VALUES));
55c4762a1bSJed Brown 
56c4762a1bSJed Brown       i = 0+block*bs; col[0] = 0+block*bs; col[1] = 1+block*bs;
57c4762a1bSJed Brown 
58c4762a1bSJed Brown       value[0]=4.0; value[1] = -1.0;
599566063dSJacob Faibussowitsch       PetscCall(MatSetValues(C,1,&i,2,col,value,INSERT_VALUES));
60c4762a1bSJed Brown     }
61c4762a1bSJed Brown     /* off-diagonal blocks */
62c4762a1bSJed Brown     value[0]=-1.0; value[1] = -0.1; value[2] = 0.0; value[3] = -1.0; /* row-oriented */
63c4762a1bSJed Brown     for (block=0; block<mbs-1; block++) {
64c4762a1bSJed Brown       for (i=0; i<bs; i++) {
65c4762a1bSJed Brown         ridx[i] = block*bs+i; cidx[i] = (block+1)*bs+i;
66c4762a1bSJed Brown       }
679566063dSJacob Faibussowitsch       PetscCall(MatSetValues(C,bs,ridx,bs,cidx,value,INSERT_VALUES));
68c4762a1bSJed Brown     }
699566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
709566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
71c4762a1bSJed Brown   }
72c4762a1bSJed Brown 
73c4762a1bSJed Brown   /* convert C to BAIJ format */
749566063dSJacob Faibussowitsch   PetscCall(MatConvert(C,MATSEQBAIJ,MAT_INITIAL_MATRIX,&B));
759566063dSJacob Faibussowitsch   PetscCall(MatMultEqual(B,C,10,&equal));
7628b400f6SJacob Faibussowitsch   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatConvert fails!");
77c4762a1bSJed Brown 
789566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
799566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
809566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
81b122ec5aSJacob Faibussowitsch   return 0;
82c4762a1bSJed Brown }
83c4762a1bSJed Brown 
84c4762a1bSJed Brown /*TEST
85c4762a1bSJed Brown 
86c4762a1bSJed Brown    test:
87c4762a1bSJed Brown      output_file: output/ex141.out
88c4762a1bSJed Brown 
89c4762a1bSJed Brown TEST*/
90