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; 162fa5cd67SKarl Rupp if (z == PETSC_DECIDE) lu->nonzerosalongdiagonaltol = 1.e-10; 172fa5cd67SKarl Rupp else lu->nonzerosalongdiagonaltol = z; 18*3ba16761SJacob 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(); 38*3ba16761SJacob 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; 71*3ba16761SJacob 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; 804dfa11a4SJacob Faibussowitsch if (!((PC_Factor *)dir)->fact) { PetscCall(MatGetFactor(pc->pmat, ((PC_Factor *)dir)->solvertype, MAT_FACTOR_LU, &((PC_Factor *)dir)->fact)); } 819566063dSJacob Faibussowitsch PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)dir)->fact, &canuseordering)); 82f73b0415SBarry Smith if (canuseordering) { 839566063dSJacob Faibussowitsch PetscCall(PCFactorSetDefaultOrdering_Factor(pc)); 849566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col)); 851baa6e33SBarry Smith if (dir->nonzerosalongdiagonal) PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, dir->nonzerosalongdiagonaltol, dir->row, dir->col)); 8603c60df9SBarry Smith } 879566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->row, dir->col, &((PC_Factor *)dir)->info)); 889566063dSJacob Faibussowitsch PetscCall(MatGetInfo(((PC_Factor *)dir)->fact, MAT_LOCAL, &info)); 893d1c1ea0SBarry Smith dir->hdr.actualfill = info.fill_ratio_needed; 909b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 91f73b0415SBarry Smith PetscBool canuseordering; 923d1c1ea0SBarry Smith if (!dir->hdr.reuseordering) { 939566063dSJacob Faibussowitsch PetscCall(MatDestroy(&((PC_Factor *)dir)->fact)); 949566063dSJacob Faibussowitsch PetscCall(MatGetFactor(pc->pmat, ((PC_Factor *)dir)->solvertype, MAT_FACTOR_LU, &((PC_Factor *)dir)->fact)); 959566063dSJacob Faibussowitsch PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)dir)->fact, &canuseordering)); 96f73b0415SBarry Smith if (canuseordering) { 979566063dSJacob Faibussowitsch if (dir->row && dir->col && dir->row != dir->col) PetscCall(ISDestroy(&dir->row)); 989566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->col)); 999566063dSJacob Faibussowitsch PetscCall(PCFactorSetDefaultOrdering_Factor(pc)); 1009566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col)); 1011baa6e33SBarry Smith if (dir->nonzerosalongdiagonal) PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, dir->nonzerosalongdiagonaltol, dir->row, dir->col)); 10203c60df9SBarry Smith } 1039b54502bSHong Zhang } 1049566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->row, dir->col, &((PC_Factor *)dir)->info)); 1059566063dSJacob Faibussowitsch PetscCall(MatGetInfo(((PC_Factor *)dir)->fact, MAT_LOCAL, &info)); 1063d1c1ea0SBarry Smith dir->hdr.actualfill = info.fill_ratio_needed; 10704545d6dSBarry Smith } else { 1089566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err)); 109160a8794SBarry Smith if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) { 1109566063dSJacob Faibussowitsch PetscCall(MatFactorClearError(((PC_Factor *)dir)->fact)); 111b8b68cfdSBarry Smith pc->failedreason = PC_NOERROR; 11204545d6dSBarry Smith } 1139b54502bSHong Zhang } 1149566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err)); 11500e125f8SBarry Smith if (err) { /* FactorSymbolic() fails */ 11600e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 117*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1188c1cd74cSHong Zhang } 1198c1cd74cSHong Zhang 1209566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(((PC_Factor *)dir)->fact, pc->pmat, &((PC_Factor *)dir)->info)); 1219566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err)); 12200e125f8SBarry Smith if (err) { /* FactorNumeric() fails */ 12300e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1248c1cd74cSHong Zhang } 1259b54502bSHong Zhang } 12600c67f3bSHong Zhang 1279566063dSJacob Faibussowitsch PetscCall(PCFactorGetMatSolverType(pc, &stype)); 12800c67f3bSHong Zhang if (!stype) { 129ea799195SBarry Smith MatSolverType solverpackage; 1309566063dSJacob Faibussowitsch PetscCall(MatFactorGetSolverType(((PC_Factor *)dir)->fact, &solverpackage)); 1319566063dSJacob Faibussowitsch PetscCall(PCFactorSetMatSolverType(pc, solverpackage)); 13200c67f3bSHong Zhang } 133*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1349b54502bSHong Zhang } 1359b54502bSHong Zhang 136d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_LU(PC pc) 137d71ae5a4SJacob Faibussowitsch { 1389b54502bSHong Zhang PC_LU *dir = (PC_LU *)pc->data; 1399b54502bSHong Zhang 1409b54502bSHong Zhang PetscFunctionBegin; 1419566063dSJacob Faibussowitsch if (!dir->hdr.inplace && ((PC_Factor *)dir)->fact) PetscCall(MatDestroy(&((PC_Factor *)dir)->fact)); 1429566063dSJacob Faibussowitsch if (dir->row && dir->col && dir->row != dir->col) PetscCall(ISDestroy(&dir->row)); 1439566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->col)); 144*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 145574deadeSBarry Smith } 146574deadeSBarry Smith 147d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_LU(PC pc) 148d71ae5a4SJacob Faibussowitsch { 149574deadeSBarry Smith PC_LU *dir = (PC_LU *)pc->data; 150574deadeSBarry Smith 151574deadeSBarry Smith PetscFunctionBegin; 1529566063dSJacob Faibussowitsch PetscCall(PCReset_LU(pc)); 1539566063dSJacob Faibussowitsch PetscCall(PetscFree(((PC_Factor *)dir)->ordering)); 1549566063dSJacob Faibussowitsch PetscCall(PetscFree(((PC_Factor *)dir)->solvertype)); 1552e956fe4SStefano Zampini PetscCall(PCFactorClearComposedFunctions(pc)); 1569566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 157*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1589b54502bSHong Zhang } 1599b54502bSHong Zhang 160d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_LU(PC pc, Vec x, Vec y) 161d71ae5a4SJacob Faibussowitsch { 1629b54502bSHong Zhang PC_LU *dir = (PC_LU *)pc->data; 1639b54502bSHong Zhang 1649b54502bSHong Zhang PetscFunctionBegin; 1653d1c1ea0SBarry Smith if (dir->hdr.inplace) { 1669566063dSJacob Faibussowitsch PetscCall(MatSolve(pc->pmat, x, y)); 1672fa5cd67SKarl Rupp } else { 1689566063dSJacob Faibussowitsch PetscCall(MatSolve(((PC_Factor *)dir)->fact, x, y)); 1692fa5cd67SKarl Rupp } 170*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1719b54502bSHong Zhang } 1729b54502bSHong Zhang 173d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_LU(PC pc, Mat X, Mat Y) 174d71ae5a4SJacob Faibussowitsch { 1757b6e2003SPierre Jolivet PC_LU *dir = (PC_LU *)pc->data; 1767b6e2003SPierre Jolivet 1777b6e2003SPierre Jolivet PetscFunctionBegin; 1787b6e2003SPierre Jolivet if (dir->hdr.inplace) { 1799566063dSJacob Faibussowitsch PetscCall(MatMatSolve(pc->pmat, X, Y)); 1807b6e2003SPierre Jolivet } else { 1819566063dSJacob Faibussowitsch PetscCall(MatMatSolve(((PC_Factor *)dir)->fact, X, Y)); 1827b6e2003SPierre Jolivet } 183*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1847b6e2003SPierre Jolivet } 1857b6e2003SPierre Jolivet 186d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_LU(PC pc, Vec x, Vec y) 187d71ae5a4SJacob Faibussowitsch { 1889b54502bSHong Zhang PC_LU *dir = (PC_LU *)pc->data; 1899b54502bSHong Zhang 1909b54502bSHong Zhang PetscFunctionBegin; 1913d1c1ea0SBarry Smith if (dir->hdr.inplace) { 1929566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(pc->pmat, x, y)); 1932fa5cd67SKarl Rupp } else { 1949566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(((PC_Factor *)dir)->fact, x, y)); 1952fa5cd67SKarl Rupp } 196*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1979b54502bSHong Zhang } 1989b54502bSHong Zhang 1999b54502bSHong Zhang /*MC 2009b54502bSHong Zhang PCLU - Uses a direct solver, based on LU factorization, as a preconditioner 2019b54502bSHong Zhang 2029b54502bSHong Zhang Options Database Keys: 203f1580f4eSBarry Smith + -pc_factor_reuse_ordering - Activate `PCFactorSetReuseOrdering()` 204f1580f4eSBarry Smith . -pc_factor_mat_solver_type - Actives `PCFactorSetMatSolverType()` to choose the direct solver, like superlu 205f1580f4eSBarry Smith . -pc_factor_reuse_fill - Activates `PCFactorSetReuseFill()` 20655ba2a51SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 2072401956bSBarry Smith . -pc_factor_in_place - Activates in-place factorization 2082401956bSBarry Smith . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 2092401956bSBarry Smith . -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase 2109b54502bSHong Zhang stability of factorization. 211f1580f4eSBarry Smith . -pc_factor_shift_type <shifttype> - Sets shift type or -1 for the default; use '-help' for a list of available types 212f1580f4eSBarry Smith . -pc_factor_shift_amount <shiftamount> - Sets shift amount or -1 for the default 21326cc229bSBarry Smith . -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the diagonal. 214f1580f4eSBarry Smith . -pc_factor_mat_solver_type <packagename> - use an external package for the solve, see `MatSolverType` for possibilities 21526cc229bSBarry Smith - -mat_solvertype_optionname - options for a specific solver package, for example -mat_mumps_cntl_1 2169b54502bSHong Zhang 2179b54502bSHong Zhang Level: beginner 2189b54502bSHong Zhang 21995452b02SPatrick Sanan Notes: 220f1580f4eSBarry Smith Not all options work for all matrix formats 221f1580f4eSBarry Smith 222f1580f4eSBarry Smith Run with -help to see additional options for particular matrix formats or factorization algorithms 223f1580f4eSBarry Smith 22495452b02SPatrick Sanan Usually this will compute an "exact" solution in one iteration and does 2259b54502bSHong Zhang not need a Krylov method (i.e. you can use -ksp_type preonly, or 226f1580f4eSBarry Smith `KSPSetType`(ksp,`KSPPREONLY`) for the Krylov method 2279b54502bSHong Zhang 228f1580f4eSBarry Smith .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `MatSolverType`, `MatGetFactor()`, `PCQR`, `PCSVD`, 229db781477SPatrick Sanan `PCILU`, `PCCHOLESKY`, `PCICC`, `PCFactorSetReuseOrdering()`, `PCFactorSetReuseFill()`, `PCFactorGetMatrix()`, 230db781477SPatrick Sanan `PCFactorSetFill()`, `PCFactorSetUseInPlace()`, `PCFactorSetMatOrderingType()`, `PCFactorSetColumnPivot()`, 231c2e3fba1SPatrick Sanan `PCFactorSetPivotInBlocks()`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()` 232db781477SPatrick Sanan `PCFactorReorderForNonzeroDiagonal()` 2339b54502bSHong Zhang M*/ 2349b54502bSHong Zhang 235d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_LU(PC pc) 236d71ae5a4SJacob Faibussowitsch { 2379b54502bSHong Zhang PC_LU *dir; 2389b54502bSHong Zhang 2399b54502bSHong Zhang PetscFunctionBegin; 2404dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&dir)); 2413d1c1ea0SBarry Smith pc->data = (void *)dir; 2429566063dSJacob Faibussowitsch PetscCall(PCFactorInitialize(pc, MAT_FACTOR_LU)); 2439b54502bSHong Zhang dir->nonzerosalongdiagonal = PETSC_FALSE; 2449b54502bSHong Zhang 245075768bcSBarry Smith ((PC_Factor *)dir)->info.fill = 5.0; 246075768bcSBarry Smith ((PC_Factor *)dir)->info.dtcol = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */ 247f4db908eSBarry Smith ((PC_Factor *)dir)->info.shifttype = (PetscReal)MAT_SHIFT_NONE; 2480a545947SLisandro Dalcin dir->col = NULL; 2490a545947SLisandro Dalcin dir->row = NULL; 2505c9eb25fSBarry Smith 251574deadeSBarry Smith pc->ops->reset = PCReset_LU; 2529b54502bSHong Zhang pc->ops->destroy = PCDestroy_LU; 2539b54502bSHong Zhang pc->ops->apply = PCApply_LU; 2547b6e2003SPierre Jolivet pc->ops->matapply = PCMatApply_LU; 2559b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_LU; 2569b54502bSHong Zhang pc->ops->setup = PCSetUp_LU; 2579b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_LU; 25892e08861SBarry Smith pc->ops->view = PCView_Factor; 2590a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 2609566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorReorderForNonzeroDiagonal_C", PCFactorReorderForNonzeroDiagonal_LU)); 261*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2629b54502bSHong Zhang } 263