1b45d2f2cSJed Brown #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/ 2679f678eSPeter Brune #include <petsc-private/dmimpl.h> /*I "petscdm.h" I*/ 36cab3a1bSJed Brown 46cab3a1bSJed Brown #undef __FUNCT__ 598120457SMatthew G Knepley #define __FUNCT__ "SNESDMComputeFunction" 698120457SMatthew G Knepley /*@C 798120457SMatthew G Knepley SNESDMComputeFunction - This is a universal function evaluation routine that 898120457SMatthew G Knepley may be used with SNESSetFunction() as long as the user context has a DM 998120457SMatthew G Knepley as its first record and the user has called DMSetLocalFunction(). 1098120457SMatthew G Knepley 1198120457SMatthew G Knepley Collective on SNES 1298120457SMatthew G Knepley 1398120457SMatthew G Knepley Input Parameters: 1498120457SMatthew G Knepley + snes - the SNES context 1598120457SMatthew G Knepley . X - input vector 1698120457SMatthew G Knepley . F - function vector 1798120457SMatthew G Knepley - ptr - pointer to a structure that must have a DM as its first entry. 1898120457SMatthew G Knepley This ptr must have been passed into SNESDMComputeFunction() as the context. 1998120457SMatthew G Knepley 2098120457SMatthew G Knepley Level: intermediate 2198120457SMatthew G Knepley 2298120457SMatthew G Knepley .seealso: DMSetLocalFunction(), DMSetLocalJacobian(), SNESSetFunction(), SNESSetJacobian() 2398120457SMatthew G Knepley @*/ 2498120457SMatthew G Knepley PetscErrorCode SNESDMComputeFunction(SNES snes, Vec X, Vec F, void *ptr) 2598120457SMatthew G Knepley { 2698120457SMatthew G Knepley DM dm = *(DM*) ptr; 2798120457SMatthew G Knepley PetscErrorCode (*lf)(DM, Vec, Vec, void *); 2898120457SMatthew G Knepley Vec localX, localF; 2998120457SMatthew G Knepley PetscInt N, n; 3098120457SMatthew G Knepley PetscErrorCode ierr; 3198120457SMatthew G Knepley 3298120457SMatthew G Knepley PetscFunctionBegin; 3398120457SMatthew G Knepley PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 3498120457SMatthew G Knepley PetscValidHeaderSpecific(X, VEC_CLASSID, 2); 3598120457SMatthew G Knepley PetscValidHeaderSpecific(F, VEC_CLASSID, 3); 3698120457SMatthew G Knepley if (!dm) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONGSTATE, "Looks like you called SNESSetFromFuntion(snes,SNESDMComputeFunction,) without the DM context"); 3798120457SMatthew G Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 4); 3898120457SMatthew G Knepley 3998120457SMatthew G Knepley /* determine whether X = localX */ 4098120457SMatthew G Knepley ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 4198120457SMatthew G Knepley ierr = DMGetLocalVector(dm, &localF);CHKERRQ(ierr); 4298120457SMatthew G Knepley ierr = VecGetSize(X, &N);CHKERRQ(ierr); 4398120457SMatthew G Knepley ierr = VecGetSize(localX, &n);CHKERRQ(ierr); 4498120457SMatthew G Knepley 4598120457SMatthew G Knepley if (n != N){ /* X != localX */ 4698120457SMatthew G Knepley /* Scatter ghost points to local vector, using the 2-step process 4798120457SMatthew G Knepley DMGlobalToLocalBegin(), DMGlobalToLocalEnd(). 4898120457SMatthew G Knepley */ 4998120457SMatthew G Knepley ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 5098120457SMatthew G Knepley ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 5198120457SMatthew G Knepley } else { 5298120457SMatthew G Knepley ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 5398120457SMatthew G Knepley localX = X; 5498120457SMatthew G Knepley } 5598120457SMatthew G Knepley ierr = DMGetLocalFunction(dm, &lf);CHKERRQ(ierr); 5698120457SMatthew G Knepley ierr = (*lf)(dm, localX, localF, ptr);CHKERRQ(ierr); 5798120457SMatthew G Knepley if (n != N){ 5898120457SMatthew G Knepley ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 5998120457SMatthew G Knepley } 6098120457SMatthew G Knepley ierr = DMLocalToGlobalBegin(dm, localF, ADD_VALUES, F);CHKERRQ(ierr); 6198120457SMatthew G Knepley ierr = DMLocalToGlobalEnd(dm, localF, ADD_VALUES, F);CHKERRQ(ierr); 6298120457SMatthew G Knepley ierr = DMRestoreLocalVector(dm, &localF);CHKERRQ(ierr); 6398120457SMatthew G Knepley PetscFunctionReturn(0); 6498120457SMatthew G Knepley } 6598120457SMatthew G Knepley 6698120457SMatthew G Knepley #undef __FUNCT__ 6798120457SMatthew G Knepley #define __FUNCT__ "SNESDMComputeJacobian" 6898120457SMatthew G Knepley /* 6998120457SMatthew G Knepley SNESDMComputeJacobian - This is a universal Jacobian evaluation routine that 7098120457SMatthew G Knepley may be used with SNESSetJacobian() as long as the user context has a DM 7198120457SMatthew G Knepley as its first record and the user has called DMSetLocalJacobian(). 7298120457SMatthew G Knepley 7398120457SMatthew G Knepley Collective on SNES 7498120457SMatthew G Knepley 7598120457SMatthew G Knepley Input Parameters: 7698120457SMatthew G Knepley + snes - the SNES context 7798120457SMatthew G Knepley . X - input vector 7898120457SMatthew G Knepley . J - Jacobian 7998120457SMatthew G Knepley . B - Jacobian used in preconditioner (usally same as J) 8098120457SMatthew G Knepley . flag - indicates if the matrix changed its structure 8198120457SMatthew G Knepley - ptr - pointer to a structure that must have a DM as its first entry. 8298120457SMatthew G Knepley This ptr must have been passed into SNESDMComputeFunction() as the context. 8398120457SMatthew G Knepley 8498120457SMatthew G Knepley Level: intermediate 8598120457SMatthew G Knepley 8698120457SMatthew G Knepley .seealso: DMSetLocalFunction(), DMSetLocalJacobian(), SNESSetFunction(), SNESSetJacobian() 8798120457SMatthew G Knepley */ 8898120457SMatthew G Knepley PetscErrorCode SNESDMComputeJacobian(SNES snes, Vec X, Mat *J, Mat *B, MatStructure *flag, void *ptr) 8998120457SMatthew G Knepley { 9098120457SMatthew G Knepley DM dm = *(DM*) ptr; 9198120457SMatthew G Knepley PetscErrorCode (*lj)(DM, Vec, Mat, Mat, void *); 9298120457SMatthew G Knepley Vec localX; 9398120457SMatthew G Knepley PetscErrorCode ierr; 9498120457SMatthew G Knepley 9598120457SMatthew G Knepley PetscFunctionBegin; 9698120457SMatthew G Knepley ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 9798120457SMatthew G Knepley ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 9898120457SMatthew G Knepley ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 9998120457SMatthew G Knepley ierr = DMGetLocalJacobian(dm, &lj);CHKERRQ(ierr); 10098120457SMatthew G Knepley ierr = (*lj)(dm, localX, *J, *B, ptr);CHKERRQ(ierr); 10198120457SMatthew G Knepley ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 10298120457SMatthew G Knepley /* Assemble true Jacobian; if it is different */ 10398120457SMatthew G Knepley if (*J != *B) { 10498120457SMatthew G Knepley ierr = MatAssemblyBegin(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 10598120457SMatthew G Knepley ierr = MatAssemblyEnd(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 10698120457SMatthew G Knepley } 10798120457SMatthew G Knepley ierr = MatSetOption(*B, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE);CHKERRQ(ierr); 10898120457SMatthew G Knepley *flag = SAME_NONZERO_PATTERN; 10998120457SMatthew G Knepley PetscFunctionReturn(0); 11098120457SMatthew G Knepley } 11198120457SMatthew G Knepley 11298120457SMatthew 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__ 15103a0fabfSPeter Brune #define __FUNCT__ "DMRefineHook_SNESDM" 15203a0fabfSPeter Brune static PetscErrorCode DMRefineHook_SNESDM(DM dm,DM dmf,void *ctx) 15303a0fabfSPeter Brune { 15403a0fabfSPeter Brune PetscErrorCode ierr; 15503a0fabfSPeter Brune 15603a0fabfSPeter Brune PetscFunctionBegin; 15703a0fabfSPeter Brune ierr = DMSNESCopyContext(dm,dmf);CHKERRQ(ierr); 15803a0fabfSPeter Brune PetscFunctionReturn(0); 15903a0fabfSPeter Brune } 16003a0fabfSPeter Brune 16103a0fabfSPeter Brune #undef __FUNCT__ 16203a0fabfSPeter Brune #define __FUNCT__ "DMInterpolateHook_SNESDM" 16303a0fabfSPeter Brune /* This could restrict auxiliary information to the coarse level. 16403a0fabfSPeter Brune */ 16503a0fabfSPeter Brune static PetscErrorCode DMInterpolateHook_SNESDM(DM dm,Mat Interp,DM dmf,void *ctx) 16603a0fabfSPeter Brune { 16703a0fabfSPeter Brune 16803a0fabfSPeter Brune PetscFunctionBegin; 16903a0fabfSPeter Brune PetscFunctionReturn(0); 17003a0fabfSPeter Brune } 17103a0fabfSPeter Brune 17203a0fabfSPeter Brune #undef __FUNCT__ 1736cab3a1bSJed Brown #define __FUNCT__ "DMSNESGetContext" 1746cab3a1bSJed Brown /*@C 1756cab3a1bSJed Brown DMSNESGetContext - get read-only private SNESDM context from a DM 1766cab3a1bSJed Brown 1776cab3a1bSJed Brown Not Collective 1786cab3a1bSJed Brown 1796cab3a1bSJed Brown Input Argument: 1806cab3a1bSJed Brown . dm - DM to be used with SNES 1816cab3a1bSJed Brown 1826cab3a1bSJed Brown Output Argument: 1836cab3a1bSJed Brown . snesdm - private SNESDM context 1846cab3a1bSJed Brown 1856cab3a1bSJed Brown Level: developer 1866cab3a1bSJed Brown 1876cab3a1bSJed Brown Notes: 1886cab3a1bSJed Brown Use DMSNESGetContextWrite() if write access is needed. The DMSNESSetXXX API should be used wherever possible. 1896cab3a1bSJed Brown 1906cab3a1bSJed Brown .seealso: DMSNESGetContextWrite() 1916cab3a1bSJed Brown @*/ 1926cab3a1bSJed Brown PetscErrorCode DMSNESGetContext(DM dm,SNESDM *snesdm) 1936cab3a1bSJed Brown { 1946cab3a1bSJed Brown PetscErrorCode ierr; 1956cab3a1bSJed Brown PetscContainer container; 1966cab3a1bSJed Brown SNESDM sdm; 1976cab3a1bSJed Brown 1986cab3a1bSJed Brown PetscFunctionBegin; 1996cab3a1bSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2006cab3a1bSJed Brown ierr = PetscObjectQuery((PetscObject)dm,"SNESDM",(PetscObject*)&container);CHKERRQ(ierr); 2016cab3a1bSJed Brown if (container) { 2026cab3a1bSJed Brown ierr = PetscContainerGetPointer(container,(void**)snesdm);CHKERRQ(ierr); 2036cab3a1bSJed Brown } else { 2046cab3a1bSJed Brown ierr = PetscInfo(dm,"Creating new SNESDM\n");CHKERRQ(ierr); 2056cab3a1bSJed Brown ierr = PetscContainerCreate(((PetscObject)dm)->comm,&container);CHKERRQ(ierr); 2066cab3a1bSJed Brown ierr = PetscNewLog(dm,struct _n_SNESDM,&sdm);CHKERRQ(ierr); 2076cab3a1bSJed Brown ierr = PetscContainerSetPointer(container,sdm);CHKERRQ(ierr); 2086cab3a1bSJed Brown ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_SNESDM);CHKERRQ(ierr); 2096cab3a1bSJed Brown ierr = PetscObjectCompose((PetscObject)dm,"SNESDM",(PetscObject)container);CHKERRQ(ierr); 210caa4e7f2SJed Brown ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_SNESDM,DMRestrictHook_SNESDM,PETSC_NULL);CHKERRQ(ierr); 2113de5108cSBarry Smith ierr = DMRefineHookAdd(dm,DMRefineHook_SNESDM,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 2126cab3a1bSJed Brown ierr = PetscContainerGetPointer(container,(void**)snesdm);CHKERRQ(ierr); 2136cab3a1bSJed Brown ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 2146cab3a1bSJed Brown } 2156cab3a1bSJed Brown PetscFunctionReturn(0); 2166cab3a1bSJed Brown } 2176cab3a1bSJed Brown 2186cab3a1bSJed Brown #undef __FUNCT__ 2196cab3a1bSJed Brown #define __FUNCT__ "DMSNESGetContextWrite" 2206cab3a1bSJed Brown /*@C 2216cab3a1bSJed Brown DMSNESGetContextWrite - get write access to private SNESDM context from a DM 2226cab3a1bSJed Brown 2236cab3a1bSJed Brown Not Collective 2246cab3a1bSJed Brown 2256cab3a1bSJed Brown Input Argument: 2266cab3a1bSJed Brown . dm - DM to be used with SNES 2276cab3a1bSJed Brown 2286cab3a1bSJed Brown Output Argument: 2296cab3a1bSJed Brown . snesdm - private SNESDM context 2306cab3a1bSJed Brown 2316cab3a1bSJed Brown Level: developer 2326cab3a1bSJed Brown 2336cab3a1bSJed Brown .seealso: DMSNESGetContext() 2346cab3a1bSJed Brown @*/ 2356cab3a1bSJed Brown PetscErrorCode DMSNESGetContextWrite(DM dm,SNESDM *snesdm) 2366cab3a1bSJed Brown { 2376cab3a1bSJed Brown PetscErrorCode ierr; 2386cab3a1bSJed Brown SNESDM sdm; 2396cab3a1bSJed Brown 2406cab3a1bSJed Brown PetscFunctionBegin; 2416cab3a1bSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2426cab3a1bSJed Brown ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 2436cab3a1bSJed Brown if (!sdm->originaldm) sdm->originaldm = dm; 2446cab3a1bSJed Brown if (sdm->originaldm != dm) { /* Copy on write */ 2456cab3a1bSJed Brown PetscContainer container; 2466cab3a1bSJed Brown SNESDM oldsdm = sdm; 2476cab3a1bSJed Brown ierr = PetscInfo(dm,"Copying SNESDM due to write\n");CHKERRQ(ierr); 2486cab3a1bSJed Brown ierr = PetscContainerCreate(((PetscObject)dm)->comm,&container);CHKERRQ(ierr); 2496cab3a1bSJed Brown ierr = PetscNewLog(dm,struct _n_SNESDM,&sdm);CHKERRQ(ierr); 2508caf3d72SBarry Smith ierr = PetscMemcpy(sdm,oldsdm,sizeof(*sdm));CHKERRQ(ierr); 2516cab3a1bSJed Brown ierr = PetscContainerSetPointer(container,sdm);CHKERRQ(ierr); 2526cab3a1bSJed Brown ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_SNESDM);CHKERRQ(ierr); 2536cab3a1bSJed Brown ierr = PetscObjectCompose((PetscObject)dm,"SNESDM",(PetscObject)container);CHKERRQ(ierr); 2546cab3a1bSJed Brown ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 255e3c0b8a2SPeter Brune /* implementation specific copy hooks */ 256*c9058dfbSJed Brown if (sdm->duplicate) {ierr = (*sdm->duplicate)(oldsdm,dm);CHKERRQ(ierr);} 2576cab3a1bSJed Brown } 2586cab3a1bSJed Brown *snesdm = sdm; 2596cab3a1bSJed Brown PetscFunctionReturn(0); 2606cab3a1bSJed Brown } 2616cab3a1bSJed Brown 2626cab3a1bSJed Brown #undef __FUNCT__ 2636cab3a1bSJed Brown #define __FUNCT__ "DMSNESCopyContext" 2646cab3a1bSJed Brown /*@C 2656cab3a1bSJed Brown DMSNESCopyContext - copies a DM context to a new DM 2666cab3a1bSJed Brown 2676cab3a1bSJed Brown Logically Collective 2686cab3a1bSJed Brown 2696cab3a1bSJed Brown Input Arguments: 2706cab3a1bSJed Brown + dmsrc - DM to obtain context from 2716cab3a1bSJed Brown - dmdest - DM to add context to 2726cab3a1bSJed Brown 2736cab3a1bSJed Brown Level: developer 2746cab3a1bSJed Brown 2756cab3a1bSJed Brown Note: 2766cab3a1bSJed Brown The context is copied by reference. This function does not ensure that a context exists. 2776cab3a1bSJed Brown 2786cab3a1bSJed Brown .seealso: DMSNESGetContext(), SNESSetDM() 2796cab3a1bSJed Brown @*/ 2806cab3a1bSJed Brown PetscErrorCode DMSNESCopyContext(DM dmsrc,DM dmdest) 2816cab3a1bSJed Brown { 2826cab3a1bSJed Brown PetscErrorCode ierr; 2836cab3a1bSJed Brown PetscContainer container; 2846cab3a1bSJed Brown 2856cab3a1bSJed Brown PetscFunctionBegin; 2866cab3a1bSJed Brown PetscValidHeaderSpecific(dmsrc,DM_CLASSID,1); 2876cab3a1bSJed Brown PetscValidHeaderSpecific(dmdest,DM_CLASSID,2); 2886cab3a1bSJed Brown ierr = PetscObjectQuery((PetscObject)dmsrc,"SNESDM",(PetscObject*)&container);CHKERRQ(ierr); 2896cab3a1bSJed Brown if (container) { 2906cab3a1bSJed Brown ierr = PetscObjectCompose((PetscObject)dmdest,"SNESDM",(PetscObject)container);CHKERRQ(ierr); 291caa4e7f2SJed Brown ierr = DMCoarsenHookAdd(dmdest,DMCoarsenHook_SNESDM,DMRestrictHook_SNESDM,PETSC_NULL);CHKERRQ(ierr); 29203a0fabfSPeter Brune ierr = DMRefineHookAdd(dmdest,DMRefineHook_SNESDM,DMInterpolateHook_SNESDM,PETSC_NULL);CHKERRQ(ierr); 2936cab3a1bSJed Brown } 2946cab3a1bSJed Brown PetscFunctionReturn(0); 2956cab3a1bSJed Brown } 2966cab3a1bSJed Brown 2976cab3a1bSJed Brown #undef __FUNCT__ 2986cab3a1bSJed Brown #define __FUNCT__ "DMSNESSetFunction" 2996cab3a1bSJed Brown /*@C 3006cab3a1bSJed Brown DMSNESSetFunction - set SNES residual evaluation function 3016cab3a1bSJed Brown 3026cab3a1bSJed Brown Not Collective 3036cab3a1bSJed Brown 3046cab3a1bSJed Brown Input Arguments: 3056cab3a1bSJed Brown + dm - DM to be used with SNES 3066cab3a1bSJed Brown . func - residual evaluation function, see SNESSetFunction() for calling sequence 3076cab3a1bSJed Brown - ctx - context for residual evaluation 3086cab3a1bSJed Brown 3096cab3a1bSJed Brown Level: advanced 3106cab3a1bSJed Brown 3116cab3a1bSJed Brown Note: 3126cab3a1bSJed Brown SNESSetFunction() is normally used, but it calls this function internally because the user context is actually 3136cab3a1bSJed Brown associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 3146cab3a1bSJed Brown not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 3156cab3a1bSJed Brown 3166cab3a1bSJed Brown .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian() 3176cab3a1bSJed Brown @*/ 3186cab3a1bSJed Brown PetscErrorCode DMSNESSetFunction(DM dm,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx) 3196cab3a1bSJed Brown { 3206cab3a1bSJed Brown PetscErrorCode ierr; 3216cab3a1bSJed Brown SNESDM sdm; 3226cab3a1bSJed Brown 3236cab3a1bSJed Brown PetscFunctionBegin; 3246cab3a1bSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 325fdaff8d6SPeter Brune if (func || ctx) { 326fdaff8d6SPeter Brune ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr); 327fdaff8d6SPeter Brune } 3286cab3a1bSJed Brown if (func) sdm->computefunction = func; 3296cab3a1bSJed Brown if (ctx) sdm->functionctx = ctx; 3306cab3a1bSJed Brown PetscFunctionReturn(0); 3316cab3a1bSJed Brown } 3326cab3a1bSJed Brown 3336cab3a1bSJed Brown #undef __FUNCT__ 3346cab3a1bSJed Brown #define __FUNCT__ "DMSNESGetFunction" 3356cab3a1bSJed Brown /*@C 3366cab3a1bSJed Brown DMSNESGetFunction - get SNES residual evaluation function 3376cab3a1bSJed Brown 3386cab3a1bSJed Brown Not Collective 3396cab3a1bSJed Brown 3406cab3a1bSJed Brown Input Argument: 3416cab3a1bSJed Brown . dm - DM to be used with SNES 3426cab3a1bSJed Brown 3436cab3a1bSJed Brown Output Arguments: 3446cab3a1bSJed Brown + func - residual evaluation function, see SNESSetFunction() for calling sequence 3456cab3a1bSJed Brown - ctx - context for residual evaluation 3466cab3a1bSJed Brown 3476cab3a1bSJed Brown Level: advanced 3486cab3a1bSJed Brown 3496cab3a1bSJed Brown Note: 3506cab3a1bSJed Brown SNESGetFunction() is normally used, but it calls this function internally because the user context is actually 3516cab3a1bSJed Brown associated with the DM. 3526cab3a1bSJed Brown 3536cab3a1bSJed Brown .seealso: DMSNESSetContext(), DMSNESSetFunction(), SNESSetFunction() 3546cab3a1bSJed Brown @*/ 3556cab3a1bSJed Brown PetscErrorCode DMSNESGetFunction(DM dm,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx) 3566cab3a1bSJed Brown { 3576cab3a1bSJed Brown PetscErrorCode ierr; 3586cab3a1bSJed Brown SNESDM sdm; 3596cab3a1bSJed Brown 3606cab3a1bSJed Brown PetscFunctionBegin; 3616cab3a1bSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3626cab3a1bSJed Brown ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 3636cab3a1bSJed Brown if (func) *func = sdm->computefunction; 3646cab3a1bSJed Brown if (ctx) *ctx = sdm->functionctx; 3656cab3a1bSJed Brown PetscFunctionReturn(0); 3666cab3a1bSJed Brown } 3676cab3a1bSJed Brown 3686cab3a1bSJed Brown #undef __FUNCT__ 3692a4ee8f2SPeter Brune #define __FUNCT__ "DMSNESSetObjective" 3702a4ee8f2SPeter Brune /*@C 371081a7dcdSPeter Brune DMSNESSetObjective - set SNES objective evaluation function 3722a4ee8f2SPeter Brune 3732a4ee8f2SPeter Brune Not Collective 3742a4ee8f2SPeter Brune 3752a4ee8f2SPeter Brune Input Arguments: 3762a4ee8f2SPeter Brune + dm - DM to be used with SNES 3772a4ee8f2SPeter Brune . func - residual evaluation function, see SNESSetObjective() for calling sequence 3782a4ee8f2SPeter Brune - ctx - context for residual evaluation 3792a4ee8f2SPeter Brune 3802a4ee8f2SPeter Brune Level: advanced 3812a4ee8f2SPeter Brune 3822a4ee8f2SPeter Brune .seealso: DMSNESSetContext(), SNESGetObjective(), DMSNESSetFunction() 3832a4ee8f2SPeter Brune @*/ 3842a4ee8f2SPeter Brune PetscErrorCode DMSNESSetObjective(DM dm,SNESObjective func,void *ctx) 3852a4ee8f2SPeter Brune { 3862a4ee8f2SPeter Brune PetscErrorCode ierr; 3872a4ee8f2SPeter Brune SNESDM sdm; 3882a4ee8f2SPeter Brune 3892a4ee8f2SPeter Brune PetscFunctionBegin; 3902a4ee8f2SPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 391fdaff8d6SPeter Brune if (func || ctx) { 392fdaff8d6SPeter Brune ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr); 393fdaff8d6SPeter Brune } 3942a4ee8f2SPeter Brune if (func) sdm->computeobjective = func; 3952a4ee8f2SPeter Brune if (ctx) sdm->objectivectx = ctx; 3962a4ee8f2SPeter Brune PetscFunctionReturn(0); 3972a4ee8f2SPeter Brune } 3982a4ee8f2SPeter Brune 3992a4ee8f2SPeter Brune #undef __FUNCT__ 4002a4ee8f2SPeter Brune #define __FUNCT__ "DMSNESGetObjective" 4012a4ee8f2SPeter Brune /*@C 4022a4ee8f2SPeter Brune DMSNESGetObjective - get SNES objective evaluation function 4032a4ee8f2SPeter Brune 4042a4ee8f2SPeter Brune Not Collective 4052a4ee8f2SPeter Brune 4062a4ee8f2SPeter Brune Input Argument: 4072a4ee8f2SPeter Brune . dm - DM to be used with SNES 4082a4ee8f2SPeter Brune 4092a4ee8f2SPeter Brune Output Arguments: 4102a4ee8f2SPeter Brune + func - residual evaluation function, see SNESSetObjective() for calling sequence 4112a4ee8f2SPeter Brune - ctx - context for residual evaluation 4122a4ee8f2SPeter Brune 4132a4ee8f2SPeter Brune Level: advanced 4142a4ee8f2SPeter Brune 4152a4ee8f2SPeter Brune Note: 4162a4ee8f2SPeter Brune SNESGetFunction() is normally used, but it calls this function internally because the user context is actually 4172a4ee8f2SPeter Brune associated with the DM. 4182a4ee8f2SPeter Brune 4192a4ee8f2SPeter Brune .seealso: DMSNESSetContext(), DMSNESSetObjective(), SNESSetFunction() 4202a4ee8f2SPeter Brune @*/ 4212a4ee8f2SPeter Brune PetscErrorCode DMSNESGetObjective(DM dm,SNESObjective *func,void **ctx) 4222a4ee8f2SPeter Brune { 4232a4ee8f2SPeter Brune PetscErrorCode ierr; 4242a4ee8f2SPeter Brune SNESDM sdm; 4252a4ee8f2SPeter Brune 4262a4ee8f2SPeter Brune PetscFunctionBegin; 4272a4ee8f2SPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 4282a4ee8f2SPeter Brune ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 4292a4ee8f2SPeter Brune if (func) *func = sdm->computeobjective; 4302a4ee8f2SPeter Brune if (ctx) *ctx = sdm->objectivectx; 4312a4ee8f2SPeter Brune PetscFunctionReturn(0); 4322a4ee8f2SPeter Brune } 4332a4ee8f2SPeter Brune 4342a4ee8f2SPeter Brune #undef __FUNCT__ 4356cab3a1bSJed Brown #define __FUNCT__ "DMSNESSetGS" 4366cab3a1bSJed Brown /*@C 4376cab3a1bSJed Brown DMSNESSetGS - set SNES Gauss-Seidel relaxation function 4386cab3a1bSJed Brown 4396cab3a1bSJed Brown Not Collective 4406cab3a1bSJed Brown 4416cab3a1bSJed Brown Input Argument: 4426cab3a1bSJed Brown + dm - DM to be used with SNES 4436cab3a1bSJed Brown . func - relaxation function, see SNESSetGS() for calling sequence 4446cab3a1bSJed Brown - ctx - context for residual evaluation 4456cab3a1bSJed Brown 4466cab3a1bSJed Brown Level: advanced 4476cab3a1bSJed Brown 4486cab3a1bSJed Brown Note: 4496cab3a1bSJed Brown SNESSetGS() is normally used, but it calls this function internally because the user context is actually 4506cab3a1bSJed Brown associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 4516cab3a1bSJed Brown not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 4526cab3a1bSJed Brown 4536cab3a1bSJed Brown .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), DMSNESSetFunction() 4546cab3a1bSJed Brown @*/ 4556cab3a1bSJed Brown PetscErrorCode DMSNESSetGS(DM dm,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx) 4566cab3a1bSJed Brown { 4576cab3a1bSJed Brown PetscErrorCode ierr; 4586cab3a1bSJed Brown SNESDM sdm; 4596cab3a1bSJed Brown 4606cab3a1bSJed Brown PetscFunctionBegin; 4616cab3a1bSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 462fdaff8d6SPeter Brune if (func || ctx) { 463fdaff8d6SPeter Brune ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr); 464fdaff8d6SPeter Brune } 4656cab3a1bSJed Brown if (func) sdm->computegs = func; 4666cab3a1bSJed Brown if (ctx) sdm->gsctx = ctx; 4676cab3a1bSJed Brown PetscFunctionReturn(0); 4686cab3a1bSJed Brown } 4696cab3a1bSJed Brown 4706cab3a1bSJed Brown #undef __FUNCT__ 4716cab3a1bSJed Brown #define __FUNCT__ "DMSNESGetGS" 4726cab3a1bSJed Brown /*@C 4736cab3a1bSJed Brown DMSNESGetGS - get SNES Gauss-Seidel relaxation function 4746cab3a1bSJed Brown 4756cab3a1bSJed Brown Not Collective 4766cab3a1bSJed Brown 4776cab3a1bSJed Brown Input Argument: 4786cab3a1bSJed Brown . dm - DM to be used with SNES 4796cab3a1bSJed Brown 4806cab3a1bSJed Brown Output Arguments: 4816cab3a1bSJed Brown + func - relaxation function, see SNESSetGS() for calling sequence 4826cab3a1bSJed Brown - ctx - context for residual evaluation 4836cab3a1bSJed Brown 4846cab3a1bSJed Brown Level: advanced 4856cab3a1bSJed Brown 4866cab3a1bSJed Brown Note: 4876cab3a1bSJed Brown SNESGetGS() is normally used, but it calls this function internally because the user context is actually 4886cab3a1bSJed Brown associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 4896cab3a1bSJed Brown not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 4906cab3a1bSJed Brown 4916cab3a1bSJed Brown .seealso: DMSNESSetContext(), SNESGetGS(), DMSNESGetJacobian(), DMSNESGetFunction() 4926cab3a1bSJed Brown @*/ 4936cab3a1bSJed Brown PetscErrorCode DMSNESGetGS(DM dm,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx) 4946cab3a1bSJed Brown { 4956cab3a1bSJed Brown PetscErrorCode ierr; 4966cab3a1bSJed Brown SNESDM sdm; 4976cab3a1bSJed Brown 4986cab3a1bSJed Brown PetscFunctionBegin; 4996cab3a1bSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 5006cab3a1bSJed Brown ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 5016cab3a1bSJed Brown if (func) *func = sdm->computegs; 5026cab3a1bSJed Brown if (ctx) *ctx = sdm->gsctx; 5036cab3a1bSJed Brown PetscFunctionReturn(0); 5046cab3a1bSJed Brown } 5056cab3a1bSJed Brown 5066cab3a1bSJed Brown #undef __FUNCT__ 5076cab3a1bSJed Brown #define __FUNCT__ "DMSNESSetJacobian" 5086cab3a1bSJed Brown /*@C 509ecfdb398SPeter Brune DMSNESSetJacobian - set SNES Jacobian evaluation function 5106cab3a1bSJed Brown 5116cab3a1bSJed Brown Not Collective 5126cab3a1bSJed Brown 5136cab3a1bSJed Brown Input Argument: 5146cab3a1bSJed Brown + dm - DM to be used with SNES 5156cab3a1bSJed Brown . func - Jacobian evaluation function, see SNESSetJacobian() for calling sequence 5166cab3a1bSJed Brown - ctx - context for residual evaluation 5176cab3a1bSJed Brown 5186cab3a1bSJed Brown Level: advanced 5196cab3a1bSJed Brown 5206cab3a1bSJed Brown Note: 5216cab3a1bSJed Brown SNESSetJacobian() is normally used, but it calls this function internally because the user context is actually 5226cab3a1bSJed Brown associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 5236cab3a1bSJed Brown not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 5246cab3a1bSJed Brown 5256cab3a1bSJed Brown .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESGetJacobian(), SNESSetJacobian() 5266cab3a1bSJed Brown @*/ 5276cab3a1bSJed Brown PetscErrorCode DMSNESSetJacobian(DM dm,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 5286cab3a1bSJed Brown { 5296cab3a1bSJed Brown PetscErrorCode ierr; 5306cab3a1bSJed Brown SNESDM sdm; 5316cab3a1bSJed Brown 5326cab3a1bSJed Brown PetscFunctionBegin; 5336cab3a1bSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 5341fdfe764SBarry Smith if (func || ctx) { 5351fdfe764SBarry Smith ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr); 5361fdfe764SBarry Smith } 5376cab3a1bSJed Brown if (func) sdm->computejacobian = func; 5386cab3a1bSJed Brown if (ctx) sdm->jacobianctx = ctx; 5396cab3a1bSJed Brown PetscFunctionReturn(0); 5406cab3a1bSJed Brown } 5416cab3a1bSJed Brown 5426cab3a1bSJed Brown #undef __FUNCT__ 5436cab3a1bSJed Brown #define __FUNCT__ "DMSNESGetJacobian" 5446cab3a1bSJed Brown /*@C 545ecfdb398SPeter Brune DMSNESGetJacobian - get SNES Jacobian evaluation function 5466cab3a1bSJed Brown 5476cab3a1bSJed Brown Not Collective 5486cab3a1bSJed Brown 5496cab3a1bSJed Brown Input Argument: 5506cab3a1bSJed Brown . dm - DM to be used with SNES 5516cab3a1bSJed Brown 5526cab3a1bSJed Brown Output Arguments: 5536cab3a1bSJed Brown + func - Jacobian evaluation function, see SNESSetJacobian() for calling sequence 5546cab3a1bSJed Brown - ctx - context for residual evaluation 5556cab3a1bSJed Brown 5566cab3a1bSJed Brown Level: advanced 5576cab3a1bSJed Brown 5586cab3a1bSJed Brown Note: 5596cab3a1bSJed Brown SNESGetJacobian() is normally used, but it calls this function internally because the user context is actually 5606cab3a1bSJed Brown associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 5616cab3a1bSJed Brown not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 5626cab3a1bSJed Brown 5636cab3a1bSJed Brown .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian() 5646cab3a1bSJed Brown @*/ 5656cab3a1bSJed Brown PetscErrorCode DMSNESGetJacobian(DM dm,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx) 5666cab3a1bSJed Brown { 5676cab3a1bSJed Brown PetscErrorCode ierr; 5686cab3a1bSJed Brown SNESDM sdm; 5696cab3a1bSJed Brown 5706cab3a1bSJed Brown PetscFunctionBegin; 5716cab3a1bSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 5726cab3a1bSJed Brown ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 5736cab3a1bSJed Brown if (func) *func = sdm->computejacobian; 5746cab3a1bSJed Brown if (ctx) *ctx = sdm->jacobianctx; 5756cab3a1bSJed Brown PetscFunctionReturn(0); 5766cab3a1bSJed Brown } 5776cab3a1bSJed Brown 5786cab3a1bSJed Brown #undef __FUNCT__ 579e03ab78fSPeter Brune #define __FUNCT__ "DMSNESSetPicard" 580e03ab78fSPeter Brune /*@C 581e03ab78fSPeter Brune DMSNESSetPicard - set SNES Picard iteration matrix and RHS evaluation functions. 582e03ab78fSPeter Brune 583e03ab78fSPeter Brune Not Collective 584e03ab78fSPeter Brune 585e03ab78fSPeter Brune Input Argument: 586e03ab78fSPeter Brune + dm - DM to be used with SNES 587e03ab78fSPeter Brune . func - RHS evaluation function, see SNESSetFunction() for calling sequence 588e03ab78fSPeter Brune . pjac - Picard matrix evaluation function, see SNESSetJacobian() for calling sequence 589e03ab78fSPeter Brune - ctx - context for residual evaluation 590e03ab78fSPeter Brune 591e03ab78fSPeter Brune Level: advanced 592e03ab78fSPeter Brune 593e03ab78fSPeter Brune .seealso: SNESSetPicard(), DMSNESSetFunction(), DMSNESSetJacobian() 594e03ab78fSPeter Brune @*/ 595e03ab78fSPeter Brune PetscErrorCode DMSNESSetPicard(DM dm,PetscErrorCode (*pfunc)(SNES,Vec,Vec,void*),PetscErrorCode (*pjac)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 596e03ab78fSPeter Brune { 597e03ab78fSPeter Brune PetscErrorCode ierr; 598e03ab78fSPeter Brune SNESDM sdm; 599e03ab78fSPeter Brune 600e03ab78fSPeter Brune PetscFunctionBegin; 601e03ab78fSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 602e03ab78fSPeter Brune ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 603e03ab78fSPeter Brune if (pfunc) sdm->computepfunction = pfunc; 604e03ab78fSPeter Brune if (pjac) sdm->computepjacobian = pjac; 605e03ab78fSPeter Brune if (ctx) sdm->pctx = ctx; 606e03ab78fSPeter Brune PetscFunctionReturn(0); 607e03ab78fSPeter Brune } 608e03ab78fSPeter Brune 6097971a8bfSPeter Brune 6107971a8bfSPeter Brune #undef __FUNCT__ 6117971a8bfSPeter Brune #define __FUNCT__ "DMSNESGetPicard" 6127971a8bfSPeter Brune /*@C 6137971a8bfSPeter Brune DMSNESGetPicard - get SNES Picard iteration evaluation functions 6147971a8bfSPeter Brune 6157971a8bfSPeter Brune Not Collective 6167971a8bfSPeter Brune 6177971a8bfSPeter Brune Input Argument: 6187971a8bfSPeter Brune . dm - DM to be used with SNES 6197971a8bfSPeter Brune 6207971a8bfSPeter Brune Output Arguments: 6217971a8bfSPeter Brune + pfunc - Jacobian evaluation function, see SNESSetJacobian() for calling sequence 6227971a8bfSPeter Brune . pjac - RHS evaluation function, see SNESSetFunction() for calling sequence 6237971a8bfSPeter Brune - ctx - context for residual evaluation 6247971a8bfSPeter Brune 6257971a8bfSPeter Brune Level: advanced 6267971a8bfSPeter Brune 6277971a8bfSPeter Brune .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian() 6287971a8bfSPeter Brune @*/ 6297971a8bfSPeter Brune PetscErrorCode DMSNESGetPicard(DM dm,PetscErrorCode (**pfunc)(SNES,Vec,Vec,void*),PetscErrorCode (**pjac)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx) 6307971a8bfSPeter Brune { 6317971a8bfSPeter Brune PetscErrorCode ierr; 6327971a8bfSPeter Brune SNESDM sdm; 6337971a8bfSPeter Brune 6347971a8bfSPeter Brune PetscFunctionBegin; 6357971a8bfSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 6367971a8bfSPeter Brune ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 6377971a8bfSPeter Brune if (pfunc) *pfunc = sdm->computepfunction; 6387971a8bfSPeter Brune if (pjac) *pjac = sdm->computepjacobian; 6397971a8bfSPeter Brune if (ctx) *ctx = sdm->pctx; 6407971a8bfSPeter Brune PetscFunctionReturn(0); 6417971a8bfSPeter Brune } 6427971a8bfSPeter Brune 6437971a8bfSPeter Brune 644e03ab78fSPeter Brune #undef __FUNCT__ 6456cab3a1bSJed Brown #define __FUNCT__ "SNESDefaultComputeFunction_DMLegacy" 6466cab3a1bSJed Brown static PetscErrorCode SNESDefaultComputeFunction_DMLegacy(SNES snes,Vec X,Vec F,void *ctx) 6476cab3a1bSJed Brown { 6486cab3a1bSJed Brown PetscErrorCode ierr; 6496cab3a1bSJed Brown DM dm; 6506cab3a1bSJed Brown 6516cab3a1bSJed Brown PetscFunctionBegin; 6526cab3a1bSJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 6536cab3a1bSJed Brown ierr = DMComputeFunction(dm,X,F);CHKERRQ(ierr); 6546cab3a1bSJed Brown PetscFunctionReturn(0); 6556cab3a1bSJed Brown } 6566cab3a1bSJed Brown 6576cab3a1bSJed Brown #undef __FUNCT__ 6586cab3a1bSJed Brown #define __FUNCT__ "SNESDefaultComputeJacobian_DMLegacy" 6596cab3a1bSJed Brown static PetscErrorCode SNESDefaultComputeJacobian_DMLegacy(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *ctx) 6606cab3a1bSJed Brown { 6616cab3a1bSJed Brown PetscErrorCode ierr; 6626cab3a1bSJed Brown DM dm; 6636cab3a1bSJed Brown 6646cab3a1bSJed Brown PetscFunctionBegin; 6656cab3a1bSJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 6666cab3a1bSJed Brown ierr = DMComputeJacobian(dm,X,*A,*B,mstr);CHKERRQ(ierr); 6676cab3a1bSJed Brown PetscFunctionReturn(0); 6686cab3a1bSJed Brown } 6696cab3a1bSJed Brown 6706cab3a1bSJed Brown #undef __FUNCT__ 6716cab3a1bSJed Brown #define __FUNCT__ "DMSNESSetUpLegacy" 6726cab3a1bSJed Brown /* Sets up calling of legacy DM routines */ 6736cab3a1bSJed Brown PetscErrorCode DMSNESSetUpLegacy(DM dm) 6746cab3a1bSJed Brown { 6756cab3a1bSJed Brown PetscErrorCode ierr; 6766cab3a1bSJed Brown SNESDM sdm; 6776cab3a1bSJed Brown 6786cab3a1bSJed Brown PetscFunctionBegin; 6796cab3a1bSJed Brown ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 6806cab3a1bSJed Brown if (!sdm->computefunction) {ierr = DMSNESSetFunction(dm,SNESDefaultComputeFunction_DMLegacy,PETSC_NULL);CHKERRQ(ierr);} 6816cab3a1bSJed Brown ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 682679f678eSPeter Brune if (!sdm->computejacobian) { 683679f678eSPeter Brune if (dm->ops->functionj) { 684679f678eSPeter Brune ierr = DMSNESSetJacobian(dm,SNESDefaultComputeJacobian_DMLegacy,PETSC_NULL);CHKERRQ(ierr); 685679f678eSPeter Brune } else { 686679f678eSPeter Brune ierr = DMSNESSetJacobian(dm,SNESDefaultComputeJacobianColor,PETSC_NULL);CHKERRQ(ierr); 687679f678eSPeter Brune } 688679f678eSPeter Brune } 6896cab3a1bSJed Brown PetscFunctionReturn(0); 6906cab3a1bSJed Brown } 6916427ad5bSPeter Brune 6926427ad5bSPeter Brune /* block functions */ 6936427ad5bSPeter Brune 6946427ad5bSPeter Brune #undef __FUNCT__ 6956427ad5bSPeter Brune #define __FUNCT__ "DMSNESSetBlockFunction" 6966427ad5bSPeter Brune /*@C 6976427ad5bSPeter Brune DMSNESSetBlockFunction - set SNES residual evaluation function 6986427ad5bSPeter Brune 6996427ad5bSPeter Brune Not Collective 7006427ad5bSPeter Brune 7016427ad5bSPeter Brune Input Arguments: 7026427ad5bSPeter Brune + dm - DM to be used with SNES 7036427ad5bSPeter Brune . func - residual evaluation function, see SNESSetFunction() for calling sequence 7046427ad5bSPeter Brune - ctx - context for residual evaluation 7056427ad5bSPeter Brune 7066427ad5bSPeter Brune Level: developer 7076427ad5bSPeter Brune 7086427ad5bSPeter Brune Note: 7096427ad5bSPeter Brune Mostly for use in DM implementations and transferred to a block function rather than being called from here. 7106427ad5bSPeter Brune 7116427ad5bSPeter Brune .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian() 7126427ad5bSPeter Brune @*/ 7136427ad5bSPeter Brune PetscErrorCode DMSNESSetBlockFunction(DM dm,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx) 7146427ad5bSPeter Brune { 7156427ad5bSPeter Brune PetscErrorCode ierr; 7166427ad5bSPeter Brune SNESDM sdm; 7176427ad5bSPeter Brune 7186427ad5bSPeter Brune PetscFunctionBegin; 7196427ad5bSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 720fdaff8d6SPeter Brune if (func || ctx) { 721fdaff8d6SPeter Brune ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr); 722fdaff8d6SPeter Brune } 7236427ad5bSPeter Brune if (func) sdm->computeblockfunction = func; 7246427ad5bSPeter Brune if (ctx) sdm->blockfunctionctx = ctx; 7256427ad5bSPeter Brune PetscFunctionReturn(0); 7266427ad5bSPeter Brune } 7276427ad5bSPeter Brune 7286427ad5bSPeter Brune #undef __FUNCT__ 7296427ad5bSPeter Brune #define __FUNCT__ "DMSNESGetBlockFunction" 7306427ad5bSPeter Brune /*@C 7316427ad5bSPeter Brune DMSNESGetBlockFunction - get SNES residual evaluation function 7326427ad5bSPeter Brune 7336427ad5bSPeter Brune Not Collective 7346427ad5bSPeter Brune 7356427ad5bSPeter Brune Input Argument: 7366427ad5bSPeter Brune . dm - DM to be used with SNES 7376427ad5bSPeter Brune 7386427ad5bSPeter Brune Output Arguments: 7396427ad5bSPeter Brune + func - residual evaluation function, see SNESSetFunction() for calling sequence 7406427ad5bSPeter Brune - ctx - context for residual evaluation 7416427ad5bSPeter Brune 7426427ad5bSPeter Brune Level: developer 7436427ad5bSPeter Brune 7446427ad5bSPeter Brune .seealso: DMSNESSetContext(), DMSNESSetFunction(), SNESSetFunction() 7456427ad5bSPeter Brune @*/ 7466427ad5bSPeter Brune PetscErrorCode DMSNESGetBlockFunction(DM dm,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx) 7476427ad5bSPeter Brune { 7486427ad5bSPeter Brune PetscErrorCode ierr; 7496427ad5bSPeter Brune SNESDM sdm; 7506427ad5bSPeter Brune 7516427ad5bSPeter Brune PetscFunctionBegin; 7526427ad5bSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 7536427ad5bSPeter Brune ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 7546427ad5bSPeter Brune if (func) *func = sdm->computeblockfunction; 7556427ad5bSPeter Brune if (ctx) *ctx = sdm->blockfunctionctx; 7566427ad5bSPeter Brune PetscFunctionReturn(0); 7576427ad5bSPeter Brune } 7586427ad5bSPeter Brune 7596427ad5bSPeter Brune 7606427ad5bSPeter Brune #undef __FUNCT__ 7616427ad5bSPeter Brune #define __FUNCT__ "DMSNESSetBlockJacobian" 7626427ad5bSPeter Brune /*@C 7636427ad5bSPeter Brune DMSNESSetJacobian - set SNES Jacobian evaluation function 7646427ad5bSPeter Brune 7656427ad5bSPeter Brune Not Collective 7666427ad5bSPeter Brune 7676427ad5bSPeter Brune Input Argument: 7686427ad5bSPeter Brune + dm - DM to be used with SNES 7696427ad5bSPeter Brune . func - Jacobian evaluation function, see SNESSetJacobian() for calling sequence 7706427ad5bSPeter Brune - ctx - context for residual evaluation 7716427ad5bSPeter Brune 7726427ad5bSPeter Brune Level: advanced 7736427ad5bSPeter Brune 7746427ad5bSPeter Brune Note: 7756427ad5bSPeter Brune Mostly for use in DM implementations and transferred to a block function rather than being called from here. 7766427ad5bSPeter Brune 7776427ad5bSPeter Brune .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESGetJacobian(), SNESSetJacobian() 7786427ad5bSPeter Brune @*/ 7796427ad5bSPeter Brune PetscErrorCode DMSNESSetBlockJacobian(DM dm,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 7806427ad5bSPeter Brune { 7816427ad5bSPeter Brune PetscErrorCode ierr; 7826427ad5bSPeter Brune SNESDM sdm; 7836427ad5bSPeter Brune 7846427ad5bSPeter Brune PetscFunctionBegin; 7856427ad5bSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 786fdaff8d6SPeter Brune if (func || ctx) { 787fdaff8d6SPeter Brune ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr); 788fdaff8d6SPeter Brune } 7896427ad5bSPeter Brune if (func) sdm->computeblockjacobian = func; 7906427ad5bSPeter Brune if (ctx) sdm->blockjacobianctx = ctx; 7916427ad5bSPeter Brune PetscFunctionReturn(0); 7926427ad5bSPeter Brune } 7936427ad5bSPeter Brune 7946427ad5bSPeter Brune #undef __FUNCT__ 7956427ad5bSPeter Brune #define __FUNCT__ "DMSNESGetBlockJacobian" 7966427ad5bSPeter Brune /*@C 7976427ad5bSPeter Brune DMSNESGetBlockJacobian - get SNES Jacobian evaluation function 7986427ad5bSPeter Brune 7996427ad5bSPeter Brune Not Collective 8006427ad5bSPeter Brune 8016427ad5bSPeter Brune Input Argument: 8026427ad5bSPeter Brune . dm - DM to be used with SNES 8036427ad5bSPeter Brune 8046427ad5bSPeter Brune Output Arguments: 8056427ad5bSPeter Brune + func - Jacobian evaluation function, see SNESSetJacobian() for calling sequence 8066427ad5bSPeter Brune - ctx - context for residual evaluation 8076427ad5bSPeter Brune 8086427ad5bSPeter Brune Level: advanced 8096427ad5bSPeter Brune 8106427ad5bSPeter Brune .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian() 8116427ad5bSPeter Brune @*/ 8126427ad5bSPeter Brune PetscErrorCode DMSNESGetBlockJacobian(DM dm,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx) 8136427ad5bSPeter Brune { 8146427ad5bSPeter Brune PetscErrorCode ierr; 8156427ad5bSPeter Brune SNESDM sdm; 8166427ad5bSPeter Brune 8176427ad5bSPeter Brune PetscFunctionBegin; 8186427ad5bSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 8196427ad5bSPeter Brune ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 8206427ad5bSPeter Brune if (func) *func = sdm->computeblockjacobian; 8216427ad5bSPeter Brune if (ctx) *ctx = sdm->blockjacobianctx; 8226427ad5bSPeter Brune PetscFunctionReturn(0); 8236427ad5bSPeter Brune } 824