1be1d678aSKris Buschelman #define PETSCMAT_DLL 214b396bbSKris Buschelman 314b396bbSKris Buschelman /* 475af56d4SHong Zhang Provides an interface to the SuperLU 3.0 sparse solver 514b396bbSKris Buschelman */ 614b396bbSKris Buschelman 714b396bbSKris Buschelman #include "src/mat/impls/aij/seq/aij.h" 814b396bbSKris Buschelman 914b396bbSKris Buschelman EXTERN_C_BEGIN 1014b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 1166f57492SHong Zhang #include "slu_zdefs.h" 1214b396bbSKris Buschelman #else 1366f57492SHong Zhang #include "slu_ddefs.h" 1414b396bbSKris Buschelman #endif 1566f57492SHong Zhang #include "slu_util.h" 1614b396bbSKris Buschelman EXTERN_C_END 1714b396bbSKris Buschelman 1814b396bbSKris Buschelman typedef struct { 1975af56d4SHong Zhang SuperMatrix A,L,U,B,X; 2075af56d4SHong Zhang superlu_options_t options; 21da7a1d00SHong Zhang PetscInt *perm_c; /* column permutation vector */ 22da7a1d00SHong Zhang PetscInt *perm_r; /* row permutations from partial pivoting */ 23da7a1d00SHong Zhang PetscInt *etree; 24da7a1d00SHong Zhang PetscReal *R, *C; 2575af56d4SHong Zhang char equed[1]; 26da7a1d00SHong Zhang PetscInt lwork; 2775af56d4SHong Zhang void *work; 28da7a1d00SHong Zhang PetscReal rpg, rcond; 2975af56d4SHong Zhang mem_usage_t mem_usage; 3075af56d4SHong Zhang MatStructure flg; 3114b396bbSKris Buschelman 3214b396bbSKris Buschelman /* A few function pointers for inheritance */ 336849ba73SBarry Smith PetscErrorCode (*MatDuplicate)(Mat,MatDuplicateOption,Mat*); 346849ba73SBarry Smith PetscErrorCode (*MatView)(Mat,PetscViewer); 356849ba73SBarry Smith PetscErrorCode (*MatAssemblyEnd)(Mat,MatAssemblyType); 366849ba73SBarry Smith PetscErrorCode (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*); 376849ba73SBarry Smith PetscErrorCode (*MatDestroy)(Mat); 3814b396bbSKris Buschelman 3914b396bbSKris Buschelman /* Flag to clean up (non-global) SuperLU objects during Destroy */ 4014b396bbSKris Buschelman PetscTruth CleanUpSuperLU; 41f0c56d0fSKris Buschelman } Mat_SuperLU; 4214b396bbSKris Buschelman 4314b396bbSKris Buschelman 44dfbe8321SBarry Smith EXTERN PetscErrorCode MatFactorInfo_SuperLU(Mat,PetscViewer); 45dfbe8321SBarry Smith EXTERN PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat,IS,IS,MatFactorInfo*,Mat*); 46b0592601SKris Buschelman 47b0592601SKris Buschelman EXTERN_C_BEGIN 4875179d2cSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SuperLU_SeqAIJ(Mat,MatType,MatReuse,Mat*); 4975179d2cSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SuperLU(Mat,MatType,MatReuse,Mat*); 50b0592601SKris Buschelman EXTERN_C_END 5114b396bbSKris Buschelman 5214b396bbSKris Buschelman #undef __FUNCT__ 53f0c56d0fSKris Buschelman #define __FUNCT__ "MatDestroy_SuperLU" 54dfbe8321SBarry Smith PetscErrorCode MatDestroy_SuperLU(Mat A) 5514b396bbSKris Buschelman { 56dfbe8321SBarry Smith PetscErrorCode ierr; 57f0c56d0fSKris Buschelman Mat_SuperLU *lu=(Mat_SuperLU*)A->spptr; 5814b396bbSKris Buschelman 5914b396bbSKris Buschelman PetscFunctionBegin; 6075af56d4SHong Zhang if (lu->CleanUpSuperLU) { /* Free the SuperLU datastructures */ 6175af56d4SHong Zhang Destroy_SuperMatrix_Store(&lu->A); 6275af56d4SHong Zhang Destroy_SuperMatrix_Store(&lu->B); 6375af56d4SHong Zhang Destroy_SuperMatrix_Store(&lu->X); 649ce50919SHong Zhang 659ce50919SHong Zhang ierr = PetscFree(lu->etree);CHKERRQ(ierr); 669ce50919SHong Zhang ierr = PetscFree(lu->perm_r);CHKERRQ(ierr); 679ce50919SHong Zhang ierr = PetscFree(lu->perm_c);CHKERRQ(ierr); 689ce50919SHong Zhang ierr = PetscFree(lu->R);CHKERRQ(ierr); 699ce50919SHong Zhang ierr = PetscFree(lu->C);CHKERRQ(ierr); 7075af56d4SHong Zhang if ( lu->lwork >= 0 ) { 71fb3e25aaSKris Buschelman Destroy_SuperNode_Matrix(&lu->L); 72fb3e25aaSKris Buschelman Destroy_CompCol_Matrix(&lu->U); 7375af56d4SHong Zhang } 74fb3e25aaSKris Buschelman } 75ceb03754SKris Buschelman ierr = MatConvert_SuperLU_SeqAIJ(A,MATSEQAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr); 76fb3e25aaSKris Buschelman ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 7714b396bbSKris Buschelman PetscFunctionReturn(0); 7814b396bbSKris Buschelman } 7914b396bbSKris Buschelman 8014b396bbSKris Buschelman #undef __FUNCT__ 81f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_SuperLU" 82dfbe8321SBarry Smith PetscErrorCode MatView_SuperLU(Mat A,PetscViewer viewer) 8314b396bbSKris Buschelman { 84dfbe8321SBarry Smith PetscErrorCode ierr; 8532077d6dSBarry Smith PetscTruth iascii; 8614b396bbSKris Buschelman PetscViewerFormat format; 87f0c56d0fSKris Buschelman Mat_SuperLU *lu=(Mat_SuperLU*)(A->spptr); 8814b396bbSKris Buschelman 8914b396bbSKris Buschelman PetscFunctionBegin; 9014b396bbSKris Buschelman ierr = (*lu->MatView)(A,viewer);CHKERRQ(ierr); 9114b396bbSKris Buschelman 9232077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 9332077d6dSBarry Smith if (iascii) { 9414b396bbSKris Buschelman ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 95*2f59403fSHong Zhang if (format == PETSC_VIEWER_ASCII_INFO) { 96f0c56d0fSKris Buschelman ierr = MatFactorInfo_SuperLU(A,viewer);CHKERRQ(ierr); 9714b396bbSKris Buschelman } 9814b396bbSKris Buschelman } 9914b396bbSKris Buschelman PetscFunctionReturn(0); 10014b396bbSKris Buschelman } 10114b396bbSKris Buschelman 10214b396bbSKris Buschelman #undef __FUNCT__ 103f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_SuperLU" 104dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_SuperLU(Mat A,MatAssemblyType mode) { 105dfbe8321SBarry Smith PetscErrorCode ierr; 106f0c56d0fSKris Buschelman Mat_SuperLU *lu=(Mat_SuperLU*)(A->spptr); 10714b396bbSKris Buschelman 10814b396bbSKris Buschelman PetscFunctionBegin; 10914b396bbSKris Buschelman ierr = (*lu->MatAssemblyEnd)(A,mode);CHKERRQ(ierr); 110b0592601SKris Buschelman lu->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 111f0c56d0fSKris Buschelman A->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU; 11214b396bbSKris Buschelman PetscFunctionReturn(0); 11314b396bbSKris Buschelman } 11414b396bbSKris Buschelman 11575af56d4SHong Zhang /* This function was written for SuperLU 2.0 by Matthew Knepley. Not tested for SuperLU 3.0! */ 11675af56d4SHong Zhang #ifdef SuperLU2 11714b396bbSKris Buschelman #include "src/mat/impls/dense/seq/dense.h" 11814b396bbSKris Buschelman #undef __FUNCT__ 119f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreateNull_SuperLU" 120dfbe8321SBarry Smith PetscErrorCode MatCreateNull_SuperLU(Mat A,Mat *nullMat) 12114b396bbSKris Buschelman { 122f0c56d0fSKris Buschelman Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 1232a4c71feSBarry Smith PetscInt numRows = A->rmap.n,numCols = A->cmap.n; 12414b396bbSKris Buschelman SCformat *Lstore; 125da7a1d00SHong Zhang PetscInt numNullCols,size; 12675af56d4SHong Zhang SuperLUStat_t stat; 12714b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 12814b396bbSKris Buschelman doublecomplex *nullVals,*workVals; 12914b396bbSKris Buschelman #else 13014b396bbSKris Buschelman PetscScalar *nullVals,*workVals; 13114b396bbSKris Buschelman #endif 132da7a1d00SHong Zhang PetscInt row,newRow,col,newCol,block,b; 1336849ba73SBarry Smith PetscErrorCode ierr; 13414b396bbSKris Buschelman 13514b396bbSKris Buschelman PetscFunctionBegin; 13614b396bbSKris Buschelman if (!A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix"); 13714b396bbSKris Buschelman numNullCols = numCols - numRows; 13814b396bbSKris Buschelman if (numNullCols < 0) SETERRQ(PETSC_ERR_ARG_WRONG,"Function only applies to underdetermined problems"); 139be5d1d56SKris Buschelman /* Create the null matrix using MATSEQDENSE explicitly */ 140f69a0ea3SMatthew Knepley ierr = MatCreate(A->comm,nullMat);CHKERRQ(ierr); 141f69a0ea3SMatthew Knepley ierr = MatSetSizes(*nullMat,numRows,numNullCols,numRows,numNullCols);CHKERRQ(ierr); 142878740d9SKris Buschelman ierr = MatSetType(*nullMat,MATSEQDENSE);CHKERRQ(ierr); 143878740d9SKris Buschelman ierr = MatSeqDenseSetPreallocation(*nullMat,PETSC_NULL);CHKERRQ(ierr); 144958c9bccSBarry Smith if (!numNullCols) { 14514b396bbSKris Buschelman ierr = MatAssemblyBegin(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14614b396bbSKris Buschelman ierr = MatAssemblyEnd(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14714b396bbSKris Buschelman PetscFunctionReturn(0); 14814b396bbSKris Buschelman } 14914b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 15014b396bbSKris Buschelman nullVals = (doublecomplex*)((Mat_SeqDense*)(*nullMat)->data)->v; 15114b396bbSKris Buschelman #else 15214b396bbSKris Buschelman nullVals = ((Mat_SeqDense*)(*nullMat)->data)->v; 15314b396bbSKris Buschelman #endif 15414b396bbSKris Buschelman /* Copy in the columns */ 15514b396bbSKris Buschelman Lstore = (SCformat*)lu->L.Store; 15614b396bbSKris Buschelman for(block = 0; block <= Lstore->nsuper; block++) { 15714b396bbSKris Buschelman newRow = Lstore->sup_to_col[block]; 15814b396bbSKris Buschelman size = Lstore->sup_to_col[block+1] - Lstore->sup_to_col[block]; 15914b396bbSKris Buschelman for(col = Lstore->rowind_colptr[newRow]; col < Lstore->rowind_colptr[newRow+1]; col++) { 16014b396bbSKris Buschelman newCol = Lstore->rowind[col]; 16114b396bbSKris Buschelman if (newCol >= numRows) { 16214b396bbSKris Buschelman for(b = 0; b < size; b++) 16314b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 16414b396bbSKris Buschelman nullVals[(newCol-numRows)*numRows+newRow+b] = ((doublecomplex*)Lstore->nzval)[Lstore->nzval_colptr[newRow+b]+col]; 16514b396bbSKris Buschelman #else 16614b396bbSKris Buschelman nullVals[(newCol-numRows)*numRows+newRow+b] = ((double*)Lstore->nzval)[Lstore->nzval_colptr[newRow+b]+col]; 16714b396bbSKris Buschelman #endif 16814b396bbSKris Buschelman } 16914b396bbSKris Buschelman } 17014b396bbSKris Buschelman } 17114b396bbSKris Buschelman /* Permute rhs to form P^T_c B */ 172da7a1d00SHong Zhang ierr = PetscMalloc(numRows*sizeof(PetscReal),&workVals);CHKERRQ(ierr); 17314b396bbSKris Buschelman for(b = 0; b < numNullCols; b++) { 17414b396bbSKris Buschelman for(row = 0; row < numRows; row++) workVals[lu->perm_c[row]] = nullVals[b*numRows+row]; 17514b396bbSKris Buschelman for(row = 0; row < numRows; row++) nullVals[b*numRows+row] = workVals[row]; 17614b396bbSKris Buschelman } 17714b396bbSKris Buschelman /* Backward solve the upper triangle A x = b */ 17814b396bbSKris Buschelman for(b = 0; b < numNullCols; b++) { 17914b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 18075af56d4SHong Zhang sp_ztrsv("L","T","U",&lu->L,&lu->U,&nullVals[b*numRows],&stat,&ierr); 18114b396bbSKris Buschelman #else 18275af56d4SHong Zhang sp_dtrsv("L","T","U",&lu->L,&lu->U,&nullVals[b*numRows],&stat,&ierr); 18314b396bbSKris Buschelman #endif 18414b396bbSKris Buschelman if (ierr < 0) 18577431f27SBarry Smith SETERRQ1(PETSC_ERR_ARG_WRONG,"The argument %D was invalid",-ierr); 18614b396bbSKris Buschelman } 18714b396bbSKris Buschelman ierr = PetscFree(workVals);CHKERRQ(ierr); 18814b396bbSKris Buschelman 18914b396bbSKris Buschelman ierr = MatAssemblyBegin(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 19014b396bbSKris Buschelman ierr = MatAssemblyEnd(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 19114b396bbSKris Buschelman PetscFunctionReturn(0); 19214b396bbSKris Buschelman } 19375af56d4SHong Zhang #endif 19414b396bbSKris Buschelman 19514b396bbSKris Buschelman #undef __FUNCT__ 196c29e884eSHong Zhang #define __FUNCT__ "MatSolve_SuperLU_Private" 197c29e884eSHong Zhang PetscErrorCode MatSolve_SuperLU_Private(Mat A,Vec b,Vec x) 19814b396bbSKris Buschelman { 199f0c56d0fSKris Buschelman Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 20075af56d4SHong Zhang PetscScalar *barray,*xarray; 201dfbe8321SBarry Smith PetscErrorCode ierr; 202da7a1d00SHong Zhang PetscInt info,i; 20375af56d4SHong Zhang SuperLUStat_t stat; 204da7a1d00SHong Zhang PetscReal ferr,berr; 20514b396bbSKris Buschelman 20614b396bbSKris Buschelman PetscFunctionBegin; 207937a6b0eSHong Zhang if ( lu->lwork == -1 ) { 208937a6b0eSHong Zhang PetscFunctionReturn(0); 209937a6b0eSHong Zhang } 21075af56d4SHong Zhang lu->B.ncol = 1; /* Set the number of right-hand side */ 21175af56d4SHong Zhang ierr = VecGetArray(b,&barray);CHKERRQ(ierr); 21275af56d4SHong Zhang ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 2135fe6150dSHong Zhang 2145fe6150dSHong Zhang #if defined(PETSC_USE_COMPLEX) 2155fe6150dSHong Zhang ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray; 2165fe6150dSHong Zhang ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray; 2175fe6150dSHong Zhang #else 2185fe6150dSHong Zhang ((DNformat*)lu->B.Store)->nzval = barray; 21975af56d4SHong Zhang ((DNformat*)lu->X.Store)->nzval = xarray; 2205fe6150dSHong Zhang #endif 22175af56d4SHong Zhang 22275af56d4SHong Zhang /* Initialize the statistics variables. */ 22375af56d4SHong Zhang StatInit(&stat); 22475af56d4SHong Zhang 22575af56d4SHong Zhang lu->options.Fact = FACTORED; /* Indicate the factored form of A is supplied. */ 2268914a3f7SHong Zhang #if defined(PETSC_USE_COMPLEX) 2278914a3f7SHong Zhang zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 2288914a3f7SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 2298914a3f7SHong Zhang &lu->mem_usage, &stat, &info); 2308914a3f7SHong Zhang #else 23175af56d4SHong Zhang dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 23275af56d4SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 23375af56d4SHong Zhang &lu->mem_usage, &stat, &info); 2348914a3f7SHong Zhang #endif 23575af56d4SHong Zhang ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr); 23675af56d4SHong Zhang ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 23775af56d4SHong Zhang 238958c9bccSBarry Smith if ( !info || info == lu->A.ncol+1 ) { 23975af56d4SHong Zhang if ( lu->options.IterRefine ) { 2408914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n"); 2418914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR"); 24275af56d4SHong Zhang for (i = 0; i < 1; ++i) 2438914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," %8d%8d%16e%16e\n", i+1, stat.RefineSteps, ferr, berr); 24475af56d4SHong Zhang } 2458914a3f7SHong Zhang } else if ( info > 0 ){ 2468914a3f7SHong Zhang if ( lu->lwork == -1 ) { 24777431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %D bytes\n", info - lu->A.ncol); 2488914a3f7SHong Zhang } else { 24977431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," Warning: gssvx() returns info %D\n",info); 2508914a3f7SHong Zhang } 2518914a3f7SHong Zhang } else if (info < 0){ 25277431f27SBarry Smith SETERRQ2(PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", info,-info); 25314b396bbSKris Buschelman } 25414b396bbSKris Buschelman 2558914a3f7SHong Zhang if ( lu->options.PrintStat ) { 2568914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n"); 2578914a3f7SHong Zhang StatPrint(&stat); 2588914a3f7SHong Zhang } 25975af56d4SHong Zhang StatFree(&stat); 26075af56d4SHong Zhang PetscFunctionReturn(0); 26175af56d4SHong Zhang } 26214b396bbSKris Buschelman 26314b396bbSKris Buschelman #undef __FUNCT__ 264c29e884eSHong Zhang #define __FUNCT__ "MatSolve_SuperLU" 265c29e884eSHong Zhang PetscErrorCode MatSolve_SuperLU(Mat A,Vec b,Vec x) 266c29e884eSHong Zhang { 267c29e884eSHong Zhang Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 268c29e884eSHong Zhang PetscErrorCode ierr; 269c29e884eSHong Zhang 270c29e884eSHong Zhang PetscFunctionBegin; 271c29e884eSHong Zhang lu->options.Trans = TRANS; 272c29e884eSHong Zhang ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr); 273c29e884eSHong Zhang PetscFunctionReturn(0); 274c29e884eSHong Zhang } 275c29e884eSHong Zhang 276c29e884eSHong Zhang #undef __FUNCT__ 277c7c1fe80SHong Zhang #define __FUNCT__ "MatSolveTranspose_SuperLU" 278c7c1fe80SHong Zhang PetscErrorCode MatSolveTranspose_SuperLU(Mat A,Vec b,Vec x) 279c7c1fe80SHong Zhang { 280c7c1fe80SHong Zhang Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 281c7c1fe80SHong Zhang PetscErrorCode ierr; 282c7c1fe80SHong Zhang 283c7c1fe80SHong Zhang PetscFunctionBegin; 284c7c1fe80SHong Zhang lu->options.Trans = NOTRANS; 285c29e884eSHong Zhang ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr); 286c7c1fe80SHong Zhang PetscFunctionReturn(0); 287c7c1fe80SHong Zhang } 288c7c1fe80SHong Zhang 289c7c1fe80SHong Zhang #undef __FUNCT__ 290f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorNumeric_SuperLU" 291af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SuperLU(Mat A,MatFactorInfo *info,Mat *F) 29214b396bbSKris Buschelman { 29314b396bbSKris Buschelman Mat_SeqAIJ *aa = (Mat_SeqAIJ*)(A)->data; 294f0c56d0fSKris Buschelman Mat_SuperLU *lu = (Mat_SuperLU*)(*F)->spptr; 295dfbe8321SBarry Smith PetscErrorCode ierr; 296da7a1d00SHong Zhang PetscInt sinfo; 29775af56d4SHong Zhang SuperLUStat_t stat; 298da7a1d00SHong Zhang PetscReal ferr, berr; 2998914a3f7SHong Zhang NCformat *Ustore; 3008914a3f7SHong Zhang SCformat *Lstore; 30114b396bbSKris Buschelman 30214b396bbSKris Buschelman PetscFunctionBegin; 3039ce50919SHong Zhang if (lu->flg == SAME_NONZERO_PATTERN){ /* successing numerical factorization */ 3049ce50919SHong Zhang lu->options.Fact = SamePattern; 30575af56d4SHong Zhang /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */ 30675af56d4SHong Zhang Destroy_SuperMatrix_Store(&lu->A); 3079ce50919SHong Zhang if ( lu->lwork >= 0 ) { 30875af56d4SHong Zhang Destroy_SuperNode_Matrix(&lu->L); 30975af56d4SHong Zhang Destroy_CompCol_Matrix(&lu->U); 31075af56d4SHong Zhang lu->options.Fact = SamePattern; 31114b396bbSKris Buschelman } 31275af56d4SHong Zhang } 31314b396bbSKris Buschelman 31475af56d4SHong Zhang /* Create the SuperMatrix for lu->A=A^T: 31575af56d4SHong Zhang Since SuperLU likes column-oriented matrices,we pass it the transpose, 31675af56d4SHong Zhang and then solve A^T X = B in MatSolve(). */ 31775af56d4SHong Zhang #if defined(PETSC_USE_COMPLEX) 3182a4c71feSBarry Smith zCreate_CompCol_Matrix(&lu->A,A->cmap.n,A->rmap.n,aa->nz,(doublecomplex*)aa->a,aa->j,aa->i, 31975af56d4SHong Zhang SLU_NC,SLU_Z,SLU_GE); 32075af56d4SHong Zhang #else 3212a4c71feSBarry Smith dCreate_CompCol_Matrix(&lu->A,A->cmap.n,A->rmap.n,aa->nz,aa->a,aa->j,aa->i, 32275af56d4SHong Zhang SLU_NC,SLU_D,SLU_GE); 32375af56d4SHong Zhang #endif 32475af56d4SHong Zhang 32575af56d4SHong Zhang /* Initialize the statistics variables. */ 32675af56d4SHong Zhang StatInit(&stat); 32775af56d4SHong Zhang 3289ce50919SHong Zhang /* Numerical factorization */ 32975af56d4SHong Zhang lu->B.ncol = 0; /* Indicate not to solve the system */ 3308914a3f7SHong Zhang #if defined(PETSC_USE_COMPLEX) 3318914a3f7SHong Zhang zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 3328914a3f7SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 333da7a1d00SHong Zhang &lu->mem_usage, &stat, &sinfo); 3348914a3f7SHong Zhang #else 33575af56d4SHong Zhang dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 33675af56d4SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 337da7a1d00SHong Zhang &lu->mem_usage, &stat, &sinfo); 3388914a3f7SHong Zhang #endif 339da7a1d00SHong Zhang if ( !sinfo || sinfo == lu->A.ncol+1 ) { 3408914a3f7SHong Zhang if ( lu->options.PivotGrowth ) 3418914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," Recip. pivot growth = %e\n", lu->rpg); 34275af56d4SHong Zhang if ( lu->options.ConditionNumber ) 3438914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," Recip. condition number = %e\n", lu->rcond); 344da7a1d00SHong Zhang } else if ( sinfo > 0 ){ 3458914a3f7SHong Zhang if ( lu->lwork == -1 ) { 346da7a1d00SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %D bytes\n", sinfo - lu->A.ncol); 3478914a3f7SHong Zhang } else { 348da7a1d00SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," Warning: gssvx() returns info %D\n",sinfo); 3498914a3f7SHong Zhang } 350da7a1d00SHong Zhang } else { /* sinfo < 0 */ 351da7a1d00SHong Zhang SETERRQ2(PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", sinfo,-sinfo); 35275af56d4SHong Zhang } 35375af56d4SHong Zhang 3548914a3f7SHong Zhang if ( lu->options.PrintStat ) { 3558914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"MatLUFactorNumeric_SuperLU():\n"); 3568914a3f7SHong Zhang StatPrint(&stat); 3578914a3f7SHong Zhang Lstore = (SCformat *) lu->L.Store; 3588914a3f7SHong Zhang Ustore = (NCformat *) lu->U.Store; 35977431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor L = %D\n", Lstore->nnz); 36077431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor U = %D\n", Ustore->nnz); 36177431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in L+U = %D\n", Lstore->nnz + Ustore->nnz - lu->A.ncol); 36277431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," L\\U MB %.3f\ttotal MB needed %.3f\texpansions %D\n", 3638914a3f7SHong Zhang lu->mem_usage.for_lu/1e6, lu->mem_usage.total_needed/1e6, 3648914a3f7SHong Zhang lu->mem_usage.expansions); 3658914a3f7SHong Zhang } 36675af56d4SHong Zhang StatFree(&stat); 36775af56d4SHong Zhang 36875af56d4SHong Zhang lu->flg = SAME_NONZERO_PATTERN; 36975af56d4SHong Zhang PetscFunctionReturn(0); 37014b396bbSKris Buschelman } 37114b396bbSKris Buschelman 37214b396bbSKris Buschelman /* 37314b396bbSKris Buschelman Note the r permutation is ignored 37414b396bbSKris Buschelman */ 37514b396bbSKris Buschelman #undef __FUNCT__ 376f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_SuperLU" 377dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) 37814b396bbSKris Buschelman { 37914b396bbSKris Buschelman Mat B; 380f0c56d0fSKris Buschelman Mat_SuperLU *lu; 3816849ba73SBarry Smith PetscErrorCode ierr; 3822a4c71feSBarry Smith PetscInt m=A->rmap.n,n=A->cmap.n,indx; 3839ce50919SHong Zhang PetscTruth flg; 3845ca28756SHong Zhang const char *colperm[]={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */ 3855ca28756SHong Zhang const char *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"}; 3865ca28756SHong Zhang const char *rowperm[]={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */ 38714b396bbSKris Buschelman 38814b396bbSKris Buschelman PetscFunctionBegin; 389f69a0ea3SMatthew Knepley ierr = MatCreate(A->comm,&B);CHKERRQ(ierr); 3902a4c71feSBarry Smith ierr = MatSetSizes(B,A->rmap.n,A->cmap.n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 391be5d1d56SKris Buschelman ierr = MatSetType(B,A->type_name);CHKERRQ(ierr); 392899d7b4fSKris Buschelman ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 393f0c56d0fSKris Buschelman 394f0c56d0fSKris Buschelman B->ops->lufactornumeric = MatLUFactorNumeric_SuperLU; 395f0c56d0fSKris Buschelman B->ops->solve = MatSolve_SuperLU; 396c7c1fe80SHong Zhang B->ops->solvetranspose = MatSolveTranspose_SuperLU; 39714b396bbSKris Buschelman B->factor = FACTOR_LU; 39894b7f48cSBarry Smith B->assembled = PETSC_TRUE; /* required by -ksp_view */ 39914b396bbSKris Buschelman 400f0c56d0fSKris Buschelman lu = (Mat_SuperLU*)(B->spptr); 4019ce50919SHong Zhang 4029ce50919SHong Zhang /* Set SuperLU options */ 4039ce50919SHong Zhang /* the default values for options argument: 4049ce50919SHong Zhang options.Fact = DOFACT; 4059ce50919SHong Zhang options.Equil = YES; 4069ce50919SHong Zhang options.ColPerm = COLAMD; 4079ce50919SHong Zhang options.DiagPivotThresh = 1.0; 4089ce50919SHong Zhang options.Trans = NOTRANS; 4099ce50919SHong Zhang options.IterRefine = NOREFINE; 4109ce50919SHong Zhang options.SymmetricMode = NO; 4119ce50919SHong Zhang options.PivotGrowth = NO; 4129ce50919SHong Zhang options.ConditionNumber = NO; 4139ce50919SHong Zhang options.PrintStat = YES; 4149ce50919SHong Zhang */ 4159ce50919SHong Zhang set_default_options(&lu->options); 4168914a3f7SHong Zhang /* equilibration causes error in solve(), thus not supported here. See dgssvx.c for possible reason. */ 4178914a3f7SHong Zhang lu->options.Equil = NO; 4189ce50919SHong Zhang lu->options.PrintStat = NO; 4198914a3f7SHong Zhang lu->lwork = 0; /* allocate space internally by system malloc */ 4209ce50919SHong Zhang 4219ce50919SHong Zhang ierr = PetscOptionsBegin(A->comm,A->prefix,"SuperLU Options","Mat");CHKERRQ(ierr); 4229ce50919SHong Zhang /* 4234dc4c822SBarry Smith ierr = PetscOptionsTruth("-mat_superlu_equil","Equil","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 4248914a3f7SHong Zhang if (flg) lu->options.Equil = YES; -- not supported by the interface !!! 4259ce50919SHong Zhang */ 4268914a3f7SHong Zhang ierr = PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);CHKERRQ(ierr); 4279ce50919SHong Zhang if (flg) {lu->options.ColPerm = (colperm_t)indx;} 4288914a3f7SHong Zhang ierr = PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);CHKERRQ(ierr); 4299ce50919SHong Zhang if (flg) { lu->options.IterRefine = (IterRefine_t)indx;} 4304dc4c822SBarry Smith ierr = PetscOptionsTruth("-mat_superlu_symmetricmode","SymmetricMode","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 4319ce50919SHong Zhang if (flg) lu->options.SymmetricMode = YES; 4328914a3f7SHong Zhang ierr = PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&lu->options.DiagPivotThresh,PETSC_NULL);CHKERRQ(ierr); 4334dc4c822SBarry Smith ierr = PetscOptionsTruth("-mat_superlu_pivotgrowth","PivotGrowth","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 4349ce50919SHong Zhang if (flg) lu->options.PivotGrowth = YES; 4354dc4c822SBarry Smith ierr = PetscOptionsTruth("-mat_superlu_conditionnumber","ConditionNumber","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 4369ce50919SHong Zhang if (flg) lu->options.ConditionNumber = YES; 4378914a3f7SHong Zhang ierr = PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[0],&indx,&flg);CHKERRQ(ierr); 4389ce50919SHong Zhang if (flg) {lu->options.RowPerm = (rowperm_t)indx;} 4394dc4c822SBarry Smith ierr = PetscOptionsTruth("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 4409ce50919SHong Zhang if (flg) lu->options.ReplaceTinyPivot = YES; 4414dc4c822SBarry Smith ierr = PetscOptionsTruth("-mat_superlu_printstat","PrintStat","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 4429ce50919SHong Zhang if (flg) lu->options.PrintStat = YES; 4438914a3f7SHong Zhang ierr = PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,PETSC_NULL);CHKERRQ(ierr); 4445fe6150dSHong Zhang if (lu->lwork > 0 ){ 4455fe6150dSHong Zhang ierr = PetscMalloc(lu->lwork,&lu->work);CHKERRQ(ierr); 4465fe6150dSHong Zhang } else if (lu->lwork != 0 && lu->lwork != -1){ 44777431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," Warning: lwork %D is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork); 4488914a3f7SHong Zhang lu->lwork = 0; 4498914a3f7SHong Zhang } 4509ce50919SHong Zhang PetscOptionsEnd(); 4519ce50919SHong Zhang 45275af56d4SHong Zhang #ifdef SUPERLU2 4532ecb5974SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatCreateNull","MatCreateNull_SuperLU", 454f0c56d0fSKris Buschelman (void(*)(void))MatCreateNull_SuperLU);CHKERRQ(ierr); 45575af56d4SHong Zhang #endif 45614b396bbSKris Buschelman 45775af56d4SHong Zhang /* Allocate spaces (notice sizes are for the transpose) */ 458da7a1d00SHong Zhang ierr = PetscMalloc(m*sizeof(PetscInt),&lu->etree);CHKERRQ(ierr); 459da7a1d00SHong Zhang ierr = PetscMalloc(n*sizeof(PetscInt),&lu->perm_r);CHKERRQ(ierr); 460da7a1d00SHong Zhang ierr = PetscMalloc(m*sizeof(PetscInt),&lu->perm_c);CHKERRQ(ierr); 461da7a1d00SHong Zhang ierr = PetscMalloc(n*sizeof(PetscInt),&lu->R);CHKERRQ(ierr); 462da7a1d00SHong Zhang ierr = PetscMalloc(m*sizeof(PetscInt),&lu->C);CHKERRQ(ierr); 46375af56d4SHong Zhang 4649ce50919SHong Zhang /* create rhs and solution x without allocate space for .Store */ 4655fe6150dSHong Zhang #if defined(PETSC_USE_COMPLEX) 466937a6b0eSHong Zhang zCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE); 467937a6b0eSHong Zhang zCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE); 4685fe6150dSHong Zhang #else 46975af56d4SHong Zhang dCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE); 47075af56d4SHong Zhang dCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE); 4715fe6150dSHong Zhang #endif 47214b396bbSKris Buschelman 47314b396bbSKris Buschelman lu->flg = DIFFERENT_NONZERO_PATTERN; 474e740cb95SKris Buschelman lu->CleanUpSuperLU = PETSC_TRUE; 475899d7b4fSKris Buschelman 476899d7b4fSKris Buschelman *F = B; 4772a4c71feSBarry Smith ierr = PetscLogObjectMemory(B,(A->rmap.n+A->cmap.n)*sizeof(PetscInt)+sizeof(Mat_SuperLU));CHKERRQ(ierr); 47814b396bbSKris Buschelman PetscFunctionReturn(0); 47914b396bbSKris Buschelman } 48014b396bbSKris Buschelman 48194b7f48cSBarry Smith /* used by -ksp_view */ 48214b396bbSKris Buschelman #undef __FUNCT__ 483f0c56d0fSKris Buschelman #define __FUNCT__ "MatFactorInfo_SuperLU" 484dfbe8321SBarry Smith PetscErrorCode MatFactorInfo_SuperLU(Mat A,PetscViewer viewer) 48514b396bbSKris Buschelman { 486f0c56d0fSKris Buschelman Mat_SuperLU *lu= (Mat_SuperLU*)A->spptr; 487da7a1d00SHong Zhang PetscErrorCode ierr; 4889ce50919SHong Zhang superlu_options_t options; 489f0c56d0fSKris Buschelman 490f0c56d0fSKris Buschelman PetscFunctionBegin; 4919ce50919SHong Zhang /* check if matrix is superlu_dist type */ 4929ce50919SHong Zhang if (A->ops->solve != MatSolve_SuperLU) PetscFunctionReturn(0); 4939ce50919SHong Zhang 4949ce50919SHong Zhang options = lu->options; 495f0c56d0fSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");CHKERRQ(ierr); 4969ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Equil: %s\n",(options.Equil != NO) ? "YES": "NO");CHKERRQ(ierr); 49777431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ColPerm: %D\n",options.ColPerm);CHKERRQ(ierr); 49877431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," IterRefine: %D\n",options.IterRefine);CHKERRQ(ierr); 4999ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," SymmetricMode: %s\n",(options.SymmetricMode != NO) ? "YES": "NO");CHKERRQ(ierr); 5009ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," DiagPivotThresh: %g\n",options.DiagPivotThresh);CHKERRQ(ierr); 5019ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," PivotGrowth: %s\n",(options.PivotGrowth != NO) ? "YES": "NO");CHKERRQ(ierr); 5029ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ConditionNumber: %s\n",(options.ConditionNumber != NO) ? "YES": "NO");CHKERRQ(ierr); 50377431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," RowPerm: %D\n",options.RowPerm);CHKERRQ(ierr); 5049ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ReplaceTinyPivot: %s\n",(options.ReplaceTinyPivot != NO) ? "YES": "NO");CHKERRQ(ierr); 5059ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," PrintStat: %s\n",(options.PrintStat != NO) ? "YES": "NO");CHKERRQ(ierr); 50677431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," lwork: %D\n",lu->lwork);CHKERRQ(ierr); 5079ce50919SHong Zhang 508f0c56d0fSKris Buschelman PetscFunctionReturn(0); 509f0c56d0fSKris Buschelman } 510f0c56d0fSKris Buschelman 511f0c56d0fSKris Buschelman #undef __FUNCT__ 512f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_SuperLU" 513dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SuperLU(Mat A, MatDuplicateOption op, Mat *M) { 514da7a1d00SHong Zhang PetscErrorCode ierr; 5158f340917SKris Buschelman Mat_SuperLU *lu=(Mat_SuperLU *)A->spptr; 5168f340917SKris Buschelman 51714b396bbSKris Buschelman PetscFunctionBegin; 5188f340917SKris Buschelman ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr); 5193f953163SKris Buschelman ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_SuperLU));CHKERRQ(ierr); 52014b396bbSKris Buschelman PetscFunctionReturn(0); 52114b396bbSKris Buschelman } 52214b396bbSKris Buschelman 52314b396bbSKris Buschelman EXTERN_C_BEGIN 52414b396bbSKris Buschelman #undef __FUNCT__ 525b0592601SKris Buschelman #define __FUNCT__ "MatConvert_SuperLU_SeqAIJ" 52675179d2cSHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SuperLU_SeqAIJ(Mat A,MatType type,MatReuse reuse,Mat *newmat) 527521d7252SBarry Smith { 528fb3e25aaSKris Buschelman /* This routine is only called to convert an unfactored PETSc-SuperLU matrix */ 529fb3e25aaSKris Buschelman /* to its base PETSc type, so we will ignore 'MatType type'. */ 530da7a1d00SHong Zhang PetscErrorCode ierr; 531b0592601SKris Buschelman Mat B=*newmat; 532f0c56d0fSKris Buschelman Mat_SuperLU *lu=(Mat_SuperLU *)A->spptr; 533b0592601SKris Buschelman 534b0592601SKris Buschelman PetscFunctionBegin; 535ceb03754SKris Buschelman if (reuse == MAT_INITIAL_MATRIX) { 536b0592601SKris Buschelman ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 537f0c56d0fSKris Buschelman } 538b0592601SKris Buschelman /* Reset the original function pointers */ 539f0c56d0fSKris Buschelman B->ops->duplicate = lu->MatDuplicate; 540b0592601SKris Buschelman B->ops->view = lu->MatView; 541b0592601SKris Buschelman B->ops->assemblyend = lu->MatAssemblyEnd; 542b0592601SKris Buschelman B->ops->lufactorsymbolic = lu->MatLUFactorSymbolic; 543b0592601SKris Buschelman B->ops->destroy = lu->MatDestroy; 544901853e0SKris Buschelman 545901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_superlu_C","",PETSC_NULL);CHKERRQ(ierr); 546901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_superlu_seqaij_C","",PETSC_NULL);CHKERRQ(ierr); 547901853e0SKris Buschelman 548b0592601SKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);CHKERRQ(ierr); 549b0592601SKris Buschelman *newmat = B; 550b0592601SKris Buschelman PetscFunctionReturn(0); 551b0592601SKris Buschelman } 552b0592601SKris Buschelman EXTERN_C_END 553b0592601SKris Buschelman 554b0592601SKris Buschelman EXTERN_C_BEGIN 555b0592601SKris Buschelman #undef __FUNCT__ 556b0592601SKris Buschelman #define __FUNCT__ "MatConvert_SeqAIJ_SuperLU" 55775179d2cSHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SuperLU(Mat A,MatType type,MatReuse reuse,Mat *newmat) 558521d7252SBarry Smith { 559fb3e25aaSKris Buschelman /* This routine is only called to convert to MATSUPERLU */ 560fb3e25aaSKris Buschelman /* from MATSEQAIJ, so we will ignore 'MatType type'. */ 561da7a1d00SHong Zhang PetscErrorCode ierr; 562b0592601SKris Buschelman Mat B=*newmat; 563f0c56d0fSKris Buschelman Mat_SuperLU *lu; 56414b396bbSKris Buschelman 56514b396bbSKris Buschelman PetscFunctionBegin; 566ceb03754SKris Buschelman if (reuse == MAT_INITIAL_MATRIX) { 567b0592601SKris Buschelman ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 568b0592601SKris Buschelman } 56914b396bbSKris Buschelman 570f0c56d0fSKris Buschelman ierr = PetscNew(Mat_SuperLU,&lu);CHKERRQ(ierr); 571f0c56d0fSKris Buschelman lu->MatDuplicate = A->ops->duplicate; 57214b396bbSKris Buschelman lu->MatView = A->ops->view; 57314b396bbSKris Buschelman lu->MatAssemblyEnd = A->ops->assemblyend; 574b0592601SKris Buschelman lu->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 57514b396bbSKris Buschelman lu->MatDestroy = A->ops->destroy; 57614b396bbSKris Buschelman lu->CleanUpSuperLU = PETSC_FALSE; 577b0592601SKris Buschelman 578b0592601SKris Buschelman B->spptr = (void*)lu; 579f0c56d0fSKris Buschelman B->ops->duplicate = MatDuplicate_SuperLU; 580f0c56d0fSKris Buschelman B->ops->view = MatView_SuperLU; 581f0c56d0fSKris Buschelman B->ops->assemblyend = MatAssemblyEnd_SuperLU; 582f0c56d0fSKris Buschelman B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU; 583e4e36384SHong Zhang B->ops->choleskyfactorsymbolic = 0; 584f0c56d0fSKris Buschelman B->ops->destroy = MatDestroy_SuperLU; 585b0592601SKris Buschelman 586b0592601SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_superlu_C", 587b0592601SKris Buschelman "MatConvert_SeqAIJ_SuperLU",MatConvert_SeqAIJ_SuperLU);CHKERRQ(ierr); 588b0592601SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_superlu_seqaij_C", 589b0592601SKris Buschelman "MatConvert_SuperLU_SeqAIJ",MatConvert_SuperLU_SeqAIJ);CHKERRQ(ierr); 590ae15b995SBarry Smith ierr = PetscInfo(0,"Using SuperLU for SeqAIJ LU factorization and solves.\n");CHKERRQ(ierr); 591fb3e25aaSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)B,MATSUPERLU);CHKERRQ(ierr); 592b0592601SKris Buschelman *newmat = B; 593b0592601SKris Buschelman PetscFunctionReturn(0); 594b0592601SKris Buschelman } 595b0592601SKris Buschelman EXTERN_C_END 596b0592601SKris Buschelman 59724b6179bSKris Buschelman /*MC 598fafad747SKris Buschelman MATSUPERLU - MATSUPERLU = "superlu" - A matrix type providing direct solvers (LU) for sequential matrices 59924b6179bSKris Buschelman via the external package SuperLU. 60024b6179bSKris Buschelman 60124b6179bSKris Buschelman If SuperLU is installed (see the manual for 60224b6179bSKris Buschelman instructions on how to declare the existence of external packages), 60324b6179bSKris Buschelman a matrix type can be constructed which invokes SuperLU solvers. 60424b6179bSKris Buschelman After calling MatCreate(...,A), simply call MatSetType(A,MATSUPERLU). 60524b6179bSKris Buschelman 60624b6179bSKris Buschelman This matrix inherits from MATSEQAIJ. As a result, MatSeqAIJSetPreallocation is 60728b08bd3SKris Buschelman supported for this matrix type. One can also call MatConvert for an inplace conversion to or from 60828b08bd3SKris Buschelman the MATSEQAIJ type without data copy. 60924b6179bSKris Buschelman 61024b6179bSKris Buschelman Options Database Keys: 6110bad9183SKris Buschelman + -mat_type superlu - sets the matrix type to "superlu" during a call to MatSetFromOptions() 61224b6179bSKris Buschelman - -mat_superlu_ordering <0,1,2,3> - 0: natural ordering, 61324b6179bSKris Buschelman 1: MMD applied to A'*A, 61424b6179bSKris Buschelman 2: MMD applied to A'+A, 61524b6179bSKris Buschelman 3: COLAMD, approximate minimum degree column ordering 61624b6179bSKris Buschelman 61724b6179bSKris Buschelman Level: beginner 61824b6179bSKris Buschelman 61924b6179bSKris Buschelman .seealso: PCLU 62024b6179bSKris Buschelman M*/ 62124b6179bSKris Buschelman 622b0592601SKris Buschelman EXTERN_C_BEGIN 623b0592601SKris Buschelman #undef __FUNCT__ 624f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_SuperLU" 625be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SuperLU(Mat A) 626dfbe8321SBarry Smith { 627dfbe8321SBarry Smith PetscErrorCode ierr; 628b0592601SKris Buschelman 629b0592601SKris Buschelman PetscFunctionBegin; 6305441df8eSKris Buschelman /* Change type name before calling MatSetType to force proper construction of SeqAIJ and SUPERLU types */ 6315441df8eSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)A,MATSUPERLU);CHKERRQ(ierr); 632b0592601SKris Buschelman ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr); 633ceb03754SKris Buschelman ierr = MatConvert_SeqAIJ_SuperLU(A,MATSUPERLU,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr); 63414b396bbSKris Buschelman PetscFunctionReturn(0); 63514b396bbSKris Buschelman } 63614b396bbSKris Buschelman EXTERN_C_END 637