xref: /petsc/src/snes/impls/shell/snesshell.c (revision a336c15037c72f93cd561f5a5e11e93175f2efd9)
1 #include <petsc/private/snesimpl.h> /*I   "petscsnes.h"   I*/
2 
3 typedef struct {
4   PetscErrorCode (*solve)(SNES, Vec);
5   PetscCtx ctx;
6 } SNES_Shell;
7 
8 /*@C
9   SNESShellSetSolve - Sets routine to apply as solver to a `SNESSHELL` `SNES` object
10 
11   Logically Collective
12 
13   Input Parameters:
14 + snes  - the `SNES` nonlinear solver context
15 - solve - the application-provided solver routine
16 
17   Calling sequence of `apply`:
18 + snes - the preconditioner, get the application context with `SNESShellGetContext()` provided with `SNESShellSetContext()`
19 - xout - solution vector
20 
21   Level: advanced
22 
23 .seealso: [](ch_snes), `SNES`, `SNESSHELL`, `SNESShellSetContext()`, `SNESShellGetContext()`
24 @*/
25 PetscErrorCode SNESShellSetSolve(SNES snes, PetscErrorCode (*solve)(SNES snes, Vec xout))
26 {
27   PetscFunctionBegin;
28   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
29   PetscTryMethod(snes, "SNESShellSetSolve_C", (SNES, PetscErrorCode (*)(SNES, Vec)), (snes, solve));
30   PetscFunctionReturn(PETSC_SUCCESS);
31 }
32 
33 static PetscErrorCode SNESDestroy_Shell(SNES snes)
34 {
35   PetscFunctionBegin;
36   PetscCall(PetscFree(snes->data));
37   PetscFunctionReturn(PETSC_SUCCESS);
38 }
39 
40 /*@
41   SNESShellGetContext - Returns the user-provided context associated with a `SNESSHELL`
42 
43   Not Collective
44 
45   Input Parameter:
46 . snes - should have been created with `SNESSetType`(snes,`SNESSHELL`);
47 
48   Output Parameter:
49 . ctx - the user provided context
50 
51   Level: advanced
52 
53   Fortran Notes:
54   This only works when the context is a Fortran derived type or a `PetscObject`. Declare `ctx` with
55 .vb
56   type(tUsertype), pointer :: ctx
57 .ve
58 
59 .seealso: [](ch_snes), `SNES`, `SNESSHELL`, `SNESCreateShell()`, `SNESShellSetContext()`
60 @*/
61 PetscErrorCode SNESShellGetContext(SNES snes, PetscCtxRt ctx)
62 {
63   PetscBool flg;
64 
65   PetscFunctionBegin;
66   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
67   PetscAssertPointer(ctx, 2);
68   PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESSHELL, &flg));
69   if (!flg) *(void **)ctx = NULL;
70   else *(void **)ctx = ((SNES_Shell *)snes->data)->ctx;
71   PetscFunctionReturn(PETSC_SUCCESS);
72 }
73 
74 /*@
75   SNESShellSetContext - sets the context for a `SNESSHELL`
76 
77   Logically Collective
78 
79   Input Parameters:
80 + snes - the `SNESSHELL`
81 - ctx  - the context
82 
83   Level: advanced
84 
85 .seealso: [](ch_snes), `SNES`, `SNESSHELL`, `SNESCreateShell()`, `SNESShellGetContext()`
86 @*/
87 PetscErrorCode SNESShellSetContext(SNES snes, PetscCtx ctx)
88 {
89   SNES_Shell *shell = (SNES_Shell *)snes->data;
90   PetscBool   flg;
91 
92   PetscFunctionBegin;
93   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
94   PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESSHELL, &flg));
95   if (flg) shell->ctx = ctx;
96   PetscFunctionReturn(PETSC_SUCCESS);
97 }
98 
99 static PetscErrorCode SNESSolve_Shell(SNES snes)
100 {
101   SNES_Shell *shell = (SNES_Shell *)snes->data;
102 
103   PetscFunctionBegin;
104   PetscCheck(shell->solve, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Must call SNESShellSetSolve() first");
105   snes->reason = SNES_CONVERGED_ITS;
106   PetscCall((*shell->solve)(snes, snes->vec_sol));
107   PetscFunctionReturn(PETSC_SUCCESS);
108 }
109 
110 static PetscErrorCode SNESShellSetSolve_Shell(SNES snes, PetscErrorCode (*solve)(SNES, Vec))
111 {
112   SNES_Shell *shell = (SNES_Shell *)snes->data;
113 
114   PetscFunctionBegin;
115   shell->solve = solve;
116   PetscFunctionReturn(PETSC_SUCCESS);
117 }
118 
119 /*MC
120   SNESSHELL - a user provided nonlinear solver
121 
122    Level: advanced
123 
124 .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESType`, `SNESShellGetContext()`, `SNESShellSetContext()`, `SNESShellSetSolve()`
125 M*/
126 
127 PETSC_EXTERN PetscErrorCode SNESCreate_Shell(SNES snes)
128 {
129   SNES_Shell *shell;
130 
131   PetscFunctionBegin;
132   snes->ops->destroy = SNESDestroy_Shell;
133   snes->ops->solve   = SNESSolve_Shell;
134 
135   snes->usesksp = PETSC_FALSE;
136   snes->usesnpc = PETSC_FALSE;
137 
138   snes->alwayscomputesfinalresidual = PETSC_FALSE;
139 
140   PetscCall(SNESParametersInitialize(snes));
141 
142   PetscCall(PetscNew(&shell));
143   snes->data = (void *)shell;
144   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESShellSetSolve_C", SNESShellSetSolve_Shell));
145   PetscFunctionReturn(PETSC_SUCCESS);
146 }
147