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/sles/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: */
 16: #undef __FUNCT__  
 18: int MGMCycle_Private(MG *mglevels,PetscTruth *converged)
 19: {
 20:   MG          mg = *mglevels,mgc;
 21:   int         cycles = mg->cycles,ierr,its;
 22:   PetscScalar zero = 0.0;

 25:   if (converged) *converged = PETSC_FALSE;

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

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

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

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

 67: */
 68: #undef __FUNCT__  
 70: static int MGCreate_Private(MPI_Comm comm,int levels,PC pc,MPI_Comm *comms,MG **result)
 71: {
 72:   MG   *mg;
 73:   int  i,ierr,size;
 74:   char *prefix;
 75:   KSP  ksp;
 76:   PC   ipc;

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

 82:   PCGetOptionsPrefix(pc,&prefix);

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

 91:     if (comms) comm = comms[i];
 92:     SLESCreate(comm,&mg[i]->smoothd);
 93:     SLESGetKSP(mg[i]->smoothd,&ksp);
 94:     KSPSetTolerances(ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);
 95:     SLESSetOptionsPrefix(mg[i]->smoothd,prefix);

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

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

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

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

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



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

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

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

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

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

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

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

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

222: #undef __FUNCT__  
224: static int PCSetFromOptions_MG(PC pc)
225: {
226:   int        ierr,m,levels = 1;
227:   PetscTruth flg;
228:   char       buff[16],*type[] = {"additive","multiplicative","full","cascade"};


232:   PetscOptionsHead("Multigrid options");
233:     if (!pc->data) {
234:       PetscOptionsInt("-pc_mg_levels","Number of Levels","MGSetLevels",levels,&levels,&flg);
235:       MGSetLevels(pc,levels,PETSC_NULL);
236:     }
237:     PetscOptionsInt("-pc_mg_cycles","1 for V cycle, 2 for W-cycle","MGSetCycles",1,&m,&flg);
238:     if (flg) {
239:       MGSetCycles(pc,m);
240:     }
241:     PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","MGSetNumberSmoothUp",1,&m,&flg);
242:     if (flg) {
243:       MGSetNumberSmoothUp(pc,m);
244:     }
245:     PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","MGSetNumberSmoothDown",1,&m,&flg);
246:     if (flg) {
247:       MGSetNumberSmoothDown(pc,m);
248:     }
249:     PetscOptionsEList("-pc_mg_type","Multigrid type","MGSetType",type,4,"multiplicative",buff,15,&flg);
250:     if (flg) {
251:       MGType     mg;
252:       PetscTruth isadd,ismult,isfull,iskask,iscasc;

254:       PetscStrcmp(buff,type[0],&isadd);
255:       PetscStrcmp(buff,type[1],&ismult);
256:       PetscStrcmp(buff,type[2],&isfull);
257:       PetscStrcmp(buff,type[3],&iscasc);
258:       PetscStrcmp(buff,"kaskade",&iskask);

260:       if      (isadd)  mg = MGADDITIVE;
261:       else if (ismult) mg = MGMULTIPLICATIVE;
262:       else if (isfull) mg = MGFULL;
263:       else if (iskask) mg = MGKASKADE;
264:       else if (iscasc) mg = MGKASKADE;
265:       else SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Unknown type: %s",buff);
266:       MGSetType(pc,mg);
267:     }
268:     PetscOptionsName("-pc_mg_log","Log times for each multigrid level","None",&flg);
269:     if (flg) {
270:       MG   *mg = (MG*)pc->data;
271:       int  i;
272:       char eventname[128];
273:       if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
274:       levels = mg[0]->levels;
275:       for (i=0; i<levels; i++) {
276:         sprintf(eventname,"MSetup Level %d",i);
277:         PetscLogEventRegister(&mg[i]->eventsetup,eventname,pc->cookie);
278:         sprintf(eventname,"MGSolve Level %d",i);
279:         PetscLogEventRegister(&mg[i]->eventsolve,eventname,pc->cookie);
280:       }
281:     }
282:   PetscOptionsTail();
283:   return(0);
284: }

286: #undef __FUNCT__  
288: static int PCView_MG(PC pc,PetscViewer viewer)
289: {
290:   MG         *mg = (MG*)pc->data;
291:   KSP        kspu,kspd;
292:   int        itu,itd,ierr,levels = mg[0]->levels,i;
293:   PetscReal  dtol,atol,rtol;
294:   char       *cstring;
295:   PetscTruth isascii;

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

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

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

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

383:   for (i=1; i<n; i++) {
384:     if (mg[i]->smoothd) {
385:       SLESGetKSP(mg[i]->smoothd,&ksp);
386:       KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
387:       if (mg[i]->eventsetup) {PetscLogEventBegin(mg[i]->eventsetup,0,0,0,0);}
388:       SLESSetUp(mg[i]->smoothd,mg[i]->b,mg[i]->x);
389:       if (mg[i]->eventsetup) {PetscLogEventEnd(mg[i]->eventsetup,0,0,0,0);}
390:     }
391:   }
392:   for (i=0; i<n; i++) {
393:     if (mg[i]->smoothu && mg[i]->smoothu != mg[i]->smoothd) {
394:       SLESGetKSP(mg[i]->smoothu,&ksp);
395:       KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
396:       if (mg[i]->eventsetup) {PetscLogEventBegin(mg[i]->eventsetup,0,0,0,0);}
397:       SLESSetUp(mg[i]->smoothu,mg[i]->b,mg[i]->x);
398:       if (mg[i]->eventsetup) {PetscLogEventEnd(mg[i]->eventsetup,0,0,0,0);}
399:     }
400:   }

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

416:   if (pc->setupcalled == 0) {
417:     if (monitor) {
418:       SLESGetKSP(mg[0]->smoothd,&ksp);
419:       PetscObjectGetComm((PetscObject)ksp,&comm);
420:       PetscViewerASCIIOpen(comm,"stdout",&ascii);
421:       PetscViewerASCIISetTab(ascii,n);
422:       KSPSetMonitor(ksp,KSPDefaultMonitor,ascii,(int(*)(void*))PetscViewerDestroy);
423:     }
424:     SLESSetFromOptions(mg[0]->smoothd);
425:   }

427:   if (mg[0]->eventsetup) {PetscLogEventBegin(mg[0]->eventsetup,0,0,0,0);}
428:   SLESSetUp(mg[0]->smoothd,mg[0]->b,mg[0]->x);
429:   if (mg[0]->eventsetup) {PetscLogEventEnd(mg[0]->eventsetup,0,0,0,0);}

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

446:   return(0);
447: }

449: /* -------------------------------------------------------------------------------------*/

451: #undef __FUNCT__  
453: /*@C
454:    MGSetLevels - Sets the number of levels to use with MG.
455:    Must be called before any other MG routine.

457:    Collective on PC

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

465:    Level: intermediate

467:    Notes:
468:      If the number of levels is one then the multigrid uses the -mg_levels prefix
469:   for setting the level options rather than the -mg_coarse prefix.

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

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


483:   if (pc->data) {
484:     SETERRQ(1,"Number levels already set for MGn
485:     make sure that you call MGSetLevels() before SLESSetFromOptions()");
486:   }
487:   ierr                     = MGCreate_Private(pc->comm,levels,pc,comms,&mg);
488:   mg[0]->am                = MGMULTIPLICATIVE;
489:   pc->data                 = (void*)mg;
490:   pc->ops->applyrichardson = PCApplyRichardson_MG;
491:   return(0);
492: }

494: #undef __FUNCT__  
496: /*@
497:    MGGetLevels - Gets the number of levels to use with MG.

499:    Not Collective

501:    Input Parameter:
502: .  pc - the preconditioner context

504:    Output parameter:
505: .  levels - the number of levels

507:    Level: advanced

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

511: .seealso: MGSetLevels()
512: @*/
513: int MGGetLevels(PC pc,int *levels)
514: {
515:   MG  *mg;


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

525: #undef __FUNCT__  
527: /*@
528:    MGSetType - Determines the form of multigrid to use:
529:    multiplicative, additive, full, or the Kaskade algorithm.

531:    Collective on PC

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

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

542:    Level: advanced

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

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

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

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

563: #undef __FUNCT__  
565: /*@
566:    MGSetCycles - Sets the type cycles to use.  Use MGSetCyclesOnLevel() for more 
567:    complicated cycling.

569:    Collective on PC

571:    Input Parameters:
572: +  mg - the multigrid context 
573: -  n - the number of cycles

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

578:    Level: advanced

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

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

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

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

601: #undef __FUNCT__  
603: /*@
604:    MGCheck - Checks that all components of the MG structure have 
605:    been set.

607:    Collective on PC

609:    Input Parameters:
610: .  mg - the MG structure

612:    Level: advanced

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

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

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

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

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


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

666:    Collective on PC

668:    Input Parameters:
669: +  mg - the multigrid context 
670: -  n - the number of smoothing steps

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

675:    Level: advanced

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

679: .seealso: MGSetNumberSmoothUp()
680: @*/
681: int MGSetNumberSmoothDown(PC pc,int n)
682: {
683:   MG  *mg;
684:   int i,levels,ierr;
685:   KSP ksp;

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

693:   for (i=0; i<levels; i++) {
694:     SLESGetKSP(mg[i]->smoothd,&ksp);
695:     KSPSetTolerances(ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);
696:     mg[i]->default_smoothd = n;
697:   }
698:   return(0);
699: }

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

708:    Collective on PC

710:    Input Parameters:
711: +  mg - the multigrid context 
712: -  n - the number of smoothing steps

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

717:    Level: advanced

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

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

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

735:   for (i=0; i<levels; i++) {
736:     SLESGetKSP(mg[i]->smoothu,&ksp);
737:     KSPSetTolerances(ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);
738:     mg[i]->default_smoothu = n;
739:   }
740:   return(0);
741: }

743: /* ----------------------------------------------------------------------------------------*/

745: EXTERN_C_BEGIN
746: #undef __FUNCT__  
748: int PCCreate_MG(PC pc)
749: {
751:   pc->ops->apply          = PCApply_MG;
752:   pc->ops->setup          = PCSetUp_MG;
753:   pc->ops->destroy        = PCDestroy_MG;
754:   pc->ops->setfromoptions = PCSetFromOptions_MG;
755:   pc->ops->view           = PCView_MG;

757:   pc->data                = (void*)0;
758:   return(0);
759: }
760: EXTERN_C_END