xref: /petsc/src/ksp/pc/impls/factor/lu/lu.c (revision 719d5645761d844e4357b7ee00a3296dffe0b787)
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__
14c7393fdbSBarry Smith #define __FUNCT__ "PCFactorSetMatSolverPackage_LU"
15c7393fdbSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatSolverPackage_LU(PC pc,const MatSolverPackage stype)
16c62fe420SBarry Smith {
17c62fe420SBarry Smith   PetscErrorCode ierr;
18c62fe420SBarry Smith   PC_LU          *lu = (PC_LU*)pc->data;
19c62fe420SBarry Smith 
20c62fe420SBarry Smith   PetscFunctionBegin;
2131a0a2fdSHong Zhang   ierr = PetscStrfree(lu->solvertype);CHKERRQ(ierr);
22c62fe420SBarry Smith   ierr = PetscStrallocpy(stype,&lu->solvertype);CHKERRQ(ierr);
23c62fe420SBarry Smith   PetscFunctionReturn(0);
24c62fe420SBarry Smith }
25c62fe420SBarry Smith EXTERN_C_END
26c62fe420SBarry Smith 
27c62fe420SBarry Smith EXTERN_C_BEGIN
28c62fe420SBarry Smith #undef __FUNCT__
29afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetZeroPivot_LU"
30dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot_LU(PC pc,PetscReal z)
31afaefe49SHong Zhang {
32c62fe420SBarry Smith   PC_LU *lu = (PC_LU*)pc->data;
33afaefe49SHong Zhang 
34afaefe49SHong Zhang   PetscFunctionBegin;
35afaefe49SHong Zhang   lu->info.zeropivot = z;
36afaefe49SHong Zhang   PetscFunctionReturn(0);
37afaefe49SHong Zhang }
38afaefe49SHong Zhang EXTERN_C_END
39afaefe49SHong Zhang 
40afaefe49SHong Zhang EXTERN_C_BEGIN
41afaefe49SHong Zhang #undef __FUNCT__
42afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetShiftNonzero_LU"
43dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero_LU(PC pc,PetscReal shift)
44afaefe49SHong Zhang {
45c62fe420SBarry Smith   PC_LU *dir = (PC_LU*)pc->data;
46afaefe49SHong Zhang 
47afaefe49SHong Zhang   PetscFunctionBegin;
48afaefe49SHong Zhang   if (shift == (PetscReal) PETSC_DECIDE) {
49afaefe49SHong Zhang     dir->info.shiftnz = 1.e-12;
50afaefe49SHong Zhang   } else {
51afaefe49SHong Zhang     dir->info.shiftnz = shift;
52afaefe49SHong Zhang   }
53afaefe49SHong Zhang   PetscFunctionReturn(0);
54afaefe49SHong Zhang }
55afaefe49SHong Zhang EXTERN_C_END
56afaefe49SHong Zhang 
57afaefe49SHong Zhang EXTERN_C_BEGIN
58afaefe49SHong Zhang #undef __FUNCT__
59afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetShiftPd_LU"
60dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd_LU(PC pc,PetscTruth shift)
61afaefe49SHong Zhang {
62c62fe420SBarry Smith   PC_LU *dir = (PC_LU*)pc->data;
63afaefe49SHong Zhang 
64afaefe49SHong Zhang   PetscFunctionBegin;
65fbf22428SSatish Balay   if (shift) {
66fbf22428SSatish Balay     dir->info.shift_fraction = 0.0;
67fbf22428SSatish Balay     dir->info.shiftpd = 1.0;
68fbf22428SSatish Balay   } else {
69fbf22428SSatish Balay     dir->info.shiftpd = 0.0;
70fbf22428SSatish Balay   }
71afaefe49SHong Zhang   PetscFunctionReturn(0);
72afaefe49SHong Zhang }
73afaefe49SHong Zhang EXTERN_C_END
74afaefe49SHong Zhang 
75afaefe49SHong Zhang EXTERN_C_BEGIN
76afaefe49SHong Zhang #undef __FUNCT__
772401956bSBarry Smith #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal_LU"
782401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal_LU(PC pc,PetscReal z)
799b54502bSHong Zhang {
809b54502bSHong Zhang   PC_LU *lu = (PC_LU*)pc->data;
819b54502bSHong Zhang 
829b54502bSHong Zhang   PetscFunctionBegin;
839b54502bSHong Zhang   lu->nonzerosalongdiagonal = PETSC_TRUE;
849b54502bSHong Zhang   if (z == PETSC_DECIDE) {
859b54502bSHong Zhang     lu->nonzerosalongdiagonaltol = 1.e-10;
869b54502bSHong Zhang   } else {
879b54502bSHong Zhang     lu->nonzerosalongdiagonaltol = z;
889b54502bSHong Zhang   }
899b54502bSHong Zhang   PetscFunctionReturn(0);
909b54502bSHong Zhang }
919b54502bSHong Zhang EXTERN_C_END
929b54502bSHong Zhang 
939b54502bSHong Zhang EXTERN_C_BEGIN
949b54502bSHong Zhang #undef __FUNCT__
952401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering_LU"
962401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering_LU(PC pc,PetscTruth flag)
979b54502bSHong Zhang {
98c62fe420SBarry Smith   PC_LU *lu = (PC_LU*)pc->data;
999b54502bSHong Zhang 
1009b54502bSHong Zhang   PetscFunctionBegin;
1019b54502bSHong Zhang   lu->reuseordering = flag;
1029b54502bSHong Zhang   PetscFunctionReturn(0);
1039b54502bSHong Zhang }
1049b54502bSHong Zhang EXTERN_C_END
1059b54502bSHong Zhang 
1069b54502bSHong Zhang EXTERN_C_BEGIN
1079b54502bSHong Zhang #undef __FUNCT__
1082401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseFill_LU"
1092401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill_LU(PC pc,PetscTruth flag)
1109b54502bSHong Zhang {
111c62fe420SBarry Smith   PC_LU *lu = (PC_LU*)pc->data;
1129b54502bSHong Zhang 
1139b54502bSHong Zhang   PetscFunctionBegin;
1149b54502bSHong Zhang   lu->reusefill = flag;
1159b54502bSHong Zhang   PetscFunctionReturn(0);
1169b54502bSHong Zhang }
1179b54502bSHong Zhang EXTERN_C_END
1189b54502bSHong Zhang 
1199b54502bSHong Zhang #undef __FUNCT__
1209b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_LU"
1219b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_LU(PC pc)
1229b54502bSHong Zhang {
1239b54502bSHong Zhang   PC_LU           *lu = (PC_LU*)pc->data;
1249b54502bSHong Zhang   PetscErrorCode  ierr;
1259b54502bSHong Zhang   PetscTruth      flg,set;
1269989ab13SBarry Smith   char            tname[256], solvertype[64];
1279b54502bSHong Zhang   PetscFList      ordlist;
1289b54502bSHong Zhang   PetscReal       tol;
1299b54502bSHong Zhang 
1309b54502bSHong Zhang   PetscFunctionBegin;
1319b54502bSHong Zhang   ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr);
1329b54502bSHong Zhang   ierr = PetscOptionsHead("LU options");CHKERRQ(ierr);
1332401956bSBarry Smith     ierr = PetscOptionsName("-pc_factor_in_place","Form LU in the same memory as the matrix","PCFactorSetUseInPlace",&flg);CHKERRQ(ierr);
1349b54502bSHong Zhang     if (flg) {
1352401956bSBarry Smith       ierr = PCFactorSetUseInPlace(pc);CHKERRQ(ierr);
1369b54502bSHong Zhang     }
13755ba2a51SBarry 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);
1389b54502bSHong Zhang 
1399f95998fSHong Zhang     ierr = PetscOptionsName("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr);
1409b54502bSHong Zhang     if (flg) {
141afaefe49SHong Zhang         ierr = PCFactorSetShiftNonzero(pc,(PetscReal) PETSC_DECIDE);CHKERRQ(ierr);
1429b54502bSHong Zhang     }
1439f95998fSHong Zhang     ierr = PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",lu->info.shiftnz,&lu->info.shiftnz,0);CHKERRQ(ierr);
1449f95998fSHong Zhang     ierr = PetscOptionsName("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",&flg);CHKERRQ(ierr);
1459b54502bSHong Zhang     if (flg) {
146afaefe49SHong Zhang       ierr = PCFactorSetShiftPd(pc,PETSC_TRUE);CHKERRQ(ierr);
1479b54502bSHong Zhang     }
148ee45ca4aSHong Zhang     ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",lu->info.zeropivot,&lu->info.zeropivot,0);CHKERRQ(ierr);
1499b54502bSHong Zhang 
1502401956bSBarry Smith     ierr = PetscOptionsName("-pc_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",&flg);CHKERRQ(ierr);
1519b54502bSHong Zhang     if (flg) {
1522401956bSBarry Smith       ierr = PCFactorSetReuseFill(pc,PETSC_TRUE);CHKERRQ(ierr);
1539b54502bSHong Zhang     }
1542401956bSBarry Smith     ierr = PetscOptionsName("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",&flg);CHKERRQ(ierr);
1559b54502bSHong Zhang     if (flg) {
1562401956bSBarry Smith       ierr = PCFactorSetReuseOrdering(pc,PETSC_TRUE);CHKERRQ(ierr);
1579b54502bSHong Zhang     }
1589b54502bSHong Zhang 
1599b54502bSHong Zhang     ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr);
160e5a9bf91SBarry Smith     ierr = PetscOptionsList("-pc_factor_mat_ordering_type","Reordering to reduce nonzeros in LU","PCFactorSetMatOrderingType",ordlist,lu->ordering,tname,256,&flg);CHKERRQ(ierr);
1619b54502bSHong Zhang     if (flg) {
162e5a9bf91SBarry Smith       ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr);
1639b54502bSHong Zhang     }
1649b54502bSHong Zhang 
1655c9eb25fSBarry Smith     /* maybe should have MatGetSolverTypes(Mat,&list) like the ordering list */
166c7393fdbSBarry Smith     ierr = PetscOptionsString("-pc_factor_mat_solver_package","Specific LU solver to use","MatGetFactor",lu->solvertype,solvertype,64,&flg);CHKERRQ(ierr);
1675c9eb25fSBarry Smith     if (flg) {
168c7393fdbSBarry Smith       ierr = PCFactorSetMatSolverPackage(pc,solvertype);CHKERRQ(ierr);
1695c9eb25fSBarry Smith     }
1705c9eb25fSBarry Smith 
1712401956bSBarry Smith     ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr);
1729b54502bSHong Zhang     if (flg) {
1739b54502bSHong Zhang       tol = PETSC_DECIDE;
1742401956bSBarry Smith       ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr);
1752401956bSBarry Smith       ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr);
1769b54502bSHong Zhang     }
1779b54502bSHong Zhang 
1782401956bSBarry Smith     ierr = PetscOptionsReal("-pc_factor_pivoting","Pivoting tolerance (used only for some factorization)","PCFactorSetPivoting",lu->info.dtcol,&lu->info.dtcol,&flg);CHKERRQ(ierr);
1799b54502bSHong Zhang 
1809b54502bSHong Zhang     flg = lu->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE;
1812401956bSBarry Smith     ierr = PetscOptionsTruth("-pc_factor_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCFactorSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr);
1829b54502bSHong Zhang     if (set) {
1832401956bSBarry Smith       ierr = PCFactorSetPivotInBlocks(pc,flg);CHKERRQ(ierr);
1849b54502bSHong Zhang     }
1859b54502bSHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
1869b54502bSHong Zhang   PetscFunctionReturn(0);
1879b54502bSHong Zhang }
1889b54502bSHong Zhang 
1899b54502bSHong Zhang #undef __FUNCT__
1909b54502bSHong Zhang #define __FUNCT__ "PCView_LU"
1919b54502bSHong Zhang static PetscErrorCode PCView_LU(PC pc,PetscViewer viewer)
1929b54502bSHong Zhang {
1939b54502bSHong Zhang   PC_LU          *lu = (PC_LU*)pc->data;
1949b54502bSHong Zhang   PetscErrorCode ierr;
1959b54502bSHong Zhang   PetscTruth     iascii,isstring;
1969b54502bSHong Zhang 
1979b54502bSHong Zhang   PetscFunctionBegin;
1989b54502bSHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
1999b54502bSHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
2009b54502bSHong Zhang   if (iascii) {
2019b54502bSHong Zhang 
2029b54502bSHong Zhang     if (lu->inplace) {ierr = PetscViewerASCIIPrintf(viewer,"  LU: in-place factorization\n");CHKERRQ(ierr);}
2039b54502bSHong Zhang     else             {ierr = PetscViewerASCIIPrintf(viewer,"  LU: out-of-place factorization\n");CHKERRQ(ierr);}
2049b54502bSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"    matrix ordering: %s\n",lu->ordering);CHKERRQ(ierr);
205a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  LU: tolerance for zero pivot %G\n",lu->info.zeropivot);CHKERRQ(ierr);
2060a29876aSHong Zhang     if (lu->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer,"  LU: using Manteuffel shift\n");CHKERRQ(ierr);}
2079b54502bSHong Zhang     if (lu->fact) {
208f3a39becSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  LU: factor fill ratio needed %G\n",lu->actualfill);CHKERRQ(ierr);
209f3a39becSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"       Factored matrix follows\n");CHKERRQ(ierr);
210f3a39becSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
211f3a39becSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
212f3a39becSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
213f3a39becSBarry Smith       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
2149b54502bSHong Zhang       ierr = MatView(lu->fact,viewer);CHKERRQ(ierr);
2159b54502bSHong Zhang       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
216f3a39becSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
217f3a39becSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
218f3a39becSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
2199b54502bSHong Zhang     }
2209b54502bSHong Zhang     if (lu->reusefill)    {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing fill from past factorization\n");CHKERRQ(ierr);}
2219b54502bSHong Zhang     if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing reordering from past factorization\n");CHKERRQ(ierr);}
2229b54502bSHong Zhang   } else if (isstring) {
2239b54502bSHong Zhang     ierr = PetscViewerStringSPrintf(viewer," order=%s",lu->ordering);CHKERRQ(ierr);CHKERRQ(ierr);
2249b54502bSHong Zhang   } else {
2259b54502bSHong Zhang     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCLU",((PetscObject)viewer)->type_name);
2269b54502bSHong Zhang   }
2279b54502bSHong Zhang   PetscFunctionReturn(0);
2289b54502bSHong Zhang }
2299b54502bSHong Zhang 
2309b54502bSHong Zhang #undef __FUNCT__
231a4fd02acSBarry Smith #define __FUNCT__ "PCFactorGetMatrix_LU"
232a4fd02acSBarry Smith static PetscErrorCode PCFactorGetMatrix_LU(PC pc,Mat *mat)
2339b54502bSHong Zhang {
2349b54502bSHong Zhang   PC_LU *dir = (PC_LU*)pc->data;
2359b54502bSHong Zhang 
2369b54502bSHong Zhang   PetscFunctionBegin;
2379b54502bSHong Zhang   if (!dir->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()");
2389b54502bSHong Zhang   *mat = dir->fact;
2399b54502bSHong Zhang   PetscFunctionReturn(0);
2409b54502bSHong Zhang }
2419b54502bSHong Zhang 
2429b54502bSHong Zhang #undef __FUNCT__
2439b54502bSHong Zhang #define __FUNCT__ "PCSetUp_LU"
2449b54502bSHong Zhang static PetscErrorCode PCSetUp_LU(PC pc)
2459b54502bSHong Zhang {
2469b54502bSHong Zhang   PetscErrorCode ierr;
2479b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
2489b54502bSHong Zhang 
2499b54502bSHong Zhang   PetscFunctionBegin;
2509b54502bSHong Zhang   if (dir->reusefill && pc->setupcalled) dir->info.fill = dir->actualfill;
2519b54502bSHong Zhang 
2529b54502bSHong Zhang   if (dir->inplace) {
2539b54502bSHong Zhang     if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);}
2549b54502bSHong Zhang     if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);}
2559b54502bSHong Zhang     ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
25603c60df9SBarry Smith     if (dir->row) {
25703c60df9SBarry Smith       ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);
25803c60df9SBarry Smith       ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr);
25903c60df9SBarry Smith     }
2609b54502bSHong Zhang     ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&dir->info);CHKERRQ(ierr);
2619b54502bSHong Zhang     dir->fact = pc->pmat;
2629b54502bSHong Zhang   } else {
2639b54502bSHong Zhang     MatInfo info;
2649b54502bSHong Zhang     if (!pc->setupcalled) {
2659b54502bSHong Zhang       ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
2669b54502bSHong Zhang       if (dir->nonzerosalongdiagonal) {
2679b54502bSHong Zhang         ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr);
2689b54502bSHong Zhang       }
26903c60df9SBarry Smith       if (dir->row) {
27003c60df9SBarry Smith         ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);
27103c60df9SBarry Smith         ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr);
27203c60df9SBarry Smith       }
2735c9eb25fSBarry Smith       ierr = MatGetFactor(pc->pmat,dir->solvertype,MAT_FACTOR_LU,&dir->fact);CHKERRQ(ierr);
274*719d5645SBarry Smith       ierr = MatLUFactorSymbolic(dir->fact,pc->pmat,dir->row,dir->col,&dir->info);CHKERRQ(ierr);
2759b54502bSHong Zhang       ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
2769b54502bSHong Zhang       dir->actualfill = info.fill_ratio_needed;
27752e6d16bSBarry Smith       ierr = PetscLogObjectParent(pc,dir->fact);CHKERRQ(ierr);
2789b54502bSHong Zhang     } else if (pc->flag != SAME_NONZERO_PATTERN) {
2799b54502bSHong Zhang       if (!dir->reuseordering) {
2809b54502bSHong Zhang         if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);}
2819b54502bSHong Zhang         if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);}
2829b54502bSHong Zhang         ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
2839b54502bSHong Zhang         if (dir->nonzerosalongdiagonal) {
2849b54502bSHong Zhang           ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr);
2859b54502bSHong Zhang         }
28603c60df9SBarry Smith         if (dir->row) {
28703c60df9SBarry Smith           ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);
28803c60df9SBarry Smith           ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr);
28903c60df9SBarry Smith         }
2909b54502bSHong Zhang       }
2919b54502bSHong Zhang       ierr = MatDestroy(dir->fact);CHKERRQ(ierr);
2925c9eb25fSBarry Smith       ierr = MatGetFactor(pc->pmat,dir->solvertype,MAT_FACTOR_LU,&dir->fact);CHKERRQ(ierr);
293*719d5645SBarry Smith       ierr = MatLUFactorSymbolic(dir->fact,pc->pmat,dir->row,dir->col,&dir->info);CHKERRQ(ierr);
2949b54502bSHong Zhang       ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
2959b54502bSHong Zhang       dir->actualfill = info.fill_ratio_needed;
29652e6d16bSBarry Smith       ierr = PetscLogObjectParent(pc,dir->fact);CHKERRQ(ierr);
2979b54502bSHong Zhang     }
298*719d5645SBarry Smith     ierr = MatLUFactorNumeric(dir->fact,pc->pmat,&dir->info);CHKERRQ(ierr);
2999b54502bSHong Zhang   }
3009b54502bSHong Zhang   PetscFunctionReturn(0);
3019b54502bSHong Zhang }
3029b54502bSHong Zhang 
3039b54502bSHong Zhang #undef __FUNCT__
3049b54502bSHong Zhang #define __FUNCT__ "PCDestroy_LU"
3059b54502bSHong Zhang static PetscErrorCode PCDestroy_LU(PC pc)
3069b54502bSHong Zhang {
3079b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
3089b54502bSHong Zhang   PetscErrorCode ierr;
3099b54502bSHong Zhang 
3109b54502bSHong Zhang   PetscFunctionBegin;
3119b54502bSHong Zhang   if (!dir->inplace && dir->fact) {ierr = MatDestroy(dir->fact);CHKERRQ(ierr);}
3129b54502bSHong Zhang   if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);}
3139b54502bSHong Zhang   if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);}
3149b54502bSHong Zhang   ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr);
3155c9eb25fSBarry Smith   ierr = PetscStrfree(dir->solvertype);CHKERRQ(ierr);
3169b54502bSHong Zhang   ierr = PetscFree(dir);CHKERRQ(ierr);
3179b54502bSHong Zhang   PetscFunctionReturn(0);
3189b54502bSHong Zhang }
3199b54502bSHong Zhang 
3209b54502bSHong Zhang #undef __FUNCT__
3219b54502bSHong Zhang #define __FUNCT__ "PCApply_LU"
3229b54502bSHong Zhang static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y)
3239b54502bSHong Zhang {
3249b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
3259b54502bSHong Zhang   PetscErrorCode ierr;
3269b54502bSHong Zhang 
3279b54502bSHong Zhang   PetscFunctionBegin;
3289b54502bSHong Zhang   if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);}
3299b54502bSHong Zhang   else              {ierr = MatSolve(dir->fact,x,y);CHKERRQ(ierr);}
3309b54502bSHong Zhang   PetscFunctionReturn(0);
3319b54502bSHong Zhang }
3329b54502bSHong Zhang 
3339b54502bSHong Zhang #undef __FUNCT__
3349b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_LU"
3359b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y)
3369b54502bSHong Zhang {
3379b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
3389b54502bSHong Zhang   PetscErrorCode ierr;
3399b54502bSHong Zhang 
3409b54502bSHong Zhang   PetscFunctionBegin;
3419b54502bSHong Zhang   if (dir->inplace) {ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);}
3429b54502bSHong Zhang   else              {ierr = MatSolveTranspose(dir->fact,x,y);CHKERRQ(ierr);}
3439b54502bSHong Zhang   PetscFunctionReturn(0);
3449b54502bSHong Zhang }
3459b54502bSHong Zhang 
3469b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/
3479b54502bSHong Zhang 
3489b54502bSHong Zhang EXTERN_C_BEGIN
3499b54502bSHong Zhang #undef __FUNCT__
35055ba2a51SBarry Smith #define __FUNCT__ "PCFactorSetFill_LU"
35155ba2a51SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill_LU(PC pc,PetscReal fill)
3529b54502bSHong Zhang {
353c62fe420SBarry Smith   PC_LU *dir = (PC_LU*)pc->data;
3549b54502bSHong Zhang 
3559b54502bSHong Zhang   PetscFunctionBegin;
3569b54502bSHong Zhang   dir->info.fill = fill;
3579b54502bSHong Zhang   PetscFunctionReturn(0);
3589b54502bSHong Zhang }
3599b54502bSHong Zhang EXTERN_C_END
3609b54502bSHong Zhang 
3619b54502bSHong Zhang EXTERN_C_BEGIN
3629b54502bSHong Zhang #undef __FUNCT__
3632401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace_LU"
3642401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_LU(PC pc)
3659b54502bSHong Zhang {
366c62fe420SBarry Smith   PC_LU *dir = (PC_LU*)pc->data;
3679b54502bSHong Zhang 
3689b54502bSHong Zhang   PetscFunctionBegin;
3699b54502bSHong Zhang   dir->inplace = PETSC_TRUE;
3709b54502bSHong Zhang   PetscFunctionReturn(0);
3719b54502bSHong Zhang }
3729b54502bSHong Zhang EXTERN_C_END
3739b54502bSHong Zhang 
3749b54502bSHong Zhang EXTERN_C_BEGIN
3759b54502bSHong Zhang #undef __FUNCT__
376e5a9bf91SBarry Smith #define __FUNCT__ "PCFactorSetMatOrderingType_LU"
377e36ff184SSatish Balay PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType_LU(PC pc,const MatOrderingType ordering)
3789b54502bSHong Zhang {
3799b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
3809b54502bSHong Zhang   PetscErrorCode ierr;
3819b54502bSHong Zhang 
3829b54502bSHong Zhang   PetscFunctionBegin;
3839b54502bSHong Zhang   ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr);
3849b54502bSHong Zhang   ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr);
3859b54502bSHong Zhang   PetscFunctionReturn(0);
3869b54502bSHong Zhang }
3879b54502bSHong Zhang EXTERN_C_END
3889b54502bSHong Zhang 
3899b54502bSHong Zhang EXTERN_C_BEGIN
3909b54502bSHong Zhang #undef __FUNCT__
3912401956bSBarry Smith #define __FUNCT__ "PCFactorSetPivoting_LU"
3922401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivoting_LU(PC pc,PetscReal dtcol)
3939b54502bSHong Zhang {
3949b54502bSHong Zhang   PC_LU *dir = (PC_LU*)pc->data;
3959b54502bSHong Zhang 
3969b54502bSHong Zhang   PetscFunctionBegin;
397a83599f4SBarry 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);
3989b54502bSHong Zhang   dir->info.dtcol = dtcol;
3999b54502bSHong Zhang   PetscFunctionReturn(0);
4009b54502bSHong Zhang }
4019b54502bSHong Zhang EXTERN_C_END
4029b54502bSHong Zhang 
4039b54502bSHong Zhang EXTERN_C_BEGIN
4049b54502bSHong Zhang #undef __FUNCT__
4052401956bSBarry Smith #define __FUNCT__ "PCFactorSetPivotInBlocks_LU"
4062401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks_LU(PC pc,PetscTruth pivot)
4079b54502bSHong Zhang {
4089b54502bSHong Zhang   PC_LU *dir = (PC_LU*)pc->data;
4099b54502bSHong Zhang 
4109b54502bSHong Zhang   PetscFunctionBegin;
4119b54502bSHong Zhang   dir->info.pivotinblocks = pivot ? 1.0 : 0.0;
4129b54502bSHong Zhang   PetscFunctionReturn(0);
4139b54502bSHong Zhang }
4149b54502bSHong Zhang EXTERN_C_END
4159b54502bSHong Zhang 
4169b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/
4179b54502bSHong Zhang 
4189b54502bSHong Zhang #undef __FUNCT__
4192401956bSBarry Smith #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal"
4209b54502bSHong Zhang /*@
4212401956bSBarry Smith    PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal
4229b54502bSHong Zhang 
4239b54502bSHong Zhang    Collective on PC
4249b54502bSHong Zhang 
4259b54502bSHong Zhang    Input Parameters:
4269b54502bSHong Zhang +  pc - the preconditioner context
4279b54502bSHong Zhang -  tol - diagonal entries smaller than this in absolute value are considered zero
4289b54502bSHong Zhang 
4299b54502bSHong Zhang    Options Database Key:
4302401956bSBarry Smith .  -pc_factor_nonzeros_along_diagonal
4319b54502bSHong Zhang 
4329b54502bSHong Zhang    Level: intermediate
4339b54502bSHong Zhang 
4349b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill
4359b54502bSHong Zhang 
43655ba2a51SBarry Smith .seealso: PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetZeroPivot(), MatReorderForNonzeroDiagonal()
4379b54502bSHong Zhang @*/
4382401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal(PC pc,PetscReal rtol)
4399b54502bSHong Zhang {
4409b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscReal);
4419b54502bSHong Zhang 
4429b54502bSHong Zhang   PetscFunctionBegin;
4439b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
4442401956bSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",(void (**)(void))&f);CHKERRQ(ierr);
4459b54502bSHong Zhang   if (f) {
4469b54502bSHong Zhang     ierr = (*f)(pc,rtol);CHKERRQ(ierr);
4479b54502bSHong Zhang   }
4489b54502bSHong Zhang   PetscFunctionReturn(0);
4499b54502bSHong Zhang }
4509b54502bSHong Zhang 
4519b54502bSHong Zhang #undef __FUNCT__
452c7393fdbSBarry Smith #define __FUNCT__ "PCFactorSetMatSolverPackage"
453c62fe420SBarry Smith /*@
454c7393fdbSBarry Smith    PCFactorSetMatSolverPackage - sets the software that is used to perform the factorization
455c62fe420SBarry Smith 
456c62fe420SBarry Smith    Collective on PC
457c62fe420SBarry Smith 
458c62fe420SBarry Smith    Input Parameters:
459c62fe420SBarry Smith +  pc - the preconditioner context
460c62fe420SBarry Smith -  stype - for example, spooles, superlu, superlu_d
461c62fe420SBarry Smith 
462c62fe420SBarry Smith    Options Database Key:
463c7393fdbSBarry Smith .  -pc_factor_mat_solver_package <stype> - spooles, petsc, superlu, superlu_dist, mumps
464c62fe420SBarry Smith 
465c62fe420SBarry Smith    Level: intermediate
466c62fe420SBarry Smith 
467c62fe420SBarry Smith    Note:
468c62fe420SBarry Smith      By default this will use the PETSc factorization if it exists
469c62fe420SBarry Smith 
470c62fe420SBarry Smith .keywords: PC, set, factorization, direct, fill
471c62fe420SBarry Smith 
472c7393fdbSBarry Smith .seealso: MatGetFactor(), MatSolverPackage
473c62fe420SBarry Smith 
474c62fe420SBarry Smith @*/
475c7393fdbSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatSolverPackage(PC pc,const MatSolverPackage stype)
476c62fe420SBarry Smith {
477c7393fdbSBarry Smith   PetscErrorCode ierr,(*f)(PC,const MatSolverPackage);
478c62fe420SBarry Smith 
479c62fe420SBarry Smith   PetscFunctionBegin;
480c62fe420SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
481c7393fdbSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",(void (**)(void))&f);CHKERRQ(ierr);
482c62fe420SBarry Smith   if (f) {
483c62fe420SBarry Smith     ierr = (*f)(pc,stype);CHKERRQ(ierr);
484c62fe420SBarry Smith   }
485c62fe420SBarry Smith   PetscFunctionReturn(0);
486c62fe420SBarry Smith }
487c62fe420SBarry Smith 
488c62fe420SBarry Smith #undef __FUNCT__
48955ba2a51SBarry Smith #define __FUNCT__ "PCFactorSetFill"
4909b54502bSHong Zhang /*@
49155ba2a51SBarry Smith    PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix,
4929b54502bSHong Zhang    fill = number nonzeros in factor/number nonzeros in original matrix.
4939b54502bSHong Zhang 
4949b54502bSHong Zhang    Collective on PC
4959b54502bSHong Zhang 
4969b54502bSHong Zhang    Input Parameters:
4979b54502bSHong Zhang +  pc - the preconditioner context
4989b54502bSHong Zhang -  fill - amount of expected fill
4999b54502bSHong Zhang 
5009b54502bSHong Zhang    Options Database Key:
50155ba2a51SBarry Smith .  -pc_factor_fill <fill> - Sets fill amount
5029b54502bSHong Zhang 
5039b54502bSHong Zhang    Level: intermediate
5049b54502bSHong Zhang 
5059b54502bSHong Zhang    Note:
5069b54502bSHong Zhang    For sparse matrix factorizations it is difficult to predict how much
5076cf91177SBarry Smith    fill to expect. By running with the option -info PETSc will print the
5089b54502bSHong Zhang    actual amount of fill used; allowing you to set the value accurately for
5099b54502bSHong Zhang    future runs. Default PETSc uses a value of 5.0
5109b54502bSHong Zhang 
5119b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill
5129b54502bSHong Zhang 
5139b54502bSHong Zhang @*/
51455ba2a51SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill(PC pc,PetscReal fill)
5159b54502bSHong Zhang {
5169b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscReal);
5179b54502bSHong Zhang 
5189b54502bSHong Zhang   PetscFunctionBegin;
5199b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
5209b54502bSHong Zhang   if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0");
52155ba2a51SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetFill_C",(void (**)(void))&f);CHKERRQ(ierr);
5229b54502bSHong Zhang   if (f) {
5239b54502bSHong Zhang     ierr = (*f)(pc,fill);CHKERRQ(ierr);
5249b54502bSHong Zhang   }
5259b54502bSHong Zhang   PetscFunctionReturn(0);
5269b54502bSHong Zhang }
5279b54502bSHong Zhang 
5289b54502bSHong Zhang #undef __FUNCT__
5292401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace"
5309b54502bSHong Zhang /*@
5312401956bSBarry Smith    PCFactorSetUseInPlace - Tells the system to do an in-place factorization.
5329b54502bSHong Zhang    For dense matrices, this enables the solution of much larger problems.
5339b54502bSHong Zhang    For sparse matrices the factorization cannot be done truly in-place
5349b54502bSHong Zhang    so this does not save memory during the factorization, but after the matrix
5359b54502bSHong Zhang    is factored, the original unfactored matrix is freed, thus recovering that
5369b54502bSHong Zhang    space.
5379b54502bSHong Zhang 
5389b54502bSHong Zhang    Collective on PC
5399b54502bSHong Zhang 
5409b54502bSHong Zhang    Input Parameters:
5419b54502bSHong Zhang .  pc - the preconditioner context
5429b54502bSHong Zhang 
5439b54502bSHong Zhang    Options Database Key:
5442401956bSBarry Smith .  -pc_factor_in_place - Activates in-place factorization
5459b54502bSHong Zhang 
5469b54502bSHong Zhang    Notes:
5472401956bSBarry Smith    PCFactorSetUseInplace() can only be used with the KSP method KSPPREONLY or when
5489b54502bSHong Zhang    a different matrix is provided for the multiply and the preconditioner in
5499b54502bSHong Zhang    a call to KSPSetOperators().
5509b54502bSHong Zhang    This is because the Krylov space methods require an application of the
5519b54502bSHong Zhang    matrix multiplication, which is not possible here because the matrix has
5529b54502bSHong Zhang    been factored in-place, replacing the original matrix.
5539b54502bSHong Zhang 
5549b54502bSHong Zhang    Level: intermediate
5559b54502bSHong Zhang 
5569b54502bSHong Zhang .keywords: PC, set, factorization, direct, inplace, in-place, LU
5579b54502bSHong Zhang 
5589b54502bSHong Zhang .seealso: PCILUSetUseInPlace()
5599b54502bSHong Zhang @*/
5602401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace(PC pc)
5619b54502bSHong Zhang {
5629b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC);
5639b54502bSHong Zhang 
5649b54502bSHong Zhang   PetscFunctionBegin;
5659b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
5662401956bSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",(void (**)(void))&f);CHKERRQ(ierr);
5679b54502bSHong Zhang   if (f) {
5689b54502bSHong Zhang     ierr = (*f)(pc);CHKERRQ(ierr);
5699b54502bSHong Zhang   }
5709b54502bSHong Zhang   PetscFunctionReturn(0);
5719b54502bSHong Zhang }
5729b54502bSHong Zhang 
5739b54502bSHong Zhang #undef __FUNCT__
574e5a9bf91SBarry Smith #define __FUNCT__ "PCFactorSetMatOrderingType"
5759b54502bSHong Zhang /*@C
576e5a9bf91SBarry Smith     PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to
5779b54502bSHong Zhang     be used in the LU factorization.
5789b54502bSHong Zhang 
5799b54502bSHong Zhang     Collective on PC
5809b54502bSHong Zhang 
5819b54502bSHong Zhang     Input Parameters:
5829b54502bSHong Zhang +   pc - the preconditioner context
5839b54502bSHong Zhang -   ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM
5849b54502bSHong Zhang 
5859b54502bSHong Zhang     Options Database Key:
5862401956bSBarry Smith .   -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
5879b54502bSHong Zhang 
5889b54502bSHong Zhang     Level: intermediate
5899b54502bSHong Zhang 
5909b54502bSHong Zhang     Notes: nested dissection is used by default
5919b54502bSHong Zhang 
5922b6de112SBarry Smith     For Cholesky and ICC and the SBAIJ format reorderings are not available,
5932b6de112SBarry Smith     since only the upper triangular part of the matrix is stored. You can use the
5942b6de112SBarry Smith     SeqAIJ format in this case to get reorderings.
5952b6de112SBarry Smith 
5969b54502bSHong Zhang @*/
597e36ff184SSatish Balay PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType(PC pc,const MatOrderingType ordering)
5989b54502bSHong Zhang {
599e36ff184SSatish Balay   PetscErrorCode ierr,(*f)(PC,const MatOrderingType);
6009b54502bSHong Zhang 
6019b54502bSHong Zhang   PetscFunctionBegin;
602e5a9bf91SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",(void (**)(void))&f);CHKERRQ(ierr);
6039b54502bSHong Zhang   if (f) {
6049b54502bSHong Zhang     ierr = (*f)(pc,ordering);CHKERRQ(ierr);
6059b54502bSHong Zhang   }
6069b54502bSHong Zhang   PetscFunctionReturn(0);
6079b54502bSHong Zhang }
6089b54502bSHong Zhang 
6099b54502bSHong Zhang #undef __FUNCT__
6102401956bSBarry Smith #define __FUNCT__ "PCFactorSetPivoting"
6119b54502bSHong Zhang /*@
6122401956bSBarry Smith     PCFactorSetPivoting - Determines when pivoting is done during LU.
6139b54502bSHong Zhang       For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices
6149b54502bSHong Zhang       it is never done. For the Matlab and SuperLU factorization this is used.
6159b54502bSHong Zhang 
6169b54502bSHong Zhang     Collective on PC
6179b54502bSHong Zhang 
6189b54502bSHong Zhang     Input Parameters:
6199b54502bSHong Zhang +   pc - the preconditioner context
6209b54502bSHong Zhang -   dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable)
6219b54502bSHong Zhang 
6229b54502bSHong Zhang     Options Database Key:
6232401956bSBarry Smith .   -pc_factor_pivoting <dtcol>
6249b54502bSHong Zhang 
6259b54502bSHong Zhang     Level: intermediate
6269b54502bSHong Zhang 
6272401956bSBarry Smith .seealso: PCILUSetMatOrdering(), PCFactorSetPivotInBlocks()
6289b54502bSHong Zhang @*/
6292401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivoting(PC pc,PetscReal dtcol)
6309b54502bSHong Zhang {
6319b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscReal);
6329b54502bSHong Zhang 
6339b54502bSHong Zhang   PetscFunctionBegin;
6342401956bSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivoting_C",(void (**)(void))&f);CHKERRQ(ierr);
6359b54502bSHong Zhang   if (f) {
6369b54502bSHong Zhang     ierr = (*f)(pc,dtcol);CHKERRQ(ierr);
6379b54502bSHong Zhang   }
6389b54502bSHong Zhang   PetscFunctionReturn(0);
6399b54502bSHong Zhang }
6409b54502bSHong Zhang 
6419b54502bSHong Zhang #undef __FUNCT__
6422401956bSBarry Smith #define __FUNCT__ "PCFactorSetPivotInBlocks"
6439b54502bSHong Zhang /*@
6442401956bSBarry Smith     PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block
6459b54502bSHong Zhang       with BAIJ or SBAIJ matrices
6469b54502bSHong Zhang 
6479b54502bSHong Zhang     Collective on PC
6489b54502bSHong Zhang 
6499b54502bSHong Zhang     Input Parameters:
6509b54502bSHong Zhang +   pc - the preconditioner context
6519b54502bSHong Zhang -   pivot - PETSC_TRUE or PETSC_FALSE
6529b54502bSHong Zhang 
6539b54502bSHong Zhang     Options Database Key:
6542401956bSBarry Smith .   -pc_factor_pivot_in_blocks <true,false>
6559b54502bSHong Zhang 
6569b54502bSHong Zhang     Level: intermediate
6579b54502bSHong Zhang 
6582401956bSBarry Smith .seealso: PCILUSetMatOrdering(), PCFactorSetPivoting()
6599b54502bSHong Zhang @*/
6602401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks(PC pc,PetscTruth pivot)
6619b54502bSHong Zhang {
6629b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscTruth);
6639b54502bSHong Zhang 
6649b54502bSHong Zhang   PetscFunctionBegin;
6652401956bSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr);
6669b54502bSHong Zhang   if (f) {
6679b54502bSHong Zhang     ierr = (*f)(pc,pivot);CHKERRQ(ierr);
6689b54502bSHong Zhang   }
6699b54502bSHong Zhang   PetscFunctionReturn(0);
6709b54502bSHong Zhang }
6719b54502bSHong Zhang 
6729b54502bSHong Zhang /* ------------------------------------------------------------------------ */
6739b54502bSHong Zhang 
6749b54502bSHong Zhang /*MC
6759b54502bSHong Zhang    PCLU - Uses a direct solver, based on LU factorization, as a preconditioner
6769b54502bSHong Zhang 
6779b54502bSHong Zhang    Options Database Keys:
6782401956bSBarry Smith +  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
679c7393fdbSBarry Smith .  -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like spooles
6802401956bSBarry Smith .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
68155ba2a51SBarry Smith .  -pc_factor_fill <fill> - Sets fill amount
6822401956bSBarry Smith .  -pc_factor_in_place - Activates in-place factorization
6832401956bSBarry Smith .  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
6842401956bSBarry Smith .  -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase
6859b54502bSHong Zhang                                          stability of factorization.
686f251bdbdSHong Zhang .  -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default
687e22d95b2SBarry Smith .  -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value
688f251bdbdSHong Zhang         is optional with PETSC_TRUE being the default
689e22d95b2SBarry Smith -   -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the
690e22d95b2SBarry Smith         diagonal.
6919b54502bSHong Zhang 
6929b54502bSHong Zhang    Notes: Not all options work for all matrix formats
6939b54502bSHong Zhang           Run with -help to see additional options for particular matrix formats or factorization
6949b54502bSHong Zhang           algorithms
6959b54502bSHong Zhang 
6969b54502bSHong Zhang    Level: beginner
6979b54502bSHong Zhang 
6989b54502bSHong Zhang    Concepts: LU factorization, direct solver
6999b54502bSHong Zhang 
7009b54502bSHong Zhang    Notes: Usually this will compute an "exact" solution in one iteration and does
7019b54502bSHong Zhang           not need a Krylov method (i.e. you can use -ksp_type preonly, or
7029b54502bSHong Zhang           KSPSetType(ksp,KSPPREONLY) for the Krylov method
7039b54502bSHong Zhang 
7049b54502bSHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
705a4fd02acSBarry Smith            PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(),
706e5a9bf91SBarry Smith            PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetPivoting(),
707e22d95b2SBarry Smith            PCFactorSetPivotingInBlocks(),PCFactorSetShiftNonzero(),PCFactorSetShiftPd(), PCFactorReorderForNonzeroDiagonal()
7089b54502bSHong Zhang M*/
7099b54502bSHong Zhang 
7109b54502bSHong Zhang EXTERN_C_BEGIN
7119b54502bSHong Zhang #undef __FUNCT__
7129b54502bSHong Zhang #define __FUNCT__ "PCCreate_LU"
713dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_LU(PC pc)
7149b54502bSHong Zhang {
7159b54502bSHong Zhang   PetscErrorCode ierr;
7169b54502bSHong Zhang   PetscMPIInt    size;
7179b54502bSHong Zhang   PC_LU          *dir;
7189b54502bSHong Zhang 
7199b54502bSHong Zhang   PetscFunctionBegin;
72038f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_LU,&dir);CHKERRQ(ierr);
7219b54502bSHong Zhang 
7229b54502bSHong Zhang   ierr = MatFactorInfoInitialize(&dir->info);CHKERRQ(ierr);
7239b54502bSHong Zhang   dir->fact                  = 0;
7249b54502bSHong Zhang   dir->inplace               = PETSC_FALSE;
7259b54502bSHong Zhang   dir->nonzerosalongdiagonal = PETSC_FALSE;
7269b54502bSHong Zhang 
7279b54502bSHong Zhang   dir->info.fill           = 5.0;
7289b54502bSHong Zhang   dir->info.dtcol          = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */
7290a29876aSHong Zhang   dir->info.shiftnz        = 0.0;
7309b54502bSHong Zhang   dir->info.zeropivot      = 1.e-12;
7319b54502bSHong Zhang   dir->info.pivotinblocks  = 1.0;
732fbf22428SSatish Balay   dir->info.shiftpd        = 0.0; /* false */
7339b54502bSHong Zhang   dir->info.shift_fraction = 0.0;
7349b54502bSHong Zhang   dir->col                 = 0;
7359b54502bSHong Zhang   dir->row                 = 0;
7365c9eb25fSBarry Smith 
73743244d56SBarry Smith   ierr = PetscStrallocpy(MAT_SOLVER_PETSC,&dir->solvertype);CHKERRQ(ierr);
7387adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
7399b54502bSHong Zhang   if (size == 1) {
7409b54502bSHong Zhang     ierr = PetscStrallocpy(MATORDERING_ND,&dir->ordering);CHKERRQ(ierr);
7419b54502bSHong Zhang   } else {
7429b54502bSHong Zhang     ierr = PetscStrallocpy(MATORDERING_NATURAL,&dir->ordering);CHKERRQ(ierr);
7439b54502bSHong Zhang   }
7449b54502bSHong Zhang   dir->reusefill        = PETSC_FALSE;
7459b54502bSHong Zhang   dir->reuseordering    = PETSC_FALSE;
7469b54502bSHong Zhang   pc->data              = (void*)dir;
7479b54502bSHong Zhang 
7489b54502bSHong Zhang   pc->ops->destroy           = PCDestroy_LU;
7499b54502bSHong Zhang   pc->ops->apply             = PCApply_LU;
7509b54502bSHong Zhang   pc->ops->applytranspose    = PCApplyTranspose_LU;
7519b54502bSHong Zhang   pc->ops->setup             = PCSetUp_LU;
7529b54502bSHong Zhang   pc->ops->setfromoptions    = PCSetFromOptions_LU;
7539b54502bSHong Zhang   pc->ops->view              = PCView_LU;
7549b54502bSHong Zhang   pc->ops->applyrichardson   = 0;
755a4fd02acSBarry Smith   pc->ops->getfactoredmatrix = PCFactorGetMatrix_LU;
7569b54502bSHong Zhang 
757c7393fdbSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_LU",
758c7393fdbSBarry Smith                     PCFactorSetMatSolverPackage_LU);CHKERRQ(ierr);
759afaefe49SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_LU",
760afaefe49SHong Zhang                     PCFactorSetZeroPivot_LU);CHKERRQ(ierr);
761afaefe49SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_LU",
762afaefe49SHong Zhang                     PCFactorSetShiftNonzero_LU);CHKERRQ(ierr);
763afaefe49SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_LU",
764afaefe49SHong Zhang                     PCFactorSetShiftPd_LU);CHKERRQ(ierr);
765afaefe49SHong Zhang 
76655ba2a51SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_LU",
76755ba2a51SBarry Smith                     PCFactorSetFill_LU);CHKERRQ(ierr);
7682401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_LU",
7692401956bSBarry Smith                     PCFactorSetUseInPlace_LU);CHKERRQ(ierr);
770e5a9bf91SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_LU",
771e5a9bf91SBarry Smith                     PCFactorSetMatOrderingType_LU);CHKERRQ(ierr);
7722401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_LU",
7732401956bSBarry Smith                     PCFactorSetReuseOrdering_LU);CHKERRQ(ierr);
7742401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_LU",
7752401956bSBarry Smith                     PCFactorSetReuseFill_LU);CHKERRQ(ierr);
7762401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivoting_C","PCFactorSetPivoting_LU",
7772401956bSBarry Smith                     PCFactorSetPivoting_LU);CHKERRQ(ierr);
7782401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_LU",
7792401956bSBarry Smith                     PCFactorSetPivotInBlocks_LU);CHKERRQ(ierr);
7802401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_LU",
7812401956bSBarry Smith                     PCFactorReorderForNonzeroDiagonal_LU);CHKERRQ(ierr);
7829b54502bSHong Zhang   PetscFunctionReturn(0);
7839b54502bSHong Zhang }
7849b54502bSHong Zhang EXTERN_C_END
785