Actual source code: mg.c

petsc-master 2017-02-19
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:       KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
105:     }
106:   }

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

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

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

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

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

152:    Logically Collective on PC

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

160:    Level: intermediate

162:    Notes:
163:      If the number of levels is one then the multigrid uses the -mg_levels prefix
164:   for setting the level options rather than the -mg_coarse prefix.

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

168: .seealso: PCMGSetType(), PCMGGetLevels()
169: @*/
170: PetscErrorCode  PCMGSetLevels(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:   PetscObjectGetComm((PetscObject)pc,&comm);
188:   if (mg->nlevels == levels) return(0);
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: }


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

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



298: extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**);
299: extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**);
300: extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**);

302: /*
303:    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
304:              or full cycle of multigrid.

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

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

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


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

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

369:   mgctype = (PCMGCycleType) mglevels[0]->cycles;
370:   PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);
371:   if (flg) {
372:     PCMGSetCycleType(pc,mgctype);
373:   }
374:   gtype = mg->galerkin;
375:   PetscOptionsEnum("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",PCMGGalerkinTypes,(PetscEnum)gtype,(PetscEnum*)&gtype,&flg);
376:   if (flg) {
377:     PCMGSetGalerkin(pc,gtype);
378:   }
379:   PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",mg->default_smoothu,&m,&flg);
380:   if (flg) {
381:     PCMGSetNumberSmoothUp(pc,m);
382:   }
383:   PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",mg->default_smoothd,&m,&flg);
384:   if (flg) {
385:     PCMGSetNumberSmoothDown(pc,m);
386:   }
387:   mgtype = mg->am;
388:   PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);
389:   if (flg) {
390:     PCMGSetType(pc,mgtype);
391:   }
392:   if (mg->am == PC_MG_MULTIPLICATIVE) {
393:     PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGMultiplicativeSetCycles",mg->cyclesperpcapply,&cycles,&flg);
394:     if (flg) {
395:       PCMGMultiplicativeSetCycles(pc,cycles);
396:     }
397:   }
398:   flg  = PETSC_FALSE;
399:   PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);
400:   if (flg) {
401:     PetscInt i;
402:     char     eventname[128];
403:     if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
404:     levels = mglevels[0]->levels;
405:     for (i=0; i<levels; i++) {
406:       sprintf(eventname,"MGSetup Level %d",(int)i);
407:       PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);
408:       sprintf(eventname,"MGSmooth Level %d",(int)i);
409:       PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);
410:       if (i) {
411:         sprintf(eventname,"MGResid Level %d",(int)i);
412:         PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);
413:         sprintf(eventname,"MGInterp Level %d",(int)i);
414:         PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);
415:       }
416:     }

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

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

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

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

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

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

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

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

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


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

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

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

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


604:   /*
605:    Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS)
606:    Skipping for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs?
607:   */
608:   if (missinginterpolate && pc->dm && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) {
609:     /* construct the interpolation from the DMs */
610:     Mat p;
611:     Vec rscale;
612:     PetscMalloc1(n,&dms);
613:     dms[n-1] = pc->dm;
614:     /* Separately create them so we do not get DMKSP interference between levels */
615:     for (i=n-2; i>-1; i--) {DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);}
616:     for (i=n-2; i>-1; i--) {
617:       DMKSP     kdm;
618:       PetscBool dmhasrestrict;
619:       KSPSetDM(mglevels[i]->smoothd,dms[i]);
620:       if (!needRestricts) {KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);}
621:       DMGetDMKSPWrite(dms[i],&kdm);
622:       /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
623:        * a bitwise OR of computing the matrix, RHS, and initial iterate. */
624:       kdm->ops->computerhs = NULL;
625:       kdm->rhsctx          = NULL;
626:       if (!mglevels[i+1]->interpolate) {
627:         DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);
628:         PCMGSetInterpolation(pc,i+1,p);
629:         if (rscale) {PCMGSetRScale(pc,i+1,rscale);}
630:         VecDestroy(&rscale);
631:         MatDestroy(&p);
632:       }
633:       DMHasCreateRestriction(dms[i],&dmhasrestrict);
634:       if (dmhasrestrict && !mglevels[i+1]->restrct){
635:         DMCreateRestriction(dms[i],dms[i+1],&p);
636:         PCMGSetRestriction(pc,i+1,p);
637:         MatDestroy(&p);
638:       }
639:     }

641:     for (i=n-2; i>-1; i--) {DMDestroy(&dms[i]);}
642:     PetscFree(dms);
643:   }

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

651:   if (mg->galerkin < PC_MG_GALERKIN_NONE) {
652:     Mat       A,B;
653:     PetscBool doA = PETSC_FALSE,doB = PETSC_FALSE;
654:     MatReuse  reuse = MAT_INITIAL_MATRIX;

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

731:   if (!pc->setupcalled) {
732:     for (i=0; i<n; i++) {
733:       KSPSetFromOptions(mglevels[i]->smoothd);
734:     }
735:     for (i=1; i<n; i++) {
736:       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
737:         KSPSetFromOptions(mglevels[i]->smoothu);
738:       }
739:     }
740:     /* insure that if either interpolation or restriction is set the other other one is set */
741:     for (i=1; i<n; i++) {
742:       PCMGGetInterpolation(pc,i,NULL);
743:       PCMGGetRestriction(pc,i,NULL);
744:     }
745:     for (i=0; i<n-1; i++) {
746:       if (!mglevels[i]->b) {
747:         Vec *vec;
748:         KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);
749:         PCMGSetRhs(pc,i,*vec);
750:         VecDestroy(vec);
751:         PetscFree(vec);
752:       }
753:       if (!mglevels[i]->r && i) {
754:         VecDuplicate(mglevels[i]->b,&tvec);
755:         PCMGSetR(pc,i,tvec);
756:         VecDestroy(&tvec);
757:       }
758:       if (!mglevels[i]->x) {
759:         VecDuplicate(mglevels[i]->b,&tvec);
760:         PCMGSetX(pc,i,tvec);
761:         VecDestroy(&tvec);
762:       }
763:     }
764:     if (n != 1 && !mglevels[n-1]->r) {
765:       /* PCMGSetR() on the finest level if user did not supply it */
766:       Vec *vec;
767:       KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);
768:       PCMGSetR(pc,n-1,*vec);
769:       VecDestroy(vec);
770:       PetscFree(vec);
771:     }
772:   }

774:   if (pc->dm) {
775:     /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
776:     for (i=0; i<n-1; i++) {
777:       if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
778:     }
779:   }

781:   for (i=1; i<n; i++) {
782:     if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1){
783:       /* if doing only down then initial guess is zero */
784:       KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);
785:     }
786:     if (mglevels[i]->eventsmoothsetup) {PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);}
787:     KSPSetUp(mglevels[i]->smoothd);
788:     if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PCSETUP_FAILED) {
789:       pc->failedreason = PC_SUBPC_ERROR;
790:     }
791:     if (mglevels[i]->eventsmoothsetup) {PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);}
792:     if (!mglevels[i]->residual) {
793:       Mat mat;
794:       KSPGetOperators(mglevels[i]->smoothd,NULL,&mat);
795:       PCMGSetResidual(pc,i,PCMGResidualDefault,mat);
796:     }
797:   }
798:   for (i=1; i<n; i++) {
799:     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
800:       Mat downmat,downpmat;

802:       /* check if operators have been set for up, if not use down operators to set them */
803:       KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);
804:       if (!opsset) {
805:         KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);
806:         KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);
807:       }

809:       KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);
810:       if (mglevels[i]->eventsmoothsetup) {PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);}
811:       KSPSetUp(mglevels[i]->smoothu);
812:       if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PCSETUP_FAILED) {
813:         pc->failedreason = PC_SUBPC_ERROR;
814:       }
815:       if (mglevels[i]->eventsmoothsetup) {PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);}
816:     }
817:   }

819:   if (mglevels[0]->eventsmoothsetup) {PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);}
820:   KSPSetUp(mglevels[0]->smoothd);
821:   if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PCSETUP_FAILED) {
822:     pc->failedreason = PC_SUBPC_ERROR;
823:   }
824:   if (mglevels[0]->eventsmoothsetup) {PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);}

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

831:    Only support one or the other at the same time.
832:   */
833: #if defined(PETSC_USE_SOCKET_VIEWER)
834:   PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);
835:   if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc));
836:   dump = PETSC_FALSE;
837: #endif
838:   PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);
839:   if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc));

841:   if (viewer) {
842:     for (i=1; i<n; i++) {
843:       MatView(mglevels[i]->restrct,viewer);
844:     }
845:     for (i=0; i<n; i++) {
846:       KSPGetPC(mglevels[i]->smoothd,&pc);
847:       MatView(pc->mat,viewer);
848:     }
849:   }
850:   return(0);
851: }

853: /* -------------------------------------------------------------------------------------*/

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

858:    Not Collective

860:    Input Parameter:
861: .  pc - the preconditioner context

863:    Output parameter:
864: .  levels - the number of levels

866:    Level: advanced

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

870: .seealso: PCMGSetLevels()
871: @*/
872: PetscErrorCode  PCMGGetLevels(PC pc,PetscInt *levels)
873: {
874:   PC_MG *mg = (PC_MG*)pc->data;

879:   *levels = mg->nlevels;
880:   return(0);
881: }

883: /*@
884:    PCMGSetType - Determines the form of multigrid to use:
885:    multiplicative, additive, full, or the Kaskade algorithm.

887:    Logically Collective on PC

889:    Input Parameters:
890: +  pc - the preconditioner context
891: -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
892:    PC_MG_FULL, PC_MG_KASKADE

894:    Options Database Key:
895: .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
896:    additive, full, kaskade

898:    Level: advanced

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

902: .seealso: PCMGSetLevels()
903: @*/
904: PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
905: {
906:   PC_MG *mg = (PC_MG*)pc->data;

911:   mg->am = form;
912:   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
913:   else pc->ops->applyrichardson = NULL;
914:   return(0);
915: }

917: /*@
918:    PCMGGetType - Determines the form of multigrid to use:
919:    multiplicative, additive, full, or the Kaskade algorithm.

921:    Logically Collective on PC

923:    Input Parameter:
924: .  pc - the preconditioner context

926:    Output Parameter:
927: .  type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE


930:    Level: advanced

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

934: .seealso: PCMGSetLevels()
935: @*/
936: PetscErrorCode  PCMGGetType(PC pc,PCMGType *type)
937: {
938:   PC_MG *mg = (PC_MG*)pc->data;

942:   *type = mg->am;
943:   return(0);
944: }

946: /*@
947:    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
948:    complicated cycling.

950:    Logically Collective on PC

952:    Input Parameters:
953: +  pc - the multigrid context
954: -  PC_MG_CYCLE_V or PC_MG_CYCLE_W

956:    Options Database Key:
957: .  -pc_mg_cycle_type <v,w>

959:    Level: advanced

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

963: .seealso: PCMGSetCycleTypeOnLevel()
964: @*/
965: PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
966: {
967:   PC_MG        *mg        = (PC_MG*)pc->data;
968:   PC_MG_Levels **mglevels = mg->levels;
969:   PetscInt     i,levels;

973:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
975:   levels = mglevels[0]->levels;

977:   for (i=0; i<levels; i++) mglevels[i]->cycles = n;
978:   return(0);
979: }

981: /*@
982:    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
983:          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used

985:    Logically Collective on PC

987:    Input Parameters:
988: +  pc - the multigrid context
989: -  n - number of cycles (default is 1)

991:    Options Database Key:
992: .  -pc_mg_multiplicative_cycles n

994:    Level: advanced

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

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

1000: .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
1001: @*/
1002: PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
1003: {
1004:   PC_MG        *mg        = (PC_MG*)pc->data;

1009:   mg->cyclesperpcapply = n;
1010:   return(0);
1011: }

1013: PetscErrorCode PCMGSetGalerkin_MG(PC pc,PCMGGalerkinType use)
1014: {
1015:   PC_MG *mg = (PC_MG*)pc->data;

1018:   mg->galerkin = use;
1019:   return(0);
1020: }

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

1026:    Logically Collective on PC

1028:    Input Parameters:
1029: +  pc - the multigrid context
1030: -  use - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, or PC_MG_GALERKIN_NONE

1032:    Options Database Key:
1033: .  -pc_mg_galerkin <both,pmat,mat,none>

1035:    Level: intermediate

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

1040: .keywords: MG, set, Galerkin

1042: .seealso: PCMGGetGalerkin(), PCMGGalerkinType

1044: @*/
1045: PetscErrorCode PCMGSetGalerkin(PC pc,PCMGGalerkinType use)
1046: {

1051:   PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PCMGGalerkinType),(pc,use));
1052:   return(0);
1053: }

1055: /*@
1056:    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
1057:       A_i-1 = r_i * A_i * p_i

1059:    Not Collective

1061:    Input Parameter:
1062: .  pc - the multigrid context

1064:    Output Parameter:
1065: .  galerkin - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, PC_MG_GALERKIN_NONE, or PC_MG_GALERKIN_EXTERNAL

1067:    Level: intermediate

1069: .keywords: MG, set, Galerkin

1071: .seealso: PCMGSetGalerkin(), PCMGGalerkinType

1073: @*/
1074: PetscErrorCode  PCMGGetGalerkin(PC pc,PCMGGalerkinType  *galerkin)
1075: {
1076:   PC_MG *mg = (PC_MG*)pc->data;

1080:   *galerkin = mg->galerkin;
1081:   return(0);
1082: }

1084: /*@
1085:    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
1086:    use on all levels. Use PCMGGetSmootherDown() to set different
1087:    pre-smoothing steps on different levels.

1089:    Logically Collective on PC

1091:    Input Parameters:
1092: +  mg - the multigrid context
1093: -  n - the number of smoothing steps

1095:    Options Database Key:
1096: .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps

1098:    Level: advanced
1099:     If the number of smoothing steps is changed in this call then the PCMGGetSmoothUp() will be called and now the
1100:    up smoother will no longer share the same KSP object as the down smoother. Use PCMGSetNumberSmooth() to set the same
1101:    number of smoothing steps for pre and post smoothing.

1103: .keywords: MG, smooth, down, pre-smoothing, steps, multigrid

1105: .seealso: PCMGSetNumberSmoothUp(), PCMGSetNumberSmooth()
1106: @*/
1107: PetscErrorCode  PCMGSetNumberSmoothDown(PC pc,PetscInt n)
1108: {
1109:   PC_MG          *mg        = (PC_MG*)pc->data;
1110:   PC_MG_Levels   **mglevels = mg->levels;
1112:   PetscInt       i,levels;

1116:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1118:   levels = mglevels[0]->levels;

1120:   for (i=1; i<levels; i++) {
1121:     PetscInt nc;
1122:     KSPGetTolerances(mglevels[i]->smoothd,NULL,NULL,NULL,&nc);
1123:     if (nc == n) continue;

1125:     /* make sure smoother up and down are different */
1126:     PCMGGetSmootherUp(pc,i,NULL);
1127:     KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);

1129:     mg->default_smoothd = n;
1130:   }
1131:   return(0);
1132: }

1134: /*@
1135:    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
1136:    on all levels. Use PCMGGetSmootherUp() to set different numbers of
1137:    post-smoothing steps on different levels.

1139:    Logically Collective on PC

1141:    Input Parameters:
1142: +  mg - the multigrid context
1143: -  n - the number of smoothing steps

1145:    Options Database Key:
1146: .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps

1148:    Level: advanced

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

1153:     If the number of smoothing steps is changed in this call then the PCMGGetSmoothUp() will be called and now the 
1154:    up smoother will no longer share the same KSP object as the down smoother. Use PCMGSetNumberSmooth() to set the same
1155:    number of smoothing steps for pre and post smoothing.


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

1160: .seealso: PCMGSetNumberSmoothDown(), PCMGSetNumberSmooth()
1161: @*/
1162: PetscErrorCode  PCMGSetNumberSmoothUp(PC pc,PetscInt n)
1163: {
1164:   PC_MG          *mg        = (PC_MG*)pc->data;
1165:   PC_MG_Levels   **mglevels = mg->levels;
1167:   PetscInt       i,levels;

1171:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1173:   levels = mglevels[0]->levels;

1175:   for (i=1; i<levels; i++) {
1176:     if (mglevels[i]->smoothu == mglevels[i]->smoothd) {
1177:       PetscInt nc;
1178:       KSPGetTolerances(mglevels[i]->smoothd,NULL,NULL,NULL,&nc);
1179:       if (nc == n) continue;
1180:     }

1182:     /* make sure smoother up and down are different */
1183:     PCMGGetSmootherUp(pc,i,NULL);
1184:     KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);

1186:     mg->default_smoothu = n;
1187:   }
1188:   return(0);
1189: }

1191: /*@
1192:    PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use
1193:    on all levels. Use PCMGSetSmoothUp() and PCMGSetSmoothDown() set different numbers of
1194:    pre ad post-smoothing steps

1196:    Logically Collective on PC

1198:    Input Parameters:
1199: +  mg - the multigrid context
1200: -  n - the number of smoothing steps

1202:    Options Database Key:
1203: +  -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps
1204: .  -pc_mg_smooth_down <n> - Sets number of pre-smoothing steps (if setting different pre and post amounts)
1205: -  -pc_mg_smooth_up <n> - Sets number of post-smoothing steps (if setting different pre and post amounts)

1207:    Level: advanced

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

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

1214: .seealso: PCMGSetNumberSmoothDown(), PCMGSetNumberSmoothUp()
1215: @*/
1216: PetscErrorCode  PCMGSetNumberSmooth(PC pc,PetscInt n)
1217: {
1218:   PC_MG          *mg        = (PC_MG*)pc->data;
1219:   PC_MG_Levels   **mglevels = mg->levels;
1221:   PetscInt       i,levels;

1225:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1227:   levels = mglevels[0]->levels;

1229:   for (i=1; i<levels; i++) {
1230:     KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);
1231:     KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);
1232:     mg->default_smoothu = n;
1233:     mg->default_smoothd = n;
1234:   }
1235:   return(0);
1236: }

1238: /* ----------------------------------------------------------------------------------------*/

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

1244:    Options Database Keys:
1245: +  -pc_mg_levels <nlevels> - number of levels including finest
1246: .  -pc_mg_cycle_type <v,w> - 
1247: .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
1248: .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
1249: .  -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
1250: .  -pc_mg_log - log information about time spent on each level of the solver
1251: .  -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
1252: .  -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
1253: .  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1254:                         to the Socket viewer for reading from MATLAB.
1255: -  -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices
1256:                         to the binary output file called binaryoutput

1258:    Notes: 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

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

1262:    Level: intermediate

1264:    Concepts: multigrid/multilevel

1266: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
1267:            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
1268:            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
1269:            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
1270:            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
1271: M*/

1273: PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc)
1274: {
1275:   PC_MG          *mg;

1279:   PetscNewLog(pc,&mg);
1280:   pc->data     = (void*)mg;
1281:   mg->nlevels  = -1;
1282:   mg->am       = PC_MG_MULTIPLICATIVE;
1283:   mg->galerkin = PC_MG_GALERKIN_NONE;

1285:   pc->useAmat = PETSC_TRUE;

1287:   pc->ops->apply          = PCApply_MG;
1288:   pc->ops->setup          = PCSetUp_MG;
1289:   pc->ops->reset          = PCReset_MG;
1290:   pc->ops->destroy        = PCDestroy_MG;
1291:   pc->ops->setfromoptions = PCSetFromOptions_MG;
1292:   pc->ops->view           = PCView_MG;

1294:   PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);
1295:   return(0);
1296: }