xref: /petsc/src/mat/tests/ex190.c (revision 6a98f8dc3f2c9149905a87dc2e9d0fedaf64e09a)
1 static char help[] = "Tests MatLoad() 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   char           file[PETSC_MAX_PATH_LEN];
10   PetscErrorCode ierr;
11   PetscBool      flg;
12   PetscMPIInt    rank;
13 
14   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
15   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
16 
17   /* Determine files from which we read the matrix */
18   ierr = PetscOptionsGetString(NULL,NULL,"-f",file,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
19   if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f");
20 
21   /* Load matrices */
22   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);CHKERRQ(ierr);
23   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
24   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
25   ierr = MatSetBlockSize(A,2);CHKERRQ(ierr);
26   if (!rank) {
27     ierr = MatSetSizes(A, 4, PETSC_DETERMINE, PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
28   } else {
29     ierr = MatSetSizes(A, 8, PETSC_DETERMINE, PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
30   }
31   ierr = MatLoad(A,fd);CHKERRQ(ierr);
32   ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr);
33   ierr = MatDestroy(&A);CHKERRQ(ierr);
34   ierr = PetscFinalize();
35   return ierr;
36 }
37 
38 
39 /*TEST
40 
41       test:
42          nsize: 2
43          args: -mat_type aij -mat_view -f ${wPETSC_DIR}/share/petsc/datafiles/matrices/ns-real-int32-float64
44          requires: double !complex !define(PETSC_USE_64BIT_INDICES)
45 
46       test:
47          suffix: 2
48          nsize: 2
49          args: -mat_type baij -mat_view -f ${wPETSC_DIR}/share/petsc/datafiles/matrices/ns-real-int32-float64
50          requires: double !complex !define(PETSC_USE_64BIT_INDICES)
51 
52       test:
53          suffix: 3
54          nsize: 2
55          args: -mat_type sbaij -mat_view -f ${wPETSC_DIR}/share/petsc/datafiles/matrices/ns-real-int32-float64
56          requires: double !complex !define(PETSC_USE_64BIT_INDICES)
57 
58 TEST*/
59