Actual source code: dmsnes.c
petsc-3.9.4 2018-09-11
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->data = kdm->data;
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;
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) sdm->originaldm = dm;
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: }
248: *snesdm = sdm;
249: return(0);
250: }
252: /*@C
253: DMCopyDMSNES - copies a DM context to a new DM
255: Logically Collective
257: Input Arguments:
258: + dmsrc - DM to obtain context from
259: - dmdest - DM to add context to
261: Level: developer
263: Note:
264: The context is copied by reference. This function does not ensure that a context exists.
266: .seealso: DMGetDMSNES(), SNESSetDM()
267: @*/
268: PetscErrorCode DMCopyDMSNES(DM dmsrc,DM dmdest)
269: {
275: if (!dmdest->dmsnes) {DMSNESCreate(PetscObjectComm((PetscObject) dmdest), (DMSNES *) &dmdest->dmsnes);}
276: DMSNESCopy((DMSNES) dmsrc->dmsnes, (DMSNES) dmdest->dmsnes);
277: DMCoarsenHookAdd(dmdest,DMCoarsenHook_DMSNES,NULL,NULL);
278: DMRefineHookAdd(dmdest,DMRefineHook_DMSNES,NULL,NULL);
279: DMSubDomainHookAdd(dmdest,DMSubDomainHook_DMSNES,DMSubDomainRestrictHook_DMSNES,NULL);
280: return(0);
281: }
283: /*@C
284: DMSNESSetFunction - set SNES residual evaluation function
286: Not Collective
288: Input Arguments:
289: + dm - DM to be used with SNES
290: . f - residual evaluation function; see SNESFunction for details
291: - ctx - context for residual evaluation
293: Level: advanced
295: Note:
296: SNESSetFunction() is normally used, but it calls this function internally because the user context is actually
297: associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or
298: not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
300: .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), SNESFunction
301: @*/
302: PetscErrorCode DMSNESSetFunction(DM dm,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
303: {
305: DMSNES sdm;
309: if (f || ctx) {
310: DMGetDMSNESWrite(dm,&sdm);
311: }
312: if (f) sdm->ops->computefunction = f;
313: if (ctx) sdm->functionctx = ctx;
314: return(0);
315: }
317: /*@C
318: DMSNESGetFunction - get SNES residual evaluation function
320: Not Collective
322: Input Argument:
323: . dm - DM to be used with SNES
325: Output Arguments:
326: + f - residual evaluation function; see SNESFunction for details
327: - ctx - context for residual evaluation
329: Level: advanced
331: Note:
332: SNESGetFunction() is normally used, but it calls this function internally because the user context is actually
333: associated with the DM.
335: .seealso: DMSNESSetContext(), DMSNESSetFunction(), SNESSetFunction(), SNESFunction
336: @*/
337: PetscErrorCode DMSNESGetFunction(DM dm,PetscErrorCode (**f)(SNES,Vec,Vec,void*),void **ctx)
338: {
340: DMSNES sdm;
344: DMGetDMSNES(dm,&sdm);
345: if (f) *f = sdm->ops->computefunction;
346: if (ctx) *ctx = sdm->functionctx;
347: return(0);
348: }
350: /*@C
351: DMSNESSetObjective - set SNES objective evaluation function
353: Not Collective
355: Input Arguments:
356: + dm - DM to be used with SNES
357: . obj - objective evaluation function; see SNESObjectiveFunction for details
358: - ctx - context for residual evaluation
360: Level: advanced
362: .seealso: DMSNESSetContext(), SNESGetObjective(), DMSNESSetFunction()
363: @*/
364: PetscErrorCode DMSNESSetObjective(DM dm,PetscErrorCode (*obj)(SNES,Vec,PetscReal*,void*),void *ctx)
365: {
367: DMSNES sdm;
371: if (obj || ctx) {
372: DMGetDMSNESWrite(dm,&sdm);
373: }
374: if (obj) sdm->ops->computeobjective = obj;
375: if (ctx) sdm->objectivectx = ctx;
376: return(0);
377: }
379: /*@C
380: DMSNESGetObjective - get SNES objective evaluation function
382: Not Collective
384: Input Argument:
385: . dm - DM to be used with SNES
387: Output Arguments:
388: + obj- residual evaluation function; see SNESObjectiveFunction for details
389: - ctx - context for residual evaluation
391: Level: advanced
393: Note:
394: SNESGetFunction() is normally used, but it calls this function internally because the user context is actually
395: associated with the DM.
397: .seealso: DMSNESSetContext(), DMSNESSetObjective(), SNESSetFunction()
398: @*/
399: PetscErrorCode DMSNESGetObjective(DM dm,PetscErrorCode (**obj)(SNES,Vec,PetscReal*,void*),void **ctx)
400: {
402: DMSNES sdm;
406: DMGetDMSNES(dm,&sdm);
407: if (obj) *obj = sdm->ops->computeobjective;
408: if (ctx) *ctx = sdm->objectivectx;
409: return(0);
410: }
412: /*@C
413: DMSNESSetNGS - set SNES Gauss-Seidel relaxation function
415: Not Collective
417: Input Argument:
418: + dm - DM to be used with SNES
419: . f - relaxation function, see SNESGSFunction
420: - ctx - context for residual evaluation
422: Level: advanced
424: Note:
425: SNESSetNGS() is normally used, but it calls this function internally because the user context is actually
426: associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or
427: not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
429: .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), DMSNESSetFunction(), SNESGSFunction
430: @*/
431: PetscErrorCode DMSNESSetNGS(DM dm,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
432: {
434: DMSNES sdm;
438: if (f || ctx) {
439: DMGetDMSNESWrite(dm,&sdm);
440: }
441: if (f) sdm->ops->computegs = f;
442: if (ctx) sdm->gsctx = ctx;
443: return(0);
444: }
446: /*@C
447: DMSNESGetNGS - get SNES Gauss-Seidel relaxation function
449: Not Collective
451: Input Argument:
452: . dm - DM to be used with SNES
454: Output Arguments:
455: + f - relaxation function which performs Gauss-Seidel sweeps, see SNESGSFunction
456: - ctx - context for residual evaluation
458: Level: advanced
460: Note:
461: SNESGetNGS() is normally used, but it calls this function internally because the user context is actually
462: associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or
463: not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
465: .seealso: DMSNESSetContext(), SNESGetNGS(), DMSNESGetJacobian(), DMSNESGetFunction(), SNESNGSFunction
466: @*/
467: PetscErrorCode DMSNESGetNGS(DM dm,PetscErrorCode (**f)(SNES,Vec,Vec,void*),void **ctx)
468: {
470: DMSNES sdm;
474: DMGetDMSNES(dm,&sdm);
475: if (f) *f = sdm->ops->computegs;
476: if (ctx) *ctx = sdm->gsctx;
477: return(0);
478: }
480: /*@C
481: DMSNESSetJacobian - set SNES Jacobian evaluation function
483: Not Collective
485: Input Argument:
486: + dm - DM to be used with SNES
487: . J - Jacobian evaluation function
488: - ctx - context for residual evaluation
490: Level: advanced
492: Note:
493: SNESSetJacobian() is normally used, but it calls this function internally because the user context is actually
494: associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or
495: not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
497: .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESGetJacobian(), SNESSetJacobian(), SNESJacobianFunction
498: @*/
499: PetscErrorCode DMSNESSetJacobian(DM dm,PetscErrorCode (*J)(SNES,Vec,Mat,Mat,void*),void *ctx)
500: {
502: DMSNES sdm;
506: if (J || ctx) {
507: DMGetDMSNESWrite(dm,&sdm);
508: }
509: if (J) sdm->ops->computejacobian = J;
510: if (ctx) sdm->jacobianctx = ctx;
511: return(0);
512: }
514: /*@C
515: DMSNESGetJacobian - get SNES Jacobian evaluation function
517: Not Collective
519: Input Argument:
520: . dm - DM to be used with SNES
522: Output Arguments:
523: + J - Jacobian evaluation function; see SNESJacobianFunction for all calling sequence
524: - ctx - context for residual evaluation
526: Level: advanced
528: Note:
529: SNESGetJacobian() is normally used, but it calls this function internally because the user context is actually
530: associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or
531: not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
533: .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), SNESJacobianFunction
534: @*/
535: PetscErrorCode DMSNESGetJacobian(DM dm,PetscErrorCode (**J)(SNES,Vec,Mat,Mat,void*),void **ctx)
536: {
538: DMSNES sdm;
542: DMGetDMSNES(dm,&sdm);
543: if (J) *J = sdm->ops->computejacobian;
544: if (ctx) *ctx = sdm->jacobianctx;
545: return(0);
546: }
548: /*@C
549: DMSNESSetPicard - set SNES Picard iteration matrix and RHS evaluation functions.
551: Not Collective
553: Input Argument:
554: + dm - DM to be used with SNES
555: . b - RHS evaluation function
556: . J - Picard matrix evaluation function
557: - ctx - context for residual evaluation
559: Level: advanced
561: .seealso: SNESSetPicard(), DMSNESSetFunction(), DMSNESSetJacobian()
562: @*/
563: PetscErrorCode DMSNESSetPicard(DM dm,PetscErrorCode (*b)(SNES,Vec,Vec,void*),PetscErrorCode (*J)(SNES,Vec,Mat,Mat,void*),void *ctx)
564: {
566: DMSNES sdm;
570: DMGetDMSNES(dm,&sdm);
571: if (b) sdm->ops->computepfunction = b;
572: if (J) sdm->ops->computepjacobian = J;
573: if (ctx) sdm->pctx = ctx;
574: return(0);
575: }
577: /*@C
578: DMSNESGetPicard - get SNES Picard iteration evaluation functions
580: Not Collective
582: Input Argument:
583: . dm - DM to be used with SNES
585: Output Arguments:
586: + b - RHS evaluation function; see SNESFunction for details
587: . J - RHS evaluation function; see SNESJacobianFunction for detailsa
588: - ctx - context for residual evaluation
590: Level: advanced
592: .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian()
593: @*/
594: PetscErrorCode DMSNESGetPicard(DM dm,PetscErrorCode (**b)(SNES,Vec,Vec,void*),PetscErrorCode (**J)(SNES,Vec,Mat,Mat,void*),void **ctx)
595: {
597: DMSNES sdm;
601: DMGetDMSNES(dm,&sdm);
602: if (b) *b = sdm->ops->computepfunction;
603: if (J) *J = sdm->ops->computepjacobian;
604: if (ctx) *ctx = sdm->pctx;
605: return(0);
606: }