xref: /petsc/src/ksp/pc/impls/factor/lu/lu.c (revision 6cf91177622bd3d3e4d9f7c593ec388633359851)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2dba47a55SKris Buschelman 
39b54502bSHong Zhang /*
49b54502bSHong Zhang    Defines a direct factorization preconditioner for any Mat implementation
59b54502bSHong Zhang    Note: this need not be consided a preconditioner since it supplies
69b54502bSHong Zhang          a direct solver.
79b54502bSHong Zhang */
8ee45ca4aSHong Zhang 
96356e834SBarry Smith #include "private/pcimpl.h"                /*I "petscpc.h" I*/
10ee45ca4aSHong Zhang #include "src/ksp/pc/impls/factor/lu/lu.h"
119b54502bSHong Zhang 
129b54502bSHong Zhang EXTERN_C_BEGIN
139b54502bSHong Zhang #undef __FUNCT__
14afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetZeroPivot_LU"
15dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot_LU(PC pc,PetscReal z)
16afaefe49SHong Zhang {
17afaefe49SHong Zhang   PC_LU *lu;
18afaefe49SHong Zhang 
19afaefe49SHong Zhang   PetscFunctionBegin;
20afaefe49SHong Zhang   lu                 = (PC_LU*)pc->data;
21afaefe49SHong Zhang   lu->info.zeropivot = z;
22afaefe49SHong Zhang   PetscFunctionReturn(0);
23afaefe49SHong Zhang }
24afaefe49SHong Zhang EXTERN_C_END
25afaefe49SHong Zhang 
26afaefe49SHong Zhang EXTERN_C_BEGIN
27afaefe49SHong Zhang #undef __FUNCT__
28afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetShiftNonzero_LU"
29dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero_LU(PC pc,PetscReal shift)
30afaefe49SHong Zhang {
31afaefe49SHong Zhang   PC_LU *dir;
32afaefe49SHong Zhang 
33afaefe49SHong Zhang   PetscFunctionBegin;
34afaefe49SHong Zhang   dir = (PC_LU*)pc->data;
35afaefe49SHong Zhang   if (shift == (PetscReal) PETSC_DECIDE) {
36afaefe49SHong Zhang     dir->info.shiftnz = 1.e-12;
37afaefe49SHong Zhang   } else {
38afaefe49SHong Zhang     dir->info.shiftnz = shift;
39afaefe49SHong Zhang   }
40afaefe49SHong Zhang   PetscFunctionReturn(0);
41afaefe49SHong Zhang }
42afaefe49SHong Zhang EXTERN_C_END
43afaefe49SHong Zhang 
44afaefe49SHong Zhang EXTERN_C_BEGIN
45afaefe49SHong Zhang #undef __FUNCT__
46afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetShiftPd_LU"
47dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd_LU(PC pc,PetscTruth shift)
48afaefe49SHong Zhang {
49afaefe49SHong Zhang   PC_LU *dir;
50afaefe49SHong Zhang 
51afaefe49SHong Zhang   PetscFunctionBegin;
52afaefe49SHong Zhang   dir = (PC_LU*)pc->data;
53fbf22428SSatish Balay   if (shift) {
54fbf22428SSatish Balay     dir->info.shift_fraction = 0.0;
55fbf22428SSatish Balay     dir->info.shiftpd = 1.0;
56fbf22428SSatish Balay   } else {
57fbf22428SSatish Balay     dir->info.shiftpd = 0.0;
58fbf22428SSatish Balay   }
59afaefe49SHong Zhang   PetscFunctionReturn(0);
60afaefe49SHong Zhang }
61afaefe49SHong Zhang EXTERN_C_END
62afaefe49SHong Zhang 
63afaefe49SHong Zhang EXTERN_C_BEGIN
64afaefe49SHong Zhang #undef __FUNCT__
652401956bSBarry Smith #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal_LU"
662401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal_LU(PC pc,PetscReal z)
679b54502bSHong Zhang {
689b54502bSHong Zhang   PC_LU *lu = (PC_LU*)pc->data;
699b54502bSHong Zhang 
709b54502bSHong Zhang   PetscFunctionBegin;
719b54502bSHong Zhang   lu->nonzerosalongdiagonal = PETSC_TRUE;
729b54502bSHong Zhang   if (z == PETSC_DECIDE) {
739b54502bSHong Zhang     lu->nonzerosalongdiagonaltol = 1.e-10;
749b54502bSHong Zhang   } else {
759b54502bSHong Zhang     lu->nonzerosalongdiagonaltol = z;
769b54502bSHong Zhang   }
779b54502bSHong Zhang   PetscFunctionReturn(0);
789b54502bSHong Zhang }
799b54502bSHong Zhang EXTERN_C_END
809b54502bSHong Zhang 
819b54502bSHong Zhang EXTERN_C_BEGIN
829b54502bSHong Zhang #undef __FUNCT__
832401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering_LU"
842401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering_LU(PC pc,PetscTruth flag)
859b54502bSHong Zhang {
869b54502bSHong Zhang   PC_LU *lu;
879b54502bSHong Zhang 
889b54502bSHong Zhang   PetscFunctionBegin;
899b54502bSHong Zhang   lu                = (PC_LU*)pc->data;
909b54502bSHong Zhang   lu->reuseordering = flag;
919b54502bSHong Zhang   PetscFunctionReturn(0);
929b54502bSHong Zhang }
939b54502bSHong Zhang EXTERN_C_END
949b54502bSHong Zhang 
959b54502bSHong Zhang EXTERN_C_BEGIN
969b54502bSHong Zhang #undef __FUNCT__
972401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseFill_LU"
982401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill_LU(PC pc,PetscTruth flag)
999b54502bSHong Zhang {
1009b54502bSHong Zhang   PC_LU *lu;
1019b54502bSHong Zhang 
1029b54502bSHong Zhang   PetscFunctionBegin;
1039b54502bSHong Zhang   lu = (PC_LU*)pc->data;
1049b54502bSHong Zhang   lu->reusefill = flag;
1059b54502bSHong Zhang   PetscFunctionReturn(0);
1069b54502bSHong Zhang }
1079b54502bSHong Zhang EXTERN_C_END
1089b54502bSHong Zhang 
1099b54502bSHong Zhang #undef __FUNCT__
1109b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_LU"
1119b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_LU(PC pc)
1129b54502bSHong Zhang {
1139b54502bSHong Zhang   PC_LU          *lu = (PC_LU*)pc->data;
1149b54502bSHong Zhang   PetscErrorCode ierr;
1159b54502bSHong Zhang   PetscTruth     flg,set;
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);
1232401956bSBarry Smith     ierr = PetscOptionsName("-pc_factor_in_place","Form LU in the same memory as the matrix","PCFactorSetUseInPlace",&flg);CHKERRQ(ierr);
1249b54502bSHong Zhang     if (flg) {
1252401956bSBarry Smith       ierr = PCFactorSetUseInPlace(pc);CHKERRQ(ierr);
1269b54502bSHong Zhang     }
12755ba2a51SBarry 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 
1402401956bSBarry Smith     ierr = PetscOptionsName("-pc_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",&flg);CHKERRQ(ierr);
1419b54502bSHong Zhang     if (flg) {
1422401956bSBarry Smith       ierr = PCFactorSetReuseFill(pc,PETSC_TRUE);CHKERRQ(ierr);
1439b54502bSHong Zhang     }
1442401956bSBarry Smith     ierr = PetscOptionsName("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",&flg);CHKERRQ(ierr);
1459b54502bSHong Zhang     if (flg) {
1462401956bSBarry Smith       ierr = PCFactorSetReuseOrdering(pc,PETSC_TRUE);CHKERRQ(ierr);
1479b54502bSHong Zhang     }
1489b54502bSHong Zhang 
1499b54502bSHong Zhang     ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr);
1502401956bSBarry Smith     ierr = PetscOptionsList("-pc_factor_mat_ordering_type","Reordering to reduce nonzeros in LU","PCFactorSetMatOrdering",ordlist,lu->ordering,tname,256,&flg);CHKERRQ(ierr);
1519b54502bSHong Zhang     if (flg) {
1522401956bSBarry Smith       ierr = PCFactorSetMatOrdering(pc,tname);CHKERRQ(ierr);
1539b54502bSHong Zhang     }
1549b54502bSHong Zhang 
1552401956bSBarry Smith     ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr);
1569b54502bSHong Zhang     if (flg) {
1579b54502bSHong Zhang       tol = PETSC_DECIDE;
1582401956bSBarry Smith       ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr);
1592401956bSBarry Smith       ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr);
1609b54502bSHong Zhang     }
1619b54502bSHong Zhang 
1622401956bSBarry Smith     ierr = PetscOptionsReal("-pc_factor_pivoting","Pivoting tolerance (used only for some factorization)","PCFactorSetPivoting",lu->info.dtcol,&lu->info.dtcol,&flg);CHKERRQ(ierr);
1639b54502bSHong Zhang 
1649b54502bSHong Zhang     flg = lu->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE;
1652401956bSBarry Smith     ierr = PetscOptionsTruth("-pc_factor_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCFactorSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr);
1669b54502bSHong Zhang     if (set) {
1672401956bSBarry Smith       ierr = PCFactorSetPivotInBlocks(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__
32655ba2a51SBarry Smith #define __FUNCT__ "PCFactorSetFill_LU"
32755ba2a51SBarry 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__
3402401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace_LU"
3412401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_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__
3542401956bSBarry Smith #define __FUNCT__ "PCFactorSetMatOrdering_LU"
3552401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrdering_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__
3692401956bSBarry Smith #define __FUNCT__ "PCFactorSetPivoting_LU"
3702401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivoting_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__
3832401956bSBarry Smith #define __FUNCT__ "PCFactorSetPivotInBlocks_LU"
3842401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks_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__
3972401956bSBarry Smith #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal"
3989b54502bSHong Zhang /*@
3992401956bSBarry Smith    PCFactorReorderForNonzeroDiagonal - 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:
4082401956bSBarry Smith .  -pc_factor_nonzeros_along_diagonal
4099b54502bSHong Zhang 
4109b54502bSHong Zhang    Level: intermediate
4119b54502bSHong Zhang 
4129b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill
4139b54502bSHong Zhang 
41455ba2a51SBarry Smith .seealso: PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetZeroPivot(), MatReorderForNonzeroDiagonal()
4159b54502bSHong Zhang @*/
4162401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal(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);
4222401956bSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_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__
43055ba2a51SBarry Smith #define __FUNCT__ "PCFactorSetFill"
4319b54502bSHong Zhang /*@
43255ba2a51SBarry Smith    PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix,
4339b54502bSHong Zhang    fill = number nonzeros in factor/number nonzeros in original matrix.
4349b54502bSHong Zhang 
4359b54502bSHong Zhang    Collective on PC
4369b54502bSHong Zhang 
4379b54502bSHong Zhang    Input Parameters:
4389b54502bSHong Zhang +  pc - the preconditioner context
4399b54502bSHong Zhang -  fill - amount of expected fill
4409b54502bSHong Zhang 
4419b54502bSHong Zhang    Options Database Key:
44255ba2a51SBarry Smith .  -pc_factor_fill <fill> - Sets fill amount
4439b54502bSHong Zhang 
4449b54502bSHong Zhang    Level: intermediate
4459b54502bSHong Zhang 
4469b54502bSHong Zhang    Note:
4479b54502bSHong Zhang    For sparse matrix factorizations it is difficult to predict how much
448*6cf91177SBarry Smith    fill to expect. By running with the option -info PETSc will print the
4499b54502bSHong Zhang    actual amount of fill used; allowing you to set the value accurately for
4509b54502bSHong Zhang    future runs. Default PETSc uses a value of 5.0
4519b54502bSHong Zhang 
4529b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill
4539b54502bSHong Zhang 
4549b54502bSHong Zhang @*/
45555ba2a51SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill(PC pc,PetscReal fill)
4569b54502bSHong Zhang {
4579b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscReal);
4589b54502bSHong Zhang 
4599b54502bSHong Zhang   PetscFunctionBegin;
4609b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
4619b54502bSHong Zhang   if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0");
46255ba2a51SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetFill_C",(void (**)(void))&f);CHKERRQ(ierr);
4639b54502bSHong Zhang   if (f) {
4649b54502bSHong Zhang     ierr = (*f)(pc,fill);CHKERRQ(ierr);
4659b54502bSHong Zhang   }
4669b54502bSHong Zhang   PetscFunctionReturn(0);
4679b54502bSHong Zhang }
4689b54502bSHong Zhang 
4699b54502bSHong Zhang #undef __FUNCT__
4702401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace"
4719b54502bSHong Zhang /*@
4722401956bSBarry Smith    PCFactorSetUseInPlace - Tells the system to do an in-place factorization.
4739b54502bSHong Zhang    For dense matrices, this enables the solution of much larger problems.
4749b54502bSHong Zhang    For sparse matrices the factorization cannot be done truly in-place
4759b54502bSHong Zhang    so this does not save memory during the factorization, but after the matrix
4769b54502bSHong Zhang    is factored, the original unfactored matrix is freed, thus recovering that
4779b54502bSHong Zhang    space.
4789b54502bSHong Zhang 
4799b54502bSHong Zhang    Collective on PC
4809b54502bSHong Zhang 
4819b54502bSHong Zhang    Input Parameters:
4829b54502bSHong Zhang .  pc - the preconditioner context
4839b54502bSHong Zhang 
4849b54502bSHong Zhang    Options Database Key:
4852401956bSBarry Smith .  -pc_factor_in_place - Activates in-place factorization
4869b54502bSHong Zhang 
4879b54502bSHong Zhang    Notes:
4882401956bSBarry Smith    PCFactorSetUseInplace() can only be used with the KSP method KSPPREONLY or when
4899b54502bSHong Zhang    a different matrix is provided for the multiply and the preconditioner in
4909b54502bSHong Zhang    a call to KSPSetOperators().
4919b54502bSHong Zhang    This is because the Krylov space methods require an application of the
4929b54502bSHong Zhang    matrix multiplication, which is not possible here because the matrix has
4939b54502bSHong Zhang    been factored in-place, replacing the original matrix.
4949b54502bSHong Zhang 
4959b54502bSHong Zhang    Level: intermediate
4969b54502bSHong Zhang 
4979b54502bSHong Zhang .keywords: PC, set, factorization, direct, inplace, in-place, LU
4989b54502bSHong Zhang 
4999b54502bSHong Zhang .seealso: PCILUSetUseInPlace()
5009b54502bSHong Zhang @*/
5012401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace(PC pc)
5029b54502bSHong Zhang {
5039b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC);
5049b54502bSHong Zhang 
5059b54502bSHong Zhang   PetscFunctionBegin;
5069b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
5072401956bSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",(void (**)(void))&f);CHKERRQ(ierr);
5089b54502bSHong Zhang   if (f) {
5099b54502bSHong Zhang     ierr = (*f)(pc);CHKERRQ(ierr);
5109b54502bSHong Zhang   }
5119b54502bSHong Zhang   PetscFunctionReturn(0);
5129b54502bSHong Zhang }
5139b54502bSHong Zhang 
5149b54502bSHong Zhang #undef __FUNCT__
5152401956bSBarry Smith #define __FUNCT__ "PCFactorSetMatOrdering"
5169b54502bSHong Zhang /*@C
5172401956bSBarry Smith     PCFactorSetMatOrdering - Sets the ordering routine (to reduce fill) to
5189b54502bSHong Zhang     be used in the LU factorization.
5199b54502bSHong Zhang 
5209b54502bSHong Zhang     Collective on PC
5219b54502bSHong Zhang 
5229b54502bSHong Zhang     Input Parameters:
5239b54502bSHong Zhang +   pc - the preconditioner context
5249b54502bSHong Zhang -   ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM
5259b54502bSHong Zhang 
5269b54502bSHong Zhang     Options Database Key:
5272401956bSBarry Smith .   -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
5289b54502bSHong Zhang 
5299b54502bSHong Zhang     Level: intermediate
5309b54502bSHong Zhang 
5319b54502bSHong Zhang     Notes: nested dissection is used by default
5329b54502bSHong Zhang 
5339b54502bSHong Zhang .seealso: PCILUSetMatOrdering()
5349b54502bSHong Zhang @*/
5352401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrdering(PC pc,MatOrderingType ordering)
5369b54502bSHong Zhang {
5379b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,MatOrderingType);
5389b54502bSHong Zhang 
5399b54502bSHong Zhang   PetscFunctionBegin;
5402401956bSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatOrdering_C",(void (**)(void))&f);CHKERRQ(ierr);
5419b54502bSHong Zhang   if (f) {
5429b54502bSHong Zhang     ierr = (*f)(pc,ordering);CHKERRQ(ierr);
5439b54502bSHong Zhang   }
5449b54502bSHong Zhang   PetscFunctionReturn(0);
5459b54502bSHong Zhang }
5469b54502bSHong Zhang 
5479b54502bSHong Zhang #undef __FUNCT__
5482401956bSBarry Smith #define __FUNCT__ "PCFactorSetPivoting"
5499b54502bSHong Zhang /*@
5502401956bSBarry Smith     PCFactorSetPivoting - Determines when pivoting is done during LU.
5519b54502bSHong Zhang       For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices
5529b54502bSHong Zhang       it is never done. For the Matlab and SuperLU factorization this is used.
5539b54502bSHong Zhang 
5549b54502bSHong Zhang     Collective on PC
5559b54502bSHong Zhang 
5569b54502bSHong Zhang     Input Parameters:
5579b54502bSHong Zhang +   pc - the preconditioner context
5589b54502bSHong Zhang -   dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable)
5599b54502bSHong Zhang 
5609b54502bSHong Zhang     Options Database Key:
5612401956bSBarry Smith .   -pc_factor_pivoting <dtcol>
5629b54502bSHong Zhang 
5639b54502bSHong Zhang     Level: intermediate
5649b54502bSHong Zhang 
5652401956bSBarry Smith .seealso: PCILUSetMatOrdering(), PCFactorSetPivotInBlocks()
5669b54502bSHong Zhang @*/
5672401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivoting(PC pc,PetscReal dtcol)
5689b54502bSHong Zhang {
5699b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscReal);
5709b54502bSHong Zhang 
5719b54502bSHong Zhang   PetscFunctionBegin;
5722401956bSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivoting_C",(void (**)(void))&f);CHKERRQ(ierr);
5739b54502bSHong Zhang   if (f) {
5749b54502bSHong Zhang     ierr = (*f)(pc,dtcol);CHKERRQ(ierr);
5759b54502bSHong Zhang   }
5769b54502bSHong Zhang   PetscFunctionReturn(0);
5779b54502bSHong Zhang }
5789b54502bSHong Zhang 
5799b54502bSHong Zhang #undef __FUNCT__
5802401956bSBarry Smith #define __FUNCT__ "PCFactorSetPivotInBlocks"
5819b54502bSHong Zhang /*@
5822401956bSBarry Smith     PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block
5839b54502bSHong Zhang       with BAIJ or SBAIJ matrices
5849b54502bSHong Zhang 
5859b54502bSHong Zhang     Collective on PC
5869b54502bSHong Zhang 
5879b54502bSHong Zhang     Input Parameters:
5889b54502bSHong Zhang +   pc - the preconditioner context
5899b54502bSHong Zhang -   pivot - PETSC_TRUE or PETSC_FALSE
5909b54502bSHong Zhang 
5919b54502bSHong Zhang     Options Database Key:
5922401956bSBarry Smith .   -pc_factor_pivot_in_blocks <true,false>
5939b54502bSHong Zhang 
5949b54502bSHong Zhang     Level: intermediate
5959b54502bSHong Zhang 
5962401956bSBarry Smith .seealso: PCILUSetMatOrdering(), PCFactorSetPivoting()
5979b54502bSHong Zhang @*/
5982401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks(PC pc,PetscTruth pivot)
5999b54502bSHong Zhang {
6009b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscTruth);
6019b54502bSHong Zhang 
6029b54502bSHong Zhang   PetscFunctionBegin;
6032401956bSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr);
6049b54502bSHong Zhang   if (f) {
6059b54502bSHong Zhang     ierr = (*f)(pc,pivot);CHKERRQ(ierr);
6069b54502bSHong Zhang   }
6079b54502bSHong Zhang   PetscFunctionReturn(0);
6089b54502bSHong Zhang }
6099b54502bSHong Zhang 
6109b54502bSHong Zhang /* ------------------------------------------------------------------------ */
6119b54502bSHong Zhang 
6129b54502bSHong Zhang /*MC
6139b54502bSHong Zhang    PCLU - Uses a direct solver, based on LU factorization, as a preconditioner
6149b54502bSHong Zhang 
6159b54502bSHong Zhang    Options Database Keys:
6162401956bSBarry Smith +  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
6172401956bSBarry Smith .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
61855ba2a51SBarry Smith .  -pc_factor_fill <fill> - Sets fill amount
6192401956bSBarry Smith .  -pc_factor_in_place - Activates in-place factorization
6202401956bSBarry Smith .  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
6212401956bSBarry Smith .  -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase
6229b54502bSHong Zhang                                          stability of factorization.
623f251bdbdSHong Zhang .  -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default
624f251bdbdSHong Zhang -  -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value
625f251bdbdSHong Zhang    is optional with PETSC_TRUE being the default
6269b54502bSHong Zhang 
6279b54502bSHong Zhang    Notes: Not all options work for all matrix formats
6289b54502bSHong Zhang           Run with -help to see additional options for particular matrix formats or factorization
6299b54502bSHong Zhang           algorithms
6309b54502bSHong Zhang 
6319b54502bSHong Zhang    Level: beginner
6329b54502bSHong Zhang 
6339b54502bSHong Zhang    Concepts: LU factorization, direct solver
6349b54502bSHong Zhang 
6359b54502bSHong Zhang    Notes: Usually this will compute an "exact" solution in one iteration and does
6369b54502bSHong Zhang           not need a Krylov method (i.e. you can use -ksp_type preonly, or
6379b54502bSHong Zhang           KSPSetType(ksp,KSPPREONLY) for the Krylov method
6389b54502bSHong Zhang 
6399b54502bSHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
6402401956bSBarry Smith            PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCGetFactoredMatrix(),
6412401956bSBarry Smith            PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrdering(), PCFactorSetPivoting(),
6422401956bSBarry Smith            PCFactorSetPivotingInBlocks(),PCFactorSetShiftNonzero(),PCFactorSetShiftPd()
6439b54502bSHong Zhang M*/
6449b54502bSHong Zhang 
6459b54502bSHong Zhang EXTERN_C_BEGIN
6469b54502bSHong Zhang #undef __FUNCT__
6479b54502bSHong Zhang #define __FUNCT__ "PCCreate_LU"
648dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_LU(PC pc)
6499b54502bSHong Zhang {
6509b54502bSHong Zhang   PetscErrorCode ierr;
6519b54502bSHong Zhang   PetscMPIInt    size;
6529b54502bSHong Zhang   PC_LU          *dir;
6539b54502bSHong Zhang 
6549b54502bSHong Zhang   PetscFunctionBegin;
6559b54502bSHong Zhang   ierr = PetscNew(PC_LU,&dir);CHKERRQ(ierr);
65652e6d16bSBarry Smith   ierr = PetscLogObjectMemory(pc,sizeof(PC_LU));CHKERRQ(ierr);
6579b54502bSHong Zhang 
6589b54502bSHong Zhang   ierr = MatFactorInfoInitialize(&dir->info);CHKERRQ(ierr);
6599b54502bSHong Zhang   dir->fact                  = 0;
6609b54502bSHong Zhang   dir->inplace               = PETSC_FALSE;
6619b54502bSHong Zhang   dir->nonzerosalongdiagonal = PETSC_FALSE;
6629b54502bSHong Zhang 
6639b54502bSHong Zhang   dir->info.fill           = 5.0;
6649b54502bSHong Zhang   dir->info.dtcol          = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */
6650a29876aSHong Zhang   dir->info.shiftnz        = 0.0;
6669b54502bSHong Zhang   dir->info.zeropivot      = 1.e-12;
6679b54502bSHong Zhang   dir->info.pivotinblocks  = 1.0;
668fbf22428SSatish Balay   dir->info.shiftpd        = 0.0; /* false */
6699b54502bSHong Zhang   dir->info.shift_fraction = 0.0;
6709b54502bSHong Zhang   dir->col                 = 0;
6719b54502bSHong Zhang   dir->row                 = 0;
6729b54502bSHong Zhang   ierr = MPI_Comm_size(pc->comm,&size);CHKERRQ(ierr);
6739b54502bSHong Zhang   if (size == 1) {
6749b54502bSHong Zhang     ierr = PetscStrallocpy(MATORDERING_ND,&dir->ordering);CHKERRQ(ierr);
6759b54502bSHong Zhang   } else {
6769b54502bSHong Zhang     ierr = PetscStrallocpy(MATORDERING_NATURAL,&dir->ordering);CHKERRQ(ierr);
6779b54502bSHong Zhang   }
6789b54502bSHong Zhang   dir->reusefill        = PETSC_FALSE;
6799b54502bSHong Zhang   dir->reuseordering    = PETSC_FALSE;
6809b54502bSHong Zhang   pc->data              = (void*)dir;
6819b54502bSHong Zhang 
6829b54502bSHong Zhang   pc->ops->destroy           = PCDestroy_LU;
6839b54502bSHong Zhang   pc->ops->apply             = PCApply_LU;
6849b54502bSHong Zhang   pc->ops->applytranspose    = PCApplyTranspose_LU;
6859b54502bSHong Zhang   pc->ops->setup             = PCSetUp_LU;
6869b54502bSHong Zhang   pc->ops->setfromoptions    = PCSetFromOptions_LU;
6879b54502bSHong Zhang   pc->ops->view              = PCView_LU;
6889b54502bSHong Zhang   pc->ops->applyrichardson   = 0;
6899b54502bSHong Zhang   pc->ops->getfactoredmatrix = PCGetFactoredMatrix_LU;
6909b54502bSHong Zhang 
691afaefe49SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_LU",
692afaefe49SHong Zhang                     PCFactorSetZeroPivot_LU);CHKERRQ(ierr);
693afaefe49SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_LU",
694afaefe49SHong Zhang                     PCFactorSetShiftNonzero_LU);CHKERRQ(ierr);
695afaefe49SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_LU",
696afaefe49SHong Zhang                     PCFactorSetShiftPd_LU);CHKERRQ(ierr);
697afaefe49SHong Zhang 
69855ba2a51SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_LU",
69955ba2a51SBarry Smith                     PCFactorSetFill_LU);CHKERRQ(ierr);
7002401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_LU",
7012401956bSBarry Smith                     PCFactorSetUseInPlace_LU);CHKERRQ(ierr);
7022401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrdering_C","PCFactorSetMatOrdering_LU",
7032401956bSBarry Smith                     PCFactorSetMatOrdering_LU);CHKERRQ(ierr);
7042401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_LU",
7052401956bSBarry Smith                     PCFactorSetReuseOrdering_LU);CHKERRQ(ierr);
7062401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_LU",
7072401956bSBarry Smith                     PCFactorSetReuseFill_LU);CHKERRQ(ierr);
7082401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivoting_C","PCFactorSetPivoting_LU",
7092401956bSBarry Smith                     PCFactorSetPivoting_LU);CHKERRQ(ierr);
7102401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_LU",
7112401956bSBarry Smith                     PCFactorSetPivotInBlocks_LU);CHKERRQ(ierr);
7122401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_LU",
7132401956bSBarry Smith                     PCFactorReorderForNonzeroDiagonal_LU);CHKERRQ(ierr);
7149b54502bSHong Zhang   PetscFunctionReturn(0);
7159b54502bSHong Zhang }
7169b54502bSHong Zhang EXTERN_C_END
717