Actual source code: lu.c

petsc-master 2016-07-27
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;
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:     }
160:     MatFactorGetError(((PC_Factor*)dir)->fact,&err);
161:     if (err) { /* FactorSymbolic() fails */
162:       pc->failedreason = (PCFailedReason)err;
163:       return(0);
164:     }

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

172:   }

174:   PCFactorGetMatSolverPackage(pc,&stype);
175:   if (!stype) {
176:     const MatSolverPackage solverpackage;
177:     MatFactorGetSolverPackage(((PC_Factor*)dir)->fact,&solverpackage);
178:     PCFactorSetMatSolverPackage(pc,solverpackage);
179:   }
180:   return(0);
181: }

185: static PetscErrorCode PCReset_LU(PC pc)
186: {
187:   PC_LU          *dir = (PC_LU*)pc->data;

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

199: static PetscErrorCode PCDestroy_LU(PC pc)
200: {
201:   PC_LU          *dir = (PC_LU*)pc->data;

205:   PCReset_LU(pc);
206:   PetscFree(((PC_Factor*)dir)->ordering);
207:   PetscFree(((PC_Factor*)dir)->solvertype);
208:   PetscFree(pc->data);
209:   return(0);
210: }

214: static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y)
215: {
216:   PC_LU          *dir = (PC_LU*)pc->data;

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

230: static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y)
231: {
232:   PC_LU          *dir = (PC_LU*)pc->data;

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

244: /* -----------------------------------------------------------------------------------*/

248: PetscErrorCode  PCFactorSetUseInPlace_LU(PC pc,PetscBool flg)
249: {
250:   PC_LU *dir = (PC_LU*)pc->data;

253:   dir->inplace = flg;
254:   return(0);
255: }

259: PetscErrorCode  PCFactorGetUseInPlace_LU(PC pc,PetscBool *flg)
260: {
261:   PC_LU *dir = (PC_LU*)pc->data;

264:   *flg = dir->inplace;
265:   return(0);
266: }

268: /* ------------------------------------------------------------------------ */

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

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

287:    Notes: Not all options work for all matrix formats
288:           Run with -help to see additional options for particular matrix formats or factorization
289:           algorithms

291:    Level: beginner

293:    Concepts: LU factorization, direct solver

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

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

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

315:   PetscNewLog(pc,&dir);

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

319:   ((PC_Factor*)dir)->fact       = NULL;
320:   ((PC_Factor*)dir)->factortype = MAT_FACTOR_LU;
321:   dir->inplace                  = PETSC_FALSE;
322:   dir->nonzerosalongdiagonal    = PETSC_FALSE;

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

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

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

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