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