114b396bbSKris Buschelman 214b396bbSKris Buschelman /* 375af56d4SHong Zhang Provides an interface to the SuperLU 3.0 sparse solver 414b396bbSKris Buschelman */ 514b396bbSKris Buschelman 614b396bbSKris Buschelman #include "src/mat/impls/aij/seq/aij.h" 714b396bbSKris Buschelman 814b396bbSKris Buschelman EXTERN_C_BEGIN 914b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 1014b396bbSKris Buschelman #include "zsp_defs.h" 1114b396bbSKris Buschelman #else 1214b396bbSKris Buschelman #include "dsp_defs.h" 1314b396bbSKris Buschelman #endif 1414b396bbSKris Buschelman #include "util.h" 1514b396bbSKris Buschelman EXTERN_C_END 1614b396bbSKris Buschelman 1714b396bbSKris Buschelman typedef struct { 1875af56d4SHong Zhang SuperMatrix A,L,U,B,X; 1975af56d4SHong Zhang superlu_options_t options; 20*da7a1d00SHong Zhang PetscInt *perm_c; /* column permutation vector */ 21*da7a1d00SHong Zhang PetscInt *perm_r; /* row permutations from partial pivoting */ 22*da7a1d00SHong Zhang PetscInt *etree; 23*da7a1d00SHong Zhang PetscReal *R, *C; 2475af56d4SHong Zhang char equed[1]; 25*da7a1d00SHong Zhang PetscInt lwork; 2675af56d4SHong Zhang void *work; 27*da7a1d00SHong Zhang PetscReal rpg, rcond; 2875af56d4SHong Zhang mem_usage_t mem_usage; 2975af56d4SHong Zhang MatStructure flg; 3014b396bbSKris Buschelman 3114b396bbSKris Buschelman /* A few function pointers for inheritance */ 326849ba73SBarry Smith PetscErrorCode (*MatDuplicate)(Mat,MatDuplicateOption,Mat*); 336849ba73SBarry Smith PetscErrorCode (*MatView)(Mat,PetscViewer); 346849ba73SBarry Smith PetscErrorCode (*MatAssemblyEnd)(Mat,MatAssemblyType); 356849ba73SBarry Smith PetscErrorCode (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*); 366849ba73SBarry Smith PetscErrorCode (*MatDestroy)(Mat); 3714b396bbSKris Buschelman 3814b396bbSKris Buschelman /* Flag to clean up (non-global) SuperLU objects during Destroy */ 3914b396bbSKris Buschelman PetscTruth CleanUpSuperLU; 40f0c56d0fSKris Buschelman } Mat_SuperLU; 4114b396bbSKris Buschelman 4214b396bbSKris Buschelman 43dfbe8321SBarry Smith EXTERN PetscErrorCode MatFactorInfo_SuperLU(Mat,PetscViewer); 44dfbe8321SBarry Smith EXTERN PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat,IS,IS,MatFactorInfo*,Mat*); 45b0592601SKris Buschelman 46b0592601SKris Buschelman EXTERN_C_BEGIN 47ceb03754SKris Buschelman EXTERN PetscErrorCode MatConvert_SuperLU_SeqAIJ(Mat,const MatType,MatReuse,Mat*); 48ceb03754SKris Buschelman EXTERN PetscErrorCode MatConvert_SeqAIJ_SuperLU(Mat,const MatType,MatReuse,Mat*); 49b0592601SKris Buschelman EXTERN_C_END 5014b396bbSKris Buschelman 5114b396bbSKris Buschelman #undef __FUNCT__ 52f0c56d0fSKris Buschelman #define __FUNCT__ "MatDestroy_SuperLU" 53dfbe8321SBarry Smith PetscErrorCode MatDestroy_SuperLU(Mat A) 5414b396bbSKris Buschelman { 55dfbe8321SBarry Smith PetscErrorCode ierr; 56f0c56d0fSKris Buschelman Mat_SuperLU *lu=(Mat_SuperLU*)A->spptr; 5714b396bbSKris Buschelman 5814b396bbSKris Buschelman PetscFunctionBegin; 5975af56d4SHong Zhang if (lu->CleanUpSuperLU) { /* Free the SuperLU datastructures */ 6075af56d4SHong Zhang Destroy_SuperMatrix_Store(&lu->A); 6175af56d4SHong Zhang Destroy_SuperMatrix_Store(&lu->B); 6275af56d4SHong Zhang Destroy_SuperMatrix_Store(&lu->X); 639ce50919SHong Zhang 649ce50919SHong Zhang ierr = PetscFree(lu->etree);CHKERRQ(ierr); 659ce50919SHong Zhang ierr = PetscFree(lu->perm_r);CHKERRQ(ierr); 669ce50919SHong Zhang ierr = PetscFree(lu->perm_c);CHKERRQ(ierr); 679ce50919SHong Zhang ierr = PetscFree(lu->R);CHKERRQ(ierr); 689ce50919SHong Zhang ierr = PetscFree(lu->C);CHKERRQ(ierr); 6975af56d4SHong Zhang if ( lu->lwork >= 0 ) { 70fb3e25aaSKris Buschelman Destroy_SuperNode_Matrix(&lu->L); 71fb3e25aaSKris Buschelman Destroy_CompCol_Matrix(&lu->U); 7275af56d4SHong Zhang } 73fb3e25aaSKris Buschelman } 74ceb03754SKris Buschelman ierr = MatConvert_SuperLU_SeqAIJ(A,MATSEQAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr); 75fb3e25aaSKris Buschelman ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 7614b396bbSKris Buschelman PetscFunctionReturn(0); 7714b396bbSKris Buschelman } 7814b396bbSKris Buschelman 7914b396bbSKris Buschelman #undef __FUNCT__ 80f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_SuperLU" 81dfbe8321SBarry Smith PetscErrorCode MatView_SuperLU(Mat A,PetscViewer viewer) 8214b396bbSKris Buschelman { 83dfbe8321SBarry Smith PetscErrorCode ierr; 8432077d6dSBarry Smith PetscTruth iascii; 8514b396bbSKris Buschelman PetscViewerFormat format; 86f0c56d0fSKris Buschelman Mat_SuperLU *lu=(Mat_SuperLU*)(A->spptr); 8714b396bbSKris Buschelman 8814b396bbSKris Buschelman PetscFunctionBegin; 8914b396bbSKris Buschelman ierr = (*lu->MatView)(A,viewer);CHKERRQ(ierr); 9014b396bbSKris Buschelman 9132077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 9232077d6dSBarry Smith if (iascii) { 9314b396bbSKris Buschelman ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 9414b396bbSKris Buschelman if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) { 95f0c56d0fSKris Buschelman ierr = MatFactorInfo_SuperLU(A,viewer);CHKERRQ(ierr); 9614b396bbSKris Buschelman } 9714b396bbSKris Buschelman } 9814b396bbSKris Buschelman PetscFunctionReturn(0); 9914b396bbSKris Buschelman } 10014b396bbSKris Buschelman 10114b396bbSKris Buschelman #undef __FUNCT__ 102f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_SuperLU" 103dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_SuperLU(Mat A,MatAssemblyType mode) { 104dfbe8321SBarry Smith PetscErrorCode ierr; 105f0c56d0fSKris Buschelman Mat_SuperLU *lu=(Mat_SuperLU*)(A->spptr); 10614b396bbSKris Buschelman 10714b396bbSKris Buschelman PetscFunctionBegin; 10814b396bbSKris Buschelman ierr = (*lu->MatAssemblyEnd)(A,mode);CHKERRQ(ierr); 109b0592601SKris Buschelman lu->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 110f0c56d0fSKris Buschelman A->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU; 11114b396bbSKris Buschelman PetscFunctionReturn(0); 11214b396bbSKris Buschelman } 11314b396bbSKris Buschelman 11475af56d4SHong Zhang /* This function was written for SuperLU 2.0 by Matthew Knepley. Not tested for SuperLU 3.0! */ 11575af56d4SHong Zhang #ifdef SuperLU2 11614b396bbSKris Buschelman #include "src/mat/impls/dense/seq/dense.h" 11714b396bbSKris Buschelman #undef __FUNCT__ 118f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreateNull_SuperLU" 119dfbe8321SBarry Smith PetscErrorCode MatCreateNull_SuperLU(Mat A,Mat *nullMat) 12014b396bbSKris Buschelman { 121f0c56d0fSKris Buschelman Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 122*da7a1d00SHong Zhang PetscInt numRows = A->m,numCols = A->n; 12314b396bbSKris Buschelman SCformat *Lstore; 124*da7a1d00SHong Zhang PetscInt numNullCols,size; 12575af56d4SHong Zhang SuperLUStat_t stat; 12614b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 12714b396bbSKris Buschelman doublecomplex *nullVals,*workVals; 12814b396bbSKris Buschelman #else 12914b396bbSKris Buschelman PetscScalar *nullVals,*workVals; 13014b396bbSKris Buschelman #endif 131*da7a1d00SHong Zhang PetscInt row,newRow,col,newCol,block,b; 1326849ba73SBarry Smith PetscErrorCode ierr; 13314b396bbSKris Buschelman 13414b396bbSKris Buschelman PetscFunctionBegin; 13514b396bbSKris Buschelman if (!A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix"); 13614b396bbSKris Buschelman numNullCols = numCols - numRows; 13714b396bbSKris Buschelman if (numNullCols < 0) SETERRQ(PETSC_ERR_ARG_WRONG,"Function only applies to underdetermined problems"); 138be5d1d56SKris Buschelman /* Create the null matrix using MATSEQDENSE explicitly */ 139878740d9SKris Buschelman ierr = MatCreate(A->comm,numRows,numNullCols,numRows,numNullCols,nullMat);CHKERRQ(ierr); 140878740d9SKris Buschelman ierr = MatSetType(*nullMat,MATSEQDENSE);CHKERRQ(ierr); 141878740d9SKris Buschelman ierr = MatSeqDenseSetPreallocation(*nullMat,PETSC_NULL);CHKERRQ(ierr); 142958c9bccSBarry Smith if (!numNullCols) { 14314b396bbSKris Buschelman ierr = MatAssemblyBegin(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14414b396bbSKris Buschelman ierr = MatAssemblyEnd(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14514b396bbSKris Buschelman PetscFunctionReturn(0); 14614b396bbSKris Buschelman } 14714b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 14814b396bbSKris Buschelman nullVals = (doublecomplex*)((Mat_SeqDense*)(*nullMat)->data)->v; 14914b396bbSKris Buschelman #else 15014b396bbSKris Buschelman nullVals = ((Mat_SeqDense*)(*nullMat)->data)->v; 15114b396bbSKris Buschelman #endif 15214b396bbSKris Buschelman /* Copy in the columns */ 15314b396bbSKris Buschelman Lstore = (SCformat*)lu->L.Store; 15414b396bbSKris Buschelman for(block = 0; block <= Lstore->nsuper; block++) { 15514b396bbSKris Buschelman newRow = Lstore->sup_to_col[block]; 15614b396bbSKris Buschelman size = Lstore->sup_to_col[block+1] - Lstore->sup_to_col[block]; 15714b396bbSKris Buschelman for(col = Lstore->rowind_colptr[newRow]; col < Lstore->rowind_colptr[newRow+1]; col++) { 15814b396bbSKris Buschelman newCol = Lstore->rowind[col]; 15914b396bbSKris Buschelman if (newCol >= numRows) { 16014b396bbSKris Buschelman for(b = 0; b < size; b++) 16114b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 16214b396bbSKris Buschelman nullVals[(newCol-numRows)*numRows+newRow+b] = ((doublecomplex*)Lstore->nzval)[Lstore->nzval_colptr[newRow+b]+col]; 16314b396bbSKris Buschelman #else 16414b396bbSKris Buschelman nullVals[(newCol-numRows)*numRows+newRow+b] = ((double*)Lstore->nzval)[Lstore->nzval_colptr[newRow+b]+col]; 16514b396bbSKris Buschelman #endif 16614b396bbSKris Buschelman } 16714b396bbSKris Buschelman } 16814b396bbSKris Buschelman } 16914b396bbSKris Buschelman /* Permute rhs to form P^T_c B */ 170*da7a1d00SHong Zhang ierr = PetscMalloc(numRows*sizeof(PetscReal),&workVals);CHKERRQ(ierr); 17114b396bbSKris Buschelman for(b = 0; b < numNullCols; b++) { 17214b396bbSKris Buschelman for(row = 0; row < numRows; row++) workVals[lu->perm_c[row]] = nullVals[b*numRows+row]; 17314b396bbSKris Buschelman for(row = 0; row < numRows; row++) nullVals[b*numRows+row] = workVals[row]; 17414b396bbSKris Buschelman } 17514b396bbSKris Buschelman /* Backward solve the upper triangle A x = b */ 17614b396bbSKris Buschelman for(b = 0; b < numNullCols; b++) { 17714b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 17875af56d4SHong Zhang sp_ztrsv("L","T","U",&lu->L,&lu->U,&nullVals[b*numRows],&stat,&ierr); 17914b396bbSKris Buschelman #else 18075af56d4SHong Zhang sp_dtrsv("L","T","U",&lu->L,&lu->U,&nullVals[b*numRows],&stat,&ierr); 18114b396bbSKris Buschelman #endif 18214b396bbSKris Buschelman if (ierr < 0) 18377431f27SBarry Smith SETERRQ1(PETSC_ERR_ARG_WRONG,"The argument %D was invalid",-ierr); 18414b396bbSKris Buschelman } 18514b396bbSKris Buschelman ierr = PetscFree(workVals);CHKERRQ(ierr); 18614b396bbSKris Buschelman 18714b396bbSKris Buschelman ierr = MatAssemblyBegin(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 18814b396bbSKris Buschelman ierr = MatAssemblyEnd(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 18914b396bbSKris Buschelman PetscFunctionReturn(0); 19014b396bbSKris Buschelman } 19175af56d4SHong Zhang #endif 19214b396bbSKris Buschelman 19314b396bbSKris Buschelman #undef __FUNCT__ 194f0c56d0fSKris Buschelman #define __FUNCT__ "MatSolve_SuperLU" 195dfbe8321SBarry Smith PetscErrorCode MatSolve_SuperLU(Mat A,Vec b,Vec x) 19614b396bbSKris Buschelman { 197f0c56d0fSKris Buschelman Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 19875af56d4SHong Zhang PetscScalar *barray,*xarray; 199dfbe8321SBarry Smith PetscErrorCode ierr; 200*da7a1d00SHong Zhang PetscInt info,i; 20175af56d4SHong Zhang SuperLUStat_t stat; 202*da7a1d00SHong Zhang PetscReal ferr,berr; 20314b396bbSKris Buschelman 20414b396bbSKris Buschelman PetscFunctionBegin; 205937a6b0eSHong Zhang if ( lu->lwork == -1 ) { 206937a6b0eSHong Zhang PetscFunctionReturn(0); 207937a6b0eSHong Zhang } 20875af56d4SHong Zhang lu->B.ncol = 1; /* Set the number of right-hand side */ 20975af56d4SHong Zhang ierr = VecGetArray(b,&barray);CHKERRQ(ierr); 21075af56d4SHong Zhang ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 2115fe6150dSHong Zhang 2125fe6150dSHong Zhang #if defined(PETSC_USE_COMPLEX) 2135fe6150dSHong Zhang ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray; 2145fe6150dSHong Zhang ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray; 2155fe6150dSHong Zhang #else 2165fe6150dSHong Zhang ((DNformat*)lu->B.Store)->nzval = barray; 21775af56d4SHong Zhang ((DNformat*)lu->X.Store)->nzval = xarray; 2185fe6150dSHong Zhang #endif 21975af56d4SHong Zhang 22075af56d4SHong Zhang /* Initialize the statistics variables. */ 22175af56d4SHong Zhang StatInit(&stat); 22275af56d4SHong Zhang 22375af56d4SHong Zhang lu->options.Fact = FACTORED; /* Indicate the factored form of A is supplied. */ 22475af56d4SHong Zhang lu->options.Trans = TRANS; 2258914a3f7SHong Zhang #if defined(PETSC_USE_COMPLEX) 2268914a3f7SHong Zhang zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 2278914a3f7SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 2288914a3f7SHong Zhang &lu->mem_usage, &stat, &info); 2298914a3f7SHong Zhang #else 23075af56d4SHong Zhang dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 23175af56d4SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 23275af56d4SHong Zhang &lu->mem_usage, &stat, &info); 2338914a3f7SHong Zhang #endif 23475af56d4SHong Zhang ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr); 23575af56d4SHong Zhang ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 23675af56d4SHong Zhang 237958c9bccSBarry Smith if ( !info || info == lu->A.ncol+1 ) { 23875af56d4SHong Zhang if ( lu->options.IterRefine ) { 2398914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n"); 2408914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR"); 24175af56d4SHong Zhang for (i = 0; i < 1; ++i) 2428914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," %8d%8d%16e%16e\n", i+1, stat.RefineSteps, ferr, berr); 24375af56d4SHong Zhang } 2448914a3f7SHong Zhang } else if ( info > 0 ){ 2458914a3f7SHong Zhang if ( lu->lwork == -1 ) { 24677431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %D bytes\n", info - lu->A.ncol); 2478914a3f7SHong Zhang } else { 24877431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," Warning: gssvx() returns info %D\n",info); 2498914a3f7SHong Zhang } 2508914a3f7SHong Zhang } else if (info < 0){ 25177431f27SBarry Smith SETERRQ2(PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", info,-info); 25214b396bbSKris Buschelman } 25314b396bbSKris Buschelman 2548914a3f7SHong Zhang if ( lu->options.PrintStat ) { 2558914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n"); 2568914a3f7SHong Zhang StatPrint(&stat); 2578914a3f7SHong Zhang } 25875af56d4SHong Zhang StatFree(&stat); 25975af56d4SHong Zhang PetscFunctionReturn(0); 26075af56d4SHong Zhang } 26114b396bbSKris Buschelman 26214b396bbSKris Buschelman #undef __FUNCT__ 263f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorNumeric_SuperLU" 264af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SuperLU(Mat A,MatFactorInfo *info,Mat *F) 26514b396bbSKris Buschelman { 26614b396bbSKris Buschelman Mat_SeqAIJ *aa = (Mat_SeqAIJ*)(A)->data; 267f0c56d0fSKris Buschelman Mat_SuperLU *lu = (Mat_SuperLU*)(*F)->spptr; 268dfbe8321SBarry Smith PetscErrorCode ierr; 269*da7a1d00SHong Zhang PetscInt sinfo; 27075af56d4SHong Zhang SuperLUStat_t stat; 271*da7a1d00SHong Zhang PetscReal ferr, berr; 2728914a3f7SHong Zhang NCformat *Ustore; 2738914a3f7SHong Zhang SCformat *Lstore; 27414b396bbSKris Buschelman 27514b396bbSKris Buschelman PetscFunctionBegin; 2769ce50919SHong Zhang if (lu->flg == SAME_NONZERO_PATTERN){ /* successing numerical factorization */ 2779ce50919SHong Zhang lu->options.Fact = SamePattern; 27875af56d4SHong Zhang /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */ 27975af56d4SHong Zhang Destroy_SuperMatrix_Store(&lu->A); 2809ce50919SHong Zhang if ( lu->lwork >= 0 ) { 28175af56d4SHong Zhang Destroy_SuperNode_Matrix(&lu->L); 28275af56d4SHong Zhang Destroy_CompCol_Matrix(&lu->U); 28375af56d4SHong Zhang lu->options.Fact = SamePattern; 28414b396bbSKris Buschelman } 28575af56d4SHong Zhang } 28614b396bbSKris Buschelman 28775af56d4SHong Zhang /* Create the SuperMatrix for lu->A=A^T: 28875af56d4SHong Zhang Since SuperLU likes column-oriented matrices,we pass it the transpose, 28975af56d4SHong Zhang and then solve A^T X = B in MatSolve(). */ 29075af56d4SHong Zhang #if defined(PETSC_USE_COMPLEX) 2915fe6150dSHong Zhang zCreate_CompCol_Matrix(&lu->A,A->n,A->m,aa->nz,(doublecomplex*)aa->a,aa->j,aa->i, 29275af56d4SHong Zhang SLU_NC,SLU_Z,SLU_GE); 29375af56d4SHong Zhang #else 29475af56d4SHong Zhang dCreate_CompCol_Matrix(&lu->A,A->n,A->m,aa->nz,aa->a,aa->j,aa->i, 29575af56d4SHong Zhang SLU_NC,SLU_D,SLU_GE); 29675af56d4SHong Zhang #endif 29775af56d4SHong Zhang 29875af56d4SHong Zhang /* Initialize the statistics variables. */ 29975af56d4SHong Zhang StatInit(&stat); 30075af56d4SHong Zhang 3019ce50919SHong Zhang /* Numerical factorization */ 30275af56d4SHong Zhang lu->B.ncol = 0; /* Indicate not to solve the system */ 3038914a3f7SHong Zhang #if defined(PETSC_USE_COMPLEX) 3048914a3f7SHong Zhang zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 3058914a3f7SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 306*da7a1d00SHong Zhang &lu->mem_usage, &stat, &sinfo); 3078914a3f7SHong Zhang #else 30875af56d4SHong Zhang dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 30975af56d4SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 310*da7a1d00SHong Zhang &lu->mem_usage, &stat, &sinfo); 3118914a3f7SHong Zhang #endif 312*da7a1d00SHong Zhang if ( !sinfo || sinfo == lu->A.ncol+1 ) { 3138914a3f7SHong Zhang if ( lu->options.PivotGrowth ) 3148914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," Recip. pivot growth = %e\n", lu->rpg); 31575af56d4SHong Zhang if ( lu->options.ConditionNumber ) 3168914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," Recip. condition number = %e\n", lu->rcond); 317*da7a1d00SHong Zhang } else if ( sinfo > 0 ){ 3188914a3f7SHong Zhang if ( lu->lwork == -1 ) { 319*da7a1d00SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %D bytes\n", sinfo - lu->A.ncol); 3208914a3f7SHong Zhang } else { 321*da7a1d00SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," Warning: gssvx() returns info %D\n",sinfo); 3228914a3f7SHong Zhang } 323*da7a1d00SHong Zhang } else { /* sinfo < 0 */ 324*da7a1d00SHong Zhang SETERRQ2(PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", sinfo,-sinfo); 32575af56d4SHong Zhang } 32675af56d4SHong Zhang 3278914a3f7SHong Zhang if ( lu->options.PrintStat ) { 3288914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"MatLUFactorNumeric_SuperLU():\n"); 3298914a3f7SHong Zhang StatPrint(&stat); 3308914a3f7SHong Zhang Lstore = (SCformat *) lu->L.Store; 3318914a3f7SHong Zhang Ustore = (NCformat *) lu->U.Store; 33277431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor L = %D\n", Lstore->nnz); 33377431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor U = %D\n", Ustore->nnz); 33477431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in L+U = %D\n", Lstore->nnz + Ustore->nnz - lu->A.ncol); 33577431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," L\\U MB %.3f\ttotal MB needed %.3f\texpansions %D\n", 3368914a3f7SHong Zhang lu->mem_usage.for_lu/1e6, lu->mem_usage.total_needed/1e6, 3378914a3f7SHong Zhang lu->mem_usage.expansions); 3388914a3f7SHong Zhang } 33975af56d4SHong Zhang StatFree(&stat); 34075af56d4SHong Zhang 34175af56d4SHong Zhang lu->flg = SAME_NONZERO_PATTERN; 34275af56d4SHong Zhang PetscFunctionReturn(0); 34314b396bbSKris Buschelman } 34414b396bbSKris Buschelman 34514b396bbSKris Buschelman /* 34614b396bbSKris Buschelman Note the r permutation is ignored 34714b396bbSKris Buschelman */ 34814b396bbSKris Buschelman #undef __FUNCT__ 349f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_SuperLU" 350dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) 35114b396bbSKris Buschelman { 35214b396bbSKris Buschelman Mat B; 353f0c56d0fSKris Buschelman Mat_SuperLU *lu; 3546849ba73SBarry Smith PetscErrorCode ierr; 355*da7a1d00SHong Zhang PetscInt m=A->m,n=A->n,indx; 3569ce50919SHong Zhang PetscTruth flg; 3575ca28756SHong Zhang const char *colperm[]={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */ 3585ca28756SHong Zhang const char *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"}; 3595ca28756SHong Zhang const char *rowperm[]={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */ 36014b396bbSKris Buschelman 36114b396bbSKris Buschelman PetscFunctionBegin; 362899d7b4fSKris Buschelman ierr = MatCreate(A->comm,A->m,A->n,PETSC_DETERMINE,PETSC_DETERMINE,&B);CHKERRQ(ierr); 363be5d1d56SKris Buschelman ierr = MatSetType(B,A->type_name);CHKERRQ(ierr); 364899d7b4fSKris Buschelman ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 365f0c56d0fSKris Buschelman 366f0c56d0fSKris Buschelman B->ops->lufactornumeric = MatLUFactorNumeric_SuperLU; 367f0c56d0fSKris Buschelman B->ops->solve = MatSolve_SuperLU; 36814b396bbSKris Buschelman B->factor = FACTOR_LU; 36994b7f48cSBarry Smith B->assembled = PETSC_TRUE; /* required by -ksp_view */ 37014b396bbSKris Buschelman 371f0c56d0fSKris Buschelman lu = (Mat_SuperLU*)(B->spptr); 3729ce50919SHong Zhang 3739ce50919SHong Zhang /* Set SuperLU options */ 3749ce50919SHong Zhang /* the default values for options argument: 3759ce50919SHong Zhang options.Fact = DOFACT; 3769ce50919SHong Zhang options.Equil = YES; 3779ce50919SHong Zhang options.ColPerm = COLAMD; 3789ce50919SHong Zhang options.DiagPivotThresh = 1.0; 3799ce50919SHong Zhang options.Trans = NOTRANS; 3809ce50919SHong Zhang options.IterRefine = NOREFINE; 3819ce50919SHong Zhang options.SymmetricMode = NO; 3829ce50919SHong Zhang options.PivotGrowth = NO; 3839ce50919SHong Zhang options.ConditionNumber = NO; 3849ce50919SHong Zhang options.PrintStat = YES; 3859ce50919SHong Zhang */ 3869ce50919SHong Zhang set_default_options(&lu->options); 3878914a3f7SHong Zhang /* equilibration causes error in solve(), thus not supported here. See dgssvx.c for possible reason. */ 3888914a3f7SHong Zhang lu->options.Equil = NO; 3899ce50919SHong Zhang lu->options.PrintStat = NO; 3908914a3f7SHong Zhang lu->lwork = 0; /* allocate space internally by system malloc */ 3919ce50919SHong Zhang 3929ce50919SHong Zhang ierr = PetscOptionsBegin(A->comm,A->prefix,"SuperLU Options","Mat");CHKERRQ(ierr); 3939ce50919SHong Zhang /* 3948914a3f7SHong Zhang ierr = PetscOptionsLogical("-mat_superlu_equil","Equil","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 3958914a3f7SHong Zhang if (flg) lu->options.Equil = YES; -- not supported by the interface !!! 3969ce50919SHong Zhang */ 3978914a3f7SHong Zhang ierr = PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);CHKERRQ(ierr); 3989ce50919SHong Zhang if (flg) {lu->options.ColPerm = (colperm_t)indx;} 3998914a3f7SHong Zhang ierr = PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);CHKERRQ(ierr); 4009ce50919SHong Zhang if (flg) { lu->options.IterRefine = (IterRefine_t)indx;} 4018914a3f7SHong Zhang ierr = PetscOptionsLogical("-mat_superlu_symmetricmode","SymmetricMode","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 4029ce50919SHong Zhang if (flg) lu->options.SymmetricMode = YES; 4038914a3f7SHong Zhang ierr = PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&lu->options.DiagPivotThresh,PETSC_NULL);CHKERRQ(ierr); 4048914a3f7SHong Zhang ierr = PetscOptionsLogical("-mat_superlu_pivotgrowth","PivotGrowth","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 4059ce50919SHong Zhang if (flg) lu->options.PivotGrowth = YES; 4068914a3f7SHong Zhang ierr = PetscOptionsLogical("-mat_superlu_conditionnumber","ConditionNumber","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 4079ce50919SHong Zhang if (flg) lu->options.ConditionNumber = YES; 4088914a3f7SHong Zhang ierr = PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[0],&indx,&flg);CHKERRQ(ierr); 4099ce50919SHong Zhang if (flg) {lu->options.RowPerm = (rowperm_t)indx;} 4108914a3f7SHong Zhang ierr = PetscOptionsLogical("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 4119ce50919SHong Zhang if (flg) lu->options.ReplaceTinyPivot = YES; 4128914a3f7SHong Zhang ierr = PetscOptionsLogical("-mat_superlu_printstat","PrintStat","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 4139ce50919SHong Zhang if (flg) lu->options.PrintStat = YES; 4148914a3f7SHong Zhang ierr = PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,PETSC_NULL);CHKERRQ(ierr); 4155fe6150dSHong Zhang if (lu->lwork > 0 ){ 4165fe6150dSHong Zhang ierr = PetscMalloc(lu->lwork,&lu->work);CHKERRQ(ierr); 4175fe6150dSHong Zhang } else if (lu->lwork != 0 && lu->lwork != -1){ 41877431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," Warning: lwork %D is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork); 4198914a3f7SHong Zhang lu->lwork = 0; 4208914a3f7SHong Zhang } 4219ce50919SHong Zhang PetscOptionsEnd(); 4229ce50919SHong Zhang 42375af56d4SHong Zhang #ifdef SUPERLU2 4242ecb5974SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatCreateNull","MatCreateNull_SuperLU", 425f0c56d0fSKris Buschelman (void(*)(void))MatCreateNull_SuperLU);CHKERRQ(ierr); 42675af56d4SHong Zhang #endif 42714b396bbSKris Buschelman 42875af56d4SHong Zhang /* Allocate spaces (notice sizes are for the transpose) */ 429*da7a1d00SHong Zhang ierr = PetscMalloc(m*sizeof(PetscInt),&lu->etree);CHKERRQ(ierr); 430*da7a1d00SHong Zhang ierr = PetscMalloc(n*sizeof(PetscInt),&lu->perm_r);CHKERRQ(ierr); 431*da7a1d00SHong Zhang ierr = PetscMalloc(m*sizeof(PetscInt),&lu->perm_c);CHKERRQ(ierr); 432*da7a1d00SHong Zhang ierr = PetscMalloc(n*sizeof(PetscInt),&lu->R);CHKERRQ(ierr); 433*da7a1d00SHong Zhang ierr = PetscMalloc(m*sizeof(PetscInt),&lu->C);CHKERRQ(ierr); 43475af56d4SHong Zhang 4359ce50919SHong Zhang /* create rhs and solution x without allocate space for .Store */ 4365fe6150dSHong Zhang #if defined(PETSC_USE_COMPLEX) 437937a6b0eSHong Zhang zCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE); 438937a6b0eSHong Zhang zCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE); 4395fe6150dSHong Zhang #else 44075af56d4SHong Zhang dCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE); 44175af56d4SHong Zhang dCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE); 4425fe6150dSHong Zhang #endif 44314b396bbSKris Buschelman 44414b396bbSKris Buschelman lu->flg = DIFFERENT_NONZERO_PATTERN; 445e740cb95SKris Buschelman lu->CleanUpSuperLU = PETSC_TRUE; 446899d7b4fSKris Buschelman 447899d7b4fSKris Buschelman *F = B; 448*da7a1d00SHong Zhang ierr = PetscLogObjectMemory(B,(A->m+A->n)*sizeof(PetscInt)+sizeof(Mat_SuperLU));CHKERRQ(ierr); 44914b396bbSKris Buschelman PetscFunctionReturn(0); 45014b396bbSKris Buschelman } 45114b396bbSKris Buschelman 45294b7f48cSBarry Smith /* used by -ksp_view */ 45314b396bbSKris Buschelman #undef __FUNCT__ 454f0c56d0fSKris Buschelman #define __FUNCT__ "MatFactorInfo_SuperLU" 455dfbe8321SBarry Smith PetscErrorCode MatFactorInfo_SuperLU(Mat A,PetscViewer viewer) 45614b396bbSKris Buschelman { 457f0c56d0fSKris Buschelman Mat_SuperLU *lu= (Mat_SuperLU*)A->spptr; 458*da7a1d00SHong Zhang PetscErrorCode ierr; 4599ce50919SHong Zhang superlu_options_t options; 460f0c56d0fSKris Buschelman 461f0c56d0fSKris Buschelman PetscFunctionBegin; 4629ce50919SHong Zhang /* check if matrix is superlu_dist type */ 4639ce50919SHong Zhang if (A->ops->solve != MatSolve_SuperLU) PetscFunctionReturn(0); 4649ce50919SHong Zhang 4659ce50919SHong Zhang options = lu->options; 466f0c56d0fSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");CHKERRQ(ierr); 4679ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Equil: %s\n",(options.Equil != NO) ? "YES": "NO");CHKERRQ(ierr); 46877431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ColPerm: %D\n",options.ColPerm);CHKERRQ(ierr); 46977431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," IterRefine: %D\n",options.IterRefine);CHKERRQ(ierr); 4709ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," SymmetricMode: %s\n",(options.SymmetricMode != NO) ? "YES": "NO");CHKERRQ(ierr); 4719ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," DiagPivotThresh: %g\n",options.DiagPivotThresh);CHKERRQ(ierr); 4729ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," PivotGrowth: %s\n",(options.PivotGrowth != NO) ? "YES": "NO");CHKERRQ(ierr); 4739ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ConditionNumber: %s\n",(options.ConditionNumber != NO) ? "YES": "NO");CHKERRQ(ierr); 47477431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," RowPerm: %D\n",options.RowPerm);CHKERRQ(ierr); 4759ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ReplaceTinyPivot: %s\n",(options.ReplaceTinyPivot != NO) ? "YES": "NO");CHKERRQ(ierr); 4769ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," PrintStat: %s\n",(options.PrintStat != NO) ? "YES": "NO");CHKERRQ(ierr); 47777431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," lwork: %D\n",lu->lwork);CHKERRQ(ierr); 4789ce50919SHong Zhang 479f0c56d0fSKris Buschelman PetscFunctionReturn(0); 480f0c56d0fSKris Buschelman } 481f0c56d0fSKris Buschelman 482f0c56d0fSKris Buschelman #undef __FUNCT__ 483f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_SuperLU" 484dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SuperLU(Mat A, MatDuplicateOption op, Mat *M) { 485*da7a1d00SHong Zhang PetscErrorCode ierr; 4868f340917SKris Buschelman Mat_SuperLU *lu=(Mat_SuperLU *)A->spptr; 4878f340917SKris Buschelman 48814b396bbSKris Buschelman PetscFunctionBegin; 4898f340917SKris Buschelman ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr); 4903f953163SKris Buschelman ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_SuperLU));CHKERRQ(ierr); 49114b396bbSKris Buschelman PetscFunctionReturn(0); 49214b396bbSKris Buschelman } 49314b396bbSKris Buschelman 49414b396bbSKris Buschelman EXTERN_C_BEGIN 49514b396bbSKris Buschelman #undef __FUNCT__ 496b0592601SKris Buschelman #define __FUNCT__ "MatConvert_SuperLU_SeqAIJ" 497ceb03754SKris Buschelman PetscErrorCode MatConvert_SuperLU_SeqAIJ(Mat A,const MatType type,MatReuse reuse,Mat *newmat) 498521d7252SBarry Smith { 499fb3e25aaSKris Buschelman /* This routine is only called to convert an unfactored PETSc-SuperLU matrix */ 500fb3e25aaSKris Buschelman /* to its base PETSc type, so we will ignore 'MatType type'. */ 501*da7a1d00SHong Zhang PetscErrorCode ierr; 502b0592601SKris Buschelman Mat B=*newmat; 503f0c56d0fSKris Buschelman Mat_SuperLU *lu=(Mat_SuperLU *)A->spptr; 504b0592601SKris Buschelman 505b0592601SKris Buschelman PetscFunctionBegin; 506ceb03754SKris Buschelman if (reuse == MAT_INITIAL_MATRIX) { 507b0592601SKris Buschelman ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 508f0c56d0fSKris Buschelman } 509b0592601SKris Buschelman /* Reset the original function pointers */ 510f0c56d0fSKris Buschelman B->ops->duplicate = lu->MatDuplicate; 511b0592601SKris Buschelman B->ops->view = lu->MatView; 512b0592601SKris Buschelman B->ops->assemblyend = lu->MatAssemblyEnd; 513b0592601SKris Buschelman B->ops->lufactorsymbolic = lu->MatLUFactorSymbolic; 514b0592601SKris Buschelman B->ops->destroy = lu->MatDestroy; 515b0592601SKris Buschelman /* lu is only a function pointer stash unless we've factored the matrix, which we haven't! */ 516b0592601SKris Buschelman ierr = PetscFree(lu);CHKERRQ(ierr); 517901853e0SKris Buschelman 518901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_superlu_C","",PETSC_NULL);CHKERRQ(ierr); 519901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_superlu_seqaij_C","",PETSC_NULL);CHKERRQ(ierr); 520901853e0SKris Buschelman 521b0592601SKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);CHKERRQ(ierr); 522b0592601SKris Buschelman *newmat = B; 523b0592601SKris Buschelman PetscFunctionReturn(0); 524b0592601SKris Buschelman } 525b0592601SKris Buschelman EXTERN_C_END 526b0592601SKris Buschelman 527b0592601SKris Buschelman EXTERN_C_BEGIN 528b0592601SKris Buschelman #undef __FUNCT__ 529b0592601SKris Buschelman #define __FUNCT__ "MatConvert_SeqAIJ_SuperLU" 530ceb03754SKris Buschelman PetscErrorCode MatConvert_SeqAIJ_SuperLU(Mat A,const MatType type,MatReuse reuse,Mat *newmat) 531521d7252SBarry Smith { 532fb3e25aaSKris Buschelman /* This routine is only called to convert to MATSUPERLU */ 533fb3e25aaSKris Buschelman /* from MATSEQAIJ, so we will ignore 'MatType type'. */ 534*da7a1d00SHong Zhang PetscErrorCode ierr; 535b0592601SKris Buschelman Mat B=*newmat; 536f0c56d0fSKris Buschelman Mat_SuperLU *lu; 53714b396bbSKris Buschelman 53814b396bbSKris Buschelman PetscFunctionBegin; 539ceb03754SKris Buschelman if (reuse == MAT_INITIAL_MATRIX) { 540b0592601SKris Buschelman ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 541b0592601SKris Buschelman } 54214b396bbSKris Buschelman 543f0c56d0fSKris Buschelman ierr = PetscNew(Mat_SuperLU,&lu);CHKERRQ(ierr); 544f0c56d0fSKris Buschelman lu->MatDuplicate = A->ops->duplicate; 54514b396bbSKris Buschelman lu->MatView = A->ops->view; 54614b396bbSKris Buschelman lu->MatAssemblyEnd = A->ops->assemblyend; 547b0592601SKris Buschelman lu->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 54814b396bbSKris Buschelman lu->MatDestroy = A->ops->destroy; 54914b396bbSKris Buschelman lu->CleanUpSuperLU = PETSC_FALSE; 550b0592601SKris Buschelman 551b0592601SKris Buschelman B->spptr = (void*)lu; 552f0c56d0fSKris Buschelman B->ops->duplicate = MatDuplicate_SuperLU; 553f0c56d0fSKris Buschelman B->ops->view = MatView_SuperLU; 554f0c56d0fSKris Buschelman B->ops->assemblyend = MatAssemblyEnd_SuperLU; 555f0c56d0fSKris Buschelman B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU; 556e4e36384SHong Zhang B->ops->choleskyfactorsymbolic = 0; 557f0c56d0fSKris Buschelman B->ops->destroy = MatDestroy_SuperLU; 558b0592601SKris Buschelman 559b0592601SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_superlu_C", 560b0592601SKris Buschelman "MatConvert_SeqAIJ_SuperLU",MatConvert_SeqAIJ_SuperLU);CHKERRQ(ierr); 561b0592601SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_superlu_seqaij_C", 562b0592601SKris Buschelman "MatConvert_SuperLU_SeqAIJ",MatConvert_SuperLU_SeqAIJ);CHKERRQ(ierr); 56352e6d16bSBarry Smith PetscLogInfo(0,"MatConvert_SeqAIJ_SuperLU:Using SuperLU for SeqAIJ LU factorization and solves."); 564fb3e25aaSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)B,MATSUPERLU);CHKERRQ(ierr); 565b0592601SKris Buschelman *newmat = B; 566b0592601SKris Buschelman PetscFunctionReturn(0); 567b0592601SKris Buschelman } 568b0592601SKris Buschelman EXTERN_C_END 569b0592601SKris Buschelman 57024b6179bSKris Buschelman /*MC 571fafad747SKris Buschelman MATSUPERLU - MATSUPERLU = "superlu" - A matrix type providing direct solvers (LU) for sequential matrices 57224b6179bSKris Buschelman via the external package SuperLU. 57324b6179bSKris Buschelman 57424b6179bSKris Buschelman If SuperLU is installed (see the manual for 57524b6179bSKris Buschelman instructions on how to declare the existence of external packages), 57624b6179bSKris Buschelman a matrix type can be constructed which invokes SuperLU solvers. 57724b6179bSKris Buschelman After calling MatCreate(...,A), simply call MatSetType(A,MATSUPERLU). 57824b6179bSKris Buschelman 57924b6179bSKris Buschelman This matrix inherits from MATSEQAIJ. As a result, MatSeqAIJSetPreallocation is 58028b08bd3SKris Buschelman supported for this matrix type. One can also call MatConvert for an inplace conversion to or from 58128b08bd3SKris Buschelman the MATSEQAIJ type without data copy. 58224b6179bSKris Buschelman 58324b6179bSKris Buschelman Options Database Keys: 5840bad9183SKris Buschelman + -mat_type superlu - sets the matrix type to "superlu" during a call to MatSetFromOptions() 58524b6179bSKris Buschelman - -mat_superlu_ordering <0,1,2,3> - 0: natural ordering, 58624b6179bSKris Buschelman 1: MMD applied to A'*A, 58724b6179bSKris Buschelman 2: MMD applied to A'+A, 58824b6179bSKris Buschelman 3: COLAMD, approximate minimum degree column ordering 58924b6179bSKris Buschelman 59024b6179bSKris Buschelman Level: beginner 59124b6179bSKris Buschelman 59224b6179bSKris Buschelman .seealso: PCLU 59324b6179bSKris Buschelman M*/ 59424b6179bSKris Buschelman 595b0592601SKris Buschelman EXTERN_C_BEGIN 596b0592601SKris Buschelman #undef __FUNCT__ 597f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_SuperLU" 598dfbe8321SBarry Smith PetscErrorCode MatCreate_SuperLU(Mat A) 599dfbe8321SBarry Smith { 600dfbe8321SBarry Smith PetscErrorCode ierr; 601b0592601SKris Buschelman 602b0592601SKris Buschelman PetscFunctionBegin; 6035441df8eSKris Buschelman /* Change type name before calling MatSetType to force proper construction of SeqAIJ and SUPERLU types */ 6045441df8eSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)A,MATSUPERLU);CHKERRQ(ierr); 605b0592601SKris Buschelman ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr); 606ceb03754SKris Buschelman ierr = MatConvert_SeqAIJ_SuperLU(A,MATSUPERLU,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr); 60714b396bbSKris Buschelman PetscFunctionReturn(0); 60814b396bbSKris Buschelman } 60914b396bbSKris Buschelman EXTERN_C_END 610