Actual source code: mg.c

petsc-master 2017-09-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:       VecDestroy(&mglevels[i+1]->rscale);
137:     }

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

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

154:    Logically Collective on PC

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

162:    Level: intermediate

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

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

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

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

206:   mg->nlevels = levels;

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

211:   PCGetOptionsPrefix(pc,&prefix);

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

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

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

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

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

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

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


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

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



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

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

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

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

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


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

359:   PetscOptionsHead(PetscOptionsObject,"Multigrid options");
360:   if (!mg->levels) {
361:     PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);
362:     if (!flg && pc->dm) {
363:       DMGetRefineLevel(pc->dm,&levels);
364:       levels++;
365:       mg->usedmfornumberoflevels = PETSC_TRUE;
366:     }
367:     PCMGSetLevels(pc,levels,NULL);
368:   }
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:   PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",mg->default_smoothu,&m,&flg);
382:   if (flg) {
383:     PCMGSetNumberSmoothUp(pc,m);
384:   }
385:   PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",mg->default_smoothd,&m,&flg);
386:   if (flg) {
387:     PCMGSetNumberSmoothDown(pc,m);
388:   }
389:   mgtype = mg->am;
390:   PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);
391:   if (flg) {
392:     PCMGSetType(pc,mgtype);
393:   }
394:   if (mg->am == PC_MG_MULTIPLICATIVE) {
395:     PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGMultiplicativeSetCycles",mg->cyclesperpcapply,&cycles,&flg);
396:     if (flg) {
397:       PCMGMultiplicativeSetCycles(pc,cycles);
398:     }
399:   }
400:   flg  = PETSC_FALSE;
401:   PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);
402:   if (flg) {
403:     PetscInt i;
404:     char     eventname[128];
405:     if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
406:     levels = mglevels[0]->levels;
407:     for (i=0; i<levels; i++) {
408:       sprintf(eventname,"MGSetup Level %d",(int)i);
409:       PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);
410:       sprintf(eventname,"MGSmooth Level %d",(int)i);
411:       PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);
412:       if (i) {
413:         sprintf(eventname,"MGResid Level %d",(int)i);
414:         PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);
415:         sprintf(eventname,"MGInterp Level %d",(int)i);
416:         PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);
417:       }
418:     }

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

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

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

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

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

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

534:  #include <petsc/private/dmimpl.h>
535:  #include <petsc/private/kspimpl.h>

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

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


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

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

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

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


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

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

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

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

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

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

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

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

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

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

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

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

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

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

855: /* -------------------------------------------------------------------------------------*/

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

860:    Not Collective

862:    Input Parameter:
863: .  pc - the preconditioner context

865:    Output parameter:
866: .  levels - the number of levels

868:    Level: advanced

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

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

881:   *levels = mg->nlevels;
882:   return(0);
883: }

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

889:    Logically Collective on PC

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

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

900:    Level: advanced

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

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

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

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

923:    Logically Collective on PC

925:    Input Parameter:
926: .  pc - the preconditioner context

928:    Output Parameter:
929: .  type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE


932:    Level: advanced

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

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

944:   *type = mg->am;
945:   return(0);
946: }

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

952:    Logically Collective on PC

954:    Input Parameters:
955: +  pc - the multigrid context
956: -  n - either PC_MG_CYCLE_V or PC_MG_CYCLE_W

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

961:    Level: advanced

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

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

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

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

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

987:    Logically Collective on PC

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

993:    Options Database Key:
994: .  -pc_mg_multiplicative_cycles n

996:    Level: advanced

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

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

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

1011:   mg->cyclesperpcapply = n;
1012:   return(0);
1013: }

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

1020:   mg->galerkin = use;
1021:   return(0);
1022: }

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

1028:    Logically Collective on PC

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

1034:    Options Database Key:
1035: .  -pc_mg_galerkin <both,pmat,mat,none>

1037:    Level: intermediate

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

1042: .keywords: MG, set, Galerkin

1044: .seealso: PCMGGetGalerkin(), PCMGGalerkinType

1046: @*/
1047: PetscErrorCode PCMGSetGalerkin(PC pc,PCMGGalerkinType use)
1048: {

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

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

1061:    Not Collective

1063:    Input Parameter:
1064: .  pc - the multigrid context

1066:    Output Parameter:
1067: .  galerkin - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, PC_MG_GALERKIN_NONE, or PC_MG_GALERKIN_EXTERNAL

1069:    Level: intermediate

1071: .keywords: MG, set, Galerkin

1073: .seealso: PCMGSetGalerkin(), PCMGGalerkinType

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

1082:   *galerkin = mg->galerkin;
1083:   return(0);
1084: }

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

1091:    Logically Collective on PC

1093:    Input Parameters:
1094: +  mg - the multigrid context
1095: -  n - the number of smoothing steps

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

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

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

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

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

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

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

1131:     mg->default_smoothd = n;
1132:   }
1133:   return(0);
1134: }

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

1141:    Logically Collective on PC

1143:    Input Parameters:
1144: +  mg - the multigrid context
1145: -  n - the number of smoothing steps

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

1150:    Level: advanced

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

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


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

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

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

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

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

1188:     mg->default_smoothu = n;
1189:   }
1190:   return(0);
1191: }

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

1198:    Logically Collective on PC

1200:    Input Parameters:
1201: +  mg - the multigrid context
1202: -  n - the number of smoothing steps

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

1209:    Level: advanced

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

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

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

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

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

1240: /* ----------------------------------------------------------------------------------------*/

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

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

1260:    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

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

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

1269:    Level: intermediate

1271:    Concepts: multigrid/multilevel

1273: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
1274:            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
1275:            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
1276:            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
1277:            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
1278: M*/

1280: PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc)
1281: {
1282:   PC_MG          *mg;

1286:   PetscNewLog(pc,&mg);
1287:   pc->data     = (void*)mg;
1288:   mg->nlevels  = -1;
1289:   mg->am       = PC_MG_MULTIPLICATIVE;
1290:   mg->galerkin = PC_MG_GALERKIN_NONE;

1292:   pc->useAmat = PETSC_TRUE;

1294:   pc->ops->apply          = PCApply_MG;
1295:   pc->ops->setup          = PCSetUp_MG;
1296:   pc->ops->reset          = PCReset_MG;
1297:   pc->ops->destroy        = PCDestroy_MG;
1298:   pc->ops->setfromoptions = PCSetFromOptions_MG;
1299:   pc->ops->view           = PCView_MG;

1301:   PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);
1302:   return(0);
1303: }