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); 30*b0592601SKris 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); 39*b0592601SKris Buschelman EXTERN int MatLUFactorSymbolic_SeqAIJ_SuperLU(Mat,IS,IS,MatFactorInfo*,Mat*); 40*b0592601SKris Buschelman 41*b0592601SKris Buschelman EXTERN_C_BEGIN 42*b0592601SKris Buschelman EXTERN int MatConvert_SuperLU_SeqAIJ(Mat,MatType,Mat*); 43*b0592601SKris Buschelman EXTERN int MatConvert_SeqAIJ_SuperLU(Mat,MatType,Mat*); 44*b0592601SKris 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 { 5014b396bbSKris Buschelman Mat_SeqAIJ_SuperLU *lu = (Mat_SeqAIJ_SuperLU*)A->spptr; 5114b396bbSKris Buschelman int ierr,(*destroy)(Mat); 5214b396bbSKris Buschelman 5314b396bbSKris Buschelman PetscFunctionBegin; 5414b396bbSKris Buschelman /* It looks like this is decreasing the reference count a second time during MatDestroy?! */ 55*b0592601SKris Buschelman /* if (--A->refct > 0)PetscFunctionReturn(0); */ 5614b396bbSKris Buschelman destroy = lu->MatDestroy; 57*b0592601SKris Buschelman ierr = MatConvert_SuperLU_SeqAIJ(A,MATSEQAIJ,&A);CHKERRQ(ierr); 5814b396bbSKris Buschelman ierr = (*destroy)(A);CHKERRQ(ierr); 5914b396bbSKris Buschelman PetscFunctionReturn(0); 6014b396bbSKris Buschelman } 6114b396bbSKris Buschelman 6214b396bbSKris Buschelman #undef __FUNCT__ 6314b396bbSKris Buschelman #define __FUNCT__ "MatView_SeqAIJ_Spooles" 6414b396bbSKris Buschelman int MatView_SeqAIJ_SuperLU(Mat A,PetscViewer viewer) 6514b396bbSKris Buschelman { 6614b396bbSKris Buschelman int ierr; 6714b396bbSKris Buschelman PetscTruth isascii; 6814b396bbSKris Buschelman PetscViewerFormat format; 6914b396bbSKris Buschelman Mat_SeqAIJ_SuperLU *lu=(Mat_SeqAIJ_SuperLU*)(A->spptr); 7014b396bbSKris Buschelman 7114b396bbSKris Buschelman PetscFunctionBegin; 7214b396bbSKris Buschelman ierr = (*lu->MatView)(A,viewer);CHKERRQ(ierr); 7314b396bbSKris Buschelman 7414b396bbSKris Buschelman ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 7514b396bbSKris Buschelman if (isascii) { 7614b396bbSKris Buschelman ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 7714b396bbSKris Buschelman if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) { 7814b396bbSKris Buschelman ierr = MatSeqAIJFactorInfo_SuperLU(A,viewer);CHKERRQ(ierr); 7914b396bbSKris Buschelman } 8014b396bbSKris Buschelman } 8114b396bbSKris Buschelman PetscFunctionReturn(0); 8214b396bbSKris Buschelman } 8314b396bbSKris Buschelman 8414b396bbSKris Buschelman #undef __FUNCT__ 8514b396bbSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_SeqAIJ_SuperLU" 8614b396bbSKris Buschelman int MatAssemblyEnd_SeqAIJ_SuperLU(Mat A,MatAssemblyType mode) { 8714b396bbSKris Buschelman int ierr; 8814b396bbSKris Buschelman Mat_SeqAIJ_SuperLU *lu=(Mat_SeqAIJ_SuperLU*)(A->spptr); 8914b396bbSKris Buschelman 9014b396bbSKris Buschelman PetscFunctionBegin; 9114b396bbSKris Buschelman ierr = (*lu->MatAssemblyEnd)(A,mode);CHKERRQ(ierr); 92*b0592601SKris Buschelman 93*b0592601SKris Buschelman lu->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 94*b0592601SKris Buschelman A->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ_SuperLU; 9514b396bbSKris Buschelman PetscFunctionReturn(0); 9614b396bbSKris Buschelman } 9714b396bbSKris Buschelman 9814b396bbSKris Buschelman #include "src/mat/impls/dense/seq/dense.h" 9914b396bbSKris Buschelman #undef __FUNCT__ 10014b396bbSKris Buschelman #define __FUNCT__ "MatCreateNull_SeqAIJ_SuperLU" 10114b396bbSKris Buschelman int MatCreateNull_SeqAIJ_SuperLU(Mat A,Mat *nullMat) 10214b396bbSKris Buschelman { 10314b396bbSKris Buschelman Mat_SeqAIJ_SuperLU *lu = (Mat_SeqAIJ_SuperLU*)A->spptr; 10414b396bbSKris Buschelman int numRows = A->m,numCols = A->n; 10514b396bbSKris Buschelman SCformat *Lstore; 10614b396bbSKris Buschelman int numNullCols,size; 10714b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 10814b396bbSKris Buschelman doublecomplex *nullVals,*workVals; 10914b396bbSKris Buschelman #else 11014b396bbSKris Buschelman PetscScalar *nullVals,*workVals; 11114b396bbSKris Buschelman #endif 11214b396bbSKris Buschelman int row,newRow,col,newCol,block,b,ierr; 11314b396bbSKris Buschelman 11414b396bbSKris Buschelman PetscFunctionBegin; 11514b396bbSKris Buschelman PetscValidHeaderSpecific(A,MAT_COOKIE); 11614b396bbSKris Buschelman PetscValidPointer(nullMat); 11714b396bbSKris Buschelman if (!A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix"); 11814b396bbSKris Buschelman numNullCols = numCols - numRows; 11914b396bbSKris Buschelman if (numNullCols < 0) SETERRQ(PETSC_ERR_ARG_WRONG,"Function only applies to underdetermined problems"); 12014b396bbSKris Buschelman /* Create the null matrix */ 12114b396bbSKris Buschelman ierr = MatCreateSeqDense(A->comm,numRows,numNullCols,PETSC_NULL,nullMat);CHKERRQ(ierr); 12214b396bbSKris Buschelman if (numNullCols == 0) { 12314b396bbSKris Buschelman ierr = MatAssemblyBegin(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 12414b396bbSKris Buschelman ierr = MatAssemblyEnd(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 12514b396bbSKris Buschelman PetscFunctionReturn(0); 12614b396bbSKris Buschelman } 12714b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 12814b396bbSKris Buschelman nullVals = (doublecomplex*)((Mat_SeqDense*)(*nullMat)->data)->v; 12914b396bbSKris Buschelman #else 13014b396bbSKris Buschelman nullVals = ((Mat_SeqDense*)(*nullMat)->data)->v; 13114b396bbSKris Buschelman #endif 13214b396bbSKris Buschelman /* Copy in the columns */ 13314b396bbSKris Buschelman Lstore = (SCformat*)lu->L.Store; 13414b396bbSKris Buschelman for(block = 0; block <= Lstore->nsuper; block++) { 13514b396bbSKris Buschelman newRow = Lstore->sup_to_col[block]; 13614b396bbSKris Buschelman size = Lstore->sup_to_col[block+1] - Lstore->sup_to_col[block]; 13714b396bbSKris Buschelman for(col = Lstore->rowind_colptr[newRow]; col < Lstore->rowind_colptr[newRow+1]; col++) { 13814b396bbSKris Buschelman newCol = Lstore->rowind[col]; 13914b396bbSKris Buschelman if (newCol >= numRows) { 14014b396bbSKris Buschelman for(b = 0; b < size; b++) 14114b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 14214b396bbSKris Buschelman nullVals[(newCol-numRows)*numRows+newRow+b] = ((doublecomplex*)Lstore->nzval)[Lstore->nzval_colptr[newRow+b]+col]; 14314b396bbSKris Buschelman #else 14414b396bbSKris Buschelman nullVals[(newCol-numRows)*numRows+newRow+b] = ((double*)Lstore->nzval)[Lstore->nzval_colptr[newRow+b]+col]; 14514b396bbSKris Buschelman #endif 14614b396bbSKris Buschelman } 14714b396bbSKris Buschelman } 14814b396bbSKris Buschelman } 14914b396bbSKris Buschelman /* Permute rhs to form P^T_c B */ 15014b396bbSKris Buschelman ierr = PetscMalloc(numRows*sizeof(double),&workVals);CHKERRQ(ierr); 15114b396bbSKris Buschelman for(b = 0; b < numNullCols; b++) { 15214b396bbSKris Buschelman for(row = 0; row < numRows; row++) workVals[lu->perm_c[row]] = nullVals[b*numRows+row]; 15314b396bbSKris Buschelman for(row = 0; row < numRows; row++) nullVals[b*numRows+row] = workVals[row]; 15414b396bbSKris Buschelman } 15514b396bbSKris Buschelman /* Backward solve the upper triangle A x = b */ 15614b396bbSKris Buschelman for(b = 0; b < numNullCols; b++) { 15714b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 15814b396bbSKris Buschelman sp_ztrsv("L","T","U",&lu->L,&lu->U,&nullVals[b*numRows],&ierr); 15914b396bbSKris Buschelman #else 16014b396bbSKris Buschelman sp_dtrsv("L","T","U",&lu->L,&lu->U,&nullVals[b*numRows],&ierr); 16114b396bbSKris Buschelman #endif 16214b396bbSKris Buschelman if (ierr < 0) 16314b396bbSKris Buschelman SETERRQ1(PETSC_ERR_ARG_WRONG,"The argument %d was invalid",-ierr); 16414b396bbSKris Buschelman } 16514b396bbSKris Buschelman ierr = PetscFree(workVals);CHKERRQ(ierr); 16614b396bbSKris Buschelman 16714b396bbSKris Buschelman ierr = MatAssemblyBegin(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 16814b396bbSKris Buschelman ierr = MatAssemblyEnd(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 16914b396bbSKris Buschelman PetscFunctionReturn(0); 17014b396bbSKris Buschelman } 17114b396bbSKris Buschelman 17214b396bbSKris Buschelman #undef __FUNCT__ 17314b396bbSKris Buschelman #define __FUNCT__ "MatSolve_SeqAIJ_SuperLU" 17414b396bbSKris Buschelman int MatSolve_SeqAIJ_SuperLU(Mat A,Vec b,Vec x) 17514b396bbSKris Buschelman { 17614b396bbSKris Buschelman Mat_SeqAIJ_SuperLU *lu = (Mat_SeqAIJ_SuperLU*)A->spptr; 17714b396bbSKris Buschelman PetscScalar *array; 17814b396bbSKris Buschelman int m,ierr; 17914b396bbSKris Buschelman 18014b396bbSKris Buschelman PetscFunctionBegin; 18114b396bbSKris Buschelman ierr = VecGetLocalSize(b,&m);CHKERRQ(ierr); 18214b396bbSKris Buschelman ierr = VecCopy(b,x);CHKERRQ(ierr); 18314b396bbSKris Buschelman ierr = VecGetArray(x,&array);CHKERRQ(ierr); 18414b396bbSKris Buschelman /* Create the Rhs */ 18514b396bbSKris Buschelman lu->B.Stype = SLU_DN; 18614b396bbSKris Buschelman lu->B.Mtype = SLU_GE; 18714b396bbSKris Buschelman lu->B.nrow = m; 18814b396bbSKris Buschelman lu->B.ncol = 1; 18914b396bbSKris Buschelman ((DNformat*)lu->B.Store)->lda = m; 19014b396bbSKris Buschelman ((DNformat*)lu->B.Store)->nzval = array; 19114b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 19214b396bbSKris Buschelman lu->B.Dtype = SLU_Z; 19314b396bbSKris Buschelman zgstrs("T",&lu->L,&lu->U,lu->perm_r,lu->perm_c,&lu->B,&ierr); 19414b396bbSKris Buschelman #else 19514b396bbSKris Buschelman lu->B.Dtype = SLU_D; 19614b396bbSKris Buschelman dgstrs("T",&lu->L,&lu->U,lu->perm_r,lu->perm_c,&lu->B,&ierr); 19714b396bbSKris Buschelman #endif 19814b396bbSKris Buschelman if (ierr < 0) SETERRQ1(PETSC_ERR_ARG_WRONG,"The diagonal element of row %d was invalid",-ierr); 19914b396bbSKris Buschelman ierr = VecRestoreArray(x,&array);CHKERRQ(ierr); 20014b396bbSKris Buschelman PetscFunctionReturn(0); 20114b396bbSKris Buschelman } 20214b396bbSKris Buschelman 20314b396bbSKris Buschelman static int StatInitCalled = 0; 20414b396bbSKris Buschelman 20514b396bbSKris Buschelman #undef __FUNCT__ 20614b396bbSKris Buschelman #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ_SuperLU" 20714b396bbSKris Buschelman int MatLUFactorNumeric_SeqAIJ_SuperLU(Mat A,Mat *F) 20814b396bbSKris Buschelman { 20914b396bbSKris Buschelman Mat_SeqAIJ *aa = (Mat_SeqAIJ*)(A)->data; 21014b396bbSKris Buschelman Mat_SeqAIJ_SuperLU *lu = (Mat_SeqAIJ_SuperLU*)(*F)->spptr; 211d54de34fSKris Buschelman int *etree,ierr; 21214b396bbSKris Buschelman PetscTruth flag; 21314b396bbSKris Buschelman 21414b396bbSKris Buschelman PetscFunctionBegin; 21514b396bbSKris Buschelman /* Create the SuperMatrix for A^T: 21614b396bbSKris Buschelman Since SuperLU only likes column-oriented matrices,we pass it the transpose, 21714b396bbSKris Buschelman and then solve A^T X = B in MatSolve(). 21814b396bbSKris Buschelman */ 21914b396bbSKris Buschelman 22014b396bbSKris Buschelman if ( lu->flg == DIFFERENT_NONZERO_PATTERN){ /* first numerical factorization */ 22114b396bbSKris Buschelman lu->A.Stype = SLU_NC; 22214b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 22314b396bbSKris Buschelman lu->A.Dtype = SLU_Z; 22414b396bbSKris Buschelman #else 22514b396bbSKris Buschelman lu->A.Dtype = SLU_D; 22614b396bbSKris Buschelman #endif 22714b396bbSKris Buschelman lu->A.Mtype = SLU_GE; 22814b396bbSKris Buschelman lu->A.nrow = A->n; 22914b396bbSKris Buschelman lu->A.ncol = A->m; 23014b396bbSKris Buschelman 23114b396bbSKris Buschelman ierr = PetscMalloc(sizeof(NCformat),&lu->store);CHKERRQ(ierr); 23214b396bbSKris Buschelman ierr = PetscMalloc(sizeof(DNformat),&lu->B.Store);CHKERRQ(ierr); 23314b396bbSKris Buschelman } 23414b396bbSKris Buschelman lu->store->nnz = aa->nz; 23514b396bbSKris Buschelman lu->store->colptr = aa->i; 23614b396bbSKris Buschelman lu->store->rowind = aa->j; 23714b396bbSKris Buschelman lu->store->nzval = aa->a; 23814b396bbSKris Buschelman lu->A.Store = lu->store; 23914b396bbSKris Buschelman 24014b396bbSKris Buschelman /* Set SuperLU options */ 24114b396bbSKris Buschelman lu->relax = sp_ienv(2); 24214b396bbSKris Buschelman lu->panel_size = sp_ienv(1); 24314b396bbSKris Buschelman /* We have to initialize global data or SuperLU crashes (sucky design) */ 24414b396bbSKris Buschelman if (!StatInitCalled) { 24514b396bbSKris Buschelman StatInit(lu->panel_size,lu->relax); 24614b396bbSKris Buschelman } 24714b396bbSKris Buschelman StatInitCalled++; 24814b396bbSKris Buschelman 24914b396bbSKris Buschelman ierr = PetscOptionsBegin(A->comm,A->prefix,"SuperLU Options","Mat");CHKERRQ(ierr); 25014b396bbSKris Buschelman /* use SuperLU mat ordeing */ 25114b396bbSKris 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); 25214b396bbSKris Buschelman if (flag) { 25314b396bbSKris Buschelman get_perm_c(lu->ispec, &lu->A, lu->perm_c); 25414b396bbSKris Buschelman lu->SuperluMatOdering = PETSC_TRUE; 25514b396bbSKris Buschelman } 25614b396bbSKris Buschelman PetscOptionsEnd(); 25714b396bbSKris Buschelman 25814b396bbSKris Buschelman /* Create the elimination tree */ 25914b396bbSKris Buschelman ierr = PetscMalloc(A->n*sizeof(int),&etree);CHKERRQ(ierr); 26014b396bbSKris Buschelman sp_preorder("N",&lu->A,lu->perm_c,etree,&lu->AC); 26114b396bbSKris Buschelman /* Factor the matrix */ 26214b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 26314b396bbSKris 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); 26414b396bbSKris Buschelman #else 26514b396bbSKris 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); 26614b396bbSKris Buschelman #endif 26714b396bbSKris Buschelman if (ierr < 0) { 26814b396bbSKris Buschelman SETERRQ1(PETSC_ERR_ARG_WRONG,"The diagonal element of row %d was invalid",-ierr); 26914b396bbSKris Buschelman } else if (ierr > 0) { 27014b396bbSKris Buschelman if (ierr <= A->m) { 27114b396bbSKris Buschelman SETERRQ1(PETSC_ERR_ARG_WRONG,"The diagonal element %d of U is exactly zero",ierr); 27214b396bbSKris Buschelman } else { 27314b396bbSKris Buschelman SETERRQ1(PETSC_ERR_ARG_WRONG,"Memory allocation failure after %d bytes were allocated",ierr-A->m); 27414b396bbSKris Buschelman } 27514b396bbSKris Buschelman } 27614b396bbSKris Buschelman 27714b396bbSKris Buschelman /* Cleanup */ 27814b396bbSKris Buschelman ierr = PetscFree(etree);CHKERRQ(ierr); 27914b396bbSKris Buschelman 28014b396bbSKris Buschelman lu->flg = SAME_NONZERO_PATTERN; 28114b396bbSKris Buschelman PetscFunctionReturn(0); 28214b396bbSKris Buschelman } 28314b396bbSKris Buschelman 28414b396bbSKris Buschelman /* 28514b396bbSKris Buschelman Note the r permutation is ignored 28614b396bbSKris Buschelman */ 28714b396bbSKris Buschelman #undef __FUNCT__ 28814b396bbSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJ_SuperLU" 28914b396bbSKris Buschelman int MatLUFactorSymbolic_SeqAIJ_SuperLU(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) 29014b396bbSKris Buschelman { 29114b396bbSKris Buschelman Mat B; 29214b396bbSKris Buschelman Mat_SeqAIJ_SuperLU *lu; 29314b396bbSKris Buschelman int ierr,*ca; 29414b396bbSKris Buschelman 29514b396bbSKris Buschelman PetscFunctionBegin; 29614b396bbSKris Buschelman 297899d7b4fSKris Buschelman ierr = MatCreate(A->comm,A->m,A->n,PETSC_DETERMINE,PETSC_DETERMINE,&B);CHKERRQ(ierr); 298899d7b4fSKris Buschelman ierr = MatSetType(B,MATSUPERLU);CHKERRQ(ierr); 299899d7b4fSKris Buschelman ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 30014b396bbSKris Buschelman B->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJ_SuperLU; 30114b396bbSKris Buschelman B->ops->solve = MatSolve_SeqAIJ_SuperLU; 30214b396bbSKris Buschelman B->factor = FACTOR_LU; 303899d7b4fSKris Buschelman B->assembled = PETSC_TRUE; /* required by -sles_view */ 30414b396bbSKris Buschelman 305899d7b4fSKris Buschelman lu = (Mat_SeqAIJ_SuperLU*)(B->spptr); 30614b396bbSKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)B,"MatCreateNull","MatCreateNull_SeqAIJ_SuperLU", 30714b396bbSKris Buschelman (void(*)(void))MatCreateNull_SeqAIJ_SuperLU);CHKERRQ(ierr); 30814b396bbSKris Buschelman 30914b396bbSKris Buschelman /* Allocate the work arrays required by SuperLU (notice sizes are for the transpose) */ 31014b396bbSKris Buschelman ierr = PetscMalloc(A->n*sizeof(int),&lu->perm_r);CHKERRQ(ierr); 31114b396bbSKris Buschelman ierr = PetscMalloc(A->m*sizeof(int),&lu->perm_c);CHKERRQ(ierr); 31214b396bbSKris Buschelman 31314b396bbSKris Buschelman /* use PETSc mat ordering */ 31414b396bbSKris Buschelman ierr = ISGetIndices(c,&ca);CHKERRQ(ierr); 31514b396bbSKris Buschelman ierr = PetscMemcpy(lu->perm_c,ca,A->m*sizeof(int));CHKERRQ(ierr); 31614b396bbSKris Buschelman ierr = ISRestoreIndices(c,&ca);CHKERRQ(ierr); 31714b396bbSKris Buschelman lu->SuperluMatOdering = PETSC_FALSE; 31814b396bbSKris Buschelman 31914b396bbSKris Buschelman lu->pivot_threshold = info->dtcol; 32014b396bbSKris Buschelman PetscLogObjectMemory(B,(A->m+A->n)*sizeof(int)+sizeof(Mat_SeqAIJ_SuperLU)); 32114b396bbSKris Buschelman 32214b396bbSKris Buschelman lu->flg = DIFFERENT_NONZERO_PATTERN; 323e740cb95SKris Buschelman lu->CleanUpSuperLU = PETSC_TRUE; 324899d7b4fSKris Buschelman 325899d7b4fSKris Buschelman *F = B; 32614b396bbSKris Buschelman PetscFunctionReturn(0); 32714b396bbSKris Buschelman } 32814b396bbSKris Buschelman 32914b396bbSKris Buschelman /* used by -sles_view */ 33014b396bbSKris Buschelman #undef __FUNCT__ 33114b396bbSKris Buschelman #define __FUNCT__ "MatSeqAIJFactorInfo_SuperLU" 33214b396bbSKris Buschelman int MatSeqAIJFactorInfo_SuperLU(Mat A,PetscViewer viewer) 33314b396bbSKris Buschelman { 33414b396bbSKris Buschelman Mat_SeqAIJ_SuperLU *lu= (Mat_SeqAIJ_SuperLU*)A->spptr; 33514b396bbSKris Buschelman int ierr; 33614b396bbSKris Buschelman PetscFunctionBegin; 33714b396bbSKris Buschelman 33814b396bbSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");CHKERRQ(ierr); 33914b396bbSKris Buschelman if(lu->SuperluMatOdering) ierr = PetscViewerASCIIPrintf(viewer," SuperLU mat ordering: %d\n",lu->ispec);CHKERRQ(ierr); 34014b396bbSKris Buschelman 34114b396bbSKris Buschelman PetscFunctionReturn(0); 34214b396bbSKris Buschelman } 34314b396bbSKris Buschelman 34414b396bbSKris Buschelman EXTERN_C_BEGIN 34514b396bbSKris Buschelman #undef __FUNCT__ 346*b0592601SKris Buschelman #define __FUNCT__ "MatConvert_SuperLU_SeqAIJ" 347*b0592601SKris Buschelman int MatConvert_SuperLU_SeqAIJ(Mat A,MatType type,Mat *newmat) { 34814b396bbSKris Buschelman int ierr; 349*b0592601SKris Buschelman Mat B=*newmat; 350*b0592601SKris Buschelman Mat_SeqAIJ_SuperLU *lu=(Mat_SeqAIJ_SuperLU *)A->spptr; 351*b0592601SKris Buschelman 352*b0592601SKris Buschelman PetscFunctionBegin; 353*b0592601SKris Buschelman if (B != A) { 354*b0592601SKris Buschelman /* This routine was inherited from SeqAIJ. */ 355*b0592601SKris Buschelman ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 356*b0592601SKris Buschelman } else { 357*b0592601SKris Buschelman if (lu->CleanUpSuperLU) { 358*b0592601SKris Buschelman /* We have to free the global data or SuperLU crashes (sucky design)*/ 359*b0592601SKris Buschelman /* Since we don't know if more solves on other matrices may be done 360*b0592601SKris Buschelman we cannot free the yucky SuperLU global data 361*b0592601SKris Buschelman StatFree(); 362*b0592601SKris Buschelman */ 363*b0592601SKris Buschelman 364*b0592601SKris Buschelman /* Free the SuperLU datastructures */ 365*b0592601SKris Buschelman Destroy_CompCol_Permuted(&lu->AC); 366*b0592601SKris Buschelman Destroy_SuperNode_Matrix(&lu->L); 367*b0592601SKris Buschelman Destroy_CompCol_Matrix(&lu->U); 368*b0592601SKris Buschelman ierr = PetscFree(lu->B.Store);CHKERRQ(ierr); 369*b0592601SKris Buschelman ierr = PetscFree(lu->A.Store);CHKERRQ(ierr); 370*b0592601SKris Buschelman ierr = PetscFree(lu->perm_r);CHKERRQ(ierr); 371*b0592601SKris Buschelman ierr = PetscFree(lu->perm_c);CHKERRQ(ierr); 372*b0592601SKris Buschelman } 373*b0592601SKris Buschelman /* Reset the original function pointers */ 374*b0592601SKris Buschelman B->ops->view = lu->MatView; 375*b0592601SKris Buschelman B->ops->assemblyend = lu->MatAssemblyEnd; 376*b0592601SKris Buschelman B->ops->lufactorsymbolic = lu->MatLUFactorSymbolic; 377*b0592601SKris Buschelman B->ops->destroy = lu->MatDestroy; 378*b0592601SKris Buschelman /* lu is only a function pointer stash unless we've factored the matrix, which we haven't! */ 379*b0592601SKris Buschelman ierr = PetscFree(lu);CHKERRQ(ierr); 380*b0592601SKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);CHKERRQ(ierr); 381*b0592601SKris Buschelman } 382*b0592601SKris Buschelman *newmat = B; 383*b0592601SKris Buschelman PetscFunctionReturn(0); 384*b0592601SKris Buschelman } 385*b0592601SKris Buschelman EXTERN_C_END 386*b0592601SKris Buschelman 387*b0592601SKris Buschelman EXTERN_C_BEGIN 388*b0592601SKris Buschelman #undef __FUNCT__ 389*b0592601SKris Buschelman #define __FUNCT__ "MatConvert_SeqAIJ_SuperLU" 390*b0592601SKris Buschelman int MatConvert_SeqAIJ_SuperLU(Mat A,MatType type,Mat *newmat) { 391*b0592601SKris Buschelman int ierr; 392*b0592601SKris Buschelman Mat B=*newmat; 39314b396bbSKris Buschelman Mat_SeqAIJ_SuperLU *lu; 39414b396bbSKris Buschelman 39514b396bbSKris Buschelman PetscFunctionBegin; 396*b0592601SKris Buschelman if (B != A) { 397*b0592601SKris Buschelman ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 398*b0592601SKris Buschelman } 39914b396bbSKris Buschelman 40014b396bbSKris Buschelman ierr = PetscNew(Mat_SeqAIJ_SuperLU,&lu);CHKERRQ(ierr); 40114b396bbSKris Buschelman lu->MatView = A->ops->view; 40214b396bbSKris Buschelman lu->MatAssemblyEnd = A->ops->assemblyend; 403*b0592601SKris Buschelman lu->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 40414b396bbSKris Buschelman lu->MatDestroy = A->ops->destroy; 40514b396bbSKris Buschelman lu->CleanUpSuperLU = PETSC_FALSE; 406*b0592601SKris Buschelman 407*b0592601SKris Buschelman B->spptr = (void*)lu; 408*b0592601SKris Buschelman B->ops->view = MatView_SeqAIJ_SuperLU; 409*b0592601SKris Buschelman B->ops->assemblyend = MatAssemblyEnd_SeqAIJ_SuperLU; 410*b0592601SKris Buschelman B->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ_SuperLU; 411*b0592601SKris Buschelman B->ops->destroy = MatDestroy_SeqAIJ_SuperLU; 412*b0592601SKris Buschelman 413*b0592601SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_superlu_C", 414*b0592601SKris Buschelman "MatConvert_SeqAIJ_SuperLU",MatConvert_SeqAIJ_SuperLU);CHKERRQ(ierr); 415*b0592601SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_superlu_seqaij_C", 416*b0592601SKris Buschelman "MatConvert_SuperLU_SeqAIJ",MatConvert_SuperLU_SeqAIJ);CHKERRQ(ierr); 417*b0592601SKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)B,type);CHKERRQ(ierr); 418*b0592601SKris Buschelman *newmat = B; 419*b0592601SKris Buschelman PetscFunctionReturn(0); 420*b0592601SKris Buschelman } 421*b0592601SKris Buschelman EXTERN_C_END 422*b0592601SKris Buschelman 423*b0592601SKris Buschelman EXTERN_C_BEGIN 424*b0592601SKris Buschelman #undef __FUNCT__ 425*b0592601SKris Buschelman #define __FUNCT__ "MatCreate_SeqAIJ_SuperLU" 426*b0592601SKris Buschelman int MatCreate_SeqAIJ_SuperLU(Mat A) { 427*b0592601SKris Buschelman int ierr; 428*b0592601SKris Buschelman 429*b0592601SKris Buschelman PetscFunctionBegin; 430*b0592601SKris Buschelman ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr); 431*b0592601SKris Buschelman ierr = MatConvert_SeqAIJ_SuperLU(A,MATSUPERLU,&A);CHKERRQ(ierr); 43214b396bbSKris Buschelman PetscFunctionReturn(0); 43314b396bbSKris Buschelman } 43414b396bbSKris Buschelman EXTERN_C_END 435