114b396bbSKris Buschelman /*$Id: superlu.c,v 1.10 2001/08/15 15:56:50 bsmith Exp $*/ 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; 2175af56d4SHong Zhang int *perm_c; /* column permutation vector */ 2275af56d4SHong Zhang int *perm_r; /* row permutations from partial pivoting */ 2375af56d4SHong Zhang int *etree; 2475af56d4SHong Zhang double *R, *C; 2575af56d4SHong Zhang char equed[1]; 2675af56d4SHong Zhang int lwork; 2775af56d4SHong Zhang void *work; 2875af56d4SHong Zhang double rpg, rcond; 2975af56d4SHong Zhang mem_usage_t mem_usage; 3075af56d4SHong Zhang MatStructure flg; 3175af56d4SHong Zhang /* 3214b396bbSKris Buschelman SuperMatrix A,B,AC,L,U; 3314b396bbSKris Buschelman int *perm_r,*perm_c,ispec,relax,panel_size; 3414b396bbSKris Buschelman double pivot_threshold; 3514b396bbSKris Buschelman NCformat *store; 3614b396bbSKris Buschelman MatStructure flg; 3714b396bbSKris Buschelman PetscTruth SuperluMatOdering; 3875af56d4SHong Zhang */ 3914b396bbSKris Buschelman 4014b396bbSKris Buschelman /* A few function pointers for inheritance */ 41f0c56d0fSKris Buschelman int (*MatDuplicate)(Mat,MatDuplicateOption,Mat*); 4214b396bbSKris Buschelman int (*MatView)(Mat,PetscViewer); 4314b396bbSKris Buschelman int (*MatAssemblyEnd)(Mat,MatAssemblyType); 44b0592601SKris Buschelman int (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*); 4514b396bbSKris Buschelman int (*MatDestroy)(Mat); 4614b396bbSKris Buschelman 4714b396bbSKris Buschelman /* Flag to clean up (non-global) SuperLU objects during Destroy */ 4814b396bbSKris Buschelman PetscTruth CleanUpSuperLU; 49f0c56d0fSKris Buschelman } Mat_SuperLU; 5014b396bbSKris Buschelman 5114b396bbSKris Buschelman 52f0c56d0fSKris Buschelman EXTERN int MatFactorInfo_SuperLU(Mat,PetscViewer); 53f0c56d0fSKris Buschelman EXTERN int MatLUFactorSymbolic_SuperLU(Mat,IS,IS,MatFactorInfo*,Mat*); 54b0592601SKris Buschelman 55b0592601SKris Buschelman EXTERN_C_BEGIN 56b0592601SKris Buschelman EXTERN int MatConvert_SuperLU_SeqAIJ(Mat,MatType,Mat*); 57b0592601SKris Buschelman EXTERN int MatConvert_SeqAIJ_SuperLU(Mat,MatType,Mat*); 58b0592601SKris Buschelman EXTERN_C_END 5914b396bbSKris Buschelman 6014b396bbSKris Buschelman #undef __FUNCT__ 61f0c56d0fSKris Buschelman #define __FUNCT__ "MatDestroy_SuperLU" 62f0c56d0fSKris Buschelman int MatDestroy_SuperLU(Mat A) 6314b396bbSKris Buschelman { 64460c4903SKris Buschelman int ierr; 65f0c56d0fSKris Buschelman Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 6614b396bbSKris Buschelman 6714b396bbSKris Buschelman PetscFunctionBegin; 6875af56d4SHong Zhang if (lu->CleanUpSuperLU) { /* Free the SuperLU datastructures */ 6975af56d4SHong Zhang /* Destroy_CompCol_Matrix(&lu->A); */ /* hangs inside memory.c! */ 7075af56d4SHong Zhang Destroy_SuperMatrix_Store(&lu->A); 7175af56d4SHong Zhang Destroy_SuperMatrix_Store(&lu->B); 7275af56d4SHong Zhang Destroy_SuperMatrix_Store(&lu->X); 73*9ce50919SHong Zhang 74*9ce50919SHong Zhang ierr = PetscFree(lu->etree);CHKERRQ(ierr); 75*9ce50919SHong Zhang ierr = PetscFree(lu->perm_r);CHKERRQ(ierr); 76*9ce50919SHong Zhang ierr = PetscFree(lu->perm_c);CHKERRQ(ierr); 77*9ce50919SHong Zhang ierr = PetscFree(lu->R);CHKERRQ(ierr); 78*9ce50919SHong Zhang ierr = PetscFree(lu->C);CHKERRQ(ierr); 7975af56d4SHong Zhang if ( lu->lwork >= 0 ) { 80fb3e25aaSKris Buschelman Destroy_SuperNode_Matrix(&lu->L); 81fb3e25aaSKris Buschelman Destroy_CompCol_Matrix(&lu->U); 8275af56d4SHong Zhang } 83fb3e25aaSKris Buschelman } 84b0592601SKris Buschelman ierr = MatConvert_SuperLU_SeqAIJ(A,MATSEQAIJ,&A);CHKERRQ(ierr); 85fb3e25aaSKris Buschelman ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 8614b396bbSKris Buschelman PetscFunctionReturn(0); 8714b396bbSKris Buschelman } 8814b396bbSKris Buschelman 8914b396bbSKris Buschelman #undef __FUNCT__ 90f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_SuperLU" 91f0c56d0fSKris Buschelman int MatView_SuperLU(Mat A,PetscViewer viewer) 9214b396bbSKris Buschelman { 9314b396bbSKris Buschelman int ierr; 9414b396bbSKris Buschelman PetscTruth isascii; 9514b396bbSKris Buschelman PetscViewerFormat format; 96f0c56d0fSKris Buschelman Mat_SuperLU *lu=(Mat_SuperLU*)(A->spptr); 9714b396bbSKris Buschelman 9814b396bbSKris Buschelman PetscFunctionBegin; 9914b396bbSKris Buschelman ierr = (*lu->MatView)(A,viewer);CHKERRQ(ierr); 10014b396bbSKris Buschelman 10114b396bbSKris Buschelman ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 10214b396bbSKris Buschelman if (isascii) { 10314b396bbSKris Buschelman ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 10414b396bbSKris Buschelman if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) { 105f0c56d0fSKris Buschelman ierr = MatFactorInfo_SuperLU(A,viewer);CHKERRQ(ierr); 10614b396bbSKris Buschelman } 10714b396bbSKris Buschelman } 10814b396bbSKris Buschelman PetscFunctionReturn(0); 10914b396bbSKris Buschelman } 11014b396bbSKris Buschelman 11114b396bbSKris Buschelman #undef __FUNCT__ 112f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_SuperLU" 113f0c56d0fSKris Buschelman int MatAssemblyEnd_SuperLU(Mat A,MatAssemblyType mode) { 11414b396bbSKris Buschelman int ierr; 115f0c56d0fSKris Buschelman Mat_SuperLU *lu=(Mat_SuperLU*)(A->spptr); 11614b396bbSKris Buschelman 11714b396bbSKris Buschelman PetscFunctionBegin; 11814b396bbSKris Buschelman ierr = (*lu->MatAssemblyEnd)(A,mode);CHKERRQ(ierr); 119b0592601SKris Buschelman 120b0592601SKris Buschelman lu->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 121f0c56d0fSKris Buschelman A->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU; 12214b396bbSKris Buschelman PetscFunctionReturn(0); 12314b396bbSKris Buschelman } 12414b396bbSKris Buschelman 12575af56d4SHong Zhang /* This function was written for SuperLU 2.0 by Matthew Knepley. Not tested for SuperLU 3.0! */ 12675af56d4SHong Zhang #ifdef SuperLU2 12714b396bbSKris Buschelman #include "src/mat/impls/dense/seq/dense.h" 12814b396bbSKris Buschelman #undef __FUNCT__ 129f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreateNull_SuperLU" 130f0c56d0fSKris Buschelman int MatCreateNull_SuperLU(Mat A,Mat *nullMat) 13114b396bbSKris Buschelman { 132f0c56d0fSKris Buschelman Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 13314b396bbSKris Buschelman int numRows = A->m,numCols = A->n; 13414b396bbSKris Buschelman SCformat *Lstore; 13514b396bbSKris Buschelman int numNullCols,size; 13675af56d4SHong Zhang SuperLUStat_t stat; 13714b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 13814b396bbSKris Buschelman doublecomplex *nullVals,*workVals; 13914b396bbSKris Buschelman #else 14014b396bbSKris Buschelman PetscScalar *nullVals,*workVals; 14114b396bbSKris Buschelman #endif 14214b396bbSKris Buschelman int row,newRow,col,newCol,block,b,ierr; 14314b396bbSKris Buschelman 14414b396bbSKris Buschelman PetscFunctionBegin; 14514b396bbSKris Buschelman PetscValidHeaderSpecific(A,MAT_COOKIE); 14614b396bbSKris Buschelman PetscValidPointer(nullMat); 14714b396bbSKris Buschelman if (!A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix"); 14814b396bbSKris Buschelman numNullCols = numCols - numRows; 14914b396bbSKris Buschelman if (numNullCols < 0) SETERRQ(PETSC_ERR_ARG_WRONG,"Function only applies to underdetermined problems"); 15014b396bbSKris Buschelman /* Create the null matrix */ 15114b396bbSKris Buschelman ierr = MatCreateSeqDense(A->comm,numRows,numNullCols,PETSC_NULL,nullMat);CHKERRQ(ierr); 15214b396bbSKris Buschelman if (numNullCols == 0) { 15314b396bbSKris Buschelman ierr = MatAssemblyBegin(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 15414b396bbSKris Buschelman ierr = MatAssemblyEnd(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 15514b396bbSKris Buschelman PetscFunctionReturn(0); 15614b396bbSKris Buschelman } 15714b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 15814b396bbSKris Buschelman nullVals = (doublecomplex*)((Mat_SeqDense*)(*nullMat)->data)->v; 15914b396bbSKris Buschelman #else 16014b396bbSKris Buschelman nullVals = ((Mat_SeqDense*)(*nullMat)->data)->v; 16114b396bbSKris Buschelman #endif 16214b396bbSKris Buschelman /* Copy in the columns */ 16314b396bbSKris Buschelman Lstore = (SCformat*)lu->L.Store; 16414b396bbSKris Buschelman for(block = 0; block <= Lstore->nsuper; block++) { 16514b396bbSKris Buschelman newRow = Lstore->sup_to_col[block]; 16614b396bbSKris Buschelman size = Lstore->sup_to_col[block+1] - Lstore->sup_to_col[block]; 16714b396bbSKris Buschelman for(col = Lstore->rowind_colptr[newRow]; col < Lstore->rowind_colptr[newRow+1]; col++) { 16814b396bbSKris Buschelman newCol = Lstore->rowind[col]; 16914b396bbSKris Buschelman if (newCol >= numRows) { 17014b396bbSKris Buschelman for(b = 0; b < size; b++) 17114b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 17214b396bbSKris Buschelman nullVals[(newCol-numRows)*numRows+newRow+b] = ((doublecomplex*)Lstore->nzval)[Lstore->nzval_colptr[newRow+b]+col]; 17314b396bbSKris Buschelman #else 17414b396bbSKris Buschelman nullVals[(newCol-numRows)*numRows+newRow+b] = ((double*)Lstore->nzval)[Lstore->nzval_colptr[newRow+b]+col]; 17514b396bbSKris Buschelman #endif 17614b396bbSKris Buschelman } 17714b396bbSKris Buschelman } 17814b396bbSKris Buschelman } 17914b396bbSKris Buschelman /* Permute rhs to form P^T_c B */ 18014b396bbSKris Buschelman ierr = PetscMalloc(numRows*sizeof(double),&workVals);CHKERRQ(ierr); 18114b396bbSKris Buschelman for(b = 0; b < numNullCols; b++) { 18214b396bbSKris Buschelman for(row = 0; row < numRows; row++) workVals[lu->perm_c[row]] = nullVals[b*numRows+row]; 18314b396bbSKris Buschelman for(row = 0; row < numRows; row++) nullVals[b*numRows+row] = workVals[row]; 18414b396bbSKris Buschelman } 18514b396bbSKris Buschelman /* Backward solve the upper triangle A x = b */ 18614b396bbSKris Buschelman for(b = 0; b < numNullCols; b++) { 18714b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 18875af56d4SHong Zhang sp_ztrsv("L","T","U",&lu->L,&lu->U,&nullVals[b*numRows],&stat,&ierr); 18914b396bbSKris Buschelman #else 19075af56d4SHong Zhang sp_dtrsv("L","T","U",&lu->L,&lu->U,&nullVals[b*numRows],&stat,&ierr); 19114b396bbSKris Buschelman #endif 19214b396bbSKris Buschelman if (ierr < 0) 19314b396bbSKris Buschelman SETERRQ1(PETSC_ERR_ARG_WRONG,"The argument %d was invalid",-ierr); 19414b396bbSKris Buschelman } 19514b396bbSKris Buschelman ierr = PetscFree(workVals);CHKERRQ(ierr); 19614b396bbSKris Buschelman 19714b396bbSKris Buschelman ierr = MatAssemblyBegin(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 19814b396bbSKris Buschelman ierr = MatAssemblyEnd(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 19914b396bbSKris Buschelman PetscFunctionReturn(0); 20014b396bbSKris Buschelman } 20175af56d4SHong Zhang #endif 20214b396bbSKris Buschelman 20314b396bbSKris Buschelman #undef __FUNCT__ 204f0c56d0fSKris Buschelman #define __FUNCT__ "MatSolve_SuperLU" 205f0c56d0fSKris Buschelman int MatSolve_SuperLU(Mat A,Vec b,Vec x) 20614b396bbSKris Buschelman { 207f0c56d0fSKris Buschelman Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 20875af56d4SHong Zhang PetscScalar *barray,*xarray; 20975af56d4SHong Zhang int m,n,ierr,lwork,info,i; 21075af56d4SHong Zhang SuperLUStat_t stat; 21175af56d4SHong Zhang double ferr,berr,*rhsb,*rhsx; 21214b396bbSKris Buschelman 21314b396bbSKris Buschelman PetscFunctionBegin; 21475af56d4SHong Zhang /* rhs vector */ 21575af56d4SHong Zhang lu->B.ncol = 1; /* Set the number of right-hand side */ 21675af56d4SHong Zhang ierr = VecGetArray(b,&barray);CHKERRQ(ierr); 21775af56d4SHong Zhang ((DNformat*)lu->B.Store)->nzval = barray; 21875af56d4SHong Zhang 21975af56d4SHong Zhang /* solution vector */ 22075af56d4SHong Zhang ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 22175af56d4SHong Zhang ((DNformat*)lu->X.Store)->nzval = xarray; 22275af56d4SHong Zhang 22375af56d4SHong Zhang /* Initialize the statistics variables. */ 22475af56d4SHong Zhang StatInit(&stat); 22575af56d4SHong Zhang 22675af56d4SHong Zhang lu->options.Fact = FACTORED; /* Indicate the factored form of A is supplied. */ 22775af56d4SHong Zhang lu->options.Trans = TRANS; 22875af56d4SHong Zhang dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 22975af56d4SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 23075af56d4SHong Zhang &lu->mem_usage, &stat, &info); 23175af56d4SHong Zhang 23275af56d4SHong Zhang ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr); 23375af56d4SHong Zhang ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 23475af56d4SHong Zhang 23575af56d4SHong Zhang if ( info == 0 || info == lu->A.ncol+1 ) { 23675af56d4SHong Zhang if ( lu->options.IterRefine ) { 23775af56d4SHong Zhang printf("Iterative Refinement:\n"); 23875af56d4SHong Zhang printf("%8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR"); 23975af56d4SHong Zhang for (i = 0; i < 1; ++i) 24075af56d4SHong Zhang printf("%8d%8d%16e%16e\n", i+1, stat.RefineSteps, ferr, berr); 24175af56d4SHong Zhang } 24275af56d4SHong Zhang fflush(stdout); 24375af56d4SHong Zhang } else if ( info > 0 && lu->lwork == -1 ) { 24475af56d4SHong Zhang printf("** Estimated memory: %d bytes\n", info - n); 24514b396bbSKris Buschelman } 24614b396bbSKris Buschelman 24775af56d4SHong Zhang if ( lu->options.PrintStat ) StatPrint(&stat); 24875af56d4SHong Zhang StatFree(&stat); 24975af56d4SHong Zhang PetscFunctionReturn(0); 25075af56d4SHong Zhang } 25114b396bbSKris Buschelman 25214b396bbSKris Buschelman #undef __FUNCT__ 253f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorNumeric_SuperLU" 254f0c56d0fSKris Buschelman int MatLUFactorNumeric_SuperLU(Mat A,Mat *F) 25514b396bbSKris Buschelman { 25614b396bbSKris Buschelman Mat_SeqAIJ *aa = (Mat_SeqAIJ*)(A)->data; 257f0c56d0fSKris Buschelman Mat_SuperLU *lu = (Mat_SuperLU*)(*F)->spptr; 25875af56d4SHong Zhang int ierr,info; 25914b396bbSKris Buschelman PetscTruth flag; 26075af56d4SHong Zhang SuperLUStat_t stat; 26175af56d4SHong Zhang double ferr, berr; 26214b396bbSKris Buschelman 26314b396bbSKris Buschelman PetscFunctionBegin; 264*9ce50919SHong Zhang if (lu->flg == SAME_NONZERO_PATTERN){ /* successing numerical factorization */ 265*9ce50919SHong Zhang lu->options.Fact = SamePattern; 26675af56d4SHong Zhang /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */ 26775af56d4SHong Zhang Destroy_SuperMatrix_Store(&lu->A); 268*9ce50919SHong Zhang if ( lu->lwork >= 0 ) { 26975af56d4SHong Zhang Destroy_SuperNode_Matrix(&lu->L); 27075af56d4SHong Zhang Destroy_CompCol_Matrix(&lu->U); 27175af56d4SHong Zhang lu->options.Fact = SamePattern; 27214b396bbSKris Buschelman } 27375af56d4SHong Zhang } 27414b396bbSKris Buschelman 27575af56d4SHong Zhang /* Create the SuperMatrix for lu->A=A^T: 27675af56d4SHong Zhang Since SuperLU likes column-oriented matrices,we pass it the transpose, 27775af56d4SHong Zhang and then solve A^T X = B in MatSolve(). */ 27875af56d4SHong Zhang #if defined(PETSC_USE_COMPLEX) 27975af56d4SHong Zhang zCreate_CompCol_Matrix(&lu->A,A->n,A->m,aa->nz,aa->a,aa->j,aa->i, 28075af56d4SHong Zhang SLU_NC,SLU_Z,SLU_GE); 28175af56d4SHong Zhang #else 28275af56d4SHong Zhang dCreate_CompCol_Matrix(&lu->A,A->n,A->m,aa->nz,aa->a,aa->j,aa->i, 28375af56d4SHong Zhang SLU_NC,SLU_D,SLU_GE); 28475af56d4SHong Zhang #endif 28575af56d4SHong Zhang 28675af56d4SHong Zhang /* Initialize the statistics variables. */ 28775af56d4SHong Zhang StatInit(&stat); 28875af56d4SHong Zhang 289*9ce50919SHong Zhang /* Numerical factorization */ 29075af56d4SHong Zhang lu->lwork = 0; /* allocate space internally by system malloc */ 29175af56d4SHong Zhang lu->B.ncol = 0; /* Indicate not to solve the system */ 29275af56d4SHong Zhang dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 29375af56d4SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 29475af56d4SHong Zhang &lu->mem_usage, &stat, &info); 29575af56d4SHong Zhang 29675af56d4SHong Zhang if ( info == 0 || info == lu->A.ncol+1 ) { 29775af56d4SHong Zhang if ( lu->options.PivotGrowth ) printf("Recip. pivot growth = %e\n", lu->rpg); 29875af56d4SHong Zhang if ( lu->options.ConditionNumber ) 29975af56d4SHong Zhang printf("Recip. condition number = %e\n", lu->rcond); 30075af56d4SHong Zhang /* 30175af56d4SHong Zhang NCformat *Ustore; 30275af56d4SHong Zhang SCformat *Lstore; 30375af56d4SHong Zhang Lstore = (SCformat *) lu->L.Store; 30475af56d4SHong Zhang Ustore = (NCformat *) lu->U.Store; 30575af56d4SHong Zhang printf("No of nonzeros in factor L = %d\n", Lstore->nnz); 30675af56d4SHong Zhang printf("No of nonzeros in factor U = %d\n", Ustore->nnz); 30775af56d4SHong Zhang printf("No of nonzeros in L+U = %d\n", Lstore->nnz + Ustore->nnz - n); 30875af56d4SHong Zhang printf("L\\U MB %.3f\ttotal MB needed %.3f\texpansions %d\n", 30975af56d4SHong Zhang mem_usage.for_lu/1e6, mem_usage.total_needed/1e6, 31075af56d4SHong Zhang mem_usage.expansions); 31175af56d4SHong Zhang */ 31275af56d4SHong Zhang fflush(stdout); 31375af56d4SHong Zhang 31475af56d4SHong Zhang } else if ( info > 0 && lu->lwork == -1 ) { 31575af56d4SHong Zhang printf("** Estimated memory: %d bytes\n", info - lu->A.ncol); 31675af56d4SHong Zhang } 31775af56d4SHong Zhang 31875af56d4SHong Zhang if ( lu->options.PrintStat ) StatPrint(&stat); 31975af56d4SHong Zhang StatFree(&stat); 32075af56d4SHong Zhang 32175af56d4SHong Zhang lu->flg = SAME_NONZERO_PATTERN; 32275af56d4SHong Zhang PetscFunctionReturn(0); 32375af56d4SHong Zhang 32475af56d4SHong Zhang #ifdef OLD 32514b396bbSKris Buschelman /* Set SuperLU options */ 32614b396bbSKris Buschelman lu->relax = sp_ienv(2); 32714b396bbSKris Buschelman lu->panel_size = sp_ienv(1); 32814b396bbSKris Buschelman /* We have to initialize global data or SuperLU crashes (sucky design) */ 32914b396bbSKris Buschelman if (!StatInitCalled) { 33014b396bbSKris Buschelman StatInit(lu->panel_size,lu->relax); 33114b396bbSKris Buschelman } 33214b396bbSKris Buschelman StatInitCalled++; 33314b396bbSKris Buschelman 33414b396bbSKris Buschelman ierr = PetscOptionsBegin(A->comm,A->prefix,"SuperLU Options","Mat");CHKERRQ(ierr); 33514b396bbSKris Buschelman /* use SuperLU mat ordeing */ 33614b396bbSKris Buschelman ierr = PetscOptionsInt("-mat_superlu_ordering","SuperLU ordering type (one of 0, 1, 2, 3)\n 0: natural ordering;\n 1: MMD applied to A'*A;\n 2: MMD applied to A'+A;\n 3: COLAMD, approximate minimum degree column ordering","None",lu->ispec,&lu->ispec,&flag);CHKERRQ(ierr); 33714b396bbSKris Buschelman if (flag) { 33814b396bbSKris Buschelman get_perm_c(lu->ispec, &lu->A, lu->perm_c); 33914b396bbSKris Buschelman lu->SuperluMatOdering = PETSC_TRUE; 34014b396bbSKris Buschelman } 34114b396bbSKris Buschelman PetscOptionsEnd(); 34214b396bbSKris Buschelman 34314b396bbSKris Buschelman /* Create the elimination tree */ 34414b396bbSKris Buschelman ierr = PetscMalloc(A->n*sizeof(int),&etree);CHKERRQ(ierr); 34514b396bbSKris Buschelman sp_preorder("N",&lu->A,lu->perm_c,etree,&lu->AC); 34614b396bbSKris Buschelman /* Factor the matrix */ 34714b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 34814b396bbSKris Buschelman zgstrf("N",&lu->AC,lu->pivot_threshold,0.0,lu->relax,lu->panel_size,etree,PETSC_NULL,0,lu->perm_r,lu->perm_c,&lu->L,&lu->U,&ierr); 34914b396bbSKris Buschelman #else 35014b396bbSKris Buschelman dgstrf("N",&lu->AC,lu->pivot_threshold,0.0,lu->relax,lu->panel_size,etree,PETSC_NULL,0,lu->perm_r,lu->perm_c,&lu->L,&lu->U,&ierr); 35114b396bbSKris Buschelman #endif 35214b396bbSKris Buschelman if (ierr < 0) { 35314b396bbSKris Buschelman SETERRQ1(PETSC_ERR_ARG_WRONG,"The diagonal element of row %d was invalid",-ierr); 35414b396bbSKris Buschelman } else if (ierr > 0) { 35514b396bbSKris Buschelman if (ierr <= A->m) { 35614b396bbSKris Buschelman SETERRQ1(PETSC_ERR_ARG_WRONG,"The diagonal element %d of U is exactly zero",ierr); 35714b396bbSKris Buschelman } else { 35814b396bbSKris Buschelman SETERRQ1(PETSC_ERR_ARG_WRONG,"Memory allocation failure after %d bytes were allocated",ierr-A->m); 35914b396bbSKris Buschelman } 36014b396bbSKris Buschelman } 36114b396bbSKris Buschelman 36214b396bbSKris Buschelman /* Cleanup */ 36314b396bbSKris Buschelman ierr = PetscFree(etree);CHKERRQ(ierr); 36475af56d4SHong Zhang #endif /* OLD */ 36514b396bbSKris Buschelman 36614b396bbSKris Buschelman } 36714b396bbSKris Buschelman 36814b396bbSKris Buschelman /* 36914b396bbSKris Buschelman Note the r permutation is ignored 37014b396bbSKris Buschelman */ 37114b396bbSKris Buschelman #undef __FUNCT__ 372f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_SuperLU" 373f0c56d0fSKris Buschelman int MatLUFactorSymbolic_SuperLU(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) 37414b396bbSKris Buschelman { 37514b396bbSKris Buschelman Mat B; 376f0c56d0fSKris Buschelman Mat_SuperLU *lu; 377*9ce50919SHong Zhang int ierr,m=A->m,n=A->n,indx; 378*9ce50919SHong Zhang PetscTruth flg; 379*9ce50919SHong Zhang char *colperm[]={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by petsc interface yet */ 380*9ce50919SHong Zhang char *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"}; 381*9ce50919SHong Zhang char *rowperm[]={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by petsc interface yet */ 38214b396bbSKris Buschelman 38314b396bbSKris Buschelman PetscFunctionBegin; 38414b396bbSKris Buschelman 385899d7b4fSKris Buschelman ierr = MatCreate(A->comm,A->m,A->n,PETSC_DETERMINE,PETSC_DETERMINE,&B);CHKERRQ(ierr); 386899d7b4fSKris Buschelman ierr = MatSetType(B,MATSUPERLU);CHKERRQ(ierr); 387899d7b4fSKris Buschelman ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 388f0c56d0fSKris Buschelman 389f0c56d0fSKris Buschelman B->ops->lufactornumeric = MatLUFactorNumeric_SuperLU; 390f0c56d0fSKris Buschelman B->ops->solve = MatSolve_SuperLU; 39114b396bbSKris Buschelman B->factor = FACTOR_LU; 392899d7b4fSKris Buschelman B->assembled = PETSC_TRUE; /* required by -sles_view */ 39314b396bbSKris Buschelman 394f0c56d0fSKris Buschelman lu = (Mat_SuperLU*)(B->spptr); 395*9ce50919SHong Zhang 396*9ce50919SHong Zhang /* Set SuperLU options */ 397*9ce50919SHong Zhang /* the default values for options argument: 398*9ce50919SHong Zhang options.Fact = DOFACT; 399*9ce50919SHong Zhang options.Equil = YES; 400*9ce50919SHong Zhang options.ColPerm = COLAMD; 401*9ce50919SHong Zhang options.DiagPivotThresh = 1.0; 402*9ce50919SHong Zhang options.Trans = NOTRANS; 403*9ce50919SHong Zhang options.IterRefine = NOREFINE; 404*9ce50919SHong Zhang options.SymmetricMode = NO; 405*9ce50919SHong Zhang options.PivotGrowth = NO; 406*9ce50919SHong Zhang options.ConditionNumber = NO; 407*9ce50919SHong Zhang options.PrintStat = YES; 408*9ce50919SHong Zhang */ 409*9ce50919SHong Zhang set_default_options(&lu->options); 410*9ce50919SHong Zhang lu->options.Equil = NO; /* equilibration causes error in solve */ 411*9ce50919SHong Zhang lu->options.PrintStat = NO; 412*9ce50919SHong Zhang 413*9ce50919SHong Zhang ierr = PetscOptionsBegin(A->comm,A->prefix,"SuperLU Options","Mat");CHKERRQ(ierr); 414*9ce50919SHong Zhang /* 415*9ce50919SHong Zhang ierr = PetscOptionsLogical("-mat_superlu_Equil","Equil","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 416*9ce50919SHong Zhang if (flg) lu->options.Equil = YES; -- do not work!!! 417*9ce50919SHong Zhang */ 418*9ce50919SHong Zhang ierr = PetscOptionsEList("-mat_superlu_ColPerm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);CHKERRQ(ierr); 419*9ce50919SHong Zhang if (flg) {lu->options.ColPerm = (colperm_t)indx;} 420*9ce50919SHong Zhang ierr = PetscOptionsEList("-mat_superlu_IterRefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);CHKERRQ(ierr); 421*9ce50919SHong Zhang if (flg) { lu->options.IterRefine = (IterRefine_t)indx;} 422*9ce50919SHong Zhang ierr = PetscOptionsLogical("-mat_superlu_SymmetricMode","SymmetricMode","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 423*9ce50919SHong Zhang if (flg) lu->options.SymmetricMode = YES; 424*9ce50919SHong Zhang ierr = PetscOptionsReal("-mat_superlu_DiagPivotThresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&lu->options.DiagPivotThresh,PETSC_NULL);CHKERRQ(ierr); 425*9ce50919SHong Zhang ierr = PetscOptionsLogical("-mat_superlu_PivotGrowth","PivotGrowth","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 426*9ce50919SHong Zhang if (flg) lu->options.PivotGrowth = YES; 427*9ce50919SHong Zhang ierr = PetscOptionsLogical("-mat_superlu_ConditionNumber","ConditionNumber","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 428*9ce50919SHong Zhang if (flg) lu->options.ConditionNumber = YES; 429*9ce50919SHong Zhang ierr = PetscOptionsEList("-mat_superlu_RowPerm","RowPerm","None",rowperm,2,rowperm[0],&indx,&flg);CHKERRQ(ierr); 430*9ce50919SHong Zhang if (flg) {lu->options.RowPerm = (rowperm_t)indx;} 431*9ce50919SHong Zhang ierr = PetscOptionsLogical("-mat_superlu_ReplaceTinyPivot","ReplaceTinyPivot","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 432*9ce50919SHong Zhang if (flg) lu->options.ReplaceTinyPivot = YES; 433*9ce50919SHong Zhang ierr = PetscOptionsLogical("-mat_superlu_PrintStat","PrintStat","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 434*9ce50919SHong Zhang if (flg) lu->options.PrintStat = YES; 435*9ce50919SHong Zhang PetscOptionsEnd(); 436*9ce50919SHong Zhang 43775af56d4SHong Zhang #ifdef SUPERLU2 4382ecb5974SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatCreateNull","MatCreateNull_SuperLU", 439f0c56d0fSKris Buschelman (void(*)(void))MatCreateNull_SuperLU);CHKERRQ(ierr); 44075af56d4SHong Zhang #endif 44114b396bbSKris Buschelman 44275af56d4SHong Zhang /* Allocate spaces (notice sizes are for the transpose) */ 443*9ce50919SHong Zhang ierr = PetscMalloc(m*sizeof(int),&lu->etree);CHKERRQ(ierr); 444*9ce50919SHong Zhang ierr = PetscMalloc(n*sizeof(int),&lu->perm_r);CHKERRQ(ierr); 445*9ce50919SHong Zhang ierr = PetscMalloc(m*sizeof(int),&lu->perm_c);CHKERRQ(ierr); 446*9ce50919SHong Zhang ierr = PetscMalloc(n*sizeof(int),&lu->R);CHKERRQ(ierr); 447*9ce50919SHong Zhang ierr = PetscMalloc(m*sizeof(int),&lu->C);CHKERRQ(ierr); 44875af56d4SHong Zhang 449*9ce50919SHong Zhang /* create rhs and solution x without allocate space for .Store */ 45075af56d4SHong Zhang dCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE); 45175af56d4SHong Zhang dCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE); 45214b396bbSKris Buschelman 45314b396bbSKris Buschelman lu->flg = DIFFERENT_NONZERO_PATTERN; 454e740cb95SKris Buschelman lu->CleanUpSuperLU = PETSC_TRUE; 455899d7b4fSKris Buschelman 456899d7b4fSKris Buschelman *F = B; 457*9ce50919SHong Zhang PetscLogObjectMemory(B,(A->m+A->n)*sizeof(int)+sizeof(Mat_SuperLU)); 45814b396bbSKris Buschelman PetscFunctionReturn(0); 45914b396bbSKris Buschelman } 46014b396bbSKris Buschelman 46114b396bbSKris Buschelman /* used by -sles_view */ 46214b396bbSKris Buschelman #undef __FUNCT__ 463f0c56d0fSKris Buschelman #define __FUNCT__ "MatFactorInfo_SuperLU" 464f0c56d0fSKris Buschelman int MatFactorInfo_SuperLU(Mat A,PetscViewer viewer) 46514b396bbSKris Buschelman { 466f0c56d0fSKris Buschelman Mat_SuperLU *lu= (Mat_SuperLU*)A->spptr; 467f0c56d0fSKris Buschelman int ierr; 468*9ce50919SHong Zhang superlu_options_t options; 469f0c56d0fSKris Buschelman 470f0c56d0fSKris Buschelman PetscFunctionBegin; 471*9ce50919SHong Zhang /* check if matrix is superlu_dist type */ 472*9ce50919SHong Zhang if (A->ops->solve != MatSolve_SuperLU) PetscFunctionReturn(0); 473*9ce50919SHong Zhang 474*9ce50919SHong Zhang options = lu->options; 475f0c56d0fSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");CHKERRQ(ierr); 476*9ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Equil: %s\n",(options.Equil != NO) ? "YES": "NO");CHKERRQ(ierr); 477*9ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ColPerm: %d\n",options.ColPerm);CHKERRQ(ierr); 478*9ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," IterRefine: %d\n",options.IterRefine);CHKERRQ(ierr); 479*9ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," SymmetricMode: %s\n",(options.SymmetricMode != NO) ? "YES": "NO");CHKERRQ(ierr); 480*9ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," DiagPivotThresh: %g\n",options.DiagPivotThresh);CHKERRQ(ierr); 481*9ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," PivotGrowth: %s\n",(options.PivotGrowth != NO) ? "YES": "NO");CHKERRQ(ierr); 482*9ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ConditionNumber: %s\n",(options.ConditionNumber != NO) ? "YES": "NO");CHKERRQ(ierr); 483*9ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RowPerm: %d\n",options.RowPerm);CHKERRQ(ierr); 484*9ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ReplaceTinyPivot: %s\n",(options.ReplaceTinyPivot != NO) ? "YES": "NO");CHKERRQ(ierr); 485*9ce50919SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," PrintStat: %s\n",(options.PrintStat != NO) ? "YES": "NO");CHKERRQ(ierr); 486*9ce50919SHong Zhang 487f0c56d0fSKris Buschelman PetscFunctionReturn(0); 488f0c56d0fSKris Buschelman } 489f0c56d0fSKris Buschelman 490f0c56d0fSKris Buschelman #undef __FUNCT__ 491f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_SuperLU" 492f0c56d0fSKris Buschelman int MatDuplicate_SuperLU(Mat A, MatDuplicateOption op, Mat *M) { 49314b396bbSKris Buschelman int ierr; 4948f340917SKris Buschelman Mat_SuperLU *lu=(Mat_SuperLU *)A->spptr; 4958f340917SKris Buschelman 49614b396bbSKris Buschelman PetscFunctionBegin; 4978f340917SKris Buschelman ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr); 498f0c56d0fSKris Buschelman ierr = MatConvert_SeqAIJ_SuperLU(*M,MATSUPERLU,M);CHKERRQ(ierr); 4993f953163SKris Buschelman ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_SuperLU));CHKERRQ(ierr); 50014b396bbSKris Buschelman PetscFunctionReturn(0); 50114b396bbSKris Buschelman } 50214b396bbSKris Buschelman 50314b396bbSKris Buschelman EXTERN_C_BEGIN 50414b396bbSKris Buschelman #undef __FUNCT__ 505b0592601SKris Buschelman #define __FUNCT__ "MatConvert_SuperLU_SeqAIJ" 506b0592601SKris Buschelman int MatConvert_SuperLU_SeqAIJ(Mat A,MatType type,Mat *newmat) { 507fb3e25aaSKris Buschelman /* This routine is only called to convert an unfactored PETSc-SuperLU matrix */ 508fb3e25aaSKris Buschelman /* to its base PETSc type, so we will ignore 'MatType type'. */ 50914b396bbSKris Buschelman int ierr; 510b0592601SKris Buschelman Mat B=*newmat; 511f0c56d0fSKris Buschelman Mat_SuperLU *lu=(Mat_SuperLU *)A->spptr; 512b0592601SKris Buschelman 513b0592601SKris Buschelman PetscFunctionBegin; 514b0592601SKris Buschelman if (B != A) { 515b0592601SKris Buschelman ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 516f0c56d0fSKris Buschelman } 517b0592601SKris Buschelman /* Reset the original function pointers */ 518f0c56d0fSKris Buschelman B->ops->duplicate = lu->MatDuplicate; 519b0592601SKris Buschelman B->ops->view = lu->MatView; 520b0592601SKris Buschelman B->ops->assemblyend = lu->MatAssemblyEnd; 521b0592601SKris Buschelman B->ops->lufactorsymbolic = lu->MatLUFactorSymbolic; 522b0592601SKris Buschelman B->ops->destroy = lu->MatDestroy; 523b0592601SKris Buschelman /* lu is only a function pointer stash unless we've factored the matrix, which we haven't! */ 524b0592601SKris Buschelman ierr = PetscFree(lu);CHKERRQ(ierr); 525b0592601SKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);CHKERRQ(ierr); 526b0592601SKris Buschelman *newmat = B; 527b0592601SKris Buschelman PetscFunctionReturn(0); 528b0592601SKris Buschelman } 529b0592601SKris Buschelman EXTERN_C_END 530b0592601SKris Buschelman 531b0592601SKris Buschelman EXTERN_C_BEGIN 532b0592601SKris Buschelman #undef __FUNCT__ 533b0592601SKris Buschelman #define __FUNCT__ "MatConvert_SeqAIJ_SuperLU" 534b0592601SKris Buschelman int MatConvert_SeqAIJ_SuperLU(Mat A,MatType type,Mat *newmat) { 535fb3e25aaSKris Buschelman /* This routine is only called to convert to MATSUPERLU */ 536fb3e25aaSKris Buschelman /* from MATSEQAIJ, so we will ignore 'MatType type'. */ 537b0592601SKris Buschelman int ierr; 538b0592601SKris Buschelman Mat B=*newmat; 539f0c56d0fSKris Buschelman Mat_SuperLU *lu; 54014b396bbSKris Buschelman 54114b396bbSKris Buschelman PetscFunctionBegin; 542b0592601SKris Buschelman if (B != A) { 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; 559f0c56d0fSKris Buschelman B->ops->destroy = MatDestroy_SuperLU; 560b0592601SKris Buschelman 561b0592601SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_superlu_C", 562b0592601SKris Buschelman "MatConvert_SeqAIJ_SuperLU",MatConvert_SeqAIJ_SuperLU);CHKERRQ(ierr); 563b0592601SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_superlu_seqaij_C", 564b0592601SKris Buschelman "MatConvert_SuperLU_SeqAIJ",MatConvert_SuperLU_SeqAIJ);CHKERRQ(ierr); 565fb3e25aaSKris Buschelman PetscLogInfo(0,"Using SuperLU for SeqAIJ LU factorization and solves."); 566fb3e25aaSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)B,MATSUPERLU);CHKERRQ(ierr); 567b0592601SKris Buschelman *newmat = B; 568b0592601SKris Buschelman PetscFunctionReturn(0); 569b0592601SKris Buschelman } 570b0592601SKris Buschelman EXTERN_C_END 571b0592601SKris Buschelman 57224b6179bSKris Buschelman /*MC 573fafad747SKris Buschelman MATSUPERLU - MATSUPERLU = "superlu" - A matrix type providing direct solvers (LU) for sequential matrices 57424b6179bSKris Buschelman via the external package SuperLU. 57524b6179bSKris Buschelman 57624b6179bSKris Buschelman If SuperLU is installed (see the manual for 57724b6179bSKris Buschelman instructions on how to declare the existence of external packages), 57824b6179bSKris Buschelman a matrix type can be constructed which invokes SuperLU solvers. 57924b6179bSKris Buschelman After calling MatCreate(...,A), simply call MatSetType(A,MATSUPERLU). 58024b6179bSKris Buschelman This matrix type is only supported for double precision real. 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" 601f0c56d0fSKris Buschelman int MatCreate_SuperLU(Mat A) { 602b0592601SKris Buschelman int ierr; 603b0592601SKris Buschelman 604b0592601SKris Buschelman PetscFunctionBegin; 6055441df8eSKris Buschelman /* Change type name before calling MatSetType to force proper construction of SeqAIJ and SUPERLU types */ 6065441df8eSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)A,MATSUPERLU);CHKERRQ(ierr); 607b0592601SKris Buschelman ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr); 608b0592601SKris Buschelman ierr = MatConvert_SeqAIJ_SuperLU(A,MATSUPERLU,&A);CHKERRQ(ierr); 60914b396bbSKris Buschelman PetscFunctionReturn(0); 61014b396bbSKris Buschelman } 61114b396bbSKris Buschelman EXTERN_C_END 612