1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Reads a PETSc matrix and computes the 2 norm of the columns\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /*T 5c4762a1bSJed Brown Concepts: Mat^loading a binary matrix; 6c4762a1bSJed Brown Processors: n 7c4762a1bSJed Brown T*/ 8c4762a1bSJed Brown 9c4762a1bSJed Brown /* 10c4762a1bSJed Brown Include "petscmat.h" so that we can use matrices. 11c4762a1bSJed Brown automatically includes: 12c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 13c4762a1bSJed Brown petscmat.h - matrices 14c4762a1bSJed Brown petscis.h - index sets petscviewer.h - viewers 15c4762a1bSJed Brown */ 16c4762a1bSJed Brown #include <petscmat.h> 17c4762a1bSJed Brown 18c4762a1bSJed Brown int main(int argc,char **args) 19c4762a1bSJed Brown { 20c4762a1bSJed Brown Mat A; /* matrix */ 21c4762a1bSJed Brown PetscViewer fd; /* viewer */ 22c4762a1bSJed Brown char file[PETSC_MAX_PATH_LEN]; /* input file name */ 23c4762a1bSJed Brown PetscErrorCode ierr; 24c4762a1bSJed Brown PetscReal *norms; 25c4762a1bSJed Brown PetscInt n,cstart,cend; 26c4762a1bSJed Brown PetscBool flg; 27c4762a1bSJed Brown PetscViewerFormat format; 28c4762a1bSJed Brown 29c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 30c4762a1bSJed Brown /* 31c4762a1bSJed Brown Determine files from which we read the matrix 32c4762a1bSJed Brown */ 33589a23caSBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg);CHKERRQ(ierr); 34c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option"); 35c4762a1bSJed Brown 36c4762a1bSJed Brown /* 37c4762a1bSJed Brown Open binary file. Note that we use FILE_MODE_READ to indicate 38c4762a1bSJed Brown reading from this file. 39c4762a1bSJed Brown */ 40c4762a1bSJed Brown ierr = PetscViewerCreate(PETSC_COMM_WORLD,&fd);CHKERRQ(ierr); 41c4762a1bSJed Brown ierr = PetscViewerSetType(fd,PETSCVIEWERBINARY);CHKERRQ(ierr); 42c4762a1bSJed Brown ierr = PetscViewerSetFromOptions(fd);CHKERRQ(ierr); 43c4762a1bSJed Brown ierr = PetscOptionsGetEnum(NULL,NULL,"-viewer_format",PetscViewerFormats,(PetscEnum*)&format,&flg);CHKERRQ(ierr); 44c4762a1bSJed Brown if (flg) {ierr = PetscViewerPushFormat(fd,format);CHKERRQ(ierr);} 45c4762a1bSJed Brown ierr = PetscViewerFileSetMode(fd,FILE_MODE_READ);CHKERRQ(ierr); 46c4762a1bSJed Brown ierr = PetscViewerFileSetName(fd,file);CHKERRQ(ierr); 47c4762a1bSJed Brown 48c4762a1bSJed Brown /* 49c4762a1bSJed Brown Load the matrix; then destroy the viewer. 50c4762a1bSJed Brown Matrix type is set automatically but you can override it by MatSetType() prior to MatLoad(). 51c4762a1bSJed Brown Do that only if you really insist on the given type. 52c4762a1bSJed Brown */ 53c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 54c4762a1bSJed Brown ierr = MatSetOptionsPrefix(A,"a_");CHKERRQ(ierr); 55c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) A,"A");CHKERRQ(ierr); 56c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 57c4762a1bSJed Brown ierr = MatLoad(A,fd);CHKERRQ(ierr); 58c4762a1bSJed Brown ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr); 59c4762a1bSJed Brown 60c4762a1bSJed Brown ierr = MatGetSize(A,NULL,&n);CHKERRQ(ierr); 61c4762a1bSJed Brown ierr = MatGetOwnershipRangeColumn(A,&cstart,&cend);CHKERRQ(ierr); 62c4762a1bSJed Brown ierr = PetscMalloc1(n,&norms);CHKERRQ(ierr); 63c4762a1bSJed Brown ierr = MatGetColumnNorms(A,NORM_2,norms);CHKERRQ(ierr); 64c4762a1bSJed Brown ierr = PetscRealView(cend-cstart,norms+cstart,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 65c4762a1bSJed Brown ierr = PetscFree(norms);CHKERRQ(ierr); 66c4762a1bSJed Brown 67c4762a1bSJed Brown ierr = PetscObjectPrintClassNamePrefixType((PetscObject)A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 68c4762a1bSJed Brown ierr = MatGetOption(A,MAT_SYMMETRIC,&flg);CHKERRQ(ierr); 69c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"MAT_SYMMETRIC: %D\n",flg);CHKERRQ(ierr); 70c4762a1bSJed Brown ierr = MatViewFromOptions(A,NULL,"-mat_view");CHKERRQ(ierr); 71c4762a1bSJed Brown 72c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 73c4762a1bSJed Brown ierr = PetscFinalize(); 74c4762a1bSJed Brown return ierr; 75c4762a1bSJed Brown } 76c4762a1bSJed Brown 77c4762a1bSJed Brown /*TEST 78c4762a1bSJed Brown 79c4762a1bSJed Brown test: 80c4762a1bSJed Brown suffix: mpiaij 81c4762a1bSJed Brown nsize: 2 82*dfd57a17SPierre Jolivet requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) 83c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/small -a_mat_type mpiaij 84c4762a1bSJed Brown args: -a_matload_symmetric 85c4762a1bSJed Brown 86c4762a1bSJed Brown test: 87c4762a1bSJed Brown suffix: mpiaij_hdf5 88c4762a1bSJed Brown nsize: 2 89*dfd57a17SPierre Jolivet requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) hdf5 defined(PETSC_HDF5_HAVE_ZLIB) 90c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/matlab/small.mat -a_mat_type mpiaij -viewer_type hdf5 -viewer_format hdf5_mat 91c4762a1bSJed Brown args: -a_matload_symmetric 92c4762a1bSJed Brown 93c4762a1bSJed Brown test: 94c4762a1bSJed Brown suffix: mpiaij_rect_hdf5 95c4762a1bSJed Brown nsize: 2 96*dfd57a17SPierre Jolivet requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) hdf5 defined(PETSC_HDF5_HAVE_ZLIB) 97c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/matlab/small_rect.mat -a_mat_type mpiaij -viewer_type hdf5 -viewer_format hdf5_mat 98c4762a1bSJed Brown 99c4762a1bSJed Brown test: 100c4762a1bSJed Brown # test for more processes than rows 101c4762a1bSJed Brown suffix: mpiaij_hdf5_tiny 102c4762a1bSJed Brown nsize: 8 103*dfd57a17SPierre Jolivet requires: double !complex !defined(PETSC_USE_64BIT_INDICES) hdf5 defined(PETSC_HDF5_HAVE_ZLIB) 104c4762a1bSJed Brown args: -f ${wPETSC_DIR}/share/petsc/datafiles/matrices/tiny_system_with_x0.mat -a_mat_type mpiaij -viewer_type hdf5 -viewer_format hdf5_mat 105c4762a1bSJed Brown args: -a_matload_symmetric 106c4762a1bSJed Brown 107c4762a1bSJed Brown test: 108c4762a1bSJed Brown # test for more processes than rows, complex 109c4762a1bSJed Brown TODO: not yet implemented for MATLAB complex format 110c4762a1bSJed Brown suffix: mpiaij_hdf5_tiny_complex 111c4762a1bSJed Brown nsize: 8 112*dfd57a17SPierre Jolivet requires: double complex !defined(PETSC_USE_64BIT_INDICES) hdf5 defined(PETSC_HDF5_HAVE_ZLIB) 113c4762a1bSJed Brown args: -f ${wPETSC_DIR}/share/petsc/datafiles/matrices/tiny_system_with_x0_complex.mat -a_mat_type mpiaij -viewer_type hdf5 -viewer_format hdf5_mat 114c4762a1bSJed Brown args: -a_matload_symmetric 115c4762a1bSJed Brown 116c4762a1bSJed Brown test: 117c4762a1bSJed Brown TODO: mpibaij not supported yet 118c4762a1bSJed Brown suffix: mpibaij_hdf5 119c4762a1bSJed Brown nsize: 2 120*dfd57a17SPierre Jolivet requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) hdf5 defined(PETSC_HDF5_HAVE_ZLIB) 121c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/matlab/small.mat -a_mat_type mpibaij -a_mat_block_size 2 -viewer_type hdf5 -viewer_format hdf5_mat 122c4762a1bSJed Brown args: -a_matload_symmetric 123c4762a1bSJed Brown 124c4762a1bSJed Brown test: 125c4762a1bSJed Brown suffix: mpidense 126c4762a1bSJed Brown nsize: 2 127*dfd57a17SPierre Jolivet requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) 128c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/small -a_mat_type mpidense 129c4762a1bSJed Brown args: -a_matload_symmetric 130c4762a1bSJed Brown 131c4762a1bSJed Brown test: 132c4762a1bSJed Brown suffix: seqaij 133*dfd57a17SPierre Jolivet requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) 134c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/small -a_mat_type seqaij 135c4762a1bSJed Brown args: -a_matload_symmetric 136c4762a1bSJed Brown 137c4762a1bSJed Brown test: 138c4762a1bSJed Brown suffix: seqaij_hdf5 139*dfd57a17SPierre Jolivet requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) hdf5 defined(PETSC_HDF5_HAVE_ZLIB) 140c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/matlab/small.mat -a_mat_type seqaij -viewer_type hdf5 -viewer_format hdf5_mat 141c4762a1bSJed Brown args: -a_matload_symmetric 142c4762a1bSJed Brown 143c4762a1bSJed Brown test: 144c4762a1bSJed Brown suffix: seqaij_rect_hdf5 145*dfd57a17SPierre Jolivet requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) hdf5 defined(PETSC_HDF5_HAVE_ZLIB) 146c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/matlab/small_rect.mat -a_mat_type seqaij -viewer_type hdf5 -viewer_format hdf5_mat 147c4762a1bSJed Brown 148c4762a1bSJed Brown test: 149c4762a1bSJed Brown suffix: seqdense 150*dfd57a17SPierre Jolivet requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) 151c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/small -a_mat_type seqdense 152c4762a1bSJed Brown args: -a_matload_symmetric 153c4762a1bSJed Brown 154c4762a1bSJed Brown test: 155c4762a1bSJed Brown suffix: seqdense_hdf5 156*dfd57a17SPierre Jolivet requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) hdf5 defined(PETSC_HDF5_HAVE_ZLIB) 157c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/matlab/small_dense.mat -a_mat_type seqdense -viewer_type hdf5 -viewer_format hdf5_mat 158c4762a1bSJed Brown args: -a_matload_symmetric 159c4762a1bSJed Brown 160c4762a1bSJed Brown test: 161c4762a1bSJed Brown suffix: seqdense_rect_hdf5 162*dfd57a17SPierre Jolivet requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) hdf5 defined(PETSC_HDF5_HAVE_ZLIB) 163c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/matlab/small_rect_dense.mat -a_mat_type seqdense -viewer_type hdf5 -viewer_format hdf5_mat 164c4762a1bSJed Brown 165c4762a1bSJed Brown test: 166c4762a1bSJed Brown suffix: mpidense_hdf5 167c4762a1bSJed Brown nsize: 2 168*dfd57a17SPierre Jolivet requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) hdf5 defined(PETSC_HDF5_HAVE_ZLIB) 169c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/matlab/small_dense.mat -a_mat_type mpidense -viewer_type hdf5 -viewer_format hdf5_mat 170c4762a1bSJed Brown args: -a_matload_symmetric 171c4762a1bSJed Brown 172c4762a1bSJed Brown test: 173c4762a1bSJed Brown suffix: mpidense_rect_hdf5 174c4762a1bSJed Brown nsize: 2 175*dfd57a17SPierre Jolivet requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) hdf5 defined(PETSC_HDF5_HAVE_ZLIB) 176c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/matlab/small_rect_dense.mat -a_mat_type mpidense -viewer_type hdf5 -viewer_format hdf5_mat 177c4762a1bSJed Brown TEST*/ 178