xref: /petsc/src/ksp/pc/impls/factor/lu/lu.c (revision 2e956fe4fc852fabc23b437482e1fb7b77fddb0d)
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;
49a7f29c02SPierre Jolivet   if (!((PetscObject)pc->pmat)->prefix) {
50f023e1d5SPierre Jolivet     PetscCall(PCGetOptionsPrefix(pc,&prefix));
51f023e1d5SPierre Jolivet     PetscCall(MatSetOptionsPrefix(pc->pmat,prefix));
52a7f29c02SPierre Jolivet   }
53c6e4fdc6SHong Zhang   pc->failedreason = PC_NOERROR;
543d1c1ea0SBarry Smith   if (dir->hdr.reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->hdr.actualfill;
559b54502bSHong Zhang 
569566063dSJacob Faibussowitsch   PetscCall(MatSetErrorIfFailure(pc->pmat,pc->erroriffailure));
573d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
589d76b4d0SMatthew G. Knepley     MatFactorType ftype;
599d76b4d0SMatthew G. Knepley 
609566063dSJacob Faibussowitsch     PetscCall(MatGetFactorType(pc->pmat, &ftype));
619d76b4d0SMatthew G. Knepley     if (ftype == MAT_FACTOR_NONE) {
629566063dSJacob Faibussowitsch       if (dir->row && dir->col && dir->row != dir->col) PetscCall(ISDestroy(&dir->row));
639566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&dir->col));
642c7c0729SBarry Smith       /* This should only get the ordering if needed, but since MatGetFactor() is not called we can't know if it is needed */
659566063dSJacob Faibussowitsch       PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
669566063dSJacob Faibussowitsch       PetscCall(MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col));
6703c60df9SBarry Smith       if (dir->row) {
689566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row));
699566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col));
7003c60df9SBarry Smith       }
719566063dSJacob Faibussowitsch       PetscCall(MatLUFactor(pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info));
729566063dSJacob Faibussowitsch       PetscCall(MatFactorGetError(pc->pmat,&err));
7300e125f8SBarry Smith       if (err) { /* Factor() fails */
7400e125f8SBarry Smith         pc->failedreason = (PCFailedReason)err;
756baea169SHong Zhang         PetscFunctionReturn(0);
766baea169SHong Zhang       }
779d76b4d0SMatthew G. Knepley     }
78075768bcSBarry Smith     ((PC_Factor*)dir)->fact = pc->pmat;
799b54502bSHong Zhang   } else {
809b54502bSHong Zhang     MatInfo info;
8100e125f8SBarry Smith 
829b54502bSHong Zhang     if (!pc->setupcalled) {
83f73b0415SBarry Smith       PetscBool canuseordering;
842c7c0729SBarry Smith       if (!((PC_Factor*)dir)->fact) {
859566063dSJacob Faibussowitsch         PetscCall(MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact));
869566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact));
872c7c0729SBarry Smith       }
889566063dSJacob Faibussowitsch       PetscCall(MatFactorGetCanUseOrdering(((PC_Factor*)dir)->fact,&canuseordering));
89f73b0415SBarry Smith       if (canuseordering) {
909566063dSJacob Faibussowitsch         PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
919566063dSJacob Faibussowitsch         PetscCall(MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col));
929b54502bSHong Zhang         if (dir->nonzerosalongdiagonal) {
939566063dSJacob Faibussowitsch           PetscCall(MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col));
949b54502bSHong Zhang         }
959566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row));
969566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col));
9703c60df9SBarry Smith       }
989566063dSJacob Faibussowitsch       PetscCall(MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info));
999566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info));
1003d1c1ea0SBarry Smith       dir->hdr.actualfill = info.fill_ratio_needed;
1019b54502bSHong Zhang     } else if (pc->flag != SAME_NONZERO_PATTERN) {
102f73b0415SBarry Smith       PetscBool canuseordering;
1033d1c1ea0SBarry Smith       if (!dir->hdr.reuseordering) {
1049566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&((PC_Factor*)dir)->fact));
1059566063dSJacob Faibussowitsch         PetscCall(MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact));
1069566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact));
1079566063dSJacob Faibussowitsch         PetscCall(MatFactorGetCanUseOrdering(((PC_Factor*)dir)->fact,&canuseordering));
108f73b0415SBarry Smith         if (canuseordering) {
1099566063dSJacob Faibussowitsch           if (dir->row && dir->col && dir->row != dir->col) PetscCall(ISDestroy(&dir->row));
1109566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&dir->col));
1119566063dSJacob Faibussowitsch           PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
1129566063dSJacob Faibussowitsch           PetscCall(MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col));
1139b54502bSHong Zhang           if (dir->nonzerosalongdiagonal) {
1149566063dSJacob Faibussowitsch             PetscCall(MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col));
1159b54502bSHong Zhang           }
1169566063dSJacob Faibussowitsch           PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row));
1179566063dSJacob Faibussowitsch           PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col));
11803c60df9SBarry Smith         }
1199b54502bSHong Zhang       }
1209566063dSJacob Faibussowitsch       PetscCall(MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info));
1219566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info));
1223d1c1ea0SBarry Smith       dir->hdr.actualfill = info.fill_ratio_needed;
12304545d6dSBarry Smith     } else {
1249566063dSJacob Faibussowitsch       PetscCall(MatFactorGetError(((PC_Factor*)dir)->fact,&err));
125160a8794SBarry Smith       if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) {
1269566063dSJacob Faibussowitsch         PetscCall(MatFactorClearError(((PC_Factor*)dir)->fact));
127b8b68cfdSBarry Smith         pc->failedreason = PC_NOERROR;
12804545d6dSBarry Smith       }
1299b54502bSHong Zhang     }
1309566063dSJacob Faibussowitsch     PetscCall(MatFactorGetError(((PC_Factor*)dir)->fact,&err));
13100e125f8SBarry Smith     if (err) { /* FactorSymbolic() fails */
13200e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1338c1cd74cSHong Zhang       PetscFunctionReturn(0);
1348c1cd74cSHong Zhang     }
1358c1cd74cSHong Zhang 
1369566063dSJacob Faibussowitsch     PetscCall(MatLUFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info));
1379566063dSJacob Faibussowitsch     PetscCall(MatFactorGetError(((PC_Factor*)dir)->fact,&err));
13800e125f8SBarry Smith     if (err) { /* FactorNumeric() fails */
13900e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1408c1cd74cSHong Zhang     }
141680c5173SHong Zhang 
1429b54502bSHong Zhang   }
14300c67f3bSHong Zhang 
1449566063dSJacob Faibussowitsch   PetscCall(PCFactorGetMatSolverType(pc,&stype));
14500c67f3bSHong Zhang   if (!stype) {
146ea799195SBarry Smith     MatSolverType solverpackage;
1479566063dSJacob Faibussowitsch     PetscCall(MatFactorGetSolverType(((PC_Factor*)dir)->fact,&solverpackage));
1489566063dSJacob Faibussowitsch     PetscCall(PCFactorSetMatSolverType(pc,solverpackage));
14900c67f3bSHong Zhang   }
1509b54502bSHong Zhang   PetscFunctionReturn(0);
1519b54502bSHong Zhang }
1529b54502bSHong Zhang 
153574deadeSBarry Smith static PetscErrorCode PCReset_LU(PC pc)
1549b54502bSHong Zhang {
1559b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
1569b54502bSHong Zhang 
1579b54502bSHong Zhang   PetscFunctionBegin;
1589566063dSJacob Faibussowitsch   if (!dir->hdr.inplace && ((PC_Factor*)dir)->fact) PetscCall(MatDestroy(&((PC_Factor*)dir)->fact));
1599566063dSJacob Faibussowitsch   if (dir->row && dir->col && dir->row != dir->col) PetscCall(ISDestroy(&dir->row));
1609566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&dir->col));
161574deadeSBarry Smith   PetscFunctionReturn(0);
162574deadeSBarry Smith }
163574deadeSBarry Smith 
164574deadeSBarry Smith static PetscErrorCode PCDestroy_LU(PC pc)
165574deadeSBarry Smith {
166574deadeSBarry Smith   PC_LU          *dir = (PC_LU*)pc->data;
167574deadeSBarry Smith 
168574deadeSBarry Smith   PetscFunctionBegin;
1699566063dSJacob Faibussowitsch   PetscCall(PCReset_LU(pc));
1709566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PC_Factor*)dir)->ordering));
1719566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PC_Factor*)dir)->solvertype));
172*2e956fe4SStefano Zampini   PetscCall(PCFactorClearComposedFunctions(pc));
1739566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
1749b54502bSHong Zhang   PetscFunctionReturn(0);
1759b54502bSHong Zhang }
1769b54502bSHong Zhang 
1779b54502bSHong Zhang static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y)
1789b54502bSHong Zhang {
1799b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
1809b54502bSHong Zhang 
1819b54502bSHong Zhang   PetscFunctionBegin;
1823d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
1839566063dSJacob Faibussowitsch     PetscCall(MatSolve(pc->pmat,x,y));
1842fa5cd67SKarl Rupp   } else {
1859566063dSJacob Faibussowitsch     PetscCall(MatSolve(((PC_Factor*)dir)->fact,x,y));
1862fa5cd67SKarl Rupp   }
1879b54502bSHong Zhang   PetscFunctionReturn(0);
1889b54502bSHong Zhang }
1899b54502bSHong Zhang 
1907b6e2003SPierre Jolivet static PetscErrorCode PCMatApply_LU(PC pc,Mat X,Mat Y)
1917b6e2003SPierre Jolivet {
1927b6e2003SPierre Jolivet   PC_LU          *dir = (PC_LU*)pc->data;
1937b6e2003SPierre Jolivet 
1947b6e2003SPierre Jolivet   PetscFunctionBegin;
1957b6e2003SPierre Jolivet   if (dir->hdr.inplace) {
1969566063dSJacob Faibussowitsch     PetscCall(MatMatSolve(pc->pmat,X,Y));
1977b6e2003SPierre Jolivet   } else {
1989566063dSJacob Faibussowitsch     PetscCall(MatMatSolve(((PC_Factor*)dir)->fact,X,Y));
1997b6e2003SPierre Jolivet   }
2007b6e2003SPierre Jolivet   PetscFunctionReturn(0);
2017b6e2003SPierre Jolivet }
2027b6e2003SPierre Jolivet 
2039b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y)
2049b54502bSHong Zhang {
2059b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
2069b54502bSHong Zhang 
2079b54502bSHong Zhang   PetscFunctionBegin;
2083d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
2099566063dSJacob Faibussowitsch     PetscCall(MatSolveTranspose(pc->pmat,x,y));
2102fa5cd67SKarl Rupp   } else {
2119566063dSJacob Faibussowitsch     PetscCall(MatSolveTranspose(((PC_Factor*)dir)->fact,x,y));
2122fa5cd67SKarl Rupp   }
2139b54502bSHong Zhang   PetscFunctionReturn(0);
2149b54502bSHong Zhang }
2159b54502bSHong Zhang 
2169b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/
2179b54502bSHong Zhang 
2189b54502bSHong Zhang /*MC
2199b54502bSHong Zhang    PCLU - Uses a direct solver, based on LU factorization, as a preconditioner
2209b54502bSHong Zhang 
2219b54502bSHong Zhang    Options Database Keys:
2222401956bSBarry Smith +  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
2233ca39a21SBarry Smith .  -pc_factor_mat_solver_type - Actives PCFactorSetMatSolverType() to choose the direct solver, like superlu
2242401956bSBarry Smith .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
22555ba2a51SBarry Smith .  -pc_factor_fill <fill> - Sets fill amount
2262401956bSBarry Smith .  -pc_factor_in_place - Activates in-place factorization
2272401956bSBarry Smith .  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
2282401956bSBarry Smith .  -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase
2299b54502bSHong Zhang                                          stability of factorization.
230145b9266SHong Zhang .  -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types
231145b9266SHong Zhang .  -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default
232e22d95b2SBarry Smith -   -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the
233e22d95b2SBarry Smith         diagonal.
2349b54502bSHong Zhang 
23595452b02SPatrick Sanan    Notes:
23695452b02SPatrick Sanan     Not all options work for all matrix formats
2379b54502bSHong Zhang           Run with -help to see additional options for particular matrix formats or factorization
2389b54502bSHong Zhang           algorithms
2399b54502bSHong Zhang 
2409b54502bSHong Zhang    Level: beginner
2419b54502bSHong Zhang 
24295452b02SPatrick Sanan    Notes:
24395452b02SPatrick Sanan     Usually this will compute an "exact" solution in one iteration and does
2449b54502bSHong Zhang           not need a Krylov method (i.e. you can use -ksp_type preonly, or
2459b54502bSHong Zhang           KSPSetType(ksp,KSPPREONLY) for the Krylov method
2469b54502bSHong Zhang 
247db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`,
248db781477SPatrick Sanan           `PCILU`, `PCCHOLESKY`, `PCICC`, `PCFactorSetReuseOrdering()`, `PCFactorSetReuseFill()`, `PCFactorGetMatrix()`,
249db781477SPatrick Sanan           `PCFactorSetFill()`, `PCFactorSetUseInPlace()`, `PCFactorSetMatOrderingType()`, `PCFactorSetColumnPivot()`,
250c2e3fba1SPatrick Sanan           `PCFactorSetPivotInBlocks()`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()`
251db781477SPatrick Sanan           `PCFactorReorderForNonzeroDiagonal()`
2529b54502bSHong Zhang M*/
2539b54502bSHong Zhang 
2548cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_LU(PC pc)
2559b54502bSHong Zhang {
2569b54502bSHong Zhang   PC_LU          *dir;
2579b54502bSHong Zhang 
2589b54502bSHong Zhang   PetscFunctionBegin;
2599566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc,&dir));
2603d1c1ea0SBarry Smith   pc->data = (void*)dir;
2619566063dSJacob Faibussowitsch   PetscCall(PCFactorInitialize(pc,MAT_FACTOR_LU));
2629b54502bSHong Zhang   dir->nonzerosalongdiagonal    = PETSC_FALSE;
2639b54502bSHong Zhang 
264075768bcSBarry Smith   ((PC_Factor*)dir)->info.fill          = 5.0;
265075768bcSBarry Smith   ((PC_Factor*)dir)->info.dtcol         = 1.e-6;  /* default to pivoting; this is only thing PETSc LU supports */
266f4db908eSBarry Smith   ((PC_Factor*)dir)->info.shifttype     = (PetscReal)MAT_SHIFT_NONE;
2670a545947SLisandro Dalcin   dir->col                              = NULL;
2680a545947SLisandro Dalcin   dir->row                              = NULL;
2695c9eb25fSBarry Smith 
270574deadeSBarry Smith   pc->ops->reset             = PCReset_LU;
2719b54502bSHong Zhang   pc->ops->destroy           = PCDestroy_LU;
2729b54502bSHong Zhang   pc->ops->apply             = PCApply_LU;
2737b6e2003SPierre Jolivet   pc->ops->matapply          = PCMatApply_LU;
2749b54502bSHong Zhang   pc->ops->applytranspose    = PCApplyTranspose_LU;
2759b54502bSHong Zhang   pc->ops->setup             = PCSetUp_LU;
2769b54502bSHong Zhang   pc->ops->setfromoptions    = PCSetFromOptions_LU;
27792e08861SBarry Smith   pc->ops->view              = PCView_Factor;
2780a545947SLisandro Dalcin   pc->ops->applyrichardson   = NULL;
2799566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",PCFactorReorderForNonzeroDiagonal_LU));
2809b54502bSHong Zhang   PetscFunctionReturn(0);
2819b54502bSHong Zhang }
282