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