xref: /petsc/src/ksp/pc/impls/factor/lu/lu.c (revision 7112b564e8d203471dda2cf2b83517e67bda591d)
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 
9075768bcSBarry Smith #include "src/ksp/pc/impls/factor/lu/lu.h"  /*I "petscpc.h" I*/
109b54502bSHong Zhang 
119b54502bSHong Zhang EXTERN_C_BEGIN
129b54502bSHong Zhang #undef __FUNCT__
132401956bSBarry Smith #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal_LU"
142401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal_LU(PC pc,PetscReal z)
159b54502bSHong Zhang {
169b54502bSHong Zhang   PC_LU *lu = (PC_LU*)pc->data;
179b54502bSHong Zhang 
189b54502bSHong Zhang   PetscFunctionBegin;
199b54502bSHong Zhang   lu->nonzerosalongdiagonal = PETSC_TRUE;
209b54502bSHong Zhang   if (z == PETSC_DECIDE) {
219b54502bSHong Zhang     lu->nonzerosalongdiagonaltol = 1.e-10;
229b54502bSHong Zhang   } else {
239b54502bSHong Zhang     lu->nonzerosalongdiagonaltol = z;
249b54502bSHong Zhang   }
259b54502bSHong Zhang   PetscFunctionReturn(0);
269b54502bSHong Zhang }
279b54502bSHong Zhang EXTERN_C_END
289b54502bSHong Zhang 
299b54502bSHong Zhang EXTERN_C_BEGIN
309b54502bSHong Zhang #undef __FUNCT__
312401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering_LU"
322401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering_LU(PC pc,PetscTruth flag)
339b54502bSHong Zhang {
34c62fe420SBarry Smith   PC_LU *lu = (PC_LU*)pc->data;
359b54502bSHong Zhang 
369b54502bSHong Zhang   PetscFunctionBegin;
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__
442401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseFill_LU"
452401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill_LU(PC pc,PetscTruth flag)
469b54502bSHong Zhang {
47c62fe420SBarry Smith   PC_LU *lu = (PC_LU*)pc->data;
489b54502bSHong Zhang 
499b54502bSHong Zhang   PetscFunctionBegin;
509b54502bSHong Zhang   lu->reusefill = flag;
519b54502bSHong Zhang   PetscFunctionReturn(0);
529b54502bSHong Zhang }
539b54502bSHong Zhang EXTERN_C_END
549b54502bSHong Zhang 
559b54502bSHong Zhang #undef __FUNCT__
569b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_LU"
579b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_LU(PC pc)
589b54502bSHong Zhang {
599b54502bSHong Zhang   PC_LU           *lu = (PC_LU*)pc->data;
609b54502bSHong Zhang   PetscErrorCode  ierr;
619b54502bSHong Zhang   PetscTruth      flg,set;
629989ab13SBarry Smith   char            tname[256], solvertype[64];
639b54502bSHong Zhang   PetscFList      ordlist;
649b54502bSHong Zhang   PetscReal       tol;
659b54502bSHong Zhang 
669b54502bSHong Zhang   PetscFunctionBegin;
679b54502bSHong Zhang   ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr);
689b54502bSHong Zhang   ierr = PetscOptionsHead("LU options");CHKERRQ(ierr);
692401956bSBarry Smith     ierr = PetscOptionsName("-pc_factor_in_place","Form LU in the same memory as the matrix","PCFactorSetUseInPlace",&flg);CHKERRQ(ierr);
709b54502bSHong Zhang     if (flg) {
712401956bSBarry Smith       ierr = PCFactorSetUseInPlace(pc);CHKERRQ(ierr);
729b54502bSHong Zhang     }
73075768bcSBarry Smith     ierr = PetscOptionsReal("-pc_factor_fill","Expected non-zeros in LU/non-zeros in matrix","PCFactorSetFill",((PC_Factor*)lu)->info.fill,&((PC_Factor*)lu)->info.fill,0);CHKERRQ(ierr);
749b54502bSHong Zhang 
759f95998fSHong Zhang     ierr = PetscOptionsName("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr);
769b54502bSHong Zhang     if (flg) {
77afaefe49SHong Zhang         ierr = PCFactorSetShiftNonzero(pc,(PetscReal) PETSC_DECIDE);CHKERRQ(ierr);
789b54502bSHong Zhang     }
79075768bcSBarry Smith     ierr = PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",((PC_Factor*)lu)->info.shiftnz,&((PC_Factor*)lu)->info.shiftnz,0);CHKERRQ(ierr);
809f95998fSHong Zhang     ierr = PetscOptionsName("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",&flg);CHKERRQ(ierr);
819b54502bSHong Zhang     if (flg) {
82afaefe49SHong Zhang       ierr = PCFactorSetShiftPd(pc,PETSC_TRUE);CHKERRQ(ierr);
839b54502bSHong Zhang     }
84075768bcSBarry Smith     ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",((PC_Factor*)lu)->info.zeropivot,&((PC_Factor*)lu)->info.zeropivot,0);CHKERRQ(ierr);
859b54502bSHong Zhang 
862401956bSBarry Smith     ierr = PetscOptionsName("-pc_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",&flg);CHKERRQ(ierr);
879b54502bSHong Zhang     if (flg) {
882401956bSBarry Smith       ierr = PCFactorSetReuseFill(pc,PETSC_TRUE);CHKERRQ(ierr);
899b54502bSHong Zhang     }
902401956bSBarry Smith     ierr = PetscOptionsName("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",&flg);CHKERRQ(ierr);
919b54502bSHong Zhang     if (flg) {
922401956bSBarry Smith       ierr = PCFactorSetReuseOrdering(pc,PETSC_TRUE);CHKERRQ(ierr);
939b54502bSHong Zhang     }
949b54502bSHong Zhang 
959b54502bSHong Zhang     ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr);
96075768bcSBarry Smith     ierr = PetscOptionsList("-pc_factor_mat_ordering_type","Reordering to reduce nonzeros in LU","PCFactorSetMatOrderingType",ordlist,((PC_Factor*)lu)->ordering,tname,256,&flg);CHKERRQ(ierr);
979b54502bSHong Zhang     if (flg) {
98e5a9bf91SBarry Smith       ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr);
999b54502bSHong Zhang     }
1009b54502bSHong Zhang 
1015c9eb25fSBarry Smith     /* maybe should have MatGetSolverTypes(Mat,&list) like the ordering list */
102075768bcSBarry Smith     ierr = PetscOptionsString("-pc_factor_mat_solver_package","Specific LU solver to use","MatGetFactor",((PC_Factor*)lu)->solvertype,solvertype,64,&flg);CHKERRQ(ierr);
1035c9eb25fSBarry Smith     if (flg) {
104c7393fdbSBarry Smith       ierr = PCFactorSetMatSolverPackage(pc,solvertype);CHKERRQ(ierr);
1055c9eb25fSBarry Smith     }
1065c9eb25fSBarry Smith 
1072401956bSBarry Smith     ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr);
1089b54502bSHong Zhang     if (flg) {
1099b54502bSHong Zhang       tol = PETSC_DECIDE;
1102401956bSBarry Smith       ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr);
1112401956bSBarry Smith       ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr);
1129b54502bSHong Zhang     }
1139b54502bSHong Zhang 
114075768bcSBarry Smith     ierr = PetscOptionsReal("-pc_factor_pivoting","Pivoting tolerance (used only for some factorization)","PCFactorSetPivoting",((PC_Factor*)lu)->info.dtcol,&((PC_Factor*)lu)->info.dtcol,&flg);CHKERRQ(ierr);
1159b54502bSHong Zhang 
116075768bcSBarry Smith     flg = ((PC_Factor*)lu)->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE;
1172401956bSBarry Smith     ierr = PetscOptionsTruth("-pc_factor_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCFactorSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr);
1189b54502bSHong Zhang     if (set) {
1192401956bSBarry Smith       ierr = PCFactorSetPivotInBlocks(pc,flg);CHKERRQ(ierr);
1209b54502bSHong Zhang     }
1219b54502bSHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
1229b54502bSHong Zhang   PetscFunctionReturn(0);
1239b54502bSHong Zhang }
1249b54502bSHong Zhang 
1259b54502bSHong Zhang #undef __FUNCT__
1269b54502bSHong Zhang #define __FUNCT__ "PCView_LU"
1279b54502bSHong Zhang static PetscErrorCode PCView_LU(PC pc,PetscViewer viewer)
1289b54502bSHong Zhang {
1299b54502bSHong Zhang   PC_LU          *lu = (PC_LU*)pc->data;
1309b54502bSHong Zhang   PetscErrorCode ierr;
1319b54502bSHong Zhang   PetscTruth     iascii,isstring;
1329b54502bSHong Zhang 
1339b54502bSHong Zhang   PetscFunctionBegin;
1349b54502bSHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
1359b54502bSHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
1369b54502bSHong Zhang   if (iascii) {
1379b54502bSHong Zhang 
1389b54502bSHong Zhang     if (lu->inplace) {ierr = PetscViewerASCIIPrintf(viewer,"  LU: in-place factorization\n");CHKERRQ(ierr);}
1399b54502bSHong Zhang     else             {ierr = PetscViewerASCIIPrintf(viewer,"  LU: out-of-place factorization\n");CHKERRQ(ierr);}
140075768bcSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    matrix ordering: %s\n",((PC_Factor*)lu)->ordering);CHKERRQ(ierr);
141075768bcSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  LU: tolerance for zero pivot %G\n",((PC_Factor*)lu)->info.zeropivot);CHKERRQ(ierr);
142075768bcSBarry Smith     if (((PC_Factor*)lu)->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer,"  LU: using Manteuffel shift\n");CHKERRQ(ierr);}
143075768bcSBarry Smith     if (((PC_Factor*)lu)->fact) {
144f3a39becSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  LU: factor fill ratio needed %G\n",lu->actualfill);CHKERRQ(ierr);
145f3a39becSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"       Factored matrix follows\n");CHKERRQ(ierr);
146f3a39becSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
147f3a39becSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
148f3a39becSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
149f3a39becSBarry Smith       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
150075768bcSBarry Smith       ierr = MatView(((PC_Factor*)lu)->fact,viewer);CHKERRQ(ierr);
1519b54502bSHong Zhang       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
152f3a39becSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
153f3a39becSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
154f3a39becSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1559b54502bSHong Zhang     }
1569b54502bSHong Zhang     if (lu->reusefill)    {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing fill from past factorization\n");CHKERRQ(ierr);}
1579b54502bSHong Zhang     if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing reordering from past factorization\n");CHKERRQ(ierr);}
1589b54502bSHong Zhang   } else if (isstring) {
159075768bcSBarry Smith     ierr = PetscViewerStringSPrintf(viewer," order=%s",((PC_Factor*)lu)->ordering);CHKERRQ(ierr);CHKERRQ(ierr);
1609b54502bSHong Zhang   } else {
1619b54502bSHong Zhang     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCLU",((PetscObject)viewer)->type_name);
1629b54502bSHong Zhang   }
1639b54502bSHong Zhang   PetscFunctionReturn(0);
1649b54502bSHong Zhang }
1659b54502bSHong Zhang 
1669b54502bSHong Zhang #undef __FUNCT__
1679b54502bSHong Zhang #define __FUNCT__ "PCSetUp_LU"
1689b54502bSHong Zhang static PetscErrorCode PCSetUp_LU(PC pc)
1699b54502bSHong Zhang {
1709b54502bSHong Zhang   PetscErrorCode ierr;
1719b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
1729b54502bSHong Zhang 
1739b54502bSHong Zhang   PetscFunctionBegin;
174075768bcSBarry Smith   if (dir->reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->actualfill;
1759b54502bSHong Zhang 
1769b54502bSHong Zhang   if (dir->inplace) {
1779b54502bSHong Zhang     if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);}
1789b54502bSHong Zhang     if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);}
179075768bcSBarry Smith     ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
18003c60df9SBarry Smith     if (dir->row) {
18103c60df9SBarry Smith       ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);
18203c60df9SBarry Smith       ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr);
18303c60df9SBarry Smith     }
184075768bcSBarry Smith     ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
185075768bcSBarry Smith     ((PC_Factor*)dir)->fact = pc->pmat;
1869b54502bSHong Zhang   } else {
1879b54502bSHong Zhang     MatInfo info;
1889b54502bSHong Zhang     if (!pc->setupcalled) {
189075768bcSBarry Smith       ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
1909b54502bSHong Zhang       if (dir->nonzerosalongdiagonal) {
1919b54502bSHong Zhang         ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr);
1929b54502bSHong Zhang       }
19303c60df9SBarry Smith       if (dir->row) {
19403c60df9SBarry Smith         ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);
19503c60df9SBarry Smith         ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr);
19603c60df9SBarry Smith       }
197075768bcSBarry Smith       ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr);
198075768bcSBarry Smith       ierr = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
199075768bcSBarry Smith       ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
2009b54502bSHong Zhang       dir->actualfill = info.fill_ratio_needed;
201075768bcSBarry Smith       ierr = PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);CHKERRQ(ierr);
2029b54502bSHong Zhang     } else if (pc->flag != SAME_NONZERO_PATTERN) {
2039b54502bSHong Zhang       if (!dir->reuseordering) {
2049b54502bSHong Zhang         if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);}
2059b54502bSHong Zhang         if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);}
206075768bcSBarry Smith         ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
2079b54502bSHong Zhang         if (dir->nonzerosalongdiagonal) {
2089b54502bSHong Zhang           ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr);
2099b54502bSHong Zhang         }
21003c60df9SBarry Smith         if (dir->row) {
21103c60df9SBarry Smith           ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);
21203c60df9SBarry Smith           ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr);
21303c60df9SBarry Smith         }
2149b54502bSHong Zhang       }
215075768bcSBarry Smith       ierr = MatDestroy(((PC_Factor*)dir)->fact);CHKERRQ(ierr);
216075768bcSBarry Smith       ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr);
217075768bcSBarry Smith       ierr = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
218075768bcSBarry Smith       ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
2199b54502bSHong Zhang       dir->actualfill = info.fill_ratio_needed;
220075768bcSBarry Smith       ierr = PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);CHKERRQ(ierr);
2219b54502bSHong Zhang     }
222075768bcSBarry Smith     ierr = MatLUFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
2239b54502bSHong Zhang   }
2249b54502bSHong Zhang   PetscFunctionReturn(0);
2259b54502bSHong Zhang }
2269b54502bSHong Zhang 
2279b54502bSHong Zhang #undef __FUNCT__
2289b54502bSHong Zhang #define __FUNCT__ "PCDestroy_LU"
2299b54502bSHong Zhang static PetscErrorCode PCDestroy_LU(PC pc)
2309b54502bSHong Zhang {
2319b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
2329b54502bSHong Zhang   PetscErrorCode ierr;
2339b54502bSHong Zhang 
2349b54502bSHong Zhang   PetscFunctionBegin;
235075768bcSBarry Smith   if (!dir->inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(((PC_Factor*)dir)->fact);CHKERRQ(ierr);}
2369b54502bSHong Zhang   if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);}
2379b54502bSHong Zhang   if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);}
238075768bcSBarry Smith   ierr = PetscStrfree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
239075768bcSBarry Smith   ierr = PetscStrfree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr);
2409b54502bSHong Zhang   ierr = PetscFree(dir);CHKERRQ(ierr);
2419b54502bSHong Zhang   PetscFunctionReturn(0);
2429b54502bSHong Zhang }
2439b54502bSHong Zhang 
2449b54502bSHong Zhang #undef __FUNCT__
2459b54502bSHong Zhang #define __FUNCT__ "PCApply_LU"
2469b54502bSHong Zhang static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y)
2479b54502bSHong Zhang {
2489b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
2499b54502bSHong Zhang   PetscErrorCode ierr;
2509b54502bSHong Zhang 
2519b54502bSHong Zhang   PetscFunctionBegin;
2529b54502bSHong Zhang   if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);}
253075768bcSBarry Smith   else              {ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);}
2549b54502bSHong Zhang   PetscFunctionReturn(0);
2559b54502bSHong Zhang }
2569b54502bSHong Zhang 
2579b54502bSHong Zhang #undef __FUNCT__
2589b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_LU"
2599b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y)
2609b54502bSHong Zhang {
2619b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
2629b54502bSHong Zhang   PetscErrorCode ierr;
2639b54502bSHong Zhang 
2649b54502bSHong Zhang   PetscFunctionBegin;
2659b54502bSHong Zhang   if (dir->inplace) {ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);}
266075768bcSBarry Smith   else              {ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);}
2679b54502bSHong Zhang   PetscFunctionReturn(0);
2689b54502bSHong Zhang }
2699b54502bSHong Zhang 
2709b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/
2719b54502bSHong Zhang 
2729b54502bSHong Zhang EXTERN_C_BEGIN
2739b54502bSHong Zhang #undef __FUNCT__
2742401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace_LU"
2752401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_LU(PC pc)
2769b54502bSHong Zhang {
277c62fe420SBarry Smith   PC_LU *dir = (PC_LU*)pc->data;
2789b54502bSHong Zhang 
2799b54502bSHong Zhang   PetscFunctionBegin;
2809b54502bSHong Zhang   dir->inplace = PETSC_TRUE;
2819b54502bSHong Zhang   PetscFunctionReturn(0);
2829b54502bSHong Zhang }
2839b54502bSHong Zhang EXTERN_C_END
2849b54502bSHong Zhang 
2859b54502bSHong Zhang /* ------------------------------------------------------------------------ */
2869b54502bSHong Zhang 
2879b54502bSHong Zhang /*MC
2889b54502bSHong Zhang    PCLU - Uses a direct solver, based on LU factorization, as a preconditioner
2899b54502bSHong Zhang 
2909b54502bSHong Zhang    Options Database Keys:
2912401956bSBarry Smith +  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
292c7393fdbSBarry Smith .  -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like spooles
2932401956bSBarry Smith .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
29455ba2a51SBarry Smith .  -pc_factor_fill <fill> - Sets fill amount
2952401956bSBarry Smith .  -pc_factor_in_place - Activates in-place factorization
2962401956bSBarry Smith .  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
2972401956bSBarry Smith .  -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase
2989b54502bSHong Zhang                                          stability of factorization.
299f251bdbdSHong Zhang .  -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default
300e22d95b2SBarry Smith .  -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value
301f251bdbdSHong Zhang         is optional with PETSC_TRUE being the default
302e22d95b2SBarry Smith -   -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the
303e22d95b2SBarry Smith         diagonal.
3049b54502bSHong Zhang 
3059b54502bSHong Zhang    Notes: Not all options work for all matrix formats
3069b54502bSHong Zhang           Run with -help to see additional options for particular matrix formats or factorization
3079b54502bSHong Zhang           algorithms
3089b54502bSHong Zhang 
3099b54502bSHong Zhang    Level: beginner
3109b54502bSHong Zhang 
3119b54502bSHong Zhang    Concepts: LU factorization, direct solver
3129b54502bSHong Zhang 
3139b54502bSHong Zhang    Notes: Usually this will compute an "exact" solution in one iteration and does
3149b54502bSHong Zhang           not need a Krylov method (i.e. you can use -ksp_type preonly, or
3159b54502bSHong Zhang           KSPSetType(ksp,KSPPREONLY) for the Krylov method
3169b54502bSHong Zhang 
3179b54502bSHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
318a4fd02acSBarry Smith            PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(),
319e5a9bf91SBarry Smith            PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetPivoting(),
320e22d95b2SBarry Smith            PCFactorSetPivotingInBlocks(),PCFactorSetShiftNonzero(),PCFactorSetShiftPd(), PCFactorReorderForNonzeroDiagonal()
3219b54502bSHong Zhang M*/
3229b54502bSHong Zhang 
3239b54502bSHong Zhang EXTERN_C_BEGIN
3249b54502bSHong Zhang #undef __FUNCT__
3259b54502bSHong Zhang #define __FUNCT__ "PCCreate_LU"
326dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_LU(PC pc)
3279b54502bSHong Zhang {
3289b54502bSHong Zhang   PetscErrorCode ierr;
3299b54502bSHong Zhang   PetscMPIInt    size;
3309b54502bSHong Zhang   PC_LU          *dir;
3319b54502bSHong Zhang 
3329b54502bSHong Zhang   PetscFunctionBegin;
33338f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_LU,&dir);CHKERRQ(ierr);
3349b54502bSHong Zhang 
335075768bcSBarry Smith   ierr = MatFactorInfoInitialize(&((PC_Factor*)dir)->info);CHKERRQ(ierr);
336075768bcSBarry Smith   ((PC_Factor*)dir)->fact                  = 0;
3379b54502bSHong Zhang   dir->inplace               = PETSC_FALSE;
3389b54502bSHong Zhang   dir->nonzerosalongdiagonal = PETSC_FALSE;
3399b54502bSHong Zhang 
340075768bcSBarry Smith   ((PC_Factor*)dir)->info.fill           = 5.0;
341075768bcSBarry Smith   ((PC_Factor*)dir)->info.dtcol          = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */
342075768bcSBarry Smith   ((PC_Factor*)dir)->info.shiftnz        = 0.0;
343075768bcSBarry Smith   ((PC_Factor*)dir)->info.zeropivot      = 1.e-12;
344075768bcSBarry Smith   ((PC_Factor*)dir)->info.pivotinblocks  = 1.0;
345075768bcSBarry Smith   ((PC_Factor*)dir)->info.shiftpd        = 0.0; /* false */
3469b54502bSHong Zhang   dir->col                 = 0;
3479b54502bSHong Zhang   dir->row                 = 0;
3485c9eb25fSBarry Smith 
349075768bcSBarry Smith   ierr = PetscStrallocpy(MAT_SOLVER_PETSC,&((PC_Factor*)dir)->solvertype);CHKERRQ(ierr);
3507adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
3519b54502bSHong Zhang   if (size == 1) {
352075768bcSBarry Smith     ierr = PetscStrallocpy(MATORDERING_ND,&((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
3539b54502bSHong Zhang   } else {
354075768bcSBarry Smith     ierr = PetscStrallocpy(MATORDERING_NATURAL,&((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
3559b54502bSHong Zhang   }
3569b54502bSHong Zhang   dir->reusefill        = PETSC_FALSE;
3579b54502bSHong Zhang   dir->reuseordering    = PETSC_FALSE;
3589b54502bSHong Zhang   pc->data              = (void*)dir;
3599b54502bSHong Zhang 
3609b54502bSHong Zhang   pc->ops->destroy           = PCDestroy_LU;
3619b54502bSHong Zhang   pc->ops->apply             = PCApply_LU;
3629b54502bSHong Zhang   pc->ops->applytranspose    = PCApplyTranspose_LU;
3639b54502bSHong Zhang   pc->ops->setup             = PCSetUp_LU;
3649b54502bSHong Zhang   pc->ops->setfromoptions    = PCSetFromOptions_LU;
3659b54502bSHong Zhang   pc->ops->view              = PCView_LU;
3669b54502bSHong Zhang   pc->ops->applyrichardson   = 0;
36785317021SBarry Smith   pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor;
3689b54502bSHong Zhang 
369*7112b564SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor",
370*7112b564SBarry Smith                     PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr);
37185317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Factor",
37285317021SBarry Smith                     PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr);
37385317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor",
37485317021SBarry Smith                     PCFactorSetZeroPivot_Factor);CHKERRQ(ierr);
37585317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_Factor",
37685317021SBarry Smith                     PCFactorSetShiftNonzero_Factor);CHKERRQ(ierr);
37785317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_Factor",
37885317021SBarry Smith                     PCFactorSetShiftPd_Factor);CHKERRQ(ierr);
379afaefe49SHong Zhang 
38085317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor",
38185317021SBarry Smith                     PCFactorSetFill_Factor);CHKERRQ(ierr);
3822401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_LU",
3832401956bSBarry Smith                     PCFactorSetUseInPlace_LU);CHKERRQ(ierr);
38485317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor",
38585317021SBarry Smith                     PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr);
3862401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_LU",
3872401956bSBarry Smith                     PCFactorSetReuseOrdering_LU);CHKERRQ(ierr);
3882401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_LU",
3892401956bSBarry Smith                     PCFactorSetReuseFill_LU);CHKERRQ(ierr);
39085317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivoting_C","PCFactorSetPivoting_Factor",
39185317021SBarry Smith                     PCFactorSetPivoting_Factor);CHKERRQ(ierr);
39285317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_Factor",
39385317021SBarry Smith                     PCFactorSetPivotInBlocks_Factor);CHKERRQ(ierr);
3942401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_LU",
3952401956bSBarry Smith                     PCFactorReorderForNonzeroDiagonal_LU);CHKERRQ(ierr);
3969b54502bSHong Zhang   PetscFunctionReturn(0);
3979b54502bSHong Zhang }
3989b54502bSHong Zhang EXTERN_C_END
399