xref: /petsc/src/ksp/pc/impls/factor/lu/lu.c (revision e36ff184ab5236b473a60bcbbae2ea4eff91e3c6)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2dba47a55SKris Buschelman 
39b54502bSHong Zhang /*
49b54502bSHong Zhang    Defines a direct factorization preconditioner for any Mat implementation
59b54502bSHong Zhang    Note: this need not be consided a preconditioner since it supplies
69b54502bSHong Zhang          a direct solver.
79b54502bSHong Zhang */
8ee45ca4aSHong Zhang 
96356e834SBarry Smith #include "private/pcimpl.h"                /*I "petscpc.h" I*/
10ee45ca4aSHong Zhang #include "src/ksp/pc/impls/factor/lu/lu.h"
119b54502bSHong Zhang 
129b54502bSHong Zhang EXTERN_C_BEGIN
139b54502bSHong Zhang #undef __FUNCT__
14afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetZeroPivot_LU"
15dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot_LU(PC pc,PetscReal z)
16afaefe49SHong Zhang {
17afaefe49SHong Zhang   PC_LU *lu;
18afaefe49SHong Zhang 
19afaefe49SHong Zhang   PetscFunctionBegin;
20afaefe49SHong Zhang   lu                 = (PC_LU*)pc->data;
21afaefe49SHong Zhang   lu->info.zeropivot = z;
22afaefe49SHong Zhang   PetscFunctionReturn(0);
23afaefe49SHong Zhang }
24afaefe49SHong Zhang EXTERN_C_END
25afaefe49SHong Zhang 
26afaefe49SHong Zhang EXTERN_C_BEGIN
27afaefe49SHong Zhang #undef __FUNCT__
28afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetShiftNonzero_LU"
29dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero_LU(PC pc,PetscReal shift)
30afaefe49SHong Zhang {
31afaefe49SHong Zhang   PC_LU *dir;
32afaefe49SHong Zhang 
33afaefe49SHong Zhang   PetscFunctionBegin;
34afaefe49SHong Zhang   dir = (PC_LU*)pc->data;
35afaefe49SHong Zhang   if (shift == (PetscReal) PETSC_DECIDE) {
36afaefe49SHong Zhang     dir->info.shiftnz = 1.e-12;
37afaefe49SHong Zhang   } else {
38afaefe49SHong Zhang     dir->info.shiftnz = shift;
39afaefe49SHong Zhang   }
40afaefe49SHong Zhang   PetscFunctionReturn(0);
41afaefe49SHong Zhang }
42afaefe49SHong Zhang EXTERN_C_END
43afaefe49SHong Zhang 
44afaefe49SHong Zhang EXTERN_C_BEGIN
45afaefe49SHong Zhang #undef __FUNCT__
46afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetShiftPd_LU"
47dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd_LU(PC pc,PetscTruth shift)
48afaefe49SHong Zhang {
49afaefe49SHong Zhang   PC_LU *dir;
50afaefe49SHong Zhang 
51afaefe49SHong Zhang   PetscFunctionBegin;
52afaefe49SHong Zhang   dir = (PC_LU*)pc->data;
53fbf22428SSatish Balay   if (shift) {
54fbf22428SSatish Balay     dir->info.shift_fraction = 0.0;
55fbf22428SSatish Balay     dir->info.shiftpd = 1.0;
56fbf22428SSatish Balay   } else {
57fbf22428SSatish Balay     dir->info.shiftpd = 0.0;
58fbf22428SSatish Balay   }
59afaefe49SHong Zhang   PetscFunctionReturn(0);
60afaefe49SHong Zhang }
61afaefe49SHong Zhang EXTERN_C_END
62afaefe49SHong Zhang 
63afaefe49SHong Zhang EXTERN_C_BEGIN
64afaefe49SHong Zhang #undef __FUNCT__
652401956bSBarry Smith #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal_LU"
662401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal_LU(PC pc,PetscReal z)
679b54502bSHong Zhang {
689b54502bSHong Zhang   PC_LU *lu = (PC_LU*)pc->data;
699b54502bSHong Zhang 
709b54502bSHong Zhang   PetscFunctionBegin;
719b54502bSHong Zhang   lu->nonzerosalongdiagonal = PETSC_TRUE;
729b54502bSHong Zhang   if (z == PETSC_DECIDE) {
739b54502bSHong Zhang     lu->nonzerosalongdiagonaltol = 1.e-10;
749b54502bSHong Zhang   } else {
759b54502bSHong Zhang     lu->nonzerosalongdiagonaltol = z;
769b54502bSHong Zhang   }
779b54502bSHong Zhang   PetscFunctionReturn(0);
789b54502bSHong Zhang }
799b54502bSHong Zhang EXTERN_C_END
809b54502bSHong Zhang 
819b54502bSHong Zhang EXTERN_C_BEGIN
829b54502bSHong Zhang #undef __FUNCT__
832401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering_LU"
842401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering_LU(PC pc,PetscTruth flag)
859b54502bSHong Zhang {
869b54502bSHong Zhang   PC_LU *lu;
879b54502bSHong Zhang 
889b54502bSHong Zhang   PetscFunctionBegin;
899b54502bSHong Zhang   lu                = (PC_LU*)pc->data;
909b54502bSHong Zhang   lu->reuseordering = flag;
919b54502bSHong Zhang   PetscFunctionReturn(0);
929b54502bSHong Zhang }
939b54502bSHong Zhang EXTERN_C_END
949b54502bSHong Zhang 
959b54502bSHong Zhang EXTERN_C_BEGIN
969b54502bSHong Zhang #undef __FUNCT__
972401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseFill_LU"
982401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill_LU(PC pc,PetscTruth flag)
999b54502bSHong Zhang {
1009b54502bSHong Zhang   PC_LU *lu;
1019b54502bSHong Zhang 
1029b54502bSHong Zhang   PetscFunctionBegin;
1039b54502bSHong Zhang   lu = (PC_LU*)pc->data;
1049b54502bSHong Zhang   lu->reusefill = flag;
1059b54502bSHong Zhang   PetscFunctionReturn(0);
1069b54502bSHong Zhang }
1079b54502bSHong Zhang EXTERN_C_END
1089b54502bSHong Zhang 
1099b54502bSHong Zhang #undef __FUNCT__
1109b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_LU"
1119b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_LU(PC pc)
1129b54502bSHong Zhang {
1139b54502bSHong Zhang   PC_LU          *lu = (PC_LU*)pc->data;
1149b54502bSHong Zhang   PetscErrorCode ierr;
1159b54502bSHong Zhang   PetscTruth     flg,set;
1169989ab13SBarry Smith   char           tname[256],solvertype[64];
1179989ab13SBarry Smith   MatSolverType  stype;
1189b54502bSHong Zhang   PetscFList     ordlist;
1199b54502bSHong Zhang   PetscReal      tol;
1209b54502bSHong Zhang 
1219b54502bSHong Zhang   PetscFunctionBegin;
1229b54502bSHong Zhang   ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr);
1239b54502bSHong Zhang   ierr = PetscOptionsHead("LU options");CHKERRQ(ierr);
1242401956bSBarry Smith     ierr = PetscOptionsName("-pc_factor_in_place","Form LU in the same memory as the matrix","PCFactorSetUseInPlace",&flg);CHKERRQ(ierr);
1259b54502bSHong Zhang     if (flg) {
1262401956bSBarry Smith       ierr = PCFactorSetUseInPlace(pc);CHKERRQ(ierr);
1279b54502bSHong Zhang     }
12855ba2a51SBarry Smith     ierr = PetscOptionsReal("-pc_factor_fill","Expected non-zeros in LU/non-zeros in matrix","PCFactorSetFill",lu->info.fill,&lu->info.fill,0);CHKERRQ(ierr);
1299b54502bSHong Zhang 
1309f95998fSHong Zhang     ierr = PetscOptionsName("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr);
1319b54502bSHong Zhang     if (flg) {
132afaefe49SHong Zhang         ierr = PCFactorSetShiftNonzero(pc,(PetscReal) PETSC_DECIDE);CHKERRQ(ierr);
1339b54502bSHong Zhang     }
1349f95998fSHong Zhang     ierr = PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",lu->info.shiftnz,&lu->info.shiftnz,0);CHKERRQ(ierr);
1359f95998fSHong Zhang     ierr = PetscOptionsName("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",&flg);CHKERRQ(ierr);
1369b54502bSHong Zhang     if (flg) {
137afaefe49SHong Zhang       ierr = PCFactorSetShiftPd(pc,PETSC_TRUE);CHKERRQ(ierr);
1389b54502bSHong Zhang     }
139ee45ca4aSHong Zhang     ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",lu->info.zeropivot,&lu->info.zeropivot,0);CHKERRQ(ierr);
1409b54502bSHong Zhang 
1412401956bSBarry Smith     ierr = PetscOptionsName("-pc_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",&flg);CHKERRQ(ierr);
1429b54502bSHong Zhang     if (flg) {
1432401956bSBarry Smith       ierr = PCFactorSetReuseFill(pc,PETSC_TRUE);CHKERRQ(ierr);
1449b54502bSHong Zhang     }
1452401956bSBarry Smith     ierr = PetscOptionsName("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",&flg);CHKERRQ(ierr);
1469b54502bSHong Zhang     if (flg) {
1472401956bSBarry Smith       ierr = PCFactorSetReuseOrdering(pc,PETSC_TRUE);CHKERRQ(ierr);
1489b54502bSHong Zhang     }
1499b54502bSHong Zhang 
1509b54502bSHong Zhang     ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr);
151e5a9bf91SBarry Smith     ierr = PetscOptionsList("-pc_factor_mat_ordering_type","Reordering to reduce nonzeros in LU","PCFactorSetMatOrderingType",ordlist,lu->ordering,tname,256,&flg);CHKERRQ(ierr);
1529b54502bSHong Zhang     if (flg) {
153e5a9bf91SBarry Smith       ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr);
1549b54502bSHong Zhang     }
1559b54502bSHong Zhang 
1562401956bSBarry Smith     ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr);
1579b54502bSHong Zhang     if (flg) {
1589b54502bSHong Zhang       tol = PETSC_DECIDE;
1592401956bSBarry Smith       ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr);
1602401956bSBarry Smith       ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr);
1619b54502bSHong Zhang     }
1629b54502bSHong Zhang 
1632401956bSBarry Smith     ierr = PetscOptionsReal("-pc_factor_pivoting","Pivoting tolerance (used only for some factorization)","PCFactorSetPivoting",lu->info.dtcol,&lu->info.dtcol,&flg);CHKERRQ(ierr);
1649b54502bSHong Zhang 
1659b54502bSHong Zhang     flg = lu->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE;
1662401956bSBarry Smith     ierr = PetscOptionsTruth("-pc_factor_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCFactorSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr);
1679b54502bSHong Zhang     if (set) {
1682401956bSBarry Smith       ierr = PCFactorSetPivotInBlocks(pc,flg);CHKERRQ(ierr);
1699b54502bSHong Zhang     }
170a542b6e8SBarry Smith     if (pc->pmat) {
1719989ab13SBarry Smith       ierr = MatGetSolverType(pc->pmat,&stype);CHKERRQ(ierr);
1729989ab13SBarry Smith       ierr = PetscOptionsString("-pc_factor_solver_type","Type of solver to use for factorization","MatSetSolverType",stype,solvertype,64,&set);CHKERRQ(ierr);
1739989ab13SBarry Smith       if (set) {
1749989ab13SBarry Smith 	ierr = MatSetSolverType(pc->pmat,solvertype);CHKERRQ(ierr);
1759989ab13SBarry Smith       }
176a542b6e8SBarry Smith     }
1779b54502bSHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
1789b54502bSHong Zhang   PetscFunctionReturn(0);
1799b54502bSHong Zhang }
1809b54502bSHong Zhang 
1819b54502bSHong Zhang #undef __FUNCT__
1829b54502bSHong Zhang #define __FUNCT__ "PCView_LU"
1839b54502bSHong Zhang static PetscErrorCode PCView_LU(PC pc,PetscViewer viewer)
1849b54502bSHong Zhang {
1859b54502bSHong Zhang   PC_LU          *lu = (PC_LU*)pc->data;
1869b54502bSHong Zhang   PetscErrorCode ierr;
1879b54502bSHong Zhang   PetscTruth     iascii,isstring;
1889b54502bSHong Zhang 
1899b54502bSHong Zhang   PetscFunctionBegin;
1909b54502bSHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
1919b54502bSHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
1929b54502bSHong Zhang   if (iascii) {
1939b54502bSHong Zhang 
1949b54502bSHong Zhang     if (lu->inplace) {ierr = PetscViewerASCIIPrintf(viewer,"  LU: in-place factorization\n");CHKERRQ(ierr);}
1959b54502bSHong Zhang     else             {ierr = PetscViewerASCIIPrintf(viewer,"  LU: out-of-place factorization\n");CHKERRQ(ierr);}
1969b54502bSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"    matrix ordering: %s\n",lu->ordering);CHKERRQ(ierr);
197a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  LU: tolerance for zero pivot %G\n",lu->info.zeropivot);CHKERRQ(ierr);
1980a29876aSHong Zhang     if (lu->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer,"  LU: using Manteuffel shift\n");CHKERRQ(ierr);}
1999b54502bSHong Zhang     if (lu->fact) {
200f3a39becSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  LU: factor fill ratio needed %G\n",lu->actualfill);CHKERRQ(ierr);
201f3a39becSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"       Factored matrix follows\n");CHKERRQ(ierr);
202f3a39becSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
203f3a39becSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
204f3a39becSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
205f3a39becSBarry Smith       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
2069b54502bSHong Zhang       ierr = MatView(lu->fact,viewer);CHKERRQ(ierr);
2079b54502bSHong Zhang       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
208f3a39becSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
209f3a39becSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
210f3a39becSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
2119b54502bSHong Zhang     }
2129b54502bSHong Zhang     if (lu->reusefill)    {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing fill from past factorization\n");CHKERRQ(ierr);}
2139b54502bSHong Zhang     if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing reordering from past factorization\n");CHKERRQ(ierr);}
2149b54502bSHong Zhang   } else if (isstring) {
2159b54502bSHong Zhang     ierr = PetscViewerStringSPrintf(viewer," order=%s",lu->ordering);CHKERRQ(ierr);CHKERRQ(ierr);
2169b54502bSHong Zhang   } else {
2179b54502bSHong Zhang     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCLU",((PetscObject)viewer)->type_name);
2189b54502bSHong Zhang   }
2199b54502bSHong Zhang   PetscFunctionReturn(0);
2209b54502bSHong Zhang }
2219b54502bSHong Zhang 
2229b54502bSHong Zhang #undef __FUNCT__
223a4fd02acSBarry Smith #define __FUNCT__ "PCFactorGetMatrix_LU"
224a4fd02acSBarry Smith static PetscErrorCode PCFactorGetMatrix_LU(PC pc,Mat *mat)
2259b54502bSHong Zhang {
2269b54502bSHong Zhang   PC_LU *dir = (PC_LU*)pc->data;
2279b54502bSHong Zhang 
2289b54502bSHong Zhang   PetscFunctionBegin;
2299b54502bSHong Zhang   if (!dir->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()");
2309b54502bSHong Zhang   *mat = dir->fact;
2319b54502bSHong Zhang   PetscFunctionReturn(0);
2329b54502bSHong Zhang }
2339b54502bSHong Zhang 
2349b54502bSHong Zhang #undef __FUNCT__
2359b54502bSHong Zhang #define __FUNCT__ "PCSetUp_LU"
2369b54502bSHong Zhang static PetscErrorCode PCSetUp_LU(PC pc)
2379b54502bSHong Zhang {
2389b54502bSHong Zhang   PetscErrorCode ierr;
2399b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
2409b54502bSHong Zhang 
2419b54502bSHong Zhang   PetscFunctionBegin;
2429b54502bSHong Zhang   if (dir->reusefill && pc->setupcalled) dir->info.fill = dir->actualfill;
2439b54502bSHong Zhang 
2449b54502bSHong Zhang   if (dir->inplace) {
2459b54502bSHong Zhang     if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);}
2469b54502bSHong Zhang     if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);}
2479b54502bSHong Zhang     ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
24803c60df9SBarry Smith     if (dir->row) {
24903c60df9SBarry Smith       ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);
25003c60df9SBarry Smith       ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr);
25103c60df9SBarry Smith     }
2529b54502bSHong Zhang     ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&dir->info);CHKERRQ(ierr);
2539b54502bSHong Zhang     dir->fact = pc->pmat;
2549b54502bSHong Zhang   } else {
2559b54502bSHong Zhang     MatInfo info;
2569b54502bSHong Zhang     if (!pc->setupcalled) {
2579b54502bSHong Zhang       ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
2589b54502bSHong Zhang       if (dir->nonzerosalongdiagonal) {
2599b54502bSHong Zhang         ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr);
2609b54502bSHong Zhang       }
26103c60df9SBarry Smith       if (dir->row) {
26203c60df9SBarry Smith         ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);
26303c60df9SBarry Smith         ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr);
26403c60df9SBarry Smith       }
2659b54502bSHong Zhang       ierr = MatLUFactorSymbolic(pc->pmat,dir->row,dir->col,&dir->info,&dir->fact);CHKERRQ(ierr);
2669b54502bSHong Zhang       ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
2679b54502bSHong Zhang       dir->actualfill = info.fill_ratio_needed;
26852e6d16bSBarry Smith       ierr = PetscLogObjectParent(pc,dir->fact);CHKERRQ(ierr);
2699b54502bSHong Zhang     } else if (pc->flag != SAME_NONZERO_PATTERN) {
2709b54502bSHong Zhang       if (!dir->reuseordering) {
2719b54502bSHong Zhang         if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);}
2729b54502bSHong Zhang         if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);}
2739b54502bSHong Zhang         ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
2749b54502bSHong Zhang         if (dir->nonzerosalongdiagonal) {
2759b54502bSHong Zhang           ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr);
2769b54502bSHong Zhang         }
27703c60df9SBarry Smith         if (dir->row) {
27803c60df9SBarry Smith           ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);
27903c60df9SBarry Smith           ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr);
28003c60df9SBarry Smith         }
2819b54502bSHong Zhang       }
2829b54502bSHong Zhang       ierr = MatDestroy(dir->fact);CHKERRQ(ierr);
2839b54502bSHong Zhang       ierr = MatLUFactorSymbolic(pc->pmat,dir->row,dir->col,&dir->info,&dir->fact);CHKERRQ(ierr);
2849b54502bSHong Zhang       ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
2859b54502bSHong Zhang       dir->actualfill = info.fill_ratio_needed;
28652e6d16bSBarry Smith       ierr = PetscLogObjectParent(pc,dir->fact);CHKERRQ(ierr);
2879b54502bSHong Zhang     }
2889b54502bSHong Zhang     ierr = MatLUFactorNumeric(pc->pmat,&dir->info,&dir->fact);CHKERRQ(ierr);
2899b54502bSHong Zhang   }
2909b54502bSHong Zhang   PetscFunctionReturn(0);
2919b54502bSHong Zhang }
2929b54502bSHong Zhang 
2939b54502bSHong Zhang #undef __FUNCT__
2949b54502bSHong Zhang #define __FUNCT__ "PCDestroy_LU"
2959b54502bSHong Zhang static PetscErrorCode PCDestroy_LU(PC pc)
2969b54502bSHong Zhang {
2979b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
2989b54502bSHong Zhang   PetscErrorCode ierr;
2999b54502bSHong Zhang 
3009b54502bSHong Zhang   PetscFunctionBegin;
3019b54502bSHong Zhang   if (!dir->inplace && dir->fact) {ierr = MatDestroy(dir->fact);CHKERRQ(ierr);}
3029b54502bSHong Zhang   if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);}
3039b54502bSHong Zhang   if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);}
30405dd20ceSSatish Balay   ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr);
3059b54502bSHong Zhang   ierr = PetscFree(dir);CHKERRQ(ierr);
3069b54502bSHong Zhang   PetscFunctionReturn(0);
3079b54502bSHong Zhang }
3089b54502bSHong Zhang 
3099b54502bSHong Zhang #undef __FUNCT__
3109b54502bSHong Zhang #define __FUNCT__ "PCApply_LU"
3119b54502bSHong Zhang static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y)
3129b54502bSHong Zhang {
3139b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
3149b54502bSHong Zhang   PetscErrorCode ierr;
3159b54502bSHong Zhang 
3169b54502bSHong Zhang   PetscFunctionBegin;
3179b54502bSHong Zhang   if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);}
3189b54502bSHong Zhang   else              {ierr = MatSolve(dir->fact,x,y);CHKERRQ(ierr);}
3199b54502bSHong Zhang   PetscFunctionReturn(0);
3209b54502bSHong Zhang }
3219b54502bSHong Zhang 
3229b54502bSHong Zhang #undef __FUNCT__
3239b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_LU"
3249b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y)
3259b54502bSHong Zhang {
3269b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
3279b54502bSHong Zhang   PetscErrorCode ierr;
3289b54502bSHong Zhang 
3299b54502bSHong Zhang   PetscFunctionBegin;
3309b54502bSHong Zhang   if (dir->inplace) {ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);}
3319b54502bSHong Zhang   else              {ierr = MatSolveTranspose(dir->fact,x,y);CHKERRQ(ierr);}
3329b54502bSHong Zhang   PetscFunctionReturn(0);
3339b54502bSHong Zhang }
3349b54502bSHong Zhang 
3359b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/
3369b54502bSHong Zhang 
3379b54502bSHong Zhang EXTERN_C_BEGIN
3389b54502bSHong Zhang #undef __FUNCT__
33955ba2a51SBarry Smith #define __FUNCT__ "PCFactorSetFill_LU"
34055ba2a51SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill_LU(PC pc,PetscReal fill)
3419b54502bSHong Zhang {
3429b54502bSHong Zhang   PC_LU *dir;
3439b54502bSHong Zhang 
3449b54502bSHong Zhang   PetscFunctionBegin;
3459b54502bSHong Zhang   dir = (PC_LU*)pc->data;
3469b54502bSHong Zhang   dir->info.fill = fill;
3479b54502bSHong Zhang   PetscFunctionReturn(0);
3489b54502bSHong Zhang }
3499b54502bSHong Zhang EXTERN_C_END
3509b54502bSHong Zhang 
3519b54502bSHong Zhang EXTERN_C_BEGIN
3529b54502bSHong Zhang #undef __FUNCT__
3532401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace_LU"
3542401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_LU(PC pc)
3559b54502bSHong Zhang {
3569b54502bSHong Zhang   PC_LU *dir;
3579b54502bSHong Zhang 
3589b54502bSHong Zhang   PetscFunctionBegin;
3599b54502bSHong Zhang   dir = (PC_LU*)pc->data;
3609b54502bSHong Zhang   dir->inplace = PETSC_TRUE;
3619b54502bSHong Zhang   PetscFunctionReturn(0);
3629b54502bSHong Zhang }
3639b54502bSHong Zhang EXTERN_C_END
3649b54502bSHong Zhang 
3659b54502bSHong Zhang EXTERN_C_BEGIN
3669b54502bSHong Zhang #undef __FUNCT__
367e5a9bf91SBarry Smith #define __FUNCT__ "PCFactorSetMatOrderingType_LU"
368*e36ff184SSatish Balay PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType_LU(PC pc,const MatOrderingType ordering)
3699b54502bSHong Zhang {
3709b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
3719b54502bSHong Zhang   PetscErrorCode ierr;
3729b54502bSHong Zhang 
3739b54502bSHong Zhang   PetscFunctionBegin;
37405dd20ceSSatish Balay   ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr);
37505dd20ceSSatish Balay   ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr);
3769b54502bSHong Zhang   PetscFunctionReturn(0);
3779b54502bSHong Zhang }
3789b54502bSHong Zhang EXTERN_C_END
3799b54502bSHong Zhang 
3809b54502bSHong Zhang EXTERN_C_BEGIN
3819b54502bSHong Zhang #undef __FUNCT__
3822401956bSBarry Smith #define __FUNCT__ "PCFactorSetPivoting_LU"
3832401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivoting_LU(PC pc,PetscReal dtcol)
3849b54502bSHong Zhang {
3859b54502bSHong Zhang   PC_LU *dir = (PC_LU*)pc->data;
3869b54502bSHong Zhang 
3879b54502bSHong Zhang   PetscFunctionBegin;
388a83599f4SBarry Smith   if (dtcol < 0.0 || dtcol > 1.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Column pivot tolerance is %G must be between 0 and 1",dtcol);
3899b54502bSHong Zhang   dir->info.dtcol = dtcol;
3909b54502bSHong Zhang   PetscFunctionReturn(0);
3919b54502bSHong Zhang }
3929b54502bSHong Zhang EXTERN_C_END
3939b54502bSHong Zhang 
3949b54502bSHong Zhang EXTERN_C_BEGIN
3959b54502bSHong Zhang #undef __FUNCT__
3962401956bSBarry Smith #define __FUNCT__ "PCFactorSetPivotInBlocks_LU"
3972401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks_LU(PC pc,PetscTruth pivot)
3989b54502bSHong Zhang {
3999b54502bSHong Zhang   PC_LU *dir = (PC_LU*)pc->data;
4009b54502bSHong Zhang 
4019b54502bSHong Zhang   PetscFunctionBegin;
4029b54502bSHong Zhang   dir->info.pivotinblocks = pivot ? 1.0 : 0.0;
4039b54502bSHong Zhang   PetscFunctionReturn(0);
4049b54502bSHong Zhang }
4059b54502bSHong Zhang EXTERN_C_END
4069b54502bSHong Zhang 
4079b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/
4089b54502bSHong Zhang 
4099b54502bSHong Zhang #undef __FUNCT__
4102401956bSBarry Smith #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal"
4119b54502bSHong Zhang /*@
4122401956bSBarry Smith    PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal
4139b54502bSHong Zhang 
4149b54502bSHong Zhang    Collective on PC
4159b54502bSHong Zhang 
4169b54502bSHong Zhang    Input Parameters:
4179b54502bSHong Zhang +  pc - the preconditioner context
4189b54502bSHong Zhang -  tol - diagonal entries smaller than this in absolute value are considered zero
4199b54502bSHong Zhang 
4209b54502bSHong Zhang    Options Database Key:
4212401956bSBarry Smith .  -pc_factor_nonzeros_along_diagonal
4229b54502bSHong Zhang 
4239b54502bSHong Zhang    Level: intermediate
4249b54502bSHong Zhang 
4259b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill
4269b54502bSHong Zhang 
42755ba2a51SBarry Smith .seealso: PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetZeroPivot(), MatReorderForNonzeroDiagonal()
4289b54502bSHong Zhang @*/
4292401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal(PC pc,PetscReal rtol)
4309b54502bSHong Zhang {
4319b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscReal);
4329b54502bSHong Zhang 
4339b54502bSHong Zhang   PetscFunctionBegin;
4349b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
4352401956bSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",(void (**)(void))&f);CHKERRQ(ierr);
4369b54502bSHong Zhang   if (f) {
4379b54502bSHong Zhang     ierr = (*f)(pc,rtol);CHKERRQ(ierr);
4389b54502bSHong Zhang   }
4399b54502bSHong Zhang   PetscFunctionReturn(0);
4409b54502bSHong Zhang }
4419b54502bSHong Zhang 
4429b54502bSHong Zhang #undef __FUNCT__
44355ba2a51SBarry Smith #define __FUNCT__ "PCFactorSetFill"
4449b54502bSHong Zhang /*@
44555ba2a51SBarry Smith    PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix,
4469b54502bSHong Zhang    fill = number nonzeros in factor/number nonzeros in original matrix.
4479b54502bSHong Zhang 
4489b54502bSHong Zhang    Collective on PC
4499b54502bSHong Zhang 
4509b54502bSHong Zhang    Input Parameters:
4519b54502bSHong Zhang +  pc - the preconditioner context
4529b54502bSHong Zhang -  fill - amount of expected fill
4539b54502bSHong Zhang 
4549b54502bSHong Zhang    Options Database Key:
45555ba2a51SBarry Smith .  -pc_factor_fill <fill> - Sets fill amount
4569b54502bSHong Zhang 
4579b54502bSHong Zhang    Level: intermediate
4589b54502bSHong Zhang 
4599b54502bSHong Zhang    Note:
4609b54502bSHong Zhang    For sparse matrix factorizations it is difficult to predict how much
4616cf91177SBarry Smith    fill to expect. By running with the option -info PETSc will print the
4629b54502bSHong Zhang    actual amount of fill used; allowing you to set the value accurately for
4639b54502bSHong Zhang    future runs. Default PETSc uses a value of 5.0
4649b54502bSHong Zhang 
4659b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill
4669b54502bSHong Zhang 
4679b54502bSHong Zhang @*/
46855ba2a51SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill(PC pc,PetscReal fill)
4699b54502bSHong Zhang {
4709b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscReal);
4719b54502bSHong Zhang 
4729b54502bSHong Zhang   PetscFunctionBegin;
4739b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
4749b54502bSHong Zhang   if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0");
47555ba2a51SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetFill_C",(void (**)(void))&f);CHKERRQ(ierr);
4769b54502bSHong Zhang   if (f) {
4779b54502bSHong Zhang     ierr = (*f)(pc,fill);CHKERRQ(ierr);
4789b54502bSHong Zhang   }
4799b54502bSHong Zhang   PetscFunctionReturn(0);
4809b54502bSHong Zhang }
4819b54502bSHong Zhang 
4829b54502bSHong Zhang #undef __FUNCT__
4832401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace"
4849b54502bSHong Zhang /*@
4852401956bSBarry Smith    PCFactorSetUseInPlace - Tells the system to do an in-place factorization.
4869b54502bSHong Zhang    For dense matrices, this enables the solution of much larger problems.
4879b54502bSHong Zhang    For sparse matrices the factorization cannot be done truly in-place
4889b54502bSHong Zhang    so this does not save memory during the factorization, but after the matrix
4899b54502bSHong Zhang    is factored, the original unfactored matrix is freed, thus recovering that
4909b54502bSHong Zhang    space.
4919b54502bSHong Zhang 
4929b54502bSHong Zhang    Collective on PC
4939b54502bSHong Zhang 
4949b54502bSHong Zhang    Input Parameters:
4959b54502bSHong Zhang .  pc - the preconditioner context
4969b54502bSHong Zhang 
4979b54502bSHong Zhang    Options Database Key:
4982401956bSBarry Smith .  -pc_factor_in_place - Activates in-place factorization
4999b54502bSHong Zhang 
5009b54502bSHong Zhang    Notes:
5012401956bSBarry Smith    PCFactorSetUseInplace() can only be used with the KSP method KSPPREONLY or when
5029b54502bSHong Zhang    a different matrix is provided for the multiply and the preconditioner in
5039b54502bSHong Zhang    a call to KSPSetOperators().
5049b54502bSHong Zhang    This is because the Krylov space methods require an application of the
5059b54502bSHong Zhang    matrix multiplication, which is not possible here because the matrix has
5069b54502bSHong Zhang    been factored in-place, replacing the original matrix.
5079b54502bSHong Zhang 
5089b54502bSHong Zhang    Level: intermediate
5099b54502bSHong Zhang 
5109b54502bSHong Zhang .keywords: PC, set, factorization, direct, inplace, in-place, LU
5119b54502bSHong Zhang 
5129b54502bSHong Zhang .seealso: PCILUSetUseInPlace()
5139b54502bSHong Zhang @*/
5142401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace(PC pc)
5159b54502bSHong Zhang {
5169b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC);
5179b54502bSHong Zhang 
5189b54502bSHong Zhang   PetscFunctionBegin;
5199b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
5202401956bSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",(void (**)(void))&f);CHKERRQ(ierr);
5219b54502bSHong Zhang   if (f) {
5229b54502bSHong Zhang     ierr = (*f)(pc);CHKERRQ(ierr);
5239b54502bSHong Zhang   }
5249b54502bSHong Zhang   PetscFunctionReturn(0);
5259b54502bSHong Zhang }
5269b54502bSHong Zhang 
5279b54502bSHong Zhang #undef __FUNCT__
528e5a9bf91SBarry Smith #define __FUNCT__ "PCFactorSetMatOrderingType"
5299b54502bSHong Zhang /*@C
530e5a9bf91SBarry Smith     PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to
5319b54502bSHong Zhang     be used in the LU factorization.
5329b54502bSHong Zhang 
5339b54502bSHong Zhang     Collective on PC
5349b54502bSHong Zhang 
5359b54502bSHong Zhang     Input Parameters:
5369b54502bSHong Zhang +   pc - the preconditioner context
5379b54502bSHong Zhang -   ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM
5389b54502bSHong Zhang 
5399b54502bSHong Zhang     Options Database Key:
5402401956bSBarry Smith .   -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
5419b54502bSHong Zhang 
5429b54502bSHong Zhang     Level: intermediate
5439b54502bSHong Zhang 
5449b54502bSHong Zhang     Notes: nested dissection is used by default
5459b54502bSHong Zhang 
5462b6de112SBarry Smith     For Cholesky and ICC and the SBAIJ format reorderings are not available,
5472b6de112SBarry Smith     since only the upper triangular part of the matrix is stored. You can use the
5482b6de112SBarry Smith     SeqAIJ format in this case to get reorderings.
5492b6de112SBarry Smith 
5509b54502bSHong Zhang @*/
551*e36ff184SSatish Balay PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType(PC pc,const MatOrderingType ordering)
5529b54502bSHong Zhang {
553*e36ff184SSatish Balay   PetscErrorCode ierr,(*f)(PC,const MatOrderingType);
5549b54502bSHong Zhang 
5559b54502bSHong Zhang   PetscFunctionBegin;
556e5a9bf91SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",(void (**)(void))&f);CHKERRQ(ierr);
5579b54502bSHong Zhang   if (f) {
5589b54502bSHong Zhang     ierr = (*f)(pc,ordering);CHKERRQ(ierr);
5599b54502bSHong Zhang   }
5609b54502bSHong Zhang   PetscFunctionReturn(0);
5619b54502bSHong Zhang }
5629b54502bSHong Zhang 
5639b54502bSHong Zhang #undef __FUNCT__
5642401956bSBarry Smith #define __FUNCT__ "PCFactorSetPivoting"
5659b54502bSHong Zhang /*@
5662401956bSBarry Smith     PCFactorSetPivoting - Determines when pivoting is done during LU.
5679b54502bSHong Zhang       For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices
5689b54502bSHong Zhang       it is never done. For the Matlab and SuperLU factorization this is used.
5699b54502bSHong Zhang 
5709b54502bSHong Zhang     Collective on PC
5719b54502bSHong Zhang 
5729b54502bSHong Zhang     Input Parameters:
5739b54502bSHong Zhang +   pc - the preconditioner context
5749b54502bSHong Zhang -   dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable)
5759b54502bSHong Zhang 
5769b54502bSHong Zhang     Options Database Key:
5772401956bSBarry Smith .   -pc_factor_pivoting <dtcol>
5789b54502bSHong Zhang 
5799b54502bSHong Zhang     Level: intermediate
5809b54502bSHong Zhang 
5812401956bSBarry Smith .seealso: PCILUSetMatOrdering(), PCFactorSetPivotInBlocks()
5829b54502bSHong Zhang @*/
5832401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivoting(PC pc,PetscReal dtcol)
5849b54502bSHong Zhang {
5859b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscReal);
5869b54502bSHong Zhang 
5879b54502bSHong Zhang   PetscFunctionBegin;
5882401956bSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivoting_C",(void (**)(void))&f);CHKERRQ(ierr);
5899b54502bSHong Zhang   if (f) {
5909b54502bSHong Zhang     ierr = (*f)(pc,dtcol);CHKERRQ(ierr);
5919b54502bSHong Zhang   }
5929b54502bSHong Zhang   PetscFunctionReturn(0);
5939b54502bSHong Zhang }
5949b54502bSHong Zhang 
5959b54502bSHong Zhang #undef __FUNCT__
5962401956bSBarry Smith #define __FUNCT__ "PCFactorSetPivotInBlocks"
5979b54502bSHong Zhang /*@
5982401956bSBarry Smith     PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block
5999b54502bSHong Zhang       with BAIJ or SBAIJ matrices
6009b54502bSHong Zhang 
6019b54502bSHong Zhang     Collective on PC
6029b54502bSHong Zhang 
6039b54502bSHong Zhang     Input Parameters:
6049b54502bSHong Zhang +   pc - the preconditioner context
6059b54502bSHong Zhang -   pivot - PETSC_TRUE or PETSC_FALSE
6069b54502bSHong Zhang 
6079b54502bSHong Zhang     Options Database Key:
6082401956bSBarry Smith .   -pc_factor_pivot_in_blocks <true,false>
6099b54502bSHong Zhang 
6109b54502bSHong Zhang     Level: intermediate
6119b54502bSHong Zhang 
6122401956bSBarry Smith .seealso: PCILUSetMatOrdering(), PCFactorSetPivoting()
6139b54502bSHong Zhang @*/
6142401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks(PC pc,PetscTruth pivot)
6159b54502bSHong Zhang {
6169b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscTruth);
6179b54502bSHong Zhang 
6189b54502bSHong Zhang   PetscFunctionBegin;
6192401956bSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr);
6209b54502bSHong Zhang   if (f) {
6219b54502bSHong Zhang     ierr = (*f)(pc,pivot);CHKERRQ(ierr);
6229b54502bSHong Zhang   }
6239b54502bSHong Zhang   PetscFunctionReturn(0);
6249b54502bSHong Zhang }
6259b54502bSHong Zhang 
6269b54502bSHong Zhang /* ------------------------------------------------------------------------ */
6279b54502bSHong Zhang 
6289b54502bSHong Zhang /*MC
6299b54502bSHong Zhang    PCLU - Uses a direct solver, based on LU factorization, as a preconditioner
6309b54502bSHong Zhang 
6319b54502bSHong Zhang    Options Database Keys:
6322401956bSBarry Smith +  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
6332401956bSBarry Smith .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
63455ba2a51SBarry Smith .  -pc_factor_fill <fill> - Sets fill amount
6352401956bSBarry Smith .  -pc_factor_in_place - Activates in-place factorization
6362401956bSBarry Smith .  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
6372401956bSBarry Smith .  -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase
6389b54502bSHong Zhang                                          stability of factorization.
639f251bdbdSHong Zhang .  -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default
640e22d95b2SBarry Smith .  -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value
641f251bdbdSHong Zhang         is optional with PETSC_TRUE being the default
642e22d95b2SBarry Smith -   -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the
643e22d95b2SBarry Smith         diagonal.
6449b54502bSHong Zhang 
6459b54502bSHong Zhang    Notes: Not all options work for all matrix formats
6469b54502bSHong Zhang           Run with -help to see additional options for particular matrix formats or factorization
6479b54502bSHong Zhang           algorithms
6489b54502bSHong Zhang 
6499b54502bSHong Zhang    Level: beginner
6509b54502bSHong Zhang 
6519b54502bSHong Zhang    Concepts: LU factorization, direct solver
6529b54502bSHong Zhang 
6539b54502bSHong Zhang    Notes: Usually this will compute an "exact" solution in one iteration and does
6549b54502bSHong Zhang           not need a Krylov method (i.e. you can use -ksp_type preonly, or
6559b54502bSHong Zhang           KSPSetType(ksp,KSPPREONLY) for the Krylov method
6569b54502bSHong Zhang 
6579b54502bSHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
658a4fd02acSBarry Smith            PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(),
659e5a9bf91SBarry Smith            PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetPivoting(),
660e22d95b2SBarry Smith            PCFactorSetPivotingInBlocks(),PCFactorSetShiftNonzero(),PCFactorSetShiftPd(), PCFactorReorderForNonzeroDiagonal()
6619b54502bSHong Zhang M*/
6629b54502bSHong Zhang 
6639b54502bSHong Zhang EXTERN_C_BEGIN
6649b54502bSHong Zhang #undef __FUNCT__
6659b54502bSHong Zhang #define __FUNCT__ "PCCreate_LU"
666dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_LU(PC pc)
6679b54502bSHong Zhang {
6689b54502bSHong Zhang   PetscErrorCode ierr;
6699b54502bSHong Zhang   PetscMPIInt    size;
6709b54502bSHong Zhang   PC_LU          *dir;
6719b54502bSHong Zhang 
6729b54502bSHong Zhang   PetscFunctionBegin;
67338f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_LU,&dir);CHKERRQ(ierr);
6749b54502bSHong Zhang 
6759b54502bSHong Zhang   ierr = MatFactorInfoInitialize(&dir->info);CHKERRQ(ierr);
6769b54502bSHong Zhang   dir->fact                  = 0;
6779b54502bSHong Zhang   dir->inplace               = PETSC_FALSE;
6789b54502bSHong Zhang   dir->nonzerosalongdiagonal = PETSC_FALSE;
6799b54502bSHong Zhang 
6809b54502bSHong Zhang   dir->info.fill           = 5.0;
6819b54502bSHong Zhang   dir->info.dtcol          = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */
6820a29876aSHong Zhang   dir->info.shiftnz        = 0.0;
6839b54502bSHong Zhang   dir->info.zeropivot      = 1.e-12;
6849b54502bSHong Zhang   dir->info.pivotinblocks  = 1.0;
685fbf22428SSatish Balay   dir->info.shiftpd        = 0.0; /* false */
6869b54502bSHong Zhang   dir->info.shift_fraction = 0.0;
6879b54502bSHong Zhang   dir->col                 = 0;
6889b54502bSHong Zhang   dir->row                 = 0;
6897adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
6909b54502bSHong Zhang   if (size == 1) {
69105dd20ceSSatish Balay     ierr = PetscStrallocpy(MATORDERING_ND,&dir->ordering);CHKERRQ(ierr);
6929b54502bSHong Zhang   } else {
69305dd20ceSSatish Balay     ierr = PetscStrallocpy(MATORDERING_NATURAL,&dir->ordering);CHKERRQ(ierr);
6949b54502bSHong Zhang   }
6959b54502bSHong Zhang   dir->reusefill        = PETSC_FALSE;
6969b54502bSHong Zhang   dir->reuseordering    = PETSC_FALSE;
6979b54502bSHong Zhang   pc->data              = (void*)dir;
6989b54502bSHong Zhang 
6999b54502bSHong Zhang   pc->ops->destroy           = PCDestroy_LU;
7009b54502bSHong Zhang   pc->ops->apply             = PCApply_LU;
7019b54502bSHong Zhang   pc->ops->applytranspose    = PCApplyTranspose_LU;
7029b54502bSHong Zhang   pc->ops->setup             = PCSetUp_LU;
7039b54502bSHong Zhang   pc->ops->setfromoptions    = PCSetFromOptions_LU;
7049b54502bSHong Zhang   pc->ops->view              = PCView_LU;
7059b54502bSHong Zhang   pc->ops->applyrichardson   = 0;
706a4fd02acSBarry Smith   pc->ops->getfactoredmatrix = PCFactorGetMatrix_LU;
7079b54502bSHong Zhang 
708afaefe49SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_LU",
709afaefe49SHong Zhang                     PCFactorSetZeroPivot_LU);CHKERRQ(ierr);
710afaefe49SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_LU",
711afaefe49SHong Zhang                     PCFactorSetShiftNonzero_LU);CHKERRQ(ierr);
712afaefe49SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_LU",
713afaefe49SHong Zhang                     PCFactorSetShiftPd_LU);CHKERRQ(ierr);
714afaefe49SHong Zhang 
71555ba2a51SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_LU",
71655ba2a51SBarry Smith                     PCFactorSetFill_LU);CHKERRQ(ierr);
7172401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_LU",
7182401956bSBarry Smith                     PCFactorSetUseInPlace_LU);CHKERRQ(ierr);
719e5a9bf91SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_LU",
720e5a9bf91SBarry Smith                     PCFactorSetMatOrderingType_LU);CHKERRQ(ierr);
7212401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_LU",
7222401956bSBarry Smith                     PCFactorSetReuseOrdering_LU);CHKERRQ(ierr);
7232401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_LU",
7242401956bSBarry Smith                     PCFactorSetReuseFill_LU);CHKERRQ(ierr);
7252401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivoting_C","PCFactorSetPivoting_LU",
7262401956bSBarry Smith                     PCFactorSetPivoting_LU);CHKERRQ(ierr);
7272401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_LU",
7282401956bSBarry Smith                     PCFactorSetPivotInBlocks_LU);CHKERRQ(ierr);
7292401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_LU",
7302401956bSBarry Smith                     PCFactorReorderForNonzeroDiagonal_LU);CHKERRQ(ierr);
7319b54502bSHong Zhang   PetscFunctionReturn(0);
7329b54502bSHong Zhang }
7339b54502bSHong Zhang EXTERN_C_END
734