Actual source code: mg.c

  1: /*$Id: mg.c,v 1.128 2001/10/03 22:12:40 balay Exp $*/
  2: /*
  3:     Defines the multigrid preconditioner interface.
  4: */
 5:  #include src/ksp/pc/impls/mg/mgimpl.h


  8: /*
  9:        MGMCycle_Private - Given an MG structure created with MGCreate() runs 
 10:                   one multiplicative cycle down through the levels and
 11:                   back up.

 13:     Input Parameter:
 14: .   mg - structure created with  MGCreate().
 15: */
 18: int MGMCycle_Private(MG *mglevels,PetscTruth *converged)
 19: {
 20:   MG          mg = *mglevels,mgc;
 21:   int         cycles = mg->cycles,ierr;
 22:   PetscScalar zero = 0.0;

 25:   if (converged) *converged = PETSC_FALSE;

 27:   if (mg->eventsolve) {PetscLogEventBegin(mg->eventsolve,0,0,0,0);}
 28:   KSPSetRhs(mg->smoothd,mg->b);
 29:   KSPSetSolution(mg->smoothd,mg->x);
 30:   KSPSolve(mg->smoothd);
 31:   if (mg->eventsolve) {PetscLogEventEnd(mg->eventsolve,0,0,0,0);}
 32:   if (mg->level) {  /* not the coarsest grid */
 33:     (*mg->residual)(mg->A,mg->b,mg->x,mg->r);

 35:     /* if on finest level and have convergence criteria set */
 36:     if (mg->level == mg->levels-1 && mg->ttol) {
 37:       PetscReal rnorm;
 38:       VecNorm(mg->r,NORM_2,&rnorm);
 39:       if (rnorm <= mg->ttol) {
 40:         *converged = PETSC_TRUE;
 41:         if (rnorm < mg->atol) {
 42:           PetscLogInfo(0,"Linear solver has converged. Residual norm %g is less than absolute tolerance %g\n",rnorm,mg->atol);
 43:         } else {
 44:           PetscLogInfo(0,"Linear solver has converged. Residual norm %g is less than relative tolerance times initial residual norm %g\n",rnorm,mg->ttol);
 45:         }
 46:         return(0);
 47:       }
 48:     }

 50:     mgc = *(mglevels - 1);
 51:     MatRestrict(mg->restrct,mg->r,mgc->b);
 52:     VecSet(&zero,mgc->x);
 53:     while (cycles--) {
 54:       MGMCycle_Private(mglevels-1,converged);
 55:     }
 56:     MatInterpolateAdd(mg->interpolate,mgc->x,mg->x,mg->x);
 57:     if (mg->eventsolve) {PetscLogEventBegin(mg->eventsolve,0,0,0,0);}
 58:     KSPSetRhs(mg->smoothu,mg->b);
 59:     KSPSetSolution(mg->smoothu,mg->x);
 60:     KSPSolve(mg->smoothu);
 61:     if (mg->eventsolve) {PetscLogEventEnd(mg->eventsolve,0,0,0,0);}
 62:   }
 63:   return(0);
 64: }

 66: /*
 67:        MGCreate_Private - Creates a MG structure for use with the
 68:                multigrid code. Level 0 is the coarsest. (But the 
 69:                finest level is stored first in the array).

 71: */
 74: static int MGCreate_Private(MPI_Comm comm,int levels,PC pc,MPI_Comm *comms,MG **result)
 75: {
 76:   MG   *mg;
 77:   int  i,ierr,size;
 78:   char *prefix;
 79:   PC   ipc;

 82:   PetscMalloc(levels*sizeof(MG),&mg);
 83:   PetscLogObjectMemory(pc,levels*(sizeof(MG)+sizeof(struct _MG)));

 85:   PCGetOptionsPrefix(pc,&prefix);

 87:   for (i=0; i<levels; i++) {
 88:     PetscNew(struct _MG,&mg[i]);
 89:     PetscMemzero(mg[i],sizeof(struct _MG));
 90:     mg[i]->level  = i;
 91:     mg[i]->levels = levels;
 92:     mg[i]->cycles = 1;

 94:     if (comms) comm = comms[i];
 95:     KSPCreate(comm,&mg[i]->smoothd);
 96:     KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);
 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:         PCRedundantGetPC(ipc,&ipc);
110:       }
111:       PCSetType(ipc,PCLU);

113:     } else {
114:       KSPAppendOptionsPrefix(mg[i]->smoothd,"mg_levels_");
115:     }
116:     PetscLogObjectParent(pc,mg[i]->smoothd);
117:     mg[i]->smoothu         = mg[i]->smoothd;
118:     mg[i]->default_smoothu = 10000;
119:     mg[i]->default_smoothd = 10000;
120:     mg[i]->rtol = 0.0;
121:     mg[i]->atol = 0.0;
122:     mg[i]->dtol = 0.0;
123:     mg[i]->ttol = 0.0;
124:     mg[i]->eventsetup = 0;
125:     mg[i]->eventsolve = 0;
126:   }
127:   *result = mg;
128:   return(0);
129: }

133: static int PCDestroy_MG(PC pc)
134: {
135:   MG  *mg = (MG*)pc->data;
136:   int i,n = mg[0]->levels,ierr;

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



152: EXTERN int MGACycle_Private(MG*);
153: EXTERN int MGFCycle_Private(MG*);
154: EXTERN int MGKCycle_Private(MG*);

156: /*
157:    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
158:              or full cycle of multigrid. 

160:   Note: 
161:   A simple wrapper which calls MGMCycle(),MGACycle(), or MGFCycle(). 
162: */
165: static int PCApply_MG(PC pc,Vec b,Vec x)
166: {
167:   MG          *mg = (MG*)pc->data;
168:   PetscScalar zero = 0.0;
169:   int         levels = mg[0]->levels,ierr;

172:   mg[levels-1]->b = b;
173:   mg[levels-1]->x = x;
174:   if (mg[0]->am == MGMULTIPLICATIVE) {
175:     VecSet(&zero,x);
176:     MGMCycle_Private(mg+levels-1,PETSC_NULL);
177:   }
178:   else if (mg[0]->am == MGADDITIVE) {
179:     MGACycle_Private(mg);
180:   }
181:   else if (mg[0]->am == MGKASKADE) {
182:     MGKCycle_Private(mg);
183:   }
184:   else {
185:     MGFCycle_Private(mg);
186:   }
187:   return(0);
188: }

192: static int PCApplyRichardson_MG(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal atol, PetscReal dtol,int its)
193: {
194:   MG         *mg = (MG*)pc->data;
195:   int        ierr,levels = mg[0]->levels;
196:   PetscTruth converged = PETSC_FALSE;

199:   mg[levels-1]->b    = b;
200:   mg[levels-1]->x    = x;

202:   mg[levels-1]->rtol = rtol;
203:   mg[levels-1]->atol = atol;
204:   mg[levels-1]->dtol = dtol;
205:   if (rtol) {
206:     /* compute initial residual norm for relative convergence test */
207:     PetscReal rnorm;
208:     (*mg[levels-1]->residual)(mg[levels-1]->A,b,x,w);
209:     VecNorm(w,NORM_2,&rnorm);
210:     mg[levels-1]->ttol = PetscMax(rtol*rnorm,atol);
211:   } else if (atol) {
212:     mg[levels-1]->ttol = atol;
213:   } else {
214:     mg[levels-1]->ttol = 0.0;
215:   }

217:   while (its-- && !converged) {
218:     MGMCycle_Private(mg+levels-1,&converged);
219:   }
220:   return(0);
221: }

225: static int PCSetFromOptions_MG(PC pc)
226: {
227:   int        ierr,indx,m,levels = 1;
228:   PetscTruth flg;
229:   const char *type[] = {"additive","multiplicative","full","cascade","kascade"};


233:   PetscOptionsHead("Multigrid options");
234:     if (!pc->data) {
235:       PetscOptionsInt("-pc_mg_levels","Number of Levels","MGSetLevels",levels,&levels,&flg);
236:       MGSetLevels(pc,levels,PETSC_NULL);
237:     }
238:     PetscOptionsInt("-pc_mg_cycles","1 for V cycle, 2 for W-cycle","MGSetCycles",1,&m,&flg);
239:     if (flg) {
240:       MGSetCycles(pc,m);
241:     }
242:     PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","MGSetNumberSmoothUp",1,&m,&flg);
243:     if (flg) {
244:       MGSetNumberSmoothUp(pc,m);
245:     }
246:     PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","MGSetNumberSmoothDown",1,&m,&flg);
247:     if (flg) {
248:       MGSetNumberSmoothDown(pc,m);
249:     }
250:     PetscOptionsEList("-pc_mg_type","Multigrid type","MGSetType",type,5,type[1],&indx,&flg);
251:     if (flg) {
252:       MGType mg = (MGType) 0;
253:       switch (indx) {
254:       case 0:
255:         mg = MGADDITIVE;
256:         break;
257:       case 1:
258:         mg = MGMULTIPLICATIVE;
259:         break;
260:       case 2:
261:         mg = MGFULL;
262:         break;
263:       case 3:
264:         mg = MGKASKADE;
265:         break;
266:       case 4:
267:         mg = MGKASKADE;
268:         break;
269:       }
270:       MGSetType(pc,mg);
271:     }
272:     PetscOptionsName("-pc_mg_log","Log times for each multigrid level","None",&flg);
273:     if (flg) {
274:       MG   *mg = (MG*)pc->data;
275:       int  i;
276:       char eventname[128];
277:       if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
278:       levels = mg[0]->levels;
279:       for (i=0; i<levels; i++) {
280:         sprintf(eventname,"MSetup Level %d",i);
281:         PetscLogEventRegister(&mg[i]->eventsetup,eventname,pc->cookie);
282:         sprintf(eventname,"MGSolve Level %d",i);
283:         PetscLogEventRegister(&mg[i]->eventsolve,eventname,pc->cookie);
284:       }
285:     }
286:   PetscOptionsTail();
287:   return(0);
288: }

292: static int PCView_MG(PC pc,PetscViewer viewer)
293: {
294:   MG         *mg = (MG*)pc->data;
295:   int        ierr,levels = mg[0]->levels,i;
296:   const char *cstring;
297:   PetscTruth isascii;

300:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
301:   if (isascii) {
302:     if (mg[0]->am == MGMULTIPLICATIVE) cstring = "multiplicative";
303:     else if (mg[0]->am == MGADDITIVE)  cstring = "additive";
304:     else if (mg[0]->am == MGFULL)      cstring = "full";
305:     else if (mg[0]->am == MGKASKADE)   cstring = "Kaskade";
306:     else cstring = "unknown";
307:     PetscViewerASCIIPrintf(viewer,"  MG: type is %s, levels=%d cycles=%d, pre-smooths=%d, post-smooths=%d\n",
308:                       cstring,levels,mg[0]->cycles,mg[0]->default_smoothd,mg[0]->default_smoothu);
309:     for (i=0; i<levels; i++) {
310:       PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %d -------------------------------\n",i);
311:       PetscViewerASCIIPushTab(viewer);
312:       KSPView(mg[i]->smoothd,viewer);
313:       PetscViewerASCIIPopTab(viewer);
314:       if (mg[i]->smoothd == mg[i]->smoothu) {
315:         PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");
316:       } else {
317:         PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %d -------------------------------\n",i);
318:         PetscViewerASCIIPushTab(viewer);
319:         KSPView(mg[i]->smoothu,viewer);
320:         PetscViewerASCIIPopTab(viewer);
321:       }
322:     }
323:   } else {
324:     SETERRQ1(1,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name);
325:   }
326:   return(0);
327: }

329: /*
330:     Calls setup for the KSP on each level
331: */
334: static int PCSetUp_MG(PC pc)
335: {
336:   MG          *mg = (MG*)pc->data;
337:   int         ierr,i,n = mg[0]->levels;
338:   PC          cpc;
339:   PetscTruth  preonly,lu,redundant,monitor = PETSC_FALSE,dump;
340:   PetscViewer ascii;
341:   MPI_Comm    comm;

344:   /*
345:      temporarily stick pc->vec into mg[0]->b and x so that 
346:    KSPSetUp is happy. Since currently those slots are empty.
347:   */
348:   mg[n-1]->x = pc->vec;
349:   mg[n-1]->b = pc->vec;

351:   if (pc->setupcalled == 0) {
352:     PetscOptionsHasName(0,"-pc_mg_monitor",&monitor);
353: 
354:     for (i=1; i<n; i++) {
355:       if (mg[i]->smoothd) {
356:         if (monitor) {
357:           PetscObjectGetComm((PetscObject)mg[i]->smoothd,&comm);
358:           PetscViewerASCIIOpen(comm,"stdout",&ascii);
359:           PetscViewerASCIISetTab(ascii,n-i);
360:           KSPSetMonitor(mg[i]->smoothd,KSPDefaultMonitor,ascii,(int(*)(void*))PetscViewerDestroy);
361:         }
362:         KSPSetFromOptions(mg[i]->smoothd);
363:       }
364:     }
365:     for (i=0; i<n; i++) {
366:       if (mg[i]->smoothu && mg[i]->smoothu != mg[i]->smoothd) {
367:         if (monitor) {
368:           PetscObjectGetComm((PetscObject)mg[i]->smoothu,&comm);
369:           PetscViewerASCIIOpen(comm,"stdout",&ascii);
370:           PetscViewerASCIISetTab(ascii,n-i);
371:           KSPSetMonitor(mg[i]->smoothu,KSPDefaultMonitor,ascii,(int(*)(void*))PetscViewerDestroy);
372:         }
373:         KSPSetFromOptions(mg[i]->smoothu);
374:       }
375:     }
376:   }

378:   for (i=1; i<n; i++) {
379:     if (mg[i]->smoothd) {
380:       KSPSetInitialGuessNonzero(mg[i]->smoothd,PETSC_TRUE);
381:       if (mg[i]->eventsetup) {PetscLogEventBegin(mg[i]->eventsetup,0,0,0,0);}
382:       KSPSetRhs(mg[i]->smoothd,mg[i]->b);
383:       KSPSetSolution(mg[i]->smoothd,mg[i]->x);
384:       KSPSetUp(mg[i]->smoothd);
385:       if (mg[i]->eventsetup) {PetscLogEventEnd(mg[i]->eventsetup,0,0,0,0);}
386:     }
387:   }
388:   for (i=0; i<n; i++) {
389:     if (mg[i]->smoothu && mg[i]->smoothu != mg[i]->smoothd) {
390:       KSPSetInitialGuessNonzero(mg[i]->smoothu,PETSC_TRUE);
391:       if (mg[i]->eventsetup) {PetscLogEventBegin(mg[i]->eventsetup,0,0,0,0);}
392:       KSPSetRhs(mg[i]->smoothu,mg[i]->b);
393:       KSPSetSolution(mg[i]->smoothu,mg[i]->x);
394:       KSPSetUp(mg[i]->smoothu);
395:       if (mg[i]->eventsetup) {PetscLogEventEnd(mg[i]->eventsetup,0,0,0,0);}
396:     }
397:   }

399:   /*
400:       If coarse solver is not direct method then DO NOT USE preonly 
401:   */
402:   PetscTypeCompare((PetscObject)mg[0]->smoothd,KSPPREONLY,&preonly);
403:   if (preonly) {
404:     KSPGetPC(mg[0]->smoothd,&cpc);
405:     PetscTypeCompare((PetscObject)cpc,PCLU,&lu);
406:     PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);
407:     if (!lu && !redundant) {
408:       KSPSetType(mg[0]->smoothd,KSPGMRES);
409:     }
410:   }

412:   if (pc->setupcalled == 0) {
413:     if (monitor) {
414:       PetscObjectGetComm((PetscObject)mg[0]->smoothd,&comm);
415:       PetscViewerASCIIOpen(comm,"stdout",&ascii);
416:       PetscViewerASCIISetTab(ascii,n);
417:       KSPSetMonitor(mg[0]->smoothd,KSPDefaultMonitor,ascii,(int(*)(void*))PetscViewerDestroy);
418:     }
419:     KSPSetFromOptions(mg[0]->smoothd);
420:   }

422:   if (mg[0]->eventsetup) {PetscLogEventBegin(mg[0]->eventsetup,0,0,0,0);}
423:   KSPSetRhs(mg[0]->smoothd,mg[0]->b);
424:   KSPSetSolution(mg[0]->smoothd,mg[0]->x);
425:   KSPSetUp(mg[0]->smoothd);
426:   if (mg[0]->eventsetup) {PetscLogEventEnd(mg[0]->eventsetup,0,0,0,0);}

428:   /*
429:      Dump the interpolation/restriction matrices to matlab plus the 
430:    Jacobian/stiffness on each level. This allows Matlab users to 
431:    easily check if the Galerkin condition A_c = R A_f R^T is satisfied */
432:   PetscOptionsHasName(pc->prefix,"-pc_mg_dump_matlab",&dump);
433:   if (dump) {
434:     for (i=1; i<n; i++) {
435:       MatView(mg[i]->restrct,PETSC_VIEWER_SOCKET_(pc->comm));
436:     }
437:     for (i=0; i<n; i++) {
438:       KSPGetPC(mg[i]->smoothd,&pc);
439:       MatView(pc->mat,PETSC_VIEWER_SOCKET_(pc->comm));
440:     }
441:   }

443:   return(0);
444: }

446: /* -------------------------------------------------------------------------------------*/

450: /*@C
451:    MGSetLevels - Sets the number of levels to use with MG.
452:    Must be called before any other MG routine.

454:    Collective on PC

456:    Input Parameters:
457: +  pc - the preconditioner context
458: .  levels - the number of levels
459: -  comms - optional communicators for each level; this is to allow solving the coarser problems
460:            on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran

462:    Level: intermediate

464:    Notes:
465:      If the number of levels is one then the multigrid uses the -mg_levels prefix
466:   for setting the level options rather than the -mg_coarse prefix.

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

470: .seealso: MGSetType(), MGGetLevels()
471: @*/
472: int MGSetLevels(PC pc,int levels,MPI_Comm *comms)
473: {
475:   MG  *mg;


480:   if (pc->data) {
481:     SETERRQ(1,"Number levels already set for MG\n\
482:     make sure that you call MGSetLevels() before KSPSetFromOptions()");
483:   }
484:   MGCreate_Private(pc->comm,levels,pc,comms,&mg);
485:   mg[0]->am                = MGMULTIPLICATIVE;
486:   pc->data                 = (void*)mg;
487:   pc->ops->applyrichardson = PCApplyRichardson_MG;
488:   return(0);
489: }

493: /*@
494:    MGGetLevels - Gets the number of levels to use with MG.

496:    Not Collective

498:    Input Parameter:
499: .  pc - the preconditioner context

501:    Output parameter:
502: .  levels - the number of levels

504:    Level: advanced

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

508: .seealso: MGSetLevels()
509: @*/
510: int MGGetLevels(PC pc,int *levels)
511: {
512:   MG  *mg;


518:   mg      = (MG*)pc->data;
519:   *levels = mg[0]->levels;
520:   return(0);
521: }

525: /*@
526:    MGSetType - Determines the form of multigrid to use:
527:    multiplicative, additive, full, or the Kaskade algorithm.

529:    Collective on PC

531:    Input Parameters:
532: +  pc - the preconditioner context
533: -  form - multigrid form, one of MGMULTIPLICATIVE, MGADDITIVE,
534:    MGFULL, MGKASKADE

536:    Options Database Key:
537: .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
538:    additive, full, kaskade   

540:    Level: advanced

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

544: .seealso: MGSetLevels()
545: @*/
546: int MGSetType(PC pc,MGType form)
547: {
548:   MG *mg;

552:   mg = (MG*)pc->data;

554:   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
555:   mg[0]->am = form;
556:   if (form == MGMULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
557:   else pc->ops->applyrichardson = 0;
558:   return(0);
559: }

563: /*@
564:    MGSetCycles - Sets the type cycles to use.  Use MGSetCyclesOnLevel() for more 
565:    complicated cycling.

567:    Collective on PC

569:    Input Parameters:
570: +  mg - the multigrid context 
571: -  n - the number of cycles

573:    Options Database Key:
574: $  -pc_mg_cycles n - 1 denotes a V-cycle; 2 denotes a W-cycle.

576:    Level: advanced

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

580: .seealso: MGSetCyclesOnLevel()
581: @*/
582: int MGSetCycles(PC pc,int n)
583: {
584:   MG  *mg;
585:   int i,levels;

589:   mg     = (MG*)pc->data;
590:   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
591:   levels = mg[0]->levels;

593:   for (i=0; i<levels; i++) {
594:     mg[i]->cycles  = n;
595:   }
596:   return(0);
597: }

601: /*@
602:    MGCheck - Checks that all components of the MG structure have 
603:    been set.

605:    Collective on PC

607:    Input Parameters:
608: .  mg - the MG structure

610:    Level: advanced

612: .keywords: MG, check, set, multigrid
613: @*/
614: int MGCheck(PC pc)
615: {
616:   MG  *mg;
617:   int i,n,count = 0;

621:   mg = (MG*)pc->data;

623:   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");

625:   n = mg[0]->levels;

627:   for (i=1; i<n; i++) {
628:     if (!mg[i]->restrct) {
629:       (*PetscErrorPrintf)("No restrict set level %d \n",n-i); count++;
630:     }
631:     if (!mg[i]->interpolate) {
632:       (*PetscErrorPrintf)("No interpolate set level %d \n",n-i); count++;
633:     }
634:     if (!mg[i]->residual) {
635:       (*PetscErrorPrintf)("No residual set level %d \n",n-i); count++;
636:     }
637:     if (!mg[i]->smoothu) {
638:       (*PetscErrorPrintf)("No smoothup set level %d \n",n-i); count++;
639:     }
640:     if (!mg[i]->smoothd) {
641:       (*PetscErrorPrintf)("No smoothdown set level %d \n",n-i); count++;
642:     }
643:     if (!mg[i]->r) {
644:       (*PetscErrorPrintf)("No r set level %d \n",n-i); count++;
645:     }
646:     if (!mg[i-1]->x) {
647:       (*PetscErrorPrintf)("No x set level %d \n",n-i); count++;
648:     }
649:     if (!mg[i-1]->b) {
650:       (*PetscErrorPrintf)("No b set level %d \n",n-i); count++;
651:     }
652:   }
653:   PetscFunctionReturn(count);
654: }


659: /*@
660:    MGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
661:    use on all levels. Use MGGetSmootherDown() to set different 
662:    pre-smoothing steps on different levels.

664:    Collective on PC

666:    Input Parameters:
667: +  mg - the multigrid context 
668: -  n - the number of smoothing steps

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

673:    Level: advanced

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

677: .seealso: MGSetNumberSmoothUp()
678: @*/
679: int MGSetNumberSmoothDown(PC pc,int n)
680: {
681:   MG  *mg;
682:   int i,levels,ierr;

686:   mg     = (MG*)pc->data;
687:   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
688:   levels = mg[0]->levels;

690:   for (i=0; i<levels; i++) {
691:     /* make sure smoother up and down are different */
692:     MGGetSmootherUp(pc,i,PETSC_NULL);
693:     KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);
694:     mg[i]->default_smoothd = n;
695:   }
696:   return(0);
697: }

701: /*@
702:    MGSetNumberSmoothUp - Sets the number of post-smoothing steps to use 
703:    on all levels. Use MGGetSmootherUp() to set different numbers of 
704:    post-smoothing steps on different levels.

706:    Collective on PC

708:    Input Parameters:
709: +  mg - the multigrid context 
710: -  n - the number of smoothing steps

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

715:    Level: advanced

717:    Note: this does not set a value on the coarsest grid, since we assume that
718:     there is no seperate smooth up on the coarsest grid. If you want to have a
719:     seperate smooth up on the coarsest grid then call MGGetSmoothUp(pc,0,&ksp);

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

723: .seealso: MGSetNumberSmoothDown()
724: @*/
725: int  MGSetNumberSmoothUp(PC pc,int n)
726: {
727:   MG  *mg;
728:   int i,levels,ierr;

732:   mg     = (MG*)pc->data;
733:   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
734:   levels = mg[0]->levels;

736:   for (i=1; i<levels; i++) {
737:     /* make sure smoother up and down are different */
738:     MGGetSmootherUp(pc,i,PETSC_NULL);
739:     KSPSetTolerances(mg[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);
740:     mg[i]->default_smoothu = n;
741:   }
742:   return(0);
743: }

745: /* ----------------------------------------------------------------------------------------*/

747: /*MC
748:    PCMG - Use geometric multigrid preconditioning. This preconditioner requires you provide additional
749:     information about the coarser grid matrices and restriction/interpolation operators.

751:    Options Database Keys:
752: +  -pc_mg_levels <nlevels> - number of levels including finest
753: .  -pc_mg_cycles 1 or 2 - for V or W-cycle
754: .  -pc_mg_smoothup <n> - number of smoothing steps before interpolation
755: .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
756: .  -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default
757: .  -pc_mg_log - log information about time spent on each level of the solver
758: .  -pc_mg_monitor - print information on the multigrid convergence
759: -  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
760:                         to the Socket viewer for reading from Matlab.

762:    Notes:

764:    Level: intermediate

766:    Concepts: multigrid

768: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 
769:            MGSetLevels(), MGGetLevels(), MGSetType(), MPSetCycles(), MGSetNumberSmoothDown(),
770:            MGSetNumberSmoothUp(), MGGetCoarseSolve(), MGSetResidual(), MGSetInterpolation(),
771:            MGSetRestriction(), MGGetSmoother(), MGGetSmootherUp(), MGGetSmootherDown(),
772:            MGSetCyclesOnLevel(), MGSetRhs(), MGSetX(), MGSetR()           
773: M*/

775: EXTERN_C_BEGIN
778: int PCCreate_MG(PC pc)
779: {
781:   pc->ops->apply          = PCApply_MG;
782:   pc->ops->setup          = PCSetUp_MG;
783:   pc->ops->destroy        = PCDestroy_MG;
784:   pc->ops->setfromoptions = PCSetFromOptions_MG;
785:   pc->ops->view           = PCView_MG;

787:   pc->data                = (void*)0;
788:   return(0);
789: }
790: EXTERN_C_END