xref: /petsc/src/mat/impls/aij/seq/aijmkl/aijmkl.c (revision 6e369cd56c19f466a165a36e147eb4df4120aa3d)
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. */
174abfa3b3SRichard Tran Mills   PetscBool sparse_optimized; /* If PETSC_TRUE, then mkl_sparse_optimize() has been called. */
18b8cbc1fbSRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_OPTIMIZE
19df555b71SRichard Tran Mills   sparse_matrix_t csrA; /* "Handle" used by SpMV2 inspector-executor routines. */
20df555b71SRichard Tran Mills   struct matrix_descr descr;
21b8cbc1fbSRichard Tran Mills #endif
224a2a386eSRichard Tran Mills } Mat_SeqAIJMKL;
234a2a386eSRichard Tran Mills 
244a2a386eSRichard Tran Mills extern PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat,MatAssemblyType);
254a2a386eSRichard Tran Mills 
264a2a386eSRichard Tran Mills PETSC_INTERN PetscErrorCode MatConvert_SeqAIJMKL_SeqAIJ(Mat A,MatType type,MatReuse reuse,Mat *newmat)
274a2a386eSRichard Tran Mills {
284a2a386eSRichard Tran Mills   /* This routine is only called to convert a MATAIJMKL to its base PETSc type, */
294a2a386eSRichard Tran Mills   /* so we will ignore 'MatType type'. */
304a2a386eSRichard Tran Mills   PetscErrorCode ierr;
314a2a386eSRichard Tran Mills   Mat            B       = *newmat;
324a2a386eSRichard Tran Mills   Mat_SeqAIJMKL  *aijmkl=(Mat_SeqAIJMKL*)A->spptr;
334a2a386eSRichard Tran Mills 
344a2a386eSRichard Tran Mills   PetscFunctionBegin;
354a2a386eSRichard Tran Mills   if (reuse == MAT_INITIAL_MATRIX) {
364a2a386eSRichard Tran Mills     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
37e9c94282SRichard Tran Mills     aijmkl = (Mat_SeqAIJMKL*)B->spptr;
384a2a386eSRichard Tran Mills   }
394a2a386eSRichard Tran Mills 
404a2a386eSRichard Tran Mills   /* Reset the original function pointers. */
4154871a98SRichard Tran Mills   B->ops->duplicate        = MatDuplicate_SeqAIJ;
424a2a386eSRichard Tran Mills   B->ops->assemblyend      = MatAssemblyEnd_SeqAIJ;
434a2a386eSRichard Tran Mills   B->ops->destroy          = MatDestroy_SeqAIJ;
4454871a98SRichard Tran Mills   B->ops->mult             = MatMult_SeqAIJ;
45ff03dc53SRichard Tran Mills   B->ops->multtranspose    = MatMultTranspose_SeqAIJ;
4654871a98SRichard Tran Mills   B->ops->multadd          = MatMultAdd_SeqAIJ;
47ff03dc53SRichard Tran Mills   B->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJ;
484a2a386eSRichard Tran Mills 
49e9c94282SRichard Tran Mills   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaijmkl_seqaij_C",NULL);CHKERRQ(ierr);
50e9c94282SRichard Tran Mills   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMult_seqdense_seqaijmkl_C",NULL);CHKERRQ(ierr);
51e9c94282SRichard Tran Mills   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMultSymbolic_seqdense_seqaijmkl_C",NULL);CHKERRQ(ierr);
52e9c94282SRichard Tran Mills   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMultNumeric_seqdense_seqaijmkl_C",NULL);CHKERRQ(ierr);
53e9c94282SRichard Tran Mills 
544abfa3b3SRichard Tran Mills   /* Free everything in the Mat_SeqAIJMKL data structure. Currently, this
55e9c94282SRichard Tran Mills    * simply involves destroying the MKL sparse matrix handle and then freeing
56e9c94282SRichard Tran Mills    * the spptr pointer. */
574abfa3b3SRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_OPTIMIZE
584abfa3b3SRichard Tran Mills   if (aijmkl->sparse_optimized) {
590632b357SRichard Tran Mills     sparse_status_t stat;
604abfa3b3SRichard Tran Mills     stat = mkl_sparse_destroy(aijmkl->csrA);
614abfa3b3SRichard Tran Mills     if (stat != SPARSE_STATUS_SUCCESS) {
624abfa3b3SRichard Tran Mills       PetscFunctionReturn(PETSC_ERR_LIB);
634abfa3b3SRichard Tran Mills     }
644abfa3b3SRichard Tran Mills   }
654abfa3b3SRichard Tran Mills #endif /* PETSC_HAVE_MKL_SPARSE_OPTIMIZE */
66e9c94282SRichard Tran Mills   ierr = PetscFree(B->spptr);CHKERRQ(ierr);
674a2a386eSRichard Tran Mills 
684a2a386eSRichard Tran Mills   /* Change the type of B to MATSEQAIJ. */
694a2a386eSRichard Tran Mills   ierr = PetscObjectChangeTypeName((PetscObject)B, MATSEQAIJ);CHKERRQ(ierr);
704a2a386eSRichard Tran Mills 
714a2a386eSRichard Tran Mills   *newmat = B;
724a2a386eSRichard Tran Mills   PetscFunctionReturn(0);
734a2a386eSRichard Tran Mills }
744a2a386eSRichard Tran Mills 
754a2a386eSRichard Tran Mills PetscErrorCode MatDestroy_SeqAIJMKL(Mat A)
764a2a386eSRichard Tran Mills {
774a2a386eSRichard Tran Mills   PetscErrorCode ierr;
784a2a386eSRichard Tran Mills   Mat_SeqAIJMKL *aijmkl = (Mat_SeqAIJMKL*) A->spptr;
794a2a386eSRichard Tran Mills 
804a2a386eSRichard Tran Mills   PetscFunctionBegin;
81e9c94282SRichard Tran Mills 
82e9c94282SRichard Tran Mills   /* If MatHeaderMerge() was used, then this SeqAIJMKL matrix will not have an
83e9c94282SRichard Tran Mills    * spptr pointer. */
84e9c94282SRichard Tran Mills   if (aijmkl) {
854a2a386eSRichard Tran Mills     /* Clean up everything in the Mat_SeqAIJMKL data structure, then free A->spptr. */
864abfa3b3SRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_OPTIMIZE
874abfa3b3SRichard Tran Mills     if (aijmkl->sparse_optimized) {
884abfa3b3SRichard Tran Mills       sparse_status_t stat = SPARSE_STATUS_SUCCESS;
894abfa3b3SRichard Tran Mills       stat = mkl_sparse_destroy(aijmkl->csrA);
904abfa3b3SRichard Tran Mills       if (stat != SPARSE_STATUS_SUCCESS) {
914abfa3b3SRichard Tran Mills         PetscFunctionReturn(PETSC_ERR_LIB);
924abfa3b3SRichard Tran Mills       }
934abfa3b3SRichard Tran Mills     }
944abfa3b3SRichard Tran Mills #endif /* PETSC_HAVE_MKL_SPARSE_OPTIMIZE */
954a2a386eSRichard Tran Mills     ierr = PetscFree(A->spptr);CHKERRQ(ierr);
96e9c94282SRichard Tran Mills   }
974a2a386eSRichard Tran Mills 
984a2a386eSRichard Tran Mills   /* Change the type of A back to SEQAIJ and use MatDestroy_SeqAIJ()
994a2a386eSRichard Tran Mills    * to destroy everything that remains. */
1004a2a386eSRichard Tran Mills   ierr = PetscObjectChangeTypeName((PetscObject)A, MATSEQAIJ);CHKERRQ(ierr);
1014a2a386eSRichard Tran Mills   /* Note that I don't call MatSetType().  I believe this is because that
1024a2a386eSRichard Tran Mills    * is only to be called when *building* a matrix.  I could be wrong, but
1034a2a386eSRichard Tran Mills    * that is how things work for the SuperLU matrix class. */
1044a2a386eSRichard Tran Mills   ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
1054a2a386eSRichard Tran Mills   PetscFunctionReturn(0);
1064a2a386eSRichard Tran Mills }
1074a2a386eSRichard Tran Mills 
108*6e369cd5SRichard Tran Mills PETSC_INTERN PetscErrorCode MatSeqAIJMKL_create_mkl_handle(Mat A)
1094a2a386eSRichard Tran Mills {
1104a2a386eSRichard Tran Mills   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
111df555b71SRichard Tran Mills   Mat_SeqAIJMKL   *aijmkl;
11258678438SRichard Tran Mills   PetscInt        m,n;
113*6e369cd5SRichard Tran Mills   MatScalar       *aa;
114df555b71SRichard Tran Mills   PetscInt        *aj,*ai;
1154a2a386eSRichard Tran Mills 
116*6e369cd5SRichard Tran Mills #ifndef PETSC_HAVE_MKL_SPARSE_OPTIMIZE
117*6e369cd5SRichard Tran Mills   /* If the MKL library does not have mkl_sparse_optimize(), then this routine
118*6e369cd5SRichard Tran Mills    * does nothing. We make it callable anyway in this case because it cuts
119*6e369cd5SRichard Tran Mills    * down on littering the code with #ifdefs. */
120*6e369cd5SRichard Tran Mills   PetscFunctionReturn(0);
121*6e369cd5SRichard Tran Mills #else
1224a2a386eSRichard Tran Mills 
123*6e369cd5SRichard Tran Mills   sparse_status_t stat;
1244a2a386eSRichard Tran Mills 
125df555b71SRichard Tran Mills   aijmkl = (Mat_SeqAIJMKL*) A->spptr;
126*6e369cd5SRichard Tran Mills 
127*6e369cd5SRichard Tran Mills   if (aijmkl->no_SpMV2) PetscFunctionReturn(0);
128*6e369cd5SRichard Tran Mills 
1290632b357SRichard Tran Mills   if (aijmkl->sparse_optimized) {
1300632b357SRichard Tran Mills     /* Matrix has been previously assembled and optimized. Must destroy old
1310632b357SRichard Tran Mills      * matrix handle before running the optimization step again. */
1320632b357SRichard Tran Mills     stat = mkl_sparse_destroy(aijmkl->csrA);
1330632b357SRichard Tran Mills     if (stat != SPARSE_STATUS_SUCCESS) {
1340632b357SRichard Tran Mills       PetscFunctionReturn(PETSC_ERR_LIB);
1350632b357SRichard Tran Mills     }
1360632b357SRichard Tran Mills   }
137*6e369cd5SRichard Tran Mills 
138c9d46305SRichard Tran Mills   /* Now perform the SpMV2 setup and matrix optimization. */
139df555b71SRichard Tran Mills   aijmkl->descr.type        = SPARSE_MATRIX_TYPE_GENERAL;
140df555b71SRichard Tran Mills   aijmkl->descr.mode        = SPARSE_FILL_MODE_LOWER;
141df555b71SRichard Tran Mills   aijmkl->descr.diag        = SPARSE_DIAG_NON_UNIT;
14258678438SRichard Tran Mills   m = A->rmap->n;
14358678438SRichard Tran Mills   n = A->cmap->n;
144df555b71SRichard Tran Mills   aj   = a->j;  /* aj[k] gives column index for element aa[k]. */
145df555b71SRichard Tran Mills   aa   = a->a;  /* Nonzero elements stored row-by-row. */
146df555b71SRichard Tran Mills   ai   = a->i;  /* ai[k] is the position in aa and aj where row k starts. */
147f36dfe3fSRichard Tran Mills   if (n>0) {
14858678438SRichard Tran Mills     stat = mkl_sparse_x_create_csr (&aijmkl->csrA,SPARSE_INDEX_BASE_ZERO,m,n,ai,ai+1,aj,aa);
149df555b71SRichard Tran Mills     stat = mkl_sparse_set_mv_hint(aijmkl->csrA,SPARSE_OPERATION_NON_TRANSPOSE,aijmkl->descr,1000);
150df555b71SRichard Tran Mills     stat = mkl_sparse_set_memory_hint(aijmkl->csrA,SPARSE_MEMORY_AGGRESSIVE);
151df555b71SRichard Tran Mills     stat = mkl_sparse_optimize(aijmkl->csrA);
152df555b71SRichard Tran Mills     if (stat != SPARSE_STATUS_SUCCESS) {
153f68ad4bdSRichard Tran Mills       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Intel MKL error: unable to create matrix handle/complete mkl_sparse_optimize");
154df555b71SRichard Tran Mills       PetscFunctionReturn(PETSC_ERR_LIB);
155df555b71SRichard Tran Mills     }
1564abfa3b3SRichard Tran Mills     aijmkl->sparse_optimized = PETSC_TRUE;
157c9d46305SRichard Tran Mills   }
158*6e369cd5SRichard Tran Mills 
159*6e369cd5SRichard Tran Mills   PetscFunctionReturn(0);
160d995685eSRichard Tran Mills #endif
161*6e369cd5SRichard Tran Mills }
162*6e369cd5SRichard Tran Mills 
163*6e369cd5SRichard Tran Mills PetscErrorCode MatDuplicate_SeqAIJMKL(Mat A, MatDuplicateOption op, Mat *M)
164*6e369cd5SRichard Tran Mills {
165*6e369cd5SRichard Tran Mills   PetscErrorCode ierr;
166*6e369cd5SRichard Tran Mills   Mat_SeqAIJMKL *aijmkl;
167*6e369cd5SRichard Tran Mills   Mat_SeqAIJMKL *aijmkl_dest;
168*6e369cd5SRichard Tran Mills 
169*6e369cd5SRichard Tran Mills   PetscFunctionBegin;
170*6e369cd5SRichard Tran Mills   ierr = MatDuplicate_SeqAIJ(A,op,M);CHKERRQ(ierr);
171*6e369cd5SRichard Tran Mills   aijmkl      = (Mat_SeqAIJMKL*) A->spptr;
172*6e369cd5SRichard Tran Mills   aijmkl_dest = (Mat_SeqAIJMKL*) (*M)->spptr;
173*6e369cd5SRichard Tran Mills   ierr = PetscMemcpy(aijmkl_dest,aijmkl,sizeof(Mat_SeqAIJMKL));CHKERRQ(ierr);
174*6e369cd5SRichard Tran Mills   aijmkl_dest->sparse_optimized = PETSC_FALSE;
175*6e369cd5SRichard Tran Mills   ierr = MatSeqAIJMKL_create_mkl_handle(A);CHKERRQ(ierr);
176*6e369cd5SRichard Tran Mills   PetscFunctionReturn(0);
177*6e369cd5SRichard Tran Mills }
178*6e369cd5SRichard Tran Mills 
179*6e369cd5SRichard Tran Mills PetscErrorCode MatAssemblyEnd_SeqAIJMKL(Mat A, MatAssemblyType mode)
180*6e369cd5SRichard Tran Mills {
181*6e369cd5SRichard Tran Mills   PetscErrorCode  ierr;
182*6e369cd5SRichard Tran Mills   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
183*6e369cd5SRichard Tran Mills 
184*6e369cd5SRichard Tran Mills   PetscFunctionBegin;
185*6e369cd5SRichard Tran Mills   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
186*6e369cd5SRichard Tran Mills 
187*6e369cd5SRichard Tran Mills   /* Since a MATSEQAIJMKL matrix is really just a MATSEQAIJ with some
188*6e369cd5SRichard Tran Mills    * extra information and some different methods, call the AssemblyEnd
189*6e369cd5SRichard Tran Mills    * routine for a MATSEQAIJ.
190*6e369cd5SRichard Tran Mills    * I'm not sure if this is the best way to do this, but it avoids
191*6e369cd5SRichard Tran Mills    * a lot of code duplication.
192*6e369cd5SRichard Tran Mills    * I also note that currently MATSEQAIJMKL doesn't know anything about
193*6e369cd5SRichard Tran Mills    * the Mat_CompressedRow data structure that SeqAIJ now uses when there
194*6e369cd5SRichard Tran Mills    * are many zero rows.  If the SeqAIJ assembly end routine decides to use
195*6e369cd5SRichard Tran Mills    * this, this may break things.  (Don't know... haven't looked at it.
196*6e369cd5SRichard Tran Mills    * Do I need to disable this somehow?) */
197*6e369cd5SRichard Tran Mills   a->inode.use = PETSC_FALSE;  /* Must disable: otherwise the MKL routines won't get used. */
198*6e369cd5SRichard Tran Mills   ierr         = MatAssemblyEnd_SeqAIJ(A, mode);CHKERRQ(ierr);
199*6e369cd5SRichard Tran Mills 
200*6e369cd5SRichard Tran Mills   /* Now create the MKL sparse handle (if needed; the function checks). */
201*6e369cd5SRichard Tran Mills   ierr = MatSeqAIJMKL_create_mkl_handle(A);CHKERRQ(ierr);
202df555b71SRichard Tran Mills 
2034a2a386eSRichard Tran Mills   PetscFunctionReturn(0);
2044a2a386eSRichard Tran Mills }
2054a2a386eSRichard Tran Mills 
2064a2a386eSRichard Tran Mills PetscErrorCode MatMult_SeqAIJMKL(Mat A,Vec xx,Vec yy)
2074a2a386eSRichard Tran Mills {
2084a2a386eSRichard Tran Mills   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
2094a2a386eSRichard Tran Mills   const PetscScalar *x;
2104a2a386eSRichard Tran Mills   PetscScalar       *y;
2114a2a386eSRichard Tran Mills   const MatScalar   *aa;
2124a2a386eSRichard Tran Mills   PetscErrorCode    ierr;
2134a2a386eSRichard Tran Mills   PetscInt          m=A->rmap->n;
2144a2a386eSRichard Tran Mills   const PetscInt    *aj,*ai;
2154a2a386eSRichard Tran Mills 
2164a2a386eSRichard Tran Mills   /* Variables not in MatMult_SeqAIJ. */
217ff03dc53SRichard Tran Mills   char transa = 'n';  /* Used to indicate to MKL that we are not computing the transpose product. */
218ff03dc53SRichard Tran Mills 
219ff03dc53SRichard Tran Mills   PetscFunctionBegin;
220ff03dc53SRichard Tran Mills   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
221ff03dc53SRichard Tran Mills   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
222ff03dc53SRichard Tran Mills   aj   = a->j;  /* aj[k] gives column index for element aa[k]. */
223ff03dc53SRichard Tran Mills   aa   = a->a;  /* Nonzero elements stored row-by-row. */
224ff03dc53SRichard Tran Mills   ai   = a->i;  /* ai[k] is the position in aa and aj where row k starts. */
225ff03dc53SRichard Tran Mills 
226ff03dc53SRichard Tran Mills   /* Call MKL sparse BLAS routine to do the MatMult. */
227ff03dc53SRichard Tran Mills   mkl_cspblas_xcsrgemv(&transa,&m,aa,ai,aj,x,y);
228ff03dc53SRichard Tran Mills 
229ff03dc53SRichard Tran Mills   ierr = PetscLogFlops(2.0*a->nz - a->nonzerorowcnt);CHKERRQ(ierr);
230ff03dc53SRichard Tran Mills   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
231ff03dc53SRichard Tran Mills   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
232ff03dc53SRichard Tran Mills   PetscFunctionReturn(0);
233ff03dc53SRichard Tran Mills }
234ff03dc53SRichard Tran Mills 
235d995685eSRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_OPTIMIZE
236df555b71SRichard Tran Mills PetscErrorCode MatMult_SeqAIJMKL_SpMV2(Mat A,Vec xx,Vec yy)
237df555b71SRichard Tran Mills {
238df555b71SRichard Tran Mills   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
239df555b71SRichard Tran Mills   Mat_SeqAIJMKL     *aijmkl=(Mat_SeqAIJMKL*)A->spptr;
240df555b71SRichard Tran Mills   const PetscScalar *x;
241df555b71SRichard Tran Mills   PetscScalar       *y;
242df555b71SRichard Tran Mills   PetscErrorCode    ierr;
243df555b71SRichard Tran Mills   sparse_status_t stat = SPARSE_STATUS_SUCCESS;
244df555b71SRichard Tran Mills 
245df555b71SRichard Tran Mills   PetscFunctionBegin;
246df555b71SRichard Tran Mills 
247f36dfe3fSRichard Tran Mills   /* If there are no rows, this is a no-op: return immediately. */
248f36dfe3fSRichard Tran Mills   if(A->cmap->n < 1) PetscFunctionReturn(0);
249f36dfe3fSRichard Tran Mills 
250df555b71SRichard Tran Mills   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
251df555b71SRichard Tran Mills   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
252df555b71SRichard Tran Mills 
253df555b71SRichard Tran Mills   /* Call MKL SpMV2 executor routine to do the MatMult. */
254df555b71SRichard Tran Mills   stat = mkl_sparse_x_mv(SPARSE_OPERATION_NON_TRANSPOSE,1.0,aijmkl->csrA,aijmkl->descr,x,0.0,y);
255df555b71SRichard Tran Mills 
256df555b71SRichard Tran Mills   ierr = PetscLogFlops(2.0*a->nz - a->nonzerorowcnt);CHKERRQ(ierr);
257df555b71SRichard Tran Mills   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
258df555b71SRichard Tran Mills   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
259df555b71SRichard Tran Mills   if (stat != SPARSE_STATUS_SUCCESS) {
260df555b71SRichard Tran Mills     PetscFunctionReturn(PETSC_ERR_LIB);
261df555b71SRichard Tran Mills   }
262df555b71SRichard Tran Mills   PetscFunctionReturn(0);
263df555b71SRichard Tran Mills }
264d995685eSRichard Tran Mills #endif /* PETSC_HAVE_MKL_SPARSE_OPTIMIZE */
265df555b71SRichard Tran Mills 
266ff03dc53SRichard Tran Mills PetscErrorCode MatMultTranspose_SeqAIJMKL(Mat A,Vec xx,Vec yy)
267ff03dc53SRichard Tran Mills {
268ff03dc53SRichard Tran Mills   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
269ff03dc53SRichard Tran Mills   const PetscScalar *x;
270ff03dc53SRichard Tran Mills   PetscScalar       *y;
271ff03dc53SRichard Tran Mills   const MatScalar   *aa;
272ff03dc53SRichard Tran Mills   PetscErrorCode    ierr;
273ff03dc53SRichard Tran Mills   PetscInt          m=A->rmap->n;
274ff03dc53SRichard Tran Mills   const PetscInt    *aj,*ai;
275ff03dc53SRichard Tran Mills 
276ff03dc53SRichard Tran Mills   /* Variables not in MatMultTranspose_SeqAIJ. */
277ff03dc53SRichard Tran Mills   char transa = 't';  /* Used to indicate to MKL that we are computing the transpose product. */
2784a2a386eSRichard Tran Mills 
2794a2a386eSRichard Tran Mills   PetscFunctionBegin;
2804a2a386eSRichard Tran Mills   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
2814a2a386eSRichard Tran Mills   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
2824a2a386eSRichard Tran Mills   aj   = a->j;  /* aj[k] gives column index for element aa[k]. */
2834a2a386eSRichard Tran Mills   aa   = a->a;  /* Nonzero elements stored row-by-row. */
2844a2a386eSRichard Tran Mills   ai   = a->i;  /* ai[k] is the position in aa and aj where row k starts. */
2854a2a386eSRichard Tran Mills 
2864a2a386eSRichard Tran Mills   /* Call MKL sparse BLAS routine to do the MatMult. */
2874a2a386eSRichard Tran Mills   mkl_cspblas_xcsrgemv(&transa,&m,aa,ai,aj,x,y);
2884a2a386eSRichard Tran Mills 
2894a2a386eSRichard Tran Mills   ierr = PetscLogFlops(2.0*a->nz - a->nonzerorowcnt);CHKERRQ(ierr);
2904a2a386eSRichard Tran Mills   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
2914a2a386eSRichard Tran Mills   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
2924a2a386eSRichard Tran Mills   PetscFunctionReturn(0);
2934a2a386eSRichard Tran Mills }
2944a2a386eSRichard Tran Mills 
295d995685eSRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_OPTIMIZE
296df555b71SRichard Tran Mills PetscErrorCode MatMultTranspose_SeqAIJMKL_SpMV2(Mat A,Vec xx,Vec yy)
297df555b71SRichard Tran Mills {
298df555b71SRichard Tran Mills   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
299df555b71SRichard Tran Mills   Mat_SeqAIJMKL     *aijmkl=(Mat_SeqAIJMKL*)A->spptr;
300df555b71SRichard Tran Mills   const PetscScalar *x;
301df555b71SRichard Tran Mills   PetscScalar       *y;
302df555b71SRichard Tran Mills   PetscErrorCode    ierr;
3030632b357SRichard Tran Mills   sparse_status_t   stat;
304df555b71SRichard Tran Mills 
305df555b71SRichard Tran Mills   PetscFunctionBegin;
306df555b71SRichard Tran Mills 
307f36dfe3fSRichard Tran Mills   /* If there are no rows, this is a no-op: return immediately. */
308f36dfe3fSRichard Tran Mills   if(A->cmap->n < 1) PetscFunctionReturn(0);
309f36dfe3fSRichard Tran Mills 
310df555b71SRichard Tran Mills   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
311df555b71SRichard Tran Mills   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
312df555b71SRichard Tran Mills 
313df555b71SRichard Tran Mills   /* Call MKL SpMV2 executor routine to do the MatMultTranspose. */
314df555b71SRichard Tran Mills   stat = mkl_sparse_x_mv(SPARSE_OPERATION_TRANSPOSE,1.0,aijmkl->csrA,aijmkl->descr,x,0.0,y);
315df555b71SRichard Tran Mills 
316df555b71SRichard Tran Mills   ierr = PetscLogFlops(2.0*a->nz - a->nonzerorowcnt);CHKERRQ(ierr);
317df555b71SRichard Tran Mills   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
318df555b71SRichard Tran Mills   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
319df555b71SRichard Tran Mills   if (stat != SPARSE_STATUS_SUCCESS) {
320df555b71SRichard Tran Mills     PetscFunctionReturn(PETSC_ERR_LIB);
321df555b71SRichard Tran Mills   }
322df555b71SRichard Tran Mills   PetscFunctionReturn(0);
323df555b71SRichard Tran Mills }
324d995685eSRichard Tran Mills #endif /* PETSC_HAVE_MKL_SPARSE_OPTIMIZE */
325df555b71SRichard Tran Mills 
3264a2a386eSRichard Tran Mills PetscErrorCode MatMultAdd_SeqAIJMKL(Mat A,Vec xx,Vec yy,Vec zz)
3274a2a386eSRichard Tran Mills {
3284a2a386eSRichard Tran Mills   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
3294a2a386eSRichard Tran Mills   const PetscScalar *x;
3304a2a386eSRichard Tran Mills   PetscScalar       *y,*z;
3314a2a386eSRichard Tran Mills   const MatScalar   *aa;
3324a2a386eSRichard Tran Mills   PetscErrorCode    ierr;
3334a2a386eSRichard Tran Mills   PetscInt          m=A->rmap->n;
3344a2a386eSRichard Tran Mills   const PetscInt    *aj,*ai;
3354a2a386eSRichard Tran Mills   PetscInt          i;
3364a2a386eSRichard Tran Mills 
337ff03dc53SRichard Tran Mills   /* Variables not in MatMultAdd_SeqAIJ. */
338ff03dc53SRichard Tran Mills   char transa = 'n';  /* Used to indicate to MKL that we are not computing the transpose product. */
339a84739b8SRichard Tran Mills   PetscScalar       alpha = 1.0;
340a84739b8SRichard Tran Mills   PetscScalar       beta = 1.0;
341a84739b8SRichard Tran Mills   char              matdescra[6];
342ff03dc53SRichard Tran Mills 
343ff03dc53SRichard Tran Mills   PetscFunctionBegin;
344a84739b8SRichard Tran Mills   matdescra[0] = 'g';  /* Indicates to MKL that we using a general CSR matrix. */
345a84739b8SRichard Tran Mills   matdescra[3] = 'c';  /* Indicates to MKL that we use C-style (0-based) indexing. */
346a84739b8SRichard Tran Mills 
347ff03dc53SRichard Tran Mills   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
348ff03dc53SRichard Tran Mills   ierr = VecGetArrayPair(yy,zz,&y,&z);CHKERRQ(ierr);
349ff03dc53SRichard Tran Mills   aj   = a->j;  /* aj[k] gives column index for element aa[k]. */
350ff03dc53SRichard Tran Mills   aa   = a->a;  /* Nonzero elements stored row-by-row. */
351ff03dc53SRichard Tran Mills   ai   = a->i;  /* ai[k] is the position in aa and aj where row k starts. */
352ff03dc53SRichard Tran Mills 
353ff03dc53SRichard Tran Mills   /* Call MKL sparse BLAS routine to do the MatMult. */
354a84739b8SRichard Tran Mills   if (zz == yy) {
355a84739b8SRichard 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. */
356a84739b8SRichard Tran Mills     mkl_xcsrmv(&transa,&m,&m,&alpha,matdescra,aa,aj,ai,ai+1,x,&beta,y);
357a84739b8SRichard Tran Mills   } else {
358a84739b8SRichard Tran Mills     /* zz and yy are different vectors, so we call mkl_cspblas_xcsrgemv(), which calculates y = A*x, and then
359a84739b8SRichard Tran Mills      * we add the contents of vector yy to the result; MKL sparse BLAS does not have a MatMultAdd equivalent. */
360ff03dc53SRichard Tran Mills     mkl_cspblas_xcsrgemv(&transa,&m,aa,ai,aj,x,z);
361ff03dc53SRichard Tran Mills     for (i=0; i<m; i++) {
362ff03dc53SRichard Tran Mills       z[i] += y[i];
363ff03dc53SRichard Tran Mills     }
364a84739b8SRichard Tran Mills   }
365ff03dc53SRichard Tran Mills 
366ff03dc53SRichard Tran Mills   ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr);
367ff03dc53SRichard Tran Mills   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
368ff03dc53SRichard Tran Mills   ierr = VecRestoreArrayPair(yy,zz,&y,&z);CHKERRQ(ierr);
369ff03dc53SRichard Tran Mills   PetscFunctionReturn(0);
370ff03dc53SRichard Tran Mills }
371ff03dc53SRichard Tran Mills 
372d995685eSRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_OPTIMIZE
373df555b71SRichard Tran Mills PetscErrorCode MatMultAdd_SeqAIJMKL_SpMV2(Mat A,Vec xx,Vec yy,Vec zz)
374df555b71SRichard Tran Mills {
375df555b71SRichard Tran Mills   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
376df555b71SRichard Tran Mills   Mat_SeqAIJMKL     *aijmkl=(Mat_SeqAIJMKL*)A->spptr;
377df555b71SRichard Tran Mills   const PetscScalar *x;
378df555b71SRichard Tran Mills   PetscScalar       *y,*z;
379df555b71SRichard Tran Mills   PetscErrorCode    ierr;
380df555b71SRichard Tran Mills   PetscInt          m=A->rmap->n;
381df555b71SRichard Tran Mills   PetscInt          i;
382df555b71SRichard Tran Mills 
383df555b71SRichard Tran Mills   /* Variables not in MatMultAdd_SeqAIJ. */
384df555b71SRichard Tran Mills   sparse_status_t stat = SPARSE_STATUS_SUCCESS;
385df555b71SRichard Tran Mills 
386df555b71SRichard Tran Mills   PetscFunctionBegin;
387df555b71SRichard Tran Mills 
388f36dfe3fSRichard Tran Mills   /* If there are no rows, this is a no-op: return immediately. */
389f36dfe3fSRichard Tran Mills   if(A->cmap->n < 1) PetscFunctionReturn(0);
390df555b71SRichard Tran Mills 
391df555b71SRichard Tran Mills   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
392df555b71SRichard Tran Mills   ierr = VecGetArrayPair(yy,zz,&y,&z);CHKERRQ(ierr);
393df555b71SRichard Tran Mills 
394df555b71SRichard Tran Mills   /* Call MKL sparse BLAS routine to do the MatMult. */
395df555b71SRichard Tran Mills   if (zz == yy) {
396df555b71SRichard 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,
397df555b71SRichard Tran Mills      * with alpha and beta both set to 1.0. */
398df555b71SRichard Tran Mills     stat = mkl_sparse_x_mv(SPARSE_OPERATION_NON_TRANSPOSE,1.0,aijmkl->csrA,aijmkl->descr,x,1.0,y);
399df555b71SRichard Tran Mills   } else {
400df555b71SRichard 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
401df555b71SRichard Tran Mills      * we add the contents of vector yy to the result; MKL sparse BLAS does not have a MatMultAdd equivalent. */
402df555b71SRichard Tran Mills     stat = mkl_sparse_x_mv(SPARSE_OPERATION_NON_TRANSPOSE,1.0,aijmkl->csrA,aijmkl->descr,x,0.0,y);
403df555b71SRichard Tran Mills     for (i=0; i<m; i++) {
404df555b71SRichard Tran Mills       z[i] += y[i];
405df555b71SRichard Tran Mills     }
406df555b71SRichard Tran Mills   }
407df555b71SRichard Tran Mills 
408df555b71SRichard Tran Mills   ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr);
409df555b71SRichard Tran Mills   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
410df555b71SRichard Tran Mills   ierr = VecRestoreArrayPair(yy,zz,&y,&z);CHKERRQ(ierr);
411df555b71SRichard Tran Mills   if (stat != SPARSE_STATUS_SUCCESS) {
412df555b71SRichard Tran Mills     PetscFunctionReturn(PETSC_ERR_LIB);
413df555b71SRichard Tran Mills   }
414df555b71SRichard Tran Mills   PetscFunctionReturn(0);
415df555b71SRichard Tran Mills }
416d995685eSRichard Tran Mills #endif /* PETSC_HAVE_MKL_SPARSE_OPTIMIZE */
417df555b71SRichard Tran Mills 
418ff03dc53SRichard Tran Mills PetscErrorCode MatMultTransposeAdd_SeqAIJMKL(Mat A,Vec xx,Vec yy,Vec zz)
419ff03dc53SRichard Tran Mills {
420ff03dc53SRichard Tran Mills   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
421ff03dc53SRichard Tran Mills   const PetscScalar *x;
422ff03dc53SRichard Tran Mills   PetscScalar       *y,*z;
423ff03dc53SRichard Tran Mills   const MatScalar   *aa;
424ff03dc53SRichard Tran Mills   PetscErrorCode    ierr;
425ff03dc53SRichard Tran Mills   PetscInt          m=A->rmap->n;
426ff03dc53SRichard Tran Mills   const PetscInt    *aj,*ai;
427ff03dc53SRichard Tran Mills   PetscInt          i;
428ff03dc53SRichard Tran Mills 
429ff03dc53SRichard Tran Mills   /* Variables not in MatMultTransposeAdd_SeqAIJ. */
430ff03dc53SRichard Tran Mills   char transa = 't';  /* Used to indicate to MKL that we are computing the transpose product. */
431a84739b8SRichard Tran Mills   PetscScalar       alpha = 1.0;
432a84739b8SRichard Tran Mills   PetscScalar       beta = 1.0;
433a84739b8SRichard Tran Mills   char              matdescra[6];
4344a2a386eSRichard Tran Mills 
4354a2a386eSRichard Tran Mills   PetscFunctionBegin;
436a84739b8SRichard Tran Mills   matdescra[0] = 'g';  /* Indicates to MKL that we using a general CSR matrix. */
437a84739b8SRichard Tran Mills   matdescra[3] = 'c';  /* Indicates to MKL that we use C-style (0-based) indexing. */
438a84739b8SRichard Tran Mills 
4394a2a386eSRichard Tran Mills   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
4404a2a386eSRichard Tran Mills   ierr = VecGetArrayPair(yy,zz,&y,&z);CHKERRQ(ierr);
4414a2a386eSRichard Tran Mills   aj   = a->j;  /* aj[k] gives column index for element aa[k]. */
4424a2a386eSRichard Tran Mills   aa   = a->a;  /* Nonzero elements stored row-by-row. */
4434a2a386eSRichard Tran Mills   ai   = a->i;  /* ai[k] is the position in aa and aj where row k starts. */
4444a2a386eSRichard Tran Mills 
4454a2a386eSRichard Tran Mills   /* Call MKL sparse BLAS routine to do the MatMult. */
446a84739b8SRichard Tran Mills   if (zz == yy) {
447a84739b8SRichard 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. */
448a84739b8SRichard Tran Mills     mkl_xcsrmv(&transa,&m,&m,&alpha,matdescra,aa,aj,ai,ai+1,x,&beta,y);
449a84739b8SRichard Tran Mills   } else {
450a84739b8SRichard Tran Mills     /* zz and yy are different vectors, so we call mkl_cspblas_xcsrgemv(), which calculates y = A*x, and then
451a84739b8SRichard Tran Mills      * we add the contents of vector yy to the result; MKL sparse BLAS does not have a MatMultAdd equivalent. */
4524a2a386eSRichard Tran Mills     mkl_cspblas_xcsrgemv(&transa,&m,aa,ai,aj,x,z);
4534a2a386eSRichard Tran Mills     for (i=0; i<m; i++) {
4544a2a386eSRichard Tran Mills       z[i] += y[i];
4554a2a386eSRichard Tran Mills     }
456a84739b8SRichard Tran Mills   }
4574a2a386eSRichard Tran Mills 
4584a2a386eSRichard Tran Mills   ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr);
4594a2a386eSRichard Tran Mills   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
4604a2a386eSRichard Tran Mills   ierr = VecRestoreArrayPair(yy,zz,&y,&z);CHKERRQ(ierr);
4614a2a386eSRichard Tran Mills   PetscFunctionReturn(0);
4624a2a386eSRichard Tran Mills }
4634a2a386eSRichard Tran Mills 
464d995685eSRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_OPTIMIZE
465df555b71SRichard Tran Mills PetscErrorCode MatMultTransposeAdd_SeqAIJMKL_SpMV2(Mat A,Vec xx,Vec yy,Vec zz)
466df555b71SRichard Tran Mills {
467df555b71SRichard Tran Mills   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
468df555b71SRichard Tran Mills   Mat_SeqAIJMKL     *aijmkl=(Mat_SeqAIJMKL*)A->spptr;
469df555b71SRichard Tran Mills   const PetscScalar *x;
470df555b71SRichard Tran Mills   PetscScalar       *y,*z;
471df555b71SRichard Tran Mills   PetscErrorCode    ierr;
472df555b71SRichard Tran Mills   PetscInt          m=A->rmap->n;
473df555b71SRichard Tran Mills   PetscInt          i;
474df555b71SRichard Tran Mills 
475df555b71SRichard Tran Mills   /* Variables not in MatMultTransposeAdd_SeqAIJ. */
476df555b71SRichard Tran Mills   sparse_status_t stat = SPARSE_STATUS_SUCCESS;
477df555b71SRichard Tran Mills 
478df555b71SRichard Tran Mills   PetscFunctionBegin;
479df555b71SRichard Tran Mills 
480f36dfe3fSRichard Tran Mills   /* If there are no rows, this is a no-op: return immediately. */
481f36dfe3fSRichard Tran Mills   if(A->cmap->n < 1) PetscFunctionReturn(0);
482f36dfe3fSRichard Tran Mills 
483df555b71SRichard Tran Mills   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
484df555b71SRichard Tran Mills   ierr = VecGetArrayPair(yy,zz,&y,&z);CHKERRQ(ierr);
485df555b71SRichard Tran Mills 
486df555b71SRichard Tran Mills   /* Call MKL sparse BLAS routine to do the MatMult. */
487df555b71SRichard Tran Mills   if (zz == yy) {
488df555b71SRichard 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,
489df555b71SRichard Tran Mills      * with alpha and beta both set to 1.0. */
490df555b71SRichard Tran Mills     stat = mkl_sparse_x_mv(SPARSE_OPERATION_TRANSPOSE,1.0,aijmkl->csrA,aijmkl->descr,x,1.0,y);
491df555b71SRichard Tran Mills   } else {
492df555b71SRichard 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
493df555b71SRichard Tran Mills      * we add the contents of vector yy to the result; MKL sparse BLAS does not have a MatMultAdd equivalent. */
494df555b71SRichard Tran Mills     stat = mkl_sparse_x_mv(SPARSE_OPERATION_TRANSPOSE,1.0,aijmkl->csrA,aijmkl->descr,x,0.0,y);
495df555b71SRichard Tran Mills     for (i=0; i<m; i++) {
496df555b71SRichard Tran Mills       z[i] += y[i];
497df555b71SRichard Tran Mills     }
498df555b71SRichard Tran Mills   }
499df555b71SRichard Tran Mills 
500df555b71SRichard Tran Mills   ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr);
501df555b71SRichard Tran Mills   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
502df555b71SRichard Tran Mills   ierr = VecRestoreArrayPair(yy,zz,&y,&z);CHKERRQ(ierr);
503df555b71SRichard Tran Mills   if (stat != SPARSE_STATUS_SUCCESS) {
504df555b71SRichard Tran Mills     PetscFunctionReturn(PETSC_ERR_LIB);
505df555b71SRichard Tran Mills   }
506df555b71SRichard Tran Mills   PetscFunctionReturn(0);
507df555b71SRichard Tran Mills }
508d995685eSRichard Tran Mills #endif /* PETSC_HAVE_MKL_SPARSE_OPTIMIZE */
509df555b71SRichard Tran Mills 
510df555b71SRichard Tran Mills 
5114a2a386eSRichard Tran Mills /* MatConvert_SeqAIJ_SeqAIJMKL converts a SeqAIJ matrix into a
5124a2a386eSRichard Tran Mills  * SeqAIJMKL matrix.  This routine is called by the MatCreate_SeqMKLAIJ()
5134a2a386eSRichard Tran Mills  * routine, but can also be used to convert an assembled SeqAIJ matrix
5144a2a386eSRichard Tran Mills  * into a SeqAIJMKL one. */
5154a2a386eSRichard Tran Mills PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJMKL(Mat A,MatType type,MatReuse reuse,Mat *newmat)
5164a2a386eSRichard Tran Mills {
5174a2a386eSRichard Tran Mills   PetscErrorCode ierr;
5184a2a386eSRichard Tran Mills   Mat            B = *newmat;
5194a2a386eSRichard Tran Mills   Mat_SeqAIJMKL  *aijmkl;
520c9d46305SRichard Tran Mills   PetscBool      set;
521e9c94282SRichard Tran Mills   PetscBool      sametype;
5224a2a386eSRichard Tran Mills 
5234a2a386eSRichard Tran Mills   PetscFunctionBegin;
5244a2a386eSRichard Tran Mills   if (reuse == MAT_INITIAL_MATRIX) {
5254a2a386eSRichard Tran Mills     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
5264a2a386eSRichard Tran Mills   }
5274a2a386eSRichard Tran Mills 
528e9c94282SRichard Tran Mills   ierr = PetscObjectTypeCompare((PetscObject)A,type,&sametype);CHKERRQ(ierr);
529e9c94282SRichard Tran Mills   if (sametype) PetscFunctionReturn(0);
530e9c94282SRichard Tran Mills 
5314a2a386eSRichard Tran Mills   ierr     = PetscNewLog(B,&aijmkl);CHKERRQ(ierr);
5324a2a386eSRichard Tran Mills   B->spptr = (void*) aijmkl;
5334a2a386eSRichard Tran Mills 
534df555b71SRichard Tran Mills   /* Set function pointers for methods that we inherit from AIJ but override.
535e9c94282SRichard Tran Mills    * We also parse some command line options below, since those determine some of the methods we point to.
536e9c94282SRichard Tran Mills    * Note: Currently the transposed operations are not being set because I encounter memory corruption
537df555b71SRichard Tran Mills    * when these are enabled.  Need to look at this with Valgrind or similar. --RTM */
5384a2a386eSRichard Tran Mills   B->ops->duplicate        = MatDuplicate_SeqAIJMKL;
5394a2a386eSRichard Tran Mills   B->ops->assemblyend      = MatAssemblyEnd_SeqAIJMKL;
5404a2a386eSRichard Tran Mills   B->ops->destroy          = MatDestroy_SeqAIJMKL;
541c9d46305SRichard Tran Mills 
5424abfa3b3SRichard Tran Mills   aijmkl->sparse_optimized = PETSC_FALSE;
543d995685eSRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_OPTIMIZE
544d995685eSRichard Tran Mills   aijmkl->no_SpMV2 = PETSC_FALSE;  /* Default to using the SpMV2 routines if our MKL supports them. */
545d995685eSRichard Tran Mills #elif
546d995685eSRichard Tran Mills   aijmkl->no_SpMV2 = PETSC_TRUE;
547d995685eSRichard Tran Mills #endif
5484abfa3b3SRichard Tran Mills 
5494abfa3b3SRichard Tran Mills   /* Parse command line options. */
550c9d46305SRichard Tran Mills   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"AIJMKL Options","Mat");CHKERRQ(ierr);
551c9d46305SRichard Tran Mills   ierr = PetscOptionsBool("-mat_aijmkl_no_spmv2","NoSPMV2","None",(PetscBool)aijmkl->no_SpMV2,(PetscBool*)&aijmkl->no_SpMV2,&set);CHKERRQ(ierr);
552c9d46305SRichard Tran Mills   ierr = PetscOptionsEnd();CHKERRQ(ierr);
553d995685eSRichard Tran Mills #ifndef PETSC_HAVE_MKL_SPARSE_OPTIMIZE
554d995685eSRichard Tran Mills   if(!aijmkl->no_SpMV2) {
555d995685eSRichard 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");
556d995685eSRichard Tran Mills     aijmkl->no_SpMV2 = PETSC_TRUE;
557d995685eSRichard Tran Mills   }
558d995685eSRichard Tran Mills #endif
559c9d46305SRichard Tran Mills 
560c9d46305SRichard Tran Mills   if(!aijmkl->no_SpMV2) {
561d995685eSRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_OPTIMIZE
562df555b71SRichard Tran Mills     B->ops->mult             = MatMult_SeqAIJMKL_SpMV2;
563df555b71SRichard Tran Mills     /* B->ops->multtranspose    = MatMultTranspose_SeqAIJMKL_SpMV2; */
564df555b71SRichard Tran Mills     B->ops->multadd          = MatMultAdd_SeqAIJMKL_SpMV2;
565df555b71SRichard Tran Mills     /* B->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJMKL_SpMV2; */
566d995685eSRichard Tran Mills #endif
567c9d46305SRichard Tran Mills   } else {
5684a2a386eSRichard Tran Mills     B->ops->mult             = MatMult_SeqAIJMKL;
569c9d46305SRichard Tran Mills     /* B->ops->multtranspose    = MatMultTranspose_SeqAIJMKL; */
5704a2a386eSRichard Tran Mills     B->ops->multadd          = MatMultAdd_SeqAIJMKL;
571c9d46305SRichard Tran Mills     /* B->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJMKL; */
572c9d46305SRichard Tran Mills   }
5734a2a386eSRichard Tran Mills 
5744a2a386eSRichard Tran Mills   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaijmkl_seqaij_C",MatConvert_SeqAIJMKL_SeqAIJ);CHKERRQ(ierr);
575e9c94282SRichard Tran Mills   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMult_seqdense_seqaijmkl_C",MatMatMult_SeqDense_SeqAIJ);CHKERRQ(ierr);
576e9c94282SRichard Tran Mills   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMultSymbolic_seqdense_seqaijmkl_C",MatMatMultSymbolic_SeqDense_SeqAIJ);CHKERRQ(ierr);
577e9c94282SRichard Tran Mills   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMultNumeric_seqdense_seqaijmkl_C",MatMatMultNumeric_SeqDense_SeqAIJ);CHKERRQ(ierr);
5784a2a386eSRichard Tran Mills 
5794a2a386eSRichard Tran Mills   ierr    = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJMKL);CHKERRQ(ierr);
5804a2a386eSRichard Tran Mills   *newmat = B;
5814a2a386eSRichard Tran Mills   PetscFunctionReturn(0);
5824a2a386eSRichard Tran Mills }
5834a2a386eSRichard Tran Mills 
5844a2a386eSRichard Tran Mills /*@C
5854a2a386eSRichard Tran Mills    MatCreateSeqAIJMKL - Creates a sparse matrix of type SEQAIJMKL.
5864a2a386eSRichard Tran Mills    This type inherits from AIJ and is largely identical, but uses sparse BLAS
5874a2a386eSRichard Tran Mills    routines from Intel MKL whenever possible.
5884a2a386eSRichard Tran Mills    Collective on MPI_Comm
5894a2a386eSRichard Tran Mills 
5904a2a386eSRichard Tran Mills    Input Parameters:
5914a2a386eSRichard Tran Mills +  comm - MPI communicator, set to PETSC_COMM_SELF
5924a2a386eSRichard Tran Mills .  m - number of rows
5934a2a386eSRichard Tran Mills .  n - number of columns
5944a2a386eSRichard Tran Mills .  nz - number of nonzeros per row (same for all rows)
5954a2a386eSRichard Tran Mills -  nnz - array containing the number of nonzeros in the various rows
5964a2a386eSRichard Tran Mills          (possibly different for each row) or NULL
5974a2a386eSRichard Tran Mills 
5984a2a386eSRichard Tran Mills    Output Parameter:
5994a2a386eSRichard Tran Mills .  A - the matrix
6004a2a386eSRichard Tran Mills 
6014a2a386eSRichard Tran Mills    Notes:
6024a2a386eSRichard Tran Mills    If nnz is given then nz is ignored
6034a2a386eSRichard Tran Mills 
6044a2a386eSRichard Tran Mills    Level: intermediate
6054a2a386eSRichard Tran Mills 
6064a2a386eSRichard Tran Mills .keywords: matrix, cray, sparse, parallel
6074a2a386eSRichard Tran Mills 
6084a2a386eSRichard Tran Mills .seealso: MatCreate(), MatCreateMPIAIJMKL(), MatSetValues()
6094a2a386eSRichard Tran Mills @*/
6104a2a386eSRichard Tran Mills PetscErrorCode  MatCreateSeqAIJMKL(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
6114a2a386eSRichard Tran Mills {
6124a2a386eSRichard Tran Mills   PetscErrorCode ierr;
6134a2a386eSRichard Tran Mills 
6144a2a386eSRichard Tran Mills   PetscFunctionBegin;
6154a2a386eSRichard Tran Mills   ierr = MatCreate(comm,A);CHKERRQ(ierr);
6164a2a386eSRichard Tran Mills   ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr);
6174a2a386eSRichard Tran Mills   ierr = MatSetType(*A,MATSEQAIJMKL);CHKERRQ(ierr);
6184a2a386eSRichard Tran Mills   ierr = MatSeqAIJSetPreallocation_SeqAIJ(*A,nz,nnz);CHKERRQ(ierr);
6194a2a386eSRichard Tran Mills   PetscFunctionReturn(0);
6204a2a386eSRichard Tran Mills }
6214a2a386eSRichard Tran Mills 
6224a2a386eSRichard Tran Mills PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJMKL(Mat A)
6234a2a386eSRichard Tran Mills {
6244a2a386eSRichard Tran Mills   PetscErrorCode ierr;
6254a2a386eSRichard Tran Mills 
6264a2a386eSRichard Tran Mills   PetscFunctionBegin;
6274a2a386eSRichard Tran Mills   ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
6284a2a386eSRichard Tran Mills   ierr = MatConvert_SeqAIJ_SeqAIJMKL(A,MATSEQAIJMKL,MAT_INPLACE_MATRIX,&A);CHKERRQ(ierr);
6294a2a386eSRichard Tran Mills   PetscFunctionReturn(0);
6304a2a386eSRichard Tran Mills }
631