xref: /petsc/src/ksp/pc/impls/factor/lu/lu.c (revision dbbe0bcd3f3a8fbab5a45420dc06f8387e5764c6)
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 
21*dbbe0bcdSBarry Smith static PetscErrorCode PCSetFromOptions_LU(PC pc,PetscOptionItems *PetscOptionsObject)
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;
28d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"LU options");
29*dbbe0bcdSBarry 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();
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;
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));
6603c60df9SBarry Smith       if (dir->row) {
679566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row));
689566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col));
6903c60df9SBarry Smith       }
709566063dSJacob Faibussowitsch       PetscCall(MatLUFactor(pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info));
719566063dSJacob Faibussowitsch       PetscCall(MatFactorGetError(pc->pmat,&err));
7200e125f8SBarry Smith       if (err) { /* Factor() fails */
7300e125f8SBarry Smith         pc->failedreason = (PCFailedReason)err;
746baea169SHong Zhang         PetscFunctionReturn(0);
756baea169SHong Zhang       }
769d76b4d0SMatthew G. Knepley     }
77075768bcSBarry Smith     ((PC_Factor*)dir)->fact = pc->pmat;
789b54502bSHong Zhang   } else {
799b54502bSHong Zhang     MatInfo info;
8000e125f8SBarry Smith 
819b54502bSHong Zhang     if (!pc->setupcalled) {
82f73b0415SBarry Smith       PetscBool canuseordering;
832c7c0729SBarry Smith       if (!((PC_Factor*)dir)->fact) {
849566063dSJacob Faibussowitsch         PetscCall(MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact));
859566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact));
862c7c0729SBarry Smith       }
879566063dSJacob Faibussowitsch       PetscCall(MatFactorGetCanUseOrdering(((PC_Factor*)dir)->fact,&canuseordering));
88f73b0415SBarry Smith       if (canuseordering) {
899566063dSJacob Faibussowitsch         PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
909566063dSJacob Faibussowitsch         PetscCall(MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col));
911baa6e33SBarry Smith         if (dir->nonzerosalongdiagonal) PetscCall(MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col));
929566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row));
939566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col));
9403c60df9SBarry Smith       }
959566063dSJacob Faibussowitsch       PetscCall(MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info));
969566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info));
973d1c1ea0SBarry Smith       dir->hdr.actualfill = info.fill_ratio_needed;
989b54502bSHong Zhang     } else if (pc->flag != SAME_NONZERO_PATTERN) {
99f73b0415SBarry Smith       PetscBool canuseordering;
1003d1c1ea0SBarry Smith       if (!dir->hdr.reuseordering) {
1019566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&((PC_Factor*)dir)->fact));
1029566063dSJacob Faibussowitsch         PetscCall(MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact));
1039566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact));
1049566063dSJacob Faibussowitsch         PetscCall(MatFactorGetCanUseOrdering(((PC_Factor*)dir)->fact,&canuseordering));
105f73b0415SBarry Smith         if (canuseordering) {
1069566063dSJacob Faibussowitsch           if (dir->row && dir->col && dir->row != dir->col) PetscCall(ISDestroy(&dir->row));
1079566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&dir->col));
1089566063dSJacob Faibussowitsch           PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
1099566063dSJacob Faibussowitsch           PetscCall(MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col));
1101baa6e33SBarry Smith           if (dir->nonzerosalongdiagonal) PetscCall(MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col));
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));
1672e956fe4SStefano Zampini   PetscCall(PCFactorClearComposedFunctions(pc));
1689566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
1699b54502bSHong Zhang   PetscFunctionReturn(0);
1709b54502bSHong Zhang }
1719b54502bSHong Zhang 
1729b54502bSHong Zhang static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y)
1739b54502bSHong Zhang {
1749b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
1759b54502bSHong Zhang 
1769b54502bSHong Zhang   PetscFunctionBegin;
1773d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
1789566063dSJacob Faibussowitsch     PetscCall(MatSolve(pc->pmat,x,y));
1792fa5cd67SKarl Rupp   } else {
1809566063dSJacob Faibussowitsch     PetscCall(MatSolve(((PC_Factor*)dir)->fact,x,y));
1812fa5cd67SKarl Rupp   }
1829b54502bSHong Zhang   PetscFunctionReturn(0);
1839b54502bSHong Zhang }
1849b54502bSHong Zhang 
1857b6e2003SPierre Jolivet static PetscErrorCode PCMatApply_LU(PC pc,Mat X,Mat Y)
1867b6e2003SPierre Jolivet {
1877b6e2003SPierre Jolivet   PC_LU          *dir = (PC_LU*)pc->data;
1887b6e2003SPierre Jolivet 
1897b6e2003SPierre Jolivet   PetscFunctionBegin;
1907b6e2003SPierre Jolivet   if (dir->hdr.inplace) {
1919566063dSJacob Faibussowitsch     PetscCall(MatMatSolve(pc->pmat,X,Y));
1927b6e2003SPierre Jolivet   } else {
1939566063dSJacob Faibussowitsch     PetscCall(MatMatSolve(((PC_Factor*)dir)->fact,X,Y));
1947b6e2003SPierre Jolivet   }
1957b6e2003SPierre Jolivet   PetscFunctionReturn(0);
1967b6e2003SPierre Jolivet }
1977b6e2003SPierre Jolivet 
1989b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y)
1999b54502bSHong Zhang {
2009b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
2019b54502bSHong Zhang 
2029b54502bSHong Zhang   PetscFunctionBegin;
2033d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
2049566063dSJacob Faibussowitsch     PetscCall(MatSolveTranspose(pc->pmat,x,y));
2052fa5cd67SKarl Rupp   } else {
2069566063dSJacob Faibussowitsch     PetscCall(MatSolveTranspose(((PC_Factor*)dir)->fact,x,y));
2072fa5cd67SKarl Rupp   }
2089b54502bSHong Zhang   PetscFunctionReturn(0);
2099b54502bSHong Zhang }
2109b54502bSHong Zhang 
2119b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/
2129b54502bSHong Zhang 
2139b54502bSHong Zhang /*MC
2149b54502bSHong Zhang    PCLU - Uses a direct solver, based on LU factorization, as a preconditioner
2159b54502bSHong Zhang 
2169b54502bSHong Zhang    Options Database Keys:
2172401956bSBarry Smith +  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
2183ca39a21SBarry Smith .  -pc_factor_mat_solver_type - Actives PCFactorSetMatSolverType() to choose the direct solver, like superlu
2192401956bSBarry Smith .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
22055ba2a51SBarry Smith .  -pc_factor_fill <fill> - Sets fill amount
2212401956bSBarry Smith .  -pc_factor_in_place - Activates in-place factorization
2222401956bSBarry Smith .  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
2232401956bSBarry Smith .  -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase
2249b54502bSHong Zhang                                          stability of factorization.
225145b9266SHong Zhang .  -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types
226145b9266SHong Zhang .  -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default
22726cc229bSBarry Smith .  -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the diagonal.
22826cc229bSBarry Smith -  -mat_solvertype_optionname - options for a specific solver package, for example -mat_mumps_cntl_1
2299b54502bSHong Zhang 
23095452b02SPatrick Sanan    Notes:
23195452b02SPatrick Sanan     Not all options work for all matrix formats
2329b54502bSHong Zhang           Run with -help to see additional options for particular matrix formats or factorization
2339b54502bSHong Zhang           algorithms
2349b54502bSHong Zhang 
2359b54502bSHong Zhang    Level: beginner
2369b54502bSHong Zhang 
23795452b02SPatrick Sanan    Notes:
23895452b02SPatrick Sanan     Usually this will compute an "exact" solution in one iteration and does
2399b54502bSHong Zhang           not need a Krylov method (i.e. you can use -ksp_type preonly, or
2409b54502bSHong Zhang           KSPSetType(ksp,KSPPREONLY) for the Krylov method
2419b54502bSHong Zhang 
242db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`,
243db781477SPatrick Sanan           `PCILU`, `PCCHOLESKY`, `PCICC`, `PCFactorSetReuseOrdering()`, `PCFactorSetReuseFill()`, `PCFactorGetMatrix()`,
244db781477SPatrick Sanan           `PCFactorSetFill()`, `PCFactorSetUseInPlace()`, `PCFactorSetMatOrderingType()`, `PCFactorSetColumnPivot()`,
245c2e3fba1SPatrick Sanan           `PCFactorSetPivotInBlocks()`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()`
246db781477SPatrick Sanan           `PCFactorReorderForNonzeroDiagonal()`
2479b54502bSHong Zhang M*/
2489b54502bSHong Zhang 
2498cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_LU(PC pc)
2509b54502bSHong Zhang {
2519b54502bSHong Zhang   PC_LU          *dir;
2529b54502bSHong Zhang 
2539b54502bSHong Zhang   PetscFunctionBegin;
2549566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc,&dir));
2553d1c1ea0SBarry Smith   pc->data = (void*)dir;
2569566063dSJacob Faibussowitsch   PetscCall(PCFactorInitialize(pc,MAT_FACTOR_LU));
2579b54502bSHong Zhang   dir->nonzerosalongdiagonal    = PETSC_FALSE;
2589b54502bSHong Zhang 
259075768bcSBarry Smith   ((PC_Factor*)dir)->info.fill          = 5.0;
260075768bcSBarry Smith   ((PC_Factor*)dir)->info.dtcol         = 1.e-6;  /* default to pivoting; this is only thing PETSc LU supports */
261f4db908eSBarry Smith   ((PC_Factor*)dir)->info.shifttype     = (PetscReal)MAT_SHIFT_NONE;
2620a545947SLisandro Dalcin   dir->col                              = NULL;
2630a545947SLisandro Dalcin   dir->row                              = NULL;
2645c9eb25fSBarry Smith 
265574deadeSBarry Smith   pc->ops->reset             = PCReset_LU;
2669b54502bSHong Zhang   pc->ops->destroy           = PCDestroy_LU;
2679b54502bSHong Zhang   pc->ops->apply             = PCApply_LU;
2687b6e2003SPierre Jolivet   pc->ops->matapply          = PCMatApply_LU;
2699b54502bSHong Zhang   pc->ops->applytranspose    = PCApplyTranspose_LU;
2709b54502bSHong Zhang   pc->ops->setup             = PCSetUp_LU;
2719b54502bSHong Zhang   pc->ops->setfromoptions    = PCSetFromOptions_LU;
27292e08861SBarry Smith   pc->ops->view              = PCView_Factor;
2730a545947SLisandro Dalcin   pc->ops->applyrichardson   = NULL;
2749566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",PCFactorReorderForNonzeroDiagonal_LU));
2759b54502bSHong Zhang   PetscFunctionReturn(0);
2769b54502bSHong Zhang }
277