Actual source code: lu.c
1: #define PETSCKSP_DLL
3: /*
4: Defines a direct factorization preconditioner for any Mat implementation
5: Note: this need not be consided a preconditioner since it supplies
6: a direct solver.
7: */
9: #include ../src/ksp/pc/impls/factor/lu/lu.h
14: PetscErrorCode PCFactorReorderForNonzeroDiagonal_LU(PC pc,PetscReal z)
15: {
16: PC_LU *lu = (PC_LU*)pc->data;
19: lu->nonzerosalongdiagonal = PETSC_TRUE;
20: if (z == PETSC_DECIDE) {
21: lu->nonzerosalongdiagonaltol = 1.e-10;
22: } else {
23: lu->nonzerosalongdiagonaltol = z;
24: }
25: return(0);
26: }
32: PetscErrorCode PCFactorSetReuseOrdering_LU(PC pc,PetscTruth flag)
33: {
34: PC_LU *lu = (PC_LU*)pc->data;
37: lu->reuseordering = flag;
38: return(0);
39: }
45: PetscErrorCode PCFactorSetReuseFill_LU(PC pc,PetscTruth flag)
46: {
47: PC_LU *lu = (PC_LU*)pc->data;
50: lu->reusefill = flag;
51: return(0);
52: }
57: static PetscErrorCode PCSetFromOptions_LU(PC pc)
58: {
59: PC_LU *lu = (PC_LU*)pc->data;
60: PetscErrorCode ierr;
61: PetscTruth flg = PETSC_FALSE;
62: PetscReal tol;
65: PetscOptionsHead("LU options");
66: PCSetFromOptions_Factor(pc);
68: PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);
69: if (flg) {
70: tol = PETSC_DECIDE;
71: PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,0);
72: PCFactorReorderForNonzeroDiagonal(pc,tol);
73: }
74: PetscOptionsTail();
75: return(0);
76: }
80: static PetscErrorCode PCView_LU(PC pc,PetscViewer viewer)
81: {
82: PC_LU *lu = (PC_LU*)pc->data;
84: PetscTruth iascii;
87: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
88: if (iascii) {
89: if (lu->inplace) {
90: PetscViewerASCIIPrintf(viewer," LU: in-place factorization\n");
91: } else {
92: PetscViewerASCIIPrintf(viewer," LU: out-of-place factorization\n");
93: }
94:
95: if (lu->reusefill) {PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");}
96: if (lu->reuseordering) {PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");}
97: }
98: PCView_Factor(pc,viewer);
99: return(0);
100: }
104: static PetscErrorCode PCSetUp_LU(PC pc)
105: {
107: PC_LU *dir = (PC_LU*)pc->data;
110: if (dir->reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->actualfill;
112: if (dir->inplace) {
113: if (dir->row && dir->col && dir->row != dir->col) {ISDestroy(dir->row);}
114: if (dir->col) {ISDestroy(dir->col);}
115: MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);
116: if (dir->row) {
117: PetscLogObjectParent(pc,dir->row);
118: PetscLogObjectParent(pc,dir->col);
119: }
120: MatLUFactor(pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);
121: ((PC_Factor*)dir)->fact = pc->pmat;
122: } else {
123: MatInfo info;
124: if (!pc->setupcalled) {
125: MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);
126: if (dir->nonzerosalongdiagonal) {
127: MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);
128: }
129: if (dir->row) {
130: PetscLogObjectParent(pc,dir->row);
131: PetscLogObjectParent(pc,dir->col);
132: }
133: MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);
134: MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);
135: MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);
136: dir->actualfill = info.fill_ratio_needed;
137: PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);
138: } else if (pc->flag != SAME_NONZERO_PATTERN) {
139: if (!dir->reuseordering) {
140: if (dir->row && dir->col && dir->row != dir->col) {ISDestroy(dir->row);}
141: if (dir->col) {ISDestroy(dir->col);}
142: MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);
143: if (dir->nonzerosalongdiagonal) {
144: MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);
145: }
146: if (dir->row) {
147: PetscLogObjectParent(pc,dir->row);
148: PetscLogObjectParent(pc,dir->col);
149: }
150: }
151: MatDestroy(((PC_Factor*)dir)->fact);
152: MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);
153: MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);
154: MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);
155: dir->actualfill = info.fill_ratio_needed;
156: PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);
157: }
158: MatLUFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);
159: }
160: return(0);
161: }
165: static PetscErrorCode PCDestroy_LU(PC pc)
166: {
167: PC_LU *dir = (PC_LU*)pc->data;
171: if (!dir->inplace && ((PC_Factor*)dir)->fact) {MatDestroy(((PC_Factor*)dir)->fact);}
172: if (dir->row && dir->col && dir->row != dir->col) {ISDestroy(dir->row);}
173: if (dir->col) {ISDestroy(dir->col);}
174: PetscStrfree(((PC_Factor*)dir)->ordering);
175: PetscStrfree(((PC_Factor*)dir)->solvertype);
176: PetscFree(dir);
177: return(0);
178: }
182: static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y)
183: {
184: PC_LU *dir = (PC_LU*)pc->data;
188: if (dir->inplace) {MatSolve(pc->pmat,x,y);}
189: else {MatSolve(((PC_Factor*)dir)->fact,x,y);}
190: return(0);
191: }
195: static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y)
196: {
197: PC_LU *dir = (PC_LU*)pc->data;
201: if (dir->inplace) {MatSolveTranspose(pc->pmat,x,y);}
202: else {MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);}
203: return(0);
204: }
206: /* -----------------------------------------------------------------------------------*/
211: PetscErrorCode PCFactorSetUseInPlace_LU(PC pc)
212: {
213: PC_LU *dir = (PC_LU*)pc->data;
216: dir->inplace = PETSC_TRUE;
217: return(0);
218: }
221: /* ------------------------------------------------------------------------ */
223: /*MC
224: PCLU - Uses a direct solver, based on LU factorization, as a preconditioner
226: Options Database Keys:
227: + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
228: . -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like spooles
229: . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
230: . -pc_factor_fill <fill> - Sets fill amount
231: . -pc_factor_in_place - Activates in-place factorization
232: . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
233: . -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase
234: stability of factorization.
235: . -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types
236: . -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default
237: - -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the
238: diagonal.
240: Notes: Not all options work for all matrix formats
241: Run with -help to see additional options for particular matrix formats or factorization
242: algorithms
244: Level: beginner
246: Concepts: LU factorization, direct solver
248: Notes: Usually this will compute an "exact" solution in one iteration and does
249: not need a Krylov method (i.e. you can use -ksp_type preonly, or
250: KSPSetType(ksp,KSPPREONLY) for the Krylov method
252: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC,
253: PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(),
254: PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetColumnPivot(),
255: PCFactorSetPivotingInBlocks(),PCFactorSetShiftType(),PCFactorSetShiftAmount()
256: PCFactorReorderForNonzeroDiagonal()
257: M*/
262: PetscErrorCode PCCreate_LU(PC pc)
263: {
265: PetscMPIInt size;
266: PC_LU *dir;
269: PetscNewLog(pc,PC_LU,&dir);
271: MatFactorInfoInitialize(&((PC_Factor*)dir)->info);
272: ((PC_Factor*)dir)->fact = PETSC_NULL;
273: ((PC_Factor*)dir)->factortype = MAT_FACTOR_LU;
274: dir->inplace = PETSC_FALSE;
275: dir->nonzerosalongdiagonal = PETSC_FALSE;
277: ((PC_Factor*)dir)->info.fill = 5.0;
278: ((PC_Factor*)dir)->info.dtcol = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */
279: ((PC_Factor*)dir)->info.shifttype = (PetscReal) MAT_SHIFT_NONE;
280: ((PC_Factor*)dir)->info.shiftamount = 0.0;
281: ((PC_Factor*)dir)->info.zeropivot = 1.e-12;
282: ((PC_Factor*)dir)->info.pivotinblocks = 1.0;
283: dir->col = 0;
284: dir->row = 0;
286: PetscStrallocpy(MAT_SOLVER_PETSC,&((PC_Factor*)dir)->solvertype);
287: MPI_Comm_size(((PetscObject)pc)->comm,&size);
288: if (size == 1) {
289: PetscStrallocpy(MATORDERING_ND,&((PC_Factor*)dir)->ordering);
290: } else {
291: PetscStrallocpy(MATORDERING_NATURAL,&((PC_Factor*)dir)->ordering);
292: }
293: dir->reusefill = PETSC_FALSE;
294: dir->reuseordering = PETSC_FALSE;
295: pc->data = (void*)dir;
297: pc->ops->destroy = PCDestroy_LU;
298: pc->ops->apply = PCApply_LU;
299: pc->ops->applytranspose = PCApplyTranspose_LU;
300: pc->ops->setup = PCSetUp_LU;
301: pc->ops->setfromoptions = PCSetFromOptions_LU;
302: pc->ops->view = PCView_LU;
303: pc->ops->applyrichardson = 0;
304: pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor;
306: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor",
307: PCFactorGetMatSolverPackage_Factor);
308: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Factor",
309: PCFactorSetMatSolverPackage_Factor);
310: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor",
311: PCFactorSetZeroPivot_Factor);
312: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftType_C","PCFactorSetShiftType_Factor",
313: PCFactorSetShiftType_Factor);
314: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftAmount_C","PCFactorSetShiftAmount_Factor",
315: PCFactorSetShiftAmount_Factor);
316: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor",
317: PCFactorSetFill_Factor);
318: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_LU",
319: PCFactorSetUseInPlace_LU);
320: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor",
321: PCFactorSetMatOrderingType_Factor);
322: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_LU",
323: PCFactorSetReuseOrdering_LU);
324: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_LU",
325: PCFactorSetReuseFill_LU);
326: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetColumnPivot_C","PCFactorSetColumnPivot_Factor",
327: PCFactorSetColumnPivot_Factor);
328: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_Factor",
329: PCFactorSetPivotInBlocks_Factor);
330: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_LU",
331: PCFactorReorderForNonzeroDiagonal_LU);
332: return(0);
333: }