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 10d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorReorderForNonzeroDiagonal_LU(PC pc, PetscReal z) 11d71ae5a4SJacob Faibussowitsch { 129b54502bSHong Zhang PC_LU *lu = (PC_LU *)pc->data; 139b54502bSHong Zhang 149b54502bSHong Zhang PetscFunctionBegin; 159b54502bSHong Zhang lu->nonzerosalongdiagonal = PETSC_TRUE; 1613bcc0bdSJacob Faibussowitsch if (z == (PetscReal)PETSC_DECIDE) lu->nonzerosalongdiagonaltol = 1.e-10; 172fa5cd67SKarl Rupp else lu->nonzerosalongdiagonaltol = z; 183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 199b54502bSHong Zhang } 209b54502bSHong Zhang 21d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_LU(PC pc, PetscOptionItems *PetscOptionsObject) 22d71ae5a4SJacob Faibussowitsch { 239b54502bSHong Zhang PC_LU *lu = (PC_LU *)pc->data; 24ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE; 259b54502bSHong Zhang PetscReal tol; 269b54502bSHong Zhang 279b54502bSHong Zhang PetscFunctionBegin; 28d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "LU options"); 29dbbe0bcdSBarry Smith PetscCall(PCSetFromOptions_Factor(pc, PetscOptionsObject)); 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 } 37d0609cedSBarry Smith PetscOptionsHeadEnd(); 383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 399b54502bSHong Zhang } 409b54502bSHong Zhang 41d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_LU(PC pc) 42d71ae5a4SJacob Faibussowitsch { 439b54502bSHong Zhang PC_LU *dir = (PC_LU *)pc->data; 44ea799195SBarry Smith MatSolverType stype; 4500e125f8SBarry Smith MatFactorError err; 46f023e1d5SPierre Jolivet const char *prefix; 473d1c1ea0SBarry Smith 489b54502bSHong Zhang PetscFunctionBegin; 49c6e4fdc6SHong Zhang pc->failedreason = PC_NOERROR; 503d1c1ea0SBarry Smith if (dir->hdr.reusefill && pc->setupcalled) ((PC_Factor *)dir)->info.fill = dir->hdr.actualfill; 519b54502bSHong Zhang 5226cc229bSBarry Smith PetscCall(PCGetOptionsPrefix(pc, &prefix)); 5326cc229bSBarry Smith PetscCall(MatSetOptionsPrefixFactor(pc->pmat, prefix)); 5426cc229bSBarry Smith 559566063dSJacob Faibussowitsch PetscCall(MatSetErrorIfFailure(pc->pmat, pc->erroriffailure)); 563d1c1ea0SBarry Smith if (dir->hdr.inplace) { 579d76b4d0SMatthew G. Knepley MatFactorType ftype; 589d76b4d0SMatthew G. Knepley 599566063dSJacob Faibussowitsch PetscCall(MatGetFactorType(pc->pmat, &ftype)); 609d76b4d0SMatthew G. Knepley if (ftype == MAT_FACTOR_NONE) { 619566063dSJacob Faibussowitsch if (dir->row && dir->col && dir->row != dir->col) PetscCall(ISDestroy(&dir->row)); 629566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->col)); 632c7c0729SBarry Smith /* This should only get the ordering if needed, but since MatGetFactor() is not called we can't know if it is needed */ 649566063dSJacob Faibussowitsch PetscCall(PCFactorSetDefaultOrdering_Factor(pc)); 659566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col)); 664dfa11a4SJacob Faibussowitsch if (dir->row) { } 679566063dSJacob Faibussowitsch PetscCall(MatLUFactor(pc->pmat, dir->row, dir->col, &((PC_Factor *)dir)->info)); 689566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(pc->pmat, &err)); 6900e125f8SBarry Smith if (err) { /* Factor() fails */ 7000e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 726baea169SHong Zhang } 739d76b4d0SMatthew G. Knepley } 74075768bcSBarry Smith ((PC_Factor *)dir)->fact = pc->pmat; 759b54502bSHong Zhang } else { 769b54502bSHong Zhang MatInfo info; 7700e125f8SBarry Smith 789b54502bSHong Zhang if (!pc->setupcalled) { 79f73b0415SBarry Smith PetscBool canuseordering; 80*03e5aca4SStefano Zampini 81*03e5aca4SStefano Zampini PetscCall(PCFactorSetUpMatSolverType(pc)); 829566063dSJacob Faibussowitsch PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)dir)->fact, &canuseordering)); 83f73b0415SBarry Smith if (canuseordering) { 849566063dSJacob Faibussowitsch PetscCall(PCFactorSetDefaultOrdering_Factor(pc)); 859566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col)); 861baa6e33SBarry Smith if (dir->nonzerosalongdiagonal) PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, dir->nonzerosalongdiagonaltol, dir->row, dir->col)); 8703c60df9SBarry Smith } 889566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->row, dir->col, &((PC_Factor *)dir)->info)); 899566063dSJacob Faibussowitsch PetscCall(MatGetInfo(((PC_Factor *)dir)->fact, MAT_LOCAL, &info)); 903d1c1ea0SBarry Smith dir->hdr.actualfill = info.fill_ratio_needed; 919b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 92f73b0415SBarry Smith PetscBool canuseordering; 93*03e5aca4SStefano Zampini 943d1c1ea0SBarry Smith if (!dir->hdr.reuseordering) { 959566063dSJacob Faibussowitsch PetscCall(MatDestroy(&((PC_Factor *)dir)->fact)); 96*03e5aca4SStefano Zampini PetscCall(PCFactorSetUpMatSolverType(pc)); 979566063dSJacob Faibussowitsch PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)dir)->fact, &canuseordering)); 98f73b0415SBarry Smith if (canuseordering) { 999566063dSJacob Faibussowitsch if (dir->row && dir->col && dir->row != dir->col) PetscCall(ISDestroy(&dir->row)); 1009566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->col)); 1019566063dSJacob Faibussowitsch PetscCall(PCFactorSetDefaultOrdering_Factor(pc)); 1029566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col)); 1031baa6e33SBarry Smith if (dir->nonzerosalongdiagonal) PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, dir->nonzerosalongdiagonaltol, dir->row, dir->col)); 10403c60df9SBarry Smith } 1059b54502bSHong Zhang } 1069566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->row, dir->col, &((PC_Factor *)dir)->info)); 1079566063dSJacob Faibussowitsch PetscCall(MatGetInfo(((PC_Factor *)dir)->fact, MAT_LOCAL, &info)); 1083d1c1ea0SBarry Smith dir->hdr.actualfill = info.fill_ratio_needed; 10904545d6dSBarry Smith } else { 1109566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err)); 111160a8794SBarry Smith if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) { 1129566063dSJacob Faibussowitsch PetscCall(MatFactorClearError(((PC_Factor *)dir)->fact)); 113b8b68cfdSBarry Smith pc->failedreason = PC_NOERROR; 11404545d6dSBarry Smith } 1159b54502bSHong Zhang } 1169566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err)); 11700e125f8SBarry Smith if (err) { /* FactorSymbolic() fails */ 11800e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1208c1cd74cSHong Zhang } 1218c1cd74cSHong Zhang 1229566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(((PC_Factor *)dir)->fact, pc->pmat, &((PC_Factor *)dir)->info)); 1239566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err)); 12400e125f8SBarry Smith if (err) { /* FactorNumeric() fails */ 12500e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1268c1cd74cSHong Zhang } 1279b54502bSHong Zhang } 12800c67f3bSHong Zhang 1299566063dSJacob Faibussowitsch PetscCall(PCFactorGetMatSolverType(pc, &stype)); 13000c67f3bSHong Zhang if (!stype) { 131ea799195SBarry Smith MatSolverType solverpackage; 1329566063dSJacob Faibussowitsch PetscCall(MatFactorGetSolverType(((PC_Factor *)dir)->fact, &solverpackage)); 1339566063dSJacob Faibussowitsch PetscCall(PCFactorSetMatSolverType(pc, solverpackage)); 13400c67f3bSHong Zhang } 1353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1369b54502bSHong Zhang } 1379b54502bSHong Zhang 138d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_LU(PC pc) 139d71ae5a4SJacob Faibussowitsch { 1409b54502bSHong Zhang PC_LU *dir = (PC_LU *)pc->data; 1419b54502bSHong Zhang 1429b54502bSHong Zhang PetscFunctionBegin; 1439566063dSJacob Faibussowitsch if (!dir->hdr.inplace && ((PC_Factor *)dir)->fact) PetscCall(MatDestroy(&((PC_Factor *)dir)->fact)); 1449566063dSJacob Faibussowitsch if (dir->row && dir->col && dir->row != dir->col) PetscCall(ISDestroy(&dir->row)); 1459566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->col)); 1463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 147574deadeSBarry Smith } 148574deadeSBarry Smith 149d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_LU(PC pc) 150d71ae5a4SJacob Faibussowitsch { 151574deadeSBarry Smith PC_LU *dir = (PC_LU *)pc->data; 152574deadeSBarry Smith 153574deadeSBarry Smith PetscFunctionBegin; 1549566063dSJacob Faibussowitsch PetscCall(PCReset_LU(pc)); 1559566063dSJacob Faibussowitsch PetscCall(PetscFree(((PC_Factor *)dir)->ordering)); 1569566063dSJacob Faibussowitsch PetscCall(PetscFree(((PC_Factor *)dir)->solvertype)); 1572e956fe4SStefano Zampini PetscCall(PCFactorClearComposedFunctions(pc)); 1589566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 1593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1609b54502bSHong Zhang } 1619b54502bSHong Zhang 162d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_LU(PC pc, Vec x, Vec y) 163d71ae5a4SJacob Faibussowitsch { 1649b54502bSHong Zhang PC_LU *dir = (PC_LU *)pc->data; 1659b54502bSHong Zhang 1669b54502bSHong Zhang PetscFunctionBegin; 1673d1c1ea0SBarry Smith if (dir->hdr.inplace) { 1689566063dSJacob Faibussowitsch PetscCall(MatSolve(pc->pmat, x, y)); 1692fa5cd67SKarl Rupp } else { 1709566063dSJacob Faibussowitsch PetscCall(MatSolve(((PC_Factor *)dir)->fact, x, y)); 1712fa5cd67SKarl Rupp } 1723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1739b54502bSHong Zhang } 1749b54502bSHong Zhang 175d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_LU(PC pc, Mat X, Mat Y) 176d71ae5a4SJacob Faibussowitsch { 1777b6e2003SPierre Jolivet PC_LU *dir = (PC_LU *)pc->data; 1787b6e2003SPierre Jolivet 1797b6e2003SPierre Jolivet PetscFunctionBegin; 1807b6e2003SPierre Jolivet if (dir->hdr.inplace) { 1819566063dSJacob Faibussowitsch PetscCall(MatMatSolve(pc->pmat, X, Y)); 1827b6e2003SPierre Jolivet } else { 1839566063dSJacob Faibussowitsch PetscCall(MatMatSolve(((PC_Factor *)dir)->fact, X, Y)); 1847b6e2003SPierre Jolivet } 1853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1867b6e2003SPierre Jolivet } 1877b6e2003SPierre Jolivet 188d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_LU(PC pc, Vec x, Vec y) 189d71ae5a4SJacob Faibussowitsch { 1909b54502bSHong Zhang PC_LU *dir = (PC_LU *)pc->data; 1919b54502bSHong Zhang 1929b54502bSHong Zhang PetscFunctionBegin; 1933d1c1ea0SBarry Smith if (dir->hdr.inplace) { 1949566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(pc->pmat, x, y)); 1952fa5cd67SKarl Rupp } else { 1969566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(((PC_Factor *)dir)->fact, x, y)); 1972fa5cd67SKarl Rupp } 1983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1999b54502bSHong Zhang } 2009b54502bSHong Zhang 2019b54502bSHong Zhang /*MC 2029b54502bSHong Zhang PCLU - Uses a direct solver, based on LU factorization, as a preconditioner 2039b54502bSHong Zhang 2049b54502bSHong Zhang Options Database Keys: 205f1580f4eSBarry Smith + -pc_factor_reuse_ordering - Activate `PCFactorSetReuseOrdering()` 206f1580f4eSBarry Smith . -pc_factor_mat_solver_type - Actives `PCFactorSetMatSolverType()` to choose the direct solver, like superlu 207f1580f4eSBarry Smith . -pc_factor_reuse_fill - Activates `PCFactorSetReuseFill()` 20855ba2a51SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 2092401956bSBarry Smith . -pc_factor_in_place - Activates in-place factorization 2102401956bSBarry Smith . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 2112401956bSBarry Smith . -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase 2129b54502bSHong Zhang stability of factorization. 213f1580f4eSBarry Smith . -pc_factor_shift_type <shifttype> - Sets shift type or -1 for the default; use '-help' for a list of available types 214f1580f4eSBarry Smith . -pc_factor_shift_amount <shiftamount> - Sets shift amount or -1 for the default 21526cc229bSBarry Smith . -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the diagonal. 216f1580f4eSBarry Smith . -pc_factor_mat_solver_type <packagename> - use an external package for the solve, see `MatSolverType` for possibilities 21726cc229bSBarry Smith - -mat_solvertype_optionname - options for a specific solver package, for example -mat_mumps_cntl_1 2189b54502bSHong Zhang 2199b54502bSHong Zhang Level: beginner 2209b54502bSHong Zhang 22195452b02SPatrick Sanan Notes: 222f1580f4eSBarry Smith Not all options work for all matrix formats 223f1580f4eSBarry Smith 224f1580f4eSBarry Smith Run with -help to see additional options for particular matrix formats or factorization algorithms 225f1580f4eSBarry Smith 22695452b02SPatrick Sanan Usually this will compute an "exact" solution in one iteration and does 2279b54502bSHong Zhang not need a Krylov method (i.e. you can use -ksp_type preonly, or 228f1580f4eSBarry Smith `KSPSetType`(ksp,`KSPPREONLY`) for the Krylov method 2299b54502bSHong Zhang 230f1580f4eSBarry Smith .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `MatSolverType`, `MatGetFactor()`, `PCQR`, `PCSVD`, 231db781477SPatrick Sanan `PCILU`, `PCCHOLESKY`, `PCICC`, `PCFactorSetReuseOrdering()`, `PCFactorSetReuseFill()`, `PCFactorGetMatrix()`, 232db781477SPatrick Sanan `PCFactorSetFill()`, `PCFactorSetUseInPlace()`, `PCFactorSetMatOrderingType()`, `PCFactorSetColumnPivot()`, 233c2e3fba1SPatrick Sanan `PCFactorSetPivotInBlocks()`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()` 234db781477SPatrick Sanan `PCFactorReorderForNonzeroDiagonal()` 2359b54502bSHong Zhang M*/ 2369b54502bSHong Zhang 237d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_LU(PC pc) 238d71ae5a4SJacob Faibussowitsch { 2399b54502bSHong Zhang PC_LU *dir; 2409b54502bSHong Zhang 2419b54502bSHong Zhang PetscFunctionBegin; 2424dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&dir)); 2433d1c1ea0SBarry Smith pc->data = (void *)dir; 2449566063dSJacob Faibussowitsch PetscCall(PCFactorInitialize(pc, MAT_FACTOR_LU)); 2459b54502bSHong Zhang dir->nonzerosalongdiagonal = PETSC_FALSE; 2469b54502bSHong Zhang 247075768bcSBarry Smith ((PC_Factor *)dir)->info.fill = 5.0; 248075768bcSBarry Smith ((PC_Factor *)dir)->info.dtcol = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */ 249f4db908eSBarry Smith ((PC_Factor *)dir)->info.shifttype = (PetscReal)MAT_SHIFT_NONE; 2500a545947SLisandro Dalcin dir->col = NULL; 2510a545947SLisandro Dalcin dir->row = NULL; 2525c9eb25fSBarry Smith 253574deadeSBarry Smith pc->ops->reset = PCReset_LU; 2549b54502bSHong Zhang pc->ops->destroy = PCDestroy_LU; 2559b54502bSHong Zhang pc->ops->apply = PCApply_LU; 2567b6e2003SPierre Jolivet pc->ops->matapply = PCMatApply_LU; 2579b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_LU; 2589b54502bSHong Zhang pc->ops->setup = PCSetUp_LU; 2599b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_LU; 26092e08861SBarry Smith pc->ops->view = PCView_Factor; 2610a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 2629566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorReorderForNonzeroDiagonal_C", PCFactorReorderForNonzeroDiagonal_LU)); 2633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2649b54502bSHong Zhang } 265