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