Actual source code: mg.c

petsc-master 2017-12-11
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:   if (mg->nlevels == levels) return(0);
190:   PetscObjectGetComm((PetscObject)pc,&comm);
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,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:   levels = PetscMax(mg->nlevels,1);
360:   PetscOptionsHead(PetscOptionsObject,"Multigrid options");
361:   PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);
362:   if (!flg && !mg->levels && pc->dm) {
363:     DMGetRefineLevel(pc->dm,&levels);
364:     levels++;
365:     mg->usedmfornumberoflevels = PETSC_TRUE;
366:   }
367:   PCMGSetLevels(pc,levels,NULL);
368:   mglevels = mg->levels;

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

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

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

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

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

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

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

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

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

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

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


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

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

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

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


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

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

652:     for (i=n-2; i>-1; i--) {DMDestroy(&dms[i]);}
653:     PetscFree(dms);
654:   }

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

662:   if (mg->galerkin < PC_MG_GALERKIN_NONE) {
663:     Mat       A,B;
664:     PetscBool doA = PETSC_FALSE,doB = PETSC_FALSE;
665:     MatReuse  reuse = MAT_INITIAL_MATRIX;

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

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

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

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

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

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

830:   if (mglevels[0]->eventsmoothsetup) {PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);}
831:   KSPSetUp(mglevels[0]->smoothd);
832:   if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PCSETUP_FAILED) {
833:     pc->failedreason = PC_SUBPC_ERROR;
834:   }
835:   if (mglevels[0]->eventsmoothsetup) {PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);}

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

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

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

864: /* -------------------------------------------------------------------------------------*/

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

869:    Not Collective

871:    Input Parameter:
872: .  pc - the preconditioner context

874:    Output parameter:
875: .  levels - the number of levels

877:    Level: advanced

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

881: .seealso: PCMGSetLevels()
882: @*/
883: PetscErrorCode  PCMGGetLevels(PC pc,PetscInt *levels)
884: {
885:   PC_MG *mg = (PC_MG*)pc->data;

890:   *levels = mg->nlevels;
891:   return(0);
892: }

894: /*@
895:    PCMGSetType - Determines the form of multigrid to use:
896:    multiplicative, additive, full, or the Kaskade algorithm.

898:    Logically Collective on PC

900:    Input Parameters:
901: +  pc - the preconditioner context
902: -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
903:    PC_MG_FULL, PC_MG_KASKADE

905:    Options Database Key:
906: .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
907:    additive, full, kaskade

909:    Level: advanced

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

913: .seealso: PCMGSetLevels()
914: @*/
915: PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
916: {
917:   PC_MG *mg = (PC_MG*)pc->data;

922:   mg->am = form;
923:   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
924:   else pc->ops->applyrichardson = NULL;
925:   return(0);
926: }

928: /*@
929:    PCMGGetType - Determines the form of multigrid to use:
930:    multiplicative, additive, full, or the Kaskade algorithm.

932:    Logically Collective on PC

934:    Input Parameter:
935: .  pc - the preconditioner context

937:    Output Parameter:
938: .  type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE


941:    Level: advanced

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

945: .seealso: PCMGSetLevels()
946: @*/
947: PetscErrorCode  PCMGGetType(PC pc,PCMGType *type)
948: {
949:   PC_MG *mg = (PC_MG*)pc->data;

953:   *type = mg->am;
954:   return(0);
955: }

957: /*@
958:    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
959:    complicated cycling.

961:    Logically Collective on PC

963:    Input Parameters:
964: +  pc - the multigrid context
965: -  n - either PC_MG_CYCLE_V or PC_MG_CYCLE_W

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

970:    Level: advanced

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

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

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

991: /*@
992:    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
993:          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used

995:    Logically Collective on PC

997:    Input Parameters:
998: +  pc - the multigrid context
999: -  n - number of cycles (default is 1)

1001:    Options Database Key:
1002: .  -pc_mg_multiplicative_cycles n

1004:    Level: advanced

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

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

1010: .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
1011: @*/
1012: PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
1013: {
1014:   PC_MG        *mg        = (PC_MG*)pc->data;

1019:   mg->cyclesperpcapply = n;
1020:   return(0);
1021: }

1023: PetscErrorCode PCMGSetGalerkin_MG(PC pc,PCMGGalerkinType use)
1024: {
1025:   PC_MG *mg = (PC_MG*)pc->data;

1028:   mg->galerkin = use;
1029:   return(0);
1030: }

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

1036:    Logically Collective on PC

1038:    Input Parameters:
1039: +  pc - the multigrid context
1040: -  use - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, or PC_MG_GALERKIN_NONE

1042:    Options Database Key:
1043: .  -pc_mg_galerkin <both,pmat,mat,none>

1045:    Level: intermediate

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

1050: .keywords: MG, set, Galerkin

1052: .seealso: PCMGGetGalerkin(), PCMGGalerkinType

1054: @*/
1055: PetscErrorCode PCMGSetGalerkin(PC pc,PCMGGalerkinType use)
1056: {

1061:   PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PCMGGalerkinType),(pc,use));
1062:   return(0);
1063: }

1065: /*@
1066:    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
1067:       A_i-1 = r_i * A_i * p_i

1069:    Not Collective

1071:    Input Parameter:
1072: .  pc - the multigrid context

1074:    Output Parameter:
1075: .  galerkin - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, PC_MG_GALERKIN_NONE, or PC_MG_GALERKIN_EXTERNAL

1077:    Level: intermediate

1079: .keywords: MG, set, Galerkin

1081: .seealso: PCMGSetGalerkin(), PCMGGalerkinType

1083: @*/
1084: PetscErrorCode  PCMGGetGalerkin(PC pc,PCMGGalerkinType  *galerkin)
1085: {
1086:   PC_MG *mg = (PC_MG*)pc->data;

1090:   *galerkin = mg->galerkin;
1091:   return(0);
1092: }

1094: /*@
1095:    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
1096:    use on all levels. Use PCMGGetSmootherDown() to set different
1097:    pre-smoothing steps on different levels.

1099:    Logically Collective on PC

1101:    Input Parameters:
1102: +  mg - the multigrid context
1103: -  n - the number of smoothing steps

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

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

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

1115: .seealso: PCMGSetNumberSmoothUp(), PCMGSetNumberSmooth()
1116: @*/
1117: PetscErrorCode  PCMGSetNumberSmoothDown(PC pc,PetscInt n)
1118: {
1119:   PC_MG          *mg        = (PC_MG*)pc->data;
1120:   PC_MG_Levels   **mglevels = mg->levels;
1122:   PetscInt       i,levels;

1127:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1128:   levels = mglevels[0]->levels;
1129:   for (i=1; i<levels; i++) {
1130:     PetscInt nc;
1131:     KSPGetTolerances(mglevels[i]->smoothd,NULL,NULL,NULL,&nc);
1132:     if (nc == n) continue;

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

1138:     mg->default_smoothd = n;
1139:   }
1140:   return(0);
1141: }

1143: /*@
1144:    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
1145:    on all levels. Use PCMGGetSmootherUp() to set different numbers of
1146:    post-smoothing steps on different levels.

1148:    Logically Collective on PC

1150:    Input Parameters:
1151: +  mg - the multigrid context
1152: -  n - the number of smoothing steps

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

1157:    Level: advanced

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

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


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

1169: .seealso: PCMGSetNumberSmoothDown(), PCMGSetNumberSmooth()
1170: @*/
1171: PetscErrorCode  PCMGSetNumberSmoothUp(PC pc,PetscInt n)
1172: {
1173:   PC_MG          *mg        = (PC_MG*)pc->data;
1174:   PC_MG_Levels   **mglevels = mg->levels;
1176:   PetscInt       i,levels;

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

1184:   for (i=1; i<levels; i++) {
1185:     if (mglevels[i]->smoothu == mglevels[i]->smoothd) {
1186:       PetscInt nc;
1187:       KSPGetTolerances(mglevels[i]->smoothd,NULL,NULL,NULL,&nc);
1188:       if (nc == n) continue;
1189:     }

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

1195:     mg->default_smoothu = n;
1196:   }
1197:   return(0);
1198: }

1200: /*@
1201:    PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use
1202:    on all levels. Use PCMGSetSmoothUp() and PCMGSetSmoothDown() set different numbers of
1203:    pre ad post-smoothing steps

1205:    Logically Collective on PC

1207:    Input Parameters:
1208: +  mg - the multigrid context
1209: -  n - the number of smoothing steps

1211:    Options Database Key:
1212: +  -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps
1213: .  -pc_mg_smooth_down <n> - Sets number of pre-smoothing steps (if setting different pre and post amounts)
1214: -  -pc_mg_smooth_up <n> - Sets number of post-smoothing steps (if setting different pre and post amounts)

1216:    Level: advanced

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

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

1223: .seealso: PCMGSetNumberSmoothDown(), PCMGSetNumberSmoothUp()
1224: @*/
1225: PetscErrorCode  PCMGSetNumberSmooth(PC pc,PetscInt n)
1226: {
1227:   PC_MG          *mg        = (PC_MG*)pc->data;
1228:   PC_MG_Levels   **mglevels = mg->levels;
1230:   PetscInt       i,levels;

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

1238:   for (i=1; i<levels; i++) {
1239:     KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);
1240:     KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);
1241:     mg->default_smoothu = n;
1242:     mg->default_smoothd = n;
1243:   }
1244:   return(0);
1245: }

1247: /* ----------------------------------------------------------------------------------------*/

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

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

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

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

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

1276:    Level: intermediate

1278:    Concepts: multigrid/multilevel

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

1287: PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc)
1288: {
1289:   PC_MG          *mg;

1293:   PetscNewLog(pc,&mg);
1294:   pc->data     = (void*)mg;
1295:   mg->nlevels  = -1;
1296:   mg->am       = PC_MG_MULTIPLICATIVE;
1297:   mg->galerkin = PC_MG_GALERKIN_NONE;

1299:   pc->useAmat = PETSC_TRUE;

1301:   pc->ops->apply          = PCApply_MG;
1302:   pc->ops->setup          = PCSetUp_MG;
1303:   pc->ops->reset          = PCReset_MG;
1304:   pc->ops->destroy        = PCDestroy_MG;
1305:   pc->ops->setfromoptions = PCSetFromOptions_MG;
1306:   pc->ops->view           = PCView_MG;

1308:   PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);
1309:   return(0);
1310: }