xref: /petsc/src/mat/tests/ex237.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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>
14792fecdfSBarry Smith #define PetscCallMKLSparse(func, args) do {               \
1556bd94f4SPierre Jolivet     sparse_status_t __ierr;                                    \
16792fecdfSBarry Smith     PetscStackPushExternal(#func);                                     \
1756bd94f4SPierre Jolivet     __ierr = func args;                                        \
1856bd94f4SPierre Jolivet     PetscStackPop;                                             \
1908401ef6SPierre Jolivet     PetscCheck(__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
22792fecdfSBarry Smith #define PetscCallMKLSparse(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 
38*327415f7SBarry Smith   PetscFunctionBeginUser;
399566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
409566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
41be096a46SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only");
429566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-f", file, PETSC_MAX_PATH_LEN, &flg));
4328b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_WORLD, PETSC_ERR_USER_INPUT, "Must indicate binary file with the -f option");
449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-trial", &trial, NULL));
459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-bs", bs, &nbs, &flg));
4656bd94f4SPierre Jolivet   if (!flg) {
4756bd94f4SPierre Jolivet     nbs = 1;
4856bd94f4SPierre Jolivet     bs[0] = 1;
4956bd94f4SPierre Jolivet   }
509566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-N", N, &nN, &flg));
5156bd94f4SPierre Jolivet   if (!flg) {
5256bd94f4SPierre Jolivet     nN = 8;
5356bd94f4SPierre Jolivet     N[0] = 1;  N[1] = 2;  N[2] = 4;  N[3] = 8;
5456bd94f4SPierre Jolivet     N[4] = 16; N[5] = 32; N[6] = 64; N[7] = 128;
5556bd94f4SPierre Jolivet   }
569566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetStringArray(NULL, NULL, "-type", type, &ntype, &flg));
5756bd94f4SPierre Jolivet   if (!flg) {
5856bd94f4SPierre Jolivet     ntype = 1;
599566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(MATSEQAIJ, &type[0]));
6056bd94f4SPierre Jolivet   }
619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-check", &check, NULL));
629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-trans", &trans, NULL));
639566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-convert_aij", &convert, NULL));
6456bd94f4SPierre Jolivet   for (PetscInt j = 0; j < nbs; ++j) {
659566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &viewer));
669566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
679566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(A));
689566063dSJacob Faibussowitsch     PetscCall(MatLoad(A, viewer));
699566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
709566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompareAny((PetscObject)A, &flg, MATSEQAIJ, MATMPIAIJ, ""));
7128b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_WORLD, PETSC_ERR_USER_INPUT, "Must indicate a MatAIJ input matrix");
729566063dSJacob Faibussowitsch     PetscCall(MatGetSize(A, &m, &M));
7356bd94f4SPierre Jolivet     if (m == M) {
7456bd94f4SPierre Jolivet       Mat oA;
759566063dSJacob Faibussowitsch       PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &oA));
769566063dSJacob Faibussowitsch       PetscCall(MatAXPY(A, 1.0, oA, DIFFERENT_NONZERO_PATTERN));
779566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&oA));
7856bd94f4SPierre Jolivet     }
799566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(A, &m, NULL));
809566063dSJacob Faibussowitsch     PetscCall(MatGetSize(A, &M, NULL));
8156bd94f4SPierre Jolivet     if (bs[j] > 1) {
8256bd94f4SPierre Jolivet       Mat               T, Tt, B;
8356bd94f4SPierre Jolivet       const PetscScalar *ptr;
8456bd94f4SPierre Jolivet       PetscScalar       *val, *Aa;
8556bd94f4SPierre Jolivet       const PetscInt    *Ai, *Aj;
8656bd94f4SPierre Jolivet       PetscInt          An, i, k;
8756bd94f4SPierre Jolivet       PetscBool         done;
8856bd94f4SPierre Jolivet 
899566063dSJacob Faibussowitsch       PetscCall(MatCreateDense(PETSC_COMM_SELF, bs[j], bs[j], bs[j], bs[j], NULL, &T));
909566063dSJacob Faibussowitsch       PetscCall(MatSetRandom(T, NULL));
919566063dSJacob Faibussowitsch       PetscCall(MatTranspose(T, MAT_INITIAL_MATRIX, &Tt));
929566063dSJacob Faibussowitsch       PetscCall(MatAXPY(T, 1.0, Tt, SAME_NONZERO_PATTERN));
939566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&Tt));
949566063dSJacob Faibussowitsch       PetscCall(MatDenseGetArrayRead(T, &ptr));
959566063dSJacob Faibussowitsch       PetscCall(MatGetRowIJ(A, 0, PETSC_FALSE, PETSC_FALSE, &An, &Ai, &Aj, &done));
96e00437b9SBarry Smith       PetscCheck(done && An == m,PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Inconsistent sizes");
979566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJGetArray(A, &Aa));
989566063dSJacob Faibussowitsch       PetscCall(MatCreate(PETSC_COMM_WORLD, &B));
999566063dSJacob Faibussowitsch       PetscCall(MatSetType(B, MATSEQBAIJ));
1009566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(B, bs[j] * An, bs[j] * An, PETSC_DECIDE, PETSC_DECIDE));
1019566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(Ai[An] * bs[j] * bs[j], &val));
10256bd94f4SPierre Jolivet       for (i = 0; i < Ai[An]; ++i)
10356bd94f4SPierre Jolivet         for (k = 0; k < bs[j] * bs[j]; ++k)
10456bd94f4SPierre Jolivet           val[i * bs[j] * bs[j] + k] = Aa[i] * ptr[k];
1059566063dSJacob Faibussowitsch       PetscCall(MatSetOption(B, MAT_ROW_ORIENTED, PETSC_FALSE));
1069566063dSJacob Faibussowitsch       PetscCall(MatSeqBAIJSetPreallocationCSR(B, bs[j], Ai, Aj, val));
1079566063dSJacob Faibussowitsch       PetscCall(PetscFree(val));
1089566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJRestoreArray(A, &Aa));
1099566063dSJacob Faibussowitsch       PetscCall(MatRestoreRowIJ(A, 0, PETSC_FALSE, PETSC_FALSE, &An, &Ai, &Aj, &done));
1109566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreArrayRead(T, &ptr));
1119566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&T));
1129566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&A));
11356bd94f4SPierre Jolivet       A    = B;
11456bd94f4SPierre Jolivet     }
11556bd94f4SPierre Jolivet     /* reconvert back to SeqAIJ before converting to the desired type later */
11656bd94f4SPierre Jolivet     if (!convert) E = A;
1179566063dSJacob Faibussowitsch     PetscCall(MatConvert(A, MATSEQAIJ, convert ? MAT_INITIAL_MATRIX : MAT_INPLACE_MATRIX, &E));
1189566063dSJacob Faibussowitsch     PetscCall(MatSetOption(E, MAT_SYMMETRIC, PETSC_TRUE));
11956bd94f4SPierre Jolivet     for (PetscInt i = 0; i < ntype; ++i) {
12056bd94f4SPierre Jolivet       char        *tmp;
12156bd94f4SPierre Jolivet       PetscInt    *ia_ptr, *ja_ptr, k;
12256bd94f4SPierre Jolivet       PetscScalar *a_ptr;
123b88df2e7SBarry Smith #if defined(PETSC_HAVE_MKL_SPARSE_OPTIMIZE)
12456bd94f4SPierre Jolivet       struct matrix_descr descr;
12556bd94f4SPierre Jolivet       sparse_matrix_t     spr;
12656bd94f4SPierre Jolivet       descr.type = SPARSE_MATRIX_TYPE_GENERAL;
12756bd94f4SPierre Jolivet       descr.diag = SPARSE_DIAG_NON_UNIT;
12856bd94f4SPierre Jolivet #endif
12956bd94f4SPierre Jolivet       if (convert) {
1309566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&A));
13156bd94f4SPierre Jolivet       }
1329566063dSJacob Faibussowitsch       PetscCall(PetscStrstr(type[i], "mkl", &tmp));
13356bd94f4SPierre Jolivet       if (tmp) {
13456bd94f4SPierre Jolivet         size_t mlen, tlen;
13556bd94f4SPierre Jolivet         char base[256];
13656bd94f4SPierre Jolivet 
13756bd94f4SPierre Jolivet         mkl  = PETSC_TRUE;
1389566063dSJacob Faibussowitsch         PetscCall(PetscStrlen(tmp, &mlen));
1399566063dSJacob Faibussowitsch         PetscCall(PetscStrlen(type[i], &tlen));
1409566063dSJacob Faibussowitsch         PetscCall(PetscStrncpy(base, type[i], tlen-mlen + 1));
1419566063dSJacob Faibussowitsch         PetscCall(MatConvert(E, base, convert ? MAT_INITIAL_MATRIX : MAT_INPLACE_MATRIX, &A));
14256bd94f4SPierre Jolivet       } else {
14356bd94f4SPierre Jolivet         mkl  = PETSC_FALSE;
1449566063dSJacob Faibussowitsch         PetscCall(PetscStrstr(type[i], "maij", &tmp));
14556bd94f4SPierre Jolivet         if (!tmp) {
1469566063dSJacob Faibussowitsch           PetscCall(MatConvert(E, type[i], convert ? MAT_INITIAL_MATRIX : MAT_INPLACE_MATRIX, &A));
14756bd94f4SPierre Jolivet         } else {
1489566063dSJacob Faibussowitsch           PetscCall(MatConvert(E, MATAIJ, convert ? MAT_INITIAL_MATRIX : MAT_INPLACE_MATRIX, &A));
14956bd94f4SPierre Jolivet           maij = PETSC_TRUE;
15056bd94f4SPierre Jolivet         }
15156bd94f4SPierre Jolivet       }
1529566063dSJacob Faibussowitsch       PetscCall(PetscObjectTypeCompareAny((PetscObject)A, &cuda, MATSEQAIJCUSPARSE, MATMPIAIJCUSPARSE, ""));
15356bd94f4SPierre Jolivet       if (mkl) {
15456bd94f4SPierre Jolivet         const PetscInt *Ai, *Aj;
15556bd94f4SPierre Jolivet         PetscInt       An;
15656bd94f4SPierre Jolivet         PetscBool      done;
15756bd94f4SPierre Jolivet 
1589566063dSJacob Faibussowitsch         PetscCall(PetscObjectTypeCompareAny((PetscObject)A, &flg, MATSEQAIJ, MATSEQBAIJ, MATSEQSBAIJ, ""));
15928b400f6SJacob Faibussowitsch         PetscCheck(flg,PETSC_COMM_WORLD, PETSC_ERR_USER_INPUT, "Not implemented");
1609566063dSJacob Faibussowitsch         PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQAIJ, &flg));
1619566063dSJacob Faibussowitsch         PetscCall(MatGetRowIJ(A, 0, PETSC_FALSE, flg ? PETSC_FALSE : PETSC_TRUE, &An, &Ai, &Aj, &done));
16228b400f6SJacob Faibussowitsch         PetscCheck(done,PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Inconsistent sizes");
1639566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(An + 1, &ia_ptr));
1649566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(Ai[An], &ja_ptr));
16556bd94f4SPierre Jolivet         if (flg) { /* SeqAIJ */
16656bd94f4SPierre Jolivet           for (k = 0; k < An + 1; ++k) ia_ptr[k] = Ai[k];
16756bd94f4SPierre Jolivet           for (k = 0; k < Ai[An]; ++k) ja_ptr[k] = Aj[k];
1689566063dSJacob Faibussowitsch           PetscCall(MatSeqAIJGetArray(A, &a_ptr));
169792fecdfSBarry Smith           PetscCallMKLSparse(mkl_sparse_d_create_csr, (&spr, SPARSE_INDEX_BASE_ZERO, An, An, ia_ptr, ia_ptr + 1, ja_ptr, a_ptr));
17056bd94f4SPierre Jolivet         } else {
1719566063dSJacob Faibussowitsch           PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQBAIJ, &flg));
17256bd94f4SPierre Jolivet           if (flg) {
17356bd94f4SPierre Jolivet             for (k = 0; k < An + 1; ++k) ia_ptr[k] = Ai[k] + 1; /* Fortran indexing to maximize cases covered by _mm routines */
17456bd94f4SPierre Jolivet             for (k = 0; k < Ai[An]; ++k) ja_ptr[k] = Aj[k] + 1; /* Fortran indexing to maximize cases covered by _mm routines */
1759566063dSJacob Faibussowitsch             PetscCall(MatSeqBAIJGetArray(A, &a_ptr));
176792fecdfSBarry 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));
17756bd94f4SPierre Jolivet           } else {
1789566063dSJacob Faibussowitsch             PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQSBAIJ, &flg));
17956bd94f4SPierre Jolivet             if (flg) {
18056bd94f4SPierre Jolivet               for (k = 0; k < An + 1; ++k) ia_ptr[k] = Ai[k] + 1; /* Fortran indexing to maximize cases covered by _mm routines */
18156bd94f4SPierre Jolivet               for (k = 0; k < Ai[An]; ++k) ja_ptr[k] = Aj[k] + 1; /* Fortran indexing to maximize cases covered by _mm routines */
1829566063dSJacob Faibussowitsch               PetscCall(MatSeqSBAIJGetArray(A, &a_ptr));
183792fecdfSBarry 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));
184b88df2e7SBarry Smith #if defined(PETSC_HAVE_MKL_SPARSE_OPTIMIZE)
18556bd94f4SPierre Jolivet               descr.type = SPARSE_MATRIX_TYPE_SYMMETRIC;
18656bd94f4SPierre Jolivet               descr.mode = SPARSE_FILL_MODE_UPPER;
18756bd94f4SPierre Jolivet               descr.diag = SPARSE_DIAG_NON_UNIT;
18856bd94f4SPierre Jolivet #endif
18956bd94f4SPierre Jolivet             }
19056bd94f4SPierre Jolivet           }
19156bd94f4SPierre Jolivet         }
1929566063dSJacob Faibussowitsch         PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQAIJ, &flg));
1939566063dSJacob Faibussowitsch         PetscCall(MatRestoreRowIJ(A, 0, PETSC_FALSE, flg ? PETSC_FALSE : PETSC_TRUE, &An, &Ai, &Aj, &done));
19456bd94f4SPierre Jolivet       }
19556bd94f4SPierre Jolivet 
1969566063dSJacob Faibussowitsch       PetscCall(MatViewFromOptions(A, NULL, "-A_view"));
19756bd94f4SPierre Jolivet 
19856bd94f4SPierre Jolivet       for (k = 0; k < nN; ++k) {
19956bd94f4SPierre Jolivet         MatType       Atype, Ctype;
20056bd94f4SPierre Jolivet         PetscInt      AM, AN, CM, CN, t;
201956f8c0dSBarry Smith #if defined(PETSC_USE_LOG)
20256bd94f4SPierre Jolivet         PetscLogStage stage, tstage;
20356bd94f4SPierre Jolivet         char          stage_s[256];
204956f8c0dSBarry Smith #endif
20556bd94f4SPierre Jolivet 
2069566063dSJacob Faibussowitsch         PetscCall(MatCreateDense(PETSC_COMM_WORLD, bs[j] * m, PETSC_DECIDE, bs[j] * M, N[k], NULL, &C));
2079566063dSJacob Faibussowitsch         PetscCall(MatCreateDense(PETSC_COMM_WORLD, bs[j] * m, PETSC_DECIDE, bs[j] * M, N[k], NULL, &D));
2089566063dSJacob Faibussowitsch         PetscCall(MatSetRandom(C, NULL));
20956bd94f4SPierre Jolivet         if (cuda) { /* convert to GPU if needed */
2109566063dSJacob Faibussowitsch           PetscCall(MatConvert(C, MATDENSECUDA, MAT_INPLACE_MATRIX, &C));
2119566063dSJacob Faibussowitsch           PetscCall(MatConvert(D, MATDENSECUDA, MAT_INPLACE_MATRIX, &D));
21256bd94f4SPierre Jolivet         }
21356bd94f4SPierre Jolivet         if (mkl) {
214792fecdfSBarry 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));
215792fecdfSBarry Smith           else          PetscCallMKLSparse(mkl_sparse_set_mv_hint, (spr, SPARSE_OPERATION_NON_TRANSPOSE, descr, 1 + trial));
216792fecdfSBarry Smith           PetscCallMKLSparse(mkl_sparse_set_memory_hint, (spr, SPARSE_MEMORY_AGGRESSIVE));
217792fecdfSBarry Smith           PetscCallMKLSparse(mkl_sparse_optimize, (spr));
21856bd94f4SPierre Jolivet         }
2199566063dSJacob Faibussowitsch         PetscCall(MatGetType(A, &Atype));
2209566063dSJacob Faibussowitsch         PetscCall(MatGetType(C, &Ctype));
2219566063dSJacob Faibussowitsch         PetscCall(MatGetSize(A, &AM, &AN));
2229566063dSJacob Faibussowitsch         PetscCall(MatGetSize(C, &CM, &CN));
22356bd94f4SPierre Jolivet 
224956f8c0dSBarry Smith #if defined(PETSC_USE_LOG)
22556bd94f4SPierre Jolivet         if (!maij || N[k] > 1) {
2269566063dSJacob Faibussowitsch           PetscCall(PetscSNPrintf(stage_s, sizeof(stage_s), "type_%s-bs_%" PetscInt_FMT "-N_%02d", type[i], bs[j], (int)N[k]));
2279566063dSJacob Faibussowitsch           PetscCall(PetscLogStageRegister(stage_s, &stage));
22856bd94f4SPierre Jolivet         }
22956bd94f4SPierre Jolivet         if (trans && N[k] > 1) {
2309566063dSJacob Faibussowitsch           PetscCall(PetscSNPrintf(stage_s, sizeof(stage_s), "trans_type_%s-bs_%" PetscInt_FMT "-N_%02d", type[i], bs[j], (int)N[k]));
2319566063dSJacob Faibussowitsch           PetscCall(PetscLogStageRegister(stage_s, &tstage));
23256bd94f4SPierre Jolivet         }
233956f8c0dSBarry Smith #endif
23456bd94f4SPierre Jolivet         /* A*B */
23556bd94f4SPierre Jolivet         if (N[k] > 1) {
23656bd94f4SPierre Jolivet           if (!maij) {
2379566063dSJacob Faibussowitsch             PetscCall(MatProductCreateWithMat(A, C, NULL, D));
2389566063dSJacob Faibussowitsch             PetscCall(MatProductSetType(D, MATPRODUCT_AB));
2399566063dSJacob Faibussowitsch             PetscCall(MatProductSetFromOptions(D));
2409566063dSJacob Faibussowitsch             PetscCall(MatProductSymbolic(D));
24156bd94f4SPierre Jolivet           }
24256bd94f4SPierre Jolivet 
24356bd94f4SPierre Jolivet           if (!mkl) {
24456bd94f4SPierre Jolivet             if (!maij) {
2459566063dSJacob Faibussowitsch               PetscCall(MatProductNumeric(D));
2469566063dSJacob 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));
2479566063dSJacob Faibussowitsch               PetscCall(PetscLogStagePush(stage));
24856bd94f4SPierre Jolivet               for (t = 0; t < trial; ++t) {
2499566063dSJacob Faibussowitsch                 PetscCall(MatProductNumeric(D));
25056bd94f4SPierre Jolivet               }
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));
26956bd94f4SPierre Jolivet               for (t = 0; t < trial; ++t) {
2709566063dSJacob Faibussowitsch                 PetscCall(MatMult(E, cC, cD));
27156bd94f4SPierre Jolivet               }
2729566063dSJacob Faibussowitsch               PetscCall(PetscLogStagePop());
2739566063dSJacob Faibussowitsch               PetscCall(VecDestroy(&cD));
2749566063dSJacob Faibussowitsch               PetscCall(VecDestroy(&cC));
2759566063dSJacob Faibussowitsch               PetscCall(MatDestroy(&E));
2769566063dSJacob Faibussowitsch               PetscCall(MatDenseRestoreArrayWrite(Dt, &d_ptr));
2779566063dSJacob Faibussowitsch               PetscCall(MatDenseRestoreArrayRead(Ct, &c_ptr));
2789566063dSJacob Faibussowitsch               PetscCall(MatTranspose(Ct, MAT_INPLACE_MATRIX, &Ct));
2799566063dSJacob Faibussowitsch               PetscCall(MatTranspose(Dt, MAT_INPLACE_MATRIX, &Dt));
28056bd94f4SPierre Jolivet             }
28156bd94f4SPierre Jolivet           } else {
28256bd94f4SPierre Jolivet             const PetscScalar *c_ptr;
28356bd94f4SPierre Jolivet             PetscScalar       *d_ptr;
28456bd94f4SPierre Jolivet 
2859566063dSJacob Faibussowitsch             PetscCall(MatDenseGetArrayRead(C, &c_ptr));
2869566063dSJacob Faibussowitsch             PetscCall(MatDenseGetArrayWrite(D, &d_ptr));
287792fecdfSBarry 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));
2889566063dSJacob 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));
2899566063dSJacob Faibussowitsch             PetscCall(PetscLogStagePush(stage));
29056bd94f4SPierre Jolivet             for (t = 0; t < trial; ++t) {
291792fecdfSBarry 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));
29256bd94f4SPierre Jolivet             }
2939566063dSJacob Faibussowitsch             PetscCall(PetscLogStagePop());
2949566063dSJacob Faibussowitsch             PetscCall(MatDenseRestoreArrayWrite(D, &d_ptr));
2959566063dSJacob Faibussowitsch             PetscCall(MatDenseRestoreArrayRead(C, &c_ptr));
29656bd94f4SPierre Jolivet           }
29756bd94f4SPierre Jolivet         } else if (maij) {
2989566063dSJacob Faibussowitsch           PetscCall(MatDestroy(&C));
2999566063dSJacob Faibussowitsch           PetscCall(MatDestroy(&D));
30056bd94f4SPierre Jolivet           continue;
30156bd94f4SPierre Jolivet         } else if (!mkl) {
30256bd94f4SPierre Jolivet           Vec cC, cD;
30356bd94f4SPierre Jolivet 
3049566063dSJacob Faibussowitsch           PetscCall(MatDenseGetColumnVecRead(C, 0, &cC));
3059566063dSJacob Faibussowitsch           PetscCall(MatDenseGetColumnVecWrite(D, 0, &cD));
3069566063dSJacob Faibussowitsch           PetscCall(MatMult(A, cC, cD));
3079566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Benchmarking MatMult: with A %s %" PetscInt_FMT "x%" PetscInt_FMT "\n", Atype, AM, AN));
3089566063dSJacob Faibussowitsch           PetscCall(PetscLogStagePush(stage));
30956bd94f4SPierre Jolivet           for (t = 0; t < trial; ++t) {
3109566063dSJacob Faibussowitsch             PetscCall(MatMult(A, cC, cD));
31156bd94f4SPierre Jolivet           }
3129566063dSJacob Faibussowitsch           PetscCall(PetscLogStagePop());
3139566063dSJacob Faibussowitsch           PetscCall(MatDenseRestoreColumnVecRead(C, 0, &cC));
3149566063dSJacob Faibussowitsch           PetscCall(MatDenseRestoreColumnVecWrite(D, 0, &cD));
31556bd94f4SPierre Jolivet         } else {
31656bd94f4SPierre Jolivet           const PetscScalar *c_ptr;
31756bd94f4SPierre Jolivet           PetscScalar       *d_ptr;
31856bd94f4SPierre Jolivet 
3199566063dSJacob Faibussowitsch           PetscCall(MatDenseGetArrayRead(C, &c_ptr));
3209566063dSJacob Faibussowitsch           PetscCall(MatDenseGetArrayWrite(D, &d_ptr));
3219566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Benchmarking mkl_sparse_d_mv: with A %s %" PetscInt_FMT "x%" PetscInt_FMT "\n", Atype, AM, AN));
322792fecdfSBarry Smith           PetscCallMKLSparse(mkl_sparse_d_mv, (SPARSE_OPERATION_NON_TRANSPOSE, 1.0, spr, descr, c_ptr, 0.0, d_ptr));
3239566063dSJacob Faibussowitsch           PetscCall(PetscLogStagePush(stage));
32456bd94f4SPierre Jolivet           for (t = 0; t < trial; ++t) {
325792fecdfSBarry Smith             PetscCallMKLSparse(mkl_sparse_d_mv, (SPARSE_OPERATION_NON_TRANSPOSE, 1.0, spr, descr, c_ptr, 0.0, d_ptr));
32656bd94f4SPierre Jolivet           }
3279566063dSJacob Faibussowitsch           PetscCall(PetscLogStagePop());
3289566063dSJacob Faibussowitsch           PetscCall(MatDenseRestoreArrayWrite(D, &d_ptr));
3299566063dSJacob Faibussowitsch           PetscCall(MatDenseRestoreArrayRead(C, &c_ptr));
33056bd94f4SPierre Jolivet         }
33156bd94f4SPierre Jolivet 
33256bd94f4SPierre Jolivet         if (check) {
3339566063dSJacob Faibussowitsch           PetscCall(MatMatMultEqual(A, C, D, 10, &flg));
33456bd94f4SPierre Jolivet           if (!flg) {
33556bd94f4SPierre Jolivet             MatType Dtype;
33656bd94f4SPierre Jolivet 
3379566063dSJacob Faibussowitsch             PetscCall(MatGetType(D, &Dtype));
3389566063dSJacob 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]));
33956bd94f4SPierre Jolivet           }
34056bd94f4SPierre Jolivet         }
34156bd94f4SPierre Jolivet 
34256bd94f4SPierre Jolivet         /* MKL implementation seems buggy for ABt */
34356bd94f4SPierre Jolivet         /* A*Bt */
34456bd94f4SPierre Jolivet         if (!mkl && trans && N[k] > 1) {
3459566063dSJacob Faibussowitsch           PetscCall(PetscObjectTypeCompareAny((PetscObject)A, &flg, MATSEQAIJ, MATMPIAIJ, ""));
34656bd94f4SPierre Jolivet           if (flg) {
3479566063dSJacob Faibussowitsch             PetscCall(MatTranspose(C, MAT_INPLACE_MATRIX, &C));
3489566063dSJacob Faibussowitsch             PetscCall(MatGetType(C, &Ctype));
34956bd94f4SPierre Jolivet             if (!mkl) {
3509566063dSJacob Faibussowitsch               PetscCall(MatProductCreateWithMat(A, C, NULL, D));
3519566063dSJacob Faibussowitsch               PetscCall(MatProductSetType(D, MATPRODUCT_ABt));
3529566063dSJacob Faibussowitsch               PetscCall(MatProductSetFromOptions(D));
3539566063dSJacob Faibussowitsch               PetscCall(MatProductSymbolic(D));
3549566063dSJacob Faibussowitsch               PetscCall(MatProductNumeric(D));
3559566063dSJacob 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));
3569566063dSJacob Faibussowitsch               PetscCall(PetscLogStagePush(tstage));
35756bd94f4SPierre Jolivet               for (t = 0; t < trial; ++t) {
3589566063dSJacob Faibussowitsch                 PetscCall(MatProductNumeric(D));
35956bd94f4SPierre Jolivet               }
3609566063dSJacob Faibussowitsch               PetscCall(PetscLogStagePop());
36156bd94f4SPierre Jolivet             } else {
36256bd94f4SPierre Jolivet               const PetscScalar *c_ptr;
36356bd94f4SPierre Jolivet               PetscScalar       *d_ptr;
36456bd94f4SPierre Jolivet 
365792fecdfSBarry Smith               PetscCallMKLSparse(mkl_sparse_set_mm_hint, (spr, SPARSE_OPERATION_NON_TRANSPOSE, descr, SPARSE_LAYOUT_ROW_MAJOR, N[k], 1 + trial));
366792fecdfSBarry Smith               PetscCallMKLSparse(mkl_sparse_optimize, (spr));
3679566063dSJacob Faibussowitsch               PetscCall(MatDenseGetArrayRead(C, &c_ptr));
3689566063dSJacob Faibussowitsch               PetscCall(MatDenseGetArrayWrite(D, &d_ptr));
3699566063dSJacob 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));
370792fecdfSBarry 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));
3719566063dSJacob Faibussowitsch               PetscCall(PetscLogStagePush(stage));
37256bd94f4SPierre Jolivet               for (t = 0; t < trial; ++t) {
373792fecdfSBarry 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));
37456bd94f4SPierre Jolivet               }
3759566063dSJacob Faibussowitsch               PetscCall(PetscLogStagePop());
3769566063dSJacob Faibussowitsch               PetscCall(MatDenseRestoreArrayWrite(D, &d_ptr));
3779566063dSJacob Faibussowitsch               PetscCall(MatDenseRestoreArrayRead(C, &c_ptr));
37856bd94f4SPierre Jolivet             }
37956bd94f4SPierre Jolivet           }
38056bd94f4SPierre Jolivet         }
38156bd94f4SPierre Jolivet 
38256bd94f4SPierre Jolivet         if (!mkl && trans && N[k] > 1 && flg && check) {
3839566063dSJacob Faibussowitsch           PetscCall(MatMatTransposeMultEqual(A, C, D, 10, &flg));
38456bd94f4SPierre Jolivet           if (!flg) {
38556bd94f4SPierre Jolivet             MatType Dtype;
3869566063dSJacob Faibussowitsch             PetscCall(MatGetType(D, &Dtype));
3879566063dSJacob 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]));
38856bd94f4SPierre Jolivet           }
38956bd94f4SPierre Jolivet         }
3909566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&C));
3919566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&D));
39256bd94f4SPierre Jolivet       }
39356bd94f4SPierre Jolivet       if (mkl) {
394792fecdfSBarry Smith         PetscCallMKLSparse(mkl_sparse_destroy, (spr));
3959566063dSJacob Faibussowitsch         PetscCall(PetscFree(ia_ptr));
3969566063dSJacob Faibussowitsch         PetscCall(PetscFree(ja_ptr));
39756bd94f4SPierre Jolivet       }
39856bd94f4SPierre Jolivet       if (cuda && i != ntype - 1) {
3999566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "AIJCUSPARSE must be last, otherwise MatConvert() to another MatType is too slow\n"));
40056bd94f4SPierre Jolivet         break;
40156bd94f4SPierre Jolivet       }
40256bd94f4SPierre Jolivet     }
4039566063dSJacob Faibussowitsch     if (E != A) PetscCall(MatDestroy(&E));
4049566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A));
40556bd94f4SPierre Jolivet   }
4069566063dSJacob Faibussowitsch   for (m = 0; m < ntype; ++m) PetscCall(PetscFree(type[m]));
4079566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
40856bd94f4SPierre Jolivet   return 0;
40956bd94f4SPierre Jolivet }
41056bd94f4SPierre Jolivet 
41156bd94f4SPierre Jolivet /*TEST
412b8582c7cSSatish Balay    build:
413dfd57a17SPierre Jolivet      requires: double !complex !defined(PETSC_USE_64BIT_INDICES)
41456bd94f4SPierre Jolivet 
41556bd94f4SPierre Jolivet    testset:
41656bd94f4SPierre Jolivet      nsize: 1
41756bd94f4SPierre Jolivet      filter: sed "/Benchmarking/d"
41856bd94f4SPierre 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}
41956bd94f4SPierre Jolivet      test:
42056bd94f4SPierre Jolivet        suffix: basic
42156bd94f4SPierre Jolivet        args: -type aij,sbaij,baij
42256bd94f4SPierre Jolivet        output_file: output/ex237.out
42356bd94f4SPierre Jolivet      test:
42456bd94f4SPierre Jolivet        suffix: maij
42556bd94f4SPierre Jolivet        args: -type aij,maij
42656bd94f4SPierre Jolivet        output_file: output/ex237.out
42756bd94f4SPierre Jolivet      test:
42856bd94f4SPierre Jolivet        suffix: cuda
42956bd94f4SPierre Jolivet        requires: cuda
43056bd94f4SPierre Jolivet        args: -type aij,aijcusparse
43156bd94f4SPierre Jolivet        output_file: output/ex237.out
43256bd94f4SPierre Jolivet      test:
43356bd94f4SPierre Jolivet        suffix: mkl
434b88df2e7SBarry Smith        requires: mkl_sparse_optimize
43556bd94f4SPierre Jolivet        args: -type aij,aijmkl,baijmkl,sbaijmkl
43656bd94f4SPierre Jolivet        output_file: output/ex237.out
43756bd94f4SPierre Jolivet 
43856bd94f4SPierre Jolivet TEST*/
439