Actual source code: fasfunc.c

petsc-3.7.3 2016-08-01
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:   PetscObjectComposedDataSetInt((PetscObject) snes, PetscMGLevelId, levels-1);
115:   for (i = levels - 1; i >= 0; i--) {
116:     if (comms) comm = comms[i];
117:     fas->level  = i;
118:     fas->levels = levels;
119:     fas->fine   = snes;
120:     fas->next   = NULL;
121:     if (i > 0) {
122:       SNESCreate(comm, &fas->next);
123:       SNESGetOptionsPrefix(fas->fine, &optionsprefix);
124:       sprintf(tprefix,"fas_levels_%d_cycle_",(int)fas->level);
125:       SNESAppendOptionsPrefix(fas->next,optionsprefix);
126:       SNESAppendOptionsPrefix(fas->next,tprefix);
127:       SNESSetType(fas->next, SNESFAS);
128:       SNESSetTolerances(fas->next, fas->next->abstol, fas->next->rtol, fas->next->stol, fas->n_cycles, fas->next->max_funcs);
129:       PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);
130:       PetscObjectComposedDataSetInt((PetscObject) fas->next, PetscMGLevelId, i-1);

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

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


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

147:    Input Parameter:
148: .  snes - the nonlinear solver context

150:    Output parameter:
151: .  levels - the number of levels

153:    Level: advanced

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

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

164:   *levels = fas->levels;
165:   return(0);
166: }


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

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

180:    Level: advanced

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

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

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

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

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

210:    Logically Collective on SNES

212:    Input Parameters:
213: +  snes - the multigrid context
214: -  n    - the number of smoothing steps

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

219:    Level: advanced

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

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

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

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

250:    Logically Collective on SNES

252:    Input Parameters:
253: +  snes - the multigrid context
254: -  n    - the number of smoothing steps

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

259:    Level: advanced

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

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

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

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


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

289:    Logically Collective on SNES

291:    Input Parameters:
292: +  snes - the multigrid context
293: -  n    - the number of smoothing steps

295:    Options Database Key:
296: .  -snes_fas_continuation - sets continuation to true

298:    Level: advanced

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

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

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

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


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

337:    Logically Collective on SNES

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

343:    Options Database Key:
344: .  -snes_fas_cycles 1 or 2

346:    Level: advanced

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

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

359:   SNESFASCycleIsFine(snes, &isFine);

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


374: /*@
375:    SNESFASSetMonitor - Sets the method-specific cycle monitoring

377:    Logically Collective on SNES

379:    Input Parameters:
380: +  snes   - the FAS context
381: .  vf     - viewer and format structure (may be NULL if flg is FALSE)
382: -  flg    - monitor or not

384:    Level: advanced

386: .keywords: FAS, monitor

388: .seealso: SNESFASSetCyclesOnLevel()
389: @*/
390: PetscErrorCode SNESFASSetMonitor(SNES snes, PetscViewerAndFormat *vf, PetscBool flg)
391: {
392:   SNES_FAS       *fas = (SNES_FAS*)snes->data;
394:   PetscBool      isFine;
395:   PetscInt       i, levels = fas->levels;
396:   SNES           levelsnes;

399:   SNESFASCycleIsFine(snes, &isFine);
400:   if (isFine) {
401:     for (i = 0; i < levels; i++) {
402:       SNESFASGetCycleSNES(snes, i, &levelsnes);
403:       fas  = (SNES_FAS*)levelsnes->data;
404:       if (flg) {
405:         /* set the monitors for the upsmoother and downsmoother */
406:         SNESMonitorCancel(levelsnes);
407:         /* Only register destroy on finest level */
408:         SNESMonitorSet(levelsnes,(PetscErrorCode (*)(SNES,PetscInt,PetscReal,void*))SNESMonitorDefault,vf,(!i ? (PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy : NULL));
409:       } else if (i != fas->levels - 1) {
410:         /* unset the monitors on the coarse levels */
411:         SNESMonitorCancel(levelsnes);
412:       }
413:     }
414:   }
415:   return(0);
416: }

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

423:    Logically Collective on SNES

425:    Input Parameters:
426: +  snes   - the FAS context
427: -  flg    - monitor or not

429:    Level: advanced

431: .keywords: FAS, logging

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

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

474: /*
475: Creates the default smoother type.

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

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

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

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

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

522:    Logically Collective on SNES

524:    Input Parameters:
525: +  snes   - the multigrid context
526: .  level  - the level to set the number of cycles on
527: -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle

529:    Level: advanced

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

533: .seealso: SNESFASSetCycles()
534: @*/
535: PetscErrorCode SNESFASCycleSetCycles(SNES snes, PetscInt cycles)
536: {
537:   SNES_FAS       *fas =  (SNES_FAS*)snes->data;

541:   fas->n_cycles = cycles;
542:   SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);
543:   return(0);
544: }


549: /*@
550:    SNESFASCycleGetSmoother - Gets the smoother on a particular cycle level.

552:    Logically Collective on SNES

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

557:    Output Parameters:
558: .  smooth - the smoother

560:    Level: advanced

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

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

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

581:    Logically Collective on SNES

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

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

589:    Notes:
590:    Returns the downsmoother if no up smoother is available.  This enables transparent
591:    default behavior in the process of the solve.

593:    Level: advanced

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

597: .seealso: SNESFASCycleGetSmoother(), SNESFASCycleGetSmootherDown()
598: @*/
599: PetscErrorCode SNESFASCycleGetSmootherUp(SNES snes, SNES *smoothu)
600: {
601:   SNES_FAS *fas;

605:   fas = (SNES_FAS*)snes->data;
606:   if (!fas->smoothu) *smoothu = fas->smoothd;
607:   else *smoothu = fas->smoothu;
608:   return(0);
609: }

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

616:    Logically Collective on SNES

618:    Input Parameters:
619: .  snes   - the multigrid context

621:    Output Parameters:
622: .  smoothd - the smoother

624:    Level: advanced

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

628: .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
629: @*/
630: PetscErrorCode SNESFASCycleGetSmootherDown(SNES snes, SNES *smoothd)
631: {
632:   SNES_FAS *fas;

636:   fas      = (SNES_FAS*)snes->data;
637:   *smoothd = fas->smoothd;
638:   return(0);
639: }


644: /*@
645:    SNESFASCycleGetCorrection - Gets the coarse correction FAS context for this level

647:    Logically Collective on SNES

649:    Input Parameters:
650: .  snes   - the multigrid context

652:    Output Parameters:
653: .  correction - the coarse correction on this level

655:    Notes:
656:    Returns NULL on the coarsest level.

658:    Level: advanced

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

662: .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
663: @*/
664: PetscErrorCode SNESFASCycleGetCorrection(SNES snes, SNES *correction)
665: {
666:   SNES_FAS *fas;

670:   fas         = (SNES_FAS*)snes->data;
671:   *correction = fas->next;
672:   return(0);
673: }

677: /*@
678:    SNESFASCycleGetInterpolation - Gets the interpolation on this level

680:    Logically Collective on SNES

682:    Input Parameters:
683: .  snes   - the multigrid context

685:    Output Parameters:
686: .  mat    - the interpolation operator on this level

688:    Level: developer

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

692: .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
693: @*/
694: PetscErrorCode SNESFASCycleGetInterpolation(SNES snes, Mat *mat)
695: {
696:   SNES_FAS *fas;

700:   fas  = (SNES_FAS*)snes->data;
701:   *mat = fas->interpolate;
702:   return(0);
703: }


708: /*@
709:    SNESFASCycleGetRestriction - Gets the restriction on this level

711:    Logically Collective on SNES

713:    Input Parameters:
714: .  snes   - the multigrid context

716:    Output Parameters:
717: .  mat    - the restriction operator on this level

719:    Level: developer

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

723: .seealso: SNESFASGetRestriction(), SNESFASCycleGetInterpolation()
724: @*/
725: PetscErrorCode SNESFASCycleGetRestriction(SNES snes, Mat *mat)
726: {
727:   SNES_FAS *fas;

731:   fas  = (SNES_FAS*)snes->data;
732:   *mat = fas->restrct;
733:   return(0);
734: }


739: /*@
740:    SNESFASCycleGetInjection - Gets the injection on this level

742:    Logically Collective on SNES

744:    Input Parameters:
745: .  snes   - the multigrid context

747:    Output Parameters:
748: .  mat    - the restriction operator on this level

750:    Level: developer

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

754: .seealso: SNESFASGetInjection(), SNESFASCycleGetRestriction()
755: @*/
756: PetscErrorCode SNESFASCycleGetInjection(SNES snes, Mat *mat)
757: {
758:   SNES_FAS *fas;

762:   fas  = (SNES_FAS*)snes->data;
763:   *mat = fas->inject;
764:   return(0);
765: }

769: /*@
770:    SNESFASCycleGetRScale - Gets the injection on this level

772:    Logically Collective on SNES

774:    Input Parameters:
775: .  snes   - the multigrid context

777:    Output Parameters:
778: .  mat    - the restriction operator on this level

780:    Level: developer

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

784: .seealso: SNESFASCycleGetRestriction(), SNESFASGetRScale()
785: @*/
786: PetscErrorCode SNESFASCycleGetRScale(SNES snes, Vec *vec)
787: {
788:   SNES_FAS *fas;

792:   fas  = (SNES_FAS*)snes->data;
793:   *vec = fas->rscale;
794:   return(0);
795: }

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

802:    Logically Collective on SNES

804:    Input Parameters:
805: .  snes   - the FAS context

807:    Output Parameters:
808: .  flg - indicates if this is the fine level or not

810:    Level: advanced

812: .keywords: SNES, FAS

814: .seealso: SNESFASSetLevels()
815: @*/
816: PetscErrorCode SNESFASCycleIsFine(SNES snes, PetscBool *flg)
817: {
818:   SNES_FAS *fas;

822:   fas = (SNES_FAS*)snes->data;
823:   if (fas->level == fas->levels - 1) *flg = PETSC_TRUE;
824:   else *flg = PETSC_FALSE;
825:   return(0);
826: }

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

832: /*@
833:    SNESFASSetInterpolation - Sets the function to be used to calculate the
834:    interpolation from l-1 to the lth level

836:    Input Parameters:
837: +  snes      - the multigrid context
838: .  mat       - the interpolation operator
839: -  level     - the level (0 is coarsest) to supply [do not supply 0]

841:    Level: advanced

843:    Notes:
844:           Usually this is the same matrix used also to set the restriction
845:     for the same level.

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

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

852: .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRScale()
853: @*/
854: PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat)
855: {
856:   SNES_FAS       *fas;
858:   SNES           levelsnes;

861:   SNESFASGetCycleSNES(snes, level, &levelsnes);
862:   fas  = (SNES_FAS*)levelsnes->data;
863:   PetscObjectReference((PetscObject)mat);
864:   MatDestroy(&fas->interpolate);

866:   fas->interpolate = mat;
867:   return(0);
868: }

872: /*@
873:    SNESFASGetInterpolation - Gets the matrix used to calculate the
874:    interpolation from l-1 to the lth level

876:    Input Parameters:
877: +  snes      - the multigrid context
878: -  level     - the level (0 is coarsest) to supply [do not supply 0]

880:    Output Parameters:
881: .  mat       - the interpolation operator

883:    Level: advanced

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

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

896:   SNESFASGetCycleSNES(snes, level, &levelsnes);
897:   fas  = (SNES_FAS*)levelsnes->data;
898:   *mat = fas->interpolate;
899:   return(0);
900: }

904: /*@
905:    SNESFASSetRestriction - Sets the function to be used to restrict the defect
906:    from level l to l-1.

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

913:    Level: advanced

915:    Notes:
916:           Usually this is the same matrix used also to set the interpolation
917:     for the same level.

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

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

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

927: .seealso: SNESFASSetInterpolation(), SNESFASSetInjection()
928: @*/
929: PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat)
930: {
931:   SNES_FAS       *fas;
933:   SNES           levelsnes;

936:   SNESFASGetCycleSNES(snes, level, &levelsnes);
937:   fas  = (SNES_FAS*)levelsnes->data;
938:   PetscObjectReference((PetscObject)mat);
939:   MatDestroy(&fas->restrct);

941:   fas->restrct = mat;
942:   return(0);
943: }

947: /*@
948:    SNESFASGetRestriction - Gets the matrix used to calculate the
949:    restriction from l to the l-1th level

951:    Input Parameters:
952: +  snes      - the multigrid context
953: -  level     - the level (0 is coarsest) to supply [do not supply 0]

955:    Output Parameters:
956: .  mat       - the interpolation operator

958:    Level: advanced

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

962: .seealso: SNESFASSetRestriction(), SNESFASGetInjection(), SNESFASGetInterpolation(), SNESFASGetRScale()
963: @*/
964: PetscErrorCode SNESFASGetRestriction(SNES snes, PetscInt level, Mat *mat)
965: {
966:   SNES_FAS       *fas;
968:   SNES           levelsnes;

971:   SNESFASGetCycleSNES(snes, level, &levelsnes);
972:   fas  = (SNES_FAS*)levelsnes->data;
973:   *mat = fas->restrct;
974:   return(0);
975: }


980: /*@
981:    SNESFASSetInjection - Sets the function to be used to inject the solution
982:    from level l to l-1.

984:    Input Parameters:
985:  +  snes  - the multigrid context
986: .  mat   - the restriction matrix
987: -  level - the level (0 is coarsest) to supply [Do not supply 0]

989:    Level: advanced

991:    Notes:
992:          If you do not set this, the restriction and rscale is used to
993:    project the solution instead.

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

997: .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction()
998: @*/
999: PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat)
1000: {
1001:   SNES_FAS       *fas;
1003:   SNES           levelsnes;

1006:   SNESFASGetCycleSNES(snes, level, &levelsnes);
1007:   fas  = (SNES_FAS*)levelsnes->data;
1008:   PetscObjectReference((PetscObject)mat);
1009:   MatDestroy(&fas->inject);

1011:   fas->inject = mat;
1012:   return(0);
1013: }


1018: /*@
1019:    SNESFASGetInjection - Gets the matrix used to calculate the
1020:    injection from l-1 to the lth level

1022:    Input Parameters:
1023: +  snes      - the multigrid context
1024: -  level     - the level (0 is coarsest) to supply [do not supply 0]

1026:    Output Parameters:
1027: .  mat       - the injection operator

1029:    Level: advanced

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

1033: .seealso: SNESFASSetInjection(), SNESFASGetRestriction(), SNESFASGetInterpolation(), SNESFASGetRScale()
1034: @*/
1035: PetscErrorCode SNESFASGetInjection(SNES snes, PetscInt level, Mat *mat)
1036: {
1037:   SNES_FAS       *fas;
1039:   SNES           levelsnes;

1042:   SNESFASGetCycleSNES(snes, level, &levelsnes);
1043:   fas  = (SNES_FAS*)levelsnes->data;
1044:   *mat = fas->inject;
1045:   return(0);
1046: }

1050: /*@
1051:    SNESFASSetRScale - Sets the scaling factor of the restriction
1052:    operator from level l to l-1.

1054:    Input Parameters:
1055: +  snes   - the multigrid context
1056: .  rscale - the restriction scaling
1057: -  level  - the level (0 is coarsest) to supply [Do not supply 0]

1059:    Level: advanced

1061:    Notes:
1062:          This is only used in the case that the injection is not set.

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

1066: .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1067: @*/
1068: PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale)
1069: {
1070:   SNES_FAS       *fas;
1072:   SNES           levelsnes;

1075:   SNESFASGetCycleSNES(snes, level, &levelsnes);
1076:   fas  = (SNES_FAS*)levelsnes->data;
1077:   PetscObjectReference((PetscObject)rscale);
1078:   VecDestroy(&fas->rscale);

1080:   fas->rscale = rscale;
1081:   return(0);
1082: }

1086: /*@
1087:    SNESFASGetSmoother - Gets the default smoother on a level.

1089:    Input Parameters:
1090: +  snes   - the multigrid context
1091: -  level  - the level (0 is coarsest) to supply

1093:    Output Parameters:
1094:    smooth  - the smoother

1096:    Level: advanced

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

1100: .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1101: @*/
1102: PetscErrorCode SNESFASGetSmoother(SNES snes, PetscInt level, SNES *smooth)
1103: {
1104:   SNES_FAS       *fas;
1106:   SNES           levelsnes;

1109:   SNESFASGetCycleSNES(snes, level, &levelsnes);
1110:   fas  = (SNES_FAS*)levelsnes->data;
1111:   if (!fas->smoothd) {
1112:     SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);
1113:   }
1114:   *smooth = fas->smoothd;
1115:   return(0);
1116: }

1120: /*@
1121:    SNESFASGetSmootherDown - Gets the downsmoother on a level.

1123:    Input Parameters:
1124: +  snes   - the multigrid context
1125: -  level  - the level (0 is coarsest) to supply

1127:    Output Parameters:
1128:    smooth  - the smoother

1130:    Level: advanced

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

1134: .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1135: @*/
1136: PetscErrorCode SNESFASGetSmootherDown(SNES snes, PetscInt level, SNES *smooth)
1137: {
1138:   SNES_FAS       *fas;
1140:   SNES           levelsnes;

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

1158: /*@
1159:    SNESFASGetSmootherUp - Gets the upsmoother on a level.

1161:    Input Parameters:
1162: +  snes   - the multigrid context
1163: -  level  - the level (0 is coarsest)

1165:    Output Parameters:
1166:    smooth  - the smoother

1168:    Level: advanced

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

1172: .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1173: @*/
1174: PetscErrorCode SNESFASGetSmootherUp(SNES snes, PetscInt level, SNES *smooth)
1175: {
1176:   SNES_FAS       *fas;
1178:   SNES           levelsnes;

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

1196: /*@
1197:   SNESFASGetCoarseSolve - Gets the coarsest solver.

1199:   Input Parameters:
1200: . snes - the multigrid context

1202:   Output Parameters:
1203: . coarse - the coarse-level solver

1205:   Level: advanced

1207: .keywords: FAS, MG, get, multigrid, solver, coarse
1208: .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1209: @*/
1210: PetscErrorCode SNESFASGetCoarseSolve(SNES snes, SNES *coarse)
1211: {
1212:   SNES_FAS       *fas;
1214:   SNES           levelsnes;

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

1229: /*@
1230:    SNESFASFullSetDownSweep - Smooth during the initial downsweep for SNESFAS

1232:    Logically Collective on SNES

1234:    Input Parameters:
1235: +  snes - the multigrid context
1236: -  swp - whether to downsweep or not

1238:    Options Database Key:
1239: .  -snes_fas_full_downsweep - Sets number of pre-smoothing steps

1241:    Level: advanced

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

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

1253:   fas->full_downsweep = swp;
1254:   if (fas->next) {
1255:     SNESFASFullSetDownSweep(fas->next,swp);
1256:   }
1257:   return(0);
1258: }