Actual source code: dmsnes.c
petsc-3.12.5 2020-03-29
1: #include <petsc/private/snesimpl.h>
2: #include <petsc/private/dmimpl.h>
4: static PetscErrorCode DMSNESDestroy(DMSNES *kdm)
5: {
9: if (!*kdm) return(0);
11: if (--((PetscObject)(*kdm))->refct > 0) {*kdm = 0; return(0);}
12: if ((*kdm)->ops->destroy) {((*kdm)->ops->destroy)(*kdm);}
13: PetscHeaderDestroy(kdm);
14: return(0);
15: }
17: PetscErrorCode DMSNESLoad(DMSNES kdm,PetscViewer viewer)
18: {
22: PetscViewerBinaryRead(viewer,&kdm->ops->computefunction,1,NULL,PETSC_FUNCTION);
23: PetscViewerBinaryRead(viewer,&kdm->ops->computejacobian,1,NULL,PETSC_FUNCTION);
24: return(0);
25: }
27: PetscErrorCode DMSNESView(DMSNES kdm,PetscViewer viewer)
28: {
30: PetscBool isascii,isbinary;
33: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
34: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
35: if (isascii) {
36: #if defined(PETSC_SERIALIZE_FUNCTIONS) && defined(PETSC_SERIALIZE_FUNCTIONS_VIEW)
37: const char *fname;
39: PetscFPTFind(kdm->ops->computefunction,&fname);
40: if (fname) {
41: PetscViewerASCIIPrintf(viewer,"Function used by SNES: %s\n",fname);
42: }
43: PetscFPTFind(kdm->ops->computejacobian,&fname);
44: if (fname) {
45: PetscViewerASCIIPrintf(viewer,"Jacobian function used by SNES: %s\n",fname);
46: }
47: #endif
48: } else if (isbinary) {
49: struct {
50: PetscErrorCode (*func)(SNES,Vec,Vec,void*);
51: } funcstruct;
52: struct {
53: PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*);
54: } jacstruct;
55: funcstruct.func = kdm->ops->computefunction;
56: jacstruct.jac = kdm->ops->computejacobian;
57: PetscViewerBinaryWrite(viewer,&funcstruct,1,PETSC_FUNCTION,PETSC_FALSE);
58: PetscViewerBinaryWrite(viewer,&jacstruct,1,PETSC_FUNCTION,PETSC_FALSE);
59: }
60: return(0);
61: }
63: static PetscErrorCode DMSNESCreate(MPI_Comm comm,DMSNES *kdm)
64: {
68: SNESInitializePackage();
69: PetscHeaderCreate(*kdm, DMSNES_CLASSID, "DMSNES", "DMSNES", "DMSNES", comm, DMSNESDestroy, DMSNESView);
70: return(0);
71: }
73: /* Attaches the DMSNES to the coarse level.
74: * Under what conditions should we copy versus duplicate?
75: */
76: static PetscErrorCode DMCoarsenHook_DMSNES(DM dm,DM dmc,void *ctx)
77: {
81: DMCopyDMSNES(dm,dmc);
82: return(0);
83: }
85: /* This could restrict auxiliary information to the coarse level.
86: */
87: static PetscErrorCode DMRestrictHook_DMSNES(DM dm,Mat Restrict,Vec rscale,Mat Inject,DM dmc,void *ctx)
88: {
91: return(0);
92: }
94: /* Attaches the DMSNES to the subdomain. */
95: static PetscErrorCode DMSubDomainHook_DMSNES(DM dm,DM subdm,void *ctx)
96: {
100: DMCopyDMSNES(dm,subdm);
101: return(0);
102: }
104: /* This could restrict auxiliary information to the coarse level.
105: */
106: static PetscErrorCode DMSubDomainRestrictHook_DMSNES(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
107: {
110: return(0);
111: }
113: static PetscErrorCode DMRefineHook_DMSNES(DM dm,DM dmf,void *ctx)
114: {
118: DMCopyDMSNES(dm,dmf);
119: return(0);
120: }
122: /* This could restrict auxiliary information to the coarse level.
123: */
124: static PetscErrorCode DMInterpolateHook_DMSNES(DM dm,Mat Interp,DM dmf,void *ctx)
125: {
128: return(0);
129: }
131: /*@C
132: DMSNESCopy - copies the information in a DMSNES to another DMSNES
134: Not Collective
136: Input Argument:
137: + kdm - Original DMSNES
138: - nkdm - DMSNES to receive the data, should have been created with DMSNESCreate()
140: Level: developer
142: .seealso: DMSNESCreate(), DMSNESDestroy()
143: @*/
144: PetscErrorCode DMSNESCopy(DMSNES kdm,DMSNES nkdm)
145: {
151: nkdm->ops->computefunction = kdm->ops->computefunction;
152: nkdm->ops->computejacobian = kdm->ops->computejacobian;
153: nkdm->ops->computegs = kdm->ops->computegs;
154: nkdm->ops->computeobjective = kdm->ops->computeobjective;
155: nkdm->ops->computepjacobian = kdm->ops->computepjacobian;
156: nkdm->ops->computepfunction = kdm->ops->computepfunction;
157: nkdm->ops->destroy = kdm->ops->destroy;
158: nkdm->ops->duplicate = kdm->ops->duplicate;
160: nkdm->functionctx = kdm->functionctx;
161: nkdm->gsctx = kdm->gsctx;
162: nkdm->pctx = kdm->pctx;
163: nkdm->jacobianctx = kdm->jacobianctx;
164: nkdm->objectivectx = kdm->objectivectx;
165: nkdm->originaldm = kdm->originaldm;
167: /*
168: nkdm->fortran_func_pointers[0] = kdm->fortran_func_pointers[0];
169: nkdm->fortran_func_pointers[1] = kdm->fortran_func_pointers[1];
170: nkdm->fortran_func_pointers[2] = kdm->fortran_func_pointers[2];
171: */
173: /* implementation specific copy hooks */
174: if (kdm->ops->duplicate) {(*kdm->ops->duplicate)(kdm,nkdm);}
175: return(0);
176: }
178: /*@C
179: DMGetDMSNES - get read-only private DMSNES context from a DM
181: Not Collective
183: Input Argument:
184: . dm - DM to be used with SNES
186: Output Argument:
187: . snesdm - private DMSNES context
189: Level: developer
191: Notes:
192: Use DMGetDMSNESWrite() if write access is needed. The DMSNESSetXXX API should be used wherever possible.
194: .seealso: DMGetDMSNESWrite()
195: @*/
196: PetscErrorCode DMGetDMSNES(DM dm,DMSNES *snesdm)
197: {
202: *snesdm = (DMSNES) dm->dmsnes;
203: if (!*snesdm) {
204: PetscInfo(dm,"Creating new DMSNES\n");
205: DMSNESCreate(PetscObjectComm((PetscObject)dm),snesdm);
207: dm->dmsnes = (PetscObject) *snesdm;
208: (*snesdm)->originaldm = dm;
209: DMCoarsenHookAdd(dm,DMCoarsenHook_DMSNES,DMRestrictHook_DMSNES,NULL);
210: DMRefineHookAdd(dm,DMRefineHook_DMSNES,DMInterpolateHook_DMSNES,NULL);
211: DMSubDomainHookAdd(dm,DMSubDomainHook_DMSNES,DMSubDomainRestrictHook_DMSNES,NULL);
212: }
213: return(0);
214: }
216: /*@C
217: DMGetDMSNESWrite - get write access to private DMSNES context from a DM
219: Not Collective
221: Input Argument:
222: . dm - DM to be used with SNES
224: Output Argument:
225: . snesdm - private DMSNES context
227: Level: developer
229: .seealso: DMGetDMSNES()
230: @*/
231: PetscErrorCode DMGetDMSNESWrite(DM dm,DMSNES *snesdm)
232: {
234: DMSNES sdm;
238: DMGetDMSNES(dm,&sdm);
239: if (!sdm->originaldm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"DMSNES has a NULL originaldm");
240: if (sdm->originaldm != dm) { /* Copy on write */
241: DMSNES oldsdm = sdm;
242: PetscInfo(dm,"Copying DMSNES due to write\n");
243: DMSNESCreate(PetscObjectComm((PetscObject)dm),&sdm);
244: DMSNESCopy(oldsdm,sdm);
245: DMSNESDestroy((DMSNES*)&dm->dmsnes);
246: dm->dmsnes = (PetscObject)sdm;
247: sdm->originaldm = dm;
248: }
249: *snesdm = sdm;
250: return(0);
251: }
253: /*@C
254: DMCopyDMSNES - copies a DM context to a new DM
256: Logically Collective
258: Input Arguments:
259: + dmsrc - DM to obtain context from
260: - dmdest - DM to add context to
262: Level: developer
264: Note:
265: The context is copied by reference. This function does not ensure that a context exists.
267: .seealso: DMGetDMSNES(), SNESSetDM()
268: @*/
269: PetscErrorCode DMCopyDMSNES(DM dmsrc,DM dmdest)
270: {
276: if (!dmdest->dmsnes) {DMSNESCreate(PetscObjectComm((PetscObject) dmdest), (DMSNES *) &dmdest->dmsnes);}
277: DMSNESCopy((DMSNES) dmsrc->dmsnes, (DMSNES) dmdest->dmsnes);
278: DMCoarsenHookAdd(dmdest,DMCoarsenHook_DMSNES,NULL,NULL);
279: DMRefineHookAdd(dmdest,DMRefineHook_DMSNES,NULL,NULL);
280: DMSubDomainHookAdd(dmdest,DMSubDomainHook_DMSNES,DMSubDomainRestrictHook_DMSNES,NULL);
281: return(0);
282: }
284: /*@C
285: DMSNESSetFunction - set SNES residual evaluation function
287: Not Collective
289: Input Arguments:
290: + dm - DM to be used with SNES
291: . f - residual evaluation function; see SNESFunction for details
292: - ctx - context for residual evaluation
294: Level: advanced
296: Note:
297: SNESSetFunction() is normally used, but it calls this function internally because the user context is actually
298: associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or
299: not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
301: .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), SNESFunction
302: @*/
303: PetscErrorCode DMSNESSetFunction(DM dm,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
304: {
306: DMSNES sdm;
310: if (f || ctx) {
311: DMGetDMSNESWrite(dm,&sdm);
312: }
313: if (f) sdm->ops->computefunction = f;
314: if (ctx) sdm->functionctx = ctx;
315: return(0);
316: }
318: /*@C
319: DMSNESGetFunction - get SNES residual evaluation function
321: Not Collective
323: Input Argument:
324: . dm - DM to be used with SNES
326: Output Arguments:
327: + f - residual evaluation function; see SNESFunction for details
328: - ctx - context for residual evaluation
330: Level: advanced
332: Note:
333: SNESGetFunction() is normally used, but it calls this function internally because the user context is actually
334: associated with the DM.
336: .seealso: DMSNESSetContext(), DMSNESSetFunction(), SNESSetFunction(), SNESFunction
337: @*/
338: PetscErrorCode DMSNESGetFunction(DM dm,PetscErrorCode (**f)(SNES,Vec,Vec,void*),void **ctx)
339: {
341: DMSNES sdm;
345: DMGetDMSNES(dm,&sdm);
346: if (f) *f = sdm->ops->computefunction;
347: if (ctx) *ctx = sdm->functionctx;
348: return(0);
349: }
351: /*@C
352: DMSNESSetObjective - set SNES objective evaluation function
354: Not Collective
356: Input Arguments:
357: + dm - DM to be used with SNES
358: . obj - objective evaluation function; see SNESObjectiveFunction for details
359: - ctx - context for residual evaluation
361: Level: advanced
363: .seealso: DMSNESSetContext(), SNESGetObjective(), DMSNESSetFunction()
364: @*/
365: PetscErrorCode DMSNESSetObjective(DM dm,PetscErrorCode (*obj)(SNES,Vec,PetscReal*,void*),void *ctx)
366: {
368: DMSNES sdm;
372: if (obj || ctx) {
373: DMGetDMSNESWrite(dm,&sdm);
374: }
375: if (obj) sdm->ops->computeobjective = obj;
376: if (ctx) sdm->objectivectx = ctx;
377: return(0);
378: }
380: /*@C
381: DMSNESGetObjective - get SNES objective evaluation function
383: Not Collective
385: Input Argument:
386: . dm - DM to be used with SNES
388: Output Arguments:
389: + obj- residual evaluation function; see SNESObjectiveFunction for details
390: - ctx - context for residual evaluation
392: Level: advanced
394: Note:
395: SNESGetFunction() is normally used, but it calls this function internally because the user context is actually
396: associated with the DM.
398: .seealso: DMSNESSetContext(), DMSNESSetObjective(), SNESSetFunction()
399: @*/
400: PetscErrorCode DMSNESGetObjective(DM dm,PetscErrorCode (**obj)(SNES,Vec,PetscReal*,void*),void **ctx)
401: {
403: DMSNES sdm;
407: DMGetDMSNES(dm,&sdm);
408: if (obj) *obj = sdm->ops->computeobjective;
409: if (ctx) *ctx = sdm->objectivectx;
410: return(0);
411: }
413: /*@C
414: DMSNESSetNGS - set SNES Gauss-Seidel relaxation function
416: Not Collective
418: Input Argument:
419: + dm - DM to be used with SNES
420: . f - relaxation function, see SNESGSFunction
421: - ctx - context for residual evaluation
423: Level: advanced
425: Note:
426: SNESSetNGS() is normally used, but it calls this function internally because the user context is actually
427: associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or
428: not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
430: .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), DMSNESSetFunction(), SNESGSFunction
431: @*/
432: PetscErrorCode DMSNESSetNGS(DM dm,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
433: {
435: DMSNES sdm;
439: if (f || ctx) {
440: DMGetDMSNESWrite(dm,&sdm);
441: }
442: if (f) sdm->ops->computegs = f;
443: if (ctx) sdm->gsctx = ctx;
444: return(0);
445: }
447: /*@C
448: DMSNESGetNGS - get SNES Gauss-Seidel relaxation function
450: Not Collective
452: Input Argument:
453: . dm - DM to be used with SNES
455: Output Arguments:
456: + f - relaxation function which performs Gauss-Seidel sweeps, see SNESGSFunction
457: - ctx - context for residual evaluation
459: Level: advanced
461: Note:
462: SNESGetNGS() is normally used, but it calls this function internally because the user context is actually
463: associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or
464: not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
466: .seealso: DMSNESSetContext(), SNESGetNGS(), DMSNESGetJacobian(), DMSNESGetFunction(), SNESNGSFunction
467: @*/
468: PetscErrorCode DMSNESGetNGS(DM dm,PetscErrorCode (**f)(SNES,Vec,Vec,void*),void **ctx)
469: {
471: DMSNES sdm;
475: DMGetDMSNES(dm,&sdm);
476: if (f) *f = sdm->ops->computegs;
477: if (ctx) *ctx = sdm->gsctx;
478: return(0);
479: }
481: /*@C
482: DMSNESSetJacobian - set SNES Jacobian evaluation function
484: Not Collective
486: Input Argument:
487: + dm - DM to be used with SNES
488: . J - Jacobian evaluation function
489: - ctx - context for residual evaluation
491: Level: advanced
493: Note:
494: SNESSetJacobian() is normally used, but it calls this function internally because the user context is actually
495: associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or
496: not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
498: .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESGetJacobian(), SNESSetJacobian(), SNESJacobianFunction
499: @*/
500: PetscErrorCode DMSNESSetJacobian(DM dm,PetscErrorCode (*J)(SNES,Vec,Mat,Mat,void*),void *ctx)
501: {
503: DMSNES sdm;
507: if (J || ctx) {
508: DMGetDMSNESWrite(dm,&sdm);
509: }
510: if (J) sdm->ops->computejacobian = J;
511: if (ctx) sdm->jacobianctx = ctx;
512: return(0);
513: }
515: /*@C
516: DMSNESGetJacobian - get SNES Jacobian evaluation function
518: Not Collective
520: Input Argument:
521: . dm - DM to be used with SNES
523: Output Arguments:
524: + J - Jacobian evaluation function; see SNESJacobianFunction for all calling sequence
525: - ctx - context for residual evaluation
527: Level: advanced
529: Note:
530: SNESGetJacobian() is normally used, but it calls this function internally because the user context is actually
531: associated with the DM. 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 Jacobian.
534: .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), SNESJacobianFunction
535: @*/
536: PetscErrorCode DMSNESGetJacobian(DM dm,PetscErrorCode (**J)(SNES,Vec,Mat,Mat,void*),void **ctx)
537: {
539: DMSNES sdm;
543: DMGetDMSNES(dm,&sdm);
544: if (J) *J = sdm->ops->computejacobian;
545: if (ctx) *ctx = sdm->jacobianctx;
546: return(0);
547: }
549: /*@C
550: DMSNESSetPicard - set SNES Picard iteration matrix and RHS evaluation functions.
552: Not Collective
554: Input Argument:
555: + dm - DM to be used with SNES
556: . b - RHS evaluation function
557: . J - Picard matrix evaluation function
558: - ctx - context for residual evaluation
560: Level: advanced
562: .seealso: SNESSetPicard(), DMSNESSetFunction(), DMSNESSetJacobian()
563: @*/
564: PetscErrorCode DMSNESSetPicard(DM dm,PetscErrorCode (*b)(SNES,Vec,Vec,void*),PetscErrorCode (*J)(SNES,Vec,Mat,Mat,void*),void *ctx)
565: {
567: DMSNES sdm;
571: DMGetDMSNES(dm,&sdm);
572: if (b) sdm->ops->computepfunction = b;
573: if (J) sdm->ops->computepjacobian = J;
574: if (ctx) sdm->pctx = ctx;
575: return(0);
576: }
578: /*@C
579: DMSNESGetPicard - get SNES Picard iteration evaluation functions
581: Not Collective
583: Input Argument:
584: . dm - DM to be used with SNES
586: Output Arguments:
587: + b - RHS evaluation function; see SNESFunction for details
588: . J - RHS evaluation function; see SNESJacobianFunction for detailsa
589: - ctx - context for residual evaluation
591: Level: advanced
593: .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian()
594: @*/
595: PetscErrorCode DMSNESGetPicard(DM dm,PetscErrorCode (**b)(SNES,Vec,Vec,void*),PetscErrorCode (**J)(SNES,Vec,Mat,Mat,void*),void **ctx)
596: {
598: DMSNES sdm;
602: DMGetDMSNES(dm,&sdm);
603: if (b) *b = sdm->ops->computepfunction;
604: if (J) *J = sdm->ops->computepjacobian;
605: if (ctx) *ctx = sdm->pctx;
606: return(0);
607: }