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