xref: /petsc/src/mat/tests/ex138.c (revision 98d129c30f3ee9fdddc40fdbc5a989b7be64f888)
1 static char help[] = "Tests MatGetColumnNorms()/Sums()/Means() for matrix read from file.";
2 
3 #include <petscmat.h>
4 
5 int main(int argc, char **args)
6 {
7   Mat          A;
8   PetscReal   *reductions_real;
9   PetscScalar *reductions_scalar;
10   char         file[PETSC_MAX_PATH_LEN];
11   PetscBool    flg;
12   PetscViewer  fd;
13   PetscInt     n;
14   PetscMPIInt  rank;
15 
16   PetscFunctionBeginUser;
17   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
18   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
19   PetscCall(PetscOptionsGetString(NULL, NULL, "-f", file, sizeof(file), &flg));
20   PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "Must indicate binary file with the -f option");
21   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &fd));
22   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
23   PetscCall(MatSetFromOptions(A));
24   PetscCall(MatLoad(A, fd));
25   PetscCall(PetscViewerDestroy(&fd));
26 
27   PetscCall(MatGetSize(A, NULL, &n));
28   PetscCall(PetscMalloc1(n, &reductions_real));
29   PetscCall(PetscMalloc1(n, &reductions_scalar));
30 
31   PetscCall(MatGetColumnNorms(A, NORM_2, reductions_real));
32   if (rank == 0) {
33     PetscCall(PetscPrintf(PETSC_COMM_SELF, "NORM_2:\n"));
34     PetscCall(PetscRealView(n, reductions_real, PETSC_VIEWER_STDOUT_SELF));
35   }
36 
37   PetscCall(MatGetColumnNorms(A, NORM_1, reductions_real));
38   if (rank == 0) {
39     PetscCall(PetscPrintf(PETSC_COMM_SELF, "NORM_1:\n"));
40     PetscCall(PetscRealView(n, reductions_real, PETSC_VIEWER_STDOUT_SELF));
41   }
42 
43   PetscCall(MatGetColumnNorms(A, NORM_INFINITY, reductions_real));
44   if (rank == 0) {
45     PetscCall(PetscPrintf(PETSC_COMM_SELF, "NORM_INFINITY:\n"));
46     PetscCall(PetscRealView(n, reductions_real, PETSC_VIEWER_STDOUT_SELF));
47   }
48 
49   PetscCall(MatGetColumnSums(A, reductions_scalar));
50   if (rank == 0) {
51     PetscCall(PetscPrintf(PETSC_COMM_SELF, "REDUCTION_SUM:\n"));
52     PetscCall(PetscScalarView(n, reductions_scalar, PETSC_VIEWER_STDOUT_SELF));
53   }
54 
55   PetscCall(MatGetColumnMeans(A, reductions_scalar));
56   if (rank == 0) {
57     PetscCall(PetscPrintf(PETSC_COMM_SELF, "REDUCTION_MEAN:\n"));
58     PetscCall(PetscScalarView(n, reductions_scalar, PETSC_VIEWER_STDOUT_SELF));
59   }
60 
61   PetscCall(PetscFree(reductions_real));
62   PetscCall(PetscFree(reductions_scalar));
63   PetscCall(MatDestroy(&A));
64   PetscCall(PetscFinalize());
65   return 0;
66 }
67 
68 /*TEST
69 
70    test:
71       suffix: 1
72       nsize: 2
73       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
74       args: -f ${DATAFILESPATH}/matrices/small -mat_type aij
75       output_file: output/ex138.out
76 
77    test:
78       suffix: 2
79       nsize: {{1 2}}
80       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
81       args: -f ${DATAFILESPATH}/matrices/small -mat_type baij -matload_block_size {{2 3}}
82       output_file: output/ex138.out
83 
84    test:
85       suffix: complex
86       nsize: 2
87       requires: datafilespath complex double !defined(PETSC_USE_64BIT_INDICES)
88       args: -f ${DATAFILESPATH}/matrices/nimrod/small_112905 -mat_type aij
89       output_file: output/ex138_complex.out
90       filter: grep -E "\ 0:|1340:|1344:"
91 
92 TEST*/
93