xref: /petsc/src/snes/utils/dmsnes.c (revision dfe153156869f5931ee6b299eecdcfd8470c0dd0)
1b45d2f2cSJed Brown #include <petsc-private/snesimpl.h>   /*I "petscsnes.h" I*/
26cab3a1bSJed Brown #include <petscdm.h>            /*I "petscdm.h" I*/
36cab3a1bSJed Brown 
46cab3a1bSJed Brown #undef __FUNCT__
56cab3a1bSJed Brown #define __FUNCT__ "PetscContainerDestroy_SNESDM"
66cab3a1bSJed Brown static PetscErrorCode PetscContainerDestroy_SNESDM(void *ctx)
76cab3a1bSJed Brown {
86cab3a1bSJed Brown   PetscErrorCode ierr;
96cab3a1bSJed Brown   SNESDM sdm = (SNESDM)ctx;
106cab3a1bSJed Brown 
116cab3a1bSJed Brown   PetscFunctionBegin;
126cab3a1bSJed Brown   if (sdm->destroy) {ierr = (*sdm->destroy)(sdm);CHKERRQ(ierr);}
136cab3a1bSJed Brown   ierr = PetscFree(sdm);CHKERRQ(ierr);
146cab3a1bSJed Brown   PetscFunctionReturn(0);
156cab3a1bSJed Brown }
166cab3a1bSJed Brown 
176cab3a1bSJed Brown #undef __FUNCT__
186cab3a1bSJed Brown #define __FUNCT__ "DMCoarsenHook_SNESDM"
196cab3a1bSJed Brown /* Attaches the SNESDM to the coarse level.
206cab3a1bSJed Brown  * Under what conditions should we copy versus duplicate?
216cab3a1bSJed Brown  */
226cab3a1bSJed Brown static PetscErrorCode DMCoarsenHook_SNESDM(DM dm,DM dmc,void *ctx)
236cab3a1bSJed Brown {
246cab3a1bSJed Brown   PetscErrorCode ierr;
256cab3a1bSJed Brown 
266cab3a1bSJed Brown   PetscFunctionBegin;
276cab3a1bSJed Brown   ierr = DMSNESCopyContext(dm,dmc);CHKERRQ(ierr);
286cab3a1bSJed Brown   PetscFunctionReturn(0);
296cab3a1bSJed Brown }
306cab3a1bSJed Brown 
316cab3a1bSJed Brown #undef __FUNCT__
32caa4e7f2SJed Brown #define __FUNCT__ "DMRestrictHook_SNESDM"
33*dfe15315SJed Brown /* This could restrict auxiliary information to the coarse level.
34caa4e7f2SJed Brown  */
35caa4e7f2SJed Brown static PetscErrorCode DMRestrictHook_SNESDM(DM dm,Mat Restrict,Vec rscale,Mat Inject,DM dmc,void *ctx)
36caa4e7f2SJed Brown {
37caa4e7f2SJed Brown 
38caa4e7f2SJed Brown   PetscFunctionBegin;
39caa4e7f2SJed Brown   PetscFunctionReturn(0);
40caa4e7f2SJed Brown }
41caa4e7f2SJed Brown 
42caa4e7f2SJed Brown #undef __FUNCT__
436cab3a1bSJed Brown #define __FUNCT__ "DMSNESGetContext"
446cab3a1bSJed Brown /*@C
456cab3a1bSJed Brown    DMSNESGetContext - get read-only private SNESDM context from a DM
466cab3a1bSJed Brown 
476cab3a1bSJed Brown    Not Collective
486cab3a1bSJed Brown 
496cab3a1bSJed Brown    Input Argument:
506cab3a1bSJed Brown .  dm - DM to be used with SNES
516cab3a1bSJed Brown 
526cab3a1bSJed Brown    Output Argument:
536cab3a1bSJed Brown .  snesdm - private SNESDM context
546cab3a1bSJed Brown 
556cab3a1bSJed Brown    Level: developer
566cab3a1bSJed Brown 
576cab3a1bSJed Brown    Notes:
586cab3a1bSJed Brown    Use DMSNESGetContextWrite() if write access is needed. The DMSNESSetXXX API should be used wherever possible.
596cab3a1bSJed Brown 
606cab3a1bSJed Brown .seealso: DMSNESGetContextWrite()
616cab3a1bSJed Brown @*/
626cab3a1bSJed Brown PetscErrorCode DMSNESGetContext(DM dm,SNESDM *snesdm)
636cab3a1bSJed Brown {
646cab3a1bSJed Brown   PetscErrorCode ierr;
656cab3a1bSJed Brown   PetscContainer container;
666cab3a1bSJed Brown   SNESDM         sdm;
676cab3a1bSJed Brown 
686cab3a1bSJed Brown   PetscFunctionBegin;
696cab3a1bSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
706cab3a1bSJed Brown   ierr = PetscObjectQuery((PetscObject)dm,"SNESDM",(PetscObject*)&container);CHKERRQ(ierr);
716cab3a1bSJed Brown   if (container) {
726cab3a1bSJed Brown     ierr = PetscContainerGetPointer(container,(void**)snesdm);CHKERRQ(ierr);
736cab3a1bSJed Brown   } else {
746cab3a1bSJed Brown     ierr = PetscInfo(dm,"Creating new SNESDM\n");CHKERRQ(ierr);
756cab3a1bSJed Brown     ierr = PetscContainerCreate(((PetscObject)dm)->comm,&container);CHKERRQ(ierr);
766cab3a1bSJed Brown     ierr = PetscNewLog(dm,struct _n_SNESDM,&sdm);CHKERRQ(ierr);
776cab3a1bSJed Brown     ierr = PetscContainerSetPointer(container,sdm);CHKERRQ(ierr);
786cab3a1bSJed Brown     ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_SNESDM);CHKERRQ(ierr);
796cab3a1bSJed Brown     ierr = PetscObjectCompose((PetscObject)dm,"SNESDM",(PetscObject)container);CHKERRQ(ierr);
80caa4e7f2SJed Brown     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_SNESDM,DMRestrictHook_SNESDM,PETSC_NULL);CHKERRQ(ierr);
816cab3a1bSJed Brown     ierr = PetscContainerGetPointer(container,(void**)snesdm);CHKERRQ(ierr);
826cab3a1bSJed Brown     ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
836cab3a1bSJed Brown   }
846cab3a1bSJed Brown   PetscFunctionReturn(0);
856cab3a1bSJed Brown }
866cab3a1bSJed Brown 
876cab3a1bSJed Brown #undef __FUNCT__
886cab3a1bSJed Brown #define __FUNCT__ "DMSNESGetContextWrite"
896cab3a1bSJed Brown /*@C
906cab3a1bSJed Brown    DMSNESGetContextWrite - get write access to private SNESDM context from a DM
916cab3a1bSJed Brown 
926cab3a1bSJed Brown    Not Collective
936cab3a1bSJed Brown 
946cab3a1bSJed Brown    Input Argument:
956cab3a1bSJed Brown .  dm - DM to be used with SNES
966cab3a1bSJed Brown 
976cab3a1bSJed Brown    Output Argument:
986cab3a1bSJed Brown .  snesdm - private SNESDM context
996cab3a1bSJed Brown 
1006cab3a1bSJed Brown    Level: developer
1016cab3a1bSJed Brown 
1026cab3a1bSJed Brown .seealso: DMSNESGetContext()
1036cab3a1bSJed Brown @*/
1046cab3a1bSJed Brown PetscErrorCode DMSNESGetContextWrite(DM dm,SNESDM *snesdm)
1056cab3a1bSJed Brown {
1066cab3a1bSJed Brown   PetscErrorCode ierr;
1076cab3a1bSJed Brown   SNESDM         sdm;
1086cab3a1bSJed Brown 
1096cab3a1bSJed Brown   PetscFunctionBegin;
1106cab3a1bSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1116cab3a1bSJed Brown   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
1126cab3a1bSJed Brown   if (!sdm->originaldm) sdm->originaldm = dm;
1136cab3a1bSJed Brown   if (sdm->originaldm != dm) {  /* Copy on write */
1146cab3a1bSJed Brown     PetscContainer container;
1156cab3a1bSJed Brown     SNESDM         oldsdm = sdm;
1166cab3a1bSJed Brown     ierr = PetscInfo(dm,"Copying SNESDM due to write\n");CHKERRQ(ierr);
1176cab3a1bSJed Brown     ierr = PetscContainerCreate(((PetscObject)dm)->comm,&container);CHKERRQ(ierr);
1186cab3a1bSJed Brown     ierr = PetscNewLog(dm,struct _n_SNESDM,&sdm);CHKERRQ(ierr);
1196cab3a1bSJed Brown     ierr = PetscMemcpy(sdm,oldsdm,sizeof *sdm);CHKERRQ(ierr);
1206cab3a1bSJed Brown     ierr = PetscContainerSetPointer(container,sdm);CHKERRQ(ierr);
1216cab3a1bSJed Brown     ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_SNESDM);CHKERRQ(ierr);
1226cab3a1bSJed Brown     ierr = PetscObjectCompose((PetscObject)dm,"SNESDM",(PetscObject)container);CHKERRQ(ierr);
1236cab3a1bSJed Brown     ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
1246cab3a1bSJed Brown   }
1256cab3a1bSJed Brown   *snesdm = sdm;
1266cab3a1bSJed Brown   PetscFunctionReturn(0);
1276cab3a1bSJed Brown }
1286cab3a1bSJed Brown 
1296cab3a1bSJed Brown #undef __FUNCT__
1306cab3a1bSJed Brown #define __FUNCT__ "DMSNESCopyContext"
1316cab3a1bSJed Brown /*@C
1326cab3a1bSJed Brown    DMSNESCopyContext - copies a DM context to a new DM
1336cab3a1bSJed Brown 
1346cab3a1bSJed Brown    Logically Collective
1356cab3a1bSJed Brown 
1366cab3a1bSJed Brown    Input Arguments:
1376cab3a1bSJed Brown +  dmsrc - DM to obtain context from
1386cab3a1bSJed Brown -  dmdest - DM to add context to
1396cab3a1bSJed Brown 
1406cab3a1bSJed Brown    Level: developer
1416cab3a1bSJed Brown 
1426cab3a1bSJed Brown    Note:
1436cab3a1bSJed Brown    The context is copied by reference. This function does not ensure that a context exists.
1446cab3a1bSJed Brown 
1456cab3a1bSJed Brown .seealso: DMSNESGetContext(), SNESSetDM()
1466cab3a1bSJed Brown @*/
1476cab3a1bSJed Brown PetscErrorCode DMSNESCopyContext(DM dmsrc,DM dmdest)
1486cab3a1bSJed Brown {
1496cab3a1bSJed Brown   PetscErrorCode ierr;
1506cab3a1bSJed Brown   PetscContainer container;
1516cab3a1bSJed Brown 
1526cab3a1bSJed Brown   PetscFunctionBegin;
1536cab3a1bSJed Brown   PetscValidHeaderSpecific(dmsrc,DM_CLASSID,1);
1546cab3a1bSJed Brown   PetscValidHeaderSpecific(dmdest,DM_CLASSID,2);
1556cab3a1bSJed Brown   ierr = PetscObjectQuery((PetscObject)dmsrc,"SNESDM",(PetscObject*)&container);CHKERRQ(ierr);
1566cab3a1bSJed Brown   if (container) {
1576cab3a1bSJed Brown     ierr = PetscObjectCompose((PetscObject)dmdest,"SNESDM",(PetscObject)container);CHKERRQ(ierr);
158caa4e7f2SJed Brown     ierr = DMCoarsenHookAdd(dmdest,DMCoarsenHook_SNESDM,DMRestrictHook_SNESDM,PETSC_NULL);CHKERRQ(ierr);
1596cab3a1bSJed Brown   }
1606cab3a1bSJed Brown   PetscFunctionReturn(0);
1616cab3a1bSJed Brown }
1626cab3a1bSJed Brown 
1636cab3a1bSJed Brown #undef __FUNCT__
1646cab3a1bSJed Brown #define __FUNCT__ "DMSNESSetFunction"
1656cab3a1bSJed Brown /*@C
1666cab3a1bSJed Brown    DMSNESSetFunction - set SNES residual evaluation function
1676cab3a1bSJed Brown 
1686cab3a1bSJed Brown    Not Collective
1696cab3a1bSJed Brown 
1706cab3a1bSJed Brown    Input Arguments:
1716cab3a1bSJed Brown +  dm - DM to be used with SNES
1726cab3a1bSJed Brown .  func - residual evaluation function, see SNESSetFunction() for calling sequence
1736cab3a1bSJed Brown -  ctx - context for residual evaluation
1746cab3a1bSJed Brown 
1756cab3a1bSJed Brown    Level: advanced
1766cab3a1bSJed Brown 
1776cab3a1bSJed Brown    Note:
1786cab3a1bSJed Brown    SNESSetFunction() is normally used, but it calls this function internally because the user context is actually
1796cab3a1bSJed Brown    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
1806cab3a1bSJed Brown    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
1816cab3a1bSJed Brown 
1826cab3a1bSJed Brown .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian()
1836cab3a1bSJed Brown @*/
1846cab3a1bSJed Brown PetscErrorCode DMSNESSetFunction(DM dm,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx)
1856cab3a1bSJed Brown {
1866cab3a1bSJed Brown   PetscErrorCode ierr;
1876cab3a1bSJed Brown   SNESDM sdm;
1886cab3a1bSJed Brown 
1896cab3a1bSJed Brown   PetscFunctionBegin;
1906cab3a1bSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1916cab3a1bSJed Brown   ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr);
1926cab3a1bSJed Brown   if (func) sdm->computefunction = func;
1936cab3a1bSJed Brown   if (ctx)  sdm->functionctx = ctx;
1946cab3a1bSJed Brown   PetscFunctionReturn(0);
1956cab3a1bSJed Brown }
1966cab3a1bSJed Brown 
1976cab3a1bSJed Brown #undef __FUNCT__
1986cab3a1bSJed Brown #define __FUNCT__ "DMSNESGetFunction"
1996cab3a1bSJed Brown /*@C
2006cab3a1bSJed Brown    DMSNESGetFunction - get SNES residual evaluation function
2016cab3a1bSJed Brown 
2026cab3a1bSJed Brown    Not Collective
2036cab3a1bSJed Brown 
2046cab3a1bSJed Brown    Input Argument:
2056cab3a1bSJed Brown .  dm - DM to be used with SNES
2066cab3a1bSJed Brown 
2076cab3a1bSJed Brown    Output Arguments:
2086cab3a1bSJed Brown +  func - residual evaluation function, see SNESSetFunction() for calling sequence
2096cab3a1bSJed Brown -  ctx - context for residual evaluation
2106cab3a1bSJed Brown 
2116cab3a1bSJed Brown    Level: advanced
2126cab3a1bSJed Brown 
2136cab3a1bSJed Brown    Note:
2146cab3a1bSJed Brown    SNESGetFunction() is normally used, but it calls this function internally because the user context is actually
2156cab3a1bSJed Brown    associated with the DM.
2166cab3a1bSJed Brown 
2176cab3a1bSJed Brown .seealso: DMSNESSetContext(), DMSNESSetFunction(), SNESSetFunction()
2186cab3a1bSJed Brown @*/
2196cab3a1bSJed Brown PetscErrorCode DMSNESGetFunction(DM dm,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx)
2206cab3a1bSJed Brown {
2216cab3a1bSJed Brown   PetscErrorCode ierr;
2226cab3a1bSJed Brown   SNESDM sdm;
2236cab3a1bSJed Brown 
2246cab3a1bSJed Brown   PetscFunctionBegin;
2256cab3a1bSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2266cab3a1bSJed Brown   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
2276cab3a1bSJed Brown   if (func) *func = sdm->computefunction;
2286cab3a1bSJed Brown   if (ctx)  *ctx = sdm->functionctx;
2296cab3a1bSJed Brown   PetscFunctionReturn(0);
2306cab3a1bSJed Brown }
2316cab3a1bSJed Brown 
2326cab3a1bSJed Brown #undef __FUNCT__
2336cab3a1bSJed Brown #define __FUNCT__ "DMSNESSetGS"
2346cab3a1bSJed Brown /*@C
2356cab3a1bSJed Brown    DMSNESSetGS - set SNES Gauss-Seidel relaxation function
2366cab3a1bSJed Brown 
2376cab3a1bSJed Brown    Not Collective
2386cab3a1bSJed Brown 
2396cab3a1bSJed Brown    Input Argument:
2406cab3a1bSJed Brown +  dm - DM to be used with SNES
2416cab3a1bSJed Brown .  func - relaxation function, see SNESSetGS() for calling sequence
2426cab3a1bSJed Brown -  ctx - context for residual evaluation
2436cab3a1bSJed Brown 
2446cab3a1bSJed Brown    Level: advanced
2456cab3a1bSJed Brown 
2466cab3a1bSJed Brown    Note:
2476cab3a1bSJed Brown    SNESSetGS() is normally used, but it calls this function internally because the user context is actually
2486cab3a1bSJed Brown    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
2496cab3a1bSJed Brown    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
2506cab3a1bSJed Brown 
2516cab3a1bSJed Brown .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), DMSNESSetFunction()
2526cab3a1bSJed Brown @*/
2536cab3a1bSJed Brown PetscErrorCode DMSNESSetGS(DM dm,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx)
2546cab3a1bSJed Brown {
2556cab3a1bSJed Brown   PetscErrorCode ierr;
2566cab3a1bSJed Brown   SNESDM sdm;
2576cab3a1bSJed Brown 
2586cab3a1bSJed Brown   PetscFunctionBegin;
2596cab3a1bSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2606cab3a1bSJed Brown   ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr);
2616cab3a1bSJed Brown   if (func) sdm->computegs = func;
2626cab3a1bSJed Brown   if (ctx)  sdm->gsctx = ctx;
2636cab3a1bSJed Brown   PetscFunctionReturn(0);
2646cab3a1bSJed Brown }
2656cab3a1bSJed Brown 
2666cab3a1bSJed Brown #undef __FUNCT__
2676cab3a1bSJed Brown #define __FUNCT__ "DMSNESGetGS"
2686cab3a1bSJed Brown /*@C
2696cab3a1bSJed Brown    DMSNESGetGS - get SNES Gauss-Seidel relaxation function
2706cab3a1bSJed Brown 
2716cab3a1bSJed Brown    Not Collective
2726cab3a1bSJed Brown 
2736cab3a1bSJed Brown    Input Argument:
2746cab3a1bSJed Brown .  dm - DM to be used with SNES
2756cab3a1bSJed Brown 
2766cab3a1bSJed Brown    Output Arguments:
2776cab3a1bSJed Brown +  func - relaxation function, see SNESSetGS() for calling sequence
2786cab3a1bSJed Brown -  ctx - context for residual evaluation
2796cab3a1bSJed Brown 
2806cab3a1bSJed Brown    Level: advanced
2816cab3a1bSJed Brown 
2826cab3a1bSJed Brown    Note:
2836cab3a1bSJed Brown    SNESGetGS() is normally used, but it calls this function internally because the user context is actually
2846cab3a1bSJed Brown    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
2856cab3a1bSJed Brown    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
2866cab3a1bSJed Brown 
2876cab3a1bSJed Brown .seealso: DMSNESSetContext(), SNESGetGS(), DMSNESGetJacobian(), DMSNESGetFunction()
2886cab3a1bSJed Brown @*/
2896cab3a1bSJed Brown PetscErrorCode DMSNESGetGS(DM dm,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx)
2906cab3a1bSJed Brown {
2916cab3a1bSJed Brown   PetscErrorCode ierr;
2926cab3a1bSJed Brown   SNESDM sdm;
2936cab3a1bSJed Brown 
2946cab3a1bSJed Brown   PetscFunctionBegin;
2956cab3a1bSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2966cab3a1bSJed Brown   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
2976cab3a1bSJed Brown   if (func) *func = sdm->computegs;
2986cab3a1bSJed Brown   if (ctx)  *ctx = sdm->gsctx;
2996cab3a1bSJed Brown   PetscFunctionReturn(0);
3006cab3a1bSJed Brown }
3016cab3a1bSJed Brown 
3026cab3a1bSJed Brown #undef __FUNCT__
3036cab3a1bSJed Brown #define __FUNCT__ "DMSNESSetJacobian"
3046cab3a1bSJed Brown /*@C
3056cab3a1bSJed Brown    DMSNESSetFunction - set SNES Jacobian evaluation function
3066cab3a1bSJed Brown 
3076cab3a1bSJed Brown    Not Collective
3086cab3a1bSJed Brown 
3096cab3a1bSJed Brown    Input Argument:
3106cab3a1bSJed Brown +  dm - DM to be used with SNES
3116cab3a1bSJed Brown .  func - Jacobian evaluation function, see SNESSetJacobian() for calling sequence
3126cab3a1bSJed Brown -  ctx - context for residual evaluation
3136cab3a1bSJed Brown 
3146cab3a1bSJed Brown    Level: advanced
3156cab3a1bSJed Brown 
3166cab3a1bSJed Brown    Note:
3176cab3a1bSJed Brown    SNESSetJacobian() is normally used, but it calls this function internally because the user context is actually
3186cab3a1bSJed Brown    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
3196cab3a1bSJed Brown    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
3206cab3a1bSJed Brown 
3216cab3a1bSJed Brown .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESGetJacobian(), SNESSetJacobian()
3226cab3a1bSJed Brown @*/
3236cab3a1bSJed Brown PetscErrorCode DMSNESSetJacobian(DM dm,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
3246cab3a1bSJed Brown {
3256cab3a1bSJed Brown   PetscErrorCode ierr;
3266cab3a1bSJed Brown   SNESDM sdm;
3276cab3a1bSJed Brown 
3286cab3a1bSJed Brown   PetscFunctionBegin;
3296cab3a1bSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3306cab3a1bSJed Brown   ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr);
3316cab3a1bSJed Brown   if (func) sdm->computejacobian = func;
3326cab3a1bSJed Brown   if (ctx)  sdm->jacobianctx = ctx;
3336cab3a1bSJed Brown   PetscFunctionReturn(0);
3346cab3a1bSJed Brown }
3356cab3a1bSJed Brown 
3366cab3a1bSJed Brown #undef __FUNCT__
3376cab3a1bSJed Brown #define __FUNCT__ "DMSNESGetJacobian"
3386cab3a1bSJed Brown /*@C
3396cab3a1bSJed Brown    DMSNESGetFunction - get SNES Jacobian evaluation function
3406cab3a1bSJed Brown 
3416cab3a1bSJed Brown    Not Collective
3426cab3a1bSJed Brown 
3436cab3a1bSJed Brown    Input Argument:
3446cab3a1bSJed Brown .  dm - DM to be used with SNES
3456cab3a1bSJed Brown 
3466cab3a1bSJed Brown    Output Arguments:
3476cab3a1bSJed Brown +  func - Jacobian evaluation function, see SNESSetJacobian() for calling sequence
3486cab3a1bSJed Brown -  ctx - context for residual evaluation
3496cab3a1bSJed Brown 
3506cab3a1bSJed Brown    Level: advanced
3516cab3a1bSJed Brown 
3526cab3a1bSJed Brown    Note:
3536cab3a1bSJed Brown    SNESGetJacobian() is normally used, but it calls this function internally because the user context is actually
3546cab3a1bSJed Brown    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
3556cab3a1bSJed Brown    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
3566cab3a1bSJed Brown 
3576cab3a1bSJed Brown .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian()
3586cab3a1bSJed Brown @*/
3596cab3a1bSJed Brown PetscErrorCode DMSNESGetJacobian(DM dm,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx)
3606cab3a1bSJed Brown {
3616cab3a1bSJed Brown   PetscErrorCode ierr;
3626cab3a1bSJed Brown   SNESDM sdm;
3636cab3a1bSJed Brown 
3646cab3a1bSJed Brown   PetscFunctionBegin;
3656cab3a1bSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3666cab3a1bSJed Brown   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
3676cab3a1bSJed Brown   if (func) *func = sdm->computejacobian;
3686cab3a1bSJed Brown   if (ctx)  *ctx = sdm->jacobianctx;
3696cab3a1bSJed Brown   PetscFunctionReturn(0);
3706cab3a1bSJed Brown }
3716cab3a1bSJed Brown 
3726cab3a1bSJed Brown #undef __FUNCT__
3736cab3a1bSJed Brown #define __FUNCT__ "SNESDefaultComputeFunction_DMLegacy"
3746cab3a1bSJed Brown static PetscErrorCode SNESDefaultComputeFunction_DMLegacy(SNES snes,Vec X,Vec F,void *ctx)
3756cab3a1bSJed Brown {
3766cab3a1bSJed Brown   PetscErrorCode ierr;
3776cab3a1bSJed Brown   DM             dm;
3786cab3a1bSJed Brown 
3796cab3a1bSJed Brown   PetscFunctionBegin;
3806cab3a1bSJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
3816cab3a1bSJed Brown   ierr = DMComputeFunction(dm,X,F);CHKERRQ(ierr);
3826cab3a1bSJed Brown   PetscFunctionReturn(0);
3836cab3a1bSJed Brown }
3846cab3a1bSJed Brown 
3856cab3a1bSJed Brown #undef __FUNCT__
3866cab3a1bSJed Brown #define __FUNCT__ "SNESDefaultComputeJacobian_DMLegacy"
3876cab3a1bSJed Brown static PetscErrorCode SNESDefaultComputeJacobian_DMLegacy(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *ctx)
3886cab3a1bSJed Brown {
3896cab3a1bSJed Brown   PetscErrorCode ierr;
3906cab3a1bSJed Brown   DM             dm;
3916cab3a1bSJed Brown 
3926cab3a1bSJed Brown   PetscFunctionBegin;
3936cab3a1bSJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
3946cab3a1bSJed Brown   ierr = DMComputeJacobian(dm,X,*A,*B,mstr);CHKERRQ(ierr);
3956cab3a1bSJed Brown   PetscFunctionReturn(0);
3966cab3a1bSJed Brown }
3976cab3a1bSJed Brown 
3986cab3a1bSJed Brown #undef __FUNCT__
3996cab3a1bSJed Brown #define __FUNCT__ "DMSNESSetUpLegacy"
4006cab3a1bSJed Brown /* Sets up calling of legacy DM routines */
4016cab3a1bSJed Brown PetscErrorCode DMSNESSetUpLegacy(DM dm)
4026cab3a1bSJed Brown {
4036cab3a1bSJed Brown   PetscErrorCode ierr;
4046cab3a1bSJed Brown   SNESDM         sdm;
4056cab3a1bSJed Brown 
4066cab3a1bSJed Brown   PetscFunctionBegin;
4076cab3a1bSJed Brown   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
4086cab3a1bSJed Brown   if (!sdm->computefunction) {ierr = DMSNESSetFunction(dm,SNESDefaultComputeFunction_DMLegacy,PETSC_NULL);CHKERRQ(ierr);}
4096cab3a1bSJed Brown   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
4106cab3a1bSJed Brown   if (!sdm->computejacobian) {ierr = DMSNESSetJacobian(dm,SNESDefaultComputeJacobian_DMLegacy,PETSC_NULL);CHKERRQ(ierr);}
4116cab3a1bSJed Brown   PetscFunctionReturn(0);
4126cab3a1bSJed Brown }
413