xref: /petsc/src/ksp/pc/impls/factor/lu/lu.c (revision 55ba2a515af5dbbcf670f4e560eb9724ec139bf1)
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__
659b54502bSHong Zhang #define __FUNCT__ "PCLUReorderForNonzeroDiagonal_LU"
66dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCLUReorderForNonzeroDiagonal_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__
839b54502bSHong Zhang #define __FUNCT__ "PCLUSetReuseOrdering_LU"
84dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCLUSetReuseOrdering_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__
979b54502bSHong Zhang #define __FUNCT__ "PCLUSetReuseFill_LU"
98dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCLUSetReuseFill_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;
1169b54502bSHong Zhang   char           tname[256];
1179b54502bSHong Zhang   PetscFList     ordlist;
1189b54502bSHong Zhang   PetscReal      tol;
1199b54502bSHong Zhang 
1209b54502bSHong Zhang   PetscFunctionBegin;
1219b54502bSHong Zhang   ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr);
1229b54502bSHong Zhang   ierr = PetscOptionsHead("LU options");CHKERRQ(ierr);
1239b54502bSHong Zhang     ierr = PetscOptionsName("-pc_lu_in_place","Form LU in the same memory as the matrix","PCLUSetUseInPlace",&flg);CHKERRQ(ierr);
1249b54502bSHong Zhang     if (flg) {
1259b54502bSHong Zhang       ierr = PCLUSetUseInPlace(pc);CHKERRQ(ierr);
1269b54502bSHong Zhang     }
127*55ba2a51SBarry 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);
1289b54502bSHong Zhang 
1299f95998fSHong Zhang     ierr = PetscOptionsName("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr);
1309b54502bSHong Zhang     if (flg) {
131afaefe49SHong Zhang         ierr = PCFactorSetShiftNonzero(pc,(PetscReal) PETSC_DECIDE);CHKERRQ(ierr);
1329b54502bSHong Zhang     }
1339f95998fSHong Zhang     ierr = PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",lu->info.shiftnz,&lu->info.shiftnz,0);CHKERRQ(ierr);
1349f95998fSHong Zhang     ierr = PetscOptionsName("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",&flg);CHKERRQ(ierr);
1359b54502bSHong Zhang     if (flg) {
136afaefe49SHong Zhang       ierr = PCFactorSetShiftPd(pc,PETSC_TRUE);CHKERRQ(ierr);
1379b54502bSHong Zhang     }
138ee45ca4aSHong Zhang     ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",lu->info.zeropivot,&lu->info.zeropivot,0);CHKERRQ(ierr);
1399b54502bSHong Zhang 
1409b54502bSHong Zhang     ierr = PetscOptionsName("-pc_lu_reuse_fill","Use fill from previous factorization","PCLUSetReuseFill",&flg);CHKERRQ(ierr);
1419b54502bSHong Zhang     if (flg) {
1429b54502bSHong Zhang       ierr = PCLUSetReuseFill(pc,PETSC_TRUE);CHKERRQ(ierr);
1439b54502bSHong Zhang     }
1449b54502bSHong Zhang     ierr = PetscOptionsName("-pc_lu_reuse_ordering","Reuse ordering from previous factorization","PCLUSetReuseOrdering",&flg);CHKERRQ(ierr);
1459b54502bSHong Zhang     if (flg) {
1469b54502bSHong Zhang       ierr = PCLUSetReuseOrdering(pc,PETSC_TRUE);CHKERRQ(ierr);
1479b54502bSHong Zhang     }
1489b54502bSHong Zhang 
1499b54502bSHong Zhang     ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr);
1509b54502bSHong Zhang     ierr = PetscOptionsList("-pc_lu_mat_ordering_type","Reordering to reduce nonzeros in LU","PCLUSetMatOrdering",ordlist,lu->ordering,tname,256,&flg);CHKERRQ(ierr);
1519b54502bSHong Zhang     if (flg) {
1529b54502bSHong Zhang       ierr = PCLUSetMatOrdering(pc,tname);CHKERRQ(ierr);
1539b54502bSHong Zhang     }
1549b54502bSHong Zhang 
1559b54502bSHong Zhang     ierr = PetscOptionsName("-pc_lu_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCLUReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr);
1569b54502bSHong Zhang     if (flg) {
1579b54502bSHong Zhang       tol = PETSC_DECIDE;
1589b54502bSHong Zhang       ierr = PetscOptionsReal("-pc_lu_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCLUReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr);
1599b54502bSHong Zhang       ierr = PCLUReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr);
1609b54502bSHong Zhang     }
1619b54502bSHong Zhang 
1629b54502bSHong Zhang     ierr = PetscOptionsReal("-pc_lu_pivoting","Pivoting tolerance (used only for some factorization)","PCLUSetPivoting",lu->info.dtcol,&lu->info.dtcol,&flg);CHKERRQ(ierr);
1639b54502bSHong Zhang 
1649b54502bSHong Zhang     flg = lu->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE;
1654dc4c822SBarry Smith     ierr = PetscOptionsTruth("-pc_lu_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCLUSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr);
1669b54502bSHong Zhang     if (set) {
1679b54502bSHong Zhang       ierr = PCLUSetPivotInBlocks(pc,flg);CHKERRQ(ierr);
1689b54502bSHong Zhang     }
1699b54502bSHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
1709b54502bSHong Zhang   PetscFunctionReturn(0);
1719b54502bSHong Zhang }
1729b54502bSHong Zhang 
1739b54502bSHong Zhang #undef __FUNCT__
1749b54502bSHong Zhang #define __FUNCT__ "PCView_LU"
1759b54502bSHong Zhang static PetscErrorCode PCView_LU(PC pc,PetscViewer viewer)
1769b54502bSHong Zhang {
1779b54502bSHong Zhang   PC_LU          *lu = (PC_LU*)pc->data;
1789b54502bSHong Zhang   PetscErrorCode ierr;
1799b54502bSHong Zhang   PetscTruth     iascii,isstring;
1809b54502bSHong Zhang 
1819b54502bSHong Zhang   PetscFunctionBegin;
1829b54502bSHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
1839b54502bSHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
1849b54502bSHong Zhang   if (iascii) {
1859b54502bSHong Zhang     MatInfo info;
1869b54502bSHong Zhang 
1879b54502bSHong Zhang     if (lu->inplace) {ierr = PetscViewerASCIIPrintf(viewer,"  LU: in-place factorization\n");CHKERRQ(ierr);}
1889b54502bSHong Zhang     else             {ierr = PetscViewerASCIIPrintf(viewer,"  LU: out-of-place factorization\n");CHKERRQ(ierr);}
1899b54502bSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"    matrix ordering: %s\n",lu->ordering);CHKERRQ(ierr);
190a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  LU: tolerance for zero pivot %G\n",lu->info.zeropivot);CHKERRQ(ierr);
1910a29876aSHong Zhang     if (lu->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer,"  LU: using Manteuffel shift\n");CHKERRQ(ierr);}
1929b54502bSHong Zhang     if (lu->fact) {
1939b54502bSHong Zhang       ierr = MatGetInfo(lu->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
194a83599f4SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    LU nonzeros %G\n",info.nz_used);CHKERRQ(ierr);
1959b54502bSHong Zhang       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_FACTOR_INFO);CHKERRQ(ierr);
1969b54502bSHong Zhang       ierr = MatView(lu->fact,viewer);CHKERRQ(ierr);
1979b54502bSHong Zhang       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
1989b54502bSHong Zhang     }
1999b54502bSHong Zhang     if (lu->reusefill)    {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing fill from past factorization\n");CHKERRQ(ierr);}
2009b54502bSHong Zhang     if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing reordering from past factorization\n");CHKERRQ(ierr);}
2019b54502bSHong Zhang   } else if (isstring) {
2029b54502bSHong Zhang     ierr = PetscViewerStringSPrintf(viewer," order=%s",lu->ordering);CHKERRQ(ierr);CHKERRQ(ierr);
2039b54502bSHong Zhang   } else {
2049b54502bSHong Zhang     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCLU",((PetscObject)viewer)->type_name);
2059b54502bSHong Zhang   }
2069b54502bSHong Zhang   PetscFunctionReturn(0);
2079b54502bSHong Zhang }
2089b54502bSHong Zhang 
2099b54502bSHong Zhang #undef __FUNCT__
2109b54502bSHong Zhang #define __FUNCT__ "PCGetFactoredMatrix_LU"
2119b54502bSHong Zhang static PetscErrorCode PCGetFactoredMatrix_LU(PC pc,Mat *mat)
2129b54502bSHong Zhang {
2139b54502bSHong Zhang   PC_LU *dir = (PC_LU*)pc->data;
2149b54502bSHong Zhang 
2159b54502bSHong Zhang   PetscFunctionBegin;
2169b54502bSHong Zhang   if (!dir->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()");
2179b54502bSHong Zhang   *mat = dir->fact;
2189b54502bSHong Zhang   PetscFunctionReturn(0);
2199b54502bSHong Zhang }
2209b54502bSHong Zhang 
2219b54502bSHong Zhang #undef __FUNCT__
2229b54502bSHong Zhang #define __FUNCT__ "PCSetUp_LU"
2239b54502bSHong Zhang static PetscErrorCode PCSetUp_LU(PC pc)
2249b54502bSHong Zhang {
2259b54502bSHong Zhang   PetscErrorCode ierr;
2269b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
2279b54502bSHong Zhang 
2289b54502bSHong Zhang   PetscFunctionBegin;
2299b54502bSHong Zhang   if (dir->reusefill && pc->setupcalled) dir->info.fill = dir->actualfill;
2309b54502bSHong Zhang 
2319b54502bSHong Zhang   if (dir->inplace) {
2329b54502bSHong Zhang     if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);}
2339b54502bSHong Zhang     if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);}
2349b54502bSHong Zhang     ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
23503c60df9SBarry Smith     if (dir->row) {
23603c60df9SBarry Smith       ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);
23703c60df9SBarry Smith       ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr);
23803c60df9SBarry Smith     }
2399b54502bSHong Zhang     ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&dir->info);CHKERRQ(ierr);
2409b54502bSHong Zhang     dir->fact = pc->pmat;
2419b54502bSHong Zhang   } else {
2429b54502bSHong Zhang     MatInfo info;
2439b54502bSHong Zhang     if (!pc->setupcalled) {
2449b54502bSHong Zhang       ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
2459b54502bSHong Zhang       if (dir->nonzerosalongdiagonal) {
2469b54502bSHong Zhang         ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr);
2479b54502bSHong Zhang       }
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 = MatLUFactorSymbolic(pc->pmat,dir->row,dir->col,&dir->info,&dir->fact);CHKERRQ(ierr);
2539b54502bSHong Zhang       ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
2549b54502bSHong Zhang       dir->actualfill = info.fill_ratio_needed;
25552e6d16bSBarry Smith       ierr = PetscLogObjectParent(pc,dir->fact);CHKERRQ(ierr);
2569b54502bSHong Zhang     } else if (pc->flag != SAME_NONZERO_PATTERN) {
2579b54502bSHong Zhang       if (!dir->reuseordering) {
2589b54502bSHong Zhang         if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);}
2599b54502bSHong Zhang         if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);}
2609b54502bSHong Zhang         ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
2619b54502bSHong Zhang         if (dir->nonzerosalongdiagonal) {
2629b54502bSHong Zhang           ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr);
2639b54502bSHong Zhang         }
26403c60df9SBarry Smith         if (dir->row) {
26503c60df9SBarry Smith           ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);
26603c60df9SBarry Smith           ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr);
26703c60df9SBarry Smith         }
2689b54502bSHong Zhang       }
2699b54502bSHong Zhang       ierr = MatDestroy(dir->fact);CHKERRQ(ierr);
2709b54502bSHong Zhang       ierr = MatLUFactorSymbolic(pc->pmat,dir->row,dir->col,&dir->info,&dir->fact);CHKERRQ(ierr);
2719b54502bSHong Zhang       ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
2729b54502bSHong Zhang       dir->actualfill = info.fill_ratio_needed;
27352e6d16bSBarry Smith       ierr = PetscLogObjectParent(pc,dir->fact);CHKERRQ(ierr);
2749b54502bSHong Zhang     }
2759b54502bSHong Zhang     ierr = MatLUFactorNumeric(pc->pmat,&dir->info,&dir->fact);CHKERRQ(ierr);
2769b54502bSHong Zhang   }
2779b54502bSHong Zhang   PetscFunctionReturn(0);
2789b54502bSHong Zhang }
2799b54502bSHong Zhang 
2809b54502bSHong Zhang #undef __FUNCT__
2819b54502bSHong Zhang #define __FUNCT__ "PCDestroy_LU"
2829b54502bSHong Zhang static PetscErrorCode PCDestroy_LU(PC pc)
2839b54502bSHong Zhang {
2849b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
2859b54502bSHong Zhang   PetscErrorCode ierr;
2869b54502bSHong Zhang 
2879b54502bSHong Zhang   PetscFunctionBegin;
2889b54502bSHong Zhang   if (!dir->inplace && dir->fact) {ierr = MatDestroy(dir->fact);CHKERRQ(ierr);}
2899b54502bSHong Zhang   if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);}
2909b54502bSHong Zhang   if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);}
2919b54502bSHong Zhang   ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr);
2929b54502bSHong Zhang   ierr = PetscFree(dir);CHKERRQ(ierr);
2939b54502bSHong Zhang   PetscFunctionReturn(0);
2949b54502bSHong Zhang }
2959b54502bSHong Zhang 
2969b54502bSHong Zhang #undef __FUNCT__
2979b54502bSHong Zhang #define __FUNCT__ "PCApply_LU"
2989b54502bSHong Zhang static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y)
2999b54502bSHong Zhang {
3009b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
3019b54502bSHong Zhang   PetscErrorCode ierr;
3029b54502bSHong Zhang 
3039b54502bSHong Zhang   PetscFunctionBegin;
3049b54502bSHong Zhang   if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);}
3059b54502bSHong Zhang   else              {ierr = MatSolve(dir->fact,x,y);CHKERRQ(ierr);}
3069b54502bSHong Zhang   PetscFunctionReturn(0);
3079b54502bSHong Zhang }
3089b54502bSHong Zhang 
3099b54502bSHong Zhang #undef __FUNCT__
3109b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_LU"
3119b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_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 = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);}
3189b54502bSHong Zhang   else              {ierr = MatSolveTranspose(dir->fact,x,y);CHKERRQ(ierr);}
3199b54502bSHong Zhang   PetscFunctionReturn(0);
3209b54502bSHong Zhang }
3219b54502bSHong Zhang 
3229b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/
3239b54502bSHong Zhang 
3249b54502bSHong Zhang EXTERN_C_BEGIN
3259b54502bSHong Zhang #undef __FUNCT__
326*55ba2a51SBarry Smith #define __FUNCT__ "PCFactorSetFill_LU"
327*55ba2a51SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill_LU(PC pc,PetscReal fill)
3289b54502bSHong Zhang {
3299b54502bSHong Zhang   PC_LU *dir;
3309b54502bSHong Zhang 
3319b54502bSHong Zhang   PetscFunctionBegin;
3329b54502bSHong Zhang   dir = (PC_LU*)pc->data;
3339b54502bSHong Zhang   dir->info.fill = fill;
3349b54502bSHong Zhang   PetscFunctionReturn(0);
3359b54502bSHong Zhang }
3369b54502bSHong Zhang EXTERN_C_END
3379b54502bSHong Zhang 
3389b54502bSHong Zhang EXTERN_C_BEGIN
3399b54502bSHong Zhang #undef __FUNCT__
3409b54502bSHong Zhang #define __FUNCT__ "PCLUSetUseInPlace_LU"
341dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCLUSetUseInPlace_LU(PC pc)
3429b54502bSHong Zhang {
3439b54502bSHong Zhang   PC_LU *dir;
3449b54502bSHong Zhang 
3459b54502bSHong Zhang   PetscFunctionBegin;
3469b54502bSHong Zhang   dir = (PC_LU*)pc->data;
3479b54502bSHong Zhang   dir->inplace = PETSC_TRUE;
3489b54502bSHong Zhang   PetscFunctionReturn(0);
3499b54502bSHong Zhang }
3509b54502bSHong Zhang EXTERN_C_END
3519b54502bSHong Zhang 
3529b54502bSHong Zhang EXTERN_C_BEGIN
3539b54502bSHong Zhang #undef __FUNCT__
3549b54502bSHong Zhang #define __FUNCT__ "PCLUSetMatOrdering_LU"
355dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCLUSetMatOrdering_LU(PC pc,MatOrderingType ordering)
3569b54502bSHong Zhang {
3579b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
3589b54502bSHong Zhang   PetscErrorCode ierr;
3599b54502bSHong Zhang 
3609b54502bSHong Zhang   PetscFunctionBegin;
3619b54502bSHong Zhang   ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr);
3629b54502bSHong Zhang   ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr);
3639b54502bSHong Zhang   PetscFunctionReturn(0);
3649b54502bSHong Zhang }
3659b54502bSHong Zhang EXTERN_C_END
3669b54502bSHong Zhang 
3679b54502bSHong Zhang EXTERN_C_BEGIN
3689b54502bSHong Zhang #undef __FUNCT__
3699b54502bSHong Zhang #define __FUNCT__ "PCLUSetPivoting_LU"
370dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCLUSetPivoting_LU(PC pc,PetscReal dtcol)
3719b54502bSHong Zhang {
3729b54502bSHong Zhang   PC_LU *dir = (PC_LU*)pc->data;
3739b54502bSHong Zhang 
3749b54502bSHong Zhang   PetscFunctionBegin;
375a83599f4SBarry 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);
3769b54502bSHong Zhang   dir->info.dtcol = dtcol;
3779b54502bSHong Zhang   PetscFunctionReturn(0);
3789b54502bSHong Zhang }
3799b54502bSHong Zhang EXTERN_C_END
3809b54502bSHong Zhang 
3819b54502bSHong Zhang EXTERN_C_BEGIN
3829b54502bSHong Zhang #undef __FUNCT__
3839b54502bSHong Zhang #define __FUNCT__ "PCLUSetPivotInBlocks_LU"
384dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCLUSetPivotInBlocks_LU(PC pc,PetscTruth pivot)
3859b54502bSHong Zhang {
3869b54502bSHong Zhang   PC_LU *dir = (PC_LU*)pc->data;
3879b54502bSHong Zhang 
3889b54502bSHong Zhang   PetscFunctionBegin;
3899b54502bSHong Zhang   dir->info.pivotinblocks = pivot ? 1.0 : 0.0;
3909b54502bSHong Zhang   PetscFunctionReturn(0);
3919b54502bSHong Zhang }
3929b54502bSHong Zhang EXTERN_C_END
3939b54502bSHong Zhang 
3949b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/
3959b54502bSHong Zhang 
3969b54502bSHong Zhang #undef __FUNCT__
3979b54502bSHong Zhang #define __FUNCT__ "PCLUReorderForNonzeroDiagonal"
3989b54502bSHong Zhang /*@
3999b54502bSHong Zhang    PCLUReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal
4009b54502bSHong Zhang 
4019b54502bSHong Zhang    Collective on PC
4029b54502bSHong Zhang 
4039b54502bSHong Zhang    Input Parameters:
4049b54502bSHong Zhang +  pc - the preconditioner context
4059b54502bSHong Zhang -  tol - diagonal entries smaller than this in absolute value are considered zero
4069b54502bSHong Zhang 
4079b54502bSHong Zhang    Options Database Key:
4089b54502bSHong Zhang .  -pc_lu_nonzeros_along_diagonal
4099b54502bSHong Zhang 
4109b54502bSHong Zhang    Level: intermediate
4119b54502bSHong Zhang 
4129b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill
4139b54502bSHong Zhang 
414*55ba2a51SBarry Smith .seealso: PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetZeroPivot(), MatReorderForNonzeroDiagonal()
4159b54502bSHong Zhang @*/
416dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCLUReorderForNonzeroDiagonal(PC pc,PetscReal rtol)
4179b54502bSHong Zhang {
4189b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscReal);
4199b54502bSHong Zhang 
4209b54502bSHong Zhang   PetscFunctionBegin;
4219b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
4229b54502bSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUReorderForNonzeroDiagonal_C",(void (**)(void))&f);CHKERRQ(ierr);
4239b54502bSHong Zhang   if (f) {
4249b54502bSHong Zhang     ierr = (*f)(pc,rtol);CHKERRQ(ierr);
4259b54502bSHong Zhang   }
4269b54502bSHong Zhang   PetscFunctionReturn(0);
4279b54502bSHong Zhang }
4289b54502bSHong Zhang 
4299b54502bSHong Zhang #undef __FUNCT__
4309b54502bSHong Zhang #define __FUNCT__ "PCLUSetReuseOrdering"
4319b54502bSHong Zhang /*@
4329b54502bSHong Zhang    PCLUSetReuseOrdering - When similar matrices are factored, this
4339b54502bSHong Zhang    causes the ordering computed in the first factor to be used for all
4349b54502bSHong Zhang    following factors; applies to both fill and drop tolerance LUs.
4359b54502bSHong Zhang 
4369b54502bSHong Zhang    Collective on PC
4379b54502bSHong Zhang 
4389b54502bSHong Zhang    Input Parameters:
4399b54502bSHong Zhang +  pc - the preconditioner context
4409b54502bSHong Zhang -  flag - PETSC_TRUE to reuse else PETSC_FALSE
4419b54502bSHong Zhang 
4429b54502bSHong Zhang    Options Database Key:
4439b54502bSHong Zhang .  -pc_lu_reuse_ordering - Activate PCLUSetReuseOrdering()
4449b54502bSHong Zhang 
4459b54502bSHong Zhang    Level: intermediate
4469b54502bSHong Zhang 
4479b54502bSHong Zhang .keywords: PC, levels, reordering, factorization, incomplete, LU
4489b54502bSHong Zhang 
4499b54502bSHong Zhang .seealso: PCLUSetReuseFill(), PCILUSetReuseOrdering(), PCILUDTSetReuseFill()
4509b54502bSHong Zhang @*/
451dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCLUSetReuseOrdering(PC pc,PetscTruth flag)
4529b54502bSHong Zhang {
4539b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscTruth);
4549b54502bSHong Zhang 
4559b54502bSHong Zhang   PetscFunctionBegin;
4569b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
4579b54502bSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetReuseOrdering_C",(void (**)(void))&f);CHKERRQ(ierr);
4589b54502bSHong Zhang   if (f) {
4599b54502bSHong Zhang     ierr = (*f)(pc,flag);CHKERRQ(ierr);
4609b54502bSHong Zhang   }
4619b54502bSHong Zhang   PetscFunctionReturn(0);
4629b54502bSHong Zhang }
4639b54502bSHong Zhang 
4649b54502bSHong Zhang #undef __FUNCT__
4659b54502bSHong Zhang #define __FUNCT__ "PCLUSetReuseFill"
4669b54502bSHong Zhang /*@
4679b54502bSHong Zhang    PCLUSetReuseFill - When matrices with same nonzero structure are LU factored,
4689b54502bSHong Zhang    this causes later ones to use the fill computed in the initial factorization.
4699b54502bSHong Zhang 
4709b54502bSHong Zhang    Collective on PC
4719b54502bSHong Zhang 
4729b54502bSHong Zhang    Input Parameters:
4739b54502bSHong Zhang +  pc - the preconditioner context
4749b54502bSHong Zhang -  flag - PETSC_TRUE to reuse else PETSC_FALSE
4759b54502bSHong Zhang 
4769b54502bSHong Zhang    Options Database Key:
4779b54502bSHong Zhang .  -pc_lu_reuse_fill - Activates PCLUSetReuseFill()
4789b54502bSHong Zhang 
4799b54502bSHong Zhang    Level: intermediate
4809b54502bSHong Zhang 
4819b54502bSHong Zhang .keywords: PC, levels, reordering, factorization, incomplete, LU
4829b54502bSHong Zhang 
4839b54502bSHong Zhang .seealso: PCILUSetReuseOrdering(), PCLUSetReuseOrdering(), PCILUDTSetReuseFill()
4849b54502bSHong Zhang @*/
485dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCLUSetReuseFill(PC pc,PetscTruth flag)
4869b54502bSHong Zhang {
4879b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscTruth);
4889b54502bSHong Zhang 
4899b54502bSHong Zhang   PetscFunctionBegin;
4909b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
4919b54502bSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetReuseFill_C",(void (**)(void))&f);CHKERRQ(ierr);
4929b54502bSHong Zhang   if (f) {
4939b54502bSHong Zhang     ierr = (*f)(pc,flag);CHKERRQ(ierr);
4949b54502bSHong Zhang   }
4959b54502bSHong Zhang   PetscFunctionReturn(0);
4969b54502bSHong Zhang }
4979b54502bSHong Zhang 
4989b54502bSHong Zhang #undef __FUNCT__
499*55ba2a51SBarry Smith #define __FUNCT__ "PCFactorSetFill"
5009b54502bSHong Zhang /*@
501*55ba2a51SBarry Smith    PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix,
5029b54502bSHong Zhang    fill = number nonzeros in factor/number nonzeros in original matrix.
5039b54502bSHong Zhang 
5049b54502bSHong Zhang    Collective on PC
5059b54502bSHong Zhang 
5069b54502bSHong Zhang    Input Parameters:
5079b54502bSHong Zhang +  pc - the preconditioner context
5089b54502bSHong Zhang -  fill - amount of expected fill
5099b54502bSHong Zhang 
5109b54502bSHong Zhang    Options Database Key:
511*55ba2a51SBarry Smith .  -pc_factor_fill <fill> - Sets fill amount
5129b54502bSHong Zhang 
5139b54502bSHong Zhang    Level: intermediate
5149b54502bSHong Zhang 
5159b54502bSHong Zhang    Note:
5169b54502bSHong Zhang    For sparse matrix factorizations it is difficult to predict how much
5179be77112SSatish Balay    fill to expect. By running with the option -verbose_info PETSc will print the
5189b54502bSHong Zhang    actual amount of fill used; allowing you to set the value accurately for
5199b54502bSHong Zhang    future runs. Default PETSc uses a value of 5.0
5209b54502bSHong Zhang 
5219b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill
5229b54502bSHong Zhang 
5239b54502bSHong Zhang @*/
524*55ba2a51SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill(PC pc,PetscReal fill)
5259b54502bSHong Zhang {
5269b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscReal);
5279b54502bSHong Zhang 
5289b54502bSHong Zhang   PetscFunctionBegin;
5299b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
5309b54502bSHong Zhang   if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0");
531*55ba2a51SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetFill_C",(void (**)(void))&f);CHKERRQ(ierr);
5329b54502bSHong Zhang   if (f) {
5339b54502bSHong Zhang     ierr = (*f)(pc,fill);CHKERRQ(ierr);
5349b54502bSHong Zhang   }
5359b54502bSHong Zhang   PetscFunctionReturn(0);
5369b54502bSHong Zhang }
5379b54502bSHong Zhang 
5389b54502bSHong Zhang #undef __FUNCT__
5399b54502bSHong Zhang #define __FUNCT__ "PCLUSetUseInPlace"
5409b54502bSHong Zhang /*@
5419b54502bSHong Zhang    PCLUSetUseInPlace - Tells the system to do an in-place factorization.
5429b54502bSHong Zhang    For dense matrices, this enables the solution of much larger problems.
5439b54502bSHong Zhang    For sparse matrices the factorization cannot be done truly in-place
5449b54502bSHong Zhang    so this does not save memory during the factorization, but after the matrix
5459b54502bSHong Zhang    is factored, the original unfactored matrix is freed, thus recovering that
5469b54502bSHong Zhang    space.
5479b54502bSHong Zhang 
5489b54502bSHong Zhang    Collective on PC
5499b54502bSHong Zhang 
5509b54502bSHong Zhang    Input Parameters:
5519b54502bSHong Zhang .  pc - the preconditioner context
5529b54502bSHong Zhang 
5539b54502bSHong Zhang    Options Database Key:
5549b54502bSHong Zhang .  -pc_lu_in_place - Activates in-place factorization
5559b54502bSHong Zhang 
5569b54502bSHong Zhang    Notes:
5579b54502bSHong Zhang    PCLUSetUseInplace() can only be used with the KSP method KSPPREONLY or when
5589b54502bSHong Zhang    a different matrix is provided for the multiply and the preconditioner in
5599b54502bSHong Zhang    a call to KSPSetOperators().
5609b54502bSHong Zhang    This is because the Krylov space methods require an application of the
5619b54502bSHong Zhang    matrix multiplication, which is not possible here because the matrix has
5629b54502bSHong Zhang    been factored in-place, replacing the original matrix.
5639b54502bSHong Zhang 
5649b54502bSHong Zhang    Level: intermediate
5659b54502bSHong Zhang 
5669b54502bSHong Zhang .keywords: PC, set, factorization, direct, inplace, in-place, LU
5679b54502bSHong Zhang 
5689b54502bSHong Zhang .seealso: PCILUSetUseInPlace()
5699b54502bSHong Zhang @*/
570dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCLUSetUseInPlace(PC pc)
5719b54502bSHong Zhang {
5729b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC);
5739b54502bSHong Zhang 
5749b54502bSHong Zhang   PetscFunctionBegin;
5759b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
5769b54502bSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetUseInPlace_C",(void (**)(void))&f);CHKERRQ(ierr);
5779b54502bSHong Zhang   if (f) {
5789b54502bSHong Zhang     ierr = (*f)(pc);CHKERRQ(ierr);
5799b54502bSHong Zhang   }
5809b54502bSHong Zhang   PetscFunctionReturn(0);
5819b54502bSHong Zhang }
5829b54502bSHong Zhang 
5839b54502bSHong Zhang #undef __FUNCT__
5849b54502bSHong Zhang #define __FUNCT__ "PCLUSetMatOrdering"
5859b54502bSHong Zhang /*@C
5869b54502bSHong Zhang     PCLUSetMatOrdering - Sets the ordering routine (to reduce fill) to
5879b54502bSHong Zhang     be used in the LU factorization.
5889b54502bSHong Zhang 
5899b54502bSHong Zhang     Collective on PC
5909b54502bSHong Zhang 
5919b54502bSHong Zhang     Input Parameters:
5929b54502bSHong Zhang +   pc - the preconditioner context
5939b54502bSHong Zhang -   ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM
5949b54502bSHong Zhang 
5959b54502bSHong Zhang     Options Database Key:
5969b54502bSHong Zhang .   -pc_lu_mat_ordering_type <nd,rcm,...> - Sets ordering routine
5979b54502bSHong Zhang 
5989b54502bSHong Zhang     Level: intermediate
5999b54502bSHong Zhang 
6009b54502bSHong Zhang     Notes: nested dissection is used by default
6019b54502bSHong Zhang 
6029b54502bSHong Zhang .seealso: PCILUSetMatOrdering()
6039b54502bSHong Zhang @*/
604dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCLUSetMatOrdering(PC pc,MatOrderingType ordering)
6059b54502bSHong Zhang {
6069b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,MatOrderingType);
6079b54502bSHong Zhang 
6089b54502bSHong Zhang   PetscFunctionBegin;
6099b54502bSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetMatOrdering_C",(void (**)(void))&f);CHKERRQ(ierr);
6109b54502bSHong Zhang   if (f) {
6119b54502bSHong Zhang     ierr = (*f)(pc,ordering);CHKERRQ(ierr);
6129b54502bSHong Zhang   }
6139b54502bSHong Zhang   PetscFunctionReturn(0);
6149b54502bSHong Zhang }
6159b54502bSHong Zhang 
6169b54502bSHong Zhang #undef __FUNCT__
6179b54502bSHong Zhang #define __FUNCT__ "PCLUSetPivoting"
6189b54502bSHong Zhang /*@
6199b54502bSHong Zhang     PCLUSetPivoting - Determines when pivoting is done during LU.
6209b54502bSHong Zhang       For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices
6219b54502bSHong Zhang       it is never done. For the Matlab and SuperLU factorization this is used.
6229b54502bSHong Zhang 
6239b54502bSHong Zhang     Collective on PC
6249b54502bSHong Zhang 
6259b54502bSHong Zhang     Input Parameters:
6269b54502bSHong Zhang +   pc - the preconditioner context
6279b54502bSHong Zhang -   dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable)
6289b54502bSHong Zhang 
6299b54502bSHong Zhang     Options Database Key:
6309b54502bSHong Zhang .   -pc_lu_pivoting <dtcol>
6319b54502bSHong Zhang 
6329b54502bSHong Zhang     Level: intermediate
6339b54502bSHong Zhang 
6349b54502bSHong Zhang .seealso: PCILUSetMatOrdering(), PCLUSetPivotInBlocks()
6359b54502bSHong Zhang @*/
636dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCLUSetPivoting(PC pc,PetscReal dtcol)
6379b54502bSHong Zhang {
6389b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscReal);
6399b54502bSHong Zhang 
6409b54502bSHong Zhang   PetscFunctionBegin;
6419b54502bSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetPivoting_C",(void (**)(void))&f);CHKERRQ(ierr);
6429b54502bSHong Zhang   if (f) {
6439b54502bSHong Zhang     ierr = (*f)(pc,dtcol);CHKERRQ(ierr);
6449b54502bSHong Zhang   }
6459b54502bSHong Zhang   PetscFunctionReturn(0);
6469b54502bSHong Zhang }
6479b54502bSHong Zhang 
6489b54502bSHong Zhang #undef __FUNCT__
6499b54502bSHong Zhang #define __FUNCT__ "PCLUSetPivotInBlocks"
6509b54502bSHong Zhang /*@
6519b54502bSHong Zhang     PCLUSetPivotInBlocks - Determines if pivoting is done while factoring each block
6529b54502bSHong Zhang       with BAIJ or SBAIJ matrices
6539b54502bSHong Zhang 
6549b54502bSHong Zhang     Collective on PC
6559b54502bSHong Zhang 
6569b54502bSHong Zhang     Input Parameters:
6579b54502bSHong Zhang +   pc - the preconditioner context
6589b54502bSHong Zhang -   pivot - PETSC_TRUE or PETSC_FALSE
6599b54502bSHong Zhang 
6609b54502bSHong Zhang     Options Database Key:
6619b54502bSHong Zhang .   -pc_lu_pivot_in_blocks <true,false>
6629b54502bSHong Zhang 
6639b54502bSHong Zhang     Level: intermediate
6649b54502bSHong Zhang 
6659b54502bSHong Zhang .seealso: PCILUSetMatOrdering(), PCLUSetPivoting()
6669b54502bSHong Zhang @*/
667dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCLUSetPivotInBlocks(PC pc,PetscTruth pivot)
6689b54502bSHong Zhang {
6699b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscTruth);
6709b54502bSHong Zhang 
6719b54502bSHong Zhang   PetscFunctionBegin;
6729b54502bSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetPivotInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr);
6739b54502bSHong Zhang   if (f) {
6749b54502bSHong Zhang     ierr = (*f)(pc,pivot);CHKERRQ(ierr);
6759b54502bSHong Zhang   }
6769b54502bSHong Zhang   PetscFunctionReturn(0);
6779b54502bSHong Zhang }
6789b54502bSHong Zhang 
6799b54502bSHong Zhang /* ------------------------------------------------------------------------ */
6809b54502bSHong Zhang 
6819b54502bSHong Zhang /*MC
6829b54502bSHong Zhang    PCLU - Uses a direct solver, based on LU factorization, as a preconditioner
6839b54502bSHong Zhang 
6849b54502bSHong Zhang    Options Database Keys:
6859b54502bSHong Zhang +  -pc_lu_reuse_ordering - Activate PCLUSetReuseOrdering()
6869b54502bSHong Zhang .  -pc_lu_reuse_fill - Activates PCLUSetReuseFill()
687*55ba2a51SBarry Smith .  -pc_factor_fill <fill> - Sets fill amount
6889b54502bSHong Zhang .  -pc_lu_in_place - Activates in-place factorization
6899b54502bSHong Zhang .  -pc_lu_mat_ordering_type <nd,rcm,...> - Sets ordering routine
690f251bdbdSHong Zhang .  -pc_lu_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase
6919b54502bSHong Zhang                                          stability of factorization.
692f251bdbdSHong Zhang .  -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default
693f251bdbdSHong Zhang -  -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value
694f251bdbdSHong Zhang    is optional with PETSC_TRUE being the default
6959b54502bSHong Zhang 
6969b54502bSHong Zhang    Notes: Not all options work for all matrix formats
6979b54502bSHong Zhang           Run with -help to see additional options for particular matrix formats or factorization
6989b54502bSHong Zhang           algorithms
6999b54502bSHong Zhang 
7009b54502bSHong Zhang    Level: beginner
7019b54502bSHong Zhang 
7029b54502bSHong Zhang    Concepts: LU factorization, direct solver
7039b54502bSHong Zhang 
7049b54502bSHong Zhang    Notes: Usually this will compute an "exact" solution in one iteration and does
7059b54502bSHong Zhang           not need a Krylov method (i.e. you can use -ksp_type preonly, or
7069b54502bSHong Zhang           KSPSetType(ksp,KSPPREONLY) for the Krylov method
7079b54502bSHong Zhang 
7089b54502bSHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
7099b54502bSHong Zhang            PCILU, PCCHOLESKY, PCICC, PCLUSetReuseOrdering(), PCLUSetReuseFill(), PCGetFactoredMatrix(),
710*55ba2a51SBarry Smith            PCFactorSetFill(), PCLUSetUseInPlace(), PCLUSetMatOrdering(), PCFactorSetPivoting(),
711f251bdbdSHong Zhang            PCLUSetPivotingInBlocks(),PCFactorSetShiftNonzero(),PCFactorSetShiftPd()
7129b54502bSHong Zhang M*/
7139b54502bSHong Zhang 
7149b54502bSHong Zhang EXTERN_C_BEGIN
7159b54502bSHong Zhang #undef __FUNCT__
7169b54502bSHong Zhang #define __FUNCT__ "PCCreate_LU"
717dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_LU(PC pc)
7189b54502bSHong Zhang {
7199b54502bSHong Zhang   PetscErrorCode ierr;
7209b54502bSHong Zhang   PetscMPIInt    size;
7219b54502bSHong Zhang   PC_LU          *dir;
7229b54502bSHong Zhang 
7239b54502bSHong Zhang   PetscFunctionBegin;
7249b54502bSHong Zhang   ierr = PetscNew(PC_LU,&dir);CHKERRQ(ierr);
72552e6d16bSBarry Smith   ierr = PetscLogObjectMemory(pc,sizeof(PC_LU));CHKERRQ(ierr);
7269b54502bSHong Zhang 
7279b54502bSHong Zhang   ierr = MatFactorInfoInitialize(&dir->info);CHKERRQ(ierr);
7289b54502bSHong Zhang   dir->fact                  = 0;
7299b54502bSHong Zhang   dir->inplace               = PETSC_FALSE;
7309b54502bSHong Zhang   dir->nonzerosalongdiagonal = PETSC_FALSE;
7319b54502bSHong Zhang 
7329b54502bSHong Zhang   dir->info.fill           = 5.0;
7339b54502bSHong Zhang   dir->info.dtcol          = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */
7340a29876aSHong Zhang   dir->info.shiftnz        = 0.0;
7359b54502bSHong Zhang   dir->info.zeropivot      = 1.e-12;
7369b54502bSHong Zhang   dir->info.pivotinblocks  = 1.0;
737fbf22428SSatish Balay   dir->info.shiftpd        = 0.0; /* false */
7389b54502bSHong Zhang   dir->info.shift_fraction = 0.0;
7399b54502bSHong Zhang   dir->col                 = 0;
7409b54502bSHong Zhang   dir->row                 = 0;
7419b54502bSHong Zhang   ierr = MPI_Comm_size(pc->comm,&size);CHKERRQ(ierr);
7429b54502bSHong Zhang   if (size == 1) {
7439b54502bSHong Zhang     ierr = PetscStrallocpy(MATORDERING_ND,&dir->ordering);CHKERRQ(ierr);
7449b54502bSHong Zhang   } else {
7459b54502bSHong Zhang     ierr = PetscStrallocpy(MATORDERING_NATURAL,&dir->ordering);CHKERRQ(ierr);
7469b54502bSHong Zhang   }
7479b54502bSHong Zhang   dir->reusefill        = PETSC_FALSE;
7489b54502bSHong Zhang   dir->reuseordering    = PETSC_FALSE;
7499b54502bSHong Zhang   pc->data              = (void*)dir;
7509b54502bSHong Zhang 
7519b54502bSHong Zhang   pc->ops->destroy           = PCDestroy_LU;
7529b54502bSHong Zhang   pc->ops->apply             = PCApply_LU;
7539b54502bSHong Zhang   pc->ops->applytranspose    = PCApplyTranspose_LU;
7549b54502bSHong Zhang   pc->ops->setup             = PCSetUp_LU;
7559b54502bSHong Zhang   pc->ops->setfromoptions    = PCSetFromOptions_LU;
7569b54502bSHong Zhang   pc->ops->view              = PCView_LU;
7579b54502bSHong Zhang   pc->ops->applyrichardson   = 0;
7589b54502bSHong Zhang   pc->ops->getfactoredmatrix = PCGetFactoredMatrix_LU;
7599b54502bSHong Zhang 
760afaefe49SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_LU",
761afaefe49SHong Zhang                     PCFactorSetZeroPivot_LU);CHKERRQ(ierr);
762afaefe49SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_LU",
763afaefe49SHong Zhang                     PCFactorSetShiftNonzero_LU);CHKERRQ(ierr);
764afaefe49SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_LU",
765afaefe49SHong Zhang                     PCFactorSetShiftPd_LU);CHKERRQ(ierr);
766afaefe49SHong Zhang 
767*55ba2a51SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_LU",
768*55ba2a51SBarry Smith                     PCFactorSetFill_LU);CHKERRQ(ierr);
7699b54502bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetUseInPlace_C","PCLUSetUseInPlace_LU",
7709b54502bSHong Zhang                     PCLUSetUseInPlace_LU);CHKERRQ(ierr);
7719b54502bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetMatOrdering_C","PCLUSetMatOrdering_LU",
7729b54502bSHong Zhang                     PCLUSetMatOrdering_LU);CHKERRQ(ierr);
7739b54502bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetReuseOrdering_C","PCLUSetReuseOrdering_LU",
7749b54502bSHong Zhang                     PCLUSetReuseOrdering_LU);CHKERRQ(ierr);
7759b54502bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetReuseFill_C","PCLUSetReuseFill_LU",
7769b54502bSHong Zhang                     PCLUSetReuseFill_LU);CHKERRQ(ierr);
7779b54502bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetPivoting_C","PCLUSetPivoting_LU",
7789b54502bSHong Zhang                     PCLUSetPivoting_LU);CHKERRQ(ierr);
7799b54502bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetPivotInBlocks_C","PCLUSetPivotInBlocks_LU",
7809b54502bSHong Zhang                     PCLUSetPivotInBlocks_LU);CHKERRQ(ierr);
7819b54502bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUReorderForNonzeroDiagonal_C","PCLUReorderForNonzeroDiagonal_LU",
7829b54502bSHong Zhang                     PCLUReorderForNonzeroDiagonal_LU);CHKERRQ(ierr);
7839b54502bSHong Zhang   PetscFunctionReturn(0);
7849b54502bSHong Zhang }
7859b54502bSHong Zhang EXTERN_C_END
786