Actual source code: cholesky.c


  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>

  9: typedef struct {
 10:   PC_Factor hdr;
 11:   IS        row, col; /* index sets used for reordering */
 12: } PC_Cholesky;

 14: static PetscErrorCode PCSetFromOptions_Cholesky(PC pc, PetscOptionItems *PetscOptionsObject)
 15: {
 16:   PetscFunctionBegin;
 17:   PetscOptionsHeadBegin(PetscOptionsObject, "Cholesky options");
 18:   PetscCall(PCSetFromOptions_Factor(pc, PetscOptionsObject));
 19:   PetscOptionsHeadEnd();
 20:   PetscFunctionReturn(PETSC_SUCCESS);
 21: }

 23: static PetscErrorCode PCSetUp_Cholesky(PC pc)
 24: {
 25:   PetscBool      flg;
 26:   PC_Cholesky   *dir = (PC_Cholesky *)pc->data;
 27:   MatSolverType  stype;
 28:   MatFactorError err;
 29:   const char    *prefix;

 31:   PetscFunctionBegin;
 32:   pc->failedreason = PC_NOERROR;
 33:   if (dir->hdr.reusefill && pc->setupcalled) ((PC_Factor *)dir)->info.fill = dir->hdr.actualfill;

 35:   PetscCall(PCGetOptionsPrefix(pc, &prefix));
 36:   PetscCall(MatSetOptionsPrefixFactor(pc->pmat, prefix));

 38:   PetscCall(MatSetErrorIfFailure(pc->pmat, pc->erroriffailure));
 39:   if (dir->hdr.inplace) {
 40:     if (dir->row && dir->col && (dir->row != dir->col)) PetscCall(ISDestroy(&dir->row));
 41:     PetscCall(ISDestroy(&dir->col));
 42:     /* should only get reordering if the factor matrix uses it but cannot determine because MatGetFactor() not called */
 43:     PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
 44:     PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col));
 45:     if (dir->col && (dir->row != dir->col)) { /* only use row ordering for SBAIJ */
 46:       PetscCall(ISDestroy(&dir->col));
 47:     }
 48:     PetscCall(MatCholeskyFactor(pc->pmat, dir->row, &((PC_Factor *)dir)->info));
 49:     PetscCall(MatFactorGetError(pc->pmat, &err));
 50:     if (err) { /* Factor() fails */
 51:       pc->failedreason = (PCFailedReason)err;
 52:       PetscFunctionReturn(PETSC_SUCCESS);
 53:     }

 55:     ((PC_Factor *)dir)->fact = pc->pmat;
 56:   } else {
 57:     MatInfo info;

 59:     if (!pc->setupcalled) {
 60:       PetscBool canuseordering;
 61:       if (!((PC_Factor *)dir)->fact) { PetscCall(MatGetFactor(pc->pmat, ((PC_Factor *)dir)->solvertype, MAT_FACTOR_CHOLESKY, &((PC_Factor *)dir)->fact)); }
 62:       PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)dir)->fact, &canuseordering));
 63:       if (canuseordering) {
 64:         PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
 65:         PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col));
 66:         /* check if dir->row == dir->col */
 67:         if (dir->row) {
 68:           PetscCall(ISEqual(dir->row, dir->col, &flg));
 69:           PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "row and column permutations must be equal");
 70:         }
 71:         PetscCall(ISDestroy(&dir->col)); /* only pass one ordering into CholeskyFactor */

 73:         flg = PETSC_FALSE;
 74:         PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_factor_nonzeros_along_diagonal", &flg, NULL));
 75:         if (flg) {
 76:           PetscReal tol = 1.e-10;
 77:           PetscCall(PetscOptionsGetReal(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_factor_nonzeros_along_diagonal", &tol, NULL));
 78:           PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, tol, dir->row, dir->row));
 79:         }
 80:       }
 81:       PetscCall(MatCholeskyFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->row, &((PC_Factor *)dir)->info));
 82:       PetscCall(MatGetInfo(((PC_Factor *)dir)->fact, MAT_LOCAL, &info));
 83:       dir->hdr.actualfill = info.fill_ratio_needed;
 84:     } else if (pc->flag != SAME_NONZERO_PATTERN) {
 85:       if (!dir->hdr.reuseordering) {
 86:         PetscBool canuseordering;
 87:         PetscCall(MatDestroy(&((PC_Factor *)dir)->fact));
 88:         PetscCall(MatGetFactor(pc->pmat, ((PC_Factor *)dir)->solvertype, MAT_FACTOR_CHOLESKY, &((PC_Factor *)dir)->fact));
 89:         PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)dir)->fact, &canuseordering));
 90:         if (canuseordering) {
 91:           PetscCall(ISDestroy(&dir->row));
 92:           PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
 93:           PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col));
 94:           PetscCall(ISDestroy(&dir->col)); /* only use dir->row ordering in CholeskyFactor */

 96:           flg = PETSC_FALSE;
 97:           PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_factor_nonzeros_along_diagonal", &flg, NULL));
 98:           if (flg) {
 99:             PetscReal tol = 1.e-10;
100:             PetscCall(PetscOptionsGetReal(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_factor_nonzeros_along_diagonal", &tol, NULL));
101:             PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, tol, dir->row, dir->row));
102:           }
103:         }
104:       }
105:       PetscCall(MatCholeskyFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->row, &((PC_Factor *)dir)->info));
106:       PetscCall(MatGetInfo(((PC_Factor *)dir)->fact, MAT_LOCAL, &info));
107:       dir->hdr.actualfill = info.fill_ratio_needed;
108:     } else {
109:       PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err));
110:       if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) {
111:         PetscCall(MatFactorClearError(((PC_Factor *)dir)->fact));
112:         pc->failedreason = PC_NOERROR;
113:       }
114:     }
115:     PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err));
116:     if (err) { /* FactorSymbolic() fails */
117:       pc->failedreason = (PCFailedReason)err;
118:       PetscFunctionReturn(PETSC_SUCCESS);
119:     }

121:     PetscCall(MatCholeskyFactorNumeric(((PC_Factor *)dir)->fact, pc->pmat, &((PC_Factor *)dir)->info));
122:     PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err));
123:     if (err) { /* FactorNumeric() fails */
124:       pc->failedreason = (PCFailedReason)err;
125:     }
126:   }

128:   PetscCall(PCFactorGetMatSolverType(pc, &stype));
129:   if (!stype) {
130:     MatSolverType solverpackage;
131:     PetscCall(MatFactorGetSolverType(((PC_Factor *)dir)->fact, &solverpackage));
132:     PetscCall(PCFactorSetMatSolverType(pc, solverpackage));
133:   }
134:   PetscFunctionReturn(PETSC_SUCCESS);
135: }

137: static PetscErrorCode PCReset_Cholesky(PC pc)
138: {
139:   PC_Cholesky *dir = (PC_Cholesky *)pc->data;

141:   PetscFunctionBegin;
142:   if (!dir->hdr.inplace && ((PC_Factor *)dir)->fact) PetscCall(MatDestroy(&((PC_Factor *)dir)->fact));
143:   PetscCall(ISDestroy(&dir->row));
144:   PetscCall(ISDestroy(&dir->col));
145:   PetscFunctionReturn(PETSC_SUCCESS);
146: }

148: static PetscErrorCode PCDestroy_Cholesky(PC pc)
149: {
150:   PC_Cholesky *dir = (PC_Cholesky *)pc->data;

152:   PetscFunctionBegin;
153:   PetscCall(PCReset_Cholesky(pc));
154:   PetscCall(PetscFree(((PC_Factor *)dir)->ordering));
155:   PetscCall(PetscFree(((PC_Factor *)dir)->solvertype));
156:   PetscCall(PCFactorClearComposedFunctions(pc));
157:   PetscCall(PetscFree(pc->data));
158:   PetscFunctionReturn(PETSC_SUCCESS);
159: }

161: static PetscErrorCode PCApply_Cholesky(PC pc, Vec x, Vec y)
162: {
163:   PC_Cholesky *dir = (PC_Cholesky *)pc->data;

165:   PetscFunctionBegin;
166:   if (dir->hdr.inplace) {
167:     PetscCall(MatSolve(pc->pmat, x, y));
168:   } else {
169:     PetscCall(MatSolve(((PC_Factor *)dir)->fact, x, y));
170:   }
171:   PetscFunctionReturn(PETSC_SUCCESS);
172: }

174: static PetscErrorCode PCMatApply_Cholesky(PC pc, Mat X, Mat Y)
175: {
176:   PC_Cholesky *dir = (PC_Cholesky *)pc->data;

178:   PetscFunctionBegin;
179:   if (dir->hdr.inplace) {
180:     PetscCall(MatMatSolve(pc->pmat, X, Y));
181:   } else {
182:     PetscCall(MatMatSolve(((PC_Factor *)dir)->fact, X, Y));
183:   }
184:   PetscFunctionReturn(PETSC_SUCCESS);
185: }

187: static PetscErrorCode PCApplySymmetricLeft_Cholesky(PC pc, Vec x, Vec y)
188: {
189:   PC_Cholesky *dir = (PC_Cholesky *)pc->data;

191:   PetscFunctionBegin;
192:   if (dir->hdr.inplace) {
193:     PetscCall(MatForwardSolve(pc->pmat, x, y));
194:   } else {
195:     PetscCall(MatForwardSolve(((PC_Factor *)dir)->fact, x, y));
196:   }
197:   PetscFunctionReturn(PETSC_SUCCESS);
198: }

200: static PetscErrorCode PCApplySymmetricRight_Cholesky(PC pc, Vec x, Vec y)
201: {
202:   PC_Cholesky *dir = (PC_Cholesky *)pc->data;

204:   PetscFunctionBegin;
205:   if (dir->hdr.inplace) {
206:     PetscCall(MatBackwardSolve(pc->pmat, x, y));
207:   } else {
208:     PetscCall(MatBackwardSolve(((PC_Factor *)dir)->fact, x, y));
209:   }
210:   PetscFunctionReturn(PETSC_SUCCESS);
211: }

213: static PetscErrorCode PCApplyTranspose_Cholesky(PC pc, Vec x, Vec y)
214: {
215:   PC_Cholesky *dir = (PC_Cholesky *)pc->data;

217:   PetscFunctionBegin;
218:   if (dir->hdr.inplace) {
219:     PetscCall(MatSolveTranspose(pc->pmat, x, y));
220:   } else {
221:     PetscCall(MatSolveTranspose(((PC_Factor *)dir)->fact, x, y));
222:   }
223:   PetscFunctionReturn(PETSC_SUCCESS);
224: }

226: /*@
227:    PCFactorSetReuseOrdering - When similar matrices are factored, this
228:    causes the ordering computed in the first factor to be used for all
229:    following factors.

231:    Logically Collective

233:    Input Parameters:
234: +  pc - the preconditioner context
235: -  flag - `PETSC_TRUE` to reuse else `PETSC_FALSE`

237:    Options Database Key:
238: .  -pc_factor_reuse_ordering - Activate `PCFactorSetReuseOrdering()`

240:    Level: intermediate

242: .seealso: `PCLU`, `PCCHOLESKY`, `PCFactorSetReuseFill()`
243: @*/
244: PetscErrorCode PCFactorSetReuseOrdering(PC pc, PetscBool flag)
245: {
246:   PetscFunctionBegin;
249:   PetscTryMethod(pc, "PCFactorSetReuseOrdering_C", (PC, PetscBool), (pc, flag));
250:   PetscFunctionReturn(PETSC_SUCCESS);
251: }

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

256:    Options Database Keys:
257: +  -pc_factor_reuse_ordering - Activate `PCFactorSetReuseOrdering()`
258: .  -pc_factor_mat_solver_type - Actives `PCFactorSetMatSolverType()` to choose the direct solver, like superlu
259: .  -pc_factor_reuse_fill - Activates `PCFactorSetReuseFill()`
260: .  -pc_factor_fill <fill> - Sets fill amount
261: .  -pc_factor_in_place - Activates in-place factorization
262: -  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine

264:    Level: beginner

266:    Notes:
267:    Not all options work for all matrix formats

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

273: .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`,
274:           `PCILU`, `PCLU`, `PCICC`, `PCFactorSetReuseOrdering()`, `PCFactorSetReuseFill()`, `PCFactorGetMatrix()`,
275:           `PCFactorSetFill()`, `PCFactorSetShiftNonzero()`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()`
276:           `PCFactorSetUseInPlace()`, `PCFactorGetUseInPlace()`, `PCFactorSetMatOrderingType()`, `PCFactorSetReuseOrdering()`
277: M*/

279: PETSC_EXTERN PetscErrorCode PCCreate_Cholesky(PC pc)
280: {
281:   PC_Cholesky *dir;

283:   PetscFunctionBegin;
284:   PetscCall(PetscNew(&dir));
285:   pc->data = (void *)dir;
286:   PetscCall(PCFactorInitialize(pc, MAT_FACTOR_CHOLESKY));

288:   ((PC_Factor *)dir)->info.fill = 5.0;

290:   pc->ops->destroy             = PCDestroy_Cholesky;
291:   pc->ops->reset               = PCReset_Cholesky;
292:   pc->ops->apply               = PCApply_Cholesky;
293:   pc->ops->matapply            = PCMatApply_Cholesky;
294:   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_Cholesky;
295:   pc->ops->applysymmetricright = PCApplySymmetricRight_Cholesky;
296:   pc->ops->applytranspose      = PCApplyTranspose_Cholesky;
297:   pc->ops->setup               = PCSetUp_Cholesky;
298:   pc->ops->setfromoptions      = PCSetFromOptions_Cholesky;
299:   pc->ops->view                = PCView_Factor;
300:   pc->ops->applyrichardson     = NULL;
301:   PetscFunctionReturn(PETSC_SUCCESS);
302: }