156bd94f4SPierre Jolivet static char help[] = "Mini-app to benchmark matrix--matrix multiplication\n\n"; 256bd94f4SPierre Jolivet 356bd94f4SPierre Jolivet /* 456bd94f4SPierre Jolivet See the paper below for more information 556bd94f4SPierre Jolivet 656bd94f4SPierre Jolivet "KSPHPDDM and PCHPDDM: Extending PETSc with Robust Overlapping Schwarz Preconditioners and Advanced Krylov Methods", 756bd94f4SPierre Jolivet P. Jolivet, J. E. Roman, and S. Zampini (2020). 856bd94f4SPierre Jolivet */ 956bd94f4SPierre Jolivet 1056bd94f4SPierre Jolivet #include <petsc.h> 1156bd94f4SPierre Jolivet 12b88df2e7SBarry Smith #if defined(PETSC_HAVE_MKL_SPARSE_OPTIMIZE) 1356bd94f4SPierre Jolivet #include <mkl.h> 149371c9d4SSatish Balay #define PetscCallMKLSparse(func, args) \ 159371c9d4SSatish Balay do { \ 1656bd94f4SPierre Jolivet sparse_status_t __ierr; \ 17792fecdfSBarry Smith PetscStackPushExternal(#func); \ 1856bd94f4SPierre Jolivet __ierr = func args; \ 1956bd94f4SPierre Jolivet PetscStackPop; \ 2008401ef6SPierre Jolivet PetscCheck(__ierr == SPARSE_STATUS_SUCCESS, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in %s(): error code %d", #func, (int)__ierr); \ 2156bd94f4SPierre Jolivet } while (0) 2256bd94f4SPierre Jolivet #else 239371c9d4SSatish Balay #define PetscCallMKLSparse(func, args) \ 249371c9d4SSatish Balay do { SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No MKL support"); } while (0) 2556bd94f4SPierre Jolivet #endif 2656bd94f4SPierre Jolivet 279371c9d4SSatish Balay int main(int argc, char **argv) { 2856bd94f4SPierre Jolivet Mat A, C, D, E; 2956bd94f4SPierre Jolivet PetscInt nbs = 10, ntype = 10, nN = 8, m, M, trial = 5; 3056bd94f4SPierre Jolivet PetscViewer viewer; 3156bd94f4SPierre Jolivet PetscInt bs[10], N[8]; 3256bd94f4SPierre Jolivet char *type[10]; 3356bd94f4SPierre Jolivet PetscMPIInt size; 3456bd94f4SPierre Jolivet PetscBool flg, cuda, maij = PETSC_FALSE, check = PETSC_FALSE, trans = PETSC_FALSE, convert = PETSC_FALSE, mkl; 3556bd94f4SPierre Jolivet char file[PETSC_MAX_PATH_LEN]; 3656bd94f4SPierre Jolivet 37327415f7SBarry Smith PetscFunctionBeginUser; 389566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 399566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 40be096a46SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only"); 419566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-f", file, PETSC_MAX_PATH_LEN, &flg)); 4228b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER_INPUT, "Must indicate binary file with the -f option"); 439566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-trial", &trial, NULL)); 449566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-bs", bs, &nbs, &flg)); 4556bd94f4SPierre Jolivet if (!flg) { 4656bd94f4SPierre Jolivet nbs = 1; 4756bd94f4SPierre Jolivet bs[0] = 1; 4856bd94f4SPierre Jolivet } 499566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-N", N, &nN, &flg)); 5056bd94f4SPierre Jolivet if (!flg) { 5156bd94f4SPierre Jolivet nN = 8; 529371c9d4SSatish Balay N[0] = 1; 539371c9d4SSatish Balay N[1] = 2; 549371c9d4SSatish Balay N[2] = 4; 559371c9d4SSatish Balay N[3] = 8; 569371c9d4SSatish Balay N[4] = 16; 579371c9d4SSatish Balay N[5] = 32; 589371c9d4SSatish Balay N[6] = 64; 599371c9d4SSatish Balay N[7] = 128; 6056bd94f4SPierre Jolivet } 619566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetStringArray(NULL, NULL, "-type", type, &ntype, &flg)); 6256bd94f4SPierre Jolivet if (!flg) { 6356bd94f4SPierre Jolivet ntype = 1; 649566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATSEQAIJ, &type[0])); 6556bd94f4SPierre Jolivet } 669566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-check", &check, NULL)); 679566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-trans", &trans, NULL)); 689566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-convert_aij", &convert, NULL)); 6956bd94f4SPierre Jolivet for (PetscInt j = 0; j < nbs; ++j) { 709566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &viewer)); 719566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 729566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 739566063dSJacob Faibussowitsch PetscCall(MatLoad(A, viewer)); 749566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 759566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompareAny((PetscObject)A, &flg, MATSEQAIJ, MATMPIAIJ, "")); 7628b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER_INPUT, "Must indicate a MatAIJ input matrix"); 779566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &m, &M)); 7856bd94f4SPierre Jolivet if (m == M) { 7956bd94f4SPierre Jolivet Mat oA; 809566063dSJacob Faibussowitsch PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &oA)); 819566063dSJacob Faibussowitsch PetscCall(MatAXPY(A, 1.0, oA, DIFFERENT_NONZERO_PATTERN)); 829566063dSJacob Faibussowitsch PetscCall(MatDestroy(&oA)); 8356bd94f4SPierre Jolivet } 849566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, NULL)); 859566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, NULL)); 8656bd94f4SPierre Jolivet if (bs[j] > 1) { 8756bd94f4SPierre Jolivet Mat T, Tt, B; 8856bd94f4SPierre Jolivet const PetscScalar *ptr; 8956bd94f4SPierre Jolivet PetscScalar *val, *Aa; 9056bd94f4SPierre Jolivet const PetscInt *Ai, *Aj; 9156bd94f4SPierre Jolivet PetscInt An, i, k; 9256bd94f4SPierre Jolivet PetscBool done; 9356bd94f4SPierre Jolivet 949566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_SELF, bs[j], bs[j], bs[j], bs[j], NULL, &T)); 959566063dSJacob Faibussowitsch PetscCall(MatSetRandom(T, NULL)); 969566063dSJacob Faibussowitsch PetscCall(MatTranspose(T, MAT_INITIAL_MATRIX, &Tt)); 979566063dSJacob Faibussowitsch PetscCall(MatAXPY(T, 1.0, Tt, SAME_NONZERO_PATTERN)); 989566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Tt)); 999566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(T, &ptr)); 1009566063dSJacob Faibussowitsch PetscCall(MatGetRowIJ(A, 0, PETSC_FALSE, PETSC_FALSE, &An, &Ai, &Aj, &done)); 101e00437b9SBarry Smith PetscCheck(done && An == m, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Inconsistent sizes"); 1029566063dSJacob Faibussowitsch PetscCall(MatSeqAIJGetArray(A, &Aa)); 1039566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &B)); 1049566063dSJacob Faibussowitsch PetscCall(MatSetType(B, MATSEQBAIJ)); 1059566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, bs[j] * An, bs[j] * An, PETSC_DECIDE, PETSC_DECIDE)); 1069566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Ai[An] * bs[j] * bs[j], &val)); 10756bd94f4SPierre Jolivet for (i = 0; i < Ai[An]; ++i) 1089371c9d4SSatish Balay for (k = 0; k < bs[j] * bs[j]; ++k) val[i * bs[j] * bs[j] + k] = Aa[i] * ptr[k]; 1099566063dSJacob Faibussowitsch PetscCall(MatSetOption(B, MAT_ROW_ORIENTED, PETSC_FALSE)); 1109566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetPreallocationCSR(B, bs[j], Ai, Aj, val)); 1119566063dSJacob Faibussowitsch PetscCall(PetscFree(val)); 1129566063dSJacob Faibussowitsch PetscCall(MatSeqAIJRestoreArray(A, &Aa)); 1139566063dSJacob Faibussowitsch PetscCall(MatRestoreRowIJ(A, 0, PETSC_FALSE, PETSC_FALSE, &An, &Ai, &Aj, &done)); 1149566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(T, &ptr)); 1159566063dSJacob Faibussowitsch PetscCall(MatDestroy(&T)); 1169566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 11756bd94f4SPierre Jolivet A = B; 11856bd94f4SPierre Jolivet } 11956bd94f4SPierre Jolivet /* reconvert back to SeqAIJ before converting to the desired type later */ 12056bd94f4SPierre Jolivet if (!convert) E = A; 1219566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATSEQAIJ, convert ? MAT_INITIAL_MATRIX : MAT_INPLACE_MATRIX, &E)); 1229566063dSJacob Faibussowitsch PetscCall(MatSetOption(E, MAT_SYMMETRIC, PETSC_TRUE)); 12356bd94f4SPierre Jolivet for (PetscInt i = 0; i < ntype; ++i) { 12456bd94f4SPierre Jolivet char *tmp; 12556bd94f4SPierre Jolivet PetscInt *ia_ptr, *ja_ptr, k; 12656bd94f4SPierre Jolivet PetscScalar *a_ptr; 127b88df2e7SBarry Smith #if defined(PETSC_HAVE_MKL_SPARSE_OPTIMIZE) 12856bd94f4SPierre Jolivet struct matrix_descr descr; 12956bd94f4SPierre Jolivet sparse_matrix_t spr; 13056bd94f4SPierre Jolivet descr.type = SPARSE_MATRIX_TYPE_GENERAL; 13156bd94f4SPierre Jolivet descr.diag = SPARSE_DIAG_NON_UNIT; 13256bd94f4SPierre Jolivet #endif 133*48a46eb9SPierre Jolivet if (convert) PetscCall(MatDestroy(&A)); 1349566063dSJacob Faibussowitsch PetscCall(PetscStrstr(type[i], "mkl", &tmp)); 13556bd94f4SPierre Jolivet if (tmp) { 13656bd94f4SPierre Jolivet size_t mlen, tlen; 13756bd94f4SPierre Jolivet char base[256]; 13856bd94f4SPierre Jolivet 13956bd94f4SPierre Jolivet mkl = PETSC_TRUE; 1409566063dSJacob Faibussowitsch PetscCall(PetscStrlen(tmp, &mlen)); 1419566063dSJacob Faibussowitsch PetscCall(PetscStrlen(type[i], &tlen)); 1429566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(base, type[i], tlen - mlen + 1)); 1439566063dSJacob Faibussowitsch PetscCall(MatConvert(E, base, convert ? MAT_INITIAL_MATRIX : MAT_INPLACE_MATRIX, &A)); 14456bd94f4SPierre Jolivet } else { 14556bd94f4SPierre Jolivet mkl = PETSC_FALSE; 1469566063dSJacob Faibussowitsch PetscCall(PetscStrstr(type[i], "maij", &tmp)); 14756bd94f4SPierre Jolivet if (!tmp) { 1489566063dSJacob Faibussowitsch PetscCall(MatConvert(E, type[i], convert ? MAT_INITIAL_MATRIX : MAT_INPLACE_MATRIX, &A)); 14956bd94f4SPierre Jolivet } else { 1509566063dSJacob Faibussowitsch PetscCall(MatConvert(E, MATAIJ, convert ? MAT_INITIAL_MATRIX : MAT_INPLACE_MATRIX, &A)); 15156bd94f4SPierre Jolivet maij = PETSC_TRUE; 15256bd94f4SPierre Jolivet } 15356bd94f4SPierre Jolivet } 1549566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompareAny((PetscObject)A, &cuda, MATSEQAIJCUSPARSE, MATMPIAIJCUSPARSE, "")); 15556bd94f4SPierre Jolivet if (mkl) { 15656bd94f4SPierre Jolivet const PetscInt *Ai, *Aj; 15756bd94f4SPierre Jolivet PetscInt An; 15856bd94f4SPierre Jolivet PetscBool done; 15956bd94f4SPierre Jolivet 1609566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompareAny((PetscObject)A, &flg, MATSEQAIJ, MATSEQBAIJ, MATSEQSBAIJ, "")); 16128b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER_INPUT, "Not implemented"); 1629566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQAIJ, &flg)); 1639566063dSJacob Faibussowitsch PetscCall(MatGetRowIJ(A, 0, PETSC_FALSE, flg ? PETSC_FALSE : PETSC_TRUE, &An, &Ai, &Aj, &done)); 16428b400f6SJacob Faibussowitsch PetscCheck(done, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Inconsistent sizes"); 1659566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(An + 1, &ia_ptr)); 1669566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Ai[An], &ja_ptr)); 16756bd94f4SPierre Jolivet if (flg) { /* SeqAIJ */ 16856bd94f4SPierre Jolivet for (k = 0; k < An + 1; ++k) ia_ptr[k] = Ai[k]; 16956bd94f4SPierre Jolivet for (k = 0; k < Ai[An]; ++k) ja_ptr[k] = Aj[k]; 1709566063dSJacob Faibussowitsch PetscCall(MatSeqAIJGetArray(A, &a_ptr)); 171792fecdfSBarry Smith PetscCallMKLSparse(mkl_sparse_d_create_csr, (&spr, SPARSE_INDEX_BASE_ZERO, An, An, ia_ptr, ia_ptr + 1, ja_ptr, a_ptr)); 17256bd94f4SPierre Jolivet } else { 1739566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQBAIJ, &flg)); 17456bd94f4SPierre Jolivet if (flg) { 17556bd94f4SPierre Jolivet for (k = 0; k < An + 1; ++k) ia_ptr[k] = Ai[k] + 1; /* Fortran indexing to maximize cases covered by _mm routines */ 17656bd94f4SPierre Jolivet for (k = 0; k < Ai[An]; ++k) ja_ptr[k] = Aj[k] + 1; /* Fortran indexing to maximize cases covered by _mm routines */ 1779566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJGetArray(A, &a_ptr)); 178792fecdfSBarry Smith PetscCallMKLSparse(mkl_sparse_d_create_bsr, (&spr, SPARSE_INDEX_BASE_ONE, SPARSE_LAYOUT_COLUMN_MAJOR, An, An, bs[j], ia_ptr, ia_ptr + 1, ja_ptr, a_ptr)); 17956bd94f4SPierre Jolivet } else { 1809566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQSBAIJ, &flg)); 18156bd94f4SPierre Jolivet if (flg) { 18256bd94f4SPierre Jolivet for (k = 0; k < An + 1; ++k) ia_ptr[k] = Ai[k] + 1; /* Fortran indexing to maximize cases covered by _mm routines */ 18356bd94f4SPierre Jolivet for (k = 0; k < Ai[An]; ++k) ja_ptr[k] = Aj[k] + 1; /* Fortran indexing to maximize cases covered by _mm routines */ 1849566063dSJacob Faibussowitsch PetscCall(MatSeqSBAIJGetArray(A, &a_ptr)); 185792fecdfSBarry Smith PetscCallMKLSparse(mkl_sparse_d_create_bsr, (&spr, SPARSE_INDEX_BASE_ONE, SPARSE_LAYOUT_COLUMN_MAJOR, An, An, bs[j], ia_ptr, ia_ptr + 1, ja_ptr, a_ptr)); 186b88df2e7SBarry Smith #if defined(PETSC_HAVE_MKL_SPARSE_OPTIMIZE) 18756bd94f4SPierre Jolivet descr.type = SPARSE_MATRIX_TYPE_SYMMETRIC; 18856bd94f4SPierre Jolivet descr.mode = SPARSE_FILL_MODE_UPPER; 18956bd94f4SPierre Jolivet descr.diag = SPARSE_DIAG_NON_UNIT; 19056bd94f4SPierre Jolivet #endif 19156bd94f4SPierre Jolivet } 19256bd94f4SPierre Jolivet } 19356bd94f4SPierre Jolivet } 1949566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQAIJ, &flg)); 1959566063dSJacob Faibussowitsch PetscCall(MatRestoreRowIJ(A, 0, PETSC_FALSE, flg ? PETSC_FALSE : PETSC_TRUE, &An, &Ai, &Aj, &done)); 19656bd94f4SPierre Jolivet } 19756bd94f4SPierre Jolivet 1989566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A, NULL, "-A_view")); 19956bd94f4SPierre Jolivet 20056bd94f4SPierre Jolivet for (k = 0; k < nN; ++k) { 20156bd94f4SPierre Jolivet MatType Atype, Ctype; 20256bd94f4SPierre Jolivet PetscInt AM, AN, CM, CN, t; 203956f8c0dSBarry Smith #if defined(PETSC_USE_LOG) 20456bd94f4SPierre Jolivet PetscLogStage stage, tstage; 20556bd94f4SPierre Jolivet char stage_s[256]; 206956f8c0dSBarry Smith #endif 20756bd94f4SPierre Jolivet 2089566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD, bs[j] * m, PETSC_DECIDE, bs[j] * M, N[k], NULL, &C)); 2099566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD, bs[j] * m, PETSC_DECIDE, bs[j] * M, N[k], NULL, &D)); 2109566063dSJacob Faibussowitsch PetscCall(MatSetRandom(C, NULL)); 21156bd94f4SPierre Jolivet if (cuda) { /* convert to GPU if needed */ 2129566063dSJacob Faibussowitsch PetscCall(MatConvert(C, MATDENSECUDA, MAT_INPLACE_MATRIX, &C)); 2139566063dSJacob Faibussowitsch PetscCall(MatConvert(D, MATDENSECUDA, MAT_INPLACE_MATRIX, &D)); 21456bd94f4SPierre Jolivet } 21556bd94f4SPierre Jolivet if (mkl) { 216792fecdfSBarry Smith if (N[k] > 1) PetscCallMKLSparse(mkl_sparse_set_mm_hint, (spr, SPARSE_OPERATION_NON_TRANSPOSE, descr, SPARSE_LAYOUT_COLUMN_MAJOR, N[k], 1 + trial)); 217792fecdfSBarry Smith else PetscCallMKLSparse(mkl_sparse_set_mv_hint, (spr, SPARSE_OPERATION_NON_TRANSPOSE, descr, 1 + trial)); 218792fecdfSBarry Smith PetscCallMKLSparse(mkl_sparse_set_memory_hint, (spr, SPARSE_MEMORY_AGGRESSIVE)); 219792fecdfSBarry Smith PetscCallMKLSparse(mkl_sparse_optimize, (spr)); 22056bd94f4SPierre Jolivet } 2219566063dSJacob Faibussowitsch PetscCall(MatGetType(A, &Atype)); 2229566063dSJacob Faibussowitsch PetscCall(MatGetType(C, &Ctype)); 2239566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &AM, &AN)); 2249566063dSJacob Faibussowitsch PetscCall(MatGetSize(C, &CM, &CN)); 22556bd94f4SPierre Jolivet 226956f8c0dSBarry Smith #if defined(PETSC_USE_LOG) 22756bd94f4SPierre Jolivet if (!maij || N[k] > 1) { 2289566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(stage_s, sizeof(stage_s), "type_%s-bs_%" PetscInt_FMT "-N_%02d", type[i], bs[j], (int)N[k])); 2299566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister(stage_s, &stage)); 23056bd94f4SPierre Jolivet } 23156bd94f4SPierre Jolivet if (trans && N[k] > 1) { 2329566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(stage_s, sizeof(stage_s), "trans_type_%s-bs_%" PetscInt_FMT "-N_%02d", type[i], bs[j], (int)N[k])); 2339566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister(stage_s, &tstage)); 23456bd94f4SPierre Jolivet } 235956f8c0dSBarry Smith #endif 23656bd94f4SPierre Jolivet /* A*B */ 23756bd94f4SPierre Jolivet if (N[k] > 1) { 23856bd94f4SPierre Jolivet if (!maij) { 2399566063dSJacob Faibussowitsch PetscCall(MatProductCreateWithMat(A, C, NULL, D)); 2409566063dSJacob Faibussowitsch PetscCall(MatProductSetType(D, MATPRODUCT_AB)); 2419566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(D)); 2429566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(D)); 24356bd94f4SPierre Jolivet } 24456bd94f4SPierre Jolivet 24556bd94f4SPierre Jolivet if (!mkl) { 24656bd94f4SPierre Jolivet if (!maij) { 2479566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(D)); 2489566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Benchmarking MatProduct %s: with A %s %" PetscInt_FMT "x%" PetscInt_FMT " and B %s %" PetscInt_FMT "x%" PetscInt_FMT "\n", MatProductTypes[MATPRODUCT_AB], Atype, AM, AN, Ctype, CM, CN)); 2499566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage)); 250*48a46eb9SPierre Jolivet for (t = 0; t < trial; ++t) PetscCall(MatProductNumeric(D)); 2519566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 25256bd94f4SPierre Jolivet } else { 25356bd94f4SPierre Jolivet Mat E, Ct, Dt; 25456bd94f4SPierre Jolivet Vec cC, cD; 25556bd94f4SPierre Jolivet const PetscScalar *c_ptr; 25656bd94f4SPierre Jolivet PetscScalar *d_ptr; 2579566063dSJacob Faibussowitsch PetscCall(MatCreateMAIJ(A, N[k], &E)); 2589566063dSJacob Faibussowitsch PetscCall(MatDenseGetLocalMatrix(C, &Ct)); 2599566063dSJacob Faibussowitsch PetscCall(MatDenseGetLocalMatrix(D, &Dt)); 2609566063dSJacob Faibussowitsch PetscCall(MatTranspose(Ct, MAT_INPLACE_MATRIX, &Ct)); 2619566063dSJacob Faibussowitsch PetscCall(MatTranspose(Dt, MAT_INPLACE_MATRIX, &Dt)); 2629566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(Ct, &c_ptr)); 2639566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(Dt, &d_ptr)); 2649566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, AM * N[k], PETSC_DECIDE, c_ptr, &cC)); 2659566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, AM * N[k], PETSC_DECIDE, d_ptr, &cD)); 2669566063dSJacob Faibussowitsch PetscCall(MatMult(E, cC, cD)); 2679566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Benchmarking MatMult: with A %s %" PetscInt_FMT "x%" PetscInt_FMT " and B %s %" PetscInt_FMT "x%" PetscInt_FMT "\n", MATMAIJ, AM, AN, VECMPI, AM * N[k], 1)); 2689566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage)); 269*48a46eb9SPierre Jolivet for (t = 0; t < trial; ++t) PetscCall(MatMult(E, cC, cD)); 2709566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 2719566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cD)); 2729566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cC)); 2739566063dSJacob Faibussowitsch PetscCall(MatDestroy(&E)); 2749566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(Dt, &d_ptr)); 2759566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(Ct, &c_ptr)); 2769566063dSJacob Faibussowitsch PetscCall(MatTranspose(Ct, MAT_INPLACE_MATRIX, &Ct)); 2779566063dSJacob Faibussowitsch PetscCall(MatTranspose(Dt, MAT_INPLACE_MATRIX, &Dt)); 27856bd94f4SPierre Jolivet } 27956bd94f4SPierre Jolivet } else { 28056bd94f4SPierre Jolivet const PetscScalar *c_ptr; 28156bd94f4SPierre Jolivet PetscScalar *d_ptr; 28256bd94f4SPierre Jolivet 2839566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(C, &c_ptr)); 2849566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(D, &d_ptr)); 285792fecdfSBarry Smith PetscCallMKLSparse(mkl_sparse_d_mm, (SPARSE_OPERATION_NON_TRANSPOSE, 1.0, spr, descr, SPARSE_LAYOUT_COLUMN_MAJOR, c_ptr, CN, CM, 0.0, d_ptr, CM)); 2869566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Benchmarking mkl_sparse_d_mm (COLUMN_MAJOR): with A %s %" PetscInt_FMT "x%" PetscInt_FMT " and B %s %" PetscInt_FMT "x%" PetscInt_FMT "\n", Atype, AM, AN, Ctype, CM, CN)); 2879566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage)); 288*48a46eb9SPierre Jolivet for (t = 0; t < trial; ++t) PetscCallMKLSparse(mkl_sparse_d_mm, (SPARSE_OPERATION_NON_TRANSPOSE, 1.0, spr, descr, SPARSE_LAYOUT_COLUMN_MAJOR, c_ptr, CN, CM, 0.0, d_ptr, CM)); 2899566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 2909566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(D, &d_ptr)); 2919566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(C, &c_ptr)); 29256bd94f4SPierre Jolivet } 29356bd94f4SPierre Jolivet } else if (maij) { 2949566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 2959566063dSJacob Faibussowitsch PetscCall(MatDestroy(&D)); 29656bd94f4SPierre Jolivet continue; 29756bd94f4SPierre Jolivet } else if (!mkl) { 29856bd94f4SPierre Jolivet Vec cC, cD; 29956bd94f4SPierre Jolivet 3009566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecRead(C, 0, &cC)); 3019566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecWrite(D, 0, &cD)); 3029566063dSJacob Faibussowitsch PetscCall(MatMult(A, cC, cD)); 3039566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Benchmarking MatMult: with A %s %" PetscInt_FMT "x%" PetscInt_FMT "\n", Atype, AM, AN)); 3049566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage)); 305*48a46eb9SPierre Jolivet for (t = 0; t < trial; ++t) PetscCall(MatMult(A, cC, cD)); 3069566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 3079566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecRead(C, 0, &cC)); 3089566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecWrite(D, 0, &cD)); 30956bd94f4SPierre Jolivet } else { 31056bd94f4SPierre Jolivet const PetscScalar *c_ptr; 31156bd94f4SPierre Jolivet PetscScalar *d_ptr; 31256bd94f4SPierre Jolivet 3139566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(C, &c_ptr)); 3149566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(D, &d_ptr)); 3159566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Benchmarking mkl_sparse_d_mv: with A %s %" PetscInt_FMT "x%" PetscInt_FMT "\n", Atype, AM, AN)); 316792fecdfSBarry Smith PetscCallMKLSparse(mkl_sparse_d_mv, (SPARSE_OPERATION_NON_TRANSPOSE, 1.0, spr, descr, c_ptr, 0.0, d_ptr)); 3179566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage)); 318*48a46eb9SPierre Jolivet for (t = 0; t < trial; ++t) PetscCallMKLSparse(mkl_sparse_d_mv, (SPARSE_OPERATION_NON_TRANSPOSE, 1.0, spr, descr, c_ptr, 0.0, d_ptr)); 3199566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 3209566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(D, &d_ptr)); 3219566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(C, &c_ptr)); 32256bd94f4SPierre Jolivet } 32356bd94f4SPierre Jolivet 32456bd94f4SPierre Jolivet if (check) { 3259566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(A, C, D, 10, &flg)); 32656bd94f4SPierre Jolivet if (!flg) { 32756bd94f4SPierre Jolivet MatType Dtype; 32856bd94f4SPierre Jolivet 3299566063dSJacob Faibussowitsch PetscCall(MatGetType(D, &Dtype)); 3309566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error with A %s%s, C %s, D %s, Nk %" PetscInt_FMT "\n", Atype, mkl ? "mkl" : "", Ctype, Dtype, N[k])); 33156bd94f4SPierre Jolivet } 33256bd94f4SPierre Jolivet } 33356bd94f4SPierre Jolivet 33456bd94f4SPierre Jolivet /* MKL implementation seems buggy for ABt */ 33556bd94f4SPierre Jolivet /* A*Bt */ 33656bd94f4SPierre Jolivet if (!mkl && trans && N[k] > 1) { 3379566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompareAny((PetscObject)A, &flg, MATSEQAIJ, MATMPIAIJ, "")); 33856bd94f4SPierre Jolivet if (flg) { 3399566063dSJacob Faibussowitsch PetscCall(MatTranspose(C, MAT_INPLACE_MATRIX, &C)); 3409566063dSJacob Faibussowitsch PetscCall(MatGetType(C, &Ctype)); 34156bd94f4SPierre Jolivet if (!mkl) { 3429566063dSJacob Faibussowitsch PetscCall(MatProductCreateWithMat(A, C, NULL, D)); 3439566063dSJacob Faibussowitsch PetscCall(MatProductSetType(D, MATPRODUCT_ABt)); 3449566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(D)); 3459566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(D)); 3469566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(D)); 3479566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Benchmarking MatProduct %s: with A %s %" PetscInt_FMT "x%" PetscInt_FMT " and Bt %s %" PetscInt_FMT "x%" PetscInt_FMT "\n", MatProductTypes[MATPRODUCT_ABt], Atype, AM, AN, Ctype, CM, CN)); 3489566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(tstage)); 349*48a46eb9SPierre Jolivet for (t = 0; t < trial; ++t) PetscCall(MatProductNumeric(D)); 3509566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 35156bd94f4SPierre Jolivet } else { 35256bd94f4SPierre Jolivet const PetscScalar *c_ptr; 35356bd94f4SPierre Jolivet PetscScalar *d_ptr; 35456bd94f4SPierre Jolivet 355792fecdfSBarry Smith PetscCallMKLSparse(mkl_sparse_set_mm_hint, (spr, SPARSE_OPERATION_NON_TRANSPOSE, descr, SPARSE_LAYOUT_ROW_MAJOR, N[k], 1 + trial)); 356792fecdfSBarry Smith PetscCallMKLSparse(mkl_sparse_optimize, (spr)); 3579566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(C, &c_ptr)); 3589566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(D, &d_ptr)); 3599566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Benchmarking mkl_sparse_d_mm (ROW_MAJOR): with A %s %" PetscInt_FMT "x%" PetscInt_FMT " and B %s %" PetscInt_FMT "x%" PetscInt_FMT "\n", Atype, AM, AN, Ctype, CM, CN)); 360792fecdfSBarry Smith PetscCallMKLSparse(mkl_sparse_d_mm, (SPARSE_OPERATION_NON_TRANSPOSE, 1.0, spr, descr, SPARSE_LAYOUT_ROW_MAJOR, c_ptr, CN, CM, 0.0, d_ptr, CM)); 3619566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage)); 362*48a46eb9SPierre Jolivet for (t = 0; t < trial; ++t) PetscCallMKLSparse(mkl_sparse_d_mm, (SPARSE_OPERATION_NON_TRANSPOSE, 1.0, spr, descr, SPARSE_LAYOUT_ROW_MAJOR, c_ptr, CN, CM, 0.0, d_ptr, CM)); 3639566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 3649566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(D, &d_ptr)); 3659566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(C, &c_ptr)); 36656bd94f4SPierre Jolivet } 36756bd94f4SPierre Jolivet } 36856bd94f4SPierre Jolivet } 36956bd94f4SPierre Jolivet 37056bd94f4SPierre Jolivet if (!mkl && trans && N[k] > 1 && flg && check) { 3719566063dSJacob Faibussowitsch PetscCall(MatMatTransposeMultEqual(A, C, D, 10, &flg)); 37256bd94f4SPierre Jolivet if (!flg) { 37356bd94f4SPierre Jolivet MatType Dtype; 3749566063dSJacob Faibussowitsch PetscCall(MatGetType(D, &Dtype)); 3759566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error with A %s%s, C %s, D %s, Nk %" PetscInt_FMT "\n", Atype, mkl ? "mkl" : "", Ctype, Dtype, N[k])); 37656bd94f4SPierre Jolivet } 37756bd94f4SPierre Jolivet } 3789566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 3799566063dSJacob Faibussowitsch PetscCall(MatDestroy(&D)); 38056bd94f4SPierre Jolivet } 38156bd94f4SPierre Jolivet if (mkl) { 382792fecdfSBarry Smith PetscCallMKLSparse(mkl_sparse_destroy, (spr)); 3839566063dSJacob Faibussowitsch PetscCall(PetscFree(ia_ptr)); 3849566063dSJacob Faibussowitsch PetscCall(PetscFree(ja_ptr)); 38556bd94f4SPierre Jolivet } 38656bd94f4SPierre Jolivet if (cuda && i != ntype - 1) { 3879566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "AIJCUSPARSE must be last, otherwise MatConvert() to another MatType is too slow\n")); 38856bd94f4SPierre Jolivet break; 38956bd94f4SPierre Jolivet } 39056bd94f4SPierre Jolivet } 3919566063dSJacob Faibussowitsch if (E != A) PetscCall(MatDestroy(&E)); 3929566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 39356bd94f4SPierre Jolivet } 3949566063dSJacob Faibussowitsch for (m = 0; m < ntype; ++m) PetscCall(PetscFree(type[m])); 3959566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 39656bd94f4SPierre Jolivet return 0; 39756bd94f4SPierre Jolivet } 39856bd94f4SPierre Jolivet 39956bd94f4SPierre Jolivet /*TEST 400b8582c7cSSatish Balay build: 401dfd57a17SPierre Jolivet requires: double !complex !defined(PETSC_USE_64BIT_INDICES) 40256bd94f4SPierre Jolivet 40356bd94f4SPierre Jolivet testset: 40456bd94f4SPierre Jolivet nsize: 1 40556bd94f4SPierre Jolivet filter: sed "/Benchmarking/d" 40656bd94f4SPierre Jolivet args: -f ${wPETSC_DIR}/share/petsc/datafiles/matrices/spd-real-int32-float64 -bs 1,2,3 -N 1,2,18 -check -trans -convert_aij {{false true}shared output} 40756bd94f4SPierre Jolivet test: 40856bd94f4SPierre Jolivet suffix: basic 40956bd94f4SPierre Jolivet args: -type aij,sbaij,baij 41056bd94f4SPierre Jolivet output_file: output/ex237.out 41156bd94f4SPierre Jolivet test: 41256bd94f4SPierre Jolivet suffix: maij 41356bd94f4SPierre Jolivet args: -type aij,maij 41456bd94f4SPierre Jolivet output_file: output/ex237.out 41556bd94f4SPierre Jolivet test: 41656bd94f4SPierre Jolivet suffix: cuda 41756bd94f4SPierre Jolivet requires: cuda 41856bd94f4SPierre Jolivet args: -type aij,aijcusparse 41956bd94f4SPierre Jolivet output_file: output/ex237.out 42056bd94f4SPierre Jolivet test: 42156bd94f4SPierre Jolivet suffix: mkl 422b88df2e7SBarry Smith requires: mkl_sparse_optimize 42356bd94f4SPierre Jolivet args: -type aij,aijmkl,baijmkl,sbaijmkl 42456bd94f4SPierre Jolivet output_file: output/ex237.out 42556bd94f4SPierre Jolivet 42656bd94f4SPierre Jolivet TEST*/ 427