Actual source code: lu.c

petsc-3.7.5 2017-01-01
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>  /*I "petscpc.h" I*/


 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;

102:   pc->failedreason = PC_NOERROR;
103:   if (dir->reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->actualfill;

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

120:     ((PC_Factor*)dir)->fact = pc->pmat;
121:   } else {
122:     MatInfo info;
123:     Mat     F;
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:       F = ((PC_Factor*)dir)->fact;
161:       if ((PCFailedReason)F->errortype == PC_FACTOR_NUMERIC_ZEROPIVOT) {
162:         F->errortype     = MAT_FACTOR_NOERROR;
163:         pc->failedreason = (PCFailedReason)F->errortype;
164:       }
165:     }
166:     F = ((PC_Factor*)dir)->fact;
167:     if (F->errortype) { /* FactorSymbolic() fails */
168:       pc->failedreason = (PCFailedReason)F->errortype;
169:       return(0);
170:     }

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

177:   }

179:   PCFactorGetMatSolverPackage(pc,&stype);
180:   if (!stype) {
181:     PCFactorSetMatSolverPackage(pc,((PC_Factor*)dir)->fact->solvertype);
182:   }
183:   return(0);
184: }

188: static PetscErrorCode PCReset_LU(PC pc)
189: {
190:   PC_LU          *dir = (PC_LU*)pc->data;

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

202: static PetscErrorCode PCDestroy_LU(PC pc)
203: {
204:   PC_LU          *dir = (PC_LU*)pc->data;

208:   PCReset_LU(pc);
209:   PetscFree(((PC_Factor*)dir)->ordering);
210:   PetscFree(((PC_Factor*)dir)->solvertype);
211:   PetscFree(pc->data);
212:   return(0);
213: }

217: static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y)
218: {
219:   PC_LU          *dir = (PC_LU*)pc->data;

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

233: static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y)
234: {
235:   PC_LU          *dir = (PC_LU*)pc->data;

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

247: /* -----------------------------------------------------------------------------------*/

251: PetscErrorCode  PCFactorSetUseInPlace_LU(PC pc,PetscBool flg)
252: {
253:   PC_LU *dir = (PC_LU*)pc->data;

256:   dir->inplace = flg;
257:   return(0);
258: }

262: PetscErrorCode  PCFactorGetUseInPlace_LU(PC pc,PetscBool *flg)
263: {
264:   PC_LU *dir = (PC_LU*)pc->data;

267:   *flg = dir->inplace;
268:   return(0);
269: }

271: /* ------------------------------------------------------------------------ */

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

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

290:    Notes: Not all options work for all matrix formats
291:           Run with -help to see additional options for particular matrix formats or factorization
292:           algorithms

294:    Level: beginner

296:    Concepts: LU factorization, direct solver

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

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

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

318:   PetscNewLog(pc,&dir);

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

322:   ((PC_Factor*)dir)->fact       = NULL;
323:   ((PC_Factor*)dir)->factortype = MAT_FACTOR_LU;
324:   dir->inplace                  = PETSC_FALSE;
325:   dir->nonzerosalongdiagonal    = PETSC_FALSE;

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

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

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

356:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",PCFactorSetZeroPivot_Factor);
357:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetZeroPivot_C",PCFactorGetZeroPivot_Factor);
358:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftType_C",PCFactorSetShiftType_Factor);
359:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetShiftType_C",PCFactorGetShiftType_Factor);
360:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftAmount_C",PCFactorSetShiftAmount_Factor);
361:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetShiftAmount_C",PCFactorGetShiftAmount_Factor);
362:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUpMatSolverPackage_C",PCFactorSetUpMatSolverPackage_Factor);
363:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",PCFactorGetMatSolverPackage_Factor);
364:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",PCFactorSetMatSolverPackage_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: }