Actual source code: dmsnes.c
petsc-3.7.3 2016-08-01
1: #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/
2: #include <petsc/private/dmimpl.h> /*I "petscdm.h" I*/
6: static PetscErrorCode DMSNESDestroy(DMSNES *kdm)
7: {
11: if (!*kdm) return(0);
13: if (--((PetscObject)(*kdm))->refct > 0) {*kdm = 0; return(0);}
14: if ((*kdm)->ops->destroy) {((*kdm)->ops->destroy)(*kdm);}
15: PetscHeaderDestroy(kdm);
16: return(0);
17: }
21: PetscErrorCode DMSNESLoad(DMSNES kdm,PetscViewer viewer)
22: {
26: PetscViewerBinaryRead(viewer,&kdm->ops->computefunction,1,NULL,PETSC_FUNCTION);
27: PetscViewerBinaryRead(viewer,&kdm->ops->computejacobian,1,NULL,PETSC_FUNCTION);
28: return(0);
29: }
33: PetscErrorCode DMSNESView(DMSNES kdm,PetscViewer viewer)
34: {
36: PetscBool isascii,isbinary;
39: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
40: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
41: if (isascii) {
42: #if defined(PETSC_SERIALIZE_FUNCTIONS) && defined(PETSC_SERIALIZE_FUNCTIONS_VIEW)
43: const char *fname;
45: PetscFPTFind(kdm->ops->computefunction,&fname);
46: if (fname) {
47: PetscViewerASCIIPrintf(viewer,"Function used by SNES: %s\n",fname);
48: }
49: PetscFPTFind(kdm->ops->computejacobian,&fname);
50: if (fname) {
51: PetscViewerASCIIPrintf(viewer,"Jacobian function used by SNES: %s\n",fname);
52: }
53: #endif
54: } else if (isbinary) {
55: struct {
56: PetscErrorCode (*func)(SNES,Vec,Vec,void*);
57: } funcstruct;
58: struct {
59: PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*);
60: } jacstruct;
61: funcstruct.func = kdm->ops->computefunction;
62: jacstruct.jac = kdm->ops->computejacobian;
63: PetscViewerBinaryWrite(viewer,&funcstruct,1,PETSC_FUNCTION,PETSC_FALSE);
64: PetscViewerBinaryWrite(viewer,&jacstruct,1,PETSC_FUNCTION,PETSC_FALSE);
65: }
66: return(0);
67: }
71: static PetscErrorCode DMSNESCreate(MPI_Comm comm,DMSNES *kdm)
72: {
76: SNESInitializePackage();
77: PetscHeaderCreate(*kdm, DMSNES_CLASSID, "DMSNES", "DMSNES", "DMSNES", comm, DMSNESDestroy, DMSNESView);
78: return(0);
79: }
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: {
91: DMCopyDMSNES(dm,dmc);
92: return(0);
93: }
97: /* This could restrict auxiliary information to the coarse level.
98: */
99: static PetscErrorCode DMRestrictHook_DMSNES(DM dm,Mat Restrict,Vec rscale,Mat Inject,DM dmc,void *ctx)
100: {
103: return(0);
104: }
108: /* Attaches the DMSNES to the subdomain. */
109: static PetscErrorCode DMSubDomainHook_DMSNES(DM dm,DM subdm,void *ctx)
110: {
114: DMCopyDMSNES(dm,subdm);
115: return(0);
116: }
120: /* This could restrict auxiliary information to the coarse level.
121: */
122: static PetscErrorCode DMSubDomainRestrictHook_DMSNES(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
123: {
126: return(0);
127: }
131: static PetscErrorCode DMRefineHook_DMSNES(DM dm,DM dmf,void *ctx)
132: {
136: DMCopyDMSNES(dm,dmf);
137: return(0);
138: }
142: /* This could restrict auxiliary information to the coarse level.
143: */
144: static PetscErrorCode DMInterpolateHook_DMSNES(DM dm,Mat Interp,DM dmf,void *ctx)
145: {
148: return(0);
149: }
153: /*@C
154: DMSNESCopy - copies the information in a DMSNES to another DMSNES
156: Not Collective
158: Input Argument:
159: + kdm - Original DMSNES
160: - nkdm - DMSNES to receive the data, should have been created with DMSNESCreate()
162: Level: developer
164: .seealso: DMSNESCreate(), DMSNESDestroy()
165: @*/
166: PetscErrorCode DMSNESCopy(DMSNES kdm,DMSNES nkdm)
167: {
173: nkdm->ops->computefunction = kdm->ops->computefunction;
174: nkdm->ops->computejacobian = kdm->ops->computejacobian;
175: nkdm->ops->computegs = kdm->ops->computegs;
176: nkdm->ops->computeobjective = kdm->ops->computeobjective;
177: nkdm->ops->computepjacobian = kdm->ops->computepjacobian;
178: nkdm->ops->computepfunction = kdm->ops->computepfunction;
179: nkdm->ops->destroy = kdm->ops->destroy;
180: nkdm->ops->duplicate = kdm->ops->duplicate;
182: nkdm->functionctx = kdm->functionctx;
183: nkdm->gsctx = kdm->gsctx;
184: nkdm->pctx = kdm->pctx;
185: nkdm->jacobianctx = kdm->jacobianctx;
186: nkdm->objectivectx = kdm->objectivectx;
187: nkdm->data = kdm->data;
189: /*
190: nkdm->fortran_func_pointers[0] = kdm->fortran_func_pointers[0];
191: nkdm->fortran_func_pointers[1] = kdm->fortran_func_pointers[1];
192: nkdm->fortran_func_pointers[2] = kdm->fortran_func_pointers[2];
193: */
195: /* implementation specific copy hooks */
196: if (kdm->ops->duplicate) {(*kdm->ops->duplicate)(kdm,nkdm);}
197: return(0);
198: }
202: /*@C
203: DMGetDMSNES - get read-only private DMSNES context from a DM
205: Not Collective
207: Input Argument:
208: . dm - DM to be used with SNES
210: Output Argument:
211: . snesdm - private DMSNES context
213: Level: developer
215: Notes:
216: Use DMGetDMSNESWrite() if write access is needed. The DMSNESSetXXX API should be used wherever possible.
218: .seealso: DMGetDMSNESWrite()
219: @*/
220: PetscErrorCode DMGetDMSNES(DM dm,DMSNES *snesdm)
221: {
226: *snesdm = (DMSNES) dm->dmsnes;
227: if (!*snesdm) {
228: PetscInfo(dm,"Creating new DMSNES\n");
229: DMSNESCreate(PetscObjectComm((PetscObject)dm),snesdm);
231: dm->dmsnes = (PetscObject) *snesdm;
233: DMCoarsenHookAdd(dm,DMCoarsenHook_DMSNES,DMRestrictHook_DMSNES,NULL);
234: DMRefineHookAdd(dm,DMRefineHook_DMSNES,DMInterpolateHook_DMSNES,NULL);
235: DMSubDomainHookAdd(dm,DMSubDomainHook_DMSNES,DMSubDomainRestrictHook_DMSNES,NULL);
236: }
237: return(0);
238: }
242: /*@C
243: DMGetDMSNESWrite - get write access to private DMSNES context from a DM
245: Not Collective
247: Input Argument:
248: . dm - DM to be used with SNES
250: Output Argument:
251: . snesdm - private DMSNES context
253: Level: developer
255: .seealso: DMGetDMSNES()
256: @*/
257: PetscErrorCode DMGetDMSNESWrite(DM dm,DMSNES *snesdm)
258: {
260: DMSNES sdm;
264: DMGetDMSNES(dm,&sdm);
265: if (!sdm->originaldm) sdm->originaldm = dm;
266: if (sdm->originaldm != dm) { /* Copy on write */
267: DMSNES oldsdm = sdm;
268: PetscInfo(dm,"Copying DMSNES due to write\n");
269: DMSNESCreate(PetscObjectComm((PetscObject)dm),&sdm);
270: DMSNESCopy(oldsdm,sdm);
271: DMSNESDestroy((DMSNES*)&dm->dmsnes);
272: dm->dmsnes = (PetscObject)sdm;
273: }
274: *snesdm = sdm;
275: return(0);
276: }
280: /*@C
281: DMCopyDMSNES - copies a DM context to a new DM
283: Logically Collective
285: Input Arguments:
286: + dmsrc - DM to obtain context from
287: - dmdest - DM to add context to
289: Level: developer
291: Note:
292: The context is copied by reference. This function does not ensure that a context exists.
294: .seealso: DMGetDMSNES(), SNESSetDM()
295: @*/
296: PetscErrorCode DMCopyDMSNES(DM dmsrc,DM dmdest)
297: {
303: if (!dmdest->dmsnes) {DMSNESCreate(PetscObjectComm((PetscObject) dmdest), (DMSNES *) &dmdest->dmsnes);}
304: DMSNESCopy((DMSNES) dmsrc->dmsnes, (DMSNES) dmdest->dmsnes);
305: DMCoarsenHookAdd(dmdest,DMCoarsenHook_DMSNES,NULL,NULL);
306: DMRefineHookAdd(dmdest,DMRefineHook_DMSNES,NULL,NULL);
307: DMSubDomainHookAdd(dmdest,DMSubDomainHook_DMSNES,DMSubDomainRestrictHook_DMSNES,NULL);
308: return(0);
309: }
313: /*@C
314: DMSNESSetFunction - set SNES residual evaluation function
316: Not Collective
318: Input Arguments:
319: + dm - DM to be used with SNES
320: . f - residual evaluation function; see SNESFunction for details
321: - ctx - context for residual evaluation
323: Level: advanced
325: Note:
326: SNESSetFunction() is normally used, but it calls this function internally because the user context is actually
327: associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or
328: not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
330: .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), SNESFunction
331: @*/
332: PetscErrorCode DMSNESSetFunction(DM dm,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
333: {
335: DMSNES sdm;
339: if (f || ctx) {
340: DMGetDMSNESWrite(dm,&sdm);
341: }
342: if (f) sdm->ops->computefunction = f;
343: if (ctx) sdm->functionctx = ctx;
344: return(0);
345: }
349: /*@C
350: DMSNESGetFunction - get SNES residual evaluation function
352: Not Collective
354: Input Argument:
355: . dm - DM to be used with SNES
357: Output Arguments:
358: + f - residual evaluation function; see SNESFunction for details
359: - ctx - context for residual evaluation
361: Level: advanced
363: Note:
364: SNESGetFunction() is normally used, but it calls this function internally because the user context is actually
365: associated with the DM.
367: .seealso: DMSNESSetContext(), DMSNESSetFunction(), SNESSetFunction(), SNESFunction
368: @*/
369: PetscErrorCode DMSNESGetFunction(DM dm,PetscErrorCode (**f)(SNES,Vec,Vec,void*),void **ctx)
370: {
372: DMSNES sdm;
376: DMGetDMSNES(dm,&sdm);
377: if (f) *f = sdm->ops->computefunction;
378: if (ctx) *ctx = sdm->functionctx;
379: return(0);
380: }
384: /*@C
385: DMSNESSetObjective - set SNES objective evaluation function
387: Not Collective
389: Input Arguments:
390: + dm - DM to be used with SNES
391: . obj - objective evaluation function; see SNESObjectiveFunction for details
392: - ctx - context for residual evaluation
394: Level: advanced
396: .seealso: DMSNESSetContext(), SNESGetObjective(), DMSNESSetFunction()
397: @*/
398: PetscErrorCode DMSNESSetObjective(DM dm,PetscErrorCode (*obj)(SNES,Vec,PetscReal*,void*),void *ctx)
399: {
401: DMSNES sdm;
405: if (obj || ctx) {
406: DMGetDMSNESWrite(dm,&sdm);
407: }
408: if (obj) sdm->ops->computeobjective = obj;
409: if (ctx) sdm->objectivectx = ctx;
410: return(0);
411: }
415: /*@C
416: DMSNESGetObjective - get SNES objective evaluation function
418: Not Collective
420: Input Argument:
421: . dm - DM to be used with SNES
423: Output Arguments:
424: + obj- residual evaluation function; see SNESObjectiveFunction for details
425: - ctx - context for residual evaluation
427: Level: advanced
429: Note:
430: SNESGetFunction() is normally used, but it calls this function internally because the user context is actually
431: associated with the DM.
433: .seealso: DMSNESSetContext(), DMSNESSetObjective(), SNESSetFunction()
434: @*/
435: PetscErrorCode DMSNESGetObjective(DM dm,PetscErrorCode (**obj)(SNES,Vec,PetscReal*,void*),void **ctx)
436: {
438: DMSNES sdm;
442: DMGetDMSNES(dm,&sdm);
443: if (obj) *obj = sdm->ops->computeobjective;
444: if (ctx) *ctx = sdm->objectivectx;
445: return(0);
446: }
450: /*@C
451: DMSNESSetNGS - set SNES Gauss-Seidel relaxation function
453: Not Collective
455: Input Argument:
456: + dm - DM to be used with SNES
457: . f - relaxation function, see SNESGSFunction
458: - ctx - context for residual evaluation
460: Level: advanced
462: Note:
463: SNESSetNGS() is normally used, but it calls this function internally because the user context is actually
464: associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or
465: not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
467: .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), DMSNESSetFunction(), SNESGSFunction
468: @*/
469: PetscErrorCode DMSNESSetNGS(DM dm,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
470: {
472: DMSNES sdm;
476: if (f || ctx) {
477: DMGetDMSNESWrite(dm,&sdm);
478: }
479: if (f) sdm->ops->computegs = f;
480: if (ctx) sdm->gsctx = ctx;
481: return(0);
482: }
486: /*@C
487: DMSNESGetNGS - get SNES Gauss-Seidel relaxation function
489: Not Collective
491: Input Argument:
492: . dm - DM to be used with SNES
494: Output Arguments:
495: + f - relaxation function which performs Gauss-Seidel sweeps, see SNESGSFunction
496: - ctx - context for residual evaluation
498: Level: advanced
500: Note:
501: SNESGetNGS() is normally used, but it calls this function internally because the user context is actually
502: associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or
503: not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
505: .seealso: DMSNESSetContext(), SNESGetNGS(), DMSNESGetJacobian(), DMSNESGetFunction(), SNESNGSFunction
506: @*/
507: PetscErrorCode DMSNESGetNGS(DM dm,PetscErrorCode (**f)(SNES,Vec,Vec,void*),void **ctx)
508: {
510: DMSNES sdm;
514: DMGetDMSNES(dm,&sdm);
515: if (f) *f = sdm->ops->computegs;
516: if (ctx) *ctx = sdm->gsctx;
517: return(0);
518: }
522: /*@C
523: DMSNESSetJacobian - set SNES Jacobian evaluation function
525: Not Collective
527: Input Argument:
528: + dm - DM to be used with SNES
529: . J - Jacobian evaluation function
530: - ctx - context for residual evaluation
532: Level: advanced
534: Note:
535: SNESSetJacobian() is normally used, but it calls this function internally because the user context is actually
536: associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or
537: not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
539: .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESGetJacobian(), SNESSetJacobian(), SNESJacobianFunction
540: @*/
541: PetscErrorCode DMSNESSetJacobian(DM dm,PetscErrorCode (*J)(SNES,Vec,Mat,Mat,void*),void *ctx)
542: {
544: DMSNES sdm;
548: if (J || ctx) {
549: DMGetDMSNESWrite(dm,&sdm);
550: }
551: if (J) sdm->ops->computejacobian = J;
552: if (ctx) sdm->jacobianctx = ctx;
553: return(0);
554: }
558: /*@C
559: DMSNESGetJacobian - get SNES Jacobian evaluation function
561: Not Collective
563: Input Argument:
564: . dm - DM to be used with SNES
566: Output Arguments:
567: + J - Jacobian evaluation function; see SNESJacobianFunction for all calling sequence
568: - ctx - context for residual evaluation
570: Level: advanced
572: Note:
573: SNESGetJacobian() is normally used, but it calls this function internally because the user context is actually
574: associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or
575: not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
577: .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), SNESJacobianFunction
578: @*/
579: PetscErrorCode DMSNESGetJacobian(DM dm,PetscErrorCode (**J)(SNES,Vec,Mat,Mat,void*),void **ctx)
580: {
582: DMSNES sdm;
586: DMGetDMSNES(dm,&sdm);
587: if (J) *J = sdm->ops->computejacobian;
588: if (ctx) *ctx = sdm->jacobianctx;
589: return(0);
590: }
594: /*@C
595: DMSNESSetPicard - set SNES Picard iteration matrix and RHS evaluation functions.
597: Not Collective
599: Input Argument:
600: + dm - DM to be used with SNES
601: . b - RHS evaluation function
602: . J - Picard matrix evaluation function
603: - ctx - context for residual evaluation
605: Level: advanced
607: .seealso: SNESSetPicard(), DMSNESSetFunction(), DMSNESSetJacobian()
608: @*/
609: PetscErrorCode DMSNESSetPicard(DM dm,PetscErrorCode (*b)(SNES,Vec,Vec,void*),PetscErrorCode (*J)(SNES,Vec,Mat,Mat,void*),void *ctx)
610: {
612: DMSNES sdm;
616: DMGetDMSNES(dm,&sdm);
617: if (b) sdm->ops->computepfunction = b;
618: if (J) sdm->ops->computepjacobian = J;
619: if (ctx) sdm->pctx = ctx;
620: return(0);
621: }
625: /*@C
626: DMSNESGetPicard - get SNES Picard iteration evaluation functions
628: Not Collective
630: Input Argument:
631: . dm - DM to be used with SNES
633: Output Arguments:
634: + b - RHS evaluation function; see SNESFunction for details
635: . J - RHS evaluation function; see SNESJacobianFunction for detailsa
636: - ctx - context for residual evaluation
638: Level: advanced
640: .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian()
641: @*/
642: PetscErrorCode DMSNESGetPicard(DM dm,PetscErrorCode (**b)(SNES,Vec,Vec,void*),PetscErrorCode (**J)(SNES,Vec,Mat,Mat,void*),void **ctx)
643: {
645: DMSNES sdm;
649: DMGetDMSNES(dm,&sdm);
650: if (b) *b = sdm->ops->computepfunction;
651: if (J) *J = sdm->ops->computepjacobian;
652: if (ctx) *ctx = sdm->pctx;
653: return(0);
654: }