Actual source code: lu.c
petsc-3.11.4 2019-09-28
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>
11: PetscErrorCode PCFactorReorderForNonzeroDiagonal_LU(PC pc,PetscReal z)
12: {
13: PC_LU *lu = (PC_LU*)pc->data;
16: lu->nonzerosalongdiagonal = PETSC_TRUE;
17: if (z == PETSC_DECIDE) lu->nonzerosalongdiagonaltol = 1.e-10;
18: else lu->nonzerosalongdiagonaltol = z;
19: return(0);
20: }
22: static PetscErrorCode PCSetFromOptions_LU(PetscOptionItems *PetscOptionsObject,PC pc)
23: {
24: PC_LU *lu = (PC_LU*)pc->data;
26: PetscBool flg = PETSC_FALSE;
27: PetscReal tol;
30: PetscOptionsHead(PetscOptionsObject,"LU options");
31: PCSetFromOptions_Factor(PetscOptionsObject,pc);
33: PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);
34: if (flg) {
35: tol = PETSC_DECIDE;
36: PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,0);
37: PCFactorReorderForNonzeroDiagonal(pc,tol);
38: }
39: PetscOptionsTail();
40: return(0);
41: }
43: static PetscErrorCode PCView_LU(PC pc,PetscViewer viewer)
44: {
48: PCView_Factor(pc,viewer);
49: return(0);
50: }
52: static PetscErrorCode PCSetUp_LU(PC pc)
53: {
54: PetscErrorCode ierr;
55: PC_LU *dir = (PC_LU*)pc->data;
56: MatSolverType stype;
57: MatFactorError err;
60: pc->failedreason = PC_NOERROR;
61: if (dir->hdr.reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->hdr.actualfill;
63: MatSetErrorIfFailure(pc->pmat,pc->erroriffailure);
64: if (dir->hdr.inplace) {
65: MatFactorType ftype;
67: MatGetFactorType(pc->pmat, &ftype);
68: if (ftype == MAT_FACTOR_NONE) {
69: if (dir->row && dir->col && dir->row != dir->col) {ISDestroy(&dir->row);}
70: ISDestroy(&dir->col);
71: MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);
72: if (dir->row) {
73: PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);
74: PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);
75: }
76: MatLUFactor(pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);
77: MatFactorGetError(pc->pmat,&err);
78: if (err) { /* Factor() fails */
79: pc->failedreason = (PCFailedReason)err;
80: return(0);
81: }
82: }
83: ((PC_Factor*)dir)->fact = pc->pmat;
84: } else {
85: MatInfo info;
87: if (!pc->setupcalled) {
88: MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);
89: if (dir->nonzerosalongdiagonal) {
90: MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);
91: }
92: if (dir->row) {
93: PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);
94: PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);
95: }
96: if (!((PC_Factor*)dir)->fact) {
97: MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);
98: }
99: MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);
100: MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);
101: dir->hdr.actualfill = info.fill_ratio_needed;
102: PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);
103: } else if (pc->flag != SAME_NONZERO_PATTERN) {
104: if (!dir->hdr.reuseordering) {
105: if (dir->row && dir->col && dir->row != dir->col) {ISDestroy(&dir->row);}
106: ISDestroy(&dir->col);
107: MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);
108: if (dir->nonzerosalongdiagonal) {
109: MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);
110: }
111: if (dir->row) {
112: PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);
113: PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);
114: }
115: }
116: MatDestroy(&((PC_Factor*)dir)->fact);
117: MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);
118: MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);
119: MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);
120: dir->hdr.actualfill = info.fill_ratio_needed;
121: PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);
122: } else {
123: MatFactorGetError(((PC_Factor*)dir)->fact,&err);
124: if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) {
125: MatFactorClearError(((PC_Factor*)dir)->fact);
126: pc->failedreason = PC_NOERROR;
127: }
128: }
129: MatFactorGetError(((PC_Factor*)dir)->fact,&err);
130: if (err) { /* FactorSymbolic() fails */
131: pc->failedreason = (PCFailedReason)err;
132: return(0);
133: }
135: MatLUFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);
136: MatFactorGetError(((PC_Factor*)dir)->fact,&err);
137: if (err) { /* FactorNumeric() fails */
138: pc->failedreason = (PCFailedReason)err;
139: }
141: }
143: PCFactorGetMatSolverType(pc,&stype);
144: if (!stype) {
145: MatSolverType solverpackage;
146: MatFactorGetSolverType(((PC_Factor*)dir)->fact,&solverpackage);
147: PCFactorSetMatSolverType(pc,solverpackage);
148: }
149: return(0);
150: }
152: static PetscErrorCode PCReset_LU(PC pc)
153: {
154: PC_LU *dir = (PC_LU*)pc->data;
158: if (!dir->hdr.inplace && ((PC_Factor*)dir)->fact) {MatDestroy(&((PC_Factor*)dir)->fact);}
159: if (dir->row && dir->col && dir->row != dir->col) {ISDestroy(&dir->row);}
160: ISDestroy(&dir->col);
161: return(0);
162: }
164: static PetscErrorCode PCDestroy_LU(PC pc)
165: {
166: PC_LU *dir = (PC_LU*)pc->data;
170: PCReset_LU(pc);
171: PetscFree(((PC_Factor*)dir)->ordering);
172: PetscFree(((PC_Factor*)dir)->solvertype);
173: PetscFree(pc->data);
174: return(0);
175: }
177: static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y)
178: {
179: PC_LU *dir = (PC_LU*)pc->data;
183: if (dir->hdr.inplace) {
184: MatSolve(pc->pmat,x,y);
185: } else {
186: MatSolve(((PC_Factor*)dir)->fact,x,y);
187: }
188: return(0);
189: }
191: static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y)
192: {
193: PC_LU *dir = (PC_LU*)pc->data;
197: if (dir->hdr.inplace) {
198: MatSolveTranspose(pc->pmat,x,y);
199: } else {
200: MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);
201: }
202: return(0);
203: }
205: /* -----------------------------------------------------------------------------------*/
207: /*MC
208: PCLU - Uses a direct solver, based on LU factorization, as a preconditioner
210: Options Database Keys:
211: + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
212: . -pc_factor_mat_solver_type - Actives PCFactorSetMatSolverType() to choose the direct solver, like superlu
213: . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
214: . -pc_factor_fill <fill> - Sets fill amount
215: . -pc_factor_in_place - Activates in-place factorization
216: . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
217: . -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase
218: stability of factorization.
219: . -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types
220: . -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default
221: - -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the
222: diagonal.
224: Notes:
225: Not all options work for all matrix formats
226: Run with -help to see additional options for particular matrix formats or factorization
227: algorithms
229: Level: beginner
231: Concepts: LU factorization, direct solver
233: Notes:
234: Usually this will compute an "exact" solution in one iteration and does
235: not need a Krylov method (i.e. you can use -ksp_type preonly, or
236: KSPSetType(ksp,KSPPREONLY) for the Krylov method
238: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC,
239: PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(),
240: PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetColumnPivot(),
241: PCFactorSetPivotingInBlocks(),PCFactorSetShiftType(),PCFactorSetShiftAmount()
242: PCFactorReorderForNonzeroDiagonal()
243: M*/
245: PETSC_EXTERN PetscErrorCode PCCreate_LU(PC pc)
246: {
248: PetscMPIInt size;
249: PC_LU *dir;
252: PetscNewLog(pc,&dir);
253: pc->data = (void*)dir;
254: PCFactorInitialize(pc);
255: ((PC_Factor*)dir)->factortype = MAT_FACTOR_LU;
256: dir->nonzerosalongdiagonal = PETSC_FALSE;
258: ((PC_Factor*)dir)->info.fill = 5.0;
259: ((PC_Factor*)dir)->info.dtcol = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */
260: ((PC_Factor*)dir)->info.shifttype = (PetscReal)MAT_SHIFT_NONE;
261: dir->col = 0;
262: dir->row = 0;
264: MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);
265: if (size == 1) {
266: PetscStrallocpy(MATORDERINGND,(char**)&((PC_Factor*)dir)->ordering);
267: } else {
268: PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)dir)->ordering);
269: }
271: pc->ops->reset = PCReset_LU;
272: pc->ops->destroy = PCDestroy_LU;
273: pc->ops->apply = PCApply_LU;
274: pc->ops->applytranspose = PCApplyTranspose_LU;
275: pc->ops->setup = PCSetUp_LU;
276: pc->ops->setfromoptions = PCSetFromOptions_LU;
277: pc->ops->view = PCView_LU;
278: pc->ops->applyrichardson = 0;
279: PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",PCFactorReorderForNonzeroDiagonal_LU);
280: return(0);
281: }