1dba47a55SKris Buschelman 29b54502bSHong Zhang /* 39b54502bSHong Zhang Defines a direct factorization preconditioner for any Mat implementation 49b54502bSHong Zhang Note: this need not be consided a preconditioner since it supplies 59b54502bSHong Zhang a direct solver. 69b54502bSHong Zhang */ 7ee45ca4aSHong Zhang 8c6db04a5SJed Brown #include <../src/ksp/pc/impls/factor/lu/lu.h> /*I "petscpc.h" I*/ 99b54502bSHong Zhang 107087cfbeSBarry Smith PetscErrorCode PCFactorReorderForNonzeroDiagonal_LU(PC pc,PetscReal z) 119b54502bSHong Zhang { 129b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 139b54502bSHong Zhang 149b54502bSHong Zhang PetscFunctionBegin; 159b54502bSHong Zhang lu->nonzerosalongdiagonal = PETSC_TRUE; 162fa5cd67SKarl Rupp if (z == PETSC_DECIDE) lu->nonzerosalongdiagonaltol = 1.e-10; 172fa5cd67SKarl Rupp else lu->nonzerosalongdiagonaltol = z; 189b54502bSHong Zhang PetscFunctionReturn(0); 199b54502bSHong Zhang } 209b54502bSHong Zhang 214416b707SBarry Smith static PetscErrorCode PCSetFromOptions_LU(PetscOptionItems *PetscOptionsObject,PC pc) 229b54502bSHong Zhang { 239b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 24ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE; 259b54502bSHong Zhang PetscReal tol; 269b54502bSHong Zhang 279b54502bSHong Zhang PetscFunctionBegin; 28*d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"LU options"); 299566063dSJacob Faibussowitsch PetscCall(PCSetFromOptions_Factor(PetscOptionsObject,pc)); 305c9eb25fSBarry Smith 319566063dSJacob Faibussowitsch PetscCall(PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg)); 329b54502bSHong Zhang if (flg) { 339b54502bSHong Zhang tol = PETSC_DECIDE; 349566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,NULL)); 359566063dSJacob Faibussowitsch PetscCall(PCFactorReorderForNonzeroDiagonal(pc,tol)); 369b54502bSHong Zhang } 37*d0609cedSBarry Smith PetscOptionsHeadEnd(); 389b54502bSHong Zhang PetscFunctionReturn(0); 399b54502bSHong Zhang } 409b54502bSHong Zhang 419b54502bSHong Zhang static PetscErrorCode PCSetUp_LU(PC pc) 429b54502bSHong Zhang { 439b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 44ea799195SBarry Smith MatSolverType stype; 4500e125f8SBarry Smith MatFactorError err; 463d1c1ea0SBarry Smith 479b54502bSHong Zhang PetscFunctionBegin; 48c6e4fdc6SHong Zhang pc->failedreason = PC_NOERROR; 493d1c1ea0SBarry Smith if (dir->hdr.reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->hdr.actualfill; 509b54502bSHong Zhang 519566063dSJacob Faibussowitsch PetscCall(MatSetErrorIfFailure(pc->pmat,pc->erroriffailure)); 523d1c1ea0SBarry Smith if (dir->hdr.inplace) { 539d76b4d0SMatthew G. Knepley MatFactorType ftype; 549d76b4d0SMatthew G. Knepley 559566063dSJacob Faibussowitsch PetscCall(MatGetFactorType(pc->pmat, &ftype)); 569d76b4d0SMatthew G. Knepley if (ftype == MAT_FACTOR_NONE) { 579566063dSJacob Faibussowitsch if (dir->row && dir->col && dir->row != dir->col) PetscCall(ISDestroy(&dir->row)); 589566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->col)); 592c7c0729SBarry Smith /* This should only get the ordering if needed, but since MatGetFactor() is not called we can't know if it is needed */ 609566063dSJacob Faibussowitsch PetscCall(PCFactorSetDefaultOrdering_Factor(pc)); 619566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col)); 6203c60df9SBarry Smith if (dir->row) { 639566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row)); 649566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col)); 6503c60df9SBarry Smith } 669566063dSJacob Faibussowitsch PetscCall(MatLUFactor(pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info)); 679566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(pc->pmat,&err)); 6800e125f8SBarry Smith if (err) { /* Factor() fails */ 6900e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 706baea169SHong Zhang PetscFunctionReturn(0); 716baea169SHong Zhang } 729d76b4d0SMatthew G. Knepley } 73075768bcSBarry Smith ((PC_Factor*)dir)->fact = pc->pmat; 749b54502bSHong Zhang } else { 759b54502bSHong Zhang MatInfo info; 7600e125f8SBarry Smith 779b54502bSHong Zhang if (!pc->setupcalled) { 78f73b0415SBarry Smith PetscBool canuseordering; 792c7c0729SBarry Smith if (!((PC_Factor*)dir)->fact) { 809566063dSJacob Faibussowitsch PetscCall(MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact)); 819566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact)); 822c7c0729SBarry Smith } 839566063dSJacob Faibussowitsch PetscCall(MatFactorGetCanUseOrdering(((PC_Factor*)dir)->fact,&canuseordering)); 84f73b0415SBarry Smith if (canuseordering) { 859566063dSJacob Faibussowitsch PetscCall(PCFactorSetDefaultOrdering_Factor(pc)); 869566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col)); 879b54502bSHong Zhang if (dir->nonzerosalongdiagonal) { 889566063dSJacob Faibussowitsch PetscCall(MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col)); 899b54502bSHong Zhang } 909566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row)); 919566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col)); 9203c60df9SBarry Smith } 939566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info)); 949566063dSJacob Faibussowitsch PetscCall(MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info)); 953d1c1ea0SBarry Smith dir->hdr.actualfill = info.fill_ratio_needed; 969b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 97f73b0415SBarry Smith PetscBool canuseordering; 983d1c1ea0SBarry Smith if (!dir->hdr.reuseordering) { 999566063dSJacob Faibussowitsch PetscCall(MatDestroy(&((PC_Factor*)dir)->fact)); 1009566063dSJacob Faibussowitsch PetscCall(MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact)); 1019566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact)); 1029566063dSJacob Faibussowitsch PetscCall(MatFactorGetCanUseOrdering(((PC_Factor*)dir)->fact,&canuseordering)); 103f73b0415SBarry Smith if (canuseordering) { 1049566063dSJacob Faibussowitsch if (dir->row && dir->col && dir->row != dir->col) PetscCall(ISDestroy(&dir->row)); 1059566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->col)); 1069566063dSJacob Faibussowitsch PetscCall(PCFactorSetDefaultOrdering_Factor(pc)); 1079566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col)); 1089b54502bSHong Zhang if (dir->nonzerosalongdiagonal) { 1099566063dSJacob Faibussowitsch PetscCall(MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col)); 1109b54502bSHong Zhang } 1119566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row)); 1129566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col)); 11303c60df9SBarry Smith } 1149b54502bSHong Zhang } 1159566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info)); 1169566063dSJacob Faibussowitsch PetscCall(MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info)); 1173d1c1ea0SBarry Smith dir->hdr.actualfill = info.fill_ratio_needed; 11804545d6dSBarry Smith } else { 1199566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor*)dir)->fact,&err)); 120160a8794SBarry Smith if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) { 1219566063dSJacob Faibussowitsch PetscCall(MatFactorClearError(((PC_Factor*)dir)->fact)); 122b8b68cfdSBarry Smith pc->failedreason = PC_NOERROR; 12304545d6dSBarry Smith } 1249b54502bSHong Zhang } 1259566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor*)dir)->fact,&err)); 12600e125f8SBarry Smith if (err) { /* FactorSymbolic() fails */ 12700e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1288c1cd74cSHong Zhang PetscFunctionReturn(0); 1298c1cd74cSHong Zhang } 1308c1cd74cSHong Zhang 1319566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info)); 1329566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor*)dir)->fact,&err)); 13300e125f8SBarry Smith if (err) { /* FactorNumeric() fails */ 13400e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1358c1cd74cSHong Zhang } 136680c5173SHong Zhang 1379b54502bSHong Zhang } 13800c67f3bSHong Zhang 1399566063dSJacob Faibussowitsch PetscCall(PCFactorGetMatSolverType(pc,&stype)); 14000c67f3bSHong Zhang if (!stype) { 141ea799195SBarry Smith MatSolverType solverpackage; 1429566063dSJacob Faibussowitsch PetscCall(MatFactorGetSolverType(((PC_Factor*)dir)->fact,&solverpackage)); 1439566063dSJacob Faibussowitsch PetscCall(PCFactorSetMatSolverType(pc,solverpackage)); 14400c67f3bSHong Zhang } 1459b54502bSHong Zhang PetscFunctionReturn(0); 1469b54502bSHong Zhang } 1479b54502bSHong Zhang 148574deadeSBarry Smith static PetscErrorCode PCReset_LU(PC pc) 1499b54502bSHong Zhang { 1509b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 1519b54502bSHong Zhang 1529b54502bSHong Zhang PetscFunctionBegin; 1539566063dSJacob Faibussowitsch if (!dir->hdr.inplace && ((PC_Factor*)dir)->fact) PetscCall(MatDestroy(&((PC_Factor*)dir)->fact)); 1549566063dSJacob Faibussowitsch if (dir->row && dir->col && dir->row != dir->col) PetscCall(ISDestroy(&dir->row)); 1559566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->col)); 156574deadeSBarry Smith PetscFunctionReturn(0); 157574deadeSBarry Smith } 158574deadeSBarry Smith 159574deadeSBarry Smith static PetscErrorCode PCDestroy_LU(PC pc) 160574deadeSBarry Smith { 161574deadeSBarry Smith PC_LU *dir = (PC_LU*)pc->data; 162574deadeSBarry Smith 163574deadeSBarry Smith PetscFunctionBegin; 1649566063dSJacob Faibussowitsch PetscCall(PCReset_LU(pc)); 1659566063dSJacob Faibussowitsch PetscCall(PetscFree(((PC_Factor*)dir)->ordering)); 1669566063dSJacob Faibussowitsch PetscCall(PetscFree(((PC_Factor*)dir)->solvertype)); 1679566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 1689b54502bSHong Zhang PetscFunctionReturn(0); 1699b54502bSHong Zhang } 1709b54502bSHong Zhang 1719b54502bSHong Zhang static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y) 1729b54502bSHong Zhang { 1739b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 1749b54502bSHong Zhang 1759b54502bSHong Zhang PetscFunctionBegin; 1763d1c1ea0SBarry Smith if (dir->hdr.inplace) { 1779566063dSJacob Faibussowitsch PetscCall(MatSolve(pc->pmat,x,y)); 1782fa5cd67SKarl Rupp } else { 1799566063dSJacob Faibussowitsch PetscCall(MatSolve(((PC_Factor*)dir)->fact,x,y)); 1802fa5cd67SKarl Rupp } 1819b54502bSHong Zhang PetscFunctionReturn(0); 1829b54502bSHong Zhang } 1839b54502bSHong Zhang 1847b6e2003SPierre Jolivet static PetscErrorCode PCMatApply_LU(PC pc,Mat X,Mat Y) 1857b6e2003SPierre Jolivet { 1867b6e2003SPierre Jolivet PC_LU *dir = (PC_LU*)pc->data; 1877b6e2003SPierre Jolivet 1887b6e2003SPierre Jolivet PetscFunctionBegin; 1897b6e2003SPierre Jolivet if (dir->hdr.inplace) { 1909566063dSJacob Faibussowitsch PetscCall(MatMatSolve(pc->pmat,X,Y)); 1917b6e2003SPierre Jolivet } else { 1929566063dSJacob Faibussowitsch PetscCall(MatMatSolve(((PC_Factor*)dir)->fact,X,Y)); 1937b6e2003SPierre Jolivet } 1947b6e2003SPierre Jolivet PetscFunctionReturn(0); 1957b6e2003SPierre Jolivet } 1967b6e2003SPierre Jolivet 1979b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y) 1989b54502bSHong Zhang { 1999b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 2009b54502bSHong Zhang 2019b54502bSHong Zhang PetscFunctionBegin; 2023d1c1ea0SBarry Smith if (dir->hdr.inplace) { 2039566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(pc->pmat,x,y)); 2042fa5cd67SKarl Rupp } else { 2059566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(((PC_Factor*)dir)->fact,x,y)); 2062fa5cd67SKarl Rupp } 2079b54502bSHong Zhang PetscFunctionReturn(0); 2089b54502bSHong Zhang } 2099b54502bSHong Zhang 2109b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 2119b54502bSHong Zhang 2129b54502bSHong Zhang /*MC 2139b54502bSHong Zhang PCLU - Uses a direct solver, based on LU factorization, as a preconditioner 2149b54502bSHong Zhang 2159b54502bSHong Zhang Options Database Keys: 2162401956bSBarry Smith + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 2173ca39a21SBarry Smith . -pc_factor_mat_solver_type - Actives PCFactorSetMatSolverType() to choose the direct solver, like superlu 2182401956bSBarry Smith . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 21955ba2a51SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 2202401956bSBarry Smith . -pc_factor_in_place - Activates in-place factorization 2212401956bSBarry Smith . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 2222401956bSBarry Smith . -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase 2239b54502bSHong Zhang stability of factorization. 224145b9266SHong Zhang . -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types 225145b9266SHong Zhang . -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default 226e22d95b2SBarry Smith - -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the 227e22d95b2SBarry Smith diagonal. 2289b54502bSHong Zhang 22995452b02SPatrick Sanan Notes: 23095452b02SPatrick Sanan Not all options work for all matrix formats 2319b54502bSHong Zhang Run with -help to see additional options for particular matrix formats or factorization 2329b54502bSHong Zhang algorithms 2339b54502bSHong Zhang 2349b54502bSHong Zhang Level: beginner 2359b54502bSHong Zhang 23695452b02SPatrick Sanan Notes: 23795452b02SPatrick Sanan Usually this will compute an "exact" solution in one iteration and does 2389b54502bSHong Zhang not need a Krylov method (i.e. you can use -ksp_type preonly, or 2399b54502bSHong Zhang KSPSetType(ksp,KSPPREONLY) for the Krylov method 2409b54502bSHong Zhang 2419b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 242a4fd02acSBarry Smith PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(), 2438ff23777SHong Zhang PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetColumnPivot(), 24452e318dbSPierre Jolivet PCFactorSetPivotInBlocks(),PCFactorSetShiftType(),PCFactorSetShiftAmount() 2458ff23777SHong Zhang PCFactorReorderForNonzeroDiagonal() 2469b54502bSHong Zhang M*/ 2479b54502bSHong Zhang 2488cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_LU(PC pc) 2499b54502bSHong Zhang { 2509b54502bSHong Zhang PC_LU *dir; 2519b54502bSHong Zhang 2529b54502bSHong Zhang PetscFunctionBegin; 2539566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc,&dir)); 2543d1c1ea0SBarry Smith pc->data = (void*)dir; 2559566063dSJacob Faibussowitsch PetscCall(PCFactorInitialize(pc,MAT_FACTOR_LU)); 2569b54502bSHong Zhang dir->nonzerosalongdiagonal = PETSC_FALSE; 2579b54502bSHong Zhang 258075768bcSBarry Smith ((PC_Factor*)dir)->info.fill = 5.0; 259075768bcSBarry Smith ((PC_Factor*)dir)->info.dtcol = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */ 260f4db908eSBarry Smith ((PC_Factor*)dir)->info.shifttype = (PetscReal)MAT_SHIFT_NONE; 2610a545947SLisandro Dalcin dir->col = NULL; 2620a545947SLisandro Dalcin dir->row = NULL; 2635c9eb25fSBarry Smith 264574deadeSBarry Smith pc->ops->reset = PCReset_LU; 2659b54502bSHong Zhang pc->ops->destroy = PCDestroy_LU; 2669b54502bSHong Zhang pc->ops->apply = PCApply_LU; 2677b6e2003SPierre Jolivet pc->ops->matapply = PCMatApply_LU; 2689b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_LU; 2699b54502bSHong Zhang pc->ops->setup = PCSetUp_LU; 2709b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_LU; 27192e08861SBarry Smith pc->ops->view = PCView_Factor; 2720a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 2739566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",PCFactorReorderForNonzeroDiagonal_LU)); 2749b54502bSHong Zhang PetscFunctionReturn(0); 2759b54502bSHong Zhang } 276