Actual source code: cholesky.c

  1: #define PETSCKSP_DLL

  3: /*
  4:    Defines a direct factorization preconditioner for any Mat implementation
  5:    Note: this need not be consided a preconditioner since it supplies
  6:          a direct solver.
  7: */
 8:  #include ../src/ksp/pc/impls/factor/factor.h

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

 22: PetscErrorCode  PCFactorSetReuseOrdering_Cholesky(PC pc,PetscTruth flag)
 23: {
 24:   PC_Cholesky *lu;
 25: 
 27:   lu               = (PC_Cholesky*)pc->data;
 28:   lu->reuseordering = flag;
 29:   return(0);
 30: }

 36: PetscErrorCode  PCFactorSetReuseFill_Cholesky(PC pc,PetscTruth flag)
 37: {
 38:   PC_Cholesky *lu;
 39: 
 41:   lu = (PC_Cholesky*)pc->data;
 42:   lu->reusefill = flag;
 43:   return(0);
 44: }

 49: static PetscErrorCode PCSetFromOptions_Cholesky(PC pc)
 50: {
 52: 
 54:   PetscOptionsHead("Cholesky options");
 55:     PCSetFromOptions_Factor(pc);
 56:   PetscOptionsTail();
 57:   return(0);
 58: }

 62: static PetscErrorCode PCView_Cholesky(PC pc,PetscViewer viewer)
 63: {
 64:   PC_Cholesky    *chol = (PC_Cholesky*)pc->data;
 66:   PetscTruth     iascii;
 67: 
 69:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
 70:   if (iascii) {
 71:     if (chol->inplace) {
 72:       PetscViewerASCIIPrintf(viewer,"  Cholesky: in-place factorization\n");
 73:     } else {
 74:       PetscViewerASCIIPrintf(viewer,"  Cholesky: out-of-place factorization\n");
 75:     }
 76: 
 77:     if (chol->reusefill)    {PetscViewerASCIIPrintf(viewer,"  Reusing fill from past factorization\n");}
 78:     if (chol->reuseordering) {PetscViewerASCIIPrintf(viewer,"  Reusing reordering from past factorization\n");}
 79:   }
 80:   PCView_Factor(pc,viewer);
 81:   return(0);
 82: }


 87: static PetscErrorCode PCSetUp_Cholesky(PC pc)
 88: {
 90:   PetscTruth     flg;
 91:   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;

 94:   if (dir->reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->actualfill;
 95: 
 96:   if (dir->inplace) {
 97:     if (dir->row && dir->col && (dir->row != dir->col)) {
 98:       ISDestroy(dir->row);
 99:       dir->row = 0;
100:     }
101:     if (dir->col) {
102:       ISDestroy(dir->col);
103:       dir->col = 0;
104:     }
105:     MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);
106:     if (dir->col && (dir->row != dir->col)) {  /* only use row ordering for SBAIJ */
107:       ISDestroy(dir->col);
108:       dir->col=0;
109:     }
110:     if (dir->row) {PetscLogObjectParent(pc,dir->row);}
111:     MatCholeskyFactor(pc->pmat,dir->row,&((PC_Factor*)dir)->info);
112:     ((PC_Factor*)dir)->fact = pc->pmat;
113:   } else {
114:     MatInfo info;
115:     if (!pc->setupcalled) {
116:       MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);
117:       /* check if dir->row == dir->col */
118:       ISEqual(dir->row,dir->col,&flg);
119:       if (!flg) SETERRQ(PETSC_ERR_ARG_INCOMP,"row and column permutations must equal");
120:       ISDestroy(dir->col); /* only pass one ordering into CholeskyFactor */
121:       dir->col=0;

123:       flg  = PETSC_FALSE;
124:       PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg,PETSC_NULL);
125:       if (flg) {
126:         PetscReal tol = 1.e-10;
127:         PetscOptionsGetReal(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,PETSC_NULL);
128:         MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);
129:       }
130:       if (dir->row) {PetscLogObjectParent(pc,dir->row);}
131:       MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact);
132:       MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info);
133:       MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);
134:       dir->actualfill = info.fill_ratio_needed;
135:       PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);
136:     } else if (pc->flag != SAME_NONZERO_PATTERN) {
137:       if (!dir->reuseordering) {
138:         if (dir->row && dir->col && (dir->row != dir->col)) {
139:           ISDestroy(dir->row);
140:           dir->row = 0;
141:         }
142:         if (dir->col) {
143:           ISDestroy(dir->col);
144:           dir->col =0;
145:         }
146:         MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);
147:         if (dir->col && (dir->row != dir->col)) {  /* only use row ordering for SBAIJ */
148:           ISDestroy(dir->col);
149:           dir->col=0;
150:         }
151:         flg  = PETSC_FALSE;
152:         PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg,PETSC_NULL);
153:         if (flg) {
154:           PetscReal tol = 1.e-10;
155:           PetscOptionsGetReal(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,PETSC_NULL);
156:           MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);
157:         }
158:         if (dir->row) {PetscLogObjectParent(pc,dir->row);}
159:       }
160:       MatDestroy(((PC_Factor*)dir)->fact);
161:       MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact);
162:       MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info);
163:       MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);
164:       dir->actualfill = info.fill_ratio_needed;
165:       PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);
166:     }
167:     MatCholeskyFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);
168:   }
169:   return(0);
170: }

174: static PetscErrorCode PCDestroy_Cholesky(PC pc)
175: {
176:   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;

180:   if (!dir->inplace && ((PC_Factor*)dir)->fact) {MatDestroy(((PC_Factor*)dir)->fact);}
181:   if (dir->row) {ISDestroy(dir->row);}
182:   if (dir->col) {ISDestroy(dir->col);}
183:   PetscStrfree(((PC_Factor*)dir)->ordering);
184:   PetscStrfree(((PC_Factor*)dir)->solvertype);
185:   PetscFree(dir);
186:   return(0);
187: }

191: static PetscErrorCode PCApply_Cholesky(PC pc,Vec x,Vec y)
192: {
193:   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
195: 
197:   if (dir->inplace) {MatSolve(pc->pmat,x,y);}
198:   else              {MatSolve(((PC_Factor*)dir)->fact,x,y);}
199:   return(0);
200: }

204: static PetscErrorCode PCApplyTranspose_Cholesky(PC pc,Vec x,Vec y)
205: {
206:   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;

210:   if (dir->inplace) {MatSolveTranspose(pc->pmat,x,y);}
211:   else              {MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);}
212:   return(0);
213: }

215: /* -----------------------------------------------------------------------------------*/

220: PetscErrorCode  PCFactorSetUseInPlace_Cholesky(PC pc)
221: {
222:   PC_Cholesky *dir;

225:   dir = (PC_Cholesky*)pc->data;
226:   dir->inplace = PETSC_TRUE;
227:   return(0);
228: }

231: /* -----------------------------------------------------------------------------------*/

235: /*@
236:    PCFactorSetReuseOrdering - When similar matrices are factored, this
237:    causes the ordering computed in the first factor to be used for all
238:    following factors.

240:    Collective on PC

242:    Input Parameters:
243: +  pc - the preconditioner context
244: -  flag - PETSC_TRUE to reuse else PETSC_FALSE

246:    Options Database Key:
247: .  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()

249:    Level: intermediate

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

253: .seealso: PCFactorSetReuseFill()
254: @*/
255: PetscErrorCode  PCFactorSetReuseOrdering(PC pc,PetscTruth flag)
256: {
257:   PetscErrorCode ierr,(*f)(PC,PetscTruth);

261:   PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",(void (**)(void))&f);
262:   if (f) {
263:     (*f)(pc,flag);
264:   }
265:   return(0);
266: }


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

272:    Options Database Keys:
273: +  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
274: .  -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like spooles
275: .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
276: .  -pc_factor_fill <fill> - Sets fill amount
277: .  -pc_factor_in_place - Activates in-place factorization
278: -  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine

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

282:    Level: beginner

284:    Concepts: Cholesky factorization, direct solver

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

290: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
291:            PCILU, PCLU, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(),
292:            PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetShiftType(), PCFactorSetShiftAmount()
293:            PCFactorSetUseInPlace(), PCFactorSetMatOrderingType()

295: M*/

300: PetscErrorCode  PCCreate_Cholesky(PC pc)
301: {
303:   PC_Cholesky    *dir;

306:   PetscNewLog(pc,PC_Cholesky,&dir);

308:   ((PC_Factor*)dir)->fact                   = 0;
309:   dir->inplace                = PETSC_FALSE;
310:   MatFactorInfoInitialize(&((PC_Factor*)dir)->info);
311:   ((PC_Factor*)dir)->factortype         = MAT_FACTOR_CHOLESKY;
312:   ((PC_Factor*)dir)->info.fill          = 5.0;
313:   ((PC_Factor*)dir)->info.shifttype     = (PetscReal) MAT_SHIFT_NONE;
314:   ((PC_Factor*)dir)->info.shiftamount   = 0.0;
315:   ((PC_Factor*)dir)->info.pivotinblocks = 1.0;
316:   dir->col                    = 0;
317:   dir->row                    = 0;
318:   PetscStrallocpy(MATORDERING_NATURAL,&((PC_Factor*)dir)->ordering);
319:   PetscStrallocpy(MAT_SOLVER_PETSC,&((PC_Factor*)dir)->solvertype);
320:   dir->reusefill        = PETSC_FALSE;
321:   dir->reuseordering    = PETSC_FALSE;
322:   pc->data              = (void*)dir;

324:   pc->ops->destroy           = PCDestroy_Cholesky;
325:   pc->ops->apply             = PCApply_Cholesky;
326:   pc->ops->applytranspose    = PCApplyTranspose_Cholesky;
327:   pc->ops->setup             = PCSetUp_Cholesky;
328:   pc->ops->setfromoptions    = PCSetFromOptions_Cholesky;
329:   pc->ops->view              = PCView_Cholesky;
330:   pc->ops->applyrichardson   = 0;
331:   pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor;

333:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Factor",
334:                     PCFactorSetMatSolverPackage_Factor);
335:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor",
336:                     PCFactorGetMatSolverPackage_Factor);
337:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor",
338:                     PCFactorSetZeroPivot_Factor);
339:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftType_C","PCFactorSetShiftType_Factor",
340:                     PCFactorSetShiftType_Factor);
341:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftAmount_C","PCFactorSetShiftAmount_Factor",
342:                     PCFactorSetShiftAmount_Factor);
343:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor",
344:                     PCFactorSetFill_Factor);
345:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_Cholesky",
346:                     PCFactorSetUseInPlace_Cholesky);
347:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor",
348:                     PCFactorSetMatOrderingType_Factor);
349:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_Cholesky",
350:                     PCFactorSetReuseOrdering_Cholesky);
351:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_Cholesky",
352:                     PCFactorSetReuseFill_Cholesky);
353:   return(0);
354: }