xref: /petsc/src/snes/utils/dmsnes.c (revision caa4e7f26594b8daaf440f4757330959cfc6cef4)
16cab3a1bSJed Brown #include <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;
12*caa4e7f2SJed Brown   if (sdm->dmloopref) {         /* Reconnect old loop before we drop our reference */
13*caa4e7f2SJed Brown     DM dminner;
14*caa4e7f2SJed Brown     ierr = PetscObjectQuery((PetscObject)sdm->vec_sol,"DM",(PetscObject*)&dminner);CHKERRQ(ierr);
15*caa4e7f2SJed Brown     if (dminner) {ierr = PetscObjectReference((PetscObject)dminner);CHKERRQ(ierr);}
16*caa4e7f2SJed Brown   }
17*caa4e7f2SJed Brown   ierr = VecDestroy(&sdm->vec_sol);CHKERRQ(ierr);
186cab3a1bSJed Brown   if (sdm->destroy) {ierr = (*sdm->destroy)(sdm);CHKERRQ(ierr);}
196cab3a1bSJed Brown   ierr = PetscFree(sdm);CHKERRQ(ierr);
206cab3a1bSJed Brown   PetscFunctionReturn(0);
216cab3a1bSJed Brown }
226cab3a1bSJed Brown 
236cab3a1bSJed Brown #undef __FUNCT__
246cab3a1bSJed Brown #define __FUNCT__ "DMCoarsenHook_SNESDM"
256cab3a1bSJed Brown /* Attaches the SNESDM to the coarse level.
266cab3a1bSJed Brown  * Under what conditions should we copy versus duplicate?
276cab3a1bSJed Brown  */
286cab3a1bSJed Brown static PetscErrorCode DMCoarsenHook_SNESDM(DM dm,DM dmc,void *ctx)
296cab3a1bSJed Brown {
306cab3a1bSJed Brown   PetscErrorCode ierr;
316cab3a1bSJed Brown 
326cab3a1bSJed Brown   PetscFunctionBegin;
336cab3a1bSJed Brown   ierr = DMSNESCopyContext(dm,dmc);CHKERRQ(ierr);
346cab3a1bSJed Brown   PetscFunctionReturn(0);
356cab3a1bSJed Brown }
366cab3a1bSJed Brown 
376cab3a1bSJed Brown #undef __FUNCT__
38*caa4e7f2SJed Brown #define __FUNCT__ "DMRestrictHook_SNESDM"
39*caa4e7f2SJed Brown /* Restrict state vector from finest level
40*caa4e7f2SJed Brown  */
41*caa4e7f2SJed Brown static PetscErrorCode DMRestrictHook_SNESDM(DM dm,Mat Restrict,Vec rscale,Mat Inject,DM dmc,void *ctx)
42*caa4e7f2SJed Brown {
43*caa4e7f2SJed Brown   PetscErrorCode ierr;
44*caa4e7f2SJed Brown   SNESDM sdm,sdmc;
45*caa4e7f2SJed Brown   Vec Xcoarse;
46*caa4e7f2SJed Brown 
47*caa4e7f2SJed Brown   PetscFunctionBegin;
48*caa4e7f2SJed Brown   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
49*caa4e7f2SJed Brown   ierr = DMSNESGetContextWrite(dmc,&sdmc);CHKERRQ(ierr);
50*caa4e7f2SJed Brown   if (sdm->vec_sol) {
51*caa4e7f2SJed Brown     ierr = DMSNESGetSolution(dmc,&Xcoarse);CHKERRQ(ierr);
52*caa4e7f2SJed Brown     ierr = MatRestrict(Restrict,sdm->vec_sol,Xcoarse);CHKERRQ(ierr);
53*caa4e7f2SJed Brown     ierr = VecPointwiseMult(Xcoarse,Xcoarse,rscale);CHKERRQ(ierr);
54*caa4e7f2SJed Brown   }
55*caa4e7f2SJed Brown   PetscFunctionReturn(0);
56*caa4e7f2SJed Brown }
57*caa4e7f2SJed Brown 
58*caa4e7f2SJed Brown #undef __FUNCT__
596cab3a1bSJed Brown #define __FUNCT__ "DMSNESGetContext"
606cab3a1bSJed Brown /*@C
616cab3a1bSJed Brown    DMSNESGetContext - get read-only private SNESDM context from a DM
626cab3a1bSJed Brown 
636cab3a1bSJed Brown    Not Collective
646cab3a1bSJed Brown 
656cab3a1bSJed Brown    Input Argument:
666cab3a1bSJed Brown .  dm - DM to be used with SNES
676cab3a1bSJed Brown 
686cab3a1bSJed Brown    Output Argument:
696cab3a1bSJed Brown .  snesdm - private SNESDM context
706cab3a1bSJed Brown 
716cab3a1bSJed Brown    Level: developer
726cab3a1bSJed Brown 
736cab3a1bSJed Brown    Notes:
746cab3a1bSJed Brown    Use DMSNESGetContextWrite() if write access is needed. The DMSNESSetXXX API should be used wherever possible.
756cab3a1bSJed Brown 
766cab3a1bSJed Brown .seealso: DMSNESGetContextWrite()
776cab3a1bSJed Brown @*/
786cab3a1bSJed Brown PetscErrorCode DMSNESGetContext(DM dm,SNESDM *snesdm)
796cab3a1bSJed Brown {
806cab3a1bSJed Brown   PetscErrorCode ierr;
816cab3a1bSJed Brown   PetscContainer container;
826cab3a1bSJed Brown   SNESDM         sdm;
836cab3a1bSJed Brown 
846cab3a1bSJed Brown   PetscFunctionBegin;
856cab3a1bSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
866cab3a1bSJed Brown   ierr = PetscObjectQuery((PetscObject)dm,"SNESDM",(PetscObject*)&container);CHKERRQ(ierr);
876cab3a1bSJed Brown   if (container) {
886cab3a1bSJed Brown     ierr = PetscContainerGetPointer(container,(void**)snesdm);CHKERRQ(ierr);
896cab3a1bSJed Brown   } else {
906cab3a1bSJed Brown     ierr = PetscInfo(dm,"Creating new SNESDM\n");CHKERRQ(ierr);
916cab3a1bSJed Brown     ierr = PetscContainerCreate(((PetscObject)dm)->comm,&container);CHKERRQ(ierr);
926cab3a1bSJed Brown     ierr = PetscNewLog(dm,struct _n_SNESDM,&sdm);CHKERRQ(ierr);
936cab3a1bSJed Brown     ierr = PetscContainerSetPointer(container,sdm);CHKERRQ(ierr);
946cab3a1bSJed Brown     ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_SNESDM);CHKERRQ(ierr);
956cab3a1bSJed Brown     ierr = PetscObjectCompose((PetscObject)dm,"SNESDM",(PetscObject)container);CHKERRQ(ierr);
96*caa4e7f2SJed Brown     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_SNESDM,DMRestrictHook_SNESDM,PETSC_NULL);CHKERRQ(ierr);
976cab3a1bSJed Brown     ierr = PetscContainerGetPointer(container,(void**)snesdm);CHKERRQ(ierr);
986cab3a1bSJed Brown     ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
996cab3a1bSJed Brown   }
1006cab3a1bSJed Brown   PetscFunctionReturn(0);
1016cab3a1bSJed Brown }
1026cab3a1bSJed Brown 
1036cab3a1bSJed Brown #undef __FUNCT__
1046cab3a1bSJed Brown #define __FUNCT__ "DMSNESGetContextWrite"
1056cab3a1bSJed Brown /*@C
1066cab3a1bSJed Brown    DMSNESGetContextWrite - get write access to private SNESDM context from a DM
1076cab3a1bSJed Brown 
1086cab3a1bSJed Brown    Not Collective
1096cab3a1bSJed Brown 
1106cab3a1bSJed Brown    Input Argument:
1116cab3a1bSJed Brown .  dm - DM to be used with SNES
1126cab3a1bSJed Brown 
1136cab3a1bSJed Brown    Output Argument:
1146cab3a1bSJed Brown .  snesdm - private SNESDM context
1156cab3a1bSJed Brown 
1166cab3a1bSJed Brown    Level: developer
1176cab3a1bSJed Brown 
1186cab3a1bSJed Brown .seealso: DMSNESGetContext()
1196cab3a1bSJed Brown @*/
1206cab3a1bSJed Brown PetscErrorCode DMSNESGetContextWrite(DM dm,SNESDM *snesdm)
1216cab3a1bSJed Brown {
1226cab3a1bSJed Brown   PetscErrorCode ierr;
1236cab3a1bSJed Brown   SNESDM         sdm;
1246cab3a1bSJed Brown 
1256cab3a1bSJed Brown   PetscFunctionBegin;
1266cab3a1bSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1276cab3a1bSJed Brown   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
1286cab3a1bSJed Brown   if (!sdm->originaldm) sdm->originaldm = dm;
1296cab3a1bSJed Brown   if (sdm->originaldm != dm) {  /* Copy on write */
1306cab3a1bSJed Brown     PetscContainer container;
1316cab3a1bSJed Brown     SNESDM         oldsdm = sdm;
1326cab3a1bSJed Brown     ierr = PetscInfo(dm,"Copying SNESDM due to write\n");CHKERRQ(ierr);
1336cab3a1bSJed Brown     ierr = PetscContainerCreate(((PetscObject)dm)->comm,&container);CHKERRQ(ierr);
1346cab3a1bSJed Brown     ierr = PetscNewLog(dm,struct _n_SNESDM,&sdm);CHKERRQ(ierr);
1356cab3a1bSJed Brown     ierr = PetscMemcpy(sdm,oldsdm,sizeof *sdm);CHKERRQ(ierr);
136*caa4e7f2SJed Brown     sdm->vec_sol = PETSC_NULL;
137*caa4e7f2SJed Brown     sdm->dmloopref = PETSC_FALSE;
1386cab3a1bSJed Brown     ierr = PetscContainerSetPointer(container,sdm);CHKERRQ(ierr);
1396cab3a1bSJed Brown     ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_SNESDM);CHKERRQ(ierr);
1406cab3a1bSJed Brown     ierr = PetscObjectCompose((PetscObject)dm,"SNESDM",(PetscObject)container);CHKERRQ(ierr);
1416cab3a1bSJed Brown     ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
1426cab3a1bSJed Brown   }
1436cab3a1bSJed Brown   *snesdm = sdm;
1446cab3a1bSJed Brown   PetscFunctionReturn(0);
1456cab3a1bSJed Brown }
1466cab3a1bSJed Brown 
1476cab3a1bSJed Brown #undef __FUNCT__
1486cab3a1bSJed Brown #define __FUNCT__ "DMSNESCopyContext"
1496cab3a1bSJed Brown /*@C
1506cab3a1bSJed Brown    DMSNESCopyContext - copies a DM context to a new DM
1516cab3a1bSJed Brown 
1526cab3a1bSJed Brown    Logically Collective
1536cab3a1bSJed Brown 
1546cab3a1bSJed Brown    Input Arguments:
1556cab3a1bSJed Brown +  dmsrc - DM to obtain context from
1566cab3a1bSJed Brown -  dmdest - DM to add context to
1576cab3a1bSJed Brown 
1586cab3a1bSJed Brown    Level: developer
1596cab3a1bSJed Brown 
1606cab3a1bSJed Brown    Note:
1616cab3a1bSJed Brown    The context is copied by reference. This function does not ensure that a context exists.
1626cab3a1bSJed Brown 
1636cab3a1bSJed Brown .seealso: DMSNESGetContext(), SNESSetDM()
1646cab3a1bSJed Brown @*/
1656cab3a1bSJed Brown PetscErrorCode DMSNESCopyContext(DM dmsrc,DM dmdest)
1666cab3a1bSJed Brown {
1676cab3a1bSJed Brown   PetscErrorCode ierr;
1686cab3a1bSJed Brown   PetscContainer container;
1696cab3a1bSJed Brown 
1706cab3a1bSJed Brown   PetscFunctionBegin;
1716cab3a1bSJed Brown   PetscValidHeaderSpecific(dmsrc,DM_CLASSID,1);
1726cab3a1bSJed Brown   PetscValidHeaderSpecific(dmdest,DM_CLASSID,2);
1736cab3a1bSJed Brown   ierr = PetscObjectQuery((PetscObject)dmsrc,"SNESDM",(PetscObject*)&container);CHKERRQ(ierr);
1746cab3a1bSJed Brown   if (container) {
1756cab3a1bSJed Brown     ierr = PetscObjectCompose((PetscObject)dmdest,"SNESDM",(PetscObject)container);CHKERRQ(ierr);
176*caa4e7f2SJed Brown     ierr = DMCoarsenHookAdd(dmdest,DMCoarsenHook_SNESDM,DMRestrictHook_SNESDM,PETSC_NULL);CHKERRQ(ierr);
1776cab3a1bSJed Brown   }
1786cab3a1bSJed Brown   PetscFunctionReturn(0);
1796cab3a1bSJed Brown }
1806cab3a1bSJed Brown 
1816cab3a1bSJed Brown #undef __FUNCT__
1826cab3a1bSJed Brown #define __FUNCT__ "DMSNESSetFunction"
1836cab3a1bSJed Brown /*@C
1846cab3a1bSJed Brown    DMSNESSetFunction - set SNES residual evaluation function
1856cab3a1bSJed Brown 
1866cab3a1bSJed Brown    Not Collective
1876cab3a1bSJed Brown 
1886cab3a1bSJed Brown    Input Arguments:
1896cab3a1bSJed Brown +  dm - DM to be used with SNES
1906cab3a1bSJed Brown .  func - residual evaluation function, see SNESSetFunction() for calling sequence
1916cab3a1bSJed Brown -  ctx - context for residual evaluation
1926cab3a1bSJed Brown 
1936cab3a1bSJed Brown    Level: advanced
1946cab3a1bSJed Brown 
1956cab3a1bSJed Brown    Note:
1966cab3a1bSJed Brown    SNESSetFunction() is normally used, but it calls this function internally because the user context is actually
1976cab3a1bSJed Brown    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
1986cab3a1bSJed Brown    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
1996cab3a1bSJed Brown 
2006cab3a1bSJed Brown .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian()
2016cab3a1bSJed Brown @*/
2026cab3a1bSJed Brown PetscErrorCode DMSNESSetFunction(DM dm,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx)
2036cab3a1bSJed Brown {
2046cab3a1bSJed Brown   PetscErrorCode ierr;
2056cab3a1bSJed Brown   SNESDM sdm;
2066cab3a1bSJed Brown 
2076cab3a1bSJed Brown   PetscFunctionBegin;
2086cab3a1bSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2096cab3a1bSJed Brown   ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr);
2106cab3a1bSJed Brown   if (func) sdm->computefunction = func;
2116cab3a1bSJed Brown   if (ctx)  sdm->functionctx = ctx;
2126cab3a1bSJed Brown   PetscFunctionReturn(0);
2136cab3a1bSJed Brown }
2146cab3a1bSJed Brown 
2156cab3a1bSJed Brown #undef __FUNCT__
2166cab3a1bSJed Brown #define __FUNCT__ "DMSNESGetFunction"
2176cab3a1bSJed Brown /*@C
2186cab3a1bSJed Brown    DMSNESGetFunction - get SNES residual evaluation function
2196cab3a1bSJed Brown 
2206cab3a1bSJed Brown    Not Collective
2216cab3a1bSJed Brown 
2226cab3a1bSJed Brown    Input Argument:
2236cab3a1bSJed Brown .  dm - DM to be used with SNES
2246cab3a1bSJed Brown 
2256cab3a1bSJed Brown    Output Arguments:
2266cab3a1bSJed Brown +  func - residual evaluation function, see SNESSetFunction() for calling sequence
2276cab3a1bSJed Brown -  ctx - context for residual evaluation
2286cab3a1bSJed Brown 
2296cab3a1bSJed Brown    Level: advanced
2306cab3a1bSJed Brown 
2316cab3a1bSJed Brown    Note:
2326cab3a1bSJed Brown    SNESGetFunction() is normally used, but it calls this function internally because the user context is actually
2336cab3a1bSJed Brown    associated with the DM.
2346cab3a1bSJed Brown 
2356cab3a1bSJed Brown .seealso: DMSNESSetContext(), DMSNESSetFunction(), SNESSetFunction()
2366cab3a1bSJed Brown @*/
2376cab3a1bSJed Brown PetscErrorCode DMSNESGetFunction(DM dm,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx)
2386cab3a1bSJed Brown {
2396cab3a1bSJed Brown   PetscErrorCode ierr;
2406cab3a1bSJed Brown   SNESDM sdm;
2416cab3a1bSJed Brown 
2426cab3a1bSJed Brown   PetscFunctionBegin;
2436cab3a1bSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2446cab3a1bSJed Brown   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
2456cab3a1bSJed Brown   if (func) *func = sdm->computefunction;
2466cab3a1bSJed Brown   if (ctx)  *ctx = sdm->functionctx;
2476cab3a1bSJed Brown   PetscFunctionReturn(0);
2486cab3a1bSJed Brown }
2496cab3a1bSJed Brown 
2506cab3a1bSJed Brown #undef __FUNCT__
2516cab3a1bSJed Brown #define __FUNCT__ "DMSNESSetGS"
2526cab3a1bSJed Brown /*@C
2536cab3a1bSJed Brown    DMSNESSetGS - set SNES Gauss-Seidel relaxation function
2546cab3a1bSJed Brown 
2556cab3a1bSJed Brown    Not Collective
2566cab3a1bSJed Brown 
2576cab3a1bSJed Brown    Input Argument:
2586cab3a1bSJed Brown +  dm - DM to be used with SNES
2596cab3a1bSJed Brown .  func - relaxation function, see SNESSetGS() for calling sequence
2606cab3a1bSJed Brown -  ctx - context for residual evaluation
2616cab3a1bSJed Brown 
2626cab3a1bSJed Brown    Level: advanced
2636cab3a1bSJed Brown 
2646cab3a1bSJed Brown    Note:
2656cab3a1bSJed Brown    SNESSetGS() is normally used, but it calls this function internally because the user context is actually
2666cab3a1bSJed Brown    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
2676cab3a1bSJed Brown    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
2686cab3a1bSJed Brown 
2696cab3a1bSJed Brown .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), DMSNESSetFunction()
2706cab3a1bSJed Brown @*/
2716cab3a1bSJed Brown PetscErrorCode DMSNESSetGS(DM dm,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx)
2726cab3a1bSJed Brown {
2736cab3a1bSJed Brown   PetscErrorCode ierr;
2746cab3a1bSJed Brown   SNESDM sdm;
2756cab3a1bSJed Brown 
2766cab3a1bSJed Brown   PetscFunctionBegin;
2776cab3a1bSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2786cab3a1bSJed Brown   ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr);
2796cab3a1bSJed Brown   if (func) sdm->computegs = func;
2806cab3a1bSJed Brown   if (ctx)  sdm->gsctx = ctx;
2816cab3a1bSJed Brown   PetscFunctionReturn(0);
2826cab3a1bSJed Brown }
2836cab3a1bSJed Brown 
2846cab3a1bSJed Brown #undef __FUNCT__
2856cab3a1bSJed Brown #define __FUNCT__ "DMSNESGetGS"
2866cab3a1bSJed Brown /*@C
2876cab3a1bSJed Brown    DMSNESGetGS - get SNES Gauss-Seidel relaxation function
2886cab3a1bSJed Brown 
2896cab3a1bSJed Brown    Not Collective
2906cab3a1bSJed Brown 
2916cab3a1bSJed Brown    Input Argument:
2926cab3a1bSJed Brown .  dm - DM to be used with SNES
2936cab3a1bSJed Brown 
2946cab3a1bSJed Brown    Output Arguments:
2956cab3a1bSJed Brown +  func - relaxation function, see SNESSetGS() for calling sequence
2966cab3a1bSJed Brown -  ctx - context for residual evaluation
2976cab3a1bSJed Brown 
2986cab3a1bSJed Brown    Level: advanced
2996cab3a1bSJed Brown 
3006cab3a1bSJed Brown    Note:
3016cab3a1bSJed Brown    SNESGetGS() is normally used, but it calls this function internally because the user context is actually
3026cab3a1bSJed Brown    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
3036cab3a1bSJed Brown    not. If DM took a more central role at some later date, this could become the primary method of setting the residual.
3046cab3a1bSJed Brown 
3056cab3a1bSJed Brown .seealso: DMSNESSetContext(), SNESGetGS(), DMSNESGetJacobian(), DMSNESGetFunction()
3066cab3a1bSJed Brown @*/
3076cab3a1bSJed Brown PetscErrorCode DMSNESGetGS(DM dm,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx)
3086cab3a1bSJed Brown {
3096cab3a1bSJed Brown   PetscErrorCode ierr;
3106cab3a1bSJed Brown   SNESDM sdm;
3116cab3a1bSJed Brown 
3126cab3a1bSJed Brown   PetscFunctionBegin;
3136cab3a1bSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3146cab3a1bSJed Brown   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
3156cab3a1bSJed Brown   if (func) *func = sdm->computegs;
3166cab3a1bSJed Brown   if (ctx)  *ctx = sdm->gsctx;
3176cab3a1bSJed Brown   PetscFunctionReturn(0);
3186cab3a1bSJed Brown }
3196cab3a1bSJed Brown 
3206cab3a1bSJed Brown #undef __FUNCT__
3216cab3a1bSJed Brown #define __FUNCT__ "DMSNESSetJacobian"
3226cab3a1bSJed Brown /*@C
3236cab3a1bSJed Brown    DMSNESSetFunction - set SNES Jacobian evaluation function
3246cab3a1bSJed Brown 
3256cab3a1bSJed Brown    Not Collective
3266cab3a1bSJed Brown 
3276cab3a1bSJed Brown    Input Argument:
3286cab3a1bSJed Brown +  dm - DM to be used with SNES
3296cab3a1bSJed Brown .  func - Jacobian evaluation function, see SNESSetJacobian() for calling sequence
3306cab3a1bSJed Brown -  ctx - context for residual evaluation
3316cab3a1bSJed Brown 
3326cab3a1bSJed Brown    Level: advanced
3336cab3a1bSJed Brown 
3346cab3a1bSJed Brown    Note:
3356cab3a1bSJed Brown    SNESSetJacobian() is normally used, but it calls this function internally because the user context is actually
3366cab3a1bSJed Brown    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
3376cab3a1bSJed Brown    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
3386cab3a1bSJed Brown 
3396cab3a1bSJed Brown .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESGetJacobian(), SNESSetJacobian()
3406cab3a1bSJed Brown @*/
3416cab3a1bSJed Brown PetscErrorCode DMSNESSetJacobian(DM dm,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
3426cab3a1bSJed Brown {
3436cab3a1bSJed Brown   PetscErrorCode ierr;
3446cab3a1bSJed Brown   SNESDM sdm;
3456cab3a1bSJed Brown 
3466cab3a1bSJed Brown   PetscFunctionBegin;
3476cab3a1bSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3486cab3a1bSJed Brown   ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr);
3496cab3a1bSJed Brown   if (func) sdm->computejacobian = func;
3506cab3a1bSJed Brown   if (ctx)  sdm->jacobianctx = ctx;
3516cab3a1bSJed Brown   PetscFunctionReturn(0);
3526cab3a1bSJed Brown }
3536cab3a1bSJed Brown 
3546cab3a1bSJed Brown #undef __FUNCT__
3556cab3a1bSJed Brown #define __FUNCT__ "DMSNESGetJacobian"
3566cab3a1bSJed Brown /*@C
3576cab3a1bSJed Brown    DMSNESGetFunction - get SNES Jacobian evaluation function
3586cab3a1bSJed Brown 
3596cab3a1bSJed Brown    Not Collective
3606cab3a1bSJed Brown 
3616cab3a1bSJed Brown    Input Argument:
3626cab3a1bSJed Brown .  dm - DM to be used with SNES
3636cab3a1bSJed Brown 
3646cab3a1bSJed Brown    Output Arguments:
3656cab3a1bSJed Brown +  func - Jacobian evaluation function, see SNESSetJacobian() for calling sequence
3666cab3a1bSJed Brown -  ctx - context for residual evaluation
3676cab3a1bSJed Brown 
3686cab3a1bSJed Brown    Level: advanced
3696cab3a1bSJed Brown 
3706cab3a1bSJed Brown    Note:
3716cab3a1bSJed Brown    SNESGetJacobian() is normally used, but it calls this function internally because the user context is actually
3726cab3a1bSJed Brown    associated with the DM.  This makes the interface consistent regardless of whether the user interacts with a DM or
3736cab3a1bSJed Brown    not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian.
3746cab3a1bSJed Brown 
3756cab3a1bSJed Brown .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian()
3766cab3a1bSJed Brown @*/
3776cab3a1bSJed Brown PetscErrorCode DMSNESGetJacobian(DM dm,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx)
3786cab3a1bSJed Brown {
3796cab3a1bSJed Brown   PetscErrorCode ierr;
3806cab3a1bSJed Brown   SNESDM sdm;
3816cab3a1bSJed Brown 
3826cab3a1bSJed Brown   PetscFunctionBegin;
3836cab3a1bSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3846cab3a1bSJed Brown   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
3856cab3a1bSJed Brown   if (func) *func = sdm->computejacobian;
3866cab3a1bSJed Brown   if (ctx)  *ctx = sdm->jacobianctx;
3876cab3a1bSJed Brown   PetscFunctionReturn(0);
3886cab3a1bSJed Brown }
3896cab3a1bSJed Brown 
3906cab3a1bSJed Brown #undef __FUNCT__
3916cab3a1bSJed Brown #define __FUNCT__ "SNESDefaultComputeFunction_DMLegacy"
3926cab3a1bSJed Brown static PetscErrorCode SNESDefaultComputeFunction_DMLegacy(SNES snes,Vec X,Vec F,void *ctx)
3936cab3a1bSJed Brown {
3946cab3a1bSJed Brown   PetscErrorCode ierr;
3956cab3a1bSJed Brown   DM             dm;
3966cab3a1bSJed Brown 
3976cab3a1bSJed Brown   PetscFunctionBegin;
3986cab3a1bSJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
3996cab3a1bSJed Brown   ierr = DMComputeFunction(dm,X,F);CHKERRQ(ierr);
4006cab3a1bSJed Brown   PetscFunctionReturn(0);
4016cab3a1bSJed Brown }
4026cab3a1bSJed Brown 
4036cab3a1bSJed Brown #undef __FUNCT__
4046cab3a1bSJed Brown #define __FUNCT__ "SNESDefaultComputeJacobian_DMLegacy"
4056cab3a1bSJed Brown static PetscErrorCode SNESDefaultComputeJacobian_DMLegacy(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *ctx)
4066cab3a1bSJed Brown {
4076cab3a1bSJed Brown   PetscErrorCode ierr;
4086cab3a1bSJed Brown   DM             dm;
4096cab3a1bSJed Brown 
4106cab3a1bSJed Brown   PetscFunctionBegin;
4116cab3a1bSJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
4126cab3a1bSJed Brown   ierr = DMComputeJacobian(dm,X,*A,*B,mstr);CHKERRQ(ierr);
4136cab3a1bSJed Brown   PetscFunctionReturn(0);
4146cab3a1bSJed Brown }
4156cab3a1bSJed Brown 
4166cab3a1bSJed Brown #undef __FUNCT__
4176cab3a1bSJed Brown #define __FUNCT__ "DMSNESSetUpLegacy"
4186cab3a1bSJed Brown /* Sets up calling of legacy DM routines */
4196cab3a1bSJed Brown PetscErrorCode DMSNESSetUpLegacy(DM dm)
4206cab3a1bSJed Brown {
4216cab3a1bSJed Brown   PetscErrorCode ierr;
4226cab3a1bSJed Brown   SNESDM         sdm;
4236cab3a1bSJed Brown 
4246cab3a1bSJed Brown   PetscFunctionBegin;
4256cab3a1bSJed Brown   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
4266cab3a1bSJed Brown   if (!sdm->computefunction) {ierr = DMSNESSetFunction(dm,SNESDefaultComputeFunction_DMLegacy,PETSC_NULL);CHKERRQ(ierr);}
4276cab3a1bSJed Brown   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
4286cab3a1bSJed Brown   if (!sdm->computejacobian) {ierr = DMSNESSetJacobian(dm,SNESDefaultComputeJacobian_DMLegacy,PETSC_NULL);CHKERRQ(ierr);}
4296cab3a1bSJed Brown   PetscFunctionReturn(0);
4306cab3a1bSJed Brown }
431*caa4e7f2SJed Brown 
432*caa4e7f2SJed Brown #undef __FUNCT__
433*caa4e7f2SJed Brown #define __FUNCT__ "DMSNESGetSolution"
434*caa4e7f2SJed Brown PetscErrorCode DMSNESGetSolution(DM dm,Vec *vec_sol)
435*caa4e7f2SJed Brown {
436*caa4e7f2SJed Brown   PetscErrorCode ierr;
437*caa4e7f2SJed Brown   SNESDM sdm;
438*caa4e7f2SJed Brown 
439*caa4e7f2SJed Brown   PetscFunctionBegin;
440*caa4e7f2SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
441*caa4e7f2SJed Brown   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
442*caa4e7f2SJed Brown   if (!sdm->vec_sol) {
443*caa4e7f2SJed Brown     ierr = DMCreateGlobalVector(dm,&sdm->vec_sol);CHKERRQ(ierr);
444*caa4e7f2SJed Brown     /* This vector holds an extra reference to the DM, decrement the DM reference count so that it remains constant */
445*caa4e7f2SJed Brown     sdm->dmloopref = PETSC_TRUE;
446*caa4e7f2SJed Brown     ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
447*caa4e7f2SJed Brown   }
448*caa4e7f2SJed Brown   *vec_sol = sdm->vec_sol;
449*caa4e7f2SJed Brown   PetscFunctionReturn(0);
450*caa4e7f2SJed Brown }
451*caa4e7f2SJed Brown 
452*caa4e7f2SJed Brown #undef __FUNCT__
453*caa4e7f2SJed Brown #define __FUNCT__ "DMSNESSetSolution"
454*caa4e7f2SJed Brown PetscErrorCode DMSNESSetSolution(DM dm,Vec vec_sol)
455*caa4e7f2SJed Brown {
456*caa4e7f2SJed Brown   PetscErrorCode ierr;
457*caa4e7f2SJed Brown   DM dminner;
458*caa4e7f2SJed Brown   SNESDM sdm;
459*caa4e7f2SJed Brown 
460*caa4e7f2SJed Brown   PetscFunctionBegin;
461*caa4e7f2SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
462*caa4e7f2SJed Brown   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
463*caa4e7f2SJed Brown   if (sdm->vec_sol == vec_sol) PetscFunctionReturn(0);
464*caa4e7f2SJed Brown   if (sdm->dmloopref) {         /* Reconnect old loop before we drop our reference */
465*caa4e7f2SJed Brown     ierr = PetscObjectQuery((PetscObject)sdm->vec_sol,"DM",(PetscObject*)&dminner);CHKERRQ(ierr);
466*caa4e7f2SJed Brown     if (dminner) {ierr = PetscObjectReference((PetscObject)dminner);CHKERRQ(ierr);}
467*caa4e7f2SJed Brown   }
468*caa4e7f2SJed Brown   sdm->dmloopref = PETSC_FALSE;
469*caa4e7f2SJed Brown 
470*caa4e7f2SJed Brown   if (vec_sol) {ierr = PetscObjectReference((PetscObject)vec_sol);CHKERRQ(ierr);}
471*caa4e7f2SJed Brown   ierr = VecDestroy(&sdm->vec_sol);CHKERRQ(ierr);
472*caa4e7f2SJed Brown   sdm->vec_sol = vec_sol;
473*caa4e7f2SJed Brown 
474*caa4e7f2SJed Brown   /* Break new loop if we just created one */
475*caa4e7f2SJed Brown   if (sdm->vec_sol) {
476*caa4e7f2SJed Brown     ierr = PetscObjectQuery((PetscObject)sdm->vec_sol,"DM",(PetscObject*)&dminner);CHKERRQ(ierr);
477*caa4e7f2SJed Brown     if (dminner == dm) {
478*caa4e7f2SJed Brown       sdm->dmloopref = PETSC_TRUE;
479*caa4e7f2SJed Brown       if (((PetscObject)dm)->refct <= 1) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONG,"The incoming Vec holds the last reference to the DM");
480*caa4e7f2SJed Brown       ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
481*caa4e7f2SJed Brown     }
482*caa4e7f2SJed Brown   }
483*caa4e7f2SJed Brown   PetscFunctionReturn(0);
484*caa4e7f2SJed Brown }
485