Actual source code: fasfunc.c

petsc-3.8.4 2018-03-24
Report Typos and Errors
  1:  #include <../src/snes/impls/fas/fasimpls.h>


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

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

  9:    Logically Collective

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

 15: Level: intermediate

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

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


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

 38: Logically Collective

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

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

 46: Level: intermediate

 48: .seealso: PCMGSetType()
 49: @*/
 50: PetscErrorCode  SNESFASGetType(SNES snes,SNESFASType *fastype)
 51: {
 52:   SNES_FAS *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: .keywords: FAS, MG, set, levels, multigrid

 79: .seealso: SNESFASGetLevels()
 80: @*/
 81: PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms)
 82: {
 84:   PetscInt       i;
 85:   const char     *optionsprefix;
 86:   char           tprefix[128];
 87:   SNES_FAS       *fas = (SNES_FAS*)snes->data;
 88:   SNES           prevsnes;
 89:   MPI_Comm       comm;

 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:   SNESReset(snes);
 98:   /* destroy any coarser levels if necessary */
 99:   if (fas->next) 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:       sprintf(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: }


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

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

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

142:    Level: advanced

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

146: .seealso: SNESFASSetLevels(), PCMGGetLevels()
147: @*/
148: PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt *levels)
149: {
150:   SNES_FAS * fas = (SNES_FAS*)snes->data;

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


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

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

167:    Level: advanced

169: .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid

171: .seealso: SNESFASSetLevels(), SNESFASGetLevels()
172: @*/
173: PetscErrorCode SNESFASGetCycleSNES(SNES snes,PetscInt level,SNES *lsnes)
174: {
175:   SNES_FAS *fas = (SNES_FAS*)snes->data;
176:   PetscInt i;

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: .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid

208: .seealso: SNESFASSetNumberSmoothDown()
209: @*/
210: PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n)
211: {
212:   SNES_FAS       *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: .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid

246: .seealso: SNESFASSetNumberSmoothUp()
247: @*/
248: PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n)
249: {
250:   SNES_FAS       *fas =  (SNES_FAS*)snes->data;
251:   PetscErrorCode 0;

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


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

270:    Logically Collective on SNES

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

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

279:    Level: advanced

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

283: .keywords: FAS, MG, smoother, continuation

285: .seealso: SNESFAS
286: @*/
287: PetscErrorCode SNESFASSetContinuation(SNES snes,PetscBool continuation)
288: {
289:   const char     *optionsprefix;
290:   char           tprefix[128];
291:   SNES_FAS       *fas =  (SNES_FAS*)snes->data;
292:   PetscErrorCode 0;

295:   SNESGetOptionsPrefix(fas->fine, &optionsprefix);
296:   if (!fas->smoothu) {
297:     SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);
298:   }
299:   sprintf(tprefix,"fas_levels_continuation_");
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: }


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

316:    Logically Collective on SNES

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

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

325:    Level: advanced

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

329: .seealso: SNESFASSetCyclesOnLevel()
330: @*/
331: PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles)
332: {
333:   SNES_FAS       *fas = (SNES_FAS*)snes->data;
335:   PetscBool      isFine;

338:   SNESFASCycleIsFine(snes, &isFine);

340:   fas->n_cycles = cycles;
341:   if (!isFine) {
342:     SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);
343:   }
344:   if (fas->next) {
345:     SNESFASSetCycles(fas->next, cycles);
346:   }
347:   return(0);
348: }


351: /*@
352:    SNESFASSetMonitor - Sets the method-specific cycle monitoring

354:    Logically Collective on SNES

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

361:    Level: advanced

363: .keywords: FAS, monitor

365: .seealso: SNESFASSetCyclesOnLevel()
366: @*/
367: PetscErrorCode SNESFASSetMonitor(SNES snes, PetscViewerAndFormat *vf, PetscBool flg)
368: {
369:   SNES_FAS       *fas = (SNES_FAS*)snes->data;
371:   PetscBool      isFine;
372:   PetscInt       i, levels = fas->levels;
373:   SNES           levelsnes;

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

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

398:    Logically Collective on SNES

400:    Input Parameters:
401: +  snes   - the FAS context
402: -  flg    - monitor or not

404:    Level: advanced

406: .keywords: FAS, logging

408: .seealso: SNESFASSetMonitor()
409: @*/
410: PetscErrorCode SNESFASSetLog(SNES snes, PetscBool flg)
411: {
412:   SNES_FAS       *fas = (SNES_FAS*)snes->data;
414:   PetscBool      isFine;
415:   PetscInt       i, levels = fas->levels;
416:   SNES           levelsnes;
417:   char           eventname[128];

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

447: /*
448: Creates the default smoother type.

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

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

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

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

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

493:    Logically Collective on SNES

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

500:    Level: advanced

502: .keywords: SNES, FAS, set, cycles, V-cycle, W-cycle, multigrid

504: .seealso: SNESFASSetCycles()
505: @*/
506: PetscErrorCode SNESFASCycleSetCycles(SNES snes, PetscInt cycles)
507: {
508:   SNES_FAS       *fas =  (SNES_FAS*)snes->data;

512:   fas->n_cycles = cycles;
513:   SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);
514:   return(0);
515: }


518: /*@
519:    SNESFASCycleGetSmoother - Gets the smoother on a particular cycle level.

521:    Logically Collective on SNES

523:    Input Parameters:
524: .  snes   - the multigrid context

526:    Output Parameters:
527: .  smooth - the smoother

529:    Level: advanced

531: .keywords: SNES, FAS, get, smoother, multigrid

533: .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmootherDown()
534: @*/
535: PetscErrorCode SNESFASCycleGetSmoother(SNES snes, SNES *smooth)
536: {
537:   SNES_FAS *fas;

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

548:    Logically Collective on SNES

550:    Input Parameters:
551: .  snes   - the multigrid context

553:    Output Parameters:
554: .  smoothu - the smoother

556:    Notes:
557:    Returns the downsmoother if no up smoother is available.  This enables transparent
558:    default behavior in the process of the solve.

560:    Level: advanced

562: .keywords: SNES, FAS, get, smoother, multigrid

564: .seealso: SNESFASCycleGetSmoother(), SNESFASCycleGetSmootherDown()
565: @*/
566: PetscErrorCode SNESFASCycleGetSmootherUp(SNES snes, SNES *smoothu)
567: {
568:   SNES_FAS *fas;

572:   fas = (SNES_FAS*)snes->data;
573:   if (!fas->smoothu) *smoothu = fas->smoothd;
574:   else *smoothu = fas->smoothu;
575:   return(0);
576: }

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

581:    Logically Collective on SNES

583:    Input Parameters:
584: .  snes   - the multigrid context

586:    Output Parameters:
587: .  smoothd - the smoother

589:    Level: advanced

591: .keywords: SNES, FAS, get, smoother, multigrid

593: .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
594: @*/
595: PetscErrorCode SNESFASCycleGetSmootherDown(SNES snes, SNES *smoothd)
596: {
597:   SNES_FAS *fas;

601:   fas      = (SNES_FAS*)snes->data;
602:   *smoothd = fas->smoothd;
603:   return(0);
604: }


607: /*@
608:    SNESFASCycleGetCorrection - Gets the coarse correction FAS context for this level

610:    Logically Collective on SNES

612:    Input Parameters:
613: .  snes   - the multigrid context

615:    Output Parameters:
616: .  correction - the coarse correction on this level

618:    Notes:
619:    Returns NULL on the coarsest level.

621:    Level: advanced

623: .keywords: SNES, FAS, get, smoother, multigrid

625: .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
626: @*/
627: PetscErrorCode SNESFASCycleGetCorrection(SNES snes, SNES *correction)
628: {
629:   SNES_FAS *fas;

633:   fas         = (SNES_FAS*)snes->data;
634:   *correction = fas->next;
635:   return(0);
636: }

638: /*@
639:    SNESFASCycleGetInterpolation - Gets the interpolation on this level

641:    Logically Collective on SNES

643:    Input Parameters:
644: .  snes   - the multigrid context

646:    Output Parameters:
647: .  mat    - the interpolation operator on this level

649:    Level: developer

651: .keywords: SNES, FAS, get, smoother, multigrid

653: .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
654: @*/
655: PetscErrorCode SNESFASCycleGetInterpolation(SNES snes, Mat *mat)
656: {
657:   SNES_FAS *fas;

661:   fas  = (SNES_FAS*)snes->data;
662:   *mat = fas->interpolate;
663:   return(0);
664: }


667: /*@
668:    SNESFASCycleGetRestriction - Gets the restriction on this level

670:    Logically Collective on SNES

672:    Input Parameters:
673: .  snes   - the multigrid context

675:    Output Parameters:
676: .  mat    - the restriction operator on this level

678:    Level: developer

680: .keywords: SNES, FAS, get, smoother, multigrid

682: .seealso: SNESFASGetRestriction(), SNESFASCycleGetInterpolation()
683: @*/
684: PetscErrorCode SNESFASCycleGetRestriction(SNES snes, Mat *mat)
685: {
686:   SNES_FAS *fas;

690:   fas  = (SNES_FAS*)snes->data;
691:   *mat = fas->restrct;
692:   return(0);
693: }


696: /*@
697:    SNESFASCycleGetInjection - Gets the injection on this level

699:    Logically Collective on SNES

701:    Input Parameters:
702: .  snes   - the multigrid context

704:    Output Parameters:
705: .  mat    - the restriction operator on this level

707:    Level: developer

709: .keywords: SNES, FAS, get, smoother, multigrid

711: .seealso: SNESFASGetInjection(), SNESFASCycleGetRestriction()
712: @*/
713: PetscErrorCode SNESFASCycleGetInjection(SNES snes, Mat *mat)
714: {
715:   SNES_FAS *fas;

719:   fas  = (SNES_FAS*)snes->data;
720:   *mat = fas->inject;
721:   return(0);
722: }

724: /*@
725:    SNESFASCycleGetRScale - Gets the injection on this level

727:    Logically Collective on SNES

729:    Input Parameters:
730: .  snes   - the multigrid context

732:    Output Parameters:
733: .  mat    - the restriction operator on this level

735:    Level: developer

737: .keywords: SNES, FAS, get, smoother, multigrid

739: .seealso: SNESFASCycleGetRestriction(), SNESFASGetRScale()
740: @*/
741: PetscErrorCode SNESFASCycleGetRScale(SNES snes, Vec *vec)
742: {
743:   SNES_FAS *fas;

747:   fas  = (SNES_FAS*)snes->data;
748:   *vec = fas->rscale;
749:   return(0);
750: }

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

755:    Logically Collective on SNES

757:    Input Parameters:
758: .  snes   - the FAS context

760:    Output Parameters:
761: .  flg - indicates if this is the fine level or not

763:    Level: advanced

765: .keywords: SNES, FAS

767: .seealso: SNESFASSetLevels()
768: @*/
769: PetscErrorCode SNESFASCycleIsFine(SNES snes, PetscBool *flg)
770: {
771:   SNES_FAS *fas;

775:   fas = (SNES_FAS*)snes->data;
776:   if (fas->level == fas->levels - 1) *flg = PETSC_TRUE;
777:   else *flg = PETSC_FALSE;
778:   return(0);
779: }

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

783: /*@
784:    SNESFASSetInterpolation - Sets the function to be used to calculate the
785:    interpolation from l-1 to the lth level

787:    Input Parameters:
788: +  snes      - the multigrid context
789: .  mat       - the interpolation operator
790: -  level     - the level (0 is coarsest) to supply [do not supply 0]

792:    Level: advanced

794:    Notes:
795:           Usually this is the same matrix used also to set the restriction
796:     for the same level.

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

801: .keywords:  FAS, multigrid, set, interpolate, level

803: .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRScale()
804: @*/
805: PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat)
806: {
807:   SNES_FAS       *fas;
809:   SNES           levelsnes;

812:   SNESFASGetCycleSNES(snes, level, &levelsnes);
813:   fas  = (SNES_FAS*)levelsnes->data;
814:   PetscObjectReference((PetscObject)mat);
815:   MatDestroy(&fas->interpolate);

817:   fas->interpolate = mat;
818:   return(0);
819: }

821: /*@
822:    SNESFASGetInterpolation - Gets the matrix used to calculate the
823:    interpolation from l-1 to the lth level

825:    Input Parameters:
826: +  snes      - the multigrid context
827: -  level     - the level (0 is coarsest) to supply [do not supply 0]

829:    Output Parameters:
830: .  mat       - the interpolation operator

832:    Level: advanced

834: .keywords:  FAS, multigrid, get, interpolate, level

836: .seealso: SNESFASSetInterpolation(), SNESFASGetInjection(), SNESFASGetRestriction(), SNESFASGetRScale()
837: @*/
838: PetscErrorCode SNESFASGetInterpolation(SNES snes, PetscInt level, Mat *mat)
839: {
840:   SNES_FAS       *fas;
842:   SNES           levelsnes;

845:   SNESFASGetCycleSNES(snes, level, &levelsnes);
846:   fas  = (SNES_FAS*)levelsnes->data;
847:   *mat = fas->interpolate;
848:   return(0);
849: }

851: /*@
852:    SNESFASSetRestriction - Sets the function to be used to restrict the defect
853:    from level l to l-1.

855:    Input Parameters:
856: +  snes  - the multigrid context
857: .  mat   - the restriction matrix
858: -  level - the level (0 is coarsest) to supply [Do not supply 0]

860:    Level: advanced

862:    Notes:
863:           Usually this is the same matrix used also to set the interpolation
864:     for the same level.

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

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

872: .keywords: FAS, MG, set, multigrid, restriction, level

874: .seealso: SNESFASSetInterpolation(), SNESFASSetInjection()
875: @*/
876: PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat)
877: {
878:   SNES_FAS       *fas;
880:   SNES           levelsnes;

883:   SNESFASGetCycleSNES(snes, level, &levelsnes);
884:   fas  = (SNES_FAS*)levelsnes->data;
885:   PetscObjectReference((PetscObject)mat);
886:   MatDestroy(&fas->restrct);

888:   fas->restrct = mat;
889:   return(0);
890: }

892: /*@
893:    SNESFASGetRestriction - Gets the matrix used to calculate the
894:    restriction from l to the l-1th level

896:    Input Parameters:
897: +  snes      - the multigrid context
898: -  level     - the level (0 is coarsest) to supply [do not supply 0]

900:    Output Parameters:
901: .  mat       - the interpolation operator

903:    Level: advanced

905: .keywords:  FAS, multigrid, get, restrict, level

907: .seealso: SNESFASSetRestriction(), SNESFASGetInjection(), SNESFASGetInterpolation(), SNESFASGetRScale()
908: @*/
909: PetscErrorCode SNESFASGetRestriction(SNES snes, PetscInt level, Mat *mat)
910: {
911:   SNES_FAS       *fas;
913:   SNES           levelsnes;

916:   SNESFASGetCycleSNES(snes, level, &levelsnes);
917:   fas  = (SNES_FAS*)levelsnes->data;
918:   *mat = fas->restrct;
919:   return(0);
920: }


923: /*@
924:    SNESFASSetInjection - Sets the function to be used to inject the solution
925:    from level l to l-1.

927:    Input Parameters:
928:  +  snes  - the multigrid context
929: .  mat   - the restriction matrix
930: -  level - the level (0 is coarsest) to supply [Do not supply 0]

932:    Level: advanced

934:    Notes:
935:          If you do not set this, the restriction and rscale is used to
936:    project the solution instead.

938: .keywords: FAS, MG, set, multigrid, restriction, level

940: .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction()
941: @*/
942: PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat)
943: {
944:   SNES_FAS       *fas;
946:   SNES           levelsnes;

949:   SNESFASGetCycleSNES(snes, level, &levelsnes);
950:   fas  = (SNES_FAS*)levelsnes->data;
951:   PetscObjectReference((PetscObject)mat);
952:   MatDestroy(&fas->inject);

954:   fas->inject = mat;
955:   return(0);
956: }


959: /*@
960:    SNESFASGetInjection - Gets the matrix used to calculate the
961:    injection from l-1 to the lth level

963:    Input Parameters:
964: +  snes      - the multigrid context
965: -  level     - the level (0 is coarsest) to supply [do not supply 0]

967:    Output Parameters:
968: .  mat       - the injection operator

970:    Level: advanced

972: .keywords:  FAS, multigrid, get, restrict, level

974: .seealso: SNESFASSetInjection(), SNESFASGetRestriction(), SNESFASGetInterpolation(), SNESFASGetRScale()
975: @*/
976: PetscErrorCode SNESFASGetInjection(SNES snes, PetscInt level, Mat *mat)
977: {
978:   SNES_FAS       *fas;
980:   SNES           levelsnes;

983:   SNESFASGetCycleSNES(snes, level, &levelsnes);
984:   fas  = (SNES_FAS*)levelsnes->data;
985:   *mat = fas->inject;
986:   return(0);
987: }

989: /*@
990:    SNESFASSetRScale - Sets the scaling factor of the restriction
991:    operator from level l to l-1.

993:    Input Parameters:
994: +  snes   - the multigrid context
995: .  rscale - the restriction scaling
996: -  level  - the level (0 is coarsest) to supply [Do not supply 0]

998:    Level: advanced

1000:    Notes:
1001:          This is only used in the case that the injection is not set.

1003: .keywords: FAS, MG, set, multigrid, restriction, level

1005: .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1006: @*/
1007: PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale)
1008: {
1009:   SNES_FAS       *fas;
1011:   SNES           levelsnes;

1014:   SNESFASGetCycleSNES(snes, level, &levelsnes);
1015:   fas  = (SNES_FAS*)levelsnes->data;
1016:   PetscObjectReference((PetscObject)rscale);
1017:   VecDestroy(&fas->rscale);

1019:   fas->rscale = rscale;
1020:   return(0);
1021: }

1023: /*@
1024:    SNESFASGetSmoother - Gets the default smoother on a level.

1026:    Input Parameters:
1027: +  snes   - the multigrid context
1028: -  level  - the level (0 is coarsest) to supply

1030:    Output Parameters:
1031:    smooth  - the smoother

1033:    Level: advanced

1035: .keywords: FAS, MG, get, multigrid, smoother, level

1037: .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1038: @*/
1039: PetscErrorCode SNESFASGetSmoother(SNES snes, PetscInt level, SNES *smooth)
1040: {
1041:   SNES_FAS       *fas;
1043:   SNES           levelsnes;

1046:   SNESFASGetCycleSNES(snes, level, &levelsnes);
1047:   fas  = (SNES_FAS*)levelsnes->data;
1048:   if (!fas->smoothd) {
1049:     SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);
1050:   }
1051:   *smooth = fas->smoothd;
1052:   return(0);
1053: }

1055: /*@
1056:    SNESFASGetSmootherDown - Gets the downsmoother on a level.

1058:    Input Parameters:
1059: +  snes   - the multigrid context
1060: -  level  - the level (0 is coarsest) to supply

1062:    Output Parameters:
1063:    smooth  - the smoother

1065:    Level: advanced

1067: .keywords: FAS, MG, get, multigrid, smoother, level

1069: .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1070: @*/
1071: PetscErrorCode SNESFASGetSmootherDown(SNES snes, PetscInt level, SNES *smooth)
1072: {
1073:   SNES_FAS       *fas;
1075:   SNES           levelsnes;

1078:   SNESFASGetCycleSNES(snes, level, &levelsnes);
1079:   fas  = (SNES_FAS*)levelsnes->data;
1080:   /* if the user chooses to differentiate smoothers, create them both at this point */
1081:   if (!fas->smoothd) {
1082:     SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);
1083:   }
1084:   if (!fas->smoothu) {
1085:     SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu);
1086:   }
1087:   *smooth = fas->smoothd;
1088:   return(0);
1089: }

1091: /*@
1092:    SNESFASGetSmootherUp - Gets the upsmoother on a level.

1094:    Input Parameters:
1095: +  snes   - the multigrid context
1096: -  level  - the level (0 is coarsest)

1098:    Output Parameters:
1099:    smooth  - the smoother

1101:    Level: advanced

1103: .keywords: FAS, MG, get, multigrid, smoother, level

1105: .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1106: @*/
1107: PetscErrorCode SNESFASGetSmootherUp(SNES snes, PetscInt level, SNES *smooth)
1108: {
1109:   SNES_FAS       *fas;
1111:   SNES           levelsnes;

1114:   SNESFASGetCycleSNES(snes, level, &levelsnes);
1115:   fas  = (SNES_FAS*)levelsnes->data;
1116:   /* if the user chooses to differentiate smoothers, create them both at this point */
1117:   if (!fas->smoothd) {
1118:     SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);
1119:   }
1120:   if (!fas->smoothu) {
1121:     SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu);
1122:   }
1123:   *smooth = fas->smoothu;
1124:   return(0);
1125: }

1127: /*@
1128:   SNESFASGetCoarseSolve - Gets the coarsest solver.

1130:   Input Parameters:
1131: . snes - the multigrid context

1133:   Output Parameters:
1134: . coarse - the coarse-level solver

1136:   Level: advanced

1138: .keywords: FAS, MG, get, multigrid, solver, coarse
1139: .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1140: @*/
1141: PetscErrorCode SNESFASGetCoarseSolve(SNES snes, SNES *coarse)
1142: {
1143:   SNES_FAS       *fas;
1145:   SNES           levelsnes;

1148:   SNESFASGetCycleSNES(snes, 0, &levelsnes);
1149:   fas  = (SNES_FAS*)levelsnes->data;
1150:   /* if the user chooses to differentiate smoothers, create them both at this point */
1151:   if (!fas->smoothd) {
1152:     SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);
1153:   }
1154:   *coarse = fas->smoothd;
1155:   return(0);
1156: }

1158: /*@
1159:    SNESFASFullSetDownSweep - Smooth during the initial downsweep for SNESFAS

1161:    Logically Collective on SNES

1163:    Input Parameters:
1164: +  snes - the multigrid context
1165: -  swp - whether to downsweep or not

1167:    Options Database Key:
1168: .  -snes_fas_full_downsweep - Sets number of pre-smoothing steps

1170:    Level: advanced

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

1174: .seealso: SNESFASSetNumberSmoothUp()
1175: @*/
1176: PetscErrorCode SNESFASFullSetDownSweep(SNES snes,PetscBool swp)
1177: {
1178:   SNES_FAS       *fas =  (SNES_FAS*)snes->data;
1179:   PetscErrorCode 0;

1182:   fas->full_downsweep = swp;
1183:   if (fas->next) {
1184:     SNESFASFullSetDownSweep(fas->next,swp);
1185:   }
1186:   return(0);
1187: }