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