Actual source code: ilu.c

petsc-3.13.6 2020-09-29
Report Typos and Errors

  2: /*
  3:    Defines a ILU factorization preconditioner for any Mat implementation
  4: */
  5:  #include <../src/ksp/pc/impls/factor/ilu/ilu.h>

  7: PetscErrorCode  PCFactorReorderForNonzeroDiagonal_ILU(PC pc,PetscReal z)
  8: {
  9:   PC_ILU *ilu = (PC_ILU*)pc->data;

 12:   ilu->nonzerosalongdiagonal = PETSC_TRUE;
 13:   if (z == PETSC_DECIDE) ilu->nonzerosalongdiagonaltol = 1.e-10;
 14:   else ilu->nonzerosalongdiagonaltol = z;
 15:   return(0);
 16: }

 18: PetscErrorCode PCReset_ILU(PC pc)
 19: {
 20:   PC_ILU         *ilu = (PC_ILU*)pc->data;

 24:   if (!ilu->hdr.inplace) {MatDestroy(&((PC_Factor*)ilu)->fact);}
 25:   if (ilu->row && ilu->col && ilu->row != ilu->col) {ISDestroy(&ilu->row);}
 26:   ISDestroy(&ilu->col);
 27:   return(0);
 28: }

 30: PetscErrorCode  PCFactorSetDropTolerance_ILU(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount)
 31: {
 32:   PC_ILU *ilu = (PC_ILU*)pc->data;

 35:   if (pc->setupcalled && (((PC_Factor*)ilu)->info.dt != dt || ((PC_Factor*)ilu)->info.dtcol != dtcol || ((PC_Factor*)ilu)->info.dtcount != dtcount)) {
 36:     SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot change drop tolerance after using PC");
 37:   }
 38:   ((PC_Factor*)ilu)->info.dt      = dt;
 39:   ((PC_Factor*)ilu)->info.dtcol   = dtcol;
 40:   ((PC_Factor*)ilu)->info.dtcount = dtcount;
 41:   ((PC_Factor*)ilu)->info.usedt   = 1.0;
 42:   return(0);
 43: }

 45: static PetscErrorCode PCSetFromOptions_ILU(PetscOptionItems *PetscOptionsObject,PC pc)
 46: {
 48:   PetscInt       itmp;
 49:   PetscBool      flg,set;
 50:   PC_ILU         *ilu = (PC_ILU*)pc->data;
 51:   PetscReal      tol;

 54:   PetscOptionsHead(PetscOptionsObject,"ILU Options");
 55:   PCSetFromOptions_Factor(PetscOptionsObject,pc);

 57:   PetscOptionsInt("-pc_factor_levels","levels of fill","PCFactorSetLevels",(PetscInt)((PC_Factor*)ilu)->info.levels,&itmp,&flg);
 58:   if (flg) ((PC_Factor*)ilu)->info.levels = itmp;

 60:   PetscOptionsBool("-pc_factor_diagonal_fill","Allow fill into empty diagonal entry","PCFactorSetAllowDiagonalFill",((PC_Factor*)ilu)->info.diagonal_fill ? PETSC_TRUE : PETSC_FALSE,&flg,&set);
 61:   if (set) ((PC_Factor*)ilu)->info.diagonal_fill = (PetscReal) flg;
 62:   PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);
 63:   if (flg) {
 64:     tol  = PETSC_DECIDE;
 65:     PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",ilu->nonzerosalongdiagonaltol,&tol,0);
 66:     PCFactorReorderForNonzeroDiagonal(pc,tol);
 67:   }

 69:   PetscOptionsTail();
 70:   return(0);
 71: }

 73: static PetscErrorCode PCSetUp_ILU(PC pc)
 74: {
 75:   PetscErrorCode         ierr;
 76:   PC_ILU                 *ilu = (PC_ILU*)pc->data;
 77:   MatInfo                info;
 78:   PetscBool              flg;
 79:   MatSolverType          stype;
 80:   MatFactorError         err;

 83:   pc->failedreason = PC_NOERROR;
 84:   /* ugly hack to change default, since it is not support by some matrix types */
 85:   if (((PC_Factor*)ilu)->info.shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
 86:     PetscObjectTypeCompare((PetscObject)pc->pmat,MATSEQAIJ,&flg);
 87:     if (!flg) {
 88:       PetscObjectTypeCompare((PetscObject)pc->pmat,MATMPIAIJ,&flg);
 89:       if (!flg) {
 90:         ((PC_Factor*)ilu)->info.shifttype = (PetscReal)MAT_SHIFT_INBLOCKS;
 91:         PetscInfo(pc,"Changing shift type from NONZERO to INBLOCKS because block matrices do not support NONZERO\n");
 92:       }
 93:     }
 94:   }

 96:   MatSetErrorIfFailure(pc->pmat,pc->erroriffailure);
 97:   if (ilu->hdr.inplace) {
 98:     if (!pc->setupcalled) {

100:       /* In-place factorization only makes sense with the natural ordering,
101:          so we only need to get the ordering once, even if nonzero structure changes */
102:       MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);
103:       if (ilu->row) {PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->row);}
104:       if (ilu->col) {PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->col);}
105:     }

107:     /* In place ILU only makes sense with fill factor of 1.0 because
108:        cannot have levels of fill */
109:     ((PC_Factor*)ilu)->info.fill          = 1.0;
110:     ((PC_Factor*)ilu)->info.diagonal_fill = 0.0;

112:     MatILUFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);
113:     MatFactorGetError(pc->pmat,&err);
114:     if (err) { /* Factor() fails */
115:       pc->failedreason = (PCFailedReason)err;
116:       return(0);
117:     }

119:     ((PC_Factor*)ilu)->fact = pc->pmat;
120:     /* must update the pc record of the matrix state or the PC will attempt to run PCSetUp() yet again */
121:     PetscObjectStateGet((PetscObject)pc->pmat,&pc->matstate);
122:   } else {
123:     if (!pc->setupcalled) {
124:       /* first time in so compute reordering and symbolic factorization */
125:       MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);
126:       if (ilu->row) {PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->row);}
127:       if (ilu->col) {PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->col);}
128:       /*  Remove zeros along diagonal?     */
129:       if (ilu->nonzerosalongdiagonal) {
130:         MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);
131:       }
132:       if (!((PC_Factor*)ilu)->fact) {
133:         MatGetFactor(pc->pmat,((PC_Factor*)ilu)->solvertype,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);
134:       }
135:       MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);
136:       MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);
137:       ilu->hdr.actualfill = info.fill_ratio_needed;

139:       PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)ilu)->fact);
140:     } else if (pc->flag != SAME_NONZERO_PATTERN) {
141:       if (!ilu->hdr.reuseordering) {
142:         /* compute a new ordering for the ILU */
143:         ISDestroy(&ilu->row);
144:         ISDestroy(&ilu->col);
145:         MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);
146:         if (ilu->row) {PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->row);}
147:         if (ilu->col) {PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->col);}
148:         /*  Remove zeros along diagonal?     */
149:         if (ilu->nonzerosalongdiagonal) {
150:           MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);
151:         }
152:       }
153:       MatDestroy(&((PC_Factor*)ilu)->fact);
154:       MatGetFactor(pc->pmat,((PC_Factor*)ilu)->solvertype,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);
155:       MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);
156:       MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);
157:       ilu->hdr.actualfill = info.fill_ratio_needed;

159:       PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)ilu)->fact);
160:     }
161:     MatFactorGetError(((PC_Factor*)ilu)->fact,&err);
162:     if (err) { /* FactorSymbolic() fails */
163:       pc->failedreason = (PCFailedReason)err;
164:       return(0);
165:     }

167:     MatLUFactorNumeric(((PC_Factor*)ilu)->fact,pc->pmat,&((PC_Factor*)ilu)->info);
168:     MatFactorGetError(((PC_Factor*)ilu)->fact,&err);
169:     if (err) { /* FactorNumeric() fails */
170:       pc->failedreason = (PCFailedReason)err;
171:     }
172:   }

174:   PCFactorGetMatSolverType(pc,&stype);
175:   if (!stype) {
176:     MatSolverType solverpackage;
177:     MatFactorGetSolverType(((PC_Factor*)ilu)->fact,&solverpackage);
178:     PCFactorSetMatSolverType(pc,solverpackage);
179:   }
180:   return(0);
181: }

183: static PetscErrorCode PCDestroy_ILU(PC pc)
184: {
185:   PC_ILU         *ilu = (PC_ILU*)pc->data;

189:   PCReset_ILU(pc);
190:   PetscFree(((PC_Factor*)ilu)->solvertype);
191:   PetscFree(((PC_Factor*)ilu)->ordering);
192:   PetscFree(pc->data);
193:   return(0);
194: }

196: static PetscErrorCode PCApply_ILU(PC pc,Vec x,Vec y)
197: {
198:   PC_ILU         *ilu = (PC_ILU*)pc->data;

202:   MatSolve(((PC_Factor*)ilu)->fact,x,y);
203:   return(0);
204: }

206: static PetscErrorCode PCApplyTranspose_ILU(PC pc,Vec x,Vec y)
207: {
208:   PC_ILU         *ilu = (PC_ILU*)pc->data;

212:   MatSolveTranspose(((PC_Factor*)ilu)->fact,x,y);
213:   return(0);
214: }

216: static PetscErrorCode PCApplySymmetricLeft_ILU(PC pc,Vec x,Vec y)
217: {
219:   PC_ILU         *icc = (PC_ILU*)pc->data;

222:   MatForwardSolve(((PC_Factor*)icc)->fact,x,y);
223:   return(0);
224: }

226: static PetscErrorCode PCApplySymmetricRight_ILU(PC pc,Vec x,Vec y)
227: {
229:   PC_ILU         *icc = (PC_ILU*)pc->data;

232:   MatBackwardSolve(((PC_Factor*)icc)->fact,x,y);
233:   return(0);
234: }

236: /*MC
237:      PCILU - Incomplete factorization preconditioners.

239:    Options Database Keys:
240: +  -pc_factor_levels <k> - number of levels of fill for ILU(k)
241: .  -pc_factor_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for
242:                       its factorization (overwrites original matrix)
243: .  -pc_factor_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill
244: .  -pc_factor_reuse_ordering - reuse ordering of factorized matrix from previous factorization
245: .  -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
246: .  -pc_factor_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal,
247:                                    this decreases the chance of getting a zero pivot
248: .  -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
249: -  -pc_factor_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger
250:                              than 1 the diagonal blocks are factored with partial pivoting (this increases the
251:                              stability of the ILU factorization

253:    Level: beginner

255:    Notes:
256:     Only implemented for some matrix formats. (for parallel see PCHYPRE for hypre's ILU)

258:           For BAIJ matrices this implements a point block ILU

260:           The "symmetric" Section 1.5 Writing Application Codes with PETSc of this preconditioner is not actually symmetric since L is not transpose(U)
261:           even when the matrix is not symmetric since the U stores the diagonals of the factorization.

263:           If you are using MATSEQAIJCUSPARSE matrices (or MATMPIAIJCUSPARESE matrices with block Jacobi), factorization 
264:           is never done on the GPU).

266:    References:
267: +  1. - T. Dupont, R. Kendall, and H. Rachford. An approximate factorization procedure for solving
268:    self adjoint elliptic difference equations. SIAM J. Numer. Anal., 5, 1968.
269: .  2. -  T.A. Oliphant. An implicit numerical method for solving two dimensional timedependent diffusion problems. Quart. Appl. Math., 19, 1961.
270: -  3. -  TONY F. CHAN AND HENK A. VAN DER VORST, APPROXIMATE AND INCOMPLETE FACTORIZATIONS, 
271:       Chapter in Parallel Numerical
272:       Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in
273:       Science and Engineering, Kluwer.


276: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType,
277:            PCFactorSetZeroPivot(), PCFactorSetShiftSetType(), PCFactorSetAmount(),
278:            PCFactorSetDropTolerance(),PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(),
279:            PCFactorSetLevels(), PCFactorSetUseInPlace(), PCFactorSetAllowDiagonalFill(), PCFactorSetPivotInBlocks(),
280:            PCFactorGetAllowDiagonalFill(), PCFactorGetUseInPlace()

282: M*/

284: PETSC_EXTERN PetscErrorCode PCCreate_ILU(PC pc)
285: {
287:   PC_ILU         *ilu;

290:   PetscNewLog(pc,&ilu);
291:   pc->data = (void*)ilu;
292:   PCFactorInitialize(pc);

294:   ((PC_Factor*)ilu)->factortype         = MAT_FACTOR_ILU;
295:   ((PC_Factor*)ilu)->info.levels        = 0.;
296:   ((PC_Factor*)ilu)->info.fill          = 1.0;
297:   ilu->col                              = 0;
298:   ilu->row                              = 0;
299:   PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)ilu)->ordering);
300:   ((PC_Factor*)ilu)->info.dt            = PETSC_DEFAULT;
301:   ((PC_Factor*)ilu)->info.dtcount       = PETSC_DEFAULT;
302:   ((PC_Factor*)ilu)->info.dtcol         = PETSC_DEFAULT;

304:   pc->ops->reset               = PCReset_ILU;
305:   pc->ops->destroy             = PCDestroy_ILU;
306:   pc->ops->apply               = PCApply_ILU;
307:   pc->ops->applytranspose      = PCApplyTranspose_ILU;
308:   pc->ops->setup               = PCSetUp_ILU;
309:   pc->ops->setfromoptions      = PCSetFromOptions_ILU;
310:   pc->ops->view                = PCView_Factor;
311:   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_ILU;
312:   pc->ops->applysymmetricright = PCApplySymmetricRight_ILU;
313:   pc->ops->applyrichardson     = 0;
314:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetDropTolerance_C",PCFactorSetDropTolerance_ILU);
315:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",PCFactorReorderForNonzeroDiagonal_ILU);
316:   return(0);
317: }