xref: /petsc/src/mat/impls/aij/seq/aijmkl/aijmkl.c (revision f8990b4a39a56ca2debaee2c8be1da38fe13b9f8)
14a2a386eSRichard Tran Mills /*
24a2a386eSRichard Tran Mills   Defines basic operations for the MATSEQAIJMKL matrix class.
34a2a386eSRichard Tran Mills   This class is derived from the MATSEQAIJ class and retains the
44a2a386eSRichard Tran Mills   compressed row storage (aka Yale sparse matrix format) but uses
54a2a386eSRichard Tran Mills   sparse BLAS operations from the Intel Math Kernel Library (MKL)
64a2a386eSRichard Tran Mills   wherever possible.
74a2a386eSRichard Tran Mills */
84a2a386eSRichard Tran Mills 
94a2a386eSRichard Tran Mills #include <../src/mat/impls/aij/seq/aij.h>
104a2a386eSRichard Tran Mills #include <../src/mat/impls/aij/seq/aijmkl/aijmkl.h>
114a2a386eSRichard Tran Mills 
124a2a386eSRichard Tran Mills /* MKL include files. */
134a2a386eSRichard Tran Mills #include <mkl_spblas.h>  /* Sparse BLAS */
144a2a386eSRichard Tran Mills 
154a2a386eSRichard Tran Mills typedef struct {
16c9d46305SRichard Tran Mills   PetscBool no_SpMV2;  /* If PETSC_TRUE, then don't use the MKL SpMV2 inspector-executor routines. */
175b49642aSRichard Tran Mills   PetscBool eager_inspection; /* If PETSC_TRUE, then call mkl_sparse_optimize() in MatDuplicate()/MatAssemblyEnd(). */
184abfa3b3SRichard Tran Mills   PetscBool sparse_optimized; /* If PETSC_TRUE, then mkl_sparse_optimize() has been called. */
19551aa5c8SRichard Tran Mills   PetscObjectState state;
20b8cbc1fbSRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_OPTIMIZE
21df555b71SRichard Tran Mills   sparse_matrix_t csrA; /* "Handle" used by SpMV2 inspector-executor routines. */
22df555b71SRichard Tran Mills   struct matrix_descr descr;
23b8cbc1fbSRichard Tran Mills #endif
244a2a386eSRichard Tran Mills } Mat_SeqAIJMKL;
254a2a386eSRichard Tran Mills 
264a2a386eSRichard Tran Mills extern PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat,MatAssemblyType);
274a2a386eSRichard Tran Mills 
284a2a386eSRichard Tran Mills PETSC_INTERN PetscErrorCode MatConvert_SeqAIJMKL_SeqAIJ(Mat A,MatType type,MatReuse reuse,Mat *newmat)
294a2a386eSRichard Tran Mills {
304a2a386eSRichard Tran Mills   /* This routine is only called to convert a MATAIJMKL to its base PETSc type, */
314a2a386eSRichard Tran Mills   /* so we will ignore 'MatType type'. */
324a2a386eSRichard Tran Mills   PetscErrorCode ierr;
334a2a386eSRichard Tran Mills   Mat            B       = *newmat;
34c1d5218aSRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_OPTIMIZE
354a2a386eSRichard Tran Mills   Mat_SeqAIJMKL  *aijmkl=(Mat_SeqAIJMKL*)A->spptr;
36c1d5218aSRichard Tran Mills #endif
374a2a386eSRichard Tran Mills 
384a2a386eSRichard Tran Mills   PetscFunctionBegin;
394a2a386eSRichard Tran Mills   if (reuse == MAT_INITIAL_MATRIX) {
404a2a386eSRichard Tran Mills     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
414a2a386eSRichard Tran Mills   }
424a2a386eSRichard Tran Mills 
434a2a386eSRichard Tran Mills   /* Reset the original function pointers. */
4454871a98SRichard Tran Mills   B->ops->duplicate        = MatDuplicate_SeqAIJ;
454a2a386eSRichard Tran Mills   B->ops->assemblyend      = MatAssemblyEnd_SeqAIJ;
464a2a386eSRichard Tran Mills   B->ops->destroy          = MatDestroy_SeqAIJ;
4754871a98SRichard Tran Mills   B->ops->mult             = MatMult_SeqAIJ;
48ff03dc53SRichard Tran Mills   B->ops->multtranspose    = MatMultTranspose_SeqAIJ;
4954871a98SRichard Tran Mills   B->ops->multadd          = MatMultAdd_SeqAIJ;
50ff03dc53SRichard Tran Mills   B->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJ;
5145fbe478SRichard Tran Mills   B->ops->matmult          = MatMatMult_SeqAIJ_SeqAIJ;
52e8be1fc7SRichard Tran Mills   B->ops->matmultnumeric   = MatMatMultNumeric_SeqAIJ_SeqAIJ;
534f53af40SRichard Tran Mills   B->ops->ptap             = MatPtAP_SeqAIJ_SeqAIJ;
544f53af40SRichard Tran Mills   B->ops->ptapnumeric      = MatPtAPNumeric_SeqAIJ_SeqAIJ;
55372ec6bbSRichard Tran Mills   B->ops->transposematmult = MatTransposeMatMult_SeqAIJ_SeqAIJ;
564a2a386eSRichard Tran Mills 
57e9c94282SRichard Tran Mills   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaijmkl_seqaij_C",NULL);CHKERRQ(ierr);
58e9c94282SRichard Tran Mills   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMult_seqdense_seqaijmkl_C",NULL);CHKERRQ(ierr);
59e9c94282SRichard Tran Mills   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMultSymbolic_seqdense_seqaijmkl_C",NULL);CHKERRQ(ierr);
60e9c94282SRichard Tran Mills   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMultNumeric_seqdense_seqaijmkl_C",NULL);CHKERRQ(ierr);
6145fbe478SRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_OPTIMIZE
624a940b00SSatish Balay   if(!aijmkl->no_SpMV2) {
6345fbe478SRichard Tran Mills     ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMult_seqaijmkl_seqaijmkl_C",NULL);CHKERRQ(ierr);
64e8be1fc7SRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_SP2M
65e8be1fc7SRichard Tran Mills     ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMultNumeric_seqaijmkl_seqaijmkl_C",NULL);CHKERRQ(ierr);
66e8be1fc7SRichard Tran Mills #endif
67372ec6bbSRichard Tran Mills     ierr = PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMult_seqaijmkl_seqaijmkl_C",NULL);CHKERRQ(ierr);
6845fbe478SRichard Tran Mills   }
69e9c94282SRichard Tran Mills 
704abfa3b3SRichard Tran Mills   /* Free everything in the Mat_SeqAIJMKL data structure. Currently, this
71e9c94282SRichard Tran Mills    * simply involves destroying the MKL sparse matrix handle and then freeing
72e9c94282SRichard Tran Mills    * the spptr pointer. */
73a8327b06SKarl Rupp   if (reuse == MAT_INITIAL_MATRIX) aijmkl = (Mat_SeqAIJMKL*)B->spptr;
74a8327b06SKarl Rupp 
754abfa3b3SRichard Tran Mills   if (aijmkl->sparse_optimized) {
760632b357SRichard Tran Mills     sparse_status_t stat;
774abfa3b3SRichard Tran Mills     stat = mkl_sparse_destroy(aijmkl->csrA);
789c46acdfSRichard Tran Mills     if (stat != SPARSE_STATUS_SUCCESS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Intel MKL error: unable to set hints/complete mkl_sparse_optimize");
794abfa3b3SRichard Tran Mills   }
804abfa3b3SRichard Tran Mills #endif /* PETSC_HAVE_MKL_SPARSE_OPTIMIZE */
81e9c94282SRichard Tran Mills   ierr = PetscFree(B->spptr);CHKERRQ(ierr);
824a2a386eSRichard Tran Mills 
834a2a386eSRichard Tran Mills   /* Change the type of B to MATSEQAIJ. */
844a2a386eSRichard Tran Mills   ierr = PetscObjectChangeTypeName((PetscObject)B, MATSEQAIJ);CHKERRQ(ierr);
854a2a386eSRichard Tran Mills 
864a2a386eSRichard Tran Mills   *newmat = B;
874a2a386eSRichard Tran Mills   PetscFunctionReturn(0);
884a2a386eSRichard Tran Mills }
894a2a386eSRichard Tran Mills 
904a2a386eSRichard Tran Mills PetscErrorCode MatDestroy_SeqAIJMKL(Mat A)
914a2a386eSRichard Tran Mills {
924a2a386eSRichard Tran Mills   PetscErrorCode ierr;
934a2a386eSRichard Tran Mills   Mat_SeqAIJMKL  *aijmkl = (Mat_SeqAIJMKL*) A->spptr;
944a2a386eSRichard Tran Mills 
954a2a386eSRichard Tran Mills   PetscFunctionBegin;
96e9c94282SRichard Tran Mills 
97e9c94282SRichard Tran Mills   /* If MatHeaderMerge() was used, then this SeqAIJMKL matrix will not have an
98e9c94282SRichard Tran Mills    * spptr pointer. */
99e9c94282SRichard Tran Mills   if (aijmkl) {
1004a2a386eSRichard Tran Mills     /* Clean up everything in the Mat_SeqAIJMKL data structure, then free A->spptr. */
1014abfa3b3SRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_OPTIMIZE
1024abfa3b3SRichard Tran Mills     if (aijmkl->sparse_optimized) {
1034abfa3b3SRichard Tran Mills       sparse_status_t stat = SPARSE_STATUS_SUCCESS;
1044abfa3b3SRichard Tran Mills       stat = mkl_sparse_destroy(aijmkl->csrA);
1059c46acdfSRichard Tran Mills       if (stat != SPARSE_STATUS_SUCCESS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Intel MKL error: error in mkl_sparse_destroy");
1064abfa3b3SRichard Tran Mills     }
1074abfa3b3SRichard Tran Mills #endif /* PETSC_HAVE_MKL_SPARSE_OPTIMIZE */
1084a2a386eSRichard Tran Mills     ierr = PetscFree(A->spptr);CHKERRQ(ierr);
109e9c94282SRichard Tran Mills   }
1104a2a386eSRichard Tran Mills 
1114a2a386eSRichard Tran Mills   /* Change the type of A back to SEQAIJ and use MatDestroy_SeqAIJ()
1124a2a386eSRichard Tran Mills    * to destroy everything that remains. */
1134a2a386eSRichard Tran Mills   ierr = PetscObjectChangeTypeName((PetscObject)A, MATSEQAIJ);CHKERRQ(ierr);
1144a2a386eSRichard Tran Mills   /* Note that I don't call MatSetType().  I believe this is because that
1154a2a386eSRichard Tran Mills    * is only to be called when *building* a matrix.  I could be wrong, but
1164a2a386eSRichard Tran Mills    * that is how things work for the SuperLU matrix class. */
1174a2a386eSRichard Tran Mills   ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
1184a2a386eSRichard Tran Mills   PetscFunctionReturn(0);
1194a2a386eSRichard Tran Mills }
1204a2a386eSRichard Tran Mills 
1215b49642aSRichard Tran Mills /* MatSeqAIJKL_create_mkl_handle(), if called with an AIJMKL matrix that has not had mkl_sparse_optimize() called for it,
1225b49642aSRichard Tran Mills  * creates an MKL sparse matrix handle from the AIJ arrays and calls mkl_sparse_optimize().
1235b49642aSRichard Tran Mills  * If called with an AIJMKL matrix for which aijmkl->sparse_optimized == PETSC_TRUE, then it destroys the old matrix
1245b49642aSRichard Tran Mills  * handle, creates a new one, and then calls mkl_sparse_optimize().
1255b49642aSRichard Tran Mills  * Although in normal MKL usage it is possible to have a valid matrix handle on which mkl_sparse_optimize() has not been
1265b49642aSRichard Tran Mills  * called, for AIJMKL the handle creation and optimization step always occur together, so we don't handle the case of
1275b49642aSRichard Tran Mills  * an unoptimized matrix handle here. */
1286e369cd5SRichard Tran Mills PETSC_INTERN PetscErrorCode MatSeqAIJMKL_create_mkl_handle(Mat A)
1294a2a386eSRichard Tran Mills {
1306e369cd5SRichard Tran Mills #ifndef PETSC_HAVE_MKL_SPARSE_OPTIMIZE
1316e369cd5SRichard Tran Mills   /* If the MKL library does not have mkl_sparse_optimize(), then this routine
1326e369cd5SRichard Tran Mills    * does nothing. We make it callable anyway in this case because it cuts
1336e369cd5SRichard Tran Mills    * down on littering the code with #ifdefs. */
13445fbe478SRichard Tran Mills   PetscFunctionBegin;
1356e369cd5SRichard Tran Mills   PetscFunctionReturn(0);
1366e369cd5SRichard Tran Mills #else
137a8327b06SKarl Rupp   Mat_SeqAIJ       *a = (Mat_SeqAIJ*)A->data;
138a8327b06SKarl Rupp   Mat_SeqAIJMKL    *aijmkl = (Mat_SeqAIJMKL*)A->spptr;
139a8327b06SKarl Rupp   PetscInt         m,n;
140a8327b06SKarl Rupp   MatScalar        *aa;
141a8327b06SKarl Rupp   PetscInt         *aj,*ai;
1426e369cd5SRichard Tran Mills   sparse_status_t  stat;
143551aa5c8SRichard Tran Mills   PetscErrorCode   ierr;
1444a2a386eSRichard Tran Mills 
145a8327b06SKarl Rupp   PetscFunctionBegin;
1466e369cd5SRichard Tran Mills   if (aijmkl->no_SpMV2) PetscFunctionReturn(0);
1476e369cd5SRichard Tran Mills 
1480632b357SRichard Tran Mills   if (aijmkl->sparse_optimized) {
1490632b357SRichard Tran Mills     /* Matrix has been previously assembled and optimized. Must destroy old
1500632b357SRichard Tran Mills      * matrix handle before running the optimization step again. */
1510632b357SRichard Tran Mills     stat = mkl_sparse_destroy(aijmkl->csrA);
1529c46acdfSRichard Tran Mills     if (stat != SPARSE_STATUS_SUCCESS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Intel MKL error: error in mkl_sparse_destroy");
1530632b357SRichard Tran Mills   }
1548d3fe1b0SRichard Tran Mills   aijmkl->sparse_optimized = PETSC_FALSE;
1556e369cd5SRichard Tran Mills 
156c9d46305SRichard Tran Mills   /* Now perform the SpMV2 setup and matrix optimization. */
157df555b71SRichard Tran Mills   aijmkl->descr.type        = SPARSE_MATRIX_TYPE_GENERAL;
158df555b71SRichard Tran Mills   aijmkl->descr.mode        = SPARSE_FILL_MODE_LOWER;
159df555b71SRichard Tran Mills   aijmkl->descr.diag        = SPARSE_DIAG_NON_UNIT;
16058678438SRichard Tran Mills   m = A->rmap->n;
16158678438SRichard Tran Mills   n = A->cmap->n;
162df555b71SRichard Tran Mills   aj   = a->j;  /* aj[k] gives column index for element aa[k]. */
163df555b71SRichard Tran Mills   aa   = a->a;  /* Nonzero elements stored row-by-row. */
164df555b71SRichard Tran Mills   ai   = a->i;  /* ai[k] is the position in aa and aj where row k starts. */
16580095d54SIrina Sokolova   if ((a->nz!=0) & !(A->structure_only)) {
1668d3fe1b0SRichard Tran Mills     /* Create a new, optimized sparse matrix handle only if the matrix has nonzero entries.
1678d3fe1b0SRichard Tran Mills      * The MKL sparse-inspector executor routines don't like being passed an empty matrix. */
16858678438SRichard Tran Mills     stat = mkl_sparse_x_create_csr(&aijmkl->csrA,SPARSE_INDEX_BASE_ZERO,m,n,ai,ai+1,aj,aa);
169e8be1fc7SRichard Tran Mills     if (stat != SPARSE_STATUS_SUCCESS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Intel MKL error: unable to create matrix handle");
170df555b71SRichard Tran Mills     stat = mkl_sparse_set_mv_hint(aijmkl->csrA,SPARSE_OPERATION_NON_TRANSPOSE,aijmkl->descr,1000);
171e8be1fc7SRichard Tran Mills     if (stat != SPARSE_STATUS_SUCCESS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Intel MKL error: unable to set mv_hint");
172df555b71SRichard Tran Mills     stat = mkl_sparse_set_memory_hint(aijmkl->csrA,SPARSE_MEMORY_AGGRESSIVE);
173e8be1fc7SRichard Tran Mills     if (stat != SPARSE_STATUS_SUCCESS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Intel MKL error: unable to set memory_hint");
174df555b71SRichard Tran Mills     stat = mkl_sparse_optimize(aijmkl->csrA);
175e8be1fc7SRichard Tran Mills     if (stat != SPARSE_STATUS_SUCCESS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Intel MKL error: unable to complete mkl_sparse_optimize");
1764abfa3b3SRichard Tran Mills     aijmkl->sparse_optimized = PETSC_TRUE;
177e995cf24SRichard Tran Mills     ierr = PetscObjectStateGet((PetscObject)A,&(aijmkl->state));CHKERRQ(ierr);
178c9d46305SRichard Tran Mills   }
1796e369cd5SRichard Tran Mills 
1806e369cd5SRichard Tran Mills   PetscFunctionReturn(0);
181d995685eSRichard Tran Mills #endif
1826e369cd5SRichard Tran Mills }
1836e369cd5SRichard Tran Mills 
18419afcda9SRichard Tran Mills /* MatSeqAIJMKL_create_from_mkl_handle() creates a sequential AIJMKL matrix from an MKL sparse matrix handle.
18519afcda9SRichard Tran Mills  * We need this to implement MatMatMult() using the MKL inspector-executor routines, which return an (unoptimized)
1866c87cf42SRichard Tran Mills  * matrix handle.
187aab60f1bSRichard Tran Mills  * Note: This routine simply destroys and replaces the original matrix if MAT_REUSE_MATRIX has been specified, as
188aab60f1bSRichard Tran Mills  * there is no good alternative. */
18919afcda9SRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_OPTIMIZE
1906c87cf42SRichard Tran Mills PETSC_INTERN PetscErrorCode MatSeqAIJMKL_create_from_mkl_handle(MPI_Comm comm,sparse_matrix_t csrA,MatReuse reuse,Mat *mat)
19119afcda9SRichard Tran Mills {
19219afcda9SRichard Tran Mills   PetscErrorCode      ierr;
19319afcda9SRichard Tran Mills   sparse_status_t     stat;
19419afcda9SRichard Tran Mills   sparse_index_base_t indexing;
19519afcda9SRichard Tran Mills   PetscInt            nrows, ncols;
19645fbe478SRichard Tran Mills   PetscInt            *aj,*ai,*dummy;
19719afcda9SRichard Tran Mills   MatScalar           *aa;
19819afcda9SRichard Tran Mills   Mat                 A;
19919afcda9SRichard Tran Mills   Mat_SeqAIJMKL       *aijmkl;
20019afcda9SRichard Tran Mills 
20145fbe478SRichard Tran Mills   /* Note: Must pass in &dummy below since MKL can't accept NULL for this output array we don't actually want. */
20245fbe478SRichard Tran Mills   stat = mkl_sparse_x_export_csr(csrA,&indexing,&nrows,&ncols,&ai,&dummy,&aj,&aa);
2039c46acdfSRichard Tran Mills   if (stat != SPARSE_STATUS_SUCCESS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Intel MKL error: unable to complete mkl_sparse_x_export_csr()");
2046c87cf42SRichard Tran Mills 
205aab60f1bSRichard Tran Mills   if (reuse == MAT_REUSE_MATRIX) {
206aab60f1bSRichard Tran Mills     ierr = MatDestroy(mat);CHKERRQ(ierr);
207aab60f1bSRichard Tran Mills   }
20819afcda9SRichard Tran Mills   ierr = MatCreate(comm,&A);CHKERRQ(ierr);
20919afcda9SRichard Tran Mills   ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
21045fbe478SRichard Tran Mills   ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,nrows,ncols);CHKERRQ(ierr);
211aab60f1bSRichard Tran Mills   /* We use MatSeqAIJSetPreallocationCSR() instead of MatCreateSeqAIJWithArrays() because we must copy the arrays exported
212aab60f1bSRichard Tran Mills    * from MKL; MKL developers tell us that modifying the arrays may cause unexpected results when using the MKL handle, and
213aab60f1bSRichard Tran Mills    * they will be destroyed when the MKL handle is destroyed.
214aab60f1bSRichard Tran Mills    * (In the interest of reducing memory consumption in future, can we figure out good ways to deal with this?) */
21519afcda9SRichard Tran Mills   ierr = MatSeqAIJSetPreallocationCSR(A,ai,aj,aa);CHKERRQ(ierr);
21619afcda9SRichard Tran Mills 
21719afcda9SRichard Tran Mills   /* We now have an assembled sequential AIJ matrix created from copies of the exported arrays from the MKL matrix handle.
21819afcda9SRichard Tran Mills    * Now turn it into a MATSEQAIJMKL. */
21919afcda9SRichard Tran Mills   ierr = MatConvert_SeqAIJ_SeqAIJMKL(A,MATSEQAIJMKL,MAT_INPLACE_MATRIX,&A);CHKERRQ(ierr);
2206c87cf42SRichard Tran Mills 
22119afcda9SRichard Tran Mills   aijmkl = (Mat_SeqAIJMKL*) A->spptr;
22219afcda9SRichard Tran Mills   aijmkl->csrA = csrA;
2236c87cf42SRichard Tran Mills 
22419afcda9SRichard Tran Mills   /* The below code duplicates much of what is in MatSeqAIJKL_create_mkl_handle(). I dislike this code duplication, but
22519afcda9SRichard Tran Mills    * MatSeqAIJMKL_create_mkl_handle() cannot be used because we don't need to create a handle -- we've already got one,
22619afcda9SRichard Tran Mills    * and just need to be able to run the MKL optimization step. */
227f3fd1758SRichard Tran Mills   aijmkl->descr.type        = SPARSE_MATRIX_TYPE_GENERAL;
228f3fd1758SRichard Tran Mills   aijmkl->descr.mode        = SPARSE_FILL_MODE_LOWER;
229f3fd1758SRichard Tran Mills   aijmkl->descr.diag        = SPARSE_DIAG_NON_UNIT;
23019afcda9SRichard Tran Mills   stat = mkl_sparse_set_mv_hint(aijmkl->csrA,SPARSE_OPERATION_NON_TRANSPOSE,aijmkl->descr,1000);
23119afcda9SRichard Tran Mills   stat = mkl_sparse_set_memory_hint(aijmkl->csrA,SPARSE_MEMORY_AGGRESSIVE);
23219afcda9SRichard Tran Mills   stat = mkl_sparse_optimize(aijmkl->csrA);
2339c46acdfSRichard Tran Mills   if (stat != SPARSE_STATUS_SUCCESS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Intel MKL error: unable to set hints/complete mkl_sparse_optimize");
23419afcda9SRichard Tran Mills   aijmkl->sparse_optimized = PETSC_TRUE;
235e995cf24SRichard Tran Mills   ierr = PetscObjectStateGet((PetscObject)A,&(aijmkl->state));CHKERRQ(ierr);
23619afcda9SRichard Tran Mills 
23719afcda9SRichard Tran Mills   *mat = A;
23819afcda9SRichard Tran Mills   PetscFunctionReturn(0);
23919afcda9SRichard Tran Mills }
24019afcda9SRichard Tran Mills #endif /* PETSC_HAVE_MKL_SPARSE_OPTIMIZE */
24119afcda9SRichard Tran Mills 
242e8be1fc7SRichard Tran Mills /* MatSeqAIJMKL_update_from_mkl_handle() updates the matrix values array from the contents of the associated MKL sparse matrix handle.
243e8be1fc7SRichard Tran Mills  * This is needed after mkl_sparse_sp2m() with SPARSE_STAGE_FINALIZE_MULT has been used to compute new values of the matrix in
244e8be1fc7SRichard Tran Mills  * MatMatMultNumeric(). */
245e8be1fc7SRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_OPTIMIZE
246e8be1fc7SRichard Tran Mills PETSC_INTERN PetscErrorCode MatSeqAIJMKL_update_from_mkl_handle(Mat A)
247e8be1fc7SRichard Tran Mills {
248e8be1fc7SRichard Tran Mills   PetscInt            i;
249e8be1fc7SRichard Tran Mills   PetscInt            nrows,ncols;
250e8be1fc7SRichard Tran Mills   PetscInt            nz;
251e8be1fc7SRichard Tran Mills   PetscInt            *ai,*aj,*dummy;
252e8be1fc7SRichard Tran Mills   PetscScalar         *aa;
253e8be1fc7SRichard Tran Mills   PetscErrorCode      ierr;
254e8be1fc7SRichard Tran Mills   Mat_SeqAIJMKL       *aijmkl;
255e8be1fc7SRichard Tran Mills   sparse_status_t     stat;
256e8be1fc7SRichard Tran Mills   sparse_index_base_t indexing;
257e8be1fc7SRichard Tran Mills 
258e8be1fc7SRichard Tran Mills   aijmkl = (Mat_SeqAIJMKL*) A->spptr;
259e8be1fc7SRichard Tran Mills 
260e8be1fc7SRichard Tran Mills   /* Note: Must pass in &dummy below since MKL can't accept NULL for this output array we don't actually want. */
261e8be1fc7SRichard Tran Mills   stat = mkl_sparse_x_export_csr(aijmkl->csrA,&indexing,&nrows,&ncols,&ai,&dummy,&aj,&aa);
262e8be1fc7SRichard Tran Mills   if (stat != SPARSE_STATUS_SUCCESS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Intel MKL error: unable to complete mkl_sparse_x_export_csr()");
263e8be1fc7SRichard Tran Mills 
264e8be1fc7SRichard Tran Mills   /* We can't just do a copy from the arrays exported by MKL to those used for the PETSc AIJ storage, because the MKL and PETSc
265e8be1fc7SRichard Tran Mills    * representations differ in small ways (e.g., more explicit nonzeros per row due to preallocation). */
266e8be1fc7SRichard Tran Mills   for (i=0; i<nrows; i++) {
267e8be1fc7SRichard Tran Mills     nz = ai[i+1] - ai[i];
268e8be1fc7SRichard Tran Mills     ierr = MatSetValues_SeqAIJ(A, 1, &i, nz, aj+ai[i], aa+ai[i], INSERT_VALUES);CHKERRQ(ierr);
269e8be1fc7SRichard Tran Mills   }
270e8be1fc7SRichard Tran Mills 
271e8be1fc7SRichard Tran Mills   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
272e8be1fc7SRichard Tran Mills   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
273e8be1fc7SRichard Tran Mills 
274e995cf24SRichard Tran Mills   ierr = PetscObjectStateGet((PetscObject)A,&(aijmkl->state));CHKERRQ(ierr);
275e995cf24SRichard Tran Mills   /* We mark our matrix as having a valid, optimized MKL handle.
276e995cf24SRichard Tran Mills    * TODO: It is valid, but I am not sure if it is optimized. Need to ask MKL developers. */
277e995cf24SRichard Tran Mills   aijmkl->sparse_optimized = PETSC_TRUE;
278e995cf24SRichard Tran Mills 
279e8be1fc7SRichard Tran Mills   PetscFunctionReturn(0);
280e8be1fc7SRichard Tran Mills }
281e8be1fc7SRichard Tran Mills #endif /* PETSC_HAVE_MKL_SPARSE_OPTIMIZE */
282e8be1fc7SRichard Tran Mills 
2836e369cd5SRichard Tran Mills PetscErrorCode MatDuplicate_SeqAIJMKL(Mat A, MatDuplicateOption op, Mat *M)
2846e369cd5SRichard Tran Mills {
2856e369cd5SRichard Tran Mills   PetscErrorCode ierr;
2866e369cd5SRichard Tran Mills   Mat_SeqAIJMKL  *aijmkl;
2876e369cd5SRichard Tran Mills   Mat_SeqAIJMKL  *aijmkl_dest;
2886e369cd5SRichard Tran Mills 
2896e369cd5SRichard Tran Mills   PetscFunctionBegin;
2906e369cd5SRichard Tran Mills   ierr = MatDuplicate_SeqAIJ(A,op,M);CHKERRQ(ierr);
2916e369cd5SRichard Tran Mills   aijmkl      = (Mat_SeqAIJMKL*) A->spptr;
2926e369cd5SRichard Tran Mills   aijmkl_dest = (Mat_SeqAIJMKL*) (*M)->spptr;
2936e369cd5SRichard Tran Mills   ierr = PetscMemcpy(aijmkl_dest,aijmkl,sizeof(Mat_SeqAIJMKL));CHKERRQ(ierr);
2946e369cd5SRichard Tran Mills   aijmkl_dest->sparse_optimized = PETSC_FALSE;
2955b49642aSRichard Tran Mills   if (aijmkl->eager_inspection) {
2966e369cd5SRichard Tran Mills     ierr = MatSeqAIJMKL_create_mkl_handle(A);CHKERRQ(ierr);
2975b49642aSRichard Tran Mills   }
2986e369cd5SRichard Tran Mills   PetscFunctionReturn(0);
2996e369cd5SRichard Tran Mills }
3006e369cd5SRichard Tran Mills 
3016e369cd5SRichard Tran Mills PetscErrorCode MatAssemblyEnd_SeqAIJMKL(Mat A, MatAssemblyType mode)
3026e369cd5SRichard Tran Mills {
3036e369cd5SRichard Tran Mills   PetscErrorCode  ierr;
3046e369cd5SRichard Tran Mills   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
3055b49642aSRichard Tran Mills   Mat_SeqAIJMKL   *aijmkl;
3066e369cd5SRichard Tran Mills 
3076e369cd5SRichard Tran Mills   PetscFunctionBegin;
3086e369cd5SRichard Tran Mills   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
3096e369cd5SRichard Tran Mills 
3106e369cd5SRichard Tran Mills   /* Since a MATSEQAIJMKL matrix is really just a MATSEQAIJ with some
3116e369cd5SRichard Tran Mills    * extra information and some different methods, call the AssemblyEnd
3126e369cd5SRichard Tran Mills    * routine for a MATSEQAIJ.
3136e369cd5SRichard Tran Mills    * I'm not sure if this is the best way to do this, but it avoids
314d96e85feSRichard Tran Mills    * a lot of code duplication. */
3156e369cd5SRichard Tran Mills   a->inode.use = PETSC_FALSE;  /* Must disable: otherwise the MKL routines won't get used. */
3166e369cd5SRichard Tran Mills   ierr = MatAssemblyEnd_SeqAIJ(A, mode);CHKERRQ(ierr);
3176e369cd5SRichard Tran Mills 
3185b49642aSRichard Tran Mills   /* If the user has requested "eager" inspection, create the optimized MKL sparse handle (if needed; the function checks).
3195b49642aSRichard Tran Mills    * (The default is to do "lazy" inspection, deferring this until something like MatMult() is called.) */
3205b49642aSRichard Tran Mills   aijmkl = (Mat_SeqAIJMKL*) A->spptr;
3215b49642aSRichard Tran Mills   if (aijmkl->eager_inspection) {
3226e369cd5SRichard Tran Mills     ierr = MatSeqAIJMKL_create_mkl_handle(A);CHKERRQ(ierr);
3235b49642aSRichard Tran Mills   }
324df555b71SRichard Tran Mills 
3254a2a386eSRichard Tran Mills   PetscFunctionReturn(0);
3264a2a386eSRichard Tran Mills }
3274a2a386eSRichard Tran Mills 
3284a2a386eSRichard Tran Mills PetscErrorCode MatMult_SeqAIJMKL(Mat A,Vec xx,Vec yy)
3294a2a386eSRichard Tran Mills {
3304a2a386eSRichard Tran Mills   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
3314a2a386eSRichard Tran Mills   const PetscScalar *x;
3324a2a386eSRichard Tran Mills   PetscScalar       *y;
3334a2a386eSRichard Tran Mills   const MatScalar   *aa;
3344a2a386eSRichard Tran Mills   PetscErrorCode    ierr;
3354a2a386eSRichard Tran Mills   PetscInt          m=A->rmap->n;
336db63039fSRichard Tran Mills   PetscInt          n=A->cmap->n;
337db63039fSRichard Tran Mills   PetscScalar       alpha = 1.0;
338db63039fSRichard Tran Mills   PetscScalar       beta = 0.0;
3394a2a386eSRichard Tran Mills   const PetscInt    *aj,*ai;
340db63039fSRichard Tran Mills   char              matdescra[6];
341db63039fSRichard Tran Mills 
3424a2a386eSRichard Tran Mills 
3434a2a386eSRichard Tran Mills   /* Variables not in MatMult_SeqAIJ. */
344ff03dc53SRichard Tran Mills   char transa = 'n';  /* Used to indicate to MKL that we are not computing the transpose product. */
345ff03dc53SRichard Tran Mills 
346ff03dc53SRichard Tran Mills   PetscFunctionBegin;
347db63039fSRichard Tran Mills   matdescra[0] = 'g';  /* Indicates to MKL that we using a general CSR matrix. */
348db63039fSRichard Tran Mills   matdescra[3] = 'c';  /* Indicates to MKL that we use C-style (0-based) indexing. */
349ff03dc53SRichard Tran Mills   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
350ff03dc53SRichard Tran Mills   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
351ff03dc53SRichard Tran Mills   aj   = a->j;  /* aj[k] gives column index for element aa[k]. */
352ff03dc53SRichard Tran Mills   aa   = a->a;  /* Nonzero elements stored row-by-row. */
353ff03dc53SRichard Tran Mills   ai   = a->i;  /* ai[k] is the position in aa and aj where row k starts. */
354ff03dc53SRichard Tran Mills 
355ff03dc53SRichard Tran Mills   /* Call MKL sparse BLAS routine to do the MatMult. */
356db63039fSRichard Tran Mills   mkl_xcsrmv(&transa,&m,&n,&alpha,matdescra,aa,aj,ai,ai+1,x,&beta,y);
357ff03dc53SRichard Tran Mills 
358ff03dc53SRichard Tran Mills   ierr = PetscLogFlops(2.0*a->nz - a->nonzerorowcnt);CHKERRQ(ierr);
359ff03dc53SRichard Tran Mills   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
360ff03dc53SRichard Tran Mills   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
361ff03dc53SRichard Tran Mills   PetscFunctionReturn(0);
362ff03dc53SRichard Tran Mills }
363ff03dc53SRichard Tran Mills 
364d995685eSRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_OPTIMIZE
365df555b71SRichard Tran Mills PetscErrorCode MatMult_SeqAIJMKL_SpMV2(Mat A,Vec xx,Vec yy)
366df555b71SRichard Tran Mills {
367df555b71SRichard Tran Mills   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
368df555b71SRichard Tran Mills   Mat_SeqAIJMKL     *aijmkl=(Mat_SeqAIJMKL*)A->spptr;
369df555b71SRichard Tran Mills   const PetscScalar *x;
370df555b71SRichard Tran Mills   PetscScalar       *y;
371df555b71SRichard Tran Mills   PetscErrorCode    ierr;
372df555b71SRichard Tran Mills   sparse_status_t   stat = SPARSE_STATUS_SUCCESS;
373551aa5c8SRichard Tran Mills   PetscObjectState  state;
374df555b71SRichard Tran Mills 
375df555b71SRichard Tran Mills   PetscFunctionBegin;
376df555b71SRichard Tran Mills 
37738987b35SRichard Tran Mills   /* If there are no nonzero entries, zero yy and return immediately. */
37838987b35SRichard Tran Mills   if(!a->nz) {
37938987b35SRichard Tran Mills     PetscInt i;
38038987b35SRichard Tran Mills     PetscInt m=A->rmap->n;
38138987b35SRichard Tran Mills     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
38238987b35SRichard Tran Mills     for (i=0; i<m; i++) {
38338987b35SRichard Tran Mills       y[i] = 0.0;
38438987b35SRichard Tran Mills     }
38538987b35SRichard Tran Mills     ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
38638987b35SRichard Tran Mills     PetscFunctionReturn(0);
38738987b35SRichard Tran Mills   }
388f36dfe3fSRichard Tran Mills 
389df555b71SRichard Tran Mills   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
390df555b71SRichard Tran Mills   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
391df555b71SRichard Tran Mills 
3923fa15762SRichard Tran Mills   /* In some cases, we get to this point without mkl_sparse_optimize() having been called, so we check and then call
3933fa15762SRichard Tran Mills    * it if needed. Eventually, when everything in PETSc is properly updating the matrix state, we should probably
3943fa15762SRichard Tran Mills    * take a "lazy" approach to creation/updating of the MKL matrix handle and plan to always do it here (when needed). */
395551aa5c8SRichard Tran Mills   ierr = PetscObjectStateGet((PetscObject)A,&state);CHKERRQ(ierr);
396551aa5c8SRichard Tran Mills   if (!aijmkl->sparse_optimized || aijmkl->state != state) {
3973fa15762SRichard Tran Mills     MatSeqAIJMKL_create_mkl_handle(A);
3983fa15762SRichard Tran Mills   }
3993fa15762SRichard Tran Mills 
400df555b71SRichard Tran Mills   /* Call MKL SpMV2 executor routine to do the MatMult. */
401df555b71SRichard Tran Mills   stat = mkl_sparse_x_mv(SPARSE_OPERATION_NON_TRANSPOSE,1.0,aijmkl->csrA,aijmkl->descr,x,0.0,y);
4029c46acdfSRichard Tran Mills   if (stat != SPARSE_STATUS_SUCCESS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Intel MKL error: error in mkl_sparse_x_mv");
403df555b71SRichard Tran Mills 
404df555b71SRichard Tran Mills   ierr = PetscLogFlops(2.0*a->nz - a->nonzerorowcnt);CHKERRQ(ierr);
405df555b71SRichard Tran Mills   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
406df555b71SRichard Tran Mills   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
407df555b71SRichard Tran Mills   PetscFunctionReturn(0);
408df555b71SRichard Tran Mills }
409d995685eSRichard Tran Mills #endif /* PETSC_HAVE_MKL_SPARSE_OPTIMIZE */
410df555b71SRichard Tran Mills 
411ff03dc53SRichard Tran Mills PetscErrorCode MatMultTranspose_SeqAIJMKL(Mat A,Vec xx,Vec yy)
412ff03dc53SRichard Tran Mills {
413ff03dc53SRichard Tran Mills   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
414ff03dc53SRichard Tran Mills   const PetscScalar *x;
415ff03dc53SRichard Tran Mills   PetscScalar       *y;
416ff03dc53SRichard Tran Mills   const MatScalar   *aa;
417ff03dc53SRichard Tran Mills   PetscErrorCode    ierr;
418ff03dc53SRichard Tran Mills   PetscInt          m=A->rmap->n;
419db63039fSRichard Tran Mills   PetscInt          n=A->cmap->n;
420db63039fSRichard Tran Mills   PetscScalar       alpha = 1.0;
421db63039fSRichard Tran Mills   PetscScalar       beta = 0.0;
422ff03dc53SRichard Tran Mills   const PetscInt    *aj,*ai;
423db63039fSRichard Tran Mills   char              matdescra[6];
424ff03dc53SRichard Tran Mills 
425ff03dc53SRichard Tran Mills   /* Variables not in MatMultTranspose_SeqAIJ. */
426ff03dc53SRichard Tran Mills   char transa = 't';  /* Used to indicate to MKL that we are computing the transpose product. */
4274a2a386eSRichard Tran Mills 
4284a2a386eSRichard Tran Mills   PetscFunctionBegin;
429969800c5SRichard Tran Mills   matdescra[0] = 'g';  /* Indicates to MKL that we using a general CSR matrix. */
430969800c5SRichard Tran Mills   matdescra[3] = 'c';  /* Indicates to MKL that we use C-style (0-based) indexing. */
4314a2a386eSRichard Tran Mills   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
4324a2a386eSRichard Tran Mills   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
4334a2a386eSRichard Tran Mills   aj   = a->j;  /* aj[k] gives column index for element aa[k]. */
4344a2a386eSRichard Tran Mills   aa   = a->a;  /* Nonzero elements stored row-by-row. */
4354a2a386eSRichard Tran Mills   ai   = a->i;  /* ai[k] is the position in aa and aj where row k starts. */
4364a2a386eSRichard Tran Mills 
4374a2a386eSRichard Tran Mills   /* Call MKL sparse BLAS routine to do the MatMult. */
438db63039fSRichard Tran Mills   mkl_xcsrmv(&transa,&m,&n,&alpha,matdescra,aa,aj,ai,ai+1,x,&beta,y);
4394a2a386eSRichard Tran Mills 
4404a2a386eSRichard Tran Mills   ierr = PetscLogFlops(2.0*a->nz - a->nonzerorowcnt);CHKERRQ(ierr);
4414a2a386eSRichard Tran Mills   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
4424a2a386eSRichard Tran Mills   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
4434a2a386eSRichard Tran Mills   PetscFunctionReturn(0);
4444a2a386eSRichard Tran Mills }
4454a2a386eSRichard Tran Mills 
446d995685eSRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_OPTIMIZE
447df555b71SRichard Tran Mills PetscErrorCode MatMultTranspose_SeqAIJMKL_SpMV2(Mat A,Vec xx,Vec yy)
448df555b71SRichard Tran Mills {
449df555b71SRichard Tran Mills   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
450df555b71SRichard Tran Mills   Mat_SeqAIJMKL     *aijmkl=(Mat_SeqAIJMKL*)A->spptr;
451df555b71SRichard Tran Mills   const PetscScalar *x;
452df555b71SRichard Tran Mills   PetscScalar       *y;
453df555b71SRichard Tran Mills   PetscErrorCode    ierr;
4540632b357SRichard Tran Mills   sparse_status_t   stat;
455551aa5c8SRichard Tran Mills   PetscObjectState  state;
456df555b71SRichard Tran Mills 
457df555b71SRichard Tran Mills   PetscFunctionBegin;
458df555b71SRichard Tran Mills 
45938987b35SRichard Tran Mills   /* If there are no nonzero entries, zero yy and return immediately. */
46038987b35SRichard Tran Mills   if(!a->nz) {
46138987b35SRichard Tran Mills     PetscInt i;
46238987b35SRichard Tran Mills     PetscInt n=A->cmap->n;
46338987b35SRichard Tran Mills     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
46438987b35SRichard Tran Mills     for (i=0; i<n; i++) {
46538987b35SRichard Tran Mills       y[i] = 0.0;
46638987b35SRichard Tran Mills     }
46738987b35SRichard Tran Mills     ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
46838987b35SRichard Tran Mills     PetscFunctionReturn(0);
46938987b35SRichard Tran Mills   }
470f36dfe3fSRichard Tran Mills 
471df555b71SRichard Tran Mills   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
472df555b71SRichard Tran Mills   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
473df555b71SRichard Tran Mills 
4743fa15762SRichard Tran Mills   /* In some cases, we get to this point without mkl_sparse_optimize() having been called, so we check and then call
4753fa15762SRichard Tran Mills    * it if needed. Eventually, when everything in PETSc is properly updating the matrix state, we should probably
4763fa15762SRichard Tran Mills    * take a "lazy" approach to creation/updating of the MKL matrix handle and plan to always do it here (when needed). */
477551aa5c8SRichard Tran Mills   ierr = PetscObjectStateGet((PetscObject)A,&state);CHKERRQ(ierr);
478551aa5c8SRichard Tran Mills   if (!aijmkl->sparse_optimized || aijmkl->state != state) {
4793fa15762SRichard Tran Mills     MatSeqAIJMKL_create_mkl_handle(A);
4803fa15762SRichard Tran Mills   }
4813fa15762SRichard Tran Mills 
482df555b71SRichard Tran Mills   /* Call MKL SpMV2 executor routine to do the MatMultTranspose. */
483df555b71SRichard Tran Mills   stat = mkl_sparse_x_mv(SPARSE_OPERATION_TRANSPOSE,1.0,aijmkl->csrA,aijmkl->descr,x,0.0,y);
4849c46acdfSRichard Tran Mills   if (stat != SPARSE_STATUS_SUCCESS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Intel MKL error: error in mkl_sparse_x_mv");
485df555b71SRichard Tran Mills 
486df555b71SRichard Tran Mills   ierr = PetscLogFlops(2.0*a->nz - a->nonzerorowcnt);CHKERRQ(ierr);
487df555b71SRichard Tran Mills   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
488df555b71SRichard Tran Mills   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
489df555b71SRichard Tran Mills   PetscFunctionReturn(0);
490df555b71SRichard Tran Mills }
491d995685eSRichard Tran Mills #endif /* PETSC_HAVE_MKL_SPARSE_OPTIMIZE */
492df555b71SRichard Tran Mills 
4934a2a386eSRichard Tran Mills PetscErrorCode MatMultAdd_SeqAIJMKL(Mat A,Vec xx,Vec yy,Vec zz)
4944a2a386eSRichard Tran Mills {
4954a2a386eSRichard Tran Mills   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
4964a2a386eSRichard Tran Mills   const PetscScalar *x;
4974a2a386eSRichard Tran Mills   PetscScalar       *y,*z;
4984a2a386eSRichard Tran Mills   const MatScalar   *aa;
4994a2a386eSRichard Tran Mills   PetscErrorCode    ierr;
5004a2a386eSRichard Tran Mills   PetscInt          m=A->rmap->n;
501db63039fSRichard Tran Mills   PetscInt          n=A->cmap->n;
5024a2a386eSRichard Tran Mills   const PetscInt    *aj,*ai;
5034a2a386eSRichard Tran Mills   PetscInt          i;
5044a2a386eSRichard Tran Mills 
505ff03dc53SRichard Tran Mills   /* Variables not in MatMultAdd_SeqAIJ. */
506ff03dc53SRichard Tran Mills   char              transa = 'n';  /* Used to indicate to MKL that we are not computing the transpose product. */
507a84739b8SRichard Tran Mills   PetscScalar       alpha = 1.0;
508db63039fSRichard Tran Mills   PetscScalar       beta;
509a84739b8SRichard Tran Mills   char              matdescra[6];
510ff03dc53SRichard Tran Mills 
511ff03dc53SRichard Tran Mills   PetscFunctionBegin;
512a84739b8SRichard Tran Mills   matdescra[0] = 'g';  /* Indicates to MKL that we using a general CSR matrix. */
513a84739b8SRichard Tran Mills   matdescra[3] = 'c';  /* Indicates to MKL that we use C-style (0-based) indexing. */
514a84739b8SRichard Tran Mills 
515ff03dc53SRichard Tran Mills   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
516ff03dc53SRichard Tran Mills   ierr = VecGetArrayPair(yy,zz,&y,&z);CHKERRQ(ierr);
517ff03dc53SRichard Tran Mills   aj   = a->j;  /* aj[k] gives column index for element aa[k]. */
518ff03dc53SRichard Tran Mills   aa   = a->a;  /* Nonzero elements stored row-by-row. */
519ff03dc53SRichard Tran Mills   ai   = a->i;  /* ai[k] is the position in aa and aj where row k starts. */
520ff03dc53SRichard Tran Mills 
521ff03dc53SRichard Tran Mills   /* Call MKL sparse BLAS routine to do the MatMult. */
522a84739b8SRichard Tran Mills   if (zz == yy) {
523a84739b8SRichard Tran Mills     /* If zz and yy are the same vector, we can use MKL's mkl_xcsrmv(), which calculates y = alpha*A*x + beta*y. */
524db63039fSRichard Tran Mills     beta = 1.0;
525db63039fSRichard Tran Mills     mkl_xcsrmv(&transa,&m,&n,&alpha,matdescra,aa,aj,ai,ai+1,x,&beta,z);
526a84739b8SRichard Tran Mills   } else {
527db63039fSRichard Tran Mills     /* zz and yy are different vectors, so call MKL's mkl_xcsrmv() with beta=0, then add the result to z.
528db63039fSRichard Tran Mills      * MKL sparse BLAS does not have a MatMultAdd equivalent. */
529db63039fSRichard Tran Mills     beta = 0.0;
530db63039fSRichard Tran Mills     mkl_xcsrmv(&transa,&m,&n,&alpha,matdescra,aa,aj,ai,ai+1,x,&beta,z);
531ff03dc53SRichard Tran Mills     for (i=0; i<m; i++) {
532ff03dc53SRichard Tran Mills       z[i] += y[i];
533ff03dc53SRichard Tran Mills     }
534a84739b8SRichard Tran Mills   }
535ff03dc53SRichard Tran Mills 
536ff03dc53SRichard Tran Mills   ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr);
537ff03dc53SRichard Tran Mills   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
538ff03dc53SRichard Tran Mills   ierr = VecRestoreArrayPair(yy,zz,&y,&z);CHKERRQ(ierr);
539ff03dc53SRichard Tran Mills   PetscFunctionReturn(0);
540ff03dc53SRichard Tran Mills }
541ff03dc53SRichard Tran Mills 
542d995685eSRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_OPTIMIZE
543df555b71SRichard Tran Mills PetscErrorCode MatMultAdd_SeqAIJMKL_SpMV2(Mat A,Vec xx,Vec yy,Vec zz)
544df555b71SRichard Tran Mills {
545df555b71SRichard Tran Mills   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
546df555b71SRichard Tran Mills   Mat_SeqAIJMKL     *aijmkl=(Mat_SeqAIJMKL*)A->spptr;
547df555b71SRichard Tran Mills   const PetscScalar *x;
548df555b71SRichard Tran Mills   PetscScalar       *y,*z;
549df555b71SRichard Tran Mills   PetscErrorCode    ierr;
550df555b71SRichard Tran Mills   PetscInt          m=A->rmap->n;
551df555b71SRichard Tran Mills   PetscInt          i;
552df555b71SRichard Tran Mills 
553df555b71SRichard Tran Mills   /* Variables not in MatMultAdd_SeqAIJ. */
554df555b71SRichard Tran Mills   sparse_status_t   stat = SPARSE_STATUS_SUCCESS;
555551aa5c8SRichard Tran Mills   PetscObjectState  state;
556df555b71SRichard Tran Mills 
557df555b71SRichard Tran Mills   PetscFunctionBegin;
558df555b71SRichard Tran Mills 
55938987b35SRichard Tran Mills   /* If there are no nonzero entries, set zz = yy and return immediately. */
56038987b35SRichard Tran Mills   if(!a->nz) {
56138987b35SRichard Tran Mills     PetscInt i;
56238987b35SRichard Tran Mills     ierr = VecGetArrayPair(yy,zz,&y,&z);CHKERRQ(ierr);
56338987b35SRichard Tran Mills     for (i=0; i<m; i++) {
56438987b35SRichard Tran Mills       z[i] = y[i];
56538987b35SRichard Tran Mills     }
56638987b35SRichard Tran Mills     ierr = VecRestoreArrayPair(yy,zz,&y,&z);CHKERRQ(ierr);
56738987b35SRichard Tran Mills     PetscFunctionReturn(0);
56838987b35SRichard Tran Mills   }
569df555b71SRichard Tran Mills 
570df555b71SRichard Tran Mills   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
571df555b71SRichard Tran Mills   ierr = VecGetArrayPair(yy,zz,&y,&z);CHKERRQ(ierr);
572df555b71SRichard Tran Mills 
5733fa15762SRichard Tran Mills   /* In some cases, we get to this point without mkl_sparse_optimize() having been called, so we check and then call
5743fa15762SRichard Tran Mills    * it if needed. Eventually, when everything in PETSc is properly updating the matrix state, we should probably
5753fa15762SRichard Tran Mills    * take a "lazy" approach to creation/updating of the MKL matrix handle and plan to always do it here (when needed). */
576551aa5c8SRichard Tran Mills   ierr = PetscObjectStateGet((PetscObject)A,&state);CHKERRQ(ierr);
577551aa5c8SRichard Tran Mills   if (!aijmkl->sparse_optimized || aijmkl->state != state) {
5783fa15762SRichard Tran Mills     MatSeqAIJMKL_create_mkl_handle(A);
5793fa15762SRichard Tran Mills   }
5803fa15762SRichard Tran Mills 
581df555b71SRichard Tran Mills   /* Call MKL sparse BLAS routine to do the MatMult. */
582df555b71SRichard Tran Mills   if (zz == yy) {
583df555b71SRichard Tran Mills     /* If zz and yy are the same vector, we can use mkl_sparse_x_mv, which calculates y = alpha*A*x + beta*y,
584df555b71SRichard Tran Mills      * with alpha and beta both set to 1.0. */
585db63039fSRichard Tran Mills     stat = mkl_sparse_x_mv(SPARSE_OPERATION_NON_TRANSPOSE,1.0,aijmkl->csrA,aijmkl->descr,x,1.0,z);
5869c46acdfSRichard Tran Mills     if (stat != SPARSE_STATUS_SUCCESS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Intel MKL error: error in mkl_sparse_x_mv");
587df555b71SRichard Tran Mills   } else {
588df555b71SRichard Tran Mills     /* zz and yy are different vectors, so we call mkl_sparse_x_mv with alpha=1.0 and beta=0.0, and then
589df555b71SRichard Tran Mills      * we add the contents of vector yy to the result; MKL sparse BLAS does not have a MatMultAdd equivalent. */
590db63039fSRichard Tran Mills     stat = mkl_sparse_x_mv(SPARSE_OPERATION_NON_TRANSPOSE,1.0,aijmkl->csrA,aijmkl->descr,x,0.0,z);
5919c46acdfSRichard Tran Mills     if (stat != SPARSE_STATUS_SUCCESS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Intel MKL error: error in mkl_sparse_x_mv");
592df555b71SRichard Tran Mills     for (i=0; i<m; i++) {
593df555b71SRichard Tran Mills       z[i] += y[i];
594df555b71SRichard Tran Mills     }
595df555b71SRichard Tran Mills   }
596df555b71SRichard Tran Mills 
597df555b71SRichard Tran Mills   ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr);
598df555b71SRichard Tran Mills   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
599df555b71SRichard Tran Mills   ierr = VecRestoreArrayPair(yy,zz,&y,&z);CHKERRQ(ierr);
600df555b71SRichard Tran Mills   PetscFunctionReturn(0);
601df555b71SRichard Tran Mills }
602d995685eSRichard Tran Mills #endif /* PETSC_HAVE_MKL_SPARSE_OPTIMIZE */
603df555b71SRichard Tran Mills 
604ff03dc53SRichard Tran Mills PetscErrorCode MatMultTransposeAdd_SeqAIJMKL(Mat A,Vec xx,Vec yy,Vec zz)
605ff03dc53SRichard Tran Mills {
606ff03dc53SRichard Tran Mills   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
607ff03dc53SRichard Tran Mills   const PetscScalar *x;
608ff03dc53SRichard Tran Mills   PetscScalar       *y,*z;
609ff03dc53SRichard Tran Mills   const MatScalar   *aa;
610ff03dc53SRichard Tran Mills   PetscErrorCode    ierr;
611ff03dc53SRichard Tran Mills   PetscInt          m=A->rmap->n;
612db63039fSRichard Tran Mills   PetscInt          n=A->cmap->n;
613ff03dc53SRichard Tran Mills   const PetscInt    *aj,*ai;
614ff03dc53SRichard Tran Mills   PetscInt          i;
615ff03dc53SRichard Tran Mills 
616ff03dc53SRichard Tran Mills   /* Variables not in MatMultTransposeAdd_SeqAIJ. */
617ff03dc53SRichard Tran Mills   char transa = 't';  /* Used to indicate to MKL that we are computing the transpose product. */
618a84739b8SRichard Tran Mills   PetscScalar       alpha = 1.0;
619db63039fSRichard Tran Mills   PetscScalar       beta;
620a84739b8SRichard Tran Mills   char              matdescra[6];
6214a2a386eSRichard Tran Mills 
6224a2a386eSRichard Tran Mills   PetscFunctionBegin;
623a84739b8SRichard Tran Mills   matdescra[0] = 'g';  /* Indicates to MKL that we using a general CSR matrix. */
624a84739b8SRichard Tran Mills   matdescra[3] = 'c';  /* Indicates to MKL that we use C-style (0-based) indexing. */
625a84739b8SRichard Tran Mills 
6264a2a386eSRichard Tran Mills   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
6274a2a386eSRichard Tran Mills   ierr = VecGetArrayPair(yy,zz,&y,&z);CHKERRQ(ierr);
6284a2a386eSRichard Tran Mills   aj   = a->j;  /* aj[k] gives column index for element aa[k]. */
6294a2a386eSRichard Tran Mills   aa   = a->a;  /* Nonzero elements stored row-by-row. */
6304a2a386eSRichard Tran Mills   ai   = a->i;  /* ai[k] is the position in aa and aj where row k starts. */
6314a2a386eSRichard Tran Mills 
6324a2a386eSRichard Tran Mills   /* Call MKL sparse BLAS routine to do the MatMult. */
633a84739b8SRichard Tran Mills   if (zz == yy) {
634a84739b8SRichard Tran Mills     /* If zz and yy are the same vector, we can use MKL's mkl_xcsrmv(), which calculates y = alpha*A*x + beta*y. */
635db63039fSRichard Tran Mills     beta = 1.0;
636969800c5SRichard Tran Mills     mkl_xcsrmv(&transa,&m,&n,&alpha,matdescra,aa,aj,ai,ai+1,x,&beta,z);
637a84739b8SRichard Tran Mills   } else {
638db63039fSRichard Tran Mills     /* zz and yy are different vectors, so call MKL's mkl_xcsrmv() with beta=0, then add the result to z.
639db63039fSRichard Tran Mills      * MKL sparse BLAS does not have a MatMultAdd equivalent. */
640db63039fSRichard Tran Mills     beta = 0.0;
641db63039fSRichard Tran Mills     mkl_xcsrmv(&transa,&m,&n,&alpha,matdescra,aa,aj,ai,ai+1,x,&beta,z);
642969800c5SRichard Tran Mills     for (i=0; i<n; i++) {
6434a2a386eSRichard Tran Mills       z[i] += y[i];
6444a2a386eSRichard Tran Mills     }
645a84739b8SRichard Tran Mills   }
6464a2a386eSRichard Tran Mills 
6474a2a386eSRichard Tran Mills   ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr);
6484a2a386eSRichard Tran Mills   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
6494a2a386eSRichard Tran Mills   ierr = VecRestoreArrayPair(yy,zz,&y,&z);CHKERRQ(ierr);
6504a2a386eSRichard Tran Mills   PetscFunctionReturn(0);
6514a2a386eSRichard Tran Mills }
6524a2a386eSRichard Tran Mills 
653d995685eSRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_OPTIMIZE
654df555b71SRichard Tran Mills PetscErrorCode MatMultTransposeAdd_SeqAIJMKL_SpMV2(Mat A,Vec xx,Vec yy,Vec zz)
655df555b71SRichard Tran Mills {
656df555b71SRichard Tran Mills   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
657df555b71SRichard Tran Mills   Mat_SeqAIJMKL     *aijmkl=(Mat_SeqAIJMKL*)A->spptr;
658df555b71SRichard Tran Mills   const PetscScalar *x;
659df555b71SRichard Tran Mills   PetscScalar       *y,*z;
660df555b71SRichard Tran Mills   PetscErrorCode    ierr;
661969800c5SRichard Tran Mills   PetscInt          n=A->cmap->n;
662df555b71SRichard Tran Mills   PetscInt          i;
663551aa5c8SRichard Tran Mills   PetscObjectState  state;
664df555b71SRichard Tran Mills 
665df555b71SRichard Tran Mills   /* Variables not in MatMultTransposeAdd_SeqAIJ. */
666df555b71SRichard Tran Mills   sparse_status_t stat = SPARSE_STATUS_SUCCESS;
667df555b71SRichard Tran Mills 
668df555b71SRichard Tran Mills   PetscFunctionBegin;
669df555b71SRichard Tran Mills 
67038987b35SRichard Tran Mills   /* If there are no nonzero entries, set zz = yy and return immediately. */
67138987b35SRichard Tran Mills   if(!a->nz) {
67238987b35SRichard Tran Mills     PetscInt i;
67338987b35SRichard Tran Mills     ierr = VecGetArrayPair(yy,zz,&y,&z);CHKERRQ(ierr);
67438987b35SRichard Tran Mills     for (i=0; i<n; i++) {
67538987b35SRichard Tran Mills       z[i] = y[i];
67638987b35SRichard Tran Mills     }
67738987b35SRichard Tran Mills     ierr = VecRestoreArrayPair(yy,zz,&y,&z);CHKERRQ(ierr);
67838987b35SRichard Tran Mills     PetscFunctionReturn(0);
67938987b35SRichard Tran Mills   }
680f36dfe3fSRichard Tran Mills 
681df555b71SRichard Tran Mills   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
682df555b71SRichard Tran Mills   ierr = VecGetArrayPair(yy,zz,&y,&z);CHKERRQ(ierr);
683df555b71SRichard Tran Mills 
6843fa15762SRichard Tran Mills   /* In some cases, we get to this point without mkl_sparse_optimize() having been called, so we check and then call
6853fa15762SRichard Tran Mills    * it if needed. Eventually, when everything in PETSc is properly updating the matrix state, we should probably
6863fa15762SRichard Tran Mills    * take a "lazy" approach to creation/updating of the MKL matrix handle and plan to always do it here (when needed). */
687551aa5c8SRichard Tran Mills   ierr = PetscObjectStateGet((PetscObject)A,&state);CHKERRQ(ierr);
688551aa5c8SRichard Tran Mills   if (!aijmkl->sparse_optimized || aijmkl->state != state) {
6893fa15762SRichard Tran Mills     MatSeqAIJMKL_create_mkl_handle(A);
6903fa15762SRichard Tran Mills   }
6913fa15762SRichard Tran Mills 
692df555b71SRichard Tran Mills   /* Call MKL sparse BLAS routine to do the MatMult. */
693df555b71SRichard Tran Mills   if (zz == yy) {
694df555b71SRichard Tran Mills     /* If zz and yy are the same vector, we can use mkl_sparse_x_mv, which calculates y = alpha*A*x + beta*y,
695df555b71SRichard Tran Mills      * with alpha and beta both set to 1.0. */
696db63039fSRichard Tran Mills     stat = mkl_sparse_x_mv(SPARSE_OPERATION_TRANSPOSE,1.0,aijmkl->csrA,aijmkl->descr,x,1.0,z);
6979c46acdfSRichard Tran Mills     if (stat != SPARSE_STATUS_SUCCESS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Intel MKL error: error in mkl_sparse_x_mv");
698df555b71SRichard Tran Mills   } else {
699df555b71SRichard Tran Mills     /* zz and yy are different vectors, so we call mkl_sparse_x_mv with alpha=1.0 and beta=0.0, and then
700df555b71SRichard Tran Mills      * we add the contents of vector yy to the result; MKL sparse BLAS does not have a MatMultAdd equivalent. */
701db63039fSRichard Tran Mills     stat = mkl_sparse_x_mv(SPARSE_OPERATION_TRANSPOSE,1.0,aijmkl->csrA,aijmkl->descr,x,0.0,z);
7029c46acdfSRichard Tran Mills     if (stat != SPARSE_STATUS_SUCCESS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Intel MKL error: error in mkl_sparse_x_mv");
703969800c5SRichard Tran Mills     for (i=0; i<n; i++) {
704df555b71SRichard Tran Mills       z[i] += y[i];
705df555b71SRichard Tran Mills     }
706df555b71SRichard Tran Mills   }
707df555b71SRichard Tran Mills 
708df555b71SRichard Tran Mills   ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr);
709df555b71SRichard Tran Mills   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
710df555b71SRichard Tran Mills   ierr = VecRestoreArrayPair(yy,zz,&y,&z);CHKERRQ(ierr);
711df555b71SRichard Tran Mills   PetscFunctionReturn(0);
712df555b71SRichard Tran Mills }
713d995685eSRichard Tran Mills #endif /* PETSC_HAVE_MKL_SPARSE_OPTIMIZE */
714df555b71SRichard Tran Mills 
71545fbe478SRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_OPTIMIZE
716aab60f1bSRichard Tran Mills /* Note that this code currently doesn't actually get used when MatMatMult() is called with MAT_REUSE_MATRIX, because
717aab60f1bSRichard Tran Mills  * the MatMatMult() interface code calls MatMatMultNumeric() in this case.
7183ecbffd0SRichard Tran Mills  * For releases of MKL prior to version 18, update 2:
719aab60f1bSRichard Tran Mills  * MKL has no notion of separately callable symbolic vs. numeric phases of sparse matrix-matrix multiply, so in the
720aab60f1bSRichard Tran Mills  * MAT_REUSE_MATRIX case, the SeqAIJ routines end up being used. Even though this means that the (hopefully more
721aab60f1bSRichard Tran Mills  * optimized) MKL routines do not get used, this probably is best because the MKL routines would waste time re-computing
722aab60f1bSRichard Tran Mills  * the symbolic portion, whereas the native PETSc SeqAIJ routines will avoid this. */
72345fbe478SRichard Tran Mills PetscErrorCode MatMatMult_SeqAIJMKL_SeqAIJMKL_SpMV2(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat*C)
72445fbe478SRichard Tran Mills {
72545fbe478SRichard Tran Mills   Mat_SeqAIJMKL    *a, *b;
72645fbe478SRichard Tran Mills   sparse_matrix_t  csrA, csrB, csrC;
72745fbe478SRichard Tran Mills   PetscErrorCode   ierr;
72845fbe478SRichard Tran Mills   sparse_status_t  stat = SPARSE_STATUS_SUCCESS;
729551aa5c8SRichard Tran Mills   PetscObjectState state;
73045fbe478SRichard Tran Mills 
73145fbe478SRichard Tran Mills   PetscFunctionBegin;
73245fbe478SRichard Tran Mills   a = (Mat_SeqAIJMKL*)A->spptr;
73345fbe478SRichard Tran Mills   b = (Mat_SeqAIJMKL*)B->spptr;
734551aa5c8SRichard Tran Mills   ierr = PetscObjectStateGet((PetscObject)A,&state);CHKERRQ(ierr);
735551aa5c8SRichard Tran Mills   if (!a->sparse_optimized || a->state != state) {
73645fbe478SRichard Tran Mills     MatSeqAIJMKL_create_mkl_handle(A);
73745fbe478SRichard Tran Mills   }
738551aa5c8SRichard Tran Mills   ierr = PetscObjectStateGet((PetscObject)B,&state);CHKERRQ(ierr);
739551aa5c8SRichard Tran Mills   if (!b->sparse_optimized || b->state != state) {
74045fbe478SRichard Tran Mills     MatSeqAIJMKL_create_mkl_handle(B);
74145fbe478SRichard Tran Mills   }
74245fbe478SRichard Tran Mills   csrA = a->csrA;
74345fbe478SRichard Tran Mills   csrB = b->csrA;
74445fbe478SRichard Tran Mills 
74545fbe478SRichard Tran Mills   stat = mkl_sparse_spmm(SPARSE_OPERATION_NON_TRANSPOSE,csrA,csrB,&csrC);
7469c46acdfSRichard Tran Mills   if (stat != SPARSE_STATUS_SUCCESS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Intel MKL error: unable to complete sparse matrix-matrix multiply");
74745fbe478SRichard Tran Mills 
7486c87cf42SRichard Tran Mills   ierr = MatSeqAIJMKL_create_from_mkl_handle(PETSC_COMM_SELF,csrC,scall,C);CHKERRQ(ierr);
74945fbe478SRichard Tran Mills 
75045fbe478SRichard Tran Mills   PetscFunctionReturn(0);
75145fbe478SRichard Tran Mills }
75245fbe478SRichard Tran Mills #endif /* PETSC_HAVE_MKL_SPARSE_OPTIMIZE */
75345fbe478SRichard Tran Mills 
754e8be1fc7SRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_SP2M
755e8be1fc7SRichard Tran Mills PetscErrorCode MatMatMultNumeric_SeqAIJMKL_SeqAIJMKL_SpMV2(Mat A,Mat B,Mat C)
756e8be1fc7SRichard Tran Mills {
757e8be1fc7SRichard Tran Mills   Mat_SeqAIJMKL       *a, *b, *c;
758e8be1fc7SRichard Tran Mills   sparse_matrix_t     csrA, csrB, csrC;
759e8be1fc7SRichard Tran Mills   PetscErrorCode      ierr;
760e8be1fc7SRichard Tran Mills   sparse_status_t     stat = SPARSE_STATUS_SUCCESS;
761e8be1fc7SRichard Tran Mills   struct matrix_descr descr_type_gen;
762e8be1fc7SRichard Tran Mills   PetscObjectState    state;
763e8be1fc7SRichard Tran Mills 
764e8be1fc7SRichard Tran Mills   PetscFunctionBegin;
765e8be1fc7SRichard Tran Mills   a = (Mat_SeqAIJMKL*)A->spptr;
766e8be1fc7SRichard Tran Mills   b = (Mat_SeqAIJMKL*)B->spptr;
767e8be1fc7SRichard Tran Mills   c = (Mat_SeqAIJMKL*)C->spptr;
768e8be1fc7SRichard Tran Mills   ierr = PetscObjectStateGet((PetscObject)A,&state);CHKERRQ(ierr);
769e8be1fc7SRichard Tran Mills   if (!a->sparse_optimized || a->state != state) {
770e8be1fc7SRichard Tran Mills     MatSeqAIJMKL_create_mkl_handle(A);
771e8be1fc7SRichard Tran Mills   }
772e8be1fc7SRichard Tran Mills   ierr = PetscObjectStateGet((PetscObject)B,&state);CHKERRQ(ierr);
773e8be1fc7SRichard Tran Mills   if (!b->sparse_optimized || b->state != state) {
774e8be1fc7SRichard Tran Mills     MatSeqAIJMKL_create_mkl_handle(B);
775e8be1fc7SRichard Tran Mills   }
776e8be1fc7SRichard Tran Mills   csrA = a->csrA;
777e8be1fc7SRichard Tran Mills   csrB = b->csrA;
778e8be1fc7SRichard Tran Mills   csrC = c->csrA;
779e8be1fc7SRichard Tran Mills   descr_type_gen.type = SPARSE_MATRIX_TYPE_GENERAL;
780e8be1fc7SRichard Tran Mills 
781e8be1fc7SRichard Tran Mills   stat = mkl_sparse_sp2m(SPARSE_OPERATION_NON_TRANSPOSE,descr_type_gen,csrA,
782e8be1fc7SRichard Tran Mills                          SPARSE_OPERATION_NON_TRANSPOSE,descr_type_gen,csrB,
783e8be1fc7SRichard Tran Mills                          SPARSE_STAGE_FINALIZE_MULT,&csrC);
784e8be1fc7SRichard Tran Mills 
785e8be1fc7SRichard Tran Mills   if (stat != SPARSE_STATUS_SUCCESS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Intel MKL error: unable to complete numerical stage of sparse matrix-matrix multiply");
786e8be1fc7SRichard Tran Mills 
787e8be1fc7SRichard Tran Mills   /* Have to update the PETSc AIJ representation for matrix C from contents of MKL handle. */
7884f53af40SRichard Tran Mills   ierr = MatSeqAIJMKL_update_from_mkl_handle(C);CHKERRQ(ierr);
789e8be1fc7SRichard Tran Mills 
790e8be1fc7SRichard Tran Mills   PetscFunctionReturn(0);
791e8be1fc7SRichard Tran Mills }
792e8be1fc7SRichard Tran Mills #endif /* PETSC_HAVE_MKL_SPARSE_SP2M */
793e8be1fc7SRichard Tran Mills 
794372ec6bbSRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_OPTIMIZE
795372ec6bbSRichard Tran Mills PetscErrorCode MatTransposeMatMult_SeqAIJMKL_SeqAIJMKL_SpMV2(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat*C)
796372ec6bbSRichard Tran Mills {
797372ec6bbSRichard Tran Mills   Mat_SeqAIJMKL    *a, *b;
798372ec6bbSRichard Tran Mills   sparse_matrix_t  csrA, csrB, csrC;
799372ec6bbSRichard Tran Mills   PetscErrorCode   ierr;
800372ec6bbSRichard Tran Mills   sparse_status_t  stat = SPARSE_STATUS_SUCCESS;
801551aa5c8SRichard Tran Mills   PetscObjectState state;
802372ec6bbSRichard Tran Mills 
803372ec6bbSRichard Tran Mills   PetscFunctionBegin;
804372ec6bbSRichard Tran Mills   a = (Mat_SeqAIJMKL*)A->spptr;
805372ec6bbSRichard Tran Mills   b = (Mat_SeqAIJMKL*)B->spptr;
806551aa5c8SRichard Tran Mills   ierr = PetscObjectStateGet((PetscObject)A,&state);CHKERRQ(ierr);
807551aa5c8SRichard Tran Mills   if (!a->sparse_optimized || a->state != state) {
808372ec6bbSRichard Tran Mills     MatSeqAIJMKL_create_mkl_handle(A);
809372ec6bbSRichard Tran Mills   }
810551aa5c8SRichard Tran Mills   ierr = PetscObjectStateGet((PetscObject)B,&state);CHKERRQ(ierr);
811551aa5c8SRichard Tran Mills   if (!b->sparse_optimized || b->state != state) {
812372ec6bbSRichard Tran Mills     MatSeqAIJMKL_create_mkl_handle(B);
813372ec6bbSRichard Tran Mills   }
814372ec6bbSRichard Tran Mills   csrA = a->csrA;
815372ec6bbSRichard Tran Mills   csrB = b->csrA;
816372ec6bbSRichard Tran Mills 
817372ec6bbSRichard Tran Mills   stat = mkl_sparse_spmm(SPARSE_OPERATION_TRANSPOSE,csrA,csrB,&csrC);
8189c46acdfSRichard Tran Mills   if (stat != SPARSE_STATUS_SUCCESS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Intel MKL error: unable to complete sparse matrix-matrix multiply");
819372ec6bbSRichard Tran Mills 
820372ec6bbSRichard Tran Mills   ierr = MatSeqAIJMKL_create_from_mkl_handle(PETSC_COMM_SELF,csrC,scall,C);CHKERRQ(ierr);
821372ec6bbSRichard Tran Mills 
822372ec6bbSRichard Tran Mills   PetscFunctionReturn(0);
823372ec6bbSRichard Tran Mills }
824372ec6bbSRichard Tran Mills #endif /* PETSC_HAVE_MKL_SPARSE_OPTIMIZE */
825372ec6bbSRichard Tran Mills 
8264f53af40SRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_SP2M
8274f53af40SRichard Tran Mills PetscErrorCode MatPtAPNumeric_SeqAIJMKL_SeqAIJMKL_SpMV2(Mat A,Mat P,Mat C)
8284f53af40SRichard Tran Mills {
8294f53af40SRichard Tran Mills   Mat_SeqAIJMKL       *a, *p, *c;
8304f53af40SRichard Tran Mills   sparse_matrix_t     csrA, csrP, csrC;
8314f53af40SRichard Tran Mills   PetscBool           set, flag;
8324f53af40SRichard Tran Mills   sparse_status_t     stat = SPARSE_STATUS_SUCCESS;
8334f53af40SRichard Tran Mills   struct matrix_descr descr_type_gen;
8344f53af40SRichard Tran Mills   PetscObjectState    state;
8354f53af40SRichard Tran Mills   PetscErrorCode      ierr;
8364f53af40SRichard Tran Mills 
8374f53af40SRichard Tran Mills   PetscFunctionBegin;
8384f53af40SRichard Tran Mills   ierr = MatIsSymmetricKnown(A,&set,&flag);
8394f53af40SRichard Tran Mills   if (!set || (set && !flag)) {
8404f53af40SRichard Tran Mills     ierr = MatPtAPNumeric_SeqAIJ_SeqAIJ(A,P,C);CHKERRQ(ierr);
8414f53af40SRichard Tran Mills     PetscFunctionReturn(0);
8424f53af40SRichard Tran Mills   }
8434f53af40SRichard Tran Mills 
8444f53af40SRichard Tran Mills   a = (Mat_SeqAIJMKL*)A->spptr;
8454f53af40SRichard Tran Mills   p = (Mat_SeqAIJMKL*)P->spptr;
8464f53af40SRichard Tran Mills   c = (Mat_SeqAIJMKL*)C->spptr;
8474f53af40SRichard Tran Mills   ierr = PetscObjectStateGet((PetscObject)A,&state);CHKERRQ(ierr);
8484f53af40SRichard Tran Mills   if (!a->sparse_optimized || a->state != state) {
8494f53af40SRichard Tran Mills     MatSeqAIJMKL_create_mkl_handle(A);
8504f53af40SRichard Tran Mills   }
8514f53af40SRichard Tran Mills   ierr = PetscObjectStateGet((PetscObject)P,&state);CHKERRQ(ierr);
8524f53af40SRichard Tran Mills   if (!p->sparse_optimized || p->state != state) {
8534f53af40SRichard Tran Mills     MatSeqAIJMKL_create_mkl_handle(P);
8544f53af40SRichard Tran Mills   }
8554f53af40SRichard Tran Mills   csrA = a->csrA;
8564f53af40SRichard Tran Mills   csrP = p->csrA;
8574f53af40SRichard Tran Mills   csrC = c->csrA;
8584f53af40SRichard Tran Mills   descr_type_gen.type = SPARSE_MATRIX_TYPE_GENERAL;
8594f53af40SRichard Tran Mills 
860*f8990b4aSRichard Tran Mills   /* Note that the call below won't work for complex matrices. (We protect this when pointers are assigned in MatConvert.) */
8614f53af40SRichard Tran Mills   stat = mkl_sparse_sypr(SPARSE_OPERATION_TRANSPOSE,csrP,csrA,descr_type_gen,&csrC,SPARSE_STAGE_FINALIZE_MULT);
8624f53af40SRichard Tran Mills   if (stat != SPARSE_STATUS_SUCCESS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Intel MKL error: unable to finalize mkl_sparse_sypr");
8634f53af40SRichard Tran Mills 
8644f53af40SRichard Tran Mills   /* Have to update the PETSc AIJ representation for matrix C from contents of MKL handle. */
8654f53af40SRichard Tran Mills   ierr = MatSeqAIJMKL_update_from_mkl_handle(C);CHKERRQ(ierr);
8664f53af40SRichard Tran Mills 
8674f53af40SRichard Tran Mills   PetscFunctionReturn(0);
8684f53af40SRichard Tran Mills }
8694f53af40SRichard Tran Mills #endif
8704f53af40SRichard Tran Mills 
8714f53af40SRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_SP2M
8724f53af40SRichard Tran Mills PetscErrorCode MatPtAP_SeqAIJMKL_SeqAIJMKL_SpMV2(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
8734f53af40SRichard Tran Mills {
8744f53af40SRichard Tran Mills   Mat_SeqAIJMKL       *a, *p;
8754f53af40SRichard Tran Mills   sparse_matrix_t     csrA, csrP, csrC;
8764f53af40SRichard Tran Mills   PetscBool           set, flag;
8774f53af40SRichard Tran Mills   sparse_status_t     stat = SPARSE_STATUS_SUCCESS;
8784f53af40SRichard Tran Mills   struct matrix_descr descr_type_gen;
8794f53af40SRichard Tran Mills   PetscObjectState    state;
8804f53af40SRichard Tran Mills   PetscErrorCode      ierr;
8814f53af40SRichard Tran Mills 
8824f53af40SRichard Tran Mills   PetscFunctionBegin;
8834f53af40SRichard Tran Mills   ierr = MatIsSymmetricKnown(A,&set,&flag);
8844f53af40SRichard Tran Mills   if (!set || (set && !flag)) {
8854f53af40SRichard Tran Mills     ierr = MatPtAP_SeqAIJ_SeqAIJ(A,P,scall,fill,C);CHKERRQ(ierr);
8864f53af40SRichard Tran Mills     PetscFunctionReturn(0);
8874f53af40SRichard Tran Mills   }
8884f53af40SRichard Tran Mills 
8894f53af40SRichard Tran Mills   if (scall == MAT_REUSE_MATRIX) {
8904f53af40SRichard Tran Mills     ierr = MatPtAPNumeric_SeqAIJMKL_SeqAIJMKL_SpMV2(A,P,*C);CHKERRQ(ierr);
8914f53af40SRichard Tran Mills     PetscFunctionReturn(0);
8924f53af40SRichard Tran Mills   }
8934f53af40SRichard Tran Mills 
8944f53af40SRichard Tran Mills   a = (Mat_SeqAIJMKL*)A->spptr;
8954f53af40SRichard Tran Mills   p = (Mat_SeqAIJMKL*)P->spptr;
8964f53af40SRichard Tran Mills   ierr = PetscObjectStateGet((PetscObject)A,&state);CHKERRQ(ierr);
8974f53af40SRichard Tran Mills   if (!a->sparse_optimized || a->state != state) {
8984f53af40SRichard Tran Mills     MatSeqAIJMKL_create_mkl_handle(A);
8994f53af40SRichard Tran Mills   }
9004f53af40SRichard Tran Mills   ierr = PetscObjectStateGet((PetscObject)P,&state);CHKERRQ(ierr);
9014f53af40SRichard Tran Mills   if (!p->sparse_optimized || p->state != state) {
9024f53af40SRichard Tran Mills     MatSeqAIJMKL_create_mkl_handle(P);
9034f53af40SRichard Tran Mills   }
9044f53af40SRichard Tran Mills   csrA = a->csrA;
9054f53af40SRichard Tran Mills   csrP = p->csrA;
9064f53af40SRichard Tran Mills   descr_type_gen.type = SPARSE_MATRIX_TYPE_GENERAL;
9074f53af40SRichard Tran Mills 
908*f8990b4aSRichard Tran Mills   /* Note that the call below won't work for complex matrices. (We protect this when pointers are assigned in MatConvert.) */
9094f53af40SRichard Tran Mills   stat = mkl_sparse_sypr(SPARSE_OPERATION_TRANSPOSE,csrP,csrA,descr_type_gen,&csrC,SPARSE_STAGE_FULL_MULT);
9104f53af40SRichard Tran Mills   if (stat != SPARSE_STATUS_SUCCESS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Intel MKL error: unable to complete full mkl_sparse_sypr");
9114f53af40SRichard Tran Mills 
9124f53af40SRichard Tran Mills   ierr = MatSeqAIJMKL_create_from_mkl_handle(PETSC_COMM_SELF,csrC,scall,C);CHKERRQ(ierr);
9134f53af40SRichard Tran Mills   ierr = MatSetOption(*C,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr);
9144f53af40SRichard Tran Mills 
9154f53af40SRichard Tran Mills   PetscFunctionReturn(0);
9164f53af40SRichard Tran Mills }
9174f53af40SRichard Tran Mills #endif
9184f53af40SRichard Tran Mills 
9194a2a386eSRichard Tran Mills /* MatConvert_SeqAIJ_SeqAIJMKL converts a SeqAIJ matrix into a
9204a2a386eSRichard Tran Mills  * SeqAIJMKL matrix.  This routine is called by the MatCreate_SeqMKLAIJ()
9214a2a386eSRichard Tran Mills  * routine, but can also be used to convert an assembled SeqAIJ matrix
9224a2a386eSRichard Tran Mills  * into a SeqAIJMKL one. */
9234a2a386eSRichard Tran Mills PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJMKL(Mat A,MatType type,MatReuse reuse,Mat *newmat)
9244a2a386eSRichard Tran Mills {
9254a2a386eSRichard Tran Mills   PetscErrorCode ierr;
9264a2a386eSRichard Tran Mills   Mat            B = *newmat;
9274a2a386eSRichard Tran Mills   Mat_SeqAIJMKL  *aijmkl;
928c9d46305SRichard Tran Mills   PetscBool      set;
929e9c94282SRichard Tran Mills   PetscBool      sametype;
9304a2a386eSRichard Tran Mills 
9314a2a386eSRichard Tran Mills   PetscFunctionBegin;
9324a2a386eSRichard Tran Mills   if (reuse == MAT_INITIAL_MATRIX) {
9334a2a386eSRichard Tran Mills     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
9344a2a386eSRichard Tran Mills   }
9354a2a386eSRichard Tran Mills 
936e9c94282SRichard Tran Mills   ierr = PetscObjectTypeCompare((PetscObject)A,type,&sametype);CHKERRQ(ierr);
937e9c94282SRichard Tran Mills   if (sametype) PetscFunctionReturn(0);
938e9c94282SRichard Tran Mills 
9394a2a386eSRichard Tran Mills   ierr     = PetscNewLog(B,&aijmkl);CHKERRQ(ierr);
9404a2a386eSRichard Tran Mills   B->spptr = (void*) aijmkl;
9414a2a386eSRichard Tran Mills 
942df555b71SRichard Tran Mills   /* Set function pointers for methods that we inherit from AIJ but override.
943969800c5SRichard Tran Mills    * We also parse some command line options below, since those determine some of the methods we point to. */
9444a2a386eSRichard Tran Mills   B->ops->duplicate        = MatDuplicate_SeqAIJMKL;
9454a2a386eSRichard Tran Mills   B->ops->assemblyend      = MatAssemblyEnd_SeqAIJMKL;
9464a2a386eSRichard Tran Mills   B->ops->destroy          = MatDestroy_SeqAIJMKL;
947c9d46305SRichard Tran Mills 
9484abfa3b3SRichard Tran Mills   aijmkl->sparse_optimized = PETSC_FALSE;
949d995685eSRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_OPTIMIZE
950d995685eSRichard Tran Mills   aijmkl->no_SpMV2 = PETSC_FALSE;  /* Default to using the SpMV2 routines if our MKL supports them. */
951a8327b06SKarl Rupp #else
952d995685eSRichard Tran Mills   aijmkl->no_SpMV2 = PETSC_TRUE;
953d995685eSRichard Tran Mills #endif
9545b49642aSRichard Tran Mills   aijmkl->eager_inspection = PETSC_FALSE;
9554abfa3b3SRichard Tran Mills 
9564abfa3b3SRichard Tran Mills   /* Parse command line options. */
957c9d46305SRichard Tran Mills   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"AIJMKL Options","Mat");CHKERRQ(ierr);
958c9d46305SRichard Tran Mills   ierr = PetscOptionsBool("-mat_aijmkl_no_spmv2","NoSPMV2","None",(PetscBool)aijmkl->no_SpMV2,(PetscBool*)&aijmkl->no_SpMV2,&set);CHKERRQ(ierr);
9595b49642aSRichard Tran Mills   ierr = PetscOptionsBool("-mat_aijmkl_eager_inspection","Eager Inspection","None",(PetscBool)aijmkl->eager_inspection,(PetscBool*)&aijmkl->eager_inspection,&set);CHKERRQ(ierr);
960c9d46305SRichard Tran Mills   ierr = PetscOptionsEnd();CHKERRQ(ierr);
961d995685eSRichard Tran Mills #ifndef PETSC_HAVE_MKL_SPARSE_OPTIMIZE
962d995685eSRichard Tran Mills   if(!aijmkl->no_SpMV2) {
963d995685eSRichard Tran Mills     ierr = PetscInfo(B,"User requested use of MKL SpMV2 routines, but MKL version does not support mkl_sparse_optimize();  defaulting to non-SpMV2 routines.\n");
964d995685eSRichard Tran Mills     aijmkl->no_SpMV2 = PETSC_TRUE;
965d995685eSRichard Tran Mills   }
966d995685eSRichard Tran Mills #endif
967c9d46305SRichard Tran Mills 
968c9d46305SRichard Tran Mills   if(!aijmkl->no_SpMV2) {
969d995685eSRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_OPTIMIZE
970df555b71SRichard Tran Mills     B->ops->mult             = MatMult_SeqAIJMKL_SpMV2;
971969800c5SRichard Tran Mills     B->ops->multtranspose    = MatMultTranspose_SeqAIJMKL_SpMV2;
972df555b71SRichard Tran Mills     B->ops->multadd          = MatMultAdd_SeqAIJMKL_SpMV2;
973969800c5SRichard Tran Mills     B->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJMKL_SpMV2;
97445fbe478SRichard Tran Mills     B->ops->matmult          = MatMatMult_SeqAIJMKL_SeqAIJMKL_SpMV2;
975e8be1fc7SRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_SP2M
976e8be1fc7SRichard Tran Mills     B->ops->matmultnumeric   = MatMatMultNumeric_SeqAIJMKL_SeqAIJMKL_SpMV2;
9774f53af40SRichard Tran Mills #ifndef PETSC_USE_COMPLEX
9784f53af40SRichard Tran Mills     B->ops->ptap             = MatPtAP_SeqAIJMKL_SeqAIJMKL_SpMV2;
9794f53af40SRichard Tran Mills     B->ops->ptapnumeric      = MatPtAPNumeric_SeqAIJMKL_SeqAIJMKL_SpMV2;
9804f53af40SRichard Tran Mills #endif
981e8be1fc7SRichard Tran Mills #endif
982a557fde5SRichard Tran Mills     B->ops->transposematmult = MatTransposeMatMult_SeqAIJMKL_SeqAIJMKL_SpMV2;
983d995685eSRichard Tran Mills #endif
984c9d46305SRichard Tran Mills   } else {
9854a2a386eSRichard Tran Mills     B->ops->mult             = MatMult_SeqAIJMKL;
986969800c5SRichard Tran Mills     B->ops->multtranspose    = MatMultTranspose_SeqAIJMKL;
9874a2a386eSRichard Tran Mills     B->ops->multadd          = MatMultAdd_SeqAIJMKL;
988969800c5SRichard Tran Mills     B->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJMKL;
989c9d46305SRichard Tran Mills   }
9904a2a386eSRichard Tran Mills 
9914a2a386eSRichard Tran Mills   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaijmkl_seqaij_C",MatConvert_SeqAIJMKL_SeqAIJ);CHKERRQ(ierr);
992e9c94282SRichard Tran Mills   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMult_seqdense_seqaijmkl_C",MatMatMult_SeqDense_SeqAIJ);CHKERRQ(ierr);
993e9c94282SRichard Tran Mills   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMultSymbolic_seqdense_seqaijmkl_C",MatMatMultSymbolic_SeqDense_SeqAIJ);CHKERRQ(ierr);
994e9c94282SRichard Tran Mills   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMultNumeric_seqdense_seqaijmkl_C",MatMatMultNumeric_SeqDense_SeqAIJ);CHKERRQ(ierr);
99545fbe478SRichard Tran Mills   if(!aijmkl->no_SpMV2) {
99645fbe478SRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_OPTIMIZE
99745fbe478SRichard Tran Mills     ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMult_seqaijmkl_seqaijmkl_C",MatMatMult_SeqAIJMKL_SeqAIJMKL_SpMV2);CHKERRQ(ierr);
998e8be1fc7SRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_SP2M
999e8be1fc7SRichard Tran Mills     ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMultNumeric_seqaijmkl_seqaijmkl_C",MatMatMultNumeric_SeqAIJMKL_SeqAIJMKL_SpMV2);CHKERRQ(ierr);
1000e8be1fc7SRichard Tran Mills #endif
1001372ec6bbSRichard Tran Mills     ierr = PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMult_seqaijmkl_seqaijmkl_C",MatTransposeMatMult_SeqAIJMKL_SeqAIJMKL_SpMV2);CHKERRQ(ierr);
100245fbe478SRichard Tran Mills #endif
100345fbe478SRichard Tran Mills   }
10044a2a386eSRichard Tran Mills 
10054a2a386eSRichard Tran Mills   ierr    = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJMKL);CHKERRQ(ierr);
10064a2a386eSRichard Tran Mills   *newmat = B;
10074a2a386eSRichard Tran Mills   PetscFunctionReturn(0);
10084a2a386eSRichard Tran Mills }
10094a2a386eSRichard Tran Mills 
10104a2a386eSRichard Tran Mills /*@C
10114a2a386eSRichard Tran Mills    MatCreateSeqAIJMKL - Creates a sparse matrix of type SEQAIJMKL.
10124a2a386eSRichard Tran Mills    This type inherits from AIJ and is largely identical, but uses sparse BLAS
10134a2a386eSRichard Tran Mills    routines from Intel MKL whenever possible.
10143af10221SRichard Tran Mills    MatMult, MatMultAdd, MatMultTranspose, MatMultTransposeAdd, MatMatMult, and MatTransposeMatMult
101590147e49SRichard Tran Mills    operations are currently supported.
101690147e49SRichard Tran Mills    If the installed version of MKL supports the "SpMV2" sparse
101790147e49SRichard Tran Mills    inspector-executor routines, then those are used by default.
101890147e49SRichard Tran Mills 
10194a2a386eSRichard Tran Mills    Collective on MPI_Comm
10204a2a386eSRichard Tran Mills 
10214a2a386eSRichard Tran Mills    Input Parameters:
10224a2a386eSRichard Tran Mills +  comm - MPI communicator, set to PETSC_COMM_SELF
10234a2a386eSRichard Tran Mills .  m - number of rows
10244a2a386eSRichard Tran Mills .  n - number of columns
10254a2a386eSRichard Tran Mills .  nz - number of nonzeros per row (same for all rows)
10264a2a386eSRichard Tran Mills -  nnz - array containing the number of nonzeros in the various rows
10274a2a386eSRichard Tran Mills          (possibly different for each row) or NULL
10284a2a386eSRichard Tran Mills 
10294a2a386eSRichard Tran Mills    Output Parameter:
10304a2a386eSRichard Tran Mills .  A - the matrix
10314a2a386eSRichard Tran Mills 
103290147e49SRichard Tran Mills    Options Database Keys:
103366b7eeb6SRichard Tran Mills +  -mat_aijmkl_no_spmv2 - disable use of the SpMV2 inspector-executor routines
103466b7eeb6SRichard Tran Mills -  -mat_aijmkl_eager_inspection - perform MKL "inspection" phase upon matrix assembly; default is to do "lazy" inspection, performing this step the first time the matrix is applied
103590147e49SRichard Tran Mills 
10364a2a386eSRichard Tran Mills    Notes:
10374a2a386eSRichard Tran Mills    If nnz is given then nz is ignored
10384a2a386eSRichard Tran Mills 
10394a2a386eSRichard Tran Mills    Level: intermediate
10404a2a386eSRichard Tran Mills 
104190147e49SRichard Tran Mills .keywords: matrix, MKL, sparse, parallel
10424a2a386eSRichard Tran Mills 
10434a2a386eSRichard Tran Mills .seealso: MatCreate(), MatCreateMPIAIJMKL(), MatSetValues()
10444a2a386eSRichard Tran Mills @*/
10454a2a386eSRichard Tran Mills PetscErrorCode  MatCreateSeqAIJMKL(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
10464a2a386eSRichard Tran Mills {
10474a2a386eSRichard Tran Mills   PetscErrorCode ierr;
10484a2a386eSRichard Tran Mills 
10494a2a386eSRichard Tran Mills   PetscFunctionBegin;
10504a2a386eSRichard Tran Mills   ierr = MatCreate(comm,A);CHKERRQ(ierr);
10514a2a386eSRichard Tran Mills   ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr);
10524a2a386eSRichard Tran Mills   ierr = MatSetType(*A,MATSEQAIJMKL);CHKERRQ(ierr);
10534a2a386eSRichard Tran Mills   ierr = MatSeqAIJSetPreallocation_SeqAIJ(*A,nz,nnz);CHKERRQ(ierr);
10544a2a386eSRichard Tran Mills   PetscFunctionReturn(0);
10554a2a386eSRichard Tran Mills }
10564a2a386eSRichard Tran Mills 
10574a2a386eSRichard Tran Mills PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJMKL(Mat A)
10584a2a386eSRichard Tran Mills {
10594a2a386eSRichard Tran Mills   PetscErrorCode ierr;
10604a2a386eSRichard Tran Mills 
10614a2a386eSRichard Tran Mills   PetscFunctionBegin;
10624a2a386eSRichard Tran Mills   ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
10634a2a386eSRichard Tran Mills   ierr = MatConvert_SeqAIJ_SeqAIJMKL(A,MATSEQAIJMKL,MAT_INPLACE_MATRIX,&A);CHKERRQ(ierr);
10644a2a386eSRichard Tran Mills   PetscFunctionReturn(0);
10654a2a386eSRichard Tran Mills }
1066