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