Actual source code: lu.c

petsc-master 2016-12-06
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>


 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: static PetscErrorCode PCSetFromOptions_LU(PetscOptionItems *PetscOptionsObject,PC pc)
 27: {
 28:   PC_LU          *lu = (PC_LU*)pc->data;
 30:   PetscBool      flg = PETSC_FALSE;
 31:   PetscReal      tol;

 34:   PetscOptionsHead(PetscOptionsObject,"LU options");
 35:   PCSetFromOptions_Factor(PetscOptionsObject,pc);

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

 49: static PetscErrorCode PCView_LU(PC pc,PetscViewer viewer)
 50: {

 54:   PCView_Factor(pc,viewer);
 55:   return(0);
 56: }

 60: static PetscErrorCode PCSetUp_LU(PC pc)
 61: {
 62:   PetscErrorCode         ierr;
 63:   PC_LU                  *dir = (PC_LU*)pc->data;
 64:   const MatSolverPackage stype;
 65:   MatFactorError         err;

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

 71:   MatSetErrorIfFailure(pc->pmat,pc->erroriffailure);
 72:   if (dir->hdr.inplace) {
 73:     if (dir->row && dir->col && dir->row != dir->col) {ISDestroy(&dir->row);}
 74:     ISDestroy(&dir->col);
 75:     MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);
 76:     if (dir->row) {
 77:       PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);
 78:       PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);
 79:     }
 80:     MatLUFactor(pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);
 81:     MatFactorGetError(pc->pmat,&err);
 82:     if (err) { /* Factor() fails */
 83:       pc->failedreason = (PCFailedReason)err;
 84:       return(0);
 85:     }

 87:     ((PC_Factor*)dir)->fact = pc->pmat;
 88:   } else {
 89:     MatInfo info;

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

139:     MatLUFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);
140:     MatFactorGetError(((PC_Factor*)dir)->fact,&err);
141:     if (err) { /* FactorNumeric() fails */
142:       pc->failedreason = (PCFailedReason)err;
143:     }

145:   }

147:   PCFactorGetMatSolverPackage(pc,&stype);
148:   if (!stype) {
149:     const MatSolverPackage solverpackage;
150:     MatFactorGetSolverPackage(((PC_Factor*)dir)->fact,&solverpackage);
151:     PCFactorSetMatSolverPackage(pc,solverpackage);
152:   }
153:   return(0);
154: }

158: static PetscErrorCode PCReset_LU(PC pc)
159: {
160:   PC_LU          *dir = (PC_LU*)pc->data;

164:   if (!dir->hdr.inplace && ((PC_Factor*)dir)->fact) {MatDestroy(&((PC_Factor*)dir)->fact);}
165:   if (dir->row && dir->col && dir->row != dir->col) {ISDestroy(&dir->row);}
166:   ISDestroy(&dir->col);
167:   return(0);
168: }

172: static PetscErrorCode PCDestroy_LU(PC pc)
173: {
174:   PC_LU          *dir = (PC_LU*)pc->data;

178:   PCReset_LU(pc);
179:   PetscFree(((PC_Factor*)dir)->ordering);
180:   PetscFree(((PC_Factor*)dir)->solvertype);
181:   PetscFree(pc->data);
182:   return(0);
183: }

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

193:   if (dir->hdr.inplace) {
194:     MatSolve(pc->pmat,x,y);
195:   } else {
196:     MatSolve(((PC_Factor*)dir)->fact,x,y);
197:   }
198:   return(0);
199: }

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

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

217: /* -----------------------------------------------------------------------------------*/

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

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

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

240:    Level: beginner

242:    Concepts: LU factorization, direct solver

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

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

257: PETSC_EXTERN PetscErrorCode PCCreate_LU(PC pc)
258: {
260:   PetscMPIInt    size;
261:   PC_LU          *dir;

264:   PetscNewLog(pc,&dir);
265:   pc->data = (void*)dir;
266:   PCFactorInitialize(pc);
267:   ((PC_Factor*)dir)->factortype = MAT_FACTOR_LU;
268:   dir->nonzerosalongdiagonal    = PETSC_FALSE;

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

276:   MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);
277:   if (size == 1) {
278:     PetscStrallocpy(MATORDERINGND,(char**)&((PC_Factor*)dir)->ordering);
279:   } else {
280:     PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)dir)->ordering);
281:   }

283:   pc->ops->reset             = PCReset_LU;
284:   pc->ops->destroy           = PCDestroy_LU;
285:   pc->ops->apply             = PCApply_LU;
286:   pc->ops->applytranspose    = PCApplyTranspose_LU;
287:   pc->ops->setup             = PCSetUp_LU;
288:   pc->ops->setfromoptions    = PCSetFromOptions_LU;
289:   pc->ops->view              = PCView_LU;
290:   pc->ops->applyrichardson   = 0;
291:   PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",PCFactorReorderForNonzeroDiagonal_LU);
292:   return(0);
293: }