Actual source code: snesmfj.c
petsc-3.13.6 2020-09-29
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;
74: MatShellGetContext(J,&j);
75: *snes = (SNES)j->ctx;
76: return(0);
77: }
79: /*
80: MatAssemblyEnd_SNESMF - Calls MatAssemblyEnd_MFFD() and then sets the
81: base from the SNES context
83: */
84: static PetscErrorCode MatAssemblyEnd_SNESMF(Mat J,MatAssemblyType mt)
85: {
87: MatMFFD j;
88: SNES snes;
89: Vec u,f;
92: MatShellGetContext(J,&j);
93: snes = (SNES)j->ctx;
94: MatAssemblyEnd_MFFD(J,mt);
96: SNESGetSolution(snes,&u);
97: if (j->func == (PetscErrorCode (*)(void*,Vec,Vec))SNESComputeFunction) {
98: SNESGetFunction(snes,&f,NULL,NULL);
99: MatMFFDSetBase_MFFD(J,u,f);
100: } else {
101: /* f value known by SNES is not correct for other differencing function */
102: MatMFFDSetBase_MFFD(J,u,NULL);
103: }
104: return(0);
105: }
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()
113: */
114: static PetscErrorCode MatAssemblyEnd_SNESMF_UseBase(Mat J,MatAssemblyType mt)
115: {
117: MatMFFD j;
118: SNES snes;
119: Vec u,f;
122: MatAssemblyEnd_MFFD(J,mt);
123: MatShellGetContext(J,&j);
124: snes = (SNES)j->ctx;
125: SNESGetSolution(snes,&u);
126: SNESGetFunction(snes,&f,NULL,NULL);
127: MatMFFDSetBase_MFFD(J,u,f);
128: return(0);
129: }
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: {
140: MatMFFDSetBase_MFFD(J,U,F);
141: J->ops->assemblyend = MatAssemblyEnd_MFFD;
142: return(0);
143: }
145: static PetscErrorCode MatSNESMFSetReuseBase_SNESMF(Mat J,PetscBool use)
146: {
148: if (use) {
149: J->ops->assemblyend = MatAssemblyEnd_SNESMF_UseBase;
150: } else {
151: J->ops->assemblyend = MatAssemblyEnd_SNESMF;
152: }
153: return(0);
154: }
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().
160: Logically Collective on Mat
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()
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
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.
176: Level: advanced
178: .seealso: MatCreateSNESMF(), MatSNESMFGetReuseBase()
179: @*/
180: PetscErrorCode MatSNESMFSetReuseBase(Mat J,PetscBool use)
181: {
186: PetscTryMethod(J,"MatSNESMFSetReuseBase_C",(Mat,PetscBool),(J,use));
187: return(0);
188: }
190: static PetscErrorCode MatSNESMFGetReuseBase_SNESMF(Mat J,PetscBool *use)
191: {
193: if (J->ops->assemblyend == MatAssemblyEnd_SNESMF_UseBase) *use = PETSC_TRUE;
194: else *use = PETSC_FALSE;
195: return(0);
196: }
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().
202: Logically Collective on Mat
204: Input Parameter:
205: . J - the MatMFFD matrix
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()
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
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.
220: Level: advanced
222: .seealso: MatCreateSNESMF(), MatSNESMFSetReuseBase()
223: @*/
224: PetscErrorCode MatSNESMFGetReuseBase(Mat J,PetscBool *use)
225: {
230: PetscUseMethod(J,"MatSNESMFGetReuseBase_C",(Mat,PetscBool*),(J,use));
231: return(0);
232: }
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.
240: Collective on SNES
242: Input Parameters:
243: . snes - the SNES context
245: Output Parameter:
246: . J - the matrix-free matrix
248: Level: advanced
251: Notes:
252: You can call SNESSetJacobian() with MatMFFDComputeJacobian() if you are using matrix and not a different
253: preconditioner matrix
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.
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().
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.
267: Warning:
268: Using a different function for the differencing will not work if you are using non-linear left preconditioning.
271: .seealso: MatDestroy(), MatMFFDSetFunction(), MatMFFDSetFunctionError(), MatMFFDDSSetUmin()
272: MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateMFFD(),
273: MatMFFDGetH(), MatMFFDRegister(), MatMFFDComputeJacobian(), MatSNESMFSetReuseBase(), MatSNESMFGetReuseBase()
275: @*/
276: PetscErrorCode MatCreateSNESMF(SNES snes,Mat *J)
277: {
279: PetscInt n,N;
280: MatMFFD mf;
283: if (snes->vec_func) {
284: VecGetLocalSize(snes->vec_func,&n);
285: VecGetSize(snes->vec_func,&N);
286: } else if (snes->dm) {
287: Vec tmp;
288: DMGetGlobalVector(snes->dm,&tmp);
289: VecGetLocalSize(tmp,&n);
290: VecGetSize(tmp,&N);
291: DMRestoreGlobalVector(snes->dm,&tmp);
292: } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() or SNESSetDM() first");
293: MatCreateMFFD(PetscObjectComm((PetscObject)snes),n,n,N,N,J);
294: MatShellGetContext(*J,&mf);
295: mf->ctx = snes;
297: if (snes->npc && snes->npcside== PC_LEFT) {
298: MatMFFDSetFunction(*J,(PetscErrorCode (*)(void*,Vec,Vec))SNESComputeFunctionDefaultNPC,snes);
299: } else {
300: MatMFFDSetFunction(*J,(PetscErrorCode (*)(void*,Vec,Vec))SNESComputeFunction,snes);
301: }
303: (*J)->ops->assemblyend = MatAssemblyEnd_SNESMF;
305: PetscObjectComposeFunction((PetscObject)*J,"MatMFFDSetBase_C",MatMFFDSetBase_SNESMF);
306: PetscObjectComposeFunction((PetscObject)*J,"MatSNESMFSetReuseBase_C",MatSNESMFSetReuseBase_SNESMF);
307: PetscObjectComposeFunction((PetscObject)*J,"MatSNESMFGetReuseBase_C",MatSNESMFGetReuseBase_SNESMF);
308: return(0);
309: }