Actual source code: mg.c

petsc-3.11.4 2019-09-28
Report Typos and Errors

  2: /*
  3:     Defines the multigrid preconditioner interface.
  4: */
  5:  #include <petsc/private/pcmgimpl.h>
  6:  #include <petscdm.h>
  7: PETSC_INTERN PetscErrorCode PCPreSolveChangeRHS(PC,PetscBool*);

  9: PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PCRichardsonConvergedReason *reason)
 10: {
 11:   PC_MG          *mg = (PC_MG*)pc->data;
 12:   PC_MG_Levels   *mgc,*mglevels = *mglevelsin;
 14:   PetscInt       cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles;

 17:   if (mglevels->eventsmoothsolve) {PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);}
 18:   KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);  /* pre-smooth */
 19:   KSPCheckSolve(mglevels->smoothd,pc,mglevels->x);
 20:   if (mglevels->eventsmoothsolve) {PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);}
 21:   if (mglevels->level) {  /* not the coarsest grid */
 22:     if (mglevels->eventresidual) {PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);}
 23:     (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);
 24:     if (mglevels->eventresidual) {PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);}

 26:     /* if on finest level and have convergence criteria set */
 27:     if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) {
 28:       PetscReal rnorm;
 29:       VecNorm(mglevels->r,NORM_2,&rnorm);
 30:       if (rnorm <= mg->ttol) {
 31:         if (rnorm < mg->abstol) {
 32:           *reason = PCRICHARDSON_CONVERGED_ATOL;
 33:           PetscInfo2(pc,"Linear solver has converged. Residual norm %g is less than absolute tolerance %g\n",(double)rnorm,(double)mg->abstol);
 34:         } else {
 35:           *reason = PCRICHARDSON_CONVERGED_RTOL;
 36:           PetscInfo2(pc,"Linear solver has converged. Residual norm %g is less than relative tolerance times initial residual norm %g\n",(double)rnorm,(double)mg->ttol);
 37:         }
 38:         return(0);
 39:       }
 40:     }

 42:     mgc = *(mglevelsin - 1);
 43:     if (mglevels->eventinterprestrict) {PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);}
 44:     MatRestrict(mglevels->restrct,mglevels->r,mgc->b);
 45:     if (mglevels->eventinterprestrict) {PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);}
 46:     VecSet(mgc->x,0.0);
 47:     while (cycles--) {
 48:       PCMGMCycle_Private(pc,mglevelsin-1,reason);
 49:     }
 50:     if (mglevels->eventinterprestrict) {PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);}
 51:     MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);
 52:     if (mglevels->eventinterprestrict) {PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);}
 53:     if (mglevels->eventsmoothsolve) {PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);}
 54:     KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);    /* post smooth */
 55:     KSPCheckSolve(mglevels->smoothu,pc,mglevels->x);
 56:     if (mglevels->eventsmoothsolve) {PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);}
 57:   }
 58:   return(0);
 59: }

 61: static PetscErrorCode PCApplyRichardson_MG(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscBool zeroguess,PetscInt *outits,PCRichardsonConvergedReason *reason)
 62: {
 63:   PC_MG          *mg        = (PC_MG*)pc->data;
 64:   PC_MG_Levels   **mglevels = mg->levels;
 66:   PC             tpc;
 67:   PetscBool      changeu,changed;
 68:   PetscInt       levels = mglevels[0]->levels,i;

 71:   /* When the DM is supplying the matrix then it will not exist until here */
 72:   for (i=0; i<levels; i++) {
 73:     if (!mglevels[i]->A) {
 74:       KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);
 75:       PetscObjectReference((PetscObject)mglevels[i]->A);
 76:     }
 77:   }

 79:   KSPGetPC(mglevels[levels-1]->smoothd,&tpc);
 80:   PCPreSolveChangeRHS(tpc,&changed);
 81:   KSPGetPC(mglevels[levels-1]->smoothu,&tpc);
 82:   PCPreSolveChangeRHS(tpc,&changeu);
 83:   if (!changed && !changeu) {
 84:     VecDestroy(&mglevels[levels-1]->b);
 85:     mglevels[levels-1]->b = b;
 86:   } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */
 87:     if (!mglevels[levels-1]->b) {
 88:       Vec *vec;

 90:       KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);
 91:       mglevels[levels-1]->b = *vec;
 92:       PetscFree(vec);
 93:     }
 94:     VecCopy(b,mglevels[levels-1]->b);
 95:   }
 96:   mglevels[levels-1]->x = x;

 98:   mg->rtol   = rtol;
 99:   mg->abstol = abstol;
100:   mg->dtol   = dtol;
101:   if (rtol) {
102:     /* compute initial residual norm for relative convergence test */
103:     PetscReal rnorm;
104:     if (zeroguess) {
105:       VecNorm(b,NORM_2,&rnorm);
106:     } else {
107:       (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);
108:       VecNorm(w,NORM_2,&rnorm);
109:     }
110:     mg->ttol = PetscMax(rtol*rnorm,abstol);
111:   } else if (abstol) mg->ttol = abstol;
112:   else mg->ttol = 0.0;

114:   /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't
115:      stop prematurely due to small residual */
116:   for (i=1; i<levels; i++) {
117:     KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
118:     if (mglevels[i]->smoothu != mglevels[i]->smoothd) {
119:       /* For Richardson the initial guess is nonzero since it is solving in each cycle the original system not just applying as a preconditioner */
120:       KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);
121:       KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
122:     }
123:   }

125:   *reason = (PCRichardsonConvergedReason)0;
126:   for (i=0; i<its; i++) {
127:     PCMGMCycle_Private(pc,mglevels+levels-1,reason);
128:     if (*reason) break;
129:   }
130:   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
131:   *outits = i;
132:   if (!changed && !changeu) mglevels[levels-1]->b = NULL;
133:   return(0);
134: }

136: PetscErrorCode PCReset_MG(PC pc)
137: {
138:   PC_MG          *mg        = (PC_MG*)pc->data;
139:   PC_MG_Levels   **mglevels = mg->levels;
141:   PetscInt       i,n;

144:   if (mglevels) {
145:     n = mglevels[0]->levels;
146:     for (i=0; i<n-1; i++) {
147:       VecDestroy(&mglevels[i+1]->r);
148:       VecDestroy(&mglevels[i]->b);
149:       VecDestroy(&mglevels[i]->x);
150:       MatDestroy(&mglevels[i+1]->restrct);
151:       MatDestroy(&mglevels[i+1]->interpolate);
152:       MatDestroy(&mglevels[i+1]->inject);
153:       VecDestroy(&mglevels[i+1]->rscale);
154:     }
155:     /* this is not null only if the smoother on the finest level
156:        changes the rhs during PreSolve */
157:     VecDestroy(&mglevels[n-1]->b);

159:     for (i=0; i<n; i++) {
160:       MatDestroy(&mglevels[i]->A);
161:       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
162:         KSPReset(mglevels[i]->smoothd);
163:       }
164:       KSPReset(mglevels[i]->smoothu);
165:     }
166:   }
167:   return(0);
168: }

170: PetscErrorCode PCMGSetLevels_MG(PC pc,PetscInt levels,MPI_Comm *comms)
171: {
173:   PC_MG          *mg        = (PC_MG*)pc->data;
174:   MPI_Comm       comm;
175:   PC_MG_Levels   **mglevels = mg->levels;
176:   PCMGType       mgtype     = mg->am;
177:   PetscInt       mgctype    = (PetscInt) PC_MG_CYCLE_V;
178:   PetscInt       i;
179:   PetscMPIInt    size;
180:   const char     *prefix;
181:   PC             ipc;
182:   PetscInt       n;

187:   if (mg->nlevels == levels) return(0);
188:   PetscObjectGetComm((PetscObject)pc,&comm);
189:   if (mglevels) {
190:     mgctype = mglevels[0]->cycles;
191:     /* changing the number of levels so free up the previous stuff */
192:     PCReset_MG(pc);
193:     n    = mglevels[0]->levels;
194:     for (i=0; i<n; i++) {
195:       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
196:         KSPDestroy(&mglevels[i]->smoothd);
197:       }
198:       KSPDestroy(&mglevels[i]->smoothu);
199:       PetscFree(mglevels[i]);
200:     }
201:     PetscFree(mg->levels);
202:   }

204:   mg->nlevels = levels;

206:   PetscMalloc1(levels,&mglevels);
207:   PetscLogObjectMemory((PetscObject)pc,levels*(sizeof(PC_MG*)));

209:   PCGetOptionsPrefix(pc,&prefix);

211:   mg->stageApply = 0;
212:   for (i=0; i<levels; i++) {
213:     PetscNewLog(pc,&mglevels[i]);

215:     mglevels[i]->level               = i;
216:     mglevels[i]->levels              = levels;
217:     mglevels[i]->cycles              = mgctype;
218:     mg->default_smoothu              = 2;
219:     mg->default_smoothd              = 2;
220:     mglevels[i]->eventsmoothsetup    = 0;
221:     mglevels[i]->eventsmoothsolve    = 0;
222:     mglevels[i]->eventresidual       = 0;
223:     mglevels[i]->eventinterprestrict = 0;

225:     if (comms) comm = comms[i];
226:     KSPCreate(comm,&mglevels[i]->smoothd);
227:     KSPSetErrorIfNotConverged(mglevels[i]->smoothd,pc->erroriffailure);
228:     PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);
229:     KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);
230:     PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->smoothd, PetscMGLevelId, mglevels[i]->level);
231:     if (i || levels == 1) {
232:       char tprefix[128];

234:       KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);
235:       KSPSetConvergenceTest(mglevels[i]->smoothd,KSPConvergedSkip,NULL,NULL);
236:       KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);
237:       KSPGetPC(mglevels[i]->smoothd,&ipc);
238:       PCSetType(ipc,PCSOR);
239:       KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd);

241:       sprintf(tprefix,"mg_levels_%d_",(int)i);
242:       KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);
243:     } else {
244:       KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");

246:       /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */
247:       KSPSetType(mglevels[0]->smoothd,KSPPREONLY);
248:       KSPGetPC(mglevels[0]->smoothd,&ipc);
249:       MPI_Comm_size(comm,&size);
250:       if (size > 1) {
251:         PCSetType(ipc,PCREDUNDANT);
252:       } else {
253:         PCSetType(ipc,PCLU);
254:       }
255:       PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS);
256:     }
257:     PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[i]->smoothd);

259:     mglevels[i]->smoothu = mglevels[i]->smoothd;
260:     mg->rtol             = 0.0;
261:     mg->abstol           = 0.0;
262:     mg->dtol             = 0.0;
263:     mg->ttol             = 0.0;
264:     mg->cyclesperpcapply = 1;
265:   }
266:   mg->levels = mglevels;
267:   PCMGSetType(pc,mgtype);
268:   return(0);
269: }

271: /*@C
272:    PCMGSetLevels - Sets the number of levels to use with MG.
273:    Must be called before any other MG routine.

275:    Logically Collective on PC

277:    Input Parameters:
278: +  pc - the preconditioner context
279: .  levels - the number of levels
280: -  comms - optional communicators for each level; this is to allow solving the coarser problems
281:            on smaller sets of processors.

283:    Level: intermediate

285:    Notes:
286:      If the number of levels is one then the multigrid uses the -mg_levels prefix
287:   for setting the level options rather than the -mg_coarse prefix.

289: .keywords: MG, set, levels, multigrid

291: .seealso: PCMGSetType(), PCMGGetLevels()
292: @*/
293: PetscErrorCode PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
294: {

300:   PetscTryMethod(pc,"PCMGSetLevels_C",(PC,PetscInt,MPI_Comm*),(pc,levels,comms));
301:   return(0);
302: }


305: PetscErrorCode PCDestroy_MG(PC pc)
306: {
308:   PC_MG          *mg        = (PC_MG*)pc->data;
309:   PC_MG_Levels   **mglevels = mg->levels;
310:   PetscInt       i,n;

313:   PCReset_MG(pc);
314:   if (mglevels) {
315:     n = mglevels[0]->levels;
316:     for (i=0; i<n; i++) {
317:       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
318:         KSPDestroy(&mglevels[i]->smoothd);
319:       }
320:       KSPDestroy(&mglevels[i]->smoothu);
321:       PetscFree(mglevels[i]);
322:     }
323:     PetscFree(mg->levels);
324:   }
325:   PetscFree(pc->data);
326:   return(0);
327: }



331: extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**);
332: extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**);
333: extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**);

335: /*
336:    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
337:              or full cycle of multigrid.

339:   Note:
340:   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
341: */
342: static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
343: {
344:   PC_MG          *mg        = (PC_MG*)pc->data;
345:   PC_MG_Levels   **mglevels = mg->levels;
347:   PC             tpc;
348:   PetscInt       levels = mglevels[0]->levels,i;
349:   PetscBool      changeu,changed;

352:   if (mg->stageApply) {PetscLogStagePush(mg->stageApply);}
353:   /* When the DM is supplying the matrix then it will not exist until here */
354:   for (i=0; i<levels; i++) {
355:     if (!mglevels[i]->A) {
356:       KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);
357:       PetscObjectReference((PetscObject)mglevels[i]->A);
358:     }
359:   }

361:   KSPGetPC(mglevels[levels-1]->smoothd,&tpc);
362:   PCPreSolveChangeRHS(tpc,&changed);
363:   KSPGetPC(mglevels[levels-1]->smoothu,&tpc);
364:   PCPreSolveChangeRHS(tpc,&changeu);
365:   if (!changeu && !changed) {
366:     VecDestroy(&mglevels[levels-1]->b);
367:     mglevels[levels-1]->b = b;
368:   } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */
369:     if (!mglevels[levels-1]->b) {
370:       Vec *vec;

372:       KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);
373:       mglevels[levels-1]->b = *vec;
374:       PetscFree(vec);
375:     }
376:     VecCopy(b,mglevels[levels-1]->b);
377:   }
378:   mglevels[levels-1]->x = x;

380:   if (mg->am == PC_MG_MULTIPLICATIVE) {
381:     VecSet(x,0.0);
382:     for (i=0; i<mg->cyclesperpcapply; i++) {
383:       PCMGMCycle_Private(pc,mglevels+levels-1,NULL);
384:     }
385:   } else if (mg->am == PC_MG_ADDITIVE) {
386:     PCMGACycle_Private(pc,mglevels);
387:   } else if (mg->am == PC_MG_KASKADE) {
388:     PCMGKCycle_Private(pc,mglevels);
389:   } else {
390:     PCMGFCycle_Private(pc,mglevels);
391:   }
392:   if (mg->stageApply) {PetscLogStagePop();}
393:   if (!changeu && !changed) mglevels[levels-1]->b = NULL;
394:   return(0);
395: }


398: PetscErrorCode PCSetFromOptions_MG(PetscOptionItems *PetscOptionsObject,PC pc)
399: {
400:   PetscErrorCode   ierr;
401:   PetscInt         levels,cycles;
402:   PetscBool        flg;
403:   PC_MG            *mg = (PC_MG*)pc->data;
404:   PC_MG_Levels     **mglevels;
405:   PCMGType         mgtype;
406:   PCMGCycleType    mgctype;
407:   PCMGGalerkinType gtype;

410:   levels = PetscMax(mg->nlevels,1);
411:   PetscOptionsHead(PetscOptionsObject,"Multigrid options");
412:   PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);
413:   if (!flg && !mg->levels && pc->dm) {
414:     DMGetRefineLevel(pc->dm,&levels);
415:     levels++;
416:     mg->usedmfornumberoflevels = PETSC_TRUE;
417:   }
418:   PCMGSetLevels(pc,levels,NULL);
419:   mglevels = mg->levels;

421:   mgctype = (PCMGCycleType) mglevels[0]->cycles;
422:   PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);
423:   if (flg) {
424:     PCMGSetCycleType(pc,mgctype);
425:   }
426:   gtype = mg->galerkin;
427:   PetscOptionsEnum("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",PCMGGalerkinTypes,(PetscEnum)gtype,(PetscEnum*)&gtype,&flg);
428:   if (flg) {
429:     PCMGSetGalerkin(pc,gtype);
430:   }
431:   flg = PETSC_FALSE;
432:   PetscOptionsBool("-pc_mg_distinct_smoothup","Create seperate smoothup KSP and append the prefix _up","PCMGSetDistinctSmoothUp",PETSC_FALSE,&flg,NULL);
433:   if (flg) {
434:     PCMGSetDistinctSmoothUp(pc);
435:   }
436:   mgtype = mg->am;
437:   PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);
438:   if (flg) {
439:     PCMGSetType(pc,mgtype);
440:   }
441:   if (mg->am == PC_MG_MULTIPLICATIVE) {
442:     PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGMultiplicativeSetCycles",mg->cyclesperpcapply,&cycles,&flg);
443:     if (flg) {
444:       PCMGMultiplicativeSetCycles(pc,cycles);
445:     }
446:   }
447:   flg  = PETSC_FALSE;
448:   PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);
449:   if (flg) {
450:     PetscInt i;
451:     char     eventname[128];

453:     levels = mglevels[0]->levels;
454:     for (i=0; i<levels; i++) {
455:       sprintf(eventname,"MGSetup Level %d",(int)i);
456:       PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);
457:       sprintf(eventname,"MGSmooth Level %d",(int)i);
458:       PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);
459:       if (i) {
460:         sprintf(eventname,"MGResid Level %d",(int)i);
461:         PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);
462:         sprintf(eventname,"MGInterp Level %d",(int)i);
463:         PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);
464:       }
465:     }

467: #if defined(PETSC_USE_LOG)
468:     {
469:       const char    *sname = "MG Apply";
470:       PetscStageLog stageLog;
471:       PetscInt      st;

473:       PetscLogGetStageLog(&stageLog);
474:       for (st = 0; st < stageLog->numStages; ++st) {
475:         PetscBool same;

477:         PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);
478:         if (same) mg->stageApply = st;
479:       }
480:       if (!mg->stageApply) {
481:         PetscLogStageRegister(sname, &mg->stageApply);
482:       }
483:     }
484: #endif
485:   }
486:   PetscOptionsTail();
487:   return(0);
488: }

490: const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
491: const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0};
492: const char *const PCMGGalerkinTypes[] = {"both","pmat","mat","none","external","PCMGGalerkinType","PC_MG_GALERKIN",0};

494:  #include <petscdraw.h>
495: PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
496: {
497:   PC_MG          *mg        = (PC_MG*)pc->data;
498:   PC_MG_Levels   **mglevels = mg->levels;
500:   PetscInt       levels = mglevels ? mglevels[0]->levels : 0,i;
501:   PetscBool      iascii,isbinary,isdraw;

504:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
505:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
506:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
507:   if (iascii) {
508:     const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown";
509:     PetscViewerASCIIPrintf(viewer,"  type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);
510:     if (mg->am == PC_MG_MULTIPLICATIVE) {
511:       PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply);
512:     }
513:     if (mg->galerkin == PC_MG_GALERKIN_BOTH) {
514:       PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");
515:     } else if (mg->galerkin == PC_MG_GALERKIN_PMAT) {
516:       PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for pmat\n");
517:     } else if (mg->galerkin == PC_MG_GALERKIN_MAT) {
518:       PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for mat\n");
519:     } else if (mg->galerkin == PC_MG_GALERKIN_EXTERNAL) {
520:       PetscViewerASCIIPrintf(viewer,"    Using externally compute Galerkin coarse grid matrices\n");
521:     } else {
522:       PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n");
523:     }
524:     if (mg->view){
525:       (*mg->view)(pc,viewer);
526:     }
527:     for (i=0; i<levels; i++) {
528:       if (!i) {
529:         PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);
530:       } else {
531:         PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);
532:       }
533:       PetscViewerASCIIPushTab(viewer);
534:       KSPView(mglevels[i]->smoothd,viewer);
535:       PetscViewerASCIIPopTab(viewer);
536:       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
537:         PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");
538:       } else if (i) {
539:         PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);
540:         PetscViewerASCIIPushTab(viewer);
541:         KSPView(mglevels[i]->smoothu,viewer);
542:         PetscViewerASCIIPopTab(viewer);
543:       }
544:     }
545:   } else if (isbinary) {
546:     for (i=levels-1; i>=0; i--) {
547:       KSPView(mglevels[i]->smoothd,viewer);
548:       if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) {
549:         KSPView(mglevels[i]->smoothu,viewer);
550:       }
551:     }
552:   } else if (isdraw) {
553:     PetscDraw draw;
554:     PetscReal x,w,y,bottom,th;
555:     PetscViewerDrawGetDraw(viewer,0,&draw);
556:     PetscDrawGetCurrentPoint(draw,&x,&y);
557:     PetscDrawStringGetSize(draw,NULL,&th);
558:     bottom = y - th;
559:     for (i=levels-1; i>=0; i--) {
560:       if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) {
561:         PetscDrawPushCurrentPoint(draw,x,bottom);
562:         KSPView(mglevels[i]->smoothd,viewer);
563:         PetscDrawPopCurrentPoint(draw);
564:       } else {
565:         w    = 0.5*PetscMin(1.0-x,x);
566:         PetscDrawPushCurrentPoint(draw,x+w,bottom);
567:         KSPView(mglevels[i]->smoothd,viewer);
568:         PetscDrawPopCurrentPoint(draw);
569:         PetscDrawPushCurrentPoint(draw,x-w,bottom);
570:         KSPView(mglevels[i]->smoothu,viewer);
571:         PetscDrawPopCurrentPoint(draw);
572:       }
573:       PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);
574:       bottom -= th;
575:     }
576:   }
577:   return(0);
578: }

580:  #include <petsc/private/dmimpl.h>
581:  #include <petsc/private/kspimpl.h>

583: /*
584:     Calls setup for the KSP on each level
585: */
586: PetscErrorCode PCSetUp_MG(PC pc)
587: {
588:   PC_MG          *mg        = (PC_MG*)pc->data;
589:   PC_MG_Levels   **mglevels = mg->levels;
591:   PetscInt       i,n;
592:   PC             cpc;
593:   PetscBool      dump = PETSC_FALSE,opsset,use_amat,missinginterpolate = PETSC_FALSE;
594:   Mat            dA,dB;
595:   Vec            tvec;
596:   DM             *dms;
597:   PetscViewer    viewer = 0;
598:   PetscBool      dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE;

601:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels with PCMGSetLevels() before setting up");
602:   n = mglevels[0]->levels;
603:   /* FIX: Move this to PCSetFromOptions_MG? */
604:   if (mg->usedmfornumberoflevels) {
605:     PetscInt levels;
606:     DMGetRefineLevel(pc->dm,&levels);
607:     levels++;
608:     if (levels > n) { /* the problem is now being solved on a finer grid */
609:       PCMGSetLevels(pc,levels,NULL);
610:       n        = levels;
611:       PCSetFromOptions(pc); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
612:       mglevels = mg->levels;
613:     }
614:   }
615:   KSPGetPC(mglevels[0]->smoothd,&cpc);


618:   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
619:   /* so use those from global PC */
620:   /* Is this what we always want? What if user wants to keep old one? */
621:   KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);
622:   if (opsset) {
623:     Mat mmat;
624:     KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);
625:     if (mmat == pc->pmat) opsset = PETSC_FALSE;
626:   }

628:   if (!opsset) {
629:     PCGetUseAmat(pc,&use_amat);
630:     if (use_amat) {
631:       PetscInfo(pc,"Using outer operators to define finest grid operator \n  because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n");
632:       KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);
633:     } else {
634:       PetscInfo(pc,"Using matrix (pmat) operators to define finest grid operator \n  because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n");
635:       KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);
636:     }
637:   }

639:   for (i=n-1; i>0; i--) {
640:     if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) {
641:       missinginterpolate = PETSC_TRUE;
642:       continue;
643:     }
644:   }

646:   KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB);
647:   if (dA == dB) dAeqdB = PETSC_TRUE;
648:   if ((mg->galerkin == PC_MG_GALERKIN_NONE) || (((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_MAT)) && !dAeqdB)) {
649:     needRestricts = PETSC_TRUE;  /* user must compute either mat, pmat, or both so must restrict x to coarser levels */
650:   }


653:   /*
654:    Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS)
655:    Skipping for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs?
656:   */
657:   if (missinginterpolate && pc->dm && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) {
658:         /* construct the interpolation from the DMs */
659:     Mat p;
660:     Vec rscale;
661:     PetscMalloc1(n,&dms);
662:     dms[n-1] = pc->dm;
663:     /* Separately create them so we do not get DMKSP interference between levels */
664:     for (i=n-2; i>-1; i--) {DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);}
665:         /*
666:            Force the mat type of coarse level operator to be AIJ because usually we want to use LU for coarse level.
667:            Notice that it can be overwritten by -mat_type because KSPSetUp() reads command line options.
668:            But it is safe to use -dm_mat_type.

670:            The mat type should not be hardcoded like this, we need to find a better way.
671:     DMSetMatType(dms[0],MATAIJ);
672:     */
673:     for (i=n-2; i>-1; i--) {
674:       DMKSP     kdm;
675:       PetscBool dmhasrestrict, dmhasinject;
676:       KSPSetDM(mglevels[i]->smoothd,dms[i]);
677:       if (!needRestricts) {KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);}
678:       DMGetDMKSPWrite(dms[i],&kdm);
679:       /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
680:        * a bitwise OR of computing the matrix, RHS, and initial iterate. */
681:       kdm->ops->computerhs = NULL;
682:       kdm->rhsctx          = NULL;
683:       if (!mglevels[i+1]->interpolate) {
684:         DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);
685:         PCMGSetInterpolation(pc,i+1,p);
686:         if (rscale) {PCMGSetRScale(pc,i+1,rscale);}
687:         VecDestroy(&rscale);
688:         MatDestroy(&p);
689:       }
690:       DMHasCreateRestriction(dms[i],&dmhasrestrict);
691:       if (dmhasrestrict && !mglevels[i+1]->restrct){
692:         DMCreateRestriction(dms[i],dms[i+1],&p);
693:         PCMGSetRestriction(pc,i+1,p);
694:         MatDestroy(&p);
695:       }
696:       DMHasCreateInjection(dms[i],&dmhasinject);
697:       if (dmhasinject && !mglevels[i+1]->inject){
698:         DMCreateInjection(dms[i],dms[i+1],&p);
699:         PCMGSetInjection(pc,i+1,p);
700:         MatDestroy(&p);
701:       }
702:     }

704:     for (i=n-2; i>-1; i--) {DMDestroy(&dms[i]);}
705:     PetscFree(dms);
706:   }

708:   if (pc->dm && !pc->setupcalled) {
709:     /* finest smoother also gets DM but it is not active, independent of whether galerkin==PC_MG_GALERKIN_EXTERNAL */
710:     KSPSetDM(mglevels[n-1]->smoothd,pc->dm);
711:     KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);
712:   }

714:   if (mg->galerkin < PC_MG_GALERKIN_NONE) {
715:     Mat       A,B;
716:     PetscBool doA = PETSC_FALSE,doB = PETSC_FALSE;
717:     MatReuse  reuse = MAT_INITIAL_MATRIX;

719:     if ((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_BOTH)) doB = PETSC_TRUE;
720:     if ((mg->galerkin == PC_MG_GALERKIN_MAT) || ((mg->galerkin == PC_MG_GALERKIN_BOTH) && (dA != dB))) doA = PETSC_TRUE;
721:     if (pc->setupcalled) reuse = MAT_REUSE_MATRIX;
722:     for (i=n-2; i>-1; i--) {
723:       if (!mglevels[i+1]->restrct && !mglevels[i+1]->interpolate) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must provide interpolation or restriction for each MG level except level 0");
724:       if (!mglevels[i+1]->interpolate) {
725:         PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);
726:       }
727:       if (!mglevels[i+1]->restrct) {
728:         PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);
729:       }
730:       if (reuse == MAT_REUSE_MATRIX) {
731:         KSPGetOperators(mglevels[i]->smoothd,&A,&B);
732:       }
733:       if (doA) {
734:         MatGalerkin(mglevels[i+1]->restrct,dA,mglevels[i+1]->interpolate,reuse,1.0,&A);
735:       }
736:       if (doB) {
737:         MatGalerkin(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,reuse,1.0,&B);
738:       }
739:       /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */
740:       if (!doA && dAeqdB) {
741:         if (reuse == MAT_INITIAL_MATRIX) {PetscObjectReference((PetscObject)B);}
742:         A = B;
743:       } else if (!doA && reuse == MAT_INITIAL_MATRIX ) {
744:         KSPGetOperators(mglevels[i]->smoothd,&A,NULL);
745:         PetscObjectReference((PetscObject)A);
746:       }
747:       if (!doB && dAeqdB) {
748:         if (reuse == MAT_INITIAL_MATRIX) {PetscObjectReference((PetscObject)A);}
749:         B = A;
750:       } else if (!doB && reuse == MAT_INITIAL_MATRIX) {
751:         KSPGetOperators(mglevels[i]->smoothd,NULL,&B);
752:         PetscObjectReference((PetscObject)B);
753:       }
754:       if (reuse == MAT_INITIAL_MATRIX) {
755:         KSPSetOperators(mglevels[i]->smoothd,A,B);
756:         PetscObjectDereference((PetscObject)A);
757:         PetscObjectDereference((PetscObject)B);
758:       }
759:       dA = A;
760:       dB = B;
761:     }
762:   }
763:   if (needRestricts && pc->dm && pc->dm->x) {
764:     /* need to restrict Jacobian location to coarser meshes for evaluation */
765:     for (i=n-2; i>-1; i--) {
766:       Mat R;
767:       Vec rscale;
768:       if (!mglevels[i]->smoothd->dm->x) {
769:         Vec *vecs;
770:         KSPCreateVecs(mglevels[i]->smoothd,1,&vecs,0,NULL);
771:         mglevels[i]->smoothd->dm->x = vecs[0];
772:         PetscFree(vecs);
773:       }
774:       PCMGGetRestriction(pc,i+1,&R);
775:       PCMGGetRScale(pc,i+1,&rscale);
776:       MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);
777:       VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);
778:     }
779:   }
780:   if (needRestricts && pc->dm) {
781:     for (i=n-2; i>=0; i--) {
782:       DM  dmfine,dmcoarse;
783:       Mat Restrict,Inject;
784:       Vec rscale;
785:       KSPGetDM(mglevels[i+1]->smoothd,&dmfine);
786:       KSPGetDM(mglevels[i]->smoothd,&dmcoarse);
787:       PCMGGetRestriction(pc,i+1,&Restrict);
788:       PCMGGetRScale(pc,i+1,&rscale);
789:       PCMGGetInjection(pc,i+1,&Inject);
790:       DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);
791:     }
792:   }

794:   if (!pc->setupcalled) {
795:     for (i=0; i<n; i++) {
796:       KSPSetFromOptions(mglevels[i]->smoothd);
797:     }
798:     for (i=1; i<n; i++) {
799:       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
800:         KSPSetFromOptions(mglevels[i]->smoothu);
801:       }
802:     }
803:     /* insure that if either interpolation or restriction is set the other other one is set */
804:     for (i=1; i<n; i++) {
805:       PCMGGetInterpolation(pc,i,NULL);
806:       PCMGGetRestriction(pc,i,NULL);
807:     }
808:     for (i=0; i<n-1; i++) {
809:       if (!mglevels[i]->b) {
810:         Vec *vec;
811:         KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);
812:         PCMGSetRhs(pc,i,*vec);
813:         VecDestroy(vec);
814:         PetscFree(vec);
815:       }
816:       if (!mglevels[i]->r && i) {
817:         VecDuplicate(mglevels[i]->b,&tvec);
818:         PCMGSetR(pc,i,tvec);
819:         VecDestroy(&tvec);
820:       }
821:       if (!mglevels[i]->x) {
822:         VecDuplicate(mglevels[i]->b,&tvec);
823:         PCMGSetX(pc,i,tvec);
824:         VecDestroy(&tvec);
825:       }
826:     }
827:     if (n != 1 && !mglevels[n-1]->r) {
828:       /* PCMGSetR() on the finest level if user did not supply it */
829:       Vec *vec;
830:       KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);
831:       PCMGSetR(pc,n-1,*vec);
832:       VecDestroy(vec);
833:       PetscFree(vec);
834:     }
835:   }

837:   if (pc->dm) {
838:     /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
839:     for (i=0; i<n-1; i++) {
840:       if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
841:     }
842:   }

844:   for (i=1; i<n; i++) {
845:     if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1){
846:       /* if doing only down then initial guess is zero */
847:       KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);
848:     }
849:     if (mglevels[i]->eventsmoothsetup) {PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);}
850:     KSPSetUp(mglevels[i]->smoothd);
851:     if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PC_FAILED) {
852:       pc->failedreason = PC_SUBPC_ERROR;
853:     }
854:     if (mglevels[i]->eventsmoothsetup) {PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);}
855:     if (!mglevels[i]->residual) {
856:       Mat mat;
857:       KSPGetOperators(mglevels[i]->smoothd,&mat,NULL);
858:       PCMGSetResidual(pc,i,PCMGResidualDefault,mat);
859:     }
860:   }
861:   for (i=1; i<n; i++) {
862:     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
863:       Mat downmat,downpmat;

865:       /* check if operators have been set for up, if not use down operators to set them */
866:       KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);
867:       if (!opsset) {
868:         KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);
869:         KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);
870:       }

872:       KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);
873:       if (mglevels[i]->eventsmoothsetup) {PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);}
874:       KSPSetUp(mglevels[i]->smoothu);
875:       if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PC_FAILED) {
876:         pc->failedreason = PC_SUBPC_ERROR;
877:       }
878:       if (mglevels[i]->eventsmoothsetup) {PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);}
879:     }
880:   }

882:   if (mglevels[0]->eventsmoothsetup) {PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);}
883:   KSPSetUp(mglevels[0]->smoothd);
884:   if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PC_FAILED) {
885:     pc->failedreason = PC_SUBPC_ERROR;
886:   }
887:   if (mglevels[0]->eventsmoothsetup) {PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);}

889:   /*
890:      Dump the interpolation/restriction matrices plus the
891:    Jacobian/stiffness on each level. This allows MATLAB users to
892:    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.

894:    Only support one or the other at the same time.
895:   */
896: #if defined(PETSC_USE_SOCKET_VIEWER)
897:   PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);
898:   if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc));
899:   dump = PETSC_FALSE;
900: #endif
901:   PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);
902:   if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc));

904:   if (viewer) {
905:     for (i=1; i<n; i++) {
906:       MatView(mglevels[i]->restrct,viewer);
907:     }
908:     for (i=0; i<n; i++) {
909:       KSPGetPC(mglevels[i]->smoothd,&pc);
910:       MatView(pc->mat,viewer);
911:     }
912:   }
913:   return(0);
914: }

916: /* -------------------------------------------------------------------------------------*/

918: PetscErrorCode PCMGGetLevels_MG(PC pc, PetscInt *levels)
919: {
920:   PC_MG *mg = (PC_MG *) pc->data;

923:   *levels = mg->nlevels;
924:   return(0);
925: }

927: /*@
928:    PCMGGetLevels - Gets the number of levels to use with MG.

930:    Not Collective

932:    Input Parameter:
933: .  pc - the preconditioner context

935:    Output parameter:
936: .  levels - the number of levels

938:    Level: advanced

940: .keywords: MG, get, levels, multigrid

942: .seealso: PCMGSetLevels()
943: @*/
944: PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels)
945: {

951:   *levels = 0;
952:   PetscTryMethod(pc,"PCMGGetLevels_C",(PC,PetscInt*),(pc,levels));
953:   return(0);
954: }

956: /*@
957:    PCMGSetType - Determines the form of multigrid to use:
958:    multiplicative, additive, full, or the Kaskade algorithm.

960:    Logically Collective on PC

962:    Input Parameters:
963: +  pc - the preconditioner context
964: -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
965:    PC_MG_FULL, PC_MG_KASKADE

967:    Options Database Key:
968: .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
969:    additive, full, kaskade

971:    Level: advanced

973: .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid

975: .seealso: PCMGSetLevels()
976: @*/
977: PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
978: {
979:   PC_MG *mg = (PC_MG*)pc->data;

984:   mg->am = form;
985:   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
986:   else pc->ops->applyrichardson = NULL;
987:   return(0);
988: }

990: /*@
991:    PCMGGetType - Determines the form of multigrid to use:
992:    multiplicative, additive, full, or the Kaskade algorithm.

994:    Logically Collective on PC

996:    Input Parameter:
997: .  pc - the preconditioner context

999:    Output Parameter:
1000: .  type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE


1003:    Level: advanced

1005: .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid

1007: .seealso: PCMGSetLevels()
1008: @*/
1009: PetscErrorCode  PCMGGetType(PC pc,PCMGType *type)
1010: {
1011:   PC_MG *mg = (PC_MG*)pc->data;

1015:   *type = mg->am;
1016:   return(0);
1017: }

1019: /*@
1020:    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
1021:    complicated cycling.

1023:    Logically Collective on PC

1025:    Input Parameters:
1026: +  pc - the multigrid context
1027: -  n - either PC_MG_CYCLE_V or PC_MG_CYCLE_W

1029:    Options Database Key:
1030: .  -pc_mg_cycle_type <v,w> - provide the cycle desired

1032:    Level: advanced

1034: .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid

1036: .seealso: PCMGSetCycleTypeOnLevel()
1037: @*/
1038: PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
1039: {
1040:   PC_MG        *mg        = (PC_MG*)pc->data;
1041:   PC_MG_Levels **mglevels = mg->levels;
1042:   PetscInt     i,levels;

1047:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1048:   levels = mglevels[0]->levels;
1049:   for (i=0; i<levels; i++) mglevels[i]->cycles = n;
1050:   return(0);
1051: }

1053: /*@
1054:    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
1055:          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used

1057:    Logically Collective on PC

1059:    Input Parameters:
1060: +  pc - the multigrid context
1061: -  n - number of cycles (default is 1)

1063:    Options Database Key:
1064: .  -pc_mg_multiplicative_cycles n

1066:    Level: advanced

1068:    Notes:
1069:     This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()

1071: .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid

1073: .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
1074: @*/
1075: PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
1076: {
1077:   PC_MG        *mg        = (PC_MG*)pc->data;

1082:   mg->cyclesperpcapply = n;
1083:   return(0);
1084: }

1086: PetscErrorCode PCMGSetGalerkin_MG(PC pc,PCMGGalerkinType use)
1087: {
1088:   PC_MG *mg = (PC_MG*)pc->data;

1091:   mg->galerkin = use;
1092:   return(0);
1093: }

1095: /*@
1096:    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
1097:       finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i

1099:    Logically Collective on PC

1101:    Input Parameters:
1102: +  pc - the multigrid context
1103: -  use - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, or PC_MG_GALERKIN_NONE

1105:    Options Database Key:
1106: .  -pc_mg_galerkin <both,pmat,mat,none>

1108:    Level: intermediate

1110:    Notes:
1111:     Some codes that use PCMG such as PCGAMG use Galerkin internally while constructing the hierarchy and thus do not
1112:      use the PCMG construction of the coarser grids.

1114: .keywords: MG, set, Galerkin

1116: .seealso: PCMGGetGalerkin(), PCMGGalerkinType

1118: @*/
1119: PetscErrorCode PCMGSetGalerkin(PC pc,PCMGGalerkinType use)
1120: {

1125:   PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PCMGGalerkinType),(pc,use));
1126:   return(0);
1127: }

1129: /*@
1130:    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
1131:       A_i-1 = r_i * A_i * p_i

1133:    Not Collective

1135:    Input Parameter:
1136: .  pc - the multigrid context

1138:    Output Parameter:
1139: .  galerkin - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, PC_MG_GALERKIN_NONE, or PC_MG_GALERKIN_EXTERNAL

1141:    Level: intermediate

1143: .keywords: MG, set, Galerkin

1145: .seealso: PCMGSetGalerkin(), PCMGGalerkinType

1147: @*/
1148: PetscErrorCode  PCMGGetGalerkin(PC pc,PCMGGalerkinType  *galerkin)
1149: {
1150:   PC_MG *mg = (PC_MG*)pc->data;

1154:   *galerkin = mg->galerkin;
1155:   return(0);
1156: }

1158: /*@
1159:    PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use
1160:    on all levels.  Use PCMGDistinctSmoothUp() to create separate up and down smoothers if you want different numbers of
1161:    pre- and post-smoothing steps.

1163:    Logically Collective on PC

1165:    Input Parameters:
1166: +  mg - the multigrid context
1167: -  n - the number of smoothing steps

1169:    Options Database Key:
1170: +  -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps

1172:    Level: advanced

1174:    Notes:
1175:     this does not set a value on the coarsest grid, since we assume that
1176:     there is no separate smooth up on the coarsest grid.

1178: .keywords: MG, smooth, up, post-smoothing, steps, multigrid

1180: .seealso: PCMGSetDistinctSmoothUp()
1181: @*/
1182: PetscErrorCode  PCMGSetNumberSmooth(PC pc,PetscInt n)
1183: {
1184:   PC_MG          *mg        = (PC_MG*)pc->data;
1185:   PC_MG_Levels   **mglevels = mg->levels;
1187:   PetscInt       i,levels;

1192:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1193:   levels = mglevels[0]->levels;

1195:   for (i=1; i<levels; i++) {
1196:     KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);
1197:     KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);
1198:     mg->default_smoothu = n;
1199:     mg->default_smoothd = n;
1200:   }
1201:   return(0);
1202: }

1204: /*@
1205:    PCMGSetDistinctSmoothUp - sets the up (post) smoother to be a seperate KSP from the down (pre) smoother on all levels
1206:        and adds the suffix _up to the options name

1208:    Logically Collective on PC

1210:    Input Parameters:
1211: .  pc - the preconditioner context

1213:    Options Database Key:
1214: .  -pc_mg_distinct_smoothup

1216:    Level: advanced

1218:    Notes:
1219:     this does not set a value on the coarsest grid, since we assume that
1220:     there is no separate smooth up on the coarsest grid.

1222: .keywords: MG, smooth, up, post-smoothing, steps, multigrid

1224: .seealso: PCMGSetNumberSmooth()
1225: @*/
1226: PetscErrorCode  PCMGSetDistinctSmoothUp(PC pc)
1227: {
1228:   PC_MG          *mg        = (PC_MG*)pc->data;
1229:   PC_MG_Levels   **mglevels = mg->levels;
1231:   PetscInt       i,levels;
1232:   KSP            subksp;

1236:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1237:   levels = mglevels[0]->levels;

1239:   for (i=1; i<levels; i++) {
1240:     const char *prefix = NULL;
1241:     /* make sure smoother up and down are different */
1242:     PCMGGetSmootherUp(pc,i,&subksp);
1243:     KSPGetOptionsPrefix(mglevels[i]->smoothd,&prefix);
1244:     KSPSetOptionsPrefix(subksp,prefix);
1245:     KSPAppendOptionsPrefix(subksp,"up_");
1246:   }
1247:   return(0);
1248: }

1250: /* ----------------------------------------------------------------------------------------*/

1252: /*MC
1253:    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
1254:     information about the coarser grid matrices and restriction/interpolation operators.

1256:    Options Database Keys:
1257: +  -pc_mg_levels <nlevels> - number of levels including finest
1258: .  -pc_mg_cycle_type <v,w> - provide the cycle desired
1259: .  -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
1260: .  -pc_mg_log - log information about time spent on each level of the solver
1261: .  -pc_mg_distinct_smoothup - configure up (after interpolation) and down (before restriction) smoothers separately (with different options prefixes)
1262: .  -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
1263: .  -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
1264: .  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1265:                         to the Socket viewer for reading from MATLAB.
1266: -  -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices
1267:                         to the binary output file called binaryoutput

1269:    Notes:
1270:     If one uses a Krylov method such GMRES or CG as the smoother then one must use KSPFGMRES, KSPGCR, or KSPRICHARDSON as the outer Krylov method

1272:        When run with a single level the smoother options are used on that level NOT the coarse grid solver options

1274:        When run with KSPRICHARDSON the convergence test changes slightly if monitor is turned on. The iteration count may change slightly. This
1275:        is because without monitoring the residual norm is computed WITHIN each multigrid cycle on the finest level after the pre-smoothing
1276:        (because the residual has just been computed for the multigrid algorithm and is hence available for free) while with monitoring the
1277:        residual is computed at the end of each cycle.

1279:    Level: intermediate

1281:    Concepts: multigrid/multilevel

1283: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
1284:            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(),
1285:            PCMGSetDistinctSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
1286:            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
1287:            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
1288: M*/

1290: PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc)
1291: {
1292:   PC_MG          *mg;

1296:   PetscNewLog(pc,&mg);
1297:   pc->data     = (void*)mg;
1298:   mg->nlevels  = -1;
1299:   mg->am       = PC_MG_MULTIPLICATIVE;
1300:   mg->galerkin = PC_MG_GALERKIN_NONE;

1302:   pc->useAmat = PETSC_TRUE;

1304:   pc->ops->apply          = PCApply_MG;
1305:   pc->ops->setup          = PCSetUp_MG;
1306:   pc->ops->reset          = PCReset_MG;
1307:   pc->ops->destroy        = PCDestroy_MG;
1308:   pc->ops->setfromoptions = PCSetFromOptions_MG;
1309:   pc->ops->view           = PCView_MG;

1311:   PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);
1312:   PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG);
1313:   PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG);
1314:   return(0);
1315: }