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