Actual source code: snesmfj.c
petsc-3.10.5 2019-03-28
2: #include <petsc/private/snesimpl.h>
3: #include <petscdm.h>
4: #include <../src/mat/impls/mffd/mffdimpl.h>
5: #include <petsc/private/matimpl.h>
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()).
12: Logically Collective on SNES
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)
22: Level: developer
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.
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.)
36: .seealso: MatMFFDGetH(), MatCreateSNESMF(), MatCreateMFFD(), MATMFFD,
37: MatMFFDSetHHistory(), MatMFFDSetFunctionError(), MatCreateMFFD(), SNESSetJacobian()
39: @*/
40: PetscErrorCode MatMFFDComputeJacobian(SNES snes,Vec x,Mat jac,Mat B,void *dummy)
41: {
45: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
46: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
47: return(0);
48: }
50: PETSC_EXTERN PetscErrorCode MatAssemblyEnd_MFFD(Mat,MatAssemblyType);
51: PETSC_EXTERN PetscErrorCode MatMFFDSetBase_MFFD(Mat,Vec,Vec);
53: /*@
54: MatSNESMFGetSNES - returns the SNES associated with a matrix created with MatCreateSNESMF()
56: Not collective
58: Input Parameter:
59: . J - the matrix
61: Output Parameter:
62: . snes - the SNES object
64: Level: advanced
66: .seealso: MatCreateSNESMF()
67: @*/
68: PetscErrorCode MatSNESMFGetSNES(Mat J,SNES *snes)
69: {
70: MatMFFD j = (MatMFFD)J->data;
73: *snes = (SNES)j->ctx;
74: return(0);
75: }
77: /*
78: MatAssemblyEnd_SNESMF - Calls MatAssemblyEnd_MFFD() and then sets the
79: base from the SNES context
81: */
82: static PetscErrorCode MatAssemblyEnd_SNESMF(Mat J,MatAssemblyType mt)
83: {
85: MatMFFD j = (MatMFFD)J->data;
86: SNES snes = (SNES)j->ctx;
87: Vec u,f;
90: MatAssemblyEnd_MFFD(J,mt);
92: SNESGetSolution(snes,&u);
93: if (j->func == (PetscErrorCode (*)(void*,Vec,Vec))SNESComputeFunction) {
94: SNESGetFunction(snes,&f,NULL,NULL);
95: MatMFFDSetBase_MFFD(J,u,f);
96: } else {
97: /* f value known by SNES is not correct for other differencing function */
98: MatMFFDSetBase_MFFD(J,u,NULL);
99: }
100: return(0);
101: }
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()
109: */
110: static PetscErrorCode MatAssemblyEnd_SNESMF_UseBase(Mat J,MatAssemblyType mt)
111: {
113: MatMFFD j = (MatMFFD)J->data;
114: SNES snes = (SNES)j->ctx;
115: Vec u,f;
118: MatAssemblyEnd_MFFD(J,mt);
120: SNESGetSolution(snes,&u);
121: SNESGetFunction(snes,&f,NULL,NULL);
122: MatMFFDSetBase_MFFD(J,u,f);
123: return(0);
124: }
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: {
135: MatMFFDSetBase_MFFD(J,U,F);
136: J->ops->assemblyend = MatAssemblyEnd_MFFD;
137: return(0);
138: }
140: static PetscErrorCode MatSNESMFSetReuseBase_SNESMF(Mat J,PetscBool use)
141: {
143: if (use) {
144: J->ops->assemblyend = MatAssemblyEnd_SNESMF_UseBase;
145: } else {
146: J->ops->assemblyend = MatAssemblyEnd_SNESMF;
147: }
148: return(0);
149: }
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().
155: Logically Collective on Mat
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()
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
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.
171: Level: advanced
173: .seealso: MatCreateSNESMF(), MatSNESMFGetReuseBase()
174: @*/
175: PetscErrorCode MatSNESMFSetReuseBase(Mat J,PetscBool use)
176: {
181: PetscTryMethod(J,"MatSNESMFSetReuseBase_C",(Mat,PetscBool),(J,use));
182: return(0);
183: }
185: static PetscErrorCode MatSNESMFGetReuseBase_SNESMF(Mat J,PetscBool *use)
186: {
188: if (J->ops->assemblyend == MatAssemblyEnd_SNESMF_UseBase) *use = PETSC_TRUE;
189: else *use = PETSC_FALSE;
190: return(0);
191: }
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().
197: Logically Collective on Mat
199: Input Parameter:
200: . J - the MatMFFD matrix
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()
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
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.
215: Level: advanced
217: .seealso: MatCreateSNESMF(), MatSNESMFSetReuseBase()
218: @*/
219: PetscErrorCode MatSNESMFGetReuseBase(Mat J,PetscBool *use)
220: {
225: PetscUseMethod(J,"MatSNESMFGetReuseBase_C",(Mat,PetscBool*),(J,use));
226: return(0);
227: }
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.
235: Collective on SNES and Vec
237: Input Parameters:
238: . snes - the SNES context
240: Output Parameter:
241: . J - the matrix-free matrix
243: Level: advanced
246: Notes:
247: You can call SNESSetJacobian() with MatMFFDComputeJacobian() if you are using matrix and not a different
248: preconditioner matrix
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.
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().
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.
262: Warning:
263: Using a different function for the differencing will not work if you are using non-linear left preconditioning.
266: .seealso: MatDestroy(), MatMFFDSetFunction(), MatMFFDSetFunctionError(), MatMFFDDSSetUmin()
267: MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateMFFD(),
268: MatMFFDGetH(), MatMFFDRegister(), MatMFFDComputeJacobian(), MatSNESMFSetReuseBase(), MatSNESMFGetReuseBase()
270: @*/
271: PetscErrorCode MatCreateSNESMF(SNES snes,Mat *J)
272: {
274: PetscInt n,N;
275: MatMFFD mf;
278: if (snes->vec_func) {
279: VecGetLocalSize(snes->vec_func,&n);
280: VecGetSize(snes->vec_func,&N);
281: } else if (snes->dm) {
282: Vec tmp;
283: DMGetGlobalVector(snes->dm,&tmp);
284: VecGetLocalSize(tmp,&n);
285: VecGetSize(tmp,&N);
286: DMRestoreGlobalVector(snes->dm,&tmp);
287: } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() or SNESSetDM() first");
288: MatCreateMFFD(PetscObjectComm((PetscObject)snes),n,n,N,N,J);
290: mf = (MatMFFD)(*J)->data;
291: mf->ctx = snes;
293: if (snes->npc && snes->npcside== PC_LEFT) {
294: MatMFFDSetFunction(*J,(PetscErrorCode (*)(void*,Vec,Vec))SNESComputeFunctionDefaultNPC,snes);
295: } else {
296: MatMFFDSetFunction(*J,(PetscErrorCode (*)(void*,Vec,Vec))SNESComputeFunction,snes);
297: }
299: (*J)->ops->assemblyend = MatAssemblyEnd_SNESMF;
301: PetscObjectComposeFunction((PetscObject)*J,"MatMFFDSetBase_C",MatMFFDSetBase_SNESMF);
302: PetscObjectComposeFunction((PetscObject)*J,"MatSNESMFSetReuseBase_C",MatSNESMFSetReuseBase_SNESMF);
303: PetscObjectComposeFunction((PetscObject)*J,"MatSNESMFGetReuseBase_C",MatSNESMFGetReuseBase_SNESMF);
304: return(0);
305: }