Actual source code: ex141.c

  2: static char help[] = "Tests converting a SBAIJ matrix to BAIJ format with MatCovert. Modified from ex55.c\n\n";

  4: #include <petscmat.h>
  5: /* Usage: ./ex141 -mat_view */

  9: int main(int argc,char **args)
 10: {
 11:   Mat            C,B;
 13:   PetscInt       i,bs=2,mbs,m,block,d_nz=6,col[3];
 14:   PetscMPIInt    size;
 15:   char           file[PETSC_MAX_PATH_LEN];
 16:   PetscViewer    fd;
 17:   PetscBool      equal,loadmat;
 18:   PetscScalar    value[4];
 19:   PetscInt       ridx[2],cidx[2];

 21:   PetscInitialize(&argc,&args,(char *)0,help);
 22:   PetscOptionsGetString(PETSC_NULL,"-f",file,PETSC_MAX_PATH_LEN,&loadmat);
 23:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 24:    if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!");

 26:   /* input matrix C */
 27:   if (loadmat){
 28:     /* Open binary file. Load a sbaij matrix, then destroy the viewer. */
 29:     PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);
 30:     MatCreate(PETSC_COMM_WORLD,&C);
 31:     MatSetType(C,MATSEQSBAIJ);
 32:     MatLoad(C,fd);
 33:     PetscViewerDestroy(&fd);
 34:   } else { /* Create a sbaij mat with bs>1  */
 35:     mbs=8;
 36:     PetscOptionsGetInt(PETSC_NULL,"-mbs",&mbs,PETSC_NULL);
 37:     m = mbs*bs;
 38:     MatCreate(PETSC_COMM_WORLD,&C);
 39:     MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m,m);
 40:     MatSetType(C,MATSBAIJ);
 41:     MatSetOption(C,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);
 42:     MatSetFromOptions(C);
 43:     MatSeqSBAIJSetPreallocation(C,bs,d_nz,PETSC_NULL);
 44: 
 45:     for (block=0; block<mbs; block++){
 46:       /* diagonal blocks */
 47:       value[0] = -1.0; value[1] = 4.0; value[2] = -1.0;
 48:       for (i=1+block*bs; i<bs-1+block*bs; i++) {
 49:         col[0] = i-1; col[1] = i; col[2] = i+1;
 50:         MatSetValues(C,1,&i,3,col,value,INSERT_VALUES);
 51:       }
 52:       i = bs - 1+block*bs; col[0] = bs - 2+block*bs; col[1] = bs - 1+block*bs;
 53:       value[0]=-1.0; value[1]=4.0;
 54:       MatSetValues(C,1,&i,2,col,value,INSERT_VALUES);

 56:       i = 0+block*bs; col[0] = 0+block*bs; col[1] = 1+block*bs;
 57:       value[0]=4.0; value[1] = -1.0;
 58:       MatSetValues(C,1,&i,2,col,value,INSERT_VALUES);
 59:     }
 60:     /* off-diagonal blocks */
 61:     value[0]=-1.0; value[1] = -0.1; value[2] = 0.0; value[3] = -1.0; /* row-oriented */
 62:     for (block=0; block<mbs-1; block++){
 63:       for (i=0; i<bs; i++){
 64:         ridx[i] = block*bs+i; cidx[i] = (block+1)*bs+i;
 65:       }
 66:       MatSetValues(C,bs,ridx,bs,cidx,value,INSERT_VALUES);
 67:     }
 68:     MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 69:     MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
 70:   }
 71: 
 72:   /* convert C to BAIJ format */
 73:   MatConvert(C,MATSEQBAIJ,MAT_INITIAL_MATRIX,&B);
 74:   MatMultEqual(B,C,10,&equal);
 75:   if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatConvert fails!");

 77:   MatDestroy(&B);
 78:   MatDestroy(&C);
 79:   PetscFinalize();
 80:   return 0;
 81: }