Actual source code: mg.c

  1: #define PETSCKSP_DLL

  3: /*
  4:     Defines the multigrid preconditioner interface.
  5: */
 6:  #include src/ksp/pc/impls/mg/mgimpl.h


 11: PetscErrorCode PCMGMCycle_Private(PC_MG **mglevels,PetscTruth *converged)
 12: {
 13:   PC_MG          *mg = *mglevels,*mgc;
 15:   PetscInt       cycles = (PetscInt) mg->cycles;

 18:   if (converged) *converged = PETSC_FALSE;

 21:   KSPSolve(mg->smoothd,mg->b,mg->x);  /* pre-smooth */
 23:   if (mg->level) {  /* not the coarsest grid */
 25:     (*mg->residual)(mg->A,mg->b,mg->x,mg->r);

 28:     /* if on finest level and have convergence criteria set */
 29:     if (mg->level == mg->levels-1 && mg->ttol) {
 30:       PetscReal rnorm;
 31:       VecNorm(mg->r,NORM_2,&rnorm);
 32:       if (rnorm <= mg->ttol) {
 33:         *converged = PETSC_TRUE;
 34:         if (rnorm < mg->abstol) {
 35:           PetscInfo2(0,"Linear solver has converged. Residual norm %G is less than absolute tolerance %G\n",rnorm,mg->abstol);
 36:         } else {
 37:           PetscInfo2(0,"Linear solver has converged. Residual norm %G is less than relative tolerance times initial residual norm %G\n",rnorm,mg->ttol);
 38:         }
 39:         return(0);
 40:       }
 41:     }

 43:     mgc = *(mglevels - 1);
 45:     MatRestrict(mg->restrct,mg->r,mgc->b);
 47:     VecSet(mgc->x,0.0);
 48:     while (cycles--) {
 49:       PCMGMCycle_Private(mglevels-1,converged);
 50:     }
 52:     MatInterpolateAdd(mg->interpolate,mgc->x,mg->x,mg->x);
 55:     KSPSolve(mg->smoothu,mg->b,mg->x);    /* post smooth */
 57:   }
 58:   return(0);
 59: }

 61: /*
 62:        PCMGCreate_Private - Creates a PC_MG structure for use with the
 63:                multigrid code. Level 0 is the coarsest. (But the 
 64:                finest level is stored first in the array).

 66: */
 69: static PetscErrorCode PCMGCreate_Private(MPI_Comm comm,PetscInt levels,PC pc,MPI_Comm *comms,PC_MG ***result)
 70: {
 71:   PC_MG          **mg;
 73:   PetscInt       i;
 74:   PetscMPIInt    size;
 75:   const char     *prefix;
 76:   PC             ipc;

 79:   PetscMalloc(levels*sizeof(PC_MG*),&mg);
 80:   PetscLogObjectMemory(pc,levels*(sizeof(PC_MG*)+sizeof(PC_MG)));

 82:   PCGetOptionsPrefix(pc,&prefix);

 84:   for (i=0; i<levels; i++) {
 85:     PetscNew(PC_MG,&mg[i]);
 86:     mg[i]->level           = i;
 87:     mg[i]->levels          = levels;
 88:     mg[i]->cycles          = PC_MG_CYCLE_V;
 89:     mg[i]->galerkin        = PETSC_FALSE;
 90:     mg[i]->galerkinused    = PETSC_FALSE;
 91:     mg[i]->default_smoothu = 1;
 92:     mg[i]->default_smoothd = 1;

 94:     if (comms) comm = comms[i];
 95:     KSPCreate(comm,&mg[i]->smoothd);
 96:     KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg[i]->default_smoothd);
 97:     KSPSetOptionsPrefix(mg[i]->smoothd,prefix);

 99:     /* do special stuff for coarse grid */
100:     if (!i && levels > 1) {
101:       KSPAppendOptionsPrefix(mg[0]->smoothd,"mg_coarse_");

103:       /* coarse solve is (redundant) LU by default */
104:       KSPSetType(mg[0]->smoothd,KSPPREONLY);
105:       KSPGetPC(mg[0]->smoothd,&ipc);
106:       MPI_Comm_size(comm,&size);
107:       if (size > 1) {
108:         PCSetType(ipc,PCREDUNDANT);
109:       } else {
110:         PCSetType(ipc,PCLU);
111:       }

113:     } else {
114:       char tprefix[128];
115:       sprintf(tprefix,"mg_levels_%d_",(int)i);
116:       KSPAppendOptionsPrefix(mg[i]->smoothd,tprefix);
117:     }
118:     PetscLogObjectParent(pc,mg[i]->smoothd);
119:     mg[i]->smoothu             = mg[i]->smoothd;
120:     mg[i]->rtol                = 0.0;
121:     mg[i]->abstol              = 0.0;
122:     mg[i]->dtol                = 0.0;
123:     mg[i]->ttol                = 0.0;
124:     mg[i]->eventsmoothsetup    = 0;
125:     mg[i]->eventsmoothsolve    = 0;
126:     mg[i]->eventresidual       = 0;
127:     mg[i]->eventinterprestrict = 0;
128:     mg[i]->cyclesperpcapply    = 1;
129:   }
130:   *result = mg;
131:   return(0);
132: }

136: static PetscErrorCode PCDestroy_MG(PC pc)
137: {
138:   PC_MG          **mg = (PC_MG**)pc->data;
140:   PetscInt       i,n = mg[0]->levels;

143:   for (i=0; i<n-1; i++) {
144:     if (mg[i+1]->r) {VecDestroy(mg[i+1]->r);}
145:     if (mg[i]->b) {VecDestroy(mg[i]->b);}
146:     if (mg[i]->x) {VecDestroy(mg[i]->x);}
147:     if (mg[i+1]->restrct) {MatDestroy(mg[i+1]->restrct);}
148:     if (mg[i+1]->interpolate) {MatDestroy(mg[i+1]->interpolate);}
149:   }

151:   for (i=0; i<n; i++) {
152:     if (mg[i]->smoothd != mg[i]->smoothu) {
153:       KSPDestroy(mg[i]->smoothd);
154:     }
155:     KSPDestroy(mg[i]->smoothu);
156:     PetscFree(mg[i]);
157:   }
158:   PetscFree(mg);
159:   return(0);
160: }



164: EXTERN PetscErrorCode PCMGACycle_Private(PC_MG**);
165: EXTERN PetscErrorCode PCMGFCycle_Private(PC_MG**);
166: EXTERN PetscErrorCode PCMGKCycle_Private(PC_MG**);

168: /*
169:    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
170:              or full cycle of multigrid. 

172:   Note: 
173:   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle(). 
174: */
177: static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
178: {
179:   PC_MG          **mg = (PC_MG**)pc->data;
181:   PetscInt       levels = mg[0]->levels,i;

184:   mg[levels-1]->b = b;
185:   mg[levels-1]->x = x;
186:   if (!mg[levels-1]->r && mg[0]->am != PC_MG_ADDITIVE && levels > 1) {
187:     Vec tvec;
188:     VecDuplicate(mg[levels-1]->b,&tvec);
189:     PCMGSetR(pc,levels-1,tvec);
190:     VecDestroy(tvec);
191:   }
192:   if (mg[0]->am == PC_MG_MULTIPLICATIVE) {
193:     VecSet(x,0.0);
194:     for (i=0; i<mg[0]->cyclesperpcapply; i++) {
195:       PCMGMCycle_Private(mg+levels-1,PETSC_NULL);
196:     }
197:   }
198:   else if (mg[0]->am == PC_MG_ADDITIVE) {
199:     PCMGACycle_Private(mg);
200:   }
201:   else if (mg[0]->am == PC_MG_KASKADE) {
202:     PCMGKCycle_Private(mg);
203:   }
204:   else {
205:     PCMGFCycle_Private(mg);
206:   }
207:   return(0);
208: }

212: static PetscErrorCode PCApplyRichardson_MG(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its)
213: {
214:   PC_MG          **mg = (PC_MG**)pc->data;
216:   PetscInt       levels = mg[0]->levels;
217:   PetscTruth     converged = PETSC_FALSE;

220:   mg[levels-1]->b    = b;
221:   mg[levels-1]->x    = x;

223:   mg[levels-1]->rtol = rtol;
224:   mg[levels-1]->abstol = abstol;
225:   mg[levels-1]->dtol = dtol;
226:   if (rtol) {
227:     /* compute initial residual norm for relative convergence test */
228:     PetscReal rnorm;
229:     (*mg[levels-1]->residual)(mg[levels-1]->A,b,x,w);
230:     VecNorm(w,NORM_2,&rnorm);
231:     mg[levels-1]->ttol = PetscMax(rtol*rnorm,abstol);
232:   } else if (abstol) {
233:     mg[levels-1]->ttol = abstol;
234:   } else {
235:     mg[levels-1]->ttol = 0.0;
236:   }

238:   while (its-- && !converged) {
239:     PCMGMCycle_Private(mg+levels-1,&converged);
240:   }
241:   return(0);
242: }

246: PetscErrorCode PCSetFromOptions_MG(PC pc)
247: {
249:   PetscInt       m,levels = 1,cycles;
250:   PetscTruth     flg;
251:   PC_MG          **mg = (PC_MG**)pc->data;
252:   PCMGType       mgtype = PC_MG_ADDITIVE;
253:   PCMGCycleType  mgctype;

256:   PetscOptionsHead("Multigrid options");
257:     if (!pc->data) {
258:       PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);
259:       PCMGSetLevels(pc,levels,PETSC_NULL);
260:       mg = (PC_MG**)pc->data;
261:     }
262:     mgctype = (PCMGCycleType) mg[0]->cycles;
263:     PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);
264:     if (flg) {
265:       PCMGSetCycleType(pc,mgctype);
266:     };
267:     PetscOptionsName("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",&flg);
268:     if (flg) {
269:       PCMGSetGalerkin(pc);
270:     }
271:     PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);
272:     if (flg) {
273:       PCMGSetNumberSmoothUp(pc,m);
274:     }
275:     PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);
276:     if (flg) {
277:       PCMGSetNumberSmoothDown(pc,m);
278:     }
279:     PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);
280:     if (flg) {
281:       PCMGSetType(pc,mgtype);
282:     }
283:     if (mg[0]->am == PC_MG_MULTIPLICATIVE) {
284:       PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg[0]->cyclesperpcapply,&cycles,&flg);
285:       if (flg) {
286:         PCMGMultiplicativeSetCycles(pc,cycles);
287:       }
288:     }
289:     PetscOptionsName("-pc_mg_log","Log times for each multigrid level","None",&flg);
290:     if (flg) {
291:       PetscInt i;
292:       char     eventname[128];
293:       if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
294:       levels = mg[0]->levels;
295:       for (i=0; i<levels; i++) {
296:         sprintf(eventname,"MGSetup Level %d",(int)i);
298:         sprintf(eventname,"MGSmooth Level %d",(int)i);
300:         if (i) {
301:           sprintf(eventname,"MGResid Level %d",(int)i);
303:           sprintf(eventname,"MGInterp Level %d",(int)i);
305:         }
306:       }
307:     }
308:   PetscOptionsTail();
309:   return(0);
310: }

312: const char *PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
313: const char *PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0};

317: static PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
318: {
319:   PC_MG          **mg = (PC_MG**)pc->data;
321:   PetscInt       levels = mg[0]->levels,i;
322:   PetscTruth     iascii;

325:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
326:   if (iascii) {
327:     PetscViewerASCIIPrintf(viewer,"  MG: type is %s, levels=%D cycles=%s, pre-smooths=%D, post-smooths=%D\n",
328:                                   PCMGTypes[mg[0]->am],levels,(mg[0]->cycles == PC_MG_CYCLE_V) ? "v" : "w",
329:                                   mg[0]->default_smoothd,mg[0]->default_smoothu);
330:     if (mg[0]->galerkin) {
331:       PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");
332:     }
333:     for (i=0; i<levels; i++) {
334:       if (!i) {
335:         PetscViewerASCIIPrintf(viewer,"Coarse gride solver -- level %D -------------------------------\n",i);
336:       } else {
337:         PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);
338:       }
339:       PetscViewerASCIIPushTab(viewer);
340:       KSPView(mg[i]->smoothd,viewer);
341:       PetscViewerASCIIPopTab(viewer);
342:       if (i && mg[i]->smoothd == mg[i]->smoothu) {
343:         PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");
344:       } else if (i){
345:         PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);
346:         PetscViewerASCIIPushTab(viewer);
347:         KSPView(mg[i]->smoothu,viewer);
348:         PetscViewerASCIIPopTab(viewer);
349:       }
350:     }
351:   } else {
352:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name);
353:   }
354:   return(0);
355: }

357: /*
358:     Calls setup for the KSP on each level
359: */
362: static PetscErrorCode PCSetUp_MG(PC pc)
363: {
364:   PC_MG                   **mg = (PC_MG**)pc->data;
365:   PetscErrorCode          ierr;
366:   PetscInt                i,n = mg[0]->levels;
367:   PC                      cpc,mpc;
368:   PetscTruth              preonly,lu,redundant,cholesky,monitor = PETSC_FALSE,dump,opsset;
369:   PetscViewerASCIIMonitor ascii;
370:   PetscViewer             viewer = PETSC_NULL;
371:   MPI_Comm                comm;
372:   Mat                     dA,dB;
373:   MatStructure            uflag;
374:   Vec                     tvec;


378:   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
379:   /* so use those from global PC */
380:   /* Is this what we always want? What if user wants to keep old one? */
381:   KSPGetOperatorsSet(mg[n-1]->smoothd,PETSC_NULL,&opsset);
382:   KSPGetPC(mg[0]->smoothd,&cpc);
383:   KSPGetPC(mg[n-1]->smoothd,&mpc);
384:   if (!opsset || ((cpc->setupcalled == 1) && (mpc->setupcalled == 2))) {
385:     PetscInfo(pc,"Using outer operators to define finest grid operator \n  because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n");
386:     KSPSetOperators(mg[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);
387:   }

389:   if (mg[0]->galerkin) {
390:     Mat B;
391:     mg[0]->galerkinused = PETSC_TRUE;
392:     /* currently only handle case where mat and pmat are the same on coarser levels */
393:     KSPGetOperators(mg[n-1]->smoothd,&dA,&dB,&uflag);
394:     if (!pc->setupcalled) {
395:       for (i=n-2; i>-1; i--) {
396:         MatPtAP(dB,mg[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);
397:         KSPSetOperators(mg[i]->smoothd,B,B,uflag);
398:         if (i != n-2) {PetscObjectDereference((PetscObject)dB);}
399:         dB   = B;
400:       }
401:       PetscObjectDereference((PetscObject)dB);
402:     } else {
403:       for (i=n-2; i>-1; i--) {
404:         KSPGetOperators(mg[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);
405:         MatPtAP(dB,mg[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);
406:         KSPSetOperators(mg[i]->smoothd,B,B,uflag);
407:         dB   = B;
408:       }
409:     }
410:   }

412:   if (!pc->setupcalled) {
413:     PetscOptionsHasName(0,"-pc_mg_monitor",&monitor);
414: 
415:     for (i=0; i<n; i++) {
416:       if (monitor) {
417:         PetscObjectGetComm((PetscObject)mg[i]->smoothd,&comm);
418:         PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);
419:         KSPMonitorSet(mg[i]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);
420:       }
421:       KSPSetFromOptions(mg[i]->smoothd);
422:     }
423:     for (i=1; i<n; i++) {
424:       if (mg[i]->smoothu && (mg[i]->smoothu != mg[i]->smoothd)) {
425:         if (monitor) {
426:           PetscObjectGetComm((PetscObject)mg[i]->smoothu,&comm);
427:           PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);
428:           KSPMonitorSet(mg[i]->smoothu,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);
429:         }
430:         KSPSetFromOptions(mg[i]->smoothu);
431:       }
432:     }
433:     for (i=1; i<n; i++) {
434:       if (!mg[i]->residual) {
435:         Mat mat;
436:         KSPGetOperators(mg[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);
437:         PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);
438:       }
439:       if (mg[i]->restrct && !mg[i]->interpolate) {
440:         PCMGSetInterpolation(pc,i,mg[i]->restrct);
441:       }
442:       if (!mg[i]->restrct && mg[i]->interpolate) {
443:         PCMGSetRestriction(pc,i,mg[i]->interpolate);
444:       }
445: #if defined(PETSC_USE_DEBUG)
446:       if (!mg[i]->restrct || !mg[i]->interpolate) {
447:         SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Need to set restriction or interpolation on level %d",(int)i);
448:       }
449: #endif
450:     }
451:     for (i=0; i<n-1; i++) {
452:       if (!mg[i]->b) {
453:         Vec *vec;
454:         KSPGetVecs(mg[i]->smoothd,1,&vec,0,PETSC_NULL);
455:         PCMGSetRhs(pc,i,*vec);
456:         PetscFree(vec);
457:       }
458:       if (!mg[i]->r && i) {
459:         VecDuplicate(mg[i]->b,&tvec);
460:         PCMGSetR(pc,i,tvec);
461:         VecDestroy(tvec);
462:       }
463:       if (!mg[i]->x) {
464:         VecDuplicate(mg[i]->b,&tvec);
465:         PCMGSetX(pc,i,tvec);
466:         VecDestroy(tvec);
467:       }
468:     }
469:   }


472:   for (i=1; i<n; i++) {
473:     if (mg[i]->smoothu == mg[i]->smoothd) {
474:       /* if doing only down then initial guess is zero */
475:       KSPSetInitialGuessNonzero(mg[i]->smoothd,PETSC_TRUE);
476:     }
478:     KSPSetUp(mg[i]->smoothd);
480:   }
481:   for (i=1; i<n; i++) {
482:     if (mg[i]->smoothu && mg[i]->smoothu != mg[i]->smoothd) {
483:       Mat          downmat,downpmat;
484:       MatStructure matflag;
485:       PetscTruth   opsset;

487:       /* check if operators have been set for up, if not use down operators to set them */
488:       KSPGetOperatorsSet(mg[i]->smoothu,&opsset,PETSC_NULL);
489:       if (!opsset) {
490:         KSPGetOperators(mg[i]->smoothd,&downmat,&downpmat,&matflag);
491:         KSPSetOperators(mg[i]->smoothu,downmat,downpmat,matflag);
492:       }

494:       KSPSetInitialGuessNonzero(mg[i]->smoothu,PETSC_TRUE);
496:       KSPSetUp(mg[i]->smoothu);
498:     }
499:   }

501:   /*
502:       If coarse solver is not direct method then DO NOT USE preonly 
503:   */
504:   PetscTypeCompare((PetscObject)mg[0]->smoothd,KSPPREONLY,&preonly);
505:   if (preonly) {
506:     PetscTypeCompare((PetscObject)cpc,PCLU,&lu);
507:     PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);
508:     PetscTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);
509:     if (!lu && !redundant && !cholesky) {
510:       KSPSetType(mg[0]->smoothd,KSPGMRES);
511:     }
512:   }

514:   if (!pc->setupcalled) {
515:     if (monitor) {
516:       PetscObjectGetComm((PetscObject)mg[0]->smoothd,&comm);
517:       PetscViewerASCIIMonitorCreate(comm,"stdout",n,&ascii);
518:       KSPMonitorSet(mg[0]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);
519:     }
520:     KSPSetFromOptions(mg[0]->smoothd);
521:   }

524:   KSPSetUp(mg[0]->smoothd);

527:   /*
528:      Dump the interpolation/restriction matrices plus the 
529:    Jacobian/stiffness on each level. This allows Matlab users to 
530:    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.

532:    Only support one or the other at the same time.
533:   */
534: #if defined(PETSC_USE_SOCKET_VIEWER)
535:   PetscOptionsHasName(pc->prefix,"-pc_mg_dump_matlab",&dump);
536:   if (dump) {
537:     viewer = PETSC_VIEWER_SOCKET_(pc->comm);
538:   }
539: #endif
540:   PetscOptionsHasName(pc->prefix,"-pc_mg_dump_binary",&dump);
541:   if (dump) {
542:     viewer = PETSC_VIEWER_BINARY_(pc->comm);
543:   }

545:   if (viewer) {
546:     for (i=1; i<n; i++) {
547:       MatView(mg[i]->restrct,viewer);
548:     }
549:     for (i=0; i<n; i++) {
550:       KSPGetPC(mg[i]->smoothd,&pc);
551:       MatView(pc->mat,viewer);
552:     }
553:   }
554:   return(0);
555: }

557: /* -------------------------------------------------------------------------------------*/

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

565:    Collective on PC

567:    Input Parameters:
568: +  pc - the preconditioner context
569: .  levels - the number of levels
570: -  comms - optional communicators for each level; this is to allow solving the coarser problems
571:            on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran

573:    Level: intermediate

575:    Notes:
576:      If the number of levels is one then the multigrid uses the -mg_levels prefix
577:   for setting the level options rather than the -mg_coarse prefix.

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

581: .seealso: PCMGSetType(), PCMGGetLevels()
582: @*/
583: PetscErrorCode  PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
584: {
586:   PC_MG          **mg=0;


591:   if (pc->data) {
592:     SETERRQ(PETSC_ERR_ORDER,"Number levels already set for MG\n\
593:     make sure that you call PCMGSetLevels() before KSPSetFromOptions()");
594:   }
595:   PCMGCreate_Private(pc->comm,levels,pc,comms,&mg);
596:   mg[0]->am                = PC_MG_MULTIPLICATIVE;
597:   pc->data                 = (void*)mg;
598:   pc->ops->applyrichardson = PCApplyRichardson_MG;
599:   return(0);
600: }

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

607:    Not Collective

609:    Input Parameter:
610: .  pc - the preconditioner context

612:    Output parameter:
613: .  levels - the number of levels

615:    Level: advanced

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

619: .seealso: PCMGSetLevels()
620: @*/
621: PetscErrorCode  PCMGGetLevels(PC pc,PetscInt *levels)
622: {
623:   PC_MG  **mg;


629:   mg      = (PC_MG**)pc->data;
630:   *levels = mg[0]->levels;
631:   return(0);
632: }

636: /*@
637:    PCMGSetType - Determines the form of multigrid to use:
638:    multiplicative, additive, full, or the Kaskade algorithm.

640:    Collective on PC

642:    Input Parameters:
643: +  pc - the preconditioner context
644: -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
645:    PC_MG_FULL, PC_MG_KASKADE

647:    Options Database Key:
648: .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
649:    additive, full, kaskade   

651:    Level: advanced

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

655: .seealso: PCMGSetLevels()
656: @*/
657: PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
658: {
659:   PC_MG **mg;

663:   mg = (PC_MG**)pc->data;

665:   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
666:   mg[0]->am = form;
667:   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
668:   else pc->ops->applyrichardson = 0;
669:   return(0);
670: }

674: /*@
675:    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more 
676:    complicated cycling.

678:    Collective on PC

680:    Input Parameters:
681: +  pc - the multigrid context 
682: -  PC_MG_CYCLE_V or PC_MG_CYCLE_W

684:    Options Database Key:
685: $  -pc_mg_cycle_type v or w

687:    Level: advanced

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

691: .seealso: PCMGSetCycleTypeOnLevel()
692: @*/
693: PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
694: {
695:   PC_MG    **mg;
696:   PetscInt i,levels;

700:   mg     = (PC_MG**)pc->data;
701:   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
702:   levels = mg[0]->levels;

704:   for (i=0; i<levels; i++) {
705:     mg[i]->cycles  = n;
706:   }
707:   return(0);
708: }

712: /*@
713:    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 
714:          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used

716:    Collective on PC

718:    Input Parameters:
719: +  pc - the multigrid context 
720: -  n - number of cycles (default is 1)

722:    Options Database Key:
723: $  -pc_mg_multiplicative_cycles n

725:    Level: advanced

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

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

731: .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
732: @*/
733: PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
734: {
735:   PC_MG    **mg;
736:   PetscInt i,levels;

740:   mg     = (PC_MG**)pc->data;
741:   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
742:   levels = mg[0]->levels;

744:   for (i=0; i<levels; i++) {
745:     mg[i]->cyclesperpcapply  = n;
746:   }
747:   return(0);
748: }

752: /*@
753:    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
754:       finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t

756:    Collective on PC

758:    Input Parameters:
759: .  pc - the multigrid context 

761:    Options Database Key:
762: $  -pc_mg_galerkin

764:    Level: intermediate

766: .keywords: MG, set, Galerkin

768: .seealso: PCMGGetGalerkin()

770: @*/
771: PetscErrorCode  PCMGSetGalerkin(PC pc)
772: {
773:   PC_MG    **mg;
774:   PetscInt i,levels;

778:   mg     = (PC_MG**)pc->data;
779:   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
780:   levels = mg[0]->levels;

782:   for (i=0; i<levels; i++) {
783:     mg[i]->galerkin = PETSC_TRUE;
784:   }
785:   return(0);
786: }

790: /*@
791:    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
792:       A_i-1 = r_i * A_i * r_i^t

794:    Not Collective

796:    Input Parameter:
797: .  pc - the multigrid context 

799:    Output Parameter:
800: .  gelerkin - PETSC_TRUE or PETSC_FALSE

802:    Options Database Key:
803: $  -pc_mg_galerkin

805:    Level: intermediate

807: .keywords: MG, set, Galerkin

809: .seealso: PCMGSetGalerkin()

811: @*/
812: PetscErrorCode  PCMGGetGalerkin(PC pc,PetscTruth *galerkin)
813: {
814:   PC_MG    **mg;

818:   mg     = (PC_MG**)pc->data;
819:   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
820:   *galerkin = mg[0]->galerkin;
821:   return(0);
822: }

826: /*@
827:    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
828:    use on all levels. Use PCMGGetSmootherDown() to set different 
829:    pre-smoothing steps on different levels.

831:    Collective on PC

833:    Input Parameters:
834: +  mg - the multigrid context 
835: -  n - the number of smoothing steps

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

840:    Level: advanced

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

844: .seealso: PCMGSetNumberSmoothUp()
845: @*/
846: PetscErrorCode  PCMGSetNumberSmoothDown(PC pc,PetscInt n)
847: {
848:   PC_MG          **mg;
850:   PetscInt       i,levels;

854:   mg     = (PC_MG**)pc->data;
855:   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
856:   levels = mg[0]->levels;

858:   for (i=1; i<levels; i++) {
859:     /* make sure smoother up and down are different */
860:     PCMGGetSmootherUp(pc,i,PETSC_NULL);
861:     KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);
862:     mg[i]->default_smoothd = n;
863:   }
864:   return(0);
865: }

869: /*@
870:    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use 
871:    on all levels. Use PCMGGetSmootherUp() to set different numbers of 
872:    post-smoothing steps on different levels.

874:    Collective on PC

876:    Input Parameters:
877: +  mg - the multigrid context 
878: -  n - the number of smoothing steps

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

883:    Level: advanced

885:    Note: this does not set a value on the coarsest grid, since we assume that
886:     there is no separate smooth up on the coarsest grid.

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

890: .seealso: PCMGSetNumberSmoothDown()
891: @*/
892: PetscErrorCode  PCMGSetNumberSmoothUp(PC pc,PetscInt n)
893: {
894:   PC_MG          **mg;
896:   PetscInt       i,levels;

900:   mg     = (PC_MG**)pc->data;
901:   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
902:   levels = mg[0]->levels;

904:   for (i=1; i<levels; i++) {
905:     /* make sure smoother up and down are different */
906:     PCMGGetSmootherUp(pc,i,PETSC_NULL);
907:     KSPSetTolerances(mg[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);
908:     mg[i]->default_smoothu = n;
909:   }
910:   return(0);
911: }

913: /* ----------------------------------------------------------------------------------------*/

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

919:    Options Database Keys:
920: +  -pc_mg_levels <nlevels> - number of levels including finest
921: .  -pc_mg_cycles v or w
922: .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
923: .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
924: .  -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default
925: .  -pc_mg_log - log information about time spent on each level of the solver
926: .  -pc_mg_monitor - print information on the multigrid convergence
927: .  -pc_mg_galerkin - use Galerkin process to compute coarser operators
928: -  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
929:                         to the Socket viewer for reading from Matlab.

931:    Notes:

933:    Level: intermediate

935:    Concepts: multigrid/multilevel

937: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 
938:            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
939:            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
940:            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
941:            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()           
942: M*/

947: PetscErrorCode  PCCreate_MG(PC pc)
948: {
950:   pc->ops->apply          = PCApply_MG;
951:   pc->ops->setup          = PCSetUp_MG;
952:   pc->ops->destroy        = PCDestroy_MG;
953:   pc->ops->setfromoptions = PCSetFromOptions_MG;
954:   pc->ops->view           = PCView_MG;

956:   pc->data                = (void*)0;
957:   return(0);
958: }