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> 1456bd94f4SPierre Jolivet #define PetscStackCallMKLSparse(func, args) do { \ 1556bd94f4SPierre Jolivet sparse_status_t __ierr; \ 1656bd94f4SPierre Jolivet PetscStackPush(#func); \ 1756bd94f4SPierre Jolivet __ierr = func args; \ 1856bd94f4SPierre Jolivet PetscStackPop; \ 192c71b3e2SJacob Faibussowitsch PetscCheckFalse(__ierr != SPARSE_STATUS_SUCCESS,PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in %s(): error code %d", #func, (int)__ierr); \ 2056bd94f4SPierre Jolivet } while (0) 2156bd94f4SPierre Jolivet #else 2256bd94f4SPierre Jolivet #define PetscStackCallMKLSparse(func, args) do { \ 2356bd94f4SPierre Jolivet SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No MKL support"); \ 2456bd94f4SPierre Jolivet } while (0) 2556bd94f4SPierre Jolivet #endif 2656bd94f4SPierre Jolivet 27956f8c0dSBarry Smith int main(int argc, char** argv) 28956f8c0dSBarry Smith { 2956bd94f4SPierre Jolivet Mat A, C, D, E; 3056bd94f4SPierre Jolivet PetscInt nbs = 10, ntype = 10, nN = 8, m, M, trial = 5; 3156bd94f4SPierre Jolivet PetscViewer viewer; 3256bd94f4SPierre Jolivet PetscInt bs[10], N[8]; 3356bd94f4SPierre Jolivet char *type[10]; 3456bd94f4SPierre Jolivet PetscMPIInt size; 3556bd94f4SPierre Jolivet PetscBool flg, cuda, maij = PETSC_FALSE, check = PETSC_FALSE, trans = PETSC_FALSE, convert = PETSC_FALSE, mkl; 3656bd94f4SPierre Jolivet char file[PETSC_MAX_PATH_LEN]; 3756bd94f4SPierre Jolivet PetscErrorCode ierr; 3856bd94f4SPierre Jolivet 392da392ccSBarry Smith ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 4056bd94f4SPierre Jolivet if (ierr) return ierr; 415f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 422c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 1,PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only"); 435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(NULL, NULL, "-f", file, PETSC_MAX_PATH_LEN, &flg)); 44*28b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD, PETSC_ERR_USER_INPUT, "Must indicate binary file with the -f option"); 455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL, NULL, "-trial", &trial, NULL)); 465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetIntArray(NULL, NULL, "-bs", bs, &nbs, &flg)); 4756bd94f4SPierre Jolivet if (!flg) { 4856bd94f4SPierre Jolivet nbs = 1; 4956bd94f4SPierre Jolivet bs[0] = 1; 5056bd94f4SPierre Jolivet } 515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetIntArray(NULL, NULL, "-N", N, &nN, &flg)); 5256bd94f4SPierre Jolivet if (!flg) { 5356bd94f4SPierre Jolivet nN = 8; 5456bd94f4SPierre Jolivet N[0] = 1; N[1] = 2; N[2] = 4; N[3] = 8; 5556bd94f4SPierre Jolivet N[4] = 16; N[5] = 32; N[6] = 64; N[7] = 128; 5656bd94f4SPierre Jolivet } 575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetStringArray(NULL, NULL, "-type", type, &ntype, &flg)); 5856bd94f4SPierre Jolivet if (!flg) { 5956bd94f4SPierre Jolivet ntype = 1; 605f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrallocpy(MATSEQAIJ, &type[0])); 6156bd94f4SPierre Jolivet } 625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL, NULL, "-check", &check, NULL)); 635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL, NULL, "-trans", &trans, NULL)); 645f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL, NULL, "-convert_aij", &convert, NULL)); 6556bd94f4SPierre Jolivet for (PetscInt j = 0; j < nbs; ++j) { 665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &viewer)); 675f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD, &A)); 685f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 695f80ce2aSJacob Faibussowitsch CHKERRQ(MatLoad(A, viewer)); 705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&viewer)); 715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompareAny((PetscObject)A, &flg, MATSEQAIJ, MATMPIAIJ, "")); 72*28b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD, PETSC_ERR_USER_INPUT, "Must indicate a MatAIJ input matrix"); 735f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(A, &m, &M)); 7456bd94f4SPierre Jolivet if (m == M) { 7556bd94f4SPierre Jolivet Mat oA; 765f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(A, MAT_INITIAL_MATRIX, &oA)); 775f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(A, 1.0, oA, DIFFERENT_NONZERO_PATTERN)); 785f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&oA)); 7956bd94f4SPierre Jolivet } 805f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(A, &m, NULL)); 815f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(A, &M, NULL)); 8256bd94f4SPierre Jolivet if (bs[j] > 1) { 8356bd94f4SPierre Jolivet Mat T, Tt, B; 8456bd94f4SPierre Jolivet const PetscScalar *ptr; 8556bd94f4SPierre Jolivet PetscScalar *val, *Aa; 8656bd94f4SPierre Jolivet const PetscInt *Ai, *Aj; 8756bd94f4SPierre Jolivet PetscInt An, i, k; 8856bd94f4SPierre Jolivet PetscBool done; 8956bd94f4SPierre Jolivet 905f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateDense(PETSC_COMM_SELF, bs[j], bs[j], bs[j], bs[j], NULL, &T)); 915f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetRandom(T, NULL)); 925f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(T, MAT_INITIAL_MATRIX, &Tt)); 935f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(T, 1.0, Tt, SAME_NONZERO_PATTERN)); 945f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Tt)); 955f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetArrayRead(T, &ptr)); 965f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowIJ(A, 0, PETSC_FALSE, PETSC_FALSE, &An, &Ai, &Aj, &done)); 972c71b3e2SJacob Faibussowitsch PetscCheckFalse(!done || An != m,PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Inconsistent sizes"); 985f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJGetArray(A, &Aa)); 995f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD, &B)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(B, MATSEQBAIJ)); 1015f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(B, bs[j] * An, bs[j] * An, PETSC_DECIDE, PETSC_DECIDE)); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(Ai[An] * bs[j] * bs[j], &val)); 10356bd94f4SPierre Jolivet for (i = 0; i < Ai[An]; ++i) 10456bd94f4SPierre Jolivet for (k = 0; k < bs[j] * bs[j]; ++k) 10556bd94f4SPierre Jolivet val[i * bs[j] * bs[j] + k] = Aa[i] * ptr[k]; 1065f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(B, MAT_ROW_ORIENTED, PETSC_FALSE)); 1075f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqBAIJSetPreallocationCSR(B, bs[j], Ai, Aj, val)); 1085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(val)); 1095f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJRestoreArray(A, &Aa)); 1105f80ce2aSJacob Faibussowitsch CHKERRQ(MatRestoreRowIJ(A, 0, PETSC_FALSE, PETSC_FALSE, &An, &Ai, &Aj, &done)); 1115f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreArrayRead(T, &ptr)); 1125f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&T)); 1135f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 11456bd94f4SPierre Jolivet A = B; 11556bd94f4SPierre Jolivet } 11656bd94f4SPierre Jolivet /* reconvert back to SeqAIJ before converting to the desired type later */ 11756bd94f4SPierre Jolivet if (!convert) E = A; 1185f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A, MATSEQAIJ, convert ? MAT_INITIAL_MATRIX : MAT_INPLACE_MATRIX, &E)); 1195f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(E, MAT_SYMMETRIC, PETSC_TRUE)); 12056bd94f4SPierre Jolivet for (PetscInt i = 0; i < ntype; ++i) { 12156bd94f4SPierre Jolivet char *tmp; 12256bd94f4SPierre Jolivet PetscInt *ia_ptr, *ja_ptr, k; 12356bd94f4SPierre Jolivet PetscScalar *a_ptr; 124b88df2e7SBarry Smith #if defined(PETSC_HAVE_MKL_SPARSE_OPTIMIZE) 12556bd94f4SPierre Jolivet struct matrix_descr descr; 12656bd94f4SPierre Jolivet sparse_matrix_t spr; 12756bd94f4SPierre Jolivet descr.type = SPARSE_MATRIX_TYPE_GENERAL; 12856bd94f4SPierre Jolivet descr.diag = SPARSE_DIAG_NON_UNIT; 12956bd94f4SPierre Jolivet #endif 13056bd94f4SPierre Jolivet if (convert) { 1315f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 13256bd94f4SPierre Jolivet } 1335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrstr(type[i], "mkl", &tmp)); 13456bd94f4SPierre Jolivet if (tmp) { 13556bd94f4SPierre Jolivet size_t mlen, tlen; 13656bd94f4SPierre Jolivet char base[256]; 13756bd94f4SPierre Jolivet 13856bd94f4SPierre Jolivet mkl = PETSC_TRUE; 1395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrlen(tmp, &mlen)); 1405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrlen(type[i], &tlen)); 1415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrncpy(base, type[i], tlen-mlen + 1)); 1425f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(E, base, convert ? MAT_INITIAL_MATRIX : MAT_INPLACE_MATRIX, &A)); 14356bd94f4SPierre Jolivet } else { 14456bd94f4SPierre Jolivet mkl = PETSC_FALSE; 1455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrstr(type[i], "maij", &tmp)); 14656bd94f4SPierre Jolivet if (!tmp) { 1475f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(E, type[i], convert ? MAT_INITIAL_MATRIX : MAT_INPLACE_MATRIX, &A)); 14856bd94f4SPierre Jolivet } else { 1495f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(E, MATAIJ, convert ? MAT_INITIAL_MATRIX : MAT_INPLACE_MATRIX, &A)); 15056bd94f4SPierre Jolivet maij = PETSC_TRUE; 15156bd94f4SPierre Jolivet } 15256bd94f4SPierre Jolivet } 1535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompareAny((PetscObject)A, &cuda, MATSEQAIJCUSPARSE, MATMPIAIJCUSPARSE, "")); 15456bd94f4SPierre Jolivet if (mkl) { 15556bd94f4SPierre Jolivet const PetscInt *Ai, *Aj; 15656bd94f4SPierre Jolivet PetscInt An; 15756bd94f4SPierre Jolivet PetscBool done; 15856bd94f4SPierre Jolivet 1595f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompareAny((PetscObject)A, &flg, MATSEQAIJ, MATSEQBAIJ, MATSEQSBAIJ, "")); 160*28b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD, PETSC_ERR_USER_INPUT, "Not implemented"); 1615f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)A, MATSEQAIJ, &flg)); 1625f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowIJ(A, 0, PETSC_FALSE, flg ? PETSC_FALSE : PETSC_TRUE, &An, &Ai, &Aj, &done)); 163*28b400f6SJacob Faibussowitsch PetscCheck(done,PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Inconsistent sizes"); 1645f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(An + 1, &ia_ptr)); 1655f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(Ai[An], &ja_ptr)); 16656bd94f4SPierre Jolivet if (flg) { /* SeqAIJ */ 16756bd94f4SPierre Jolivet for (k = 0; k < An + 1; ++k) ia_ptr[k] = Ai[k]; 16856bd94f4SPierre Jolivet for (k = 0; k < Ai[An]; ++k) ja_ptr[k] = Aj[k]; 1695f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJGetArray(A, &a_ptr)); 17056bd94f4SPierre Jolivet PetscStackCallMKLSparse(mkl_sparse_d_create_csr, (&spr, SPARSE_INDEX_BASE_ZERO, An, An, ia_ptr, ia_ptr + 1, ja_ptr, a_ptr)); 17156bd94f4SPierre Jolivet } else { 1725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)A, MATSEQBAIJ, &flg)); 17356bd94f4SPierre Jolivet if (flg) { 17456bd94f4SPierre Jolivet for (k = 0; k < An + 1; ++k) ia_ptr[k] = Ai[k] + 1; /* Fortran indexing to maximize cases covered by _mm routines */ 17556bd94f4SPierre Jolivet for (k = 0; k < Ai[An]; ++k) ja_ptr[k] = Aj[k] + 1; /* Fortran indexing to maximize cases covered by _mm routines */ 1765f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqBAIJGetArray(A, &a_ptr)); 17756bd94f4SPierre Jolivet PetscStackCallMKLSparse(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)); 17856bd94f4SPierre Jolivet } else { 1795f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)A, MATSEQSBAIJ, &flg)); 18056bd94f4SPierre Jolivet if (flg) { 18156bd94f4SPierre Jolivet for (k = 0; k < An + 1; ++k) ia_ptr[k] = Ai[k] + 1; /* Fortran indexing to maximize cases covered by _mm routines */ 18256bd94f4SPierre Jolivet for (k = 0; k < Ai[An]; ++k) ja_ptr[k] = Aj[k] + 1; /* Fortran indexing to maximize cases covered by _mm routines */ 1835f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqSBAIJGetArray(A, &a_ptr)); 18456bd94f4SPierre Jolivet PetscStackCallMKLSparse(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)); 185b88df2e7SBarry Smith #if defined(PETSC_HAVE_MKL_SPARSE_OPTIMIZE) 18656bd94f4SPierre Jolivet descr.type = SPARSE_MATRIX_TYPE_SYMMETRIC; 18756bd94f4SPierre Jolivet descr.mode = SPARSE_FILL_MODE_UPPER; 18856bd94f4SPierre Jolivet descr.diag = SPARSE_DIAG_NON_UNIT; 18956bd94f4SPierre Jolivet #endif 19056bd94f4SPierre Jolivet } 19156bd94f4SPierre Jolivet } 19256bd94f4SPierre Jolivet } 1935f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)A, MATSEQAIJ, &flg)); 1945f80ce2aSJacob Faibussowitsch CHKERRQ(MatRestoreRowIJ(A, 0, PETSC_FALSE, flg ? PETSC_FALSE : PETSC_TRUE, &An, &Ai, &Aj, &done)); 19556bd94f4SPierre Jolivet } 19656bd94f4SPierre Jolivet 1975f80ce2aSJacob Faibussowitsch CHKERRQ(MatViewFromOptions(A, NULL, "-A_view")); 19856bd94f4SPierre Jolivet 19956bd94f4SPierre Jolivet for (k = 0; k < nN; ++k) { 20056bd94f4SPierre Jolivet MatType Atype, Ctype; 20156bd94f4SPierre Jolivet PetscInt AM, AN, CM, CN, t; 202956f8c0dSBarry Smith #if defined(PETSC_USE_LOG) 20356bd94f4SPierre Jolivet PetscLogStage stage, tstage; 20456bd94f4SPierre Jolivet char stage_s[256]; 205956f8c0dSBarry Smith #endif 20656bd94f4SPierre Jolivet 2075f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateDense(PETSC_COMM_WORLD, bs[j] * m, PETSC_DECIDE, bs[j] * M, N[k], NULL, &C)); 2085f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateDense(PETSC_COMM_WORLD, bs[j] * m, PETSC_DECIDE, bs[j] * M, N[k], NULL, &D)); 2095f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetRandom(C, NULL)); 21056bd94f4SPierre Jolivet if (cuda) { /* convert to GPU if needed */ 2115f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(C, MATDENSECUDA, MAT_INPLACE_MATRIX, &C)); 2125f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(D, MATDENSECUDA, MAT_INPLACE_MATRIX, &D)); 21356bd94f4SPierre Jolivet } 21456bd94f4SPierre Jolivet if (mkl) { 21556bd94f4SPierre Jolivet if (N[k] > 1) PetscStackCallMKLSparse(mkl_sparse_set_mm_hint, (spr, SPARSE_OPERATION_NON_TRANSPOSE, descr, SPARSE_LAYOUT_COLUMN_MAJOR, N[k], 1 + trial)); 21656bd94f4SPierre Jolivet else PetscStackCallMKLSparse(mkl_sparse_set_mv_hint, (spr, SPARSE_OPERATION_NON_TRANSPOSE, descr, 1 + trial)); 21756bd94f4SPierre Jolivet PetscStackCallMKLSparse(mkl_sparse_set_memory_hint, (spr, SPARSE_MEMORY_AGGRESSIVE)); 21856bd94f4SPierre Jolivet PetscStackCallMKLSparse(mkl_sparse_optimize, (spr)); 21956bd94f4SPierre Jolivet } 2205f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetType(A, &Atype)); 2215f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetType(C, &Ctype)); 2225f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(A, &AM, &AN)); 2235f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(C, &CM, &CN)); 22456bd94f4SPierre Jolivet 225956f8c0dSBarry Smith #if defined(PETSC_USE_LOG) 22656bd94f4SPierre Jolivet if (!maij || N[k] > 1) { 2275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSNPrintf(stage_s, sizeof(stage_s), "type_%s-bs_%" PetscInt_FMT "-N_%02d", type[i], bs[j], (int)N[k])); 2285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStageRegister(stage_s, &stage)); 22956bd94f4SPierre Jolivet } 23056bd94f4SPierre Jolivet if (trans && N[k] > 1) { 2315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSNPrintf(stage_s, sizeof(stage_s), "trans_type_%s-bs_%" PetscInt_FMT "-N_%02d", type[i], bs[j], (int)N[k])); 2325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStageRegister(stage_s, &tstage)); 23356bd94f4SPierre Jolivet } 234956f8c0dSBarry Smith #endif 23556bd94f4SPierre Jolivet /* A*B */ 23656bd94f4SPierre Jolivet if (N[k] > 1) { 23756bd94f4SPierre Jolivet if (!maij) { 2385f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductCreateWithMat(A, C, NULL, D)); 2395f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSetType(D, MATPRODUCT_AB)); 2405f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSetFromOptions(D)); 2415f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSymbolic(D)); 24256bd94f4SPierre Jolivet } 24356bd94f4SPierre Jolivet 24456bd94f4SPierre Jolivet if (!mkl) { 24556bd94f4SPierre Jolivet if (!maij) { 2465f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductNumeric(D)); 2475f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 2485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePush(stage)); 24956bd94f4SPierre Jolivet for (t = 0; t < trial; ++t) { 2505f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductNumeric(D)); 25156bd94f4SPierre Jolivet } 2525f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePop()); 25356bd94f4SPierre Jolivet } else { 25456bd94f4SPierre Jolivet Mat E, Ct, Dt; 25556bd94f4SPierre Jolivet Vec cC, cD; 25656bd94f4SPierre Jolivet const PetscScalar *c_ptr; 25756bd94f4SPierre Jolivet PetscScalar *d_ptr; 2585f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateMAIJ(A, N[k], &E)); 2595f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetLocalMatrix(C, &Ct)); 2605f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetLocalMatrix(D, &Dt)); 2615f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(Ct, MAT_INPLACE_MATRIX, &Ct)); 2625f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(Dt, MAT_INPLACE_MATRIX, &Dt)); 2635f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetArrayRead(Ct, &c_ptr)); 2645f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetArrayWrite(Dt, &d_ptr)); 2655f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, AM * N[k], PETSC_DECIDE, c_ptr, &cC)); 2665f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, AM * N[k], PETSC_DECIDE, d_ptr, &cD)); 2675f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(E, cC, cD)); 2685f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 2695f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePush(stage)); 27056bd94f4SPierre Jolivet for (t = 0; t < trial; ++t) { 2715f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(E, cC, cD)); 27256bd94f4SPierre Jolivet } 2735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePop()); 2745f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&cD)); 2755f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&cC)); 2765f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&E)); 2775f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreArrayWrite(Dt, &d_ptr)); 2785f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreArrayRead(Ct, &c_ptr)); 2795f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(Ct, MAT_INPLACE_MATRIX, &Ct)); 2805f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(Dt, MAT_INPLACE_MATRIX, &Dt)); 28156bd94f4SPierre Jolivet } 28256bd94f4SPierre Jolivet } else { 28356bd94f4SPierre Jolivet const PetscScalar *c_ptr; 28456bd94f4SPierre Jolivet PetscScalar *d_ptr; 28556bd94f4SPierre Jolivet 2865f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetArrayRead(C, &c_ptr)); 2875f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetArrayWrite(D, &d_ptr)); 28856bd94f4SPierre Jolivet PetscStackCallMKLSparse(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)); 2895f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 2905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePush(stage)); 29156bd94f4SPierre Jolivet for (t = 0; t < trial; ++t) { 29256bd94f4SPierre Jolivet PetscStackCallMKLSparse(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)); 29356bd94f4SPierre Jolivet } 2945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePop()); 2955f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreArrayWrite(D, &d_ptr)); 2965f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreArrayRead(C, &c_ptr)); 29756bd94f4SPierre Jolivet } 29856bd94f4SPierre Jolivet } else if (maij) { 2995f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 3005f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&D)); 30156bd94f4SPierre Jolivet continue; 30256bd94f4SPierre Jolivet } else if (!mkl) { 30356bd94f4SPierre Jolivet Vec cC, cD; 30456bd94f4SPierre Jolivet 3055f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetColumnVecRead(C, 0, &cC)); 3065f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetColumnVecWrite(D, 0, &cD)); 3075f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A, cC, cD)); 3085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "Benchmarking MatMult: with A %s %" PetscInt_FMT "x%" PetscInt_FMT "\n", Atype, AM, AN)); 3095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePush(stage)); 31056bd94f4SPierre Jolivet for (t = 0; t < trial; ++t) { 3115f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A, cC, cD)); 31256bd94f4SPierre Jolivet } 3135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePop()); 3145f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreColumnVecRead(C, 0, &cC)); 3155f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreColumnVecWrite(D, 0, &cD)); 31656bd94f4SPierre Jolivet } else { 31756bd94f4SPierre Jolivet const PetscScalar *c_ptr; 31856bd94f4SPierre Jolivet PetscScalar *d_ptr; 31956bd94f4SPierre Jolivet 3205f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetArrayRead(C, &c_ptr)); 3215f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetArrayWrite(D, &d_ptr)); 3225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "Benchmarking mkl_sparse_d_mv: with A %s %" PetscInt_FMT "x%" PetscInt_FMT "\n", Atype, AM, AN)); 32356bd94f4SPierre Jolivet PetscStackCallMKLSparse(mkl_sparse_d_mv, (SPARSE_OPERATION_NON_TRANSPOSE, 1.0, spr, descr, c_ptr, 0.0, d_ptr)); 3245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePush(stage)); 32556bd94f4SPierre Jolivet for (t = 0; t < trial; ++t) { 32656bd94f4SPierre Jolivet PetscStackCallMKLSparse(mkl_sparse_d_mv, (SPARSE_OPERATION_NON_TRANSPOSE, 1.0, spr, descr, c_ptr, 0.0, d_ptr)); 32756bd94f4SPierre Jolivet } 3285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePop()); 3295f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreArrayWrite(D, &d_ptr)); 3305f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreArrayRead(C, &c_ptr)); 33156bd94f4SPierre Jolivet } 33256bd94f4SPierre Jolivet 33356bd94f4SPierre Jolivet if (check) { 3345f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMultEqual(A, C, D, 10, &flg)); 33556bd94f4SPierre Jolivet if (!flg) { 33656bd94f4SPierre Jolivet MatType Dtype; 33756bd94f4SPierre Jolivet 3385f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetType(D, &Dtype)); 3395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "Error with A %s%s, C %s, D %s, Nk %" PetscInt_FMT "\n", Atype, mkl ? "mkl" : "", Ctype, Dtype, N[k])); 34056bd94f4SPierre Jolivet } 34156bd94f4SPierre Jolivet } 34256bd94f4SPierre Jolivet 34356bd94f4SPierre Jolivet /* MKL implementation seems buggy for ABt */ 34456bd94f4SPierre Jolivet /* A*Bt */ 34556bd94f4SPierre Jolivet if (!mkl && trans && N[k] > 1) { 3465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompareAny((PetscObject)A, &flg, MATSEQAIJ, MATMPIAIJ, "")); 34756bd94f4SPierre Jolivet if (flg) { 3485f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(C, MAT_INPLACE_MATRIX, &C)); 3495f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetType(C, &Ctype)); 35056bd94f4SPierre Jolivet if (!mkl) { 3515f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductCreateWithMat(A, C, NULL, D)); 3525f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSetType(D, MATPRODUCT_ABt)); 3535f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSetFromOptions(D)); 3545f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSymbolic(D)); 3555f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductNumeric(D)); 3565f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 3575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePush(tstage)); 35856bd94f4SPierre Jolivet for (t = 0; t < trial; ++t) { 3595f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductNumeric(D)); 36056bd94f4SPierre Jolivet } 3615f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePop()); 36256bd94f4SPierre Jolivet } else { 36356bd94f4SPierre Jolivet const PetscScalar *c_ptr; 36456bd94f4SPierre Jolivet PetscScalar *d_ptr; 36556bd94f4SPierre Jolivet 36656bd94f4SPierre Jolivet PetscStackCallMKLSparse(mkl_sparse_set_mm_hint, (spr, SPARSE_OPERATION_NON_TRANSPOSE, descr, SPARSE_LAYOUT_ROW_MAJOR, N[k], 1 + trial)); 36756bd94f4SPierre Jolivet PetscStackCallMKLSparse(mkl_sparse_optimize, (spr)); 3685f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetArrayRead(C, &c_ptr)); 3695f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetArrayWrite(D, &d_ptr)); 3705f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 37156bd94f4SPierre Jolivet PetscStackCallMKLSparse(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)); 3725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePush(stage)); 37356bd94f4SPierre Jolivet for (t = 0; t < trial; ++t) { 37456bd94f4SPierre Jolivet PetscStackCallMKLSparse(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)); 37556bd94f4SPierre Jolivet } 3765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePop()); 3775f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreArrayWrite(D, &d_ptr)); 3785f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreArrayRead(C, &c_ptr)); 37956bd94f4SPierre Jolivet } 38056bd94f4SPierre Jolivet } 38156bd94f4SPierre Jolivet } 38256bd94f4SPierre Jolivet 38356bd94f4SPierre Jolivet if (!mkl && trans && N[k] > 1 && flg && check) { 3845f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatTransposeMultEqual(A, C, D, 10, &flg)); 38556bd94f4SPierre Jolivet if (!flg) { 38656bd94f4SPierre Jolivet MatType Dtype; 3875f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetType(D, &Dtype)); 3885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "Error with A %s%s, C %s, D %s, Nk %" PetscInt_FMT "\n", Atype, mkl ? "mkl" : "", Ctype, Dtype, N[k])); 38956bd94f4SPierre Jolivet } 39056bd94f4SPierre Jolivet } 3915f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 3925f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&D)); 39356bd94f4SPierre Jolivet } 39456bd94f4SPierre Jolivet if (mkl) { 39556bd94f4SPierre Jolivet PetscStackCallMKLSparse(mkl_sparse_destroy, (spr)); 3965f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(ia_ptr)); 3975f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(ja_ptr)); 39856bd94f4SPierre Jolivet } 39956bd94f4SPierre Jolivet if (cuda && i != ntype - 1) { 4005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "AIJCUSPARSE must be last, otherwise MatConvert() to another MatType is too slow\n")); 40156bd94f4SPierre Jolivet break; 40256bd94f4SPierre Jolivet } 40356bd94f4SPierre Jolivet } 40456bd94f4SPierre Jolivet if (E != A) { 4055f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&E)); 40656bd94f4SPierre Jolivet } 4075f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 40856bd94f4SPierre Jolivet } 40956bd94f4SPierre Jolivet for (m = 0; m < ntype; ++m) { 4105f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(type[m])); 41156bd94f4SPierre Jolivet } 41256bd94f4SPierre Jolivet ierr = PetscFinalize(); 41356bd94f4SPierre Jolivet return 0; 41456bd94f4SPierre Jolivet } 41556bd94f4SPierre Jolivet 41656bd94f4SPierre Jolivet /*TEST 417b8582c7cSSatish Balay build: 418dfd57a17SPierre Jolivet requires: double !complex !defined(PETSC_USE_64BIT_INDICES) 41956bd94f4SPierre Jolivet 42056bd94f4SPierre Jolivet testset: 42156bd94f4SPierre Jolivet nsize: 1 42256bd94f4SPierre Jolivet filter: sed "/Benchmarking/d" 42356bd94f4SPierre 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} 42456bd94f4SPierre Jolivet test: 42556bd94f4SPierre Jolivet suffix: basic 42656bd94f4SPierre Jolivet args: -type aij,sbaij,baij 42756bd94f4SPierre Jolivet output_file: output/ex237.out 42856bd94f4SPierre Jolivet test: 42956bd94f4SPierre Jolivet suffix: maij 43056bd94f4SPierre Jolivet args: -type aij,maij 43156bd94f4SPierre Jolivet output_file: output/ex237.out 43256bd94f4SPierre Jolivet test: 43356bd94f4SPierre Jolivet suffix: cuda 43456bd94f4SPierre Jolivet requires: cuda 43556bd94f4SPierre Jolivet args: -type aij,aijcusparse 43656bd94f4SPierre Jolivet output_file: output/ex237.out 43756bd94f4SPierre Jolivet test: 43856bd94f4SPierre Jolivet suffix: mkl 439b88df2e7SBarry Smith requires: mkl_sparse_optimize 44056bd94f4SPierre Jolivet args: -type aij,aijmkl,baijmkl,sbaijmkl 44156bd94f4SPierre Jolivet output_file: output/ex237.out 44256bd94f4SPierre Jolivet 44356bd94f4SPierre Jolivet TEST*/ 444