xref: /petsc/src/ksp/pc/impls/factor/lu/lu.c (revision f8260c8f9a86972e69213fdea9e81a0caf515f42)
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"
147087cfbeSBarry Smith PetscErrorCode  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"
327087cfbeSBarry Smith PetscErrorCode  PCFactorSetReuseOrdering_LU(PC pc,PetscBool  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"
457087cfbeSBarry Smith PetscErrorCode  PCFactorSetReuseFill_LU(PC pc,PetscBool  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;
61ace3abfcSBarry Smith   PetscBool       flg = PETSC_FALSE;
629b54502bSHong Zhang   PetscReal       tol;
639b54502bSHong Zhang 
649b54502bSHong Zhang   PetscFunctionBegin;
659b54502bSHong Zhang   ierr = PetscOptionsHead("LU options");CHKERRQ(ierr);
668ff23777SHong 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;
84ace3abfcSBarry Smith   PetscBool      iascii;
859b54502bSHong Zhang 
869b54502bSHong Zhang   PetscFunctionBegin;
872692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
889b54502bSHong Zhang   if (iascii) {
89914a5d51SHong Zhang     if (lu->inplace) {
90914a5d51SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  LU: in-place factorization\n");CHKERRQ(ierr);
91914a5d51SHong Zhang     } else {
92914a5d51SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  LU: out-of-place factorization\n");CHKERRQ(ierr);
93145b9266SHong Zhang     }
94145b9266SHong Zhang 
959b54502bSHong Zhang     if (lu->reusefill)    {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing fill from past factorization\n");CHKERRQ(ierr);}
969b54502bSHong Zhang     if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing reordering from past factorization\n");CHKERRQ(ierr);}
979b54502bSHong Zhang   }
98914a5d51SHong Zhang   ierr = PCView_Factor(pc,viewer);CHKERRQ(ierr);
999b54502bSHong Zhang   PetscFunctionReturn(0);
1009b54502bSHong Zhang }
1019b54502bSHong Zhang 
1029b54502bSHong Zhang #undef __FUNCT__
1039b54502bSHong Zhang #define __FUNCT__ "PCSetUp_LU"
1049b54502bSHong Zhang static PetscErrorCode PCSetUp_LU(PC pc)
1059b54502bSHong Zhang {
1069b54502bSHong Zhang   PetscErrorCode ierr;
1079b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
1089b54502bSHong Zhang 
1099b54502bSHong Zhang   PetscFunctionBegin;
110075768bcSBarry Smith   if (dir->reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->actualfill;
1119b54502bSHong Zhang 
1129b54502bSHong Zhang   if (dir->inplace) {
1139b54502bSHong Zhang     if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);}
1149b54502bSHong Zhang     if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);}
115075768bcSBarry Smith     ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
11603c60df9SBarry Smith     if (dir->row) {
11703c60df9SBarry Smith       ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);
11803c60df9SBarry Smith       ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr);
11903c60df9SBarry Smith     }
120075768bcSBarry Smith     ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
121075768bcSBarry Smith     ((PC_Factor*)dir)->fact = pc->pmat;
1229b54502bSHong Zhang   } else {
1239b54502bSHong Zhang     MatInfo info;
1249b54502bSHong Zhang     if (!pc->setupcalled) {
125075768bcSBarry Smith       ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
1269b54502bSHong Zhang       if (dir->nonzerosalongdiagonal) {
1279b54502bSHong Zhang         ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr);
1289b54502bSHong Zhang       }
12903c60df9SBarry Smith       if (dir->row) {
13003c60df9SBarry Smith         ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);
13103c60df9SBarry Smith         ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr);
13203c60df9SBarry Smith       }
133d09a07f4SBarry Smith       if (!((PC_Factor*)dir)->fact){
134075768bcSBarry Smith         ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr);
135a1f19f5aSHong Zhang       }
136075768bcSBarry Smith       ierr = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
137075768bcSBarry Smith       ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
1389b54502bSHong Zhang       dir->actualfill = info.fill_ratio_needed;
139075768bcSBarry Smith       ierr = PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);CHKERRQ(ierr);
1409b54502bSHong Zhang     } else if (pc->flag != SAME_NONZERO_PATTERN) {
1419b54502bSHong Zhang       if (!dir->reuseordering) {
1429b54502bSHong Zhang         if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);}
1439b54502bSHong Zhang         if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);}
144075768bcSBarry Smith         ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
1459b54502bSHong Zhang         if (dir->nonzerosalongdiagonal) {
1469b54502bSHong Zhang           ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr);
1479b54502bSHong Zhang         }
14803c60df9SBarry Smith         if (dir->row) {
14903c60df9SBarry Smith           ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);
15003c60df9SBarry Smith           ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr);
15103c60df9SBarry Smith         }
1529b54502bSHong Zhang       }
153075768bcSBarry Smith       ierr = MatDestroy(((PC_Factor*)dir)->fact);CHKERRQ(ierr);
154075768bcSBarry Smith       ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr);
155075768bcSBarry Smith       ierr = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
156075768bcSBarry Smith       ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
1579b54502bSHong Zhang       dir->actualfill = info.fill_ratio_needed;
158075768bcSBarry Smith       ierr = PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);CHKERRQ(ierr);
1599b54502bSHong Zhang     }
160075768bcSBarry Smith     ierr = MatLUFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
1619b54502bSHong Zhang   }
1629b54502bSHong Zhang   PetscFunctionReturn(0);
1639b54502bSHong Zhang }
1649b54502bSHong Zhang 
1659b54502bSHong Zhang #undef __FUNCT__
1669b54502bSHong Zhang #define __FUNCT__ "PCDestroy_LU"
1679b54502bSHong Zhang static PetscErrorCode PCDestroy_LU(PC pc)
1689b54502bSHong Zhang {
1699b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
1709b54502bSHong Zhang   PetscErrorCode ierr;
1719b54502bSHong Zhang 
1729b54502bSHong Zhang   PetscFunctionBegin;
173075768bcSBarry Smith   if (!dir->inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(((PC_Factor*)dir)->fact);CHKERRQ(ierr);}
1749b54502bSHong Zhang   if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);}
1759b54502bSHong Zhang   if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);}
176503cfb0cSBarry Smith   ierr = PetscFree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
177503cfb0cSBarry Smith   ierr = PetscFree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr);
1789b54502bSHong Zhang   ierr = PetscFree(dir);CHKERRQ(ierr);
1799b54502bSHong Zhang   PetscFunctionReturn(0);
1809b54502bSHong Zhang }
1819b54502bSHong Zhang 
1829b54502bSHong Zhang #undef __FUNCT__
1839b54502bSHong Zhang #define __FUNCT__ "PCApply_LU"
1849b54502bSHong Zhang static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y)
1859b54502bSHong Zhang {
1869b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
1879b54502bSHong Zhang   PetscErrorCode ierr;
1889b54502bSHong Zhang 
1899b54502bSHong Zhang   PetscFunctionBegin;
1909b54502bSHong Zhang   if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);}
191075768bcSBarry Smith   else              {ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);}
1929b54502bSHong Zhang   PetscFunctionReturn(0);
1939b54502bSHong Zhang }
1949b54502bSHong Zhang 
1959b54502bSHong Zhang #undef __FUNCT__
1969b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_LU"
1979b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y)
1989b54502bSHong Zhang {
1999b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
2009b54502bSHong Zhang   PetscErrorCode ierr;
2019b54502bSHong Zhang 
2029b54502bSHong Zhang   PetscFunctionBegin;
2039b54502bSHong Zhang   if (dir->inplace) {ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);}
204075768bcSBarry Smith   else              {ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);}
2059b54502bSHong Zhang   PetscFunctionReturn(0);
2069b54502bSHong Zhang }
2079b54502bSHong Zhang 
2089b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/
2099b54502bSHong Zhang 
2109b54502bSHong Zhang EXTERN_C_BEGIN
2119b54502bSHong Zhang #undef __FUNCT__
2122401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace_LU"
2137087cfbeSBarry Smith PetscErrorCode  PCFactorSetUseInPlace_LU(PC pc)
2149b54502bSHong Zhang {
215c62fe420SBarry Smith   PC_LU *dir = (PC_LU*)pc->data;
2169b54502bSHong Zhang 
2179b54502bSHong Zhang   PetscFunctionBegin;
2189b54502bSHong Zhang   dir->inplace = PETSC_TRUE;
2199b54502bSHong Zhang   PetscFunctionReturn(0);
2209b54502bSHong Zhang }
2219b54502bSHong Zhang EXTERN_C_END
2229b54502bSHong Zhang 
2239b54502bSHong Zhang /* ------------------------------------------------------------------------ */
2249b54502bSHong Zhang 
2259b54502bSHong Zhang /*MC
2269b54502bSHong Zhang    PCLU - Uses a direct solver, based on LU factorization, as a preconditioner
2279b54502bSHong Zhang 
2289b54502bSHong Zhang    Options Database Keys:
2292401956bSBarry Smith +  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
230c7393fdbSBarry Smith .  -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like spooles
2312401956bSBarry Smith .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
23255ba2a51SBarry Smith .  -pc_factor_fill <fill> - Sets fill amount
2332401956bSBarry Smith .  -pc_factor_in_place - Activates in-place factorization
2342401956bSBarry Smith .  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
2352401956bSBarry Smith .  -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase
2369b54502bSHong Zhang                                          stability of factorization.
237145b9266SHong Zhang .  -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types
238145b9266SHong Zhang .  -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default
239e22d95b2SBarry Smith -   -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the
240e22d95b2SBarry Smith         diagonal.
2419b54502bSHong Zhang 
2429b54502bSHong Zhang    Notes: Not all options work for all matrix formats
2439b54502bSHong Zhang           Run with -help to see additional options for particular matrix formats or factorization
2449b54502bSHong Zhang           algorithms
2459b54502bSHong Zhang 
2469b54502bSHong Zhang    Level: beginner
2479b54502bSHong Zhang 
2489b54502bSHong Zhang    Concepts: LU factorization, direct solver
2499b54502bSHong Zhang 
2509b54502bSHong Zhang    Notes: Usually this will compute an "exact" solution in one iteration and does
2519b54502bSHong Zhang           not need a Krylov method (i.e. you can use -ksp_type preonly, or
2529b54502bSHong Zhang           KSPSetType(ksp,KSPPREONLY) for the Krylov method
2539b54502bSHong Zhang 
2549b54502bSHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
255a4fd02acSBarry Smith            PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(),
2568ff23777SHong Zhang            PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetColumnPivot(),
257145b9266SHong Zhang            PCFactorSetPivotingInBlocks(),PCFactorSetShiftType(),PCFactorSetShiftAmount()
2588ff23777SHong Zhang            PCFactorReorderForNonzeroDiagonal()
2599b54502bSHong Zhang M*/
2609b54502bSHong Zhang 
2619b54502bSHong Zhang EXTERN_C_BEGIN
2629b54502bSHong Zhang #undef __FUNCT__
2639b54502bSHong Zhang #define __FUNCT__ "PCCreate_LU"
2647087cfbeSBarry Smith PetscErrorCode  PCCreate_LU(PC pc)
2659b54502bSHong Zhang {
2669b54502bSHong Zhang   PetscErrorCode ierr;
2679b54502bSHong Zhang   PetscMPIInt    size;
2689b54502bSHong Zhang   PC_LU          *dir;
2699b54502bSHong Zhang 
2709b54502bSHong Zhang   PetscFunctionBegin;
27138f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_LU,&dir);CHKERRQ(ierr);
2729b54502bSHong Zhang 
273075768bcSBarry Smith   ierr = MatFactorInfoInitialize(&((PC_Factor*)dir)->info);CHKERRQ(ierr);
274879e8a4dSBarry Smith   ((PC_Factor*)dir)->fact       = PETSC_NULL;
275d5f3da31SBarry Smith   ((PC_Factor*)dir)->factortype = MAT_FACTOR_LU;
2769b54502bSHong Zhang   dir->inplace                  = PETSC_FALSE;
2779b54502bSHong Zhang   dir->nonzerosalongdiagonal    = PETSC_FALSE;
2789b54502bSHong Zhang 
279075768bcSBarry Smith   ((PC_Factor*)dir)->info.fill           = 5.0;
280075768bcSBarry Smith   ((PC_Factor*)dir)->info.dtcol          = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */
281f4db908eSBarry Smith   ((PC_Factor*)dir)->info.shifttype      = (PetscReal) MAT_SHIFT_NONE;
282915743fcSHong Zhang   ((PC_Factor*)dir)->info.shiftamount    = 0.0;
283075768bcSBarry Smith   ((PC_Factor*)dir)->info.zeropivot      = 1.e-12;
284075768bcSBarry Smith   ((PC_Factor*)dir)->info.pivotinblocks  = 1.0;
2859b54502bSHong Zhang   dir->col                 = 0;
2869b54502bSHong Zhang   dir->row                 = 0;
2875c9eb25fSBarry Smith 
288a1f19f5aSHong Zhang   ierr = PetscStrallocpy(MATSOLVERPETSC,&((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); /* default SolverPackage */
2897adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
2909b54502bSHong Zhang   if (size == 1) {
2912692d6eeSBarry Smith     ierr = PetscStrallocpy(MATORDERINGND,&((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
2929b54502bSHong Zhang   } else {
2932692d6eeSBarry Smith     ierr = PetscStrallocpy(MATORDERINGNATURAL,&((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
2949b54502bSHong Zhang   }
2959b54502bSHong Zhang   dir->reusefill        = PETSC_FALSE;
2969b54502bSHong Zhang   dir->reuseordering    = PETSC_FALSE;
2979b54502bSHong Zhang   pc->data              = (void*)dir;
2989b54502bSHong Zhang 
2999b54502bSHong Zhang   pc->ops->destroy           = PCDestroy_LU;
3009b54502bSHong Zhang   pc->ops->apply             = PCApply_LU;
3019b54502bSHong Zhang   pc->ops->applytranspose    = PCApplyTranspose_LU;
3029b54502bSHong Zhang   pc->ops->setup             = PCSetUp_LU;
3039b54502bSHong Zhang   pc->ops->setfromoptions    = PCSetFromOptions_LU;
3049b54502bSHong Zhang   pc->ops->view              = PCView_LU;
3059b54502bSHong Zhang   pc->ops->applyrichardson   = 0;
30685317021SBarry Smith   pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor;
3079b54502bSHong Zhang 
308*f8260c8fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUpMatSolverPackage_C","PCFactorSetUpMatSolverPackage_Factor",
309*f8260c8fSBarry Smith                     PCFactorSetUpMatSolverPackage_Factor);CHKERRQ(ierr);
3107112b564SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor",
3117112b564SBarry Smith                     PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr);
31285317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Factor",
31385317021SBarry Smith                     PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr);
31485317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor",
31585317021SBarry Smith                     PCFactorSetZeroPivot_Factor);CHKERRQ(ierr);
316d90ac83dSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftType_C","PCFactorSetShiftType_Factor",
317d90ac83dSHong Zhang                     PCFactorSetShiftType_Factor);CHKERRQ(ierr);
318d90ac83dSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftAmount_C","PCFactorSetShiftAmount_Factor",
319d90ac83dSHong Zhang                     PCFactorSetShiftAmount_Factor);CHKERRQ(ierr);
32085317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor",
32185317021SBarry Smith                     PCFactorSetFill_Factor);CHKERRQ(ierr);
3222401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_LU",
3232401956bSBarry Smith                     PCFactorSetUseInPlace_LU);CHKERRQ(ierr);
32485317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor",
32585317021SBarry Smith                     PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr);
3262401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_LU",
3272401956bSBarry Smith                     PCFactorSetReuseOrdering_LU);CHKERRQ(ierr);
3282401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_LU",
3292401956bSBarry Smith                     PCFactorSetReuseFill_LU);CHKERRQ(ierr);
3308ff23777SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetColumnPivot_C","PCFactorSetColumnPivot_Factor",
3318ff23777SHong Zhang                     PCFactorSetColumnPivot_Factor);CHKERRQ(ierr);
33285317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_Factor",
33385317021SBarry Smith                     PCFactorSetPivotInBlocks_Factor);CHKERRQ(ierr);
3342401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_LU",
3352401956bSBarry Smith                     PCFactorReorderForNonzeroDiagonal_LU);CHKERRQ(ierr);
3369b54502bSHong Zhang   PetscFunctionReturn(0);
3379b54502bSHong Zhang }
3389b54502bSHong Zhang EXTERN_C_END
339