xref: /petsc/src/ksp/pc/impls/factor/lu/lu.c (revision 35bd34fa3344fb53ae05f23bf99eddd2f9bb15e6)
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;
21*35bd34faSBarry Smith   if (lu->fact) {
22*35bd34faSBarry Smith     const MatSolverPackage ltype;
23*35bd34faSBarry Smith     PetscTruth             flg;
24*35bd34faSBarry Smith     ierr = MatFactorGetSolverPackage(lu->fact,&ltype);CHKERRQ(ierr);
25*35bd34faSBarry Smith     ierr = PetscStrcmp(stype,ltype,&flg);CHKERRQ(ierr);
26*35bd34faSBarry Smith     if (!flg) {
27*35bd34faSBarry Smith       SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot change solver matrix package after PC has been setup or used");
28*35bd34faSBarry Smith     }
29*35bd34faSBarry Smith   }
3031a0a2fdSHong Zhang   ierr = PetscStrfree(lu->solvertype);CHKERRQ(ierr);
31c62fe420SBarry Smith   ierr = PetscStrallocpy(stype,&lu->solvertype);CHKERRQ(ierr);
32c62fe420SBarry Smith   PetscFunctionReturn(0);
33c62fe420SBarry Smith }
34c62fe420SBarry Smith EXTERN_C_END
35c62fe420SBarry Smith 
36c62fe420SBarry Smith EXTERN_C_BEGIN
37c62fe420SBarry Smith #undef __FUNCT__
38afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetZeroPivot_LU"
39dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot_LU(PC pc,PetscReal z)
40afaefe49SHong Zhang {
41c62fe420SBarry Smith   PC_LU *lu = (PC_LU*)pc->data;
42afaefe49SHong Zhang 
43afaefe49SHong Zhang   PetscFunctionBegin;
44afaefe49SHong Zhang   lu->info.zeropivot = z;
45afaefe49SHong Zhang   PetscFunctionReturn(0);
46afaefe49SHong Zhang }
47afaefe49SHong Zhang EXTERN_C_END
48afaefe49SHong Zhang 
49afaefe49SHong Zhang EXTERN_C_BEGIN
50afaefe49SHong Zhang #undef __FUNCT__
51afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetShiftNonzero_LU"
52dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero_LU(PC pc,PetscReal shift)
53afaefe49SHong Zhang {
54c62fe420SBarry Smith   PC_LU *dir = (PC_LU*)pc->data;
55afaefe49SHong Zhang 
56afaefe49SHong Zhang   PetscFunctionBegin;
57afaefe49SHong Zhang   if (shift == (PetscReal) PETSC_DECIDE) {
58afaefe49SHong Zhang     dir->info.shiftnz = 1.e-12;
59afaefe49SHong Zhang   } else {
60afaefe49SHong Zhang     dir->info.shiftnz = shift;
61afaefe49SHong Zhang   }
62afaefe49SHong Zhang   PetscFunctionReturn(0);
63afaefe49SHong Zhang }
64afaefe49SHong Zhang EXTERN_C_END
65afaefe49SHong Zhang 
66afaefe49SHong Zhang EXTERN_C_BEGIN
67afaefe49SHong Zhang #undef __FUNCT__
68afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetShiftPd_LU"
69dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd_LU(PC pc,PetscTruth shift)
70afaefe49SHong Zhang {
71c62fe420SBarry Smith   PC_LU *dir = (PC_LU*)pc->data;
72afaefe49SHong Zhang 
73afaefe49SHong Zhang   PetscFunctionBegin;
74fbf22428SSatish Balay   if (shift) {
75fbf22428SSatish Balay     dir->info.shiftpd = 1.0;
76fbf22428SSatish Balay   } else {
77fbf22428SSatish Balay     dir->info.shiftpd = 0.0;
78fbf22428SSatish Balay   }
79afaefe49SHong Zhang   PetscFunctionReturn(0);
80afaefe49SHong Zhang }
81afaefe49SHong Zhang EXTERN_C_END
82afaefe49SHong Zhang 
83afaefe49SHong Zhang EXTERN_C_BEGIN
84afaefe49SHong Zhang #undef __FUNCT__
852401956bSBarry Smith #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal_LU"
862401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal_LU(PC pc,PetscReal z)
879b54502bSHong Zhang {
889b54502bSHong Zhang   PC_LU *lu = (PC_LU*)pc->data;
899b54502bSHong Zhang 
909b54502bSHong Zhang   PetscFunctionBegin;
919b54502bSHong Zhang   lu->nonzerosalongdiagonal = PETSC_TRUE;
929b54502bSHong Zhang   if (z == PETSC_DECIDE) {
939b54502bSHong Zhang     lu->nonzerosalongdiagonaltol = 1.e-10;
949b54502bSHong Zhang   } else {
959b54502bSHong Zhang     lu->nonzerosalongdiagonaltol = z;
969b54502bSHong Zhang   }
979b54502bSHong Zhang   PetscFunctionReturn(0);
989b54502bSHong Zhang }
999b54502bSHong Zhang EXTERN_C_END
1009b54502bSHong Zhang 
1019b54502bSHong Zhang EXTERN_C_BEGIN
1029b54502bSHong Zhang #undef __FUNCT__
1032401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering_LU"
1042401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering_LU(PC pc,PetscTruth flag)
1059b54502bSHong Zhang {
106c62fe420SBarry Smith   PC_LU *lu = (PC_LU*)pc->data;
1079b54502bSHong Zhang 
1089b54502bSHong Zhang   PetscFunctionBegin;
1099b54502bSHong Zhang   lu->reuseordering = flag;
1109b54502bSHong Zhang   PetscFunctionReturn(0);
1119b54502bSHong Zhang }
1129b54502bSHong Zhang EXTERN_C_END
1139b54502bSHong Zhang 
1149b54502bSHong Zhang EXTERN_C_BEGIN
1159b54502bSHong Zhang #undef __FUNCT__
1162401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseFill_LU"
1172401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill_LU(PC pc,PetscTruth flag)
1189b54502bSHong Zhang {
119c62fe420SBarry Smith   PC_LU *lu = (PC_LU*)pc->data;
1209b54502bSHong Zhang 
1219b54502bSHong Zhang   PetscFunctionBegin;
1229b54502bSHong Zhang   lu->reusefill = flag;
1239b54502bSHong Zhang   PetscFunctionReturn(0);
1249b54502bSHong Zhang }
1259b54502bSHong Zhang EXTERN_C_END
1269b54502bSHong Zhang 
1279b54502bSHong Zhang #undef __FUNCT__
1289b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_LU"
1299b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_LU(PC pc)
1309b54502bSHong Zhang {
1319b54502bSHong Zhang   PC_LU           *lu = (PC_LU*)pc->data;
1329b54502bSHong Zhang   PetscErrorCode  ierr;
1339b54502bSHong Zhang   PetscTruth      flg,set;
1349989ab13SBarry Smith   char            tname[256], solvertype[64];
1359b54502bSHong Zhang   PetscFList      ordlist;
1369b54502bSHong Zhang   PetscReal       tol;
1379b54502bSHong Zhang 
1389b54502bSHong Zhang   PetscFunctionBegin;
1399b54502bSHong Zhang   ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr);
1409b54502bSHong Zhang   ierr = PetscOptionsHead("LU options");CHKERRQ(ierr);
1412401956bSBarry Smith     ierr = PetscOptionsName("-pc_factor_in_place","Form LU in the same memory as the matrix","PCFactorSetUseInPlace",&flg);CHKERRQ(ierr);
1429b54502bSHong Zhang     if (flg) {
1432401956bSBarry Smith       ierr = PCFactorSetUseInPlace(pc);CHKERRQ(ierr);
1449b54502bSHong Zhang     }
14555ba2a51SBarry 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);
1469b54502bSHong Zhang 
1479f95998fSHong Zhang     ierr = PetscOptionsName("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr);
1489b54502bSHong Zhang     if (flg) {
149afaefe49SHong Zhang         ierr = PCFactorSetShiftNonzero(pc,(PetscReal) PETSC_DECIDE);CHKERRQ(ierr);
1509b54502bSHong Zhang     }
1519f95998fSHong Zhang     ierr = PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",lu->info.shiftnz,&lu->info.shiftnz,0);CHKERRQ(ierr);
1529f95998fSHong Zhang     ierr = PetscOptionsName("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",&flg);CHKERRQ(ierr);
1539b54502bSHong Zhang     if (flg) {
154afaefe49SHong Zhang       ierr = PCFactorSetShiftPd(pc,PETSC_TRUE);CHKERRQ(ierr);
1559b54502bSHong Zhang     }
156ee45ca4aSHong Zhang     ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",lu->info.zeropivot,&lu->info.zeropivot,0);CHKERRQ(ierr);
1579b54502bSHong Zhang 
1582401956bSBarry Smith     ierr = PetscOptionsName("-pc_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",&flg);CHKERRQ(ierr);
1599b54502bSHong Zhang     if (flg) {
1602401956bSBarry Smith       ierr = PCFactorSetReuseFill(pc,PETSC_TRUE);CHKERRQ(ierr);
1619b54502bSHong Zhang     }
1622401956bSBarry Smith     ierr = PetscOptionsName("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",&flg);CHKERRQ(ierr);
1639b54502bSHong Zhang     if (flg) {
1642401956bSBarry Smith       ierr = PCFactorSetReuseOrdering(pc,PETSC_TRUE);CHKERRQ(ierr);
1659b54502bSHong Zhang     }
1669b54502bSHong Zhang 
1679b54502bSHong Zhang     ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr);
168e5a9bf91SBarry Smith     ierr = PetscOptionsList("-pc_factor_mat_ordering_type","Reordering to reduce nonzeros in LU","PCFactorSetMatOrderingType",ordlist,lu->ordering,tname,256,&flg);CHKERRQ(ierr);
1699b54502bSHong Zhang     if (flg) {
170e5a9bf91SBarry Smith       ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr);
1719b54502bSHong Zhang     }
1729b54502bSHong Zhang 
1735c9eb25fSBarry Smith     /* maybe should have MatGetSolverTypes(Mat,&list) like the ordering list */
174c7393fdbSBarry Smith     ierr = PetscOptionsString("-pc_factor_mat_solver_package","Specific LU solver to use","MatGetFactor",lu->solvertype,solvertype,64,&flg);CHKERRQ(ierr);
1755c9eb25fSBarry Smith     if (flg) {
176c7393fdbSBarry Smith       ierr = PCFactorSetMatSolverPackage(pc,solvertype);CHKERRQ(ierr);
1775c9eb25fSBarry Smith     }
1785c9eb25fSBarry Smith 
1792401956bSBarry Smith     ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr);
1809b54502bSHong Zhang     if (flg) {
1819b54502bSHong Zhang       tol = PETSC_DECIDE;
1822401956bSBarry Smith       ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr);
1832401956bSBarry Smith       ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr);
1849b54502bSHong Zhang     }
1859b54502bSHong Zhang 
1862401956bSBarry Smith     ierr = PetscOptionsReal("-pc_factor_pivoting","Pivoting tolerance (used only for some factorization)","PCFactorSetPivoting",lu->info.dtcol,&lu->info.dtcol,&flg);CHKERRQ(ierr);
1879b54502bSHong Zhang 
1889b54502bSHong Zhang     flg = lu->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE;
1892401956bSBarry Smith     ierr = PetscOptionsTruth("-pc_factor_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCFactorSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr);
1909b54502bSHong Zhang     if (set) {
1912401956bSBarry Smith       ierr = PCFactorSetPivotInBlocks(pc,flg);CHKERRQ(ierr);
1929b54502bSHong Zhang     }
1939b54502bSHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
1949b54502bSHong Zhang   PetscFunctionReturn(0);
1959b54502bSHong Zhang }
1969b54502bSHong Zhang 
1979b54502bSHong Zhang #undef __FUNCT__
1989b54502bSHong Zhang #define __FUNCT__ "PCView_LU"
1999b54502bSHong Zhang static PetscErrorCode PCView_LU(PC pc,PetscViewer viewer)
2009b54502bSHong Zhang {
2019b54502bSHong Zhang   PC_LU          *lu = (PC_LU*)pc->data;
2029b54502bSHong Zhang   PetscErrorCode ierr;
2039b54502bSHong Zhang   PetscTruth     iascii,isstring;
2049b54502bSHong Zhang 
2059b54502bSHong Zhang   PetscFunctionBegin;
2069b54502bSHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
2079b54502bSHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
2089b54502bSHong Zhang   if (iascii) {
2099b54502bSHong Zhang 
2109b54502bSHong Zhang     if (lu->inplace) {ierr = PetscViewerASCIIPrintf(viewer,"  LU: in-place factorization\n");CHKERRQ(ierr);}
2119b54502bSHong Zhang     else             {ierr = PetscViewerASCIIPrintf(viewer,"  LU: out-of-place factorization\n");CHKERRQ(ierr);}
2129b54502bSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"    matrix ordering: %s\n",lu->ordering);CHKERRQ(ierr);
213a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  LU: tolerance for zero pivot %G\n",lu->info.zeropivot);CHKERRQ(ierr);
2140a29876aSHong Zhang     if (lu->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer,"  LU: using Manteuffel shift\n");CHKERRQ(ierr);}
2159b54502bSHong Zhang     if (lu->fact) {
216f3a39becSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  LU: factor fill ratio needed %G\n",lu->actualfill);CHKERRQ(ierr);
217f3a39becSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"       Factored matrix follows\n");CHKERRQ(ierr);
218f3a39becSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
219f3a39becSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
220f3a39becSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
221f3a39becSBarry Smith       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
2229b54502bSHong Zhang       ierr = MatView(lu->fact,viewer);CHKERRQ(ierr);
2239b54502bSHong Zhang       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
224f3a39becSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
225f3a39becSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
226f3a39becSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
2279b54502bSHong Zhang     }
2289b54502bSHong Zhang     if (lu->reusefill)    {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing fill from past factorization\n");CHKERRQ(ierr);}
2299b54502bSHong Zhang     if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing reordering from past factorization\n");CHKERRQ(ierr);}
2309b54502bSHong Zhang   } else if (isstring) {
2319b54502bSHong Zhang     ierr = PetscViewerStringSPrintf(viewer," order=%s",lu->ordering);CHKERRQ(ierr);CHKERRQ(ierr);
2329b54502bSHong Zhang   } else {
2339b54502bSHong Zhang     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCLU",((PetscObject)viewer)->type_name);
2349b54502bSHong Zhang   }
2359b54502bSHong Zhang   PetscFunctionReturn(0);
2369b54502bSHong Zhang }
2379b54502bSHong Zhang 
2389b54502bSHong Zhang #undef __FUNCT__
239a4fd02acSBarry Smith #define __FUNCT__ "PCFactorGetMatrix_LU"
240a4fd02acSBarry Smith static PetscErrorCode PCFactorGetMatrix_LU(PC pc,Mat *mat)
2419b54502bSHong Zhang {
2429b54502bSHong Zhang   PC_LU *dir = (PC_LU*)pc->data;
2439b54502bSHong Zhang 
2449b54502bSHong Zhang   PetscFunctionBegin;
2459b54502bSHong Zhang   if (!dir->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()");
2469b54502bSHong Zhang   *mat = dir->fact;
2479b54502bSHong Zhang   PetscFunctionReturn(0);
2489b54502bSHong Zhang }
2499b54502bSHong Zhang 
2509b54502bSHong Zhang #undef __FUNCT__
2519b54502bSHong Zhang #define __FUNCT__ "PCSetUp_LU"
2529b54502bSHong Zhang static PetscErrorCode PCSetUp_LU(PC pc)
2539b54502bSHong Zhang {
2549b54502bSHong Zhang   PetscErrorCode ierr;
2559b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
2569b54502bSHong Zhang 
2579b54502bSHong Zhang   PetscFunctionBegin;
2589b54502bSHong Zhang   if (dir->reusefill && pc->setupcalled) dir->info.fill = dir->actualfill;
2599b54502bSHong Zhang 
2609b54502bSHong Zhang   if (dir->inplace) {
2619b54502bSHong Zhang     if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);}
2629b54502bSHong Zhang     if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);}
2639b54502bSHong Zhang     ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
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     ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&dir->info);CHKERRQ(ierr);
2699b54502bSHong Zhang     dir->fact = pc->pmat;
2709b54502bSHong Zhang   } else {
2719b54502bSHong Zhang     MatInfo info;
2729b54502bSHong Zhang     if (!pc->setupcalled) {
2739b54502bSHong Zhang       ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
2749b54502bSHong Zhang       if (dir->nonzerosalongdiagonal) {
2759b54502bSHong Zhang         ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr);
2769b54502bSHong Zhang       }
27703c60df9SBarry Smith       if (dir->row) {
27803c60df9SBarry Smith         ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);
27903c60df9SBarry Smith         ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr);
28003c60df9SBarry Smith       }
2815c9eb25fSBarry Smith       ierr = MatGetFactor(pc->pmat,dir->solvertype,MAT_FACTOR_LU,&dir->fact);CHKERRQ(ierr);
282719d5645SBarry Smith       ierr = MatLUFactorSymbolic(dir->fact,pc->pmat,dir->row,dir->col,&dir->info);CHKERRQ(ierr);
2839b54502bSHong Zhang       ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
2849b54502bSHong Zhang       dir->actualfill = info.fill_ratio_needed;
28552e6d16bSBarry Smith       ierr = PetscLogObjectParent(pc,dir->fact);CHKERRQ(ierr);
2869b54502bSHong Zhang     } else if (pc->flag != SAME_NONZERO_PATTERN) {
2879b54502bSHong Zhang       if (!dir->reuseordering) {
2889b54502bSHong Zhang         if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);}
2899b54502bSHong Zhang         if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);}
2909b54502bSHong Zhang         ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
2919b54502bSHong Zhang         if (dir->nonzerosalongdiagonal) {
2929b54502bSHong Zhang           ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr);
2939b54502bSHong Zhang         }
29403c60df9SBarry Smith         if (dir->row) {
29503c60df9SBarry Smith           ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);
29603c60df9SBarry Smith           ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr);
29703c60df9SBarry Smith         }
2989b54502bSHong Zhang       }
2999b54502bSHong Zhang       ierr = MatDestroy(dir->fact);CHKERRQ(ierr);
3005c9eb25fSBarry Smith       ierr = MatGetFactor(pc->pmat,dir->solvertype,MAT_FACTOR_LU,&dir->fact);CHKERRQ(ierr);
301719d5645SBarry Smith       ierr = MatLUFactorSymbolic(dir->fact,pc->pmat,dir->row,dir->col,&dir->info);CHKERRQ(ierr);
3029b54502bSHong Zhang       ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
3039b54502bSHong Zhang       dir->actualfill = info.fill_ratio_needed;
30452e6d16bSBarry Smith       ierr = PetscLogObjectParent(pc,dir->fact);CHKERRQ(ierr);
3059b54502bSHong Zhang     }
306719d5645SBarry Smith     ierr = MatLUFactorNumeric(dir->fact,pc->pmat,&dir->info);CHKERRQ(ierr);
3079b54502bSHong Zhang   }
3089b54502bSHong Zhang   PetscFunctionReturn(0);
3099b54502bSHong Zhang }
3109b54502bSHong Zhang 
3119b54502bSHong Zhang #undef __FUNCT__
3129b54502bSHong Zhang #define __FUNCT__ "PCDestroy_LU"
3139b54502bSHong Zhang static PetscErrorCode PCDestroy_LU(PC pc)
3149b54502bSHong Zhang {
3159b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
3169b54502bSHong Zhang   PetscErrorCode ierr;
3179b54502bSHong Zhang 
3189b54502bSHong Zhang   PetscFunctionBegin;
3199b54502bSHong Zhang   if (!dir->inplace && dir->fact) {ierr = MatDestroy(dir->fact);CHKERRQ(ierr);}
3209b54502bSHong Zhang   if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);}
3219b54502bSHong Zhang   if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);}
3229b54502bSHong Zhang   ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr);
3235c9eb25fSBarry Smith   ierr = PetscStrfree(dir->solvertype);CHKERRQ(ierr);
3249b54502bSHong Zhang   ierr = PetscFree(dir);CHKERRQ(ierr);
3259b54502bSHong Zhang   PetscFunctionReturn(0);
3269b54502bSHong Zhang }
3279b54502bSHong Zhang 
3289b54502bSHong Zhang #undef __FUNCT__
3299b54502bSHong Zhang #define __FUNCT__ "PCApply_LU"
3309b54502bSHong Zhang static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y)
3319b54502bSHong Zhang {
3329b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
3339b54502bSHong Zhang   PetscErrorCode ierr;
3349b54502bSHong Zhang 
3359b54502bSHong Zhang   PetscFunctionBegin;
3369b54502bSHong Zhang   if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);}
3379b54502bSHong Zhang   else              {ierr = MatSolve(dir->fact,x,y);CHKERRQ(ierr);}
3389b54502bSHong Zhang   PetscFunctionReturn(0);
3399b54502bSHong Zhang }
3409b54502bSHong Zhang 
3419b54502bSHong Zhang #undef __FUNCT__
3429b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_LU"
3439b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y)
3449b54502bSHong Zhang {
3459b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
3469b54502bSHong Zhang   PetscErrorCode ierr;
3479b54502bSHong Zhang 
3489b54502bSHong Zhang   PetscFunctionBegin;
3499b54502bSHong Zhang   if (dir->inplace) {ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);}
3509b54502bSHong Zhang   else              {ierr = MatSolveTranspose(dir->fact,x,y);CHKERRQ(ierr);}
3519b54502bSHong Zhang   PetscFunctionReturn(0);
3529b54502bSHong Zhang }
3539b54502bSHong Zhang 
3549b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/
3559b54502bSHong Zhang 
3569b54502bSHong Zhang EXTERN_C_BEGIN
3579b54502bSHong Zhang #undef __FUNCT__
35855ba2a51SBarry Smith #define __FUNCT__ "PCFactorSetFill_LU"
35955ba2a51SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill_LU(PC pc,PetscReal fill)
3609b54502bSHong Zhang {
361c62fe420SBarry Smith   PC_LU *dir = (PC_LU*)pc->data;
3629b54502bSHong Zhang 
3639b54502bSHong Zhang   PetscFunctionBegin;
3649b54502bSHong Zhang   dir->info.fill = fill;
3659b54502bSHong Zhang   PetscFunctionReturn(0);
3669b54502bSHong Zhang }
3679b54502bSHong Zhang EXTERN_C_END
3689b54502bSHong Zhang 
3699b54502bSHong Zhang EXTERN_C_BEGIN
3709b54502bSHong Zhang #undef __FUNCT__
3712401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace_LU"
3722401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_LU(PC pc)
3739b54502bSHong Zhang {
374c62fe420SBarry Smith   PC_LU *dir = (PC_LU*)pc->data;
3759b54502bSHong Zhang 
3769b54502bSHong Zhang   PetscFunctionBegin;
3779b54502bSHong Zhang   dir->inplace = PETSC_TRUE;
3789b54502bSHong Zhang   PetscFunctionReturn(0);
3799b54502bSHong Zhang }
3809b54502bSHong Zhang EXTERN_C_END
3819b54502bSHong Zhang 
3829b54502bSHong Zhang EXTERN_C_BEGIN
3839b54502bSHong Zhang #undef __FUNCT__
384e5a9bf91SBarry Smith #define __FUNCT__ "PCFactorSetMatOrderingType_LU"
385e36ff184SSatish Balay PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType_LU(PC pc,const MatOrderingType ordering)
3869b54502bSHong Zhang {
3879b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
3889b54502bSHong Zhang   PetscErrorCode ierr;
3899b54502bSHong Zhang 
3909b54502bSHong Zhang   PetscFunctionBegin;
3919b54502bSHong Zhang   ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr);
3929b54502bSHong Zhang   ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr);
3939b54502bSHong Zhang   PetscFunctionReturn(0);
3949b54502bSHong Zhang }
3959b54502bSHong Zhang EXTERN_C_END
3969b54502bSHong Zhang 
3979b54502bSHong Zhang EXTERN_C_BEGIN
3989b54502bSHong Zhang #undef __FUNCT__
3992401956bSBarry Smith #define __FUNCT__ "PCFactorSetPivoting_LU"
4002401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivoting_LU(PC pc,PetscReal dtcol)
4019b54502bSHong Zhang {
4029b54502bSHong Zhang   PC_LU *dir = (PC_LU*)pc->data;
4039b54502bSHong Zhang 
4049b54502bSHong Zhang   PetscFunctionBegin;
405a83599f4SBarry 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);
4069b54502bSHong Zhang   dir->info.dtcol = dtcol;
4079b54502bSHong Zhang   PetscFunctionReturn(0);
4089b54502bSHong Zhang }
4099b54502bSHong Zhang EXTERN_C_END
4109b54502bSHong Zhang 
4119b54502bSHong Zhang EXTERN_C_BEGIN
4129b54502bSHong Zhang #undef __FUNCT__
4132401956bSBarry Smith #define __FUNCT__ "PCFactorSetPivotInBlocks_LU"
4142401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks_LU(PC pc,PetscTruth pivot)
4159b54502bSHong Zhang {
4169b54502bSHong Zhang   PC_LU *dir = (PC_LU*)pc->data;
4179b54502bSHong Zhang 
4189b54502bSHong Zhang   PetscFunctionBegin;
4199b54502bSHong Zhang   dir->info.pivotinblocks = pivot ? 1.0 : 0.0;
4209b54502bSHong Zhang   PetscFunctionReturn(0);
4219b54502bSHong Zhang }
4229b54502bSHong Zhang EXTERN_C_END
4239b54502bSHong Zhang 
4249b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/
4259b54502bSHong Zhang 
4269b54502bSHong Zhang #undef __FUNCT__
4272401956bSBarry Smith #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal"
4289b54502bSHong Zhang /*@
4292401956bSBarry Smith    PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal
4309b54502bSHong Zhang 
4319b54502bSHong Zhang    Collective on PC
4329b54502bSHong Zhang 
4339b54502bSHong Zhang    Input Parameters:
4349b54502bSHong Zhang +  pc - the preconditioner context
4359b54502bSHong Zhang -  tol - diagonal entries smaller than this in absolute value are considered zero
4369b54502bSHong Zhang 
4379b54502bSHong Zhang    Options Database Key:
4382401956bSBarry Smith .  -pc_factor_nonzeros_along_diagonal
4399b54502bSHong Zhang 
4409b54502bSHong Zhang    Level: intermediate
4419b54502bSHong Zhang 
4429b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill
4439b54502bSHong Zhang 
44455ba2a51SBarry Smith .seealso: PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetZeroPivot(), MatReorderForNonzeroDiagonal()
4459b54502bSHong Zhang @*/
4462401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal(PC pc,PetscReal rtol)
4479b54502bSHong Zhang {
4489b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscReal);
4499b54502bSHong Zhang 
4509b54502bSHong Zhang   PetscFunctionBegin;
4519b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
4522401956bSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",(void (**)(void))&f);CHKERRQ(ierr);
4539b54502bSHong Zhang   if (f) {
4549b54502bSHong Zhang     ierr = (*f)(pc,rtol);CHKERRQ(ierr);
4559b54502bSHong Zhang   }
4569b54502bSHong Zhang   PetscFunctionReturn(0);
4579b54502bSHong Zhang }
4589b54502bSHong Zhang 
4599b54502bSHong Zhang #undef __FUNCT__
460c7393fdbSBarry Smith #define __FUNCT__ "PCFactorSetMatSolverPackage"
461c62fe420SBarry Smith /*@
462c7393fdbSBarry Smith    PCFactorSetMatSolverPackage - sets the software that is used to perform the factorization
463c62fe420SBarry Smith 
464c62fe420SBarry Smith    Collective on PC
465c62fe420SBarry Smith 
466c62fe420SBarry Smith    Input Parameters:
467c62fe420SBarry Smith +  pc - the preconditioner context
468*35bd34faSBarry Smith -  stype - for example, spooles, superlu, superlu_dist
469c62fe420SBarry Smith 
470c62fe420SBarry Smith    Options Database Key:
471c7393fdbSBarry Smith .  -pc_factor_mat_solver_package <stype> - spooles, petsc, superlu, superlu_dist, mumps
472c62fe420SBarry Smith 
473c62fe420SBarry Smith    Level: intermediate
474c62fe420SBarry Smith 
475c62fe420SBarry Smith    Note:
476c62fe420SBarry Smith      By default this will use the PETSc factorization if it exists
477c62fe420SBarry Smith 
478*35bd34faSBarry Smith      Use PCFactorGetMatrix(pc,&mat); MatFactorGetPackage(mat,&package); to see what package is being used
479*35bd34faSBarry Smith 
480*35bd34faSBarry Smith 
481c62fe420SBarry Smith .keywords: PC, set, factorization, direct, fill
482c62fe420SBarry Smith 
483c7393fdbSBarry Smith .seealso: MatGetFactor(), MatSolverPackage
484c62fe420SBarry Smith 
485c62fe420SBarry Smith @*/
486c7393fdbSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatSolverPackage(PC pc,const MatSolverPackage stype)
487c62fe420SBarry Smith {
488c7393fdbSBarry Smith   PetscErrorCode ierr,(*f)(PC,const MatSolverPackage);
489c62fe420SBarry Smith 
490c62fe420SBarry Smith   PetscFunctionBegin;
491c62fe420SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
492c7393fdbSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",(void (**)(void))&f);CHKERRQ(ierr);
493c62fe420SBarry Smith   if (f) {
494c62fe420SBarry Smith     ierr = (*f)(pc,stype);CHKERRQ(ierr);
495c62fe420SBarry Smith   }
496c62fe420SBarry Smith   PetscFunctionReturn(0);
497c62fe420SBarry Smith }
498c62fe420SBarry Smith 
499c62fe420SBarry Smith #undef __FUNCT__
50055ba2a51SBarry Smith #define __FUNCT__ "PCFactorSetFill"
5019b54502bSHong Zhang /*@
50255ba2a51SBarry Smith    PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix,
5039b54502bSHong Zhang    fill = number nonzeros in factor/number nonzeros in original matrix.
5049b54502bSHong Zhang 
5059b54502bSHong Zhang    Collective on PC
5069b54502bSHong Zhang 
5079b54502bSHong Zhang    Input Parameters:
5089b54502bSHong Zhang +  pc - the preconditioner context
5099b54502bSHong Zhang -  fill - amount of expected fill
5109b54502bSHong Zhang 
5119b54502bSHong Zhang    Options Database Key:
51255ba2a51SBarry Smith .  -pc_factor_fill <fill> - Sets fill amount
5139b54502bSHong Zhang 
5149b54502bSHong Zhang    Level: intermediate
5159b54502bSHong Zhang 
5169b54502bSHong Zhang    Note:
5179b54502bSHong Zhang    For sparse matrix factorizations it is difficult to predict how much
5186cf91177SBarry Smith    fill to expect. By running with the option -info PETSc will print the
5199b54502bSHong Zhang    actual amount of fill used; allowing you to set the value accurately for
5209b54502bSHong Zhang    future runs. Default PETSc uses a value of 5.0
5219b54502bSHong Zhang 
5229b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill
5239b54502bSHong Zhang 
5249b54502bSHong Zhang @*/
52555ba2a51SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill(PC pc,PetscReal fill)
5269b54502bSHong Zhang {
5279b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscReal);
5289b54502bSHong Zhang 
5299b54502bSHong Zhang   PetscFunctionBegin;
5309b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
5319b54502bSHong Zhang   if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0");
53255ba2a51SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetFill_C",(void (**)(void))&f);CHKERRQ(ierr);
5339b54502bSHong Zhang   if (f) {
5349b54502bSHong Zhang     ierr = (*f)(pc,fill);CHKERRQ(ierr);
5359b54502bSHong Zhang   }
5369b54502bSHong Zhang   PetscFunctionReturn(0);
5379b54502bSHong Zhang }
5389b54502bSHong Zhang 
5399b54502bSHong Zhang #undef __FUNCT__
5402401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace"
5419b54502bSHong Zhang /*@
5422401956bSBarry Smith    PCFactorSetUseInPlace - Tells the system to do an in-place factorization.
5439b54502bSHong Zhang    For dense matrices, this enables the solution of much larger problems.
5449b54502bSHong Zhang    For sparse matrices the factorization cannot be done truly in-place
5459b54502bSHong Zhang    so this does not save memory during the factorization, but after the matrix
5469b54502bSHong Zhang    is factored, the original unfactored matrix is freed, thus recovering that
5479b54502bSHong Zhang    space.
5489b54502bSHong Zhang 
5499b54502bSHong Zhang    Collective on PC
5509b54502bSHong Zhang 
5519b54502bSHong Zhang    Input Parameters:
5529b54502bSHong Zhang .  pc - the preconditioner context
5539b54502bSHong Zhang 
5549b54502bSHong Zhang    Options Database Key:
5552401956bSBarry Smith .  -pc_factor_in_place - Activates in-place factorization
5569b54502bSHong Zhang 
5579b54502bSHong Zhang    Notes:
5582401956bSBarry Smith    PCFactorSetUseInplace() can only be used with the KSP method KSPPREONLY or when
5599b54502bSHong Zhang    a different matrix is provided for the multiply and the preconditioner in
5609b54502bSHong Zhang    a call to KSPSetOperators().
5619b54502bSHong Zhang    This is because the Krylov space methods require an application of the
5629b54502bSHong Zhang    matrix multiplication, which is not possible here because the matrix has
5639b54502bSHong Zhang    been factored in-place, replacing the original matrix.
5649b54502bSHong Zhang 
5659b54502bSHong Zhang    Level: intermediate
5669b54502bSHong Zhang 
5679b54502bSHong Zhang .keywords: PC, set, factorization, direct, inplace, in-place, LU
5689b54502bSHong Zhang 
5699b54502bSHong Zhang .seealso: PCILUSetUseInPlace()
5709b54502bSHong Zhang @*/
5712401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace(PC pc)
5729b54502bSHong Zhang {
5739b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC);
5749b54502bSHong Zhang 
5759b54502bSHong Zhang   PetscFunctionBegin;
5769b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
5772401956bSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",(void (**)(void))&f);CHKERRQ(ierr);
5789b54502bSHong Zhang   if (f) {
5799b54502bSHong Zhang     ierr = (*f)(pc);CHKERRQ(ierr);
5809b54502bSHong Zhang   }
5819b54502bSHong Zhang   PetscFunctionReturn(0);
5829b54502bSHong Zhang }
5839b54502bSHong Zhang 
5849b54502bSHong Zhang #undef __FUNCT__
585e5a9bf91SBarry Smith #define __FUNCT__ "PCFactorSetMatOrderingType"
5869b54502bSHong Zhang /*@C
587e5a9bf91SBarry Smith     PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to
5889b54502bSHong Zhang     be used in the LU factorization.
5899b54502bSHong Zhang 
5909b54502bSHong Zhang     Collective on PC
5919b54502bSHong Zhang 
5929b54502bSHong Zhang     Input Parameters:
5939b54502bSHong Zhang +   pc - the preconditioner context
5949b54502bSHong Zhang -   ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM
5959b54502bSHong Zhang 
5969b54502bSHong Zhang     Options Database Key:
5972401956bSBarry Smith .   -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
5989b54502bSHong Zhang 
5999b54502bSHong Zhang     Level: intermediate
6009b54502bSHong Zhang 
6019b54502bSHong Zhang     Notes: nested dissection is used by default
6029b54502bSHong Zhang 
6032b6de112SBarry Smith     For Cholesky and ICC and the SBAIJ format reorderings are not available,
6042b6de112SBarry Smith     since only the upper triangular part of the matrix is stored. You can use the
6052b6de112SBarry Smith     SeqAIJ format in this case to get reorderings.
6062b6de112SBarry Smith 
6079b54502bSHong Zhang @*/
608e36ff184SSatish Balay PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType(PC pc,const MatOrderingType ordering)
6099b54502bSHong Zhang {
610e36ff184SSatish Balay   PetscErrorCode ierr,(*f)(PC,const MatOrderingType);
6119b54502bSHong Zhang 
6129b54502bSHong Zhang   PetscFunctionBegin;
613e5a9bf91SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",(void (**)(void))&f);CHKERRQ(ierr);
6149b54502bSHong Zhang   if (f) {
6159b54502bSHong Zhang     ierr = (*f)(pc,ordering);CHKERRQ(ierr);
6169b54502bSHong Zhang   }
6179b54502bSHong Zhang   PetscFunctionReturn(0);
6189b54502bSHong Zhang }
6199b54502bSHong Zhang 
6209b54502bSHong Zhang #undef __FUNCT__
6212401956bSBarry Smith #define __FUNCT__ "PCFactorSetPivoting"
6229b54502bSHong Zhang /*@
6232401956bSBarry Smith     PCFactorSetPivoting - Determines when pivoting is done during LU.
6249b54502bSHong Zhang       For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices
6259b54502bSHong Zhang       it is never done. For the Matlab and SuperLU factorization this is used.
6269b54502bSHong Zhang 
6279b54502bSHong Zhang     Collective on PC
6289b54502bSHong Zhang 
6299b54502bSHong Zhang     Input Parameters:
6309b54502bSHong Zhang +   pc - the preconditioner context
6319b54502bSHong Zhang -   dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable)
6329b54502bSHong Zhang 
6339b54502bSHong Zhang     Options Database Key:
6342401956bSBarry Smith .   -pc_factor_pivoting <dtcol>
6359b54502bSHong Zhang 
6369b54502bSHong Zhang     Level: intermediate
6379b54502bSHong Zhang 
6382401956bSBarry Smith .seealso: PCILUSetMatOrdering(), PCFactorSetPivotInBlocks()
6399b54502bSHong Zhang @*/
6402401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivoting(PC pc,PetscReal dtcol)
6419b54502bSHong Zhang {
6429b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscReal);
6439b54502bSHong Zhang 
6449b54502bSHong Zhang   PetscFunctionBegin;
6452401956bSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivoting_C",(void (**)(void))&f);CHKERRQ(ierr);
6469b54502bSHong Zhang   if (f) {
6479b54502bSHong Zhang     ierr = (*f)(pc,dtcol);CHKERRQ(ierr);
6489b54502bSHong Zhang   }
6499b54502bSHong Zhang   PetscFunctionReturn(0);
6509b54502bSHong Zhang }
6519b54502bSHong Zhang 
6529b54502bSHong Zhang #undef __FUNCT__
6532401956bSBarry Smith #define __FUNCT__ "PCFactorSetPivotInBlocks"
6549b54502bSHong Zhang /*@
6552401956bSBarry Smith     PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block
6569b54502bSHong Zhang       with BAIJ or SBAIJ matrices
6579b54502bSHong Zhang 
6589b54502bSHong Zhang     Collective on PC
6599b54502bSHong Zhang 
6609b54502bSHong Zhang     Input Parameters:
6619b54502bSHong Zhang +   pc - the preconditioner context
6629b54502bSHong Zhang -   pivot - PETSC_TRUE or PETSC_FALSE
6639b54502bSHong Zhang 
6649b54502bSHong Zhang     Options Database Key:
6652401956bSBarry Smith .   -pc_factor_pivot_in_blocks <true,false>
6669b54502bSHong Zhang 
6679b54502bSHong Zhang     Level: intermediate
6689b54502bSHong Zhang 
6692401956bSBarry Smith .seealso: PCILUSetMatOrdering(), PCFactorSetPivoting()
6709b54502bSHong Zhang @*/
6712401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks(PC pc,PetscTruth pivot)
6729b54502bSHong Zhang {
6739b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscTruth);
6749b54502bSHong Zhang 
6759b54502bSHong Zhang   PetscFunctionBegin;
6762401956bSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr);
6779b54502bSHong Zhang   if (f) {
6789b54502bSHong Zhang     ierr = (*f)(pc,pivot);CHKERRQ(ierr);
6799b54502bSHong Zhang   }
6809b54502bSHong Zhang   PetscFunctionReturn(0);
6819b54502bSHong Zhang }
6829b54502bSHong Zhang 
6839b54502bSHong Zhang /* ------------------------------------------------------------------------ */
6849b54502bSHong Zhang 
6859b54502bSHong Zhang /*MC
6869b54502bSHong Zhang    PCLU - Uses a direct solver, based on LU factorization, as a preconditioner
6879b54502bSHong Zhang 
6889b54502bSHong Zhang    Options Database Keys:
6892401956bSBarry Smith +  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
690c7393fdbSBarry Smith .  -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like spooles
6912401956bSBarry Smith .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
69255ba2a51SBarry Smith .  -pc_factor_fill <fill> - Sets fill amount
6932401956bSBarry Smith .  -pc_factor_in_place - Activates in-place factorization
6942401956bSBarry Smith .  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
6952401956bSBarry Smith .  -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase
6969b54502bSHong Zhang                                          stability of factorization.
697f251bdbdSHong Zhang .  -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default
698e22d95b2SBarry Smith .  -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value
699f251bdbdSHong Zhang         is optional with PETSC_TRUE being the default
700e22d95b2SBarry Smith -   -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the
701e22d95b2SBarry Smith         diagonal.
7029b54502bSHong Zhang 
7039b54502bSHong Zhang    Notes: Not all options work for all matrix formats
7049b54502bSHong Zhang           Run with -help to see additional options for particular matrix formats or factorization
7059b54502bSHong Zhang           algorithms
7069b54502bSHong Zhang 
7079b54502bSHong Zhang    Level: beginner
7089b54502bSHong Zhang 
7099b54502bSHong Zhang    Concepts: LU factorization, direct solver
7109b54502bSHong Zhang 
7119b54502bSHong Zhang    Notes: Usually this will compute an "exact" solution in one iteration and does
7129b54502bSHong Zhang           not need a Krylov method (i.e. you can use -ksp_type preonly, or
7139b54502bSHong Zhang           KSPSetType(ksp,KSPPREONLY) for the Krylov method
7149b54502bSHong Zhang 
7159b54502bSHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
716a4fd02acSBarry Smith            PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(),
717e5a9bf91SBarry Smith            PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetPivoting(),
718e22d95b2SBarry Smith            PCFactorSetPivotingInBlocks(),PCFactorSetShiftNonzero(),PCFactorSetShiftPd(), PCFactorReorderForNonzeroDiagonal()
7199b54502bSHong Zhang M*/
7209b54502bSHong Zhang 
7219b54502bSHong Zhang EXTERN_C_BEGIN
7229b54502bSHong Zhang #undef __FUNCT__
7239b54502bSHong Zhang #define __FUNCT__ "PCCreate_LU"
724dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_LU(PC pc)
7259b54502bSHong Zhang {
7269b54502bSHong Zhang   PetscErrorCode ierr;
7279b54502bSHong Zhang   PetscMPIInt    size;
7289b54502bSHong Zhang   PC_LU          *dir;
7299b54502bSHong Zhang 
7309b54502bSHong Zhang   PetscFunctionBegin;
73138f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_LU,&dir);CHKERRQ(ierr);
7329b54502bSHong Zhang 
7339b54502bSHong Zhang   ierr = MatFactorInfoInitialize(&dir->info);CHKERRQ(ierr);
7349b54502bSHong Zhang   dir->fact                  = 0;
7359b54502bSHong Zhang   dir->inplace               = PETSC_FALSE;
7369b54502bSHong Zhang   dir->nonzerosalongdiagonal = PETSC_FALSE;
7379b54502bSHong Zhang 
7389b54502bSHong Zhang   dir->info.fill           = 5.0;
7399b54502bSHong Zhang   dir->info.dtcol          = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */
7400a29876aSHong Zhang   dir->info.shiftnz        = 0.0;
7419b54502bSHong Zhang   dir->info.zeropivot      = 1.e-12;
7429b54502bSHong Zhang   dir->info.pivotinblocks  = 1.0;
743fbf22428SSatish Balay   dir->info.shiftpd        = 0.0; /* false */
7449b54502bSHong Zhang   dir->col                 = 0;
7459b54502bSHong Zhang   dir->row                 = 0;
7465c9eb25fSBarry Smith 
74743244d56SBarry Smith   ierr = PetscStrallocpy(MAT_SOLVER_PETSC,&dir->solvertype);CHKERRQ(ierr);
7487adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
7499b54502bSHong Zhang   if (size == 1) {
7509b54502bSHong Zhang     ierr = PetscStrallocpy(MATORDERING_ND,&dir->ordering);CHKERRQ(ierr);
7519b54502bSHong Zhang   } else {
7529b54502bSHong Zhang     ierr = PetscStrallocpy(MATORDERING_NATURAL,&dir->ordering);CHKERRQ(ierr);
7539b54502bSHong Zhang   }
7549b54502bSHong Zhang   dir->reusefill        = PETSC_FALSE;
7559b54502bSHong Zhang   dir->reuseordering    = PETSC_FALSE;
7569b54502bSHong Zhang   pc->data              = (void*)dir;
7579b54502bSHong Zhang 
7589b54502bSHong Zhang   pc->ops->destroy           = PCDestroy_LU;
7599b54502bSHong Zhang   pc->ops->apply             = PCApply_LU;
7609b54502bSHong Zhang   pc->ops->applytranspose    = PCApplyTranspose_LU;
7619b54502bSHong Zhang   pc->ops->setup             = PCSetUp_LU;
7629b54502bSHong Zhang   pc->ops->setfromoptions    = PCSetFromOptions_LU;
7639b54502bSHong Zhang   pc->ops->view              = PCView_LU;
7649b54502bSHong Zhang   pc->ops->applyrichardson   = 0;
765a4fd02acSBarry Smith   pc->ops->getfactoredmatrix = PCFactorGetMatrix_LU;
7669b54502bSHong Zhang 
767c7393fdbSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_LU",
768c7393fdbSBarry Smith                     PCFactorSetMatSolverPackage_LU);CHKERRQ(ierr);
769afaefe49SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_LU",
770afaefe49SHong Zhang                     PCFactorSetZeroPivot_LU);CHKERRQ(ierr);
771afaefe49SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_LU",
772afaefe49SHong Zhang                     PCFactorSetShiftNonzero_LU);CHKERRQ(ierr);
773afaefe49SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_LU",
774afaefe49SHong Zhang                     PCFactorSetShiftPd_LU);CHKERRQ(ierr);
775afaefe49SHong Zhang 
77655ba2a51SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_LU",
77755ba2a51SBarry Smith                     PCFactorSetFill_LU);CHKERRQ(ierr);
7782401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_LU",
7792401956bSBarry Smith                     PCFactorSetUseInPlace_LU);CHKERRQ(ierr);
780e5a9bf91SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_LU",
781e5a9bf91SBarry Smith                     PCFactorSetMatOrderingType_LU);CHKERRQ(ierr);
7822401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_LU",
7832401956bSBarry Smith                     PCFactorSetReuseOrdering_LU);CHKERRQ(ierr);
7842401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_LU",
7852401956bSBarry Smith                     PCFactorSetReuseFill_LU);CHKERRQ(ierr);
7862401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivoting_C","PCFactorSetPivoting_LU",
7872401956bSBarry Smith                     PCFactorSetPivoting_LU);CHKERRQ(ierr);
7882401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_LU",
7892401956bSBarry Smith                     PCFactorSetPivotInBlocks_LU);CHKERRQ(ierr);
7902401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_LU",
7912401956bSBarry Smith                     PCFactorReorderForNonzeroDiagonal_LU);CHKERRQ(ierr);
7929b54502bSHong Zhang   PetscFunctionReturn(0);
7939b54502bSHong Zhang }
7949b54502bSHong Zhang EXTERN_C_END
795