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; 2075af56d4SHong Zhang int *perm_c; /* column permutation vector */ 2175af56d4SHong Zhang int *perm_r; /* row permutations from partial pivoting */ 2275af56d4SHong Zhang int *etree; 2375af56d4SHong Zhang double *R, *C; 2475af56d4SHong Zhang char equed[1]; 2575af56d4SHong Zhang int lwork; 2675af56d4SHong Zhang void *work; 2775af56d4SHong Zhang double rpg, rcond; 2875af56d4SHong Zhang mem_usage_t mem_usage; 2975af56d4SHong Zhang MatStructure flg; 3014b396bbSKris Buschelman 3114b396bbSKris Buschelman /* A few function pointers for inheritance */ 32f0c56d0fSKris Buschelman int (*MatDuplicate)(Mat,MatDuplicateOption,Mat*); 3314b396bbSKris Buschelman int (*MatView)(Mat,PetscViewer); 3414b396bbSKris Buschelman int (*MatAssemblyEnd)(Mat,MatAssemblyType); 35b0592601SKris Buschelman int (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*); 3614b396bbSKris Buschelman int (*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 43f0c56d0fSKris Buschelman EXTERN int MatFactorInfo_SuperLU(Mat,PetscViewer); 44f0c56d0fSKris Buschelman EXTERN int MatLUFactorSymbolic_SuperLU(Mat,IS,IS,MatFactorInfo*,Mat*); 45b0592601SKris Buschelman 46b0592601SKris Buschelman EXTERN_C_BEGIN 475ca28756SHong Zhang EXTERN int MatConvert_SuperLU_SeqAIJ(Mat,const MatType,Mat*); 485ca28756SHong Zhang EXTERN int MatConvert_SeqAIJ_SuperLU(Mat,const MatType,Mat*); 49b0592601SKris Buschelman EXTERN_C_END 5014b396bbSKris Buschelman 5114b396bbSKris Buschelman #undef __FUNCT__ 52f0c56d0fSKris Buschelman #define __FUNCT__ "MatDestroy_SuperLU" 53f0c56d0fSKris Buschelman int MatDestroy_SuperLU(Mat A) 5414b396bbSKris Buschelman { 55460c4903SKris Buschelman int 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 } 74b0592601SKris Buschelman ierr = MatConvert_SuperLU_SeqAIJ(A,MATSEQAIJ,&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" 81f0c56d0fSKris Buschelman int MatView_SuperLU(Mat A,PetscViewer viewer) 8214b396bbSKris Buschelman { 8314b396bbSKris Buschelman int 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" 103f0c56d0fSKris Buschelman int MatAssemblyEnd_SuperLU(Mat A,MatAssemblyType mode) { 10414b396bbSKris Buschelman int 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 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" 120f0c56d0fSKris Buschelman int MatCreateNull_SuperLU(Mat A,Mat *nullMat) 12114b396bbSKris Buschelman { 122f0c56d0fSKris Buschelman Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 12314b396bbSKris Buschelman int numRows = A->m,numCols = A->n; 12414b396bbSKris Buschelman SCformat *Lstore; 12514b396bbSKris Buschelman int 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 13214b396bbSKris Buschelman int row,newRow,col,newCol,block,b,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); 14214b396bbSKris Buschelman if (numNullCols == 0) { 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 */ 17014b396bbSKris Buschelman ierr = PetscMalloc(numRows*sizeof(double),&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) 18314b396bbSKris Buschelman 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" 195f0c56d0fSKris Buschelman int MatSolve_SuperLU(Mat A,Vec b,Vec x) 19614b396bbSKris Buschelman { 197f0c56d0fSKris Buschelman Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 19875af56d4SHong Zhang PetscScalar *barray,*xarray; 1998914a3f7SHong Zhang int ierr,info,i; 20075af56d4SHong Zhang SuperLUStat_t stat; 2018914a3f7SHong Zhang double ferr,berr; 20214b396bbSKris Buschelman 20314b396bbSKris Buschelman PetscFunctionBegin; 204937a6b0eSHong Zhang if ( lu->lwork == -1 ) { 205937a6b0eSHong Zhang PetscFunctionReturn(0); 206937a6b0eSHong Zhang } 20775af56d4SHong Zhang lu->B.ncol = 1; /* Set the number of right-hand side */ 20875af56d4SHong Zhang ierr = VecGetArray(b,&barray);CHKERRQ(ierr); 20975af56d4SHong Zhang ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 2105fe6150dSHong Zhang 2115fe6150dSHong Zhang #if defined(PETSC_USE_COMPLEX) 2125fe6150dSHong Zhang ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray; 2135fe6150dSHong Zhang ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray; 2145fe6150dSHong Zhang #else 2155fe6150dSHong Zhang ((DNformat*)lu->B.Store)->nzval = barray; 21675af56d4SHong Zhang ((DNformat*)lu->X.Store)->nzval = xarray; 2175fe6150dSHong Zhang #endif 21875af56d4SHong Zhang 21975af56d4SHong Zhang /* Initialize the statistics variables. */ 22075af56d4SHong Zhang StatInit(&stat); 22175af56d4SHong Zhang 22275af56d4SHong Zhang lu->options.Fact = FACTORED; /* Indicate the factored form of A is supplied. */ 22375af56d4SHong Zhang lu->options.Trans = TRANS; 2248914a3f7SHong Zhang #if defined(PETSC_USE_COMPLEX) 2258914a3f7SHong Zhang zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 2268914a3f7SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 2278914a3f7SHong Zhang &lu->mem_usage, &stat, &info); 2288914a3f7SHong Zhang #else 22975af56d4SHong Zhang dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 23075af56d4SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 23175af56d4SHong Zhang &lu->mem_usage, &stat, &info); 2328914a3f7SHong Zhang #endif 23375af56d4SHong Zhang ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr); 23475af56d4SHong Zhang ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 23575af56d4SHong Zhang 23675af56d4SHong Zhang if ( info == 0 || info == lu->A.ncol+1 ) { 23775af56d4SHong Zhang if ( lu->options.IterRefine ) { 2388914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n"); 2398914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR"); 24075af56d4SHong Zhang for (i = 0; i < 1; ++i) 2418914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," %8d%8d%16e%16e\n", i+1, stat.RefineSteps, ferr, berr); 24275af56d4SHong Zhang } 2438914a3f7SHong Zhang } else if ( info > 0 ){ 2448914a3f7SHong Zhang if ( lu->lwork == -1 ) { 2458914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %d bytes\n", info - lu->A.ncol); 2468914a3f7SHong Zhang } else { 2478914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," Warning: gssvx() returns info %d\n",info); 2488914a3f7SHong Zhang } 2498914a3f7SHong Zhang } else if (info < 0){ 2508914a3f7SHong Zhang SETERRQ2(1, "info = %d, the %d-th argument in gssvx() had an illegal value", info,-info); 25114b396bbSKris Buschelman } 25214b396bbSKris Buschelman 2538914a3f7SHong Zhang if ( lu->options.PrintStat ) { 2548914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n"); 2558914a3f7SHong Zhang StatPrint(&stat); 2568914a3f7SHong Zhang } 25775af56d4SHong Zhang StatFree(&stat); 25875af56d4SHong Zhang PetscFunctionReturn(0); 25975af56d4SHong Zhang } 26014b396bbSKris Buschelman 26114b396bbSKris Buschelman #undef __FUNCT__ 262f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorNumeric_SuperLU" 263f0c56d0fSKris Buschelman int MatLUFactorNumeric_SuperLU(Mat A,Mat *F) 26414b396bbSKris Buschelman { 26514b396bbSKris Buschelman Mat_SeqAIJ *aa = (Mat_SeqAIJ*)(A)->data; 266f0c56d0fSKris Buschelman Mat_SuperLU *lu = (Mat_SuperLU*)(*F)->spptr; 26775af56d4SHong Zhang int ierr,info; 26875af56d4SHong Zhang SuperLUStat_t stat; 26975af56d4SHong Zhang double ferr, berr; 2708914a3f7SHong Zhang NCformat *Ustore; 2718914a3f7SHong Zhang SCformat *Lstore; 27214b396bbSKris Buschelman 27314b396bbSKris Buschelman PetscFunctionBegin; 2749ce50919SHong Zhang if (lu->flg == SAME_NONZERO_PATTERN){ /* successing numerical factorization */ 2759ce50919SHong Zhang lu->options.Fact = SamePattern; 27675af56d4SHong Zhang /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */ 27775af56d4SHong Zhang Destroy_SuperMatrix_Store(&lu->A); 2789ce50919SHong Zhang if ( lu->lwork >= 0 ) { 27975af56d4SHong Zhang Destroy_SuperNode_Matrix(&lu->L); 28075af56d4SHong Zhang Destroy_CompCol_Matrix(&lu->U); 28175af56d4SHong Zhang lu->options.Fact = SamePattern; 28214b396bbSKris Buschelman } 28375af56d4SHong Zhang } 28414b396bbSKris Buschelman 28575af56d4SHong Zhang /* Create the SuperMatrix for lu->A=A^T: 28675af56d4SHong Zhang Since SuperLU likes column-oriented matrices,we pass it the transpose, 28775af56d4SHong Zhang and then solve A^T X = B in MatSolve(). */ 28875af56d4SHong Zhang #if defined(PETSC_USE_COMPLEX) 2895fe6150dSHong Zhang zCreate_CompCol_Matrix(&lu->A,A->n,A->m,aa->nz,(doublecomplex*)aa->a,aa->j,aa->i, 29075af56d4SHong Zhang SLU_NC,SLU_Z,SLU_GE); 29175af56d4SHong Zhang #else 29275af56d4SHong Zhang dCreate_CompCol_Matrix(&lu->A,A->n,A->m,aa->nz,aa->a,aa->j,aa->i, 29375af56d4SHong Zhang SLU_NC,SLU_D,SLU_GE); 29475af56d4SHong Zhang #endif 29575af56d4SHong Zhang 29675af56d4SHong Zhang /* Initialize the statistics variables. */ 29775af56d4SHong Zhang StatInit(&stat); 29875af56d4SHong Zhang 2999ce50919SHong Zhang /* Numerical factorization */ 30075af56d4SHong Zhang lu->B.ncol = 0; /* Indicate not to solve the system */ 3018914a3f7SHong Zhang #if defined(PETSC_USE_COMPLEX) 3028914a3f7SHong Zhang zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 3038914a3f7SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 3048914a3f7SHong Zhang &lu->mem_usage, &stat, &info); 3058914a3f7SHong Zhang #else 30675af56d4SHong Zhang dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 30775af56d4SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 30875af56d4SHong Zhang &lu->mem_usage, &stat, &info); 3098914a3f7SHong Zhang #endif 31075af56d4SHong Zhang if ( info == 0 || info == lu->A.ncol+1 ) { 3118914a3f7SHong Zhang if ( lu->options.PivotGrowth ) 3128914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," Recip. pivot growth = %e\n", lu->rpg); 31375af56d4SHong Zhang if ( lu->options.ConditionNumber ) 3148914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," Recip. condition number = %e\n", lu->rcond); 3158914a3f7SHong Zhang } else if ( info > 0 ){ 3168914a3f7SHong Zhang if ( lu->lwork == -1 ) { 3178914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %d bytes\n", info - lu->A.ncol); 3188914a3f7SHong Zhang } else { 3198914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," Warning: gssvx() returns info %d\n",info); 3208914a3f7SHong Zhang } 321937a6b0eSHong Zhang } else { /* info < 0 */ 3228914a3f7SHong Zhang SETERRQ2(1, "info = %d, the %d-th argument in gssvx() had an illegal value", info,-info); 32375af56d4SHong Zhang } 32475af56d4SHong Zhang 3258914a3f7SHong Zhang if ( lu->options.PrintStat ) { 3268914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"MatLUFactorNumeric_SuperLU():\n"); 3278914a3f7SHong Zhang StatPrint(&stat); 3288914a3f7SHong Zhang Lstore = (SCformat *) lu->L.Store; 3298914a3f7SHong Zhang Ustore = (NCformat *) lu->U.Store; 3308914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor L = %d\n", Lstore->nnz); 3318914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor U = %d\n", Ustore->nnz); 3328914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in L+U = %d\n", Lstore->nnz + Ustore->nnz - lu->A.ncol); 3338914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," L\\U MB %.3f\ttotal MB needed %.3f\texpansions %d\n", 3348914a3f7SHong Zhang lu->mem_usage.for_lu/1e6, lu->mem_usage.total_needed/1e6, 3358914a3f7SHong Zhang lu->mem_usage.expansions); 3368914a3f7SHong Zhang } 33775af56d4SHong Zhang StatFree(&stat); 33875af56d4SHong Zhang 33975af56d4SHong Zhang lu->flg = SAME_NONZERO_PATTERN; 34075af56d4SHong Zhang PetscFunctionReturn(0); 34114b396bbSKris Buschelman } 34214b396bbSKris Buschelman 34314b396bbSKris Buschelman /* 34414b396bbSKris Buschelman Note the r permutation is ignored 34514b396bbSKris Buschelman */ 34614b396bbSKris Buschelman #undef __FUNCT__ 347f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_SuperLU" 348f0c56d0fSKris Buschelman int MatLUFactorSymbolic_SuperLU(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) 34914b396bbSKris Buschelman { 35014b396bbSKris Buschelman Mat B; 351f0c56d0fSKris Buschelman Mat_SuperLU *lu; 3529ce50919SHong Zhang int ierr,m=A->m,n=A->n,indx; 3539ce50919SHong Zhang PetscTruth flg; 3545ca28756SHong Zhang const char *colperm[]={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */ 3555ca28756SHong Zhang const char *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"}; 3565ca28756SHong Zhang const char *rowperm[]={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */ 35714b396bbSKris Buschelman 35814b396bbSKris Buschelman PetscFunctionBegin; 35914b396bbSKris Buschelman 360899d7b4fSKris Buschelman ierr = MatCreate(A->comm,A->m,A->n,PETSC_DETERMINE,PETSC_DETERMINE,&B);CHKERRQ(ierr); 361be5d1d56SKris Buschelman ierr = MatSetType(B,A->type_name);CHKERRQ(ierr); 362899d7b4fSKris Buschelman ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 363f0c56d0fSKris Buschelman 364f0c56d0fSKris Buschelman B->ops->lufactornumeric = MatLUFactorNumeric_SuperLU; 365f0c56d0fSKris Buschelman B->ops->solve = MatSolve_SuperLU; 36614b396bbSKris Buschelman B->factor = FACTOR_LU; 36794b7f48cSBarry Smith B->assembled = PETSC_TRUE; /* required by -ksp_view */ 36814b396bbSKris Buschelman 369f0c56d0fSKris Buschelman lu = (Mat_SuperLU*)(B->spptr); 3709ce50919SHong Zhang 3719ce50919SHong Zhang /* Set SuperLU options */ 3729ce50919SHong Zhang /* the default values for options argument: 3739ce50919SHong Zhang options.Fact = DOFACT; 3749ce50919SHong Zhang options.Equil = YES; 3759ce50919SHong Zhang options.ColPerm = COLAMD; 3769ce50919SHong Zhang options.DiagPivotThresh = 1.0; 3779ce50919SHong Zhang options.Trans = NOTRANS; 3789ce50919SHong Zhang options.IterRefine = NOREFINE; 3799ce50919SHong Zhang options.SymmetricMode = NO; 3809ce50919SHong Zhang options.PivotGrowth = NO; 3819ce50919SHong Zhang options.ConditionNumber = NO; 3829ce50919SHong Zhang options.PrintStat = YES; 3839ce50919SHong Zhang */ 3849ce50919SHong Zhang set_default_options(&lu->options); 3858914a3f7SHong Zhang /* equilibration causes error in solve(), thus not supported here. See dgssvx.c for possible reason. */ 3868914a3f7SHong Zhang lu->options.Equil = NO; 3879ce50919SHong Zhang lu->options.PrintStat = NO; 3888914a3f7SHong Zhang lu->lwork = 0; /* allocate space internally by system malloc */ 3899ce50919SHong Zhang 3909ce50919SHong Zhang ierr = PetscOptionsBegin(A->comm,A->prefix,"SuperLU Options","Mat");CHKERRQ(ierr); 3919ce50919SHong Zhang /* 3928914a3f7SHong Zhang ierr = PetscOptionsLogical("-mat_superlu_equil","Equil","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 3938914a3f7SHong Zhang if (flg) lu->options.Equil = YES; -- not supported by the interface !!! 3949ce50919SHong Zhang */ 3958914a3f7SHong Zhang ierr = PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);CHKERRQ(ierr); 3969ce50919SHong Zhang if (flg) {lu->options.ColPerm = (colperm_t)indx;} 3978914a3f7SHong Zhang ierr = PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);CHKERRQ(ierr); 3989ce50919SHong Zhang if (flg) { lu->options.IterRefine = (IterRefine_t)indx;} 3998914a3f7SHong Zhang ierr = PetscOptionsLogical("-mat_superlu_symmetricmode","SymmetricMode","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 4009ce50919SHong Zhang if (flg) lu->options.SymmetricMode = YES; 4018914a3f7SHong Zhang ierr = PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&lu->options.DiagPivotThresh,PETSC_NULL);CHKERRQ(ierr); 4028914a3f7SHong Zhang ierr = PetscOptionsLogical("-mat_superlu_pivotgrowth","PivotGrowth","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 4039ce50919SHong Zhang if (flg) lu->options.PivotGrowth = YES; 4048914a3f7SHong Zhang ierr = PetscOptionsLogical("-mat_superlu_conditionnumber","ConditionNumber","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 4059ce50919SHong Zhang if (flg) lu->options.ConditionNumber = YES; 4068914a3f7SHong Zhang ierr = PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[0],&indx,&flg);CHKERRQ(ierr); 4079ce50919SHong Zhang if (flg) {lu->options.RowPerm = (rowperm_t)indx;} 4088914a3f7SHong Zhang ierr = PetscOptionsLogical("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 4099ce50919SHong Zhang if (flg) lu->options.ReplaceTinyPivot = YES; 4108914a3f7SHong Zhang ierr = PetscOptionsLogical("-mat_superlu_printstat","PrintStat","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 4119ce50919SHong Zhang if (flg) lu->options.PrintStat = YES; 4128914a3f7SHong Zhang ierr = PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,PETSC_NULL);CHKERRQ(ierr); 4135fe6150dSHong Zhang if (lu->lwork > 0 ){ 4145fe6150dSHong Zhang ierr = PetscMalloc(lu->lwork,&lu->work);CHKERRQ(ierr); 4155fe6150dSHong Zhang } else if (lu->lwork != 0 && lu->lwork != -1){ 416937a6b0eSHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," Warning: lwork %d is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork); 4178914a3f7SHong Zhang lu->lwork = 0; 4188914a3f7SHong Zhang } 4199ce50919SHong Zhang PetscOptionsEnd(); 4209ce50919SHong Zhang 42175af56d4SHong Zhang #ifdef SUPERLU2 4222ecb5974SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatCreateNull","MatCreateNull_SuperLU", 423f0c56d0fSKris Buschelman (void(*)(void))MatCreateNull_SuperLU);CHKERRQ(ierr); 42475af56d4SHong Zhang #endif 42514b396bbSKris Buschelman 42675af56d4SHong Zhang /* Allocate spaces (notice sizes are for the transpose) */ 4279ce50919SHong Zhang ierr = PetscMalloc(m*sizeof(int),&lu->etree);CHKERRQ(ierr); 4289ce50919SHong Zhang ierr = PetscMalloc(n*sizeof(int),&lu->perm_r);CHKERRQ(ierr); 4299ce50919SHong Zhang ierr = PetscMalloc(m*sizeof(int),&lu->perm_c);CHKERRQ(ierr); 4309ce50919SHong Zhang ierr = PetscMalloc(n*sizeof(int),&lu->R);CHKERRQ(ierr); 4319ce50919SHong Zhang ierr = PetscMalloc(m*sizeof(int),&lu->C);CHKERRQ(ierr); 43275af56d4SHong Zhang 4339ce50919SHong Zhang /* create rhs and solution x without allocate space for .Store */ 4345fe6150dSHong Zhang #if defined(PETSC_USE_COMPLEX) 435937a6b0eSHong Zhang zCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE); 436937a6b0eSHong Zhang zCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE); 4375fe6150dSHong Zhang #else 43875af56d4SHong Zhang dCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE); 43975af56d4SHong Zhang dCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE); 4405fe6150dSHong Zhang #endif 44114b396bbSKris Buschelman 44214b396bbSKris Buschelman lu->flg = DIFFERENT_NONZERO_PATTERN; 443e740cb95SKris Buschelman lu->CleanUpSuperLU = PETSC_TRUE; 444899d7b4fSKris Buschelman 445899d7b4fSKris Buschelman *F = B; 4469ce50919SHong Zhang PetscLogObjectMemory(B,(A->m+A->n)*sizeof(int)+sizeof(Mat_SuperLU)); 44714b396bbSKris Buschelman PetscFunctionReturn(0); 44814b396bbSKris Buschelman } 44914b396bbSKris Buschelman 45094b7f48cSBarry Smith /* used by -ksp_view */ 45114b396bbSKris Buschelman #undef __FUNCT__ 452f0c56d0fSKris Buschelman #define __FUNCT__ "MatFactorInfo_SuperLU" 453f0c56d0fSKris Buschelman int MatFactorInfo_SuperLU(Mat A,PetscViewer viewer) 45414b396bbSKris Buschelman { 455f0c56d0fSKris Buschelman Mat_SuperLU *lu= (Mat_SuperLU*)A->spptr; 456f0c56d0fSKris Buschelman int ierr; 4579ce50919SHong Zhang superlu_options_t options; 458f0c56d0fSKris Buschelman 459f0c56d0fSKris Buschelman PetscFunctionBegin; 4609ce50919SHong Zhang /* check if matrix is superlu_dist type */ 4619ce50919SHong Zhang if (A->ops->solve != MatSolve_SuperLU) PetscFunctionReturn(0); 4629ce50919SHong Zhang 4639ce50919SHong Zhang options = lu->options; 464f0c56d0fSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");CHKERRQ(ierr); 4659ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Equil: %s\n",(options.Equil != NO) ? "YES": "NO");CHKERRQ(ierr); 4669ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ColPerm: %d\n",options.ColPerm);CHKERRQ(ierr); 4679ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," IterRefine: %d\n",options.IterRefine);CHKERRQ(ierr); 4689ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," SymmetricMode: %s\n",(options.SymmetricMode != NO) ? "YES": "NO");CHKERRQ(ierr); 4699ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," DiagPivotThresh: %g\n",options.DiagPivotThresh);CHKERRQ(ierr); 4709ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," PivotGrowth: %s\n",(options.PivotGrowth != NO) ? "YES": "NO");CHKERRQ(ierr); 4719ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ConditionNumber: %s\n",(options.ConditionNumber != NO) ? "YES": "NO");CHKERRQ(ierr); 4729ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RowPerm: %d\n",options.RowPerm);CHKERRQ(ierr); 4739ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ReplaceTinyPivot: %s\n",(options.ReplaceTinyPivot != NO) ? "YES": "NO");CHKERRQ(ierr); 4749ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," PrintStat: %s\n",(options.PrintStat != NO) ? "YES": "NO");CHKERRQ(ierr); 4758914a3f7SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," lwork: %d\n",lu->lwork);CHKERRQ(ierr); 4769ce50919SHong Zhang 477f0c56d0fSKris Buschelman PetscFunctionReturn(0); 478f0c56d0fSKris Buschelman } 479f0c56d0fSKris Buschelman 480f0c56d0fSKris Buschelman #undef __FUNCT__ 481f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_SuperLU" 482f0c56d0fSKris Buschelman int MatDuplicate_SuperLU(Mat A, MatDuplicateOption op, Mat *M) { 48314b396bbSKris Buschelman int ierr; 4848f340917SKris Buschelman Mat_SuperLU *lu=(Mat_SuperLU *)A->spptr; 4858f340917SKris Buschelman 48614b396bbSKris Buschelman PetscFunctionBegin; 4878f340917SKris Buschelman ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr); 4883f953163SKris Buschelman ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_SuperLU));CHKERRQ(ierr); 48914b396bbSKris Buschelman PetscFunctionReturn(0); 49014b396bbSKris Buschelman } 49114b396bbSKris Buschelman 49214b396bbSKris Buschelman EXTERN_C_BEGIN 49314b396bbSKris Buschelman #undef __FUNCT__ 494b0592601SKris Buschelman #define __FUNCT__ "MatConvert_SuperLU_SeqAIJ" 4958e9aea5cSBarry Smith int MatConvert_SuperLU_SeqAIJ(Mat A,const MatType type,Mat *newmat) { 496fb3e25aaSKris Buschelman /* This routine is only called to convert an unfactored PETSc-SuperLU matrix */ 497fb3e25aaSKris Buschelman /* to its base PETSc type, so we will ignore 'MatType type'. */ 49814b396bbSKris Buschelman int ierr; 499b0592601SKris Buschelman Mat B=*newmat; 500f0c56d0fSKris Buschelman Mat_SuperLU *lu=(Mat_SuperLU *)A->spptr; 501b0592601SKris Buschelman 502b0592601SKris Buschelman PetscFunctionBegin; 503b0592601SKris Buschelman if (B != A) { 504b0592601SKris Buschelman ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 505f0c56d0fSKris Buschelman } 506b0592601SKris Buschelman /* Reset the original function pointers */ 507f0c56d0fSKris Buschelman B->ops->duplicate = lu->MatDuplicate; 508b0592601SKris Buschelman B->ops->view = lu->MatView; 509b0592601SKris Buschelman B->ops->assemblyend = lu->MatAssemblyEnd; 510b0592601SKris Buschelman B->ops->lufactorsymbolic = lu->MatLUFactorSymbolic; 511b0592601SKris Buschelman B->ops->destroy = lu->MatDestroy; 512b0592601SKris Buschelman /* lu is only a function pointer stash unless we've factored the matrix, which we haven't! */ 513b0592601SKris Buschelman ierr = PetscFree(lu);CHKERRQ(ierr); 514b0592601SKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);CHKERRQ(ierr); 515b0592601SKris Buschelman *newmat = B; 516b0592601SKris Buschelman PetscFunctionReturn(0); 517b0592601SKris Buschelman } 518b0592601SKris Buschelman EXTERN_C_END 519b0592601SKris Buschelman 520b0592601SKris Buschelman EXTERN_C_BEGIN 521b0592601SKris Buschelman #undef __FUNCT__ 522b0592601SKris Buschelman #define __FUNCT__ "MatConvert_SeqAIJ_SuperLU" 5238e9aea5cSBarry Smith int MatConvert_SeqAIJ_SuperLU(Mat A,const MatType type,Mat *newmat) { 524fb3e25aaSKris Buschelman /* This routine is only called to convert to MATSUPERLU */ 525fb3e25aaSKris Buschelman /* from MATSEQAIJ, so we will ignore 'MatType type'. */ 526b0592601SKris Buschelman int ierr; 527b0592601SKris Buschelman Mat B=*newmat; 528f0c56d0fSKris Buschelman Mat_SuperLU *lu; 52914b396bbSKris Buschelman 53014b396bbSKris Buschelman PetscFunctionBegin; 531b0592601SKris Buschelman if (B != A) { 532b0592601SKris Buschelman ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 533b0592601SKris Buschelman } 53414b396bbSKris Buschelman 535f0c56d0fSKris Buschelman ierr = PetscNew(Mat_SuperLU,&lu);CHKERRQ(ierr); 536f0c56d0fSKris Buschelman lu->MatDuplicate = A->ops->duplicate; 53714b396bbSKris Buschelman lu->MatView = A->ops->view; 53814b396bbSKris Buschelman lu->MatAssemblyEnd = A->ops->assemblyend; 539b0592601SKris Buschelman lu->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 54014b396bbSKris Buschelman lu->MatDestroy = A->ops->destroy; 54114b396bbSKris Buschelman lu->CleanUpSuperLU = PETSC_FALSE; 542b0592601SKris Buschelman 543b0592601SKris Buschelman B->spptr = (void*)lu; 544f0c56d0fSKris Buschelman B->ops->duplicate = MatDuplicate_SuperLU; 545f0c56d0fSKris Buschelman B->ops->view = MatView_SuperLU; 546f0c56d0fSKris Buschelman B->ops->assemblyend = MatAssemblyEnd_SuperLU; 547f0c56d0fSKris Buschelman B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU; 548*e4e36384SHong Zhang B->ops->choleskyfactorsymbolic = 0; 549f0c56d0fSKris Buschelman B->ops->destroy = MatDestroy_SuperLU; 550b0592601SKris Buschelman 551b0592601SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_superlu_C", 552b0592601SKris Buschelman "MatConvert_SeqAIJ_SuperLU",MatConvert_SeqAIJ_SuperLU);CHKERRQ(ierr); 553b0592601SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_superlu_seqaij_C", 554b0592601SKris Buschelman "MatConvert_SuperLU_SeqAIJ",MatConvert_SuperLU_SeqAIJ);CHKERRQ(ierr); 555fb3e25aaSKris Buschelman PetscLogInfo(0,"Using SuperLU for SeqAIJ LU factorization and solves."); 556fb3e25aaSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)B,MATSUPERLU);CHKERRQ(ierr); 557b0592601SKris Buschelman *newmat = B; 558b0592601SKris Buschelman PetscFunctionReturn(0); 559b0592601SKris Buschelman } 560b0592601SKris Buschelman EXTERN_C_END 561b0592601SKris Buschelman 56224b6179bSKris Buschelman /*MC 563fafad747SKris Buschelman MATSUPERLU - MATSUPERLU = "superlu" - A matrix type providing direct solvers (LU) for sequential matrices 56424b6179bSKris Buschelman via the external package SuperLU. 56524b6179bSKris Buschelman 56624b6179bSKris Buschelman If SuperLU is installed (see the manual for 56724b6179bSKris Buschelman instructions on how to declare the existence of external packages), 56824b6179bSKris Buschelman a matrix type can be constructed which invokes SuperLU solvers. 56924b6179bSKris Buschelman After calling MatCreate(...,A), simply call MatSetType(A,MATSUPERLU). 57024b6179bSKris Buschelman This matrix type is only supported for double precision real. 57124b6179bSKris Buschelman 57224b6179bSKris Buschelman This matrix inherits from MATSEQAIJ. As a result, MatSeqAIJSetPreallocation is 57328b08bd3SKris Buschelman supported for this matrix type. One can also call MatConvert for an inplace conversion to or from 57428b08bd3SKris Buschelman the MATSEQAIJ type without data copy. 57524b6179bSKris Buschelman 57624b6179bSKris Buschelman Options Database Keys: 5770bad9183SKris Buschelman + -mat_type superlu - sets the matrix type to "superlu" during a call to MatSetFromOptions() 57824b6179bSKris Buschelman - -mat_superlu_ordering <0,1,2,3> - 0: natural ordering, 57924b6179bSKris Buschelman 1: MMD applied to A'*A, 58024b6179bSKris Buschelman 2: MMD applied to A'+A, 58124b6179bSKris Buschelman 3: COLAMD, approximate minimum degree column ordering 58224b6179bSKris Buschelman 58324b6179bSKris Buschelman Level: beginner 58424b6179bSKris Buschelman 58524b6179bSKris Buschelman .seealso: PCLU 58624b6179bSKris Buschelman M*/ 58724b6179bSKris Buschelman 588b0592601SKris Buschelman EXTERN_C_BEGIN 589b0592601SKris Buschelman #undef __FUNCT__ 590f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_SuperLU" 591f0c56d0fSKris Buschelman int MatCreate_SuperLU(Mat A) { 592b0592601SKris Buschelman int ierr; 593b0592601SKris Buschelman 594b0592601SKris Buschelman PetscFunctionBegin; 5955441df8eSKris Buschelman /* Change type name before calling MatSetType to force proper construction of SeqAIJ and SUPERLU types */ 5965441df8eSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)A,MATSUPERLU);CHKERRQ(ierr); 597b0592601SKris Buschelman ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr); 598b0592601SKris Buschelman ierr = MatConvert_SeqAIJ_SuperLU(A,MATSUPERLU,&A);CHKERRQ(ierr); 59914b396bbSKris Buschelman PetscFunctionReturn(0); 60014b396bbSKris Buschelman } 60114b396bbSKris Buschelman EXTERN_C_END 602