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