Actual source code: fasfunc.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  1: #include <../src/snes/impls/fas/fasimpls.h> /*I  "petscsnes.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, SNES_FAS_MULTIPLICATIVE, SNES_FAS_FULL, or SNES_FAS_KASKADE

 19: Level: intermediate

 21: .seealso: PCMGSetType()
 22: @*/
 23: PetscErrorCode  SNESFASSetType(SNES snes,SNESFASType fastype)
 24: {
 25:   SNES_FAS       *fas = (SNES_FAS*)snes->data;

 31:   fas->fastype = fastype;
 32:   if (fas->next) {
 33:     SNESFASSetType(fas->next, fastype);
 34:   }
 35:   return(0);
 36: }


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

 44: Logically Collective

 46: Input Parameters:
 47: . snes - FAS context

 49: Output Parameters:
 50: . fastype - SNES_FAS_ADDITIVE or SNES_FAS_MULTIPLICATIVE

 52: Level: intermediate

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

 63:   *fastype = fas->fastype;
 64:   return(0);
 65: }

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

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

 80:    Level: intermediate

 82:    Notes:
 83:      If the number of levels is one then the multigrid uses the -fas_levels prefix
 84:   for setting the level options rather than the -fas_coarse prefix.

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

 88: .seealso: SNESFASGetLevels()
 89: @*/
 90: PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms)
 91: {
 93:   PetscInt       i;
 94:   const char     *optionsprefix;
 95:   char           tprefix[128];
 96:   SNES_FAS       *fas = (SNES_FAS*)snes->data;
 97:   SNES           prevsnes;
 98:   MPI_Comm       comm;

101:   PetscObjectGetComm((PetscObject)snes,&comm);
102:   if (levels == fas->levels) {
103:     if (!comms) return(0);
104:   }
105:   /* user has changed the number of levels; reset */
106:   SNESReset(snes);
107:   /* destroy any coarser levels if necessary */
108:   if (fas->next) SNESDestroy(&fas->next);
109:   fas->next     = NULL;
110:   fas->previous = NULL;
111:   prevsnes      = snes;
112:   /* setup the finest level */
113:   SNESGetOptionsPrefix(snes, &optionsprefix);
114:   for (i = levels - 1; i >= 0; i--) {
115:     if (comms) comm = comms[i];
116:     fas->level  = i;
117:     fas->levels = levels;
118:     fas->fine   = snes;
119:     fas->next   = NULL;
120:     if (i > 0) {
121:       SNESCreate(comm, &fas->next);
122:       SNESGetOptionsPrefix(fas->fine, &optionsprefix);
123:       sprintf(tprefix,"fas_levels_%d_cycle_",(int)fas->level);
124:       SNESAppendOptionsPrefix(fas->next,optionsprefix);
125:       SNESAppendOptionsPrefix(fas->next,tprefix);
126:       SNESSetType(fas->next, SNESFAS);
127:       SNESSetTolerances(fas->next, fas->next->abstol, fas->next->rtol, fas->next->stol, fas->n_cycles, fas->next->max_funcs);
128:       PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);

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

132:       prevsnes = fas->next;
133:       fas      = (SNES_FAS*)prevsnes->data;
134:     }
135:   }
136:   return(0);
137: }


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

145:    Input Parameter:
146: .  snes - the nonlinear solver context

148:    Output parameter:
149: .  levels - the number of levels

151:    Level: advanced

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

155: .seealso: SNESFASSetLevels(), PCMGGetLevels()
156: @*/
157: PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt *levels)
158: {
159:   SNES_FAS * fas = (SNES_FAS*)snes->data;

162:   *levels = fas->levels;
163:   return(0);
164: }


169: /*@
170:    SNESFASGetCycleSNES - Gets the SNES corresponding to a particular
171:    level of the FAS hierarchy.

173:    Input Parameters:
174: +  snes    - the multigrid context
175:    level   - the level to get
176: -  lsnes   - whether to use the nonlinear smoother or not

178:    Level: advanced

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

182: .seealso: SNESFASSetLevels(), SNESFASGetLevels()
183: @*/
184: PetscErrorCode SNESFASGetCycleSNES(SNES snes,PetscInt level,SNES *lsnes)
185: {
186:   SNES_FAS *fas = (SNES_FAS*)snes->data;
187:   PetscInt i;

190:   if (level > fas->levels-1) SETERRQ2(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Requested level %D from SNESFAS containing %D levels",level,fas->levels);
191:   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);

193:   *lsnes = snes;
194:   for (i = fas->level; i > level; i--) {
195:     *lsnes = fas->next;
196:     fas    = (SNES_FAS*)(*lsnes)->data;
197:   }
198:   if (fas->level != level) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt");
199:   return(0);
200: }

204: /*@
205:    SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to
206:    use on all levels.

208:    Logically Collective on SNES

210:    Input Parameters:
211: +  snes - the multigrid context
212: -  n    - the number of smoothing steps

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

217:    Level: advanced

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

221: .seealso: SNESFASSetNumberSmoothDown()
222: @*/
223: PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n)
224: {
225:   SNES_FAS       *fas =  (SNES_FAS*)snes->data;

229:   fas->max_up_it = n;
230:   if (!fas->smoothu && fas->level != 0) {
231:     SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);
232:   }
233:   if (fas->smoothu) {
234:     SNESSetTolerances(fas->smoothu, fas->smoothu->abstol, fas->smoothu->rtol, fas->smoothu->stol, n, fas->smoothu->max_funcs);
235:   }
236:   if (fas->next) {
237:     SNESFASSetNumberSmoothUp(fas->next, n);
238:   }
239:   return(0);
240: }

244: /*@
245:    SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to
246:    use on all levels.

248:    Logically Collective on SNES

250:    Input Parameters:
251: +  snes - the multigrid context
252: -  n    - the number of smoothing steps

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

257:    Level: advanced

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

261: .seealso: SNESFASSetNumberSmoothUp()
262: @*/
263: PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n)
264: {
265:   SNES_FAS       *fas =  (SNES_FAS*)snes->data;
266:   PetscErrorCode 0;

269:   if (!fas->smoothd) {
270:     SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);
271:   }
272:   SNESSetTolerances(fas->smoothd, fas->smoothd->abstol, fas->smoothd->rtol, fas->smoothd->stol, n, fas->smoothd->max_funcs);

274:   fas->max_down_it = n;
275:   if (fas->next) {
276:     SNESFASSetNumberSmoothDown(fas->next, n);
277:   }
278:   return(0);
279: }


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

287:    Logically Collective on SNES

289:    Input Parameters:
290: +  snes - the multigrid context
291: -  n    - the number of smoothing steps

293:    Options Database Key:
294: .  -snes_fas_continuation - sets continuation to true

296:    Level: advanced

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

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

302: .seealso: SNESFAS
303: @*/
304: PetscErrorCode SNESFASSetContinuation(SNES snes,PetscBool continuation)
305: {
306:   const char     *optionsprefix;
307:   char           tprefix[128];
308:   SNES_FAS       *fas =  (SNES_FAS*)snes->data;
309:   PetscErrorCode 0;

312:   SNESGetOptionsPrefix(fas->fine, &optionsprefix);
313:   if (!fas->smoothu) {
314:     SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);
315:   }
316:   sprintf(tprefix,"fas_levels_continuation_");
317:   SNESSetOptionsPrefix(fas->smoothu, optionsprefix);
318:   SNESAppendOptionsPrefix(fas->smoothu, tprefix);
319:   SNESSetType(fas->smoothu,SNESNEWTONLS);
320:   SNESSetTolerances(fas->smoothu,fas->fine->abstol,fas->fine->rtol,fas->fine->stol,50,100);
321:   fas->continuation = continuation;
322:   if (fas->next) {
323:     SNESFASSetContinuation(fas->next,continuation);
324:   }
325:   return(0);
326: }


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

335:    Logically Collective on SNES

337:    Input Parameters:
338: +  snes   - the multigrid context
339: -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle

341:    Options Database Key:
342: .  -snes_fas_cycles 1 or 2

344:    Level: advanced

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

348: .seealso: SNESFASSetCyclesOnLevel()
349: @*/
350: PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles)
351: {
352:   SNES_FAS       *fas = (SNES_FAS*)snes->data;
354:   PetscBool      isFine;

357:   SNESFASCycleIsFine(snes, &isFine);

359:   fas->n_cycles = cycles;
360:   if (!isFine) {
361:     SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);
362:   }
363:   if (fas->next) {
364:     SNESFASSetCycles(fas->next, cycles);
365:   }
366:   return(0);
367: }


372: /*@
373:    SNESFASSetMonitor - Sets the method-specific cycle monitoring

375:    Logically Collective on SNES

377:    Input Parameters:
378: +  snes   - the FAS context
379: -  flg    - monitor or not

381:    Level: advanced

383: .keywords: FAS, monitor

385: .seealso: SNESFASSetCyclesOnLevel()
386: @*/
387: PetscErrorCode SNESFASSetMonitor(SNES snes, PetscBool flg)
388: {
389:   SNES_FAS       *fas = (SNES_FAS*)snes->data;
391:   PetscBool      isFine;
392:   PetscInt       i, levels = fas->levels;
393:   SNES           levelsnes;

396:   SNESFASCycleIsFine(snes, &isFine);
397:   if (isFine) {
398:     for (i = 0; i < levels; i++) {
399:       SNESFASGetCycleSNES(snes, i, &levelsnes);
400:       fas  = (SNES_FAS*)levelsnes->data;
401:       if (flg) {
402:         fas->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)levelsnes));
403:         /* set the monitors for the upsmoother and downsmoother */
404:         SNESMonitorCancel(levelsnes);
405:         SNESMonitorSet(levelsnes,SNESMonitorDefault,NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);
406:       } else if (i != fas->levels - 1) {
407:         /* unset the monitors on the coarse levels */
408:         SNESMonitorCancel(levelsnes);
409:       }
410:     }
411:   }
412:   return(0);
413: }

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

420:    Logically Collective on SNES

422:    Input Parameters:
423: +  snes   - the FAS context
424: -  flg    - monitor or not

426:    Level: advanced

428: .keywords: FAS, logging

430: .seealso: SNESFASSetMonitor()
431: @*/
432: PetscErrorCode SNESFASSetLog(SNES snes, PetscBool flg)
433: {
434:   SNES_FAS       *fas = (SNES_FAS*)snes->data;
436:   PetscBool      isFine;
437:   PetscInt       i, levels = fas->levels;
438:   SNES           levelsnes;
439:   char           eventname[128];

442:   SNESFASCycleIsFine(snes, &isFine);
443:   if (isFine) {
444:     for (i = 0; i < levels; i++) {
445:       SNESFASGetCycleSNES(snes, i, &levelsnes);
446:       fas  = (SNES_FAS*)levelsnes->data;
447:       if (flg) {
448:         sprintf(eventname,"FASSetup %d",(int)i);
449:         PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventsmoothsetup);
450:         sprintf(eventname,"FASSmooth %d",(int)i);
451:         PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventsmoothsolve);
452:         sprintf(eventname,"FASResid %d",(int)i);
453:         PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventresidual);
454:         if (i) {
455:           sprintf(eventname,"FASInterp %d",(int)i);
456:           PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventinterprestrict);
457:         }
458:       } else {
459:         fas->eventsmoothsetup    = 0;
460:         fas->eventsmoothsolve    = 0;
461:         fas->eventresidual       = 0;
462:         fas->eventinterprestrict = 0;
463:       }
464:     }
465:   }
466:   return(0);
467: }

471: /*
472: Creates the default smoother type.

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

476:  */
477: PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES snes, SNES *smooth)
478: {
479:   SNES_FAS       *fas;
480:   const char     *optionsprefix;
481:   char           tprefix[128];
483:   SNES           nsmooth;

487:   fas  = (SNES_FAS*)snes->data;
488:   SNESGetOptionsPrefix(fas->fine, &optionsprefix);
489:   /* create the default smoother */
490:   SNESCreate(PetscObjectComm((PetscObject)snes), &nsmooth);
491:   if (fas->level == 0) {
492:     sprintf(tprefix,"fas_coarse_");
493:     SNESAppendOptionsPrefix(nsmooth, optionsprefix);
494:     SNESAppendOptionsPrefix(nsmooth, tprefix);
495:     SNESSetType(nsmooth, SNESNEWTONLS);
496:     SNESSetTolerances(nsmooth, nsmooth->abstol, nsmooth->rtol, nsmooth->stol, nsmooth->max_its, nsmooth->max_funcs);
497:   } else {
498:     sprintf(tprefix,"fas_levels_%d_",(int)fas->level);
499:     SNESAppendOptionsPrefix(nsmooth, optionsprefix);
500:     SNESAppendOptionsPrefix(nsmooth, tprefix);
501:     SNESSetType(nsmooth, SNESNRICHARDSON);
502:     SNESSetTolerances(nsmooth, 0.0, 0.0, 0.0, fas->max_down_it, nsmooth->max_funcs);
503:   }
504:   PetscObjectIncrementTabLevel((PetscObject)nsmooth, (PetscObject)snes, 1);
505:   PetscLogObjectParent((PetscObject)snes,(PetscObject)nsmooth);
506:   PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)nsmooth);
507:   *smooth = nsmooth;
508:   return(0);
509: }

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

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

518:    Logically Collective on SNES

520:    Input Parameters:
521: +  snes   - the multigrid context
522: .  level  - the level to set the number of cycles on
523: -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle

525:    Level: advanced

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

529: .seealso: SNESFASSetCycles()
530: @*/
531: PetscErrorCode SNESFASCycleSetCycles(SNES snes, PetscInt cycles)
532: {
533:   SNES_FAS       *fas =  (SNES_FAS*)snes->data;

537:   fas->n_cycles = cycles;
538:   SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);
539:   return(0);
540: }


545: /*@
546:    SNESFASCycleGetSmoother - Gets the smoother on a particular cycle level.

548:    Logically Collective on SNES

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

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

556:    Level: advanced

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

560: .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmootherDown()
561: @*/
562: PetscErrorCode SNESFASCycleGetSmoother(SNES snes, SNES *smooth)
563: {
564:   SNES_FAS *fas;

568:   fas     = (SNES_FAS*)snes->data;
569:   *smooth = fas->smoothd;
570:   return(0);
571: }
574: /*@
575:    SNESFASCycleGetSmootherUp - Gets the up smoother on a particular cycle level.

577:    Logically Collective on SNES

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

582:    Output Parameters:
583: .  smoothu - the smoother

585:    Notes:
586:    Returns the downsmoother if no up smoother is available.  This enables transparent
587:    default behavior in the process of the solve.

589:    Level: advanced

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

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

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

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

612:    Logically Collective on SNES

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

617:    Output Parameters:
618: .  smoothd - the smoother

620:    Level: advanced

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

624: .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
625: @*/
626: PetscErrorCode SNESFASCycleGetSmootherDown(SNES snes, SNES *smoothd)
627: {
628:   SNES_FAS *fas;

632:   fas      = (SNES_FAS*)snes->data;
633:   *smoothd = fas->smoothd;
634:   return(0);
635: }


640: /*@
641:    SNESFASCycleGetCorrection - Gets the coarse correction FAS context for this level

643:    Logically Collective on SNES

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

648:    Output Parameters:
649: .  correction - the coarse correction on this level

651:    Notes:
652:    Returns NULL on the coarsest level.

654:    Level: advanced

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

658: .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
659: @*/
660: PetscErrorCode SNESFASCycleGetCorrection(SNES snes, SNES *correction)
661: {
662:   SNES_FAS *fas;

666:   fas         = (SNES_FAS*)snes->data;
667:   *correction = fas->next;
668:   return(0);
669: }

673: /*@
674:    SNESFASCycleGetInterpolation - Gets the interpolation on this level

676:    Logically Collective on SNES

678:    Input Parameters:
679: .  snes   - the multigrid context

681:    Output Parameters:
682: .  mat    - the interpolation operator on this level

684:    Level: developer

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

688: .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
689: @*/
690: PetscErrorCode SNESFASCycleGetInterpolation(SNES snes, Mat *mat)
691: {
692:   SNES_FAS *fas;

696:   fas  = (SNES_FAS*)snes->data;
697:   *mat = fas->interpolate;
698:   return(0);
699: }


704: /*@
705:    SNESFASCycleGetRestriction - Gets the restriction on this level

707:    Logically Collective on SNES

709:    Input Parameters:
710: .  snes   - the multigrid context

712:    Output Parameters:
713: .  mat    - the restriction operator on this level

715:    Level: developer

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

719: .seealso: SNESFASGetRestriction(), SNESFASCycleGetInterpolation()
720: @*/
721: PetscErrorCode SNESFASCycleGetRestriction(SNES snes, Mat *mat)
722: {
723:   SNES_FAS *fas;

727:   fas  = (SNES_FAS*)snes->data;
728:   *mat = fas->restrct;
729:   return(0);
730: }


735: /*@
736:    SNESFASCycleGetInjection - Gets the injection on this level

738:    Logically Collective on SNES

740:    Input Parameters:
741: .  snes   - the multigrid context

743:    Output Parameters:
744: .  mat    - the restriction operator on this level

746:    Level: developer

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

750: .seealso: SNESFASGetInjection(), SNESFASCycleGetRestriction()
751: @*/
752: PetscErrorCode SNESFASCycleGetInjection(SNES snes, Mat *mat)
753: {
754:   SNES_FAS *fas;

758:   fas  = (SNES_FAS*)snes->data;
759:   *mat = fas->inject;
760:   return(0);
761: }

765: /*@
766:    SNESFASCycleGetRScale - Gets the injection on this level

768:    Logically Collective on SNES

770:    Input Parameters:
771: .  snes   - the multigrid context

773:    Output Parameters:
774: .  mat    - the restriction operator on this level

776:    Level: developer

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

780: .seealso: SNESFASCycleGetRestriction(), SNESFASGetRScale()
781: @*/
782: PetscErrorCode SNESFASCycleGetRScale(SNES snes, Vec *vec)
783: {
784:   SNES_FAS *fas;

788:   fas  = (SNES_FAS*)snes->data;
789:   *vec = fas->rscale;
790:   return(0);
791: }

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

798:    Logically Collective on SNES

800:    Input Parameters:
801: .  snes   - the FAS context

803:    Output Parameters:
804: .  flg - indicates if this is the fine level or not

806:    Level: advanced

808: .keywords: SNES, FAS

810: .seealso: SNESFASSetLevels()
811: @*/
812: PetscErrorCode SNESFASCycleIsFine(SNES snes, PetscBool *flg)
813: {
814:   SNES_FAS *fas;

818:   fas = (SNES_FAS*)snes->data;
819:   if (fas->level == fas->levels - 1) *flg = PETSC_TRUE;
820:   else *flg = PETSC_FALSE;
821:   return(0);
822: }

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

828: /*@
829:    SNESFASSetInterpolation - Sets the function to be used to calculate the
830:    interpolation from l-1 to the lth level

832:    Input Parameters:
833: +  snes      - the multigrid context
834: .  mat       - the interpolation operator
835: -  level     - the level (0 is coarsest) to supply [do not supply 0]

837:    Level: advanced

839:    Notes:
840:           Usually this is the same matrix used also to set the restriction
841:     for the same level.

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

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

848: .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRScale()
849: @*/
850: PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat)
851: {
852:   SNES_FAS       *fas;
854:   SNES           levelsnes;

857:   SNESFASGetCycleSNES(snes, level, &levelsnes);
858:   fas  = (SNES_FAS*)levelsnes->data;
859:   PetscObjectReference((PetscObject)mat);
860:   MatDestroy(&fas->interpolate);

862:   fas->interpolate = mat;
863:   return(0);
864: }

868: /*@
869:    SNESFASGetInterpolation - Gets the matrix used to calculate the
870:    interpolation from l-1 to the lth level

872:    Input Parameters:
873: +  snes      - the multigrid context
874: -  level     - the level (0 is coarsest) to supply [do not supply 0]

876:    Output Parameters:
877: .  mat       - the interpolation operator

879:    Level: advanced

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

883: .seealso: SNESFASSetInterpolation(), SNESFASGetInjection(), SNESFASGetRestriction(), SNESFASGetRScale()
884: @*/
885: PetscErrorCode SNESFASGetInterpolation(SNES snes, PetscInt level, Mat *mat)
886: {
887:   SNES_FAS       *fas;
889:   SNES           levelsnes;

892:   SNESFASGetCycleSNES(snes, level, &levelsnes);
893:   fas  = (SNES_FAS*)levelsnes->data;
894:   *mat = fas->interpolate;
895:   return(0);
896: }

900: /*@
901:    SNESFASSetRestriction - Sets the function to be used to restrict the defect
902:    from level l to l-1.

904:    Input Parameters:
905: +  snes  - the multigrid context
906: .  mat   - the restriction matrix
907: -  level - the level (0 is coarsest) to supply [Do not supply 0]

909:    Level: advanced

911:    Notes:
912:           Usually this is the same matrix used also to set the interpolation
913:     for the same level.

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

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

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

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

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

937:   fas->restrct = mat;
938:   return(0);
939: }

943: /*@
944:    SNESFASGetRestriction - Gets the matrix used to calculate the
945:    restriction from l to the l-1th level

947:    Input Parameters:
948: +  snes      - the multigrid context
949: -  level     - the level (0 is coarsest) to supply [do not supply 0]

951:    Output Parameters:
952: .  mat       - the interpolation operator

954:    Level: advanced

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

958: .seealso: SNESFASSetRestriction(), SNESFASGetInjection(), SNESFASGetInterpolation(), SNESFASGetRScale()
959: @*/
960: PetscErrorCode SNESFASGetRestriction(SNES snes, PetscInt level, Mat *mat)
961: {
962:   SNES_FAS       *fas;
964:   SNES           levelsnes;

967:   SNESFASGetCycleSNES(snes, level, &levelsnes);
968:   fas  = (SNES_FAS*)levelsnes->data;
969:   *mat = fas->restrct;
970:   return(0);
971: }


976: /*@
977:    SNESFASSetInjection - Sets the function to be used to inject the solution
978:    from level l to l-1.

980:    Input Parameters:
981:  +  snes  - the multigrid context
982: .  mat   - the restriction matrix
983: -  level - the level (0 is coarsest) to supply [Do not supply 0]

985:    Level: advanced

987:    Notes:
988:          If you do not set this, the restriction and rscale is used to
989:    project the solution instead.

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

993: .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction()
994: @*/
995: PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat)
996: {
997:   SNES_FAS       *fas;
999:   SNES           levelsnes;

1002:   SNESFASGetCycleSNES(snes, level, &levelsnes);
1003:   fas  = (SNES_FAS*)levelsnes->data;
1004:   PetscObjectReference((PetscObject)mat);
1005:   MatDestroy(&fas->inject);

1007:   fas->inject = mat;
1008:   return(0);
1009: }


1014: /*@
1015:    SNESFASGetInjection - Gets the matrix used to calculate the
1016:    injection from l-1 to the lth level

1018:    Input Parameters:
1019: +  snes      - the multigrid context
1020: -  level     - the level (0 is coarsest) to supply [do not supply 0]

1022:    Output Parameters:
1023: .  mat       - the injection operator

1025:    Level: advanced

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

1029: .seealso: SNESFASSetInjection(), SNESFASGetRestriction(), SNESFASGetInterpolation(), SNESFASGetRScale()
1030: @*/
1031: PetscErrorCode SNESFASGetInjection(SNES snes, PetscInt level, Mat *mat)
1032: {
1033:   SNES_FAS       *fas;
1035:   SNES           levelsnes;

1038:   SNESFASGetCycleSNES(snes, level, &levelsnes);
1039:   fas  = (SNES_FAS*)levelsnes->data;
1040:   *mat = fas->inject;
1041:   return(0);
1042: }

1046: /*@
1047:    SNESFASSetRScale - Sets the scaling factor of the restriction
1048:    operator from level l to l-1.

1050:    Input Parameters:
1051: +  snes   - the multigrid context
1052: .  rscale - the restriction scaling
1053: -  level  - the level (0 is coarsest) to supply [Do not supply 0]

1055:    Level: advanced

1057:    Notes:
1058:          This is only used in the case that the injection is not set.

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

1062: .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1063: @*/
1064: PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale)
1065: {
1066:   SNES_FAS       *fas;
1068:   SNES           levelsnes;

1071:   SNESFASGetCycleSNES(snes, level, &levelsnes);
1072:   fas  = (SNES_FAS*)levelsnes->data;
1073:   PetscObjectReference((PetscObject)rscale);
1074:   VecDestroy(&fas->rscale);

1076:   fas->rscale = rscale;
1077:   return(0);
1078: }

1082: /*@
1083:    SNESFASGetSmoother - Gets the default smoother on a level.

1085:    Input Parameters:
1086: +  snes   - the multigrid context
1087: -  level  - the level (0 is coarsest) to supply

1089:    Output Parameters:
1090:    smooth  - the smoother

1092:    Level: advanced

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

1096: .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1097: @*/
1098: PetscErrorCode SNESFASGetSmoother(SNES snes, PetscInt level, SNES *smooth)
1099: {
1100:   SNES_FAS       *fas;
1102:   SNES           levelsnes;

1105:   SNESFASGetCycleSNES(snes, level, &levelsnes);
1106:   fas  = (SNES_FAS*)levelsnes->data;
1107:   if (!fas->smoothd) {
1108:     SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);
1109:   }
1110:   *smooth = fas->smoothd;
1111:   return(0);
1112: }

1116: /*@
1117:    SNESFASGetSmootherDown - Gets the downsmoother on a level.

1119:    Input Parameters:
1120: +  snes   - the multigrid context
1121: -  level  - the level (0 is coarsest) to supply

1123:    Output Parameters:
1124:    smooth  - the smoother

1126:    Level: advanced

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

1130: .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1131: @*/
1132: PetscErrorCode SNESFASGetSmootherDown(SNES snes, PetscInt level, SNES *smooth)
1133: {
1134:   SNES_FAS       *fas;
1136:   SNES           levelsnes;

1139:   SNESFASGetCycleSNES(snes, level, &levelsnes);
1140:   fas  = (SNES_FAS*)levelsnes->data;
1141:   /* if the user chooses to differentiate smoothers, create them both at this point */
1142:   if (!fas->smoothd) {
1143:     SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);
1144:   }
1145:   if (!fas->smoothu) {
1146:     SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu);
1147:   }
1148:   *smooth = fas->smoothd;
1149:   return(0);
1150: }

1154: /*@
1155:    SNESFASGetSmootherUp - Gets the upsmoother on a level.

1157:    Input Parameters:
1158: +  snes   - the multigrid context
1159: -  level  - the level (0 is coarsest)

1161:    Output Parameters:
1162:    smooth  - the smoother

1164:    Level: advanced

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

1168: .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1169: @*/
1170: PetscErrorCode SNESFASGetSmootherUp(SNES snes, PetscInt level, SNES *smooth)
1171: {
1172:   SNES_FAS       *fas;
1174:   SNES           levelsnes;

1177:   SNESFASGetCycleSNES(snes, level, &levelsnes);
1178:   fas  = (SNES_FAS*)levelsnes->data;
1179:   /* if the user chooses to differentiate smoothers, create them both at this point */
1180:   if (!fas->smoothd) {
1181:     SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);
1182:   }
1183:   if (!fas->smoothu) {
1184:     SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu);
1185:   }
1186:   *smooth = fas->smoothu;
1187:   return(0);
1188: }

1192: /*@
1193:    SNESFASGetCoarseSolve - Gets the coarsest solver.

1195:    Input Parameters:
1196: +  snes   - the multigrid context

1198:    Output Parameters:
1199:    solve  - the coarse-level solver

1201:    Level: advanced

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

1205: .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1206: @*/
1207: PetscErrorCode SNESFASGetCoarseSolve(SNES snes, SNES *smooth)
1208: {
1209:   SNES_FAS       *fas;
1211:   SNES           levelsnes;

1214:   SNESFASGetCycleSNES(snes, 0, &levelsnes);
1215:   fas  = (SNES_FAS*)levelsnes->data;
1216:   /* if the user chooses to differentiate smoothers, create them both at this point */
1217:   if (!fas->smoothd) {
1218:     SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);
1219:   }
1220:   *smooth = fas->smoothd;
1221:   return(0);
1222: }

1226: /*@
1227:    SNESFASFullSetDownSweep - Smooth during the initial downsweep for SNESFAS

1229:    Logically Collective on SNES

1231:    Input Parameters:
1232: +  snes - the multigrid context
1233: -  swp - whether to downsweep or not

1235:    Options Database Key:
1236: .  -snes_fas_full_downsweep - Sets number of pre-smoothing steps

1238:    Level: advanced

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

1242: .seealso: SNESFASSetNumberSmoothUp()
1243: @*/
1244: PetscErrorCode SNESFASFullSetDownSweep(SNES snes,PetscBool swp)
1245: {
1246:   SNES_FAS       *fas =  (SNES_FAS*)snes->data;
1247:   PetscErrorCode 0;

1250:   fas->full_downsweep = swp;
1251:   if (fas->next) {
1252:     SNESFASFullSetDownSweep(fas->next,swp);
1253:   }
1254:   return(0);
1255: }