xref: /petsc/src/ksp/pc/impls/factor/lu/lu.c (revision c6e4fdc6cd3657a39cd296b432bf7685bf018768)
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__
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 
939b54502bSHong Zhang #undef __FUNCT__
949b54502bSHong Zhang #define __FUNCT__ "PCSetUp_LU"
959b54502bSHong Zhang static PetscErrorCode PCSetUp_LU(PC pc)
969b54502bSHong Zhang {
979b54502bSHong Zhang   PetscErrorCode         ierr;
989b54502bSHong Zhang   PC_LU                  *dir = (PC_LU*)pc->data;
9900c67f3bSHong Zhang   const MatSolverPackage stype;
1009b54502bSHong Zhang 
1019b54502bSHong Zhang   PetscFunctionBegin;
102*c6e4fdc6SHong Zhang   pc->failedreason = PC_NOERROR;
103075768bcSBarry Smith   if (dir->reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->actualfill;
1049b54502bSHong Zhang 
10584d44b13SHong Zhang   ierr = MatSetErrorIfFailure(pc->pmat,pc->erroriffailure);CHKERRQ(ierr);
1069b54502bSHong Zhang   if (dir->inplace) {
107fcfd50ebSBarry Smith     if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);}
108fcfd50ebSBarry Smith     ierr = ISDestroy(&dir->col);CHKERRQ(ierr);
109075768bcSBarry Smith     ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
11003c60df9SBarry Smith     if (dir->row) {
1113bb1ff40SBarry Smith       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);
1123bb1ff40SBarry Smith       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);CHKERRQ(ierr);
11303c60df9SBarry Smith     }
114075768bcSBarry Smith     ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
115680c5173SHong Zhang     if (pc->pmat->errortype) { /* Factor() fails */
116680c5173SHong Zhang       pc->failedreason = (PCFailedReason)pc->pmat->errortype;
1176baea169SHong Zhang       PetscFunctionReturn(0);
1186baea169SHong Zhang     }
1196baea169SHong Zhang 
120075768bcSBarry Smith     ((PC_Factor*)dir)->fact = pc->pmat;
1219b54502bSHong Zhang   } else {
1229b54502bSHong Zhang     MatInfo info;
123680c5173SHong Zhang     Mat     F;
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) {
1303bb1ff40SBarry Smith         ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);
1313bb1ff40SBarry Smith         ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)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;
1393bb1ff40SBarry Smith       ierr            = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr);
1409b54502bSHong Zhang     } else if (pc->flag != SAME_NONZERO_PATTERN) {
1419b54502bSHong Zhang       if (!dir->reuseordering) {
142fcfd50ebSBarry Smith         if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);}
143fcfd50ebSBarry Smith         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) {
1493bb1ff40SBarry Smith           ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);
1503bb1ff40SBarry Smith           ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);CHKERRQ(ierr);
15103c60df9SBarry Smith         }
1529b54502bSHong Zhang       }
1536bf464f9SBarry 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;
1583bb1ff40SBarry Smith       ierr            = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr);
15904545d6dSBarry Smith     } else {
16004545d6dSBarry Smith       F = ((PC_Factor*)dir)->fact;
16104545d6dSBarry Smith       if ((PCFailedReason)F->errortype == PC_FACTOR_NUMERIC_ZEROPIVOT) {
16204545d6dSBarry Smith         F->errortype     = MAT_FACTOR_NOERROR;
16304545d6dSBarry Smith         pc->failedreason = (PCFailedReason)F->errortype;
16404545d6dSBarry Smith       }
1659b54502bSHong Zhang     }
166680c5173SHong Zhang     F = ((PC_Factor*)dir)->fact;
167680c5173SHong Zhang     if (F->errortype) { /* FactorSymbolic() fails */
168680c5173SHong Zhang       pc->failedreason = (PCFailedReason)F->errortype;
1698c1cd74cSHong Zhang       PetscFunctionReturn(0);
1708c1cd74cSHong Zhang     }
1718c1cd74cSHong Zhang 
172075768bcSBarry Smith     ierr = MatLUFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
173680c5173SHong Zhang     if (F->errortype) { /* FactorNumeric() fails */
174680c5173SHong Zhang       pc->failedreason = (PCFailedReason)F->errortype;
1758c1cd74cSHong Zhang     }
176680c5173SHong Zhang 
1779b54502bSHong Zhang   }
17800c67f3bSHong Zhang 
17900c67f3bSHong Zhang   ierr = PCFactorGetMatSolverPackage(pc,&stype);CHKERRQ(ierr);
18000c67f3bSHong Zhang   if (!stype) {
18100c67f3bSHong Zhang     ierr = PCFactorSetMatSolverPackage(pc,((PC_Factor*)dir)->fact->solvertype);CHKERRQ(ierr);
18200c67f3bSHong Zhang   }
1839b54502bSHong Zhang   PetscFunctionReturn(0);
1849b54502bSHong Zhang }
1859b54502bSHong Zhang 
1869b54502bSHong Zhang #undef __FUNCT__
187574deadeSBarry Smith #define __FUNCT__ "PCReset_LU"
188574deadeSBarry Smith static PetscErrorCode PCReset_LU(PC pc)
1899b54502bSHong Zhang {
1909b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
1919b54502bSHong Zhang   PetscErrorCode ierr;
1929b54502bSHong Zhang 
1939b54502bSHong Zhang   PetscFunctionBegin;
1946bf464f9SBarry Smith   if (!dir->inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr);}
195fcfd50ebSBarry Smith   if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);}
196fcfd50ebSBarry Smith   ierr = ISDestroy(&dir->col);CHKERRQ(ierr);
197574deadeSBarry Smith   PetscFunctionReturn(0);
198574deadeSBarry Smith }
199574deadeSBarry Smith 
200574deadeSBarry Smith #undef __FUNCT__
201574deadeSBarry Smith #define __FUNCT__ "PCDestroy_LU"
202574deadeSBarry Smith static PetscErrorCode PCDestroy_LU(PC pc)
203574deadeSBarry Smith {
204574deadeSBarry Smith   PC_LU          *dir = (PC_LU*)pc->data;
205574deadeSBarry Smith   PetscErrorCode ierr;
206574deadeSBarry Smith 
207574deadeSBarry Smith   PetscFunctionBegin;
208574deadeSBarry Smith   ierr = PCReset_LU(pc);CHKERRQ(ierr);
209503cfb0cSBarry Smith   ierr = PetscFree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
210503cfb0cSBarry Smith   ierr = PetscFree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr);
211c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
2129b54502bSHong Zhang   PetscFunctionReturn(0);
2139b54502bSHong Zhang }
2149b54502bSHong Zhang 
2159b54502bSHong Zhang #undef __FUNCT__
2169b54502bSHong Zhang #define __FUNCT__ "PCApply_LU"
2179b54502bSHong Zhang static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y)
2189b54502bSHong Zhang {
2199b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
2209b54502bSHong Zhang   PetscErrorCode ierr;
2219b54502bSHong Zhang 
2229b54502bSHong Zhang   PetscFunctionBegin;
2232fa5cd67SKarl Rupp   if (dir->inplace) {
2242fa5cd67SKarl Rupp     ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);
2252fa5cd67SKarl Rupp   } else {
2262fa5cd67SKarl Rupp     ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);
2272fa5cd67SKarl Rupp   }
2289b54502bSHong Zhang   PetscFunctionReturn(0);
2299b54502bSHong Zhang }
2309b54502bSHong Zhang 
2319b54502bSHong Zhang #undef __FUNCT__
2329b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_LU"
2339b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y)
2349b54502bSHong Zhang {
2359b54502bSHong Zhang   PC_LU          *dir = (PC_LU*)pc->data;
2369b54502bSHong Zhang   PetscErrorCode ierr;
2379b54502bSHong Zhang 
2389b54502bSHong Zhang   PetscFunctionBegin;
2392fa5cd67SKarl Rupp   if (dir->inplace) {
2402fa5cd67SKarl Rupp     ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);
2412fa5cd67SKarl Rupp   } else {
2422fa5cd67SKarl Rupp     ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);
2432fa5cd67SKarl Rupp   }
2449b54502bSHong Zhang   PetscFunctionReturn(0);
2459b54502bSHong Zhang }
2469b54502bSHong Zhang 
2479b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/
2489b54502bSHong Zhang 
2499b54502bSHong Zhang #undef __FUNCT__
2502401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace_LU"
2518e37d05fSBarry Smith PetscErrorCode  PCFactorSetUseInPlace_LU(PC pc,PetscBool flg)
2529b54502bSHong Zhang {
253c62fe420SBarry Smith   PC_LU *dir = (PC_LU*)pc->data;
2549b54502bSHong Zhang 
2559b54502bSHong Zhang   PetscFunctionBegin;
2568e37d05fSBarry Smith   dir->inplace = flg;
2578e37d05fSBarry Smith   PetscFunctionReturn(0);
2588e37d05fSBarry Smith }
2598e37d05fSBarry Smith 
2608e37d05fSBarry Smith #undef __FUNCT__
2618e37d05fSBarry Smith #define __FUNCT__ "PCFactorGetUseInPlace_LU"
2628e37d05fSBarry Smith PetscErrorCode  PCFactorGetUseInPlace_LU(PC pc,PetscBool *flg)
2638e37d05fSBarry Smith {
2648e37d05fSBarry Smith   PC_LU *dir = (PC_LU*)pc->data;
2658e37d05fSBarry Smith 
2668e37d05fSBarry Smith   PetscFunctionBegin;
2678e37d05fSBarry Smith   *flg = dir->inplace;
2689b54502bSHong Zhang   PetscFunctionReturn(0);
2699b54502bSHong Zhang }
2709b54502bSHong Zhang 
2719b54502bSHong Zhang /* ------------------------------------------------------------------------ */
2729b54502bSHong Zhang 
2739b54502bSHong Zhang /*MC
2749b54502bSHong Zhang    PCLU - Uses a direct solver, based on LU factorization, as a preconditioner
2759b54502bSHong Zhang 
2769b54502bSHong Zhang    Options Database Keys:
2772401956bSBarry Smith +  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
278d45987f3SHong Zhang .  -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like superlu
2792401956bSBarry Smith .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
28055ba2a51SBarry Smith .  -pc_factor_fill <fill> - Sets fill amount
2812401956bSBarry Smith .  -pc_factor_in_place - Activates in-place factorization
2822401956bSBarry Smith .  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
2832401956bSBarry Smith .  -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase
2849b54502bSHong Zhang                                          stability of factorization.
285145b9266SHong Zhang .  -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types
286145b9266SHong Zhang .  -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default
287e22d95b2SBarry Smith -   -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the
288e22d95b2SBarry Smith         diagonal.
2899b54502bSHong Zhang 
2909b54502bSHong Zhang    Notes: Not all options work for all matrix formats
2919b54502bSHong Zhang           Run with -help to see additional options for particular matrix formats or factorization
2929b54502bSHong Zhang           algorithms
2939b54502bSHong Zhang 
2949b54502bSHong Zhang    Level: beginner
2959b54502bSHong Zhang 
2969b54502bSHong Zhang    Concepts: LU factorization, direct solver
2979b54502bSHong Zhang 
2989b54502bSHong Zhang    Notes: Usually this will compute an "exact" solution in one iteration and does
2999b54502bSHong Zhang           not need a Krylov method (i.e. you can use -ksp_type preonly, or
3009b54502bSHong Zhang           KSPSetType(ksp,KSPPREONLY) for the Krylov method
3019b54502bSHong Zhang 
3029b54502bSHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
303a4fd02acSBarry Smith            PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(),
3048ff23777SHong Zhang            PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetColumnPivot(),
305145b9266SHong Zhang            PCFactorSetPivotingInBlocks(),PCFactorSetShiftType(),PCFactorSetShiftAmount()
3068ff23777SHong Zhang            PCFactorReorderForNonzeroDiagonal()
3079b54502bSHong Zhang M*/
3089b54502bSHong Zhang 
3099b54502bSHong Zhang #undef __FUNCT__
3109b54502bSHong Zhang #define __FUNCT__ "PCCreate_LU"
3118cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_LU(PC pc)
3129b54502bSHong Zhang {
3139b54502bSHong Zhang   PetscErrorCode ierr;
3149b54502bSHong Zhang   PetscMPIInt    size;
3159b54502bSHong Zhang   PC_LU          *dir;
3169b54502bSHong Zhang 
3179b54502bSHong Zhang   PetscFunctionBegin;
318b00a9115SJed Brown   ierr = PetscNewLog(pc,&dir);CHKERRQ(ierr);
3199b54502bSHong Zhang 
320075768bcSBarry Smith   ierr = MatFactorInfoInitialize(&((PC_Factor*)dir)->info);CHKERRQ(ierr);
3212fa5cd67SKarl Rupp 
3220298fd71SBarry Smith   ((PC_Factor*)dir)->fact       = NULL;
323d5f3da31SBarry Smith   ((PC_Factor*)dir)->factortype = MAT_FACTOR_LU;
3249b54502bSHong Zhang   dir->inplace                  = PETSC_FALSE;
3259b54502bSHong Zhang   dir->nonzerosalongdiagonal    = PETSC_FALSE;
3269b54502bSHong Zhang 
327075768bcSBarry Smith   ((PC_Factor*)dir)->info.fill          = 5.0;
328075768bcSBarry Smith   ((PC_Factor*)dir)->info.dtcol         = 1.e-6;  /* default to pivoting; this is only thing PETSc LU supports */
329f4db908eSBarry Smith   ((PC_Factor*)dir)->info.shifttype     = (PetscReal)MAT_SHIFT_NONE;
330915743fcSHong Zhang   ((PC_Factor*)dir)->info.shiftamount   = 0.0;
3310ed735ceSHong Zhang   ((PC_Factor*)dir)->info.zeropivot     = 100.0*PETSC_MACHINE_EPSILON;
332075768bcSBarry Smith   ((PC_Factor*)dir)->info.pivotinblocks = 1.0;
3339b54502bSHong Zhang   dir->col                              = 0;
3349b54502bSHong Zhang   dir->row                              = 0;
3355c9eb25fSBarry Smith 
336ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
3379b54502bSHong Zhang   if (size == 1) {
33819fd82e9SBarry Smith     ierr = PetscStrallocpy(MATORDERINGND,(char**)&((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
3399b54502bSHong Zhang   } else {
34019fd82e9SBarry Smith     ierr = PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
3419b54502bSHong Zhang   }
3429b54502bSHong Zhang   dir->reusefill     = PETSC_FALSE;
3439b54502bSHong Zhang   dir->reuseordering = PETSC_FALSE;
3449b54502bSHong Zhang   pc->data           = (void*)dir;
3459b54502bSHong Zhang 
346574deadeSBarry Smith   pc->ops->reset             = PCReset_LU;
3479b54502bSHong Zhang   pc->ops->destroy           = PCDestroy_LU;
3489b54502bSHong Zhang   pc->ops->apply             = PCApply_LU;
3499b54502bSHong Zhang   pc->ops->applytranspose    = PCApplyTranspose_LU;
3509b54502bSHong Zhang   pc->ops->setup             = PCSetUp_LU;
3519b54502bSHong Zhang   pc->ops->setfromoptions    = PCSetFromOptions_LU;
3529b54502bSHong Zhang   pc->ops->view              = PCView_LU;
3539b54502bSHong Zhang   pc->ops->applyrichardson   = 0;
35485317021SBarry Smith   pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor;
3559b54502bSHong Zhang 
356c7f610a1SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",PCFactorSetZeroPivot_Factor);CHKERRQ(ierr);
357c7f610a1SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetZeroPivot_C",PCFactorGetZeroPivot_Factor);CHKERRQ(ierr);
358c7f610a1SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftType_C",PCFactorSetShiftType_Factor);CHKERRQ(ierr);
359c7f610a1SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetShiftType_C",PCFactorGetShiftType_Factor);CHKERRQ(ierr);
360c7f610a1SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftAmount_C",PCFactorSetShiftAmount_Factor);CHKERRQ(ierr);
361c7f610a1SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetShiftAmount_C",PCFactorGetShiftAmount_Factor);CHKERRQ(ierr);
362bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUpMatSolverPackage_C",PCFactorSetUpMatSolverPackage_Factor);CHKERRQ(ierr);
363bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr);
364bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr);
365bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetFill_C",PCFactorSetFill_Factor);CHKERRQ(ierr);
366bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",PCFactorSetUseInPlace_LU);CHKERRQ(ierr);
3678e37d05fSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetUseInPlace_C",PCFactorGetUseInPlace_LU);CHKERRQ(ierr);
368bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr);
369bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",PCFactorSetReuseOrdering_LU);CHKERRQ(ierr);
370bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseFill_C",PCFactorSetReuseFill_LU);CHKERRQ(ierr);
371bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetColumnPivot_C",PCFactorSetColumnPivot_Factor);CHKERRQ(ierr);
372bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",PCFactorSetPivotInBlocks_Factor);CHKERRQ(ierr);
373bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",PCFactorReorderForNonzeroDiagonal_LU);CHKERRQ(ierr);
3749b54502bSHong Zhang   PetscFunctionReturn(0);
3759b54502bSHong Zhang }
376