xref: /petsc/src/snes/mf/snesmfj.c (revision a4e35b1925eceef64945ea472b84f2bf06a67b5e)
1 
2 #include <petsc/private/snesimpl.h> /*I  "petscsnes.h" I*/
3 #include <petscdm.h>                /*I  "petscdm.h"   I*/
4 #include <../src/mat/impls/mffd/mffdimpl.h>
5 #include <petsc/private/matimpl.h>
6 
7 /*@C
8   MatMFFDComputeJacobian - Tells the matrix-free Jacobian object the new location at which
9   Jacobian matrix vector products will be computed at, i.e. J(x) * a. The x is obtained
10   from the `SNES` object (using `SNESGetSolution()`).
11 
12   Collective
13 
14   Input Parameters:
15 + snes  - the nonlinear solver context
16 . x     - the point at which the Jacobian vector products will be performed
17 . jac   - the matrix-free Jacobian object of `MatType` `MATMFFD`, likely obtained with `MatCreateSNESMF()`
18 . B     - either the same as jac or another matrix type (ignored)
19 - dummy - the user context (ignored)
20 
21   Options Database Key:
22 . -snes_mf - use the matrix created with `MatSNESMFCreate()` to setup the Jacobian for each new solution in the Newton process
23 
24   Level: developer
25 
26   Notes:
27   If `MatMFFDSetBase()` is ever called on jac then this routine will NO longer get
28   the x from the `SNES` object and `MatMFFDSetBase()` must from that point on be used to
29   change the base vector x.
30 
31   This can be passed into `SNESSetJacobian()` as the Jacobian evaluation function argument
32   when using a completely matrix-free solver,
33   that is the B matrix is also the same matrix operator. This is used when you select
34   -snes_mf but rarely used directly by users. (All this routine does is call `MatAssemblyBegin/End()` on
35   the `Mat` jac.)
36 
37 .seealso: `MatMFFDGetH()`, `MatCreateSNESMF()`, `MatMFFDSetBase()`, `MatCreateMFFD()`, `MATMFFD`,
38           `MatMFFDSetHHistory()`, `MatMFFDSetFunctionError()`, `SNESSetJacobian()`
39 @*/
40 PetscErrorCode MatMFFDComputeJacobian(SNES snes, Vec x, Mat jac, Mat B, void *dummy)
41 {
42   PetscFunctionBegin;
43   PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
44   PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
45   PetscFunctionReturn(PETSC_SUCCESS);
46 }
47 
48 PETSC_EXTERN PetscErrorCode MatAssemblyEnd_MFFD(Mat, MatAssemblyType);
49 PETSC_EXTERN PetscErrorCode MatMFFDSetBase_MFFD(Mat, Vec, Vec);
50 
51 /*@
52   MatSNESMFGetSNES - returns the `SNES` associated with a matrix created with `MatCreateSNESMF()`
53 
54   Not Collective
55 
56   Input Parameter:
57 . J - the matrix
58 
59   Output Parameter:
60 . snes - the `SNES` object
61 
62   Level: advanced
63 
64 .seealso: `Mat`, `SNES`, `MatCreateSNESMF()`
65 @*/
66 PetscErrorCode MatSNESMFGetSNES(Mat J, SNES *snes)
67 {
68   MatMFFD j;
69 
70   PetscFunctionBegin;
71   PetscCall(MatShellGetContext(J, &j));
72   *snes = (SNES)j->ctx;
73   PetscFunctionReturn(PETSC_SUCCESS);
74 }
75 
76 /*
77    MatAssemblyEnd_SNESMF - Calls MatAssemblyEnd_MFFD() and then sets the
78     base from the SNES context
79 
80 */
81 static PetscErrorCode MatAssemblyEnd_SNESMF(Mat J, MatAssemblyType mt)
82 {
83   MatMFFD j;
84   SNES    snes;
85   Vec     u, f;
86   DM      dm;
87   DMSNES  dms;
88 
89   PetscFunctionBegin;
90   PetscCall(MatShellGetContext(J, &j));
91   snes = (SNES)j->ctx;
92   PetscCall(MatAssemblyEnd_MFFD(J, mt));
93 
94   PetscCall(SNESGetSolution(snes, &u));
95   PetscCall(SNESGetDM(snes, &dm));
96   PetscCall(DMGetDMSNES(dm, &dms));
97   if ((j->func == (PetscErrorCode(*)(void *, Vec, Vec))SNESComputeFunction) && !dms->ops->computemffunction) {
98     PetscCall(SNESGetFunction(snes, &f, NULL, NULL));
99     PetscCall(MatMFFDSetBase_MFFD(J, u, f));
100   } else {
101     /* f value known by SNES is not correct for other differencing function */
102     PetscCall(MatMFFDSetBase_MFFD(J, u, NULL));
103   }
104   PetscFunctionReturn(PETSC_SUCCESS);
105 }
106 
107 /*
108    MatAssemblyEnd_SNESMF_UseBase - Calls MatAssemblyEnd_MFFD() and then sets the
109     base from the SNES context. This version will cause the base to be used for differencing
110     even if the func is not SNESComputeFunction. See: MatSNESMFUseBase()
111 
112 */
113 static PetscErrorCode MatAssemblyEnd_SNESMF_UseBase(Mat J, MatAssemblyType mt)
114 {
115   MatMFFD j;
116   SNES    snes;
117   Vec     u, f;
118 
119   PetscFunctionBegin;
120   PetscCall(MatAssemblyEnd_MFFD(J, mt));
121   PetscCall(MatShellGetContext(J, &j));
122   snes = (SNES)j->ctx;
123   PetscCall(SNESGetSolution(snes, &u));
124   PetscCall(SNESGetFunction(snes, &f, NULL, NULL));
125   PetscCall(MatMFFDSetBase_MFFD(J, u, f));
126   PetscFunctionReturn(PETSC_SUCCESS);
127 }
128 
129 /*
130     This routine resets the MatAssemblyEnd() for the MatMFFD created from MatCreateSNESMF() so that it NO longer
131   uses the solution in the SNES object to update the base. See the warning in MatCreateSNESMF().
132 */
133 static PetscErrorCode MatMFFDSetBase_SNESMF(Mat J, Vec U, Vec F)
134 {
135   PetscFunctionBegin;
136   PetscCall(MatMFFDSetBase_MFFD(J, U, F));
137   J->ops->assemblyend = MatAssemblyEnd_MFFD;
138   PetscFunctionReturn(PETSC_SUCCESS);
139 }
140 
141 static PetscErrorCode MatSNESMFSetReuseBase_SNESMF(Mat J, PetscBool use)
142 {
143   PetscFunctionBegin;
144   if (use) {
145     J->ops->assemblyend = MatAssemblyEnd_SNESMF_UseBase;
146   } else {
147     J->ops->assemblyend = MatAssemblyEnd_SNESMF;
148   }
149   PetscFunctionReturn(PETSC_SUCCESS);
150 }
151 
152 /*@
153   MatSNESMFSetReuseBase - Causes the base vector to be used for differencing even if the function provided to `SNESSetFunction()` is not the
154   same as that provided to `MatMFFDSetFunction()`.
155 
156   Logically Collective
157 
158   Input Parameters:
159 + J   - the `MATMFFD` matrix
160 - use - if true always reuse the base vector instead of recomputing f(u) even if the function in the `MATMFFD` is
161           not `SNESComputeFunction()`
162 
163   Level: advanced
164 
165   Note:
166   Care must be taken when using this routine to insure that the function provided to `MatMFFDSetFunction()`, call it F_MF() is compatible with
167   with that provided to `SNESSetFunction()`, call it F_SNES(). That is, (F_MF(u + h*d) - F_SNES(u))/h has to approximate the derivative
168 
169   Developer Notes:
170   This was provided for the MOOSE team who desired to have a `SNESSetFunction()` function that could change configurations (similar to variable
171   switching) to contacts while the function provided to `MatMFFDSetFunction()` cannot. Except for the possibility of changing the configuration
172   both functions compute the same mathematical function so the differencing makes sense.
173 
174 .seealso: `MATMFFD`, `MatMFFDSetFunction()`, `SNESSetFunction()`, `MatCreateSNESMF()`, `MatSNESMFGetReuseBase()`
175 @*/
176 PetscErrorCode MatSNESMFSetReuseBase(Mat J, PetscBool use)
177 {
178   PetscFunctionBegin;
179   PetscValidHeaderSpecific(J, MAT_CLASSID, 1);
180   PetscTryMethod(J, "MatSNESMFSetReuseBase_C", (Mat, PetscBool), (J, use));
181   PetscFunctionReturn(PETSC_SUCCESS);
182 }
183 
184 static PetscErrorCode MatSNESMFGetReuseBase_SNESMF(Mat J, PetscBool *use)
185 {
186   PetscFunctionBegin;
187   if (J->ops->assemblyend == MatAssemblyEnd_SNESMF_UseBase) *use = PETSC_TRUE;
188   else *use = PETSC_FALSE;
189   PetscFunctionReturn(PETSC_SUCCESS);
190 }
191 
192 /*@
193   MatSNESMFGetReuseBase - Determines if the base vector is to be used for differencing even if the function provided to `SNESSetFunction()` is not the
194   same as that provided to `MatMFFDSetFunction()`.
195 
196   Logically Collective
197 
198   Input Parameter:
199 . J - the `MATMFFD` matrix
200 
201   Output Parameter:
202 . use - if true always reuse the base vector instead of recomputing f(u) even if the function in the `MATMFFD` is
203           not `SNESComputeFunction()`
204 
205   Level: advanced
206 
207   Note:
208   See `MatSNESMFSetReuseBase()`
209 
210 .seealso: `MatSNESMFSetReuseBase()`, `MatCreateSNESMF()`
211 @*/
212 PetscErrorCode MatSNESMFGetReuseBase(Mat J, PetscBool *use)
213 {
214   PetscFunctionBegin;
215   PetscValidHeaderSpecific(J, MAT_CLASSID, 1);
216   PetscUseMethod(J, "MatSNESMFGetReuseBase_C", (Mat, PetscBool *), (J, use));
217   PetscFunctionReturn(PETSC_SUCCESS);
218 }
219 
220 /*@
221   MatCreateSNESMF - Creates a matrix-free matrix context for use with
222   a `SNES` solver.  This matrix can be used as the Jacobian argument for
223   the routine `SNESSetJacobian()`. See `MatCreateMFFD()` for details on how
224   the finite difference computation is done.
225 
226   Collective
227 
228   Input Parameters:
229 . snes - the `SNES` context
230 
231   Output Parameter:
232 . J - the matrix-free matrix which is of type `MATMFFD`
233 
234   Level: advanced
235 
236   Notes:
237   You can call `SNESSetJacobian()` with `MatMFFDComputeJacobian()` if you are using matrix and not a different
238   preconditioner matrix
239 
240   If you wish to provide a different function to do differencing on to compute the matrix-free operator than
241   that provided to `SNESSetFunction()` then call `MatMFFDSetFunction()` with your function after this call.
242 
243   The difference between this routine and `MatCreateMFFD()` is that this matrix
244   automatically gets the current base vector from the `SNES` object and not from an
245   explicit call to `MatMFFDSetBase()`.
246 
247   If `MatMFFDSetBase()` is ever called on jac then this routine will NO longer get
248   the x from the `SNES` object and `MatMFFDSetBase()` must from that point on be used to
249   change the base vector x.
250 
251   Using a different function for the differencing will not work if you are using non-linear left preconditioning.
252 
253 .seealso: `MATMFFD`, `MatDestroy()`, `MatMFFDSetFunction()`, `MatMFFDSetFunctionError()`, `MatMFFDDSSetUmin()`
254           `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `MatCreateMFFD()`,
255           `MatMFFDGetH()`, `MatMFFDRegister()`, `MatMFFDComputeJacobian()`, `MatSNESMFSetReuseBase()`, `MatSNESMFGetReuseBase()`
256 @*/
257 PetscErrorCode MatCreateSNESMF(SNES snes, Mat *J)
258 {
259   PetscInt n, N;
260   MatMFFD  mf;
261 
262   PetscFunctionBegin;
263   if (snes->vec_func) {
264     PetscCall(VecGetLocalSize(snes->vec_func, &n));
265     PetscCall(VecGetSize(snes->vec_func, &N));
266   } else if (snes->dm) {
267     Vec tmp;
268     PetscCall(DMGetGlobalVector(snes->dm, &tmp));
269     PetscCall(VecGetLocalSize(tmp, &n));
270     PetscCall(VecGetSize(tmp, &N));
271     PetscCall(DMRestoreGlobalVector(snes->dm, &tmp));
272   } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() or SNESSetDM() first");
273   PetscCall(MatCreateMFFD(PetscObjectComm((PetscObject)snes), n, n, N, N, J));
274   PetscCall(MatShellGetContext(*J, &mf));
275   mf->ctx = snes;
276 
277   if (snes->npc && snes->npcside == PC_LEFT) {
278     PetscCall(MatMFFDSetFunction(*J, (PetscErrorCode(*)(void *, Vec, Vec))SNESComputeFunctionDefaultNPC, snes));
279   } else {
280     DM     dm;
281     DMSNES dms;
282 
283     PetscCall(SNESGetDM(snes, &dm));
284     PetscCall(DMGetDMSNES(dm, &dms));
285     PetscCall(MatMFFDSetFunction(*J, (PetscErrorCode(*)(void *, Vec, Vec))(dms->ops->computemffunction ? SNESComputeMFFunction : SNESComputeFunction), snes));
286   }
287   (*J)->ops->assemblyend = MatAssemblyEnd_SNESMF;
288 
289   PetscCall(PetscObjectComposeFunction((PetscObject)*J, "MatMFFDSetBase_C", MatMFFDSetBase_SNESMF));
290   PetscCall(PetscObjectComposeFunction((PetscObject)*J, "MatSNESMFSetReuseBase_C", MatSNESMFSetReuseBase_SNESMF));
291   PetscCall(PetscObjectComposeFunction((PetscObject)*J, "MatSNESMFGetReuseBase_C", MatSNESMFGetReuseBase_SNESMF));
292   PetscFunctionReturn(PETSC_SUCCESS);
293 }
294