Actual source code: mg.c

petsc-master 2018-05-23
Report Typos and Errors

  2: /*
  3:     Defines the multigrid preconditioner interface.
  4: */
  5:  #include <petsc/private/pcmgimpl.h>
  6:  #include <petscdm.h>

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

 18:   if (mglevels->eventsmoothsolve) {PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);}
 19:   KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);  /* pre-smooth */
 20:   KSPGetPC(mglevels->smoothd,&subpc);
 21:   PCGetSetUpFailedReason(subpc,&pcreason);
 22:   if (pcreason) {
 23:     pc->failedreason = PC_SUBPC_ERROR;
 24:   }
 25:   if (mglevels->eventsmoothsolve) {PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);}
 26:   if (mglevels->level) {  /* not the coarsest grid */
 27:     if (mglevels->eventresidual) {PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);}
 28:     (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);
 29:     if (mglevels->eventresidual) {PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);}

 31:     /* if on finest level and have convergence criteria set */
 32:     if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) {
 33:       PetscReal rnorm;
 34:       VecNorm(mglevels->r,NORM_2,&rnorm);
 35:       if (rnorm <= mg->ttol) {
 36:         if (rnorm < mg->abstol) {
 37:           *reason = PCRICHARDSON_CONVERGED_ATOL;
 38:           PetscInfo2(pc,"Linear solver has converged. Residual norm %g is less than absolute tolerance %g\n",(double)rnorm,(double)mg->abstol);
 39:         } else {
 40:           *reason = PCRICHARDSON_CONVERGED_RTOL;
 41:           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);
 42:         }
 43:         return(0);
 44:       }
 45:     }

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

 65: 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)
 66: {
 67:   PC_MG          *mg        = (PC_MG*)pc->data;
 68:   PC_MG_Levels   **mglevels = mg->levels;
 70:   PetscInt       levels = mglevels[0]->levels,i;

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

 83:   mg->rtol   = rtol;
 84:   mg->abstol = abstol;
 85:   mg->dtol   = dtol;
 86:   if (rtol) {
 87:     /* compute initial residual norm for relative convergence test */
 88:     PetscReal rnorm;
 89:     if (zeroguess) {
 90:       VecNorm(b,NORM_2,&rnorm);
 91:     } else {
 92:       (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);
 93:       VecNorm(w,NORM_2,&rnorm);
 94:     }
 95:     mg->ttol = PetscMax(rtol*rnorm,abstol);
 96:   } else if (abstol) mg->ttol = abstol;
 97:   else mg->ttol = 0.0;

 99:   /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't
100:      stop prematurely due to small residual */
101:   for (i=1; i<levels; i++) {
102:     KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
103:     if (mglevels[i]->smoothu != mglevels[i]->smoothd) {
104:       /* For Richardson the initial guess is nonzero since it is solving in each cycle the original system not just applying as a preconditioner */
105:       KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);
106:       KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
107:     }
108:   }

110:   *reason = (PCRichardsonConvergedReason)0;
111:   for (i=0; i<its; i++) {
112:     PCMGMCycle_Private(pc,mglevels+levels-1,reason);
113:     if (*reason) break;
114:   }
115:   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
116:   *outits = i;
117:   return(0);
118: }

120: PetscErrorCode PCReset_MG(PC pc)
121: {
122:   PC_MG          *mg        = (PC_MG*)pc->data;
123:   PC_MG_Levels   **mglevels = mg->levels;
125:   PetscInt       i,n;

128:   if (mglevels) {
129:     n = mglevels[0]->levels;
130:     for (i=0; i<n-1; i++) {
131:       VecDestroy(&mglevels[i+1]->r);
132:       VecDestroy(&mglevels[i]->b);
133:       VecDestroy(&mglevels[i]->x);
134:       MatDestroy(&mglevels[i+1]->restrct);
135:       MatDestroy(&mglevels[i+1]->interpolate);
136:       MatDestroy(&mglevels[i+1]->inject);
137:       VecDestroy(&mglevels[i+1]->rscale);
138:     }

140:     for (i=0; i<n; i++) {
141:       MatDestroy(&mglevels[i]->A);
142:       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
143:         KSPReset(mglevels[i]->smoothd);
144:       }
145:       KSPReset(mglevels[i]->smoothu);
146:     }
147:   }
148:   return(0);
149: }

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

155:    Logically Collective on PC

157:    Input Parameters:
158: +  pc - the preconditioner context
159: .  levels - the number of levels
160: -  comms - optional communicators for each level; this is to allow solving the coarser problems
161:            on smaller sets of processors.

163:    Level: intermediate

165:    Notes:
166:      If the number of levels is one then the multigrid uses the -mg_levels prefix
167:   for setting the level options rather than the -mg_coarse prefix.

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

171: .seealso: PCMGSetType(), PCMGGetLevels()
172: @*/
173: PetscErrorCode  PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
174: {
176:   PC_MG          *mg        = (PC_MG*)pc->data;
177:   MPI_Comm       comm;
178:   PC_MG_Levels   **mglevels = mg->levels;
179:   PCMGType       mgtype     = mg->am;
180:   PetscInt       mgctype    = (PetscInt) PC_MG_CYCLE_V;
181:   PetscInt       i;
182:   PetscMPIInt    size;
183:   const char     *prefix;
184:   PC             ipc;
185:   PetscInt       n;

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

207:   mg->nlevels = levels;

209:   PetscMalloc1(levels,&mglevels);
210:   PetscLogObjectMemory((PetscObject)pc,levels*(sizeof(PC_MG*)));

212:   PCGetOptionsPrefix(pc,&prefix);

214:   mg->stageApply = 0;
215:   for (i=0; i<levels; i++) {
216:     PetscNewLog(pc,&mglevels[i]);

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

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

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

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

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

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


275: PetscErrorCode PCDestroy_MG(PC pc)
276: {
278:   PC_MG          *mg        = (PC_MG*)pc->data;
279:   PC_MG_Levels   **mglevels = mg->levels;
280:   PetscInt       i,n;

283:   PCReset_MG(pc);
284:   if (mglevels) {
285:     n = mglevels[0]->levels;
286:     for (i=0; i<n; i++) {
287:       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
288:         KSPDestroy(&mglevels[i]->smoothd);
289:       }
290:       KSPDestroy(&mglevels[i]->smoothu);
291:       PetscFree(mglevels[i]);
292:     }
293:     PetscFree(mg->levels);
294:   }
295:   PetscFree(pc->data);
296:   return(0);
297: }



301: extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**);
302: extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**);
303: extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**);

305: /*
306:    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
307:              or full cycle of multigrid.

309:   Note:
310:   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
311: */
312: static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
313: {
314:   PC_MG          *mg        = (PC_MG*)pc->data;
315:   PC_MG_Levels   **mglevels = mg->levels;
317:   PetscInt       levels = mglevels[0]->levels,i;

320:   if (mg->stageApply) {PetscLogStagePush(mg->stageApply);}
321:   /* When the DM is supplying the matrix then it will not exist until here */
322:   for (i=0; i<levels; i++) {
323:     if (!mglevels[i]->A) {
324:       KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);
325:       PetscObjectReference((PetscObject)mglevels[i]->A);
326:     }
327:   }

329:   mglevels[levels-1]->b = b;
330:   mglevels[levels-1]->x = x;
331:   if (mg->am == PC_MG_MULTIPLICATIVE) {
332:     VecSet(x,0.0);
333:     for (i=0; i<mg->cyclesperpcapply; i++) {
334:       PCMGMCycle_Private(pc,mglevels+levels-1,NULL);
335:     }
336:   } else if (mg->am == PC_MG_ADDITIVE) {
337:     PCMGACycle_Private(pc,mglevels);
338:   } else if (mg->am == PC_MG_KASKADE) {
339:     PCMGKCycle_Private(pc,mglevels);
340:   } else {
341:     PCMGFCycle_Private(pc,mglevels);
342:   }
343:   if (mg->stageApply) {PetscLogStagePop();}
344:   return(0);
345: }


348: PetscErrorCode PCSetFromOptions_MG(PetscOptionItems *PetscOptionsObject,PC pc)
349: {
350:   PetscErrorCode   ierr;
351:   PetscInt         levels,cycles;
352:   PetscBool        flg;
353:   PC_MG            *mg = (PC_MG*)pc->data;
354:   PC_MG_Levels     **mglevels;
355:   PCMGType         mgtype;
356:   PCMGCycleType    mgctype;
357:   PCMGGalerkinType gtype;

360:   levels = PetscMax(mg->nlevels,1);
361:   PetscOptionsHead(PetscOptionsObject,"Multigrid options");
362:   PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);
363:   if (!flg && !mg->levels && pc->dm) {
364:     DMGetRefineLevel(pc->dm,&levels);
365:     levels++;
366:     mg->usedmfornumberoflevels = PETSC_TRUE;
367:   }
368:   PCMGSetLevels(pc,levels,NULL);
369:   mglevels = mg->levels;

371:   mgctype = (PCMGCycleType) mglevels[0]->cycles;
372:   PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);
373:   if (flg) {
374:     PCMGSetCycleType(pc,mgctype);
375:   }
376:   gtype = mg->galerkin;
377:   PetscOptionsEnum("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",PCMGGalerkinTypes,(PetscEnum)gtype,(PetscEnum*)&gtype,&flg);
378:   if (flg) {
379:     PCMGSetGalerkin(pc,gtype);
380:   }
381:   flg = PETSC_FALSE;
382:   PetscOptionsBool("-pc_mg_distinct_smoothup","Create seperate smoothup KSP and append the prefix _up","PCMGSetDistinctSmoothUp",PETSC_FALSE,&flg,NULL);
383:   if (flg) {
384:     PCMGSetDistinctSmoothUp(pc);
385:   }
386:   mgtype = mg->am;
387:   PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);
388:   if (flg) {
389:     PCMGSetType(pc,mgtype);
390:   }
391:   if (mg->am == PC_MG_MULTIPLICATIVE) {
392:     PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGMultiplicativeSetCycles",mg->cyclesperpcapply,&cycles,&flg);
393:     if (flg) {
394:       PCMGMultiplicativeSetCycles(pc,cycles);
395:     }
396:   }
397:   flg  = PETSC_FALSE;
398:   PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);
399:   if (flg) {
400:     PetscInt i;
401:     char     eventname[128];

403:     levels = mglevels[0]->levels;
404:     for (i=0; i<levels; i++) {
405:       sprintf(eventname,"MGSetup Level %d",(int)i);
406:       PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);
407:       sprintf(eventname,"MGSmooth Level %d",(int)i);
408:       PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);
409:       if (i) {
410:         sprintf(eventname,"MGResid Level %d",(int)i);
411:         PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);
412:         sprintf(eventname,"MGInterp Level %d",(int)i);
413:         PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);
414:       }
415:     }

417: #if defined(PETSC_USE_LOG)
418:     {
419:       const char    *sname = "MG Apply";
420:       PetscStageLog stageLog;
421:       PetscInt      st;

424:       PetscLogGetStageLog(&stageLog);
425:       for (st = 0; st < stageLog->numStages; ++st) {
426:         PetscBool same;

428:         PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);
429:         if (same) mg->stageApply = st;
430:       }
431:       if (!mg->stageApply) {
432:         PetscLogStageRegister(sname, &mg->stageApply);
433:       }
434:     }
435: #endif
436:   }
437:   PetscOptionsTail();
438:   return(0);
439: }

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

445:  #include <petscdraw.h>
446: PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
447: {
448:   PC_MG          *mg        = (PC_MG*)pc->data;
449:   PC_MG_Levels   **mglevels = mg->levels;
451:   PetscInt       levels = mglevels ? mglevels[0]->levels : 0,i;
452:   PetscBool      iascii,isbinary,isdraw;

455:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
456:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
457:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
458:   if (iascii) {
459:     const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown";
460:     PetscViewerASCIIPrintf(viewer,"  type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);
461:     if (mg->am == PC_MG_MULTIPLICATIVE) {
462:       PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply);
463:     }
464:     if (mg->galerkin == PC_MG_GALERKIN_BOTH) {
465:       PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");
466:     } else if (mg->galerkin == PC_MG_GALERKIN_PMAT) {
467:       PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for pmat\n");
468:     } else if (mg->galerkin == PC_MG_GALERKIN_MAT) {
469:       PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for mat\n");
470:     } else if (mg->galerkin == PC_MG_GALERKIN_EXTERNAL) {
471:       PetscViewerASCIIPrintf(viewer,"    Using externally compute Galerkin coarse grid matrices\n");
472:     } else {
473:       PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n");
474:     }
475:     if (mg->view){
476:       (*mg->view)(pc,viewer);
477:     }
478:     for (i=0; i<levels; i++) {
479:       if (!i) {
480:         PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);
481:       } else {
482:         PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);
483:       }
484:       PetscViewerASCIIPushTab(viewer);
485:       KSPView(mglevels[i]->smoothd,viewer);
486:       PetscViewerASCIIPopTab(viewer);
487:       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
488:         PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");
489:       } else if (i) {
490:         PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);
491:         PetscViewerASCIIPushTab(viewer);
492:         KSPView(mglevels[i]->smoothu,viewer);
493:         PetscViewerASCIIPopTab(viewer);
494:       }
495:     }
496:   } else if (isbinary) {
497:     for (i=levels-1; i>=0; i--) {
498:       KSPView(mglevels[i]->smoothd,viewer);
499:       if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) {
500:         KSPView(mglevels[i]->smoothu,viewer);
501:       }
502:     }
503:   } else if (isdraw) {
504:     PetscDraw draw;
505:     PetscReal x,w,y,bottom,th;
506:     PetscViewerDrawGetDraw(viewer,0,&draw);
507:     PetscDrawGetCurrentPoint(draw,&x,&y);
508:     PetscDrawStringGetSize(draw,NULL,&th);
509:     bottom = y - th;
510:     for (i=levels-1; i>=0; i--) {
511:       if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) {
512:         PetscDrawPushCurrentPoint(draw,x,bottom);
513:         KSPView(mglevels[i]->smoothd,viewer);
514:         PetscDrawPopCurrentPoint(draw);
515:       } else {
516:         w    = 0.5*PetscMin(1.0-x,x);
517:         PetscDrawPushCurrentPoint(draw,x+w,bottom);
518:         KSPView(mglevels[i]->smoothd,viewer);
519:         PetscDrawPopCurrentPoint(draw);
520:         PetscDrawPushCurrentPoint(draw,x-w,bottom);
521:         KSPView(mglevels[i]->smoothu,viewer);
522:         PetscDrawPopCurrentPoint(draw);
523:       }
524:       PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);
525:       bottom -= th;
526:     }
527:   }
528:   return(0);
529: }

531:  #include <petsc/private/dmimpl.h>
532:  #include <petsc/private/kspimpl.h>

534: /*
535:     Calls setup for the KSP on each level
536: */
537: PetscErrorCode PCSetUp_MG(PC pc)
538: {
539:   PC_MG          *mg        = (PC_MG*)pc->data;
540:   PC_MG_Levels   **mglevels = mg->levels;
542:   PetscInt       i,n;
543:   PC             cpc;
544:   PetscBool      dump = PETSC_FALSE,opsset,use_amat,missinginterpolate = PETSC_FALSE;
545:   Mat            dA,dB;
546:   Vec            tvec;
547:   DM             *dms;
548:   PetscViewer    viewer = 0;
549:   PetscBool      dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE;

552:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels with PCMGSetLevels() before setting up");
553:   n = mglevels[0]->levels;
554:   /* FIX: Move this to PCSetFromOptions_MG? */
555:   if (mg->usedmfornumberoflevels) {
556:     PetscInt levels;
557:     DMGetRefineLevel(pc->dm,&levels);
558:     levels++;
559:     if (levels > n) { /* the problem is now being solved on a finer grid */
560:       PCMGSetLevels(pc,levels,NULL);
561:       n        = levels;
562:       PCSetFromOptions(pc); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
563:       mglevels =  mg->levels;
564:     }
565:   }
566:   KSPGetPC(mglevels[0]->smoothd,&cpc);


569:   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
570:   /* so use those from global PC */
571:   /* Is this what we always want? What if user wants to keep old one? */
572:   KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);
573:   if (opsset) {
574:     Mat mmat;
575:     KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);
576:     if (mmat == pc->pmat) opsset = PETSC_FALSE;
577:   }

579:   if (!opsset) {
580:     PCGetUseAmat(pc,&use_amat);
581:     if(use_amat){
582:       PetscInfo(pc,"Using outer operators to define finest grid operator \n  because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n");
583:       KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);
584:     }
585:     else {
586:       PetscInfo(pc,"Using matrix (pmat) operators to define finest grid operator \n  because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n");
587:       KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);
588:     }
589:   }

591:   for (i=n-1; i>0; i--) {
592:     if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) {
593:       missinginterpolate = PETSC_TRUE;
594:       continue;
595:     }
596:   }

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


605:   /*
606:    Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS)
607:    Skipping for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs?
608:   */
609:   if (missinginterpolate && pc->dm && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) {
610:         /* construct the interpolation from the DMs */
611:     Mat p;
612:     Vec rscale;
613:     PetscMalloc1(n,&dms);
614:     dms[n-1] = pc->dm;
615:     /* Separately create them so we do not get DMKSP interference between levels */
616:     for (i=n-2; i>-1; i--) {DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);}
617:         /*
618:            Force the mat type of coarse level operator to be AIJ because usually we want to use LU for coarse level.
619:            Notice that it can be overwritten by -mat_type because KSPSetUp() reads command line options.
620:            But it is safe to use -dm_mat_type.

622:            The mat type should not be hardcoded like this, we need to find a better way.
623:     DMSetMatType(dms[0],MATAIJ);
624:     */
625:     for (i=n-2; i>-1; i--) {
626:       DMKSP     kdm;
627:       PetscBool dmhasrestrict, dmhasinject;
628:       KSPSetDM(mglevels[i]->smoothd,dms[i]);
629:       if (!needRestricts) {KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);}
630:       DMGetDMKSPWrite(dms[i],&kdm);
631:       /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
632:        * a bitwise OR of computing the matrix, RHS, and initial iterate. */
633:       kdm->ops->computerhs = NULL;
634:       kdm->rhsctx          = NULL;
635:       if (!mglevels[i+1]->interpolate) {
636:         DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);
637:         PCMGSetInterpolation(pc,i+1,p);
638:         if (rscale) {PCMGSetRScale(pc,i+1,rscale);}
639:         VecDestroy(&rscale);
640:         MatDestroy(&p);
641:       }
642:       DMHasCreateRestriction(dms[i],&dmhasrestrict);
643:       if (dmhasrestrict && !mglevels[i+1]->restrct){
644:         DMCreateRestriction(dms[i],dms[i+1],&p);
645:         PCMGSetRestriction(pc,i+1,p);
646:         MatDestroy(&p);
647:       }
648:       DMHasCreateInjection(dms[i],&dmhasinject);
649:       if (dmhasinject && !mglevels[i+1]->inject){
650:         DMCreateInjection(dms[i],dms[i+1],&p);
651:         PCMGSetInjection(pc,i+1,p);
652:         MatDestroy(&p);
653:       }
654:     }

656:     for (i=n-2; i>-1; i--) {DMDestroy(&dms[i]);}
657:     PetscFree(dms);
658:   }

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

666:   if (mg->galerkin < PC_MG_GALERKIN_NONE) {
667:     Mat       A,B;
668:     PetscBool doA = PETSC_FALSE,doB = PETSC_FALSE;
669:     MatReuse  reuse = MAT_INITIAL_MATRIX;

671:     if ((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_BOTH)) doB = PETSC_TRUE;
672:     if ((mg->galerkin == PC_MG_GALERKIN_MAT) || ((mg->galerkin == PC_MG_GALERKIN_BOTH) && (dA != dB))) doA = PETSC_TRUE;
673:     if (pc->setupcalled) reuse = MAT_REUSE_MATRIX;
674:     for (i=n-2; i>-1; i--) {
675:       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");
676:       if (!mglevels[i+1]->interpolate) {
677:         PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);
678:       }
679:       if (!mglevels[i+1]->restrct) {
680:         PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);
681:       }
682:       if (reuse == MAT_REUSE_MATRIX) {
683:         KSPGetOperators(mglevels[i]->smoothd,&A,&B);
684:       }
685:       if (doA) {
686:         MatGalerkin(mglevels[i+1]->restrct,dA,mglevels[i+1]->interpolate,reuse,1.0,&A);
687:       }
688:       if (doB) {
689:         MatGalerkin(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,reuse,1.0,&B);
690:       }
691:       /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */
692:       if (!doA && dAeqdB) {
693:         if (reuse == MAT_INITIAL_MATRIX) {PetscObjectReference((PetscObject)B);}
694:         A = B;
695:       } else if (!doA && reuse == MAT_INITIAL_MATRIX ) {
696:         KSPGetOperators(mglevels[i]->smoothd,&A,NULL);
697:         PetscObjectReference((PetscObject)A);
698:       }
699:       if (!doB && dAeqdB) {
700:         if (reuse == MAT_INITIAL_MATRIX) {PetscObjectReference((PetscObject)A);}
701:         B = A;
702:       } else if (!doB && reuse == MAT_INITIAL_MATRIX) {
703:         KSPGetOperators(mglevels[i]->smoothd,NULL,&B);
704:         PetscObjectReference((PetscObject)B);
705:       }
706:       if (reuse == MAT_INITIAL_MATRIX) {
707:         KSPSetOperators(mglevels[i]->smoothd,A,B);
708:         PetscObjectDereference((PetscObject)A);
709:         PetscObjectDereference((PetscObject)B);
710:       }
711:       dA = A;
712:       dB = B;
713:     }
714:   }
715:   if (needRestricts && pc->dm && pc->dm->x) {
716:     /* need to restrict Jacobian location to coarser meshes for evaluation */
717:     for (i=n-2; i>-1; i--) {
718:       Mat R;
719:       Vec rscale;
720:       if (!mglevels[i]->smoothd->dm->x) {
721:         Vec *vecs;
722:         KSPCreateVecs(mglevels[i]->smoothd,1,&vecs,0,NULL);
723:         mglevels[i]->smoothd->dm->x = vecs[0];
724:         PetscFree(vecs);
725:       }
726:       PCMGGetRestriction(pc,i+1,&R);
727:       PCMGGetRScale(pc,i+1,&rscale);
728:       MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);
729:       VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);
730:     }
731:   }
732:   if (needRestricts && pc->dm) {
733:     for (i=n-2; i>=0; i--) {
734:       DM  dmfine,dmcoarse;
735:       Mat Restrict,Inject;
736:       Vec rscale;
737:       KSPGetDM(mglevels[i+1]->smoothd,&dmfine);
738:       KSPGetDM(mglevels[i]->smoothd,&dmcoarse);
739:       PCMGGetRestriction(pc,i+1,&Restrict);
740:       PCMGGetRScale(pc,i+1,&rscale);
741:       PCMGGetInjection(pc,i+1,&Inject);
742:       DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);
743:     }
744:   }

746:   if (!pc->setupcalled) {
747:     for (i=0; i<n; i++) {
748:       KSPSetFromOptions(mglevels[i]->smoothd);
749:     }
750:     for (i=1; i<n; i++) {
751:       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
752:         KSPSetFromOptions(mglevels[i]->smoothu);
753:       }
754:     }
755:     /* insure that if either interpolation or restriction is set the other other one is set */
756:     for (i=1; i<n; i++) {
757:       PCMGGetInterpolation(pc,i,NULL);
758:       PCMGGetRestriction(pc,i,NULL);
759:     }
760:     for (i=0; i<n-1; i++) {
761:       if (!mglevels[i]->b) {
762:         Vec *vec;
763:         KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);
764:         PCMGSetRhs(pc,i,*vec);
765:         VecDestroy(vec);
766:         PetscFree(vec);
767:       }
768:       if (!mglevels[i]->r && i) {
769:         VecDuplicate(mglevels[i]->b,&tvec);
770:         PCMGSetR(pc,i,tvec);
771:         VecDestroy(&tvec);
772:       }
773:       if (!mglevels[i]->x) {
774:         VecDuplicate(mglevels[i]->b,&tvec);
775:         PCMGSetX(pc,i,tvec);
776:         VecDestroy(&tvec);
777:       }
778:     }
779:     if (n != 1 && !mglevels[n-1]->r) {
780:       /* PCMGSetR() on the finest level if user did not supply it */
781:       Vec *vec;
782:       KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);
783:       PCMGSetR(pc,n-1,*vec);
784:       VecDestroy(vec);
785:       PetscFree(vec);
786:     }
787:   }

789:   if (pc->dm) {
790:     /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
791:     for (i=0; i<n-1; i++) {
792:       if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
793:     }
794:   }

796:   for (i=1; i<n; i++) {
797:     if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1){
798:       /* if doing only down then initial guess is zero */
799:       KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);
800:     }
801:     if (mglevels[i]->eventsmoothsetup) {PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);}
802:     KSPSetUp(mglevels[i]->smoothd);
803:     if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PCSETUP_FAILED) {
804:       pc->failedreason = PC_SUBPC_ERROR;
805:     }
806:     if (mglevels[i]->eventsmoothsetup) {PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);}
807:     if (!mglevels[i]->residual) {
808:       Mat mat;
809:       KSPGetOperators(mglevels[i]->smoothd,&mat,NULL);
810:       PCMGSetResidual(pc,i,PCMGResidualDefault,mat);
811:     }
812:   }
813:   for (i=1; i<n; i++) {
814:     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
815:       Mat downmat,downpmat;

817:       /* check if operators have been set for up, if not use down operators to set them */
818:       KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);
819:       if (!opsset) {
820:         KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);
821:         KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);
822:       }

824:       KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);
825:       if (mglevels[i]->eventsmoothsetup) {PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);}
826:       KSPSetUp(mglevels[i]->smoothu);
827:       if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PCSETUP_FAILED) {
828:         pc->failedreason = PC_SUBPC_ERROR;
829:       }
830:       if (mglevels[i]->eventsmoothsetup) {PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);}
831:     }
832:   }

834:   if (mglevels[0]->eventsmoothsetup) {PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);}
835:   KSPSetUp(mglevels[0]->smoothd);
836:   if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PCSETUP_FAILED) {
837:     pc->failedreason = PC_SUBPC_ERROR;
838:   }
839:   if (mglevels[0]->eventsmoothsetup) {PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);}

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

846:    Only support one or the other at the same time.
847:   */
848: #if defined(PETSC_USE_SOCKET_VIEWER)
849:   PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);
850:   if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc));
851:   dump = PETSC_FALSE;
852: #endif
853:   PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);
854:   if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc));

856:   if (viewer) {
857:     for (i=1; i<n; i++) {
858:       MatView(mglevels[i]->restrct,viewer);
859:     }
860:     for (i=0; i<n; i++) {
861:       KSPGetPC(mglevels[i]->smoothd,&pc);
862:       MatView(pc->mat,viewer);
863:     }
864:   }
865:   return(0);
866: }

868: /* -------------------------------------------------------------------------------------*/

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

873:    Not Collective

875:    Input Parameter:
876: .  pc - the preconditioner context

878:    Output parameter:
879: .  levels - the number of levels

881:    Level: advanced

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

885: .seealso: PCMGSetLevels()
886: @*/
887: PetscErrorCode  PCMGGetLevels(PC pc,PetscInt *levels)
888: {
889:   PC_MG *mg = (PC_MG*)pc->data;

894:   *levels = mg->nlevels;
895:   return(0);
896: }

898: /*@
899:    PCMGSetType - Determines the form of multigrid to use:
900:    multiplicative, additive, full, or the Kaskade algorithm.

902:    Logically Collective on PC

904:    Input Parameters:
905: +  pc - the preconditioner context
906: -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
907:    PC_MG_FULL, PC_MG_KASKADE

909:    Options Database Key:
910: .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
911:    additive, full, kaskade

913:    Level: advanced

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

917: .seealso: PCMGSetLevels()
918: @*/
919: PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
920: {
921:   PC_MG *mg = (PC_MG*)pc->data;

926:   mg->am = form;
927:   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
928:   else pc->ops->applyrichardson = NULL;
929:   return(0);
930: }

932: /*@
933:    PCMGGetType - Determines the form of multigrid to use:
934:    multiplicative, additive, full, or the Kaskade algorithm.

936:    Logically Collective on PC

938:    Input Parameter:
939: .  pc - the preconditioner context

941:    Output Parameter:
942: .  type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE


945:    Level: advanced

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

949: .seealso: PCMGSetLevels()
950: @*/
951: PetscErrorCode  PCMGGetType(PC pc,PCMGType *type)
952: {
953:   PC_MG *mg = (PC_MG*)pc->data;

957:   *type = mg->am;
958:   return(0);
959: }

961: /*@
962:    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
963:    complicated cycling.

965:    Logically Collective on PC

967:    Input Parameters:
968: +  pc - the multigrid context
969: -  n - either PC_MG_CYCLE_V or PC_MG_CYCLE_W

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

974:    Level: advanced

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

978: .seealso: PCMGSetCycleTypeOnLevel()
979: @*/
980: PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
981: {
982:   PC_MG        *mg        = (PC_MG*)pc->data;
983:   PC_MG_Levels **mglevels = mg->levels;
984:   PetscInt     i,levels;

989:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
990:   levels = mglevels[0]->levels;
991:   for (i=0; i<levels; i++) mglevels[i]->cycles = n;
992:   return(0);
993: }

995: /*@
996:    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
997:          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used

999:    Logically Collective on PC

1001:    Input Parameters:
1002: +  pc - the multigrid context
1003: -  n - number of cycles (default is 1)

1005:    Options Database Key:
1006: .  -pc_mg_multiplicative_cycles n

1008:    Level: advanced

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

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

1015: .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
1016: @*/
1017: PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
1018: {
1019:   PC_MG        *mg        = (PC_MG*)pc->data;

1024:   mg->cyclesperpcapply = n;
1025:   return(0);
1026: }

1028: PetscErrorCode PCMGSetGalerkin_MG(PC pc,PCMGGalerkinType use)
1029: {
1030:   PC_MG *mg = (PC_MG*)pc->data;

1033:   mg->galerkin = use;
1034:   return(0);
1035: }

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

1041:    Logically Collective on PC

1043:    Input Parameters:
1044: +  pc - the multigrid context
1045: -  use - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, or PC_MG_GALERKIN_NONE

1047:    Options Database Key:
1048: .  -pc_mg_galerkin <both,pmat,mat,none>

1050:    Level: intermediate

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

1056: .keywords: MG, set, Galerkin

1058: .seealso: PCMGGetGalerkin(), PCMGGalerkinType

1060: @*/
1061: PetscErrorCode PCMGSetGalerkin(PC pc,PCMGGalerkinType use)
1062: {

1067:   PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PCMGGalerkinType),(pc,use));
1068:   return(0);
1069: }

1071: /*@
1072:    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
1073:       A_i-1 = r_i * A_i * p_i

1075:    Not Collective

1077:    Input Parameter:
1078: .  pc - the multigrid context

1080:    Output Parameter:
1081: .  galerkin - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, PC_MG_GALERKIN_NONE, or PC_MG_GALERKIN_EXTERNAL

1083:    Level: intermediate

1085: .keywords: MG, set, Galerkin

1087: .seealso: PCMGSetGalerkin(), PCMGGalerkinType

1089: @*/
1090: PetscErrorCode  PCMGGetGalerkin(PC pc,PCMGGalerkinType  *galerkin)
1091: {
1092:   PC_MG *mg = (PC_MG*)pc->data;

1096:   *galerkin = mg->galerkin;
1097:   return(0);
1098: }

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

1105:    Logically Collective on PC

1107:    Input Parameters:
1108: +  mg - the multigrid context
1109: -  n - the number of smoothing steps

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

1114:    Level: advanced

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

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

1122: .seealso: PCMGSetDistinctSmoothUp()
1123: @*/
1124: PetscErrorCode  PCMGSetNumberSmooth(PC pc,PetscInt n)
1125: {
1126:   PC_MG          *mg        = (PC_MG*)pc->data;
1127:   PC_MG_Levels   **mglevels = mg->levels;
1129:   PetscInt       i,levels;

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

1137:   for (i=1; i<levels; i++) {
1138:     KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);
1139:     KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);
1140:     mg->default_smoothu = n;
1141:     mg->default_smoothd = n;
1142:   }
1143:   return(0);
1144: }

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

1150:    Logically Collective on PC

1152:    Input Parameters:
1153: .  pc - the preconditioner context

1155:    Options Database Key:
1156: .  -pc_mg_distinct_smoothup

1158:    Level: advanced

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

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

1166: .seealso: PCMGSetNumberSmooth()
1167: @*/
1168: PetscErrorCode  PCMGSetDistinctSmoothUp(PC pc)
1169: {
1170:   PC_MG          *mg        = (PC_MG*)pc->data;
1171:   PC_MG_Levels   **mglevels = mg->levels;
1173:   PetscInt       i,levels;
1174:   KSP            subksp;

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

1181:   for (i=1; i<levels; i++) {
1182:     const char *prefix = NULL;
1183:     /* make sure smoother up and down are different */
1184:     PCMGGetSmootherUp(pc,i,&subksp);
1185:     KSPGetOptionsPrefix(mglevels[i]->smoothd,&prefix);
1186:     KSPSetOptionsPrefix(subksp,prefix);
1187:     KSPAppendOptionsPrefix(subksp,"up_");
1188:   }
1189:   return(0);
1190: }

1192: /* ----------------------------------------------------------------------------------------*/

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

1198:    Options Database Keys:
1199: +  -pc_mg_levels <nlevels> - number of levels including finest
1200: .  -pc_mg_cycle_type <v,w> - 
1201: .  -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
1202: .  -pc_mg_log - log information about time spent on each level of the solver
1203: .  -pc_mg_distinct_smoothup - configure up (after interpolation) and down (before restriction) smoothers separately (with different options prefixes)
1204: .  -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
1205: .  -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
1206: .  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1207:                         to the Socket viewer for reading from MATLAB.
1208: -  -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices
1209:                         to the binary output file called binaryoutput

1211:    Notes:
1212:     If one uses a Krylov method such GMRES or CG as the smoother than one must use KSPFGMRES, KSPGCG, or KSPRICHARDSON as the outer Krylov method

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

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

1221:    Level: intermediate

1223:    Concepts: multigrid/multilevel

1225: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
1226:            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(),
1227:            PCMGSetDistinctSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
1228:            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
1229:            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
1230: M*/

1232: PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc)
1233: {
1234:   PC_MG          *mg;

1238:   PetscNewLog(pc,&mg);
1239:   pc->data     = (void*)mg;
1240:   mg->nlevels  = -1;
1241:   mg->am       = PC_MG_MULTIPLICATIVE;
1242:   mg->galerkin = PC_MG_GALERKIN_NONE;

1244:   pc->useAmat = PETSC_TRUE;

1246:   pc->ops->apply          = PCApply_MG;
1247:   pc->ops->setup          = PCSetUp_MG;
1248:   pc->ops->reset          = PCReset_MG;
1249:   pc->ops->destroy        = PCDestroy_MG;
1250:   pc->ops->setfromoptions = PCSetFromOptions_MG;
1251:   pc->ops->view           = PCView_MG;

1253:   PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);
1254:   return(0);
1255: }