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) 1114b396bbSKris Buschelman #include "zsp_defs.h" 1214b396bbSKris Buschelman #else 1314b396bbSKris Buschelman #include "dsp_defs.h" 1414b396bbSKris Buschelman #endif 1514b396bbSKris Buschelman #include "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 48be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SuperLU_SeqAIJ(Mat,const MatType,MatReuse,Mat*); 49be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SuperLU(Mat,const 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); 9514b396bbSKris Buschelman if (format == PETSC_VIEWER_ASCII_FACTOR_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; 123da7a1d00SHong Zhang PetscInt numRows = A->m,numCols = A->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 */ 140*f69a0ea3SMatthew Knepley ierr = MatCreate(A->comm,nullMat);CHKERRQ(ierr); 141*f69a0ea3SMatthew 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__ 196f0c56d0fSKris Buschelman #define __FUNCT__ "MatSolve_SuperLU" 197dfbe8321SBarry Smith PetscErrorCode MatSolve_SuperLU(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. */ 22675af56d4SHong Zhang lu->options.Trans = TRANS; 2278914a3f7SHong Zhang #if defined(PETSC_USE_COMPLEX) 2288914a3f7SHong Zhang zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 2298914a3f7SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 2308914a3f7SHong Zhang &lu->mem_usage, &stat, &info); 2318914a3f7SHong Zhang #else 23275af56d4SHong Zhang dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 23375af56d4SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 23475af56d4SHong Zhang &lu->mem_usage, &stat, &info); 2358914a3f7SHong Zhang #endif 23675af56d4SHong Zhang ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr); 23775af56d4SHong Zhang ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 23875af56d4SHong Zhang 239958c9bccSBarry Smith if ( !info || info == lu->A.ncol+1 ) { 24075af56d4SHong Zhang if ( lu->options.IterRefine ) { 2418914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n"); 2428914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR"); 24375af56d4SHong Zhang for (i = 0; i < 1; ++i) 2448914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," %8d%8d%16e%16e\n", i+1, stat.RefineSteps, ferr, berr); 24575af56d4SHong Zhang } 2468914a3f7SHong Zhang } else if ( info > 0 ){ 2478914a3f7SHong Zhang if ( lu->lwork == -1 ) { 24877431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %D bytes\n", info - lu->A.ncol); 2498914a3f7SHong Zhang } else { 25077431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," Warning: gssvx() returns info %D\n",info); 2518914a3f7SHong Zhang } 2528914a3f7SHong Zhang } else if (info < 0){ 25377431f27SBarry Smith SETERRQ2(PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", info,-info); 25414b396bbSKris Buschelman } 25514b396bbSKris Buschelman 2568914a3f7SHong Zhang if ( lu->options.PrintStat ) { 2578914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n"); 2588914a3f7SHong Zhang StatPrint(&stat); 2598914a3f7SHong Zhang } 26075af56d4SHong Zhang StatFree(&stat); 26175af56d4SHong Zhang PetscFunctionReturn(0); 26275af56d4SHong Zhang } 26314b396bbSKris Buschelman 26414b396bbSKris Buschelman #undef __FUNCT__ 265f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorNumeric_SuperLU" 266af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SuperLU(Mat A,MatFactorInfo *info,Mat *F) 26714b396bbSKris Buschelman { 26814b396bbSKris Buschelman Mat_SeqAIJ *aa = (Mat_SeqAIJ*)(A)->data; 269f0c56d0fSKris Buschelman Mat_SuperLU *lu = (Mat_SuperLU*)(*F)->spptr; 270dfbe8321SBarry Smith PetscErrorCode ierr; 271da7a1d00SHong Zhang PetscInt sinfo; 27275af56d4SHong Zhang SuperLUStat_t stat; 273da7a1d00SHong Zhang PetscReal ferr, berr; 2748914a3f7SHong Zhang NCformat *Ustore; 2758914a3f7SHong Zhang SCformat *Lstore; 27614b396bbSKris Buschelman 27714b396bbSKris Buschelman PetscFunctionBegin; 2789ce50919SHong Zhang if (lu->flg == SAME_NONZERO_PATTERN){ /* successing numerical factorization */ 2799ce50919SHong Zhang lu->options.Fact = SamePattern; 28075af56d4SHong Zhang /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */ 28175af56d4SHong Zhang Destroy_SuperMatrix_Store(&lu->A); 2829ce50919SHong Zhang if ( lu->lwork >= 0 ) { 28375af56d4SHong Zhang Destroy_SuperNode_Matrix(&lu->L); 28475af56d4SHong Zhang Destroy_CompCol_Matrix(&lu->U); 28575af56d4SHong Zhang lu->options.Fact = SamePattern; 28614b396bbSKris Buschelman } 28775af56d4SHong Zhang } 28814b396bbSKris Buschelman 28975af56d4SHong Zhang /* Create the SuperMatrix for lu->A=A^T: 29075af56d4SHong Zhang Since SuperLU likes column-oriented matrices,we pass it the transpose, 29175af56d4SHong Zhang and then solve A^T X = B in MatSolve(). */ 29275af56d4SHong Zhang #if defined(PETSC_USE_COMPLEX) 2935fe6150dSHong Zhang zCreate_CompCol_Matrix(&lu->A,A->n,A->m,aa->nz,(doublecomplex*)aa->a,aa->j,aa->i, 29475af56d4SHong Zhang SLU_NC,SLU_Z,SLU_GE); 29575af56d4SHong Zhang #else 29675af56d4SHong Zhang dCreate_CompCol_Matrix(&lu->A,A->n,A->m,aa->nz,aa->a,aa->j,aa->i, 29775af56d4SHong Zhang SLU_NC,SLU_D,SLU_GE); 29875af56d4SHong Zhang #endif 29975af56d4SHong Zhang 30075af56d4SHong Zhang /* Initialize the statistics variables. */ 30175af56d4SHong Zhang StatInit(&stat); 30275af56d4SHong Zhang 3039ce50919SHong Zhang /* Numerical factorization */ 30475af56d4SHong Zhang lu->B.ncol = 0; /* Indicate not to solve the system */ 3058914a3f7SHong Zhang #if defined(PETSC_USE_COMPLEX) 3068914a3f7SHong Zhang zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 3078914a3f7SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 308da7a1d00SHong Zhang &lu->mem_usage, &stat, &sinfo); 3098914a3f7SHong Zhang #else 31075af56d4SHong Zhang dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 31175af56d4SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 312da7a1d00SHong Zhang &lu->mem_usage, &stat, &sinfo); 3138914a3f7SHong Zhang #endif 314da7a1d00SHong Zhang if ( !sinfo || sinfo == lu->A.ncol+1 ) { 3158914a3f7SHong Zhang if ( lu->options.PivotGrowth ) 3168914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," Recip. pivot growth = %e\n", lu->rpg); 31775af56d4SHong Zhang if ( lu->options.ConditionNumber ) 3188914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," Recip. condition number = %e\n", lu->rcond); 319da7a1d00SHong Zhang } else if ( sinfo > 0 ){ 3208914a3f7SHong Zhang if ( lu->lwork == -1 ) { 321da7a1d00SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %D bytes\n", sinfo - lu->A.ncol); 3228914a3f7SHong Zhang } else { 323da7a1d00SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," Warning: gssvx() returns info %D\n",sinfo); 3248914a3f7SHong Zhang } 325da7a1d00SHong Zhang } else { /* sinfo < 0 */ 326da7a1d00SHong Zhang SETERRQ2(PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", sinfo,-sinfo); 32775af56d4SHong Zhang } 32875af56d4SHong Zhang 3298914a3f7SHong Zhang if ( lu->options.PrintStat ) { 3308914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"MatLUFactorNumeric_SuperLU():\n"); 3318914a3f7SHong Zhang StatPrint(&stat); 3328914a3f7SHong Zhang Lstore = (SCformat *) lu->L.Store; 3338914a3f7SHong Zhang Ustore = (NCformat *) lu->U.Store; 33477431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor L = %D\n", Lstore->nnz); 33577431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor U = %D\n", Ustore->nnz); 33677431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in L+U = %D\n", Lstore->nnz + Ustore->nnz - lu->A.ncol); 33777431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," L\\U MB %.3f\ttotal MB needed %.3f\texpansions %D\n", 3388914a3f7SHong Zhang lu->mem_usage.for_lu/1e6, lu->mem_usage.total_needed/1e6, 3398914a3f7SHong Zhang lu->mem_usage.expansions); 3408914a3f7SHong Zhang } 34175af56d4SHong Zhang StatFree(&stat); 34275af56d4SHong Zhang 34375af56d4SHong Zhang lu->flg = SAME_NONZERO_PATTERN; 34475af56d4SHong Zhang PetscFunctionReturn(0); 34514b396bbSKris Buschelman } 34614b396bbSKris Buschelman 34714b396bbSKris Buschelman /* 34814b396bbSKris Buschelman Note the r permutation is ignored 34914b396bbSKris Buschelman */ 35014b396bbSKris Buschelman #undef __FUNCT__ 351f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_SuperLU" 352dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) 35314b396bbSKris Buschelman { 35414b396bbSKris Buschelman Mat B; 355f0c56d0fSKris Buschelman Mat_SuperLU *lu; 3566849ba73SBarry Smith PetscErrorCode ierr; 357da7a1d00SHong Zhang PetscInt m=A->m,n=A->n,indx; 3589ce50919SHong Zhang PetscTruth flg; 3595ca28756SHong Zhang const char *colperm[]={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */ 3605ca28756SHong Zhang const char *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"}; 3615ca28756SHong Zhang const char *rowperm[]={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */ 36214b396bbSKris Buschelman 36314b396bbSKris Buschelman PetscFunctionBegin; 364*f69a0ea3SMatthew Knepley ierr = MatCreate(A->comm,&B);CHKERRQ(ierr); 365*f69a0ea3SMatthew Knepley ierr = MatSetSizes(B,A->m,A->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 366be5d1d56SKris Buschelman ierr = MatSetType(B,A->type_name);CHKERRQ(ierr); 367899d7b4fSKris Buschelman ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 368f0c56d0fSKris Buschelman 369f0c56d0fSKris Buschelman B->ops->lufactornumeric = MatLUFactorNumeric_SuperLU; 370f0c56d0fSKris Buschelman B->ops->solve = MatSolve_SuperLU; 37114b396bbSKris Buschelman B->factor = FACTOR_LU; 37294b7f48cSBarry Smith B->assembled = PETSC_TRUE; /* required by -ksp_view */ 37314b396bbSKris Buschelman 374f0c56d0fSKris Buschelman lu = (Mat_SuperLU*)(B->spptr); 3759ce50919SHong Zhang 3769ce50919SHong Zhang /* Set SuperLU options */ 3779ce50919SHong Zhang /* the default values for options argument: 3789ce50919SHong Zhang options.Fact = DOFACT; 3799ce50919SHong Zhang options.Equil = YES; 3809ce50919SHong Zhang options.ColPerm = COLAMD; 3819ce50919SHong Zhang options.DiagPivotThresh = 1.0; 3829ce50919SHong Zhang options.Trans = NOTRANS; 3839ce50919SHong Zhang options.IterRefine = NOREFINE; 3849ce50919SHong Zhang options.SymmetricMode = NO; 3859ce50919SHong Zhang options.PivotGrowth = NO; 3869ce50919SHong Zhang options.ConditionNumber = NO; 3879ce50919SHong Zhang options.PrintStat = YES; 3889ce50919SHong Zhang */ 3899ce50919SHong Zhang set_default_options(&lu->options); 3908914a3f7SHong Zhang /* equilibration causes error in solve(), thus not supported here. See dgssvx.c for possible reason. */ 3918914a3f7SHong Zhang lu->options.Equil = NO; 3929ce50919SHong Zhang lu->options.PrintStat = NO; 3938914a3f7SHong Zhang lu->lwork = 0; /* allocate space internally by system malloc */ 3949ce50919SHong Zhang 3959ce50919SHong Zhang ierr = PetscOptionsBegin(A->comm,A->prefix,"SuperLU Options","Mat");CHKERRQ(ierr); 3969ce50919SHong Zhang /* 3978914a3f7SHong Zhang ierr = PetscOptionsLogical("-mat_superlu_equil","Equil","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 3988914a3f7SHong Zhang if (flg) lu->options.Equil = YES; -- not supported by the interface !!! 3999ce50919SHong Zhang */ 4008914a3f7SHong Zhang ierr = PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);CHKERRQ(ierr); 4019ce50919SHong Zhang if (flg) {lu->options.ColPerm = (colperm_t)indx;} 4028914a3f7SHong Zhang ierr = PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);CHKERRQ(ierr); 4039ce50919SHong Zhang if (flg) { lu->options.IterRefine = (IterRefine_t)indx;} 4048914a3f7SHong Zhang ierr = PetscOptionsLogical("-mat_superlu_symmetricmode","SymmetricMode","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 4059ce50919SHong Zhang if (flg) lu->options.SymmetricMode = YES; 4068914a3f7SHong Zhang ierr = PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&lu->options.DiagPivotThresh,PETSC_NULL);CHKERRQ(ierr); 4078914a3f7SHong Zhang ierr = PetscOptionsLogical("-mat_superlu_pivotgrowth","PivotGrowth","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 4089ce50919SHong Zhang if (flg) lu->options.PivotGrowth = YES; 4098914a3f7SHong Zhang ierr = PetscOptionsLogical("-mat_superlu_conditionnumber","ConditionNumber","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 4109ce50919SHong Zhang if (flg) lu->options.ConditionNumber = YES; 4118914a3f7SHong Zhang ierr = PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[0],&indx,&flg);CHKERRQ(ierr); 4129ce50919SHong Zhang if (flg) {lu->options.RowPerm = (rowperm_t)indx;} 4138914a3f7SHong Zhang ierr = PetscOptionsLogical("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 4149ce50919SHong Zhang if (flg) lu->options.ReplaceTinyPivot = YES; 4158914a3f7SHong Zhang ierr = PetscOptionsLogical("-mat_superlu_printstat","PrintStat","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 4169ce50919SHong Zhang if (flg) lu->options.PrintStat = YES; 4178914a3f7SHong Zhang ierr = PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,PETSC_NULL);CHKERRQ(ierr); 4185fe6150dSHong Zhang if (lu->lwork > 0 ){ 4195fe6150dSHong Zhang ierr = PetscMalloc(lu->lwork,&lu->work);CHKERRQ(ierr); 4205fe6150dSHong Zhang } else if (lu->lwork != 0 && lu->lwork != -1){ 42177431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," Warning: lwork %D is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork); 4228914a3f7SHong Zhang lu->lwork = 0; 4238914a3f7SHong Zhang } 4249ce50919SHong Zhang PetscOptionsEnd(); 4259ce50919SHong Zhang 42675af56d4SHong Zhang #ifdef SUPERLU2 4272ecb5974SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatCreateNull","MatCreateNull_SuperLU", 428f0c56d0fSKris Buschelman (void(*)(void))MatCreateNull_SuperLU);CHKERRQ(ierr); 42975af56d4SHong Zhang #endif 43014b396bbSKris Buschelman 43175af56d4SHong Zhang /* Allocate spaces (notice sizes are for the transpose) */ 432da7a1d00SHong Zhang ierr = PetscMalloc(m*sizeof(PetscInt),&lu->etree);CHKERRQ(ierr); 433da7a1d00SHong Zhang ierr = PetscMalloc(n*sizeof(PetscInt),&lu->perm_r);CHKERRQ(ierr); 434da7a1d00SHong Zhang ierr = PetscMalloc(m*sizeof(PetscInt),&lu->perm_c);CHKERRQ(ierr); 435da7a1d00SHong Zhang ierr = PetscMalloc(n*sizeof(PetscInt),&lu->R);CHKERRQ(ierr); 436da7a1d00SHong Zhang ierr = PetscMalloc(m*sizeof(PetscInt),&lu->C);CHKERRQ(ierr); 43775af56d4SHong Zhang 4389ce50919SHong Zhang /* create rhs and solution x without allocate space for .Store */ 4395fe6150dSHong Zhang #if defined(PETSC_USE_COMPLEX) 440937a6b0eSHong Zhang zCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE); 441937a6b0eSHong Zhang zCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE); 4425fe6150dSHong Zhang #else 44375af56d4SHong Zhang dCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE); 44475af56d4SHong Zhang dCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE); 4455fe6150dSHong Zhang #endif 44614b396bbSKris Buschelman 44714b396bbSKris Buschelman lu->flg = DIFFERENT_NONZERO_PATTERN; 448e740cb95SKris Buschelman lu->CleanUpSuperLU = PETSC_TRUE; 449899d7b4fSKris Buschelman 450899d7b4fSKris Buschelman *F = B; 451da7a1d00SHong Zhang ierr = PetscLogObjectMemory(B,(A->m+A->n)*sizeof(PetscInt)+sizeof(Mat_SuperLU));CHKERRQ(ierr); 45214b396bbSKris Buschelman PetscFunctionReturn(0); 45314b396bbSKris Buschelman } 45414b396bbSKris Buschelman 45594b7f48cSBarry Smith /* used by -ksp_view */ 45614b396bbSKris Buschelman #undef __FUNCT__ 457f0c56d0fSKris Buschelman #define __FUNCT__ "MatFactorInfo_SuperLU" 458dfbe8321SBarry Smith PetscErrorCode MatFactorInfo_SuperLU(Mat A,PetscViewer viewer) 45914b396bbSKris Buschelman { 460f0c56d0fSKris Buschelman Mat_SuperLU *lu= (Mat_SuperLU*)A->spptr; 461da7a1d00SHong Zhang PetscErrorCode ierr; 4629ce50919SHong Zhang superlu_options_t options; 463f0c56d0fSKris Buschelman 464f0c56d0fSKris Buschelman PetscFunctionBegin; 4659ce50919SHong Zhang /* check if matrix is superlu_dist type */ 4669ce50919SHong Zhang if (A->ops->solve != MatSolve_SuperLU) PetscFunctionReturn(0); 4679ce50919SHong Zhang 4689ce50919SHong Zhang options = lu->options; 469f0c56d0fSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");CHKERRQ(ierr); 4709ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Equil: %s\n",(options.Equil != NO) ? "YES": "NO");CHKERRQ(ierr); 47177431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ColPerm: %D\n",options.ColPerm);CHKERRQ(ierr); 47277431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," IterRefine: %D\n",options.IterRefine);CHKERRQ(ierr); 4739ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," SymmetricMode: %s\n",(options.SymmetricMode != NO) ? "YES": "NO");CHKERRQ(ierr); 4749ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," DiagPivotThresh: %g\n",options.DiagPivotThresh);CHKERRQ(ierr); 4759ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," PivotGrowth: %s\n",(options.PivotGrowth != NO) ? "YES": "NO");CHKERRQ(ierr); 4769ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ConditionNumber: %s\n",(options.ConditionNumber != NO) ? "YES": "NO");CHKERRQ(ierr); 47777431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," RowPerm: %D\n",options.RowPerm);CHKERRQ(ierr); 4789ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ReplaceTinyPivot: %s\n",(options.ReplaceTinyPivot != NO) ? "YES": "NO");CHKERRQ(ierr); 4799ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," PrintStat: %s\n",(options.PrintStat != NO) ? "YES": "NO");CHKERRQ(ierr); 48077431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," lwork: %D\n",lu->lwork);CHKERRQ(ierr); 4819ce50919SHong Zhang 482f0c56d0fSKris Buschelman PetscFunctionReturn(0); 483f0c56d0fSKris Buschelman } 484f0c56d0fSKris Buschelman 485f0c56d0fSKris Buschelman #undef __FUNCT__ 486f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_SuperLU" 487dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SuperLU(Mat A, MatDuplicateOption op, Mat *M) { 488da7a1d00SHong Zhang PetscErrorCode ierr; 4898f340917SKris Buschelman Mat_SuperLU *lu=(Mat_SuperLU *)A->spptr; 4908f340917SKris Buschelman 49114b396bbSKris Buschelman PetscFunctionBegin; 4928f340917SKris Buschelman ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr); 4933f953163SKris Buschelman ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_SuperLU));CHKERRQ(ierr); 49414b396bbSKris Buschelman PetscFunctionReturn(0); 49514b396bbSKris Buschelman } 49614b396bbSKris Buschelman 49714b396bbSKris Buschelman EXTERN_C_BEGIN 49814b396bbSKris Buschelman #undef __FUNCT__ 499b0592601SKris Buschelman #define __FUNCT__ "MatConvert_SuperLU_SeqAIJ" 500be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SuperLU_SeqAIJ(Mat A,const MatType type,MatReuse reuse,Mat *newmat) 501521d7252SBarry Smith { 502fb3e25aaSKris Buschelman /* This routine is only called to convert an unfactored PETSc-SuperLU matrix */ 503fb3e25aaSKris Buschelman /* to its base PETSc type, so we will ignore 'MatType type'. */ 504da7a1d00SHong Zhang PetscErrorCode ierr; 505b0592601SKris Buschelman Mat B=*newmat; 506f0c56d0fSKris Buschelman Mat_SuperLU *lu=(Mat_SuperLU *)A->spptr; 507b0592601SKris Buschelman 508b0592601SKris Buschelman PetscFunctionBegin; 509ceb03754SKris Buschelman if (reuse == MAT_INITIAL_MATRIX) { 510b0592601SKris Buschelman ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 511f0c56d0fSKris Buschelman } 512b0592601SKris Buschelman /* Reset the original function pointers */ 513f0c56d0fSKris Buschelman B->ops->duplicate = lu->MatDuplicate; 514b0592601SKris Buschelman B->ops->view = lu->MatView; 515b0592601SKris Buschelman B->ops->assemblyend = lu->MatAssemblyEnd; 516b0592601SKris Buschelman B->ops->lufactorsymbolic = lu->MatLUFactorSymbolic; 517b0592601SKris Buschelman B->ops->destroy = lu->MatDestroy; 518b0592601SKris Buschelman /* lu is only a function pointer stash unless we've factored the matrix, which we haven't! */ 519b0592601SKris Buschelman ierr = PetscFree(lu);CHKERRQ(ierr); 520901853e0SKris Buschelman 521901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_superlu_C","",PETSC_NULL);CHKERRQ(ierr); 522901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_superlu_seqaij_C","",PETSC_NULL);CHKERRQ(ierr); 523901853e0SKris Buschelman 524b0592601SKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);CHKERRQ(ierr); 525b0592601SKris Buschelman *newmat = B; 526b0592601SKris Buschelman PetscFunctionReturn(0); 527b0592601SKris Buschelman } 528b0592601SKris Buschelman EXTERN_C_END 529b0592601SKris Buschelman 530b0592601SKris Buschelman EXTERN_C_BEGIN 531b0592601SKris Buschelman #undef __FUNCT__ 532b0592601SKris Buschelman #define __FUNCT__ "MatConvert_SeqAIJ_SuperLU" 533be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SuperLU(Mat A,const MatType type,MatReuse reuse,Mat *newmat) 534521d7252SBarry Smith { 535fb3e25aaSKris Buschelman /* This routine is only called to convert to MATSUPERLU */ 536fb3e25aaSKris Buschelman /* from MATSEQAIJ, so we will ignore 'MatType type'. */ 537da7a1d00SHong Zhang PetscErrorCode ierr; 538b0592601SKris Buschelman Mat B=*newmat; 539f0c56d0fSKris Buschelman Mat_SuperLU *lu; 54014b396bbSKris Buschelman 54114b396bbSKris Buschelman PetscFunctionBegin; 542ceb03754SKris Buschelman if (reuse == MAT_INITIAL_MATRIX) { 543b0592601SKris Buschelman ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 544b0592601SKris Buschelman } 54514b396bbSKris Buschelman 546f0c56d0fSKris Buschelman ierr = PetscNew(Mat_SuperLU,&lu);CHKERRQ(ierr); 547f0c56d0fSKris Buschelman lu->MatDuplicate = A->ops->duplicate; 54814b396bbSKris Buschelman lu->MatView = A->ops->view; 54914b396bbSKris Buschelman lu->MatAssemblyEnd = A->ops->assemblyend; 550b0592601SKris Buschelman lu->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 55114b396bbSKris Buschelman lu->MatDestroy = A->ops->destroy; 55214b396bbSKris Buschelman lu->CleanUpSuperLU = PETSC_FALSE; 553b0592601SKris Buschelman 554b0592601SKris Buschelman B->spptr = (void*)lu; 555f0c56d0fSKris Buschelman B->ops->duplicate = MatDuplicate_SuperLU; 556f0c56d0fSKris Buschelman B->ops->view = MatView_SuperLU; 557f0c56d0fSKris Buschelman B->ops->assemblyend = MatAssemblyEnd_SuperLU; 558f0c56d0fSKris Buschelman B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU; 559e4e36384SHong Zhang B->ops->choleskyfactorsymbolic = 0; 560f0c56d0fSKris Buschelman B->ops->destroy = MatDestroy_SuperLU; 561b0592601SKris Buschelman 562b0592601SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_superlu_C", 563b0592601SKris Buschelman "MatConvert_SeqAIJ_SuperLU",MatConvert_SeqAIJ_SuperLU);CHKERRQ(ierr); 564b0592601SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_superlu_seqaij_C", 565b0592601SKris Buschelman "MatConvert_SuperLU_SeqAIJ",MatConvert_SuperLU_SeqAIJ);CHKERRQ(ierr); 56663ba0a88SBarry Smith ierr = PetscLogInfo((0,"MatConvert_SeqAIJ_SuperLU:Using SuperLU for SeqAIJ LU factorization and solves.\n"));CHKERRQ(ierr); 567fb3e25aaSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)B,MATSUPERLU);CHKERRQ(ierr); 568b0592601SKris Buschelman *newmat = B; 569b0592601SKris Buschelman PetscFunctionReturn(0); 570b0592601SKris Buschelman } 571b0592601SKris Buschelman EXTERN_C_END 572b0592601SKris Buschelman 57324b6179bSKris Buschelman /*MC 574fafad747SKris Buschelman MATSUPERLU - MATSUPERLU = "superlu" - A matrix type providing direct solvers (LU) for sequential matrices 57524b6179bSKris Buschelman via the external package SuperLU. 57624b6179bSKris Buschelman 57724b6179bSKris Buschelman If SuperLU is installed (see the manual for 57824b6179bSKris Buschelman instructions on how to declare the existence of external packages), 57924b6179bSKris Buschelman a matrix type can be constructed which invokes SuperLU solvers. 58024b6179bSKris Buschelman After calling MatCreate(...,A), simply call MatSetType(A,MATSUPERLU). 58124b6179bSKris Buschelman 58224b6179bSKris Buschelman This matrix inherits from MATSEQAIJ. As a result, MatSeqAIJSetPreallocation is 58328b08bd3SKris Buschelman supported for this matrix type. One can also call MatConvert for an inplace conversion to or from 58428b08bd3SKris Buschelman the MATSEQAIJ type without data copy. 58524b6179bSKris Buschelman 58624b6179bSKris Buschelman Options Database Keys: 5870bad9183SKris Buschelman + -mat_type superlu - sets the matrix type to "superlu" during a call to MatSetFromOptions() 58824b6179bSKris Buschelman - -mat_superlu_ordering <0,1,2,3> - 0: natural ordering, 58924b6179bSKris Buschelman 1: MMD applied to A'*A, 59024b6179bSKris Buschelman 2: MMD applied to A'+A, 59124b6179bSKris Buschelman 3: COLAMD, approximate minimum degree column ordering 59224b6179bSKris Buschelman 59324b6179bSKris Buschelman Level: beginner 59424b6179bSKris Buschelman 59524b6179bSKris Buschelman .seealso: PCLU 59624b6179bSKris Buschelman M*/ 59724b6179bSKris Buschelman 598b0592601SKris Buschelman EXTERN_C_BEGIN 599b0592601SKris Buschelman #undef __FUNCT__ 600f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_SuperLU" 601be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SuperLU(Mat A) 602dfbe8321SBarry Smith { 603dfbe8321SBarry Smith PetscErrorCode ierr; 604b0592601SKris Buschelman 605b0592601SKris Buschelman PetscFunctionBegin; 6065441df8eSKris Buschelman /* Change type name before calling MatSetType to force proper construction of SeqAIJ and SUPERLU types */ 6075441df8eSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)A,MATSUPERLU);CHKERRQ(ierr); 608b0592601SKris Buschelman ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr); 609ceb03754SKris Buschelman ierr = MatConvert_SeqAIJ_SuperLU(A,MATSUPERLU,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr); 61014b396bbSKris Buschelman PetscFunctionReturn(0); 61114b396bbSKris Buschelman } 61214b396bbSKris Buschelman EXTERN_C_END 613