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 */ 28f0c56d0fSKris Buschelman int (*MatDuplicate)(Mat,MatDuplicateOption,Mat*); 2914b396bbSKris Buschelman int (*MatView)(Mat,PetscViewer); 3014b396bbSKris Buschelman int (*MatAssemblyEnd)(Mat,MatAssemblyType); 31b0592601SKris Buschelman int (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*); 3214b396bbSKris Buschelman int (*MatDestroy)(Mat); 3314b396bbSKris Buschelman 3414b396bbSKris Buschelman /* Flag to clean up (non-global) SuperLU objects during Destroy */ 3514b396bbSKris Buschelman PetscTruth CleanUpSuperLU; 36f0c56d0fSKris Buschelman } Mat_SuperLU; 3714b396bbSKris Buschelman 3814b396bbSKris Buschelman 39f0c56d0fSKris Buschelman EXTERN int MatFactorInfo_SuperLU(Mat,PetscViewer); 40f0c56d0fSKris Buschelman EXTERN int MatLUFactorSymbolic_SuperLU(Mat,IS,IS,MatFactorInfo*,Mat*); 41b0592601SKris Buschelman 42b0592601SKris Buschelman EXTERN_C_BEGIN 43b0592601SKris Buschelman EXTERN int MatConvert_SuperLU_SeqAIJ(Mat,MatType,Mat*); 44b0592601SKris Buschelman EXTERN int MatConvert_SeqAIJ_SuperLU(Mat,MatType,Mat*); 45b0592601SKris Buschelman EXTERN_C_END 4614b396bbSKris Buschelman 4714b396bbSKris Buschelman #undef __FUNCT__ 48f0c56d0fSKris Buschelman #define __FUNCT__ "MatDestroy_SuperLU" 49f0c56d0fSKris Buschelman int MatDestroy_SuperLU(Mat A) 5014b396bbSKris Buschelman { 51460c4903SKris Buschelman int ierr; 52f0c56d0fSKris Buschelman Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 5314b396bbSKris Buschelman 5414b396bbSKris Buschelman PetscFunctionBegin; 55fb3e25aaSKris Buschelman if (lu->CleanUpSuperLU) { 56fb3e25aaSKris Buschelman /* We have to free the global data or SuperLU crashes (sucky design)*/ 57fb3e25aaSKris Buschelman /* Since we don't know if more solves on other matrices may be done 58fb3e25aaSKris Buschelman we cannot free the yucky SuperLU global data 59fb3e25aaSKris Buschelman StatFree(); 60fb3e25aaSKris Buschelman */ 61fb3e25aaSKris Buschelman 62fb3e25aaSKris Buschelman /* Free the SuperLU datastructures */ 63fb3e25aaSKris Buschelman Destroy_CompCol_Permuted(&lu->AC); 64fb3e25aaSKris Buschelman Destroy_SuperNode_Matrix(&lu->L); 65fb3e25aaSKris Buschelman Destroy_CompCol_Matrix(&lu->U); 66fb3e25aaSKris Buschelman ierr = PetscFree(lu->B.Store);CHKERRQ(ierr); 67fb3e25aaSKris Buschelman ierr = PetscFree(lu->A.Store);CHKERRQ(ierr); 68fb3e25aaSKris Buschelman ierr = PetscFree(lu->perm_r);CHKERRQ(ierr); 69fb3e25aaSKris Buschelman ierr = PetscFree(lu->perm_c);CHKERRQ(ierr); 70fb3e25aaSKris Buschelman } 71b0592601SKris Buschelman ierr = MatConvert_SuperLU_SeqAIJ(A,MATSEQAIJ,&A);CHKERRQ(ierr); 72fb3e25aaSKris Buschelman ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 7314b396bbSKris Buschelman PetscFunctionReturn(0); 7414b396bbSKris Buschelman } 7514b396bbSKris Buschelman 7614b396bbSKris Buschelman #undef __FUNCT__ 77f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_SuperLU" 78f0c56d0fSKris Buschelman int MatView_SuperLU(Mat A,PetscViewer viewer) 7914b396bbSKris Buschelman { 8014b396bbSKris Buschelman int ierr; 8114b396bbSKris Buschelman PetscTruth isascii; 8214b396bbSKris Buschelman PetscViewerFormat format; 83f0c56d0fSKris Buschelman Mat_SuperLU *lu=(Mat_SuperLU*)(A->spptr); 8414b396bbSKris Buschelman 8514b396bbSKris Buschelman PetscFunctionBegin; 8614b396bbSKris Buschelman ierr = (*lu->MatView)(A,viewer);CHKERRQ(ierr); 8714b396bbSKris Buschelman 8814b396bbSKris Buschelman ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 8914b396bbSKris Buschelman if (isascii) { 9014b396bbSKris Buschelman ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 9114b396bbSKris Buschelman if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) { 92f0c56d0fSKris Buschelman ierr = MatFactorInfo_SuperLU(A,viewer);CHKERRQ(ierr); 9314b396bbSKris Buschelman } 9414b396bbSKris Buschelman } 9514b396bbSKris Buschelman PetscFunctionReturn(0); 9614b396bbSKris Buschelman } 9714b396bbSKris Buschelman 9814b396bbSKris Buschelman #undef __FUNCT__ 99f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_SuperLU" 100f0c56d0fSKris Buschelman int MatAssemblyEnd_SuperLU(Mat A,MatAssemblyType mode) { 10114b396bbSKris Buschelman int ierr; 102f0c56d0fSKris Buschelman Mat_SuperLU *lu=(Mat_SuperLU*)(A->spptr); 10314b396bbSKris Buschelman 10414b396bbSKris Buschelman PetscFunctionBegin; 10514b396bbSKris Buschelman ierr = (*lu->MatAssemblyEnd)(A,mode);CHKERRQ(ierr); 106b0592601SKris Buschelman 107b0592601SKris Buschelman lu->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 108f0c56d0fSKris Buschelman A->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU; 10914b396bbSKris Buschelman PetscFunctionReturn(0); 11014b396bbSKris Buschelman } 11114b396bbSKris Buschelman 11214b396bbSKris Buschelman #include "src/mat/impls/dense/seq/dense.h" 11314b396bbSKris Buschelman #undef __FUNCT__ 114f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreateNull_SuperLU" 115f0c56d0fSKris Buschelman int MatCreateNull_SuperLU(Mat A,Mat *nullMat) 11614b396bbSKris Buschelman { 117f0c56d0fSKris Buschelman Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 11814b396bbSKris Buschelman int numRows = A->m,numCols = A->n; 11914b396bbSKris Buschelman SCformat *Lstore; 12014b396bbSKris Buschelman int numNullCols,size; 12114b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 12214b396bbSKris Buschelman doublecomplex *nullVals,*workVals; 12314b396bbSKris Buschelman #else 12414b396bbSKris Buschelman PetscScalar *nullVals,*workVals; 12514b396bbSKris Buschelman #endif 12614b396bbSKris Buschelman int row,newRow,col,newCol,block,b,ierr; 12714b396bbSKris Buschelman 12814b396bbSKris Buschelman PetscFunctionBegin; 12914b396bbSKris Buschelman PetscValidHeaderSpecific(A,MAT_COOKIE); 13014b396bbSKris Buschelman PetscValidPointer(nullMat); 13114b396bbSKris Buschelman if (!A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix"); 13214b396bbSKris Buschelman numNullCols = numCols - numRows; 13314b396bbSKris Buschelman if (numNullCols < 0) SETERRQ(PETSC_ERR_ARG_WRONG,"Function only applies to underdetermined problems"); 13414b396bbSKris Buschelman /* Create the null matrix */ 13514b396bbSKris Buschelman ierr = MatCreateSeqDense(A->comm,numRows,numNullCols,PETSC_NULL,nullMat);CHKERRQ(ierr); 13614b396bbSKris Buschelman if (numNullCols == 0) { 13714b396bbSKris Buschelman ierr = MatAssemblyBegin(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13814b396bbSKris Buschelman ierr = MatAssemblyEnd(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13914b396bbSKris Buschelman PetscFunctionReturn(0); 14014b396bbSKris Buschelman } 14114b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 14214b396bbSKris Buschelman nullVals = (doublecomplex*)((Mat_SeqDense*)(*nullMat)->data)->v; 14314b396bbSKris Buschelman #else 14414b396bbSKris Buschelman nullVals = ((Mat_SeqDense*)(*nullMat)->data)->v; 14514b396bbSKris Buschelman #endif 14614b396bbSKris Buschelman /* Copy in the columns */ 14714b396bbSKris Buschelman Lstore = (SCformat*)lu->L.Store; 14814b396bbSKris Buschelman for(block = 0; block <= Lstore->nsuper; block++) { 14914b396bbSKris Buschelman newRow = Lstore->sup_to_col[block]; 15014b396bbSKris Buschelman size = Lstore->sup_to_col[block+1] - Lstore->sup_to_col[block]; 15114b396bbSKris Buschelman for(col = Lstore->rowind_colptr[newRow]; col < Lstore->rowind_colptr[newRow+1]; col++) { 15214b396bbSKris Buschelman newCol = Lstore->rowind[col]; 15314b396bbSKris Buschelman if (newCol >= numRows) { 15414b396bbSKris Buschelman for(b = 0; b < size; b++) 15514b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 15614b396bbSKris Buschelman nullVals[(newCol-numRows)*numRows+newRow+b] = ((doublecomplex*)Lstore->nzval)[Lstore->nzval_colptr[newRow+b]+col]; 15714b396bbSKris Buschelman #else 15814b396bbSKris Buschelman nullVals[(newCol-numRows)*numRows+newRow+b] = ((double*)Lstore->nzval)[Lstore->nzval_colptr[newRow+b]+col]; 15914b396bbSKris Buschelman #endif 16014b396bbSKris Buschelman } 16114b396bbSKris Buschelman } 16214b396bbSKris Buschelman } 16314b396bbSKris Buschelman /* Permute rhs to form P^T_c B */ 16414b396bbSKris Buschelman ierr = PetscMalloc(numRows*sizeof(double),&workVals);CHKERRQ(ierr); 16514b396bbSKris Buschelman for(b = 0; b < numNullCols; b++) { 16614b396bbSKris Buschelman for(row = 0; row < numRows; row++) workVals[lu->perm_c[row]] = nullVals[b*numRows+row]; 16714b396bbSKris Buschelman for(row = 0; row < numRows; row++) nullVals[b*numRows+row] = workVals[row]; 16814b396bbSKris Buschelman } 16914b396bbSKris Buschelman /* Backward solve the upper triangle A x = b */ 17014b396bbSKris Buschelman for(b = 0; b < numNullCols; b++) { 17114b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 17214b396bbSKris Buschelman sp_ztrsv("L","T","U",&lu->L,&lu->U,&nullVals[b*numRows],&ierr); 17314b396bbSKris Buschelman #else 17414b396bbSKris Buschelman sp_dtrsv("L","T","U",&lu->L,&lu->U,&nullVals[b*numRows],&ierr); 17514b396bbSKris Buschelman #endif 17614b396bbSKris Buschelman if (ierr < 0) 17714b396bbSKris Buschelman SETERRQ1(PETSC_ERR_ARG_WRONG,"The argument %d was invalid",-ierr); 17814b396bbSKris Buschelman } 17914b396bbSKris Buschelman ierr = PetscFree(workVals);CHKERRQ(ierr); 18014b396bbSKris Buschelman 18114b396bbSKris Buschelman ierr = MatAssemblyBegin(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 18214b396bbSKris Buschelman ierr = MatAssemblyEnd(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 18314b396bbSKris Buschelman PetscFunctionReturn(0); 18414b396bbSKris Buschelman } 18514b396bbSKris Buschelman 18614b396bbSKris Buschelman #undef __FUNCT__ 187f0c56d0fSKris Buschelman #define __FUNCT__ "MatSolve_SuperLU" 188f0c56d0fSKris Buschelman int MatSolve_SuperLU(Mat A,Vec b,Vec x) 18914b396bbSKris Buschelman { 190f0c56d0fSKris Buschelman Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 19114b396bbSKris Buschelman PetscScalar *array; 19214b396bbSKris Buschelman int m,ierr; 19314b396bbSKris Buschelman 19414b396bbSKris Buschelman PetscFunctionBegin; 19514b396bbSKris Buschelman ierr = VecGetLocalSize(b,&m);CHKERRQ(ierr); 19614b396bbSKris Buschelman ierr = VecCopy(b,x);CHKERRQ(ierr); 19714b396bbSKris Buschelman ierr = VecGetArray(x,&array);CHKERRQ(ierr); 19814b396bbSKris Buschelman /* Create the Rhs */ 19914b396bbSKris Buschelman lu->B.Stype = SLU_DN; 20014b396bbSKris Buschelman lu->B.Mtype = SLU_GE; 20114b396bbSKris Buschelman lu->B.nrow = m; 20214b396bbSKris Buschelman lu->B.ncol = 1; 20314b396bbSKris Buschelman ((DNformat*)lu->B.Store)->lda = m; 20414b396bbSKris Buschelman ((DNformat*)lu->B.Store)->nzval = array; 20514b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 20614b396bbSKris Buschelman lu->B.Dtype = SLU_Z; 20714b396bbSKris Buschelman zgstrs("T",&lu->L,&lu->U,lu->perm_r,lu->perm_c,&lu->B,&ierr); 20814b396bbSKris Buschelman #else 20914b396bbSKris Buschelman lu->B.Dtype = SLU_D; 21014b396bbSKris Buschelman dgstrs("T",&lu->L,&lu->U,lu->perm_r,lu->perm_c,&lu->B,&ierr); 21114b396bbSKris Buschelman #endif 21214b396bbSKris Buschelman if (ierr < 0) SETERRQ1(PETSC_ERR_ARG_WRONG,"The diagonal element of row %d was invalid",-ierr); 21314b396bbSKris Buschelman ierr = VecRestoreArray(x,&array);CHKERRQ(ierr); 21414b396bbSKris Buschelman PetscFunctionReturn(0); 21514b396bbSKris Buschelman } 21614b396bbSKris Buschelman 21714b396bbSKris Buschelman static int StatInitCalled = 0; 21814b396bbSKris Buschelman 21914b396bbSKris Buschelman #undef __FUNCT__ 220f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorNumeric_SuperLU" 221f0c56d0fSKris Buschelman int MatLUFactorNumeric_SuperLU(Mat A,Mat *F) 22214b396bbSKris Buschelman { 22314b396bbSKris Buschelman Mat_SeqAIJ *aa = (Mat_SeqAIJ*)(A)->data; 224f0c56d0fSKris Buschelman Mat_SuperLU *lu = (Mat_SuperLU*)(*F)->spptr; 225d54de34fSKris Buschelman int *etree,ierr; 22614b396bbSKris Buschelman PetscTruth flag; 22714b396bbSKris Buschelman 22814b396bbSKris Buschelman PetscFunctionBegin; 22914b396bbSKris Buschelman /* Create the SuperMatrix for A^T: 23014b396bbSKris Buschelman Since SuperLU only likes column-oriented matrices,we pass it the transpose, 23114b396bbSKris Buschelman and then solve A^T X = B in MatSolve(). 23214b396bbSKris Buschelman */ 23314b396bbSKris Buschelman 23414b396bbSKris Buschelman if ( lu->flg == DIFFERENT_NONZERO_PATTERN){ /* first numerical factorization */ 23514b396bbSKris Buschelman lu->A.Stype = SLU_NC; 23614b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 23714b396bbSKris Buschelman lu->A.Dtype = SLU_Z; 23814b396bbSKris Buschelman #else 23914b396bbSKris Buschelman lu->A.Dtype = SLU_D; 24014b396bbSKris Buschelman #endif 24114b396bbSKris Buschelman lu->A.Mtype = SLU_GE; 24214b396bbSKris Buschelman lu->A.nrow = A->n; 24314b396bbSKris Buschelman lu->A.ncol = A->m; 24414b396bbSKris Buschelman 24514b396bbSKris Buschelman ierr = PetscMalloc(sizeof(NCformat),&lu->store);CHKERRQ(ierr); 24614b396bbSKris Buschelman ierr = PetscMalloc(sizeof(DNformat),&lu->B.Store);CHKERRQ(ierr); 24714b396bbSKris Buschelman } 24814b396bbSKris Buschelman lu->store->nnz = aa->nz; 24914b396bbSKris Buschelman lu->store->colptr = aa->i; 25014b396bbSKris Buschelman lu->store->rowind = aa->j; 25114b396bbSKris Buschelman lu->store->nzval = aa->a; 25214b396bbSKris Buschelman lu->A.Store = lu->store; 25314b396bbSKris Buschelman 25414b396bbSKris Buschelman /* Set SuperLU options */ 25514b396bbSKris Buschelman lu->relax = sp_ienv(2); 25614b396bbSKris Buschelman lu->panel_size = sp_ienv(1); 25714b396bbSKris Buschelman /* We have to initialize global data or SuperLU crashes (sucky design) */ 25814b396bbSKris Buschelman if (!StatInitCalled) { 25914b396bbSKris Buschelman StatInit(lu->panel_size,lu->relax); 26014b396bbSKris Buschelman } 26114b396bbSKris Buschelman StatInitCalled++; 26214b396bbSKris Buschelman 26314b396bbSKris Buschelman ierr = PetscOptionsBegin(A->comm,A->prefix,"SuperLU Options","Mat");CHKERRQ(ierr); 26414b396bbSKris Buschelman /* use SuperLU mat ordeing */ 26514b396bbSKris 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); 26614b396bbSKris Buschelman if (flag) { 26714b396bbSKris Buschelman get_perm_c(lu->ispec, &lu->A, lu->perm_c); 26814b396bbSKris Buschelman lu->SuperluMatOdering = PETSC_TRUE; 26914b396bbSKris Buschelman } 27014b396bbSKris Buschelman PetscOptionsEnd(); 27114b396bbSKris Buschelman 27214b396bbSKris Buschelman /* Create the elimination tree */ 27314b396bbSKris Buschelman ierr = PetscMalloc(A->n*sizeof(int),&etree);CHKERRQ(ierr); 27414b396bbSKris Buschelman sp_preorder("N",&lu->A,lu->perm_c,etree,&lu->AC); 27514b396bbSKris Buschelman /* Factor the matrix */ 27614b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 27714b396bbSKris 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); 27814b396bbSKris Buschelman #else 27914b396bbSKris 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); 28014b396bbSKris Buschelman #endif 28114b396bbSKris Buschelman if (ierr < 0) { 28214b396bbSKris Buschelman SETERRQ1(PETSC_ERR_ARG_WRONG,"The diagonal element of row %d was invalid",-ierr); 28314b396bbSKris Buschelman } else if (ierr > 0) { 28414b396bbSKris Buschelman if (ierr <= A->m) { 28514b396bbSKris Buschelman SETERRQ1(PETSC_ERR_ARG_WRONG,"The diagonal element %d of U is exactly zero",ierr); 28614b396bbSKris Buschelman } else { 28714b396bbSKris Buschelman SETERRQ1(PETSC_ERR_ARG_WRONG,"Memory allocation failure after %d bytes were allocated",ierr-A->m); 28814b396bbSKris Buschelman } 28914b396bbSKris Buschelman } 29014b396bbSKris Buschelman 29114b396bbSKris Buschelman /* Cleanup */ 29214b396bbSKris Buschelman ierr = PetscFree(etree);CHKERRQ(ierr); 29314b396bbSKris Buschelman 29414b396bbSKris Buschelman lu->flg = SAME_NONZERO_PATTERN; 29514b396bbSKris Buschelman PetscFunctionReturn(0); 29614b396bbSKris Buschelman } 29714b396bbSKris Buschelman 29814b396bbSKris Buschelman /* 29914b396bbSKris Buschelman Note the r permutation is ignored 30014b396bbSKris Buschelman */ 30114b396bbSKris Buschelman #undef __FUNCT__ 302f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_SuperLU" 303f0c56d0fSKris Buschelman int MatLUFactorSymbolic_SuperLU(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) 30414b396bbSKris Buschelman { 30514b396bbSKris Buschelman Mat B; 306f0c56d0fSKris Buschelman Mat_SuperLU *lu; 30714b396bbSKris Buschelman int ierr,*ca; 30814b396bbSKris Buschelman 30914b396bbSKris Buschelman PetscFunctionBegin; 31014b396bbSKris Buschelman 311899d7b4fSKris Buschelman ierr = MatCreate(A->comm,A->m,A->n,PETSC_DETERMINE,PETSC_DETERMINE,&B);CHKERRQ(ierr); 312899d7b4fSKris Buschelman ierr = MatSetType(B,MATSUPERLU);CHKERRQ(ierr); 313899d7b4fSKris Buschelman ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 314f0c56d0fSKris Buschelman 315f0c56d0fSKris Buschelman B->ops->lufactornumeric = MatLUFactorNumeric_SuperLU; 316f0c56d0fSKris Buschelman B->ops->solve = MatSolve_SuperLU; 31714b396bbSKris Buschelman B->factor = FACTOR_LU; 318*94b7f48cSBarry Smith B->assembled = PETSC_TRUE; /* required by -ksp_view */ 31914b396bbSKris Buschelman 320f0c56d0fSKris Buschelman lu = (Mat_SuperLU*)(B->spptr); 3212ecb5974SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatCreateNull","MatCreateNull_SuperLU", 322f0c56d0fSKris Buschelman (void(*)(void))MatCreateNull_SuperLU);CHKERRQ(ierr); 32314b396bbSKris Buschelman 32414b396bbSKris Buschelman /* Allocate the work arrays required by SuperLU (notice sizes are for the transpose) */ 32514b396bbSKris Buschelman ierr = PetscMalloc(A->n*sizeof(int),&lu->perm_r);CHKERRQ(ierr); 32614b396bbSKris Buschelman ierr = PetscMalloc(A->m*sizeof(int),&lu->perm_c);CHKERRQ(ierr); 32714b396bbSKris Buschelman 32814b396bbSKris Buschelman /* use PETSc mat ordering */ 32914b396bbSKris Buschelman ierr = ISGetIndices(c,&ca);CHKERRQ(ierr); 33014b396bbSKris Buschelman ierr = PetscMemcpy(lu->perm_c,ca,A->m*sizeof(int));CHKERRQ(ierr); 33114b396bbSKris Buschelman ierr = ISRestoreIndices(c,&ca);CHKERRQ(ierr); 33214b396bbSKris Buschelman lu->SuperluMatOdering = PETSC_FALSE; 33314b396bbSKris Buschelman 33414b396bbSKris Buschelman lu->pivot_threshold = info->dtcol; 335f0c56d0fSKris Buschelman PetscLogObjectMemory(B,(A->m+A->n)*sizeof(int)+sizeof(Mat_SuperLU)); 33614b396bbSKris Buschelman 33714b396bbSKris Buschelman lu->flg = DIFFERENT_NONZERO_PATTERN; 338e740cb95SKris Buschelman lu->CleanUpSuperLU = PETSC_TRUE; 339899d7b4fSKris Buschelman 340899d7b4fSKris Buschelman *F = B; 34114b396bbSKris Buschelman PetscFunctionReturn(0); 34214b396bbSKris Buschelman } 34314b396bbSKris Buschelman 344*94b7f48cSBarry Smith /* used by -ksp_view */ 34514b396bbSKris Buschelman #undef __FUNCT__ 346f0c56d0fSKris Buschelman #define __FUNCT__ "MatFactorInfo_SuperLU" 347f0c56d0fSKris Buschelman int MatFactorInfo_SuperLU(Mat A,PetscViewer viewer) 34814b396bbSKris Buschelman { 349f0c56d0fSKris Buschelman Mat_SuperLU *lu= (Mat_SuperLU*)A->spptr; 350f0c56d0fSKris Buschelman int ierr; 351f0c56d0fSKris Buschelman 352f0c56d0fSKris Buschelman PetscFunctionBegin; 353f0c56d0fSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");CHKERRQ(ierr); 354f0c56d0fSKris Buschelman if(lu->SuperluMatOdering) { 355f0c56d0fSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," SuperLU mat ordering: %d\n",lu->ispec);CHKERRQ(ierr); 356f0c56d0fSKris Buschelman } 357f0c56d0fSKris Buschelman PetscFunctionReturn(0); 358f0c56d0fSKris Buschelman } 359f0c56d0fSKris Buschelman 360f0c56d0fSKris Buschelman #undef __FUNCT__ 361f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_SuperLU" 362f0c56d0fSKris Buschelman int MatDuplicate_SuperLU(Mat A, MatDuplicateOption op, Mat *M) { 36314b396bbSKris Buschelman int ierr; 3648f340917SKris Buschelman Mat_SuperLU *lu=(Mat_SuperLU *)A->spptr; 3658f340917SKris Buschelman 36614b396bbSKris Buschelman PetscFunctionBegin; 3678f340917SKris Buschelman ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr); 368f0c56d0fSKris Buschelman ierr = MatConvert_SeqAIJ_SuperLU(*M,MATSUPERLU,M);CHKERRQ(ierr); 3693f953163SKris Buschelman ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_SuperLU));CHKERRQ(ierr); 37014b396bbSKris Buschelman PetscFunctionReturn(0); 37114b396bbSKris Buschelman } 37214b396bbSKris Buschelman 37314b396bbSKris Buschelman EXTERN_C_BEGIN 37414b396bbSKris Buschelman #undef __FUNCT__ 375b0592601SKris Buschelman #define __FUNCT__ "MatConvert_SuperLU_SeqAIJ" 376b0592601SKris Buschelman int MatConvert_SuperLU_SeqAIJ(Mat A,MatType type,Mat *newmat) { 377fb3e25aaSKris Buschelman /* This routine is only called to convert an unfactored PETSc-SuperLU matrix */ 378fb3e25aaSKris Buschelman /* to its base PETSc type, so we will ignore 'MatType type'. */ 37914b396bbSKris Buschelman int ierr; 380b0592601SKris Buschelman Mat B=*newmat; 381f0c56d0fSKris Buschelman Mat_SuperLU *lu=(Mat_SuperLU *)A->spptr; 382b0592601SKris Buschelman 383b0592601SKris Buschelman PetscFunctionBegin; 384b0592601SKris Buschelman if (B != A) { 385b0592601SKris Buschelman ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 386f0c56d0fSKris Buschelman } 387b0592601SKris Buschelman /* Reset the original function pointers */ 388f0c56d0fSKris Buschelman B->ops->duplicate = lu->MatDuplicate; 389b0592601SKris Buschelman B->ops->view = lu->MatView; 390b0592601SKris Buschelman B->ops->assemblyend = lu->MatAssemblyEnd; 391b0592601SKris Buschelman B->ops->lufactorsymbolic = lu->MatLUFactorSymbolic; 392b0592601SKris Buschelman B->ops->destroy = lu->MatDestroy; 393b0592601SKris Buschelman /* lu is only a function pointer stash unless we've factored the matrix, which we haven't! */ 394b0592601SKris Buschelman ierr = PetscFree(lu);CHKERRQ(ierr); 395b0592601SKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);CHKERRQ(ierr); 396b0592601SKris Buschelman *newmat = B; 397b0592601SKris Buschelman PetscFunctionReturn(0); 398b0592601SKris Buschelman } 399b0592601SKris Buschelman EXTERN_C_END 400b0592601SKris Buschelman 401b0592601SKris Buschelman EXTERN_C_BEGIN 402b0592601SKris Buschelman #undef __FUNCT__ 403b0592601SKris Buschelman #define __FUNCT__ "MatConvert_SeqAIJ_SuperLU" 404b0592601SKris Buschelman int MatConvert_SeqAIJ_SuperLU(Mat A,MatType type,Mat *newmat) { 405fb3e25aaSKris Buschelman /* This routine is only called to convert to MATSUPERLU */ 406fb3e25aaSKris Buschelman /* from MATSEQAIJ, so we will ignore 'MatType type'. */ 407b0592601SKris Buschelman int ierr; 408b0592601SKris Buschelman Mat B=*newmat; 409f0c56d0fSKris Buschelman Mat_SuperLU *lu; 41014b396bbSKris Buschelman 41114b396bbSKris Buschelman PetscFunctionBegin; 412b0592601SKris Buschelman if (B != A) { 413b0592601SKris Buschelman ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 414b0592601SKris Buschelman } 41514b396bbSKris Buschelman 416f0c56d0fSKris Buschelman ierr = PetscNew(Mat_SuperLU,&lu);CHKERRQ(ierr); 417f0c56d0fSKris Buschelman lu->MatDuplicate = A->ops->duplicate; 41814b396bbSKris Buschelman lu->MatView = A->ops->view; 41914b396bbSKris Buschelman lu->MatAssemblyEnd = A->ops->assemblyend; 420b0592601SKris Buschelman lu->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 42114b396bbSKris Buschelman lu->MatDestroy = A->ops->destroy; 42214b396bbSKris Buschelman lu->CleanUpSuperLU = PETSC_FALSE; 423b0592601SKris Buschelman 424b0592601SKris Buschelman B->spptr = (void*)lu; 425f0c56d0fSKris Buschelman B->ops->duplicate = MatDuplicate_SuperLU; 426f0c56d0fSKris Buschelman B->ops->view = MatView_SuperLU; 427f0c56d0fSKris Buschelman B->ops->assemblyend = MatAssemblyEnd_SuperLU; 428f0c56d0fSKris Buschelman B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU; 429f0c56d0fSKris Buschelman B->ops->destroy = MatDestroy_SuperLU; 430b0592601SKris Buschelman 431b0592601SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_superlu_C", 432b0592601SKris Buschelman "MatConvert_SeqAIJ_SuperLU",MatConvert_SeqAIJ_SuperLU);CHKERRQ(ierr); 433b0592601SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_superlu_seqaij_C", 434b0592601SKris Buschelman "MatConvert_SuperLU_SeqAIJ",MatConvert_SuperLU_SeqAIJ);CHKERRQ(ierr); 435fb3e25aaSKris Buschelman PetscLogInfo(0,"Using SuperLU for SeqAIJ LU factorization and solves."); 436fb3e25aaSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)B,MATSUPERLU);CHKERRQ(ierr); 437b0592601SKris Buschelman *newmat = B; 438b0592601SKris Buschelman PetscFunctionReturn(0); 439b0592601SKris Buschelman } 440b0592601SKris Buschelman EXTERN_C_END 441b0592601SKris Buschelman 44224b6179bSKris Buschelman /*MC 443fafad747SKris Buschelman MATSUPERLU - MATSUPERLU = "superlu" - A matrix type providing direct solvers (LU) for sequential matrices 44424b6179bSKris Buschelman via the external package SuperLU. 44524b6179bSKris Buschelman 44624b6179bSKris Buschelman If SuperLU is installed (see the manual for 44724b6179bSKris Buschelman instructions on how to declare the existence of external packages), 44824b6179bSKris Buschelman a matrix type can be constructed which invokes SuperLU solvers. 44924b6179bSKris Buschelman After calling MatCreate(...,A), simply call MatSetType(A,MATSUPERLU). 45024b6179bSKris Buschelman This matrix type is only supported for double precision real. 45124b6179bSKris Buschelman 45224b6179bSKris Buschelman This matrix inherits from MATSEQAIJ. As a result, MatSeqAIJSetPreallocation is 45328b08bd3SKris Buschelman supported for this matrix type. One can also call MatConvert for an inplace conversion to or from 45428b08bd3SKris Buschelman the MATSEQAIJ type without data copy. 45524b6179bSKris Buschelman 45624b6179bSKris Buschelman Options Database Keys: 4570bad9183SKris Buschelman + -mat_type superlu - sets the matrix type to "superlu" during a call to MatSetFromOptions() 45824b6179bSKris Buschelman - -mat_superlu_ordering <0,1,2,3> - 0: natural ordering, 45924b6179bSKris Buschelman 1: MMD applied to A'*A, 46024b6179bSKris Buschelman 2: MMD applied to A'+A, 46124b6179bSKris Buschelman 3: COLAMD, approximate minimum degree column ordering 46224b6179bSKris Buschelman 46324b6179bSKris Buschelman Level: beginner 46424b6179bSKris Buschelman 46524b6179bSKris Buschelman .seealso: PCLU 46624b6179bSKris Buschelman M*/ 46724b6179bSKris Buschelman 468b0592601SKris Buschelman EXTERN_C_BEGIN 469b0592601SKris Buschelman #undef __FUNCT__ 470f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_SuperLU" 471f0c56d0fSKris Buschelman int MatCreate_SuperLU(Mat A) { 472b0592601SKris Buschelman int ierr; 473b0592601SKris Buschelman 474b0592601SKris Buschelman PetscFunctionBegin; 4755441df8eSKris Buschelman /* Change type name before calling MatSetType to force proper construction of SeqAIJ and SUPERLU types */ 4765441df8eSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)A,MATSUPERLU);CHKERRQ(ierr); 477b0592601SKris Buschelman ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr); 478b0592601SKris Buschelman ierr = MatConvert_SeqAIJ_SuperLU(A,MATSUPERLU,&A);CHKERRQ(ierr); 47914b396bbSKris Buschelman PetscFunctionReturn(0); 48014b396bbSKris Buschelman } 48114b396bbSKris Buschelman EXTERN_C_END 482