Actual source code: lu.c
petsc-3.8.4 2018-03-24
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: const MatSolverPackage 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: if (dir->row && dir->col && dir->row != dir->col) {ISDestroy(&dir->row);}
66: ISDestroy(&dir->col);
67: MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);
68: if (dir->row) {
69: PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);
70: PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);
71: }
72: MatLUFactor(pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);
73: MatFactorGetError(pc->pmat,&err);
74: if (err) { /* Factor() fails */
75: pc->failedreason = (PCFailedReason)err;
76: return(0);
77: }
79: ((PC_Factor*)dir)->fact = pc->pmat;
80: } else {
81: MatInfo info;
83: if (!pc->setupcalled) {
84: MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);
85: if (dir->nonzerosalongdiagonal) {
86: MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);
87: }
88: if (dir->row) {
89: PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);
90: PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);
91: }
92: if (!((PC_Factor*)dir)->fact) {
93: MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);
94: }
95: MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);
96: MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);
97: dir->hdr.actualfill = info.fill_ratio_needed;
98: PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);
99: } else if (pc->flag != SAME_NONZERO_PATTERN) {
100: if (!dir->hdr.reuseordering) {
101: if (dir->row && dir->col && dir->row != dir->col) {ISDestroy(&dir->row);}
102: ISDestroy(&dir->col);
103: MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);
104: if (dir->nonzerosalongdiagonal) {
105: MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);
106: }
107: if (dir->row) {
108: PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);
109: PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);
110: }
111: }
112: MatDestroy(&((PC_Factor*)dir)->fact);
113: MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);
114: MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);
115: MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);
116: dir->hdr.actualfill = info.fill_ratio_needed;
117: PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);
118: } else {
119: MatFactorGetError(((PC_Factor*)dir)->fact,&err);
120: if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) {
121: MatFactorClearError(((PC_Factor*)dir)->fact);
122: pc->failedreason = PC_NOERROR;
123: }
124: }
125: MatFactorGetError(((PC_Factor*)dir)->fact,&err);
126: if (err) { /* FactorSymbolic() fails */
127: pc->failedreason = (PCFailedReason)err;
128: return(0);
129: }
131: MatLUFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);
132: MatFactorGetError(((PC_Factor*)dir)->fact,&err);
133: if (err) { /* FactorNumeric() fails */
134: pc->failedreason = (PCFailedReason)err;
135: }
137: }
139: PCFactorGetMatSolverPackage(pc,&stype);
140: if (!stype) {
141: const MatSolverPackage solverpackage;
142: MatFactorGetSolverPackage(((PC_Factor*)dir)->fact,&solverpackage);
143: PCFactorSetMatSolverPackage(pc,solverpackage);
144: }
145: return(0);
146: }
148: static PetscErrorCode PCReset_LU(PC pc)
149: {
150: PC_LU *dir = (PC_LU*)pc->data;
154: if (!dir->hdr.inplace && ((PC_Factor*)dir)->fact) {MatDestroy(&((PC_Factor*)dir)->fact);}
155: if (dir->row && dir->col && dir->row != dir->col) {ISDestroy(&dir->row);}
156: ISDestroy(&dir->col);
157: return(0);
158: }
160: static PetscErrorCode PCDestroy_LU(PC pc)
161: {
162: PC_LU *dir = (PC_LU*)pc->data;
166: PCReset_LU(pc);
167: PetscFree(((PC_Factor*)dir)->ordering);
168: PetscFree(((PC_Factor*)dir)->solvertype);
169: PetscFree(pc->data);
170: return(0);
171: }
173: static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y)
174: {
175: PC_LU *dir = (PC_LU*)pc->data;
179: if (dir->hdr.inplace) {
180: MatSolve(pc->pmat,x,y);
181: } else {
182: MatSolve(((PC_Factor*)dir)->fact,x,y);
183: }
184: return(0);
185: }
187: static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y)
188: {
189: PC_LU *dir = (PC_LU*)pc->data;
193: if (dir->hdr.inplace) {
194: MatSolveTranspose(pc->pmat,x,y);
195: } else {
196: MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);
197: }
198: return(0);
199: }
201: /* -----------------------------------------------------------------------------------*/
203: /*MC
204: PCLU - Uses a direct solver, based on LU factorization, as a preconditioner
206: Options Database Keys:
207: + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
208: . -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like superlu
209: . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
210: . -pc_factor_fill <fill> - Sets fill amount
211: . -pc_factor_in_place - Activates in-place factorization
212: . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
213: . -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase
214: stability of factorization.
215: . -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types
216: . -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default
217: - -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the
218: diagonal.
220: Notes: Not all options work for all matrix formats
221: Run with -help to see additional options for particular matrix formats or factorization
222: algorithms
224: Level: beginner
226: Concepts: LU factorization, direct solver
228: Notes: Usually this will compute an "exact" solution in one iteration and does
229: not need a Krylov method (i.e. you can use -ksp_type preonly, or
230: KSPSetType(ksp,KSPPREONLY) for the Krylov method
232: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC,
233: PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(),
234: PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetColumnPivot(),
235: PCFactorSetPivotingInBlocks(),PCFactorSetShiftType(),PCFactorSetShiftAmount()
236: PCFactorReorderForNonzeroDiagonal()
237: M*/
239: PETSC_EXTERN PetscErrorCode PCCreate_LU(PC pc)
240: {
242: PetscMPIInt size;
243: PC_LU *dir;
246: PetscNewLog(pc,&dir);
247: pc->data = (void*)dir;
248: PCFactorInitialize(pc);
249: ((PC_Factor*)dir)->factortype = MAT_FACTOR_LU;
250: dir->nonzerosalongdiagonal = PETSC_FALSE;
252: ((PC_Factor*)dir)->info.fill = 5.0;
253: ((PC_Factor*)dir)->info.dtcol = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */
254: ((PC_Factor*)dir)->info.shifttype = (PetscReal)MAT_SHIFT_NONE;
255: dir->col = 0;
256: dir->row = 0;
258: MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);
259: if (size == 1) {
260: PetscStrallocpy(MATORDERINGND,(char**)&((PC_Factor*)dir)->ordering);
261: } else {
262: PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)dir)->ordering);
263: }
265: pc->ops->reset = PCReset_LU;
266: pc->ops->destroy = PCDestroy_LU;
267: pc->ops->apply = PCApply_LU;
268: pc->ops->applytranspose = PCApplyTranspose_LU;
269: pc->ops->setup = PCSetUp_LU;
270: pc->ops->setfromoptions = PCSetFromOptions_LU;
271: pc->ops->view = PCView_LU;
272: pc->ops->applyrichardson = 0;
273: PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",PCFactorReorderForNonzeroDiagonal_LU);
274: return(0);
275: }