xref: /petsc/src/ksp/pc/impls/factor/lu/lu.c (revision c2e3fba1fe1cda7e6350bbca19c4ed35ce95940a)
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;
28d0609cedSBarry 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   }
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;
49f023e1d5SPierre Jolivet   PetscCall(PCGetOptionsPrefix(pc,&prefix));
50f023e1d5SPierre Jolivet   PetscCall(MatSetOptionsPrefix(pc->pmat,prefix));
51c6e4fdc6SHong Zhang   pc->failedreason = PC_NOERROR;
523d1c1ea0SBarry Smith   if (dir->hdr.reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->hdr.actualfill;
539b54502bSHong Zhang 
549566063dSJacob Faibussowitsch   PetscCall(MatSetErrorIfFailure(pc->pmat,pc->erroriffailure));
553d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
569d76b4d0SMatthew G. Knepley     MatFactorType ftype;
579d76b4d0SMatthew G. Knepley 
589566063dSJacob Faibussowitsch     PetscCall(MatGetFactorType(pc->pmat, &ftype));
599d76b4d0SMatthew G. Knepley     if (ftype == MAT_FACTOR_NONE) {
609566063dSJacob Faibussowitsch       if (dir->row && dir->col && dir->row != dir->col) PetscCall(ISDestroy(&dir->row));
619566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&dir->col));
622c7c0729SBarry Smith       /* This should only get the ordering if needed, but since MatGetFactor() is not called we can't know if it is needed */
639566063dSJacob Faibussowitsch       PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
649566063dSJacob Faibussowitsch       PetscCall(MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col));
6503c60df9SBarry Smith       if (dir->row) {
669566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row));
679566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col));
6803c60df9SBarry Smith       }
699566063dSJacob Faibussowitsch       PetscCall(MatLUFactor(pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info));
709566063dSJacob Faibussowitsch       PetscCall(MatFactorGetError(pc->pmat,&err));
7100e125f8SBarry Smith       if (err) { /* Factor() fails */
7200e125f8SBarry Smith         pc->failedreason = (PCFailedReason)err;
736baea169SHong Zhang         PetscFunctionReturn(0);
746baea169SHong Zhang       }
759d76b4d0SMatthew G. Knepley     }
76075768bcSBarry Smith     ((PC_Factor*)dir)->fact = pc->pmat;
779b54502bSHong Zhang   } else {
789b54502bSHong Zhang     MatInfo info;
7900e125f8SBarry Smith 
809b54502bSHong Zhang     if (!pc->setupcalled) {
81f73b0415SBarry Smith       PetscBool canuseordering;
822c7c0729SBarry Smith       if (!((PC_Factor*)dir)->fact) {
839566063dSJacob Faibussowitsch         PetscCall(MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact));
849566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact));
852c7c0729SBarry Smith       }
869566063dSJacob Faibussowitsch       PetscCall(MatFactorGetCanUseOrdering(((PC_Factor*)dir)->fact,&canuseordering));
87f73b0415SBarry Smith       if (canuseordering) {
889566063dSJacob Faibussowitsch         PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
899566063dSJacob Faibussowitsch         PetscCall(MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col));
909b54502bSHong Zhang         if (dir->nonzerosalongdiagonal) {
919566063dSJacob Faibussowitsch           PetscCall(MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col));
929b54502bSHong Zhang         }
939566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row));
949566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col));
9503c60df9SBarry Smith       }
969566063dSJacob Faibussowitsch       PetscCall(MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info));
979566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info));
983d1c1ea0SBarry Smith       dir->hdr.actualfill = info.fill_ratio_needed;
999b54502bSHong Zhang     } else if (pc->flag != SAME_NONZERO_PATTERN) {
100f73b0415SBarry Smith       PetscBool canuseordering;
1013d1c1ea0SBarry Smith       if (!dir->hdr.reuseordering) {
1029566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&((PC_Factor*)dir)->fact));
1039566063dSJacob Faibussowitsch         PetscCall(MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact));
1049566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact));
1059566063dSJacob Faibussowitsch         PetscCall(MatFactorGetCanUseOrdering(((PC_Factor*)dir)->fact,&canuseordering));
106f73b0415SBarry Smith         if (canuseordering) {
1079566063dSJacob Faibussowitsch           if (dir->row && dir->col && dir->row != dir->col) PetscCall(ISDestroy(&dir->row));
1089566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&dir->col));
1099566063dSJacob Faibussowitsch           PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
1109566063dSJacob Faibussowitsch           PetscCall(MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col));
1119b54502bSHong Zhang           if (dir->nonzerosalongdiagonal) {
1129566063dSJacob Faibussowitsch             PetscCall(MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col));
1139b54502bSHong Zhang           }
1149566063dSJacob Faibussowitsch           PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row));
1159566063dSJacob Faibussowitsch           PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col));
11603c60df9SBarry Smith         }
1179b54502bSHong Zhang       }
1189566063dSJacob Faibussowitsch       PetscCall(MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info));
1199566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info));
1203d1c1ea0SBarry Smith       dir->hdr.actualfill = info.fill_ratio_needed;
12104545d6dSBarry Smith     } else {
1229566063dSJacob Faibussowitsch       PetscCall(MatFactorGetError(((PC_Factor*)dir)->fact,&err));
123160a8794SBarry Smith       if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) {
1249566063dSJacob Faibussowitsch         PetscCall(MatFactorClearError(((PC_Factor*)dir)->fact));
125b8b68cfdSBarry Smith         pc->failedreason = PC_NOERROR;
12604545d6dSBarry Smith       }
1279b54502bSHong Zhang     }
1289566063dSJacob Faibussowitsch     PetscCall(MatFactorGetError(((PC_Factor*)dir)->fact,&err));
12900e125f8SBarry Smith     if (err) { /* FactorSymbolic() fails */
13000e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1318c1cd74cSHong Zhang       PetscFunctionReturn(0);
1328c1cd74cSHong Zhang     }
1338c1cd74cSHong Zhang 
1349566063dSJacob Faibussowitsch     PetscCall(MatLUFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info));
1359566063dSJacob Faibussowitsch     PetscCall(MatFactorGetError(((PC_Factor*)dir)->fact,&err));
13600e125f8SBarry Smith     if (err) { /* FactorNumeric() fails */
13700e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1388c1cd74cSHong Zhang     }
139680c5173SHong Zhang 
1409b54502bSHong Zhang   }
14100c67f3bSHong Zhang 
1429566063dSJacob Faibussowitsch   PetscCall(PCFactorGetMatSolverType(pc,&stype));
14300c67f3bSHong Zhang   if (!stype) {
144ea799195SBarry Smith     MatSolverType solverpackage;
1459566063dSJacob Faibussowitsch     PetscCall(MatFactorGetSolverType(((PC_Factor*)dir)->fact,&solverpackage));
1469566063dSJacob Faibussowitsch     PetscCall(PCFactorSetMatSolverType(pc,solverpackage));
14700c67f3bSHong Zhang   }
1489b54502bSHong Zhang   PetscFunctionReturn(0);
1499b54502bSHong Zhang }
1509b54502bSHong Zhang 
151574deadeSBarry Smith static PetscErrorCode PCReset_LU(PC pc)
1529b54502bSHong Zhang {
1539b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
1549b54502bSHong Zhang 
1559b54502bSHong Zhang   PetscFunctionBegin;
1569566063dSJacob Faibussowitsch   if (!dir->hdr.inplace && ((PC_Factor*)dir)->fact) PetscCall(MatDestroy(&((PC_Factor*)dir)->fact));
1579566063dSJacob Faibussowitsch   if (dir->row && dir->col && dir->row != dir->col) PetscCall(ISDestroy(&dir->row));
1589566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&dir->col));
159574deadeSBarry Smith   PetscFunctionReturn(0);
160574deadeSBarry Smith }
161574deadeSBarry Smith 
162574deadeSBarry Smith static PetscErrorCode PCDestroy_LU(PC pc)
163574deadeSBarry Smith {
164574deadeSBarry Smith   PC_LU          *dir = (PC_LU*)pc->data;
165574deadeSBarry Smith 
166574deadeSBarry Smith   PetscFunctionBegin;
1679566063dSJacob Faibussowitsch   PetscCall(PCReset_LU(pc));
1689566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PC_Factor*)dir)->ordering));
1699566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PC_Factor*)dir)->solvertype));
1709566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
1719b54502bSHong Zhang   PetscFunctionReturn(0);
1729b54502bSHong Zhang }
1739b54502bSHong Zhang 
1749b54502bSHong Zhang static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y)
1759b54502bSHong Zhang {
1769b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
1779b54502bSHong Zhang 
1789b54502bSHong Zhang   PetscFunctionBegin;
1793d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
1809566063dSJacob Faibussowitsch     PetscCall(MatSolve(pc->pmat,x,y));
1812fa5cd67SKarl Rupp   } else {
1829566063dSJacob Faibussowitsch     PetscCall(MatSolve(((PC_Factor*)dir)->fact,x,y));
1832fa5cd67SKarl Rupp   }
1849b54502bSHong Zhang   PetscFunctionReturn(0);
1859b54502bSHong Zhang }
1869b54502bSHong Zhang 
1877b6e2003SPierre Jolivet static PetscErrorCode PCMatApply_LU(PC pc,Mat X,Mat Y)
1887b6e2003SPierre Jolivet {
1897b6e2003SPierre Jolivet   PC_LU          *dir = (PC_LU*)pc->data;
1907b6e2003SPierre Jolivet 
1917b6e2003SPierre Jolivet   PetscFunctionBegin;
1927b6e2003SPierre Jolivet   if (dir->hdr.inplace) {
1939566063dSJacob Faibussowitsch     PetscCall(MatMatSolve(pc->pmat,X,Y));
1947b6e2003SPierre Jolivet   } else {
1959566063dSJacob Faibussowitsch     PetscCall(MatMatSolve(((PC_Factor*)dir)->fact,X,Y));
1967b6e2003SPierre Jolivet   }
1977b6e2003SPierre Jolivet   PetscFunctionReturn(0);
1987b6e2003SPierre Jolivet }
1997b6e2003SPierre Jolivet 
2009b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y)
2019b54502bSHong Zhang {
2029b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
2039b54502bSHong Zhang 
2049b54502bSHong Zhang   PetscFunctionBegin;
2053d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
2069566063dSJacob Faibussowitsch     PetscCall(MatSolveTranspose(pc->pmat,x,y));
2072fa5cd67SKarl Rupp   } else {
2089566063dSJacob Faibussowitsch     PetscCall(MatSolveTranspose(((PC_Factor*)dir)->fact,x,y));
2092fa5cd67SKarl Rupp   }
2109b54502bSHong Zhang   PetscFunctionReturn(0);
2119b54502bSHong Zhang }
2129b54502bSHong Zhang 
2139b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/
2149b54502bSHong Zhang 
2159b54502bSHong Zhang /*MC
2169b54502bSHong Zhang    PCLU - Uses a direct solver, based on LU factorization, as a preconditioner
2179b54502bSHong Zhang 
2189b54502bSHong Zhang    Options Database Keys:
2192401956bSBarry Smith +  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
2203ca39a21SBarry Smith .  -pc_factor_mat_solver_type - Actives PCFactorSetMatSolverType() to choose the direct solver, like superlu
2212401956bSBarry Smith .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
22255ba2a51SBarry Smith .  -pc_factor_fill <fill> - Sets fill amount
2232401956bSBarry Smith .  -pc_factor_in_place - Activates in-place factorization
2242401956bSBarry Smith .  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
2252401956bSBarry Smith .  -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase
2269b54502bSHong Zhang                                          stability of factorization.
227145b9266SHong Zhang .  -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types
228145b9266SHong Zhang .  -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default
229e22d95b2SBarry Smith -   -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the
230e22d95b2SBarry Smith         diagonal.
2319b54502bSHong Zhang 
23295452b02SPatrick Sanan    Notes:
23395452b02SPatrick Sanan     Not all options work for all matrix formats
2349b54502bSHong Zhang           Run with -help to see additional options for particular matrix formats or factorization
2359b54502bSHong Zhang           algorithms
2369b54502bSHong Zhang 
2379b54502bSHong Zhang    Level: beginner
2389b54502bSHong Zhang 
23995452b02SPatrick Sanan    Notes:
24095452b02SPatrick Sanan     Usually this will compute an "exact" solution in one iteration and does
2419b54502bSHong Zhang           not need a Krylov method (i.e. you can use -ksp_type preonly, or
2429b54502bSHong Zhang           KSPSetType(ksp,KSPPREONLY) for the Krylov method
2439b54502bSHong Zhang 
244db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`,
245db781477SPatrick Sanan           `PCILU`, `PCCHOLESKY`, `PCICC`, `PCFactorSetReuseOrdering()`, `PCFactorSetReuseFill()`, `PCFactorGetMatrix()`,
246db781477SPatrick Sanan           `PCFactorSetFill()`, `PCFactorSetUseInPlace()`, `PCFactorSetMatOrderingType()`, `PCFactorSetColumnPivot()`,
247*c2e3fba1SPatrick Sanan           `PCFactorSetPivotInBlocks()`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()`
248db781477SPatrick Sanan           `PCFactorReorderForNonzeroDiagonal()`
2499b54502bSHong Zhang M*/
2509b54502bSHong Zhang 
2518cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_LU(PC pc)
2529b54502bSHong Zhang {
2539b54502bSHong Zhang   PC_LU          *dir;
2549b54502bSHong Zhang 
2559b54502bSHong Zhang   PetscFunctionBegin;
2569566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc,&dir));
2573d1c1ea0SBarry Smith   pc->data = (void*)dir;
2589566063dSJacob Faibussowitsch   PetscCall(PCFactorInitialize(pc,MAT_FACTOR_LU));
2599b54502bSHong Zhang   dir->nonzerosalongdiagonal    = PETSC_FALSE;
2609b54502bSHong Zhang 
261075768bcSBarry Smith   ((PC_Factor*)dir)->info.fill          = 5.0;
262075768bcSBarry Smith   ((PC_Factor*)dir)->info.dtcol         = 1.e-6;  /* default to pivoting; this is only thing PETSc LU supports */
263f4db908eSBarry Smith   ((PC_Factor*)dir)->info.shifttype     = (PetscReal)MAT_SHIFT_NONE;
2640a545947SLisandro Dalcin   dir->col                              = NULL;
2650a545947SLisandro Dalcin   dir->row                              = NULL;
2665c9eb25fSBarry Smith 
267574deadeSBarry Smith   pc->ops->reset             = PCReset_LU;
2689b54502bSHong Zhang   pc->ops->destroy           = PCDestroy_LU;
2699b54502bSHong Zhang   pc->ops->apply             = PCApply_LU;
2707b6e2003SPierre Jolivet   pc->ops->matapply          = PCMatApply_LU;
2719b54502bSHong Zhang   pc->ops->applytranspose    = PCApplyTranspose_LU;
2729b54502bSHong Zhang   pc->ops->setup             = PCSetUp_LU;
2739b54502bSHong Zhang   pc->ops->setfromoptions    = PCSetFromOptions_LU;
27492e08861SBarry Smith   pc->ops->view              = PCView_Factor;
2750a545947SLisandro Dalcin   pc->ops->applyrichardson   = NULL;
2769566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",PCFactorReorderForNonzeroDiagonal_LU));
2779b54502bSHong Zhang   PetscFunctionReturn(0);
2789b54502bSHong Zhang }
279