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