1f478cf82SStefano Zampini #include <petscmat.h> 2f478cf82SStefano Zampini 3308ba433SVaclav Hapla #define NNORMS 6 4308ba433SVaclav Hapla 5308ba433SVaclav Hapla static PetscErrorCode MatLoadComputeNorms(Mat data_mat, PetscViewer inp_viewer, PetscReal norms[]) 6f478cf82SStefano Zampini { 7308ba433SVaclav Hapla Mat corr_mat; 8308ba433SVaclav Hapla PetscInt M,N; 9f478cf82SStefano Zampini 10308ba433SVaclav Hapla PetscFunctionBegin; 119566063dSJacob Faibussowitsch PetscCall(MatLoad(data_mat, inp_viewer)); 129566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(data_mat, MAT_FINAL_ASSEMBLY)); 139566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(data_mat, MAT_FINAL_ASSEMBLY)); 149566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(data_mat, NULL, "-view_mat")); 15f478cf82SStefano Zampini 169566063dSJacob Faibussowitsch PetscCall(MatGetSize(data_mat, &M, &N)); 179566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Data matrix size: %" PetscInt_FMT " %" PetscInt_FMT "\n", M,N)); 18308ba433SVaclav Hapla 19308ba433SVaclav Hapla /* compute matrix norms */ 209566063dSJacob Faibussowitsch PetscCall(MatNorm(data_mat, NORM_1, &norms[0])); 219566063dSJacob Faibussowitsch PetscCall(MatNorm(data_mat, NORM_INFINITY, &norms[1])); 229566063dSJacob Faibussowitsch PetscCall(MatNorm(data_mat, NORM_FROBENIUS, &norms[2])); 239566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Data matrix norms: %g %g %g\n", (double)norms[0],(double)norms[1],(double)norms[2])); 24f478cf82SStefano Zampini 25f478cf82SStefano Zampini /* compute autocorrelation matrix */ 269566063dSJacob Faibussowitsch PetscCall(MatMatTransposeMult(data_mat, data_mat, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &corr_mat)); 27f478cf82SStefano Zampini 28308ba433SVaclav Hapla /* compute autocorrelation matrix norms */ 299566063dSJacob Faibussowitsch PetscCall(MatNorm(corr_mat, NORM_1, &norms[3])); 309566063dSJacob Faibussowitsch PetscCall(MatNorm(corr_mat, NORM_INFINITY, &norms[4])); 319566063dSJacob Faibussowitsch PetscCall(MatNorm(corr_mat, NORM_FROBENIUS, &norms[5])); 329566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Autocorrelation matrix norms: %g %g %g\n", (double)norms[3],(double)norms[4],(double)norms[5])); 33f478cf82SStefano Zampini 349566063dSJacob Faibussowitsch PetscCall(MatDestroy(&corr_mat)); 35308ba433SVaclav Hapla PetscFunctionReturn(0); 36f478cf82SStefano Zampini } 37308ba433SVaclav Hapla 38308ba433SVaclav Hapla static PetscErrorCode GetReader(MPI_Comm comm, const char option[], PetscViewer *r, PetscViewerFormat *fmt) 39308ba433SVaclav Hapla { 40308ba433SVaclav Hapla PetscBool flg; 41308ba433SVaclav Hapla 42308ba433SVaclav Hapla PetscFunctionBegin; 439566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetViewer(PETSC_COMM_SELF, NULL, NULL, option, r, fmt, &flg)); 44308ba433SVaclav Hapla if (flg) { 45308ba433SVaclav Hapla PetscFileMode mode; 469566063dSJacob Faibussowitsch PetscCall(PetscViewerFileGetMode(*r, &mode)); 47308ba433SVaclav Hapla flg = (PetscBool) (mode == FILE_MODE_READ); 48308ba433SVaclav Hapla } 4928b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER_INPUT,"Need to specify %s viewer_type:file:format:read", option); 50308ba433SVaclav Hapla PetscFunctionReturn(0); 51308ba433SVaclav Hapla } 52308ba433SVaclav Hapla 53308ba433SVaclav Hapla int main(int argc, char **argv) 54308ba433SVaclav Hapla { 55308ba433SVaclav Hapla PetscInt i; 56308ba433SVaclav Hapla PetscReal norms0[NNORMS], norms1[NNORMS]; 57308ba433SVaclav Hapla PetscViewer inp_viewer; 58308ba433SVaclav Hapla PetscViewerFormat fmt; 59308ba433SVaclav Hapla Mat data_mat; 60308ba433SVaclav Hapla char mat_name[PETSC_MAX_PATH_LEN]="dmatrix"; 61308ba433SVaclav Hapla 62*327415f7SBarry Smith PetscFunctionBeginUser; 639566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, NULL)); 649566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-mat_name",mat_name,sizeof(mat_name),NULL)); 65308ba433SVaclav Hapla 66308ba433SVaclav Hapla /* load matrix sequentially */ 679566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &data_mat)); 689566063dSJacob Faibussowitsch PetscCall(MatSetType(data_mat,MATDENSE)); 699566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)data_mat, mat_name)); 709566063dSJacob Faibussowitsch PetscCall(GetReader(PETSC_COMM_SELF, "-serial_reader", &inp_viewer, &fmt)); 719566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(inp_viewer, fmt)); 729566063dSJacob Faibussowitsch PetscCall(MatLoadComputeNorms(data_mat, inp_viewer, norms0)); 739566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(inp_viewer)); 749566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&inp_viewer)); 759566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(data_mat, NULL, "-view_serial_mat")); 769566063dSJacob Faibussowitsch PetscCall(MatDestroy(&data_mat)); 77308ba433SVaclav Hapla 78308ba433SVaclav Hapla /* load matrix in parallel */ 799566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &data_mat)); 809566063dSJacob Faibussowitsch PetscCall(MatSetType(data_mat,MATDENSE)); 819566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)data_mat, mat_name)); 829566063dSJacob Faibussowitsch PetscCall(GetReader(PETSC_COMM_WORLD, "-parallel_reader", &inp_viewer, &fmt)); 839566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(inp_viewer, fmt)); 849566063dSJacob Faibussowitsch PetscCall(MatLoadComputeNorms(data_mat, inp_viewer, norms1)); 859566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(inp_viewer)); 869566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&inp_viewer)); 879566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(data_mat, NULL, "-view_parallel_mat")); 889566063dSJacob Faibussowitsch PetscCall(MatDestroy(&data_mat)); 89308ba433SVaclav Hapla 90308ba433SVaclav Hapla for (i=0; i<NNORMS; i++) { 9108401ef6SPierre Jolivet 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); 92308ba433SVaclav Hapla } 93308ba433SVaclav Hapla 949566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 95b122ec5aSJacob Faibussowitsch return 0; 96308ba433SVaclav Hapla } 97308ba433SVaclav Hapla 98308ba433SVaclav Hapla /*TEST 99308ba433SVaclav Hapla 100308ba433SVaclav Hapla test: 101308ba433SVaclav Hapla suffix: 1 102308ba433SVaclav Hapla requires: hdf5 datafilespath complex 103308ba433SVaclav Hapla args: -serial_reader hdf5:${DATAFILESPATH}/matrices/hdf5/sample_data.h5::read -parallel_reader hdf5:${DATAFILESPATH}/matrices/hdf5/sample_data.h5::read 104308ba433SVaclav Hapla nsize: {{1 2 4}} 105308ba433SVaclav Hapla 106308ba433SVaclav Hapla test: 107308ba433SVaclav Hapla requires: hdf5 datafilespath 108308ba433SVaclav Hapla args: -serial_reader hdf5:${DATAFILESPATH}/matrices/hdf5/tiny_rectangular_mat.h5::read -parallel_reader hdf5:${DATAFILESPATH}/matrices/hdf5/tiny_rectangular_mat.h5::read 109308ba433SVaclav Hapla nsize: {{1 2}} 110308ba433SVaclav Hapla test: 111308ba433SVaclav Hapla suffix: 2-complex 112308ba433SVaclav Hapla requires: complex 113308ba433SVaclav Hapla args: -mat_name ComplexMat 114308ba433SVaclav Hapla test: 115308ba433SVaclav Hapla suffix: 2-real 116308ba433SVaclav Hapla requires: !complex 117308ba433SVaclav Hapla args: -mat_name RealMat 118308ba433SVaclav Hapla 119308ba433SVaclav Hapla TEST*/ 120