1a7b5fb5fSBarry Smith #include <../src/snes/impls/fas/fasimpls.h> /*I "petscsnes.h" I*/
2ab8d36c9SPeter Brune
3ab8d36c9SPeter Brune /*@
4420bcc1bSBarry Smith SNESFASSetType - Sets the update and correction type used for `SNESFAS`.
5ab8d36c9SPeter Brune
6ab8d36c9SPeter Brune Logically Collective
7ab8d36c9SPeter Brune
8ab8d36c9SPeter Brune Input Parameters:
9420bcc1bSBarry Smith + snes - `SNESFAS` context
10f6dfbefdSBarry Smith - fastype - `SNES_FAS_ADDITIVE`, `SNES_FAS_MULTIPLICATIVE`, `SNES_FAS_FULL`, or `SNES_FAS_KASKADE`
11583a1250SSatish Balay
12583a1250SSatish Balay Level: intermediate
13ab8d36c9SPeter Brune
14420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESFAS`, `PCMGSetType()`, `SNESFASGetType()`, `SNES_FAS_ADDITIVE`, `SNES_FAS_MULTIPLICATIVE`, `SNES_FAS_FULL`, `SNES_FAS_KASKADE`
15ab8d36c9SPeter Brune @*/
SNESFASSetType(SNES snes,SNESFASType fastype)16d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetType(SNES snes, SNESFASType fastype)
17d71ae5a4SJacob Faibussowitsch {
18f833ba53SLisandro Dalcin SNES_FAS *fas;
1922d28d08SBarry Smith
20ab8d36c9SPeter Brune PetscFunctionBegin;
21f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
22ab8d36c9SPeter Brune PetscValidLogicalCollectiveEnum(snes, fastype, 2);
23f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data;
24ab8d36c9SPeter Brune fas->fastype = fastype;
251baa6e33SBarry Smith if (fas->next) PetscCall(SNESFASSetType(fas->next, fastype));
263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
27ab8d36c9SPeter Brune }
28ab8d36c9SPeter Brune
29ab8d36c9SPeter Brune /*@
30420bcc1bSBarry Smith SNESFASGetType - Gets the update and correction type used for `SNESFAS`.
31ab8d36c9SPeter Brune
32ab8d36c9SPeter Brune Logically Collective
33ab8d36c9SPeter Brune
34f6dfbefdSBarry Smith Input Parameter:
35f6dfbefdSBarry Smith . snes - `SNESFAS` context
36ab8d36c9SPeter Brune
37f6dfbefdSBarry Smith Output Parameter:
38420bcc1bSBarry Smith . fastype - `SNES_FAS_ADDITIVE`, `SNES_FAS_MULTIPLICATIVE`, `SNES_FAS_FULL`, or `SNES_FAS_KASKADE`
39ab8d36c9SPeter Brune
40583a1250SSatish Balay Level: intermediate
41583a1250SSatish Balay
42420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESFAS`, `PCMGSetType()`, `SNESFASSetType()`, `SNES_FAS_ADDITIVE`, `SNES_FAS_MULTIPLICATIVE`, `SNES_FAS_FULL`, `SNES_FAS_KASKADE`
43ab8d36c9SPeter Brune @*/
SNESFASGetType(SNES snes,SNESFASType * fastype)44d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetType(SNES snes, SNESFASType *fastype)
45d71ae5a4SJacob Faibussowitsch {
46f833ba53SLisandro Dalcin SNES_FAS *fas;
47ab8d36c9SPeter Brune
48ab8d36c9SPeter Brune PetscFunctionBegin;
49f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
504f572ea9SToby Isaac PetscAssertPointer(fastype, 2);
51f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data;
52ab8d36c9SPeter Brune *fastype = fas->fastype;
533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
54ab8d36c9SPeter Brune }
55ab8d36c9SPeter Brune
56ab8d36c9SPeter Brune /*@C
57f6dfbefdSBarry Smith SNESFASSetLevels - Sets the number of levels to use with `SNESFAS`.
58420bcc1bSBarry Smith Must be called before any other `SNESFAS` routine.
59ab8d36c9SPeter Brune
60ab8d36c9SPeter Brune Input Parameters:
61420bcc1bSBarry Smith + snes - the `SNES` context of `SNESType` `SNESFAS`
62ab8d36c9SPeter Brune . levels - the number of levels
63ab8d36c9SPeter Brune - comms - optional communicators for each level; this is to allow solving the coarser
642bf68e3eSBarry Smith problems on smaller sets of processors.
65ab8d36c9SPeter Brune
66ab8d36c9SPeter Brune Level: intermediate
67ab8d36c9SPeter Brune
68f6dfbefdSBarry Smith Note:
69420bcc1bSBarry Smith If the number of levels is one then the solver uses the `-fas_levels` prefix
702fe279fdSBarry Smith for setting the level options rather than the `-fas_coarse` prefix.
71ab8d36c9SPeter Brune
72420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASGetLevels()`
73ab8d36c9SPeter Brune @*/
SNESFASSetLevels(SNES snes,PetscInt levels,MPI_Comm * comms)74d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm *comms)
75d71ae5a4SJacob Faibussowitsch {
76ab8d36c9SPeter Brune PetscInt i;
77ab8d36c9SPeter Brune const char *optionsprefix;
78ab8d36c9SPeter Brune char tprefix[128];
79f833ba53SLisandro Dalcin SNES_FAS *fas;
80ab8d36c9SPeter Brune SNES prevsnes;
81ab8d36c9SPeter Brune MPI_Comm comm;
8222d28d08SBarry Smith
83ab8d36c9SPeter Brune PetscFunctionBegin;
84f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
85f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data;
869566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)snes, &comm));
87ab8d36c9SPeter Brune if (levels == fas->levels) {
883ba16761SJacob Faibussowitsch if (!comms) PetscFunctionReturn(PETSC_SUCCESS);
89ab8d36c9SPeter Brune }
90ab8d36c9SPeter Brune /* user has changed the number of levels; reset */
91dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, reset);
92ab8d36c9SPeter Brune /* destroy any coarser levels if necessary */
939566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&fas->next));
940298fd71SBarry Smith fas->next = NULL;
950298fd71SBarry Smith fas->previous = NULL;
96ab8d36c9SPeter Brune prevsnes = snes;
97ab8d36c9SPeter Brune /* setup the finest level */
989566063dSJacob Faibussowitsch PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix));
999566063dSJacob Faibussowitsch PetscCall(PetscObjectComposedDataSetInt((PetscObject)snes, PetscMGLevelId, levels - 1));
100ab8d36c9SPeter Brune for (i = levels - 1; i >= 0; i--) {
101ab8d36c9SPeter Brune if (comms) comm = comms[i];
102ab8d36c9SPeter Brune fas->level = i;
103ab8d36c9SPeter Brune fas->levels = levels;
104ab8d36c9SPeter Brune fas->fine = snes;
1050298fd71SBarry Smith fas->next = NULL;
106ab8d36c9SPeter Brune if (i > 0) {
1079566063dSJacob Faibussowitsch PetscCall(SNESCreate(comm, &fas->next));
1089566063dSJacob Faibussowitsch PetscCall(SNESGetOptionsPrefix(fas->fine, &optionsprefix));
109*835f2295SStefano Zampini PetscCall(PetscSNPrintf(tprefix, sizeof(tprefix), "fas_levels_%" PetscInt_FMT "_cycle_", fas->level));
1109566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(fas->next, optionsprefix));
1119566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(fas->next, tprefix));
1129566063dSJacob Faibussowitsch PetscCall(SNESSetType(fas->next, SNESFAS));
1139566063dSJacob Faibussowitsch PetscCall(SNESSetTolerances(fas->next, fas->next->abstol, fas->next->rtol, fas->next->stol, fas->n_cycles, fas->next->max_funcs));
1149566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i));
1159566063dSJacob Faibussowitsch PetscCall(PetscObjectComposedDataSetInt((PetscObject)fas->next, PetscMGLevelId, i - 1));
1161aa26658SKarl Rupp
117ab8d36c9SPeter Brune ((SNES_FAS *)fas->next->data)->previous = prevsnes;
1181aa26658SKarl Rupp
119ab8d36c9SPeter Brune prevsnes = fas->next;
120ab8d36c9SPeter Brune fas = (SNES_FAS *)prevsnes->data;
121ab8d36c9SPeter Brune }
122ab8d36c9SPeter Brune }
1233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
124ab8d36c9SPeter Brune }
125ab8d36c9SPeter Brune
126ab8d36c9SPeter Brune /*@
127f6dfbefdSBarry Smith SNESFASGetLevels - Gets the number of levels in a `SNESFAS`, including fine and coarse grids
128ab8d36c9SPeter Brune
129ab8d36c9SPeter Brune Input Parameter:
130420bcc1bSBarry Smith . snes - the `SNES` nonlinear solver context of `SNESType` `SNESFAS`
131ab8d36c9SPeter Brune
132e4094ef1SJacob Faibussowitsch Output Parameter:
133ab8d36c9SPeter Brune . levels - the number of levels
134ab8d36c9SPeter Brune
135ab8d36c9SPeter Brune Level: advanced
136ab8d36c9SPeter Brune
137420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetLevels()`, `PCMGGetLevels()`
138ab8d36c9SPeter Brune @*/
SNESFASGetLevels(SNES snes,PetscInt * levels)139d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt *levels)
140d71ae5a4SJacob Faibussowitsch {
141f833ba53SLisandro Dalcin SNES_FAS *fas;
1425fd66863SKarl Rupp
143ab8d36c9SPeter Brune PetscFunctionBegin;
144f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
1454f572ea9SToby Isaac PetscAssertPointer(levels, 2);
146f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data;
147ab8d36c9SPeter Brune *levels = fas->levels;
1483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
149ab8d36c9SPeter Brune }
150ab8d36c9SPeter Brune
151ab8d36c9SPeter Brune /*@
152ceaaa498SBarry Smith SNESFASGetCycleSNES - Gets the `SNES` corresponding to a particular level of the `SNESFAS` hierarchy
153ab8d36c9SPeter Brune
154ab8d36c9SPeter Brune Input Parameters:
155f6dfbefdSBarry Smith + snes - the `SNES` nonlinear multigrid context
156f6dfbefdSBarry Smith - level - the level to get
157f6dfbefdSBarry Smith
158f6dfbefdSBarry Smith Output Parameter:
159f6dfbefdSBarry Smith . lsnes - the `SNES` for the requested level
160ab8d36c9SPeter Brune
161ab8d36c9SPeter Brune Level: advanced
162ab8d36c9SPeter Brune
163420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetLevels()`, `SNESFASGetLevels()`
164ab8d36c9SPeter Brune @*/
SNESFASGetCycleSNES(SNES snes,PetscInt level,SNES * lsnes)165d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetCycleSNES(SNES snes, PetscInt level, SNES *lsnes)
166d71ae5a4SJacob Faibussowitsch {
167f833ba53SLisandro Dalcin SNES_FAS *fas;
168ab8d36c9SPeter Brune PetscInt i;
169ab8d36c9SPeter Brune
170ab8d36c9SPeter Brune PetscFunctionBegin;
171f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
1724f572ea9SToby Isaac PetscAssertPointer(lsnes, 3);
173f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data;
17463a3b9bcSJacob Faibussowitsch 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);
1750b121fc5SBarry Smith PetscCheck(fas->level == fas->levels - 1, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetCycleSNES may only be called on the finest-level SNES");
176ab8d36c9SPeter Brune
177ab8d36c9SPeter Brune *lsnes = snes;
178ab8d36c9SPeter Brune for (i = fas->level; i > level; i--) {
179ab8d36c9SPeter Brune *lsnes = fas->next;
180ab8d36c9SPeter Brune fas = (SNES_FAS *)(*lsnes)->data;
181ab8d36c9SPeter Brune }
18208401ef6SPierre Jolivet PetscCheck(fas->level == level, PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "SNESFAS level hierarchy corrupt");
1833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
184ab8d36c9SPeter Brune }
185ab8d36c9SPeter Brune
186ab8d36c9SPeter Brune /*@
187ab8d36c9SPeter Brune SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to
188ab8d36c9SPeter Brune use on all levels.
189ab8d36c9SPeter Brune
190c3339decSBarry Smith Logically Collective
191ab8d36c9SPeter Brune
192ab8d36c9SPeter Brune Input Parameters:
193f6dfbefdSBarry Smith + snes - the `SNES` nonlinear multigrid context
194420bcc1bSBarry Smith - n - the number of smoothing steps to use
195ab8d36c9SPeter Brune
196ab8d36c9SPeter Brune Options Database Key:
197ab8d36c9SPeter Brune . -snes_fas_smoothup <n> - Sets number of pre-smoothing steps
198ab8d36c9SPeter Brune
199ab8d36c9SPeter Brune Level: advanced
200ab8d36c9SPeter Brune
201420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetNumberSmoothDown()`
202ab8d36c9SPeter Brune @*/
SNESFASSetNumberSmoothUp(SNES snes,PetscInt n)203d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n)
204d71ae5a4SJacob Faibussowitsch {
205f833ba53SLisandro Dalcin SNES_FAS *fas;
20622d28d08SBarry Smith
207ab8d36c9SPeter Brune PetscFunctionBegin;
208f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
209f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data;
210ab8d36c9SPeter Brune fas->max_up_it = n;
21148a46eb9SPierre Jolivet if (!fas->smoothu && fas->level != 0) PetscCall(SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu));
2121baa6e33SBarry Smith if (fas->smoothu) PetscCall(SNESSetTolerances(fas->smoothu, fas->smoothu->abstol, fas->smoothu->rtol, fas->smoothu->stol, n, fas->smoothu->max_funcs));
2131baa6e33SBarry Smith if (fas->next) PetscCall(SNESFASSetNumberSmoothUp(fas->next, n));
2143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
215ab8d36c9SPeter Brune }
216ab8d36c9SPeter Brune
217ab8d36c9SPeter Brune /*@
218ab8d36c9SPeter Brune SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to
219ab8d36c9SPeter Brune use on all levels.
220ab8d36c9SPeter Brune
221c3339decSBarry Smith Logically Collective
222ab8d36c9SPeter Brune
223ab8d36c9SPeter Brune Input Parameters:
224f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context
225420bcc1bSBarry Smith - n - the number of smoothing steps to use
226ab8d36c9SPeter Brune
227ab8d36c9SPeter Brune Options Database Key:
228ab8d36c9SPeter Brune . -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps
229ab8d36c9SPeter Brune
230ab8d36c9SPeter Brune Level: advanced
231ab8d36c9SPeter Brune
232420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetNumberSmoothUp()`
233ab8d36c9SPeter Brune @*/
SNESFASSetNumberSmoothDown(SNES snes,PetscInt n)234d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n)
235d71ae5a4SJacob Faibussowitsch {
236f833ba53SLisandro Dalcin SNES_FAS *fas;
23722d28d08SBarry Smith
238ab8d36c9SPeter Brune PetscFunctionBegin;
239f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
240f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data;
24148a46eb9SPierre Jolivet if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd));
2429566063dSJacob Faibussowitsch PetscCall(SNESSetTolerances(fas->smoothd, fas->smoothd->abstol, fas->smoothd->rtol, fas->smoothd->stol, n, fas->smoothd->max_funcs));
2431aa26658SKarl Rupp
244ab8d36c9SPeter Brune fas->max_down_it = n;
2451baa6e33SBarry Smith if (fas->next) PetscCall(SNESFASSetNumberSmoothDown(fas->next, n));
2463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
247ab8d36c9SPeter Brune }
248ab8d36c9SPeter Brune
24987f44e3fSPeter Brune /*@
250420bcc1bSBarry Smith SNESFASSetContinuation - Sets the `SNESFAS` cycle to default to using exact Newton solves on the upsweep
25187f44e3fSPeter Brune
252c3339decSBarry Smith Logically Collective
25387f44e3fSPeter Brune
25487f44e3fSPeter Brune Input Parameters:
255f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context
256420bcc1bSBarry Smith - continuation - whether to use continuation
25787f44e3fSPeter Brune
25887f44e3fSPeter Brune Options Database Key:
25987f44e3fSPeter Brune . -snes_fas_continuation - sets continuation to true
26087f44e3fSPeter Brune
26187f44e3fSPeter Brune Level: advanced
26287f44e3fSPeter Brune
263f6dfbefdSBarry Smith Note:
26495452b02SPatrick Sanan This sets the prefix on the upsweep smoothers to -fas_continuation
26587f44e3fSPeter Brune
266420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetNumberSmoothUp()`
26787f44e3fSPeter Brune @*/
SNESFASSetContinuation(SNES snes,PetscBool continuation)268d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetContinuation(SNES snes, PetscBool continuation)
269d71ae5a4SJacob Faibussowitsch {
27087f44e3fSPeter Brune const char *optionsprefix;
27187f44e3fSPeter Brune char tprefix[128];
272f833ba53SLisandro Dalcin SNES_FAS *fas;
27387f44e3fSPeter Brune
27487f44e3fSPeter Brune PetscFunctionBegin;
275f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
276f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data;
2779566063dSJacob Faibussowitsch PetscCall(SNESGetOptionsPrefix(fas->fine, &optionsprefix));
27848a46eb9SPierre Jolivet if (!fas->smoothu) PetscCall(SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu));
2799566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(tprefix, "fas_levels_continuation_", sizeof(tprefix)));
2809566063dSJacob Faibussowitsch PetscCall(SNESSetOptionsPrefix(fas->smoothu, optionsprefix));
2819566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(fas->smoothu, tprefix));
2829566063dSJacob Faibussowitsch PetscCall(SNESSetType(fas->smoothu, SNESNEWTONLS));
2839566063dSJacob Faibussowitsch PetscCall(SNESSetTolerances(fas->smoothu, fas->fine->abstol, fas->fine->rtol, fas->fine->stol, 50, 100));
28487f44e3fSPeter Brune fas->continuation = continuation;
2851baa6e33SBarry Smith if (fas->next) PetscCall(SNESFASSetContinuation(fas->next, continuation));
2863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
28787f44e3fSPeter Brune }
28887f44e3fSPeter Brune
289ab8d36c9SPeter Brune /*@
290420bcc1bSBarry Smith SNESFASSetCycles - Sets the number of `SNESFAS` multigrid cycles to use each time a grid is visited. Use `SNESFASSetCyclesOnLevel()` for more
291ab8d36c9SPeter Brune complicated cycling.
292ab8d36c9SPeter Brune
293c3339decSBarry Smith Logically Collective
294ab8d36c9SPeter Brune
295ab8d36c9SPeter Brune Input Parameters:
296f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context
297ab8d36c9SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
298ab8d36c9SPeter Brune
299ab8d36c9SPeter Brune Options Database Key:
30067b8a455SSatish Balay . -snes_fas_cycles <1,2> - 1 for V-cycle, 2 for W-cycle
301ab8d36c9SPeter Brune
302ab8d36c9SPeter Brune Level: advanced
303ab8d36c9SPeter Brune
304420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetCyclesOnLevel()`
305ab8d36c9SPeter Brune @*/
SNESFASSetCycles(SNES snes,PetscInt cycles)306d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles)
307d71ae5a4SJacob Faibussowitsch {
308f833ba53SLisandro Dalcin SNES_FAS *fas;
309ab8d36c9SPeter Brune PetscBool isFine;
31022d28d08SBarry Smith
311ab8d36c9SPeter Brune PetscFunctionBegin;
312f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
3139566063dSJacob Faibussowitsch PetscCall(SNESFASCycleIsFine(snes, &isFine));
314f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data;
315ab8d36c9SPeter Brune fas->n_cycles = cycles;
31648a46eb9SPierre Jolivet if (!isFine) PetscCall(SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs));
3171baa6e33SBarry Smith if (fas->next) PetscCall(SNESFASSetCycles(fas->next, cycles));
3183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
319ab8d36c9SPeter Brune }
320ab8d36c9SPeter Brune
3215d83a8b1SBarry Smith /*@C
322c8c899caSPeter Brune SNESFASSetMonitor - Sets the method-specific cycle monitoring
323c8c899caSPeter Brune
324c3339decSBarry Smith Logically Collective
325c8c899caSPeter Brune
326c8c899caSPeter Brune Input Parameters:
327f6dfbefdSBarry Smith + snes - the `SNESFAS` context
328420bcc1bSBarry Smith . vf - viewer and format structure (may be `NULL` if `flg` is `PETSC_FALSE`)
329c8c899caSPeter Brune - flg - monitor or not
330c8c899caSPeter Brune
331c8c899caSPeter Brune Level: advanced
332c8c899caSPeter Brune
333420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESMonitorSet()`, `SNESFASSetCyclesOnLevel()`
334c8c899caSPeter Brune @*/
SNESFASSetMonitor(SNES snes,PetscViewerAndFormat * vf,PetscBool flg)335d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetMonitor(SNES snes, PetscViewerAndFormat *vf, PetscBool flg)
336d71ae5a4SJacob Faibussowitsch {
337f833ba53SLisandro Dalcin SNES_FAS *fas;
338c8c899caSPeter Brune PetscBool isFine;
339f833ba53SLisandro Dalcin PetscInt i, levels;
340c8c899caSPeter Brune SNES levelsnes;
34122d28d08SBarry Smith
342c8c899caSPeter Brune PetscFunctionBegin;
343f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
3449566063dSJacob Faibussowitsch PetscCall(SNESFASCycleIsFine(snes, &isFine));
345f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data;
346f833ba53SLisandro Dalcin levels = fas->levels;
347c8c899caSPeter Brune if (isFine) {
348c8c899caSPeter Brune for (i = 0; i < levels; i++) {
3499566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, i, &levelsnes));
350c8c899caSPeter Brune fas = (SNES_FAS *)levelsnes->data;
351c8c899caSPeter Brune if (flg) {
352c8c899caSPeter Brune /* set the monitors for the upsmoother and downsmoother */
3539566063dSJacob Faibussowitsch PetscCall(SNESMonitorCancel(levelsnes));
354d142ab34SLawrence Mitchell /* Only register destroy on finest level */
35549abdd8aSBarry Smith PetscCall(SNESMonitorSet(levelsnes, (PetscErrorCode (*)(SNES, PetscInt, PetscReal, void *))SNESMonitorDefault, vf, !i ? (PetscCtxDestroyFn *)PetscViewerAndFormatDestroy : NULL));
3561aa26658SKarl Rupp } else if (i != fas->levels - 1) {
357c8c899caSPeter Brune /* unset the monitors on the coarse levels */
3589566063dSJacob Faibussowitsch PetscCall(SNESMonitorCancel(levelsnes));
359c8c899caSPeter Brune }
360c8c899caSPeter Brune }
361c8c899caSPeter Brune }
3623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
363c8c899caSPeter Brune }
364c8c899caSPeter Brune
3650dd27c6cSPeter Brune /*@
366f6dfbefdSBarry Smith SNESFASSetLog - Sets or unsets time logging for various `SNESFAS` stages on all levels
3670dd27c6cSPeter Brune
368c3339decSBarry Smith Logically Collective
3690dd27c6cSPeter Brune
3700dd27c6cSPeter Brune Input Parameters:
371f6dfbefdSBarry Smith + snes - the `SNESFAS` context
372420bcc1bSBarry Smith - flg - whether to log or not
3730dd27c6cSPeter Brune
3740dd27c6cSPeter Brune Level: advanced
3750dd27c6cSPeter Brune
376420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetMonitor()`
3770dd27c6cSPeter Brune @*/
SNESFASSetLog(SNES snes,PetscBool flg)378d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetLog(SNES snes, PetscBool flg)
379d71ae5a4SJacob Faibussowitsch {
380f833ba53SLisandro Dalcin SNES_FAS *fas;
3810dd27c6cSPeter Brune PetscBool isFine;
382f833ba53SLisandro Dalcin PetscInt i, levels;
3830dd27c6cSPeter Brune SNES levelsnes;
3840dd27c6cSPeter Brune char eventname[128];
3850dd27c6cSPeter Brune
3860dd27c6cSPeter Brune PetscFunctionBegin;
387f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
3889566063dSJacob Faibussowitsch PetscCall(SNESFASCycleIsFine(snes, &isFine));
389f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data;
390f833ba53SLisandro Dalcin levels = fas->levels;
3910dd27c6cSPeter Brune if (isFine) {
3920dd27c6cSPeter Brune for (i = 0; i < levels; i++) {
3939566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, i, &levelsnes));
3940dd27c6cSPeter Brune fas = (SNES_FAS *)levelsnes->data;
3950dd27c6cSPeter Brune if (flg) {
396*835f2295SStefano Zampini PetscCall(PetscSNPrintf(eventname, sizeof(eventname), "FASSetup %" PetscInt_FMT, i));
3979566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(eventname, ((PetscObject)snes)->classid, &fas->eventsmoothsetup));
398*835f2295SStefano Zampini PetscCall(PetscSNPrintf(eventname, sizeof(eventname), "FASSmooth %" PetscInt_FMT, i));
3999566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(eventname, ((PetscObject)snes)->classid, &fas->eventsmoothsolve));
400*835f2295SStefano Zampini PetscCall(PetscSNPrintf(eventname, sizeof(eventname), "FASResid %" PetscInt_FMT, i));
4019566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(eventname, ((PetscObject)snes)->classid, &fas->eventresidual));
402*835f2295SStefano Zampini PetscCall(PetscSNPrintf(eventname, sizeof(eventname), "FASInterp %" PetscInt_FMT, i));
4039566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(eventname, ((PetscObject)snes)->classid, &fas->eventinterprestrict));
4040dd27c6cSPeter Brune } else {
4050298fd71SBarry Smith fas->eventsmoothsetup = 0;
4060298fd71SBarry Smith fas->eventsmoothsolve = 0;
4070298fd71SBarry Smith fas->eventresidual = 0;
4080298fd71SBarry Smith fas->eventinterprestrict = 0;
4090dd27c6cSPeter Brune }
4100dd27c6cSPeter Brune }
4110dd27c6cSPeter Brune }
4123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4130dd27c6cSPeter Brune }
4140dd27c6cSPeter Brune
415ab8d36c9SPeter Brune /*
416ab8d36c9SPeter Brune Creates the default smoother type.
417ab8d36c9SPeter Brune
41804d7464bSBarry Smith This is SNESNRICHARDSON on each fine level and SNESNEWTONLS on the coarse level.
419ab8d36c9SPeter Brune */
SNESFASCycleCreateSmoother_Private(SNES snes,SNES * smooth)420d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES snes, SNES *smooth)
421d71ae5a4SJacob Faibussowitsch {
422ab8d36c9SPeter Brune SNES_FAS *fas;
423ab8d36c9SPeter Brune const char *optionsprefix;
424ab8d36c9SPeter Brune char tprefix[128];
425ab8d36c9SPeter Brune SNES nsmooth;
42622d28d08SBarry Smith
427ab8d36c9SPeter Brune PetscFunctionBegin;
428f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
4294f572ea9SToby Isaac PetscAssertPointer(smooth, 2);
430ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data;
4319566063dSJacob Faibussowitsch PetscCall(SNESGetOptionsPrefix(fas->fine, &optionsprefix));
432ab8d36c9SPeter Brune /* create the default smoother */
4339566063dSJacob Faibussowitsch PetscCall(SNESCreate(PetscObjectComm((PetscObject)snes), &nsmooth));
434ab8d36c9SPeter Brune if (fas->level == 0) {
4359566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(tprefix, "fas_coarse_", sizeof(tprefix)));
4369566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(nsmooth, optionsprefix));
4379566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(nsmooth, tprefix));
4389566063dSJacob Faibussowitsch PetscCall(SNESSetType(nsmooth, SNESNEWTONLS));
4399566063dSJacob Faibussowitsch PetscCall(SNESSetTolerances(nsmooth, nsmooth->abstol, nsmooth->rtol, nsmooth->stol, nsmooth->max_its, nsmooth->max_funcs));
440ab8d36c9SPeter Brune } else {
441*835f2295SStefano Zampini PetscCall(PetscSNPrintf(tprefix, sizeof(tprefix), "fas_levels_%" PetscInt_FMT "_", fas->level));
4429566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(nsmooth, optionsprefix));
4439566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(nsmooth, tprefix));
4449566063dSJacob Faibussowitsch PetscCall(SNESSetType(nsmooth, SNESNRICHARDSON));
4459566063dSJacob Faibussowitsch PetscCall(SNESSetTolerances(nsmooth, 0.0, 0.0, 0.0, fas->max_down_it, nsmooth->max_funcs));
446ab8d36c9SPeter Brune }
4479566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)nsmooth, (PetscObject)snes, 1));
4489566063dSJacob Faibussowitsch PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)nsmooth));
4499566063dSJacob Faibussowitsch PetscCall(PetscObjectComposedDataSetInt((PetscObject)nsmooth, PetscMGLevelId, fas->level));
450ab8d36c9SPeter Brune *smooth = nsmooth;
4513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
452ab8d36c9SPeter Brune }
453ab8d36c9SPeter Brune
454ab8d36c9SPeter Brune /* ------------- Functions called on a particular level ----------------- */
455ab8d36c9SPeter Brune
456ab8d36c9SPeter Brune /*@
457ceaaa498SBarry Smith SNESFASCycleSetCycles - Sets the number of cycles for all levels in a `SNESFAS`
458ab8d36c9SPeter Brune
459c3339decSBarry Smith Logically Collective
460ab8d36c9SPeter Brune
461ab8d36c9SPeter Brune Input Parameters:
462f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context
463ab8d36c9SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
464ab8d36c9SPeter Brune
465ab8d36c9SPeter Brune Level: advanced
466ab8d36c9SPeter Brune
467420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetCycles()`
468ab8d36c9SPeter Brune @*/
SNESFASCycleSetCycles(SNES snes,PetscInt cycles)469d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleSetCycles(SNES snes, PetscInt cycles)
470d71ae5a4SJacob Faibussowitsch {
471f833ba53SLisandro Dalcin SNES_FAS *fas;
47222d28d08SBarry Smith
473ab8d36c9SPeter Brune PetscFunctionBegin;
474f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
475f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data;
476ab8d36c9SPeter Brune fas->n_cycles = cycles;
4779566063dSJacob Faibussowitsch PetscCall(SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs));
4783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
479ab8d36c9SPeter Brune }
480ab8d36c9SPeter Brune
481ab8d36c9SPeter Brune /*@
482ab8d36c9SPeter Brune SNESFASCycleGetSmoother - Gets the smoother on a particular cycle level.
483ab8d36c9SPeter Brune
484c3339decSBarry Smith Logically Collective
485ab8d36c9SPeter Brune
486f6dfbefdSBarry Smith Input Parameter:
487420bcc1bSBarry Smith . snes - the `SNESFAS` obtained with `SNESFASGetCycleSNES()`
488ab8d36c9SPeter Brune
489f6dfbefdSBarry Smith Output Parameter:
490ab8d36c9SPeter Brune . smooth - the smoother
491ab8d36c9SPeter Brune
492ab8d36c9SPeter Brune Level: advanced
493ab8d36c9SPeter Brune
494420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASCycleGetSmootherUp()`, `SNESFASCycleGetSmootherDown()`, `SNESFASGetCycleSNES()`
495ab8d36c9SPeter Brune @*/
SNESFASCycleGetSmoother(SNES snes,SNES * smooth)496d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleGetSmoother(SNES snes, SNES *smooth)
497d71ae5a4SJacob Faibussowitsch {
498ab8d36c9SPeter Brune SNES_FAS *fas;
4995fd66863SKarl Rupp
500ab8d36c9SPeter Brune PetscFunctionBegin;
501f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
5024f572ea9SToby Isaac PetscAssertPointer(smooth, 2);
503ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data;
504ab8d36c9SPeter Brune *smooth = fas->smoothd;
5053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
506ab8d36c9SPeter Brune }
507ab8d36c9SPeter Brune /*@
508ab8d36c9SPeter Brune SNESFASCycleGetSmootherUp - Gets the up smoother on a particular cycle level.
509ab8d36c9SPeter Brune
510c3339decSBarry Smith Logically Collective
511ab8d36c9SPeter Brune
512f6dfbefdSBarry Smith Input Parameter:
513420bcc1bSBarry Smith . snes - the `SNESFAS` obtained with `SNESFASGetCycleSNES()`
514ab8d36c9SPeter Brune
515f6dfbefdSBarry Smith Output Parameter:
516ab8d36c9SPeter Brune . smoothu - the smoother
517ab8d36c9SPeter Brune
518ceaaa498SBarry Smith Level: advanced
519ceaaa498SBarry Smith
520f6dfbefdSBarry Smith Note:
521ab8d36c9SPeter Brune Returns the downsmoother if no up smoother is available. This enables transparent
522ab8d36c9SPeter Brune default behavior in the process of the solve.
523ab8d36c9SPeter Brune
524420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASCycleGetSmoother()`, `SNESFASCycleGetSmootherDown()`, `SNESFASGetCycleSNES()`
525ab8d36c9SPeter Brune @*/
SNESFASCycleGetSmootherUp(SNES snes,SNES * smoothu)526d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleGetSmootherUp(SNES snes, SNES *smoothu)
527d71ae5a4SJacob Faibussowitsch {
528ab8d36c9SPeter Brune SNES_FAS *fas;
5295fd66863SKarl Rupp
530ab8d36c9SPeter Brune PetscFunctionBegin;
531f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
5324f572ea9SToby Isaac PetscAssertPointer(smoothu, 2);
533ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data;
5341aa26658SKarl Rupp if (!fas->smoothu) *smoothu = fas->smoothd;
5351aa26658SKarl Rupp else *smoothu = fas->smoothu;
5363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
537ab8d36c9SPeter Brune }
538ab8d36c9SPeter Brune
539ab8d36c9SPeter Brune /*@
540ab8d36c9SPeter Brune SNESFASCycleGetSmootherDown - Gets the down smoother on a particular cycle level.
541ab8d36c9SPeter Brune
542c3339decSBarry Smith Logically Collective
543ab8d36c9SPeter Brune
544f6dfbefdSBarry Smith Input Parameter:
545420bcc1bSBarry Smith . snes - `SNESFAS` obtained with `SNESFASGetCycleSNES()`
546ab8d36c9SPeter Brune
547f6dfbefdSBarry Smith Output Parameter:
548ab8d36c9SPeter Brune . smoothd - the smoother
549ab8d36c9SPeter Brune
550ab8d36c9SPeter Brune Level: advanced
551ab8d36c9SPeter Brune
552420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASCycleGetSmootherUp()`, `SNESFASCycleGetSmoother()`, `SNESFASGetCycleSNES()`
553ab8d36c9SPeter Brune @*/
SNESFASCycleGetSmootherDown(SNES snes,SNES * smoothd)554d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleGetSmootherDown(SNES snes, SNES *smoothd)
555d71ae5a4SJacob Faibussowitsch {
556ab8d36c9SPeter Brune SNES_FAS *fas;
5575fd66863SKarl Rupp
558ab8d36c9SPeter Brune PetscFunctionBegin;
559f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
5604f572ea9SToby Isaac PetscAssertPointer(smoothd, 2);
561ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data;
562ab8d36c9SPeter Brune *smoothd = fas->smoothd;
5633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
564ab8d36c9SPeter Brune }
565ab8d36c9SPeter Brune
566ab8d36c9SPeter Brune /*@
567420bcc1bSBarry Smith SNESFASCycleGetCorrection - Gets the coarse correction `SNESFAS` context for this level
568ab8d36c9SPeter Brune
569c3339decSBarry Smith Logically Collective
570ab8d36c9SPeter Brune
571f6dfbefdSBarry Smith Input Parameter:
572420bcc1bSBarry Smith . snes - the `SNESFAS` obtained with `SNESFASGetCycleSNES()`
573ab8d36c9SPeter Brune
574f6dfbefdSBarry Smith Output Parameter:
575f6dfbefdSBarry Smith . correction - the coarse correction solve on this level
576ab8d36c9SPeter Brune
577ab8d36c9SPeter Brune Level: advanced
578ab8d36c9SPeter Brune
579420bcc1bSBarry Smith Note:
580420bcc1bSBarry Smith Returns `NULL` on the coarsest level.
581420bcc1bSBarry Smith
582420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESFAS` `SNESFASCycleGetSmootherUp()`, `SNESFASCycleGetSmoother()`
583ab8d36c9SPeter Brune @*/
SNESFASCycleGetCorrection(SNES snes,SNES * correction)584d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleGetCorrection(SNES snes, SNES *correction)
585d71ae5a4SJacob Faibussowitsch {
586ab8d36c9SPeter Brune SNES_FAS *fas;
5875fd66863SKarl Rupp
588ab8d36c9SPeter Brune PetscFunctionBegin;
589f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
5904f572ea9SToby Isaac PetscAssertPointer(correction, 2);
591ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data;
592ab8d36c9SPeter Brune *correction = fas->next;
5933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
594ab8d36c9SPeter Brune }
595ab8d36c9SPeter Brune
596ab8d36c9SPeter Brune /*@
597420bcc1bSBarry Smith SNESFASCycleGetInterpolation - Gets the interpolation on a level
598ab8d36c9SPeter Brune
599c3339decSBarry Smith Logically Collective
600ab8d36c9SPeter Brune
601f6dfbefdSBarry Smith Input Parameter:
602420bcc1bSBarry Smith . snes - the `SNESFAS` obtained with `SNESFASGetCycleSNES()`
603ab8d36c9SPeter Brune
604f6dfbefdSBarry Smith Output Parameter:
605ab8d36c9SPeter Brune . mat - the interpolation operator on this level
606ab8d36c9SPeter Brune
607f6dfbefdSBarry Smith Level: advanced
608ab8d36c9SPeter Brune
609420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASCycleGetSmootherUp()`, `SNESFASCycleGetSmoother()`
610ab8d36c9SPeter Brune @*/
SNESFASCycleGetInterpolation(SNES snes,Mat * mat)611d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleGetInterpolation(SNES snes, Mat *mat)
612d71ae5a4SJacob Faibussowitsch {
613ab8d36c9SPeter Brune SNES_FAS *fas;
6145fd66863SKarl Rupp
615ab8d36c9SPeter Brune PetscFunctionBegin;
616f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
6174f572ea9SToby Isaac PetscAssertPointer(mat, 2);
618ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data;
619ab8d36c9SPeter Brune *mat = fas->interpolate;
6203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
621ab8d36c9SPeter Brune }
622ab8d36c9SPeter Brune
623ab8d36c9SPeter Brune /*@
624420bcc1bSBarry Smith SNESFASCycleGetRestriction - Gets the restriction on a level
625ab8d36c9SPeter Brune
626c3339decSBarry Smith Logically Collective
627ab8d36c9SPeter Brune
628f6dfbefdSBarry Smith Input Parameter:
629420bcc1bSBarry Smith . snes - the `SNESFAS` obtained with `SNESFASGetCycleSNES()`
630ab8d36c9SPeter Brune
631f6dfbefdSBarry Smith Output Parameter:
632ab8d36c9SPeter Brune . mat - the restriction operator on this level
633ab8d36c9SPeter Brune
634f6dfbefdSBarry Smith Level: advanced
635ab8d36c9SPeter Brune
636420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASGetRestriction()`, `SNESFASCycleGetInterpolation()`
637ab8d36c9SPeter Brune @*/
SNESFASCycleGetRestriction(SNES snes,Mat * mat)638d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleGetRestriction(SNES snes, Mat *mat)
639d71ae5a4SJacob Faibussowitsch {
640ab8d36c9SPeter Brune SNES_FAS *fas;
6415fd66863SKarl Rupp
642ab8d36c9SPeter Brune PetscFunctionBegin;
643f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
6444f572ea9SToby Isaac PetscAssertPointer(mat, 2);
645ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data;
646ab8d36c9SPeter Brune *mat = fas->restrct;
6473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
648ab8d36c9SPeter Brune }
649ab8d36c9SPeter Brune
650ab8d36c9SPeter Brune /*@
651420bcc1bSBarry Smith SNESFASCycleGetInjection - Gets the injection on a level
652ab8d36c9SPeter Brune
653c3339decSBarry Smith Logically Collective
654ab8d36c9SPeter Brune
655f6dfbefdSBarry Smith Input Parameter:
656420bcc1bSBarry Smith . snes - the `SNESFAS` obtained with `SNESFASGetCycleSNES()`
657ab8d36c9SPeter Brune
658f6dfbefdSBarry Smith Output Parameter:
659ab8d36c9SPeter Brune . mat - the restriction operator on this level
660ab8d36c9SPeter Brune
661f6dfbefdSBarry Smith Level: advanced
662ab8d36c9SPeter Brune
663420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASGetInjection()`, `SNESFASCycleGetRestriction()`
664ab8d36c9SPeter Brune @*/
SNESFASCycleGetInjection(SNES snes,Mat * mat)665d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleGetInjection(SNES snes, Mat *mat)
666d71ae5a4SJacob Faibussowitsch {
667ab8d36c9SPeter Brune SNES_FAS *fas;
6685fd66863SKarl Rupp
669ab8d36c9SPeter Brune PetscFunctionBegin;
670f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
6714f572ea9SToby Isaac PetscAssertPointer(mat, 2);
672ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data;
673ab8d36c9SPeter Brune *mat = fas->inject;
6743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
675ab8d36c9SPeter Brune }
676ab8d36c9SPeter Brune
677ab8d36c9SPeter Brune /*@
678420bcc1bSBarry Smith SNESFASCycleGetRScale - Gets the injection scale-factor on a level
679ab8d36c9SPeter Brune
680c3339decSBarry Smith Logically Collective
681ab8d36c9SPeter Brune
682f6dfbefdSBarry Smith Input Parameter:
683420bcc1bSBarry Smith . snes - the `SNESFAS` obtained with `SNESFASGetCycleSNES()`
684ab8d36c9SPeter Brune
685f6dfbefdSBarry Smith Output Parameter:
686e4094ef1SJacob Faibussowitsch . vec - the restriction operator on this level
687ab8d36c9SPeter Brune
688f6dfbefdSBarry Smith Level: advanced
689ab8d36c9SPeter Brune
690420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASCycleGetRestriction()`, `SNESFASGetRScale()`
691ab8d36c9SPeter Brune @*/
SNESFASCycleGetRScale(SNES snes,Vec * vec)692d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleGetRScale(SNES snes, Vec *vec)
693d71ae5a4SJacob Faibussowitsch {
694ab8d36c9SPeter Brune SNES_FAS *fas;
6955fd66863SKarl Rupp
696ab8d36c9SPeter Brune PetscFunctionBegin;
697f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
6984f572ea9SToby Isaac PetscAssertPointer(vec, 2);
699ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data;
700ab8d36c9SPeter Brune *vec = fas->rscale;
7013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
702ab8d36c9SPeter Brune }
703ab8d36c9SPeter Brune
704ab8d36c9SPeter Brune /*@
705420bcc1bSBarry Smith SNESFASCycleIsFine - Determines if a given `SNES` is the finest level in a `SNESFAS`
706ab8d36c9SPeter Brune
707c3339decSBarry Smith Logically Collective
708ab8d36c9SPeter Brune
709f6dfbefdSBarry Smith Input Parameter:
710420bcc1bSBarry Smith . snes - the `SNESFAS` context obtained with `SNESFASGetCycleSNES()`
711ab8d36c9SPeter Brune
712f6dfbefdSBarry Smith Output Parameter:
713ab8d36c9SPeter Brune . flg - indicates if this is the fine level or not
714ab8d36c9SPeter Brune
715ab8d36c9SPeter Brune Level: advanced
716ab8d36c9SPeter Brune
717420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetLevels()`
718ab8d36c9SPeter Brune @*/
SNESFASCycleIsFine(SNES snes,PetscBool * flg)719d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleIsFine(SNES snes, PetscBool *flg)
720d71ae5a4SJacob Faibussowitsch {
721ab8d36c9SPeter Brune SNES_FAS *fas;
7225fd66863SKarl Rupp
723ab8d36c9SPeter Brune PetscFunctionBegin;
724f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
7254f572ea9SToby Isaac PetscAssertPointer(flg, 2);
726ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data;
7271aa26658SKarl Rupp if (fas->level == fas->levels - 1) *flg = PETSC_TRUE;
7281aa26658SKarl Rupp else *flg = PETSC_FALSE;
7293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
730ab8d36c9SPeter Brune }
731ab8d36c9SPeter Brune
732f6dfbefdSBarry Smith /* functions called on the finest level that return level-specific information */
733ab8d36c9SPeter Brune
734ab8d36c9SPeter Brune /*@
735f6dfbefdSBarry Smith SNESFASSetInterpolation - Sets the `Mat` to be used to apply the
736ab8d36c9SPeter Brune interpolation from l-1 to the lth level
737ab8d36c9SPeter Brune
738ab8d36c9SPeter Brune Input Parameters:
739f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context
740ab8d36c9SPeter Brune . mat - the interpolation operator
741ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0]
742ab8d36c9SPeter Brune
743ab8d36c9SPeter Brune Level: advanced
744ab8d36c9SPeter Brune
745ab8d36c9SPeter Brune Notes:
746ab8d36c9SPeter Brune Usually this is the same matrix used also to set the restriction
747ab8d36c9SPeter Brune for the same level.
748ab8d36c9SPeter Brune
749ab8d36c9SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures
750420bcc1bSBarry Smith out from the matrix dimensions which one it is.
751ab8d36c9SPeter Brune
752420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()`, `SNESFASSetRScale()`
753ab8d36c9SPeter Brune @*/
SNESFASSetInterpolation(SNES snes,PetscInt level,Mat mat)754d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat)
755d71ae5a4SJacob Faibussowitsch {
75622d28d08SBarry Smith SNES_FAS *fas;
757ab8d36c9SPeter Brune SNES levelsnes;
75822d28d08SBarry Smith
759ab8d36c9SPeter Brune PetscFunctionBegin;
760f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
761f833ba53SLisandro Dalcin if (mat) PetscValidHeaderSpecific(mat, MAT_CLASSID, 3);
7629566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
763ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data;
7649566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat));
7659566063dSJacob Faibussowitsch PetscCall(MatDestroy(&fas->interpolate));
766ab8d36c9SPeter Brune fas->interpolate = mat;
7673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
768ab8d36c9SPeter Brune }
769ab8d36c9SPeter Brune
770ab8d36c9SPeter Brune /*@
771ab8d36c9SPeter Brune SNESFASGetInterpolation - Gets the matrix used to calculate the
772ab8d36c9SPeter Brune interpolation from l-1 to the lth level
773ab8d36c9SPeter Brune
774ab8d36c9SPeter Brune Input Parameters:
775f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context
776ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0]
777ab8d36c9SPeter Brune
778f6dfbefdSBarry Smith Output Parameter:
779ab8d36c9SPeter Brune . mat - the interpolation operator
780ab8d36c9SPeter Brune
781ab8d36c9SPeter Brune Level: advanced
782ab8d36c9SPeter Brune
783420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetInterpolation()`, `SNESFASGetInjection()`, `SNESFASGetRestriction()`, `SNESFASGetRScale()`
784ab8d36c9SPeter Brune @*/
SNESFASGetInterpolation(SNES snes,PetscInt level,Mat * mat)785d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetInterpolation(SNES snes, PetscInt level, Mat *mat)
786d71ae5a4SJacob Faibussowitsch {
78722d28d08SBarry Smith SNES_FAS *fas;
788ab8d36c9SPeter Brune SNES levelsnes;
78922d28d08SBarry Smith
790ab8d36c9SPeter Brune PetscFunctionBegin;
791f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
7924f572ea9SToby Isaac PetscAssertPointer(mat, 3);
7939566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
794ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data;
795ab8d36c9SPeter Brune *mat = fas->interpolate;
7963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
797ab8d36c9SPeter Brune }
798ab8d36c9SPeter Brune
799ab8d36c9SPeter Brune /*@
800f6dfbefdSBarry Smith SNESFASSetRestriction - Sets the matrix to be used to restrict the defect
801ab8d36c9SPeter Brune from level l to l-1.
802ab8d36c9SPeter Brune
803ab8d36c9SPeter Brune Input Parameters:
804f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context
805ab8d36c9SPeter Brune . mat - the restriction matrix
806ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0]
807ab8d36c9SPeter Brune
808ab8d36c9SPeter Brune Level: advanced
809ab8d36c9SPeter Brune
810ab8d36c9SPeter Brune Notes:
811ab8d36c9SPeter Brune Usually this is the same matrix used also to set the interpolation
812ab8d36c9SPeter Brune for the same level.
813ab8d36c9SPeter Brune
814ab8d36c9SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures
815420bcc1bSBarry Smith out from the matrix dimensions which one it is.
816ab8d36c9SPeter Brune
817420bcc1bSBarry Smith If you do not set this, the transpose of the `Mat` set with SNESFASSetInterpolation()
818ab8d36c9SPeter Brune is used.
819ab8d36c9SPeter Brune
820420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetInterpolation()`, `SNESFASSetInjection()`
821ab8d36c9SPeter Brune @*/
SNESFASSetRestriction(SNES snes,PetscInt level,Mat mat)822d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat)
823d71ae5a4SJacob Faibussowitsch {
82422d28d08SBarry Smith SNES_FAS *fas;
825ab8d36c9SPeter Brune SNES levelsnes;
82622d28d08SBarry Smith
827ab8d36c9SPeter Brune PetscFunctionBegin;
828f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
829f833ba53SLisandro Dalcin if (mat) PetscValidHeaderSpecific(mat, MAT_CLASSID, 3);
8309566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
831ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data;
8329566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat));
8339566063dSJacob Faibussowitsch PetscCall(MatDestroy(&fas->restrct));
834ab8d36c9SPeter Brune fas->restrct = mat;
8353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
836ab8d36c9SPeter Brune }
837ab8d36c9SPeter Brune
838ab8d36c9SPeter Brune /*@
839ab8d36c9SPeter Brune SNESFASGetRestriction - Gets the matrix used to calculate the
840ab8d36c9SPeter Brune restriction from l to the l-1th level
841ab8d36c9SPeter Brune
842ab8d36c9SPeter Brune Input Parameters:
843f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context
844ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0]
845ab8d36c9SPeter Brune
846f6dfbefdSBarry Smith Output Parameter:
847ab8d36c9SPeter Brune . mat - the interpolation operator
848ab8d36c9SPeter Brune
849ab8d36c9SPeter Brune Level: advanced
850ab8d36c9SPeter Brune
851420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetRestriction()`, `SNESFASGetInjection()`, `SNESFASGetInterpolation()`, `SNESFASGetRScale()`
852ab8d36c9SPeter Brune @*/
SNESFASGetRestriction(SNES snes,PetscInt level,Mat * mat)853d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetRestriction(SNES snes, PetscInt level, Mat *mat)
854d71ae5a4SJacob Faibussowitsch {
85522d28d08SBarry Smith SNES_FAS *fas;
856ab8d36c9SPeter Brune SNES levelsnes;
85722d28d08SBarry Smith
858ab8d36c9SPeter Brune PetscFunctionBegin;
859f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
8604f572ea9SToby Isaac PetscAssertPointer(mat, 3);
8619566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
862ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data;
863ab8d36c9SPeter Brune *mat = fas->restrct;
8643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
865ab8d36c9SPeter Brune }
866ab8d36c9SPeter Brune
867ab8d36c9SPeter Brune /*@
868420bcc1bSBarry Smith SNESFASSetInjection - Sets the matrix to be used to inject the solution
869420bcc1bSBarry Smith from `level` to `level-1`.
870ab8d36c9SPeter Brune
871ab8d36c9SPeter Brune Input Parameters:
872f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context
873420bcc1bSBarry Smith . mat - the injection matrix
874ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0]
875ab8d36c9SPeter Brune
876ab8d36c9SPeter Brune Level: advanced
877ab8d36c9SPeter Brune
878f6dfbefdSBarry Smith Note:
879ab8d36c9SPeter Brune If you do not set this, the restriction and rscale is used to
880ab8d36c9SPeter Brune project the solution instead.
881ab8d36c9SPeter Brune
882420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetInterpolation()`, `SNESFASSetRestriction()`
883ab8d36c9SPeter Brune @*/
SNESFASSetInjection(SNES snes,PetscInt level,Mat mat)884d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat)
885d71ae5a4SJacob Faibussowitsch {
88622d28d08SBarry Smith SNES_FAS *fas;
887ab8d36c9SPeter Brune SNES levelsnes;
88822d28d08SBarry Smith
889ab8d36c9SPeter Brune PetscFunctionBegin;
890f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
891f833ba53SLisandro Dalcin if (mat) PetscValidHeaderSpecific(mat, MAT_CLASSID, 3);
8929566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
893ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data;
8949566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat));
8959566063dSJacob Faibussowitsch PetscCall(MatDestroy(&fas->inject));
8961aa26658SKarl Rupp
897ab8d36c9SPeter Brune fas->inject = mat;
8983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
899ab8d36c9SPeter Brune }
900ab8d36c9SPeter Brune
901ab8d36c9SPeter Brune /*@
902ab8d36c9SPeter Brune SNESFASGetInjection - Gets the matrix used to calculate the
903ab8d36c9SPeter Brune injection from l-1 to the lth level
904ab8d36c9SPeter Brune
905ab8d36c9SPeter Brune Input Parameters:
906f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context
907ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0]
908ab8d36c9SPeter Brune
909f6dfbefdSBarry Smith Output Parameter:
910ab8d36c9SPeter Brune . mat - the injection operator
911ab8d36c9SPeter Brune
912ab8d36c9SPeter Brune Level: advanced
913ab8d36c9SPeter Brune
914420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetInjection()`, `SNESFASGetRestriction()`, `SNESFASGetInterpolation()`, `SNESFASGetRScale()`
915ab8d36c9SPeter Brune @*/
SNESFASGetInjection(SNES snes,PetscInt level,Mat * mat)916d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetInjection(SNES snes, PetscInt level, Mat *mat)
917d71ae5a4SJacob Faibussowitsch {
91822d28d08SBarry Smith SNES_FAS *fas;
919ab8d36c9SPeter Brune SNES levelsnes;
92022d28d08SBarry Smith
921ab8d36c9SPeter Brune PetscFunctionBegin;
922f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
9234f572ea9SToby Isaac PetscAssertPointer(mat, 3);
9249566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
925ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data;
926ab8d36c9SPeter Brune *mat = fas->inject;
9273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
928ab8d36c9SPeter Brune }
929ab8d36c9SPeter Brune
930ab8d36c9SPeter Brune /*@
931ab8d36c9SPeter Brune SNESFASSetRScale - Sets the scaling factor of the restriction
932ab8d36c9SPeter Brune operator from level l to l-1.
933ab8d36c9SPeter Brune
934ab8d36c9SPeter Brune Input Parameters:
935f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context
936ab8d36c9SPeter Brune . rscale - the restriction scaling
937ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0]
938ab8d36c9SPeter Brune
939ab8d36c9SPeter Brune Level: advanced
940ab8d36c9SPeter Brune
941f6dfbefdSBarry Smith Note:
942420bcc1bSBarry Smith This is only used when the injection is not set.
943ab8d36c9SPeter Brune
944420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()`
945ab8d36c9SPeter Brune @*/
SNESFASSetRScale(SNES snes,PetscInt level,Vec rscale)946d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale)
947d71ae5a4SJacob Faibussowitsch {
948ab8d36c9SPeter Brune SNES_FAS *fas;
949ab8d36c9SPeter Brune SNES levelsnes;
95022d28d08SBarry Smith
951ab8d36c9SPeter Brune PetscFunctionBegin;
952f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
953f833ba53SLisandro Dalcin if (rscale) PetscValidHeaderSpecific(rscale, VEC_CLASSID, 3);
9549566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
955ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data;
9569566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)rscale));
9579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&fas->rscale));
958ab8d36c9SPeter Brune fas->rscale = rscale;
9593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
960ab8d36c9SPeter Brune }
961ab8d36c9SPeter Brune
962ab8d36c9SPeter Brune /*@
963ab8d36c9SPeter Brune SNESFASGetSmoother - Gets the default smoother on a level.
964ab8d36c9SPeter Brune
965ab8d36c9SPeter Brune Input Parameters:
966f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context
967ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply
968ab8d36c9SPeter Brune
969f6dfbefdSBarry Smith Output Parameter:
970ceaaa498SBarry Smith . smooth - the smoother
971ab8d36c9SPeter Brune
972ab8d36c9SPeter Brune Level: advanced
973ab8d36c9SPeter Brune
974420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()`
975ab8d36c9SPeter Brune @*/
SNESFASGetSmoother(SNES snes,PetscInt level,SNES * smooth)976d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetSmoother(SNES snes, PetscInt level, SNES *smooth)
977d71ae5a4SJacob Faibussowitsch {
978ab8d36c9SPeter Brune SNES_FAS *fas;
979ab8d36c9SPeter Brune SNES levelsnes;
98022d28d08SBarry Smith
981ab8d36c9SPeter Brune PetscFunctionBegin;
982f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
9834f572ea9SToby Isaac PetscAssertPointer(smooth, 3);
9849566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
985ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data;
98648a46eb9SPierre Jolivet if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd));
987ab8d36c9SPeter Brune *smooth = fas->smoothd;
9883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
989ab8d36c9SPeter Brune }
990ab8d36c9SPeter Brune
991ab8d36c9SPeter Brune /*@
992ab8d36c9SPeter Brune SNESFASGetSmootherDown - Gets the downsmoother on a level.
993ab8d36c9SPeter Brune
994ab8d36c9SPeter Brune Input Parameters:
995f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context
996ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply
997ab8d36c9SPeter Brune
998f6dfbefdSBarry Smith Output Parameter:
999ceaaa498SBarry Smith . smooth - the smoother
1000ab8d36c9SPeter Brune
1001ab8d36c9SPeter Brune Level: advanced
1002ab8d36c9SPeter Brune
1003420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()`
1004ab8d36c9SPeter Brune @*/
SNESFASGetSmootherDown(SNES snes,PetscInt level,SNES * smooth)1005d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetSmootherDown(SNES snes, PetscInt level, SNES *smooth)
1006d71ae5a4SJacob Faibussowitsch {
1007ab8d36c9SPeter Brune SNES_FAS *fas;
1008ab8d36c9SPeter Brune SNES levelsnes;
100922d28d08SBarry Smith
1010ab8d36c9SPeter Brune PetscFunctionBegin;
1011f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
10124f572ea9SToby Isaac PetscAssertPointer(smooth, 3);
10139566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
1014ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data;
1015ab8d36c9SPeter Brune /* if the user chooses to differentiate smoothers, create them both at this point */
101648a46eb9SPierre Jolivet if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd));
101748a46eb9SPierre Jolivet if (!fas->smoothu) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu));
1018ab8d36c9SPeter Brune *smooth = fas->smoothd;
10193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1020ab8d36c9SPeter Brune }
1021ab8d36c9SPeter Brune
1022ab8d36c9SPeter Brune /*@
1023ab8d36c9SPeter Brune SNESFASGetSmootherUp - Gets the upsmoother on a level.
1024ab8d36c9SPeter Brune
1025ab8d36c9SPeter Brune Input Parameters:
1026f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context
1027ab8d36c9SPeter Brune - level - the level (0 is coarsest)
1028ab8d36c9SPeter Brune
1029f6dfbefdSBarry Smith Output Parameter:
1030ceaaa498SBarry Smith . smooth - the smoother
1031ab8d36c9SPeter Brune
1032ab8d36c9SPeter Brune Level: advanced
1033ab8d36c9SPeter Brune
1034420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()`
1035ab8d36c9SPeter Brune @*/
SNESFASGetSmootherUp(SNES snes,PetscInt level,SNES * smooth)1036d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetSmootherUp(SNES snes, PetscInt level, SNES *smooth)
1037d71ae5a4SJacob Faibussowitsch {
1038ab8d36c9SPeter Brune SNES_FAS *fas;
1039ab8d36c9SPeter Brune SNES levelsnes;
104022d28d08SBarry Smith
1041ab8d36c9SPeter Brune PetscFunctionBegin;
1042f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
10434f572ea9SToby Isaac PetscAssertPointer(smooth, 3);
10449566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
1045ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data;
1046ab8d36c9SPeter Brune /* if the user chooses to differentiate smoothers, create them both at this point */
104748a46eb9SPierre Jolivet if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd));
104848a46eb9SPierre Jolivet if (!fas->smoothu) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu));
1049ab8d36c9SPeter Brune *smooth = fas->smoothu;
10503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1051ab8d36c9SPeter Brune }
1052d6ad1212SPeter Brune
1053d6ad1212SPeter Brune /*@
1054420bcc1bSBarry Smith SNESFASGetCoarseSolve - Gets the coarsest level solver.
1055d6ad1212SPeter Brune
1056f6dfbefdSBarry Smith Input Parameter:
1057f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context
1058d6ad1212SPeter Brune
1059f6dfbefdSBarry Smith Output Parameter:
1060a3a80b83SMatthew G. Knepley . coarse - the coarse-level solver
1061d6ad1212SPeter Brune
1062d6ad1212SPeter Brune Level: advanced
1063d6ad1212SPeter Brune
1064420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()`
1065d6ad1212SPeter Brune @*/
SNESFASGetCoarseSolve(SNES snes,SNES * coarse)1066d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetCoarseSolve(SNES snes, SNES *coarse)
1067d71ae5a4SJacob Faibussowitsch {
1068d6ad1212SPeter Brune SNES_FAS *fas;
1069d6ad1212SPeter Brune SNES levelsnes;
107022d28d08SBarry Smith
1071d6ad1212SPeter Brune PetscFunctionBegin;
1072f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
10734f572ea9SToby Isaac PetscAssertPointer(coarse, 2);
10749566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, 0, &levelsnes));
1075d6ad1212SPeter Brune fas = (SNES_FAS *)levelsnes->data;
1076d6ad1212SPeter Brune /* if the user chooses to differentiate smoothers, create them both at this point */
107748a46eb9SPierre Jolivet if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd));
1078a3a80b83SMatthew G. Knepley *coarse = fas->smoothd;
10793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1080d6ad1212SPeter Brune }
1081928e959bSPeter Brune
1082928e959bSPeter Brune /*@
1083f6dfbefdSBarry Smith SNESFASFullSetDownSweep - Smooth during the initial downsweep for `SNESFAS`
1084928e959bSPeter Brune
1085c3339decSBarry Smith Logically Collective
1086928e959bSPeter Brune
1087928e959bSPeter Brune Input Parameters:
1088f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context
1089928e959bSPeter Brune - swp - whether to downsweep or not
1090928e959bSPeter Brune
1091928e959bSPeter Brune Options Database Key:
1092420bcc1bSBarry Smith . -snes_fas_full_downsweep - Sets whether to smooth on the initial downsweep
1093928e959bSPeter Brune
1094928e959bSPeter Brune Level: advanced
1095928e959bSPeter Brune
1096420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetNumberSmoothUp()`
1097928e959bSPeter Brune @*/
SNESFASFullSetDownSweep(SNES snes,PetscBool swp)1098d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASFullSetDownSweep(SNES snes, PetscBool swp)
1099d71ae5a4SJacob Faibussowitsch {
1100f833ba53SLisandro Dalcin SNES_FAS *fas;
1101928e959bSPeter Brune
1102928e959bSPeter Brune PetscFunctionBegin;
1103f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
1104f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data;
1105928e959bSPeter Brune fas->full_downsweep = swp;
11061baa6e33SBarry Smith if (fas->next) PetscCall(SNESFASFullSetDownSweep(fas->next, swp));
11073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1108928e959bSPeter Brune }
11096dfbafefSToby Isaac
11106dfbafefSToby Isaac /*@
1111420bcc1bSBarry Smith SNESFASFullSetTotal - Use total residual restriction and total interpolation on the initial down and up sweep of full `SNESFAS` cycles
11126dfbafefSToby Isaac
1113c3339decSBarry Smith Logically Collective
11146dfbafefSToby Isaac
11156dfbafefSToby Isaac Input Parameters:
1116f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context
11176dfbafefSToby Isaac - total - whether to use total restriction / interpolatiaon or not (the alternative is defect restriction and correction interpolation)
11186dfbafefSToby Isaac
11196dfbafefSToby Isaac Options Database Key:
1120420bcc1bSBarry Smith . -snes_fas_full_total - Use total restriction and interpolation on the initial down and up sweeps for the full `SNESFAS` cycle
11216dfbafefSToby Isaac
11226dfbafefSToby Isaac Level: advanced
11236dfbafefSToby Isaac
1124f6dfbefdSBarry Smith Note:
1125f6dfbefdSBarry Smith This option is only significant if the interpolation of a coarse correction (`MatInterpolate()`) is significantly different from total
1126f6dfbefdSBarry Smith solution interpolation (`DMInterpolateSolution()`).
11276dfbafefSToby Isaac
1128420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetNumberSmoothUp()`, `DMInterpolateSolution()`
11296dfbafefSToby Isaac @*/
SNESFASFullSetTotal(SNES snes,PetscBool total)1130d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASFullSetTotal(SNES snes, PetscBool total)
1131d71ae5a4SJacob Faibussowitsch {
11326dfbafefSToby Isaac SNES_FAS *fas;
11336dfbafefSToby Isaac
11346dfbafefSToby Isaac PetscFunctionBegin;
11356dfbafefSToby Isaac PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
11366dfbafefSToby Isaac fas = (SNES_FAS *)snes->data;
11376dfbafefSToby Isaac fas->full_total = total;
11381baa6e33SBarry Smith if (fas->next) PetscCall(SNESFASFullSetTotal(fas->next, total));
11393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
11406dfbafefSToby Isaac }
11416dfbafefSToby Isaac
11426dfbafefSToby Isaac /*@
11436dfbafefSToby Isaac SNESFASFullGetTotal - Use total residual restriction and total interpolation on the initial down and up sweep of full FAS cycles
11446dfbafefSToby Isaac
1145c3339decSBarry Smith Logically Collective
11466dfbafefSToby Isaac
1147f6dfbefdSBarry Smith Input Parameter:
1148f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context
11496dfbafefSToby Isaac
1150ceaaa498SBarry Smith Output Parameter:
11516dfbafefSToby Isaac . total - whether to use total restriction / interpolatiaon or not (the alternative is defect restriction and correction interpolation)
11526dfbafefSToby Isaac
11536dfbafefSToby Isaac Level: advanced
11546dfbafefSToby Isaac
1155420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetNumberSmoothUp()`, `DMInterpolateSolution()`, `SNESFullSetTotal()`
11566dfbafefSToby Isaac @*/
SNESFASFullGetTotal(SNES snes,PetscBool * total)1157d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASFullGetTotal(SNES snes, PetscBool *total)
1158d71ae5a4SJacob Faibussowitsch {
11596dfbafefSToby Isaac SNES_FAS *fas;
11606dfbafefSToby Isaac
11616dfbafefSToby Isaac PetscFunctionBegin;
11626dfbafefSToby Isaac PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
11636dfbafefSToby Isaac fas = (SNES_FAS *)snes->data;
11646dfbafefSToby Isaac *total = fas->full_total;
11653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
11666dfbafefSToby Isaac }
1167