1f478cf82SStefano Zampini #include <petscmat.h> 2f478cf82SStefano Zampini 3*308ba433SVaclav Hapla #define NNORMS 6 4*308ba433SVaclav Hapla 5*308ba433SVaclav Hapla static PetscErrorCode MatLoadComputeNorms(Mat data_mat, PetscViewer inp_viewer, PetscReal norms[]) 6f478cf82SStefano Zampini { 7*308ba433SVaclav Hapla Mat corr_mat; 8*308ba433SVaclav Hapla PetscInt M,N; 9f478cf82SStefano Zampini PetscErrorCode ierr; 10f478cf82SStefano Zampini 11*308ba433SVaclav Hapla PetscFunctionBegin; 12f478cf82SStefano Zampini ierr = MatLoad(data_mat, inp_viewer);CHKERRQ(ierr); 13f478cf82SStefano Zampini ierr = MatAssemblyBegin(data_mat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14f478cf82SStefano Zampini ierr = MatAssemblyEnd(data_mat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 15f478cf82SStefano Zampini ierr = MatViewFromOptions(data_mat, NULL, "-view_mat");CHKERRQ(ierr); 16f478cf82SStefano Zampini 17*308ba433SVaclav Hapla ierr = MatGetSize(data_mat, &M, &N);CHKERRQ(ierr); 18*308ba433SVaclav Hapla ierr = PetscPrintf(PETSC_COMM_WORLD, "Data matrix size: %D %D\n", M,N);CHKERRQ(ierr); 19*308ba433SVaclav Hapla 20*308ba433SVaclav Hapla /* compute matrix norms */ 21*308ba433SVaclav Hapla ierr = MatNorm(data_mat, NORM_1, &norms[0]);CHKERRQ(ierr); 22*308ba433SVaclav Hapla ierr = MatNorm(data_mat, NORM_INFINITY, &norms[1]);CHKERRQ(ierr); 23*308ba433SVaclav Hapla ierr = MatNorm(data_mat, NORM_FROBENIUS, &norms[2]);CHKERRQ(ierr); 24*308ba433SVaclav Hapla ierr = PetscPrintf(PETSC_COMM_WORLD, "Data matrix norms: %g %g %g\n", (double)norms[0],(double)norms[1],(double)norms[2]);CHKERRQ(ierr); 25f478cf82SStefano Zampini 26f478cf82SStefano Zampini /* compute autocorrelation matrix */ 27f478cf82SStefano Zampini ierr = MatMatTransposeMult(data_mat, data_mat, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &corr_mat);CHKERRQ(ierr); 28f478cf82SStefano Zampini 29*308ba433SVaclav Hapla /* compute autocorrelation matrix norms */ 30*308ba433SVaclav Hapla ierr = MatNorm(corr_mat, NORM_1, &norms[3]);CHKERRQ(ierr); 31*308ba433SVaclav Hapla ierr = MatNorm(corr_mat, NORM_INFINITY, &norms[4]);CHKERRQ(ierr); 32*308ba433SVaclav Hapla ierr = MatNorm(corr_mat, NORM_FROBENIUS, &norms[5]);CHKERRQ(ierr); 33*308ba433SVaclav Hapla ierr = PetscPrintf(PETSC_COMM_WORLD, "Autocorrelation matrix norms: %g %g %g\n", (double)norms[3],(double)norms[4],(double)norms[5]);CHKERRQ(ierr); 34f478cf82SStefano Zampini 35f478cf82SStefano Zampini ierr = MatDestroy(&corr_mat);CHKERRQ(ierr); 36*308ba433SVaclav Hapla PetscFunctionReturn(0); 37f478cf82SStefano Zampini } 38*308ba433SVaclav Hapla 39*308ba433SVaclav Hapla static PetscErrorCode GetReader(MPI_Comm comm, const char option[], PetscViewer *r, PetscViewerFormat *fmt) 40*308ba433SVaclav Hapla { 41*308ba433SVaclav Hapla PetscBool flg; 42*308ba433SVaclav Hapla PetscErrorCode ierr; 43*308ba433SVaclav Hapla 44*308ba433SVaclav Hapla PetscFunctionBegin; 45*308ba433SVaclav Hapla ierr = PetscOptionsGetViewer(PETSC_COMM_SELF, NULL, NULL, option, r, fmt, &flg);CHKERRQ(ierr); 46*308ba433SVaclav Hapla if (flg) { 47*308ba433SVaclav Hapla PetscFileMode mode; 48*308ba433SVaclav Hapla ierr = PetscViewerFileGetMode(*r, &mode);CHKERRQ(ierr); 49*308ba433SVaclav Hapla flg = (PetscBool) (mode == FILE_MODE_READ); 50*308ba433SVaclav Hapla } 51*308ba433SVaclav Hapla if (!flg) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_USER_INPUT,"Need to specify %s viewer_type:file:format:read", option); 52*308ba433SVaclav Hapla PetscFunctionReturn(0); 53*308ba433SVaclav Hapla } 54*308ba433SVaclav Hapla 55*308ba433SVaclav Hapla int main(int argc, char **argv) 56*308ba433SVaclav Hapla { 57*308ba433SVaclav Hapla PetscErrorCode ierr; 58*308ba433SVaclav Hapla PetscInt i; 59*308ba433SVaclav Hapla PetscReal norms0[NNORMS], norms1[NNORMS]; 60*308ba433SVaclav Hapla PetscViewer inp_viewer; 61*308ba433SVaclav Hapla PetscViewerFormat fmt; 62*308ba433SVaclav Hapla Mat data_mat; 63*308ba433SVaclav Hapla char mat_name[PETSC_MAX_PATH_LEN]="dmatrix"; 64*308ba433SVaclav Hapla 65*308ba433SVaclav Hapla ierr = PetscInitialize(&argc, &argv, NULL, NULL);if (ierr) return ierr; 66*308ba433SVaclav Hapla ierr = PetscOptionsGetString(NULL,NULL,"-mat_name",mat_name,sizeof(mat_name),NULL);CHKERRQ(ierr); 67*308ba433SVaclav Hapla 68*308ba433SVaclav Hapla /* load matrix sequentially */ 69*308ba433SVaclav Hapla ierr = MatCreate(PETSC_COMM_SELF, &data_mat);CHKERRQ(ierr); 70*308ba433SVaclav Hapla ierr = MatSetType(data_mat,MATDENSE);CHKERRQ(ierr); 71*308ba433SVaclav Hapla ierr = PetscObjectSetName((PetscObject)data_mat, mat_name);CHKERRQ(ierr); 72*308ba433SVaclav Hapla ierr = GetReader(PETSC_COMM_SELF, "-serial_reader", &inp_viewer, &fmt);CHKERRQ(ierr); 73*308ba433SVaclav Hapla ierr = PetscViewerPushFormat(inp_viewer, fmt);CHKERRQ(ierr); 74*308ba433SVaclav Hapla ierr = MatLoadComputeNorms(data_mat, inp_viewer, norms0);CHKERRQ(ierr); 75*308ba433SVaclav Hapla ierr = PetscViewerPopFormat(inp_viewer);CHKERRQ(ierr); 76*308ba433SVaclav Hapla ierr = PetscViewerDestroy(&inp_viewer);CHKERRQ(ierr); 77*308ba433SVaclav Hapla ierr = MatViewFromOptions(data_mat, NULL, "-view_serial_mat");CHKERRQ(ierr); 78*308ba433SVaclav Hapla ierr = MatDestroy(&data_mat);CHKERRQ(ierr); 79*308ba433SVaclav Hapla 80*308ba433SVaclav Hapla /* load matrix in parallel */ 81*308ba433SVaclav Hapla ierr = MatCreate(PETSC_COMM_WORLD, &data_mat);CHKERRQ(ierr); 82*308ba433SVaclav Hapla ierr = MatSetType(data_mat,MATDENSE);CHKERRQ(ierr); 83*308ba433SVaclav Hapla ierr = PetscObjectSetName((PetscObject)data_mat, mat_name);CHKERRQ(ierr); 84*308ba433SVaclav Hapla ierr = GetReader(PETSC_COMM_WORLD, "-parallel_reader", &inp_viewer, &fmt);CHKERRQ(ierr); 85*308ba433SVaclav Hapla ierr = PetscViewerPushFormat(inp_viewer, fmt);CHKERRQ(ierr); 86*308ba433SVaclav Hapla ierr = MatLoadComputeNorms(data_mat, inp_viewer, norms1);CHKERRQ(ierr); 87*308ba433SVaclav Hapla ierr = PetscViewerPopFormat(inp_viewer);CHKERRQ(ierr); 88*308ba433SVaclav Hapla ierr = PetscViewerDestroy(&inp_viewer);CHKERRQ(ierr); 89*308ba433SVaclav Hapla ierr = MatViewFromOptions(data_mat, NULL, "-view_parallel_mat");CHKERRQ(ierr); 90*308ba433SVaclav Hapla ierr = MatDestroy(&data_mat);CHKERRQ(ierr); 91*308ba433SVaclav Hapla 92*308ba433SVaclav Hapla for (i=0; i<NNORMS; i++) { 93*308ba433SVaclav Hapla if (PetscAbs(norms0[i] - norms1[i]) > PETSC_SMALL) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_PLIB, "norm0[%D] = %g != %g = norms1[%D]", i, (double)norms0[i], (double)norms1[i], i); 94*308ba433SVaclav Hapla } 95*308ba433SVaclav Hapla 96*308ba433SVaclav Hapla ierr = PetscFinalize(); 97*308ba433SVaclav Hapla return ierr; 98*308ba433SVaclav Hapla } 99*308ba433SVaclav Hapla 100*308ba433SVaclav Hapla /*TEST 101*308ba433SVaclav Hapla 102*308ba433SVaclav Hapla test: 103*308ba433SVaclav Hapla suffix: 1 104*308ba433SVaclav Hapla requires: hdf5 datafilespath complex 105*308ba433SVaclav Hapla args: -serial_reader hdf5:${DATAFILESPATH}/matrices/hdf5/sample_data.h5::read -parallel_reader hdf5:${DATAFILESPATH}/matrices/hdf5/sample_data.h5::read 106*308ba433SVaclav Hapla nsize: {{1 2 4}} 107*308ba433SVaclav Hapla 108*308ba433SVaclav Hapla test: 109*308ba433SVaclav Hapla requires: hdf5 datafilespath 110*308ba433SVaclav Hapla args: -serial_reader hdf5:${DATAFILESPATH}/matrices/hdf5/tiny_rectangular_mat.h5::read -parallel_reader hdf5:${DATAFILESPATH}/matrices/hdf5/tiny_rectangular_mat.h5::read 111*308ba433SVaclav Hapla nsize: {{1 2}} 112*308ba433SVaclav Hapla test: 113*308ba433SVaclav Hapla suffix: 2-complex 114*308ba433SVaclav Hapla requires: complex 115*308ba433SVaclav Hapla args: -mat_name ComplexMat 116*308ba433SVaclav Hapla test: 117*308ba433SVaclav Hapla suffix: 2-real 118*308ba433SVaclav Hapla requires: !complex 119*308ba433SVaclav Hapla args: -mat_name RealMat 120*308ba433SVaclav Hapla 121*308ba433SVaclav Hapla TEST*/ 122