xref: /petsc/src/ksp/pc/impls/factor/lu/lu.c (revision 680c5173fc5ccf4352d0ebdebcd816da6a3f83ed)
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 
10*680c5173SHong 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__
252401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering_LU"
267087cfbeSBarry Smith PetscErrorCode  PCFactorSetReuseOrdering_LU(PC pc,PetscBool flag)
279b54502bSHong Zhang {
28c62fe420SBarry Smith   PC_LU *lu = (PC_LU*)pc->data;
299b54502bSHong Zhang 
309b54502bSHong Zhang   PetscFunctionBegin;
319b54502bSHong Zhang   lu->reuseordering = flag;
329b54502bSHong Zhang   PetscFunctionReturn(0);
339b54502bSHong Zhang }
349b54502bSHong Zhang 
359b54502bSHong Zhang #undef __FUNCT__
362401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseFill_LU"
377087cfbeSBarry Smith PetscErrorCode  PCFactorSetReuseFill_LU(PC pc,PetscBool flag)
389b54502bSHong Zhang {
39c62fe420SBarry Smith   PC_LU *lu = (PC_LU*)pc->data;
409b54502bSHong Zhang 
419b54502bSHong Zhang   PetscFunctionBegin;
429b54502bSHong Zhang   lu->reusefill = flag;
439b54502bSHong Zhang   PetscFunctionReturn(0);
449b54502bSHong Zhang }
459b54502bSHong Zhang 
469b54502bSHong Zhang #undef __FUNCT__
479b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_LU"
484416b707SBarry Smith static PetscErrorCode PCSetFromOptions_LU(PetscOptionItems *PetscOptionsObject,PC pc)
499b54502bSHong Zhang {
509b54502bSHong Zhang   PC_LU          *lu = (PC_LU*)pc->data;
519b54502bSHong Zhang   PetscErrorCode ierr;
52ace3abfcSBarry Smith   PetscBool      flg = PETSC_FALSE;
539b54502bSHong Zhang   PetscReal      tol;
549b54502bSHong Zhang 
559b54502bSHong Zhang   PetscFunctionBegin;
56e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"LU options");CHKERRQ(ierr);
57e55864a3SBarry Smith   ierr = PCSetFromOptions_Factor(PetscOptionsObject,pc);CHKERRQ(ierr);
585c9eb25fSBarry Smith 
592401956bSBarry Smith   ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr);
609b54502bSHong Zhang   if (flg) {
619b54502bSHong Zhang     tol  = PETSC_DECIDE;
622401956bSBarry Smith     ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr);
632401956bSBarry Smith     ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr);
649b54502bSHong Zhang   }
659b54502bSHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
669b54502bSHong Zhang   PetscFunctionReturn(0);
679b54502bSHong Zhang }
689b54502bSHong Zhang 
699b54502bSHong Zhang #undef __FUNCT__
709b54502bSHong Zhang #define __FUNCT__ "PCView_LU"
719b54502bSHong Zhang static PetscErrorCode PCView_LU(PC pc,PetscViewer viewer)
729b54502bSHong Zhang {
739b54502bSHong Zhang   PC_LU          *lu = (PC_LU*)pc->data;
749b54502bSHong Zhang   PetscErrorCode ierr;
75ace3abfcSBarry Smith   PetscBool      iascii;
769b54502bSHong Zhang 
779b54502bSHong Zhang   PetscFunctionBegin;
78251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
799b54502bSHong Zhang   if (iascii) {
80914a5d51SHong Zhang     if (lu->inplace) {
81914a5d51SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  LU: in-place factorization\n");CHKERRQ(ierr);
82914a5d51SHong Zhang     } else {
83914a5d51SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  LU: out-of-place factorization\n");CHKERRQ(ierr);
84145b9266SHong Zhang     }
85145b9266SHong Zhang 
869b54502bSHong Zhang     if (lu->reusefill)    {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing fill from past factorization\n");CHKERRQ(ierr);}
879b54502bSHong Zhang     if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing reordering from past factorization\n");CHKERRQ(ierr);}
889b54502bSHong Zhang   }
89914a5d51SHong Zhang   ierr = PCView_Factor(pc,viewer);CHKERRQ(ierr);
909b54502bSHong Zhang   PetscFunctionReturn(0);
919b54502bSHong Zhang }
929b54502bSHong Zhang 
93*680c5173SHong Zhang #include <petsc/private/matimpl.h>
949b54502bSHong Zhang #undef __FUNCT__
959b54502bSHong Zhang #define __FUNCT__ "PCSetUp_LU"
969b54502bSHong Zhang static PetscErrorCode PCSetUp_LU(PC pc)
979b54502bSHong Zhang {
989b54502bSHong Zhang   PetscErrorCode ierr;
999b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
1009b54502bSHong Zhang 
1019b54502bSHong Zhang   PetscFunctionBegin;
102075768bcSBarry Smith   if (dir->reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->actualfill;
1039b54502bSHong Zhang 
10484d44b13SHong Zhang   ierr = MatSetErrorIfFailure(pc->pmat,pc->erroriffailure);CHKERRQ(ierr);
1059b54502bSHong Zhang   if (dir->inplace) {
106fcfd50ebSBarry Smith     if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);}
107fcfd50ebSBarry Smith     ierr = ISDestroy(&dir->col);CHKERRQ(ierr);
108075768bcSBarry Smith     ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
10903c60df9SBarry Smith     if (dir->row) {
1103bb1ff40SBarry Smith       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);
1113bb1ff40SBarry Smith       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);CHKERRQ(ierr);
11203c60df9SBarry Smith     }
113075768bcSBarry Smith     ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
114*680c5173SHong Zhang     if (pc->pmat->errortype) { /* Factor() fails */
115*680c5173SHong Zhang       pc->failedreason = (PCFailedReason)pc->pmat->errortype;
1166baea169SHong Zhang       PetscFunctionReturn(0);
1176baea169SHong Zhang     }
1186baea169SHong Zhang 
119075768bcSBarry Smith     ((PC_Factor*)dir)->fact = pc->pmat;
1209b54502bSHong Zhang   } else {
1219b54502bSHong Zhang     MatInfo info;
122*680c5173SHong Zhang     Mat     F;
1239b54502bSHong Zhang     if (!pc->setupcalled) {
124075768bcSBarry Smith       ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
1259b54502bSHong Zhang       if (dir->nonzerosalongdiagonal) {
1269b54502bSHong Zhang         ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr);
1279b54502bSHong Zhang       }
12803c60df9SBarry Smith       if (dir->row) {
1293bb1ff40SBarry Smith         ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);
1303bb1ff40SBarry Smith         ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);CHKERRQ(ierr);
13103c60df9SBarry Smith       }
132d09a07f4SBarry Smith       if (!((PC_Factor*)dir)->fact) {
133075768bcSBarry Smith         ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr);
134a1f19f5aSHong Zhang       }
135075768bcSBarry Smith       ierr            = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
136075768bcSBarry Smith       ierr            = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
1379b54502bSHong Zhang       dir->actualfill = info.fill_ratio_needed;
1383bb1ff40SBarry Smith       ierr            = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr);
1399b54502bSHong Zhang     } else if (pc->flag != SAME_NONZERO_PATTERN) {
1409b54502bSHong Zhang       if (!dir->reuseordering) {
141fcfd50ebSBarry Smith         if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);}
142fcfd50ebSBarry Smith         ierr = ISDestroy(&dir->col);CHKERRQ(ierr);
143075768bcSBarry Smith         ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
1449b54502bSHong Zhang         if (dir->nonzerosalongdiagonal) {
1459b54502bSHong Zhang           ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr);
1469b54502bSHong Zhang         }
14703c60df9SBarry Smith         if (dir->row) {
1483bb1ff40SBarry Smith           ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);
1493bb1ff40SBarry Smith           ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);CHKERRQ(ierr);
15003c60df9SBarry Smith         }
1519b54502bSHong Zhang       }
1526bf464f9SBarry Smith       ierr            = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr);
153075768bcSBarry Smith       ierr            = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr);
154075768bcSBarry Smith       ierr            = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
155075768bcSBarry Smith       ierr            = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
1569b54502bSHong Zhang       dir->actualfill = info.fill_ratio_needed;
1573bb1ff40SBarry Smith       ierr            = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr);
1589b54502bSHong Zhang     }
159*680c5173SHong Zhang     F = ((PC_Factor*)dir)->fact;
160*680c5173SHong Zhang     if (F->errortype) { /* FactorSymbolic() fails */
161*680c5173SHong Zhang       pc->failedreason = (PCFailedReason)F->errortype;
1628c1cd74cSHong Zhang       PetscFunctionReturn(0);
1638c1cd74cSHong Zhang     }
1648c1cd74cSHong Zhang 
165075768bcSBarry Smith     ierr = MatLUFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
166*680c5173SHong Zhang     if (F->errortype) { /* FactorNumeric() fails */
167*680c5173SHong Zhang       pc->failedreason = (PCFailedReason)F->errortype;
1688c1cd74cSHong Zhang     }
169*680c5173SHong Zhang 
1709b54502bSHong Zhang   }
1719b54502bSHong Zhang   PetscFunctionReturn(0);
1729b54502bSHong Zhang }
1739b54502bSHong Zhang 
1749b54502bSHong Zhang #undef __FUNCT__
175574deadeSBarry Smith #define __FUNCT__ "PCReset_LU"
176574deadeSBarry Smith static PetscErrorCode PCReset_LU(PC pc)
1779b54502bSHong Zhang {
1789b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
1799b54502bSHong Zhang   PetscErrorCode ierr;
1809b54502bSHong Zhang 
1819b54502bSHong Zhang   PetscFunctionBegin;
1826bf464f9SBarry Smith   if (!dir->inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr);}
183fcfd50ebSBarry Smith   if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);}
184fcfd50ebSBarry Smith   ierr = ISDestroy(&dir->col);CHKERRQ(ierr);
185574deadeSBarry Smith   PetscFunctionReturn(0);
186574deadeSBarry Smith }
187574deadeSBarry Smith 
188574deadeSBarry Smith #undef __FUNCT__
189574deadeSBarry Smith #define __FUNCT__ "PCDestroy_LU"
190574deadeSBarry Smith static PetscErrorCode PCDestroy_LU(PC pc)
191574deadeSBarry Smith {
192574deadeSBarry Smith   PC_LU          *dir = (PC_LU*)pc->data;
193574deadeSBarry Smith   PetscErrorCode ierr;
194574deadeSBarry Smith 
195574deadeSBarry Smith   PetscFunctionBegin;
196574deadeSBarry Smith   ierr = PCReset_LU(pc);CHKERRQ(ierr);
197503cfb0cSBarry Smith   ierr = PetscFree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
198503cfb0cSBarry Smith   ierr = PetscFree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr);
199c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
2009b54502bSHong Zhang   PetscFunctionReturn(0);
2019b54502bSHong Zhang }
2029b54502bSHong Zhang 
2039b54502bSHong Zhang #undef __FUNCT__
2049b54502bSHong Zhang #define __FUNCT__ "PCApply_LU"
2059b54502bSHong Zhang static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y)
2069b54502bSHong Zhang {
2079b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
2089b54502bSHong Zhang   PetscErrorCode ierr;
2099b54502bSHong Zhang 
2109b54502bSHong Zhang   PetscFunctionBegin;
2112fa5cd67SKarl Rupp   if (dir->inplace) {
2122fa5cd67SKarl Rupp     ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);
2132fa5cd67SKarl Rupp   } else {
2142fa5cd67SKarl Rupp     ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);
2152fa5cd67SKarl Rupp   }
2169b54502bSHong Zhang   PetscFunctionReturn(0);
2179b54502bSHong Zhang }
2189b54502bSHong Zhang 
2199b54502bSHong Zhang #undef __FUNCT__
2209b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_LU"
2219b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y)
2229b54502bSHong Zhang {
2239b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
2249b54502bSHong Zhang   PetscErrorCode ierr;
2259b54502bSHong Zhang 
2269b54502bSHong Zhang   PetscFunctionBegin;
2272fa5cd67SKarl Rupp   if (dir->inplace) {
2282fa5cd67SKarl Rupp     ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);
2292fa5cd67SKarl Rupp   } else {
2302fa5cd67SKarl Rupp     ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);
2312fa5cd67SKarl Rupp   }
2329b54502bSHong Zhang   PetscFunctionReturn(0);
2339b54502bSHong Zhang }
2349b54502bSHong Zhang 
2359b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/
2369b54502bSHong Zhang 
2379b54502bSHong Zhang #undef __FUNCT__
2382401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace_LU"
2398e37d05fSBarry Smith PetscErrorCode  PCFactorSetUseInPlace_LU(PC pc,PetscBool flg)
2409b54502bSHong Zhang {
241c62fe420SBarry Smith   PC_LU *dir = (PC_LU*)pc->data;
2429b54502bSHong Zhang 
2439b54502bSHong Zhang   PetscFunctionBegin;
2448e37d05fSBarry Smith   dir->inplace = flg;
2458e37d05fSBarry Smith   PetscFunctionReturn(0);
2468e37d05fSBarry Smith }
2478e37d05fSBarry Smith 
2488e37d05fSBarry Smith #undef __FUNCT__
2498e37d05fSBarry Smith #define __FUNCT__ "PCFactorGetUseInPlace_LU"
2508e37d05fSBarry Smith PetscErrorCode  PCFactorGetUseInPlace_LU(PC pc,PetscBool *flg)
2518e37d05fSBarry Smith {
2528e37d05fSBarry Smith   PC_LU *dir = (PC_LU*)pc->data;
2538e37d05fSBarry Smith 
2548e37d05fSBarry Smith   PetscFunctionBegin;
2558e37d05fSBarry Smith   *flg = dir->inplace;
2569b54502bSHong Zhang   PetscFunctionReturn(0);
2579b54502bSHong Zhang }
2589b54502bSHong Zhang 
2599b54502bSHong Zhang /* ------------------------------------------------------------------------ */
2609b54502bSHong Zhang 
2619b54502bSHong Zhang /*MC
2629b54502bSHong Zhang    PCLU - Uses a direct solver, based on LU factorization, as a preconditioner
2639b54502bSHong Zhang 
2649b54502bSHong Zhang    Options Database Keys:
2652401956bSBarry Smith +  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
266d45987f3SHong Zhang .  -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like superlu
2672401956bSBarry Smith .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
26855ba2a51SBarry Smith .  -pc_factor_fill <fill> - Sets fill amount
2692401956bSBarry Smith .  -pc_factor_in_place - Activates in-place factorization
2702401956bSBarry Smith .  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
2712401956bSBarry Smith .  -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase
2729b54502bSHong Zhang                                          stability of factorization.
273145b9266SHong Zhang .  -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types
274145b9266SHong Zhang .  -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default
275e22d95b2SBarry Smith -   -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the
276e22d95b2SBarry Smith         diagonal.
2779b54502bSHong Zhang 
2789b54502bSHong Zhang    Notes: Not all options work for all matrix formats
2799b54502bSHong Zhang           Run with -help to see additional options for particular matrix formats or factorization
2809b54502bSHong Zhang           algorithms
2819b54502bSHong Zhang 
2829b54502bSHong Zhang    Level: beginner
2839b54502bSHong Zhang 
2849b54502bSHong Zhang    Concepts: LU factorization, direct solver
2859b54502bSHong Zhang 
2869b54502bSHong Zhang    Notes: Usually this will compute an "exact" solution in one iteration and does
2879b54502bSHong Zhang           not need a Krylov method (i.e. you can use -ksp_type preonly, or
2889b54502bSHong Zhang           KSPSetType(ksp,KSPPREONLY) for the Krylov method
2899b54502bSHong Zhang 
2909b54502bSHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
291a4fd02acSBarry Smith            PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(),
2928ff23777SHong Zhang            PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetColumnPivot(),
293145b9266SHong Zhang            PCFactorSetPivotingInBlocks(),PCFactorSetShiftType(),PCFactorSetShiftAmount()
2948ff23777SHong Zhang            PCFactorReorderForNonzeroDiagonal()
2959b54502bSHong Zhang M*/
2969b54502bSHong Zhang 
2979b54502bSHong Zhang #undef __FUNCT__
2989b54502bSHong Zhang #define __FUNCT__ "PCCreate_LU"
2998cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_LU(PC pc)
3009b54502bSHong Zhang {
3019b54502bSHong Zhang   PetscErrorCode ierr;
3029b54502bSHong Zhang   PetscMPIInt    size;
3039b54502bSHong Zhang   PC_LU          *dir;
3049b54502bSHong Zhang 
3059b54502bSHong Zhang   PetscFunctionBegin;
306b00a9115SJed Brown   ierr = PetscNewLog(pc,&dir);CHKERRQ(ierr);
3079b54502bSHong Zhang 
308075768bcSBarry Smith   ierr = MatFactorInfoInitialize(&((PC_Factor*)dir)->info);CHKERRQ(ierr);
3092fa5cd67SKarl Rupp 
3100298fd71SBarry Smith   ((PC_Factor*)dir)->fact       = NULL;
311d5f3da31SBarry Smith   ((PC_Factor*)dir)->factortype = MAT_FACTOR_LU;
3129b54502bSHong Zhang   dir->inplace                  = PETSC_FALSE;
3139b54502bSHong Zhang   dir->nonzerosalongdiagonal    = PETSC_FALSE;
3149b54502bSHong Zhang 
315075768bcSBarry Smith   ((PC_Factor*)dir)->info.fill          = 5.0;
316075768bcSBarry Smith   ((PC_Factor*)dir)->info.dtcol         = 1.e-6;  /* default to pivoting; this is only thing PETSc LU supports */
317f4db908eSBarry Smith   ((PC_Factor*)dir)->info.shifttype     = (PetscReal)MAT_SHIFT_NONE;
318915743fcSHong Zhang   ((PC_Factor*)dir)->info.shiftamount   = 0.0;
3190ed735ceSHong Zhang   ((PC_Factor*)dir)->info.zeropivot     = 100.0*PETSC_MACHINE_EPSILON;
320075768bcSBarry Smith   ((PC_Factor*)dir)->info.pivotinblocks = 1.0;
3219b54502bSHong Zhang   dir->col                              = 0;
3229b54502bSHong Zhang   dir->row                              = 0;
3235c9eb25fSBarry Smith 
324a1f19f5aSHong Zhang   ierr = PetscStrallocpy(MATSOLVERPETSC,&((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); /* default SolverPackage */
325ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
3269b54502bSHong Zhang   if (size == 1) {
32719fd82e9SBarry Smith     ierr = PetscStrallocpy(MATORDERINGND,(char**)&((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
3289b54502bSHong Zhang   } else {
32919fd82e9SBarry Smith     ierr = PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
3309b54502bSHong Zhang   }
3319b54502bSHong Zhang   dir->reusefill     = PETSC_FALSE;
3329b54502bSHong Zhang   dir->reuseordering = PETSC_FALSE;
3339b54502bSHong Zhang   pc->data           = (void*)dir;
3349b54502bSHong Zhang 
335574deadeSBarry Smith   pc->ops->reset             = PCReset_LU;
3369b54502bSHong Zhang   pc->ops->destroy           = PCDestroy_LU;
3379b54502bSHong Zhang   pc->ops->apply             = PCApply_LU;
3389b54502bSHong Zhang   pc->ops->applytranspose    = PCApplyTranspose_LU;
3399b54502bSHong Zhang   pc->ops->setup             = PCSetUp_LU;
3409b54502bSHong Zhang   pc->ops->setfromoptions    = PCSetFromOptions_LU;
3419b54502bSHong Zhang   pc->ops->view              = PCView_LU;
3429b54502bSHong Zhang   pc->ops->applyrichardson   = 0;
34385317021SBarry Smith   pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor;
3449b54502bSHong Zhang 
345bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUpMatSolverPackage_C",PCFactorSetUpMatSolverPackage_Factor);CHKERRQ(ierr);
346bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr);
347bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr);
348bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",PCFactorSetZeroPivot_Factor);CHKERRQ(ierr);
349bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftType_C",PCFactorSetShiftType_Factor);CHKERRQ(ierr);
350bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftAmount_C",PCFactorSetShiftAmount_Factor);CHKERRQ(ierr);
351bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetFill_C",PCFactorSetFill_Factor);CHKERRQ(ierr);
352bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",PCFactorSetUseInPlace_LU);CHKERRQ(ierr);
3538e37d05fSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetUseInPlace_C",PCFactorGetUseInPlace_LU);CHKERRQ(ierr);
354bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr);
355bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",PCFactorSetReuseOrdering_LU);CHKERRQ(ierr);
356bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseFill_C",PCFactorSetReuseFill_LU);CHKERRQ(ierr);
357bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetColumnPivot_C",PCFactorSetColumnPivot_Factor);CHKERRQ(ierr);
358bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",PCFactorSetPivotInBlocks_Factor);CHKERRQ(ierr);
359bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",PCFactorReorderForNonzeroDiagonal_LU);CHKERRQ(ierr);
3609b54502bSHong Zhang   PetscFunctionReturn(0);
3619b54502bSHong Zhang }
362