Actual source code: fasfunc.c
1: #include <../src/snes/impls/fas/fasimpls.h>
3: /*@
4: SNESFASSetType - Sets the update and correction type used for `SNESFAS`.
6: Logically Collective
8: Input Parameters:
9: + snes - `SNESFAS` context
10: - fastype - `SNES_FAS_ADDITIVE`, `SNES_FAS_MULTIPLICATIVE`, `SNES_FAS_FULL`, or `SNES_FAS_KASKADE`
12: Level: intermediate
14: .seealso: [](ch_snes), `SNES`, `SNESFAS`, `PCMGSetType()`, `SNESFASGetType()`, `SNES_FAS_ADDITIVE`, `SNES_FAS_MULTIPLICATIVE`, `SNES_FAS_FULL`, `SNES_FAS_KASKADE`
15: @*/
16: PetscErrorCode SNESFASSetType(SNES snes, SNESFASType fastype)
17: {
18: SNES_FAS *fas;
20: PetscFunctionBegin;
23: fas = (SNES_FAS *)snes->data;
24: fas->fastype = fastype;
25: if (fas->next) PetscCall(SNESFASSetType(fas->next, fastype));
26: PetscFunctionReturn(PETSC_SUCCESS);
27: }
29: /*@
30: SNESFASGetType - Gets the update and correction type used for `SNESFAS`.
32: Logically Collective
34: Input Parameter:
35: . snes - `SNESFAS` context
37: Output Parameter:
38: . fastype - `SNES_FAS_ADDITIVE`, `SNES_FAS_MULTIPLICATIVE`, `SNES_FAS_FULL`, or `SNES_FAS_KASKADE`
40: Level: intermediate
42: .seealso: [](ch_snes), `SNES`, `SNESFAS`, `PCMGSetType()`, `SNESFASSetType()`, `SNES_FAS_ADDITIVE`, `SNES_FAS_MULTIPLICATIVE`, `SNES_FAS_FULL`, `SNES_FAS_KASKADE`
43: @*/
44: PetscErrorCode SNESFASGetType(SNES snes, SNESFASType *fastype)
45: {
46: SNES_FAS *fas;
48: PetscFunctionBegin;
50: PetscAssertPointer(fastype, 2);
51: fas = (SNES_FAS *)snes->data;
52: *fastype = fas->fastype;
53: PetscFunctionReturn(PETSC_SUCCESS);
54: }
56: /*@C
57: SNESFASSetLevels - Sets the number of levels to use with `SNESFAS`.
58: Must be called before any other `SNESFAS` routine.
60: Input Parameters:
61: + snes - the `SNES` context of `SNESType` `SNESFAS`
62: . levels - the number of levels
63: - comms - optional communicators for each level; this is to allow solving the coarser
64: problems on smaller sets of processors.
66: Level: intermediate
68: Note:
69: If the number of levels is one then the solver uses the `-fas_levels` prefix
70: for setting the level options rather than the `-fas_coarse` prefix.
72: .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASGetLevels()`
73: @*/
74: PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm *comms)
75: {
76: PetscInt i;
77: const char *optionsprefix;
78: char tprefix[128];
79: SNES_FAS *fas;
80: SNES prevsnes;
81: MPI_Comm comm;
83: PetscFunctionBegin;
85: fas = (SNES_FAS *)snes->data;
86: PetscCall(PetscObjectGetComm((PetscObject)snes, &comm));
87: if (levels == fas->levels) {
88: if (!comms) PetscFunctionReturn(PETSC_SUCCESS);
89: }
90: /* user has changed the number of levels; reset */
91: PetscUseTypeMethod(snes, reset);
92: /* destroy any coarser levels if necessary */
93: PetscCall(SNESDestroy(&fas->next));
94: fas->next = NULL;
95: fas->previous = NULL;
96: prevsnes = snes;
97: /* setup the finest level */
98: PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix));
99: PetscCall(PetscObjectComposedDataSetInt((PetscObject)snes, PetscMGLevelId, levels - 1));
100: for (i = levels - 1; i >= 0; i--) {
101: if (comms) comm = comms[i];
102: fas->level = i;
103: fas->levels = levels;
104: fas->fine = snes;
105: fas->next = NULL;
106: if (i > 0) {
107: PetscCall(SNESCreate(comm, &fas->next));
108: PetscCall(SNESGetOptionsPrefix(fas->fine, &optionsprefix));
109: PetscCall(PetscSNPrintf(tprefix, sizeof(tprefix), "fas_levels_%" PetscInt_FMT "_cycle_", fas->level));
110: PetscCall(SNESAppendOptionsPrefix(fas->next, optionsprefix));
111: PetscCall(SNESAppendOptionsPrefix(fas->next, tprefix));
112: PetscCall(SNESSetType(fas->next, SNESFAS));
113: PetscCall(SNESSetTolerances(fas->next, fas->next->abstol, fas->next->rtol, fas->next->stol, fas->n_cycles, fas->next->max_funcs));
114: PetscCall(PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i));
115: PetscCall(PetscObjectComposedDataSetInt((PetscObject)fas->next, PetscMGLevelId, i - 1));
117: ((SNES_FAS *)fas->next->data)->previous = prevsnes;
119: prevsnes = fas->next;
120: fas = (SNES_FAS *)prevsnes->data;
121: }
122: }
123: PetscFunctionReturn(PETSC_SUCCESS);
124: }
126: /*@
127: SNESFASGetLevels - Gets the number of levels in a `SNESFAS`, including fine and coarse grids
129: Input Parameter:
130: . snes - the `SNES` nonlinear solver context of `SNESType` `SNESFAS`
132: Output Parameter:
133: . levels - the number of levels
135: Level: advanced
137: .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetLevels()`, `PCMGGetLevels()`
138: @*/
139: PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt *levels)
140: {
141: SNES_FAS *fas;
143: PetscFunctionBegin;
145: PetscAssertPointer(levels, 2);
146: fas = (SNES_FAS *)snes->data;
147: *levels = fas->levels;
148: PetscFunctionReturn(PETSC_SUCCESS);
149: }
151: /*@
152: SNESFASGetCycleSNES - Gets the `SNES` corresponding to a particular level of the `SNESFAS` hierarchy
154: Input Parameters:
155: + snes - the `SNES` nonlinear multigrid context
156: - level - the level to get
158: Output Parameter:
159: . lsnes - the `SNES` for the requested level
161: Level: advanced
163: .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetLevels()`, `SNESFASGetLevels()`
164: @*/
165: PetscErrorCode SNESFASGetCycleSNES(SNES snes, PetscInt level, SNES *lsnes)
166: {
167: SNES_FAS *fas;
168: PetscInt i;
170: PetscFunctionBegin;
172: PetscAssertPointer(lsnes, 3);
173: fas = (SNES_FAS *)snes->data;
174: PetscCheck(level <= fas->levels - 1, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "Requested level %" PetscInt_FMT " from SNESFAS containing %" PetscInt_FMT " levels", level, fas->levels);
175: PetscCheck(fas->level == fas->levels - 1, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetCycleSNES may only be called on the finest-level SNES");
177: *lsnes = snes;
178: for (i = fas->level; i > level; i--) {
179: *lsnes = fas->next;
180: fas = (SNES_FAS *)(*lsnes)->data;
181: }
182: PetscCheck(fas->level == level, PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "SNESFAS level hierarchy corrupt");
183: PetscFunctionReturn(PETSC_SUCCESS);
184: }
186: /*@
187: SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to
188: use on all levels.
190: Logically Collective
192: Input Parameters:
193: + snes - the `SNES` nonlinear multigrid context
194: - n - the number of smoothing steps to use
196: Options Database Key:
197: . -snes_fas_smoothup <n> - Sets number of pre-smoothing steps
199: Level: advanced
201: .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetNumberSmoothDown()`
202: @*/
203: PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n)
204: {
205: SNES_FAS *fas;
207: PetscFunctionBegin;
209: fas = (SNES_FAS *)snes->data;
210: fas->max_up_it = n;
211: if (!fas->smoothu && fas->level != 0) PetscCall(SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu));
212: if (fas->smoothu) PetscCall(SNESSetTolerances(fas->smoothu, fas->smoothu->abstol, fas->smoothu->rtol, fas->smoothu->stol, n, fas->smoothu->max_funcs));
213: if (fas->next) PetscCall(SNESFASSetNumberSmoothUp(fas->next, n));
214: PetscFunctionReturn(PETSC_SUCCESS);
215: }
217: /*@
218: SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to
219: use on all levels.
221: Logically Collective
223: Input Parameters:
224: + snes - the `SNESFAS` nonlinear multigrid context
225: - n - the number of smoothing steps to use
227: Options Database Key:
228: . -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps
230: Level: advanced
232: .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetNumberSmoothUp()`
233: @*/
234: PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n)
235: {
236: SNES_FAS *fas;
238: PetscFunctionBegin;
240: fas = (SNES_FAS *)snes->data;
241: if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd));
242: PetscCall(SNESSetTolerances(fas->smoothd, fas->smoothd->abstol, fas->smoothd->rtol, fas->smoothd->stol, n, fas->smoothd->max_funcs));
244: fas->max_down_it = n;
245: if (fas->next) PetscCall(SNESFASSetNumberSmoothDown(fas->next, n));
246: PetscFunctionReturn(PETSC_SUCCESS);
247: }
249: /*@
250: SNESFASSetContinuation - Sets the `SNESFAS` cycle to default to using exact Newton solves on the upsweep
252: Logically Collective
254: Input Parameters:
255: + snes - the `SNESFAS` nonlinear multigrid context
256: - continuation - whether to use continuation
258: Options Database Key:
259: . -snes_fas_continuation - sets continuation to true
261: Level: advanced
263: Note:
264: This sets the prefix on the upsweep smoothers to -fas_continuation
266: .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetNumberSmoothUp()`
267: @*/
268: PetscErrorCode SNESFASSetContinuation(SNES snes, PetscBool continuation)
269: {
270: const char *optionsprefix;
271: char tprefix[128];
272: SNES_FAS *fas;
274: PetscFunctionBegin;
276: fas = (SNES_FAS *)snes->data;
277: PetscCall(SNESGetOptionsPrefix(fas->fine, &optionsprefix));
278: if (!fas->smoothu) PetscCall(SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu));
279: PetscCall(PetscStrncpy(tprefix, "fas_levels_continuation_", sizeof(tprefix)));
280: PetscCall(SNESSetOptionsPrefix(fas->smoothu, optionsprefix));
281: PetscCall(SNESAppendOptionsPrefix(fas->smoothu, tprefix));
282: PetscCall(SNESSetType(fas->smoothu, SNESNEWTONLS));
283: PetscCall(SNESSetTolerances(fas->smoothu, fas->fine->abstol, fas->fine->rtol, fas->fine->stol, 50, 100));
284: fas->continuation = continuation;
285: if (fas->next) PetscCall(SNESFASSetContinuation(fas->next, continuation));
286: PetscFunctionReturn(PETSC_SUCCESS);
287: }
289: /*@
290: SNESFASSetCycles - Sets the number of `SNESFAS` multigrid cycles to use each time a grid is visited. Use `SNESFASSetCyclesOnLevel()` for more
291: complicated cycling.
293: Logically Collective
295: Input Parameters:
296: + snes - the `SNESFAS` nonlinear multigrid context
297: - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
299: Options Database Key:
300: . -snes_fas_cycles <1,2> - 1 for V-cycle, 2 for W-cycle
302: Level: advanced
304: .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetCyclesOnLevel()`
305: @*/
306: PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles)
307: {
308: SNES_FAS *fas;
309: PetscBool isFine;
311: PetscFunctionBegin;
313: PetscCall(SNESFASCycleIsFine(snes, &isFine));
314: fas = (SNES_FAS *)snes->data;
315: fas->n_cycles = cycles;
316: if (!isFine) PetscCall(SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs));
317: if (fas->next) PetscCall(SNESFASSetCycles(fas->next, cycles));
318: PetscFunctionReturn(PETSC_SUCCESS);
319: }
321: /*@C
322: SNESFASSetMonitor - Sets the method-specific cycle monitoring
324: Logically Collective
326: Input Parameters:
327: + snes - the `SNESFAS` context
328: . vf - viewer and format structure (may be `NULL` if `flg` is `PETSC_FALSE`)
329: - flg - monitor or not
331: Level: advanced
333: .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESMonitorSet()`, `SNESFASSetCyclesOnLevel()`
334: @*/
335: PetscErrorCode SNESFASSetMonitor(SNES snes, PetscViewerAndFormat *vf, PetscBool flg)
336: {
337: SNES_FAS *fas;
338: PetscBool isFine;
339: PetscInt i, levels;
340: SNES levelsnes;
342: PetscFunctionBegin;
344: PetscCall(SNESFASCycleIsFine(snes, &isFine));
345: fas = (SNES_FAS *)snes->data;
346: levels = fas->levels;
347: if (isFine) {
348: for (i = 0; i < levels; i++) {
349: PetscCall(SNESFASGetCycleSNES(snes, i, &levelsnes));
350: fas = (SNES_FAS *)levelsnes->data;
351: if (flg) {
352: /* set the monitors for the upsmoother and downsmoother */
353: PetscCall(SNESMonitorCancel(levelsnes));
354: /* Only register destroy on finest level */
355: PetscCall(SNESMonitorSet(levelsnes, (PetscErrorCode (*)(SNES, PetscInt, PetscReal, void *))SNESMonitorDefault, vf, !i ? (PetscCtxDestroyFn *)PetscViewerAndFormatDestroy : NULL));
356: } else if (i != fas->levels - 1) {
357: /* unset the monitors on the coarse levels */
358: PetscCall(SNESMonitorCancel(levelsnes));
359: }
360: }
361: }
362: PetscFunctionReturn(PETSC_SUCCESS);
363: }
365: /*@
366: SNESFASSetLog - Sets or unsets time logging for various `SNESFAS` stages on all levels
368: Logically Collective
370: Input Parameters:
371: + snes - the `SNESFAS` context
372: - flg - whether to log or not
374: Level: advanced
376: .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetMonitor()`
377: @*/
378: PetscErrorCode SNESFASSetLog(SNES snes, PetscBool flg)
379: {
380: SNES_FAS *fas;
381: PetscBool isFine;
382: PetscInt i, levels;
383: SNES levelsnes;
384: char eventname[128];
386: PetscFunctionBegin;
388: PetscCall(SNESFASCycleIsFine(snes, &isFine));
389: fas = (SNES_FAS *)snes->data;
390: levels = fas->levels;
391: if (isFine) {
392: for (i = 0; i < levels; i++) {
393: PetscCall(SNESFASGetCycleSNES(snes, i, &levelsnes));
394: fas = (SNES_FAS *)levelsnes->data;
395: if (flg) {
396: PetscCall(PetscSNPrintf(eventname, sizeof(eventname), "FASSetup %" PetscInt_FMT, i));
397: PetscCall(PetscLogEventRegister(eventname, ((PetscObject)snes)->classid, &fas->eventsmoothsetup));
398: PetscCall(PetscSNPrintf(eventname, sizeof(eventname), "FASSmooth %" PetscInt_FMT, i));
399: PetscCall(PetscLogEventRegister(eventname, ((PetscObject)snes)->classid, &fas->eventsmoothsolve));
400: PetscCall(PetscSNPrintf(eventname, sizeof(eventname), "FASResid %" PetscInt_FMT, i));
401: PetscCall(PetscLogEventRegister(eventname, ((PetscObject)snes)->classid, &fas->eventresidual));
402: PetscCall(PetscSNPrintf(eventname, sizeof(eventname), "FASInterp %" PetscInt_FMT, i));
403: PetscCall(PetscLogEventRegister(eventname, ((PetscObject)snes)->classid, &fas->eventinterprestrict));
404: } else {
405: fas->eventsmoothsetup = 0;
406: fas->eventsmoothsolve = 0;
407: fas->eventresidual = 0;
408: fas->eventinterprestrict = 0;
409: }
410: }
411: }
412: PetscFunctionReturn(PETSC_SUCCESS);
413: }
415: /*
416: Creates the default smoother type.
418: This is SNESNRICHARDSON on each fine level and SNESNEWTONLS on the coarse level.
419: */
420: PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES snes, SNES *smooth)
421: {
422: SNES_FAS *fas;
423: const char *optionsprefix;
424: char tprefix[128];
425: SNES nsmooth;
427: PetscFunctionBegin;
429: PetscAssertPointer(smooth, 2);
430: fas = (SNES_FAS *)snes->data;
431: PetscCall(SNESGetOptionsPrefix(fas->fine, &optionsprefix));
432: /* create the default smoother */
433: PetscCall(SNESCreate(PetscObjectComm((PetscObject)snes), &nsmooth));
434: if (fas->level == 0) {
435: PetscCall(PetscStrncpy(tprefix, "fas_coarse_", sizeof(tprefix)));
436: PetscCall(SNESAppendOptionsPrefix(nsmooth, optionsprefix));
437: PetscCall(SNESAppendOptionsPrefix(nsmooth, tprefix));
438: PetscCall(SNESSetType(nsmooth, SNESNEWTONLS));
439: PetscCall(SNESSetTolerances(nsmooth, nsmooth->abstol, nsmooth->rtol, nsmooth->stol, nsmooth->max_its, nsmooth->max_funcs));
440: } else {
441: PetscCall(PetscSNPrintf(tprefix, sizeof(tprefix), "fas_levels_%" PetscInt_FMT "_", fas->level));
442: PetscCall(SNESAppendOptionsPrefix(nsmooth, optionsprefix));
443: PetscCall(SNESAppendOptionsPrefix(nsmooth, tprefix));
444: PetscCall(SNESSetType(nsmooth, SNESNRICHARDSON));
445: PetscCall(SNESSetTolerances(nsmooth, 0.0, 0.0, 0.0, fas->max_down_it, nsmooth->max_funcs));
446: }
447: PetscCall(PetscObjectIncrementTabLevel((PetscObject)nsmooth, (PetscObject)snes, 1));
448: PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)nsmooth));
449: PetscCall(PetscObjectComposedDataSetInt((PetscObject)nsmooth, PetscMGLevelId, fas->level));
450: *smooth = nsmooth;
451: PetscFunctionReturn(PETSC_SUCCESS);
452: }
454: /* ------------- Functions called on a particular level ----------------- */
456: /*@
457: SNESFASCycleSetCycles - Sets the number of cycles for all levels in a `SNESFAS`
459: Logically Collective
461: Input Parameters:
462: + snes - the `SNESFAS` nonlinear multigrid context
463: - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
465: Level: advanced
467: .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetCycles()`
468: @*/
469: PetscErrorCode SNESFASCycleSetCycles(SNES snes, PetscInt cycles)
470: {
471: SNES_FAS *fas;
473: PetscFunctionBegin;
475: fas = (SNES_FAS *)snes->data;
476: fas->n_cycles = cycles;
477: PetscCall(SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs));
478: PetscFunctionReturn(PETSC_SUCCESS);
479: }
481: /*@
482: SNESFASCycleGetSmoother - Gets the smoother on a particular cycle level.
484: Logically Collective
486: Input Parameter:
487: . snes - the `SNESFAS` obtained with `SNESFASGetCycleSNES()`
489: Output Parameter:
490: . smooth - the smoother
492: Level: advanced
494: .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASCycleGetSmootherUp()`, `SNESFASCycleGetSmootherDown()`, `SNESFASGetCycleSNES()`
495: @*/
496: PetscErrorCode SNESFASCycleGetSmoother(SNES snes, SNES *smooth)
497: {
498: SNES_FAS *fas;
500: PetscFunctionBegin;
502: PetscAssertPointer(smooth, 2);
503: fas = (SNES_FAS *)snes->data;
504: *smooth = fas->smoothd;
505: PetscFunctionReturn(PETSC_SUCCESS);
506: }
507: /*@
508: SNESFASCycleGetSmootherUp - Gets the up smoother on a particular cycle level.
510: Logically Collective
512: Input Parameter:
513: . snes - the `SNESFAS` obtained with `SNESFASGetCycleSNES()`
515: Output Parameter:
516: . smoothu - the smoother
518: Level: advanced
520: Note:
521: Returns the downsmoother if no up smoother is available. This enables transparent
522: default behavior in the process of the solve.
524: .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASCycleGetSmoother()`, `SNESFASCycleGetSmootherDown()`, `SNESFASGetCycleSNES()`
525: @*/
526: PetscErrorCode SNESFASCycleGetSmootherUp(SNES snes, SNES *smoothu)
527: {
528: SNES_FAS *fas;
530: PetscFunctionBegin;
532: PetscAssertPointer(smoothu, 2);
533: fas = (SNES_FAS *)snes->data;
534: if (!fas->smoothu) *smoothu = fas->smoothd;
535: else *smoothu = fas->smoothu;
536: PetscFunctionReturn(PETSC_SUCCESS);
537: }
539: /*@
540: SNESFASCycleGetSmootherDown - Gets the down smoother on a particular cycle level.
542: Logically Collective
544: Input Parameter:
545: . snes - `SNESFAS` obtained with `SNESFASGetCycleSNES()`
547: Output Parameter:
548: . smoothd - the smoother
550: Level: advanced
552: .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASCycleGetSmootherUp()`, `SNESFASCycleGetSmoother()`, `SNESFASGetCycleSNES()`
553: @*/
554: PetscErrorCode SNESFASCycleGetSmootherDown(SNES snes, SNES *smoothd)
555: {
556: SNES_FAS *fas;
558: PetscFunctionBegin;
560: PetscAssertPointer(smoothd, 2);
561: fas = (SNES_FAS *)snes->data;
562: *smoothd = fas->smoothd;
563: PetscFunctionReturn(PETSC_SUCCESS);
564: }
566: /*@
567: SNESFASCycleGetCorrection - Gets the coarse correction `SNESFAS` context for this level
569: Logically Collective
571: Input Parameter:
572: . snes - the `SNESFAS` obtained with `SNESFASGetCycleSNES()`
574: Output Parameter:
575: . correction - the coarse correction solve on this level
577: Level: advanced
579: Note:
580: Returns `NULL` on the coarsest level.
582: .seealso: [](ch_snes), `SNES`, `SNESFAS` `SNESFASCycleGetSmootherUp()`, `SNESFASCycleGetSmoother()`
583: @*/
584: PetscErrorCode SNESFASCycleGetCorrection(SNES snes, SNES *correction)
585: {
586: SNES_FAS *fas;
588: PetscFunctionBegin;
590: PetscAssertPointer(correction, 2);
591: fas = (SNES_FAS *)snes->data;
592: *correction = fas->next;
593: PetscFunctionReturn(PETSC_SUCCESS);
594: }
596: /*@
597: SNESFASCycleGetInterpolation - Gets the interpolation on a level
599: Logically Collective
601: Input Parameter:
602: . snes - the `SNESFAS` obtained with `SNESFASGetCycleSNES()`
604: Output Parameter:
605: . mat - the interpolation operator on this level
607: Level: advanced
609: .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASCycleGetSmootherUp()`, `SNESFASCycleGetSmoother()`
610: @*/
611: PetscErrorCode SNESFASCycleGetInterpolation(SNES snes, Mat *mat)
612: {
613: SNES_FAS *fas;
615: PetscFunctionBegin;
617: PetscAssertPointer(mat, 2);
618: fas = (SNES_FAS *)snes->data;
619: *mat = fas->interpolate;
620: PetscFunctionReturn(PETSC_SUCCESS);
621: }
623: /*@
624: SNESFASCycleGetRestriction - Gets the restriction on a level
626: Logically Collective
628: Input Parameter:
629: . snes - the `SNESFAS` obtained with `SNESFASGetCycleSNES()`
631: Output Parameter:
632: . mat - the restriction operator on this level
634: Level: advanced
636: .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASGetRestriction()`, `SNESFASCycleGetInterpolation()`
637: @*/
638: PetscErrorCode SNESFASCycleGetRestriction(SNES snes, Mat *mat)
639: {
640: SNES_FAS *fas;
642: PetscFunctionBegin;
644: PetscAssertPointer(mat, 2);
645: fas = (SNES_FAS *)snes->data;
646: *mat = fas->restrct;
647: PetscFunctionReturn(PETSC_SUCCESS);
648: }
650: /*@
651: SNESFASCycleGetInjection - Gets the injection on a level
653: Logically Collective
655: Input Parameter:
656: . snes - the `SNESFAS` obtained with `SNESFASGetCycleSNES()`
658: Output Parameter:
659: . mat - the restriction operator on this level
661: Level: advanced
663: .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASGetInjection()`, `SNESFASCycleGetRestriction()`
664: @*/
665: PetscErrorCode SNESFASCycleGetInjection(SNES snes, Mat *mat)
666: {
667: SNES_FAS *fas;
669: PetscFunctionBegin;
671: PetscAssertPointer(mat, 2);
672: fas = (SNES_FAS *)snes->data;
673: *mat = fas->inject;
674: PetscFunctionReturn(PETSC_SUCCESS);
675: }
677: /*@
678: SNESFASCycleGetRScale - Gets the injection scale-factor on a level
680: Logically Collective
682: Input Parameter:
683: . snes - the `SNESFAS` obtained with `SNESFASGetCycleSNES()`
685: Output Parameter:
686: . vec - the restriction operator on this level
688: Level: advanced
690: .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASCycleGetRestriction()`, `SNESFASGetRScale()`
691: @*/
692: PetscErrorCode SNESFASCycleGetRScale(SNES snes, Vec *vec)
693: {
694: SNES_FAS *fas;
696: PetscFunctionBegin;
698: PetscAssertPointer(vec, 2);
699: fas = (SNES_FAS *)snes->data;
700: *vec = fas->rscale;
701: PetscFunctionReturn(PETSC_SUCCESS);
702: }
704: /*@
705: SNESFASCycleIsFine - Determines if a given `SNES` is the finest level in a `SNESFAS`
707: Logically Collective
709: Input Parameter:
710: . snes - the `SNESFAS` context obtained with `SNESFASGetCycleSNES()`
712: Output Parameter:
713: . flg - indicates if this is the fine level or not
715: Level: advanced
717: .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetLevels()`
718: @*/
719: PetscErrorCode SNESFASCycleIsFine(SNES snes, PetscBool *flg)
720: {
721: SNES_FAS *fas;
723: PetscFunctionBegin;
725: PetscAssertPointer(flg, 2);
726: fas = (SNES_FAS *)snes->data;
727: if (fas->level == fas->levels - 1) *flg = PETSC_TRUE;
728: else *flg = PETSC_FALSE;
729: PetscFunctionReturn(PETSC_SUCCESS);
730: }
732: /* functions called on the finest level that return level-specific information */
734: /*@
735: SNESFASSetInterpolation - Sets the `Mat` to be used to apply the
736: interpolation from l-1 to the lth level
738: Input Parameters:
739: + snes - the `SNESFAS` nonlinear multigrid context
740: . mat - the interpolation operator
741: - level - the level (0 is coarsest) to supply [do not supply 0]
743: Level: advanced
745: Notes:
746: Usually this is the same matrix used also to set the restriction
747: for the same level.
749: One can pass in the interpolation matrix or its transpose; PETSc figures
750: out from the matrix dimensions which one it is.
752: .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()`, `SNESFASSetRScale()`
753: @*/
754: PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat)
755: {
756: SNES_FAS *fas;
757: SNES levelsnes;
759: PetscFunctionBegin;
762: PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
763: fas = (SNES_FAS *)levelsnes->data;
764: PetscCall(PetscObjectReference((PetscObject)mat));
765: PetscCall(MatDestroy(&fas->interpolate));
766: fas->interpolate = mat;
767: PetscFunctionReturn(PETSC_SUCCESS);
768: }
770: /*@
771: SNESFASGetInterpolation - Gets the matrix used to calculate the
772: interpolation from l-1 to the lth level
774: Input Parameters:
775: + snes - the `SNESFAS` nonlinear multigrid context
776: - level - the level (0 is coarsest) to supply [do not supply 0]
778: Output Parameter:
779: . mat - the interpolation operator
781: Level: advanced
783: .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetInterpolation()`, `SNESFASGetInjection()`, `SNESFASGetRestriction()`, `SNESFASGetRScale()`
784: @*/
785: PetscErrorCode SNESFASGetInterpolation(SNES snes, PetscInt level, Mat *mat)
786: {
787: SNES_FAS *fas;
788: SNES levelsnes;
790: PetscFunctionBegin;
792: PetscAssertPointer(mat, 3);
793: PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
794: fas = (SNES_FAS *)levelsnes->data;
795: *mat = fas->interpolate;
796: PetscFunctionReturn(PETSC_SUCCESS);
797: }
799: /*@
800: SNESFASSetRestriction - Sets the matrix to be used to restrict the defect
801: from level l to l-1.
803: Input Parameters:
804: + snes - the `SNESFAS` nonlinear multigrid context
805: . mat - the restriction matrix
806: - level - the level (0 is coarsest) to supply [Do not supply 0]
808: Level: advanced
810: Notes:
811: Usually this is the same matrix used also to set the interpolation
812: for the same level.
814: One can pass in the interpolation matrix or its transpose; PETSc figures
815: out from the matrix dimensions which one it is.
817: If you do not set this, the transpose of the `Mat` set with SNESFASSetInterpolation()
818: is used.
820: .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetInterpolation()`, `SNESFASSetInjection()`
821: @*/
822: PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat)
823: {
824: SNES_FAS *fas;
825: SNES levelsnes;
827: PetscFunctionBegin;
830: PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
831: fas = (SNES_FAS *)levelsnes->data;
832: PetscCall(PetscObjectReference((PetscObject)mat));
833: PetscCall(MatDestroy(&fas->restrct));
834: fas->restrct = mat;
835: PetscFunctionReturn(PETSC_SUCCESS);
836: }
838: /*@
839: SNESFASGetRestriction - Gets the matrix used to calculate the
840: restriction from l to the l-1th level
842: Input Parameters:
843: + snes - the `SNESFAS` nonlinear multigrid context
844: - level - the level (0 is coarsest) to supply [do not supply 0]
846: Output Parameter:
847: . mat - the interpolation operator
849: Level: advanced
851: .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetRestriction()`, `SNESFASGetInjection()`, `SNESFASGetInterpolation()`, `SNESFASGetRScale()`
852: @*/
853: PetscErrorCode SNESFASGetRestriction(SNES snes, PetscInt level, Mat *mat)
854: {
855: SNES_FAS *fas;
856: SNES levelsnes;
858: PetscFunctionBegin;
860: PetscAssertPointer(mat, 3);
861: PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
862: fas = (SNES_FAS *)levelsnes->data;
863: *mat = fas->restrct;
864: PetscFunctionReturn(PETSC_SUCCESS);
865: }
867: /*@
868: SNESFASSetInjection - Sets the matrix to be used to inject the solution
869: from `level` to `level-1`.
871: Input Parameters:
872: + snes - the `SNESFAS` nonlinear multigrid context
873: . mat - the injection matrix
874: - level - the level (0 is coarsest) to supply [Do not supply 0]
876: Level: advanced
878: Note:
879: If you do not set this, the restriction and rscale is used to
880: project the solution instead.
882: .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetInterpolation()`, `SNESFASSetRestriction()`
883: @*/
884: PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat)
885: {
886: SNES_FAS *fas;
887: SNES levelsnes;
889: PetscFunctionBegin;
892: PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
893: fas = (SNES_FAS *)levelsnes->data;
894: PetscCall(PetscObjectReference((PetscObject)mat));
895: PetscCall(MatDestroy(&fas->inject));
897: fas->inject = mat;
898: PetscFunctionReturn(PETSC_SUCCESS);
899: }
901: /*@
902: SNESFASGetInjection - Gets the matrix used to calculate the
903: injection from l-1 to the lth level
905: Input Parameters:
906: + snes - the `SNESFAS` nonlinear multigrid context
907: - level - the level (0 is coarsest) to supply [do not supply 0]
909: Output Parameter:
910: . mat - the injection operator
912: Level: advanced
914: .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetInjection()`, `SNESFASGetRestriction()`, `SNESFASGetInterpolation()`, `SNESFASGetRScale()`
915: @*/
916: PetscErrorCode SNESFASGetInjection(SNES snes, PetscInt level, Mat *mat)
917: {
918: SNES_FAS *fas;
919: SNES levelsnes;
921: PetscFunctionBegin;
923: PetscAssertPointer(mat, 3);
924: PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
925: fas = (SNES_FAS *)levelsnes->data;
926: *mat = fas->inject;
927: PetscFunctionReturn(PETSC_SUCCESS);
928: }
930: /*@
931: SNESFASSetRScale - Sets the scaling factor of the restriction
932: operator from level l to l-1.
934: Input Parameters:
935: + snes - the `SNESFAS` nonlinear multigrid context
936: . rscale - the restriction scaling
937: - level - the level (0 is coarsest) to supply [Do not supply 0]
939: Level: advanced
941: Note:
942: This is only used when the injection is not set.
944: .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()`
945: @*/
946: PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale)
947: {
948: SNES_FAS *fas;
949: SNES levelsnes;
951: PetscFunctionBegin;
954: PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
955: fas = (SNES_FAS *)levelsnes->data;
956: PetscCall(PetscObjectReference((PetscObject)rscale));
957: PetscCall(VecDestroy(&fas->rscale));
958: fas->rscale = rscale;
959: PetscFunctionReturn(PETSC_SUCCESS);
960: }
962: /*@
963: SNESFASGetSmoother - Gets the default smoother on a level.
965: Input Parameters:
966: + snes - the `SNESFAS` nonlinear multigrid context
967: - level - the level (0 is coarsest) to supply
969: Output Parameter:
970: . smooth - the smoother
972: Level: advanced
974: .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()`
975: @*/
976: PetscErrorCode SNESFASGetSmoother(SNES snes, PetscInt level, SNES *smooth)
977: {
978: SNES_FAS *fas;
979: SNES levelsnes;
981: PetscFunctionBegin;
983: PetscAssertPointer(smooth, 3);
984: PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
985: fas = (SNES_FAS *)levelsnes->data;
986: if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd));
987: *smooth = fas->smoothd;
988: PetscFunctionReturn(PETSC_SUCCESS);
989: }
991: /*@
992: SNESFASGetSmootherDown - Gets the downsmoother on a level.
994: Input Parameters:
995: + snes - the `SNESFAS` nonlinear multigrid context
996: - level - the level (0 is coarsest) to supply
998: Output Parameter:
999: . smooth - the smoother
1001: Level: advanced
1003: .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()`
1004: @*/
1005: PetscErrorCode SNESFASGetSmootherDown(SNES snes, PetscInt level, SNES *smooth)
1006: {
1007: SNES_FAS *fas;
1008: SNES levelsnes;
1010: PetscFunctionBegin;
1012: PetscAssertPointer(smooth, 3);
1013: PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
1014: fas = (SNES_FAS *)levelsnes->data;
1015: /* if the user chooses to differentiate smoothers, create them both at this point */
1016: if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd));
1017: if (!fas->smoothu) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu));
1018: *smooth = fas->smoothd;
1019: PetscFunctionReturn(PETSC_SUCCESS);
1020: }
1022: /*@
1023: SNESFASGetSmootherUp - Gets the upsmoother on a level.
1025: Input Parameters:
1026: + snes - the `SNESFAS` nonlinear multigrid context
1027: - level - the level (0 is coarsest)
1029: Output Parameter:
1030: . smooth - the smoother
1032: Level: advanced
1034: .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()`
1035: @*/
1036: PetscErrorCode SNESFASGetSmootherUp(SNES snes, PetscInt level, SNES *smooth)
1037: {
1038: SNES_FAS *fas;
1039: SNES levelsnes;
1041: PetscFunctionBegin;
1043: PetscAssertPointer(smooth, 3);
1044: PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
1045: fas = (SNES_FAS *)levelsnes->data;
1046: /* if the user chooses to differentiate smoothers, create them both at this point */
1047: if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd));
1048: if (!fas->smoothu) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu));
1049: *smooth = fas->smoothu;
1050: PetscFunctionReturn(PETSC_SUCCESS);
1051: }
1053: /*@
1054: SNESFASGetCoarseSolve - Gets the coarsest level solver.
1056: Input Parameter:
1057: . snes - the `SNESFAS` nonlinear multigrid context
1059: Output Parameter:
1060: . coarse - the coarse-level solver
1062: Level: advanced
1064: .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()`
1065: @*/
1066: PetscErrorCode SNESFASGetCoarseSolve(SNES snes, SNES *coarse)
1067: {
1068: SNES_FAS *fas;
1069: SNES levelsnes;
1071: PetscFunctionBegin;
1073: PetscAssertPointer(coarse, 2);
1074: PetscCall(SNESFASGetCycleSNES(snes, 0, &levelsnes));
1075: fas = (SNES_FAS *)levelsnes->data;
1076: /* if the user chooses to differentiate smoothers, create them both at this point */
1077: if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd));
1078: *coarse = fas->smoothd;
1079: PetscFunctionReturn(PETSC_SUCCESS);
1080: }
1082: /*@
1083: SNESFASFullSetDownSweep - Smooth during the initial downsweep for `SNESFAS`
1085: Logically Collective
1087: Input Parameters:
1088: + snes - the `SNESFAS` nonlinear multigrid context
1089: - swp - whether to downsweep or not
1091: Options Database Key:
1092: . -snes_fas_full_downsweep - Sets whether to smooth on the initial downsweep
1094: Level: advanced
1096: .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetNumberSmoothUp()`
1097: @*/
1098: PetscErrorCode SNESFASFullSetDownSweep(SNES snes, PetscBool swp)
1099: {
1100: SNES_FAS *fas;
1102: PetscFunctionBegin;
1104: fas = (SNES_FAS *)snes->data;
1105: fas->full_downsweep = swp;
1106: if (fas->next) PetscCall(SNESFASFullSetDownSweep(fas->next, swp));
1107: PetscFunctionReturn(PETSC_SUCCESS);
1108: }
1110: /*@
1111: SNESFASFullSetTotal - Use total residual restriction and total interpolation on the initial down and up sweep of full `SNESFAS` cycles
1113: Logically Collective
1115: Input Parameters:
1116: + snes - the `SNESFAS` nonlinear multigrid context
1117: - total - whether to use total restriction / interpolatiaon or not (the alternative is defect restriction and correction interpolation)
1119: Options Database Key:
1120: . -snes_fas_full_total - Use total restriction and interpolation on the initial down and up sweeps for the full `SNESFAS` cycle
1122: Level: advanced
1124: Note:
1125: This option is only significant if the interpolation of a coarse correction (`MatInterpolate()`) is significantly different from total
1126: solution interpolation (`DMInterpolateSolution()`).
1128: .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetNumberSmoothUp()`, `DMInterpolateSolution()`
1129: @*/
1130: PetscErrorCode SNESFASFullSetTotal(SNES snes, PetscBool total)
1131: {
1132: SNES_FAS *fas;
1134: PetscFunctionBegin;
1136: fas = (SNES_FAS *)snes->data;
1137: fas->full_total = total;
1138: if (fas->next) PetscCall(SNESFASFullSetTotal(fas->next, total));
1139: PetscFunctionReturn(PETSC_SUCCESS);
1140: }
1142: /*@
1143: SNESFASFullGetTotal - Use total residual restriction and total interpolation on the initial down and up sweep of full FAS cycles
1145: Logically Collective
1147: Input Parameter:
1148: . snes - the `SNESFAS` nonlinear multigrid context
1150: Output Parameter:
1151: . total - whether to use total restriction / interpolatiaon or not (the alternative is defect restriction and correction interpolation)
1153: Level: advanced
1155: .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetNumberSmoothUp()`, `DMInterpolateSolution()`, `SNESFullSetTotal()`
1156: @*/
1157: PetscErrorCode SNESFASFullGetTotal(SNES snes, PetscBool *total)
1158: {
1159: SNES_FAS *fas;
1161: PetscFunctionBegin;
1163: fas = (SNES_FAS *)snes->data;
1164: *total = fas->full_total;
1165: PetscFunctionReturn(PETSC_SUCCESS);
1166: }