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