xref: /petsc/src/snes/impls/shell/snesshell.c (revision 8a6b6cad2f7b2cdc69b9bd79694a63724703a50a)
1 #include <petsc/private/snesimpl.h> /*I   "petscsnes.h"   I*/
2 
3 typedef struct {
4   PetscErrorCode (*solve)(SNES, Vec);
5   void *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 .seealso: [](ch_snes), `SNES`, `SNESSHELL`, `SNESCreateShell()`, `SNESShellSetContext()`
54 @*/
55 PetscErrorCode SNESShellGetContext(SNES snes, void *ctx)
56 {
57   PetscBool flg;
58 
59   PetscFunctionBegin;
60   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
61   PetscAssertPointer(ctx, 2);
62   PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESSHELL, &flg));
63   if (!flg) *(void **)ctx = NULL;
64   else *(void **)ctx = ((SNES_Shell *)snes->data)->ctx;
65   PetscFunctionReturn(PETSC_SUCCESS);
66 }
67 
68 /*@
69   SNESShellSetContext - sets the context for a `SNESSHELL`
70 
71   Logically Collective
72 
73   Input Parameters:
74 + snes - the `SNESSHELL`
75 - ctx  - the context
76 
77   Level: advanced
78 
79   Fortran Note:
80   The context can only be an integer or a `PetscObject` it cannot be a Fortran array or derived type.
81 
82 .seealso: [](ch_snes), `SNES`, `SNESSHELL`, `SNESCreateShell()`, `SNESShellGetContext()`
83 @*/
84 PetscErrorCode SNESShellSetContext(SNES snes, void *ctx)
85 {
86   SNES_Shell *shell = (SNES_Shell *)snes->data;
87   PetscBool   flg;
88 
89   PetscFunctionBegin;
90   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
91   PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESSHELL, &flg));
92   if (flg) shell->ctx = ctx;
93   PetscFunctionReturn(PETSC_SUCCESS);
94 }
95 
96 static PetscErrorCode SNESSolve_Shell(SNES snes)
97 {
98   SNES_Shell *shell = (SNES_Shell *)snes->data;
99 
100   PetscFunctionBegin;
101   PetscCheck(shell->solve, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Must call SNESShellSetSolve() first");
102   snes->reason = SNES_CONVERGED_ITS;
103   PetscCall((*shell->solve)(snes, snes->vec_sol));
104   PetscFunctionReturn(PETSC_SUCCESS);
105 }
106 
107 static PetscErrorCode SNESShellSetSolve_Shell(SNES snes, PetscErrorCode (*solve)(SNES, Vec))
108 {
109   SNES_Shell *shell = (SNES_Shell *)snes->data;
110 
111   PetscFunctionBegin;
112   shell->solve = solve;
113   PetscFunctionReturn(PETSC_SUCCESS);
114 }
115 
116 /*MC
117   SNESSHELL - a user provided nonlinear solver
118 
119    Level: advanced
120 
121 .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESType`, `SNESShellGetContext()`, `SNESShellSetContext()`, `SNESShellSetSolve()`
122 M*/
123 
124 PETSC_EXTERN PetscErrorCode SNESCreate_Shell(SNES snes)
125 {
126   SNES_Shell *shell;
127 
128   PetscFunctionBegin;
129   snes->ops->destroy = SNESDestroy_Shell;
130   snes->ops->solve   = SNESSolve_Shell;
131 
132   snes->usesksp = PETSC_FALSE;
133   snes->usesnpc = PETSC_FALSE;
134 
135   snes->alwayscomputesfinalresidual = PETSC_FALSE;
136 
137   PetscCall(SNESParametersInitialize(snes));
138 
139   PetscCall(PetscNew(&shell));
140   snes->data = (void *)shell;
141   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESShellSetSolve_C", SNESShellSetSolve_Shell));
142   PetscFunctionReturn(PETSC_SUCCESS);
143 }
144