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__ 5*98120457SMatthew G Knepley #define __FUNCT__ "SNESDMComputeFunction" 6*98120457SMatthew G Knepley /*@C 7*98120457SMatthew G Knepley SNESDMComputeFunction - This is a universal function evaluation routine that 8*98120457SMatthew G Knepley may be used with SNESSetFunction() as long as the user context has a DM 9*98120457SMatthew G Knepley as its first record and the user has called DMSetLocalFunction(). 10*98120457SMatthew G Knepley 11*98120457SMatthew G Knepley Collective on SNES 12*98120457SMatthew G Knepley 13*98120457SMatthew G Knepley Input Parameters: 14*98120457SMatthew G Knepley + snes - the SNES context 15*98120457SMatthew G Knepley . X - input vector 16*98120457SMatthew G Knepley . F - function vector 17*98120457SMatthew G Knepley - ptr - pointer to a structure that must have a DM as its first entry. 18*98120457SMatthew G Knepley This ptr must have been passed into SNESDMComputeFunction() as the context. 19*98120457SMatthew G Knepley 20*98120457SMatthew G Knepley Level: intermediate 21*98120457SMatthew G Knepley 22*98120457SMatthew G Knepley .seealso: DMSetLocalFunction(), DMSetLocalJacobian(), SNESSetFunction(), SNESSetJacobian() 23*98120457SMatthew G Knepley @*/ 24*98120457SMatthew G Knepley PetscErrorCode SNESDMComputeFunction(SNES snes, Vec X, Vec F, void *ptr) 25*98120457SMatthew G Knepley { 26*98120457SMatthew G Knepley DM dm = *(DM*) ptr; 27*98120457SMatthew G Knepley PetscErrorCode (*lf)(DM, Vec, Vec, void *); 28*98120457SMatthew G Knepley Vec localX, localF; 29*98120457SMatthew G Knepley PetscInt N, n; 30*98120457SMatthew G Knepley PetscErrorCode ierr; 31*98120457SMatthew G Knepley 32*98120457SMatthew G Knepley PetscFunctionBegin; 33*98120457SMatthew G Knepley PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 34*98120457SMatthew G Knepley PetscValidHeaderSpecific(X, VEC_CLASSID, 2); 35*98120457SMatthew G Knepley PetscValidHeaderSpecific(F, VEC_CLASSID, 3); 36*98120457SMatthew G Knepley if (!dm) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONGSTATE, "Looks like you called SNESSetFromFuntion(snes,SNESDMComputeFunction,) without the DM context"); 37*98120457SMatthew G Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 4); 38*98120457SMatthew G Knepley 39*98120457SMatthew G Knepley /* determine whether X = localX */ 40*98120457SMatthew G Knepley ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 41*98120457SMatthew G Knepley ierr = DMGetLocalVector(dm, &localF);CHKERRQ(ierr); 42*98120457SMatthew G Knepley ierr = VecGetSize(X, &N);CHKERRQ(ierr); 43*98120457SMatthew G Knepley ierr = VecGetSize(localX, &n);CHKERRQ(ierr); 44*98120457SMatthew G Knepley 45*98120457SMatthew G Knepley if (n != N){ /* X != localX */ 46*98120457SMatthew G Knepley /* Scatter ghost points to local vector, using the 2-step process 47*98120457SMatthew G Knepley DMGlobalToLocalBegin(), DMGlobalToLocalEnd(). 48*98120457SMatthew G Knepley */ 49*98120457SMatthew G Knepley ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 50*98120457SMatthew G Knepley ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 51*98120457SMatthew G Knepley } else { 52*98120457SMatthew G Knepley ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 53*98120457SMatthew G Knepley localX = X; 54*98120457SMatthew G Knepley } 55*98120457SMatthew G Knepley ierr = DMGetLocalFunction(dm, &lf);CHKERRQ(ierr); 56*98120457SMatthew G Knepley ierr = (*lf)(dm, localX, localF, ptr);CHKERRQ(ierr); 57*98120457SMatthew G Knepley if (n != N){ 58*98120457SMatthew G Knepley ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 59*98120457SMatthew G Knepley } 60*98120457SMatthew G Knepley ierr = DMLocalToGlobalBegin(dm, localF, ADD_VALUES, F);CHKERRQ(ierr); 61*98120457SMatthew G Knepley ierr = DMLocalToGlobalEnd(dm, localF, ADD_VALUES, F);CHKERRQ(ierr); 62*98120457SMatthew G Knepley ierr = DMRestoreLocalVector(dm, &localF);CHKERRQ(ierr); 63*98120457SMatthew G Knepley PetscFunctionReturn(0); 64*98120457SMatthew G Knepley } 65*98120457SMatthew G Knepley 66*98120457SMatthew G Knepley #undef __FUNCT__ 67*98120457SMatthew G Knepley #define __FUNCT__ "SNESDMComputeJacobian" 68*98120457SMatthew G Knepley /* 69*98120457SMatthew G Knepley SNESDMComputeJacobian - This is a universal Jacobian evaluation routine that 70*98120457SMatthew G Knepley may be used with SNESSetJacobian() as long as the user context has a DM 71*98120457SMatthew G Knepley as its first record and the user has called DMSetLocalJacobian(). 72*98120457SMatthew G Knepley 73*98120457SMatthew G Knepley Collective on SNES 74*98120457SMatthew G Knepley 75*98120457SMatthew G Knepley Input Parameters: 76*98120457SMatthew G Knepley + snes - the SNES context 77*98120457SMatthew G Knepley . X - input vector 78*98120457SMatthew G Knepley . J - Jacobian 79*98120457SMatthew G Knepley . B - Jacobian used in preconditioner (usally same as J) 80*98120457SMatthew G Knepley . flag - indicates if the matrix changed its structure 81*98120457SMatthew G Knepley - ptr - pointer to a structure that must have a DM as its first entry. 82*98120457SMatthew G Knepley This ptr must have been passed into SNESDMComputeFunction() as the context. 83*98120457SMatthew G Knepley 84*98120457SMatthew G Knepley Level: intermediate 85*98120457SMatthew G Knepley 86*98120457SMatthew G Knepley .seealso: DMSetLocalFunction(), DMSetLocalJacobian(), SNESSetFunction(), SNESSetJacobian() 87*98120457SMatthew G Knepley */ 88*98120457SMatthew G Knepley PetscErrorCode SNESDMComputeJacobian(SNES snes, Vec X, Mat *J, Mat *B, MatStructure *flag, void *ptr) 89*98120457SMatthew G Knepley { 90*98120457SMatthew G Knepley DM dm = *(DM*) ptr; 91*98120457SMatthew G Knepley PetscErrorCode (*lj)(DM, Vec, Mat, Mat, void *); 92*98120457SMatthew G Knepley Vec localX; 93*98120457SMatthew G Knepley PetscErrorCode ierr; 94*98120457SMatthew G Knepley 95*98120457SMatthew G Knepley PetscFunctionBegin; 96*98120457SMatthew G Knepley ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 97*98120457SMatthew G Knepley ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 98*98120457SMatthew G Knepley ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 99*98120457SMatthew G Knepley ierr = DMGetLocalJacobian(dm, &lj);CHKERRQ(ierr); 100*98120457SMatthew G Knepley ierr = (*lj)(dm, localX, *J, *B, ptr);CHKERRQ(ierr); 101*98120457SMatthew G Knepley ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 102*98120457SMatthew G Knepley /* Assemble true Jacobian; if it is different */ 103*98120457SMatthew G Knepley if (*J != *B) { 104*98120457SMatthew G Knepley ierr = MatAssemblyBegin(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 105*98120457SMatthew G Knepley ierr = MatAssemblyEnd(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 106*98120457SMatthew G Knepley } 107*98120457SMatthew G Knepley ierr = MatSetOption(*B, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE);CHKERRQ(ierr); 108*98120457SMatthew G Knepley *flag = SAME_NONZERO_PATTERN; 109*98120457SMatthew G Knepley PetscFunctionReturn(0); 110*98120457SMatthew G Knepley } 111*98120457SMatthew G Knepley 112*98120457SMatthew G Knepley #undef __FUNCT__ 1136cab3a1bSJed Brown #define __FUNCT__ "PetscContainerDestroy_SNESDM" 1146cab3a1bSJed Brown static PetscErrorCode PetscContainerDestroy_SNESDM(void *ctx) 1156cab3a1bSJed Brown { 1166cab3a1bSJed Brown PetscErrorCode ierr; 1176cab3a1bSJed Brown SNESDM sdm = (SNESDM)ctx; 1186cab3a1bSJed Brown 1196cab3a1bSJed Brown PetscFunctionBegin; 1206cab3a1bSJed Brown if (sdm->destroy) {ierr = (*sdm->destroy)(sdm);CHKERRQ(ierr);} 1216cab3a1bSJed Brown ierr = PetscFree(sdm);CHKERRQ(ierr); 1226cab3a1bSJed Brown PetscFunctionReturn(0); 1236cab3a1bSJed Brown } 1246cab3a1bSJed Brown 1256cab3a1bSJed Brown #undef __FUNCT__ 1266cab3a1bSJed Brown #define __FUNCT__ "DMCoarsenHook_SNESDM" 1276cab3a1bSJed Brown /* Attaches the SNESDM to the coarse level. 1286cab3a1bSJed Brown * Under what conditions should we copy versus duplicate? 1296cab3a1bSJed Brown */ 1306cab3a1bSJed Brown static PetscErrorCode DMCoarsenHook_SNESDM(DM dm,DM dmc,void *ctx) 1316cab3a1bSJed Brown { 1326cab3a1bSJed Brown PetscErrorCode ierr; 1336cab3a1bSJed Brown 1346cab3a1bSJed Brown PetscFunctionBegin; 1356cab3a1bSJed Brown ierr = DMSNESCopyContext(dm,dmc);CHKERRQ(ierr); 1366cab3a1bSJed Brown PetscFunctionReturn(0); 1376cab3a1bSJed Brown } 1386cab3a1bSJed Brown 1396cab3a1bSJed Brown #undef __FUNCT__ 140caa4e7f2SJed Brown #define __FUNCT__ "DMRestrictHook_SNESDM" 141dfe15315SJed Brown /* This could restrict auxiliary information to the coarse level. 142caa4e7f2SJed Brown */ 143caa4e7f2SJed Brown static PetscErrorCode DMRestrictHook_SNESDM(DM dm,Mat Restrict,Vec rscale,Mat Inject,DM dmc,void *ctx) 144caa4e7f2SJed Brown { 145caa4e7f2SJed Brown 146caa4e7f2SJed Brown PetscFunctionBegin; 147caa4e7f2SJed Brown PetscFunctionReturn(0); 148caa4e7f2SJed Brown } 149caa4e7f2SJed Brown 150caa4e7f2SJed Brown #undef __FUNCT__ 1516cab3a1bSJed Brown #define __FUNCT__ "DMSNESGetContext" 1526cab3a1bSJed Brown /*@C 1536cab3a1bSJed Brown DMSNESGetContext - get read-only private SNESDM context from a DM 1546cab3a1bSJed Brown 1556cab3a1bSJed Brown Not Collective 1566cab3a1bSJed Brown 1576cab3a1bSJed Brown Input Argument: 1586cab3a1bSJed Brown . dm - DM to be used with SNES 1596cab3a1bSJed Brown 1606cab3a1bSJed Brown Output Argument: 1616cab3a1bSJed Brown . snesdm - private SNESDM context 1626cab3a1bSJed Brown 1636cab3a1bSJed Brown Level: developer 1646cab3a1bSJed Brown 1656cab3a1bSJed Brown Notes: 1666cab3a1bSJed Brown Use DMSNESGetContextWrite() if write access is needed. The DMSNESSetXXX API should be used wherever possible. 1676cab3a1bSJed Brown 1686cab3a1bSJed Brown .seealso: DMSNESGetContextWrite() 1696cab3a1bSJed Brown @*/ 1706cab3a1bSJed Brown PetscErrorCode DMSNESGetContext(DM dm,SNESDM *snesdm) 1716cab3a1bSJed Brown { 1726cab3a1bSJed Brown PetscErrorCode ierr; 1736cab3a1bSJed Brown PetscContainer container; 1746cab3a1bSJed Brown SNESDM sdm; 1756cab3a1bSJed Brown 1766cab3a1bSJed Brown PetscFunctionBegin; 1776cab3a1bSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1786cab3a1bSJed Brown ierr = PetscObjectQuery((PetscObject)dm,"SNESDM",(PetscObject*)&container);CHKERRQ(ierr); 1796cab3a1bSJed Brown if (container) { 1806cab3a1bSJed Brown ierr = PetscContainerGetPointer(container,(void**)snesdm);CHKERRQ(ierr); 1816cab3a1bSJed Brown } else { 1826cab3a1bSJed Brown ierr = PetscInfo(dm,"Creating new SNESDM\n");CHKERRQ(ierr); 1836cab3a1bSJed Brown ierr = PetscContainerCreate(((PetscObject)dm)->comm,&container);CHKERRQ(ierr); 1846cab3a1bSJed Brown ierr = PetscNewLog(dm,struct _n_SNESDM,&sdm);CHKERRQ(ierr); 1856cab3a1bSJed Brown ierr = PetscContainerSetPointer(container,sdm);CHKERRQ(ierr); 1866cab3a1bSJed Brown ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_SNESDM);CHKERRQ(ierr); 1876cab3a1bSJed Brown ierr = PetscObjectCompose((PetscObject)dm,"SNESDM",(PetscObject)container);CHKERRQ(ierr); 188caa4e7f2SJed Brown ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_SNESDM,DMRestrictHook_SNESDM,PETSC_NULL);CHKERRQ(ierr); 1896cab3a1bSJed Brown ierr = PetscContainerGetPointer(container,(void**)snesdm);CHKERRQ(ierr); 1906cab3a1bSJed Brown ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 1916cab3a1bSJed Brown } 1926cab3a1bSJed Brown PetscFunctionReturn(0); 1936cab3a1bSJed Brown } 1946cab3a1bSJed Brown 1956cab3a1bSJed Brown #undef __FUNCT__ 1966cab3a1bSJed Brown #define __FUNCT__ "DMSNESGetContextWrite" 1976cab3a1bSJed Brown /*@C 1986cab3a1bSJed Brown DMSNESGetContextWrite - get write access to private SNESDM context from a DM 1996cab3a1bSJed Brown 2006cab3a1bSJed Brown Not Collective 2016cab3a1bSJed Brown 2026cab3a1bSJed Brown Input Argument: 2036cab3a1bSJed Brown . dm - DM to be used with SNES 2046cab3a1bSJed Brown 2056cab3a1bSJed Brown Output Argument: 2066cab3a1bSJed Brown . snesdm - private SNESDM context 2076cab3a1bSJed Brown 2086cab3a1bSJed Brown Level: developer 2096cab3a1bSJed Brown 2106cab3a1bSJed Brown .seealso: DMSNESGetContext() 2116cab3a1bSJed Brown @*/ 2126cab3a1bSJed Brown PetscErrorCode DMSNESGetContextWrite(DM dm,SNESDM *snesdm) 2136cab3a1bSJed Brown { 2146cab3a1bSJed Brown PetscErrorCode ierr; 2156cab3a1bSJed Brown SNESDM sdm; 2166cab3a1bSJed Brown 2176cab3a1bSJed Brown PetscFunctionBegin; 2186cab3a1bSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2196cab3a1bSJed Brown ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 2206cab3a1bSJed Brown if (!sdm->originaldm) sdm->originaldm = dm; 2216cab3a1bSJed Brown if (sdm->originaldm != dm) { /* Copy on write */ 2226cab3a1bSJed Brown PetscContainer container; 2236cab3a1bSJed Brown SNESDM oldsdm = sdm; 2246cab3a1bSJed Brown ierr = PetscInfo(dm,"Copying SNESDM due to write\n");CHKERRQ(ierr); 2256cab3a1bSJed Brown ierr = PetscContainerCreate(((PetscObject)dm)->comm,&container);CHKERRQ(ierr); 2266cab3a1bSJed Brown ierr = PetscNewLog(dm,struct _n_SNESDM,&sdm);CHKERRQ(ierr); 2276cab3a1bSJed Brown ierr = PetscMemcpy(sdm,oldsdm,sizeof *sdm);CHKERRQ(ierr); 2286cab3a1bSJed Brown ierr = PetscContainerSetPointer(container,sdm);CHKERRQ(ierr); 2296cab3a1bSJed Brown ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_SNESDM);CHKERRQ(ierr); 2306cab3a1bSJed Brown ierr = PetscObjectCompose((PetscObject)dm,"SNESDM",(PetscObject)container);CHKERRQ(ierr); 2316cab3a1bSJed Brown ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 2326cab3a1bSJed Brown } 2336cab3a1bSJed Brown *snesdm = sdm; 2346cab3a1bSJed Brown PetscFunctionReturn(0); 2356cab3a1bSJed Brown } 2366cab3a1bSJed Brown 2376cab3a1bSJed Brown #undef __FUNCT__ 2386cab3a1bSJed Brown #define __FUNCT__ "DMSNESCopyContext" 2396cab3a1bSJed Brown /*@C 2406cab3a1bSJed Brown DMSNESCopyContext - copies a DM context to a new DM 2416cab3a1bSJed Brown 2426cab3a1bSJed Brown Logically Collective 2436cab3a1bSJed Brown 2446cab3a1bSJed Brown Input Arguments: 2456cab3a1bSJed Brown + dmsrc - DM to obtain context from 2466cab3a1bSJed Brown - dmdest - DM to add context to 2476cab3a1bSJed Brown 2486cab3a1bSJed Brown Level: developer 2496cab3a1bSJed Brown 2506cab3a1bSJed Brown Note: 2516cab3a1bSJed Brown The context is copied by reference. This function does not ensure that a context exists. 2526cab3a1bSJed Brown 2536cab3a1bSJed Brown .seealso: DMSNESGetContext(), SNESSetDM() 2546cab3a1bSJed Brown @*/ 2556cab3a1bSJed Brown PetscErrorCode DMSNESCopyContext(DM dmsrc,DM dmdest) 2566cab3a1bSJed Brown { 2576cab3a1bSJed Brown PetscErrorCode ierr; 2586cab3a1bSJed Brown PetscContainer container; 2596cab3a1bSJed Brown 2606cab3a1bSJed Brown PetscFunctionBegin; 2616cab3a1bSJed Brown PetscValidHeaderSpecific(dmsrc,DM_CLASSID,1); 2626cab3a1bSJed Brown PetscValidHeaderSpecific(dmdest,DM_CLASSID,2); 2636cab3a1bSJed Brown ierr = PetscObjectQuery((PetscObject)dmsrc,"SNESDM",(PetscObject*)&container);CHKERRQ(ierr); 2646cab3a1bSJed Brown if (container) { 2656cab3a1bSJed Brown ierr = PetscObjectCompose((PetscObject)dmdest,"SNESDM",(PetscObject)container);CHKERRQ(ierr); 266caa4e7f2SJed Brown ierr = DMCoarsenHookAdd(dmdest,DMCoarsenHook_SNESDM,DMRestrictHook_SNESDM,PETSC_NULL);CHKERRQ(ierr); 2676cab3a1bSJed Brown } 2686cab3a1bSJed Brown PetscFunctionReturn(0); 2696cab3a1bSJed Brown } 2706cab3a1bSJed Brown 2716cab3a1bSJed Brown #undef __FUNCT__ 2726cab3a1bSJed Brown #define __FUNCT__ "DMSNESSetFunction" 2736cab3a1bSJed Brown /*@C 2746cab3a1bSJed Brown DMSNESSetFunction - set SNES residual evaluation function 2756cab3a1bSJed Brown 2766cab3a1bSJed Brown Not Collective 2776cab3a1bSJed Brown 2786cab3a1bSJed Brown Input Arguments: 2796cab3a1bSJed Brown + dm - DM to be used with SNES 2806cab3a1bSJed Brown . func - residual evaluation function, see SNESSetFunction() for calling sequence 2816cab3a1bSJed Brown - ctx - context for residual evaluation 2826cab3a1bSJed Brown 2836cab3a1bSJed Brown Level: advanced 2846cab3a1bSJed Brown 2856cab3a1bSJed Brown Note: 2866cab3a1bSJed Brown SNESSetFunction() is normally used, but it calls this function internally because the user context is actually 2876cab3a1bSJed Brown associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 2886cab3a1bSJed Brown not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 2896cab3a1bSJed Brown 2906cab3a1bSJed Brown .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian() 2916cab3a1bSJed Brown @*/ 2926cab3a1bSJed Brown PetscErrorCode DMSNESSetFunction(DM dm,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx) 2936cab3a1bSJed Brown { 2946cab3a1bSJed Brown PetscErrorCode ierr; 2956cab3a1bSJed Brown SNESDM sdm; 2966cab3a1bSJed Brown 2976cab3a1bSJed Brown PetscFunctionBegin; 2986cab3a1bSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2996cab3a1bSJed Brown ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr); 3006cab3a1bSJed Brown if (func) sdm->computefunction = func; 3016cab3a1bSJed Brown if (ctx) sdm->functionctx = ctx; 3026cab3a1bSJed Brown PetscFunctionReturn(0); 3036cab3a1bSJed Brown } 3046cab3a1bSJed Brown 3056cab3a1bSJed Brown #undef __FUNCT__ 3066cab3a1bSJed Brown #define __FUNCT__ "DMSNESGetFunction" 3076cab3a1bSJed Brown /*@C 3086cab3a1bSJed Brown DMSNESGetFunction - get SNES residual evaluation function 3096cab3a1bSJed Brown 3106cab3a1bSJed Brown Not Collective 3116cab3a1bSJed Brown 3126cab3a1bSJed Brown Input Argument: 3136cab3a1bSJed Brown . dm - DM to be used with SNES 3146cab3a1bSJed Brown 3156cab3a1bSJed Brown Output Arguments: 3166cab3a1bSJed Brown + func - residual evaluation function, see SNESSetFunction() for calling sequence 3176cab3a1bSJed Brown - ctx - context for residual evaluation 3186cab3a1bSJed Brown 3196cab3a1bSJed Brown Level: advanced 3206cab3a1bSJed Brown 3216cab3a1bSJed Brown Note: 3226cab3a1bSJed Brown SNESGetFunction() is normally used, but it calls this function internally because the user context is actually 3236cab3a1bSJed Brown associated with the DM. 3246cab3a1bSJed Brown 3256cab3a1bSJed Brown .seealso: DMSNESSetContext(), DMSNESSetFunction(), SNESSetFunction() 3266cab3a1bSJed Brown @*/ 3276cab3a1bSJed Brown PetscErrorCode DMSNESGetFunction(DM dm,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx) 3286cab3a1bSJed Brown { 3296cab3a1bSJed Brown PetscErrorCode ierr; 3306cab3a1bSJed Brown SNESDM sdm; 3316cab3a1bSJed Brown 3326cab3a1bSJed Brown PetscFunctionBegin; 3336cab3a1bSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3346cab3a1bSJed Brown ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 3356cab3a1bSJed Brown if (func) *func = sdm->computefunction; 3366cab3a1bSJed Brown if (ctx) *ctx = sdm->functionctx; 3376cab3a1bSJed Brown PetscFunctionReturn(0); 3386cab3a1bSJed Brown } 3396cab3a1bSJed Brown 3406cab3a1bSJed Brown #undef __FUNCT__ 3416cab3a1bSJed Brown #define __FUNCT__ "DMSNESSetGS" 3426cab3a1bSJed Brown /*@C 3436cab3a1bSJed Brown DMSNESSetGS - set SNES Gauss-Seidel relaxation function 3446cab3a1bSJed Brown 3456cab3a1bSJed Brown Not Collective 3466cab3a1bSJed Brown 3476cab3a1bSJed Brown Input Argument: 3486cab3a1bSJed Brown + dm - DM to be used with SNES 3496cab3a1bSJed Brown . func - relaxation function, see SNESSetGS() for calling sequence 3506cab3a1bSJed Brown - ctx - context for residual evaluation 3516cab3a1bSJed Brown 3526cab3a1bSJed Brown Level: advanced 3536cab3a1bSJed Brown 3546cab3a1bSJed Brown Note: 3556cab3a1bSJed Brown SNESSetGS() is normally used, but it calls this function internally because the user context is actually 3566cab3a1bSJed Brown associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 3576cab3a1bSJed Brown not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 3586cab3a1bSJed Brown 3596cab3a1bSJed Brown .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), DMSNESSetFunction() 3606cab3a1bSJed Brown @*/ 3616cab3a1bSJed Brown PetscErrorCode DMSNESSetGS(DM dm,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx) 3626cab3a1bSJed Brown { 3636cab3a1bSJed Brown PetscErrorCode ierr; 3646cab3a1bSJed Brown SNESDM sdm; 3656cab3a1bSJed Brown 3666cab3a1bSJed Brown PetscFunctionBegin; 3676cab3a1bSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3686cab3a1bSJed Brown ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr); 3696cab3a1bSJed Brown if (func) sdm->computegs = func; 3706cab3a1bSJed Brown if (ctx) sdm->gsctx = ctx; 3716cab3a1bSJed Brown PetscFunctionReturn(0); 3726cab3a1bSJed Brown } 3736cab3a1bSJed Brown 3746cab3a1bSJed Brown #undef __FUNCT__ 3756cab3a1bSJed Brown #define __FUNCT__ "DMSNESGetGS" 3766cab3a1bSJed Brown /*@C 3776cab3a1bSJed Brown DMSNESGetGS - get SNES Gauss-Seidel relaxation function 3786cab3a1bSJed Brown 3796cab3a1bSJed Brown Not Collective 3806cab3a1bSJed Brown 3816cab3a1bSJed Brown Input Argument: 3826cab3a1bSJed Brown . dm - DM to be used with SNES 3836cab3a1bSJed Brown 3846cab3a1bSJed Brown Output Arguments: 3856cab3a1bSJed Brown + func - relaxation function, see SNESSetGS() for calling sequence 3866cab3a1bSJed Brown - ctx - context for residual evaluation 3876cab3a1bSJed Brown 3886cab3a1bSJed Brown Level: advanced 3896cab3a1bSJed Brown 3906cab3a1bSJed Brown Note: 3916cab3a1bSJed Brown SNESGetGS() is normally used, but it calls this function internally because the user context is actually 3926cab3a1bSJed Brown associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 3936cab3a1bSJed Brown not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 3946cab3a1bSJed Brown 3956cab3a1bSJed Brown .seealso: DMSNESSetContext(), SNESGetGS(), DMSNESGetJacobian(), DMSNESGetFunction() 3966cab3a1bSJed Brown @*/ 3976cab3a1bSJed Brown PetscErrorCode DMSNESGetGS(DM dm,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx) 3986cab3a1bSJed Brown { 3996cab3a1bSJed Brown PetscErrorCode ierr; 4006cab3a1bSJed Brown SNESDM sdm; 4016cab3a1bSJed Brown 4026cab3a1bSJed Brown PetscFunctionBegin; 4036cab3a1bSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 4046cab3a1bSJed Brown ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 4056cab3a1bSJed Brown if (func) *func = sdm->computegs; 4066cab3a1bSJed Brown if (ctx) *ctx = sdm->gsctx; 4076cab3a1bSJed Brown PetscFunctionReturn(0); 4086cab3a1bSJed Brown } 4096cab3a1bSJed Brown 4106cab3a1bSJed Brown #undef __FUNCT__ 4116cab3a1bSJed Brown #define __FUNCT__ "DMSNESSetJacobian" 4126cab3a1bSJed Brown /*@C 4136cab3a1bSJed Brown DMSNESSetFunction - set SNES Jacobian evaluation function 4146cab3a1bSJed Brown 4156cab3a1bSJed Brown Not Collective 4166cab3a1bSJed Brown 4176cab3a1bSJed Brown Input Argument: 4186cab3a1bSJed Brown + dm - DM to be used with SNES 4196cab3a1bSJed Brown . func - Jacobian evaluation function, see SNESSetJacobian() for calling sequence 4206cab3a1bSJed Brown - ctx - context for residual evaluation 4216cab3a1bSJed Brown 4226cab3a1bSJed Brown Level: advanced 4236cab3a1bSJed Brown 4246cab3a1bSJed Brown Note: 4256cab3a1bSJed Brown SNESSetJacobian() is normally used, but it calls this function internally because the user context is actually 4266cab3a1bSJed Brown associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 4276cab3a1bSJed Brown not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 4286cab3a1bSJed Brown 4296cab3a1bSJed Brown .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESGetJacobian(), SNESSetJacobian() 4306cab3a1bSJed Brown @*/ 4316cab3a1bSJed Brown PetscErrorCode DMSNESSetJacobian(DM dm,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 4326cab3a1bSJed Brown { 4336cab3a1bSJed Brown PetscErrorCode ierr; 4346cab3a1bSJed Brown SNESDM sdm; 4356cab3a1bSJed Brown 4366cab3a1bSJed Brown PetscFunctionBegin; 4376cab3a1bSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 4386cab3a1bSJed Brown ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr); 4396cab3a1bSJed Brown if (func) sdm->computejacobian = func; 4406cab3a1bSJed Brown if (ctx) sdm->jacobianctx = ctx; 4416cab3a1bSJed Brown PetscFunctionReturn(0); 4426cab3a1bSJed Brown } 4436cab3a1bSJed Brown 4446cab3a1bSJed Brown #undef __FUNCT__ 4456cab3a1bSJed Brown #define __FUNCT__ "DMSNESGetJacobian" 4466cab3a1bSJed Brown /*@C 4476cab3a1bSJed Brown DMSNESGetFunction - get SNES Jacobian evaluation function 4486cab3a1bSJed Brown 4496cab3a1bSJed Brown Not Collective 4506cab3a1bSJed Brown 4516cab3a1bSJed Brown Input Argument: 4526cab3a1bSJed Brown . dm - DM to be used with SNES 4536cab3a1bSJed Brown 4546cab3a1bSJed Brown Output Arguments: 4556cab3a1bSJed Brown + func - Jacobian evaluation function, see SNESSetJacobian() for calling sequence 4566cab3a1bSJed Brown - ctx - context for residual evaluation 4576cab3a1bSJed Brown 4586cab3a1bSJed Brown Level: advanced 4596cab3a1bSJed Brown 4606cab3a1bSJed Brown Note: 4616cab3a1bSJed Brown SNESGetJacobian() is normally used, but it calls this function internally because the user context is actually 4626cab3a1bSJed Brown associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 4636cab3a1bSJed Brown not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 4646cab3a1bSJed Brown 4656cab3a1bSJed Brown .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian() 4666cab3a1bSJed Brown @*/ 4676cab3a1bSJed Brown PetscErrorCode DMSNESGetJacobian(DM dm,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx) 4686cab3a1bSJed Brown { 4696cab3a1bSJed Brown PetscErrorCode ierr; 4706cab3a1bSJed Brown SNESDM sdm; 4716cab3a1bSJed Brown 4726cab3a1bSJed Brown PetscFunctionBegin; 4736cab3a1bSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 4746cab3a1bSJed Brown ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 4756cab3a1bSJed Brown if (func) *func = sdm->computejacobian; 4766cab3a1bSJed Brown if (ctx) *ctx = sdm->jacobianctx; 4776cab3a1bSJed Brown PetscFunctionReturn(0); 4786cab3a1bSJed Brown } 4796cab3a1bSJed Brown 4806cab3a1bSJed Brown #undef __FUNCT__ 4816cab3a1bSJed Brown #define __FUNCT__ "SNESDefaultComputeFunction_DMLegacy" 4826cab3a1bSJed Brown static PetscErrorCode SNESDefaultComputeFunction_DMLegacy(SNES snes,Vec X,Vec F,void *ctx) 4836cab3a1bSJed Brown { 4846cab3a1bSJed Brown PetscErrorCode ierr; 4856cab3a1bSJed Brown DM dm; 4866cab3a1bSJed Brown 4876cab3a1bSJed Brown PetscFunctionBegin; 4886cab3a1bSJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 4896cab3a1bSJed Brown ierr = DMComputeFunction(dm,X,F);CHKERRQ(ierr); 4906cab3a1bSJed Brown PetscFunctionReturn(0); 4916cab3a1bSJed Brown } 4926cab3a1bSJed Brown 4936cab3a1bSJed Brown #undef __FUNCT__ 4946cab3a1bSJed Brown #define __FUNCT__ "SNESDefaultComputeJacobian_DMLegacy" 4956cab3a1bSJed Brown static PetscErrorCode SNESDefaultComputeJacobian_DMLegacy(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *ctx) 4966cab3a1bSJed Brown { 4976cab3a1bSJed Brown PetscErrorCode ierr; 4986cab3a1bSJed Brown DM dm; 4996cab3a1bSJed Brown 5006cab3a1bSJed Brown PetscFunctionBegin; 5016cab3a1bSJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 5026cab3a1bSJed Brown ierr = DMComputeJacobian(dm,X,*A,*B,mstr);CHKERRQ(ierr); 5036cab3a1bSJed Brown PetscFunctionReturn(0); 5046cab3a1bSJed Brown } 5056cab3a1bSJed Brown 5066cab3a1bSJed Brown #undef __FUNCT__ 5076cab3a1bSJed Brown #define __FUNCT__ "DMSNESSetUpLegacy" 5086cab3a1bSJed Brown /* Sets up calling of legacy DM routines */ 5096cab3a1bSJed Brown PetscErrorCode DMSNESSetUpLegacy(DM dm) 5106cab3a1bSJed Brown { 5116cab3a1bSJed Brown PetscErrorCode ierr; 5126cab3a1bSJed Brown SNESDM sdm; 5136cab3a1bSJed Brown 5146cab3a1bSJed Brown PetscFunctionBegin; 5156cab3a1bSJed Brown ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 5166cab3a1bSJed Brown if (!sdm->computefunction) {ierr = DMSNESSetFunction(dm,SNESDefaultComputeFunction_DMLegacy,PETSC_NULL);CHKERRQ(ierr);} 5176cab3a1bSJed Brown ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 5186cab3a1bSJed Brown if (!sdm->computejacobian) {ierr = DMSNESSetJacobian(dm,SNESDefaultComputeJacobian_DMLegacy,PETSC_NULL);CHKERRQ(ierr);} 5196cab3a1bSJed Brown PetscFunctionReturn(0); 5206cab3a1bSJed Brown } 521