xref: /petsc/src/mat/tests/ex190.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown static char help[] = "Tests MatLoad() 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   char           file[PETSC_MAX_PATH_LEN];
10c4762a1bSJed Brown   PetscBool      flg;
11c4762a1bSJed Brown   PetscMPIInt    rank;
12c4762a1bSJed Brown 
13*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
145f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
15c4762a1bSJed Brown 
16c4762a1bSJed Brown   /* Determine files from which we read the matrix */
175f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg));
1828b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f");
19c4762a1bSJed Brown 
20c4762a1bSJed Brown   /* Load matrices */
215f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd));
225f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
235f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
245f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetBlockSize(A,2));
25dd400576SPatrick Sanan   if (rank == 0) {
265f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetSizes(A, 4, PETSC_DETERMINE, PETSC_DETERMINE,PETSC_DETERMINE));
27c4762a1bSJed Brown   } else {
285f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetSizes(A, 8, PETSC_DETERMINE, PETSC_DETERMINE,PETSC_DETERMINE));
29c4762a1bSJed Brown   }
305f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(A,fd));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&fd));
325f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
33*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
34*b122ec5aSJacob Faibussowitsch   return 0;
35c4762a1bSJed Brown }
36c4762a1bSJed Brown 
37c4762a1bSJed Brown /*TEST
38c4762a1bSJed Brown 
39c4762a1bSJed Brown       test:
40c4762a1bSJed Brown          nsize: 2
41c4762a1bSJed Brown          args: -mat_type aij -mat_view -f ${wPETSC_DIR}/share/petsc/datafiles/matrices/ns-real-int32-float64
42dfd57a17SPierre Jolivet          requires: double !complex !defined(PETSC_USE_64BIT_INDICES)
43c4762a1bSJed Brown 
44c4762a1bSJed Brown       test:
45c4762a1bSJed Brown          suffix: 2
46c4762a1bSJed Brown          nsize: 2
47c4762a1bSJed Brown          args: -mat_type baij -mat_view -f ${wPETSC_DIR}/share/petsc/datafiles/matrices/ns-real-int32-float64
48dfd57a17SPierre Jolivet          requires: double !complex !defined(PETSC_USE_64BIT_INDICES)
49c4762a1bSJed Brown 
50c4762a1bSJed Brown       test:
51c4762a1bSJed Brown          suffix: 3
52c4762a1bSJed Brown          nsize: 2
53c4762a1bSJed Brown          args: -mat_type sbaij -mat_view -f ${wPETSC_DIR}/share/petsc/datafiles/matrices/ns-real-int32-float64
54dfd57a17SPierre Jolivet          requires: double !complex !defined(PETSC_USE_64BIT_INDICES)
55c4762a1bSJed Brown 
56c4762a1bSJed Brown TEST*/
57