Actual source code: lu.c

petsc-main 2021-04-20
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>


 11: PetscErrorCode PCFactorReorderForNonzeroDiagonal_LU(PC pc,PetscReal z)
 12: {
 13:   PC_LU *lu = (PC_LU*)pc->data;

 16:   lu->nonzerosalongdiagonal = PETSC_TRUE;
 17:   if (z == PETSC_DECIDE) lu->nonzerosalongdiagonaltol = 1.e-10;
 18:   else lu->nonzerosalongdiagonaltol = z;
 19:   return(0);
 20: }

 22: static PetscErrorCode PCSetFromOptions_LU(PetscOptionItems *PetscOptionsObject,PC pc)
 23: {
 24:   PC_LU          *lu = (PC_LU*)pc->data;
 26:   PetscBool      flg = PETSC_FALSE;
 27:   PetscReal      tol;

 30:   PetscOptionsHead(PetscOptionsObject,"LU options");
 31:   PCSetFromOptions_Factor(PetscOptionsObject,pc);

 33:   PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);
 34:   if (flg) {
 35:     tol  = PETSC_DECIDE;
 36:     PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,NULL);
 37:     PCFactorReorderForNonzeroDiagonal(pc,tol);
 38:   }
 39:   PetscOptionsTail();
 40:   return(0);
 41: }

 43: static PetscErrorCode PCSetUp_LU(PC pc)
 44: {
 45:   PetscErrorCode         ierr;
 46:   PC_LU                  *dir = (PC_LU*)pc->data;
 47:   MatSolverType          stype;
 48:   MatFactorError         err;

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

 54:   MatSetErrorIfFailure(pc->pmat,pc->erroriffailure);
 55:   if (dir->hdr.inplace) {
 56:     MatFactorType ftype;

 58:     MatGetFactorType(pc->pmat, &ftype);
 59:     if (ftype == MAT_FACTOR_NONE) {
 60:       if (dir->row && dir->col && dir->row != dir->col) {ISDestroy(&dir->row);}
 61:       ISDestroy(&dir->col);
 62:       /* This should only get the ordering if needed, but since MatGetFactor() is not called we can't know if it is needed */
 63:       PCFactorSetDefaultOrdering_Factor(pc);
 64:       MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);
 65:       if (dir->row) {
 66:         PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);
 67:         PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);
 68:       }
 69:       MatLUFactor(pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);
 70:       MatFactorGetError(pc->pmat,&err);
 71:       if (err) { /* Factor() fails */
 72:         pc->failedreason = (PCFailedReason)err;
 73:         return(0);
 74:       }
 75:     }
 76:     ((PC_Factor*)dir)->fact = pc->pmat;
 77:   } else {
 78:     MatInfo info;

 80:     if (!pc->setupcalled) {
 81:       PetscBool canuseordering;
 82:       if (!((PC_Factor*)dir)->fact) {
 83:         MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);
 84:         PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);
 85:       }
 86:       MatFactorGetCanUseOrdering(((PC_Factor*)dir)->fact,&canuseordering);
 87:       if (canuseordering) {
 88:         PCFactorSetDefaultOrdering_Factor(pc);
 89:         MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);
 90:         if (dir->nonzerosalongdiagonal) {
 91:           MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);
 92:         }
 93:         PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);
 94:         PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);
 95:       }
 96:       MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);
 97:       MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);
 98:       dir->hdr.actualfill = info.fill_ratio_needed;
 99:     } else if (pc->flag != SAME_NONZERO_PATTERN) {
100:       PetscBool canuseordering;
101:       if (!dir->hdr.reuseordering) {
102:         MatDestroy(&((PC_Factor*)dir)->fact);
103:         MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);
104:         PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);
105:         MatFactorGetCanUseOrdering(((PC_Factor*)dir)->fact,&canuseordering);
106:         if (canuseordering) {
107:           if (dir->row && dir->col && dir->row != dir->col) {ISDestroy(&dir->row);}
108:           ISDestroy(&dir->col);
109:           PCFactorSetDefaultOrdering_Factor(pc);
110:           MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);
111:           if (dir->nonzerosalongdiagonal) {
112:             MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);
113:           }
114:           PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);
115:           PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);
116:         }
117:       }
118:       MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);
119:       MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);
120:       dir->hdr.actualfill = info.fill_ratio_needed;
121:     } else {
122:       MatFactorGetError(((PC_Factor*)dir)->fact,&err);
123:       if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) {
124:         MatFactorClearError(((PC_Factor*)dir)->fact);
125:         pc->failedreason = PC_NOERROR;
126:       }
127:     }
128:     MatFactorGetError(((PC_Factor*)dir)->fact,&err);
129:     if (err) { /* FactorSymbolic() fails */
130:       pc->failedreason = (PCFailedReason)err;
131:       return(0);
132:     }

134:     MatLUFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);
135:     MatFactorGetError(((PC_Factor*)dir)->fact,&err);
136:     if (err) { /* FactorNumeric() fails */
137:       pc->failedreason = (PCFailedReason)err;
138:     }

140:   }

142:   PCFactorGetMatSolverType(pc,&stype);
143:   if (!stype) {
144:     MatSolverType solverpackage;
145:     MatFactorGetSolverType(((PC_Factor*)dir)->fact,&solverpackage);
146:     PCFactorSetMatSolverType(pc,solverpackage);
147:   }
148:   return(0);
149: }

151: static PetscErrorCode PCReset_LU(PC pc)
152: {
153:   PC_LU          *dir = (PC_LU*)pc->data;

157:   if (!dir->hdr.inplace && ((PC_Factor*)dir)->fact) {MatDestroy(&((PC_Factor*)dir)->fact);}
158:   if (dir->row && dir->col && dir->row != dir->col) {ISDestroy(&dir->row);}
159:   ISDestroy(&dir->col);
160:   return(0);
161: }

163: static PetscErrorCode PCDestroy_LU(PC pc)
164: {
165:   PC_LU          *dir = (PC_LU*)pc->data;

169:   PCReset_LU(pc);
170:   PetscFree(((PC_Factor*)dir)->ordering);
171:   PetscFree(((PC_Factor*)dir)->solvertype);
172:   PetscFree(pc->data);
173:   return(0);
174: }

176: static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y)
177: {
178:   PC_LU          *dir = (PC_LU*)pc->data;

182:   if (dir->hdr.inplace) {
183:     MatSolve(pc->pmat,x,y);
184:   } else {
185:     MatSolve(((PC_Factor*)dir)->fact,x,y);
186:   }
187:   return(0);
188: }

190: static PetscErrorCode PCMatApply_LU(PC pc,Mat X,Mat Y)
191: {
192:   PC_LU          *dir = (PC_LU*)pc->data;

196:   if (dir->hdr.inplace) {
197:     MatMatSolve(pc->pmat,X,Y);
198:   } else {
199:     MatMatSolve(((PC_Factor*)dir)->fact,X,Y);
200:   }
201:   return(0);
202: }

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

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

218: /* -----------------------------------------------------------------------------------*/

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

223:    Options Database Keys:
224: +  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
225: .  -pc_factor_mat_solver_type - Actives PCFactorSetMatSolverType() to choose the direct solver, like superlu
226: .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
227: .  -pc_factor_fill <fill> - Sets fill amount
228: .  -pc_factor_in_place - Activates in-place factorization
229: .  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
230: .  -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase
231:                                          stability of factorization.
232: .  -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types
233: .  -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default
234: -   -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the
235:         diagonal.

237:    Notes:
238:     Not all options work for all matrix formats
239:           Run with -help to see additional options for particular matrix formats or factorization
240:           algorithms

242:    Level: beginner

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

249: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
250:            PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(),
251:            PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetColumnPivot(),
252:            PCFactorSetPivotInBlocks(),PCFactorSetShiftType(),PCFactorSetShiftAmount()
253:            PCFactorReorderForNonzeroDiagonal()
254: M*/

256: PETSC_EXTERN PetscErrorCode PCCreate_LU(PC pc)
257: {
259:   PC_LU          *dir;

262:   PetscNewLog(pc,&dir);
263:   pc->data = (void*)dir;
264:   PCFactorInitialize(pc,MAT_FACTOR_LU);
265:   dir->nonzerosalongdiagonal    = PETSC_FALSE;

267:   ((PC_Factor*)dir)->info.fill          = 5.0;
268:   ((PC_Factor*)dir)->info.dtcol         = 1.e-6;  /* default to pivoting; this is only thing PETSc LU supports */
269:   ((PC_Factor*)dir)->info.shifttype     = (PetscReal)MAT_SHIFT_NONE;
270:   dir->col                              = NULL;
271:   dir->row                              = NULL;

273:   pc->ops->reset             = PCReset_LU;
274:   pc->ops->destroy           = PCDestroy_LU;
275:   pc->ops->apply             = PCApply_LU;
276:   pc->ops->matapply          = PCMatApply_LU;
277:   pc->ops->applytranspose    = PCApplyTranspose_LU;
278:   pc->ops->setup             = PCSetUp_LU;
279:   pc->ops->setfromoptions    = PCSetFromOptions_LU;
280:   pc->ops->view              = PCView_Factor;
281:   pc->ops->applyrichardson   = NULL;
282:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",PCFactorReorderForNonzeroDiagonal_LU);
283:   return(0);
284: }