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) \ 24d71ae5a4SJacob Faibussowitsch do { \ 25d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No MKL support"); \ 26d71ae5a4SJacob Faibussowitsch } while (0) 2756bd94f4SPierre Jolivet #endif 2856bd94f4SPierre Jolivet 29d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 30d71ae5a4SJacob Faibussowitsch { 3156bd94f4SPierre Jolivet Mat A, C, D, E; 3256bd94f4SPierre Jolivet PetscInt nbs = 10, ntype = 10, nN = 8, m, M, trial = 5; 3356bd94f4SPierre Jolivet PetscViewer viewer; 3456bd94f4SPierre Jolivet PetscInt bs[10], N[8]; 3556bd94f4SPierre Jolivet char *type[10]; 3656bd94f4SPierre Jolivet PetscMPIInt size; 3756bd94f4SPierre Jolivet PetscBool flg, cuda, maij = PETSC_FALSE, check = PETSC_FALSE, trans = PETSC_FALSE, convert = PETSC_FALSE, mkl; 3856bd94f4SPierre Jolivet char file[PETSC_MAX_PATH_LEN]; 3956bd94f4SPierre Jolivet 40327415f7SBarry Smith PetscFunctionBeginUser; 419566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 429566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 43be096a46SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only"); 449566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-f", file, PETSC_MAX_PATH_LEN, &flg)); 4528b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER_INPUT, "Must indicate binary file with the -f option"); 469566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-trial", &trial, NULL)); 479566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-bs", bs, &nbs, &flg)); 4856bd94f4SPierre Jolivet if (!flg) { 4956bd94f4SPierre Jolivet nbs = 1; 5056bd94f4SPierre Jolivet bs[0] = 1; 5156bd94f4SPierre Jolivet } 529566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-N", N, &nN, &flg)); 5356bd94f4SPierre Jolivet if (!flg) { 5456bd94f4SPierre Jolivet nN = 8; 559371c9d4SSatish Balay N[0] = 1; 569371c9d4SSatish Balay N[1] = 2; 579371c9d4SSatish Balay N[2] = 4; 589371c9d4SSatish Balay N[3] = 8; 599371c9d4SSatish Balay N[4] = 16; 609371c9d4SSatish Balay N[5] = 32; 619371c9d4SSatish Balay N[6] = 64; 629371c9d4SSatish Balay N[7] = 128; 6356bd94f4SPierre Jolivet } 649566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetStringArray(NULL, NULL, "-type", type, &ntype, &flg)); 6556bd94f4SPierre Jolivet if (!flg) { 6656bd94f4SPierre Jolivet ntype = 1; 679566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATSEQAIJ, &type[0])); 6856bd94f4SPierre Jolivet } 699566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-check", &check, NULL)); 709566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-trans", &trans, NULL)); 719566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-convert_aij", &convert, NULL)); 7256bd94f4SPierre Jolivet for (PetscInt j = 0; j < nbs; ++j) { 739566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &viewer)); 749566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 759566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 769566063dSJacob Faibussowitsch PetscCall(MatLoad(A, viewer)); 779566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 789566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompareAny((PetscObject)A, &flg, MATSEQAIJ, MATMPIAIJ, "")); 7928b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER_INPUT, "Must indicate a MatAIJ input matrix"); 809566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &m, &M)); 8156bd94f4SPierre Jolivet if (m == M) { 8256bd94f4SPierre Jolivet Mat oA; 839566063dSJacob Faibussowitsch PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &oA)); 849566063dSJacob Faibussowitsch PetscCall(MatAXPY(A, 1.0, oA, DIFFERENT_NONZERO_PATTERN)); 859566063dSJacob Faibussowitsch PetscCall(MatDestroy(&oA)); 8656bd94f4SPierre Jolivet } 879566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, NULL)); 889566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, NULL)); 8956bd94f4SPierre Jolivet if (bs[j] > 1) { 9056bd94f4SPierre Jolivet Mat T, Tt, B; 9156bd94f4SPierre Jolivet const PetscScalar *ptr; 9256bd94f4SPierre Jolivet PetscScalar *val, *Aa; 9356bd94f4SPierre Jolivet const PetscInt *Ai, *Aj; 9456bd94f4SPierre Jolivet PetscInt An, i, k; 9556bd94f4SPierre Jolivet PetscBool done; 9656bd94f4SPierre Jolivet 979566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_SELF, bs[j], bs[j], bs[j], bs[j], NULL, &T)); 989566063dSJacob Faibussowitsch PetscCall(MatSetRandom(T, NULL)); 999566063dSJacob Faibussowitsch PetscCall(MatTranspose(T, MAT_INITIAL_MATRIX, &Tt)); 1009566063dSJacob Faibussowitsch PetscCall(MatAXPY(T, 1.0, Tt, SAME_NONZERO_PATTERN)); 1019566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Tt)); 1029566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(T, &ptr)); 1039566063dSJacob Faibussowitsch PetscCall(MatGetRowIJ(A, 0, PETSC_FALSE, PETSC_FALSE, &An, &Ai, &Aj, &done)); 104e00437b9SBarry Smith PetscCheck(done && An == m, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Inconsistent sizes"); 1059566063dSJacob Faibussowitsch PetscCall(MatSeqAIJGetArray(A, &Aa)); 1069566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &B)); 1079566063dSJacob Faibussowitsch PetscCall(MatSetType(B, MATSEQBAIJ)); 1089566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, bs[j] * An, bs[j] * An, PETSC_DECIDE, PETSC_DECIDE)); 1099566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Ai[An] * bs[j] * bs[j], &val)); 11056bd94f4SPierre Jolivet for (i = 0; i < Ai[An]; ++i) 1119371c9d4SSatish Balay for (k = 0; k < bs[j] * bs[j]; ++k) val[i * bs[j] * bs[j] + k] = Aa[i] * ptr[k]; 1129566063dSJacob Faibussowitsch PetscCall(MatSetOption(B, MAT_ROW_ORIENTED, PETSC_FALSE)); 1139566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetPreallocationCSR(B, bs[j], Ai, Aj, val)); 1149566063dSJacob Faibussowitsch PetscCall(PetscFree(val)); 1159566063dSJacob Faibussowitsch PetscCall(MatSeqAIJRestoreArray(A, &Aa)); 1169566063dSJacob Faibussowitsch PetscCall(MatRestoreRowIJ(A, 0, PETSC_FALSE, PETSC_FALSE, &An, &Ai, &Aj, &done)); 1179566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(T, &ptr)); 1189566063dSJacob Faibussowitsch PetscCall(MatDestroy(&T)); 1199566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 12056bd94f4SPierre Jolivet A = B; 12156bd94f4SPierre Jolivet } 12256bd94f4SPierre Jolivet /* reconvert back to SeqAIJ before converting to the desired type later */ 12356bd94f4SPierre Jolivet if (!convert) E = A; 1249566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATSEQAIJ, convert ? MAT_INITIAL_MATRIX : MAT_INPLACE_MATRIX, &E)); 1259566063dSJacob Faibussowitsch PetscCall(MatSetOption(E, MAT_SYMMETRIC, PETSC_TRUE)); 12656bd94f4SPierre Jolivet for (PetscInt i = 0; i < ntype; ++i) { 127*7864358aSSatish Balay char *tmp = NULL; 12856bd94f4SPierre Jolivet PetscInt *ia_ptr, *ja_ptr, k; 12956bd94f4SPierre Jolivet PetscScalar *a_ptr; 130b88df2e7SBarry Smith #if defined(PETSC_HAVE_MKL_SPARSE_OPTIMIZE) 13156bd94f4SPierre Jolivet struct matrix_descr descr; 13256bd94f4SPierre Jolivet sparse_matrix_t spr; 13356bd94f4SPierre Jolivet descr.type = SPARSE_MATRIX_TYPE_GENERAL; 13456bd94f4SPierre Jolivet descr.diag = SPARSE_DIAG_NON_UNIT; 13556bd94f4SPierre Jolivet #endif 13648a46eb9SPierre Jolivet if (convert) PetscCall(MatDestroy(&A)); 1379566063dSJacob Faibussowitsch PetscCall(PetscStrstr(type[i], "mkl", &tmp)); 13856bd94f4SPierre Jolivet if (tmp) { 13956bd94f4SPierre Jolivet size_t mlen, tlen; 14056bd94f4SPierre Jolivet char base[256]; 14156bd94f4SPierre Jolivet 14256bd94f4SPierre Jolivet mkl = PETSC_TRUE; 1439566063dSJacob Faibussowitsch PetscCall(PetscStrlen(tmp, &mlen)); 1449566063dSJacob Faibussowitsch PetscCall(PetscStrlen(type[i], &tlen)); 1459566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(base, type[i], tlen - mlen + 1)); 1469566063dSJacob Faibussowitsch PetscCall(MatConvert(E, base, convert ? MAT_INITIAL_MATRIX : MAT_INPLACE_MATRIX, &A)); 14756bd94f4SPierre Jolivet } else { 14856bd94f4SPierre Jolivet mkl = PETSC_FALSE; 1499566063dSJacob Faibussowitsch PetscCall(PetscStrstr(type[i], "maij", &tmp)); 15056bd94f4SPierre Jolivet if (!tmp) { 1519566063dSJacob Faibussowitsch PetscCall(MatConvert(E, type[i], convert ? MAT_INITIAL_MATRIX : MAT_INPLACE_MATRIX, &A)); 15256bd94f4SPierre Jolivet } else { 1539566063dSJacob Faibussowitsch PetscCall(MatConvert(E, MATAIJ, convert ? MAT_INITIAL_MATRIX : MAT_INPLACE_MATRIX, &A)); 15456bd94f4SPierre Jolivet maij = PETSC_TRUE; 15556bd94f4SPierre Jolivet } 15656bd94f4SPierre Jolivet } 1579566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompareAny((PetscObject)A, &cuda, MATSEQAIJCUSPARSE, MATMPIAIJCUSPARSE, "")); 15856bd94f4SPierre Jolivet if (mkl) { 15956bd94f4SPierre Jolivet const PetscInt *Ai, *Aj; 16056bd94f4SPierre Jolivet PetscInt An; 16156bd94f4SPierre Jolivet PetscBool done; 16256bd94f4SPierre Jolivet 1639566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompareAny((PetscObject)A, &flg, MATSEQAIJ, MATSEQBAIJ, MATSEQSBAIJ, "")); 16428b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER_INPUT, "Not implemented"); 1659566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQAIJ, &flg)); 1669566063dSJacob Faibussowitsch PetscCall(MatGetRowIJ(A, 0, PETSC_FALSE, flg ? PETSC_FALSE : PETSC_TRUE, &An, &Ai, &Aj, &done)); 16728b400f6SJacob Faibussowitsch PetscCheck(done, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Inconsistent sizes"); 1689566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(An + 1, &ia_ptr)); 1699566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Ai[An], &ja_ptr)); 17056bd94f4SPierre Jolivet if (flg) { /* SeqAIJ */ 17156bd94f4SPierre Jolivet for (k = 0; k < An + 1; ++k) ia_ptr[k] = Ai[k]; 17256bd94f4SPierre Jolivet for (k = 0; k < Ai[An]; ++k) ja_ptr[k] = Aj[k]; 1739566063dSJacob Faibussowitsch PetscCall(MatSeqAIJGetArray(A, &a_ptr)); 174792fecdfSBarry Smith PetscCallMKLSparse(mkl_sparse_d_create_csr, (&spr, SPARSE_INDEX_BASE_ZERO, An, An, ia_ptr, ia_ptr + 1, ja_ptr, a_ptr)); 17556bd94f4SPierre Jolivet } else { 1769566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQBAIJ, &flg)); 17756bd94f4SPierre Jolivet if (flg) { 17856bd94f4SPierre Jolivet for (k = 0; k < An + 1; ++k) ia_ptr[k] = Ai[k] + 1; /* Fortran indexing to maximize cases covered by _mm routines */ 17956bd94f4SPierre Jolivet for (k = 0; k < Ai[An]; ++k) ja_ptr[k] = Aj[k] + 1; /* Fortran indexing to maximize cases covered by _mm routines */ 1809566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJGetArray(A, &a_ptr)); 181792fecdfSBarry 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)); 18256bd94f4SPierre Jolivet } else { 1839566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQSBAIJ, &flg)); 18456bd94f4SPierre Jolivet if (flg) { 18556bd94f4SPierre Jolivet for (k = 0; k < An + 1; ++k) ia_ptr[k] = Ai[k] + 1; /* Fortran indexing to maximize cases covered by _mm routines */ 18656bd94f4SPierre Jolivet for (k = 0; k < Ai[An]; ++k) ja_ptr[k] = Aj[k] + 1; /* Fortran indexing to maximize cases covered by _mm routines */ 1879566063dSJacob Faibussowitsch PetscCall(MatSeqSBAIJGetArray(A, &a_ptr)); 188792fecdfSBarry 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)); 189b88df2e7SBarry Smith #if defined(PETSC_HAVE_MKL_SPARSE_OPTIMIZE) 19056bd94f4SPierre Jolivet descr.type = SPARSE_MATRIX_TYPE_SYMMETRIC; 19156bd94f4SPierre Jolivet descr.mode = SPARSE_FILL_MODE_UPPER; 19256bd94f4SPierre Jolivet descr.diag = SPARSE_DIAG_NON_UNIT; 19356bd94f4SPierre Jolivet #endif 19456bd94f4SPierre Jolivet } 19556bd94f4SPierre Jolivet } 19656bd94f4SPierre Jolivet } 1979566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQAIJ, &flg)); 1989566063dSJacob Faibussowitsch PetscCall(MatRestoreRowIJ(A, 0, PETSC_FALSE, flg ? PETSC_FALSE : PETSC_TRUE, &An, &Ai, &Aj, &done)); 19956bd94f4SPierre Jolivet } 20056bd94f4SPierre Jolivet 2019566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A, NULL, "-A_view")); 20256bd94f4SPierre Jolivet 20356bd94f4SPierre Jolivet for (k = 0; k < nN; ++k) { 20456bd94f4SPierre Jolivet MatType Atype, Ctype; 20556bd94f4SPierre Jolivet PetscInt AM, AN, CM, CN, t; 206956f8c0dSBarry Smith #if defined(PETSC_USE_LOG) 20756bd94f4SPierre Jolivet PetscLogStage stage, tstage; 20856bd94f4SPierre Jolivet char stage_s[256]; 209956f8c0dSBarry Smith #endif 21056bd94f4SPierre Jolivet 2119566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD, bs[j] * m, PETSC_DECIDE, bs[j] * M, N[k], NULL, &C)); 2129566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD, bs[j] * m, PETSC_DECIDE, bs[j] * M, N[k], NULL, &D)); 2139566063dSJacob Faibussowitsch PetscCall(MatSetRandom(C, NULL)); 21456bd94f4SPierre Jolivet if (cuda) { /* convert to GPU if needed */ 2159566063dSJacob Faibussowitsch PetscCall(MatConvert(C, MATDENSECUDA, MAT_INPLACE_MATRIX, &C)); 2169566063dSJacob Faibussowitsch PetscCall(MatConvert(D, MATDENSECUDA, MAT_INPLACE_MATRIX, &D)); 21756bd94f4SPierre Jolivet } 21856bd94f4SPierre Jolivet if (mkl) { 219792fecdfSBarry 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)); 220792fecdfSBarry Smith else PetscCallMKLSparse(mkl_sparse_set_mv_hint, (spr, SPARSE_OPERATION_NON_TRANSPOSE, descr, 1 + trial)); 221792fecdfSBarry Smith PetscCallMKLSparse(mkl_sparse_set_memory_hint, (spr, SPARSE_MEMORY_AGGRESSIVE)); 222792fecdfSBarry Smith PetscCallMKLSparse(mkl_sparse_optimize, (spr)); 22356bd94f4SPierre Jolivet } 2249566063dSJacob Faibussowitsch PetscCall(MatGetType(A, &Atype)); 2259566063dSJacob Faibussowitsch PetscCall(MatGetType(C, &Ctype)); 2269566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &AM, &AN)); 2279566063dSJacob Faibussowitsch PetscCall(MatGetSize(C, &CM, &CN)); 22856bd94f4SPierre Jolivet 229956f8c0dSBarry Smith #if defined(PETSC_USE_LOG) 23056bd94f4SPierre Jolivet if (!maij || N[k] > 1) { 2319566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(stage_s, sizeof(stage_s), "type_%s-bs_%" PetscInt_FMT "-N_%02d", type[i], bs[j], (int)N[k])); 2329566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister(stage_s, &stage)); 23356bd94f4SPierre Jolivet } 23456bd94f4SPierre Jolivet if (trans && N[k] > 1) { 2359566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(stage_s, sizeof(stage_s), "trans_type_%s-bs_%" PetscInt_FMT "-N_%02d", type[i], bs[j], (int)N[k])); 2369566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister(stage_s, &tstage)); 23756bd94f4SPierre Jolivet } 238956f8c0dSBarry Smith #endif 23956bd94f4SPierre Jolivet /* A*B */ 24056bd94f4SPierre Jolivet if (N[k] > 1) { 24156bd94f4SPierre Jolivet if (!maij) { 2429566063dSJacob Faibussowitsch PetscCall(MatProductCreateWithMat(A, C, NULL, D)); 2439566063dSJacob Faibussowitsch PetscCall(MatProductSetType(D, MATPRODUCT_AB)); 2449566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(D)); 2459566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(D)); 24656bd94f4SPierre Jolivet } 24756bd94f4SPierre Jolivet 24856bd94f4SPierre Jolivet if (!mkl) { 24956bd94f4SPierre Jolivet if (!maij) { 2509566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(D)); 2519566063dSJacob 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)); 2529566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage)); 25348a46eb9SPierre Jolivet for (t = 0; t < trial; ++t) PetscCall(MatProductNumeric(D)); 2549566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 25556bd94f4SPierre Jolivet } else { 25656bd94f4SPierre Jolivet Mat E, Ct, Dt; 25756bd94f4SPierre Jolivet Vec cC, cD; 25856bd94f4SPierre Jolivet const PetscScalar *c_ptr; 25956bd94f4SPierre Jolivet PetscScalar *d_ptr; 2609566063dSJacob Faibussowitsch PetscCall(MatCreateMAIJ(A, N[k], &E)); 2619566063dSJacob Faibussowitsch PetscCall(MatDenseGetLocalMatrix(C, &Ct)); 2629566063dSJacob Faibussowitsch PetscCall(MatDenseGetLocalMatrix(D, &Dt)); 2639566063dSJacob Faibussowitsch PetscCall(MatTranspose(Ct, MAT_INPLACE_MATRIX, &Ct)); 2649566063dSJacob Faibussowitsch PetscCall(MatTranspose(Dt, MAT_INPLACE_MATRIX, &Dt)); 2659566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(Ct, &c_ptr)); 2669566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(Dt, &d_ptr)); 2679566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, AM * N[k], PETSC_DECIDE, c_ptr, &cC)); 2689566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, AM * N[k], PETSC_DECIDE, d_ptr, &cD)); 2699566063dSJacob Faibussowitsch PetscCall(MatMult(E, cC, cD)); 2709566063dSJacob 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)); 2719566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage)); 27248a46eb9SPierre Jolivet for (t = 0; t < trial; ++t) PetscCall(MatMult(E, cC, cD)); 2739566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 2749566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cD)); 2759566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cC)); 2769566063dSJacob Faibussowitsch PetscCall(MatDestroy(&E)); 2779566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(Dt, &d_ptr)); 2789566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(Ct, &c_ptr)); 2799566063dSJacob Faibussowitsch PetscCall(MatTranspose(Ct, MAT_INPLACE_MATRIX, &Ct)); 2809566063dSJacob Faibussowitsch PetscCall(MatTranspose(Dt, MAT_INPLACE_MATRIX, &Dt)); 28156bd94f4SPierre Jolivet } 28256bd94f4SPierre Jolivet } else { 28356bd94f4SPierre Jolivet const PetscScalar *c_ptr; 28456bd94f4SPierre Jolivet PetscScalar *d_ptr; 28556bd94f4SPierre Jolivet 2869566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(C, &c_ptr)); 2879566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(D, &d_ptr)); 288792fecdfSBarry 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)); 2899566063dSJacob 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)); 2909566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage)); 29148a46eb9SPierre 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)); 2929566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 2939566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(D, &d_ptr)); 2949566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(C, &c_ptr)); 29556bd94f4SPierre Jolivet } 29656bd94f4SPierre Jolivet } else if (maij) { 2979566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 2989566063dSJacob Faibussowitsch PetscCall(MatDestroy(&D)); 29956bd94f4SPierre Jolivet continue; 30056bd94f4SPierre Jolivet } else if (!mkl) { 30156bd94f4SPierre Jolivet Vec cC, cD; 30256bd94f4SPierre Jolivet 3039566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecRead(C, 0, &cC)); 3049566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecWrite(D, 0, &cD)); 3059566063dSJacob Faibussowitsch PetscCall(MatMult(A, cC, cD)); 3069566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Benchmarking MatMult: with A %s %" PetscInt_FMT "x%" PetscInt_FMT "\n", Atype, AM, AN)); 3079566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage)); 30848a46eb9SPierre Jolivet for (t = 0; t < trial; ++t) PetscCall(MatMult(A, cC, cD)); 3099566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 3109566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecRead(C, 0, &cC)); 3119566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecWrite(D, 0, &cD)); 31256bd94f4SPierre Jolivet } else { 31356bd94f4SPierre Jolivet const PetscScalar *c_ptr; 31456bd94f4SPierre Jolivet PetscScalar *d_ptr; 31556bd94f4SPierre Jolivet 3169566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(C, &c_ptr)); 3179566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(D, &d_ptr)); 3189566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Benchmarking mkl_sparse_d_mv: with A %s %" PetscInt_FMT "x%" PetscInt_FMT "\n", Atype, AM, AN)); 319792fecdfSBarry Smith PetscCallMKLSparse(mkl_sparse_d_mv, (SPARSE_OPERATION_NON_TRANSPOSE, 1.0, spr, descr, c_ptr, 0.0, d_ptr)); 3209566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage)); 32148a46eb9SPierre 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)); 3229566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 3239566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(D, &d_ptr)); 3249566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(C, &c_ptr)); 32556bd94f4SPierre Jolivet } 32656bd94f4SPierre Jolivet 32756bd94f4SPierre Jolivet if (check) { 3289566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(A, C, D, 10, &flg)); 32956bd94f4SPierre Jolivet if (!flg) { 33056bd94f4SPierre Jolivet MatType Dtype; 33156bd94f4SPierre Jolivet 3329566063dSJacob Faibussowitsch PetscCall(MatGetType(D, &Dtype)); 3339566063dSJacob 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])); 33456bd94f4SPierre Jolivet } 33556bd94f4SPierre Jolivet } 33656bd94f4SPierre Jolivet 33756bd94f4SPierre Jolivet /* MKL implementation seems buggy for ABt */ 33856bd94f4SPierre Jolivet /* A*Bt */ 33956bd94f4SPierre Jolivet if (!mkl && trans && N[k] > 1) { 3409566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompareAny((PetscObject)A, &flg, MATSEQAIJ, MATMPIAIJ, "")); 34156bd94f4SPierre Jolivet if (flg) { 3429566063dSJacob Faibussowitsch PetscCall(MatTranspose(C, MAT_INPLACE_MATRIX, &C)); 3439566063dSJacob Faibussowitsch PetscCall(MatGetType(C, &Ctype)); 34456bd94f4SPierre Jolivet if (!mkl) { 3459566063dSJacob Faibussowitsch PetscCall(MatProductCreateWithMat(A, C, NULL, D)); 3469566063dSJacob Faibussowitsch PetscCall(MatProductSetType(D, MATPRODUCT_ABt)); 3479566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(D)); 3489566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(D)); 3499566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(D)); 3509566063dSJacob 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)); 3519566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(tstage)); 35248a46eb9SPierre Jolivet for (t = 0; t < trial; ++t) PetscCall(MatProductNumeric(D)); 3539566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 35456bd94f4SPierre Jolivet } else { 35556bd94f4SPierre Jolivet const PetscScalar *c_ptr; 35656bd94f4SPierre Jolivet PetscScalar *d_ptr; 35756bd94f4SPierre Jolivet 358792fecdfSBarry Smith PetscCallMKLSparse(mkl_sparse_set_mm_hint, (spr, SPARSE_OPERATION_NON_TRANSPOSE, descr, SPARSE_LAYOUT_ROW_MAJOR, N[k], 1 + trial)); 359792fecdfSBarry Smith PetscCallMKLSparse(mkl_sparse_optimize, (spr)); 3609566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(C, &c_ptr)); 3619566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(D, &d_ptr)); 3629566063dSJacob 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)); 363792fecdfSBarry 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)); 3649566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage)); 36548a46eb9SPierre 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)); 3669566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 3679566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(D, &d_ptr)); 3689566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(C, &c_ptr)); 36956bd94f4SPierre Jolivet } 37056bd94f4SPierre Jolivet } 37156bd94f4SPierre Jolivet } 37256bd94f4SPierre Jolivet 37356bd94f4SPierre Jolivet if (!mkl && trans && N[k] > 1 && flg && check) { 3749566063dSJacob Faibussowitsch PetscCall(MatMatTransposeMultEqual(A, C, D, 10, &flg)); 37556bd94f4SPierre Jolivet if (!flg) { 37656bd94f4SPierre Jolivet MatType Dtype; 3779566063dSJacob Faibussowitsch PetscCall(MatGetType(D, &Dtype)); 3789566063dSJacob 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])); 37956bd94f4SPierre Jolivet } 38056bd94f4SPierre Jolivet } 3819566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 3829566063dSJacob Faibussowitsch PetscCall(MatDestroy(&D)); 38356bd94f4SPierre Jolivet } 38456bd94f4SPierre Jolivet if (mkl) { 385792fecdfSBarry Smith PetscCallMKLSparse(mkl_sparse_destroy, (spr)); 3869566063dSJacob Faibussowitsch PetscCall(PetscFree(ia_ptr)); 3879566063dSJacob Faibussowitsch PetscCall(PetscFree(ja_ptr)); 38856bd94f4SPierre Jolivet } 38956bd94f4SPierre Jolivet if (cuda && i != ntype - 1) { 3909566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "AIJCUSPARSE must be last, otherwise MatConvert() to another MatType is too slow\n")); 39156bd94f4SPierre Jolivet break; 39256bd94f4SPierre Jolivet } 39356bd94f4SPierre Jolivet } 3949566063dSJacob Faibussowitsch if (E != A) PetscCall(MatDestroy(&E)); 3959566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 39656bd94f4SPierre Jolivet } 3979566063dSJacob Faibussowitsch for (m = 0; m < ntype; ++m) PetscCall(PetscFree(type[m])); 3989566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 39956bd94f4SPierre Jolivet return 0; 40056bd94f4SPierre Jolivet } 40156bd94f4SPierre Jolivet 40256bd94f4SPierre Jolivet /*TEST 403b8582c7cSSatish Balay build: 404dfd57a17SPierre Jolivet requires: double !complex !defined(PETSC_USE_64BIT_INDICES) 40556bd94f4SPierre Jolivet 40656bd94f4SPierre Jolivet testset: 40756bd94f4SPierre Jolivet nsize: 1 40856bd94f4SPierre Jolivet filter: sed "/Benchmarking/d" 40956bd94f4SPierre 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} 41056bd94f4SPierre Jolivet test: 41156bd94f4SPierre Jolivet suffix: basic 41256bd94f4SPierre Jolivet args: -type aij,sbaij,baij 41356bd94f4SPierre Jolivet output_file: output/ex237.out 41456bd94f4SPierre Jolivet test: 41556bd94f4SPierre Jolivet suffix: maij 41656bd94f4SPierre Jolivet args: -type aij,maij 41756bd94f4SPierre Jolivet output_file: output/ex237.out 41856bd94f4SPierre Jolivet test: 41956bd94f4SPierre Jolivet suffix: cuda 42056bd94f4SPierre Jolivet requires: cuda 42156bd94f4SPierre Jolivet args: -type aij,aijcusparse 42256bd94f4SPierre Jolivet output_file: output/ex237.out 42356bd94f4SPierre Jolivet test: 42456bd94f4SPierre Jolivet suffix: mkl 425b88df2e7SBarry Smith requires: mkl_sparse_optimize 42656bd94f4SPierre Jolivet args: -type aij,aijmkl,baijmkl,sbaijmkl 42756bd94f4SPierre Jolivet output_file: output/ex237.out 42856bd94f4SPierre Jolivet 42956bd94f4SPierre Jolivet TEST*/ 430