xref: /petsc/src/mat/tests/ex191.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4762a1bSJed Brown static char help[] = "Tests MatLoad() for dense matrix with uneven dimensions set in program\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown int main(int argc,char **args)
6c4762a1bSJed Brown {
7c4762a1bSJed Brown   Mat            A;
8c4762a1bSJed Brown   PetscViewer    fd;
9c4762a1bSJed Brown   PetscMPIInt    rank;
10c4762a1bSJed Brown   PetscScalar    *Av;
11c4762a1bSJed Brown   PetscInt       i;
12c4762a1bSJed Brown 
13*327415f7SBarry Smith   PetscFunctionBeginUser;
149566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
159566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
16c4762a1bSJed Brown 
179566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PETSC_COMM_WORLD,6,6,12,12,NULL,&A));
189566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(A,&Av));
19c4762a1bSJed Brown   for (i=0; i<6*12; i++) Av[i] = (PetscScalar) i;
209566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(A,&Av));
21c4762a1bSJed Brown 
22c4762a1bSJed Brown   /* Load matrices */
239566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"ex191matrix",FILE_MODE_WRITE,&fd));
249566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(fd,PETSC_VIEWER_NATIVE));
259566063dSJacob Faibussowitsch   PetscCall(MatView(A,fd));
269566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
279566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(fd));
289566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&fd));
29c4762a1bSJed Brown 
309566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
319566063dSJacob Faibussowitsch   PetscCall(MatSetType(A,MATDENSE));
32dd400576SPatrick Sanan   if (rank == 0) {
339566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(A, 4, PETSC_DETERMINE, PETSC_DETERMINE,PETSC_DETERMINE));
34c4762a1bSJed Brown   } else {
359566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(A, 8, PETSC_DETERMINE, PETSC_DETERMINE,PETSC_DETERMINE));
36c4762a1bSJed Brown   }
379566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"ex191matrix",FILE_MODE_READ,&fd));
389566063dSJacob Faibussowitsch   PetscCall(MatLoad(A,fd));
399566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&fd));
409566063dSJacob Faibussowitsch   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
419566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
429566063dSJacob Faibussowitsch   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
439566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
449566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
459566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
46b122ec5aSJacob Faibussowitsch   return 0;
47c4762a1bSJed Brown }
48c4762a1bSJed Brown 
49c4762a1bSJed Brown /*TEST
50c4762a1bSJed Brown 
51c4762a1bSJed Brown    test:
52c4762a1bSJed Brown       nsize: 2
53c4762a1bSJed Brown       filter: grep -v alloced
54c4762a1bSJed Brown 
55c4762a1bSJed Brown TEST*/
56