Actual source code: lu.c

petsc-master 2016-08-28
Report Typos and Errors
  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>


 13: PetscErrorCode PCFactorReorderForNonzeroDiagonal_LU(PC pc,PetscReal z)
 14: {
 15:   PC_LU *lu = (PC_LU*)pc->data;

 18:   lu->nonzerosalongdiagonal = PETSC_TRUE;
 19:   if (z == PETSC_DECIDE) lu->nonzerosalongdiagonaltol = 1.e-10;
 20:   else lu->nonzerosalongdiagonaltol = z;
 21:   return(0);
 22: }

 26: PetscErrorCode PCFactorSetReuseOrdering_LU(PC pc,PetscBool flag)
 27: {
 28:   PC_LU *lu = (PC_LU*)pc->data;

 31:   lu->reuseordering = flag;
 32:   return(0);
 33: }

 37: PetscErrorCode PCFactorSetReuseFill_LU(PC pc,PetscBool flag)
 38: {
 39:   PC_LU *lu = (PC_LU*)pc->data;

 42:   lu->reusefill = flag;
 43:   return(0);
 44: }

 48: static PetscErrorCode PCSetFromOptions_LU(PetscOptionItems *PetscOptionsObject,PC pc)
 49: {
 50:   PC_LU          *lu = (PC_LU*)pc->data;
 52:   PetscBool      flg = PETSC_FALSE;
 53:   PetscReal      tol;

 56:   PetscOptionsHead(PetscOptionsObject,"LU options");
 57:   PCSetFromOptions_Factor(PetscOptionsObject,pc);

 59:   PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);
 60:   if (flg) {
 61:     tol  = PETSC_DECIDE;
 62:     PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,0);
 63:     PCFactorReorderForNonzeroDiagonal(pc,tol);
 64:   }
 65:   PetscOptionsTail();
 66:   return(0);
 67: }

 71: static PetscErrorCode PCView_LU(PC pc,PetscViewer viewer)
 72: {
 73:   PC_LU          *lu = (PC_LU*)pc->data;
 75:   PetscBool      iascii;

 78:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
 79:   if (iascii) {
 80:     if (lu->inplace) {
 81:       PetscViewerASCIIPrintf(viewer,"  LU: in-place factorization\n");
 82:     } else {
 83:       PetscViewerASCIIPrintf(viewer,"  LU: out-of-place factorization\n");
 84:     }

 86:     if (lu->reusefill)    {PetscViewerASCIIPrintf(viewer,"       Reusing fill from past factorization\n");}
 87:     if (lu->reuseordering) {PetscViewerASCIIPrintf(viewer,"       Reusing reordering from past factorization\n");}
 88:   }
 89:   PCView_Factor(pc,viewer);
 90:   return(0);
 91: }

 95: static PetscErrorCode PCSetUp_LU(PC pc)
 96: {
 97:   PetscErrorCode         ierr;
 98:   PC_LU                  *dir = (PC_LU*)pc->data;
 99:   const MatSolverPackage stype;
100:   MatFactorError         err;
102:   if (dir->reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->actualfill;

104:   MatSetErrorIfFailure(pc->pmat,pc->erroriffailure);
105:   if (dir->inplace) {
106:     if (dir->row && dir->col && dir->row != dir->col) {ISDestroy(&dir->row);}
107:     ISDestroy(&dir->col);
108:     MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);
109:     if (dir->row) {
110:       PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);
111:       PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);
112:     }
113:     MatLUFactor(pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);
114:     MatFactorGetError(pc->pmat,&err);
115:     if (err) { /* Factor() fails */
116:       pc->failedreason = (PCFailedReason)err;
117:       return(0);
118:     }

120:     ((PC_Factor*)dir)->fact = pc->pmat;
121:   } else {
122:     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((PetscObject)pc,(PetscObject)dir->row);
131:         PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);
132:       }
133:       if (!((PC_Factor*)dir)->fact) {
134:         MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);
135:       }
136:       MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);
137:       MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);
138:       dir->actualfill = info.fill_ratio_needed;
139:       PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);
140:     } else if (pc->flag != SAME_NONZERO_PATTERN) {
141:       if (!dir->reuseordering) {
142:         if (dir->row && dir->col && dir->row != dir->col) {ISDestroy(&dir->row);}
143:         ISDestroy(&dir->col);
144:         MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);
145:         if (dir->nonzerosalongdiagonal) {
146:           MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);
147:         }
148:         if (dir->row) {
149:           PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);
150:           PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);
151:         }
152:       }
153:       MatDestroy(&((PC_Factor*)dir)->fact);
154:       MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);
155:       MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);
156:       MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);
157:       dir->actualfill = info.fill_ratio_needed;
158:       PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);
159:     } else {
160:       MatFactorGetError(((PC_Factor*)dir)->fact,&err);
161:       if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) {
162:         MatFactorClearError(((PC_Factor*)dir)->fact);
163:         pc->failedreason = PC_NOERROR;
164:       }
165:     }
166:     MatFactorGetError(((PC_Factor*)dir)->fact,&err);
167:     if (err) { /* FactorSymbolic() fails */
168:       pc->failedreason = (PCFailedReason)err;
169:       return(0);
170:     }

172:     MatLUFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);
173:     MatFactorGetError(((PC_Factor*)dir)->fact,&err);
174:     if (err) { /* FactorNumeric() fails */
175:       pc->failedreason = (PCFailedReason)err;
176:     }

178:   }

180:   PCFactorGetMatSolverPackage(pc,&stype);
181:   if (!stype) {
182:     const MatSolverPackage solverpackage;
183:     MatFactorGetSolverPackage(((PC_Factor*)dir)->fact,&solverpackage);
184:     PCFactorSetMatSolverPackage(pc,solverpackage);
185:   }
186:   return(0);
187: }

191: static PetscErrorCode PCReset_LU(PC pc)
192: {
193:   PC_LU          *dir = (PC_LU*)pc->data;

197:   if (!dir->inplace && ((PC_Factor*)dir)->fact) {MatDestroy(&((PC_Factor*)dir)->fact);}
198:   if (dir->row && dir->col && dir->row != dir->col) {ISDestroy(&dir->row);}
199:   ISDestroy(&dir->col);
200:   return(0);
201: }

205: static PetscErrorCode PCDestroy_LU(PC pc)
206: {
207:   PC_LU          *dir = (PC_LU*)pc->data;

211:   PCReset_LU(pc);
212:   PetscFree(((PC_Factor*)dir)->ordering);
213:   PetscFree(((PC_Factor*)dir)->solvertype);
214:   PetscFree(pc->data);
215:   return(0);
216: }

220: static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y)
221: {
222:   PC_LU          *dir = (PC_LU*)pc->data;

226:   if (dir->inplace) {
227:     MatSolve(pc->pmat,x,y);
228:   } else {
229:     MatSolve(((PC_Factor*)dir)->fact,x,y);
230:   }
231:   return(0);
232: }

236: static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y)
237: {
238:   PC_LU          *dir = (PC_LU*)pc->data;

242:   if (dir->inplace) {
243:     MatSolveTranspose(pc->pmat,x,y);
244:   } else {
245:     MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);
246:   }
247:   return(0);
248: }

250: /* -----------------------------------------------------------------------------------*/

254: PetscErrorCode  PCFactorSetUseInPlace_LU(PC pc,PetscBool flg)
255: {
256:   PC_LU *dir = (PC_LU*)pc->data;

259:   dir->inplace = flg;
260:   return(0);
261: }

265: PetscErrorCode  PCFactorGetUseInPlace_LU(PC pc,PetscBool *flg)
266: {
267:   PC_LU *dir = (PC_LU*)pc->data;

270:   *flg = dir->inplace;
271:   return(0);
272: }

274: /* ------------------------------------------------------------------------ */

276: /*MC
277:    PCLU - Uses a direct solver, based on LU factorization, as a preconditioner

279:    Options Database Keys:
280: +  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
281: .  -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like superlu
282: .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
283: .  -pc_factor_fill <fill> - Sets fill amount
284: .  -pc_factor_in_place - Activates in-place factorization
285: .  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
286: .  -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase
287:                                          stability of factorization.
288: .  -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types
289: .  -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default
290: -   -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the
291:         diagonal.

293:    Notes: Not all options work for all matrix formats
294:           Run with -help to see additional options for particular matrix formats or factorization
295:           algorithms

297:    Level: beginner

299:    Concepts: LU factorization, direct solver

301:    Notes: Usually this will compute an "exact" solution in one iteration and does
302:           not need a Krylov method (i.e. you can use -ksp_type preonly, or
303:           KSPSetType(ksp,KSPPREONLY) for the Krylov method

305: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
306:            PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(),
307:            PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetColumnPivot(),
308:            PCFactorSetPivotingInBlocks(),PCFactorSetShiftType(),PCFactorSetShiftAmount()
309:            PCFactorReorderForNonzeroDiagonal()
310: M*/

314: PETSC_EXTERN PetscErrorCode PCCreate_LU(PC pc)
315: {
317:   PetscMPIInt    size;
318:   PC_LU          *dir;

321:   PetscNewLog(pc,&dir);

323:   MatFactorInfoInitialize(&((PC_Factor*)dir)->info);

325:   ((PC_Factor*)dir)->fact       = NULL;
326:   ((PC_Factor*)dir)->factortype = MAT_FACTOR_LU;
327:   dir->inplace                  = PETSC_FALSE;
328:   dir->nonzerosalongdiagonal    = PETSC_FALSE;

330:   ((PC_Factor*)dir)->info.fill          = 5.0;
331:   ((PC_Factor*)dir)->info.dtcol         = 1.e-6;  /* default to pivoting; this is only thing PETSc LU supports */
332:   ((PC_Factor*)dir)->info.shifttype     = (PetscReal)MAT_SHIFT_NONE;
333:   ((PC_Factor*)dir)->info.shiftamount   = 0.0;
334:   ((PC_Factor*)dir)->info.zeropivot     = 100.0*PETSC_MACHINE_EPSILON;
335:   ((PC_Factor*)dir)->info.pivotinblocks = 1.0;
336:   dir->col                              = 0;
337:   dir->row                              = 0;

339:   MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);
340:   if (size == 1) {
341:     PetscStrallocpy(MATORDERINGND,(char**)&((PC_Factor*)dir)->ordering);
342:   } else {
343:     PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)dir)->ordering);
344:   }
345:   dir->reusefill     = PETSC_FALSE;
346:   dir->reuseordering = PETSC_FALSE;
347:   pc->data           = (void*)dir;

349:   pc->ops->reset             = PCReset_LU;
350:   pc->ops->destroy           = PCDestroy_LU;
351:   pc->ops->apply             = PCApply_LU;
352:   pc->ops->applytranspose    = PCApplyTranspose_LU;
353:   pc->ops->setup             = PCSetUp_LU;
354:   pc->ops->setfromoptions    = PCSetFromOptions_LU;
355:   pc->ops->view              = PCView_LU;
356:   pc->ops->applyrichardson   = 0;
357:   pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor;

359:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUpMatSolverPackage_C",PCFactorSetUpMatSolverPackage_Factor);
360:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",PCFactorGetMatSolverPackage_Factor);
361:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",PCFactorSetMatSolverPackage_Factor);
362:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",PCFactorSetZeroPivot_Factor);
363:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftType_C",PCFactorSetShiftType_Factor);
364:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftAmount_C",PCFactorSetShiftAmount_Factor);
365:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetFill_C",PCFactorSetFill_Factor);
366:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",PCFactorSetUseInPlace_LU);
367:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetUseInPlace_C",PCFactorGetUseInPlace_LU);
368:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",PCFactorSetMatOrderingType_Factor);
369:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",PCFactorSetReuseOrdering_LU);
370:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseFill_C",PCFactorSetReuseFill_LU);
371:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetColumnPivot_C",PCFactorSetColumnPivot_Factor);
372:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",PCFactorSetPivotInBlocks_Factor);
373:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",PCFactorReorderForNonzeroDiagonal_LU);
374:   return(0);
375: }