Actual source code: cholesky.c

petsc-master 2016-07-22
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: */
  7: #include <../src/ksp/pc/impls/factor/factor.h>         /*I "petscpc.h" I*/

  9: typedef struct {
 10:   PC_Factor hdr;
 11:   PetscReal actualfill;              /* actual fill in factor */
 12:   PetscBool inplace;                 /* flag indicating in-place factorization */
 13:   IS        row,col;                 /* index sets used for reordering */
 14:   PetscBool reuseordering;           /* reuses previous reordering computed */
 15:   PetscBool reusefill;               /* reuse fill from previous Cholesky */
 16: } PC_Cholesky;

 20: static PetscErrorCode  PCFactorSetReuseOrdering_Cholesky(PC pc,PetscBool flag)
 21: {
 22:   PC_Cholesky *lu = (PC_Cholesky*)pc->data;

 25:   lu->reuseordering = flag;
 26:   return(0);
 27: }

 31: static PetscErrorCode  PCFactorSetReuseFill_Cholesky(PC pc,PetscBool flag)
 32: {
 33:   PC_Cholesky *lu = (PC_Cholesky*)pc->data;

 36:   lu->reusefill = flag;
 37:   return(0);
 38: }

 42: static PetscErrorCode PCSetFromOptions_Cholesky(PetscOptionItems *PetscOptionsObject,PC pc)
 43: {

 47:   PetscOptionsHead(PetscOptionsObject,"Cholesky options");
 48:   PCSetFromOptions_Factor(PetscOptionsObject,pc);
 49:   PetscOptionsTail();
 50:   return(0);
 51: }

 55: static PetscErrorCode PCView_Cholesky(PC pc,PetscViewer viewer)
 56: {
 57:   PC_Cholesky    *chol = (PC_Cholesky*)pc->data;
 59:   PetscBool      iascii;

 62:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
 63:   if (iascii) {
 64:     if (chol->inplace) {
 65:       PetscViewerASCIIPrintf(viewer,"  Cholesky: in-place factorization\n");
 66:     } else {
 67:       PetscViewerASCIIPrintf(viewer,"  Cholesky: out-of-place factorization\n");
 68:     }

 70:     if (chol->reusefill)    {PetscViewerASCIIPrintf(viewer,"  Reusing fill from past factorization\n");}
 71:     if (chol->reuseordering) {PetscViewerASCIIPrintf(viewer,"  Reusing reordering from past factorization\n");}
 72:   }
 73:   PCView_Factor(pc,viewer);
 74:   return(0);
 75: }

 79: static PetscErrorCode PCSetUp_Cholesky(PC pc)
 80: {
 81:   PetscErrorCode         ierr;
 82:   PetscBool              flg;
 83:   PC_Cholesky            *dir = (PC_Cholesky*)pc->data;
 84:   const MatSolverPackage stype;
 85:   MatFactorError         err;

 88:   if (dir->reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->actualfill;

 90:   MatSetErrorIfFailure(pc->pmat,pc->erroriffailure);
 91:   if (dir->inplace) {
 92:     if (dir->row && dir->col && (dir->row != dir->col)) {
 93:       ISDestroy(&dir->row);
 94:     }
 95:     ISDestroy(&dir->col);
 96:     MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);
 97:     if (dir->col && (dir->row != dir->col)) {  /* only use row ordering for SBAIJ */
 98:       ISDestroy(&dir->col);
 99:     }
100:     if (dir->row) {PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);}
101:     MatCholeskyFactor(pc->pmat,dir->row,&((PC_Factor*)dir)->info);
102:     MatFactorGetError(pc->pmat,&err);
103:     if (err) { /* Factor() fails */
104:       pc->failedreason = (PCFailedReason)err;
105:       return(0);
106:     }

108:     ((PC_Factor*)dir)->fact = pc->pmat;
109:   } else {
110:     MatInfo info;

112:     if (!pc->setupcalled) {
113:       MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);
114:       /* check if dir->row == dir->col */
115:       ISEqual(dir->row,dir->col,&flg);
116:       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"row and column permutations must equal");
117:       ISDestroy(&dir->col); /* only pass one ordering into CholeskyFactor */

119:       flg  = PETSC_FALSE;
120:       PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg,NULL);
121:       if (flg) {
122:         PetscReal tol = 1.e-10;
123:         PetscOptionsGetReal(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,NULL);
124:         MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);
125:       }
126:       if (dir->row) {PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);}
127:       if (!((PC_Factor*)dir)->fact) {
128:         MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact);
129:       }
130:       MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info);
131:       MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);
132:       dir->actualfill = info.fill_ratio_needed;
133:       PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);
134:     } else if (pc->flag != SAME_NONZERO_PATTERN) {
135:       if (!dir->reuseordering) {
136:         ISDestroy(&dir->row);
137:         MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);
138:         ISDestroy(&dir->col); /* only use dir->row ordering in CholeskyFactor */

140:         flg  = PETSC_FALSE;
141:         PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg,NULL);
142:         if (flg) {
143:           PetscReal tol = 1.e-10;
144:           PetscOptionsGetReal(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,NULL);
145:           MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);
146:         }
147:         if (dir->row) {PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);}
148:       }
149:       MatDestroy(&((PC_Factor*)dir)->fact);
150:       MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact);
151:       MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info);
152:       MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);
153:       dir->actualfill = info.fill_ratio_needed;
154:       PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);
155:     }
156:     MatFactorGetError(((PC_Factor*)dir)->fact,&err);
157:     if (err) { /* FactorSymbolic() fails */
158:       pc->failedreason = (PCFailedReason)err;
159:       return(0);
160:     }

162:     MatCholeskyFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);
163:     MatFactorGetError(((PC_Factor*)dir)->fact,&err);
164:     if (err) { /* FactorNumeric() fails */
165:       pc->failedreason = (PCFailedReason)err;
166:     }
167:   }

169:   PCFactorGetMatSolverPackage(pc,&stype);
170:   if (!stype) {
171:     const MatSolverPackage solverpackage;
172:     MatFactorGetSolverPackage(((PC_Factor*)dir)->fact,&solverpackage);
173:     PCFactorSetMatSolverPackage(pc,solverpackage);
174:   }
175:   return(0);
176: }

180: static PetscErrorCode PCReset_Cholesky(PC pc)
181: {
182:   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;

186:   if (!dir->inplace && ((PC_Factor*)dir)->fact) {MatDestroy(&((PC_Factor*)dir)->fact);}
187:   ISDestroy(&dir->row);
188:   ISDestroy(&dir->col);
189:   return(0);
190: }

194: static PetscErrorCode PCDestroy_Cholesky(PC pc)
195: {
196:   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;

200:   PCReset_Cholesky(pc);
201:   PetscFree(((PC_Factor*)dir)->ordering);
202:   PetscFree(((PC_Factor*)dir)->solvertype);
203:   PetscFree(pc->data);
204:   return(0);
205: }

209: static PetscErrorCode PCApply_Cholesky(PC pc,Vec x,Vec y)
210: {
211:   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;

215:   if (dir->inplace) {
216:     MatSolve(pc->pmat,x,y);
217:   } else {
218:     MatSolve(((PC_Factor*)dir)->fact,x,y);
219:   }
220:   return(0);
221: }

225: static PetscErrorCode PCApplyTranspose_Cholesky(PC pc,Vec x,Vec y)
226: {
227:   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;

231:   if (dir->inplace) {
232:     MatSolveTranspose(pc->pmat,x,y);
233:   } else {
234:     MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);
235:   }
236:   return(0);
237: }

239: /* -----------------------------------------------------------------------------------*/

243: static PetscErrorCode  PCFactorSetUseInPlace_Cholesky(PC pc,PetscBool flg)
244: {
245:   PC_Cholesky *dir = (PC_Cholesky*)pc->data;

248:   dir->inplace = flg;
249:   return(0);
250: }

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

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

263: /* -----------------------------------------------------------------------------------*/

267: /*@
268:    PCFactorSetReuseOrdering - When similar matrices are factored, this
269:    causes the ordering computed in the first factor to be used for all
270:    following factors.

272:    Logically Collective on PC

274:    Input Parameters:
275: +  pc - the preconditioner context
276: -  flag - PETSC_TRUE to reuse else PETSC_FALSE

278:    Options Database Key:
279: .  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()

281:    Level: intermediate

283: .keywords: PC, levels, reordering, factorization, incomplete, LU

285: .seealso: PCFactorSetReuseFill()
286: @*/
287: PetscErrorCode  PCFactorSetReuseOrdering(PC pc,PetscBool flag)
288: {

294:   PetscTryMethod(pc,"PCFactorSetReuseOrdering_C",(PC,PetscBool),(pc,flag));
295:   return(0);
296: }

298: /*MC
299:    PCCHOLESKY - Uses a direct solver, based on Cholesky factorization, as a preconditioner

301:    Options Database Keys:
302: +  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
303: .  -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like superlu
304: .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
305: .  -pc_factor_fill <fill> - Sets fill amount
306: .  -pc_factor_in_place - Activates in-place factorization
307: -  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine

309:    Notes: Not all options work for all matrix formats

311:    Level: beginner

313:    Concepts: Cholesky factorization, direct solver

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

319: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
320:            PCILU, PCLU, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(),
321:            PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetShiftType(), PCFactorSetShiftAmount()
322:            PCFactorSetUseInPlace(), PCFactorGetUseInPlace(), PCFactorSetMatOrderingType()

324: M*/

328: PETSC_EXTERN PetscErrorCode PCCreate_Cholesky(PC pc)
329: {
331:   PC_Cholesky    *dir;

334:   PetscNewLog(pc,&dir);

336:   ((PC_Factor*)dir)->fact = 0;
337:   dir->inplace            = PETSC_FALSE;

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

341:   ((PC_Factor*)dir)->factortype         = MAT_FACTOR_CHOLESKY;
342:   ((PC_Factor*)dir)->info.fill          = 5.0;
343:   ((PC_Factor*)dir)->info.shifttype     = (PetscReal) MAT_SHIFT_NONE;
344:   ((PC_Factor*)dir)->info.shiftamount   = 0.0;
345:   ((PC_Factor*)dir)->info.zeropivot     = 100.0*PETSC_MACHINE_EPSILON;
346:   ((PC_Factor*)dir)->info.pivotinblocks = 1.0;

348:   dir->col = 0;
349:   dir->row = 0;

351:   PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)dir)->ordering);
352:   dir->reusefill        = PETSC_FALSE;
353:   dir->reuseordering    = PETSC_FALSE;
354:   pc->data              = (void*)dir;

356:   pc->ops->destroy           = PCDestroy_Cholesky;
357:   pc->ops->reset             = PCReset_Cholesky;
358:   pc->ops->apply             = PCApply_Cholesky;
359:   pc->ops->applytranspose    = PCApplyTranspose_Cholesky;
360:   pc->ops->setup             = PCSetUp_Cholesky;
361:   pc->ops->setfromoptions    = PCSetFromOptions_Cholesky;
362:   pc->ops->view              = PCView_Cholesky;
363:   pc->ops->applyrichardson   = 0;
364:   pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor;

366:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUpMatSolverPackage_C",PCFactorSetUpMatSolverPackage_Factor);
367:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",PCFactorSetMatSolverPackage_Factor);
368:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",PCFactorGetMatSolverPackage_Factor);
369:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",PCFactorSetZeroPivot_Factor);
370:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftType_C",PCFactorSetShiftType_Factor);
371:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftAmount_C",PCFactorSetShiftAmount_Factor);
372:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetFill_C",PCFactorSetFill_Factor);
373:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",PCFactorSetUseInPlace_Cholesky);
374:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetUseInPlace_C",PCFactorGetUseInPlace_Cholesky);
375:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",PCFactorSetMatOrderingType_Factor);
376:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",PCFactorSetReuseOrdering_Cholesky);
377:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseFill_C",PCFactorSetReuseFill_Cholesky);
378:   return(0);
379: }