16cab3a1bSJed Brown #include <petscdmda.h> /*I "petscdmda.h" I*/ 2af0996ceSBarry Smith #include <petsc/private/dmimpl.h> 3af0996ceSBarry Smith #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/ 46cab3a1bSJed Brown 56cab3a1bSJed Brown /* This structure holds the user-provided DMDA callbacks */ 66cab3a1bSJed Brown typedef struct { 70bf52853SStefano Zampini /* array versions for vector data */ 88434afd1SBarry Smith DMDASNESFunctionFn *residuallocal; 98434afd1SBarry Smith DMDASNESJacobianFn *jacobianlocal; 108434afd1SBarry Smith DMDASNESObjectiveFn *objectivelocal; 115505f7afSJunchao Zhang 120bf52853SStefano Zampini /* Vec version for vector data */ 138434afd1SBarry Smith DMDASNESFunctionVecFn *residuallocalvec; 148434afd1SBarry Smith DMDASNESJacobianVecFn *jacobianlocalvec; 158434afd1SBarry Smith DMDASNESObjectiveVecFn *objectivelocalvec; 160bf52853SStefano Zampini 170bf52853SStefano Zampini /* user contexts */ 186cab3a1bSJed Brown void *residuallocalctx; 196cab3a1bSJed Brown void *jacobianlocalctx; 202a4ee8f2SPeter Brune void *objectivelocalctx; 21709ac0d2SJed Brown InsertMode residuallocalimode; 22dcad5f1cSBarry Smith 23dcad5f1cSBarry Smith /* For Picard iteration defined locally */ 24dcad5f1cSBarry Smith PetscErrorCode (*rhsplocal)(DMDALocalInfo *, void *, void *, void *); 25d1e9a80fSBarry Smith PetscErrorCode (*jacobianplocal)(DMDALocalInfo *, void *, Mat, Mat, void *); 26dcad5f1cSBarry Smith void *picardlocalctx; 27942e3340SBarry Smith } DMSNES_DA; 286cab3a1bSJed Brown 29d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSNESDestroy_DMDA(DMSNES sdm) 30d71ae5a4SJacob Faibussowitsch { 316cab3a1bSJed Brown PetscFunctionBegin; 329566063dSJacob Faibussowitsch PetscCall(PetscFree(sdm->data)); 333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 346cab3a1bSJed Brown } 356cab3a1bSJed Brown 36d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSNESDuplicate_DMDA(DMSNES oldsdm, DMSNES sdm) 37d71ae5a4SJacob Faibussowitsch { 38e3c0b8a2SPeter Brune PetscFunctionBegin; 394dfa11a4SJacob Faibussowitsch PetscCall(PetscNew((DMSNES_DA **)&sdm->data)); 401baa6e33SBarry Smith if (oldsdm->data) PetscCall(PetscMemcpy(sdm->data, oldsdm->data, sizeof(DMSNES_DA))); 413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 42e3c0b8a2SPeter Brune } 43e3c0b8a2SPeter Brune 44d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDASNESGetContext(DM dm, DMSNES sdm, DMSNES_DA **dmdasnes) 45d71ae5a4SJacob Faibussowitsch { 466cab3a1bSJed Brown PetscFunctionBegin; 470298fd71SBarry Smith *dmdasnes = NULL; 486cab3a1bSJed Brown if (!sdm->data) { 494dfa11a4SJacob Faibussowitsch PetscCall(PetscNew((DMSNES_DA **)&sdm->data)); 5022c6f798SBarry Smith sdm->ops->destroy = DMSNESDestroy_DMDA; 5122c6f798SBarry Smith sdm->ops->duplicate = DMSNESDuplicate_DMDA; 526cab3a1bSJed Brown } 53942e3340SBarry Smith *dmdasnes = (DMSNES_DA *)sdm->data; 543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 556cab3a1bSJed Brown } 566cab3a1bSJed Brown 57d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESComputeFunction_DMDA(SNES snes, Vec X, Vec F, void *ctx) 58d71ae5a4SJacob Faibussowitsch { 596cab3a1bSJed Brown DM dm; 60942e3340SBarry Smith DMSNES_DA *dmdasnes = (DMSNES_DA *)ctx; 616cab3a1bSJed Brown DMDALocalInfo info; 626cab3a1bSJed Brown Vec Xloc; 630bf52853SStefano Zampini void *x, *f, *rctx; 646cab3a1bSJed Brown 656cab3a1bSJed Brown PetscFunctionBegin; 662961e153SJed Brown PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 672961e153SJed Brown PetscValidHeaderSpecific(X, VEC_CLASSID, 2); 682961e153SJed Brown PetscValidHeaderSpecific(F, VEC_CLASSID, 3); 695505f7afSJunchao Zhang PetscCheck(dmdasnes->residuallocal || dmdasnes->residuallocalvec, PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Corrupt context"); 709566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 719566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &Xloc)); 729566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc)); 739566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc)); 749566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(dm, &info)); 75*49abdd8aSBarry Smith rctx = dmdasnes->residuallocalctx ? dmdasnes->residuallocalctx : snes->ctx; 76709ac0d2SJed Brown switch (dmdasnes->residuallocalimode) { 77709ac0d2SJed Brown case INSERT_VALUES: { 789566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(SNES_FunctionEval, snes, X, F, 0)); 790bf52853SStefano Zampini if (dmdasnes->residuallocalvec) PetscCallBack("SNES DMDA local callback function", (*dmdasnes->residuallocalvec)(&info, Xloc, F, rctx)); 80ef1023bdSBarry Smith else { 815505f7afSJunchao Zhang PetscCall(DMDAVecGetArray(dm, Xloc, &x)); 825505f7afSJunchao Zhang PetscCall(DMDAVecGetArray(dm, F, &f)); 830bf52853SStefano Zampini PetscCallBack("SNES DMDA local callback function", (*dmdasnes->residuallocal)(&info, x, f, rctx)); 845505f7afSJunchao Zhang PetscCall(DMDAVecRestoreArray(dm, Xloc, &x)); 855505f7afSJunchao Zhang PetscCall(DMDAVecRestoreArray(dm, F, &f)); 865505f7afSJunchao Zhang } 879566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(SNES_FunctionEval, snes, X, F, 0)); 88709ac0d2SJed Brown } break; 89709ac0d2SJed Brown case ADD_VALUES: { 90709ac0d2SJed Brown Vec Floc; 919566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &Floc)); 929566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Floc)); 939566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(SNES_FunctionEval, snes, X, F, 0)); 940bf52853SStefano Zampini if (dmdasnes->residuallocalvec) PetscCallBack("SNES DMDA local callback function", (*dmdasnes->residuallocalvec)(&info, Xloc, Floc, rctx)); 95ef1023bdSBarry Smith else { 965505f7afSJunchao Zhang PetscCall(DMDAVecGetArray(dm, Xloc, &x)); 975505f7afSJunchao Zhang PetscCall(DMDAVecGetArray(dm, Floc, &f)); 980bf52853SStefano Zampini PetscCallBack("SNES DMDA local callback function", (*dmdasnes->residuallocal)(&info, x, f, rctx)); 995505f7afSJunchao Zhang PetscCall(DMDAVecRestoreArray(dm, Xloc, &x)); 1005505f7afSJunchao Zhang PetscCall(DMDAVecRestoreArray(dm, Floc, &f)); 1015505f7afSJunchao Zhang } 1029566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(SNES_FunctionEval, snes, X, F, 0)); 1039566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 1049566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(dm, Floc, ADD_VALUES, F)); 1059566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(dm, Floc, ADD_VALUES, F)); 1069566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &Floc)); 107709ac0d2SJed Brown } break; 108d71ae5a4SJacob Faibussowitsch default: 109d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_INCOMP, "Cannot use imode=%d", (int)dmdasnes->residuallocalimode); 110709ac0d2SJed Brown } 1119566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &Xloc)); 112f480ea8aSBarry Smith PetscCall(VecFlag(F, snes->domainerror)); 1133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1146cab3a1bSJed Brown } 1156cab3a1bSJed Brown 116d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESComputeObjective_DMDA(SNES snes, Vec X, PetscReal *ob, void *ctx) 117d71ae5a4SJacob Faibussowitsch { 1182a4ee8f2SPeter Brune DM dm; 119942e3340SBarry Smith DMSNES_DA *dmdasnes = (DMSNES_DA *)ctx; 1202a4ee8f2SPeter Brune DMDALocalInfo info; 1212a4ee8f2SPeter Brune Vec Xloc; 1220bf52853SStefano Zampini void *x, *octx; 1232a4ee8f2SPeter Brune 1242a4ee8f2SPeter Brune PetscFunctionBegin; 1252a4ee8f2SPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 1262a4ee8f2SPeter Brune PetscValidHeaderSpecific(X, VEC_CLASSID, 2); 1274f572ea9SToby Isaac PetscAssertPointer(ob, 3); 1285505f7afSJunchao Zhang PetscCheck(dmdasnes->objectivelocal || dmdasnes->objectivelocalvec, PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Corrupt context"); 1299566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 1309566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &Xloc)); 1319566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc)); 1329566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc)); 1339566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(dm, &info)); 134*49abdd8aSBarry Smith octx = dmdasnes->objectivelocalctx ? dmdasnes->objectivelocalctx : snes->ctx; 1350bf52853SStefano Zampini if (dmdasnes->objectivelocalvec) PetscCallBack("SNES DMDA local callback objective", (*dmdasnes->objectivelocalvec)(&info, Xloc, ob, octx)); 136ef1023bdSBarry Smith else { 1379566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(dm, Xloc, &x)); 1380bf52853SStefano Zampini PetscCallBack("SNES DMDA local callback objective", (*dmdasnes->objectivelocal)(&info, x, ob, octx)); 1399566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dm, Xloc, &x)); 1405505f7afSJunchao Zhang } 1419566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &Xloc)); 142462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, ob, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)snes))); 1433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1442a4ee8f2SPeter Brune } 1452a4ee8f2SPeter Brune 146cc2e6a90SBarry Smith /* Routine is called by example, hence must be labeled PETSC_EXTERN */ 147d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESComputeJacobian_DMDA(SNES snes, Vec X, Mat A, Mat B, void *ctx) 148d71ae5a4SJacob Faibussowitsch { 1492961e153SJed Brown DM dm; 150942e3340SBarry Smith DMSNES_DA *dmdasnes = (DMSNES_DA *)ctx; 1512961e153SJed Brown DMDALocalInfo info; 1522961e153SJed Brown Vec Xloc; 1530bf52853SStefano Zampini void *x, *jctx; 1542961e153SJed Brown 1552961e153SJed Brown PetscFunctionBegin; 1565505f7afSJunchao Zhang PetscCheck(dmdasnes->residuallocal || dmdasnes->residuallocalvec, PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Corrupt context"); 1579566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 158*49abdd8aSBarry Smith jctx = dmdasnes->jacobianlocalctx ? dmdasnes->jacobianlocalctx : snes->ctx; 1590bf52853SStefano Zampini if (dmdasnes->jacobianlocal || dmdasnes->jacobianlocalvec) { 1609566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &Xloc)); 1619566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc)); 1629566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc)); 1639566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(dm, &info)); 1640bf52853SStefano Zampini if (dmdasnes->jacobianlocalvec) PetscCallBack("SNES DMDA local callback Jacobian", (*dmdasnes->jacobianlocalvec)(&info, Xloc, A, B, jctx)); 165ef1023bdSBarry Smith else { 1669566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(dm, Xloc, &x)); 1670bf52853SStefano Zampini PetscCallBack("SNES DMDA local callback Jacobian", (*dmdasnes->jacobianlocal)(&info, x, A, B, jctx)); 1689566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dm, Xloc, &x)); 1695505f7afSJunchao Zhang } 1709566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &Xloc)); 1712961e153SJed Brown } else { 1722961e153SJed Brown MatFDColoring fdcoloring; 1739566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)dm, "DMDASNES_FDCOLORING", (PetscObject *)&fdcoloring)); 1742961e153SJed Brown if (!fdcoloring) { 1752961e153SJed Brown ISColoring coloring; 1762961e153SJed Brown 1779566063dSJacob Faibussowitsch PetscCall(DMCreateColoring(dm, dm->coloringtype, &coloring)); 1789566063dSJacob Faibussowitsch PetscCall(MatFDColoringCreate(B, coloring, &fdcoloring)); 1792961e153SJed Brown switch (dm->coloringtype) { 180d71ae5a4SJacob Faibussowitsch case IS_COLORING_GLOBAL: 181d71ae5a4SJacob Faibussowitsch PetscCall(MatFDColoringSetFunction(fdcoloring, (PetscErrorCode (*)(void))SNESComputeFunction_DMDA, dmdasnes)); 182d71ae5a4SJacob Faibussowitsch break; 183d71ae5a4SJacob Faibussowitsch default: 184d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "No support for coloring type '%s'", ISColoringTypes[dm->coloringtype]); 1852961e153SJed Brown } 1869566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fdcoloring, ((PetscObject)dm)->prefix)); 1879566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetFromOptions(fdcoloring)); 1889566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetUp(B, coloring, fdcoloring)); 1899566063dSJacob Faibussowitsch PetscCall(ISColoringDestroy(&coloring)); 1909566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)dm, "DMDASNES_FDCOLORING", (PetscObject)fdcoloring)); 1919566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)fdcoloring)); 1922961e153SJed Brown 1932961e153SJed Brown /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call 1942961e153SJed Brown * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the 1952961e153SJed Brown * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually 1962961e153SJed Brown * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be 197140e18c1SBarry Smith * taken when the PetscObjectList for the Vec inside MatFDColoring is destroyed. 1982961e153SJed Brown */ 1999566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)dm)); 2002961e153SJed Brown } 2019566063dSJacob Faibussowitsch PetscCall(MatFDColoringApply(B, fdcoloring, X, snes)); 2022961e153SJed Brown } 2032961e153SJed Brown /* This will be redundant if the user called both, but it's too common to forget. */ 20494ab13aaSBarry Smith if (A != B) { 2059566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 2069566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 2072961e153SJed Brown } 2083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2092961e153SJed Brown } 2102961e153SJed Brown 2116cab3a1bSJed Brown /*@C 212f6dfbefdSBarry Smith DMDASNESSetFunctionLocal - set a local residual evaluation function for use with `DMDA` 2136cab3a1bSJed Brown 2146cab3a1bSJed Brown Logically Collective 2156cab3a1bSJed Brown 2164165533cSJose E. Roman Input Parameters: 217f6dfbefdSBarry Smith + dm - `DM` to associate callback with 218f6dfbefdSBarry Smith . imode - `INSERT_VALUES` if local function computes owned part, `ADD_VALUES` if it contributes to ghosted part 2196cab3a1bSJed Brown . func - local residual evaluation 2206cab3a1bSJed Brown - ctx - optional context for local residual evaluation 2216cab3a1bSJed Brown 22220f4b53cSBarry Smith Calling sequence of `func`: 223f6dfbefdSBarry Smith + info - `DMDALocalInfo` defining the subdomain to evaluate the residual on 2240038884cSBarry Smith . x - dimensional pointer to state at which to evaluate residual (e.g. PetscScalar *x or **x or ***x) 2250038884cSBarry Smith . f - dimensional pointer to residual, write the residual here (e.g. PetscScalar *f or **f or ***f) 2266cab3a1bSJed Brown - ctx - optional context passed above 2276cab3a1bSJed Brown 2286cab3a1bSJed Brown Level: beginner 2296cab3a1bSJed Brown 230420bcc1bSBarry Smith .seealso: [](ch_snes), `DMDA`, `DMDASNESSetJacobianLocal()`, `DMSNESSetFunction()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()` 2316cab3a1bSJed Brown @*/ 2320b4db180SJacob Faibussowitsch PetscErrorCode DMDASNESSetFunctionLocal(DM dm, InsertMode imode, PetscErrorCode (*func)(DMDALocalInfo *info, void *x, void *f, void *ctx), void *ctx) 233d71ae5a4SJacob Faibussowitsch { 234942e3340SBarry Smith DMSNES sdm; 235942e3340SBarry Smith DMSNES_DA *dmdasnes; 2366cab3a1bSJed Brown 2376cab3a1bSJed Brown PetscFunctionBegin; 2386cab3a1bSJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2399566063dSJacob Faibussowitsch PetscCall(DMGetDMSNESWrite(dm, &sdm)); 2409566063dSJacob Faibussowitsch PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes)); 2411aa26658SKarl Rupp 242709ac0d2SJed Brown dmdasnes->residuallocalimode = imode; 2436cab3a1bSJed Brown dmdasnes->residuallocal = func; 2446cab3a1bSJed Brown dmdasnes->residuallocalctx = ctx; 2451aa26658SKarl Rupp 2469566063dSJacob Faibussowitsch PetscCall(DMSNESSetFunction(dm, SNESComputeFunction_DMDA, dmdasnes)); 24722c6f798SBarry Smith if (!sdm->ops->computejacobian) { /* Call us for the Jacobian too, can be overridden by the user. */ 2489566063dSJacob Faibussowitsch PetscCall(DMSNESSetJacobian(dm, SNESComputeJacobian_DMDA, dmdasnes)); 2492961e153SJed Brown } 2503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2512961e153SJed Brown } 2522961e153SJed Brown 2532961e153SJed Brown /*@C 254f6dfbefdSBarry Smith DMDASNESSetFunctionLocalVec - set a local residual evaluation function that operates on a local vector for `DMDA` 2555505f7afSJunchao Zhang 2565505f7afSJunchao Zhang Logically Collective 2575505f7afSJunchao Zhang 2585505f7afSJunchao Zhang Input Parameters: 259f6dfbefdSBarry Smith + dm - `DM` to associate callback with 260f6dfbefdSBarry Smith . imode - `INSERT_VALUES` if local function computes owned part, `ADD_VALUES` if it contributes to ghosted part 2615505f7afSJunchao Zhang . func - local residual evaluation 2625505f7afSJunchao Zhang - ctx - optional context for local residual evaluation 2635505f7afSJunchao Zhang 26420f4b53cSBarry Smith Calling sequence of `func`: 265f6dfbefdSBarry Smith + info - `DMDALocalInfo` defining the subdomain to evaluate the residual on 2665505f7afSJunchao Zhang . x - state vector at which to evaluate residual 2675505f7afSJunchao Zhang . f - residual vector 2685505f7afSJunchao Zhang - ctx - optional context passed above 2695505f7afSJunchao Zhang 2705505f7afSJunchao Zhang Level: beginner 2715505f7afSJunchao Zhang 272420bcc1bSBarry Smith .seealso: [](ch_snes), `DMDA`, `DMDASNESSetFunctionLocal()`, `DMDASNESSetJacobianLocalVec()`, `DMSNESSetFunction()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()` 2735505f7afSJunchao Zhang @*/ 2740b4db180SJacob Faibussowitsch PetscErrorCode DMDASNESSetFunctionLocalVec(DM dm, InsertMode imode, PetscErrorCode (*func)(DMDALocalInfo *info, Vec x, Vec f, void *ctx), void *ctx) 275d71ae5a4SJacob Faibussowitsch { 2765505f7afSJunchao Zhang DMSNES sdm; 2775505f7afSJunchao Zhang DMSNES_DA *dmdasnes; 2785505f7afSJunchao Zhang 2795505f7afSJunchao Zhang PetscFunctionBegin; 2805505f7afSJunchao Zhang PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2815505f7afSJunchao Zhang PetscCall(DMGetDMSNESWrite(dm, &sdm)); 2825505f7afSJunchao Zhang PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes)); 2835505f7afSJunchao Zhang 2845505f7afSJunchao Zhang dmdasnes->residuallocalimode = imode; 2855505f7afSJunchao Zhang dmdasnes->residuallocalvec = func; 2865505f7afSJunchao Zhang dmdasnes->residuallocalctx = ctx; 2875505f7afSJunchao Zhang 2885505f7afSJunchao Zhang PetscCall(DMSNESSetFunction(dm, SNESComputeFunction_DMDA, dmdasnes)); 2895505f7afSJunchao Zhang if (!sdm->ops->computejacobian) { /* Call us for the Jacobian too, can be overridden by the user. */ 2905505f7afSJunchao Zhang PetscCall(DMSNESSetJacobian(dm, SNESComputeJacobian_DMDA, dmdasnes)); 2915505f7afSJunchao Zhang } 2923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2935505f7afSJunchao Zhang } 2945505f7afSJunchao Zhang 2955505f7afSJunchao Zhang /*@C 296f6dfbefdSBarry Smith DMDASNESSetJacobianLocal - set a local Jacobian evaluation function for use with `DMDA` 2972961e153SJed Brown 2982961e153SJed Brown Logically Collective 2992961e153SJed Brown 3004165533cSJose E. Roman Input Parameters: 301f6dfbefdSBarry Smith + dm - `DM` to associate callback with 302c118d280SBarry Smith . func - local Jacobian evaluation function 3036b7eb5ceSMatthew G. Knepley - ctx - optional context for local Jacobian evaluation 3042961e153SJed Brown 30520f4b53cSBarry Smith Calling sequence of `func`: 306f6dfbefdSBarry Smith + info - `DMDALocalInfo` defining the subdomain to evaluate the Jacobian at 3070038884cSBarry Smith . x - dimensional pointer to state at which to evaluate Jacobian (e.g. PetscScalar *x or **x or ***x) 308c118d280SBarry Smith . J - `Mat` object for the Jacobian 309c118d280SBarry Smith . M - `Mat` object for the Jacobian preconditioner matrix, often `J` 3102961e153SJed Brown - ctx - optional context passed above 3112961e153SJed Brown 3122961e153SJed Brown Level: beginner 3132961e153SJed Brown 314c118d280SBarry Smith Note: 315c118d280SBarry Smith The `J` and `M` matrices are created internally by `DMCreateMatrix()` 316c118d280SBarry Smith 317420bcc1bSBarry Smith .seealso: [](ch_snes), `DMDA`, `DMDASNESSetFunctionLocal()`, `DMSNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()` 3182961e153SJed Brown @*/ 3190b4db180SJacob Faibussowitsch PetscErrorCode DMDASNESSetJacobianLocal(DM dm, PetscErrorCode (*func)(DMDALocalInfo *info, void *x, Mat J, Mat M, void *ctx), void *ctx) 320d71ae5a4SJacob Faibussowitsch { 321942e3340SBarry Smith DMSNES sdm; 322942e3340SBarry Smith DMSNES_DA *dmdasnes; 3232961e153SJed Brown 3242961e153SJed Brown PetscFunctionBegin; 3252961e153SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3269566063dSJacob Faibussowitsch PetscCall(DMGetDMSNESWrite(dm, &sdm)); 3279566063dSJacob Faibussowitsch PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes)); 3281aa26658SKarl Rupp 3292961e153SJed Brown dmdasnes->jacobianlocal = func; 3302961e153SJed Brown dmdasnes->jacobianlocalctx = ctx; 3311aa26658SKarl Rupp 3329566063dSJacob Faibussowitsch PetscCall(DMSNESSetJacobian(dm, SNESComputeJacobian_DMDA, dmdasnes)); 3333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3346cab3a1bSJed Brown } 3352a4ee8f2SPeter Brune 3362a4ee8f2SPeter Brune /*@C 337f6dfbefdSBarry Smith DMDASNESSetJacobianLocalVec - set a local Jacobian evaluation function that operates on a local vector with `DMDA` 3385505f7afSJunchao Zhang 3395505f7afSJunchao Zhang Logically Collective 3405505f7afSJunchao Zhang 3415505f7afSJunchao Zhang Input Parameters: 342f6dfbefdSBarry Smith + dm - `DM` to associate callback with 3435505f7afSJunchao Zhang . func - local Jacobian evaluation 3445505f7afSJunchao Zhang - ctx - optional context for local Jacobian evaluation 3455505f7afSJunchao Zhang 34620f4b53cSBarry Smith Calling sequence of `func`: 347f6dfbefdSBarry Smith + info - `DMDALocalInfo` defining the subdomain to evaluate the Jacobian at 3485505f7afSJunchao Zhang . x - state vector at which to evaluate Jacobian 3490b4db180SJacob Faibussowitsch . J - the Jacobian 3500b4db180SJacob Faibussowitsch . M - approximate Jacobian from which the preconditioner will be computed, often `J` 3515505f7afSJunchao Zhang - ctx - optional context passed above 3525505f7afSJunchao Zhang 3535505f7afSJunchao Zhang Level: beginner 3545505f7afSJunchao Zhang 355420bcc1bSBarry Smith .seealso: [](ch_snes), `DMDA`, `DMDASNESSetJacobianLocal()`, `DMDASNESSetFunctionLocalVec()`, `DMSNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()` 3565505f7afSJunchao Zhang @*/ 3570b4db180SJacob Faibussowitsch PetscErrorCode DMDASNESSetJacobianLocalVec(DM dm, PetscErrorCode (*func)(DMDALocalInfo *info, Vec x, Mat J, Mat M, void *), void *ctx) 358d71ae5a4SJacob Faibussowitsch { 3595505f7afSJunchao Zhang DMSNES sdm; 3605505f7afSJunchao Zhang DMSNES_DA *dmdasnes; 3615505f7afSJunchao Zhang 3625505f7afSJunchao Zhang PetscFunctionBegin; 3635505f7afSJunchao Zhang PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3645505f7afSJunchao Zhang PetscCall(DMGetDMSNESWrite(dm, &sdm)); 3655505f7afSJunchao Zhang PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes)); 3665505f7afSJunchao Zhang 3675505f7afSJunchao Zhang dmdasnes->jacobianlocalvec = func; 3685505f7afSJunchao Zhang dmdasnes->jacobianlocalctx = ctx; 3695505f7afSJunchao Zhang 3705505f7afSJunchao Zhang PetscCall(DMSNESSetJacobian(dm, SNESComputeJacobian_DMDA, dmdasnes)); 3713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3725505f7afSJunchao Zhang } 3735505f7afSJunchao Zhang 3745505f7afSJunchao Zhang /*@C 375f6dfbefdSBarry Smith DMDASNESSetObjectiveLocal - set a local residual evaluation function to used with a `DMDA` 3762a4ee8f2SPeter Brune 3772a4ee8f2SPeter Brune Logically Collective 3782a4ee8f2SPeter Brune 3794165533cSJose E. Roman Input Parameters: 380f6dfbefdSBarry Smith + dm - `DM` to associate callback with 381420bcc1bSBarry Smith . func - local objective evaluation, see `DMDASNESSetObjectiveLocal` for the calling sequence 3822a4ee8f2SPeter Brune - ctx - optional context for local residual evaluation 3832a4ee8f2SPeter Brune 3840bf52853SStefano Zampini Calling sequence of `func`: 3850bf52853SStefano Zampini + info - `DMDALocalInfo` defining the subdomain to evaluate the Jacobian at 3860bf52853SStefano Zampini . x - dimensional pointer to state at which to evaluate the objective (e.g. PetscScalar *x or **x or ***x) 3870bf52853SStefano Zampini . obj - returned objective value for the local subdomain 3880bf52853SStefano Zampini - ctx - optional context passed above 3890bf52853SStefano Zampini 3902a4ee8f2SPeter Brune Level: beginner 3912a4ee8f2SPeter Brune 3928434afd1SBarry Smith .seealso: [](ch_snes), `DMDA`, `DMSNESSetFunction()`, `DMDASNESSetJacobianLocal()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`, `DMDASNESObjectiveFn` 3932a4ee8f2SPeter Brune @*/ 3940bf52853SStefano Zampini PetscErrorCode DMDASNESSetObjectiveLocal(DM dm, PetscErrorCode (*func)(DMDALocalInfo *info, void *x, PetscReal *obj, void *), void *ctx) 395d71ae5a4SJacob Faibussowitsch { 396942e3340SBarry Smith DMSNES sdm; 397942e3340SBarry Smith DMSNES_DA *dmdasnes; 3982a4ee8f2SPeter Brune 3992a4ee8f2SPeter Brune PetscFunctionBegin; 4002a4ee8f2SPeter Brune PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4019566063dSJacob Faibussowitsch PetscCall(DMGetDMSNESWrite(dm, &sdm)); 4029566063dSJacob Faibussowitsch PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes)); 4031aa26658SKarl Rupp 4042a4ee8f2SPeter Brune dmdasnes->objectivelocal = func; 4052a4ee8f2SPeter Brune dmdasnes->objectivelocalctx = ctx; 4061aa26658SKarl Rupp 4079566063dSJacob Faibussowitsch PetscCall(DMSNESSetObjective(dm, SNESComputeObjective_DMDA, dmdasnes)); 4083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4092a4ee8f2SPeter Brune } 410dcad5f1cSBarry Smith 4115505f7afSJunchao Zhang /*@C 412e4094ef1SJacob Faibussowitsch DMDASNESSetObjectiveLocalVec - set a local residual evaluation function that operates on a local vector with `DMDA` 4135505f7afSJunchao Zhang 4145505f7afSJunchao Zhang Logically Collective 4155505f7afSJunchao Zhang 4165505f7afSJunchao Zhang Input Parameters: 417f6dfbefdSBarry Smith + dm - `DM` to associate callback with 418420bcc1bSBarry Smith . func - local objective evaluation, see `DMDASNESSetObjectiveLocalVec` for the calling sequence 4195505f7afSJunchao Zhang - ctx - optional context for local residual evaluation 4205505f7afSJunchao Zhang 4210bf52853SStefano Zampini Calling sequence of `func`: 4220bf52853SStefano Zampini + info - `DMDALocalInfo` defining the subdomain to evaluate the Jacobian at 4230bf52853SStefano Zampini . x - state vector at which to evaluate the objective 4240bf52853SStefano Zampini . obj - returned objective value for the local subdomain 4250bf52853SStefano Zampini - ctx - optional context passed above 4260bf52853SStefano Zampini 4275505f7afSJunchao Zhang Level: beginner 4285505f7afSJunchao Zhang 4298434afd1SBarry Smith .seealso: [](ch_snes), `DMDA`, `DMDASNESSetObjectiveLocal()`, `DMSNESSetFunction()`, `DMDASNESSetJacobianLocalVec()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`, `DMDASNESObjectiveVecFn` 4305505f7afSJunchao Zhang @*/ 4310bf52853SStefano Zampini PetscErrorCode DMDASNESSetObjectiveLocalVec(DM dm, PetscErrorCode (*func)(DMDALocalInfo *info, Vec x, PetscReal *obj, void *), void *ctx) 432d71ae5a4SJacob Faibussowitsch { 4335505f7afSJunchao Zhang DMSNES sdm; 4345505f7afSJunchao Zhang DMSNES_DA *dmdasnes; 4355505f7afSJunchao Zhang 4365505f7afSJunchao Zhang PetscFunctionBegin; 4375505f7afSJunchao Zhang PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4385505f7afSJunchao Zhang PetscCall(DMGetDMSNESWrite(dm, &sdm)); 4395505f7afSJunchao Zhang PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes)); 4405505f7afSJunchao Zhang 4415505f7afSJunchao Zhang dmdasnes->objectivelocalvec = func; 4425505f7afSJunchao Zhang dmdasnes->objectivelocalctx = ctx; 4435505f7afSJunchao Zhang 4445505f7afSJunchao Zhang PetscCall(DMSNESSetObjective(dm, SNESComputeObjective_DMDA, dmdasnes)); 4453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4465505f7afSJunchao Zhang } 4475505f7afSJunchao Zhang 448d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESComputePicard_DMDA(SNES snes, Vec X, Vec F, void *ctx) 449d71ae5a4SJacob Faibussowitsch { 450dcad5f1cSBarry Smith DM dm; 451942e3340SBarry Smith DMSNES_DA *dmdasnes = (DMSNES_DA *)ctx; 452dcad5f1cSBarry Smith DMDALocalInfo info; 453dcad5f1cSBarry Smith Vec Xloc; 454dcad5f1cSBarry Smith void *x, *f; 455dcad5f1cSBarry Smith 456dcad5f1cSBarry Smith PetscFunctionBegin; 457dcad5f1cSBarry Smith PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 458dcad5f1cSBarry Smith PetscValidHeaderSpecific(X, VEC_CLASSID, 2); 459dcad5f1cSBarry Smith PetscValidHeaderSpecific(F, VEC_CLASSID, 3); 46028b400f6SJacob Faibussowitsch PetscCheck(dmdasnes->rhsplocal, PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Corrupt context"); 4619566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 4629566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &Xloc)); 4639566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc)); 4649566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc)); 4659566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(dm, &info)); 4669566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(dm, Xloc, &x)); 467dcad5f1cSBarry Smith switch (dmdasnes->residuallocalimode) { 468dcad5f1cSBarry Smith case INSERT_VALUES: { 4699566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(dm, F, &f)); 470792fecdfSBarry Smith PetscCallBack("SNES Picard DMDA local callback function", (*dmdasnes->rhsplocal)(&info, x, f, dmdasnes->picardlocalctx)); 4719566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dm, F, &f)); 472dcad5f1cSBarry Smith } break; 473dcad5f1cSBarry Smith case ADD_VALUES: { 474dcad5f1cSBarry Smith Vec Floc; 4759566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &Floc)); 4769566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Floc)); 4779566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(dm, Floc, &f)); 478792fecdfSBarry Smith PetscCallBack("SNES Picard DMDA local callback function", (*dmdasnes->rhsplocal)(&info, x, f, dmdasnes->picardlocalctx)); 4799566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dm, Floc, &f)); 4809566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 4819566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(dm, Floc, ADD_VALUES, F)); 4829566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(dm, Floc, ADD_VALUES, F)); 4839566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &Floc)); 484dcad5f1cSBarry Smith } break; 485d71ae5a4SJacob Faibussowitsch default: 486d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_INCOMP, "Cannot use imode=%d", (int)dmdasnes->residuallocalimode); 487dcad5f1cSBarry Smith } 4889566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dm, Xloc, &x)); 4899566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &Xloc)); 4903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 491dcad5f1cSBarry Smith } 492dcad5f1cSBarry Smith 493d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESComputePicardJacobian_DMDA(SNES snes, Vec X, Mat A, Mat B, void *ctx) 494d71ae5a4SJacob Faibussowitsch { 495dcad5f1cSBarry Smith DM dm; 496942e3340SBarry Smith DMSNES_DA *dmdasnes = (DMSNES_DA *)ctx; 497dcad5f1cSBarry Smith DMDALocalInfo info; 498dcad5f1cSBarry Smith Vec Xloc; 499dcad5f1cSBarry Smith void *x; 500dcad5f1cSBarry Smith 501dcad5f1cSBarry Smith PetscFunctionBegin; 50228b400f6SJacob Faibussowitsch PetscCheck(dmdasnes->jacobianplocal, PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Corrupt context"); 5039566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 504dcad5f1cSBarry Smith 5059566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &Xloc)); 5069566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc)); 5079566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc)); 5089566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(dm, &info)); 5099566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(dm, Xloc, &x)); 510792fecdfSBarry Smith PetscCallBack("SNES Picard DMDA local callback Jacobian", (*dmdasnes->jacobianplocal)(&info, x, A, B, dmdasnes->picardlocalctx)); 5119566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dm, Xloc, &x)); 5129566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &Xloc)); 51394ab13aaSBarry Smith if (A != B) { 5149566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 5159566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 516dcad5f1cSBarry Smith } 5173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 518dcad5f1cSBarry Smith } 519dcad5f1cSBarry Smith 520dcad5f1cSBarry Smith /*@C 521dd8e379bSPierre Jolivet DMDASNESSetPicardLocal - set a local right-hand side and matrix evaluation function for Picard iteration with `DMDA` 522dcad5f1cSBarry Smith 523dcad5f1cSBarry Smith Logically Collective 524dcad5f1cSBarry Smith 5254165533cSJose E. Roman Input Parameters: 526f6dfbefdSBarry Smith + dm - `DM` to associate callback with 527f6dfbefdSBarry Smith . imode - `INSERT_VALUES` if local function computes owned part, `ADD_VALUES` if it contributes to ghosted part 528dcad5f1cSBarry Smith . func - local residual evaluation 529ceaaa498SBarry Smith . jac - function to compute Jacobian 530dcad5f1cSBarry Smith - ctx - optional context for local residual evaluation 531dcad5f1cSBarry Smith 53220f4b53cSBarry Smith Calling sequence of `func`: 53320f4b53cSBarry Smith + info - defines the subdomain to evaluate the residual on 534dcad5f1cSBarry Smith . x - dimensional pointer to state at which to evaluate residual 5350b4db180SJacob Faibussowitsch . f - dimensional pointer to residual, write the residual here 5360b4db180SJacob Faibussowitsch - ctx - optional context passed above 5370b4db180SJacob Faibussowitsch 5380b4db180SJacob Faibussowitsch Calling sequence of `jac`: 5390b4db180SJacob Faibussowitsch + info - defines the subdomain to evaluate the residual on 5400b4db180SJacob Faibussowitsch . x - dimensional pointer to state at which to evaluate residual 5410b4db180SJacob Faibussowitsch . jac - the Jacobian 5420b4db180SJacob Faibussowitsch . Jp - approximation to the Jacobian used to compute the preconditioner, often `J` 543dcad5f1cSBarry Smith - ctx - optional context passed above 544dcad5f1cSBarry Smith 545dcad5f1cSBarry Smith Level: beginner 546dcad5f1cSBarry Smith 54720f4b53cSBarry Smith Note: 548420bcc1bSBarry Smith The user must use `SNESSetFunction`(`snes`,`NULL`,`SNESPicardComputeFunction`,&user)); 54920f4b53cSBarry Smith in their code before calling this routine. 55020f4b53cSBarry Smith 551420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `DMDA`, `DMSNESSetFunction()`, `DMDASNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()` 552dcad5f1cSBarry Smith @*/ 5530b4db180SJacob Faibussowitsch PetscErrorCode DMDASNESSetPicardLocal(DM dm, InsertMode imode, PetscErrorCode (*func)(DMDALocalInfo *info, void *x, void *f, void *ctx), PetscErrorCode (*jac)(DMDALocalInfo *info, void *x, Mat jac, Mat Jp, void *ctx), void *ctx) 554d71ae5a4SJacob Faibussowitsch { 555942e3340SBarry Smith DMSNES sdm; 556942e3340SBarry Smith DMSNES_DA *dmdasnes; 557dcad5f1cSBarry Smith 558dcad5f1cSBarry Smith PetscFunctionBegin; 559dcad5f1cSBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5609566063dSJacob Faibussowitsch PetscCall(DMGetDMSNESWrite(dm, &sdm)); 5619566063dSJacob Faibussowitsch PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes)); 5621aa26658SKarl Rupp 563dcad5f1cSBarry Smith dmdasnes->residuallocalimode = imode; 564dcad5f1cSBarry Smith dmdasnes->rhsplocal = func; 565dcad5f1cSBarry Smith dmdasnes->jacobianplocal = jac; 566dcad5f1cSBarry Smith dmdasnes->picardlocalctx = ctx; 5671aa26658SKarl Rupp 5689566063dSJacob Faibussowitsch PetscCall(DMSNESSetPicard(dm, SNESComputePicard_DMDA, SNESComputePicardJacobian_DMDA, dmdasnes)); 5699566063dSJacob Faibussowitsch PetscCall(DMSNESSetMFFunction(dm, SNESComputeFunction_DMDA, dmdasnes)); 5703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 571dcad5f1cSBarry Smith } 572