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); 3014b396bbSKris Buschelman int (*MatDestroy)(Mat); 3114b396bbSKris Buschelman 3214b396bbSKris Buschelman /* Flag to clean up (non-global) SuperLU objects during Destroy */ 3314b396bbSKris Buschelman PetscTruth CleanUpSuperLU; 3414b396bbSKris Buschelman } Mat_SeqAIJ_SuperLU; 3514b396bbSKris Buschelman 3614b396bbSKris Buschelman 3714b396bbSKris Buschelman extern int MatDestroy_SeqAIJ(Mat); 3814b396bbSKris Buschelman 3914b396bbSKris Buschelman #undef __FUNCT__ 4014b396bbSKris Buschelman #define __FUNCT__ "MatDestroy_SeqAIJ_SuperLU" 4114b396bbSKris Buschelman int MatDestroy_SeqAIJ_SuperLU(Mat A) 4214b396bbSKris Buschelman { 4314b396bbSKris Buschelman Mat_SeqAIJ_SuperLU *lu = (Mat_SeqAIJ_SuperLU*)A->spptr; 4414b396bbSKris Buschelman 4514b396bbSKris Buschelman int ierr,(*destroy)(Mat); 4614b396bbSKris Buschelman 4714b396bbSKris Buschelman PetscFunctionBegin; 4814b396bbSKris Buschelman /* It looks like this is decreasing the reference count a second time during MatDestroy?! */ 4914b396bbSKris Buschelman if (--A->refct > 0)PetscFunctionReturn(0); 5014b396bbSKris Buschelman /* We have to free the global data or SuperLU crashes (sucky design)*/ 5114b396bbSKris Buschelman /* Since we don't know if more solves on other matrices may be done 5214b396bbSKris Buschelman we cannot free the yucky SuperLU global data 5314b396bbSKris Buschelman StatFree(); 5414b396bbSKris Buschelman */ 5514b396bbSKris Buschelman 5614b396bbSKris Buschelman /* Free the SuperLU datastructures */ 5714b396bbSKris Buschelman if (lu->CleanUpSuperLU) { 5814b396bbSKris Buschelman Destroy_CompCol_Permuted(&lu->AC); 5914b396bbSKris Buschelman Destroy_SuperNode_Matrix(&lu->L); 6014b396bbSKris Buschelman Destroy_CompCol_Matrix(&lu->U); 6114b396bbSKris Buschelman ierr = PetscFree(lu->B.Store);CHKERRQ(ierr); 6214b396bbSKris Buschelman ierr = PetscFree(lu->A.Store);CHKERRQ(ierr); 6314b396bbSKris Buschelman ierr = PetscFree(lu->perm_r);CHKERRQ(ierr); 6414b396bbSKris Buschelman ierr = PetscFree(lu->perm_c);CHKERRQ(ierr); 6514b396bbSKris Buschelman } 6614b396bbSKris Buschelman 6714b396bbSKris Buschelman destroy = lu->MatDestroy; 6814b396bbSKris Buschelman ierr = PetscFree(lu);CHKERRQ(ierr); 6914b396bbSKris Buschelman ierr = (*destroy)(A);CHKERRQ(ierr); 7014b396bbSKris Buschelman PetscFunctionReturn(0); 7114b396bbSKris Buschelman } 7214b396bbSKris Buschelman 7314b396bbSKris Buschelman #undef __FUNCT__ 7414b396bbSKris Buschelman #define __FUNCT__ "MatView_SeqAIJ_Spooles" 7514b396bbSKris Buschelman int MatView_SeqAIJ_SuperLU(Mat A,PetscViewer viewer) 7614b396bbSKris Buschelman { 7714b396bbSKris Buschelman int ierr; 7814b396bbSKris Buschelman PetscTruth isascii; 7914b396bbSKris Buschelman PetscViewerFormat format; 8014b396bbSKris Buschelman Mat_SeqAIJ_SuperLU *lu=(Mat_SeqAIJ_SuperLU*)(A->spptr); 8114b396bbSKris Buschelman 8214b396bbSKris Buschelman PetscFunctionBegin; 8314b396bbSKris Buschelman ierr = (*lu->MatView)(A,viewer);CHKERRQ(ierr); 8414b396bbSKris Buschelman 8514b396bbSKris Buschelman ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 8614b396bbSKris Buschelman if (isascii) { 8714b396bbSKris Buschelman ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 8814b396bbSKris Buschelman if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) { 8914b396bbSKris Buschelman ierr = MatSeqAIJFactorInfo_SuperLU(A,viewer);CHKERRQ(ierr); 9014b396bbSKris Buschelman } 9114b396bbSKris Buschelman } 9214b396bbSKris Buschelman PetscFunctionReturn(0); 9314b396bbSKris Buschelman } 9414b396bbSKris Buschelman 9514b396bbSKris Buschelman #undef __FUNCT__ 9614b396bbSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_SeqAIJ_SuperLU" 9714b396bbSKris Buschelman int MatAssemblyEnd_SeqAIJ_SuperLU(Mat A,MatAssemblyType mode) { 9814b396bbSKris Buschelman int ierr; 9914b396bbSKris Buschelman Mat_SeqAIJ_SuperLU *lu=(Mat_SeqAIJ_SuperLU*)(A->spptr); 10014b396bbSKris Buschelman 10114b396bbSKris Buschelman PetscFunctionBegin; 10214b396bbSKris Buschelman ierr = (*lu->MatAssemblyEnd)(A,mode);CHKERRQ(ierr); 10314b396bbSKris Buschelman ierr = MatUseSuperLU_SeqAIJ(A);CHKERRQ(ierr); 10414b396bbSKris Buschelman PetscFunctionReturn(0); 10514b396bbSKris Buschelman } 10614b396bbSKris Buschelman 10714b396bbSKris Buschelman #include "src/mat/impls/dense/seq/dense.h" 10814b396bbSKris Buschelman #undef __FUNCT__ 10914b396bbSKris Buschelman #define __FUNCT__ "MatCreateNull_SeqAIJ_SuperLU" 11014b396bbSKris Buschelman int MatCreateNull_SeqAIJ_SuperLU(Mat A,Mat *nullMat) 11114b396bbSKris Buschelman { 11214b396bbSKris Buschelman Mat_SeqAIJ_SuperLU *lu = (Mat_SeqAIJ_SuperLU*)A->spptr; 11314b396bbSKris Buschelman int numRows = A->m,numCols = A->n; 11414b396bbSKris Buschelman SCformat *Lstore; 11514b396bbSKris Buschelman int numNullCols,size; 11614b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 11714b396bbSKris Buschelman doublecomplex *nullVals,*workVals; 11814b396bbSKris Buschelman #else 11914b396bbSKris Buschelman PetscScalar *nullVals,*workVals; 12014b396bbSKris Buschelman #endif 12114b396bbSKris Buschelman int row,newRow,col,newCol,block,b,ierr; 12214b396bbSKris Buschelman 12314b396bbSKris Buschelman PetscFunctionBegin; 12414b396bbSKris Buschelman PetscValidHeaderSpecific(A,MAT_COOKIE); 12514b396bbSKris Buschelman PetscValidPointer(nullMat); 12614b396bbSKris Buschelman if (!A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix"); 12714b396bbSKris Buschelman numNullCols = numCols - numRows; 12814b396bbSKris Buschelman if (numNullCols < 0) SETERRQ(PETSC_ERR_ARG_WRONG,"Function only applies to underdetermined problems"); 12914b396bbSKris Buschelman /* Create the null matrix */ 13014b396bbSKris Buschelman ierr = MatCreateSeqDense(A->comm,numRows,numNullCols,PETSC_NULL,nullMat);CHKERRQ(ierr); 13114b396bbSKris Buschelman if (numNullCols == 0) { 13214b396bbSKris Buschelman ierr = MatAssemblyBegin(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13314b396bbSKris Buschelman ierr = MatAssemblyEnd(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13414b396bbSKris Buschelman PetscFunctionReturn(0); 13514b396bbSKris Buschelman } 13614b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 13714b396bbSKris Buschelman nullVals = (doublecomplex*)((Mat_SeqDense*)(*nullMat)->data)->v; 13814b396bbSKris Buschelman #else 13914b396bbSKris Buschelman nullVals = ((Mat_SeqDense*)(*nullMat)->data)->v; 14014b396bbSKris Buschelman #endif 14114b396bbSKris Buschelman /* Copy in the columns */ 14214b396bbSKris Buschelman Lstore = (SCformat*)lu->L.Store; 14314b396bbSKris Buschelman for(block = 0; block <= Lstore->nsuper; block++) { 14414b396bbSKris Buschelman newRow = Lstore->sup_to_col[block]; 14514b396bbSKris Buschelman size = Lstore->sup_to_col[block+1] - Lstore->sup_to_col[block]; 14614b396bbSKris Buschelman for(col = Lstore->rowind_colptr[newRow]; col < Lstore->rowind_colptr[newRow+1]; col++) { 14714b396bbSKris Buschelman newCol = Lstore->rowind[col]; 14814b396bbSKris Buschelman if (newCol >= numRows) { 14914b396bbSKris Buschelman for(b = 0; b < size; b++) 15014b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 15114b396bbSKris Buschelman nullVals[(newCol-numRows)*numRows+newRow+b] = ((doublecomplex*)Lstore->nzval)[Lstore->nzval_colptr[newRow+b]+col]; 15214b396bbSKris Buschelman #else 15314b396bbSKris Buschelman nullVals[(newCol-numRows)*numRows+newRow+b] = ((double*)Lstore->nzval)[Lstore->nzval_colptr[newRow+b]+col]; 15414b396bbSKris Buschelman #endif 15514b396bbSKris Buschelman } 15614b396bbSKris Buschelman } 15714b396bbSKris Buschelman } 15814b396bbSKris Buschelman /* Permute rhs to form P^T_c B */ 15914b396bbSKris Buschelman ierr = PetscMalloc(numRows*sizeof(double),&workVals);CHKERRQ(ierr); 16014b396bbSKris Buschelman for(b = 0; b < numNullCols; b++) { 16114b396bbSKris Buschelman for(row = 0; row < numRows; row++) workVals[lu->perm_c[row]] = nullVals[b*numRows+row]; 16214b396bbSKris Buschelman for(row = 0; row < numRows; row++) nullVals[b*numRows+row] = workVals[row]; 16314b396bbSKris Buschelman } 16414b396bbSKris Buschelman /* Backward solve the upper triangle A x = b */ 16514b396bbSKris Buschelman for(b = 0; b < numNullCols; b++) { 16614b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 16714b396bbSKris Buschelman sp_ztrsv("L","T","U",&lu->L,&lu->U,&nullVals[b*numRows],&ierr); 16814b396bbSKris Buschelman #else 16914b396bbSKris Buschelman sp_dtrsv("L","T","U",&lu->L,&lu->U,&nullVals[b*numRows],&ierr); 17014b396bbSKris Buschelman #endif 17114b396bbSKris Buschelman if (ierr < 0) 17214b396bbSKris Buschelman SETERRQ1(PETSC_ERR_ARG_WRONG,"The argument %d was invalid",-ierr); 17314b396bbSKris Buschelman } 17414b396bbSKris Buschelman ierr = PetscFree(workVals);CHKERRQ(ierr); 17514b396bbSKris Buschelman 17614b396bbSKris Buschelman ierr = MatAssemblyBegin(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 17714b396bbSKris Buschelman ierr = MatAssemblyEnd(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 17814b396bbSKris Buschelman PetscFunctionReturn(0); 17914b396bbSKris Buschelman } 18014b396bbSKris Buschelman 18114b396bbSKris Buschelman #undef __FUNCT__ 18214b396bbSKris Buschelman #define __FUNCT__ "MatSolve_SeqAIJ_SuperLU" 18314b396bbSKris Buschelman int MatSolve_SeqAIJ_SuperLU(Mat A,Vec b,Vec x) 18414b396bbSKris Buschelman { 18514b396bbSKris Buschelman Mat_SeqAIJ_SuperLU *lu = (Mat_SeqAIJ_SuperLU*)A->spptr; 18614b396bbSKris Buschelman PetscScalar *array; 18714b396bbSKris Buschelman int m,ierr; 18814b396bbSKris Buschelman 18914b396bbSKris Buschelman PetscFunctionBegin; 19014b396bbSKris Buschelman ierr = VecGetLocalSize(b,&m);CHKERRQ(ierr); 19114b396bbSKris Buschelman ierr = VecCopy(b,x);CHKERRQ(ierr); 19214b396bbSKris Buschelman ierr = VecGetArray(x,&array);CHKERRQ(ierr); 19314b396bbSKris Buschelman /* Create the Rhs */ 19414b396bbSKris Buschelman lu->B.Stype = SLU_DN; 19514b396bbSKris Buschelman lu->B.Mtype = SLU_GE; 19614b396bbSKris Buschelman lu->B.nrow = m; 19714b396bbSKris Buschelman lu->B.ncol = 1; 19814b396bbSKris Buschelman ((DNformat*)lu->B.Store)->lda = m; 19914b396bbSKris Buschelman ((DNformat*)lu->B.Store)->nzval = array; 20014b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 20114b396bbSKris Buschelman lu->B.Dtype = SLU_Z; 20214b396bbSKris Buschelman zgstrs("T",&lu->L,&lu->U,lu->perm_r,lu->perm_c,&lu->B,&ierr); 20314b396bbSKris Buschelman #else 20414b396bbSKris Buschelman lu->B.Dtype = SLU_D; 20514b396bbSKris Buschelman dgstrs("T",&lu->L,&lu->U,lu->perm_r,lu->perm_c,&lu->B,&ierr); 20614b396bbSKris Buschelman #endif 20714b396bbSKris Buschelman if (ierr < 0) SETERRQ1(PETSC_ERR_ARG_WRONG,"The diagonal element of row %d was invalid",-ierr); 20814b396bbSKris Buschelman ierr = VecRestoreArray(x,&array);CHKERRQ(ierr); 20914b396bbSKris Buschelman PetscFunctionReturn(0); 21014b396bbSKris Buschelman } 21114b396bbSKris Buschelman 21214b396bbSKris Buschelman static int StatInitCalled = 0; 21314b396bbSKris Buschelman 21414b396bbSKris Buschelman #undef __FUNCT__ 21514b396bbSKris Buschelman #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ_SuperLU" 21614b396bbSKris Buschelman int MatLUFactorNumeric_SeqAIJ_SuperLU(Mat A,Mat *F) 21714b396bbSKris Buschelman { 21814b396bbSKris Buschelman Mat_SeqAIJ *aa = (Mat_SeqAIJ*)(A)->data; 21914b396bbSKris Buschelman Mat_SeqAIJ_SuperLU *lu = (Mat_SeqAIJ_SuperLU*)(*F)->spptr; 22014b396bbSKris Buschelman int *etree,i,ierr; 22114b396bbSKris Buschelman PetscTruth flag; 22214b396bbSKris Buschelman 22314b396bbSKris Buschelman PetscFunctionBegin; 22414b396bbSKris Buschelman /* Create the SuperMatrix for A^T: 22514b396bbSKris Buschelman Since SuperLU only likes column-oriented matrices,we pass it the transpose, 22614b396bbSKris Buschelman and then solve A^T X = B in MatSolve(). 22714b396bbSKris Buschelman */ 22814b396bbSKris Buschelman 22914b396bbSKris Buschelman if ( lu->flg == DIFFERENT_NONZERO_PATTERN){ /* first numerical factorization */ 23014b396bbSKris Buschelman lu->A.Stype = SLU_NC; 23114b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 23214b396bbSKris Buschelman lu->A.Dtype = SLU_Z; 23314b396bbSKris Buschelman #else 23414b396bbSKris Buschelman lu->A.Dtype = SLU_D; 23514b396bbSKris Buschelman #endif 23614b396bbSKris Buschelman lu->A.Mtype = SLU_GE; 23714b396bbSKris Buschelman lu->A.nrow = A->n; 23814b396bbSKris Buschelman lu->A.ncol = A->m; 23914b396bbSKris Buschelman 24014b396bbSKris Buschelman ierr = PetscMalloc(sizeof(NCformat),&lu->store);CHKERRQ(ierr); 24114b396bbSKris Buschelman ierr = PetscMalloc(sizeof(DNformat),&lu->B.Store);CHKERRQ(ierr); 24214b396bbSKris Buschelman } 24314b396bbSKris Buschelman lu->store->nnz = aa->nz; 24414b396bbSKris Buschelman lu->store->colptr = aa->i; 24514b396bbSKris Buschelman lu->store->rowind = aa->j; 24614b396bbSKris Buschelman lu->store->nzval = aa->a; 24714b396bbSKris Buschelman lu->A.Store = lu->store; 24814b396bbSKris Buschelman 24914b396bbSKris Buschelman /* Set SuperLU options */ 25014b396bbSKris Buschelman lu->relax = sp_ienv(2); 25114b396bbSKris Buschelman lu->panel_size = sp_ienv(1); 25214b396bbSKris Buschelman /* We have to initialize global data or SuperLU crashes (sucky design) */ 25314b396bbSKris Buschelman if (!StatInitCalled) { 25414b396bbSKris Buschelman StatInit(lu->panel_size,lu->relax); 25514b396bbSKris Buschelman } 25614b396bbSKris Buschelman StatInitCalled++; 25714b396bbSKris Buschelman 25814b396bbSKris Buschelman ierr = PetscOptionsBegin(A->comm,A->prefix,"SuperLU Options","Mat");CHKERRQ(ierr); 25914b396bbSKris Buschelman /* use SuperLU mat ordeing */ 26014b396bbSKris 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); 26114b396bbSKris Buschelman if (flag) { 26214b396bbSKris Buschelman get_perm_c(lu->ispec, &lu->A, lu->perm_c); 26314b396bbSKris Buschelman lu->SuperluMatOdering = PETSC_TRUE; 26414b396bbSKris Buschelman } 26514b396bbSKris Buschelman PetscOptionsEnd(); 26614b396bbSKris Buschelman 26714b396bbSKris Buschelman /* Create the elimination tree */ 26814b396bbSKris Buschelman ierr = PetscMalloc(A->n*sizeof(int),&etree);CHKERRQ(ierr); 26914b396bbSKris Buschelman sp_preorder("N",&lu->A,lu->perm_c,etree,&lu->AC); 27014b396bbSKris Buschelman /* Factor the matrix */ 27114b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 27214b396bbSKris 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); 27314b396bbSKris Buschelman #else 27414b396bbSKris 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); 27514b396bbSKris Buschelman #endif 27614b396bbSKris Buschelman if (ierr < 0) { 27714b396bbSKris Buschelman SETERRQ1(PETSC_ERR_ARG_WRONG,"The diagonal element of row %d was invalid",-ierr); 27814b396bbSKris Buschelman } else if (ierr > 0) { 27914b396bbSKris Buschelman if (ierr <= A->m) { 28014b396bbSKris Buschelman SETERRQ1(PETSC_ERR_ARG_WRONG,"The diagonal element %d of U is exactly zero",ierr); 28114b396bbSKris Buschelman } else { 28214b396bbSKris Buschelman SETERRQ1(PETSC_ERR_ARG_WRONG,"Memory allocation failure after %d bytes were allocated",ierr-A->m); 28314b396bbSKris Buschelman } 28414b396bbSKris Buschelman } 28514b396bbSKris Buschelman 28614b396bbSKris Buschelman /* Cleanup */ 28714b396bbSKris Buschelman ierr = PetscFree(etree);CHKERRQ(ierr); 28814b396bbSKris Buschelman 28914b396bbSKris Buschelman lu->flg = SAME_NONZERO_PATTERN; 29014b396bbSKris Buschelman PetscFunctionReturn(0); 29114b396bbSKris Buschelman } 29214b396bbSKris Buschelman 29314b396bbSKris Buschelman /* 29414b396bbSKris Buschelman Note the r permutation is ignored 29514b396bbSKris Buschelman */ 29614b396bbSKris Buschelman #undef __FUNCT__ 29714b396bbSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJ_SuperLU" 29814b396bbSKris Buschelman int MatLUFactorSymbolic_SeqAIJ_SuperLU(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) 29914b396bbSKris Buschelman { 30014b396bbSKris Buschelman Mat B; 30114b396bbSKris Buschelman Mat_SeqAIJ_SuperLU *lu; 30214b396bbSKris Buschelman int ierr,*ca; 30314b396bbSKris Buschelman 30414b396bbSKris Buschelman PetscFunctionBegin; 30514b396bbSKris Buschelman 30614b396bbSKris Buschelman ierr = MatCreateSeqAIJ(A->comm,A->m,A->n,0,PETSC_NULL,F);CHKERRQ(ierr); 30714b396bbSKris Buschelman B = *F; 30814b396bbSKris Buschelman B->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJ_SuperLU; 30914b396bbSKris Buschelman B->ops->solve = MatSolve_SeqAIJ_SuperLU; 31014b396bbSKris Buschelman B->ops->destroy = MatDestroy_SeqAIJ_SuperLU; 31114b396bbSKris Buschelman B->factor = FACTOR_LU; 31214b396bbSKris Buschelman (*F)->assembled = PETSC_TRUE; /* required by -sles_view */ 31314b396bbSKris Buschelman 31414b396bbSKris Buschelman ierr = PetscNew(Mat_SeqAIJ_SuperLU,&lu);CHKERRQ(ierr); 31514b396bbSKris Buschelman B->spptr = (void*)lu; 31614b396bbSKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)B,"MatCreateNull","MatCreateNull_SeqAIJ_SuperLU", 31714b396bbSKris Buschelman (void(*)(void))MatCreateNull_SeqAIJ_SuperLU);CHKERRQ(ierr); 31814b396bbSKris Buschelman 31914b396bbSKris Buschelman /* Allocate the work arrays required by SuperLU (notice sizes are for the transpose) */ 32014b396bbSKris Buschelman ierr = PetscMalloc(A->n*sizeof(int),&lu->perm_r);CHKERRQ(ierr); 32114b396bbSKris Buschelman ierr = PetscMalloc(A->m*sizeof(int),&lu->perm_c);CHKERRQ(ierr); 32214b396bbSKris Buschelman 32314b396bbSKris Buschelman /* use PETSc mat ordering */ 32414b396bbSKris Buschelman ierr = ISGetIndices(c,&ca);CHKERRQ(ierr); 32514b396bbSKris Buschelman ierr = PetscMemcpy(lu->perm_c,ca,A->m*sizeof(int));CHKERRQ(ierr); 32614b396bbSKris Buschelman ierr = ISRestoreIndices(c,&ca);CHKERRQ(ierr); 32714b396bbSKris Buschelman lu->SuperluMatOdering = PETSC_FALSE; 32814b396bbSKris Buschelman 32914b396bbSKris Buschelman lu->pivot_threshold = info->dtcol; 33014b396bbSKris Buschelman PetscLogObjectMemory(B,(A->m+A->n)*sizeof(int)+sizeof(Mat_SeqAIJ_SuperLU)); 33114b396bbSKris Buschelman 33214b396bbSKris Buschelman lu->flg = DIFFERENT_NONZERO_PATTERN; 333*e740cb95SKris Buschelman lu->CleanUpSuperLU = PETSC_TRUE; 33414b396bbSKris Buschelman PetscFunctionReturn(0); 33514b396bbSKris Buschelman } 33614b396bbSKris Buschelman 33714b396bbSKris Buschelman #undef __FUNCT__ 33814b396bbSKris Buschelman #define __FUNCT__ "MatUseSuperLU_SeqAIJ" 33914b396bbSKris Buschelman int MatUseSuperLU_SeqAIJ(Mat A) 34014b396bbSKris Buschelman { 34114b396bbSKris Buschelman PetscTruth flg; 34214b396bbSKris Buschelman int ierr; 34314b396bbSKris Buschelman 34414b396bbSKris Buschelman PetscFunctionBegin; 34514b396bbSKris Buschelman PetscValidHeaderSpecific(A,MAT_COOKIE); 34614b396bbSKris Buschelman ierr = PetscTypeCompare((PetscObject)A,MATSUPERLU,&flg); 34714b396bbSKris Buschelman if (!flg) PetscFunctionReturn(0); 34814b396bbSKris Buschelman 34914b396bbSKris Buschelman A->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ_SuperLU; 35014b396bbSKris Buschelman 35114b396bbSKris Buschelman PetscFunctionReturn(0); 35214b396bbSKris Buschelman } 35314b396bbSKris Buschelman 35414b396bbSKris Buschelman /* used by -sles_view */ 35514b396bbSKris Buschelman #undef __FUNCT__ 35614b396bbSKris Buschelman #define __FUNCT__ "MatSeqAIJFactorInfo_SuperLU" 35714b396bbSKris Buschelman int MatSeqAIJFactorInfo_SuperLU(Mat A,PetscViewer viewer) 35814b396bbSKris Buschelman { 35914b396bbSKris Buschelman Mat_SeqAIJ_SuperLU *lu= (Mat_SeqAIJ_SuperLU*)A->spptr; 36014b396bbSKris Buschelman int ierr; 36114b396bbSKris Buschelman PetscFunctionBegin; 36214b396bbSKris Buschelman /* check if matrix is SuperLU type */ 36314b396bbSKris Buschelman if (A->ops->solve != MatSolve_SeqAIJ_SuperLU) PetscFunctionReturn(0); 36414b396bbSKris Buschelman 36514b396bbSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");CHKERRQ(ierr); 36614b396bbSKris Buschelman if(lu->SuperluMatOdering) ierr = PetscViewerASCIIPrintf(viewer," SuperLU mat ordering: %d\n",lu->ispec);CHKERRQ(ierr); 36714b396bbSKris Buschelman 36814b396bbSKris Buschelman PetscFunctionReturn(0); 36914b396bbSKris Buschelman } 37014b396bbSKris Buschelman 37114b396bbSKris Buschelman EXTERN_C_BEGIN 37214b396bbSKris Buschelman #undef __FUNCT__ 37314b396bbSKris Buschelman #define __FUNCT__ "MatCreate_SeqAIJ_SuperLU" 37414b396bbSKris Buschelman int MatCreate_SeqAIJ_SuperLU(Mat A) { 37514b396bbSKris Buschelman int ierr; 37614b396bbSKris Buschelman Mat_SeqAIJ_SuperLU *lu; 37714b396bbSKris Buschelman 37814b396bbSKris Buschelman PetscFunctionBegin; 37914b396bbSKris Buschelman ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr); 38014b396bbSKris Buschelman ierr = MatUseSuperLU_SeqAIJ(A);CHKERRQ(ierr); 38114b396bbSKris Buschelman 38214b396bbSKris Buschelman ierr = PetscNew(Mat_SeqAIJ_SuperLU,&lu);CHKERRQ(ierr); 38314b396bbSKris Buschelman lu->MatView = A->ops->view; 38414b396bbSKris Buschelman lu->MatAssemblyEnd = A->ops->assemblyend; 38514b396bbSKris Buschelman lu->MatDestroy = A->ops->destroy; 38614b396bbSKris Buschelman lu->CleanUpSuperLU = PETSC_FALSE; 38714b396bbSKris Buschelman A->spptr = (void*)lu; 38814b396bbSKris Buschelman A->ops->view = MatView_SeqAIJ_SuperLU; 38914b396bbSKris Buschelman A->ops->assemblyend = MatAssemblyEnd_SeqAIJ_SuperLU; 39014b396bbSKris Buschelman A->ops->destroy = MatDestroy_SeqAIJ_SuperLU; 39114b396bbSKris Buschelman PetscFunctionReturn(0); 39214b396bbSKris Buschelman } 39314b396bbSKris Buschelman EXTERN_C_END 394