xref: /petsc/src/ksp/pc/impls/factor/lu/lu.c (revision 8ff23777322cfdc407588772804a6ce977f33084)
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 
97c4f633dSBarry 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;
61*8ff23777SHong Zhang   PetscTruth      flg = PETSC_FALSE;
629b54502bSHong Zhang   PetscReal       tol;
639b54502bSHong Zhang 
649b54502bSHong Zhang   PetscFunctionBegin;
659b54502bSHong Zhang   ierr = PetscOptionsHead("LU options");CHKERRQ(ierr);
66*8ff23777SHong Zhang     ierr = PCSetFromOptions_Factor(pc);CHKERRQ(ierr);
675c9eb25fSBarry Smith 
682401956bSBarry Smith     ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr);
699b54502bSHong Zhang     if (flg) {
709b54502bSHong Zhang       tol = PETSC_DECIDE;
712401956bSBarry Smith       ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr);
722401956bSBarry Smith       ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr);
739b54502bSHong Zhang     }
749b54502bSHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
759b54502bSHong Zhang   PetscFunctionReturn(0);
769b54502bSHong Zhang }
779b54502bSHong Zhang 
789b54502bSHong Zhang #undef __FUNCT__
799b54502bSHong Zhang #define __FUNCT__ "PCView_LU"
809b54502bSHong Zhang static PetscErrorCode PCView_LU(PC pc,PetscViewer viewer)
819b54502bSHong Zhang {
829b54502bSHong Zhang   PC_LU          *lu = (PC_LU*)pc->data;
839b54502bSHong Zhang   PetscErrorCode ierr;
849b54502bSHong Zhang   PetscTruth     iascii,isstring;
859b54502bSHong Zhang 
869b54502bSHong Zhang   PetscFunctionBegin;
879b54502bSHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
889b54502bSHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
899b54502bSHong Zhang   if (iascii) {
909b54502bSHong Zhang 
919b54502bSHong Zhang     if (lu->inplace) {ierr = PetscViewerASCIIPrintf(viewer,"  LU: in-place factorization\n");CHKERRQ(ierr);}
929b54502bSHong Zhang     else             {ierr = PetscViewerASCIIPrintf(viewer,"  LU: out-of-place factorization\n");CHKERRQ(ierr);}
93075768bcSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    matrix ordering: %s\n",((PC_Factor*)lu)->ordering);CHKERRQ(ierr);
94075768bcSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  LU: tolerance for zero pivot %G\n",((PC_Factor*)lu)->info.zeropivot);CHKERRQ(ierr);
95075768bcSBarry Smith     if (((PC_Factor*)lu)->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer,"  LU: using Manteuffel shift\n");CHKERRQ(ierr);}
96075768bcSBarry Smith     if (((PC_Factor*)lu)->fact) {
97f3a39becSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  LU: factor fill ratio needed %G\n",lu->actualfill);CHKERRQ(ierr);
98f3a39becSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"       Factored matrix follows\n");CHKERRQ(ierr);
99f3a39becSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
100f3a39becSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
101f3a39becSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
102f3a39becSBarry Smith       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
103075768bcSBarry Smith       ierr = MatView(((PC_Factor*)lu)->fact,viewer);CHKERRQ(ierr);
1049b54502bSHong Zhang       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
105f3a39becSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
106f3a39becSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
107f3a39becSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1089b54502bSHong Zhang     }
1099b54502bSHong Zhang     if (lu->reusefill)    {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing fill from past factorization\n");CHKERRQ(ierr);}
1109b54502bSHong Zhang     if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing reordering from past factorization\n");CHKERRQ(ierr);}
1119b54502bSHong Zhang   } else if (isstring) {
112075768bcSBarry Smith     ierr = PetscViewerStringSPrintf(viewer," order=%s",((PC_Factor*)lu)->ordering);CHKERRQ(ierr);CHKERRQ(ierr);
1139b54502bSHong Zhang   } else {
1149b54502bSHong Zhang     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCLU",((PetscObject)viewer)->type_name);
1159b54502bSHong Zhang   }
1169b54502bSHong Zhang   PetscFunctionReturn(0);
1179b54502bSHong Zhang }
1189b54502bSHong Zhang 
1199b54502bSHong Zhang #undef __FUNCT__
1209b54502bSHong Zhang #define __FUNCT__ "PCSetUp_LU"
1219b54502bSHong Zhang static PetscErrorCode PCSetUp_LU(PC pc)
1229b54502bSHong Zhang {
1239b54502bSHong Zhang   PetscErrorCode ierr;
1249b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
1259b54502bSHong Zhang 
1269b54502bSHong Zhang   PetscFunctionBegin;
127075768bcSBarry Smith   if (dir->reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->actualfill;
1289b54502bSHong Zhang 
1299b54502bSHong Zhang   if (dir->inplace) {
1309b54502bSHong Zhang     if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);}
1319b54502bSHong Zhang     if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);}
132075768bcSBarry Smith     ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
13303c60df9SBarry Smith     if (dir->row) {
13403c60df9SBarry Smith       ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);
13503c60df9SBarry Smith       ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr);
13603c60df9SBarry Smith     }
137075768bcSBarry Smith     ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
138075768bcSBarry Smith     ((PC_Factor*)dir)->fact = pc->pmat;
1399b54502bSHong Zhang   } else {
1409b54502bSHong Zhang     MatInfo info;
1419b54502bSHong Zhang     if (!pc->setupcalled) {
142075768bcSBarry Smith       ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
1439b54502bSHong Zhang       if (dir->nonzerosalongdiagonal) {
1449b54502bSHong Zhang         ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr);
1459b54502bSHong Zhang       }
14603c60df9SBarry Smith       if (dir->row) {
14703c60df9SBarry Smith         ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);
14803c60df9SBarry Smith         ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr);
14903c60df9SBarry Smith       }
150075768bcSBarry Smith       ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr);
151075768bcSBarry Smith       ierr = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
152075768bcSBarry Smith       ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
1539b54502bSHong Zhang       dir->actualfill = info.fill_ratio_needed;
154075768bcSBarry Smith       ierr = PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);CHKERRQ(ierr);
1559b54502bSHong Zhang     } else if (pc->flag != SAME_NONZERO_PATTERN) {
1569b54502bSHong Zhang       if (!dir->reuseordering) {
1579b54502bSHong Zhang         if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);}
1589b54502bSHong Zhang         if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);}
159075768bcSBarry Smith         ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
1609b54502bSHong Zhang         if (dir->nonzerosalongdiagonal) {
1619b54502bSHong Zhang           ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr);
1629b54502bSHong Zhang         }
16303c60df9SBarry Smith         if (dir->row) {
16403c60df9SBarry Smith           ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);
16503c60df9SBarry Smith           ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr);
16603c60df9SBarry Smith         }
1679b54502bSHong Zhang       }
168075768bcSBarry Smith       ierr = MatDestroy(((PC_Factor*)dir)->fact);CHKERRQ(ierr);
169075768bcSBarry Smith       ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr);
170075768bcSBarry Smith       ierr = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
171075768bcSBarry Smith       ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
1729b54502bSHong Zhang       dir->actualfill = info.fill_ratio_needed;
173075768bcSBarry Smith       ierr = PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);CHKERRQ(ierr);
1749b54502bSHong Zhang     }
175075768bcSBarry Smith     ierr = MatLUFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
1769b54502bSHong Zhang   }
1779b54502bSHong Zhang   PetscFunctionReturn(0);
1789b54502bSHong Zhang }
1799b54502bSHong Zhang 
1809b54502bSHong Zhang #undef __FUNCT__
1819b54502bSHong Zhang #define __FUNCT__ "PCDestroy_LU"
1829b54502bSHong Zhang static PetscErrorCode PCDestroy_LU(PC pc)
1839b54502bSHong Zhang {
1849b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
1859b54502bSHong Zhang   PetscErrorCode ierr;
1869b54502bSHong Zhang 
1879b54502bSHong Zhang   PetscFunctionBegin;
188075768bcSBarry Smith   if (!dir->inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(((PC_Factor*)dir)->fact);CHKERRQ(ierr);}
1899b54502bSHong Zhang   if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);}
1909b54502bSHong Zhang   if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);}
191075768bcSBarry Smith   ierr = PetscStrfree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
192075768bcSBarry Smith   ierr = PetscStrfree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr);
1939b54502bSHong Zhang   ierr = PetscFree(dir);CHKERRQ(ierr);
1949b54502bSHong Zhang   PetscFunctionReturn(0);
1959b54502bSHong Zhang }
1969b54502bSHong Zhang 
1979b54502bSHong Zhang #undef __FUNCT__
1989b54502bSHong Zhang #define __FUNCT__ "PCApply_LU"
1999b54502bSHong Zhang static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y)
2009b54502bSHong Zhang {
2019b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
2029b54502bSHong Zhang   PetscErrorCode ierr;
2039b54502bSHong Zhang 
2049b54502bSHong Zhang   PetscFunctionBegin;
2059b54502bSHong Zhang   if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);}
206075768bcSBarry Smith   else              {ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);}
2079b54502bSHong Zhang   PetscFunctionReturn(0);
2089b54502bSHong Zhang }
2099b54502bSHong Zhang 
2109b54502bSHong Zhang #undef __FUNCT__
2119b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_LU"
2129b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y)
2139b54502bSHong Zhang {
2149b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
2159b54502bSHong Zhang   PetscErrorCode ierr;
2169b54502bSHong Zhang 
2179b54502bSHong Zhang   PetscFunctionBegin;
2189b54502bSHong Zhang   if (dir->inplace) {ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);}
219075768bcSBarry Smith   else              {ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);}
2209b54502bSHong Zhang   PetscFunctionReturn(0);
2219b54502bSHong Zhang }
2229b54502bSHong Zhang 
2239b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/
2249b54502bSHong Zhang 
2259b54502bSHong Zhang EXTERN_C_BEGIN
2269b54502bSHong Zhang #undef __FUNCT__
2272401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace_LU"
2282401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_LU(PC pc)
2299b54502bSHong Zhang {
230c62fe420SBarry Smith   PC_LU *dir = (PC_LU*)pc->data;
2319b54502bSHong Zhang 
2329b54502bSHong Zhang   PetscFunctionBegin;
2339b54502bSHong Zhang   dir->inplace = PETSC_TRUE;
2349b54502bSHong Zhang   PetscFunctionReturn(0);
2359b54502bSHong Zhang }
2369b54502bSHong Zhang EXTERN_C_END
2379b54502bSHong Zhang 
2389b54502bSHong Zhang /* ------------------------------------------------------------------------ */
2399b54502bSHong Zhang 
2409b54502bSHong Zhang /*MC
2419b54502bSHong Zhang    PCLU - Uses a direct solver, based on LU factorization, as a preconditioner
2429b54502bSHong Zhang 
2439b54502bSHong Zhang    Options Database Keys:
2442401956bSBarry Smith +  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
245c7393fdbSBarry Smith .  -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like spooles
2462401956bSBarry Smith .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
24755ba2a51SBarry Smith .  -pc_factor_fill <fill> - Sets fill amount
2482401956bSBarry Smith .  -pc_factor_in_place - Activates in-place factorization
2492401956bSBarry Smith .  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
2502401956bSBarry Smith .  -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase
2519b54502bSHong Zhang                                          stability of factorization.
252f251bdbdSHong Zhang .  -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default
253e22d95b2SBarry Smith .  -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value
254f251bdbdSHong Zhang         is optional with PETSC_TRUE being the default
255e22d95b2SBarry Smith -   -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the
256e22d95b2SBarry Smith         diagonal.
2579b54502bSHong Zhang 
2589b54502bSHong Zhang    Notes: Not all options work for all matrix formats
2599b54502bSHong Zhang           Run with -help to see additional options for particular matrix formats or factorization
2609b54502bSHong Zhang           algorithms
2619b54502bSHong Zhang 
2629b54502bSHong Zhang    Level: beginner
2639b54502bSHong Zhang 
2649b54502bSHong Zhang    Concepts: LU factorization, direct solver
2659b54502bSHong Zhang 
2669b54502bSHong Zhang    Notes: Usually this will compute an "exact" solution in one iteration and does
2679b54502bSHong Zhang           not need a Krylov method (i.e. you can use -ksp_type preonly, or
2689b54502bSHong Zhang           KSPSetType(ksp,KSPPREONLY) for the Krylov method
2699b54502bSHong Zhang 
2709b54502bSHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
271a4fd02acSBarry Smith            PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(),
272*8ff23777SHong Zhang            PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetColumnPivot(),
273*8ff23777SHong Zhang            PCFactorSetPivotingInBlocks(),PCFactorSetShiftNonzero(),PCFactorSetShiftPd(),PCFactorSetShiftInBlocks()
274*8ff23777SHong Zhang            PCFactorReorderForNonzeroDiagonal()
2759b54502bSHong Zhang M*/
2769b54502bSHong Zhang 
2779b54502bSHong Zhang EXTERN_C_BEGIN
2789b54502bSHong Zhang #undef __FUNCT__
2799b54502bSHong Zhang #define __FUNCT__ "PCCreate_LU"
280dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_LU(PC pc)
2819b54502bSHong Zhang {
2829b54502bSHong Zhang   PetscErrorCode ierr;
2839b54502bSHong Zhang   PetscMPIInt    size;
2849b54502bSHong Zhang   PC_LU          *dir;
2859b54502bSHong Zhang 
2869b54502bSHong Zhang   PetscFunctionBegin;
28738f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_LU,&dir);CHKERRQ(ierr);
2889b54502bSHong Zhang 
289075768bcSBarry Smith   ierr = MatFactorInfoInitialize(&((PC_Factor*)dir)->info);CHKERRQ(ierr);
290075768bcSBarry Smith   ((PC_Factor*)dir)->fact                  = 0;
2919b54502bSHong Zhang   dir->inplace               = PETSC_FALSE;
2929b54502bSHong Zhang   dir->nonzerosalongdiagonal = PETSC_FALSE;
2939b54502bSHong Zhang 
294075768bcSBarry Smith   ((PC_Factor*)dir)->info.fill           = 5.0;
295075768bcSBarry Smith   ((PC_Factor*)dir)->info.dtcol          = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */
296075768bcSBarry Smith   ((PC_Factor*)dir)->info.shiftnz        = 0.0;
297075768bcSBarry Smith   ((PC_Factor*)dir)->info.zeropivot      = 1.e-12;
298075768bcSBarry Smith   ((PC_Factor*)dir)->info.pivotinblocks  = 1.0;
299075768bcSBarry Smith   ((PC_Factor*)dir)->info.shiftpd        = 0.0; /* false */
3009b54502bSHong Zhang   dir->col                 = 0;
3019b54502bSHong Zhang   dir->row                 = 0;
3025c9eb25fSBarry Smith 
303075768bcSBarry Smith   ierr = PetscStrallocpy(MAT_SOLVER_PETSC,&((PC_Factor*)dir)->solvertype);CHKERRQ(ierr);
3047adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
3059b54502bSHong Zhang   if (size == 1) {
306075768bcSBarry Smith     ierr = PetscStrallocpy(MATORDERING_ND,&((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
3079b54502bSHong Zhang   } else {
308075768bcSBarry Smith     ierr = PetscStrallocpy(MATORDERING_NATURAL,&((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
3099b54502bSHong Zhang   }
3109b54502bSHong Zhang   dir->reusefill        = PETSC_FALSE;
3119b54502bSHong Zhang   dir->reuseordering    = PETSC_FALSE;
3129b54502bSHong Zhang   pc->data              = (void*)dir;
3139b54502bSHong Zhang 
3149b54502bSHong Zhang   pc->ops->destroy           = PCDestroy_LU;
3159b54502bSHong Zhang   pc->ops->apply             = PCApply_LU;
3169b54502bSHong Zhang   pc->ops->applytranspose    = PCApplyTranspose_LU;
3179b54502bSHong Zhang   pc->ops->setup             = PCSetUp_LU;
3189b54502bSHong Zhang   pc->ops->setfromoptions    = PCSetFromOptions_LU;
3199b54502bSHong Zhang   pc->ops->view              = PCView_LU;
3209b54502bSHong Zhang   pc->ops->applyrichardson   = 0;
32185317021SBarry Smith   pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor;
3229b54502bSHong Zhang 
3237112b564SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor",
3247112b564SBarry Smith                     PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr);
32585317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Factor",
32685317021SBarry Smith                     PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr);
32785317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor",
32885317021SBarry Smith                     PCFactorSetZeroPivot_Factor);CHKERRQ(ierr);
32985317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_Factor",
33085317021SBarry Smith                     PCFactorSetShiftNonzero_Factor);CHKERRQ(ierr);
33185317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_Factor",
33285317021SBarry Smith                     PCFactorSetShiftPd_Factor);CHKERRQ(ierr);
333*8ff23777SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftInBlocks_C","PCFactorSetShiftInBlocks_Factor",
334*8ff23777SHong Zhang                     PCFactorSetShiftInBlocks_Factor);CHKERRQ(ierr);
335afaefe49SHong Zhang 
33685317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor",
33785317021SBarry Smith                     PCFactorSetFill_Factor);CHKERRQ(ierr);
3382401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_LU",
3392401956bSBarry Smith                     PCFactorSetUseInPlace_LU);CHKERRQ(ierr);
34085317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor",
34185317021SBarry Smith                     PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr);
3422401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_LU",
3432401956bSBarry Smith                     PCFactorSetReuseOrdering_LU);CHKERRQ(ierr);
3442401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_LU",
3452401956bSBarry Smith                     PCFactorSetReuseFill_LU);CHKERRQ(ierr);
346*8ff23777SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetColumnPivot_C","PCFactorSetColumnPivot_Factor",
347*8ff23777SHong Zhang                     PCFactorSetColumnPivot_Factor);CHKERRQ(ierr);
34885317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_Factor",
34985317021SBarry Smith                     PCFactorSetPivotInBlocks_Factor);CHKERRQ(ierr);
3502401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_LU",
3512401956bSBarry Smith                     PCFactorReorderForNonzeroDiagonal_LU);CHKERRQ(ierr);
3529b54502bSHong Zhang   PetscFunctionReturn(0);
3539b54502bSHong Zhang }
3549b54502bSHong Zhang EXTERN_C_END
355