xref: /petsc/src/mat/tests/ex84.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1f478cf82SStefano Zampini #include <petscmat.h>
2f478cf82SStefano Zampini 
3308ba433SVaclav Hapla #define NNORMS 6
4308ba433SVaclav Hapla 
5*9371c9d4SSatish Balay static PetscErrorCode MatLoadComputeNorms(Mat data_mat, PetscViewer inp_viewer, PetscReal norms[]) {
6308ba433SVaclav Hapla   Mat      corr_mat;
7308ba433SVaclav Hapla   PetscInt M, N;
8f478cf82SStefano Zampini 
9308ba433SVaclav Hapla   PetscFunctionBegin;
109566063dSJacob Faibussowitsch   PetscCall(MatLoad(data_mat, inp_viewer));
119566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(data_mat, MAT_FINAL_ASSEMBLY));
129566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(data_mat, MAT_FINAL_ASSEMBLY));
139566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(data_mat, NULL, "-view_mat"));
14f478cf82SStefano Zampini 
159566063dSJacob Faibussowitsch   PetscCall(MatGetSize(data_mat, &M, &N));
169566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Data matrix size: %" PetscInt_FMT " %" PetscInt_FMT "\n", M, N));
17308ba433SVaclav Hapla 
18308ba433SVaclav Hapla   /* compute matrix norms */
199566063dSJacob Faibussowitsch   PetscCall(MatNorm(data_mat, NORM_1, &norms[0]));
209566063dSJacob Faibussowitsch   PetscCall(MatNorm(data_mat, NORM_INFINITY, &norms[1]));
219566063dSJacob Faibussowitsch   PetscCall(MatNorm(data_mat, NORM_FROBENIUS, &norms[2]));
229566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Data matrix norms: %g %g %g\n", (double)norms[0], (double)norms[1], (double)norms[2]));
23f478cf82SStefano Zampini 
24f478cf82SStefano Zampini   /* compute autocorrelation matrix */
259566063dSJacob Faibussowitsch   PetscCall(MatMatTransposeMult(data_mat, data_mat, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &corr_mat));
26f478cf82SStefano Zampini 
27308ba433SVaclav Hapla   /* compute autocorrelation matrix norms */
289566063dSJacob Faibussowitsch   PetscCall(MatNorm(corr_mat, NORM_1, &norms[3]));
299566063dSJacob Faibussowitsch   PetscCall(MatNorm(corr_mat, NORM_INFINITY, &norms[4]));
309566063dSJacob Faibussowitsch   PetscCall(MatNorm(corr_mat, NORM_FROBENIUS, &norms[5]));
319566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Autocorrelation matrix norms: %g %g %g\n", (double)norms[3], (double)norms[4], (double)norms[5]));
32f478cf82SStefano Zampini 
339566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&corr_mat));
34308ba433SVaclav Hapla   PetscFunctionReturn(0);
35f478cf82SStefano Zampini }
36308ba433SVaclav Hapla 
37*9371c9d4SSatish Balay static PetscErrorCode GetReader(MPI_Comm comm, const char option[], PetscViewer *r, PetscViewerFormat *fmt) {
38308ba433SVaclav Hapla   PetscBool flg;
39308ba433SVaclav Hapla 
40308ba433SVaclav Hapla   PetscFunctionBegin;
419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetViewer(PETSC_COMM_SELF, NULL, NULL, option, r, fmt, &flg));
42308ba433SVaclav Hapla   if (flg) {
43308ba433SVaclav Hapla     PetscFileMode mode;
449566063dSJacob Faibussowitsch     PetscCall(PetscViewerFileGetMode(*r, &mode));
45308ba433SVaclav Hapla     flg = (PetscBool)(mode == FILE_MODE_READ);
46308ba433SVaclav Hapla   }
4728b400f6SJacob Faibussowitsch   PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER_INPUT, "Need to specify %s viewer_type:file:format:read", option);
48308ba433SVaclav Hapla   PetscFunctionReturn(0);
49308ba433SVaclav Hapla }
50308ba433SVaclav Hapla 
51*9371c9d4SSatish Balay int main(int argc, char **argv) {
52308ba433SVaclav Hapla   PetscInt          i;
53308ba433SVaclav Hapla   PetscReal         norms0[NNORMS], norms1[NNORMS];
54308ba433SVaclav Hapla   PetscViewer       inp_viewer;
55308ba433SVaclav Hapla   PetscViewerFormat fmt;
56308ba433SVaclav Hapla   Mat               data_mat;
57308ba433SVaclav Hapla   char              mat_name[PETSC_MAX_PATH_LEN] = "dmatrix";
58308ba433SVaclav Hapla 
59327415f7SBarry Smith   PetscFunctionBeginUser;
609566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, NULL));
619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-mat_name", mat_name, sizeof(mat_name), NULL));
62308ba433SVaclav Hapla 
63308ba433SVaclav Hapla   /* load matrix sequentially */
649566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &data_mat));
659566063dSJacob Faibussowitsch   PetscCall(MatSetType(data_mat, MATDENSE));
669566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)data_mat, mat_name));
679566063dSJacob Faibussowitsch   PetscCall(GetReader(PETSC_COMM_SELF, "-serial_reader", &inp_viewer, &fmt));
689566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(inp_viewer, fmt));
699566063dSJacob Faibussowitsch   PetscCall(MatLoadComputeNorms(data_mat, inp_viewer, norms0));
709566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(inp_viewer));
719566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&inp_viewer));
729566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(data_mat, NULL, "-view_serial_mat"));
739566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&data_mat));
74308ba433SVaclav Hapla 
75308ba433SVaclav Hapla   /* load matrix in parallel */
769566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &data_mat));
779566063dSJacob Faibussowitsch   PetscCall(MatSetType(data_mat, MATDENSE));
789566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)data_mat, mat_name));
799566063dSJacob Faibussowitsch   PetscCall(GetReader(PETSC_COMM_WORLD, "-parallel_reader", &inp_viewer, &fmt));
809566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(inp_viewer, fmt));
819566063dSJacob Faibussowitsch   PetscCall(MatLoadComputeNorms(data_mat, inp_viewer, norms1));
829566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(inp_viewer));
839566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&inp_viewer));
849566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(data_mat, NULL, "-view_parallel_mat"));
859566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&data_mat));
86308ba433SVaclav Hapla 
87*9371c9d4SSatish Balay   for (i = 0; i < NNORMS; i++) { PetscCheck(PetscAbs(norms0[i] - norms1[i]) <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "norm0[%" PetscInt_FMT "] = %g != %g = norms1[%" PetscInt_FMT "]", i, (double)norms0[i], (double)norms1[i], i); }
88308ba433SVaclav Hapla 
899566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
90b122ec5aSJacob Faibussowitsch   return 0;
91308ba433SVaclav Hapla }
92308ba433SVaclav Hapla 
93308ba433SVaclav Hapla /*TEST
94308ba433SVaclav Hapla 
95308ba433SVaclav Hapla   test:
96308ba433SVaclav Hapla     suffix: 1
97308ba433SVaclav Hapla     requires: hdf5 datafilespath complex
98308ba433SVaclav Hapla     args:  -serial_reader hdf5:${DATAFILESPATH}/matrices/hdf5/sample_data.h5::read -parallel_reader hdf5:${DATAFILESPATH}/matrices/hdf5/sample_data.h5::read
99308ba433SVaclav Hapla     nsize: {{1 2 4}}
100308ba433SVaclav Hapla 
101308ba433SVaclav Hapla   test:
102308ba433SVaclav Hapla     requires: hdf5 datafilespath
103308ba433SVaclav Hapla     args:  -serial_reader hdf5:${DATAFILESPATH}/matrices/hdf5/tiny_rectangular_mat.h5::read -parallel_reader hdf5:${DATAFILESPATH}/matrices/hdf5/tiny_rectangular_mat.h5::read
104308ba433SVaclav Hapla     nsize: {{1 2}}
105308ba433SVaclav Hapla     test:
106308ba433SVaclav Hapla       suffix: 2-complex
107308ba433SVaclav Hapla       requires: complex
108308ba433SVaclav Hapla       args: -mat_name ComplexMat
109308ba433SVaclav Hapla     test:
110308ba433SVaclav Hapla       suffix: 2-real
111308ba433SVaclav Hapla       requires: !complex
112308ba433SVaclav Hapla       args: -mat_name RealMat
113308ba433SVaclav Hapla 
114308ba433SVaclav Hapla TEST*/
115