xref: /petsc/src/mat/tests/ex191.c (revision 9566063d113dddea24716c546802770db7481bc0)
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*9566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
14*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
15c4762a1bSJed Brown 
16*9566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PETSC_COMM_WORLD,6,6,12,12,NULL,&A));
17*9566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(A,&Av));
18c4762a1bSJed Brown   for (i=0; i<6*12; i++) Av[i] = (PetscScalar) i;
19*9566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(A,&Av));
20c4762a1bSJed Brown 
21c4762a1bSJed Brown   /* Load matrices */
22*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"ex191matrix",FILE_MODE_WRITE,&fd));
23*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(fd,PETSC_VIEWER_NATIVE));
24*9566063dSJacob Faibussowitsch   PetscCall(MatView(A,fd));
25*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
26*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(fd));
27*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&fd));
28c4762a1bSJed Brown 
29*9566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
30*9566063dSJacob Faibussowitsch   PetscCall(MatSetType(A,MATDENSE));
31dd400576SPatrick Sanan   if (rank == 0) {
32*9566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(A, 4, PETSC_DETERMINE, PETSC_DETERMINE,PETSC_DETERMINE));
33c4762a1bSJed Brown   } else {
34*9566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(A, 8, PETSC_DETERMINE, PETSC_DETERMINE,PETSC_DETERMINE));
35c4762a1bSJed Brown   }
36*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"ex191matrix",FILE_MODE_READ,&fd));
37*9566063dSJacob Faibussowitsch   PetscCall(MatLoad(A,fd));
38*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&fd));
39*9566063dSJacob Faibussowitsch   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
40*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
41*9566063dSJacob Faibussowitsch   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
42*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
43*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
44*9566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
45b122ec5aSJacob Faibussowitsch   return 0;
46c4762a1bSJed Brown }
47c4762a1bSJed Brown 
48c4762a1bSJed Brown /*TEST
49c4762a1bSJed Brown 
50c4762a1bSJed Brown    test:
51c4762a1bSJed Brown       nsize: 2
52c4762a1bSJed Brown       filter: grep -v alloced
53c4762a1bSJed Brown 
54c4762a1bSJed Brown TEST*/
55