xref: /petsc/src/mat/tests/ex191.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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   PetscErrorCode ierr;
10c4762a1bSJed Brown   PetscMPIInt    rank;
11c4762a1bSJed Brown   PetscScalar    *Av;
12c4762a1bSJed Brown   PetscInt       i;
13c4762a1bSJed Brown 
14c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
15*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
16c4762a1bSJed Brown 
17*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateDense(PETSC_COMM_WORLD,6,6,12,12,NULL,&A));
18*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDenseGetArray(A,&Av));
19c4762a1bSJed Brown   for (i=0; i<6*12; i++) Av[i] = (PetscScalar) i;
20*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDenseRestoreArray(A,&Av));
21c4762a1bSJed Brown 
22c4762a1bSJed Brown   /* Load matrices */
23*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"ex191matrix",FILE_MODE_WRITE,&fd));
24*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerPushFormat(fd,PETSC_VIEWER_NATIVE));
25*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(A,fd));
26*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
27*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerPopFormat(fd));
28*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&fd));
29c4762a1bSJed Brown 
30*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
31*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A,MATDENSE));
32dd400576SPatrick Sanan   if (rank == 0) {
33*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetSizes(A, 4, PETSC_DETERMINE, PETSC_DETERMINE,PETSC_DETERMINE));
34c4762a1bSJed Brown   } else {
35*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetSizes(A, 8, PETSC_DETERMINE, PETSC_DETERMINE,PETSC_DETERMINE));
36c4762a1bSJed Brown   }
37*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"ex191matrix",FILE_MODE_READ,&fd));
38*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(A,fd));
39*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&fd));
40*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
41*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
42*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
43*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
44*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
45c4762a1bSJed Brown   ierr = PetscFinalize();
46c4762a1bSJed Brown   return ierr;
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