Actual source code: mg.c

petsc-master 2016-08-26
Report Typos and Errors
  2: /*
  3:     Defines the multigrid preconditioner interface.
  4: */
 5:  #include <petsc/private/pcmgimpl.h>
 6:  #include <petscdm.h>

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

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

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

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

 69: 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)
 70: {
 71:   PC_MG          *mg        = (PC_MG*)pc->data;
 72:   PC_MG_Levels   **mglevels = mg->levels;
 74:   PetscInt       levels = mglevels[0]->levels,i;

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

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

103:   /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't
104:      stop prematurely due to small residual */
105:   for (i=1; i<levels; i++) {
106:     KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
107:     if (mglevels[i]->smoothu != mglevels[i]->smoothd) {
108:       KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
109:     }
110:   }

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

124: PetscErrorCode PCReset_MG(PC pc)
125: {
126:   PC_MG          *mg        = (PC_MG*)pc->data;
127:   PC_MG_Levels   **mglevels = mg->levels;
129:   PetscInt       i,n;

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

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

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

160:    Logically Collective on PC

162:    Input Parameters:
163: +  pc - the preconditioner context
164: .  levels - the number of levels
165: -  comms - optional communicators for each level; this is to allow solving the coarser problems
166:            on smaller sets of processors. Use NULL_OBJECT for default in Fortran

168:    Level: intermediate

170:    Notes:
171:      If the number of levels is one then the multigrid uses the -mg_levels prefix
172:   for setting the level options rather than the -mg_coarse prefix.

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

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

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

212:   mg->nlevels = levels;

214:   PetscMalloc1(levels,&mglevels);
215:   PetscLogObjectMemory((PetscObject)pc,levels*(sizeof(PC_MG*)));

217:   PCGetOptionsPrefix(pc,&prefix);

219:   mg->stageApply = 0;
220:   for (i=0; i<levels; i++) {
221:     PetscNewLog(pc,&mglevels[i]);

223:     mglevels[i]->level               = i;
224:     mglevels[i]->levels              = levels;
225:     mglevels[i]->cycles              = mgctype;
226:     mg->default_smoothu              = 2;
227:     mg->default_smoothd              = 2;
228:     mglevels[i]->eventsmoothsetup    = 0;
229:     mglevels[i]->eventsmoothsolve    = 0;
230:     mglevels[i]->eventresidual       = 0;
231:     mglevels[i]->eventinterprestrict = 0;

233:     if (comms) comm = comms[i];
234:     KSPCreate(comm,&mglevels[i]->smoothd);
235:     KSPSetErrorIfNotConverged(mglevels[i]->smoothd,pc->erroriffailure);
236:     PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);
237:     KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);
238:     PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->smoothd, PetscMGLevelId, mglevels[i]->level);
239:     if (i || levels == 1) {
240:       char tprefix[128];

242:       KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);
243:       KSPSetConvergenceTest(mglevels[i]->smoothd,KSPConvergedSkip,NULL,NULL);
244:       KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);
245:       KSPGetPC(mglevels[i]->smoothd,&ipc);
246:       PCSetType(ipc,PCSOR);
247:       KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd);

249:       sprintf(tprefix,"mg_levels_%d_",(int)i);
250:       KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);
251:     } else {
252:       KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");

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

267:     mglevels[i]->smoothu = mglevels[i]->smoothd;
268:     mg->rtol             = 0.0;
269:     mg->abstol           = 0.0;
270:     mg->dtol             = 0.0;
271:     mg->ttol             = 0.0;
272:     mg->cyclesperpcapply = 1;
273:   }
274:   mg->levels = mglevels;
275:   PCMGSetType(pc,mgtype);
276:   return(0);
277: }


282: PetscErrorCode PCDestroy_MG(PC pc)
283: {
285:   PC_MG          *mg        = (PC_MG*)pc->data;
286:   PC_MG_Levels   **mglevels = mg->levels;
287:   PetscInt       i,n;

290:   PCReset_MG(pc);
291:   if (mglevels) {
292:     n = mglevels[0]->levels;
293:     for (i=0; i<n; i++) {
294:       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
295:         KSPDestroy(&mglevels[i]->smoothd);
296:       }
297:       KSPDestroy(&mglevels[i]->smoothu);
298:       PetscFree(mglevels[i]);
299:     }
300:     PetscFree(mg->levels);
301:   }
302:   PetscFree(pc->data);
303:   return(0);
304: }



308: extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**);
309: extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**);
310: extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**);

312: /*
313:    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
314:              or full cycle of multigrid.

316:   Note:
317:   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
318: */
321: static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
322: {
323:   PC_MG          *mg        = (PC_MG*)pc->data;
324:   PC_MG_Levels   **mglevels = mg->levels;
326:   PetscInt       levels = mglevels[0]->levels,i;

329:   if (mg->stageApply) {PetscLogStagePush(mg->stageApply);}
330:   /* When the DM is supplying the matrix then it will not exist until here */
331:   for (i=0; i<levels; i++) {
332:     if (!mglevels[i]->A) {
333:       KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);
334:       PetscObjectReference((PetscObject)mglevels[i]->A);
335:     }
336:   }

338:   mglevels[levels-1]->b = b;
339:   mglevels[levels-1]->x = x;
340:   if (mg->am == PC_MG_MULTIPLICATIVE) {
341:     VecSet(x,0.0);
342:     for (i=0; i<mg->cyclesperpcapply; i++) {
343:       PCMGMCycle_Private(pc,mglevels+levels-1,NULL);
344:     }
345:   } else if (mg->am == PC_MG_ADDITIVE) {
346:     PCMGACycle_Private(pc,mglevels);
347:   } else if (mg->am == PC_MG_KASKADE) {
348:     PCMGKCycle_Private(pc,mglevels);
349:   } else {
350:     PCMGFCycle_Private(pc,mglevels);
351:   }
352:   if (mg->stageApply) {PetscLogStagePop();}
353:   return(0);
354: }


359: PetscErrorCode PCSetFromOptions_MG(PetscOptionItems *PetscOptionsObject,PC pc)
360: {
361:   PetscErrorCode   ierr;
362:   PetscInt         m,levels = 1,cycles;
363:   PetscBool        flg;
364:   PC_MG            *mg        = (PC_MG*)pc->data;
365:   PC_MG_Levels     **mglevels;
366:   PCMGType         mgtype;
367:   PCMGCycleType    mgctype;
368:   PCMGGalerkinType gtype;

371:   PetscOptionsHead(PetscOptionsObject,"Multigrid options");
372:   if (!mg->levels) {
373:     PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);
374:     if (!flg && pc->dm) {
375:       DMGetRefineLevel(pc->dm,&levels);
376:       levels++;
377:       mg->usedmfornumberoflevels = PETSC_TRUE;
378:     }
379:     PCMGSetLevels(pc,levels,NULL);
380:   }
381:   mglevels = mg->levels;

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

432: #if defined(PETSC_USE_LOG)
433:     {
434:       const char    *sname = "MG Apply";
435:       PetscStageLog stageLog;
436:       PetscInt      st;

439:       PetscLogGetStageLog(&stageLog);
440:       for (st = 0; st < stageLog->numStages; ++st) {
441:         PetscBool same;

443:         PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);
444:         if (same) mg->stageApply = st;
445:       }
446:       if (!mg->stageApply) {
447:         PetscLogStageRegister(sname, &mg->stageApply);
448:       }
449:     }
450: #endif
451:   }
452:   PetscOptionsTail();
453:   return(0);
454: }

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

460:  #include <petscdraw.h>
463: PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
464: {
465:   PC_MG          *mg        = (PC_MG*)pc->data;
466:   PC_MG_Levels   **mglevels = mg->levels;
468:   PetscInt       levels = mglevels ? mglevels[0]->levels : 0,i;
469:   PetscBool      iascii,isbinary,isdraw;

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

548:  #include <petsc/private/dmimpl.h>
549:  #include <petsc/private/kspimpl.h>

551: /*
552:     Calls setup for the KSP on each level
553: */
556: PetscErrorCode PCSetUp_MG(PC pc)
557: {
558:   PC_MG          *mg        = (PC_MG*)pc->data;
559:   PC_MG_Levels   **mglevels = mg->levels;
561:   PetscInt       i,n = mglevels[0]->levels;
562:   PC             cpc;
563:   PetscBool      dump = PETSC_FALSE,opsset,use_amat,missinginterpolate = PETSC_FALSE;
564:   Mat            dA,dB;
565:   Vec            tvec;
566:   DM             *dms;
567:   PetscViewer    viewer = 0;
568:   PetscBool      dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE;

571:   /* FIX: Move this to PCSetFromOptions_MG? */
572:   if (mg->usedmfornumberoflevels) {
573:     PetscInt levels;
574:     DMGetRefineLevel(pc->dm,&levels);
575:     levels++;
576:     if (levels > n) { /* the problem is now being solved on a finer grid */
577:       PCMGSetLevels(pc,levels,NULL);
578:       n        = levels;
579:       PCSetFromOptions(pc); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
580:       mglevels =  mg->levels;
581:     }
582:   }
583:   KSPGetPC(mglevels[0]->smoothd,&cpc);


586:   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
587:   /* so use those from global PC */
588:   /* Is this what we always want? What if user wants to keep old one? */
589:   KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);
590:   if (opsset) {
591:     Mat mmat;
592:     KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);
593:     if (mmat == pc->pmat) opsset = PETSC_FALSE;
594:   }

596:   if (!opsset) {
597:     PCGetUseAmat(pc,&use_amat);
598:     if(use_amat){
599:       PetscInfo(pc,"Using outer operators to define finest grid operator \n  because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n");
600:       KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);
601:     }
602:     else {
603:       PetscInfo(pc,"Using matrix (pmat) operators to define finest grid operator \n  because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n");
604:       KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);
605:     }
606:   }

608:   for (i=n-1; i>0; i--) {
609:     if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) {
610:       missinginterpolate = PETSC_TRUE;
611:       continue;
612:     }
613:   }

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


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

659:     for (i=n-2; i>-1; i--) {DMDestroy(&dms[i]);}
660:     PetscFree(dms);
661:   }

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

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

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

758:   if (!pc->setupcalled) {
759:     for (i=0; i<n; i++) {
760:       KSPSetFromOptions(mglevels[i]->smoothd);
761:     }
762:     for (i=1; i<n; i++) {
763:       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
764:         KSPSetFromOptions(mglevels[i]->smoothu);
765:       }
766:     }
767:     /* insure that if either interpolation or restriction is set the other other one is set */
768:     for (i=1; i<n; i++) {
769:       PCMGGetInterpolation(pc,i,NULL);
770:       PCMGGetRestriction(pc,i,NULL);
771:     }
772:     for (i=0; i<n-1; i++) {
773:       if (!mglevels[i]->b) {
774:         Vec *vec;
775:         KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);
776:         PCMGSetRhs(pc,i,*vec);
777:         VecDestroy(vec);
778:         PetscFree(vec);
779:       }
780:       if (!mglevels[i]->r && i) {
781:         VecDuplicate(mglevels[i]->b,&tvec);
782:         PCMGSetR(pc,i,tvec);
783:         VecDestroy(&tvec);
784:       }
785:       if (!mglevels[i]->x) {
786:         VecDuplicate(mglevels[i]->b,&tvec);
787:         PCMGSetX(pc,i,tvec);
788:         VecDestroy(&tvec);
789:       }
790:     }
791:     if (n != 1 && !mglevels[n-1]->r) {
792:       /* PCMGSetR() on the finest level if user did not supply it */
793:       Vec *vec;
794:       KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);
795:       PCMGSetR(pc,n-1,*vec);
796:       VecDestroy(vec);
797:       PetscFree(vec);
798:     }
799:   }

801:   if (pc->dm) {
802:     /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
803:     for (i=0; i<n-1; i++) {
804:       if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
805:     }
806:   }

808:   for (i=1; i<n; i++) {
809:     if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1){
810:       /* if doing only down then initial guess is zero */
811:       KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);
812:     }
813:     if (mglevels[i]->eventsmoothsetup) {PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);}
814:     KSPSetUp(mglevels[i]->smoothd);
815:     if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PCSETUP_FAILED) {
816:       pc->failedreason = PC_SUBPC_ERROR;
817:     }
818:     if (mglevels[i]->eventsmoothsetup) {PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);}
819:     if (!mglevels[i]->residual) {
820:       Mat mat;
821:       KSPGetOperators(mglevels[i]->smoothd,NULL,&mat);
822:       PCMGSetResidual(pc,i,PCMGResidualDefault,mat);
823:     }
824:   }
825:   for (i=1; i<n; i++) {
826:     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
827:       Mat downmat,downpmat;

829:       /* check if operators have been set for up, if not use down operators to set them */
830:       KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);
831:       if (!opsset) {
832:         KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);
833:         KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);
834:       }

836:       KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);
837:       if (mglevels[i]->eventsmoothsetup) {PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);}
838:       KSPSetUp(mglevels[i]->smoothu);
839:       if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PCSETUP_FAILED) {
840:         pc->failedreason = PC_SUBPC_ERROR;
841:       }
842:       if (mglevels[i]->eventsmoothsetup) {PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);}
843:     }
844:   }

846:   if (mglevels[0]->eventsmoothsetup) {PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);}
847:   KSPSetUp(mglevels[0]->smoothd);
848:   if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PCSETUP_FAILED) {
849:     pc->failedreason = PC_SUBPC_ERROR;
850:   }
851:   if (mglevels[0]->eventsmoothsetup) {PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);}

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

858:    Only support one or the other at the same time.
859:   */
860: #if defined(PETSC_USE_SOCKET_VIEWER)
861:   PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);
862:   if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc));
863:   dump = PETSC_FALSE;
864: #endif
865:   PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);
866:   if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc));

868:   if (viewer) {
869:     for (i=1; i<n; i++) {
870:       MatView(mglevels[i]->restrct,viewer);
871:     }
872:     for (i=0; i<n; i++) {
873:       KSPGetPC(mglevels[i]->smoothd,&pc);
874:       MatView(pc->mat,viewer);
875:     }
876:   }
877:   return(0);
878: }

880: /* -------------------------------------------------------------------------------------*/

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

887:    Not Collective

889:    Input Parameter:
890: .  pc - the preconditioner context

892:    Output parameter:
893: .  levels - the number of levels

895:    Level: advanced

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

899: .seealso: PCMGSetLevels()
900: @*/
901: PetscErrorCode  PCMGGetLevels(PC pc,PetscInt *levels)
902: {
903:   PC_MG *mg = (PC_MG*)pc->data;

908:   *levels = mg->nlevels;
909:   return(0);
910: }

914: /*@
915:    PCMGSetType - Determines the form of multigrid to use:
916:    multiplicative, additive, full, or the Kaskade algorithm.

918:    Logically Collective on PC

920:    Input Parameters:
921: +  pc - the preconditioner context
922: -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
923:    PC_MG_FULL, PC_MG_KASKADE

925:    Options Database Key:
926: .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
927:    additive, full, kaskade

929:    Level: advanced

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

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

942:   mg->am = form;
943:   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
944:   else pc->ops->applyrichardson = NULL;
945:   return(0);
946: }

950: /*@
951:    PCMGGetType - Determines the form of multigrid to use:
952:    multiplicative, additive, full, or the Kaskade algorithm.

954:    Logically Collective on PC

956:    Input Parameter:
957: .  pc - the preconditioner context

959:    Output Parameter:
960: .  type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE


963:    Level: advanced

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

967: .seealso: PCMGSetLevels()
968: @*/
969: PetscErrorCode  PCMGGetType(PC pc,PCMGType *type)
970: {
971:   PC_MG *mg = (PC_MG*)pc->data;

975:   *type = mg->am;
976:   return(0);
977: }

981: /*@
982:    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
983:    complicated cycling.

985:    Logically Collective on PC

987:    Input Parameters:
988: +  pc - the multigrid context
989: -  PC_MG_CYCLE_V or PC_MG_CYCLE_W

991:    Options Database Key:
992: .  -pc_mg_cycle_type <v,w>

994:    Level: advanced

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

998: .seealso: PCMGSetCycleTypeOnLevel()
999: @*/
1000: PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
1001: {
1002:   PC_MG        *mg        = (PC_MG*)pc->data;
1003:   PC_MG_Levels **mglevels = mg->levels;
1004:   PetscInt     i,levels;

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

1012:   for (i=0; i<levels; i++) mglevels[i]->cycles = n;
1013:   return(0);
1014: }

1018: /*@
1019:    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
1020:          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used

1022:    Logically Collective on PC

1024:    Input Parameters:
1025: +  pc - the multigrid context
1026: -  n - number of cycles (default is 1)

1028:    Options Database Key:
1029: .  -pc_mg_multiplicative_cycles n

1031:    Level: advanced

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

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

1037: .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
1038: @*/
1039: PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
1040: {
1041:   PC_MG        *mg        = (PC_MG*)pc->data;

1046:   mg->cyclesperpcapply = n;
1047:   return(0);
1048: }

1052: PetscErrorCode PCMGSetGalerkin_MG(PC pc,PCMGGalerkinType use)
1053: {
1054:   PC_MG *mg = (PC_MG*)pc->data;

1057:   mg->galerkin = use;
1058:   return(0);
1059: }

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

1067:    Logically Collective on PC

1069:    Input Parameters:
1070: +  pc - the multigrid context
1071: -  use - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, or PC_MG_GALERKIN_NONE

1073:    Options Database Key:
1074: .  -pc_mg_galerkin <both,pmat,mat,none>

1076:    Level: intermediate

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

1081: .keywords: MG, set, Galerkin

1083: .seealso: PCMGGetGalerkin(), PCMGGalerkinType

1085: @*/
1086: PetscErrorCode PCMGSetGalerkin(PC pc,PCMGGalerkinType use)
1087: {

1092:   PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PCMGGalerkinType),(pc,use));
1093:   return(0);
1094: }

1098: /*@
1099:    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
1100:       A_i-1 = r_i * A_i * p_i

1102:    Not Collective

1104:    Input Parameter:
1105: .  pc - the multigrid context

1107:    Output Parameter:
1108: .  galerkin - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, PC_MG_GALERKIN_NONE, or PC_MG_GALERKIN_EXTERNAL

1110:    Level: intermediate

1112: .keywords: MG, set, Galerkin

1114: .seealso: PCMGSetGalerkin(), PCMGGalerkinType

1116: @*/
1117: PetscErrorCode  PCMGGetGalerkin(PC pc,PCMGGalerkinType  *galerkin)
1118: {
1119:   PC_MG *mg = (PC_MG*)pc->data;

1123:   *galerkin = mg->galerkin;
1124:   return(0);
1125: }

1129: /*@
1130:    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
1131:    use on all levels. Use PCMGGetSmootherDown() to set different
1132:    pre-smoothing steps on different levels.

1134:    Logically Collective on PC

1136:    Input Parameters:
1137: +  mg - the multigrid context
1138: -  n - the number of smoothing steps

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

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

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

1150: .seealso: PCMGSetNumberSmoothUp(), PCMGSetNumberSmooth()
1151: @*/
1152: PetscErrorCode  PCMGSetNumberSmoothDown(PC pc,PetscInt n)
1153: {
1154:   PC_MG          *mg        = (PC_MG*)pc->data;
1155:   PC_MG_Levels   **mglevels = mg->levels;
1157:   PetscInt       i,levels;

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

1165:   for (i=1; i<levels; i++) {
1166:     PetscInt nc;
1167:     KSPGetTolerances(mglevels[i]->smoothd,NULL,NULL,NULL,&nc);
1168:     if (nc == n) continue;

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

1174:     mg->default_smoothd = n;
1175:   }
1176:   return(0);
1177: }

1181: /*@
1182:    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
1183:    on all levels. Use PCMGGetSmootherUp() to set different numbers of
1184:    post-smoothing steps on different levels.

1186:    Logically Collective on PC

1188:    Input Parameters:
1189: +  mg - the multigrid context
1190: -  n - the number of smoothing steps

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

1195:    Level: advanced

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

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


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

1207: .seealso: PCMGSetNumberSmoothDown(), PCMGSetNumberSmooth()
1208: @*/
1209: PetscErrorCode  PCMGSetNumberSmoothUp(PC pc,PetscInt n)
1210: {
1211:   PC_MG          *mg        = (PC_MG*)pc->data;
1212:   PC_MG_Levels   **mglevels = mg->levels;
1214:   PetscInt       i,levels;

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

1222:   for (i=1; i<levels; i++) {
1223:     if (mglevels[i]->smoothu == mglevels[i]->smoothd) {
1224:       PetscInt nc;
1225:       KSPGetTolerances(mglevels[i]->smoothd,NULL,NULL,NULL,&nc);
1226:       if (nc == n) continue;
1227:     }

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

1233:     mg->default_smoothu = n;
1234:   }
1235:   return(0);
1236: }

1240: /*@
1241:    PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use
1242:    on all levels. Use PCMGSetSmoothUp() and PCMGSetSmoothDown() set different numbers of
1243:    pre ad post-smoothing steps

1245:    Logically Collective on PC

1247:    Input Parameters:
1248: +  mg - the multigrid context
1249: -  n - the number of smoothing steps

1251:    Options Database Key:
1252: +  -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps
1253: .  -pc_mg_smooth_down <n> - Sets number of pre-smoothing steps (if setting different pre and post amounts)
1254: -  -pc_mg_smooth_up <n> - Sets number of post-smoothing steps (if setting different pre and post amounts)

1256:    Level: advanced

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

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

1263: .seealso: PCMGSetNumberSmoothDown(), PCMGSetNumberSmoothUp()
1264: @*/
1265: PetscErrorCode  PCMGSetNumberSmooth(PC pc,PetscInt n)
1266: {
1267:   PC_MG          *mg        = (PC_MG*)pc->data;
1268:   PC_MG_Levels   **mglevels = mg->levels;
1270:   PetscInt       i,levels;

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

1278:   for (i=1; i<levels; i++) {
1279:     KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);
1280:     KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);
1281:     mg->default_smoothu = n;
1282:     mg->default_smoothd = n;
1283:   }
1284:   return(0);
1285: }

1287: /* ----------------------------------------------------------------------------------------*/

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

1293:    Options Database Keys:
1294: +  -pc_mg_levels <nlevels> - number of levels including finest
1295: .  -pc_mg_cycle_type <v,w> - 
1296: .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
1297: .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
1298: .  -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
1299: .  -pc_mg_log - log information about time spent on each level of the solver
1300: .  -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
1301: .  -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
1302: .  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1303:                         to the Socket viewer for reading from MATLAB.
1304: -  -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices
1305:                         to the binary output file called binaryoutput

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

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

1311:    Level: intermediate

1313:    Concepts: multigrid/multilevel

1315: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
1316:            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
1317:            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
1318:            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
1319:            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
1320: M*/

1324: PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc)
1325: {
1326:   PC_MG          *mg;

1330:   PetscNewLog(pc,&mg);
1331:   pc->data     = (void*)mg;
1332:   mg->nlevels  = -1;
1333:   mg->am       = PC_MG_MULTIPLICATIVE;
1334:   mg->galerkin = PC_MG_GALERKIN_NONE;

1336:   pc->useAmat = PETSC_TRUE;

1338:   pc->ops->apply          = PCApply_MG;
1339:   pc->ops->setup          = PCSetUp_MG;
1340:   pc->ops->reset          = PCReset_MG;
1341:   pc->ops->destroy        = PCDestroy_MG;
1342:   pc->ops->setfromoptions = PCSetFromOptions_MG;
1343:   pc->ops->view           = PCView_MG;

1345:   PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);
1346:   return(0);
1347: }