xref: /petsc/src/mat/impls/aij/seq/aijmkl/aijmkl.c (revision ff03dc5341f8e08c034bb25cd05c5e9387d3fbf2)
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 {
164a2a386eSRichard Tran Mills   /* "Handle" used by SpMV2 inspector-executor routines. */
174a2a386eSRichard Tran Mills   sparse_matrix_t csrA;
184a2a386eSRichard Tran Mills } Mat_SeqAIJMKL;
194a2a386eSRichard Tran Mills 
204a2a386eSRichard Tran Mills extern PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat,MatAssemblyType);
214a2a386eSRichard Tran Mills 
224a2a386eSRichard Tran Mills #undef __FUNCT__
234a2a386eSRichard Tran Mills #define __FUNCT__ "MatConvert_SeqAIJMKL_SeqAIJ"
244a2a386eSRichard Tran Mills PETSC_INTERN PetscErrorCode MatConvert_SeqAIJMKL_SeqAIJ(Mat A,MatType type,MatReuse reuse,Mat *newmat)
254a2a386eSRichard Tran Mills {
264a2a386eSRichard Tran Mills   /* This routine is only called to convert a MATAIJMKL to its base PETSc type, */
274a2a386eSRichard Tran Mills   /* so we will ignore 'MatType type'. */
284a2a386eSRichard Tran Mills   PetscErrorCode ierr;
294a2a386eSRichard Tran Mills   Mat            B       = *newmat;
304a2a386eSRichard Tran Mills   Mat_SeqAIJMKL *aijmkl=(Mat_SeqAIJMKL*)A->spptr;
314a2a386eSRichard Tran Mills 
324a2a386eSRichard Tran Mills   PetscFunctionBegin;
334a2a386eSRichard Tran Mills   if (reuse == MAT_INITIAL_MATRIX) {
344a2a386eSRichard Tran Mills     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
354a2a386eSRichard Tran Mills   }
364a2a386eSRichard Tran Mills 
374a2a386eSRichard Tran Mills   /* Reset the original function pointers. */
3854871a98SRichard Tran Mills   B->ops->duplicate        = MatDuplicate_SeqAIJ;
394a2a386eSRichard Tran Mills   B->ops->assemblyend      = MatAssemblyEnd_SeqAIJ;
404a2a386eSRichard Tran Mills   B->ops->destroy          = MatDestroy_SeqAIJ;
4154871a98SRichard Tran Mills   B->ops->mult             = MatMult_SeqAIJ;
42*ff03dc53SRichard Tran Mills   B->ops->multtranspose    = MatMultTranspose_SeqAIJ;
4354871a98SRichard Tran Mills   B->ops->multadd          = MatMultAdd_SeqAIJ;
44*ff03dc53SRichard Tran Mills   B->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJ;
454a2a386eSRichard Tran Mills 
464a2a386eSRichard Tran Mills   /* Free everything in the Mat_SeqAIJMKL data structure.
474a2a386eSRichard Tran Mills    * We don't free the Mat_SeqAIJMKL struct itself, as this will
484a2a386eSRichard Tran Mills    * cause problems later when MatDestroy() tries to free it. */
494a2a386eSRichard Tran Mills   /* Actually there is nothing to do here right now.
504a2a386eSRichard Tran Mills    * When I've added use of the MKL SpMV2 inspector-executor routines, I should
514a2a386eSRichard Tran Mills    * see if there is some way to clean up the "handle" used by SpMV2. */
524a2a386eSRichard Tran Mills 
534a2a386eSRichard Tran Mills   /* Change the type of B to MATSEQAIJ. */
544a2a386eSRichard Tran Mills   ierr = PetscObjectChangeTypeName((PetscObject)B, MATSEQAIJ);CHKERRQ(ierr);
554a2a386eSRichard Tran Mills 
564a2a386eSRichard Tran Mills   *newmat = B;
574a2a386eSRichard Tran Mills   PetscFunctionReturn(0);
584a2a386eSRichard Tran Mills }
594a2a386eSRichard Tran Mills 
604a2a386eSRichard Tran Mills #undef __FUNCT__
614a2a386eSRichard Tran Mills #define __FUNCT__ "MatDestroy_SeqAIJMKL"
624a2a386eSRichard Tran Mills PetscErrorCode MatDestroy_SeqAIJMKL(Mat A)
634a2a386eSRichard Tran Mills {
644a2a386eSRichard Tran Mills   PetscErrorCode ierr;
654a2a386eSRichard Tran Mills   Mat_SeqAIJMKL *aijmkl = (Mat_SeqAIJMKL*) A->spptr;
664a2a386eSRichard Tran Mills 
674a2a386eSRichard Tran Mills   PetscFunctionBegin;
684a2a386eSRichard Tran Mills   /* Clean up everything in the Mat_SeqAIJMKL data structure, then free A->spptr. */
694a2a386eSRichard Tran Mills   mkl_sparse_destroy(aijmkl->csrA);
704a2a386eSRichard Tran Mills   ierr = PetscFree(A->spptr);CHKERRQ(ierr);
714a2a386eSRichard Tran Mills 
724a2a386eSRichard Tran Mills   /* Change the type of A back to SEQAIJ and use MatDestroy_SeqAIJ()
734a2a386eSRichard Tran Mills    * to destroy everything that remains. */
744a2a386eSRichard Tran Mills   ierr = PetscObjectChangeTypeName((PetscObject)A, MATSEQAIJ);CHKERRQ(ierr);
754a2a386eSRichard Tran Mills   /* Note that I don't call MatSetType().  I believe this is because that
764a2a386eSRichard Tran Mills    * is only to be called when *building* a matrix.  I could be wrong, but
774a2a386eSRichard Tran Mills    * that is how things work for the SuperLU matrix class. */
784a2a386eSRichard Tran Mills   ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
794a2a386eSRichard Tran Mills   PetscFunctionReturn(0);
804a2a386eSRichard Tran Mills }
814a2a386eSRichard Tran Mills 
824a2a386eSRichard Tran Mills #undef __FUNCT__
834a2a386eSRichard Tran Mills #define __FUNCT__ "MatDuplicate_SeqAIJMKL"
844a2a386eSRichard Tran Mills PetscErrorCode MatDuplicate_SeqAIJMKL(Mat A, MatDuplicateOption op, Mat *M)
854a2a386eSRichard Tran Mills {
864a2a386eSRichard Tran Mills   PetscErrorCode ierr;
874a2a386eSRichard Tran Mills   Mat_SeqAIJMKL *aijmkl      = (Mat_SeqAIJMKL*) A->spptr;
884a2a386eSRichard Tran Mills   Mat_SeqAIJMKL *aijmkl_dest = (Mat_SeqAIJMKL*) (*M)->spptr;
894a2a386eSRichard Tran Mills 
904a2a386eSRichard Tran Mills   PetscFunctionBegin;
914a2a386eSRichard Tran Mills   ierr = MatDuplicate_SeqAIJ(A,op,M);CHKERRQ(ierr);
924a2a386eSRichard Tran Mills   ierr = PetscMemcpy((*M)->spptr,aijmkl,sizeof(Mat_SeqAIJMKL));CHKERRQ(ierr);
934a2a386eSRichard Tran Mills   PetscFunctionReturn(0);
944a2a386eSRichard Tran Mills }
954a2a386eSRichard Tran Mills 
964a2a386eSRichard Tran Mills #undef __FUNCT__
974a2a386eSRichard Tran Mills #define __FUNCT__ "MatAssemblyEnd_SeqAIJMKL"
984a2a386eSRichard Tran Mills PetscErrorCode MatAssemblyEnd_SeqAIJMKL(Mat A, MatAssemblyType mode)
994a2a386eSRichard Tran Mills {
1004a2a386eSRichard Tran Mills   PetscErrorCode ierr;
1014a2a386eSRichard Tran Mills   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1024a2a386eSRichard Tran Mills 
1034a2a386eSRichard Tran Mills   PetscFunctionBegin;
1044a2a386eSRichard Tran Mills   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
1054a2a386eSRichard Tran Mills 
1064a2a386eSRichard Tran Mills   /* Since a MATSEQAIJMKL matrix is really just a MATSEQAIJ with some
1074a2a386eSRichard Tran Mills    * extra information and some different methods, call the AssemblyEnd
1084a2a386eSRichard Tran Mills    * routine for a MATSEQAIJ.
1094a2a386eSRichard Tran Mills    * I'm not sure if this is the best way to do this, but it avoids
1104a2a386eSRichard Tran Mills    * a lot of code duplication.
1114a2a386eSRichard Tran Mills    * I also note that currently MATSEQAIJMKL doesn't know anything about
1124a2a386eSRichard Tran Mills    * the Mat_CompressedRow data structure that SeqAIJ now uses when there
1134a2a386eSRichard Tran Mills    * are many zero rows.  If the SeqAIJ assembly end routine decides to use
1144a2a386eSRichard Tran Mills    * this, this may break things.  (Don't know... haven't looked at it.
1154a2a386eSRichard Tran Mills    * Do I need to disable this somehow?) */
1164a2a386eSRichard Tran Mills   a->inode.use = PETSC_FALSE;  /* Must disable: otherwise the MKL routines won't get used. */
1174a2a386eSRichard Tran Mills   ierr         = MatAssemblyEnd_SeqAIJ(A, mode);CHKERRQ(ierr);
1184a2a386eSRichard Tran Mills 
1194a2a386eSRichard Tran Mills   PetscFunctionReturn(0);
1204a2a386eSRichard Tran Mills }
1214a2a386eSRichard Tran Mills 
1224a2a386eSRichard Tran Mills #undef __FUNCT__
1234a2a386eSRichard Tran Mills #define __FUNCT__ "MatMult_SeqAIJMKL"
1244a2a386eSRichard Tran Mills PetscErrorCode MatMult_SeqAIJMKL(Mat A,Vec xx,Vec yy)
1254a2a386eSRichard Tran Mills {
1264a2a386eSRichard Tran Mills   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
1274a2a386eSRichard Tran Mills   const PetscScalar *x;
1284a2a386eSRichard Tran Mills   PetscScalar       *y;
1294a2a386eSRichard Tran Mills   const MatScalar   *aa;
1304a2a386eSRichard Tran Mills   PetscErrorCode    ierr;
1314a2a386eSRichard Tran Mills   PetscInt          m=A->rmap->n;
1324a2a386eSRichard Tran Mills   const PetscInt    *aj,*ai;
1334a2a386eSRichard Tran Mills   PetscInt          i;
1344a2a386eSRichard Tran Mills 
1354a2a386eSRichard Tran Mills   /* Variables not in MatMult_SeqAIJ. */
136*ff03dc53SRichard Tran Mills   char transa = 'n';  /* Used to indicate to MKL that we are not computing the transpose product. */
137*ff03dc53SRichard Tran Mills 
138*ff03dc53SRichard Tran Mills   PetscFunctionBegin;
139*ff03dc53SRichard Tran Mills   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
140*ff03dc53SRichard Tran Mills   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
141*ff03dc53SRichard Tran Mills   aj   = a->j;  /* aj[k] gives column index for element aa[k]. */
142*ff03dc53SRichard Tran Mills   aa   = a->a;  /* Nonzero elements stored row-by-row. */
143*ff03dc53SRichard Tran Mills   ai   = a->i;  /* ai[k] is the position in aa and aj where row k starts. */
144*ff03dc53SRichard Tran Mills 
145*ff03dc53SRichard Tran Mills   /* Call MKL sparse BLAS routine to do the MatMult. */
146*ff03dc53SRichard Tran Mills   mkl_cspblas_xcsrgemv(&transa,&m,aa,ai,aj,x,y);
147*ff03dc53SRichard Tran Mills 
148*ff03dc53SRichard Tran Mills   ierr = PetscLogFlops(2.0*a->nz - a->nonzerorowcnt);CHKERRQ(ierr);
149*ff03dc53SRichard Tran Mills   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
150*ff03dc53SRichard Tran Mills   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
151*ff03dc53SRichard Tran Mills   PetscFunctionReturn(0);
152*ff03dc53SRichard Tran Mills }
153*ff03dc53SRichard Tran Mills 
154*ff03dc53SRichard Tran Mills #undef __FUNCT__
155*ff03dc53SRichard Tran Mills #define __FUNCT__ "MatMultTranspose_SeqAIJMKL"
156*ff03dc53SRichard Tran Mills PetscErrorCode MatMultTranspose_SeqAIJMKL(Mat A,Vec xx,Vec yy)
157*ff03dc53SRichard Tran Mills {
158*ff03dc53SRichard Tran Mills   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
159*ff03dc53SRichard Tran Mills   const PetscScalar *x;
160*ff03dc53SRichard Tran Mills   PetscScalar       *y;
161*ff03dc53SRichard Tran Mills   const MatScalar   *aa;
162*ff03dc53SRichard Tran Mills   PetscErrorCode    ierr;
163*ff03dc53SRichard Tran Mills   PetscInt          m=A->rmap->n;
164*ff03dc53SRichard Tran Mills   const PetscInt    *aj,*ai;
165*ff03dc53SRichard Tran Mills   PetscInt          i;
166*ff03dc53SRichard Tran Mills 
167*ff03dc53SRichard Tran Mills   /* Variables not in MatMultTranspose_SeqAIJ. */
168*ff03dc53SRichard Tran Mills   char transa = 't';  /* Used to indicate to MKL that we are computing the transpose product. */
1694a2a386eSRichard Tran Mills 
1704a2a386eSRichard Tran Mills   PetscFunctionBegin;
1714a2a386eSRichard Tran Mills   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
1724a2a386eSRichard Tran Mills   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
1734a2a386eSRichard Tran Mills   aj   = a->j;  /* aj[k] gives column index for element aa[k]. */
1744a2a386eSRichard Tran Mills   aa   = a->a;  /* Nonzero elements stored row-by-row. */
1754a2a386eSRichard Tran Mills   ai   = a->i;  /* ai[k] is the position in aa and aj where row k starts. */
1764a2a386eSRichard Tran Mills 
1774a2a386eSRichard Tran Mills   /* Call MKL sparse BLAS routine to do the MatMult. */
1784a2a386eSRichard Tran Mills   mkl_cspblas_xcsrgemv(&transa,&m,aa,ai,aj,x,y);
1794a2a386eSRichard Tran Mills 
1804a2a386eSRichard Tran Mills   ierr = PetscLogFlops(2.0*a->nz - a->nonzerorowcnt);CHKERRQ(ierr);
1814a2a386eSRichard Tran Mills   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1824a2a386eSRichard Tran Mills   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
1834a2a386eSRichard Tran Mills   PetscFunctionReturn(0);
1844a2a386eSRichard Tran Mills }
1854a2a386eSRichard Tran Mills 
1864a2a386eSRichard Tran Mills #undef __FUNCT__
1874a2a386eSRichard Tran Mills #define __FUNCT__ "MatMultAdd_SeqAIJMKL"
1884a2a386eSRichard Tran Mills PetscErrorCode MatMultAdd_SeqAIJMKL(Mat A,Vec xx,Vec yy,Vec zz)
1894a2a386eSRichard Tran Mills {
1904a2a386eSRichard Tran Mills   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
1914a2a386eSRichard Tran Mills   const PetscScalar *x;
1924a2a386eSRichard Tran Mills   PetscScalar       *y,*z;
1934a2a386eSRichard Tran Mills   const MatScalar   *aa;
1944a2a386eSRichard Tran Mills   PetscErrorCode    ierr;
1954a2a386eSRichard Tran Mills   PetscInt          m=A->rmap->n;
1964a2a386eSRichard Tran Mills   const PetscInt    *aj,*ai;
1974a2a386eSRichard Tran Mills   PetscInt          i;
1984a2a386eSRichard Tran Mills 
199*ff03dc53SRichard Tran Mills   /* Variables not in MatMultAdd_SeqAIJ. */
200*ff03dc53SRichard Tran Mills   char transa = 'n';  /* Used to indicate to MKL that we are not computing the transpose product. */
201*ff03dc53SRichard Tran Mills 
202*ff03dc53SRichard Tran Mills   PetscFunctionBegin;
203*ff03dc53SRichard Tran Mills   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
204*ff03dc53SRichard Tran Mills   ierr = VecGetArrayPair(yy,zz,&y,&z);CHKERRQ(ierr);
205*ff03dc53SRichard Tran Mills   aj   = a->j;  /* aj[k] gives column index for element aa[k]. */
206*ff03dc53SRichard Tran Mills   aa   = a->a;  /* Nonzero elements stored row-by-row. */
207*ff03dc53SRichard Tran Mills   ai   = a->i;  /* ai[k] is the position in aa and aj where row k starts. */
208*ff03dc53SRichard Tran Mills 
209*ff03dc53SRichard Tran Mills   /* Call MKL sparse BLAS routine to do the MatMult. */
210*ff03dc53SRichard Tran Mills   mkl_cspblas_xcsrgemv(&transa,&m,aa,ai,aj,x,z);
211*ff03dc53SRichard Tran Mills 
212*ff03dc53SRichard Tran Mills   /* Now add the vector y; MKL sparse BLAS does not have a MatMultAdd equivalent. */
213*ff03dc53SRichard Tran Mills   for (i=0; i<m; i++) {
214*ff03dc53SRichard Tran Mills     z[i] += y[i];
215*ff03dc53SRichard Tran Mills   }
216*ff03dc53SRichard Tran Mills 
217*ff03dc53SRichard Tran Mills   ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr);
218*ff03dc53SRichard Tran Mills   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
219*ff03dc53SRichard Tran Mills   ierr = VecRestoreArrayPair(yy,zz,&y,&z);CHKERRQ(ierr);
220*ff03dc53SRichard Tran Mills   PetscFunctionReturn(0);
221*ff03dc53SRichard Tran Mills }
222*ff03dc53SRichard Tran Mills 
223*ff03dc53SRichard Tran Mills #undef __FUNCT__
224*ff03dc53SRichard Tran Mills #define __FUNCT__ "MatMultTransposeAdd_SeqAIJMKL"
225*ff03dc53SRichard Tran Mills PetscErrorCode MatMultTransposeAdd_SeqAIJMKL(Mat A,Vec xx,Vec yy,Vec zz)
226*ff03dc53SRichard Tran Mills {
227*ff03dc53SRichard Tran Mills   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
228*ff03dc53SRichard Tran Mills   const PetscScalar *x;
229*ff03dc53SRichard Tran Mills   PetscScalar       *y,*z;
230*ff03dc53SRichard Tran Mills   const MatScalar   *aa;
231*ff03dc53SRichard Tran Mills   PetscErrorCode    ierr;
232*ff03dc53SRichard Tran Mills   PetscInt          m=A->rmap->n;
233*ff03dc53SRichard Tran Mills   const PetscInt    *aj,*ai;
234*ff03dc53SRichard Tran Mills   PetscInt          i;
235*ff03dc53SRichard Tran Mills 
236*ff03dc53SRichard Tran Mills   /* Variables not in MatMultTransposeAdd_SeqAIJ. */
237*ff03dc53SRichard Tran Mills   char transa = 't';  /* Used to indicate to MKL that we are computing the transpose product. */
2384a2a386eSRichard Tran Mills 
2394a2a386eSRichard Tran Mills   PetscFunctionBegin;
2404a2a386eSRichard Tran Mills   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
2414a2a386eSRichard Tran Mills   ierr = VecGetArrayPair(yy,zz,&y,&z);CHKERRQ(ierr);
2424a2a386eSRichard Tran Mills   aj   = a->j;  /* aj[k] gives column index for element aa[k]. */
2434a2a386eSRichard Tran Mills   aa   = a->a;  /* Nonzero elements stored row-by-row. */
2444a2a386eSRichard Tran Mills   ai   = a->i;  /* ai[k] is the position in aa and aj where row k starts. */
2454a2a386eSRichard Tran Mills 
2464a2a386eSRichard Tran Mills   /* Call MKL sparse BLAS routine to do the MatMult. */
2474a2a386eSRichard Tran Mills   mkl_cspblas_xcsrgemv(&transa,&m,aa,ai,aj,x,z);
2484a2a386eSRichard Tran Mills 
2494a2a386eSRichard Tran Mills   /* Now add the vector y; MKL sparse BLAS does not have a MatMultAdd equivalent. */
2504a2a386eSRichard Tran Mills   for (i=0; i<m; i++) {
2514a2a386eSRichard Tran Mills     z[i] += y[i];
2524a2a386eSRichard Tran Mills   }
2534a2a386eSRichard Tran Mills 
2544a2a386eSRichard Tran Mills   ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr);
2554a2a386eSRichard Tran Mills   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
2564a2a386eSRichard Tran Mills   ierr = VecRestoreArrayPair(yy,zz,&y,&z);CHKERRQ(ierr);
2574a2a386eSRichard Tran Mills   PetscFunctionReturn(0);
2584a2a386eSRichard Tran Mills }
2594a2a386eSRichard Tran Mills 
2604a2a386eSRichard Tran Mills /* MatConvert_SeqAIJ_SeqAIJMKL converts a SeqAIJ matrix into a
2614a2a386eSRichard Tran Mills  * SeqAIJMKL matrix.  This routine is called by the MatCreate_SeqMKLAIJ()
2624a2a386eSRichard Tran Mills  * routine, but can also be used to convert an assembled SeqAIJ matrix
2634a2a386eSRichard Tran Mills  * into a SeqAIJMKL one. */
2644a2a386eSRichard Tran Mills #undef __FUNCT__
2654a2a386eSRichard Tran Mills #define __FUNCT__ "MatConvert_SeqAIJ_SeqAIJMKL"
2664a2a386eSRichard Tran Mills PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJMKL(Mat A,MatType type,MatReuse reuse,Mat *newmat)
2674a2a386eSRichard Tran Mills {
2684a2a386eSRichard Tran Mills   PetscErrorCode ierr;
2694a2a386eSRichard Tran Mills   Mat            B = *newmat;
2704a2a386eSRichard Tran Mills   Mat_SeqAIJMKL *aijmkl;
2714a2a386eSRichard Tran Mills 
2724a2a386eSRichard Tran Mills   PetscFunctionBegin;
2734a2a386eSRichard Tran Mills   if (reuse == MAT_INITIAL_MATRIX) {
2744a2a386eSRichard Tran Mills     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
2754a2a386eSRichard Tran Mills   }
2764a2a386eSRichard Tran Mills 
2774a2a386eSRichard Tran Mills   ierr     = PetscNewLog(B,&aijmkl);CHKERRQ(ierr);
2784a2a386eSRichard Tran Mills   B->spptr = (void*) aijmkl;
2794a2a386eSRichard Tran Mills 
2804a2a386eSRichard Tran Mills   /* Set function pointers for methods that we inherit from AIJ but override. */
2814a2a386eSRichard Tran Mills   B->ops->duplicate        = MatDuplicate_SeqAIJMKL;
2824a2a386eSRichard Tran Mills   B->ops->assemblyend      = MatAssemblyEnd_SeqAIJMKL;
2834a2a386eSRichard Tran Mills   B->ops->destroy          = MatDestroy_SeqAIJMKL;
2844a2a386eSRichard Tran Mills   B->ops->mult             = MatMult_SeqAIJMKL;
285*ff03dc53SRichard Tran Mills   B->ops->multtranspose    = MatMultTranspose_SeqAIJMKL;
2864a2a386eSRichard Tran Mills   B->ops->multadd          = MatMultAdd_SeqAIJMKL;
287*ff03dc53SRichard Tran Mills   B->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJMKL;
2884a2a386eSRichard Tran Mills 
2894a2a386eSRichard Tran Mills   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaijmkl_seqaij_C",MatConvert_SeqAIJMKL_SeqAIJ);CHKERRQ(ierr);
2904a2a386eSRichard Tran Mills 
2914a2a386eSRichard Tran Mills   ierr    = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJMKL);CHKERRQ(ierr);
2924a2a386eSRichard Tran Mills   *newmat = B;
2934a2a386eSRichard Tran Mills   PetscFunctionReturn(0);
2944a2a386eSRichard Tran Mills }
2954a2a386eSRichard Tran Mills 
2964a2a386eSRichard Tran Mills #undef __FUNCT__
2974a2a386eSRichard Tran Mills #define __FUNCT__ "MatCreateSeqAIJMKL"
2984a2a386eSRichard Tran Mills /*@C
2994a2a386eSRichard Tran Mills    MatCreateSeqAIJMKL - Creates a sparse matrix of type SEQAIJMKL.
3004a2a386eSRichard Tran Mills    This type inherits from AIJ and is largely identical, but uses sparse BLAS
3014a2a386eSRichard Tran Mills    routines from Intel MKL whenever possible.
3024a2a386eSRichard Tran Mills    Collective on MPI_Comm
3034a2a386eSRichard Tran Mills 
3044a2a386eSRichard Tran Mills    Input Parameters:
3054a2a386eSRichard Tran Mills +  comm - MPI communicator, set to PETSC_COMM_SELF
3064a2a386eSRichard Tran Mills .  m - number of rows
3074a2a386eSRichard Tran Mills .  n - number of columns
3084a2a386eSRichard Tran Mills .  nz - number of nonzeros per row (same for all rows)
3094a2a386eSRichard Tran Mills -  nnz - array containing the number of nonzeros in the various rows
3104a2a386eSRichard Tran Mills          (possibly different for each row) or NULL
3114a2a386eSRichard Tran Mills 
3124a2a386eSRichard Tran Mills    Output Parameter:
3134a2a386eSRichard Tran Mills .  A - the matrix
3144a2a386eSRichard Tran Mills 
3154a2a386eSRichard Tran Mills    Notes:
3164a2a386eSRichard Tran Mills    If nnz is given then nz is ignored
3174a2a386eSRichard Tran Mills 
3184a2a386eSRichard Tran Mills    Level: intermediate
3194a2a386eSRichard Tran Mills 
3204a2a386eSRichard Tran Mills .keywords: matrix, cray, sparse, parallel
3214a2a386eSRichard Tran Mills 
3224a2a386eSRichard Tran Mills .seealso: MatCreate(), MatCreateMPIAIJMKL(), MatSetValues()
3234a2a386eSRichard Tran Mills @*/
3244a2a386eSRichard Tran Mills PetscErrorCode  MatCreateSeqAIJMKL(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
3254a2a386eSRichard Tran Mills {
3264a2a386eSRichard Tran Mills   PetscErrorCode ierr;
3274a2a386eSRichard Tran Mills 
3284a2a386eSRichard Tran Mills   PetscFunctionBegin;
3294a2a386eSRichard Tran Mills   ierr = MatCreate(comm,A);CHKERRQ(ierr);
3304a2a386eSRichard Tran Mills   ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr);
3314a2a386eSRichard Tran Mills   ierr = MatSetType(*A,MATSEQAIJMKL);CHKERRQ(ierr);
3324a2a386eSRichard Tran Mills   ierr = MatSeqAIJSetPreallocation_SeqAIJ(*A,nz,nnz);CHKERRQ(ierr);
3334a2a386eSRichard Tran Mills   PetscFunctionReturn(0);
3344a2a386eSRichard Tran Mills }
3354a2a386eSRichard Tran Mills 
3364a2a386eSRichard Tran Mills #undef __FUNCT__
3374a2a386eSRichard Tran Mills #define __FUNCT__ "MatCreate_SeqAIJMKL"
3384a2a386eSRichard Tran Mills PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJMKL(Mat A)
3394a2a386eSRichard Tran Mills {
3404a2a386eSRichard Tran Mills   PetscErrorCode ierr;
3414a2a386eSRichard Tran Mills 
3424a2a386eSRichard Tran Mills   PetscFunctionBegin;
3434a2a386eSRichard Tran Mills   ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
3444a2a386eSRichard Tran Mills   ierr = MatConvert_SeqAIJ_SeqAIJMKL(A,MATSEQAIJMKL,MAT_INPLACE_MATRIX,&A);CHKERRQ(ierr);
3454a2a386eSRichard Tran Mills   PetscFunctionReturn(0);
3464a2a386eSRichard Tran Mills }
3474a2a386eSRichard Tran Mills 
3484a2a386eSRichard Tran Mills 
349