xref: /petsc/src/mat/tests/ex191.c (revision 40badf4fbc550ac1f60bd080eaff6de6d55b946d)
1 static char help[] = "Tests MatLoad() for dense matrix with uneven dimensions set in program\n\n";
2 
3 #include <petscmat.h>
4 
5 int main(int argc,char **args)
6 {
7   Mat            A;
8   PetscViewer    fd;
9   PetscErrorCode ierr;
10   PetscMPIInt    rank;
11   PetscScalar    *Av;
12   PetscInt       i;
13 
14   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
15   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
16 
17   CHKERRQ(MatCreateDense(PETSC_COMM_WORLD,6,6,12,12,NULL,&A));
18   CHKERRQ(MatDenseGetArray(A,&Av));
19   for (i=0; i<6*12; i++) Av[i] = (PetscScalar) i;
20   CHKERRQ(MatDenseRestoreArray(A,&Av));
21 
22   /* Load matrices */
23   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"ex191matrix",FILE_MODE_WRITE,&fd));
24   CHKERRQ(PetscViewerPushFormat(fd,PETSC_VIEWER_NATIVE));
25   CHKERRQ(MatView(A,fd));
26   CHKERRQ(MatDestroy(&A));
27   CHKERRQ(PetscViewerPopFormat(fd));
28   CHKERRQ(PetscViewerDestroy(&fd));
29 
30   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
31   CHKERRQ(MatSetType(A,MATDENSE));
32   if (rank == 0) {
33     CHKERRQ(MatSetSizes(A, 4, PETSC_DETERMINE, PETSC_DETERMINE,PETSC_DETERMINE));
34   } else {
35     CHKERRQ(MatSetSizes(A, 8, PETSC_DETERMINE, PETSC_DETERMINE,PETSC_DETERMINE));
36   }
37   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"ex191matrix",FILE_MODE_READ,&fd));
38   CHKERRQ(MatLoad(A,fd));
39   CHKERRQ(PetscViewerDestroy(&fd));
40   CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
41   CHKERRQ(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
42   CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
43   CHKERRQ(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
44   CHKERRQ(MatDestroy(&A));
45   ierr = PetscFinalize();
46   return ierr;
47 }
48 
49 /*TEST
50 
51    test:
52       nsize: 2
53       filter: grep -v alloced
54 
55 TEST*/
56