xref: /petsc/src/ksp/pc/impls/factor/lu/lu.c (revision ee45ca4afdde1af4d1deda7cd4dc1a4a63a3ea97)
19b54502bSHong Zhang /*
29b54502bSHong Zhang    Defines a direct factorization preconditioner for any Mat implementation
39b54502bSHong Zhang    Note: this need not be consided a preconditioner since it supplies
49b54502bSHong Zhang          a direct solver.
59b54502bSHong Zhang */
6*ee45ca4aSHong Zhang 
79b54502bSHong Zhang #include "src/ksp/pc/pcimpl.h"                /*I "petscpc.h" I*/
8*ee45ca4aSHong Zhang #include "src/ksp/pc/impls/factor/lu/lu.h"
99b54502bSHong Zhang 
109b54502bSHong Zhang EXTERN_C_BEGIN
119b54502bSHong Zhang #undef __FUNCT__
129b54502bSHong Zhang #define __FUNCT__ "PCLUReorderForNonzeroDiagonal_LU"
139b54502bSHong Zhang PetscErrorCode PCLUReorderForNonzeroDiagonal_LU(PC pc,PetscReal z)
149b54502bSHong Zhang {
159b54502bSHong Zhang   PC_LU *lu = (PC_LU*)pc->data;
169b54502bSHong Zhang 
179b54502bSHong Zhang   PetscFunctionBegin;
189b54502bSHong Zhang   lu->nonzerosalongdiagonal = PETSC_TRUE;
199b54502bSHong Zhang   if (z == PETSC_DECIDE) {
209b54502bSHong Zhang     lu->nonzerosalongdiagonaltol = 1.e-10;
219b54502bSHong Zhang   } else {
229b54502bSHong Zhang     lu->nonzerosalongdiagonaltol = z;
239b54502bSHong Zhang   }
249b54502bSHong Zhang   PetscFunctionReturn(0);
259b54502bSHong Zhang }
269b54502bSHong Zhang EXTERN_C_END
279b54502bSHong Zhang 
289b54502bSHong Zhang EXTERN_C_BEGIN
299b54502bSHong Zhang #undef __FUNCT__
309b54502bSHong Zhang #define __FUNCT__ "PCLUSetReuseOrdering_LU"
319b54502bSHong Zhang PetscErrorCode PCLUSetReuseOrdering_LU(PC pc,PetscTruth flag)
329b54502bSHong Zhang {
339b54502bSHong Zhang   PC_LU *lu;
349b54502bSHong Zhang 
359b54502bSHong Zhang   PetscFunctionBegin;
369b54502bSHong Zhang   lu                = (PC_LU*)pc->data;
379b54502bSHong Zhang   lu->reuseordering = flag;
389b54502bSHong Zhang   PetscFunctionReturn(0);
399b54502bSHong Zhang }
409b54502bSHong Zhang EXTERN_C_END
419b54502bSHong Zhang 
429b54502bSHong Zhang EXTERN_C_BEGIN
439b54502bSHong Zhang #undef __FUNCT__
449b54502bSHong Zhang #define __FUNCT__ "PCLUSetReuseFill_LU"
459b54502bSHong Zhang PetscErrorCode PCLUSetReuseFill_LU(PC pc,PetscTruth flag)
469b54502bSHong Zhang {
479b54502bSHong Zhang   PC_LU *lu;
489b54502bSHong Zhang 
499b54502bSHong Zhang   PetscFunctionBegin;
509b54502bSHong Zhang   lu = (PC_LU*)pc->data;
519b54502bSHong Zhang   lu->reusefill = flag;
529b54502bSHong Zhang   PetscFunctionReturn(0);
539b54502bSHong Zhang }
549b54502bSHong Zhang EXTERN_C_END
559b54502bSHong Zhang 
569b54502bSHong Zhang #undef __FUNCT__
579b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_LU"
589b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_LU(PC pc)
599b54502bSHong Zhang {
609b54502bSHong Zhang   PC_LU          *lu = (PC_LU*)pc->data;
619b54502bSHong Zhang   PetscErrorCode ierr;
629b54502bSHong Zhang   PetscTruth     flg,set;
639b54502bSHong Zhang   char           tname[256];
649b54502bSHong Zhang   PetscFList     ordlist;
659b54502bSHong Zhang   PetscReal      tol;
669b54502bSHong Zhang 
679b54502bSHong Zhang   PetscFunctionBegin;
689b54502bSHong Zhang   ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr);
699b54502bSHong Zhang   ierr = PetscOptionsHead("LU options");CHKERRQ(ierr);
709b54502bSHong Zhang     ierr = PetscOptionsName("-pc_lu_in_place","Form LU in the same memory as the matrix","PCLUSetUseInPlace",&flg);CHKERRQ(ierr);
719b54502bSHong Zhang     if (flg) {
729b54502bSHong Zhang       ierr = PCLUSetUseInPlace(pc);CHKERRQ(ierr);
739b54502bSHong Zhang     }
749b54502bSHong Zhang     ierr = PetscOptionsReal("-pc_lu_fill","Expected non-zeros in LU/non-zeros in matrix","PCLUSetFill",lu->info.fill,&lu->info.fill,0);CHKERRQ(ierr);
759b54502bSHong Zhang 
76*ee45ca4aSHong Zhang     ierr = PetscOptionsName("-pc_factor_shiftnonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr);
779b54502bSHong Zhang     if (flg) {
78*ee45ca4aSHong Zhang         ierr = PCFactorSetShiftNonzero((PetscReal) PETSC_DECIDE,&lu->info);CHKERRQ(ierr);
799b54502bSHong Zhang     }
800a29876aSHong Zhang     ierr = PetscOptionsReal("-pc_factor_shiftnonzero","Shift added to diagonal","PCFactorSetShiftNonzero",lu->info.shiftnz,&lu->info.shiftnz,0);CHKERRQ(ierr);
81*ee45ca4aSHong Zhang     ierr = PetscOptionsName("-pc_factor_shiftpd","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",&flg);CHKERRQ(ierr);
829b54502bSHong Zhang     if (flg) {
83*ee45ca4aSHong Zhang       ierr = PCFactorSetShiftPd(PETSC_TRUE,&lu->info);CHKERRQ(ierr);
849b54502bSHong Zhang     }
85*ee45ca4aSHong Zhang     ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",lu->info.zeropivot,&lu->info.zeropivot,0);CHKERRQ(ierr);
869b54502bSHong Zhang 
879b54502bSHong Zhang     ierr = PetscOptionsName("-pc_lu_reuse_fill","Use fill from previous factorization","PCLUSetReuseFill",&flg);CHKERRQ(ierr);
889b54502bSHong Zhang     if (flg) {
899b54502bSHong Zhang       ierr = PCLUSetReuseFill(pc,PETSC_TRUE);CHKERRQ(ierr);
909b54502bSHong Zhang     }
919b54502bSHong Zhang     ierr = PetscOptionsName("-pc_lu_reuse_ordering","Reuse ordering from previous factorization","PCLUSetReuseOrdering",&flg);CHKERRQ(ierr);
929b54502bSHong Zhang     if (flg) {
939b54502bSHong Zhang       ierr = PCLUSetReuseOrdering(pc,PETSC_TRUE);CHKERRQ(ierr);
949b54502bSHong Zhang     }
959b54502bSHong Zhang 
969b54502bSHong Zhang     ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr);
979b54502bSHong Zhang     ierr = PetscOptionsList("-pc_lu_mat_ordering_type","Reordering to reduce nonzeros in LU","PCLUSetMatOrdering",ordlist,lu->ordering,tname,256,&flg);CHKERRQ(ierr);
989b54502bSHong Zhang     if (flg) {
999b54502bSHong Zhang       ierr = PCLUSetMatOrdering(pc,tname);CHKERRQ(ierr);
1009b54502bSHong Zhang     }
1019b54502bSHong Zhang 
1029b54502bSHong Zhang     ierr = PetscOptionsName("-pc_lu_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCLUReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr);
1039b54502bSHong Zhang     if (flg) {
1049b54502bSHong Zhang       tol = PETSC_DECIDE;
1059b54502bSHong Zhang       ierr = PetscOptionsReal("-pc_lu_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCLUReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr);
1069b54502bSHong Zhang       ierr = PCLUReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr);
1079b54502bSHong Zhang     }
1089b54502bSHong Zhang 
1099b54502bSHong Zhang     ierr = PetscOptionsReal("-pc_lu_pivoting","Pivoting tolerance (used only for some factorization)","PCLUSetPivoting",lu->info.dtcol,&lu->info.dtcol,&flg);CHKERRQ(ierr);
1109b54502bSHong Zhang 
1119b54502bSHong Zhang     flg = lu->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE;
1129b54502bSHong Zhang     ierr = PetscOptionsLogical("-pc_lu_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCLUSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr);
1139b54502bSHong Zhang     if (set) {
1149b54502bSHong Zhang       ierr = PCLUSetPivotInBlocks(pc,flg);CHKERRQ(ierr);
1159b54502bSHong Zhang     }
1169b54502bSHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
1179b54502bSHong Zhang   PetscFunctionReturn(0);
1189b54502bSHong Zhang }
1199b54502bSHong Zhang 
1209b54502bSHong Zhang #undef __FUNCT__
1219b54502bSHong Zhang #define __FUNCT__ "PCView_LU"
1229b54502bSHong Zhang static PetscErrorCode PCView_LU(PC pc,PetscViewer viewer)
1239b54502bSHong Zhang {
1249b54502bSHong Zhang   PC_LU          *lu = (PC_LU*)pc->data;
1259b54502bSHong Zhang   PetscErrorCode ierr;
1269b54502bSHong Zhang   PetscTruth     iascii,isstring;
1279b54502bSHong Zhang 
1289b54502bSHong Zhang   PetscFunctionBegin;
1299b54502bSHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
1309b54502bSHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
1319b54502bSHong Zhang   if (iascii) {
1329b54502bSHong Zhang     MatInfo info;
1339b54502bSHong Zhang 
1349b54502bSHong Zhang     if (lu->inplace) {ierr = PetscViewerASCIIPrintf(viewer,"  LU: in-place factorization\n");CHKERRQ(ierr);}
1359b54502bSHong Zhang     else             {ierr = PetscViewerASCIIPrintf(viewer,"  LU: out-of-place factorization\n");CHKERRQ(ierr);}
1369b54502bSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"    matrix ordering: %s\n",lu->ordering);CHKERRQ(ierr);
1379b54502bSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  LU: tolerance for zero pivot %g\n",lu->info.zeropivot);CHKERRQ(ierr);
1380a29876aSHong Zhang     if (lu->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer,"  LU: using Manteuffel shift\n");CHKERRQ(ierr);}
1399b54502bSHong Zhang     if (lu->fact) {
1409b54502bSHong Zhang       ierr = MatGetInfo(lu->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
1419b54502bSHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"    LU nonzeros %g\n",info.nz_used);CHKERRQ(ierr);
1429b54502bSHong Zhang       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_FACTOR_INFO);CHKERRQ(ierr);
1439b54502bSHong Zhang       ierr = MatView(lu->fact,viewer);CHKERRQ(ierr);
1449b54502bSHong Zhang       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
1459b54502bSHong Zhang     }
1469b54502bSHong Zhang     if (lu->reusefill)    {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing fill from past factorization\n");CHKERRQ(ierr);}
1479b54502bSHong Zhang     if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing reordering from past factorization\n");CHKERRQ(ierr);}
1489b54502bSHong Zhang   } else if (isstring) {
1499b54502bSHong Zhang     ierr = PetscViewerStringSPrintf(viewer," order=%s",lu->ordering);CHKERRQ(ierr);CHKERRQ(ierr);
1509b54502bSHong Zhang   } else {
1519b54502bSHong Zhang     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCLU",((PetscObject)viewer)->type_name);
1529b54502bSHong Zhang   }
1539b54502bSHong Zhang   PetscFunctionReturn(0);
1549b54502bSHong Zhang }
1559b54502bSHong Zhang 
1569b54502bSHong Zhang #undef __FUNCT__
1579b54502bSHong Zhang #define __FUNCT__ "PCGetFactoredMatrix_LU"
1589b54502bSHong Zhang static PetscErrorCode PCGetFactoredMatrix_LU(PC pc,Mat *mat)
1599b54502bSHong Zhang {
1609b54502bSHong Zhang   PC_LU *dir = (PC_LU*)pc->data;
1619b54502bSHong Zhang 
1629b54502bSHong Zhang   PetscFunctionBegin;
1639b54502bSHong Zhang   if (!dir->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()");
1649b54502bSHong Zhang   *mat = dir->fact;
1659b54502bSHong Zhang   PetscFunctionReturn(0);
1669b54502bSHong Zhang }
1679b54502bSHong Zhang 
1689b54502bSHong Zhang #undef __FUNCT__
1699b54502bSHong Zhang #define __FUNCT__ "PCSetUp_LU"
1709b54502bSHong Zhang static PetscErrorCode PCSetUp_LU(PC pc)
1719b54502bSHong Zhang {
1729b54502bSHong Zhang   PetscErrorCode ierr;
1739b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
1749b54502bSHong Zhang 
1759b54502bSHong Zhang   PetscFunctionBegin;
1769b54502bSHong Zhang   if (dir->reusefill && pc->setupcalled) dir->info.fill = dir->actualfill;
1779b54502bSHong Zhang 
1789b54502bSHong Zhang   if (dir->inplace) {
1799b54502bSHong Zhang     if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);}
1809b54502bSHong Zhang     if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);}
1819b54502bSHong Zhang     ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
1829b54502bSHong Zhang     if (dir->row) {PetscLogObjectParent(pc,dir->row); PetscLogObjectParent(pc,dir->col);}
1839b54502bSHong Zhang     ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&dir->info);CHKERRQ(ierr);
1849b54502bSHong Zhang     dir->fact = pc->pmat;
1859b54502bSHong Zhang   } else {
1869b54502bSHong Zhang     MatInfo info;
1879b54502bSHong Zhang     if (!pc->setupcalled) {
1889b54502bSHong Zhang       ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
1899b54502bSHong Zhang       if (dir->nonzerosalongdiagonal) {
1909b54502bSHong Zhang         ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr);
1919b54502bSHong Zhang       }
1929b54502bSHong Zhang       if (dir->row) {PetscLogObjectParent(pc,dir->row); PetscLogObjectParent(pc,dir->col);}
1939b54502bSHong Zhang       ierr = MatLUFactorSymbolic(pc->pmat,dir->row,dir->col,&dir->info,&dir->fact);CHKERRQ(ierr);
1949b54502bSHong Zhang       ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
1959b54502bSHong Zhang       dir->actualfill = info.fill_ratio_needed;
1969b54502bSHong Zhang       PetscLogObjectParent(pc,dir->fact);
1979b54502bSHong Zhang     } else if (pc->flag != SAME_NONZERO_PATTERN) {
1989b54502bSHong Zhang       if (!dir->reuseordering) {
1999b54502bSHong Zhang         if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);}
2009b54502bSHong Zhang         if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);}
2019b54502bSHong Zhang         ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
2029b54502bSHong Zhang         if (dir->nonzerosalongdiagonal) {
2039b54502bSHong Zhang          ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr);
2049b54502bSHong Zhang         }
2059b54502bSHong Zhang         if (dir->row) {PetscLogObjectParent(pc,dir->row); PetscLogObjectParent(pc,dir->col);}
2069b54502bSHong Zhang       }
2079b54502bSHong Zhang       ierr = MatDestroy(dir->fact);CHKERRQ(ierr);
2089b54502bSHong Zhang       ierr = MatLUFactorSymbolic(pc->pmat,dir->row,dir->col,&dir->info,&dir->fact);CHKERRQ(ierr);
2099b54502bSHong Zhang       ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
2109b54502bSHong Zhang       dir->actualfill = info.fill_ratio_needed;
2119b54502bSHong Zhang       PetscLogObjectParent(pc,dir->fact);
2129b54502bSHong Zhang     }
2139b54502bSHong Zhang     ierr = MatLUFactorNumeric(pc->pmat,&dir->info,&dir->fact);CHKERRQ(ierr);
2149b54502bSHong Zhang   }
2159b54502bSHong Zhang   PetscFunctionReturn(0);
2169b54502bSHong Zhang }
2179b54502bSHong Zhang 
2189b54502bSHong Zhang #undef __FUNCT__
2199b54502bSHong Zhang #define __FUNCT__ "PCDestroy_LU"
2209b54502bSHong Zhang static PetscErrorCode PCDestroy_LU(PC pc)
2219b54502bSHong Zhang {
2229b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
2239b54502bSHong Zhang   PetscErrorCode ierr;
2249b54502bSHong Zhang 
2259b54502bSHong Zhang   PetscFunctionBegin;
2269b54502bSHong Zhang   if (!dir->inplace && dir->fact) {ierr = MatDestroy(dir->fact);CHKERRQ(ierr);}
2279b54502bSHong Zhang   if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);}
2289b54502bSHong Zhang   if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);}
2299b54502bSHong Zhang   ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr);
2309b54502bSHong Zhang   ierr = PetscFree(dir);CHKERRQ(ierr);
2319b54502bSHong Zhang   PetscFunctionReturn(0);
2329b54502bSHong Zhang }
2339b54502bSHong Zhang 
2349b54502bSHong Zhang #undef __FUNCT__
2359b54502bSHong Zhang #define __FUNCT__ "PCApply_LU"
2369b54502bSHong Zhang static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y)
2379b54502bSHong Zhang {
2389b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
2399b54502bSHong Zhang   PetscErrorCode ierr;
2409b54502bSHong Zhang 
2419b54502bSHong Zhang   PetscFunctionBegin;
2429b54502bSHong Zhang   if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);}
2439b54502bSHong Zhang   else              {ierr = MatSolve(dir->fact,x,y);CHKERRQ(ierr);}
2449b54502bSHong Zhang   PetscFunctionReturn(0);
2459b54502bSHong Zhang }
2469b54502bSHong Zhang 
2479b54502bSHong Zhang #undef __FUNCT__
2489b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_LU"
2499b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y)
2509b54502bSHong Zhang {
2519b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
2529b54502bSHong Zhang   PetscErrorCode ierr;
2539b54502bSHong Zhang 
2549b54502bSHong Zhang   PetscFunctionBegin;
2559b54502bSHong Zhang   if (dir->inplace) {ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);}
2569b54502bSHong Zhang   else              {ierr = MatSolveTranspose(dir->fact,x,y);CHKERRQ(ierr);}
2579b54502bSHong Zhang   PetscFunctionReturn(0);
2589b54502bSHong Zhang }
2599b54502bSHong Zhang 
2609b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/
2619b54502bSHong Zhang 
2629b54502bSHong Zhang EXTERN_C_BEGIN
2639b54502bSHong Zhang #undef __FUNCT__
2649b54502bSHong Zhang #define __FUNCT__ "PCLUSetFill_LU"
2659b54502bSHong Zhang PetscErrorCode PCLUSetFill_LU(PC pc,PetscReal fill)
2669b54502bSHong Zhang {
2679b54502bSHong Zhang   PC_LU *dir;
2689b54502bSHong Zhang 
2699b54502bSHong Zhang   PetscFunctionBegin;
2709b54502bSHong Zhang   dir = (PC_LU*)pc->data;
2719b54502bSHong Zhang   dir->info.fill = fill;
2729b54502bSHong Zhang   PetscFunctionReturn(0);
2739b54502bSHong Zhang }
2749b54502bSHong Zhang EXTERN_C_END
2759b54502bSHong Zhang 
2769b54502bSHong Zhang EXTERN_C_BEGIN
2779b54502bSHong Zhang #undef __FUNCT__
2789b54502bSHong Zhang #define __FUNCT__ "PCLUSetUseInPlace_LU"
2799b54502bSHong Zhang PetscErrorCode PCLUSetUseInPlace_LU(PC pc)
2809b54502bSHong Zhang {
2819b54502bSHong Zhang   PC_LU *dir;
2829b54502bSHong Zhang 
2839b54502bSHong Zhang   PetscFunctionBegin;
2849b54502bSHong Zhang   dir = (PC_LU*)pc->data;
2859b54502bSHong Zhang   dir->inplace = PETSC_TRUE;
2869b54502bSHong Zhang   PetscFunctionReturn(0);
2879b54502bSHong Zhang }
2889b54502bSHong Zhang EXTERN_C_END
2899b54502bSHong Zhang 
2909b54502bSHong Zhang EXTERN_C_BEGIN
2919b54502bSHong Zhang #undef __FUNCT__
2929b54502bSHong Zhang #define __FUNCT__ "PCLUSetMatOrdering_LU"
2939b54502bSHong Zhang PetscErrorCode PCLUSetMatOrdering_LU(PC pc,MatOrderingType ordering)
2949b54502bSHong Zhang {
2959b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
2969b54502bSHong Zhang   PetscErrorCode ierr;
2979b54502bSHong Zhang 
2989b54502bSHong Zhang   PetscFunctionBegin;
2999b54502bSHong Zhang   ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr);
3009b54502bSHong Zhang   ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr);
3019b54502bSHong Zhang   PetscFunctionReturn(0);
3029b54502bSHong Zhang }
3039b54502bSHong Zhang EXTERN_C_END
3049b54502bSHong Zhang 
3059b54502bSHong Zhang EXTERN_C_BEGIN
3069b54502bSHong Zhang #undef __FUNCT__
3079b54502bSHong Zhang #define __FUNCT__ "PCLUSetPivoting_LU"
3089b54502bSHong Zhang PetscErrorCode PCLUSetPivoting_LU(PC pc,PetscReal dtcol)
3099b54502bSHong Zhang {
3109b54502bSHong Zhang   PC_LU *dir = (PC_LU*)pc->data;
3119b54502bSHong Zhang 
3129b54502bSHong Zhang   PetscFunctionBegin;
3139b54502bSHong Zhang   if (dtcol < 0.0 || dtcol > 1.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Column pivot tolerance is %g must be between 0 and 1",dtcol);
3149b54502bSHong Zhang   dir->info.dtcol = dtcol;
3159b54502bSHong Zhang   PetscFunctionReturn(0);
3169b54502bSHong Zhang }
3179b54502bSHong Zhang EXTERN_C_END
3189b54502bSHong Zhang 
3199b54502bSHong Zhang EXTERN_C_BEGIN
3209b54502bSHong Zhang #undef __FUNCT__
3219b54502bSHong Zhang #define __FUNCT__ "PCLUSetPivotInBlocks_LU"
3229b54502bSHong Zhang PetscErrorCode PCLUSetPivotInBlocks_LU(PC pc,PetscTruth pivot)
3239b54502bSHong Zhang {
3249b54502bSHong Zhang   PC_LU *dir = (PC_LU*)pc->data;
3259b54502bSHong Zhang 
3269b54502bSHong Zhang   PetscFunctionBegin;
3279b54502bSHong Zhang   dir->info.pivotinblocks = pivot ? 1.0 : 0.0;
3289b54502bSHong Zhang   PetscFunctionReturn(0);
3299b54502bSHong Zhang }
3309b54502bSHong Zhang EXTERN_C_END
3319b54502bSHong Zhang 
3329b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/
3339b54502bSHong Zhang 
3349b54502bSHong Zhang #undef __FUNCT__
3359b54502bSHong Zhang #define __FUNCT__ "PCLUReorderForNonzeroDiagonal"
3369b54502bSHong Zhang /*@
3379b54502bSHong Zhang    PCLUReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal
3389b54502bSHong Zhang 
3399b54502bSHong Zhang    Collective on PC
3409b54502bSHong Zhang 
3419b54502bSHong Zhang    Input Parameters:
3429b54502bSHong Zhang +  pc - the preconditioner context
3439b54502bSHong Zhang -  tol - diagonal entries smaller than this in absolute value are considered zero
3449b54502bSHong Zhang 
3459b54502bSHong Zhang    Options Database Key:
3469b54502bSHong Zhang .  -pc_lu_nonzeros_along_diagonal
3479b54502bSHong Zhang 
3489b54502bSHong Zhang    Level: intermediate
3499b54502bSHong Zhang 
3509b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill
3519b54502bSHong Zhang 
352*ee45ca4aSHong Zhang .seealso: PCLUSetFill(), PCFactorSetShiftNonzero(), PCFactorSetZeroPivot(), MatReorderForNonzeroDiagonal()
3539b54502bSHong Zhang @*/
3549b54502bSHong Zhang PetscErrorCode PCLUReorderForNonzeroDiagonal(PC pc,PetscReal rtol)
3559b54502bSHong Zhang {
3569b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscReal);
3579b54502bSHong Zhang 
3589b54502bSHong Zhang   PetscFunctionBegin;
3599b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
3609b54502bSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUReorderForNonzeroDiagonal_C",(void (**)(void))&f);CHKERRQ(ierr);
3619b54502bSHong Zhang   if (f) {
3629b54502bSHong Zhang     ierr = (*f)(pc,rtol);CHKERRQ(ierr);
3639b54502bSHong Zhang   }
3649b54502bSHong Zhang   PetscFunctionReturn(0);
3659b54502bSHong Zhang }
3669b54502bSHong Zhang 
3679b54502bSHong Zhang #undef __FUNCT__
3689b54502bSHong Zhang #define __FUNCT__ "PCLUSetReuseOrdering"
3699b54502bSHong Zhang /*@
3709b54502bSHong Zhang    PCLUSetReuseOrdering - When similar matrices are factored, this
3719b54502bSHong Zhang    causes the ordering computed in the first factor to be used for all
3729b54502bSHong Zhang    following factors; applies to both fill and drop tolerance LUs.
3739b54502bSHong Zhang 
3749b54502bSHong Zhang    Collective on PC
3759b54502bSHong Zhang 
3769b54502bSHong Zhang    Input Parameters:
3779b54502bSHong Zhang +  pc - the preconditioner context
3789b54502bSHong Zhang -  flag - PETSC_TRUE to reuse else PETSC_FALSE
3799b54502bSHong Zhang 
3809b54502bSHong Zhang    Options Database Key:
3819b54502bSHong Zhang .  -pc_lu_reuse_ordering - Activate PCLUSetReuseOrdering()
3829b54502bSHong Zhang 
3839b54502bSHong Zhang    Level: intermediate
3849b54502bSHong Zhang 
3859b54502bSHong Zhang .keywords: PC, levels, reordering, factorization, incomplete, LU
3869b54502bSHong Zhang 
3879b54502bSHong Zhang .seealso: PCLUSetReuseFill(), PCILUSetReuseOrdering(), PCILUDTSetReuseFill()
3889b54502bSHong Zhang @*/
3899b54502bSHong Zhang PetscErrorCode PCLUSetReuseOrdering(PC pc,PetscTruth flag)
3909b54502bSHong Zhang {
3919b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscTruth);
3929b54502bSHong Zhang 
3939b54502bSHong Zhang   PetscFunctionBegin;
3949b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
3959b54502bSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetReuseOrdering_C",(void (**)(void))&f);CHKERRQ(ierr);
3969b54502bSHong Zhang   if (f) {
3979b54502bSHong Zhang     ierr = (*f)(pc,flag);CHKERRQ(ierr);
3989b54502bSHong Zhang   }
3999b54502bSHong Zhang   PetscFunctionReturn(0);
4009b54502bSHong Zhang }
4019b54502bSHong Zhang 
4029b54502bSHong Zhang #undef __FUNCT__
4039b54502bSHong Zhang #define __FUNCT__ "PCLUSetReuseFill"
4049b54502bSHong Zhang /*@
4059b54502bSHong Zhang    PCLUSetReuseFill - When matrices with same nonzero structure are LU factored,
4069b54502bSHong Zhang    this causes later ones to use the fill computed in the initial factorization.
4079b54502bSHong Zhang 
4089b54502bSHong Zhang    Collective on PC
4099b54502bSHong Zhang 
4109b54502bSHong Zhang    Input Parameters:
4119b54502bSHong Zhang +  pc - the preconditioner context
4129b54502bSHong Zhang -  flag - PETSC_TRUE to reuse else PETSC_FALSE
4139b54502bSHong Zhang 
4149b54502bSHong Zhang    Options Database Key:
4159b54502bSHong Zhang .  -pc_lu_reuse_fill - Activates PCLUSetReuseFill()
4169b54502bSHong Zhang 
4179b54502bSHong Zhang    Level: intermediate
4189b54502bSHong Zhang 
4199b54502bSHong Zhang .keywords: PC, levels, reordering, factorization, incomplete, LU
4209b54502bSHong Zhang 
4219b54502bSHong Zhang .seealso: PCILUSetReuseOrdering(), PCLUSetReuseOrdering(), PCILUDTSetReuseFill()
4229b54502bSHong Zhang @*/
4239b54502bSHong Zhang PetscErrorCode PCLUSetReuseFill(PC pc,PetscTruth flag)
4249b54502bSHong Zhang {
4259b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscTruth);
4269b54502bSHong Zhang 
4279b54502bSHong Zhang   PetscFunctionBegin;
4289b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
4299b54502bSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetReuseFill_C",(void (**)(void))&f);CHKERRQ(ierr);
4309b54502bSHong Zhang   if (f) {
4319b54502bSHong Zhang     ierr = (*f)(pc,flag);CHKERRQ(ierr);
4329b54502bSHong Zhang   }
4339b54502bSHong Zhang   PetscFunctionReturn(0);
4349b54502bSHong Zhang }
4359b54502bSHong Zhang 
4369b54502bSHong Zhang #undef __FUNCT__
4379b54502bSHong Zhang #define __FUNCT__ "PCLUSetFill"
4389b54502bSHong Zhang /*@
4399b54502bSHong Zhang    PCLUSetFill - Indicate the amount of fill you expect in the factored matrix,
4409b54502bSHong Zhang    fill = number nonzeros in factor/number nonzeros in original matrix.
4419b54502bSHong Zhang 
4429b54502bSHong Zhang    Collective on PC
4439b54502bSHong Zhang 
4449b54502bSHong Zhang    Input Parameters:
4459b54502bSHong Zhang +  pc - the preconditioner context
4469b54502bSHong Zhang -  fill - amount of expected fill
4479b54502bSHong Zhang 
4489b54502bSHong Zhang    Options Database Key:
4499b54502bSHong Zhang .  -pc_lu_fill <fill> - Sets fill amount
4509b54502bSHong Zhang 
4519b54502bSHong Zhang    Level: intermediate
4529b54502bSHong Zhang 
4539b54502bSHong Zhang    Note:
4549b54502bSHong Zhang    For sparse matrix factorizations it is difficult to predict how much
4559b54502bSHong Zhang    fill to expect. By running with the option -log_info PETSc will print the
4569b54502bSHong Zhang    actual amount of fill used; allowing you to set the value accurately for
4579b54502bSHong Zhang    future runs. Default PETSc uses a value of 5.0
4589b54502bSHong Zhang 
4599b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill
4609b54502bSHong Zhang 
4619b54502bSHong Zhang .seealso: PCILUSetFill()
4629b54502bSHong Zhang @*/
4639b54502bSHong Zhang PetscErrorCode PCLUSetFill(PC pc,PetscReal fill)
4649b54502bSHong Zhang {
4659b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscReal);
4669b54502bSHong Zhang 
4679b54502bSHong Zhang   PetscFunctionBegin;
4689b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
4699b54502bSHong Zhang   if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0");
4709b54502bSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetFill_C",(void (**)(void))&f);CHKERRQ(ierr);
4719b54502bSHong Zhang   if (f) {
4729b54502bSHong Zhang     ierr = (*f)(pc,fill);CHKERRQ(ierr);
4739b54502bSHong Zhang   }
4749b54502bSHong Zhang   PetscFunctionReturn(0);
4759b54502bSHong Zhang }
4769b54502bSHong Zhang 
4779b54502bSHong Zhang #undef __FUNCT__
4789b54502bSHong Zhang #define __FUNCT__ "PCLUSetUseInPlace"
4799b54502bSHong Zhang /*@
4809b54502bSHong Zhang    PCLUSetUseInPlace - Tells the system to do an in-place factorization.
4819b54502bSHong Zhang    For dense matrices, this enables the solution of much larger problems.
4829b54502bSHong Zhang    For sparse matrices the factorization cannot be done truly in-place
4839b54502bSHong Zhang    so this does not save memory during the factorization, but after the matrix
4849b54502bSHong Zhang    is factored, the original unfactored matrix is freed, thus recovering that
4859b54502bSHong Zhang    space.
4869b54502bSHong Zhang 
4879b54502bSHong Zhang    Collective on PC
4889b54502bSHong Zhang 
4899b54502bSHong Zhang    Input Parameters:
4909b54502bSHong Zhang .  pc - the preconditioner context
4919b54502bSHong Zhang 
4929b54502bSHong Zhang    Options Database Key:
4939b54502bSHong Zhang .  -pc_lu_in_place - Activates in-place factorization
4949b54502bSHong Zhang 
4959b54502bSHong Zhang    Notes:
4969b54502bSHong Zhang    PCLUSetUseInplace() can only be used with the KSP method KSPPREONLY or when
4979b54502bSHong Zhang    a different matrix is provided for the multiply and the preconditioner in
4989b54502bSHong Zhang    a call to KSPSetOperators().
4999b54502bSHong Zhang    This is because the Krylov space methods require an application of the
5009b54502bSHong Zhang    matrix multiplication, which is not possible here because the matrix has
5019b54502bSHong Zhang    been factored in-place, replacing the original matrix.
5029b54502bSHong Zhang 
5039b54502bSHong Zhang    Level: intermediate
5049b54502bSHong Zhang 
5059b54502bSHong Zhang .keywords: PC, set, factorization, direct, inplace, in-place, LU
5069b54502bSHong Zhang 
5079b54502bSHong Zhang .seealso: PCILUSetUseInPlace()
5089b54502bSHong Zhang @*/
5099b54502bSHong Zhang PetscErrorCode PCLUSetUseInPlace(PC pc)
5109b54502bSHong Zhang {
5119b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC);
5129b54502bSHong Zhang 
5139b54502bSHong Zhang   PetscFunctionBegin;
5149b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
5159b54502bSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetUseInPlace_C",(void (**)(void))&f);CHKERRQ(ierr);
5169b54502bSHong Zhang   if (f) {
5179b54502bSHong Zhang     ierr = (*f)(pc);CHKERRQ(ierr);
5189b54502bSHong Zhang   }
5199b54502bSHong Zhang   PetscFunctionReturn(0);
5209b54502bSHong Zhang }
5219b54502bSHong Zhang 
5229b54502bSHong Zhang #undef __FUNCT__
5239b54502bSHong Zhang #define __FUNCT__ "PCLUSetMatOrdering"
5249b54502bSHong Zhang /*@C
5259b54502bSHong Zhang     PCLUSetMatOrdering - Sets the ordering routine (to reduce fill) to
5269b54502bSHong Zhang     be used in the LU factorization.
5279b54502bSHong Zhang 
5289b54502bSHong Zhang     Collective on PC
5299b54502bSHong Zhang 
5309b54502bSHong Zhang     Input Parameters:
5319b54502bSHong Zhang +   pc - the preconditioner context
5329b54502bSHong Zhang -   ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM
5339b54502bSHong Zhang 
5349b54502bSHong Zhang     Options Database Key:
5359b54502bSHong Zhang .   -pc_lu_mat_ordering_type <nd,rcm,...> - Sets ordering routine
5369b54502bSHong Zhang 
5379b54502bSHong Zhang     Level: intermediate
5389b54502bSHong Zhang 
5399b54502bSHong Zhang     Notes: nested dissection is used by default
5409b54502bSHong Zhang 
5419b54502bSHong Zhang .seealso: PCILUSetMatOrdering()
5429b54502bSHong Zhang @*/
5439b54502bSHong Zhang PetscErrorCode PCLUSetMatOrdering(PC pc,MatOrderingType ordering)
5449b54502bSHong Zhang {
5459b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,MatOrderingType);
5469b54502bSHong Zhang 
5479b54502bSHong Zhang   PetscFunctionBegin;
5489b54502bSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetMatOrdering_C",(void (**)(void))&f);CHKERRQ(ierr);
5499b54502bSHong Zhang   if (f) {
5509b54502bSHong Zhang     ierr = (*f)(pc,ordering);CHKERRQ(ierr);
5519b54502bSHong Zhang   }
5529b54502bSHong Zhang   PetscFunctionReturn(0);
5539b54502bSHong Zhang }
5549b54502bSHong Zhang 
5559b54502bSHong Zhang #undef __FUNCT__
5569b54502bSHong Zhang #define __FUNCT__ "PCLUSetPivoting"
5579b54502bSHong Zhang /*@
5589b54502bSHong Zhang     PCLUSetPivoting - Determines when pivoting is done during LU.
5599b54502bSHong Zhang       For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices
5609b54502bSHong Zhang       it is never done. For the Matlab and SuperLU factorization this is used.
5619b54502bSHong Zhang 
5629b54502bSHong Zhang     Collective on PC
5639b54502bSHong Zhang 
5649b54502bSHong Zhang     Input Parameters:
5659b54502bSHong Zhang +   pc - the preconditioner context
5669b54502bSHong Zhang -   dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable)
5679b54502bSHong Zhang 
5689b54502bSHong Zhang     Options Database Key:
5699b54502bSHong Zhang .   -pc_lu_pivoting <dtcol>
5709b54502bSHong Zhang 
5719b54502bSHong Zhang     Level: intermediate
5729b54502bSHong Zhang 
5739b54502bSHong Zhang .seealso: PCILUSetMatOrdering(), PCLUSetPivotInBlocks()
5749b54502bSHong Zhang @*/
5759b54502bSHong Zhang PetscErrorCode PCLUSetPivoting(PC pc,PetscReal dtcol)
5769b54502bSHong Zhang {
5779b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscReal);
5789b54502bSHong Zhang 
5799b54502bSHong Zhang   PetscFunctionBegin;
5809b54502bSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetPivoting_C",(void (**)(void))&f);CHKERRQ(ierr);
5819b54502bSHong Zhang   if (f) {
5829b54502bSHong Zhang     ierr = (*f)(pc,dtcol);CHKERRQ(ierr);
5839b54502bSHong Zhang   }
5849b54502bSHong Zhang   PetscFunctionReturn(0);
5859b54502bSHong Zhang }
5869b54502bSHong Zhang 
5879b54502bSHong Zhang #undef __FUNCT__
5889b54502bSHong Zhang #define __FUNCT__ "PCLUSetPivotInBlocks"
5899b54502bSHong Zhang /*@
5909b54502bSHong Zhang     PCLUSetPivotInBlocks - Determines if pivoting is done while factoring each block
5919b54502bSHong Zhang       with BAIJ or SBAIJ matrices
5929b54502bSHong Zhang 
5939b54502bSHong Zhang     Collective on PC
5949b54502bSHong Zhang 
5959b54502bSHong Zhang     Input Parameters:
5969b54502bSHong Zhang +   pc - the preconditioner context
5979b54502bSHong Zhang -   pivot - PETSC_TRUE or PETSC_FALSE
5989b54502bSHong Zhang 
5999b54502bSHong Zhang     Options Database Key:
6009b54502bSHong Zhang .   -pc_lu_pivot_in_blocks <true,false>
6019b54502bSHong Zhang 
6029b54502bSHong Zhang     Level: intermediate
6039b54502bSHong Zhang 
6049b54502bSHong Zhang .seealso: PCILUSetMatOrdering(), PCLUSetPivoting()
6059b54502bSHong Zhang @*/
6069b54502bSHong Zhang PetscErrorCode PCLUSetPivotInBlocks(PC pc,PetscTruth pivot)
6079b54502bSHong Zhang {
6089b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscTruth);
6099b54502bSHong Zhang 
6109b54502bSHong Zhang   PetscFunctionBegin;
6119b54502bSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetPivotInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr);
6129b54502bSHong Zhang   if (f) {
6139b54502bSHong Zhang     ierr = (*f)(pc,pivot);CHKERRQ(ierr);
6149b54502bSHong Zhang   }
6159b54502bSHong Zhang   PetscFunctionReturn(0);
6169b54502bSHong Zhang }
6179b54502bSHong Zhang 
6189b54502bSHong Zhang /* ------------------------------------------------------------------------ */
6199b54502bSHong Zhang 
6209b54502bSHong Zhang /*MC
6219b54502bSHong Zhang    PCLU - Uses a direct solver, based on LU factorization, as a preconditioner
6229b54502bSHong Zhang 
6239b54502bSHong Zhang    Options Database Keys:
6249b54502bSHong Zhang +  -pc_lu_reuse_ordering - Activate PCLUSetReuseOrdering()
6259b54502bSHong Zhang .  -pc_lu_reuse_fill - Activates PCLUSetReuseFill()
6269b54502bSHong Zhang .  -pc_lu_fill <fill> - Sets fill amount
6279b54502bSHong Zhang .  -pc_lu_in_place - Activates in-place factorization
6289b54502bSHong Zhang .  -pc_lu_mat_ordering_type <nd,rcm,...> - Sets ordering routine
6299b54502bSHong Zhang -  -pc_lu_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase
6309b54502bSHong Zhang                                          stability of factorization.
6319b54502bSHong Zhang 
6329b54502bSHong Zhang    Notes: Not all options work for all matrix formats
6339b54502bSHong Zhang           Run with -help to see additional options for particular matrix formats or factorization
6349b54502bSHong Zhang           algorithms
6359b54502bSHong Zhang 
6369b54502bSHong Zhang    Level: beginner
6379b54502bSHong Zhang 
6389b54502bSHong Zhang    Concepts: LU factorization, direct solver
6399b54502bSHong Zhang 
6409b54502bSHong Zhang    Notes: Usually this will compute an "exact" solution in one iteration and does
6419b54502bSHong Zhang           not need a Krylov method (i.e. you can use -ksp_type preonly, or
6429b54502bSHong Zhang           KSPSetType(ksp,KSPPREONLY) for the Krylov method
6439b54502bSHong Zhang 
6449b54502bSHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
6459b54502bSHong Zhang            PCILU, PCCHOLESKY, PCICC, PCLUSetReuseOrdering(), PCLUSetReuseFill(), PCGetFactoredMatrix(),
646*ee45ca4aSHong Zhang            PCLUSetFill(), PCLUSetUseInPlace(), PCLUSetMatOrdering(), PCFactorSetPivoting(),
6479b54502bSHong Zhang            PCLUSetPivotingInBlocks()
6489b54502bSHong Zhang M*/
6499b54502bSHong Zhang 
6509b54502bSHong Zhang EXTERN_C_BEGIN
6519b54502bSHong Zhang #undef __FUNCT__
6529b54502bSHong Zhang #define __FUNCT__ "PCCreate_LU"
6539b54502bSHong Zhang PetscErrorCode PCCreate_LU(PC pc)
6549b54502bSHong Zhang {
6559b54502bSHong Zhang   PetscErrorCode ierr;
6569b54502bSHong Zhang   PetscMPIInt    size;
6579b54502bSHong Zhang   PC_LU          *dir;
6589b54502bSHong Zhang 
6599b54502bSHong Zhang   PetscFunctionBegin;
6609b54502bSHong Zhang   ierr = PetscNew(PC_LU,&dir);CHKERRQ(ierr);
6619b54502bSHong Zhang   PetscLogObjectMemory(pc,sizeof(PC_LU));
6629b54502bSHong Zhang 
6639b54502bSHong Zhang   ierr = MatFactorInfoInitialize(&dir->info);CHKERRQ(ierr);
6649b54502bSHong Zhang   dir->fact                  = 0;
6659b54502bSHong Zhang   dir->inplace               = PETSC_FALSE;
6669b54502bSHong Zhang   dir->nonzerosalongdiagonal = PETSC_FALSE;
6679b54502bSHong Zhang 
6689b54502bSHong Zhang   dir->info.fill           = 5.0;
6699b54502bSHong Zhang   dir->info.dtcol          = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */
6700a29876aSHong Zhang   dir->info.shiftnz        = 0.0;
6719b54502bSHong Zhang   dir->info.zeropivot      = 1.e-12;
6729b54502bSHong Zhang   dir->info.pivotinblocks  = 1.0;
6730a29876aSHong Zhang   dir->info.shiftpd        = PETSC_FALSE;
6749b54502bSHong Zhang   dir->info.shift_fraction = 0.0;
6759b54502bSHong Zhang   dir->col                 = 0;
6769b54502bSHong Zhang   dir->row                 = 0;
6779b54502bSHong Zhang   ierr = MPI_Comm_size(pc->comm,&size);CHKERRQ(ierr);
6789b54502bSHong Zhang   if (size == 1) {
6799b54502bSHong Zhang     ierr = PetscStrallocpy(MATORDERING_ND,&dir->ordering);CHKERRQ(ierr);
6809b54502bSHong Zhang   } else {
6819b54502bSHong Zhang     ierr = PetscStrallocpy(MATORDERING_NATURAL,&dir->ordering);CHKERRQ(ierr);
6829b54502bSHong Zhang   }
6839b54502bSHong Zhang   dir->reusefill        = PETSC_FALSE;
6849b54502bSHong Zhang   dir->reuseordering    = PETSC_FALSE;
6859b54502bSHong Zhang   pc->data              = (void*)dir;
6869b54502bSHong Zhang 
6879b54502bSHong Zhang   pc->ops->destroy           = PCDestroy_LU;
6889b54502bSHong Zhang   pc->ops->apply             = PCApply_LU;
6899b54502bSHong Zhang   pc->ops->applytranspose    = PCApplyTranspose_LU;
6909b54502bSHong Zhang   pc->ops->setup             = PCSetUp_LU;
6919b54502bSHong Zhang   pc->ops->setfromoptions    = PCSetFromOptions_LU;
6929b54502bSHong Zhang   pc->ops->view              = PCView_LU;
6939b54502bSHong Zhang   pc->ops->applyrichardson   = 0;
6949b54502bSHong Zhang   pc->ops->getfactoredmatrix = PCGetFactoredMatrix_LU;
6959b54502bSHong Zhang 
6969b54502bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetFill_C","PCLUSetFill_LU",
6979b54502bSHong Zhang                     PCLUSetFill_LU);CHKERRQ(ierr);
6989b54502bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetUseInPlace_C","PCLUSetUseInPlace_LU",
6999b54502bSHong Zhang                     PCLUSetUseInPlace_LU);CHKERRQ(ierr);
7009b54502bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetMatOrdering_C","PCLUSetMatOrdering_LU",
7019b54502bSHong Zhang                     PCLUSetMatOrdering_LU);CHKERRQ(ierr);
7029b54502bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetReuseOrdering_C","PCLUSetReuseOrdering_LU",
7039b54502bSHong Zhang                     PCLUSetReuseOrdering_LU);CHKERRQ(ierr);
7049b54502bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetReuseFill_C","PCLUSetReuseFill_LU",
7059b54502bSHong Zhang                     PCLUSetReuseFill_LU);CHKERRQ(ierr);
7069b54502bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetPivoting_C","PCLUSetPivoting_LU",
7079b54502bSHong Zhang                     PCLUSetPivoting_LU);CHKERRQ(ierr);
7089b54502bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetPivotInBlocks_C","PCLUSetPivotInBlocks_LU",
7099b54502bSHong Zhang                     PCLUSetPivotInBlocks_LU);CHKERRQ(ierr);
7109b54502bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUReorderForNonzeroDiagonal_C","PCLUReorderForNonzeroDiagonal_LU",
7119b54502bSHong Zhang                     PCLUReorderForNonzeroDiagonal_LU);CHKERRQ(ierr);
7129b54502bSHong Zhang   PetscFunctionReturn(0);
7139b54502bSHong Zhang }
7149b54502bSHong Zhang EXTERN_C_END
715