Actual source code: ex141.c

petsc-3.10.5 2019-03-28
Report Typos and Errors

  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 */

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

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

 24:   /* input matrix C */
 25:   if (loadmat) {
 26:     /* Open binary file. Load a sbaij matrix, then destroy the viewer. */
 27:     PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);
 28:     MatCreate(PETSC_COMM_WORLD,&C);
 29:     MatSetType(C,MATSEQSBAIJ);
 30:     MatLoad(C,fd);
 31:     PetscViewerDestroy(&fd);
 32:   } else { /* Create a sbaij mat with bs>1  */
 33:     mbs  =8;
 34:     PetscOptionsGetInt(NULL,NULL,"-mbs",&mbs,NULL);
 35:     m    = mbs*bs;
 36:     MatCreate(PETSC_COMM_WORLD,&C);
 37:     MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m,m);
 38:     MatSetType(C,MATSBAIJ);
 39:     MatSetFromOptions(C);
 40:     MatSeqSBAIJSetPreallocation(C,bs,d_nz,NULL);
 41:     MatSetUp(C);
 42:     MatSetOption(C,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);

 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:         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;

 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;

 58:       value[0]=4.0; value[1] = -1.0;
 59:       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:       MatSetValues(C,bs,ridx,bs,cidx,value,INSERT_VALUES);
 68:     }
 69:     MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 70:     MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
 71:   }

 73:   /* convert C to BAIJ format */
 74:   MatConvert(C,MATSEQBAIJ,MAT_INITIAL_MATRIX,&B);
 75:   MatMultEqual(B,C,10,&equal);
 76:   if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatConvert fails!");

 78:   MatDestroy(&B);
 79:   MatDestroy(&C);
 80:   PetscFinalize();
 81:   return ierr;
 82: }

 84: /*TEST

 86:    test:
 87:      output_file: output/ex141.out

 89: TEST*/