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 #undef __FUNCT__ 274a2a386eSRichard Tran Mills #define __FUNCT__ "MatConvert_SeqAIJMKL_SeqAIJ" 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; 344a2a386eSRichard Tran Mills Mat_SeqAIJMKL *aijmkl=(Mat_SeqAIJMKL*)A->spptr; 354a2a386eSRichard Tran Mills 364a2a386eSRichard Tran Mills PetscFunctionBegin; 374a2a386eSRichard Tran Mills if (reuse == MAT_INITIAL_MATRIX) { 384a2a386eSRichard Tran Mills ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 39e9c94282SRichard Tran Mills aijmkl = (Mat_SeqAIJMKL*)B->spptr; 404a2a386eSRichard Tran Mills } 414a2a386eSRichard Tran Mills 424a2a386eSRichard Tran Mills /* Reset the original function pointers. */ 4354871a98SRichard Tran Mills B->ops->duplicate = MatDuplicate_SeqAIJ; 444a2a386eSRichard Tran Mills B->ops->assemblyend = MatAssemblyEnd_SeqAIJ; 454a2a386eSRichard Tran Mills B->ops->destroy = MatDestroy_SeqAIJ; 4654871a98SRichard Tran Mills B->ops->mult = MatMult_SeqAIJ; 47ff03dc53SRichard Tran Mills B->ops->multtranspose = MatMultTranspose_SeqAIJ; 4854871a98SRichard Tran Mills B->ops->multadd = MatMultAdd_SeqAIJ; 49ff03dc53SRichard Tran Mills B->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJ; 504a2a386eSRichard Tran Mills 51e9c94282SRichard Tran Mills ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaijmkl_seqaij_C",NULL);CHKERRQ(ierr); 52e9c94282SRichard Tran Mills ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMult_seqdense_seqaijmkl_C",NULL);CHKERRQ(ierr); 53e9c94282SRichard Tran Mills ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMultSymbolic_seqdense_seqaijmkl_C",NULL);CHKERRQ(ierr); 54e9c94282SRichard Tran Mills ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMultNumeric_seqdense_seqaijmkl_C",NULL);CHKERRQ(ierr); 55e9c94282SRichard Tran Mills 564abfa3b3SRichard Tran Mills /* Free everything in the Mat_SeqAIJMKL data structure. Currently, this 57e9c94282SRichard Tran Mills * simply involves destroying the MKL sparse matrix handle and then freeing 58e9c94282SRichard Tran Mills * the spptr pointer. */ 594abfa3b3SRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_OPTIMIZE 604abfa3b3SRichard Tran Mills if (aijmkl->sparse_optimized) { 610632b357SRichard Tran Mills sparse_status_t stat; 624abfa3b3SRichard Tran Mills stat = mkl_sparse_destroy(aijmkl->csrA); 634abfa3b3SRichard Tran Mills if (stat != SPARSE_STATUS_SUCCESS) { 644abfa3b3SRichard Tran Mills PetscFunctionReturn(PETSC_ERR_LIB); 654abfa3b3SRichard Tran Mills } 664abfa3b3SRichard Tran Mills } 674abfa3b3SRichard Tran Mills #endif /* PETSC_HAVE_MKL_SPARSE_OPTIMIZE */ 68e9c94282SRichard Tran Mills ierr = PetscFree(B->spptr);CHKERRQ(ierr); 694a2a386eSRichard Tran Mills 704a2a386eSRichard Tran Mills /* Change the type of B to MATSEQAIJ. */ 714a2a386eSRichard Tran Mills ierr = PetscObjectChangeTypeName((PetscObject)B, MATSEQAIJ);CHKERRQ(ierr); 724a2a386eSRichard Tran Mills 734a2a386eSRichard Tran Mills *newmat = B; 744a2a386eSRichard Tran Mills PetscFunctionReturn(0); 754a2a386eSRichard Tran Mills } 764a2a386eSRichard Tran Mills 774a2a386eSRichard Tran Mills #undef __FUNCT__ 784a2a386eSRichard Tran Mills #define __FUNCT__ "MatDestroy_SeqAIJMKL" 794a2a386eSRichard Tran Mills PetscErrorCode MatDestroy_SeqAIJMKL(Mat A) 804a2a386eSRichard Tran Mills { 814a2a386eSRichard Tran Mills PetscErrorCode ierr; 824a2a386eSRichard Tran Mills Mat_SeqAIJMKL *aijmkl = (Mat_SeqAIJMKL*) A->spptr; 834a2a386eSRichard Tran Mills 844a2a386eSRichard Tran Mills PetscFunctionBegin; 85e9c94282SRichard Tran Mills 86e9c94282SRichard Tran Mills /* If MatHeaderMerge() was used, then this SeqAIJMKL matrix will not have an 87e9c94282SRichard Tran Mills * spptr pointer. */ 88e9c94282SRichard Tran Mills if (aijmkl) { 894a2a386eSRichard Tran Mills /* Clean up everything in the Mat_SeqAIJMKL data structure, then free A->spptr. */ 904abfa3b3SRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_OPTIMIZE 914abfa3b3SRichard Tran Mills if (aijmkl->sparse_optimized) { 924abfa3b3SRichard Tran Mills sparse_status_t stat = SPARSE_STATUS_SUCCESS; 934abfa3b3SRichard Tran Mills stat = mkl_sparse_destroy(aijmkl->csrA); 944abfa3b3SRichard Tran Mills if (stat != SPARSE_STATUS_SUCCESS) { 954abfa3b3SRichard Tran Mills PetscFunctionReturn(PETSC_ERR_LIB); 964abfa3b3SRichard Tran Mills } 974abfa3b3SRichard Tran Mills } 984abfa3b3SRichard Tran Mills #endif /* PETSC_HAVE_MKL_SPARSE_OPTIMIZE */ 994a2a386eSRichard Tran Mills ierr = PetscFree(A->spptr);CHKERRQ(ierr); 100e9c94282SRichard Tran Mills } 1014a2a386eSRichard Tran Mills 1024a2a386eSRichard Tran Mills /* Change the type of A back to SEQAIJ and use MatDestroy_SeqAIJ() 1034a2a386eSRichard Tran Mills * to destroy everything that remains. */ 1044a2a386eSRichard Tran Mills ierr = PetscObjectChangeTypeName((PetscObject)A, MATSEQAIJ);CHKERRQ(ierr); 1054a2a386eSRichard Tran Mills /* Note that I don't call MatSetType(). I believe this is because that 1064a2a386eSRichard Tran Mills * is only to be called when *building* a matrix. I could be wrong, but 1074a2a386eSRichard Tran Mills * that is how things work for the SuperLU matrix class. */ 1084a2a386eSRichard Tran Mills ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr); 1094a2a386eSRichard Tran Mills PetscFunctionReturn(0); 1104a2a386eSRichard Tran Mills } 1114a2a386eSRichard Tran Mills 1124a2a386eSRichard Tran Mills #undef __FUNCT__ 1134a2a386eSRichard Tran Mills #define __FUNCT__ "MatDuplicate_SeqAIJMKL" 1144a2a386eSRichard Tran Mills PetscErrorCode MatDuplicate_SeqAIJMKL(Mat A, MatDuplicateOption op, Mat *M) 1154a2a386eSRichard Tran Mills { 1164a2a386eSRichard Tran Mills PetscErrorCode ierr; 1170632b357SRichard Tran Mills Mat_SeqAIJMKL *aijmkl; 1180632b357SRichard Tran Mills Mat_SeqAIJMKL *aijmkl_dest; 1194a2a386eSRichard Tran Mills 1204a2a386eSRichard Tran Mills PetscFunctionBegin; 1214a2a386eSRichard Tran Mills ierr = MatDuplicate_SeqAIJ(A,op,M);CHKERRQ(ierr); 1220632b357SRichard Tran Mills aijmkl = (Mat_SeqAIJMKL*) A->spptr; 1230632b357SRichard Tran Mills aijmkl_dest = (Mat_SeqAIJMKL*) (*M)->spptr; 124a9041576SRichard Tran Mills ierr = PetscMemcpy(aijmkl_dest,aijmkl,sizeof(Mat_SeqAIJMKL));CHKERRQ(ierr); 1250632b357SRichard Tran Mills aijmkl_dest->sparse_optimized = PETSC_FALSE; 1260632b357SRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_OPTIMIZE 1270632b357SRichard Tran Mills aijmkl_dest->csrA = NULL; 128*f36dfe3fSRichard Tran Mills if (!aijmkl->no_SpMV2 && A->cmap->n > 0) { 1290632b357SRichard Tran Mills sparse_status_t stat; 1300632b357SRichard Tran Mills stat = mkl_sparse_copy(aijmkl->csrA,aijmkl->descr,&aijmkl_dest->csrA); 131f68ad4bdSRichard Tran Mills if (stat != SPARSE_STATUS_SUCCESS) { 132f68ad4bdSRichard Tran Mills SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Intel MKL error: unable to complete mkl_sparse_copy"); 133f68ad4bdSRichard Tran Mills PetscFunctionReturn(PETSC_ERR_LIB); 134f68ad4bdSRichard Tran Mills } 1350632b357SRichard Tran Mills stat = mkl_sparse_optimize(aijmkl_dest->csrA); 1360632b357SRichard Tran Mills if (stat != SPARSE_STATUS_SUCCESS) { 137f68ad4bdSRichard Tran Mills SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Intel MKL error: unable to complete mkl_sparse_optimize"); 1380632b357SRichard Tran Mills PetscFunctionReturn(PETSC_ERR_LIB); 1390632b357SRichard Tran Mills } 1400632b357SRichard Tran Mills aijmkl_dest->sparse_optimized = PETSC_TRUE; 1410632b357SRichard Tran Mills } 1420632b357SRichard Tran Mills #endif /* PETSC_HAVE_MKL_SPARSE_OPTIMIZE */ 1434a2a386eSRichard Tran Mills PetscFunctionReturn(0); 1444a2a386eSRichard Tran Mills } 1454a2a386eSRichard Tran Mills 1464a2a386eSRichard Tran Mills #undef __FUNCT__ 1474a2a386eSRichard Tran Mills #define __FUNCT__ "MatAssemblyEnd_SeqAIJMKL" 1484a2a386eSRichard Tran Mills PetscErrorCode MatAssemblyEnd_SeqAIJMKL(Mat A, MatAssemblyType mode) 1494a2a386eSRichard Tran Mills { 1504a2a386eSRichard Tran Mills PetscErrorCode ierr; 1514a2a386eSRichard Tran Mills Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 152df555b71SRichard Tran Mills Mat_SeqAIJMKL *aijmkl; 153df555b71SRichard Tran Mills 154df555b71SRichard Tran Mills MatScalar *aa; 15558678438SRichard Tran Mills PetscInt m,n; 156df555b71SRichard Tran Mills PetscInt *aj,*ai; 1574a2a386eSRichard Tran Mills 1584a2a386eSRichard Tran Mills PetscFunctionBegin; 1594a2a386eSRichard Tran Mills if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 1604a2a386eSRichard Tran Mills 1614a2a386eSRichard Tran Mills /* Since a MATSEQAIJMKL matrix is really just a MATSEQAIJ with some 1624a2a386eSRichard Tran Mills * extra information and some different methods, call the AssemblyEnd 1634a2a386eSRichard Tran Mills * routine for a MATSEQAIJ. 1644a2a386eSRichard Tran Mills * I'm not sure if this is the best way to do this, but it avoids 1654a2a386eSRichard Tran Mills * a lot of code duplication. 1664a2a386eSRichard Tran Mills * I also note that currently MATSEQAIJMKL doesn't know anything about 1674a2a386eSRichard Tran Mills * the Mat_CompressedRow data structure that SeqAIJ now uses when there 1684a2a386eSRichard Tran Mills * are many zero rows. If the SeqAIJ assembly end routine decides to use 1694a2a386eSRichard Tran Mills * this, this may break things. (Don't know... haven't looked at it. 1704a2a386eSRichard Tran Mills * Do I need to disable this somehow?) */ 1714a2a386eSRichard Tran Mills a->inode.use = PETSC_FALSE; /* Must disable: otherwise the MKL routines won't get used. */ 1724a2a386eSRichard Tran Mills ierr = MatAssemblyEnd_SeqAIJ(A, mode);CHKERRQ(ierr); 1734a2a386eSRichard Tran Mills 174df555b71SRichard Tran Mills aijmkl = (Mat_SeqAIJMKL*) A->spptr; 175d995685eSRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_OPTIMIZE 176c9d46305SRichard Tran Mills if (!aijmkl->no_SpMV2) { 1770632b357SRichard Tran Mills sparse_status_t stat; 1780632b357SRichard Tran Mills if (aijmkl->sparse_optimized) { 1790632b357SRichard Tran Mills /* Matrix has been previously assembled and optimized. Must destroy old 1800632b357SRichard Tran Mills * matrix handle before running the optimization step again. */ 1810632b357SRichard Tran Mills sparse_status_t stat; 1820632b357SRichard Tran Mills stat = mkl_sparse_destroy(aijmkl->csrA); 1830632b357SRichard Tran Mills if (stat != SPARSE_STATUS_SUCCESS) { 1840632b357SRichard Tran Mills PetscFunctionReturn(PETSC_ERR_LIB); 1850632b357SRichard Tran Mills } 1860632b357SRichard Tran Mills } 187c9d46305SRichard Tran Mills /* Now perform the SpMV2 setup and matrix optimization. */ 188df555b71SRichard Tran Mills aijmkl->descr.type = SPARSE_MATRIX_TYPE_GENERAL; 189df555b71SRichard Tran Mills aijmkl->descr.mode = SPARSE_FILL_MODE_LOWER; 190df555b71SRichard Tran Mills aijmkl->descr.diag = SPARSE_DIAG_NON_UNIT; 19158678438SRichard Tran Mills m = A->rmap->n; 19258678438SRichard Tran Mills n = A->cmap->n; 193df555b71SRichard Tran Mills aj = a->j; /* aj[k] gives column index for element aa[k]. */ 194df555b71SRichard Tran Mills aa = a->a; /* Nonzero elements stored row-by-row. */ 195df555b71SRichard Tran Mills ai = a->i; /* ai[k] is the position in aa and aj where row k starts. */ 196*f36dfe3fSRichard Tran Mills if (n>0) { 19758678438SRichard Tran Mills stat = mkl_sparse_x_create_csr (&aijmkl->csrA,SPARSE_INDEX_BASE_ZERO,m,n,ai,ai+1,aj,aa); 198df555b71SRichard Tran Mills stat = mkl_sparse_set_mv_hint(aijmkl->csrA,SPARSE_OPERATION_NON_TRANSPOSE,aijmkl->descr,1000); 199df555b71SRichard Tran Mills stat = mkl_sparse_set_memory_hint(aijmkl->csrA,SPARSE_MEMORY_AGGRESSIVE); 200df555b71SRichard Tran Mills stat = mkl_sparse_optimize(aijmkl->csrA); 201df555b71SRichard Tran Mills if (stat != SPARSE_STATUS_SUCCESS) { 202f68ad4bdSRichard Tran Mills SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Intel MKL error: unable to create matrix handle/complete mkl_sparse_optimize"); 203df555b71SRichard Tran Mills PetscFunctionReturn(PETSC_ERR_LIB); 204df555b71SRichard Tran Mills } 2054abfa3b3SRichard Tran Mills aijmkl->sparse_optimized = PETSC_TRUE; 206c9d46305SRichard Tran Mills } 207*f36dfe3fSRichard Tran Mills } 208d995685eSRichard Tran Mills #endif 209df555b71SRichard Tran Mills 2104a2a386eSRichard Tran Mills PetscFunctionReturn(0); 2114a2a386eSRichard Tran Mills } 2124a2a386eSRichard Tran Mills 2134a2a386eSRichard Tran Mills #undef __FUNCT__ 2144a2a386eSRichard Tran Mills #define __FUNCT__ "MatMult_SeqAIJMKL" 2154a2a386eSRichard Tran Mills PetscErrorCode MatMult_SeqAIJMKL(Mat A,Vec xx,Vec yy) 2164a2a386eSRichard Tran Mills { 2174a2a386eSRichard Tran Mills Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 2184a2a386eSRichard Tran Mills const PetscScalar *x; 2194a2a386eSRichard Tran Mills PetscScalar *y; 2204a2a386eSRichard Tran Mills const MatScalar *aa; 2214a2a386eSRichard Tran Mills PetscErrorCode ierr; 2224a2a386eSRichard Tran Mills PetscInt m=A->rmap->n; 2234a2a386eSRichard Tran Mills const PetscInt *aj,*ai; 2244a2a386eSRichard Tran Mills 2254a2a386eSRichard Tran Mills /* Variables not in MatMult_SeqAIJ. */ 226ff03dc53SRichard Tran Mills char transa = 'n'; /* Used to indicate to MKL that we are not computing the transpose product. */ 227ff03dc53SRichard Tran Mills 228ff03dc53SRichard Tran Mills PetscFunctionBegin; 229ff03dc53SRichard Tran Mills ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 230ff03dc53SRichard Tran Mills ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 231ff03dc53SRichard Tran Mills aj = a->j; /* aj[k] gives column index for element aa[k]. */ 232ff03dc53SRichard Tran Mills aa = a->a; /* Nonzero elements stored row-by-row. */ 233ff03dc53SRichard Tran Mills ai = a->i; /* ai[k] is the position in aa and aj where row k starts. */ 234ff03dc53SRichard Tran Mills 235ff03dc53SRichard Tran Mills /* Call MKL sparse BLAS routine to do the MatMult. */ 236ff03dc53SRichard Tran Mills mkl_cspblas_xcsrgemv(&transa,&m,aa,ai,aj,x,y); 237ff03dc53SRichard Tran Mills 238ff03dc53SRichard Tran Mills ierr = PetscLogFlops(2.0*a->nz - a->nonzerorowcnt);CHKERRQ(ierr); 239ff03dc53SRichard Tran Mills ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 240ff03dc53SRichard Tran Mills ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 241ff03dc53SRichard Tran Mills PetscFunctionReturn(0); 242ff03dc53SRichard Tran Mills } 243ff03dc53SRichard Tran Mills 244d995685eSRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_OPTIMIZE 245ff03dc53SRichard Tran Mills #undef __FUNCT__ 246df555b71SRichard Tran Mills #define __FUNCT__ "MatMult_SeqAIJMKL_SpMV2" 247df555b71SRichard Tran Mills PetscErrorCode MatMult_SeqAIJMKL_SpMV2(Mat A,Vec xx,Vec yy) 248df555b71SRichard Tran Mills { 249df555b71SRichard Tran Mills Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 250df555b71SRichard Tran Mills Mat_SeqAIJMKL *aijmkl=(Mat_SeqAIJMKL*)A->spptr; 251df555b71SRichard Tran Mills const PetscScalar *x; 252df555b71SRichard Tran Mills PetscScalar *y; 253df555b71SRichard Tran Mills PetscErrorCode ierr; 254df555b71SRichard Tran Mills sparse_status_t stat = SPARSE_STATUS_SUCCESS; 255df555b71SRichard Tran Mills 256df555b71SRichard Tran Mills PetscFunctionBegin; 257df555b71SRichard Tran Mills 258*f36dfe3fSRichard Tran Mills /* If there are no rows, this is a no-op: return immediately. */ 259*f36dfe3fSRichard Tran Mills if(A->cmap->n < 1) PetscFunctionReturn(0); 260*f36dfe3fSRichard Tran Mills 261df555b71SRichard Tran Mills ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 262df555b71SRichard Tran Mills ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 263df555b71SRichard Tran Mills 264df555b71SRichard Tran Mills /* Call MKL SpMV2 executor routine to do the MatMult. */ 265df555b71SRichard Tran Mills stat = mkl_sparse_x_mv(SPARSE_OPERATION_NON_TRANSPOSE,1.0,aijmkl->csrA,aijmkl->descr,x,0.0,y); 266df555b71SRichard Tran Mills 267df555b71SRichard Tran Mills ierr = PetscLogFlops(2.0*a->nz - a->nonzerorowcnt);CHKERRQ(ierr); 268df555b71SRichard Tran Mills ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 269df555b71SRichard Tran Mills ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 270df555b71SRichard Tran Mills if (stat != SPARSE_STATUS_SUCCESS) { 271df555b71SRichard Tran Mills PetscFunctionReturn(PETSC_ERR_LIB); 272df555b71SRichard Tran Mills } 273df555b71SRichard Tran Mills PetscFunctionReturn(0); 274df555b71SRichard Tran Mills } 275d995685eSRichard Tran Mills #endif /* PETSC_HAVE_MKL_SPARSE_OPTIMIZE */ 276df555b71SRichard Tran Mills 277df555b71SRichard Tran Mills #undef __FUNCT__ 278ff03dc53SRichard Tran Mills #define __FUNCT__ "MatMultTranspose_SeqAIJMKL" 279ff03dc53SRichard Tran Mills PetscErrorCode MatMultTranspose_SeqAIJMKL(Mat A,Vec xx,Vec yy) 280ff03dc53SRichard Tran Mills { 281ff03dc53SRichard Tran Mills Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 282ff03dc53SRichard Tran Mills const PetscScalar *x; 283ff03dc53SRichard Tran Mills PetscScalar *y; 284ff03dc53SRichard Tran Mills const MatScalar *aa; 285ff03dc53SRichard Tran Mills PetscErrorCode ierr; 286ff03dc53SRichard Tran Mills PetscInt m=A->rmap->n; 287ff03dc53SRichard Tran Mills const PetscInt *aj,*ai; 288ff03dc53SRichard Tran Mills 289ff03dc53SRichard Tran Mills /* Variables not in MatMultTranspose_SeqAIJ. */ 290ff03dc53SRichard Tran Mills char transa = 't'; /* Used to indicate to MKL that we are computing the transpose product. */ 2914a2a386eSRichard Tran Mills 2924a2a386eSRichard Tran Mills PetscFunctionBegin; 2934a2a386eSRichard Tran Mills ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 2944a2a386eSRichard Tran Mills ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 2954a2a386eSRichard Tran Mills aj = a->j; /* aj[k] gives column index for element aa[k]. */ 2964a2a386eSRichard Tran Mills aa = a->a; /* Nonzero elements stored row-by-row. */ 2974a2a386eSRichard Tran Mills ai = a->i; /* ai[k] is the position in aa and aj where row k starts. */ 2984a2a386eSRichard Tran Mills 2994a2a386eSRichard Tran Mills /* Call MKL sparse BLAS routine to do the MatMult. */ 3004a2a386eSRichard Tran Mills mkl_cspblas_xcsrgemv(&transa,&m,aa,ai,aj,x,y); 3014a2a386eSRichard Tran Mills 3024a2a386eSRichard Tran Mills ierr = PetscLogFlops(2.0*a->nz - a->nonzerorowcnt);CHKERRQ(ierr); 3034a2a386eSRichard Tran Mills ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 3044a2a386eSRichard Tran Mills ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 3054a2a386eSRichard Tran Mills PetscFunctionReturn(0); 3064a2a386eSRichard Tran Mills } 3074a2a386eSRichard Tran Mills 308d995685eSRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_OPTIMIZE 3094a2a386eSRichard Tran Mills #undef __FUNCT__ 310df555b71SRichard Tran Mills #define __FUNCT__ "MatMultTranspose_SeqAIJMKL_SpMV2" 311df555b71SRichard Tran Mills PetscErrorCode MatMultTranspose_SeqAIJMKL_SpMV2(Mat A,Vec xx,Vec yy) 312df555b71SRichard Tran Mills { 313df555b71SRichard Tran Mills Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 314df555b71SRichard Tran Mills Mat_SeqAIJMKL *aijmkl=(Mat_SeqAIJMKL*)A->spptr; 315df555b71SRichard Tran Mills const PetscScalar *x; 316df555b71SRichard Tran Mills PetscScalar *y; 317df555b71SRichard Tran Mills PetscErrorCode ierr; 3180632b357SRichard Tran Mills sparse_status_t stat; 319df555b71SRichard Tran Mills 320df555b71SRichard Tran Mills PetscFunctionBegin; 321df555b71SRichard Tran Mills 322*f36dfe3fSRichard Tran Mills /* If there are no rows, this is a no-op: return immediately. */ 323*f36dfe3fSRichard Tran Mills if(A->cmap->n < 1) PetscFunctionReturn(0); 324*f36dfe3fSRichard Tran Mills 325df555b71SRichard Tran Mills ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 326df555b71SRichard Tran Mills ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 327df555b71SRichard Tran Mills 328df555b71SRichard Tran Mills /* Call MKL SpMV2 executor routine to do the MatMultTranspose. */ 329df555b71SRichard Tran Mills stat = mkl_sparse_x_mv(SPARSE_OPERATION_TRANSPOSE,1.0,aijmkl->csrA,aijmkl->descr,x,0.0,y); 330df555b71SRichard Tran Mills 331df555b71SRichard Tran Mills ierr = PetscLogFlops(2.0*a->nz - a->nonzerorowcnt);CHKERRQ(ierr); 332df555b71SRichard Tran Mills ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 333df555b71SRichard Tran Mills ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 334df555b71SRichard Tran Mills if (stat != SPARSE_STATUS_SUCCESS) { 335df555b71SRichard Tran Mills PetscFunctionReturn(PETSC_ERR_LIB); 336df555b71SRichard Tran Mills } 337df555b71SRichard Tran Mills PetscFunctionReturn(0); 338df555b71SRichard Tran Mills } 339d995685eSRichard Tran Mills #endif /* PETSC_HAVE_MKL_SPARSE_OPTIMIZE */ 340df555b71SRichard Tran Mills 341df555b71SRichard Tran Mills #undef __FUNCT__ 3424a2a386eSRichard Tran Mills #define __FUNCT__ "MatMultAdd_SeqAIJMKL" 3434a2a386eSRichard Tran Mills PetscErrorCode MatMultAdd_SeqAIJMKL(Mat A,Vec xx,Vec yy,Vec zz) 3444a2a386eSRichard Tran Mills { 3454a2a386eSRichard Tran Mills Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 3464a2a386eSRichard Tran Mills const PetscScalar *x; 3474a2a386eSRichard Tran Mills PetscScalar *y,*z; 3484a2a386eSRichard Tran Mills const MatScalar *aa; 3494a2a386eSRichard Tran Mills PetscErrorCode ierr; 3504a2a386eSRichard Tran Mills PetscInt m=A->rmap->n; 3514a2a386eSRichard Tran Mills const PetscInt *aj,*ai; 3524a2a386eSRichard Tran Mills PetscInt i; 3534a2a386eSRichard Tran Mills 354ff03dc53SRichard Tran Mills /* Variables not in MatMultAdd_SeqAIJ. */ 355ff03dc53SRichard Tran Mills char transa = 'n'; /* Used to indicate to MKL that we are not computing the transpose product. */ 356a84739b8SRichard Tran Mills PetscScalar alpha = 1.0; 357a84739b8SRichard Tran Mills PetscScalar beta = 1.0; 358a84739b8SRichard Tran Mills char matdescra[6]; 359ff03dc53SRichard Tran Mills 360ff03dc53SRichard Tran Mills PetscFunctionBegin; 361a84739b8SRichard Tran Mills matdescra[0] = 'g'; /* Indicates to MKL that we using a general CSR matrix. */ 362a84739b8SRichard Tran Mills matdescra[3] = 'c'; /* Indicates to MKL that we use C-style (0-based) indexing. */ 363a84739b8SRichard Tran Mills 364ff03dc53SRichard Tran Mills ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 365ff03dc53SRichard Tran Mills ierr = VecGetArrayPair(yy,zz,&y,&z);CHKERRQ(ierr); 366ff03dc53SRichard Tran Mills aj = a->j; /* aj[k] gives column index for element aa[k]. */ 367ff03dc53SRichard Tran Mills aa = a->a; /* Nonzero elements stored row-by-row. */ 368ff03dc53SRichard Tran Mills ai = a->i; /* ai[k] is the position in aa and aj where row k starts. */ 369ff03dc53SRichard Tran Mills 370ff03dc53SRichard Tran Mills /* Call MKL sparse BLAS routine to do the MatMult. */ 371a84739b8SRichard Tran Mills if (zz == yy) { 372a84739b8SRichard 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. */ 373a84739b8SRichard Tran Mills mkl_xcsrmv(&transa,&m,&m,&alpha,matdescra,aa,aj,ai,ai+1,x,&beta,y); 374a84739b8SRichard Tran Mills } else { 375a84739b8SRichard Tran Mills /* zz and yy are different vectors, so we call mkl_cspblas_xcsrgemv(), which calculates y = A*x, and then 376a84739b8SRichard Tran Mills * we add the contents of vector yy to the result; MKL sparse BLAS does not have a MatMultAdd equivalent. */ 377ff03dc53SRichard Tran Mills mkl_cspblas_xcsrgemv(&transa,&m,aa,ai,aj,x,z); 378ff03dc53SRichard Tran Mills for (i=0; i<m; i++) { 379ff03dc53SRichard Tran Mills z[i] += y[i]; 380ff03dc53SRichard Tran Mills } 381a84739b8SRichard Tran Mills } 382ff03dc53SRichard Tran Mills 383ff03dc53SRichard Tran Mills ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr); 384ff03dc53SRichard Tran Mills ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 385ff03dc53SRichard Tran Mills ierr = VecRestoreArrayPair(yy,zz,&y,&z);CHKERRQ(ierr); 386ff03dc53SRichard Tran Mills PetscFunctionReturn(0); 387ff03dc53SRichard Tran Mills } 388ff03dc53SRichard Tran Mills 389d995685eSRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_OPTIMIZE 390ff03dc53SRichard Tran Mills #undef __FUNCT__ 391df555b71SRichard Tran Mills #define __FUNCT__ "MatMultAdd_SeqAIJMKL_SpMV2" 392df555b71SRichard Tran Mills PetscErrorCode MatMultAdd_SeqAIJMKL_SpMV2(Mat A,Vec xx,Vec yy,Vec zz) 393df555b71SRichard Tran Mills { 394df555b71SRichard Tran Mills Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 395df555b71SRichard Tran Mills Mat_SeqAIJMKL *aijmkl=(Mat_SeqAIJMKL*)A->spptr; 396df555b71SRichard Tran Mills const PetscScalar *x; 397df555b71SRichard Tran Mills PetscScalar *y,*z; 398df555b71SRichard Tran Mills PetscErrorCode ierr; 399df555b71SRichard Tran Mills PetscInt m=A->rmap->n; 400df555b71SRichard Tran Mills PetscInt i; 401df555b71SRichard Tran Mills 402df555b71SRichard Tran Mills /* Variables not in MatMultAdd_SeqAIJ. */ 403df555b71SRichard Tran Mills sparse_status_t stat = SPARSE_STATUS_SUCCESS; 404df555b71SRichard Tran Mills 405df555b71SRichard Tran Mills PetscFunctionBegin; 406df555b71SRichard Tran Mills 407*f36dfe3fSRichard Tran Mills /* If there are no rows, this is a no-op: return immediately. */ 408*f36dfe3fSRichard Tran Mills if(A->cmap->n < 1) PetscFunctionReturn(0); 409df555b71SRichard Tran Mills 410df555b71SRichard Tran Mills ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 411df555b71SRichard Tran Mills ierr = VecGetArrayPair(yy,zz,&y,&z);CHKERRQ(ierr); 412df555b71SRichard Tran Mills 413df555b71SRichard Tran Mills /* Call MKL sparse BLAS routine to do the MatMult. */ 414df555b71SRichard Tran Mills if (zz == yy) { 415df555b71SRichard 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, 416df555b71SRichard Tran Mills * with alpha and beta both set to 1.0. */ 417df555b71SRichard Tran Mills stat = mkl_sparse_x_mv(SPARSE_OPERATION_NON_TRANSPOSE,1.0,aijmkl->csrA,aijmkl->descr,x,1.0,y); 418df555b71SRichard Tran Mills } else { 419df555b71SRichard 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 420df555b71SRichard Tran Mills * we add the contents of vector yy to the result; MKL sparse BLAS does not have a MatMultAdd equivalent. */ 421df555b71SRichard Tran Mills stat = mkl_sparse_x_mv(SPARSE_OPERATION_NON_TRANSPOSE,1.0,aijmkl->csrA,aijmkl->descr,x,0.0,y); 422df555b71SRichard Tran Mills for (i=0; i<m; i++) { 423df555b71SRichard Tran Mills z[i] += y[i]; 424df555b71SRichard Tran Mills } 425df555b71SRichard Tran Mills } 426df555b71SRichard Tran Mills 427df555b71SRichard Tran Mills ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr); 428df555b71SRichard Tran Mills ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 429df555b71SRichard Tran Mills ierr = VecRestoreArrayPair(yy,zz,&y,&z);CHKERRQ(ierr); 430df555b71SRichard Tran Mills if (stat != SPARSE_STATUS_SUCCESS) { 431df555b71SRichard Tran Mills PetscFunctionReturn(PETSC_ERR_LIB); 432df555b71SRichard Tran Mills } 433df555b71SRichard Tran Mills PetscFunctionReturn(0); 434df555b71SRichard Tran Mills } 435d995685eSRichard Tran Mills #endif /* PETSC_HAVE_MKL_SPARSE_OPTIMIZE */ 436df555b71SRichard Tran Mills 437df555b71SRichard Tran Mills #undef __FUNCT__ 438ff03dc53SRichard Tran Mills #define __FUNCT__ "MatMultTransposeAdd_SeqAIJMKL" 439ff03dc53SRichard Tran Mills PetscErrorCode MatMultTransposeAdd_SeqAIJMKL(Mat A,Vec xx,Vec yy,Vec zz) 440ff03dc53SRichard Tran Mills { 441ff03dc53SRichard Tran Mills Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 442ff03dc53SRichard Tran Mills const PetscScalar *x; 443ff03dc53SRichard Tran Mills PetscScalar *y,*z; 444ff03dc53SRichard Tran Mills const MatScalar *aa; 445ff03dc53SRichard Tran Mills PetscErrorCode ierr; 446ff03dc53SRichard Tran Mills PetscInt m=A->rmap->n; 447ff03dc53SRichard Tran Mills const PetscInt *aj,*ai; 448ff03dc53SRichard Tran Mills PetscInt i; 449ff03dc53SRichard Tran Mills 450ff03dc53SRichard Tran Mills /* Variables not in MatMultTransposeAdd_SeqAIJ. */ 451ff03dc53SRichard Tran Mills char transa = 't'; /* Used to indicate to MKL that we are computing the transpose product. */ 452a84739b8SRichard Tran Mills PetscScalar alpha = 1.0; 453a84739b8SRichard Tran Mills PetscScalar beta = 1.0; 454a84739b8SRichard Tran Mills char matdescra[6]; 4554a2a386eSRichard Tran Mills 4564a2a386eSRichard Tran Mills PetscFunctionBegin; 457a84739b8SRichard Tran Mills matdescra[0] = 'g'; /* Indicates to MKL that we using a general CSR matrix. */ 458a84739b8SRichard Tran Mills matdescra[3] = 'c'; /* Indicates to MKL that we use C-style (0-based) indexing. */ 459a84739b8SRichard Tran Mills 4604a2a386eSRichard Tran Mills ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 4614a2a386eSRichard Tran Mills ierr = VecGetArrayPair(yy,zz,&y,&z);CHKERRQ(ierr); 4624a2a386eSRichard Tran Mills aj = a->j; /* aj[k] gives column index for element aa[k]. */ 4634a2a386eSRichard Tran Mills aa = a->a; /* Nonzero elements stored row-by-row. */ 4644a2a386eSRichard Tran Mills ai = a->i; /* ai[k] is the position in aa and aj where row k starts. */ 4654a2a386eSRichard Tran Mills 4664a2a386eSRichard Tran Mills /* Call MKL sparse BLAS routine to do the MatMult. */ 467a84739b8SRichard Tran Mills if (zz == yy) { 468a84739b8SRichard 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. */ 469a84739b8SRichard Tran Mills mkl_xcsrmv(&transa,&m,&m,&alpha,matdescra,aa,aj,ai,ai+1,x,&beta,y); 470a84739b8SRichard Tran Mills } else { 471a84739b8SRichard Tran Mills /* zz and yy are different vectors, so we call mkl_cspblas_xcsrgemv(), which calculates y = A*x, and then 472a84739b8SRichard Tran Mills * we add the contents of vector yy to the result; MKL sparse BLAS does not have a MatMultAdd equivalent. */ 4734a2a386eSRichard Tran Mills mkl_cspblas_xcsrgemv(&transa,&m,aa,ai,aj,x,z); 4744a2a386eSRichard Tran Mills for (i=0; i<m; i++) { 4754a2a386eSRichard Tran Mills z[i] += y[i]; 4764a2a386eSRichard Tran Mills } 477a84739b8SRichard Tran Mills } 4784a2a386eSRichard Tran Mills 4794a2a386eSRichard Tran Mills ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr); 4804a2a386eSRichard Tran Mills ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 4814a2a386eSRichard Tran Mills ierr = VecRestoreArrayPair(yy,zz,&y,&z);CHKERRQ(ierr); 4824a2a386eSRichard Tran Mills PetscFunctionReturn(0); 4834a2a386eSRichard Tran Mills } 4844a2a386eSRichard Tran Mills 485d995685eSRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_OPTIMIZE 486df555b71SRichard Tran Mills #undef __FUNCT__ 487df555b71SRichard Tran Mills #define __FUNCT__ "MatMultTransposeAdd_SeqAIJMKL_SpMV2" 488df555b71SRichard Tran Mills PetscErrorCode MatMultTransposeAdd_SeqAIJMKL_SpMV2(Mat A,Vec xx,Vec yy,Vec zz) 489df555b71SRichard Tran Mills { 490df555b71SRichard Tran Mills Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 491df555b71SRichard Tran Mills Mat_SeqAIJMKL *aijmkl=(Mat_SeqAIJMKL*)A->spptr; 492df555b71SRichard Tran Mills const PetscScalar *x; 493df555b71SRichard Tran Mills PetscScalar *y,*z; 494df555b71SRichard Tran Mills PetscErrorCode ierr; 495df555b71SRichard Tran Mills PetscInt m=A->rmap->n; 496df555b71SRichard Tran Mills PetscInt i; 497df555b71SRichard Tran Mills 498df555b71SRichard Tran Mills /* Variables not in MatMultTransposeAdd_SeqAIJ. */ 499df555b71SRichard Tran Mills sparse_status_t stat = SPARSE_STATUS_SUCCESS; 500df555b71SRichard Tran Mills 501df555b71SRichard Tran Mills PetscFunctionBegin; 502df555b71SRichard Tran Mills 503*f36dfe3fSRichard Tran Mills /* If there are no rows, this is a no-op: return immediately. */ 504*f36dfe3fSRichard Tran Mills if(A->cmap->n < 1) PetscFunctionReturn(0); 505*f36dfe3fSRichard Tran Mills 506df555b71SRichard Tran Mills ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 507df555b71SRichard Tran Mills ierr = VecGetArrayPair(yy,zz,&y,&z);CHKERRQ(ierr); 508df555b71SRichard Tran Mills 509df555b71SRichard Tran Mills /* Call MKL sparse BLAS routine to do the MatMult. */ 510df555b71SRichard Tran Mills if (zz == yy) { 511df555b71SRichard 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, 512df555b71SRichard Tran Mills * with alpha and beta both set to 1.0. */ 513df555b71SRichard Tran Mills stat = mkl_sparse_x_mv(SPARSE_OPERATION_TRANSPOSE,1.0,aijmkl->csrA,aijmkl->descr,x,1.0,y); 514df555b71SRichard Tran Mills } else { 515df555b71SRichard 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 516df555b71SRichard Tran Mills * we add the contents of vector yy to the result; MKL sparse BLAS does not have a MatMultAdd equivalent. */ 517df555b71SRichard Tran Mills stat = mkl_sparse_x_mv(SPARSE_OPERATION_TRANSPOSE,1.0,aijmkl->csrA,aijmkl->descr,x,0.0,y); 518df555b71SRichard Tran Mills for (i=0; i<m; i++) { 519df555b71SRichard Tran Mills z[i] += y[i]; 520df555b71SRichard Tran Mills } 521df555b71SRichard Tran Mills } 522df555b71SRichard Tran Mills 523df555b71SRichard Tran Mills ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr); 524df555b71SRichard Tran Mills ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 525df555b71SRichard Tran Mills ierr = VecRestoreArrayPair(yy,zz,&y,&z);CHKERRQ(ierr); 526df555b71SRichard Tran Mills if (stat != SPARSE_STATUS_SUCCESS) { 527df555b71SRichard Tran Mills PetscFunctionReturn(PETSC_ERR_LIB); 528df555b71SRichard Tran Mills } 529df555b71SRichard Tran Mills PetscFunctionReturn(0); 530df555b71SRichard Tran Mills } 531d995685eSRichard Tran Mills #endif /* PETSC_HAVE_MKL_SPARSE_OPTIMIZE */ 532df555b71SRichard Tran Mills 533df555b71SRichard Tran Mills 5344a2a386eSRichard Tran Mills /* MatConvert_SeqAIJ_SeqAIJMKL converts a SeqAIJ matrix into a 5354a2a386eSRichard Tran Mills * SeqAIJMKL matrix. This routine is called by the MatCreate_SeqMKLAIJ() 5364a2a386eSRichard Tran Mills * routine, but can also be used to convert an assembled SeqAIJ matrix 5374a2a386eSRichard Tran Mills * into a SeqAIJMKL one. */ 5384a2a386eSRichard Tran Mills #undef __FUNCT__ 5394a2a386eSRichard Tran Mills #define __FUNCT__ "MatConvert_SeqAIJ_SeqAIJMKL" 5404a2a386eSRichard Tran Mills PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJMKL(Mat A,MatType type,MatReuse reuse,Mat *newmat) 5414a2a386eSRichard Tran Mills { 5424a2a386eSRichard Tran Mills PetscErrorCode ierr; 5434a2a386eSRichard Tran Mills Mat B = *newmat; 5444a2a386eSRichard Tran Mills Mat_SeqAIJMKL *aijmkl; 545c9d46305SRichard Tran Mills PetscBool set; 546e9c94282SRichard Tran Mills PetscBool sametype; 5474a2a386eSRichard Tran Mills 5484a2a386eSRichard Tran Mills PetscFunctionBegin; 5494a2a386eSRichard Tran Mills if (reuse == MAT_INITIAL_MATRIX) { 5504a2a386eSRichard Tran Mills ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 5514a2a386eSRichard Tran Mills } 5524a2a386eSRichard Tran Mills 553e9c94282SRichard Tran Mills ierr = PetscObjectTypeCompare((PetscObject)A,type,&sametype);CHKERRQ(ierr); 554e9c94282SRichard Tran Mills if (sametype) PetscFunctionReturn(0); 555e9c94282SRichard Tran Mills 5564a2a386eSRichard Tran Mills ierr = PetscNewLog(B,&aijmkl);CHKERRQ(ierr); 5574a2a386eSRichard Tran Mills B->spptr = (void*) aijmkl; 5584a2a386eSRichard Tran Mills 559df555b71SRichard Tran Mills /* Set function pointers for methods that we inherit from AIJ but override. 560e9c94282SRichard Tran Mills * We also parse some command line options below, since those determine some of the methods we point to. 561e9c94282SRichard Tran Mills * Note: Currently the transposed operations are not being set because I encounter memory corruption 562df555b71SRichard Tran Mills * when these are enabled. Need to look at this with Valgrind or similar. --RTM */ 5634a2a386eSRichard Tran Mills B->ops->duplicate = MatDuplicate_SeqAIJMKL; 5644a2a386eSRichard Tran Mills B->ops->assemblyend = MatAssemblyEnd_SeqAIJMKL; 5654a2a386eSRichard Tran Mills B->ops->destroy = MatDestroy_SeqAIJMKL; 566c9d46305SRichard Tran Mills 5674abfa3b3SRichard Tran Mills aijmkl->sparse_optimized = PETSC_FALSE; 568d995685eSRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_OPTIMIZE 569d995685eSRichard Tran Mills aijmkl->no_SpMV2 = PETSC_FALSE; /* Default to using the SpMV2 routines if our MKL supports them. */ 570d995685eSRichard Tran Mills #elif 571d995685eSRichard Tran Mills aijmkl->no_SpMV2 = PETSC_TRUE; 572d995685eSRichard Tran Mills #endif 5734abfa3b3SRichard Tran Mills 5744abfa3b3SRichard Tran Mills /* Parse command line options. */ 575c9d46305SRichard Tran Mills ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"AIJMKL Options","Mat");CHKERRQ(ierr); 576c9d46305SRichard Tran Mills ierr = PetscOptionsBool("-mat_aijmkl_no_spmv2","NoSPMV2","None",(PetscBool)aijmkl->no_SpMV2,(PetscBool*)&aijmkl->no_SpMV2,&set);CHKERRQ(ierr); 577c9d46305SRichard Tran Mills ierr = PetscOptionsEnd();CHKERRQ(ierr); 578d995685eSRichard Tran Mills #ifndef PETSC_HAVE_MKL_SPARSE_OPTIMIZE 579d995685eSRichard Tran Mills if(!aijmkl->no_SpMV2) { 580d995685eSRichard 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"); 581d995685eSRichard Tran Mills aijmkl->no_SpMV2 = PETSC_TRUE; 582d995685eSRichard Tran Mills } 583d995685eSRichard Tran Mills #endif 584c9d46305SRichard Tran Mills 585c9d46305SRichard Tran Mills if(!aijmkl->no_SpMV2) { 586d995685eSRichard Tran Mills #ifdef PETSC_HAVE_MKL_SPARSE_OPTIMIZE 587df555b71SRichard Tran Mills B->ops->mult = MatMult_SeqAIJMKL_SpMV2; 588df555b71SRichard Tran Mills /* B->ops->multtranspose = MatMultTranspose_SeqAIJMKL_SpMV2; */ 589df555b71SRichard Tran Mills B->ops->multadd = MatMultAdd_SeqAIJMKL_SpMV2; 590df555b71SRichard Tran Mills /* B->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJMKL_SpMV2; */ 591d995685eSRichard Tran Mills #endif 592c9d46305SRichard Tran Mills } else { 5934a2a386eSRichard Tran Mills B->ops->mult = MatMult_SeqAIJMKL; 594c9d46305SRichard Tran Mills /* B->ops->multtranspose = MatMultTranspose_SeqAIJMKL; */ 5954a2a386eSRichard Tran Mills B->ops->multadd = MatMultAdd_SeqAIJMKL; 596c9d46305SRichard Tran Mills /* B->ops->multtransposeadd = MatMultTransposeAdd_SeqAIJMKL; */ 597c9d46305SRichard Tran Mills } 5984a2a386eSRichard Tran Mills 5994a2a386eSRichard Tran Mills ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaijmkl_seqaij_C",MatConvert_SeqAIJMKL_SeqAIJ);CHKERRQ(ierr); 600e9c94282SRichard Tran Mills ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMult_seqdense_seqaijmkl_C",MatMatMult_SeqDense_SeqAIJ);CHKERRQ(ierr); 601e9c94282SRichard Tran Mills ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMultSymbolic_seqdense_seqaijmkl_C",MatMatMultSymbolic_SeqDense_SeqAIJ);CHKERRQ(ierr); 602e9c94282SRichard Tran Mills ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMultNumeric_seqdense_seqaijmkl_C",MatMatMultNumeric_SeqDense_SeqAIJ);CHKERRQ(ierr); 6034a2a386eSRichard Tran Mills 6044a2a386eSRichard Tran Mills ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJMKL);CHKERRQ(ierr); 6054a2a386eSRichard Tran Mills *newmat = B; 6064a2a386eSRichard Tran Mills PetscFunctionReturn(0); 6074a2a386eSRichard Tran Mills } 6084a2a386eSRichard Tran Mills 6094a2a386eSRichard Tran Mills #undef __FUNCT__ 6104a2a386eSRichard Tran Mills #define __FUNCT__ "MatCreateSeqAIJMKL" 6114a2a386eSRichard Tran Mills /*@C 6124a2a386eSRichard Tran Mills MatCreateSeqAIJMKL - Creates a sparse matrix of type SEQAIJMKL. 6134a2a386eSRichard Tran Mills This type inherits from AIJ and is largely identical, but uses sparse BLAS 6144a2a386eSRichard Tran Mills routines from Intel MKL whenever possible. 6154a2a386eSRichard Tran Mills Collective on MPI_Comm 6164a2a386eSRichard Tran Mills 6174a2a386eSRichard Tran Mills Input Parameters: 6184a2a386eSRichard Tran Mills + comm - MPI communicator, set to PETSC_COMM_SELF 6194a2a386eSRichard Tran Mills . m - number of rows 6204a2a386eSRichard Tran Mills . n - number of columns 6214a2a386eSRichard Tran Mills . nz - number of nonzeros per row (same for all rows) 6224a2a386eSRichard Tran Mills - nnz - array containing the number of nonzeros in the various rows 6234a2a386eSRichard Tran Mills (possibly different for each row) or NULL 6244a2a386eSRichard Tran Mills 6254a2a386eSRichard Tran Mills Output Parameter: 6264a2a386eSRichard Tran Mills . A - the matrix 6274a2a386eSRichard Tran Mills 6284a2a386eSRichard Tran Mills Notes: 6294a2a386eSRichard Tran Mills If nnz is given then nz is ignored 6304a2a386eSRichard Tran Mills 6314a2a386eSRichard Tran Mills Level: intermediate 6324a2a386eSRichard Tran Mills 6334a2a386eSRichard Tran Mills .keywords: matrix, cray, sparse, parallel 6344a2a386eSRichard Tran Mills 6354a2a386eSRichard Tran Mills .seealso: MatCreate(), MatCreateMPIAIJMKL(), MatSetValues() 6364a2a386eSRichard Tran Mills @*/ 6374a2a386eSRichard Tran Mills PetscErrorCode MatCreateSeqAIJMKL(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A) 6384a2a386eSRichard Tran Mills { 6394a2a386eSRichard Tran Mills PetscErrorCode ierr; 6404a2a386eSRichard Tran Mills 6414a2a386eSRichard Tran Mills PetscFunctionBegin; 6424a2a386eSRichard Tran Mills ierr = MatCreate(comm,A);CHKERRQ(ierr); 6434a2a386eSRichard Tran Mills ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr); 6444a2a386eSRichard Tran Mills ierr = MatSetType(*A,MATSEQAIJMKL);CHKERRQ(ierr); 6454a2a386eSRichard Tran Mills ierr = MatSeqAIJSetPreallocation_SeqAIJ(*A,nz,nnz);CHKERRQ(ierr); 6464a2a386eSRichard Tran Mills PetscFunctionReturn(0); 6474a2a386eSRichard Tran Mills } 6484a2a386eSRichard Tran Mills 6494a2a386eSRichard Tran Mills #undef __FUNCT__ 6504a2a386eSRichard Tran Mills #define __FUNCT__ "MatCreate_SeqAIJMKL" 6514a2a386eSRichard Tran Mills PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJMKL(Mat A) 6524a2a386eSRichard Tran Mills { 6534a2a386eSRichard Tran Mills PetscErrorCode ierr; 6544a2a386eSRichard Tran Mills 6554a2a386eSRichard Tran Mills PetscFunctionBegin; 6564a2a386eSRichard Tran Mills ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr); 6574a2a386eSRichard Tran Mills ierr = MatConvert_SeqAIJ_SeqAIJMKL(A,MATSEQAIJMKL,MAT_INPLACE_MATRIX,&A);CHKERRQ(ierr); 6584a2a386eSRichard Tran Mills PetscFunctionReturn(0); 6594a2a386eSRichard Tran Mills } 660