Actual source code: lu.c
2: /*
3: Defines a direct factorization preconditioner for any Mat implementation
4: Note: this need not be consided a preconditioner since it supplies
5: a direct solver.
6: */
8: #include <../src/ksp/pc/impls/factor/lu/lu.h>
10: PetscErrorCode PCFactorReorderForNonzeroDiagonal_LU(PC pc,PetscReal z)
11: {
12: PC_LU *lu = (PC_LU*)pc->data;
14: lu->nonzerosalongdiagonal = PETSC_TRUE;
15: if (z == PETSC_DECIDE) lu->nonzerosalongdiagonaltol = 1.e-10;
16: else lu->nonzerosalongdiagonaltol = z;
17: return 0;
18: }
20: static PetscErrorCode PCSetFromOptions_LU(PetscOptionItems *PetscOptionsObject,PC pc)
21: {
22: PC_LU *lu = (PC_LU*)pc->data;
23: PetscBool flg = PETSC_FALSE;
24: PetscReal tol;
26: PetscOptionsHead(PetscOptionsObject,"LU options");
27: PCSetFromOptions_Factor(PetscOptionsObject,pc);
29: PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);
30: if (flg) {
31: tol = PETSC_DECIDE;
32: PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,NULL);
33: PCFactorReorderForNonzeroDiagonal(pc,tol);
34: }
35: PetscOptionsTail();
36: return 0;
37: }
39: static PetscErrorCode PCSetUp_LU(PC pc)
40: {
41: PC_LU *dir = (PC_LU*)pc->data;
42: MatSolverType stype;
43: MatFactorError err;
45: pc->failedreason = PC_NOERROR;
46: if (dir->hdr.reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->hdr.actualfill;
48: MatSetErrorIfFailure(pc->pmat,pc->erroriffailure);
49: if (dir->hdr.inplace) {
50: MatFactorType ftype;
52: MatGetFactorType(pc->pmat, &ftype);
53: if (ftype == MAT_FACTOR_NONE) {
54: if (dir->row && dir->col && dir->row != dir->col) ISDestroy(&dir->row);
55: ISDestroy(&dir->col);
56: /* This should only get the ordering if needed, but since MatGetFactor() is not called we can't know if it is needed */
57: PCFactorSetDefaultOrdering_Factor(pc);
58: MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);
59: if (dir->row) {
60: PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);
61: PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);
62: }
63: MatLUFactor(pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);
64: MatFactorGetError(pc->pmat,&err);
65: if (err) { /* Factor() fails */
66: pc->failedreason = (PCFailedReason)err;
67: return 0;
68: }
69: }
70: ((PC_Factor*)dir)->fact = pc->pmat;
71: } else {
72: MatInfo info;
74: if (!pc->setupcalled) {
75: PetscBool canuseordering;
76: if (!((PC_Factor*)dir)->fact) {
77: MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);
78: PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);
79: }
80: MatFactorGetCanUseOrdering(((PC_Factor*)dir)->fact,&canuseordering);
81: if (canuseordering) {
82: PCFactorSetDefaultOrdering_Factor(pc);
83: MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);
84: if (dir->nonzerosalongdiagonal) {
85: MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);
86: }
87: PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);
88: PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);
89: }
90: MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);
91: MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);
92: dir->hdr.actualfill = info.fill_ratio_needed;
93: } else if (pc->flag != SAME_NONZERO_PATTERN) {
94: PetscBool canuseordering;
95: if (!dir->hdr.reuseordering) {
96: MatDestroy(&((PC_Factor*)dir)->fact);
97: MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);
98: PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);
99: MatFactorGetCanUseOrdering(((PC_Factor*)dir)->fact,&canuseordering);
100: if (canuseordering) {
101: if (dir->row && dir->col && dir->row != dir->col) ISDestroy(&dir->row);
102: ISDestroy(&dir->col);
103: PCFactorSetDefaultOrdering_Factor(pc);
104: MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);
105: if (dir->nonzerosalongdiagonal) {
106: MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);
107: }
108: PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);
109: PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);
110: }
111: }
112: MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);
113: MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);
114: dir->hdr.actualfill = info.fill_ratio_needed;
115: } else {
116: MatFactorGetError(((PC_Factor*)dir)->fact,&err);
117: if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) {
118: MatFactorClearError(((PC_Factor*)dir)->fact);
119: pc->failedreason = PC_NOERROR;
120: }
121: }
122: MatFactorGetError(((PC_Factor*)dir)->fact,&err);
123: if (err) { /* FactorSymbolic() fails */
124: pc->failedreason = (PCFailedReason)err;
125: return 0;
126: }
128: MatLUFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);
129: MatFactorGetError(((PC_Factor*)dir)->fact,&err);
130: if (err) { /* FactorNumeric() fails */
131: pc->failedreason = (PCFailedReason)err;
132: }
134: }
136: PCFactorGetMatSolverType(pc,&stype);
137: if (!stype) {
138: MatSolverType solverpackage;
139: MatFactorGetSolverType(((PC_Factor*)dir)->fact,&solverpackage);
140: PCFactorSetMatSolverType(pc,solverpackage);
141: }
142: return 0;
143: }
145: static PetscErrorCode PCReset_LU(PC pc)
146: {
147: PC_LU *dir = (PC_LU*)pc->data;
149: if (!dir->hdr.inplace && ((PC_Factor*)dir)->fact) MatDestroy(&((PC_Factor*)dir)->fact);
150: if (dir->row && dir->col && dir->row != dir->col) ISDestroy(&dir->row);
151: ISDestroy(&dir->col);
152: return 0;
153: }
155: static PetscErrorCode PCDestroy_LU(PC pc)
156: {
157: PC_LU *dir = (PC_LU*)pc->data;
159: PCReset_LU(pc);
160: PetscFree(((PC_Factor*)dir)->ordering);
161: PetscFree(((PC_Factor*)dir)->solvertype);
162: PetscFree(pc->data);
163: return 0;
164: }
166: static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y)
167: {
168: PC_LU *dir = (PC_LU*)pc->data;
170: if (dir->hdr.inplace) {
171: MatSolve(pc->pmat,x,y);
172: } else {
173: MatSolve(((PC_Factor*)dir)->fact,x,y);
174: }
175: return 0;
176: }
178: static PetscErrorCode PCMatApply_LU(PC pc,Mat X,Mat Y)
179: {
180: PC_LU *dir = (PC_LU*)pc->data;
182: if (dir->hdr.inplace) {
183: MatMatSolve(pc->pmat,X,Y);
184: } else {
185: MatMatSolve(((PC_Factor*)dir)->fact,X,Y);
186: }
187: return 0;
188: }
190: static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y)
191: {
192: PC_LU *dir = (PC_LU*)pc->data;
194: if (dir->hdr.inplace) {
195: MatSolveTranspose(pc->pmat,x,y);
196: } else {
197: MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);
198: }
199: return 0;
200: }
202: /* -----------------------------------------------------------------------------------*/
204: /*MC
205: PCLU - Uses a direct solver, based on LU factorization, as a preconditioner
207: Options Database Keys:
208: + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
209: . -pc_factor_mat_solver_type - Actives PCFactorSetMatSolverType() to choose the direct solver, like superlu
210: . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
211: . -pc_factor_fill <fill> - Sets fill amount
212: . -pc_factor_in_place - Activates in-place factorization
213: . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
214: . -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase
215: stability of factorization.
216: . -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types
217: . -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default
218: - -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the
219: diagonal.
221: Notes:
222: Not all options work for all matrix formats
223: Run with -help to see additional options for particular matrix formats or factorization
224: algorithms
226: Level: beginner
228: Notes:
229: Usually this will compute an "exact" solution in one iteration and does
230: not need a Krylov method (i.e. you can use -ksp_type preonly, or
231: KSPSetType(ksp,KSPPREONLY) for the Krylov method
233: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC,
234: PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(),
235: PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetColumnPivot(),
236: PCFactorSetPivotInBlocks(),PCFactorSetShiftType(),PCFactorSetShiftAmount()
237: PCFactorReorderForNonzeroDiagonal()
238: M*/
240: PETSC_EXTERN PetscErrorCode PCCreate_LU(PC pc)
241: {
242: PC_LU *dir;
244: PetscNewLog(pc,&dir);
245: pc->data = (void*)dir;
246: PCFactorInitialize(pc,MAT_FACTOR_LU);
247: dir->nonzerosalongdiagonal = PETSC_FALSE;
249: ((PC_Factor*)dir)->info.fill = 5.0;
250: ((PC_Factor*)dir)->info.dtcol = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */
251: ((PC_Factor*)dir)->info.shifttype = (PetscReal)MAT_SHIFT_NONE;
252: dir->col = NULL;
253: dir->row = NULL;
255: pc->ops->reset = PCReset_LU;
256: pc->ops->destroy = PCDestroy_LU;
257: pc->ops->apply = PCApply_LU;
258: pc->ops->matapply = PCMatApply_LU;
259: pc->ops->applytranspose = PCApplyTranspose_LU;
260: pc->ops->setup = PCSetUp_LU;
261: pc->ops->setfromoptions = PCSetFromOptions_LU;
262: pc->ops->view = PCView_Factor;
263: pc->ops->applyrichardson = NULL;
264: PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",PCFactorReorderForNonzeroDiagonal_LU);
265: return 0;
266: }