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