Actual source code: fasfunc.c

  1: #include <../src/snes/impls/fas/fasimpls.h>

  3: /* -------------- functions called on the fine level -------------- */

  5: /*@
  6:     SNESFASSetType - Sets the update and correction type used for FAS.

  8:    Logically Collective

 10: Input Parameters:
 11: + snes  - FAS context
 12: - fastype  - SNES_FAS_ADDITIVE, SNES_FAS_MULTIPLICATIVE, SNES_FAS_FULL, or SNES_FAS_KASKADE

 14: Level: intermediate

 16: .seealso: PCMGSetType()
 17: @*/
 18: PetscErrorCode  SNESFASSetType(SNES snes,SNESFASType fastype)
 19: {
 20:   SNES_FAS       *fas;

 26:   fas = (SNES_FAS*)snes->data;
 27:   fas->fastype = fastype;
 28:   if (fas->next) {
 29:     SNESFASSetType(fas->next, fastype);
 30:   }
 31:   return(0);
 32: }

 34: /*@
 35: SNESFASGetType - Sets the update and correction type used for FAS.

 37: Logically Collective

 39: Input Parameters:
 40: . snes - FAS context

 42: Output Parameters:
 43: . fastype - SNES_FAS_ADDITIVE or SNES_FAS_MULTIPLICATIVE

 45: Level: intermediate

 47: .seealso: PCMGSetType()
 48: @*/
 49: PetscErrorCode  SNESFASGetType(SNES snes,SNESFASType *fastype)
 50: {
 51:   SNES_FAS *fas;

 56:   fas = (SNES_FAS*)snes->data;
 57:   *fastype = fas->fastype;
 58:   return(0);
 59: }

 61: /*@C
 62:    SNESFASSetLevels - Sets the number of levels to use with FAS.
 63:    Must be called before any other FAS routine.

 65:    Input Parameters:
 66: +  snes   - the snes context
 67: .  levels - the number of levels
 68: -  comms  - optional communicators for each level; this is to allow solving the coarser
 69:             problems on smaller sets of processors.

 71:    Level: intermediate

 73:    Notes:
 74:      If the number of levels is one then the multigrid uses the -fas_levels prefix
 75:   for setting the level options rather than the -fas_coarse prefix.

 77: .seealso: SNESFASGetLevels()
 78: @*/
 79: PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm *comms)
 80: {
 82:   PetscInt       i;
 83:   const char     *optionsprefix;
 84:   char           tprefix[128];
 85:   SNES_FAS       *fas;
 86:   SNES           prevsnes;
 87:   MPI_Comm       comm;

 91:   fas = (SNES_FAS*)snes->data;
 92:   PetscObjectGetComm((PetscObject)snes,&comm);
 93:   if (levels == fas->levels) {
 94:     if (!comms) return(0);
 95:   }
 96:   /* user has changed the number of levels; reset */
 97:   (*snes->ops->reset)(snes);
 98:   /* destroy any coarser levels if necessary */
 99:   SNESDestroy(&fas->next);
100:   fas->next     = NULL;
101:   fas->previous = NULL;
102:   prevsnes      = snes;
103:   /* setup the finest level */
104:   SNESGetOptionsPrefix(snes, &optionsprefix);
105:   PetscObjectComposedDataSetInt((PetscObject) snes, PetscMGLevelId, levels-1);
106:   for (i = levels - 1; i >= 0; i--) {
107:     if (comms) comm = comms[i];
108:     fas->level  = i;
109:     fas->levels = levels;
110:     fas->fine   = snes;
111:     fas->next   = NULL;
112:     if (i > 0) {
113:       SNESCreate(comm, &fas->next);
114:       SNESGetOptionsPrefix(fas->fine, &optionsprefix);
115:       PetscSNPrintf(tprefix,sizeof(tprefix),"fas_levels_%d_cycle_",(int)fas->level);
116:       SNESAppendOptionsPrefix(fas->next,optionsprefix);
117:       SNESAppendOptionsPrefix(fas->next,tprefix);
118:       SNESSetType(fas->next, SNESFAS);
119:       SNESSetTolerances(fas->next, fas->next->abstol, fas->next->rtol, fas->next->stol, fas->n_cycles, fas->next->max_funcs);
120:       PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);
121:       PetscObjectComposedDataSetInt((PetscObject) fas->next, PetscMGLevelId, i-1);

123:       ((SNES_FAS*)fas->next->data)->previous = prevsnes;

125:       prevsnes = fas->next;
126:       fas      = (SNES_FAS*)prevsnes->data;
127:     }
128:   }
129:   return(0);
130: }

132: /*@
133:    SNESFASGetLevels - Gets the number of levels in a FAS, including fine and coarse grids

135:    Input Parameter:
136: .  snes - the nonlinear solver context

138:    Output parameter:
139: .  levels - the number of levels

141:    Level: advanced

143: .seealso: SNESFASSetLevels(), PCMGGetLevels()
144: @*/
145: PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt *levels)
146: {
147:   SNES_FAS *fas;

152:   fas = (SNES_FAS*)snes->data;
153:   *levels = fas->levels;
154:   return(0);
155: }

157: /*@
158:    SNESFASGetCycleSNES - Gets the SNES corresponding to a particular
159:    level of the FAS hierarchy.

161:    Input Parameters:
162: +  snes    - the multigrid context
163:    level   - the level to get
164: -  lsnes   - whether to use the nonlinear smoother or not

166:    Level: advanced

168: .seealso: SNESFASSetLevels(), SNESFASGetLevels()
169: @*/
170: PetscErrorCode SNESFASGetCycleSNES(SNES snes,PetscInt level,SNES *lsnes)
171: {
172:   SNES_FAS *fas;
173:   PetscInt i;

178:   fas = (SNES_FAS*)snes->data;
179:   if (level > fas->levels-1) SETERRQ2(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Requested level %D from SNESFAS containing %D levels",level,fas->levels);
180:   if (fas->level !=  fas->levels - 1) SETERRQ2(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"SNESFASGetCycleSNES may only be called on the finest-level SNES.",level,fas->level);

182:   *lsnes = snes;
183:   for (i = fas->level; i > level; i--) {
184:     *lsnes = fas->next;
185:     fas    = (SNES_FAS*)(*lsnes)->data;
186:   }
187:   if (fas->level != level) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt");
188:   return(0);
189: }

191: /*@
192:    SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to
193:    use on all levels.

195:    Logically Collective on SNES

197:    Input Parameters:
198: +  snes - the multigrid context
199: -  n    - the number of smoothing steps

201:    Options Database Key:
202: .  -snes_fas_smoothup <n> - Sets number of pre-smoothing steps

204:    Level: advanced

206: .seealso: SNESFASSetNumberSmoothDown()
207: @*/
208: PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n)
209: {
210:   SNES_FAS       *fas;

215:   fas =  (SNES_FAS*)snes->data;
216:   fas->max_up_it = n;
217:   if (!fas->smoothu && fas->level != 0) {
218:     SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);
219:   }
220:   if (fas->smoothu) {
221:     SNESSetTolerances(fas->smoothu, fas->smoothu->abstol, fas->smoothu->rtol, fas->smoothu->stol, n, fas->smoothu->max_funcs);
222:   }
223:   if (fas->next) {
224:     SNESFASSetNumberSmoothUp(fas->next, n);
225:   }
226:   return(0);
227: }

229: /*@
230:    SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to
231:    use on all levels.

233:    Logically Collective on SNES

235:    Input Parameters:
236: +  snes - the multigrid context
237: -  n    - the number of smoothing steps

239:    Options Database Key:
240: .  -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps

242:    Level: advanced

244: .seealso: SNESFASSetNumberSmoothUp()
245: @*/
246: PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n)
247: {
248:   SNES_FAS       *fas;

253:   fas = (SNES_FAS*)snes->data;
254:   if (!fas->smoothd) {
255:     SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);
256:   }
257:   SNESSetTolerances(fas->smoothd, fas->smoothd->abstol, fas->smoothd->rtol, fas->smoothd->stol, n, fas->smoothd->max_funcs);

259:   fas->max_down_it = n;
260:   if (fas->next) {
261:     SNESFASSetNumberSmoothDown(fas->next, n);
262:   }
263:   return(0);
264: }

266: /*@
267:    SNESFASSetContinuation - Sets the FAS cycle to default to exact Newton solves on the upsweep

269:    Logically Collective on SNES

271:    Input Parameters:
272: +  snes - the multigrid context
273: -  n    - the number of smoothing steps

275:    Options Database Key:
276: .  -snes_fas_continuation - sets continuation to true

278:    Level: advanced

280:    Notes:
281:     This sets the prefix on the upsweep smoothers to -fas_continuation

283: .seealso: SNESFAS
284: @*/
285: PetscErrorCode SNESFASSetContinuation(SNES snes,PetscBool continuation)
286: {
287:   const char     *optionsprefix;
288:   char           tprefix[128];
289:   SNES_FAS       *fas;

294:   fas  = (SNES_FAS*)snes->data;
295:   SNESGetOptionsPrefix(fas->fine, &optionsprefix);
296:   if (!fas->smoothu) {
297:     SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);
298:   }
299:   PetscStrncpy(tprefix,"fas_levels_continuation_",sizeof(tprefix));
300:   SNESSetOptionsPrefix(fas->smoothu, optionsprefix);
301:   SNESAppendOptionsPrefix(fas->smoothu, tprefix);
302:   SNESSetType(fas->smoothu,SNESNEWTONLS);
303:   SNESSetTolerances(fas->smoothu,fas->fine->abstol,fas->fine->rtol,fas->fine->stol,50,100);
304:   fas->continuation = continuation;
305:   if (fas->next) {
306:     SNESFASSetContinuation(fas->next,continuation);
307:   }
308:   return(0);
309: }

311: /*@
312:    SNESFASSetCycles - Sets the number of FAS multigrid cycles to use each time a grid is visited.  Use SNESFASSetCyclesOnLevel() for more
313:    complicated cycling.

315:    Logically Collective on SNES

317:    Input Parameters:
318: +  snes   - the multigrid context
319: -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle

321:    Options Database Key:
322: .  -snes_fas_cycles 1 or 2

324:    Level: advanced

326: .seealso: SNESFASSetCyclesOnLevel()
327: @*/
328: PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles)
329: {
330:   SNES_FAS       *fas;
332:   PetscBool      isFine;

336:   SNESFASCycleIsFine(snes, &isFine);
337:   fas = (SNES_FAS*)snes->data;
338:   fas->n_cycles = cycles;
339:   if (!isFine) {
340:     SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);
341:   }
342:   if (fas->next) {
343:     SNESFASSetCycles(fas->next, cycles);
344:   }
345:   return(0);
346: }

348: /*@
349:    SNESFASSetMonitor - Sets the method-specific cycle monitoring

351:    Logically Collective on SNES

353:    Input Parameters:
354: +  snes   - the FAS context
355: .  vf     - viewer and format structure (may be NULL if flg is FALSE)
356: -  flg    - monitor or not

358:    Level: advanced

360: .seealso: SNESFASSetCyclesOnLevel()
361: @*/
362: PetscErrorCode SNESFASSetMonitor(SNES snes, PetscViewerAndFormat *vf, PetscBool flg)
363: {
364:   SNES_FAS       *fas;
366:   PetscBool      isFine;
367:   PetscInt       i, levels;
368:   SNES           levelsnes;

372:   SNESFASCycleIsFine(snes, &isFine);
373:   fas = (SNES_FAS*)snes->data;
374:   levels = fas->levels;
375:   if (isFine) {
376:     for (i = 0; i < levels; i++) {
377:       SNESFASGetCycleSNES(snes, i, &levelsnes);
378:       fas  = (SNES_FAS*)levelsnes->data;
379:       if (flg) {
380:         /* set the monitors for the upsmoother and downsmoother */
381:         SNESMonitorCancel(levelsnes);
382:         /* Only register destroy on finest level */
383:         SNESMonitorSet(levelsnes,(PetscErrorCode (*)(SNES,PetscInt,PetscReal,void*))SNESMonitorDefault,vf,(!i ? (PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy : NULL));
384:       } else if (i != fas->levels - 1) {
385:         /* unset the monitors on the coarse levels */
386:         SNESMonitorCancel(levelsnes);
387:       }
388:     }
389:   }
390:   return(0);
391: }

393: /*@
394:    SNESFASSetLog - Sets or unsets time logging for various FAS stages on all levels

396:    Logically Collective on SNES

398:    Input Parameters:
399: +  snes   - the FAS context
400: -  flg    - monitor or not

402:    Level: advanced

404: .seealso: SNESFASSetMonitor()
405: @*/
406: PetscErrorCode SNESFASSetLog(SNES snes, PetscBool flg)
407: {
408:   SNES_FAS       *fas;
410:   PetscBool      isFine;
411:   PetscInt       i, levels;
412:   SNES           levelsnes;
413:   char           eventname[128];

417:   SNESFASCycleIsFine(snes, &isFine);
418:   fas = (SNES_FAS*)snes->data;
419:   levels = fas->levels;
420:   if (isFine) {
421:     for (i = 0; i < levels; i++) {
422:       SNESFASGetCycleSNES(snes, i, &levelsnes);
423:       fas  = (SNES_FAS*)levelsnes->data;
424:       if (flg) {
425:         PetscSNPrintf(eventname,sizeof(eventname),"FASSetup  %d",(int)i);
426:         PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventsmoothsetup);
427:         PetscSNPrintf(eventname,sizeof(eventname),"FASSmooth %d",(int)i);
428:         PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventsmoothsolve);
429:         PetscSNPrintf(eventname,sizeof(eventname),"FASResid  %d",(int)i);
430:         PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventresidual);
431:         PetscSNPrintf(eventname,sizeof(eventname),"FASInterp %d",(int)i);
432:         PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventinterprestrict);
433:       } else {
434:         fas->eventsmoothsetup    = 0;
435:         fas->eventsmoothsolve    = 0;
436:         fas->eventresidual       = 0;
437:         fas->eventinterprestrict = 0;
438:       }
439:     }
440:   }
441:   return(0);
442: }

444: /*
445: Creates the default smoother type.

447: This is SNESNRICHARDSON on each fine level and SNESNEWTONLS on the coarse level.

449:  */
450: PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES snes, SNES *smooth)
451: {
452:   SNES_FAS       *fas;
453:   const char     *optionsprefix;
454:   char           tprefix[128];
456:   SNES           nsmooth;

461:   fas  = (SNES_FAS*)snes->data;
462:   SNESGetOptionsPrefix(fas->fine, &optionsprefix);
463:   /* create the default smoother */
464:   SNESCreate(PetscObjectComm((PetscObject)snes), &nsmooth);
465:   if (fas->level == 0) {
466:     PetscStrncpy(tprefix,"fas_coarse_",sizeof(tprefix));
467:     SNESAppendOptionsPrefix(nsmooth, optionsprefix);
468:     SNESAppendOptionsPrefix(nsmooth, tprefix);
469:     SNESSetType(nsmooth, SNESNEWTONLS);
470:     SNESSetTolerances(nsmooth, nsmooth->abstol, nsmooth->rtol, nsmooth->stol, nsmooth->max_its, nsmooth->max_funcs);
471:   } else {
472:     PetscSNPrintf(tprefix,sizeof(tprefix),"fas_levels_%d_",(int)fas->level);
473:     SNESAppendOptionsPrefix(nsmooth, optionsprefix);
474:     SNESAppendOptionsPrefix(nsmooth, tprefix);
475:     SNESSetType(nsmooth, SNESNRICHARDSON);
476:     SNESSetTolerances(nsmooth, 0.0, 0.0, 0.0, fas->max_down_it, nsmooth->max_funcs);
477:   }
478:   PetscObjectIncrementTabLevel((PetscObject)nsmooth, (PetscObject)snes, 1);
479:   PetscLogObjectParent((PetscObject)snes,(PetscObject)nsmooth);
480:   PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)nsmooth);
481:   PetscObjectComposedDataSetInt((PetscObject) nsmooth, PetscMGLevelId, fas->level);
482:   *smooth = nsmooth;
483:   return(0);
484: }

486: /* ------------- Functions called on a particular level ----------------- */

488: /*@
489:    SNESFASCycleSetCycles - Sets the number of cycles on a particular level.

491:    Logically Collective on SNES

493:    Input Parameters:
494: +  snes   - the multigrid context
495: .  level  - the level to set the number of cycles on
496: -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle

498:    Level: advanced

500: .seealso: SNESFASSetCycles()
501: @*/
502: PetscErrorCode SNESFASCycleSetCycles(SNES snes, PetscInt cycles)
503: {
504:   SNES_FAS       *fas;

509:   fas = (SNES_FAS*)snes->data;
510:   fas->n_cycles = cycles;
511:   SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);
512:   return(0);
513: }

515: /*@
516:    SNESFASCycleGetSmoother - Gets the smoother on a particular cycle level.

518:    Logically Collective on SNES

520:    Input Parameters:
521: .  snes   - the multigrid context

523:    Output Parameters:
524: .  smooth - the smoother

526:    Level: advanced

528: .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmootherDown()
529: @*/
530: PetscErrorCode SNESFASCycleGetSmoother(SNES snes, SNES *smooth)
531: {
532:   SNES_FAS *fas;

537:   fas     = (SNES_FAS*)snes->data;
538:   *smooth = fas->smoothd;
539:   return(0);
540: }
541: /*@
542:    SNESFASCycleGetSmootherUp - Gets the up smoother on a particular cycle level.

544:    Logically Collective on SNES

546:    Input Parameters:
547: .  snes   - the multigrid context

549:    Output Parameters:
550: .  smoothu - the smoother

552:    Notes:
553:    Returns the downsmoother if no up smoother is available.  This enables transparent
554:    default behavior in the process of the solve.

556:    Level: advanced

558: .seealso: SNESFASCycleGetSmoother(), SNESFASCycleGetSmootherDown()
559: @*/
560: PetscErrorCode SNESFASCycleGetSmootherUp(SNES snes, SNES *smoothu)
561: {
562:   SNES_FAS *fas;

567:   fas = (SNES_FAS*)snes->data;
568:   if (!fas->smoothu) *smoothu = fas->smoothd;
569:   else *smoothu = fas->smoothu;
570:   return(0);
571: }

573: /*@
574:    SNESFASCycleGetSmootherDown - Gets the down smoother on a particular cycle level.

576:    Logically Collective on SNES

578:    Input Parameters:
579: .  snes   - the multigrid context

581:    Output Parameters:
582: .  smoothd - the smoother

584:    Level: advanced

586: .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
587: @*/
588: PetscErrorCode SNESFASCycleGetSmootherDown(SNES snes, SNES *smoothd)
589: {
590:   SNES_FAS *fas;

595:   fas = (SNES_FAS*)snes->data;
596:   *smoothd = fas->smoothd;
597:   return(0);
598: }

600: /*@
601:    SNESFASCycleGetCorrection - Gets the coarse correction FAS context for this level

603:    Logically Collective on SNES

605:    Input Parameters:
606: .  snes   - the multigrid context

608:    Output Parameters:
609: .  correction - the coarse correction on this level

611:    Notes:
612:    Returns NULL on the coarsest level.

614:    Level: advanced

616: .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
617: @*/
618: PetscErrorCode SNESFASCycleGetCorrection(SNES snes, SNES *correction)
619: {
620:   SNES_FAS *fas;

625:   fas = (SNES_FAS*)snes->data;
626:   *correction = fas->next;
627:   return(0);
628: }

630: /*@
631:    SNESFASCycleGetInterpolation - Gets the interpolation on this level

633:    Logically Collective on SNES

635:    Input Parameters:
636: .  snes   - the multigrid context

638:    Output Parameters:
639: .  mat    - the interpolation operator on this level

641:    Level: developer

643: .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
644: @*/
645: PetscErrorCode SNESFASCycleGetInterpolation(SNES snes, Mat *mat)
646: {
647:   SNES_FAS *fas;

652:   fas = (SNES_FAS*)snes->data;
653:   *mat = fas->interpolate;
654:   return(0);
655: }

657: /*@
658:    SNESFASCycleGetRestriction - Gets the restriction on this level

660:    Logically Collective on SNES

662:    Input Parameters:
663: .  snes   - the multigrid context

665:    Output Parameters:
666: .  mat    - the restriction operator on this level

668:    Level: developer

670: .seealso: SNESFASGetRestriction(), SNESFASCycleGetInterpolation()
671: @*/
672: PetscErrorCode SNESFASCycleGetRestriction(SNES snes, Mat *mat)
673: {
674:   SNES_FAS *fas;

679:   fas = (SNES_FAS*)snes->data;
680:   *mat = fas->restrct;
681:   return(0);
682: }

684: /*@
685:    SNESFASCycleGetInjection - Gets the injection on this level

687:    Logically Collective on SNES

689:    Input Parameters:
690: .  snes   - the multigrid context

692:    Output Parameters:
693: .  mat    - the restriction operator on this level

695:    Level: developer

697: .seealso: SNESFASGetInjection(), SNESFASCycleGetRestriction()
698: @*/
699: PetscErrorCode SNESFASCycleGetInjection(SNES snes, Mat *mat)
700: {
701:   SNES_FAS *fas;

706:   fas = (SNES_FAS*)snes->data;
707:   *mat = fas->inject;
708:   return(0);
709: }

711: /*@
712:    SNESFASCycleGetRScale - Gets the injection on this level

714:    Logically Collective on SNES

716:    Input Parameters:
717: .  snes   - the multigrid context

719:    Output Parameters:
720: .  mat    - the restriction operator on this level

722:    Level: developer

724: .seealso: SNESFASCycleGetRestriction(), SNESFASGetRScale()
725: @*/
726: PetscErrorCode SNESFASCycleGetRScale(SNES snes, Vec *vec)
727: {
728:   SNES_FAS *fas;

733:   fas  = (SNES_FAS*)snes->data;
734:   *vec = fas->rscale;
735:   return(0);
736: }

738: /*@
739:    SNESFASCycleIsFine - Determines if a given cycle is the fine level.

741:    Logically Collective on SNES

743:    Input Parameters:
744: .  snes   - the FAS context

746:    Output Parameters:
747: .  flg - indicates if this is the fine level or not

749:    Level: advanced

751: .seealso: SNESFASSetLevels()
752: @*/
753: PetscErrorCode SNESFASCycleIsFine(SNES snes, PetscBool *flg)
754: {
755:   SNES_FAS *fas;

760:   fas = (SNES_FAS*)snes->data;
761:   if (fas->level == fas->levels - 1) *flg = PETSC_TRUE;
762:   else *flg = PETSC_FALSE;
763:   return(0);
764: }

766: /* ---------- functions called on the finest level that return level-specific information ---------- */

768: /*@
769:    SNESFASSetInterpolation - Sets the function to be used to calculate the
770:    interpolation from l-1 to the lth level

772:    Input Parameters:
773: +  snes      - the multigrid context
774: .  mat       - the interpolation operator
775: -  level     - the level (0 is coarsest) to supply [do not supply 0]

777:    Level: advanced

779:    Notes:
780:           Usually this is the same matrix used also to set the restriction
781:     for the same level.

783:           One can pass in the interpolation matrix or its transpose; PETSc figures
784:     out from the matrix size which one it is.

786: .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRScale()
787: @*/
788: PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat)
789: {
790:   SNES_FAS       *fas;
792:   SNES           levelsnes;

797:   SNESFASGetCycleSNES(snes, level, &levelsnes);
798:   fas  = (SNES_FAS*)levelsnes->data;
799:   PetscObjectReference((PetscObject)mat);
800:   MatDestroy(&fas->interpolate);
801:   fas->interpolate = mat;
802:   return(0);
803: }

805: /*@
806:    SNESFASGetInterpolation - Gets the matrix used to calculate the
807:    interpolation from l-1 to the lth level

809:    Input Parameters:
810: +  snes      - the multigrid context
811: -  level     - the level (0 is coarsest) to supply [do not supply 0]

813:    Output Parameters:
814: .  mat       - the interpolation operator

816:    Level: advanced

818: .seealso: SNESFASSetInterpolation(), SNESFASGetInjection(), SNESFASGetRestriction(), SNESFASGetRScale()
819: @*/
820: PetscErrorCode SNESFASGetInterpolation(SNES snes, PetscInt level, Mat *mat)
821: {
822:   SNES_FAS       *fas;
824:   SNES           levelsnes;

829:   SNESFASGetCycleSNES(snes, level, &levelsnes);
830:   fas  = (SNES_FAS*)levelsnes->data;
831:   *mat = fas->interpolate;
832:   return(0);
833: }

835: /*@
836:    SNESFASSetRestriction - Sets the function to be used to restrict the defect
837:    from level l to l-1.

839:    Input Parameters:
840: +  snes  - the multigrid context
841: .  mat   - the restriction matrix
842: -  level - the level (0 is coarsest) to supply [Do not supply 0]

844:    Level: advanced

846:    Notes:
847:           Usually this is the same matrix used also to set the interpolation
848:     for the same level.

850:           One can pass in the interpolation matrix or its transpose; PETSc figures
851:     out from the matrix size which one it is.

853:          If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation()
854:     is used.

856: .seealso: SNESFASSetInterpolation(), SNESFASSetInjection()
857: @*/
858: PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat)
859: {
860:   SNES_FAS       *fas;
862:   SNES           levelsnes;

867:   SNESFASGetCycleSNES(snes, level, &levelsnes);
868:   fas  = (SNES_FAS*)levelsnes->data;
869:   PetscObjectReference((PetscObject)mat);
870:   MatDestroy(&fas->restrct);
871:   fas->restrct = mat;
872:   return(0);
873: }

875: /*@
876:    SNESFASGetRestriction - Gets the matrix used to calculate the
877:    restriction from l to the l-1th level

879:    Input Parameters:
880: +  snes      - the multigrid context
881: -  level     - the level (0 is coarsest) to supply [do not supply 0]

883:    Output Parameters:
884: .  mat       - the interpolation operator

886:    Level: advanced

888: .seealso: SNESFASSetRestriction(), SNESFASGetInjection(), SNESFASGetInterpolation(), SNESFASGetRScale()
889: @*/
890: PetscErrorCode SNESFASGetRestriction(SNES snes, PetscInt level, Mat *mat)
891: {
892:   SNES_FAS       *fas;
894:   SNES           levelsnes;

899:   SNESFASGetCycleSNES(snes, level, &levelsnes);
900:   fas  = (SNES_FAS*)levelsnes->data;
901:   *mat = fas->restrct;
902:   return(0);
903: }

905: /*@
906:    SNESFASSetInjection - Sets the function to be used to inject the solution
907:    from level l to l-1.

909:    Input Parameters:
910:  +  snes  - the multigrid context
911: .  mat   - the restriction matrix
912: -  level - the level (0 is coarsest) to supply [Do not supply 0]

914:    Level: advanced

916:    Notes:
917:          If you do not set this, the restriction and rscale is used to
918:    project the solution instead.

920: .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction()
921: @*/
922: PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat)
923: {
924:   SNES_FAS       *fas;
926:   SNES           levelsnes;

931:   SNESFASGetCycleSNES(snes, level, &levelsnes);
932:   fas  = (SNES_FAS*)levelsnes->data;
933:   PetscObjectReference((PetscObject)mat);
934:   MatDestroy(&fas->inject);

936:   fas->inject = mat;
937:   return(0);
938: }

940: /*@
941:    SNESFASGetInjection - Gets the matrix used to calculate the
942:    injection from l-1 to the lth level

944:    Input Parameters:
945: +  snes      - the multigrid context
946: -  level     - the level (0 is coarsest) to supply [do not supply 0]

948:    Output Parameters:
949: .  mat       - the injection operator

951:    Level: advanced

953: .seealso: SNESFASSetInjection(), SNESFASGetRestriction(), SNESFASGetInterpolation(), SNESFASGetRScale()
954: @*/
955: PetscErrorCode SNESFASGetInjection(SNES snes, PetscInt level, Mat *mat)
956: {
957:   SNES_FAS       *fas;
959:   SNES           levelsnes;

964:   SNESFASGetCycleSNES(snes, level, &levelsnes);
965:   fas  = (SNES_FAS*)levelsnes->data;
966:   *mat = fas->inject;
967:   return(0);
968: }

970: /*@
971:    SNESFASSetRScale - Sets the scaling factor of the restriction
972:    operator from level l to l-1.

974:    Input Parameters:
975: +  snes   - the multigrid context
976: .  rscale - the restriction scaling
977: -  level  - the level (0 is coarsest) to supply [Do not supply 0]

979:    Level: advanced

981:    Notes:
982:          This is only used in the case that the injection is not set.

984: .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
985: @*/
986: PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale)
987: {
988:   SNES_FAS       *fas;
990:   SNES           levelsnes;

995:   SNESFASGetCycleSNES(snes, level, &levelsnes);
996:   fas  = (SNES_FAS*)levelsnes->data;
997:   PetscObjectReference((PetscObject)rscale);
998:   VecDestroy(&fas->rscale);
999:   fas->rscale = rscale;
1000:   return(0);
1001: }

1003: /*@
1004:    SNESFASGetSmoother - Gets the default smoother on a level.

1006:    Input Parameters:
1007: +  snes   - the multigrid context
1008: -  level  - the level (0 is coarsest) to supply

1010:    Output Parameters:
1011:    smooth  - the smoother

1013:    Level: advanced

1015: .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1016: @*/
1017: PetscErrorCode SNESFASGetSmoother(SNES snes, PetscInt level, SNES *smooth)
1018: {
1019:   SNES_FAS       *fas;
1021:   SNES           levelsnes;

1026:   SNESFASGetCycleSNES(snes, level, &levelsnes);
1027:   fas  = (SNES_FAS*)levelsnes->data;
1028:   if (!fas->smoothd) {
1029:     SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);
1030:   }
1031:   *smooth = fas->smoothd;
1032:   return(0);
1033: }

1035: /*@
1036:    SNESFASGetSmootherDown - Gets the downsmoother on a level.

1038:    Input Parameters:
1039: +  snes   - the multigrid context
1040: -  level  - the level (0 is coarsest) to supply

1042:    Output Parameters:
1043:    smooth  - the smoother

1045:    Level: advanced

1047: .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1048: @*/
1049: PetscErrorCode SNESFASGetSmootherDown(SNES snes, PetscInt level, SNES *smooth)
1050: {
1051:   SNES_FAS       *fas;
1053:   SNES           levelsnes;

1058:   SNESFASGetCycleSNES(snes, level, &levelsnes);
1059:   fas  = (SNES_FAS*)levelsnes->data;
1060:   /* if the user chooses to differentiate smoothers, create them both at this point */
1061:   if (!fas->smoothd) {
1062:     SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);
1063:   }
1064:   if (!fas->smoothu) {
1065:     SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu);
1066:   }
1067:   *smooth = fas->smoothd;
1068:   return(0);
1069: }

1071: /*@
1072:    SNESFASGetSmootherUp - Gets the upsmoother on a level.

1074:    Input Parameters:
1075: +  snes   - the multigrid context
1076: -  level  - the level (0 is coarsest)

1078:    Output Parameters:
1079:    smooth  - the smoother

1081:    Level: advanced

1083: .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1084: @*/
1085: PetscErrorCode SNESFASGetSmootherUp(SNES snes, PetscInt level, SNES *smooth)
1086: {
1087:   SNES_FAS       *fas;
1089:   SNES           levelsnes;

1094:   SNESFASGetCycleSNES(snes, level, &levelsnes);
1095:   fas  = (SNES_FAS*)levelsnes->data;
1096:   /* if the user chooses to differentiate smoothers, create them both at this point */
1097:   if (!fas->smoothd) {
1098:     SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);
1099:   }
1100:   if (!fas->smoothu) {
1101:     SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu);
1102:   }
1103:   *smooth = fas->smoothu;
1104:   return(0);
1105: }

1107: /*@
1108:   SNESFASGetCoarseSolve - Gets the coarsest solver.

1110:   Input Parameters:
1111: . snes - the multigrid context

1113:   Output Parameters:
1114: . coarse - the coarse-level solver

1116:   Level: advanced

1118: .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1119: @*/
1120: PetscErrorCode SNESFASGetCoarseSolve(SNES snes, SNES *coarse)
1121: {
1122:   SNES_FAS       *fas;
1124:   SNES           levelsnes;

1129:   SNESFASGetCycleSNES(snes, 0, &levelsnes);
1130:   fas  = (SNES_FAS*)levelsnes->data;
1131:   /* if the user chooses to differentiate smoothers, create them both at this point */
1132:   if (!fas->smoothd) {
1133:     SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);
1134:   }
1135:   *coarse = fas->smoothd;
1136:   return(0);
1137: }

1139: /*@
1140:    SNESFASFullSetDownSweep - Smooth during the initial downsweep for SNESFAS

1142:    Logically Collective on SNES

1144:    Input Parameters:
1145: +  snes - the multigrid context
1146: -  swp - whether to downsweep or not

1148:    Options Database Key:
1149: .  -snes_fas_full_downsweep - Sets number of pre-smoothing steps

1151:    Level: advanced

1153: .seealso: SNESFASSetNumberSmoothUp()
1154: @*/
1155: PetscErrorCode SNESFASFullSetDownSweep(SNES snes,PetscBool swp)
1156: {
1157:   SNES_FAS       *fas;

1162:   fas = (SNES_FAS*)snes->data;
1163:   fas->full_downsweep = swp;
1164:   if (fas->next) {
1165:     SNESFASFullSetDownSweep(fas->next,swp);
1166:   }
1167:   return(0);
1168: }

1170: /*@
1171:    SNESFASFullSetTotal - Use total residual restriction and total interpolation on the initial down and up sweep of full FAS cycles

1173:    Logically Collective on SNES

1175:    Input Parameters:
1176: +  snes - the multigrid context
1177: -  total - whether to use total restriction / interpolatiaon or not (the alternative is defect restriction and correction interpolation)

1179:    Options Database Key:
1180: .  -snes_fas_full_total

1182:    Level: advanced

1184:    Note: this option is only significant if the interpolation of a coarse correction (MatInterpolate()) is significantly different from total solution interpolation (DMInterpolateSolution()).

1186: .seealso: SNESFASSetNumberSmoothUp(), DMInterpolateSolution()
1187: @*/
1188: PetscErrorCode SNESFASFullSetTotal(SNES snes,PetscBool total)
1189: {
1190:   SNES_FAS       *fas;

1195:   fas = (SNES_FAS*)snes->data;
1196:   fas->full_total = total;
1197:   if (fas->next) {
1198:     SNESFASFullSetTotal(fas->next,total);
1199:   }
1200:   return(0);
1201: }

1203: /*@
1204:    SNESFASFullGetTotal - Use total residual restriction and total interpolation on the initial down and up sweep of full FAS cycles

1206:    Logically Collective on SNES

1208:    Input Parameters:
1209: .  snes - the multigrid context

1211:    Output:
1212: .  total - whether to use total restriction / interpolatiaon or not (the alternative is defect restriction and correction interpolation)

1214:    Options Database Key:
1215: .  -snes_fas_full_total

1217:    Level: advanced

1219: .seealso: SNESFASSetNumberSmoothUp(), DMInterpolateSolution(), SNESFullSetTotal()
1220: @*/
1221: PetscErrorCode SNESFASFullGetTotal(SNES snes,PetscBool *total)
1222: {
1223:   SNES_FAS       *fas;

1227:   fas = (SNES_FAS*)snes->data;
1228:   *total = fas->full_total;
1229:   return(0);
1230: }