Actual source code: dmsnes.c

  1: #include <petsc/private/snesimpl.h>
  2: #include <petsc/private/dmimpl.h>

  4: static PetscErrorCode DMSNESUnsetFunctionContext_DMSNES(DMSNES sdm)
  5: {
  6:   PetscFunctionBegin;
  7:   PetscCall(PetscObjectCompose((PetscObject)sdm, "function ctx", NULL));
  8:   sdm->functionctxcontainer = NULL;
  9:   PetscFunctionReturn(PETSC_SUCCESS);
 10: }

 12: static PetscErrorCode DMSNESUnsetJacobianContext_DMSNES(DMSNES sdm)
 13: {
 14:   PetscFunctionBegin;
 15:   PetscCall(PetscObjectCompose((PetscObject)sdm, "jacobian ctx", NULL));
 16:   sdm->jacobianctxcontainer = NULL;
 17:   PetscFunctionReturn(PETSC_SUCCESS);
 18: }

 20: static PetscErrorCode DMSNESDestroy(DMSNES *kdm)
 21: {
 22:   PetscFunctionBegin;
 23:   if (!*kdm) PetscFunctionReturn(PETSC_SUCCESS);
 25:   if (--((PetscObject)*kdm)->refct > 0) {
 26:     *kdm = NULL;
 27:     PetscFunctionReturn(PETSC_SUCCESS);
 28:   }
 29:   PetscCall(DMSNESUnsetFunctionContext_DMSNES(*kdm));
 30:   PetscCall(DMSNESUnsetJacobianContext_DMSNES(*kdm));
 31:   PetscTryTypeMethod(*kdm, destroy);
 32:   PetscCall(PetscHeaderDestroy(kdm));
 33:   PetscFunctionReturn(PETSC_SUCCESS);
 34: }

 36: PetscErrorCode DMSNESLoad(DMSNES kdm, PetscViewer viewer)
 37: {
 38:   PetscFunctionBegin;
 39:   PetscCall(PetscViewerBinaryRead(viewer, &kdm->ops->computefunction, 1, NULL, PETSC_FUNCTION));
 40:   PetscCall(PetscViewerBinaryRead(viewer, &kdm->ops->computejacobian, 1, NULL, PETSC_FUNCTION));
 41:   PetscFunctionReturn(PETSC_SUCCESS);
 42: }

 44: PetscErrorCode DMSNESView(DMSNES kdm, PetscViewer viewer)
 45: {
 46:   PetscBool isascii, isbinary;

 48:   PetscFunctionBegin;
 49:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
 50:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
 51:   if (isascii) {
 52: #if defined(PETSC_SERIALIZE_FUNCTIONS)
 53:     const char *fname;

 55:     PetscCall(PetscFPTFind(kdm->ops->computefunction, &fname));
 56:     if (fname) PetscCall(PetscViewerASCIIPrintf(viewer, "Function used by SNES: %s\n", fname));
 57:     PetscCall(PetscFPTFind(kdm->ops->computejacobian, &fname));
 58:     if (fname) PetscCall(PetscViewerASCIIPrintf(viewer, "Jacobian function used by SNES: %s\n", fname));
 59: #endif
 60:   } else if (isbinary) {
 61:     struct {
 62:       SNESFunctionFn *func;
 63:     } funcstruct;
 64:     struct {
 65:       SNESJacobianFn *jac;
 66:     } jacstruct;
 67:     funcstruct.func = kdm->ops->computefunction;
 68:     jacstruct.jac   = kdm->ops->computejacobian;
 69:     PetscCall(PetscViewerBinaryWrite(viewer, &funcstruct, 1, PETSC_FUNCTION));
 70:     PetscCall(PetscViewerBinaryWrite(viewer, &jacstruct, 1, PETSC_FUNCTION));
 71:   }
 72:   PetscFunctionReturn(PETSC_SUCCESS);
 73: }

 75: static PetscErrorCode DMSNESCreate(MPI_Comm comm, DMSNES *kdm)
 76: {
 77:   PetscFunctionBegin;
 78:   PetscCall(SNESInitializePackage());
 79:   PetscCall(PetscHeaderCreate(*kdm, DMSNES_CLASSID, "DMSNES", "DMSNES", "DMSNES", comm, DMSNESDestroy, DMSNESView));
 80:   PetscFunctionReturn(PETSC_SUCCESS);
 81: }

 83: /* Attaches the DMSNES to the coarse level.
 84:  * Under what conditions should we copy versus duplicate?
 85:  */
 86: static PetscErrorCode DMCoarsenHook_DMSNES(DM dm, DM dmc, void *ctx)
 87: {
 88:   PetscFunctionBegin;
 89:   PetscCall(DMCopyDMSNES(dm, dmc));
 90:   PetscFunctionReturn(PETSC_SUCCESS);
 91: }

 93: /* This could restrict auxiliary information to the coarse level.
 94:  */
 95: static PetscErrorCode DMRestrictHook_DMSNES(DM dm, Mat Restrict, Vec rscale, Mat Inject, DM dmc, void *ctx)
 96: {
 97:   PetscFunctionBegin;
 98:   PetscFunctionReturn(PETSC_SUCCESS);
 99: }

101: /* Attaches the DMSNES to the subdomain. */
102: static PetscErrorCode DMSubDomainHook_DMSNES(DM dm, DM subdm, void *ctx)
103: {
104:   PetscFunctionBegin;
105:   PetscCall(DMCopyDMSNES(dm, subdm));
106:   PetscFunctionReturn(PETSC_SUCCESS);
107: }

109: /* This could restrict auxiliary information to the coarse level.
110:  */
111: static PetscErrorCode DMSubDomainRestrictHook_DMSNES(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx)
112: {
113:   PetscFunctionBegin;
114:   PetscFunctionReturn(PETSC_SUCCESS);
115: }

117: static PetscErrorCode DMRefineHook_DMSNES(DM dm, DM dmf, void *ctx)
118: {
119:   PetscFunctionBegin;
120:   PetscCall(DMCopyDMSNES(dm, dmf));
121:   PetscFunctionReturn(PETSC_SUCCESS);
122: }

124: /* This could restrict auxiliary information to the coarse level.
125:  */
126: static PetscErrorCode DMInterpolateHook_DMSNES(DM dm, Mat Interp, DM dmf, void *ctx)
127: {
128:   PetscFunctionBegin;
129:   PetscFunctionReturn(PETSC_SUCCESS);
130: }

132: /*
133:   DMSNESCopy - copies the information in a `DMSNES` to another `DMSNES`

135:   Not Collective

137:   Input Parameters:
138: + kdm  - Original `DMSNES`
139: - nkdm - `DMSNES` to receive the data, should have been created with `DMSNESCreate()`

141:   Level: developer

143: .seealso: [](ch_snes), `DMSNES`, `DMSNESCreate()`, `DMSNESDestroy()`
144: */
145: static PetscErrorCode DMSNESCopy(DMSNES kdm, DMSNES nkdm)
146: {
147:   PetscFunctionBegin;
150:   nkdm->ops->computefunction  = kdm->ops->computefunction;
151:   nkdm->ops->computejacobian  = kdm->ops->computejacobian;
152:   nkdm->ops->computegs        = kdm->ops->computegs;
153:   nkdm->ops->computeobjective = kdm->ops->computeobjective;
154:   nkdm->ops->computepjacobian = kdm->ops->computepjacobian;
155:   nkdm->ops->computepfunction = kdm->ops->computepfunction;
156:   nkdm->ops->destroy          = kdm->ops->destroy;
157:   nkdm->ops->duplicate        = kdm->ops->duplicate;

159:   nkdm->gsctx                = kdm->gsctx;
160:   nkdm->pctx                 = kdm->pctx;
161:   nkdm->objectivectx         = kdm->objectivectx;
162:   nkdm->originaldm           = kdm->originaldm;
163:   nkdm->functionctxcontainer = kdm->functionctxcontainer;
164:   nkdm->jacobianctxcontainer = kdm->jacobianctxcontainer;
165:   if (nkdm->functionctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm, "function ctx", (PetscObject)nkdm->functionctxcontainer));
166:   if (nkdm->jacobianctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm, "jacobian ctx", (PetscObject)nkdm->jacobianctxcontainer));

168:   /*
169:   nkdm->fortran_func_pointers[0] = kdm->fortran_func_pointers[0];
170:   nkdm->fortran_func_pointers[1] = kdm->fortran_func_pointers[1];
171:   nkdm->fortran_func_pointers[2] = kdm->fortran_func_pointers[2];
172:   */

174:   /* implementation specific copy hooks */
175:   PetscTryTypeMethod(kdm, duplicate, nkdm);
176:   PetscFunctionReturn(PETSC_SUCCESS);
177: }

179: /*@C
180:   DMGetDMSNES - get read-only private `DMSNES` context from a `DM`

182:   Not Collective

184:   Input Parameter:
185: . dm - `DM` to be used with `SNES`

187:   Output Parameter:
188: . snesdm - private `DMSNES` context

190:   Level: developer

192:   Note:
193:   Use `DMGetDMSNESWrite()` if write access is needed. The DMSNESSetXXX API should be used wherever possible.

195: .seealso: [](ch_snes), `DMSNES`, `DMGetDMSNESWrite()`
196: @*/
197: PetscErrorCode DMGetDMSNES(DM dm, DMSNES *snesdm)
198: {
199:   PetscFunctionBegin;
201:   *snesdm = (DMSNES)dm->dmsnes;
202:   if (!*snesdm) {
203:     PetscCall(PetscInfo(dm, "Creating new DMSNES\n"));
204:     PetscCall(DMSNESCreate(PetscObjectComm((PetscObject)dm), snesdm));

206:     dm->dmsnes            = (PetscObject)*snesdm;
207:     (*snesdm)->originaldm = dm;
208:     PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_DMSNES, DMRestrictHook_DMSNES, NULL));
209:     PetscCall(DMRefineHookAdd(dm, DMRefineHook_DMSNES, DMInterpolateHook_DMSNES, NULL));
210:     PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_DMSNES, DMSubDomainRestrictHook_DMSNES, NULL));
211:   }
212:   PetscFunctionReturn(PETSC_SUCCESS);
213: }

215: /*@C
216:   DMGetDMSNESWrite - get write access to private `DMSNES` context from a `DM`

218:   Not Collective

220:   Input Parameter:
221: . dm - `DM` to be used with `SNES`

223:   Output Parameter:
224: . snesdm - private `DMSNES` context

226:   Level: developer

228: .seealso: [](ch_snes), `DMSNES`, `DMGetDMSNES()`
229: @*/
230: PetscErrorCode DMGetDMSNESWrite(DM dm, DMSNES *snesdm)
231: {
232:   DMSNES sdm;

234:   PetscFunctionBegin;
236:   PetscCall(DMGetDMSNES(dm, &sdm));
237:   PetscCheck(sdm->originaldm, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "DMSNES has a NULL originaldm");
238:   if (sdm->originaldm != dm) { /* Copy on write */
239:     DMSNES oldsdm = sdm;
240:     PetscCall(PetscInfo(dm, "Copying DMSNES due to write\n"));
241:     PetscCall(DMSNESCreate(PetscObjectComm((PetscObject)dm), &sdm));
242:     PetscCall(DMSNESCopy(oldsdm, sdm));
243:     PetscCall(DMSNESDestroy((DMSNES *)&dm->dmsnes));
244:     dm->dmsnes      = (PetscObject)sdm;
245:     sdm->originaldm = dm;
246:   }
247:   *snesdm = sdm;
248:   PetscFunctionReturn(PETSC_SUCCESS);
249: }

251: /*@C
252:   DMCopyDMSNES - copies a `DMSNES` context to a new `DM`

254:   Logically Collective

256:   Input Parameters:
257: + dmsrc  - `DM` to obtain context from
258: - dmdest - `DM` to add context to

260:   Level: developer

262:   Note:
263:   The context is copied by reference. This function does not ensure that a context exists.

265: .seealso: [](ch_snes), `DMSNES`, `DMGetDMSNES()`, `SNESSetDM()`
266: @*/
267: PetscErrorCode DMCopyDMSNES(DM dmsrc, DM dmdest)
268: {
269:   PetscFunctionBegin;
272:   if (!dmdest->dmsnes) PetscCall(DMSNESCreate(PetscObjectComm((PetscObject)dmdest), (DMSNES *)&dmdest->dmsnes));
273:   PetscCall(DMSNESCopy((DMSNES)dmsrc->dmsnes, (DMSNES)dmdest->dmsnes));
274:   PetscCall(DMCoarsenHookAdd(dmdest, DMCoarsenHook_DMSNES, NULL, NULL));
275:   PetscCall(DMRefineHookAdd(dmdest, DMRefineHook_DMSNES, NULL, NULL));
276:   PetscCall(DMSubDomainHookAdd(dmdest, DMSubDomainHook_DMSNES, DMSubDomainRestrictHook_DMSNES, NULL));
277:   PetscFunctionReturn(PETSC_SUCCESS);
278: }

280: /*@C
281:   DMSNESSetFunction - set `SNES` residual evaluation function

283:   Not Collective

285:   Input Parameters:
286: + dm  - DM to be used with `SNES`
287: . f   - residual evaluation function; see `SNESFunctionFn` for calling sequence
288: - ctx - context for residual evaluation

290:   Level: developer

292:   Note:
293:   `SNESSetFunction()` is normally used, but it calls this function internally because the user context is actually
294:   associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
295:   not.

297:   Developer Note:
298:   If `DM` took a more central role at some later date, this could become the primary method of setting the residual.

300: .seealso: [](ch_snes), `DMSNES`, `DMSNESSetContext()`, `SNESSetFunction()`, `DMSNESSetJacobian()`, `SNESFunctionFn`
301: @*/
302: PetscErrorCode DMSNESSetFunction(DM dm, SNESFunctionFn *f, void *ctx)
303: {
304:   DMSNES sdm;

306:   PetscFunctionBegin;
308:   PetscCall(DMGetDMSNESWrite(dm, &sdm));
309:   if (f) sdm->ops->computefunction = f;
310:   if (ctx) {
311:     PetscContainer ctxcontainer;
312:     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)sdm), &ctxcontainer));
313:     PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
314:     PetscCall(PetscObjectCompose((PetscObject)sdm, "function ctx", (PetscObject)ctxcontainer));
315:     sdm->functionctxcontainer = ctxcontainer;
316:     PetscCall(PetscContainerDestroy(&ctxcontainer));
317:   }
318:   PetscFunctionReturn(PETSC_SUCCESS);
319: }

321: /*@C
322:   DMSNESSetFunctionContextDestroy - set `SNES` residual evaluation context destroy function

324:   Not Collective

326:   Input Parameters:
327: + dm - `DM` to be used with `SNES`
328: - f  - residual evaluation context destroy function

330:   Level: developer

332: .seealso: [](ch_snes), `DMSNES`, `DMSNESSetFunction()`, `SNESSetFunction()`
333: @*/
334: PetscErrorCode DMSNESSetFunctionContextDestroy(DM dm, PetscErrorCode (*f)(void *))
335: {
336:   DMSNES sdm;

338:   PetscFunctionBegin;
340:   PetscCall(DMGetDMSNESWrite(dm, &sdm));
341:   if (sdm->functionctxcontainer) PetscCall(PetscContainerSetUserDestroy(sdm->functionctxcontainer, f));
342:   PetscFunctionReturn(PETSC_SUCCESS);
343: }

345: PetscErrorCode DMSNESUnsetFunctionContext_Internal(DM dm)
346: {
347:   DMSNES sdm;

349:   PetscFunctionBegin;
351:   PetscCall(DMGetDMSNESWrite(dm, &sdm));
352:   PetscCall(DMSNESUnsetFunctionContext_DMSNES(sdm));
353:   PetscFunctionReturn(PETSC_SUCCESS);
354: }

356: /*@C
357:   DMSNESSetMFFunction - set `SNES` residual evaluation function used in applying the matrix-free Jacobian with `-snes_mf_operator`

359:   Logically Collective

361:   Input Parameters:
362: + dm   - `DM` to be used with `SNES`
363: . func - residual evaluation function; see `SNESFunctionFn` for calling sequence
364: - ctx  - optional function context

366:   Level: developer

368:   Note:
369:   If not provided then the function provided with `SNESSetFunction()` is used

371: .seealso: [](ch_snes), `DMSNES`, `DMSNESSetContext()`, `SNESSetFunction()`, `DMSNESSetJacobian()`, `SNESFunction`, `DMSNESSetFunction()`, `SNESFunctionFn`
372: @*/
373: PetscErrorCode DMSNESSetMFFunction(DM dm, SNESFunctionFn *func, void *ctx)
374: {
375:   DMSNES sdm;

377:   PetscFunctionBegin;
379:   if (func || ctx) PetscCall(DMGetDMSNESWrite(dm, &sdm));
380:   if (func) sdm->ops->computemffunction = func;
381:   if (ctx) sdm->mffunctionctx = ctx;
382:   PetscFunctionReturn(PETSC_SUCCESS);
383: }

385: /*@C
386:   DMSNESGetFunction - get `SNES` residual evaluation function from a `DMSNES` object

388:   Not Collective

390:   Input Parameter:
391: . dm - `DM` to be used with `SNES`

393:   Output Parameters:
394: + f   - residual evaluation function; see `SNESFunctionFn` for calling sequence
395: - ctx - context for residual evaluation

397:   Level: developer

399:   Note:
400:   `SNESGetFunction()` is normally used, but it calls this function internally because the user context is actually
401:   associated with the `DM`.

403: .seealso: [](ch_snes), `DMSNES`, `DMSNESSetContext()`, `DMSNESSetFunction()`, `SNESSetFunction()`, `SNESFunctionFn`
404: @*/
405: PetscErrorCode DMSNESGetFunction(DM dm, SNESFunctionFn **f, void **ctx)
406: {
407:   DMSNES sdm;

409:   PetscFunctionBegin;
411:   PetscCall(DMGetDMSNES(dm, &sdm));
412:   if (f) *f = sdm->ops->computefunction;
413:   if (ctx) {
414:     if (sdm->functionctxcontainer) PetscCall(PetscContainerGetPointer(sdm->functionctxcontainer, ctx));
415:     else *ctx = NULL;
416:   }
417:   PetscFunctionReturn(PETSC_SUCCESS);
418: }

420: /*@C
421:   DMSNESSetObjective - Sets the objective function minimized by some of the `SNES` linesearch methods into a `DMSNES` object, used instead of the 2-norm of the residual

423:   Not Collective

425:   Input Parameters:
426: + dm  - `DM` to be used with `SNES`
427: . obj - objective evaluation routine; see `SNESObjectiveFn` for the calling sequence
428: - ctx - [optional] user-defined context for private data for the objective evaluation routine (may be `NULL`)

430:   Level: developer

432: .seealso: [](ch_snes), `DMSNES`, `DMSNESSetContext()`, `SNESGetObjective()`, `DMSNESSetFunction()`, `SNESObjectiveFn`
433: @*/
434: PetscErrorCode DMSNESSetObjective(DM dm, SNESObjectiveFn *obj, void *ctx)
435: {
436:   DMSNES sdm;

438:   PetscFunctionBegin;
440:   if (obj || ctx) PetscCall(DMGetDMSNESWrite(dm, &sdm));
441:   if (obj) sdm->ops->computeobjective = obj;
442:   if (ctx) sdm->objectivectx = ctx;
443:   PetscFunctionReturn(PETSC_SUCCESS);
444: }

446: /*@C
447:   DMSNESGetObjective - Returns the objective function set with `DMSNESSetObjective()`

449:   Not Collective

451:   Input Parameter:
452: . dm - `DM` to be used with `SNES`

454:   Output Parameters:
455: + obj - objective evaluation routine (or `NULL`); see `SNESObjectiveFn` for the calling sequence
456: - ctx - the function context (or `NULL`)

458:   Level: developer

460:   Note:
461:   `SNESGetFunction()` is normally used, but it calls this function internally because the user context is actually
462:   associated with the `DM`.

464: .seealso: [](ch_snes), `DMSNES`, `DMSNESSetContext()`, `DMSNESSetObjective()`, `SNESSetFunction()`, `SNESObjectiveFn`
465: @*/
466: PetscErrorCode DMSNESGetObjective(DM dm, SNESObjectiveFn **obj, void **ctx)
467: {
468:   DMSNES sdm;

470:   PetscFunctionBegin;
472:   PetscCall(DMGetDMSNES(dm, &sdm));
473:   if (obj) *obj = sdm->ops->computeobjective;
474:   if (ctx) *ctx = sdm->objectivectx;
475:   PetscFunctionReturn(PETSC_SUCCESS);
476: }

478: /*@C
479:   DMSNESSetNGS - set `SNES` Gauss-Seidel relaxation function into a `DMSNES` object

481:   Not Collective

483:   Input Parameters:
484: + dm  - `DM` to be used with `SNES`
485: . f   - relaxation function, see `SNESGSFunction`
486: - ctx - context for residual evaluation

488:   Level: developer

490:   Note:
491:   `SNESSetNGS()` is normally used, but it calls this function internally because the user context is actually
492:   associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
493:   not.

495:   Developer Note:
496:   If `DM` took a more central role at some later date, this could become the primary method of supplying the smoother

498: .seealso: [](ch_snes), `DMSNES`, `DMSNESSetContext()`, `SNESSetFunction()`, `DMSNESSetJacobian()`, `DMSNESSetFunction()`, `SNESGSFunction`
499: @*/
500: PetscErrorCode DMSNESSetNGS(DM dm, PetscErrorCode (*f)(SNES, Vec, Vec, void *), void *ctx)
501: {
502:   DMSNES sdm;

504:   PetscFunctionBegin;
506:   if (f || ctx) PetscCall(DMGetDMSNESWrite(dm, &sdm));
507:   if (f) sdm->ops->computegs = f;
508:   if (ctx) sdm->gsctx = ctx;
509:   PetscFunctionReturn(PETSC_SUCCESS);
510: }

512: /*@C
513:   DMSNESGetNGS - get `SNES` Gauss-Seidel relaxation function from a `DMSNES` object

515:   Not Collective

517:   Input Parameter:
518: . dm - `DM` to be used with `SNES`

520:   Output Parameters:
521: + f   - relaxation function which performs Gauss-Seidel sweeps, see `SNESSetNGS()`
522: - ctx - context for residual evaluation

524:   Level: developer

526:   Note:
527:   `SNESGetNGS()` is normally used, but it calls this function internally because the user context is actually
528:   associated with the `DM`.

530:   Developer Note:
531:   This makes the interface consistent regardless of whether the user interacts with a `DM` or
532:   not. If `DM` took a more central role at some later date, this could become the primary method of setting the residual.

534: .seealso: [](ch_snes), `DMSNES`, `DMSNESSetContext()`, `SNESGetNGS()`, `DMSNESGetJacobian()`, `DMSNESGetFunction()`
535: @*/
536: PetscErrorCode DMSNESGetNGS(DM dm, PetscErrorCode (**f)(SNES, Vec, Vec, void *), void **ctx)
537: {
538:   DMSNES sdm;

540:   PetscFunctionBegin;
542:   PetscCall(DMGetDMSNES(dm, &sdm));
543:   if (f) *f = sdm->ops->computegs;
544:   if (ctx) *ctx = sdm->gsctx;
545:   PetscFunctionReturn(PETSC_SUCCESS);
546: }

548: /*@C
549:   DMSNESSetJacobian - set `SNES` Jacobian evaluation function into a `DMSNES` object

551:   Not Collective

553:   Input Parameters:
554: + dm  - `DM` to be used with `SNES`
555: . J   - Jacobian evaluation function, see `SNESJacobianFn`
556: - ctx - context for Jacobian evaluation

558:   Level: developer

560:   Note:
561:   `SNESSetJacobian()` is normally used, but it calls this function internally because the user context is actually
562:   associated with the `DM`.

564:   Developer Note:
565:   This makes the interface consistent regardless of whether the user interacts with a `DM` or
566:   not. If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian.

568: .seealso: [](ch_snes), `DMSNES`, `DMSNESSetContext()`, `SNESSetFunction()`, `DMSNESGetJacobian()`, `SNESSetJacobian()`, `SNESJacobianFn`
569: @*/
570: PetscErrorCode DMSNESSetJacobian(DM dm, SNESJacobianFn *J, void *ctx)
571: {
572:   DMSNES sdm;

574:   PetscFunctionBegin;
576:   if (J || ctx) PetscCall(DMGetDMSNESWrite(dm, &sdm));
577:   if (J) sdm->ops->computejacobian = J;
578:   if (ctx) {
579:     PetscContainer ctxcontainer;
580:     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)sdm), &ctxcontainer));
581:     PetscCall(PetscContainerSetPointer(ctxcontainer, ctx));
582:     PetscCall(PetscObjectCompose((PetscObject)sdm, "jacobian ctx", (PetscObject)ctxcontainer));
583:     sdm->jacobianctxcontainer = ctxcontainer;
584:     PetscCall(PetscContainerDestroy(&ctxcontainer));
585:   }
586:   PetscFunctionReturn(PETSC_SUCCESS);
587: }

589: /*@C
590:   DMSNESSetJacobianContextDestroy - set `SNES` Jacobian evaluation context destroy function into a `DMSNES` object

592:   Not Collective

594:   Input Parameters:
595: + dm - `DM` to be used with `SNES`
596: - f  - Jacobian evaluation context destroy function

598:   Level: developer

600: .seealso: [](ch_snes), `DMSNES`, `DMSNESSetJacobian()`
601: @*/
602: PetscErrorCode DMSNESSetJacobianContextDestroy(DM dm, PetscErrorCode (*f)(void *))
603: {
604:   DMSNES sdm;

606:   PetscFunctionBegin;
608:   PetscCall(DMGetDMSNESWrite(dm, &sdm));
609:   if (sdm->jacobianctxcontainer) PetscCall(PetscContainerSetUserDestroy(sdm->jacobianctxcontainer, f));
610:   PetscFunctionReturn(PETSC_SUCCESS);
611: }

613: PetscErrorCode DMSNESUnsetJacobianContext_Internal(DM dm)
614: {
615:   DMSNES sdm;

617:   PetscFunctionBegin;
619:   PetscCall(DMGetDMSNESWrite(dm, &sdm));
620:   PetscCall(DMSNESUnsetJacobianContext_DMSNES(sdm));
621:   PetscFunctionReturn(PETSC_SUCCESS);
622: }

624: /*@C
625:   DMSNESGetJacobian - get `SNES` Jacobian evaluation function from a `DMSNES` object

627:   Not Collective

629:   Input Parameter:
630: . dm - `DM` to be used with `SNES`

632:   Output Parameters:
633: + J   - Jacobian evaluation function; for all calling sequence see `SNESJacobianFn`
634: - ctx - context for residual evaluation

636:   Level: developer

638:   Note:
639:   `SNESGetJacobian()` is normally used, but it calls this function internally because the user context is actually
640:   associated with the `DM`.  This makes the interface consistent regardless of whether the user interacts with a `DM` or
641:   not.

643:   Developer Note:
644:   If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian.

646: .seealso: [](ch_snes), `DMSNES`, `DMSNESSetContext()`, `SNESSetFunction()`, `DMSNESSetJacobian()`, `SNESJacobianFn`
647: @*/
648: PetscErrorCode DMSNESGetJacobian(DM dm, SNESJacobianFn **J, void **ctx)
649: {
650:   DMSNES sdm;

652:   PetscFunctionBegin;
654:   PetscCall(DMGetDMSNES(dm, &sdm));
655:   if (J) *J = sdm->ops->computejacobian;
656:   if (ctx) {
657:     if (sdm->jacobianctxcontainer) PetscCall(PetscContainerGetPointer(sdm->jacobianctxcontainer, ctx));
658:     else *ctx = NULL;
659:   }
660:   PetscFunctionReturn(PETSC_SUCCESS);
661: }

663: /*@C
664:   DMSNESSetPicard - set SNES Picard iteration matrix and RHS evaluation functions into a `DMSNES` object

666:   Not Collective

668:   Input Parameters:
669: + dm  - `DM` to be used with `SNES`
670: . b   - RHS evaluation function; see `SNESFunctionFn` for calling sequence
671: . J   - Picard matrix evaluation function; see `SNESJacobianFn` for calling sequence
672: - ctx - context for residual and matrix evaluation

674:   Level: developer

676: .seealso: [](ch_snes), `DMSNES`, `SNESSetPicard()`, `DMSNESSetFunction()`, `DMSNESSetJacobian()`, `SNESFunctionFn`, `SNESJacobianFn`
677: @*/
678: PetscErrorCode DMSNESSetPicard(DM dm, SNESFunctionFn *b, SNESJacobianFn *J, void *ctx)
679: {
680:   DMSNES sdm;

682:   PetscFunctionBegin;
684:   PetscCall(DMGetDMSNES(dm, &sdm));
685:   if (b) sdm->ops->computepfunction = b;
686:   if (J) sdm->ops->computepjacobian = J;
687:   if (ctx) sdm->pctx = ctx;
688:   PetscFunctionReturn(PETSC_SUCCESS);
689: }

691: /*@C
692:   DMSNESGetPicard - get `SNES` Picard iteration evaluation functions from a `DMSNES` object

694:   Not Collective

696:   Input Parameter:
697: . dm - `DM` to be used with `SNES`

699:   Output Parameters:
700: + b   - RHS evaluation function; see `SNESFunctionFn` for calling sequence
701: . J   - Jacobian evaluation function; see `SNESJacobianFn` for calling sequence
702: - ctx - context for residual and matrix evaluation

704:   Level: developer

706: .seealso: [](ch_snes), `DMSNES`, `DMSNESSetContext()`, `SNESSetFunction()`, `DMSNESSetJacobian()`, `SNESFunctionFn`, `SNESJacobianFn`
707: @*/
708: PetscErrorCode DMSNESGetPicard(DM dm, SNESFunctionFn **b, SNESJacobianFn **J, void **ctx)
709: {
710:   DMSNES sdm;

712:   PetscFunctionBegin;
714:   PetscCall(DMGetDMSNES(dm, &sdm));
715:   if (b) *b = sdm->ops->computepfunction;
716:   if (J) *J = sdm->ops->computepjacobian;
717:   if (ctx) *ctx = sdm->pctx;
718:   PetscFunctionReturn(PETSC_SUCCESS);
719: }