xref: /petsc/src/mat/tests/ex138.c (revision dfd57a172ac9fa6c7b5fe6de6ab5df85cefc2996)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatGetColumnNorms() 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;
10c4762a1bSJed Brown   PetscReal      *norms;
11c4762a1bSJed Brown   char           file[PETSC_MAX_PATH_LEN];
12c4762a1bSJed Brown   PetscBool      flg;
13c4762a1bSJed Brown   PetscViewer    fd;
14c4762a1bSJed Brown   PetscInt       n;
15c4762a1bSJed Brown   PetscMPIInt    rank;
16c4762a1bSJed Brown 
17c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
18ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr);
19589a23caSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg);CHKERRQ(ierr);
20c4762a1bSJed Brown   if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option");
21c4762a1bSJed Brown   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);CHKERRQ(ierr);
22c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
23c4762a1bSJed Brown   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
24c4762a1bSJed Brown   ierr = MatLoad(A,fd);CHKERRQ(ierr);
25c4762a1bSJed Brown   ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr);
26c4762a1bSJed Brown 
27c4762a1bSJed Brown   ierr = MatGetSize(A,NULL,&n);CHKERRQ(ierr);
28c4762a1bSJed Brown   ierr = PetscMalloc1(n,&norms);CHKERRQ(ierr);
29c4762a1bSJed Brown 
30c4762a1bSJed Brown   ierr = MatGetColumnNorms(A,NORM_2,norms);CHKERRQ(ierr);
31c4762a1bSJed Brown   if (!rank) {
32c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF,"NORM_2:\n");CHKERRQ(ierr);
33c4762a1bSJed Brown     ierr = PetscRealView(n,norms,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
34c4762a1bSJed Brown   }
35c4762a1bSJed Brown 
36c4762a1bSJed Brown   ierr = MatGetColumnNorms(A,NORM_1,norms);CHKERRQ(ierr);
37c4762a1bSJed Brown   if (!rank) {
38c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF,"NORM_1:\n");CHKERRQ(ierr);
39c4762a1bSJed Brown     ierr = PetscRealView(n,norms,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
40c4762a1bSJed Brown   }
41c4762a1bSJed Brown 
42c4762a1bSJed Brown   ierr = MatGetColumnNorms(A,NORM_INFINITY,norms);CHKERRQ(ierr);
43c4762a1bSJed Brown   if (!rank) {
44c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF,"NORM_INFINITY:\n");CHKERRQ(ierr);
45c4762a1bSJed Brown     ierr = PetscRealView(n,norms,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
46c4762a1bSJed Brown   }
47c4762a1bSJed Brown 
48c4762a1bSJed Brown   ierr = PetscFree(norms);CHKERRQ(ierr);
49c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
50c4762a1bSJed Brown   ierr = PetscFinalize();
51c4762a1bSJed Brown   return ierr;
52c4762a1bSJed Brown }
53c4762a1bSJed Brown 
54c4762a1bSJed Brown /*TEST
55c4762a1bSJed Brown 
56c4762a1bSJed Brown    test:
57c4762a1bSJed Brown       suffix: 1
58c4762a1bSJed Brown       nsize: 2
59*dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
60c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small -mat_type aij
61c4762a1bSJed Brown       output_file: output/ex138.out
62c4762a1bSJed Brown 
63c4762a1bSJed Brown    test:
64c4762a1bSJed Brown       suffix: 2
65c4762a1bSJed Brown       nsize: 2
66*dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
67c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small -mat_type baij -matload_block_size {{2 3}}
68c4762a1bSJed Brown       output_file: output/ex138.out
69c4762a1bSJed Brown 
70c4762a1bSJed Brown TEST*/
71