xref: /petsc/src/mat/tests/ex141.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
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   PetscErrorCode ierr;
11c4762a1bSJed Brown   PetscInt       i,bs=2,mbs,m,block,d_nz=6,col[3];
12c4762a1bSJed Brown   PetscMPIInt    size;
13c4762a1bSJed Brown   char           file[PETSC_MAX_PATH_LEN];
14c4762a1bSJed Brown   PetscViewer    fd;
15c4762a1bSJed Brown   PetscBool      equal,loadmat;
16c4762a1bSJed Brown   PetscScalar    value[4];
17c4762a1bSJed Brown   PetscInt       ridx[2],cidx[2];
18c4762a1bSJed Brown 
19c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
205f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&loadmat));
215f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
222c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"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. */
275f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd));
285f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C));
295f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetType(C,MATSEQSBAIJ));
305f80ce2aSJacob Faibussowitsch     CHKERRQ(MatLoad(C,fd));
315f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerDestroy(&fd));
32c4762a1bSJed Brown   } else { /* Create a sbaij mat with bs>1  */
33c4762a1bSJed Brown     mbs  =8;
345f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mbs",&mbs,NULL));
35c4762a1bSJed Brown     m    = mbs*bs;
365f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C));
375f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m,m));
385f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetType(C,MATSBAIJ));
395f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetFromOptions(C));
405f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSeqSBAIJSetPreallocation(C,bs,d_nz,NULL));
415f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetUp(C));
425f80ce2aSJacob Faibussowitsch     CHKERRQ(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;
495f80ce2aSJacob Faibussowitsch         CHKERRQ(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;
545f80ce2aSJacob Faibussowitsch       CHKERRQ(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;
595f80ce2aSJacob Faibussowitsch       CHKERRQ(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       }
675f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(C,bs,ridx,bs,cidx,value,INSERT_VALUES));
68c4762a1bSJed Brown     }
695f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
705f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
71c4762a1bSJed Brown   }
72c4762a1bSJed Brown 
73c4762a1bSJed Brown   /* convert C to BAIJ format */
745f80ce2aSJacob Faibussowitsch   CHKERRQ(MatConvert(C,MATSEQBAIJ,MAT_INITIAL_MATRIX,&B));
755f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultEqual(B,C,10,&equal));
76*28b400f6SJacob Faibussowitsch   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatConvert fails!");
77c4762a1bSJed Brown 
785f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
80c4762a1bSJed Brown   ierr = PetscFinalize();
81c4762a1bSJed Brown   return ierr;
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