Actual source code: fasfunc.c

petsc-3.3-p7 2013-05-11
  1: #include "../src/snes/impls/fas/fasimpls.h" /*I  "petscsnesfas.h"  I*/


  4: extern PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES, SNES *);

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

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

 13: Logically Collective

 15: Input Parameters:
 16: + snes  - FAS context
 17: - fastype  - SNES_FAS_ADDITIVE or SNES_FAS_MULTIPLICATIVE

 19: Level: intermediate

 21: .seealso: PCMGSetType()
 22: @*/
 23: PetscErrorCode  SNESFASSetType(SNES snes,SNESFASType fastype)
 24: {
 25:   SNES_FAS                   *fas = (SNES_FAS*)snes->data;
 26:   PetscErrorCode             ierr;
 30:   fas->fastype = fastype;
 31:   if (fas->next) {SNESFASSetType(fas->next, fastype);}
 32:   return(0);
 33: }


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

 41: Logically Collective

 43: Input Parameters:
 44: . snes - FAS context

 46: Output Parameters:
 47: . fastype - SNES_FAS_ADDITIVE or SNES_FAS_MULTIPLICATIVE

 49: Level: intermediate

 51: .seealso: PCMGSetType()
 52: @*/
 53: PetscErrorCode  SNESFASGetType(SNES snes,SNESFASType *fastype)
 54: {
 55:   SNES_FAS                   *fas = (SNES_FAS*)snes->data;

 60:   *fastype = fas->fastype;
 61:   return(0);
 62: }

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

 70:    Input Parameters:
 71: +  snes   - the snes context
 72: .  levels - the number of levels
 73: -  comms  - optional communicators for each level; this is to allow solving the coarser
 74:             problems on smaller sets of processors. Use PETSC_NULL_OBJECT for default in
 75:             Fortran.

 77:    Level: intermediate

 79:    Notes:
 80:      If the number of levels is one then the multigrid uses the -fas_levels prefix
 81:   for setting the level options rather than the -fas_coarse prefix.

 83: .keywords: FAS, MG, set, levels, multigrid

 85: .seealso: SNESFASGetLevels()
 86: @*/
 87: PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) {
 89:   PetscInt i;
 90:   const char     *optionsprefix;
 91:   char           tprefix[128];
 92:   SNES_FAS * fas = (SNES_FAS *)snes->data;
 93:   SNES prevsnes;
 94:   MPI_Comm comm;
 96:   comm = ((PetscObject)snes)->comm;
 97:   if (levels == fas->levels) {
 98:     if (!comms)
 99:       return(0);
100:   }
101:   /* user has changed the number of levels; reset */
102:   SNESReset(snes);
103:   /* destroy any coarser levels if necessary */
104:   if (fas->next) SNESDestroy(&fas->next);
105:   fas->next = PETSC_NULL;
106:   fas->previous = PETSC_NULL;
107:   prevsnes = snes;
108:   /* setup the finest level */
109:   SNESGetOptionsPrefix(snes, &optionsprefix);
110:   for (i = levels - 1; i >= 0; i--) {
111:     if (comms) comm = comms[i];
112:     fas->level = i;
113:     fas->levels = levels;
114:     fas->fine = snes;
115:     fas->next = PETSC_NULL;
116:     if (i > 0) {
117:       SNESCreate(comm, &fas->next);
118:       sprintf(tprefix,"fas_levels_%d_cycle_",(int)fas->level);
119:       SNESAppendOptionsPrefix(fas->next,tprefix);
120:       SNESAppendOptionsPrefix(fas->next,optionsprefix);
121:       SNESSetType(fas->next, SNESFAS);
122:       SNESSetTolerances(fas->next, fas->next->abstol, fas->next->rtol, fas->next->stol, fas->n_cycles, fas->next->max_funcs);
123:       PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);
124:       ((SNES_FAS *)fas->next->data)->previous = prevsnes;
125:       prevsnes = fas->next;
126:       fas = (SNES_FAS *)prevsnes->data;
127:     }
128:   }
129:   return(0);
130: }


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

138:    Input Parameter:
139: .  snes - the nonlinear solver context

141:    Output parameter:
142: .  levels - the number of levels

144:    Level: advanced

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

148: .seealso: SNESFASSetLevels(), PCMGGetLevels()
149: @*/
150: PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) {
151:   SNES_FAS * fas = (SNES_FAS *)snes->data;
153:   *levels = fas->levels;
154:   return(0);
155: }


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

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

169:    Level: advanced

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

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

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

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

195: /*@
196:    SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to
197:    use on all levels.

199:    Logically Collective on SNES

201:    Input Parameters:
202: +  snes - the multigrid context
203: -  n    - the number of smoothing steps

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

208:    Level: advanced

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

212: .seealso: SNESFASSetNumberSmoothDown()
213: @*/
214: PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) {
215:   SNES_FAS * fas =  (SNES_FAS *)snes->data;
216:   PetscErrorCode 0;
218:   fas->max_up_it = n;
219:   if (!fas->smoothu) {
220:     SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);
221:   }
222:   SNESSetTolerances(fas->smoothu, fas->smoothu->abstol, fas->smoothu->rtol, fas->smoothu->stol, n, fas->smoothu->max_funcs);
223:   if (fas->next) {
224:     SNESFASSetNumberSmoothUp(fas->next, n);
225:   }
226:   return(0);
227: }

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

235:    Logically Collective on SNES

237:    Input Parameters:
238: +  snes - the multigrid context
239: -  n    - the number of smoothing steps

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

244:    Level: advanced

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

248: .seealso: SNESFASSetNumberSmoothUp()
249: @*/
250: PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) {
251:   SNES_FAS * fas =  (SNES_FAS *)snes->data;
252:   PetscErrorCode 0;
254:   if (!fas->smoothu) {
255:     SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);
256:   }
257:   if (!fas->smoothd) {
258:     SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);
259:   }
260:   SNESSetTolerances(fas->smoothd, fas->smoothd->abstol, fas->smoothd->rtol, fas->smoothd->stol, n, fas->smoothd->max_funcs);
261:   fas->max_down_it = n;
262:   if (fas->next) {
263:     SNESFASSetNumberSmoothDown(fas->next, n);
264:   }
265:   return(0);
266: }


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

275:    Logically Collective on SNES

277:    Input Parameters:
278: +  snes   - the multigrid context
279: -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle

281:    Options Database Key:
282: $  -snes_fas_cycles 1 or 2

284:    Level: advanced

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

288: .seealso: SNESFASSetCyclesOnLevel()
289: @*/
290: PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) {
291:   SNES_FAS * fas = (SNES_FAS *)snes->data;
293:   PetscBool isFine;
295:   SNESFASCycleIsFine(snes, &isFine);
296:   fas->n_cycles = cycles;
297:   if (!isFine)
298:     SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);
299:   if (fas->next) {
300:     SNESFASSetCycles(fas->next, cycles);
301:   }
302:   return(0);
303: }


308: /*@
309:    SNESFASSetMonitor - Sets the method-specific cycle monitoring

311:    Logically Collective on SNES

313:    Input Parameters:
314: +  snes   - the FAS context
315: -  flg    - monitor or not

317:    Level: advanced

319: .keywords: FAS, monitor

321: .seealso: SNESFASSetCyclesOnLevel()
322: @*/
323: PetscErrorCode SNESFASSetMonitor(SNES snes, PetscBool flg) {
324:   SNES_FAS       *fas = (SNES_FAS *)snes->data;
326:   PetscBool      isFine;
327:   PetscInt       i, levels = fas->levels;
328:   SNES           levelsnes;
330:   SNESFASCycleIsFine(snes, &isFine);
331:   if (isFine) {
332:     for (i = 0; i < levels; i++) {
333:       SNESFASGetCycleSNES(snes, i, &levelsnes);
334:       fas = (SNES_FAS *)levelsnes->data;
335:       if (flg) {
336:         fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)levelsnes)->comm);
337:         /* set the monitors for the upsmoother and downsmoother */
338:         SNESMonitorCancel(levelsnes);
339:         SNESMonitorSet(levelsnes,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);
340:       } else {
341:         /* unset the monitors on the coarse levels */
342:         if (i != fas->levels - 1) {
343:           SNESMonitorCancel(levelsnes);
344:         }
345:       }
346:     }
347:   }
348:   return(0);
349: }

353: /*
354: Creates the default smoother type.

356: This is SNESNRICHARDSON on each fine level and SNESLS on the coarse level.

358:  */
359: PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES snes, SNES *smooth) {
360:   SNES_FAS       *fas;
361:   const char     *optionsprefix;
362:   char           tprefix[128];
364:   SNES           nsmooth;
367:   fas = (SNES_FAS*)snes->data;
368:   SNESGetOptionsPrefix(fas->fine, &optionsprefix);
369:   /* create the default smoother */
370:   SNESCreate(((PetscObject)snes)->comm, &nsmooth);
371:   if (fas->level == 0) {
372:     sprintf(tprefix,"fas_coarse_");
373:     SNESAppendOptionsPrefix(nsmooth, optionsprefix);
374:     SNESAppendOptionsPrefix(nsmooth, tprefix);
375:     SNESSetType(nsmooth, SNESLS);
376:     SNESSetTolerances(nsmooth, nsmooth->abstol, nsmooth->rtol, nsmooth->stol, nsmooth->max_its, nsmooth->max_funcs);
377:   } else {
378:     sprintf(tprefix,"fas_levels_%d_",(int)fas->level);
379:     SNESAppendOptionsPrefix(nsmooth, optionsprefix);
380:     SNESAppendOptionsPrefix(nsmooth, tprefix);
381:     SNESSetType(nsmooth, SNESNRICHARDSON);
382:     SNESSetTolerances(nsmooth, 0.0, 0.0, 0.0, fas->max_down_it, nsmooth->max_funcs);
383:   }
384:   PetscObjectIncrementTabLevel((PetscObject)nsmooth, (PetscObject)snes, 1);
385:   PetscLogObjectParent(snes,nsmooth);
386:   PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)nsmooth);
387:   *smooth = nsmooth;
388:   return(0);
389: }

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

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

398:    Logically Collective on SNES

400:    Input Parameters:
401: +  snes   - the multigrid context
402: .  level  - the level to set the number of cycles on
403: -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle

405:    Level: advanced

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

409: .seealso: SNESFASSetCycles()
410: @*/
411: PetscErrorCode SNESFASCycleSetCycles(SNES snes, PetscInt cycles) {
412:   SNES_FAS * fas =  (SNES_FAS *)snes->data;
415:   fas->n_cycles = cycles;
416:   SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);
417:   return(0);
418: }


423: /*@
424:    SNESFASCycleGetSmoother - Gets the smoother on a particular cycle level.

426:    Logically Collective on SNES

428:    Input Parameters:
429: .  snes   - the multigrid context

431:    Output Parameters:
432: .  smooth - the smoother

434:    Level: advanced

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

438: .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmootherDown()
439: @*/
440: PetscErrorCode SNESFASCycleGetSmoother(SNES snes, SNES *smooth)
441: {
442:   SNES_FAS       *fas;
445:   fas = (SNES_FAS*)snes->data;
446:   *smooth = fas->smoothd;
447:   return(0);
448: }
451: /*@
452:    SNESFASCycleGetSmootherUp - Gets the up smoother on a particular cycle level.

454:    Logically Collective on SNES

456:    Input Parameters:
457: .  snes   - the multigrid context

459:    Output Parameters:
460: .  smoothu - the smoother

462:    Notes:
463:    Returns the downsmoother if no up smoother is available.  This enables transparent
464:    default behavior in the process of the solve.

466:    Level: advanced

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

470: .seealso: SNESFASCycleGetSmoother(), SNESFASCycleGetSmootherDown()
471: @*/
472: PetscErrorCode SNESFASCycleGetSmootherUp(SNES snes, SNES *smoothu)
473: {
474:   SNES_FAS       *fas;
477:   fas = (SNES_FAS*)snes->data;
478:   if (!fas->smoothu) {
479:     *smoothu = fas->smoothd;
480:   } else {
481:     *smoothu = fas->smoothu;
482:   }
483:   return(0);
484: }

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

491:    Logically Collective on SNES

493:    Input Parameters:
494: .  snes   - the multigrid context

496:    Output Parameters:
497: .  smoothd - the smoother

499:    Level: advanced

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

503: .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
504: @*/
505: PetscErrorCode SNESFASCycleGetSmootherDown(SNES snes, SNES *smoothd)
506: {
507:   SNES_FAS       *fas;
510:   fas = (SNES_FAS*)snes->data;
511:   *smoothd = fas->smoothd;
512:   return(0);
513: }


518: /*@
519:    SNESFASCycleGetCorrection - Gets the coarse correction FAS context for this level

521:    Logically Collective on SNES

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

526:    Output Parameters:
527: .  correction - the coarse correction on this level

529:    Notes:
530:    Returns PETSC_NULL on the coarsest level.

532:    Level: advanced

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

536: .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
537: @*/
538: PetscErrorCode SNESFASCycleGetCorrection(SNES snes, SNES *correction)
539: {
540:   SNES_FAS       *fas;
543:   fas = (SNES_FAS*)snes->data;
544:   *correction = fas->next;
545:   return(0);
546: }

550: /*@
551:    SNESFASCycleGetInterpolation - Gets the interpolation on this level

553:    Logically Collective on SNES

555:    Input Parameters:
556: .  snes   - the multigrid context

558:    Output Parameters:
559: .  mat    - the interpolation operator on this level

561:    Level: developer

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

565: .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
566: @*/
567: PetscErrorCode SNESFASCycleGetInterpolation(SNES snes, Mat *mat)
568: {
569:   SNES_FAS       *fas;
572:   fas = (SNES_FAS*)snes->data;
573:   *mat = fas->interpolate;
574:   return(0);
575: }


580: /*@
581:    SNESFASCycleGetRestriction - Gets the restriction on this level

583:    Logically Collective on SNES

585:    Input Parameters:
586: .  snes   - the multigrid context

588:    Output Parameters:
589: .  mat    - the restriction operator on this level

591:    Level: developer

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

595: .seealso: SNESFASGetRestriction(), SNESFASCycleGetInterpolation()
596: @*/
597: PetscErrorCode SNESFASCycleGetRestriction(SNES snes, Mat *mat)
598: {
599:   SNES_FAS       *fas;
602:   fas = (SNES_FAS*)snes->data;
603:   *mat = fas->restrct;
604:   return(0);
605: }


610: /*@
611:    SNESFASCycleGetInjection - Gets the injection on this level

613:    Logically Collective on SNES

615:    Input Parameters:
616: .  snes   - the multigrid context

618:    Output Parameters:
619: .  mat    - the restriction operator on this level

621:    Level: developer

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

625: .seealso: SNESFASGetInjection(), SNESFASCycleGetRestriction()
626: @*/
627: PetscErrorCode SNESFASCycleGetInjection(SNES snes, Mat *mat)
628: {
629:   SNES_FAS       *fas;
632:   fas = (SNES_FAS*)snes->data;
633:   *mat = fas->inject;
634:   return(0);
635: }

639: /*@
640:    SNESFASCycleGetRScale - Gets the injection on this level

642:    Logically Collective on SNES

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

647:    Output Parameters:
648: .  mat    - the restriction operator on this level

650:    Level: developer

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

654: .seealso: SNESFASCycleGetRestriction(), SNESFASGetRScale()
655: @*/
656: PetscErrorCode SNESFASCycleGetRScale(SNES snes, Vec *vec)
657: {
658:   SNES_FAS       *fas;
661:   fas = (SNES_FAS*)snes->data;
662:   *vec = fas->rscale;
663:   return(0);
664: }

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

671:    Logically Collective on SNES

673:    Input Parameters:
674: .  snes   - the FAS context

676:    Output Parameters:
677: .  flg - indicates if this is the fine level or not

679:    Level: advanced

681: .keywords: SNES, FAS

683: .seealso: SNESFASSetLevels()
684: @*/
685: PetscErrorCode SNESFASCycleIsFine(SNES snes, PetscBool *flg)
686: {
687:   SNES_FAS       *fas;
690:   fas = (SNES_FAS*)snes->data;
691:   if (fas->level == fas->levels - 1) {
692:     *flg = PETSC_TRUE;
693:   } else {
694:     *flg = PETSC_FALSE;
695:   }
696:   return(0);
697: }

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

703: /*@
704:    SNESFASSetInterpolation - Sets the function to be used to calculate the
705:    interpolation from l-1 to the lth level

707:    Input Parameters:
708: +  snes      - the multigrid context
709: .  mat       - the interpolation operator
710: -  level     - the level (0 is coarsest) to supply [do not supply 0]

712:    Level: advanced

714:    Notes:
715:           Usually this is the same matrix used also to set the restriction
716:     for the same level.

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

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

723: .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRScale()
724: @*/
725: PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) {
726:   SNES_FAS       *fas =  (SNES_FAS *)snes->data;
728:   SNES           levelsnes;
730:   SNESFASGetCycleSNES(snes, level, &levelsnes);
731:   fas = (SNES_FAS *)levelsnes->data;
732:   PetscObjectReference((PetscObject)mat);
733:   MatDestroy(&fas->interpolate);
734:   fas->interpolate = mat;
735:   return(0);
736: }

740: /*@
741:    SNESFASGetInterpolation - Gets the matrix used to calculate the
742:    interpolation from l-1 to the lth level

744:    Input Parameters:
745: +  snes      - the multigrid context
746: -  level     - the level (0 is coarsest) to supply [do not supply 0]

748:    Output Parameters:
749: .  mat       - the interpolation operator

751:    Level: advanced

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

755: .seealso: SNESFASSetInterpolation(), SNESFASGetInjection(), SNESFASGetRestriction(), SNESFASGetRScale()
756: @*/
757: PetscErrorCode SNESFASGetInterpolation(SNES snes, PetscInt level, Mat *mat) {
758:   SNES_FAS       *fas =  (SNES_FAS *)snes->data;
760:   SNES           levelsnes;
762:   SNESFASGetCycleSNES(snes, level, &levelsnes);
763:   fas = (SNES_FAS *)levelsnes->data;
764:   *mat = fas->interpolate;
765:   return(0);
766: }

770: /*@
771:    SNESFASSetRestriction - Sets the function to be used to restrict the defect
772:    from level l to l-1.

774:    Input Parameters:
775: +  snes  - the multigrid context
776: .  mat   - the restriction matrix
777: -  level - the level (0 is coarsest) to supply [Do not supply 0]

779:    Level: advanced

781:    Notes:
782:           Usually this is the same matrix used also to set the interpolation
783:     for the same level.

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

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

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

793: .seealso: SNESFASSetInterpolation(), SNESFASSetInjection()
794: @*/
795: PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) {
796:   SNES_FAS * fas =  (SNES_FAS *)snes->data;
798:   SNES levelsnes;
800:   SNESFASGetCycleSNES(snes, level, &levelsnes);
801:   fas = (SNES_FAS *)levelsnes->data;
802:   PetscObjectReference((PetscObject)mat);
803:   MatDestroy(&fas->restrct);
804:   fas->restrct = mat;
805:   return(0);
806: }

810: /*@
811:    SNESFASGetRestriction - Gets the matrix used to calculate the
812:    restriction from l to the l-1th level

814:    Input Parameters:
815: +  snes      - the multigrid context
816: -  level     - the level (0 is coarsest) to supply [do not supply 0]

818:    Output Parameters:
819: .  mat       - the interpolation operator

821:    Level: advanced

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

825: .seealso: SNESFASSetRestriction(), SNESFASGetInjection(), SNESFASGetInterpolation(), SNESFASGetRScale()
826: @*/
827: PetscErrorCode SNESFASGetRestriction(SNES snes, PetscInt level, Mat *mat) {
828:   SNES_FAS       *fas =  (SNES_FAS *)snes->data;
830:   SNES           levelsnes;
832:   SNESFASGetCycleSNES(snes, level, &levelsnes);
833:   fas = (SNES_FAS *)levelsnes->data;
834:   *mat = fas->restrct;
835:   return(0);
836: }


841: /*@
842:    SNESFASSetInjection - Sets the function to be used to inject the solution
843:    from level l to l-1.

845:    Input Parameters:
846:  +  snes  - the multigrid context
847: .  mat   - the restriction matrix
848: -  level - the level (0 is coarsest) to supply [Do not supply 0]

850:    Level: advanced

852:    Notes:
853:          If you do not set this, the restriction and rscale is used to
854:    project the solution instead.

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

858: .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction()
859: @*/
860: PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) {
861:   SNES_FAS       *fas =  (SNES_FAS *)snes->data;
863:   SNES           levelsnes;
865:   SNESFASGetCycleSNES(snes, level, &levelsnes);
866:   fas = (SNES_FAS *)levelsnes->data;
867:   PetscObjectReference((PetscObject)mat);
868:   MatDestroy(&fas->inject);
869:   fas->inject = mat;
870:   return(0);
871: }


876: /*@
877:    SNESFASGetInjection - Gets the matrix used to calculate the
878:    injection from l-1 to the lth level

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

884:    Output Parameters:
885: .  mat       - the injection operator

887:    Level: advanced

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

891: .seealso: SNESFASSetInjection(), SNESFASGetRestriction(), SNESFASGetInterpolation(), SNESFASGetRScale()
892: @*/
893: PetscErrorCode SNESFASGetInjection(SNES snes, PetscInt level, Mat *mat) {
894:   SNES_FAS       *fas =  (SNES_FAS *)snes->data;
896:   SNES           levelsnes;
898:   SNESFASGetCycleSNES(snes, level, &levelsnes);
899:   fas = (SNES_FAS *)levelsnes->data;
900:   *mat = fas->inject;
901:   return(0);
902: }

906: /*@
907:    SNESFASSetRScale - Sets the scaling factor of the restriction
908:    operator from level l to l-1.

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

915:    Level: advanced

917:    Notes:
918:          This is only used in the case that the injection is not set.

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

922: .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
923: @*/
924: PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) {
925:   SNES_FAS * fas;
927:   SNES levelsnes;
929:   SNESFASGetCycleSNES(snes, level, &levelsnes);
930:   fas = (SNES_FAS *)levelsnes->data;
931:   PetscObjectReference((PetscObject)rscale);
932:   VecDestroy(&fas->rscale);
933:   fas->rscale = rscale;
934:   return(0);
935: }

939: /*@
940:    SNESFASGetSmoother - Gets the default smoother on a level.

942:    Input Parameters:
943: +  snes   - the multigrid context
944: -  level  - the level (0 is coarsest) to supply

946:    Output Parameters:
947:    smooth  - the smoother

949:    Level: advanced

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

953: .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
954: @*/
955: PetscErrorCode SNESFASGetSmoother(SNES snes, PetscInt level, SNES *smooth) {
956:   SNES_FAS * fas;
958:   SNES levelsnes;
960:   SNESFASGetCycleSNES(snes, level, &levelsnes);
961:   fas = (SNES_FAS *)levelsnes->data;
962:   if (!fas->smoothd) {
963:     SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);
964:   }
965:   *smooth = fas->smoothd;
966:   return(0);
967: }

971: /*@
972:    SNESFASGetSmootherDown - Gets the downsmoother on a level.

974:    Input Parameters:
975: +  snes   - the multigrid context
976: -  level  - the level (0 is coarsest) to supply

978:    Output Parameters:
979:    smooth  - the smoother

981:    Level: advanced

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

985: .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
986: @*/
987: PetscErrorCode SNESFASGetSmootherDown(SNES snes, PetscInt level, SNES *smooth) {
988:   SNES_FAS * fas;
990:   SNES levelsnes;
992:   SNESFASGetCycleSNES(snes, level, &levelsnes);
993:   fas = (SNES_FAS *)levelsnes->data;
994:   /* if the user chooses to differentiate smoothers, create them both at this point */
995:   if (!fas->smoothd) {
996:     SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);
997:   }
998:   if (!fas->smoothu) {
999:     SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);
1000:   }
1001:   *smooth = fas->smoothd;
1002:   return(0);
1003: }

1007: /*@
1008:    SNESFASGetSmootherUp - Gets the upsmoother on a level.

1010:    Input Parameters:
1011: +  snes   - the multigrid context
1012: -  level  - the level (0 is coarsest)

1014:    Output Parameters:
1015:    smooth  - the smoother

1017:    Level: advanced

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

1021: .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1022: @*/
1023: PetscErrorCode SNESFASGetSmootherUp(SNES snes, PetscInt level, SNES *smooth) {
1024:   SNES_FAS * fas;
1026:   SNES levelsnes;
1028:   SNESFASGetCycleSNES(snes, level, &levelsnes);
1029:   fas = (SNES_FAS *)levelsnes->data;
1030:   /* if the user chooses to differentiate smoothers, create them both at this point */
1031:   if (!fas->smoothd) {
1032:     SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);
1033:   }
1034:   if (!fas->smoothu) {
1035:     SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);
1036:   }
1037:   *smooth = fas->smoothu;
1038:   return(0);
1039: }

1043: /*@
1044:    SNESFASGetCoarseSolve - Gets the coarsest solver.

1046:    Input Parameters:
1047: +  snes   - the multigrid context

1049:    Output Parameters:
1050:    solve  - the coarse-level solver

1052:    Level: advanced

1054: .keywords: FAS, MG, get, multigrid, solver, coarse

1056: .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1057: @*/
1058: PetscErrorCode SNESFASGetCoarseSolve(SNES snes, SNES *smooth) {
1059:   SNES_FAS * fas;
1061:   SNES levelsnes;
1063:   SNESFASGetCycleSNES(snes, 0, &levelsnes);
1064:   fas = (SNES_FAS *)levelsnes->data;
1065:   /* if the user chooses to differentiate smoothers, create them both at this point */
1066:   if (!fas->smoothd) {
1067:     SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);
1068:   }
1069:   *smooth = fas->smoothd;
1070:   return(0);
1071: }