xref: /petsc/src/mat/tests/ex138.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
18   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
19   CHKERRQ(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   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd));
22   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
23   CHKERRQ(MatSetFromOptions(A));
24   CHKERRQ(MatLoad(A,fd));
25   CHKERRQ(PetscViewerDestroy(&fd));
26 
27   CHKERRQ(MatGetSize(A,NULL,&n));
28   CHKERRQ(PetscMalloc1(n,&reductions_real));
29   CHKERRQ(PetscMalloc1(n,&reductions_scalar));
30 
31   CHKERRQ(MatGetColumnNorms(A,NORM_2,reductions_real));
32   if (rank == 0) {
33     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"NORM_2:\n"));
34     CHKERRQ(PetscRealView(n,reductions_real,PETSC_VIEWER_STDOUT_SELF));
35   }
36 
37   CHKERRQ(MatGetColumnNorms(A,NORM_1,reductions_real));
38   if (rank == 0) {
39     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"NORM_1:\n"));
40     CHKERRQ(PetscRealView(n,reductions_real,PETSC_VIEWER_STDOUT_SELF));
41   }
42 
43   CHKERRQ(MatGetColumnNorms(A,NORM_INFINITY,reductions_real));
44   if (rank == 0) {
45     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"NORM_INFINITY:\n"));
46     CHKERRQ(PetscRealView(n,reductions_real,PETSC_VIEWER_STDOUT_SELF));
47   }
48 
49   CHKERRQ(MatGetColumnSums(A,reductions_scalar));
50   if (!rank) {
51     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"REDUCTION_SUM:\n"));
52     CHKERRQ(PetscScalarView(n,reductions_scalar,PETSC_VIEWER_STDOUT_SELF));
53   }
54 
55   CHKERRQ(MatGetColumnMeans(A,reductions_scalar));
56   if (!rank) {
57     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"REDUCTION_MEAN:\n"));
58     CHKERRQ(PetscScalarView(n,reductions_scalar,PETSC_VIEWER_STDOUT_SELF));
59   }
60 
61   CHKERRQ(PetscFree(reductions_real));
62   CHKERRQ(PetscFree(reductions_scalar));
63   CHKERRQ(MatDestroy(&A));
64   CHKERRQ(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