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