114b396bbSKris Buschelman /*$Id: superlu.c,v 1.10 2001/08/15 15:56:50 bsmith Exp $*/ 214b396bbSKris Buschelman 314b396bbSKris Buschelman /* 414b396bbSKris Buschelman Provides an interface to the SuperLU sparse solver 514b396bbSKris Buschelman Modified for SuperLU 2.0 by Matthew Knepley 614b396bbSKris Buschelman */ 714b396bbSKris Buschelman 814b396bbSKris Buschelman #include "src/mat/impls/aij/seq/aij.h" 914b396bbSKris Buschelman 1014b396bbSKris Buschelman EXTERN_C_BEGIN 1114b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 1214b396bbSKris Buschelman #include "zsp_defs.h" 1314b396bbSKris Buschelman #else 1414b396bbSKris Buschelman #include "dsp_defs.h" 1514b396bbSKris Buschelman #endif 1614b396bbSKris Buschelman #include "util.h" 1714b396bbSKris Buschelman EXTERN_C_END 1814b396bbSKris Buschelman 1914b396bbSKris Buschelman typedef struct { 2014b396bbSKris Buschelman SuperMatrix A,B,AC,L,U; 2114b396bbSKris Buschelman int *perm_r,*perm_c,ispec,relax,panel_size; 2214b396bbSKris Buschelman double pivot_threshold; 2314b396bbSKris Buschelman NCformat *store; 2414b396bbSKris Buschelman MatStructure flg; 2514b396bbSKris Buschelman PetscTruth SuperluMatOdering; 2614b396bbSKris Buschelman 2714b396bbSKris Buschelman /* A few function pointers for inheritance */ 2814b396bbSKris Buschelman int (*MatView)(Mat,PetscViewer); 2914b396bbSKris Buschelman int (*MatAssemblyEnd)(Mat,MatAssemblyType); 30b0592601SKris Buschelman int (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*); 3114b396bbSKris Buschelman int (*MatDestroy)(Mat); 3214b396bbSKris Buschelman 3314b396bbSKris Buschelman /* Flag to clean up (non-global) SuperLU objects during Destroy */ 3414b396bbSKris Buschelman PetscTruth CleanUpSuperLU; 3514b396bbSKris Buschelman } Mat_SeqAIJ_SuperLU; 3614b396bbSKris Buschelman 3714b396bbSKris Buschelman 380e3434eeSKris Buschelman EXTERN int MatSeqAIJFactorInfo_SuperLU(Mat,PetscViewer); 39b0592601SKris Buschelman EXTERN int MatLUFactorSymbolic_SeqAIJ_SuperLU(Mat,IS,IS,MatFactorInfo*,Mat*); 40b0592601SKris Buschelman 41b0592601SKris Buschelman EXTERN_C_BEGIN 42b0592601SKris Buschelman EXTERN int MatConvert_SuperLU_SeqAIJ(Mat,MatType,Mat*); 43b0592601SKris Buschelman EXTERN int MatConvert_SeqAIJ_SuperLU(Mat,MatType,Mat*); 44b0592601SKris Buschelman EXTERN_C_END 4514b396bbSKris Buschelman 4614b396bbSKris Buschelman #undef __FUNCT__ 4714b396bbSKris Buschelman #define __FUNCT__ "MatDestroy_SeqAIJ_SuperLU" 4814b396bbSKris Buschelman int MatDestroy_SeqAIJ_SuperLU(Mat A) 4914b396bbSKris Buschelman { 50460c4903SKris Buschelman int ierr; 5114b396bbSKris Buschelman Mat_SeqAIJ_SuperLU *lu = (Mat_SeqAIJ_SuperLU*)A->spptr; 5214b396bbSKris Buschelman 5314b396bbSKris Buschelman PetscFunctionBegin; 54fb3e25aaSKris Buschelman if (lu->CleanUpSuperLU) { 55fb3e25aaSKris Buschelman /* We have to free the global data or SuperLU crashes (sucky design)*/ 56fb3e25aaSKris Buschelman /* Since we don't know if more solves on other matrices may be done 57fb3e25aaSKris Buschelman we cannot free the yucky SuperLU global data 58fb3e25aaSKris Buschelman StatFree(); 59fb3e25aaSKris Buschelman */ 60fb3e25aaSKris Buschelman 61fb3e25aaSKris Buschelman /* Free the SuperLU datastructures */ 62fb3e25aaSKris Buschelman Destroy_CompCol_Permuted(&lu->AC); 63fb3e25aaSKris Buschelman Destroy_SuperNode_Matrix(&lu->L); 64fb3e25aaSKris Buschelman Destroy_CompCol_Matrix(&lu->U); 65fb3e25aaSKris Buschelman ierr = PetscFree(lu->B.Store);CHKERRQ(ierr); 66fb3e25aaSKris Buschelman ierr = PetscFree(lu->A.Store);CHKERRQ(ierr); 67fb3e25aaSKris Buschelman ierr = PetscFree(lu->perm_r);CHKERRQ(ierr); 68fb3e25aaSKris Buschelman ierr = PetscFree(lu->perm_c);CHKERRQ(ierr); 69fb3e25aaSKris Buschelman } 70b0592601SKris Buschelman ierr = MatConvert_SuperLU_SeqAIJ(A,MATSEQAIJ,&A);CHKERRQ(ierr); 71fb3e25aaSKris Buschelman ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 7214b396bbSKris Buschelman PetscFunctionReturn(0); 7314b396bbSKris Buschelman } 7414b396bbSKris Buschelman 7514b396bbSKris Buschelman #undef __FUNCT__ 7614b396bbSKris Buschelman #define __FUNCT__ "MatView_SeqAIJ_Spooles" 7714b396bbSKris Buschelman int MatView_SeqAIJ_SuperLU(Mat A,PetscViewer viewer) 7814b396bbSKris Buschelman { 7914b396bbSKris Buschelman int ierr; 8014b396bbSKris Buschelman PetscTruth isascii; 8114b396bbSKris Buschelman PetscViewerFormat format; 8214b396bbSKris Buschelman Mat_SeqAIJ_SuperLU *lu=(Mat_SeqAIJ_SuperLU*)(A->spptr); 8314b396bbSKris Buschelman 8414b396bbSKris Buschelman PetscFunctionBegin; 8514b396bbSKris Buschelman ierr = (*lu->MatView)(A,viewer);CHKERRQ(ierr); 8614b396bbSKris Buschelman 8714b396bbSKris Buschelman ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 8814b396bbSKris Buschelman if (isascii) { 8914b396bbSKris Buschelman ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 9014b396bbSKris Buschelman if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) { 9114b396bbSKris Buschelman ierr = MatSeqAIJFactorInfo_SuperLU(A,viewer);CHKERRQ(ierr); 9214b396bbSKris Buschelman } 9314b396bbSKris Buschelman } 9414b396bbSKris Buschelman PetscFunctionReturn(0); 9514b396bbSKris Buschelman } 9614b396bbSKris Buschelman 9714b396bbSKris Buschelman #undef __FUNCT__ 9814b396bbSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_SeqAIJ_SuperLU" 9914b396bbSKris Buschelman int MatAssemblyEnd_SeqAIJ_SuperLU(Mat A,MatAssemblyType mode) { 10014b396bbSKris Buschelman int ierr; 10114b396bbSKris Buschelman Mat_SeqAIJ_SuperLU *lu=(Mat_SeqAIJ_SuperLU*)(A->spptr); 10214b396bbSKris Buschelman 10314b396bbSKris Buschelman PetscFunctionBegin; 10414b396bbSKris Buschelman ierr = (*lu->MatAssemblyEnd)(A,mode);CHKERRQ(ierr); 105b0592601SKris Buschelman 106b0592601SKris Buschelman lu->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 107b0592601SKris Buschelman A->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ_SuperLU; 10814b396bbSKris Buschelman PetscFunctionReturn(0); 10914b396bbSKris Buschelman } 11014b396bbSKris Buschelman 11114b396bbSKris Buschelman #include "src/mat/impls/dense/seq/dense.h" 11214b396bbSKris Buschelman #undef __FUNCT__ 11314b396bbSKris Buschelman #define __FUNCT__ "MatCreateNull_SeqAIJ_SuperLU" 11414b396bbSKris Buschelman int MatCreateNull_SeqAIJ_SuperLU(Mat A,Mat *nullMat) 11514b396bbSKris Buschelman { 11614b396bbSKris Buschelman Mat_SeqAIJ_SuperLU *lu = (Mat_SeqAIJ_SuperLU*)A->spptr; 11714b396bbSKris Buschelman int numRows = A->m,numCols = A->n; 11814b396bbSKris Buschelman SCformat *Lstore; 11914b396bbSKris Buschelman int numNullCols,size; 12014b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 12114b396bbSKris Buschelman doublecomplex *nullVals,*workVals; 12214b396bbSKris Buschelman #else 12314b396bbSKris Buschelman PetscScalar *nullVals,*workVals; 12414b396bbSKris Buschelman #endif 12514b396bbSKris Buschelman int row,newRow,col,newCol,block,b,ierr; 12614b396bbSKris Buschelman 12714b396bbSKris Buschelman PetscFunctionBegin; 12814b396bbSKris Buschelman PetscValidHeaderSpecific(A,MAT_COOKIE); 12914b396bbSKris Buschelman PetscValidPointer(nullMat); 13014b396bbSKris Buschelman if (!A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix"); 13114b396bbSKris Buschelman numNullCols = numCols - numRows; 13214b396bbSKris Buschelman if (numNullCols < 0) SETERRQ(PETSC_ERR_ARG_WRONG,"Function only applies to underdetermined problems"); 13314b396bbSKris Buschelman /* Create the null matrix */ 13414b396bbSKris Buschelman ierr = MatCreateSeqDense(A->comm,numRows,numNullCols,PETSC_NULL,nullMat);CHKERRQ(ierr); 13514b396bbSKris Buschelman if (numNullCols == 0) { 13614b396bbSKris Buschelman ierr = MatAssemblyBegin(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13714b396bbSKris Buschelman ierr = MatAssemblyEnd(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13814b396bbSKris Buschelman PetscFunctionReturn(0); 13914b396bbSKris Buschelman } 14014b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 14114b396bbSKris Buschelman nullVals = (doublecomplex*)((Mat_SeqDense*)(*nullMat)->data)->v; 14214b396bbSKris Buschelman #else 14314b396bbSKris Buschelman nullVals = ((Mat_SeqDense*)(*nullMat)->data)->v; 14414b396bbSKris Buschelman #endif 14514b396bbSKris Buschelman /* Copy in the columns */ 14614b396bbSKris Buschelman Lstore = (SCformat*)lu->L.Store; 14714b396bbSKris Buschelman for(block = 0; block <= Lstore->nsuper; block++) { 14814b396bbSKris Buschelman newRow = Lstore->sup_to_col[block]; 14914b396bbSKris Buschelman size = Lstore->sup_to_col[block+1] - Lstore->sup_to_col[block]; 15014b396bbSKris Buschelman for(col = Lstore->rowind_colptr[newRow]; col < Lstore->rowind_colptr[newRow+1]; col++) { 15114b396bbSKris Buschelman newCol = Lstore->rowind[col]; 15214b396bbSKris Buschelman if (newCol >= numRows) { 15314b396bbSKris Buschelman for(b = 0; b < size; b++) 15414b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 15514b396bbSKris Buschelman nullVals[(newCol-numRows)*numRows+newRow+b] = ((doublecomplex*)Lstore->nzval)[Lstore->nzval_colptr[newRow+b]+col]; 15614b396bbSKris Buschelman #else 15714b396bbSKris Buschelman nullVals[(newCol-numRows)*numRows+newRow+b] = ((double*)Lstore->nzval)[Lstore->nzval_colptr[newRow+b]+col]; 15814b396bbSKris Buschelman #endif 15914b396bbSKris Buschelman } 16014b396bbSKris Buschelman } 16114b396bbSKris Buschelman } 16214b396bbSKris Buschelman /* Permute rhs to form P^T_c B */ 16314b396bbSKris Buschelman ierr = PetscMalloc(numRows*sizeof(double),&workVals);CHKERRQ(ierr); 16414b396bbSKris Buschelman for(b = 0; b < numNullCols; b++) { 16514b396bbSKris Buschelman for(row = 0; row < numRows; row++) workVals[lu->perm_c[row]] = nullVals[b*numRows+row]; 16614b396bbSKris Buschelman for(row = 0; row < numRows; row++) nullVals[b*numRows+row] = workVals[row]; 16714b396bbSKris Buschelman } 16814b396bbSKris Buschelman /* Backward solve the upper triangle A x = b */ 16914b396bbSKris Buschelman for(b = 0; b < numNullCols; b++) { 17014b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 17114b396bbSKris Buschelman sp_ztrsv("L","T","U",&lu->L,&lu->U,&nullVals[b*numRows],&ierr); 17214b396bbSKris Buschelman #else 17314b396bbSKris Buschelman sp_dtrsv("L","T","U",&lu->L,&lu->U,&nullVals[b*numRows],&ierr); 17414b396bbSKris Buschelman #endif 17514b396bbSKris Buschelman if (ierr < 0) 17614b396bbSKris Buschelman SETERRQ1(PETSC_ERR_ARG_WRONG,"The argument %d was invalid",-ierr); 17714b396bbSKris Buschelman } 17814b396bbSKris Buschelman ierr = PetscFree(workVals);CHKERRQ(ierr); 17914b396bbSKris Buschelman 18014b396bbSKris Buschelman ierr = MatAssemblyBegin(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 18114b396bbSKris Buschelman ierr = MatAssemblyEnd(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 18214b396bbSKris Buschelman PetscFunctionReturn(0); 18314b396bbSKris Buschelman } 18414b396bbSKris Buschelman 18514b396bbSKris Buschelman #undef __FUNCT__ 18614b396bbSKris Buschelman #define __FUNCT__ "MatSolve_SeqAIJ_SuperLU" 18714b396bbSKris Buschelman int MatSolve_SeqAIJ_SuperLU(Mat A,Vec b,Vec x) 18814b396bbSKris Buschelman { 18914b396bbSKris Buschelman Mat_SeqAIJ_SuperLU *lu = (Mat_SeqAIJ_SuperLU*)A->spptr; 19014b396bbSKris Buschelman PetscScalar *array; 19114b396bbSKris Buschelman int m,ierr; 19214b396bbSKris Buschelman 19314b396bbSKris Buschelman PetscFunctionBegin; 19414b396bbSKris Buschelman ierr = VecGetLocalSize(b,&m);CHKERRQ(ierr); 19514b396bbSKris Buschelman ierr = VecCopy(b,x);CHKERRQ(ierr); 19614b396bbSKris Buschelman ierr = VecGetArray(x,&array);CHKERRQ(ierr); 19714b396bbSKris Buschelman /* Create the Rhs */ 19814b396bbSKris Buschelman lu->B.Stype = SLU_DN; 19914b396bbSKris Buschelman lu->B.Mtype = SLU_GE; 20014b396bbSKris Buschelman lu->B.nrow = m; 20114b396bbSKris Buschelman lu->B.ncol = 1; 20214b396bbSKris Buschelman ((DNformat*)lu->B.Store)->lda = m; 20314b396bbSKris Buschelman ((DNformat*)lu->B.Store)->nzval = array; 20414b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 20514b396bbSKris Buschelman lu->B.Dtype = SLU_Z; 20614b396bbSKris Buschelman zgstrs("T",&lu->L,&lu->U,lu->perm_r,lu->perm_c,&lu->B,&ierr); 20714b396bbSKris Buschelman #else 20814b396bbSKris Buschelman lu->B.Dtype = SLU_D; 20914b396bbSKris Buschelman dgstrs("T",&lu->L,&lu->U,lu->perm_r,lu->perm_c,&lu->B,&ierr); 21014b396bbSKris Buschelman #endif 21114b396bbSKris Buschelman if (ierr < 0) SETERRQ1(PETSC_ERR_ARG_WRONG,"The diagonal element of row %d was invalid",-ierr); 21214b396bbSKris Buschelman ierr = VecRestoreArray(x,&array);CHKERRQ(ierr); 21314b396bbSKris Buschelman PetscFunctionReturn(0); 21414b396bbSKris Buschelman } 21514b396bbSKris Buschelman 21614b396bbSKris Buschelman static int StatInitCalled = 0; 21714b396bbSKris Buschelman 21814b396bbSKris Buschelman #undef __FUNCT__ 21914b396bbSKris Buschelman #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ_SuperLU" 22014b396bbSKris Buschelman int MatLUFactorNumeric_SeqAIJ_SuperLU(Mat A,Mat *F) 22114b396bbSKris Buschelman { 22214b396bbSKris Buschelman Mat_SeqAIJ *aa = (Mat_SeqAIJ*)(A)->data; 22314b396bbSKris Buschelman Mat_SeqAIJ_SuperLU *lu = (Mat_SeqAIJ_SuperLU*)(*F)->spptr; 224d54de34fSKris Buschelman int *etree,ierr; 22514b396bbSKris Buschelman PetscTruth flag; 22614b396bbSKris Buschelman 22714b396bbSKris Buschelman PetscFunctionBegin; 22814b396bbSKris Buschelman /* Create the SuperMatrix for A^T: 22914b396bbSKris Buschelman Since SuperLU only likes column-oriented matrices,we pass it the transpose, 23014b396bbSKris Buschelman and then solve A^T X = B in MatSolve(). 23114b396bbSKris Buschelman */ 23214b396bbSKris Buschelman 23314b396bbSKris Buschelman if ( lu->flg == DIFFERENT_NONZERO_PATTERN){ /* first numerical factorization */ 23414b396bbSKris Buschelman lu->A.Stype = SLU_NC; 23514b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 23614b396bbSKris Buschelman lu->A.Dtype = SLU_Z; 23714b396bbSKris Buschelman #else 23814b396bbSKris Buschelman lu->A.Dtype = SLU_D; 23914b396bbSKris Buschelman #endif 24014b396bbSKris Buschelman lu->A.Mtype = SLU_GE; 24114b396bbSKris Buschelman lu->A.nrow = A->n; 24214b396bbSKris Buschelman lu->A.ncol = A->m; 24314b396bbSKris Buschelman 24414b396bbSKris Buschelman ierr = PetscMalloc(sizeof(NCformat),&lu->store);CHKERRQ(ierr); 24514b396bbSKris Buschelman ierr = PetscMalloc(sizeof(DNformat),&lu->B.Store);CHKERRQ(ierr); 24614b396bbSKris Buschelman } 24714b396bbSKris Buschelman lu->store->nnz = aa->nz; 24814b396bbSKris Buschelman lu->store->colptr = aa->i; 24914b396bbSKris Buschelman lu->store->rowind = aa->j; 25014b396bbSKris Buschelman lu->store->nzval = aa->a; 25114b396bbSKris Buschelman lu->A.Store = lu->store; 25214b396bbSKris Buschelman 25314b396bbSKris Buschelman /* Set SuperLU options */ 25414b396bbSKris Buschelman lu->relax = sp_ienv(2); 25514b396bbSKris Buschelman lu->panel_size = sp_ienv(1); 25614b396bbSKris Buschelman /* We have to initialize global data or SuperLU crashes (sucky design) */ 25714b396bbSKris Buschelman if (!StatInitCalled) { 25814b396bbSKris Buschelman StatInit(lu->panel_size,lu->relax); 25914b396bbSKris Buschelman } 26014b396bbSKris Buschelman StatInitCalled++; 26114b396bbSKris Buschelman 26214b396bbSKris Buschelman ierr = PetscOptionsBegin(A->comm,A->prefix,"SuperLU Options","Mat");CHKERRQ(ierr); 26314b396bbSKris Buschelman /* use SuperLU mat ordeing */ 26414b396bbSKris Buschelman ierr = PetscOptionsInt("-mat_superlu_ordering","SuperLU ordering type (one of 0, 1, 2, 3)\n 0: natural ordering;\n 1: MMD applied to A'*A;\n 2: MMD applied to A'+A;\n 3: COLAMD, approximate minimum degree column ordering","None",lu->ispec,&lu->ispec,&flag);CHKERRQ(ierr); 26514b396bbSKris Buschelman if (flag) { 26614b396bbSKris Buschelman get_perm_c(lu->ispec, &lu->A, lu->perm_c); 26714b396bbSKris Buschelman lu->SuperluMatOdering = PETSC_TRUE; 26814b396bbSKris Buschelman } 26914b396bbSKris Buschelman PetscOptionsEnd(); 27014b396bbSKris Buschelman 27114b396bbSKris Buschelman /* Create the elimination tree */ 27214b396bbSKris Buschelman ierr = PetscMalloc(A->n*sizeof(int),&etree);CHKERRQ(ierr); 27314b396bbSKris Buschelman sp_preorder("N",&lu->A,lu->perm_c,etree,&lu->AC); 27414b396bbSKris Buschelman /* Factor the matrix */ 27514b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 27614b396bbSKris Buschelman zgstrf("N",&lu->AC,lu->pivot_threshold,0.0,lu->relax,lu->panel_size,etree,PETSC_NULL,0,lu->perm_r,lu->perm_c,&lu->L,&lu->U,&ierr); 27714b396bbSKris Buschelman #else 27814b396bbSKris Buschelman dgstrf("N",&lu->AC,lu->pivot_threshold,0.0,lu->relax,lu->panel_size,etree,PETSC_NULL,0,lu->perm_r,lu->perm_c,&lu->L,&lu->U,&ierr); 27914b396bbSKris Buschelman #endif 28014b396bbSKris Buschelman if (ierr < 0) { 28114b396bbSKris Buschelman SETERRQ1(PETSC_ERR_ARG_WRONG,"The diagonal element of row %d was invalid",-ierr); 28214b396bbSKris Buschelman } else if (ierr > 0) { 28314b396bbSKris Buschelman if (ierr <= A->m) { 28414b396bbSKris Buschelman SETERRQ1(PETSC_ERR_ARG_WRONG,"The diagonal element %d of U is exactly zero",ierr); 28514b396bbSKris Buschelman } else { 28614b396bbSKris Buschelman SETERRQ1(PETSC_ERR_ARG_WRONG,"Memory allocation failure after %d bytes were allocated",ierr-A->m); 28714b396bbSKris Buschelman } 28814b396bbSKris Buschelman } 28914b396bbSKris Buschelman 29014b396bbSKris Buschelman /* Cleanup */ 29114b396bbSKris Buschelman ierr = PetscFree(etree);CHKERRQ(ierr); 29214b396bbSKris Buschelman 29314b396bbSKris Buschelman lu->flg = SAME_NONZERO_PATTERN; 29414b396bbSKris Buschelman PetscFunctionReturn(0); 29514b396bbSKris Buschelman } 29614b396bbSKris Buschelman 29714b396bbSKris Buschelman /* 29814b396bbSKris Buschelman Note the r permutation is ignored 29914b396bbSKris Buschelman */ 30014b396bbSKris Buschelman #undef __FUNCT__ 30114b396bbSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJ_SuperLU" 30214b396bbSKris Buschelman int MatLUFactorSymbolic_SeqAIJ_SuperLU(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) 30314b396bbSKris Buschelman { 30414b396bbSKris Buschelman Mat B; 30514b396bbSKris Buschelman Mat_SeqAIJ_SuperLU *lu; 30614b396bbSKris Buschelman int ierr,*ca; 30714b396bbSKris Buschelman 30814b396bbSKris Buschelman PetscFunctionBegin; 30914b396bbSKris Buschelman 310899d7b4fSKris Buschelman ierr = MatCreate(A->comm,A->m,A->n,PETSC_DETERMINE,PETSC_DETERMINE,&B);CHKERRQ(ierr); 311899d7b4fSKris Buschelman ierr = MatSetType(B,MATSUPERLU);CHKERRQ(ierr); 312899d7b4fSKris Buschelman ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 31314b396bbSKris Buschelman B->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJ_SuperLU; 31414b396bbSKris Buschelman B->ops->solve = MatSolve_SeqAIJ_SuperLU; 31514b396bbSKris Buschelman B->factor = FACTOR_LU; 316899d7b4fSKris Buschelman B->assembled = PETSC_TRUE; /* required by -sles_view */ 31714b396bbSKris Buschelman 318899d7b4fSKris Buschelman lu = (Mat_SeqAIJ_SuperLU*)(B->spptr); 31914b396bbSKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)B,"MatCreateNull","MatCreateNull_SeqAIJ_SuperLU", 32014b396bbSKris Buschelman (void(*)(void))MatCreateNull_SeqAIJ_SuperLU);CHKERRQ(ierr); 32114b396bbSKris Buschelman 32214b396bbSKris Buschelman /* Allocate the work arrays required by SuperLU (notice sizes are for the transpose) */ 32314b396bbSKris Buschelman ierr = PetscMalloc(A->n*sizeof(int),&lu->perm_r);CHKERRQ(ierr); 32414b396bbSKris Buschelman ierr = PetscMalloc(A->m*sizeof(int),&lu->perm_c);CHKERRQ(ierr); 32514b396bbSKris Buschelman 32614b396bbSKris Buschelman /* use PETSc mat ordering */ 32714b396bbSKris Buschelman ierr = ISGetIndices(c,&ca);CHKERRQ(ierr); 32814b396bbSKris Buschelman ierr = PetscMemcpy(lu->perm_c,ca,A->m*sizeof(int));CHKERRQ(ierr); 32914b396bbSKris Buschelman ierr = ISRestoreIndices(c,&ca);CHKERRQ(ierr); 33014b396bbSKris Buschelman lu->SuperluMatOdering = PETSC_FALSE; 33114b396bbSKris Buschelman 33214b396bbSKris Buschelman lu->pivot_threshold = info->dtcol; 33314b396bbSKris Buschelman PetscLogObjectMemory(B,(A->m+A->n)*sizeof(int)+sizeof(Mat_SeqAIJ_SuperLU)); 33414b396bbSKris Buschelman 33514b396bbSKris Buschelman lu->flg = DIFFERENT_NONZERO_PATTERN; 336e740cb95SKris Buschelman lu->CleanUpSuperLU = PETSC_TRUE; 337899d7b4fSKris Buschelman 338899d7b4fSKris Buschelman *F = B; 33914b396bbSKris Buschelman PetscFunctionReturn(0); 34014b396bbSKris Buschelman } 34114b396bbSKris Buschelman 34214b396bbSKris Buschelman /* used by -sles_view */ 34314b396bbSKris Buschelman #undef __FUNCT__ 34414b396bbSKris Buschelman #define __FUNCT__ "MatSeqAIJFactorInfo_SuperLU" 34514b396bbSKris Buschelman int MatSeqAIJFactorInfo_SuperLU(Mat A,PetscViewer viewer) 34614b396bbSKris Buschelman { 34714b396bbSKris Buschelman Mat_SeqAIJ_SuperLU *lu= (Mat_SeqAIJ_SuperLU*)A->spptr; 34814b396bbSKris Buschelman int ierr; 34914b396bbSKris Buschelman PetscFunctionBegin; 35014b396bbSKris Buschelman 35114b396bbSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");CHKERRQ(ierr); 35214b396bbSKris Buschelman if(lu->SuperluMatOdering) ierr = PetscViewerASCIIPrintf(viewer," SuperLU mat ordering: %d\n",lu->ispec);CHKERRQ(ierr); 35314b396bbSKris Buschelman 35414b396bbSKris Buschelman PetscFunctionReturn(0); 35514b396bbSKris Buschelman } 35614b396bbSKris Buschelman 35714b396bbSKris Buschelman EXTERN_C_BEGIN 35814b396bbSKris Buschelman #undef __FUNCT__ 359b0592601SKris Buschelman #define __FUNCT__ "MatConvert_SuperLU_SeqAIJ" 360b0592601SKris Buschelman int MatConvert_SuperLU_SeqAIJ(Mat A,MatType type,Mat *newmat) { 361fb3e25aaSKris Buschelman /* This routine is only called to convert an unfactored PETSc-SuperLU matrix */ 362fb3e25aaSKris Buschelman /* to its base PETSc type, so we will ignore 'MatType type'. */ 36314b396bbSKris Buschelman int ierr; 364b0592601SKris Buschelman Mat B=*newmat; 365b0592601SKris Buschelman Mat_SeqAIJ_SuperLU *lu=(Mat_SeqAIJ_SuperLU *)A->spptr; 366b0592601SKris Buschelman 367b0592601SKris Buschelman PetscFunctionBegin; 368b0592601SKris Buschelman if (B != A) { 369b0592601SKris Buschelman /* This routine was inherited from SeqAIJ. */ 370b0592601SKris Buschelman ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 371b0592601SKris Buschelman } else { 372b0592601SKris Buschelman /* Reset the original function pointers */ 373b0592601SKris Buschelman B->ops->view = lu->MatView; 374b0592601SKris Buschelman B->ops->assemblyend = lu->MatAssemblyEnd; 375b0592601SKris Buschelman B->ops->lufactorsymbolic = lu->MatLUFactorSymbolic; 376b0592601SKris Buschelman B->ops->destroy = lu->MatDestroy; 377b0592601SKris Buschelman /* lu is only a function pointer stash unless we've factored the matrix, which we haven't! */ 378b0592601SKris Buschelman ierr = PetscFree(lu);CHKERRQ(ierr); 379b0592601SKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);CHKERRQ(ierr); 380b0592601SKris Buschelman } 381b0592601SKris Buschelman *newmat = B; 382b0592601SKris Buschelman PetscFunctionReturn(0); 383b0592601SKris Buschelman } 384b0592601SKris Buschelman EXTERN_C_END 385b0592601SKris Buschelman 386b0592601SKris Buschelman EXTERN_C_BEGIN 387b0592601SKris Buschelman #undef __FUNCT__ 388b0592601SKris Buschelman #define __FUNCT__ "MatConvert_SeqAIJ_SuperLU" 389b0592601SKris Buschelman int MatConvert_SeqAIJ_SuperLU(Mat A,MatType type,Mat *newmat) { 390fb3e25aaSKris Buschelman /* This routine is only called to convert to MATSUPERLU */ 391fb3e25aaSKris Buschelman /* from MATSEQAIJ, so we will ignore 'MatType type'. */ 392b0592601SKris Buschelman int ierr; 393b0592601SKris Buschelman Mat B=*newmat; 39414b396bbSKris Buschelman Mat_SeqAIJ_SuperLU *lu; 39514b396bbSKris Buschelman 39614b396bbSKris Buschelman PetscFunctionBegin; 397b0592601SKris Buschelman if (B != A) { 398b0592601SKris Buschelman ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 399b0592601SKris Buschelman } 40014b396bbSKris Buschelman 40114b396bbSKris Buschelman ierr = PetscNew(Mat_SeqAIJ_SuperLU,&lu);CHKERRQ(ierr); 40214b396bbSKris Buschelman lu->MatView = A->ops->view; 40314b396bbSKris Buschelman lu->MatAssemblyEnd = A->ops->assemblyend; 404b0592601SKris Buschelman lu->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 40514b396bbSKris Buschelman lu->MatDestroy = A->ops->destroy; 40614b396bbSKris Buschelman lu->CleanUpSuperLU = PETSC_FALSE; 407b0592601SKris Buschelman 408b0592601SKris Buschelman B->spptr = (void*)lu; 409b0592601SKris Buschelman B->ops->view = MatView_SeqAIJ_SuperLU; 410b0592601SKris Buschelman B->ops->assemblyend = MatAssemblyEnd_SeqAIJ_SuperLU; 411b0592601SKris Buschelman B->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ_SuperLU; 412b0592601SKris Buschelman B->ops->destroy = MatDestroy_SeqAIJ_SuperLU; 413b0592601SKris Buschelman 414b0592601SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_superlu_C", 415b0592601SKris Buschelman "MatConvert_SeqAIJ_SuperLU",MatConvert_SeqAIJ_SuperLU);CHKERRQ(ierr); 416b0592601SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_superlu_seqaij_C", 417b0592601SKris Buschelman "MatConvert_SuperLU_SeqAIJ",MatConvert_SuperLU_SeqAIJ);CHKERRQ(ierr); 418fb3e25aaSKris Buschelman PetscLogInfo(0,"Using SuperLU for SeqAIJ LU factorization and solves."); 419fb3e25aaSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)B,MATSUPERLU);CHKERRQ(ierr); 420b0592601SKris Buschelman *newmat = B; 421b0592601SKris Buschelman PetscFunctionReturn(0); 422b0592601SKris Buschelman } 423b0592601SKris Buschelman EXTERN_C_END 424b0592601SKris Buschelman 42524b6179bSKris Buschelman /*MC 42624b6179bSKris Buschelman MATSUPERLU - a matrix type providing direct solvers (LU) for sequential matrices 42724b6179bSKris Buschelman via the external package SuperLU. 42824b6179bSKris Buschelman 42924b6179bSKris Buschelman If SuperLU is installed (see the manual for 43024b6179bSKris Buschelman instructions on how to declare the existence of external packages), 43124b6179bSKris Buschelman a matrix type can be constructed which invokes SuperLU solvers. 43224b6179bSKris Buschelman After calling MatCreate(...,A), simply call MatSetType(A,MATSUPERLU). 43324b6179bSKris Buschelman This matrix type is only supported for double precision real. 43424b6179bSKris Buschelman 43524b6179bSKris Buschelman This matrix inherits from MATSEQAIJ. As a result, MatSeqAIJSetPreallocation is 43624b6179bSKris Buschelman supported for this matrix type. 43724b6179bSKris Buschelman 43824b6179bSKris Buschelman Options Database Keys: 439*2f71e704SKris Buschelman + -mat_type superlu - sets the matrix type to superlu during a call to MatSetFromOptions() 44024b6179bSKris Buschelman - -mat_superlu_ordering <0,1,2,3> - 0: natural ordering, 44124b6179bSKris Buschelman 1: MMD applied to A'*A, 44224b6179bSKris Buschelman 2: MMD applied to A'+A, 44324b6179bSKris Buschelman 3: COLAMD, approximate minimum degree column ordering 44424b6179bSKris Buschelman 44524b6179bSKris Buschelman Level: beginner 44624b6179bSKris Buschelman 44724b6179bSKris Buschelman .seealso: PCLU 44824b6179bSKris Buschelman M*/ 44924b6179bSKris Buschelman 450b0592601SKris Buschelman EXTERN_C_BEGIN 451b0592601SKris Buschelman #undef __FUNCT__ 452b0592601SKris Buschelman #define __FUNCT__ "MatCreate_SeqAIJ_SuperLU" 453b0592601SKris Buschelman int MatCreate_SeqAIJ_SuperLU(Mat A) { 454b0592601SKris Buschelman int ierr; 455b0592601SKris Buschelman 456b0592601SKris Buschelman PetscFunctionBegin; 457b0592601SKris Buschelman ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr); 458b0592601SKris Buschelman ierr = MatConvert_SeqAIJ_SuperLU(A,MATSUPERLU,&A);CHKERRQ(ierr); 45914b396bbSKris Buschelman PetscFunctionReturn(0); 46014b396bbSKris Buschelman } 46114b396bbSKris Buschelman EXTERN_C_END 462