Actual source code: ex191.c

petsc-3.6.1 2015-08-06
Report Typos and Errors
  1: static char help[] = "Tests MatLoad() for dense matrix with uneven dimensions set in program\n\n";

  3: #include <petscmat.h>

  7: int main(int argc,char **args)
  8: {
  9:   Mat            A;
 10:   PetscViewer    fd;
 12:   PetscMPIInt    rank;
 13:   PetscScalar    *Av;
 14:   PetscInt       i;

 16:   PetscInitialize(&argc,&args,(char*)0,help);
 17:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);

 19:   MatCreateDense(PETSC_COMM_WORLD,6,6,12,12,NULL,&A);
 20:   MatDenseGetArray(A,&Av);
 21:   for (i=0; i<6*12; i++) Av[i] = (PetscScalar) i;
 22:   MatDenseRestoreArray(A,&Av);

 24:   /* Load matrices */
 25:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"ex191matrix",FILE_MODE_WRITE,&fd);
 26:   PetscViewerSetFormat(fd,PETSC_VIEWER_NATIVE);
 27:   MatView(A,fd);
 28:   MatDestroy(&A);
 29:   PetscViewerDestroy(&fd);

 31:   MatCreate(PETSC_COMM_WORLD,&A);
 32:   MatSetType(A,MATDENSE);
 33:   if (!rank) {
 34:     MatSetSizes(A, 4, PETSC_DETERMINE, PETSC_DETERMINE,PETSC_DETERMINE);
 35:   } else {
 36:     MatSetSizes(A, 8, PETSC_DETERMINE, PETSC_DETERMINE,PETSC_DETERMINE);
 37:   }
 38:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"ex191matrix",FILE_MODE_READ,&fd);
 39:   MatLoad(A,fd);
 40:   PetscViewerDestroy(&fd);
 41:   MatView(A,PETSC_VIEWER_STDOUT_WORLD);
 42:   PetscViewerSetFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
 43:   MatView(A,PETSC_VIEWER_STDOUT_WORLD);
 44:   MatDestroy(&A);
 45:   PetscFinalize();
 46:   return 0;
 47: }