xref: /petsc/src/ksp/pc/impls/factor/lu/lu.c (revision 3d1c1ea0033fce25879716d8a26212643c58f132)
1dba47a55SKris Buschelman 
29b54502bSHong Zhang /*
39b54502bSHong Zhang    Defines a direct factorization preconditioner for any Mat implementation
49b54502bSHong Zhang    Note: this need not be consided a preconditioner since it supplies
59b54502bSHong Zhang          a direct solver.
69b54502bSHong Zhang */
7ee45ca4aSHong Zhang 
8c6db04a5SJed Brown #include <../src/ksp/pc/impls/factor/lu/lu.h>  /*I "petscpc.h" I*/
99b54502bSHong Zhang 
10680c5173SHong Zhang 
119b54502bSHong Zhang #undef __FUNCT__
122401956bSBarry Smith #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal_LU"
137087cfbeSBarry Smith PetscErrorCode PCFactorReorderForNonzeroDiagonal_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;
192fa5cd67SKarl Rupp   if (z == PETSC_DECIDE) lu->nonzerosalongdiagonaltol = 1.e-10;
202fa5cd67SKarl Rupp   else lu->nonzerosalongdiagonaltol = z;
219b54502bSHong Zhang   PetscFunctionReturn(0);
229b54502bSHong Zhang }
239b54502bSHong Zhang 
249b54502bSHong Zhang #undef __FUNCT__
259b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_LU"
264416b707SBarry Smith static PetscErrorCode PCSetFromOptions_LU(PetscOptionItems *PetscOptionsObject,PC pc)
279b54502bSHong Zhang {
289b54502bSHong Zhang   PC_LU          *lu = (PC_LU*)pc->data;
299b54502bSHong Zhang   PetscErrorCode ierr;
30ace3abfcSBarry Smith   PetscBool      flg = PETSC_FALSE;
319b54502bSHong Zhang   PetscReal      tol;
329b54502bSHong Zhang 
339b54502bSHong Zhang   PetscFunctionBegin;
34e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"LU options");CHKERRQ(ierr);
35e55864a3SBarry Smith   ierr = PCSetFromOptions_Factor(PetscOptionsObject,pc);CHKERRQ(ierr);
365c9eb25fSBarry Smith 
372401956bSBarry Smith   ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr);
389b54502bSHong Zhang   if (flg) {
399b54502bSHong Zhang     tol  = PETSC_DECIDE;
402401956bSBarry Smith     ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr);
412401956bSBarry Smith     ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr);
429b54502bSHong Zhang   }
439b54502bSHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
449b54502bSHong Zhang   PetscFunctionReturn(0);
459b54502bSHong Zhang }
469b54502bSHong Zhang 
479b54502bSHong Zhang #undef __FUNCT__
489b54502bSHong Zhang #define __FUNCT__ "PCView_LU"
499b54502bSHong Zhang static PetscErrorCode PCView_LU(PC pc,PetscViewer viewer)
509b54502bSHong Zhang {
519b54502bSHong Zhang   PetscErrorCode ierr;
529b54502bSHong Zhang 
539b54502bSHong Zhang   PetscFunctionBegin;
54914a5d51SHong Zhang   ierr = PCView_Factor(pc,viewer);CHKERRQ(ierr);
559b54502bSHong Zhang   PetscFunctionReturn(0);
569b54502bSHong Zhang }
579b54502bSHong Zhang 
589b54502bSHong Zhang #undef __FUNCT__
599b54502bSHong Zhang #define __FUNCT__ "PCSetUp_LU"
609b54502bSHong Zhang static PetscErrorCode PCSetUp_LU(PC pc)
619b54502bSHong Zhang {
629b54502bSHong Zhang   PetscErrorCode         ierr;
639b54502bSHong Zhang   PC_LU                  *dir = (PC_LU*)pc->data;
6400c67f3bSHong Zhang   const MatSolverPackage stype;
6500e125f8SBarry Smith   MatFactorError         err;
66*3d1c1ea0SBarry Smith 
679b54502bSHong Zhang   PetscFunctionBegin;
68*3d1c1ea0SBarry Smith   if (dir->hdr.reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->hdr.actualfill;
699b54502bSHong Zhang 
7084d44b13SHong Zhang   ierr = MatSetErrorIfFailure(pc->pmat,pc->erroriffailure);CHKERRQ(ierr);
71*3d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
72fcfd50ebSBarry Smith     if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);}
73fcfd50ebSBarry Smith     ierr = ISDestroy(&dir->col);CHKERRQ(ierr);
74075768bcSBarry Smith     ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
7503c60df9SBarry Smith     if (dir->row) {
763bb1ff40SBarry Smith       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);
773bb1ff40SBarry Smith       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);CHKERRQ(ierr);
7803c60df9SBarry Smith     }
79075768bcSBarry Smith     ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
8000e125f8SBarry Smith     ierr = MatFactorGetError(pc->pmat,&err);CHKERRQ(ierr);
8100e125f8SBarry Smith     if (err) { /* Factor() fails */
8200e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
836baea169SHong Zhang       PetscFunctionReturn(0);
846baea169SHong Zhang     }
856baea169SHong Zhang 
86075768bcSBarry Smith     ((PC_Factor*)dir)->fact = pc->pmat;
879b54502bSHong Zhang   } else {
889b54502bSHong Zhang     MatInfo info;
8900e125f8SBarry Smith 
909b54502bSHong Zhang     if (!pc->setupcalled) {
91075768bcSBarry Smith       ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
929b54502bSHong Zhang       if (dir->nonzerosalongdiagonal) {
939b54502bSHong Zhang         ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr);
949b54502bSHong Zhang       }
9503c60df9SBarry Smith       if (dir->row) {
963bb1ff40SBarry Smith         ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);
973bb1ff40SBarry Smith         ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);CHKERRQ(ierr);
9803c60df9SBarry Smith       }
99d09a07f4SBarry Smith       if (!((PC_Factor*)dir)->fact) {
100075768bcSBarry Smith         ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr);
101a1f19f5aSHong Zhang       }
102075768bcSBarry Smith       ierr                = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
103075768bcSBarry Smith       ierr                = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
104*3d1c1ea0SBarry Smith       dir->hdr.actualfill = info.fill_ratio_needed;
1053bb1ff40SBarry Smith       ierr                = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr);
1069b54502bSHong Zhang     } else if (pc->flag != SAME_NONZERO_PATTERN) {
107*3d1c1ea0SBarry Smith       if (!dir->hdr.reuseordering) {
108fcfd50ebSBarry Smith         if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);}
109fcfd50ebSBarry Smith         ierr = ISDestroy(&dir->col);CHKERRQ(ierr);
110075768bcSBarry Smith         ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
1119b54502bSHong Zhang         if (dir->nonzerosalongdiagonal) {
1129b54502bSHong Zhang           ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr);
1139b54502bSHong Zhang         }
11403c60df9SBarry Smith         if (dir->row) {
1153bb1ff40SBarry Smith           ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);
1163bb1ff40SBarry Smith           ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);CHKERRQ(ierr);
11703c60df9SBarry Smith         }
1189b54502bSHong Zhang       }
1196bf464f9SBarry Smith       ierr                = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr);
120075768bcSBarry Smith       ierr                = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr);
121075768bcSBarry Smith       ierr                = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
122075768bcSBarry Smith       ierr                = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
123*3d1c1ea0SBarry Smith       dir->hdr.actualfill = info.fill_ratio_needed;
1243bb1ff40SBarry Smith       ierr                = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr);
12504545d6dSBarry Smith     } else {
126b8b68cfdSBarry Smith       ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr);
127160a8794SBarry Smith       if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) {
128b8b68cfdSBarry Smith         ierr = MatFactorClearError(((PC_Factor*)dir)->fact);CHKERRQ(ierr);
129b8b68cfdSBarry Smith         pc->failedreason = PC_NOERROR;
13004545d6dSBarry Smith       }
1319b54502bSHong Zhang     }
13200e125f8SBarry Smith     ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr);
13300e125f8SBarry Smith     if (err) { /* FactorSymbolic() fails */
13400e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1358c1cd74cSHong Zhang       PetscFunctionReturn(0);
1368c1cd74cSHong Zhang     }
1378c1cd74cSHong Zhang 
138075768bcSBarry Smith     ierr = MatLUFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
13900e125f8SBarry Smith     ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr);
14000e125f8SBarry Smith     if (err) { /* FactorNumeric() fails */
14100e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1428c1cd74cSHong Zhang     }
143680c5173SHong Zhang 
1449b54502bSHong Zhang   }
14500c67f3bSHong Zhang 
14600c67f3bSHong Zhang   ierr = PCFactorGetMatSolverPackage(pc,&stype);CHKERRQ(ierr);
14700c67f3bSHong Zhang   if (!stype) {
14800e125f8SBarry Smith     const MatSolverPackage solverpackage;
14900e125f8SBarry Smith     ierr = MatFactorGetSolverPackage(((PC_Factor*)dir)->fact,&solverpackage);CHKERRQ(ierr);
15000e125f8SBarry Smith     ierr = PCFactorSetMatSolverPackage(pc,solverpackage);CHKERRQ(ierr);
15100c67f3bSHong Zhang   }
1529b54502bSHong Zhang   PetscFunctionReturn(0);
1539b54502bSHong Zhang }
1549b54502bSHong Zhang 
1559b54502bSHong Zhang #undef __FUNCT__
156574deadeSBarry Smith #define __FUNCT__ "PCReset_LU"
157574deadeSBarry Smith static PetscErrorCode PCReset_LU(PC pc)
1589b54502bSHong Zhang {
1599b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
1609b54502bSHong Zhang   PetscErrorCode ierr;
1619b54502bSHong Zhang 
1629b54502bSHong Zhang   PetscFunctionBegin;
163*3d1c1ea0SBarry Smith   if (!dir->hdr.inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr);}
164fcfd50ebSBarry Smith   if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);}
165fcfd50ebSBarry Smith   ierr = ISDestroy(&dir->col);CHKERRQ(ierr);
166574deadeSBarry Smith   PetscFunctionReturn(0);
167574deadeSBarry Smith }
168574deadeSBarry Smith 
169574deadeSBarry Smith #undef __FUNCT__
170574deadeSBarry Smith #define __FUNCT__ "PCDestroy_LU"
171574deadeSBarry Smith static PetscErrorCode PCDestroy_LU(PC pc)
172574deadeSBarry Smith {
173574deadeSBarry Smith   PC_LU          *dir = (PC_LU*)pc->data;
174574deadeSBarry Smith   PetscErrorCode ierr;
175574deadeSBarry Smith 
176574deadeSBarry Smith   PetscFunctionBegin;
177574deadeSBarry Smith   ierr = PCReset_LU(pc);CHKERRQ(ierr);
178503cfb0cSBarry Smith   ierr = PetscFree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
179503cfb0cSBarry Smith   ierr = PetscFree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr);
180c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1819b54502bSHong Zhang   PetscFunctionReturn(0);
1829b54502bSHong Zhang }
1839b54502bSHong Zhang 
1849b54502bSHong Zhang #undef __FUNCT__
1859b54502bSHong Zhang #define __FUNCT__ "PCApply_LU"
1869b54502bSHong Zhang static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y)
1879b54502bSHong Zhang {
1889b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
1899b54502bSHong Zhang   PetscErrorCode ierr;
1909b54502bSHong Zhang 
1919b54502bSHong Zhang   PetscFunctionBegin;
192*3d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
1932fa5cd67SKarl Rupp     ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);
1942fa5cd67SKarl Rupp   } else {
1952fa5cd67SKarl Rupp     ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);
1962fa5cd67SKarl Rupp   }
1979b54502bSHong Zhang   PetscFunctionReturn(0);
1989b54502bSHong Zhang }
1999b54502bSHong Zhang 
2009b54502bSHong Zhang #undef __FUNCT__
2019b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_LU"
2029b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y)
2039b54502bSHong Zhang {
2049b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
2059b54502bSHong Zhang   PetscErrorCode ierr;
2069b54502bSHong Zhang 
2079b54502bSHong Zhang   PetscFunctionBegin;
208*3d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
2092fa5cd67SKarl Rupp     ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);
2102fa5cd67SKarl Rupp   } else {
2112fa5cd67SKarl Rupp     ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);
2122fa5cd67SKarl Rupp   }
2139b54502bSHong Zhang   PetscFunctionReturn(0);
2149b54502bSHong Zhang }
2159b54502bSHong Zhang 
2169b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/
2179b54502bSHong Zhang 
2189b54502bSHong Zhang /*MC
2199b54502bSHong Zhang    PCLU - Uses a direct solver, based on LU factorization, as a preconditioner
2209b54502bSHong Zhang 
2219b54502bSHong Zhang    Options Database Keys:
2222401956bSBarry Smith +  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
223d45987f3SHong Zhang .  -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like superlu
2242401956bSBarry Smith .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
22555ba2a51SBarry Smith .  -pc_factor_fill <fill> - Sets fill amount
2262401956bSBarry Smith .  -pc_factor_in_place - Activates in-place factorization
2272401956bSBarry Smith .  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
2282401956bSBarry Smith .  -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase
2299b54502bSHong Zhang                                          stability of factorization.
230145b9266SHong Zhang .  -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types
231145b9266SHong Zhang .  -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default
232e22d95b2SBarry Smith -   -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the
233e22d95b2SBarry Smith         diagonal.
2349b54502bSHong Zhang 
2359b54502bSHong Zhang    Notes: Not all options work for all matrix formats
2369b54502bSHong Zhang           Run with -help to see additional options for particular matrix formats or factorization
2379b54502bSHong Zhang           algorithms
2389b54502bSHong Zhang 
2399b54502bSHong Zhang    Level: beginner
2409b54502bSHong Zhang 
2419b54502bSHong Zhang    Concepts: LU factorization, direct solver
2429b54502bSHong Zhang 
2439b54502bSHong Zhang    Notes: Usually this will compute an "exact" solution in one iteration and does
2449b54502bSHong Zhang           not need a Krylov method (i.e. you can use -ksp_type preonly, or
2459b54502bSHong Zhang           KSPSetType(ksp,KSPPREONLY) for the Krylov method
2469b54502bSHong Zhang 
2479b54502bSHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
248a4fd02acSBarry Smith            PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(),
2498ff23777SHong Zhang            PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetColumnPivot(),
250145b9266SHong Zhang            PCFactorSetPivotingInBlocks(),PCFactorSetShiftType(),PCFactorSetShiftAmount()
2518ff23777SHong Zhang            PCFactorReorderForNonzeroDiagonal()
2529b54502bSHong Zhang M*/
2539b54502bSHong Zhang 
2549b54502bSHong Zhang #undef __FUNCT__
2559b54502bSHong Zhang #define __FUNCT__ "PCCreate_LU"
2568cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_LU(PC pc)
2579b54502bSHong Zhang {
2589b54502bSHong Zhang   PetscErrorCode ierr;
2599b54502bSHong Zhang   PetscMPIInt    size;
2609b54502bSHong Zhang   PC_LU          *dir;
2619b54502bSHong Zhang 
2629b54502bSHong Zhang   PetscFunctionBegin;
263b00a9115SJed Brown   ierr     = PetscNewLog(pc,&dir);CHKERRQ(ierr);
264*3d1c1ea0SBarry Smith   pc->data = (void*)dir;
265*3d1c1ea0SBarry Smith   ierr     = PCFactorInitialize(pc);CHKERRQ(ierr);
266d5f3da31SBarry Smith   ((PC_Factor*)dir)->factortype = MAT_FACTOR_LU;
2679b54502bSHong Zhang   dir->nonzerosalongdiagonal    = PETSC_FALSE;
2689b54502bSHong Zhang 
269075768bcSBarry Smith   ((PC_Factor*)dir)->info.fill          = 5.0;
270075768bcSBarry Smith   ((PC_Factor*)dir)->info.dtcol         = 1.e-6;  /* default to pivoting; this is only thing PETSc LU supports */
271f4db908eSBarry Smith   ((PC_Factor*)dir)->info.shifttype     = (PetscReal)MAT_SHIFT_NONE;
2729b54502bSHong Zhang   dir->col                              = 0;
2739b54502bSHong Zhang   dir->row                              = 0;
2745c9eb25fSBarry Smith 
275ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
2769b54502bSHong Zhang   if (size == 1) {
27719fd82e9SBarry Smith     ierr = PetscStrallocpy(MATORDERINGND,(char**)&((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
2789b54502bSHong Zhang   } else {
27919fd82e9SBarry Smith     ierr = PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
2809b54502bSHong Zhang   }
2819b54502bSHong Zhang 
282574deadeSBarry Smith   pc->ops->reset             = PCReset_LU;
2839b54502bSHong Zhang   pc->ops->destroy           = PCDestroy_LU;
2849b54502bSHong Zhang   pc->ops->apply             = PCApply_LU;
2859b54502bSHong Zhang   pc->ops->applytranspose    = PCApplyTranspose_LU;
2869b54502bSHong Zhang   pc->ops->setup             = PCSetUp_LU;
2879b54502bSHong Zhang   pc->ops->setfromoptions    = PCSetFromOptions_LU;
2889b54502bSHong Zhang   pc->ops->view              = PCView_LU;
2899b54502bSHong Zhang   pc->ops->applyrichardson   = 0;
290bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",PCFactorReorderForNonzeroDiagonal_LU);CHKERRQ(ierr);
2919b54502bSHong Zhang   PetscFunctionReturn(0);
2929b54502bSHong Zhang }
293