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