xref: /petsc/src/snes/utils/dmdasnes.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
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 {
76cab3a1bSJed Brown   PetscErrorCode (*residuallocal)(DMDALocalInfo *, void *, void *, void *);
8d1e9a80fSBarry Smith   PetscErrorCode (*jacobianlocal)(DMDALocalInfo *, void *, Mat, Mat, void *);
92a4ee8f2SPeter Brune   PetscErrorCode (*objectivelocal)(DMDALocalInfo *, void *, PetscReal *, void *);
105505f7afSJunchao Zhang 
115505f7afSJunchao Zhang   PetscErrorCode (*residuallocalvec)(DMDALocalInfo *, Vec, Vec, void *);
125505f7afSJunchao Zhang   PetscErrorCode (*jacobianlocalvec)(DMDALocalInfo *, Vec, Mat, Mat, void *);
135505f7afSJunchao Zhang   PetscErrorCode (*objectivelocalvec)(DMDALocalInfo *, Vec, PetscReal *, void *);
146cab3a1bSJed Brown   void      *residuallocalctx;
156cab3a1bSJed Brown   void      *jacobianlocalctx;
162a4ee8f2SPeter Brune   void      *objectivelocalctx;
17709ac0d2SJed Brown   InsertMode residuallocalimode;
18dcad5f1cSBarry Smith 
19dcad5f1cSBarry Smith   /*   For Picard iteration defined locally */
20dcad5f1cSBarry Smith   PetscErrorCode (*rhsplocal)(DMDALocalInfo *, void *, void *, void *);
21d1e9a80fSBarry Smith   PetscErrorCode (*jacobianplocal)(DMDALocalInfo *, void *, Mat, Mat, void *);
22dcad5f1cSBarry Smith   void *picardlocalctx;
23942e3340SBarry Smith } DMSNES_DA;
246cab3a1bSJed Brown 
259371c9d4SSatish Balay static PetscErrorCode DMSNESDestroy_DMDA(DMSNES sdm) {
266cab3a1bSJed Brown   PetscFunctionBegin;
279566063dSJacob Faibussowitsch   PetscCall(PetscFree(sdm->data));
286cab3a1bSJed Brown   PetscFunctionReturn(0);
296cab3a1bSJed Brown }
306cab3a1bSJed Brown 
319371c9d4SSatish Balay static PetscErrorCode DMSNESDuplicate_DMDA(DMSNES oldsdm, DMSNES sdm) {
32e3c0b8a2SPeter Brune   PetscFunctionBegin;
33*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew((DMSNES_DA **)&sdm->data));
341baa6e33SBarry Smith   if (oldsdm->data) PetscCall(PetscMemcpy(sdm->data, oldsdm->data, sizeof(DMSNES_DA)));
35e3c0b8a2SPeter Brune   PetscFunctionReturn(0);
36e3c0b8a2SPeter Brune }
37e3c0b8a2SPeter Brune 
389371c9d4SSatish Balay static PetscErrorCode DMDASNESGetContext(DM dm, DMSNES sdm, DMSNES_DA **dmdasnes) {
396cab3a1bSJed Brown   PetscFunctionBegin;
400298fd71SBarry Smith   *dmdasnes = NULL;
416cab3a1bSJed Brown   if (!sdm->data) {
42*4dfa11a4SJacob Faibussowitsch     PetscCall(PetscNew((DMSNES_DA **)&sdm->data));
4322c6f798SBarry Smith     sdm->ops->destroy   = DMSNESDestroy_DMDA;
4422c6f798SBarry Smith     sdm->ops->duplicate = DMSNESDuplicate_DMDA;
456cab3a1bSJed Brown   }
46942e3340SBarry Smith   *dmdasnes = (DMSNES_DA *)sdm->data;
476cab3a1bSJed Brown   PetscFunctionReturn(0);
486cab3a1bSJed Brown }
496cab3a1bSJed Brown 
509371c9d4SSatish Balay static PetscErrorCode SNESComputeFunction_DMDA(SNES snes, Vec X, Vec F, void *ctx) {
516cab3a1bSJed Brown   DM            dm;
52942e3340SBarry Smith   DMSNES_DA    *dmdasnes = (DMSNES_DA *)ctx;
536cab3a1bSJed Brown   DMDALocalInfo info;
546cab3a1bSJed Brown   Vec           Xloc;
556cab3a1bSJed Brown   void         *x, *f;
566cab3a1bSJed Brown 
576cab3a1bSJed Brown   PetscFunctionBegin;
582961e153SJed Brown   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
592961e153SJed Brown   PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
602961e153SJed Brown   PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
615505f7afSJunchao Zhang   PetscCheck(dmdasnes->residuallocal || dmdasnes->residuallocalvec, PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Corrupt context");
629566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
639566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &Xloc));
649566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
659566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
669566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm, &info));
67709ac0d2SJed Brown   switch (dmdasnes->residuallocalimode) {
68709ac0d2SJed Brown   case INSERT_VALUES: {
699566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(SNES_FunctionEval, snes, X, F, 0));
70792fecdfSBarry Smith     if (dmdasnes->residuallocalvec) PetscCallBack("SNES DMDA local callback function", (*dmdasnes->residuallocalvec)(&info, Xloc, F, dmdasnes->residuallocalctx));
71ef1023bdSBarry Smith     else {
725505f7afSJunchao Zhang       PetscCall(DMDAVecGetArray(dm, Xloc, &x));
735505f7afSJunchao Zhang       PetscCall(DMDAVecGetArray(dm, F, &f));
74792fecdfSBarry Smith       PetscCallBack("SNES DMDA local callback function", (*dmdasnes->residuallocal)(&info, x, f, dmdasnes->residuallocalctx));
755505f7afSJunchao Zhang       PetscCall(DMDAVecRestoreArray(dm, Xloc, &x));
765505f7afSJunchao Zhang       PetscCall(DMDAVecRestoreArray(dm, F, &f));
775505f7afSJunchao Zhang     }
789566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(SNES_FunctionEval, snes, X, F, 0));
79709ac0d2SJed Brown   } break;
80709ac0d2SJed Brown   case ADD_VALUES: {
81709ac0d2SJed Brown     Vec Floc;
829566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm, &Floc));
839566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(Floc));
849566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(SNES_FunctionEval, snes, X, F, 0));
85792fecdfSBarry Smith     if (dmdasnes->residuallocalvec) PetscCallBack("SNES DMDA local callback function", (*dmdasnes->residuallocalvec)(&info, Xloc, Floc, dmdasnes->residuallocalctx));
86ef1023bdSBarry Smith     else {
875505f7afSJunchao Zhang       PetscCall(DMDAVecGetArray(dm, Xloc, &x));
885505f7afSJunchao Zhang       PetscCall(DMDAVecGetArray(dm, Floc, &f));
89792fecdfSBarry Smith       PetscCallBack("SNES DMDA local callback function", (*dmdasnes->residuallocal)(&info, x, f, dmdasnes->residuallocalctx));
905505f7afSJunchao Zhang       PetscCall(DMDAVecRestoreArray(dm, Xloc, &x));
915505f7afSJunchao Zhang       PetscCall(DMDAVecRestoreArray(dm, Floc, &f));
925505f7afSJunchao Zhang     }
939566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(SNES_FunctionEval, snes, X, F, 0));
949566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(F));
959566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalBegin(dm, Floc, ADD_VALUES, F));
969566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalEnd(dm, Floc, ADD_VALUES, F));
979566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm, &Floc));
98709ac0d2SJed Brown   } break;
9998921bdaSJacob Faibussowitsch   default: SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_INCOMP, "Cannot use imode=%d", (int)dmdasnes->residuallocalimode);
100709ac0d2SJed Brown   }
1019566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &Xloc));
1021baa6e33SBarry Smith   if (snes->domainerror) PetscCall(VecSetInf(F));
1036cab3a1bSJed Brown   PetscFunctionReturn(0);
1046cab3a1bSJed Brown }
1056cab3a1bSJed Brown 
1069371c9d4SSatish Balay static PetscErrorCode SNESComputeObjective_DMDA(SNES snes, Vec X, PetscReal *ob, void *ctx) {
1072a4ee8f2SPeter Brune   DM            dm;
108942e3340SBarry Smith   DMSNES_DA    *dmdasnes = (DMSNES_DA *)ctx;
1092a4ee8f2SPeter Brune   DMDALocalInfo info;
1102a4ee8f2SPeter Brune   Vec           Xloc;
1112a4ee8f2SPeter Brune   void         *x;
1122a4ee8f2SPeter Brune 
1132a4ee8f2SPeter Brune   PetscFunctionBegin;
1142a4ee8f2SPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1152a4ee8f2SPeter Brune   PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
116dadcf809SJacob Faibussowitsch   PetscValidRealPointer(ob, 3);
1175505f7afSJunchao Zhang   PetscCheck(dmdasnes->objectivelocal || dmdasnes->objectivelocalvec, PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Corrupt context");
1189566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
1199566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &Xloc));
1209566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
1219566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
1229566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm, &info));
123792fecdfSBarry Smith   if (dmdasnes->objectivelocalvec) PetscCallBack("SNES DMDA local callback objective", (*dmdasnes->objectivelocalvec)(&info, Xloc, ob, dmdasnes->objectivelocalctx));
124ef1023bdSBarry Smith   else {
1259566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dm, Xloc, &x));
126792fecdfSBarry Smith     PetscCallBack("SNES DMDA local callback objective", (*dmdasnes->objectivelocal)(&info, x, ob, dmdasnes->objectivelocalctx));
1279566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dm, Xloc, &x));
1285505f7afSJunchao Zhang   }
1299566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &Xloc));
1302a4ee8f2SPeter Brune   PetscFunctionReturn(0);
1312a4ee8f2SPeter Brune }
1322a4ee8f2SPeter Brune 
133cc2e6a90SBarry Smith /* Routine is called by example, hence must be labeled PETSC_EXTERN */
1349371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode SNESComputeJacobian_DMDA(SNES snes, Vec X, Mat A, Mat B, void *ctx) {
1352961e153SJed Brown   DM            dm;
136942e3340SBarry Smith   DMSNES_DA    *dmdasnes = (DMSNES_DA *)ctx;
1372961e153SJed Brown   DMDALocalInfo info;
1382961e153SJed Brown   Vec           Xloc;
1392961e153SJed Brown   void         *x;
1402961e153SJed Brown 
1412961e153SJed Brown   PetscFunctionBegin;
1425505f7afSJunchao Zhang   PetscCheck(dmdasnes->residuallocal || dmdasnes->residuallocalvec, PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Corrupt context");
1439566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
1442961e153SJed Brown 
1455505f7afSJunchao Zhang   if (dmdasnes->jacobianlocal || dmdasnes->jacobianlocalctx) {
1469566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm, &Xloc));
1479566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
1489566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
1499566063dSJacob Faibussowitsch     PetscCall(DMDAGetLocalInfo(dm, &info));
150792fecdfSBarry Smith     if (dmdasnes->jacobianlocalvec) PetscCallBack("SNES DMDA local callback Jacobian", (*dmdasnes->jacobianlocalvec)(&info, Xloc, A, B, dmdasnes->jacobianlocalctx));
151ef1023bdSBarry Smith     else {
1529566063dSJacob Faibussowitsch       PetscCall(DMDAVecGetArray(dm, Xloc, &x));
153792fecdfSBarry Smith       PetscCallBack("SNES DMDA local callback Jacobian", (*dmdasnes->jacobianlocal)(&info, x, A, B, dmdasnes->jacobianlocalctx));
1549566063dSJacob Faibussowitsch       PetscCall(DMDAVecRestoreArray(dm, Xloc, &x));
1555505f7afSJunchao Zhang     }
1569566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm, &Xloc));
1572961e153SJed Brown   } else {
1582961e153SJed Brown     MatFDColoring fdcoloring;
1599566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)dm, "DMDASNES_FDCOLORING", (PetscObject *)&fdcoloring));
1602961e153SJed Brown     if (!fdcoloring) {
1612961e153SJed Brown       ISColoring coloring;
1622961e153SJed Brown 
1639566063dSJacob Faibussowitsch       PetscCall(DMCreateColoring(dm, dm->coloringtype, &coloring));
1649566063dSJacob Faibussowitsch       PetscCall(MatFDColoringCreate(B, coloring, &fdcoloring));
1652961e153SJed Brown       switch (dm->coloringtype) {
1669371c9d4SSatish Balay       case IS_COLORING_GLOBAL: PetscCall(MatFDColoringSetFunction(fdcoloring, (PetscErrorCode(*)(void))SNESComputeFunction_DMDA, dmdasnes)); break;
16798921bdaSJacob Faibussowitsch       default: SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "No support for coloring type '%s'", ISColoringTypes[dm->coloringtype]);
1682961e153SJed Brown       }
1699566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fdcoloring, ((PetscObject)dm)->prefix));
1709566063dSJacob Faibussowitsch       PetscCall(MatFDColoringSetFromOptions(fdcoloring));
1719566063dSJacob Faibussowitsch       PetscCall(MatFDColoringSetUp(B, coloring, fdcoloring));
1729566063dSJacob Faibussowitsch       PetscCall(ISColoringDestroy(&coloring));
1739566063dSJacob Faibussowitsch       PetscCall(PetscObjectCompose((PetscObject)dm, "DMDASNES_FDCOLORING", (PetscObject)fdcoloring));
1749566063dSJacob Faibussowitsch       PetscCall(PetscObjectDereference((PetscObject)fdcoloring));
1752961e153SJed Brown 
1762961e153SJed Brown       /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call
1772961e153SJed Brown        * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the
1782961e153SJed Brown        * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually
1792961e153SJed Brown        * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be
180140e18c1SBarry Smith        * taken when the PetscObjectList for the Vec inside MatFDColoring is destroyed.
1812961e153SJed Brown        */
1829566063dSJacob Faibussowitsch       PetscCall(PetscObjectDereference((PetscObject)dm));
1832961e153SJed Brown     }
1849566063dSJacob Faibussowitsch     PetscCall(MatFDColoringApply(B, fdcoloring, X, snes));
1852961e153SJed Brown   }
1862961e153SJed Brown   /* This will be redundant if the user called both, but it's too common to forget. */
18794ab13aaSBarry Smith   if (A != B) {
1889566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1899566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
1902961e153SJed Brown   }
1912961e153SJed Brown   PetscFunctionReturn(0);
1922961e153SJed Brown }
1932961e153SJed Brown 
1946cab3a1bSJed Brown /*@C
195f6dfbefdSBarry Smith    DMDASNESSetFunctionLocal - set a local residual evaluation function for use with `DMDA`
1966cab3a1bSJed Brown 
1976cab3a1bSJed Brown    Logically Collective
1986cab3a1bSJed Brown 
1994165533cSJose E. Roman    Input Parameters:
200f6dfbefdSBarry Smith +  dm - `DM` to associate callback with
201f6dfbefdSBarry Smith .  imode - `INSERT_VALUES` if local function computes owned part, `ADD_VALUES` if it contributes to ghosted part
2026cab3a1bSJed Brown .  func - local residual evaluation
2036cab3a1bSJed Brown -  ctx - optional context for local residual evaluation
2046cab3a1bSJed Brown 
2050038884cSBarry Smith    Calling sequence:
2060038884cSBarry Smith    For PetscErrorCode (*func)(DMDALocalInfo *info,void *x, void *f, void *ctx),
207f6dfbefdSBarry Smith +  info - `DMDALocalInfo` defining the subdomain to evaluate the residual on
2080038884cSBarry Smith .  x - dimensional pointer to state at which to evaluate residual (e.g. PetscScalar *x or **x or ***x)
2090038884cSBarry Smith .  f - dimensional pointer to residual, write the residual here (e.g. PetscScalar *f or **f or ***f)
2106cab3a1bSJed Brown -  ctx - optional context passed above
2116cab3a1bSJed Brown 
2126cab3a1bSJed Brown    Level: beginner
2136cab3a1bSJed Brown 
214f6dfbefdSBarry Smith .seealso: `DMDA`, `DMDASNESSetJacobianLocal()`, `DMSNESSetFunction()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
2156cab3a1bSJed Brown @*/
2169371c9d4SSatish Balay PetscErrorCode DMDASNESSetFunctionLocal(DM dm, InsertMode imode, PetscErrorCode (*func)(DMDALocalInfo *, void *, void *, void *), void *ctx) {
217942e3340SBarry Smith   DMSNES     sdm;
218942e3340SBarry Smith   DMSNES_DA *dmdasnes;
2196cab3a1bSJed Brown 
2206cab3a1bSJed Brown   PetscFunctionBegin;
2216cab3a1bSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2229566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm, &sdm));
2239566063dSJacob Faibussowitsch   PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes));
2241aa26658SKarl Rupp 
225709ac0d2SJed Brown   dmdasnes->residuallocalimode = imode;
2266cab3a1bSJed Brown   dmdasnes->residuallocal      = func;
2276cab3a1bSJed Brown   dmdasnes->residuallocalctx   = ctx;
2281aa26658SKarl Rupp 
2299566063dSJacob Faibussowitsch   PetscCall(DMSNESSetFunction(dm, SNESComputeFunction_DMDA, dmdasnes));
23022c6f798SBarry Smith   if (!sdm->ops->computejacobian) { /* Call us for the Jacobian too, can be overridden by the user. */
2319566063dSJacob Faibussowitsch     PetscCall(DMSNESSetJacobian(dm, SNESComputeJacobian_DMDA, dmdasnes));
2322961e153SJed Brown   }
2332961e153SJed Brown   PetscFunctionReturn(0);
2342961e153SJed Brown }
2352961e153SJed Brown 
2362961e153SJed Brown /*@C
237f6dfbefdSBarry Smith    DMDASNESSetFunctionLocalVec - set a local residual evaluation function that operates on a local vector for `DMDA`
2385505f7afSJunchao Zhang 
2395505f7afSJunchao Zhang    Logically Collective
2405505f7afSJunchao Zhang 
2415505f7afSJunchao Zhang    Input Parameters:
242f6dfbefdSBarry Smith +  dm - `DM` to associate callback with
243f6dfbefdSBarry Smith .  imode - `INSERT_VALUES` if local function computes owned part, `ADD_VALUES` if it contributes to ghosted part
2445505f7afSJunchao Zhang .  func - local residual evaluation
2455505f7afSJunchao Zhang -  ctx - optional context for local residual evaluation
2465505f7afSJunchao Zhang 
2475505f7afSJunchao Zhang    Calling sequence:
2485505f7afSJunchao Zhang    For PetscErrorCode (*func)(DMDALocalInfo *info,Vec x, Vec f, void *ctx),
249f6dfbefdSBarry Smith +  info - `DMDALocalInfo` defining the subdomain to evaluate the residual on
2505505f7afSJunchao Zhang .  x - state vector at which to evaluate residual
2515505f7afSJunchao Zhang .  f - residual vector
2525505f7afSJunchao Zhang -  ctx - optional context passed above
2535505f7afSJunchao Zhang 
2545505f7afSJunchao Zhang    Level: beginner
2555505f7afSJunchao Zhang 
256f6dfbefdSBarry Smith .seealso: `DMDA`, `DMDASNESSetFunctionLocal()`, `DMDASNESSetJacobianLocalVec()`, `DMSNESSetFunction()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
2575505f7afSJunchao Zhang @*/
2589371c9d4SSatish Balay PetscErrorCode DMDASNESSetFunctionLocalVec(DM dm, InsertMode imode, PetscErrorCode (*func)(DMDALocalInfo *, Vec, Vec, void *), void *ctx) {
2595505f7afSJunchao Zhang   DMSNES     sdm;
2605505f7afSJunchao Zhang   DMSNES_DA *dmdasnes;
2615505f7afSJunchao Zhang 
2625505f7afSJunchao Zhang   PetscFunctionBegin;
2635505f7afSJunchao Zhang   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2645505f7afSJunchao Zhang   PetscCall(DMGetDMSNESWrite(dm, &sdm));
2655505f7afSJunchao Zhang   PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes));
2665505f7afSJunchao Zhang 
2675505f7afSJunchao Zhang   dmdasnes->residuallocalimode = imode;
2685505f7afSJunchao Zhang   dmdasnes->residuallocalvec   = func;
2695505f7afSJunchao Zhang   dmdasnes->residuallocalctx   = ctx;
2705505f7afSJunchao Zhang 
2715505f7afSJunchao Zhang   PetscCall(DMSNESSetFunction(dm, SNESComputeFunction_DMDA, dmdasnes));
2725505f7afSJunchao Zhang   if (!sdm->ops->computejacobian) { /* Call us for the Jacobian too, can be overridden by the user. */
2735505f7afSJunchao Zhang     PetscCall(DMSNESSetJacobian(dm, SNESComputeJacobian_DMDA, dmdasnes));
2745505f7afSJunchao Zhang   }
2755505f7afSJunchao Zhang   PetscFunctionReturn(0);
2765505f7afSJunchao Zhang }
2775505f7afSJunchao Zhang 
2785505f7afSJunchao Zhang /*@C
279f6dfbefdSBarry Smith    DMDASNESSetJacobianLocal - set a local Jacobian evaluation function for use with `DMDA`
2802961e153SJed Brown 
2812961e153SJed Brown    Logically Collective
2822961e153SJed Brown 
2834165533cSJose E. Roman    Input Parameters:
284f6dfbefdSBarry Smith +  dm - `DM` to associate callback with
2856b7eb5ceSMatthew G. Knepley .  func - local Jacobian evaluation
2866b7eb5ceSMatthew G. Knepley -  ctx - optional context for local Jacobian evaluation
2872961e153SJed Brown 
2880038884cSBarry Smith    Calling sequence:
2890038884cSBarry Smith    For PetscErrorCode (*func)(DMDALocalInfo *info,void *x,Mat J,Mat M,void *ctx),
290f6dfbefdSBarry Smith +  info - `DMDALocalInfo` defining the subdomain to evaluate the Jacobian at
2910038884cSBarry Smith .  x - dimensional pointer to state at which to evaluate Jacobian (e.g. PetscScalar *x or **x or ***x)
2926b7eb5ceSMatthew G. Knepley .  J - Mat object for the Jacobian
2936b7eb5ceSMatthew G. Knepley .  M - Mat object for the Jacobian preconditioner matrix
2942961e153SJed Brown -  ctx - optional context passed above
2952961e153SJed Brown 
2962961e153SJed Brown    Level: beginner
2972961e153SJed Brown 
298f6dfbefdSBarry Smith .seealso: `DMDA`, `DMDASNESSetFunctionLocal()`, `DMSNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
2992961e153SJed Brown @*/
3009371c9d4SSatish Balay PetscErrorCode DMDASNESSetJacobianLocal(DM dm, PetscErrorCode (*func)(DMDALocalInfo *, void *, Mat, Mat, void *), void *ctx) {
301942e3340SBarry Smith   DMSNES     sdm;
302942e3340SBarry Smith   DMSNES_DA *dmdasnes;
3032961e153SJed Brown 
3042961e153SJed Brown   PetscFunctionBegin;
3052961e153SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3069566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm, &sdm));
3079566063dSJacob Faibussowitsch   PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes));
3081aa26658SKarl Rupp 
3092961e153SJed Brown   dmdasnes->jacobianlocal    = func;
3102961e153SJed Brown   dmdasnes->jacobianlocalctx = ctx;
3111aa26658SKarl Rupp 
3129566063dSJacob Faibussowitsch   PetscCall(DMSNESSetJacobian(dm, SNESComputeJacobian_DMDA, dmdasnes));
3136cab3a1bSJed Brown   PetscFunctionReturn(0);
3146cab3a1bSJed Brown }
3152a4ee8f2SPeter Brune 
3162a4ee8f2SPeter Brune /*@C
317f6dfbefdSBarry Smith    DMDASNESSetJacobianLocalVec - set a local Jacobian evaluation function that operates on a local vector with `DMDA`
3185505f7afSJunchao Zhang 
3195505f7afSJunchao Zhang    Logically Collective
3205505f7afSJunchao Zhang 
3215505f7afSJunchao Zhang    Input Parameters:
322f6dfbefdSBarry Smith +  dm - `DM` to associate callback with
3235505f7afSJunchao Zhang .  func - local Jacobian evaluation
3245505f7afSJunchao Zhang -  ctx - optional context for local Jacobian evaluation
3255505f7afSJunchao Zhang 
3265505f7afSJunchao Zhang    Calling sequence:
3275505f7afSJunchao Zhang    For PetscErrorCode (*func)(DMDALocalInfo *info,Vec x,Mat J,Mat M,void *ctx),
328f6dfbefdSBarry Smith +  info - `DMDALocalInfo` defining the subdomain to evaluate the Jacobian at
3295505f7afSJunchao Zhang .  x - state vector at which to evaluate Jacobian
3305505f7afSJunchao Zhang .  J - Mat object for the Jacobian
3315505f7afSJunchao Zhang .  M - Mat object for the Jacobian preconditioner matrix
3325505f7afSJunchao Zhang -  ctx - optional context passed above
3335505f7afSJunchao Zhang 
3345505f7afSJunchao Zhang    Level: beginner
3355505f7afSJunchao Zhang 
336f6dfbefdSBarry Smith .seealso: `DMDA`, `DMDASNESSetJacobianLocal()`, `DMDASNESSetFunctionLocalVec()`, `DMSNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
3375505f7afSJunchao Zhang @*/
3389371c9d4SSatish Balay PetscErrorCode DMDASNESSetJacobianLocalVec(DM dm, PetscErrorCode (*func)(DMDALocalInfo *, Vec, Mat, Mat, void *), void *ctx) {
3395505f7afSJunchao Zhang   DMSNES     sdm;
3405505f7afSJunchao Zhang   DMSNES_DA *dmdasnes;
3415505f7afSJunchao Zhang 
3425505f7afSJunchao Zhang   PetscFunctionBegin;
3435505f7afSJunchao Zhang   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3445505f7afSJunchao Zhang   PetscCall(DMGetDMSNESWrite(dm, &sdm));
3455505f7afSJunchao Zhang   PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes));
3465505f7afSJunchao Zhang 
3475505f7afSJunchao Zhang   dmdasnes->jacobianlocalvec = func;
3485505f7afSJunchao Zhang   dmdasnes->jacobianlocalctx = ctx;
3495505f7afSJunchao Zhang 
3505505f7afSJunchao Zhang   PetscCall(DMSNESSetJacobian(dm, SNESComputeJacobian_DMDA, dmdasnes));
3515505f7afSJunchao Zhang   PetscFunctionReturn(0);
3525505f7afSJunchao Zhang }
3535505f7afSJunchao Zhang 
3545505f7afSJunchao Zhang /*@C
355f6dfbefdSBarry Smith    DMDASNESSetObjectiveLocal - set a local residual evaluation function to used with a `DMDA`
3562a4ee8f2SPeter Brune 
3572a4ee8f2SPeter Brune    Logically Collective
3582a4ee8f2SPeter Brune 
3594165533cSJose E. Roman    Input Parameters:
360f6dfbefdSBarry Smith +  dm - `DM` to associate callback with
3612a4ee8f2SPeter Brune .  func - local objective evaluation
3622a4ee8f2SPeter Brune -  ctx - optional context for local residual evaluation
3632a4ee8f2SPeter Brune 
3642a4ee8f2SPeter Brune    Calling sequence for func:
365f6dfbefdSBarry Smith +  info - `DMDALocalInfo` defining the subdomain to evaluate the residual on
3662a4ee8f2SPeter Brune .  x - dimensional pointer to state at which to evaluate residual
3672a4ee8f2SPeter Brune .  ob - eventual objective value
3682a4ee8f2SPeter Brune -  ctx - optional context passed above
3692a4ee8f2SPeter Brune 
3702a4ee8f2SPeter Brune    Level: beginner
3712a4ee8f2SPeter Brune 
372f6dfbefdSBarry Smith .seealso: `DMDA`, `DMSNESSetFunction()`, `DMDASNESSetJacobianLocal()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
3732a4ee8f2SPeter Brune @*/
3749371c9d4SSatish Balay PetscErrorCode DMDASNESSetObjectiveLocal(DM dm, DMDASNESObjective func, void *ctx) {
375942e3340SBarry Smith   DMSNES     sdm;
376942e3340SBarry Smith   DMSNES_DA *dmdasnes;
3772a4ee8f2SPeter Brune 
3782a4ee8f2SPeter Brune   PetscFunctionBegin;
3792a4ee8f2SPeter Brune   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3809566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm, &sdm));
3819566063dSJacob Faibussowitsch   PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes));
3821aa26658SKarl Rupp 
3832a4ee8f2SPeter Brune   dmdasnes->objectivelocal    = func;
3842a4ee8f2SPeter Brune   dmdasnes->objectivelocalctx = ctx;
3851aa26658SKarl Rupp 
3869566063dSJacob Faibussowitsch   PetscCall(DMSNESSetObjective(dm, SNESComputeObjective_DMDA, dmdasnes));
3872a4ee8f2SPeter Brune   PetscFunctionReturn(0);
3882a4ee8f2SPeter Brune }
389dcad5f1cSBarry Smith 
3905505f7afSJunchao Zhang /*@C
391f6dfbefdSBarry Smith    DMDASNESSetObjectiveLocal - set a local residual evaluation function that operates on a local vector with `DMDA`
3925505f7afSJunchao Zhang 
3935505f7afSJunchao Zhang    Logically Collective
3945505f7afSJunchao Zhang 
3955505f7afSJunchao Zhang    Input Parameters:
396f6dfbefdSBarry Smith +  dm - `DM` to associate callback with
3975505f7afSJunchao Zhang .  func - local objective evaluation
3985505f7afSJunchao Zhang -  ctx - optional context for local residual evaluation
3995505f7afSJunchao Zhang 
4005505f7afSJunchao Zhang    Calling sequence
4015505f7afSJunchao Zhang    For PetscErrorCode (*func)(DMDALocalInfo *info,Vec x,PetscReal *ob,void *ctx);
402f6dfbefdSBarry Smith +  info - `DMDALocalInfo` defining the subdomain to evaluate the residual on
4035505f7afSJunchao Zhang .  x - state vector at which to evaluate residual
4045505f7afSJunchao Zhang .  ob - eventual objective value
4055505f7afSJunchao Zhang -  ctx - optional context passed above
4065505f7afSJunchao Zhang 
4075505f7afSJunchao Zhang    Level: beginner
4085505f7afSJunchao Zhang 
409f6dfbefdSBarry Smith .seealso: `DMDA`, `DMDASNESSetObjectiveLocal()`, `DMSNESSetFunction()`, `DMDASNESSetJacobianLocalVec()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
4105505f7afSJunchao Zhang @*/
4119371c9d4SSatish Balay PetscErrorCode DMDASNESSetObjectiveLocalVec(DM dm, DMDASNESObjectiveVec func, void *ctx) {
4125505f7afSJunchao Zhang   DMSNES     sdm;
4135505f7afSJunchao Zhang   DMSNES_DA *dmdasnes;
4145505f7afSJunchao Zhang 
4155505f7afSJunchao Zhang   PetscFunctionBegin;
4165505f7afSJunchao Zhang   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4175505f7afSJunchao Zhang   PetscCall(DMGetDMSNESWrite(dm, &sdm));
4185505f7afSJunchao Zhang   PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes));
4195505f7afSJunchao Zhang 
4205505f7afSJunchao Zhang   dmdasnes->objectivelocalvec = func;
4215505f7afSJunchao Zhang   dmdasnes->objectivelocalctx = ctx;
4225505f7afSJunchao Zhang 
4235505f7afSJunchao Zhang   PetscCall(DMSNESSetObjective(dm, SNESComputeObjective_DMDA, dmdasnes));
4245505f7afSJunchao Zhang   PetscFunctionReturn(0);
4255505f7afSJunchao Zhang }
4265505f7afSJunchao Zhang 
4279371c9d4SSatish Balay static PetscErrorCode SNESComputePicard_DMDA(SNES snes, Vec X, Vec F, void *ctx) {
428dcad5f1cSBarry Smith   DM            dm;
429942e3340SBarry Smith   DMSNES_DA    *dmdasnes = (DMSNES_DA *)ctx;
430dcad5f1cSBarry Smith   DMDALocalInfo info;
431dcad5f1cSBarry Smith   Vec           Xloc;
432dcad5f1cSBarry Smith   void         *x, *f;
433dcad5f1cSBarry Smith 
434dcad5f1cSBarry Smith   PetscFunctionBegin;
435dcad5f1cSBarry Smith   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
436dcad5f1cSBarry Smith   PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
437dcad5f1cSBarry Smith   PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
43828b400f6SJacob Faibussowitsch   PetscCheck(dmdasnes->rhsplocal, PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Corrupt context");
4399566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
4409566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &Xloc));
4419566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
4429566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
4439566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm, &info));
4449566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(dm, Xloc, &x));
445dcad5f1cSBarry Smith   switch (dmdasnes->residuallocalimode) {
446dcad5f1cSBarry Smith   case INSERT_VALUES: {
4479566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dm, F, &f));
448792fecdfSBarry Smith     PetscCallBack("SNES Picard DMDA local callback function", (*dmdasnes->rhsplocal)(&info, x, f, dmdasnes->picardlocalctx));
4499566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dm, F, &f));
450dcad5f1cSBarry Smith   } break;
451dcad5f1cSBarry Smith   case ADD_VALUES: {
452dcad5f1cSBarry Smith     Vec Floc;
4539566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm, &Floc));
4549566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(Floc));
4559566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dm, Floc, &f));
456792fecdfSBarry Smith     PetscCallBack("SNES Picard DMDA local callback function", (*dmdasnes->rhsplocal)(&info, x, f, dmdasnes->picardlocalctx));
4579566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dm, Floc, &f));
4589566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(F));
4599566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalBegin(dm, Floc, ADD_VALUES, F));
4609566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalEnd(dm, Floc, ADD_VALUES, F));
4619566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm, &Floc));
462dcad5f1cSBarry Smith   } break;
46398921bdaSJacob Faibussowitsch   default: SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_INCOMP, "Cannot use imode=%d", (int)dmdasnes->residuallocalimode);
464dcad5f1cSBarry Smith   }
4659566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(dm, Xloc, &x));
4669566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &Xloc));
467dcad5f1cSBarry Smith   PetscFunctionReturn(0);
468dcad5f1cSBarry Smith }
469dcad5f1cSBarry Smith 
4709371c9d4SSatish Balay static PetscErrorCode SNESComputePicardJacobian_DMDA(SNES snes, Vec X, Mat A, Mat B, void *ctx) {
471dcad5f1cSBarry Smith   DM            dm;
472942e3340SBarry Smith   DMSNES_DA    *dmdasnes = (DMSNES_DA *)ctx;
473dcad5f1cSBarry Smith   DMDALocalInfo info;
474dcad5f1cSBarry Smith   Vec           Xloc;
475dcad5f1cSBarry Smith   void         *x;
476dcad5f1cSBarry Smith 
477dcad5f1cSBarry Smith   PetscFunctionBegin;
47828b400f6SJacob Faibussowitsch   PetscCheck(dmdasnes->jacobianplocal, PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "Corrupt context");
4799566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
480dcad5f1cSBarry Smith 
4819566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &Xloc));
4829566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc));
4839566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc));
4849566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm, &info));
4859566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(dm, Xloc, &x));
486792fecdfSBarry Smith   PetscCallBack("SNES Picard DMDA local callback Jacobian", (*dmdasnes->jacobianplocal)(&info, x, A, B, dmdasnes->picardlocalctx));
4879566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(dm, Xloc, &x));
4889566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &Xloc));
48994ab13aaSBarry Smith   if (A != B) {
4909566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
4919566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
492dcad5f1cSBarry Smith   }
493dcad5f1cSBarry Smith   PetscFunctionReturn(0);
494dcad5f1cSBarry Smith }
495dcad5f1cSBarry Smith 
496dcad5f1cSBarry Smith /*@C
497f6dfbefdSBarry Smith    DMDASNESSetPicardLocal - set a local right hand side and matrix evaluation function for Picard iteration with `DMDA`
498dcad5f1cSBarry Smith 
499dcad5f1cSBarry Smith    Logically Collective
500dcad5f1cSBarry Smith 
5014165533cSJose E. Roman    Input Parameters:
502f6dfbefdSBarry Smith +  dm - `DM` to associate callback with
503f6dfbefdSBarry Smith .  imode - `INSERT_VALUES` if local function computes owned part, `ADD_VALUES` if it contributes to ghosted part
504dcad5f1cSBarry Smith .  func - local residual evaluation
505dcad5f1cSBarry Smith -  ctx - optional context for local residual evaluation
506dcad5f1cSBarry Smith 
507dcad5f1cSBarry Smith    Calling sequence for func:
508f6dfbefdSBarry Smith +  info - `DMDALocalInfo` defining the subdomain to evaluate the residual on
509dcad5f1cSBarry Smith .  x - dimensional pointer to state at which to evaluate residual
510dcad5f1cSBarry Smith .  f - dimensional pointer to residual, write the residual here
511dcad5f1cSBarry Smith -  ctx - optional context passed above
512dcad5f1cSBarry Smith 
513f6dfbefdSBarry Smith    Note:
514f6dfbefdSBarry Smith     The user must use `SNESSetFunction`(snes,NULL,`SNESPicardComputeFunction`,&user));
515dcad5f1cSBarry Smith     in their code before calling this routine.
516dcad5f1cSBarry Smith 
517dcad5f1cSBarry Smith    Level: beginner
518dcad5f1cSBarry Smith 
519f6dfbefdSBarry Smith .seealso: `DMDA`, `DMSNESSetFunction()`, `DMDASNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`
520dcad5f1cSBarry Smith @*/
5219371c9d4SSatish Balay PetscErrorCode DMDASNESSetPicardLocal(DM dm, InsertMode imode, PetscErrorCode (*func)(DMDALocalInfo *, void *, void *, void *), PetscErrorCode (*jac)(DMDALocalInfo *, void *, Mat, Mat, void *), void *ctx) {
522942e3340SBarry Smith   DMSNES     sdm;
523942e3340SBarry Smith   DMSNES_DA *dmdasnes;
524dcad5f1cSBarry Smith 
525dcad5f1cSBarry Smith   PetscFunctionBegin;
526dcad5f1cSBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5279566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNESWrite(dm, &sdm));
5289566063dSJacob Faibussowitsch   PetscCall(DMDASNESGetContext(dm, sdm, &dmdasnes));
5291aa26658SKarl Rupp 
530dcad5f1cSBarry Smith   dmdasnes->residuallocalimode = imode;
531dcad5f1cSBarry Smith   dmdasnes->rhsplocal          = func;
532dcad5f1cSBarry Smith   dmdasnes->jacobianplocal     = jac;
533dcad5f1cSBarry Smith   dmdasnes->picardlocalctx     = ctx;
5341aa26658SKarl Rupp 
5359566063dSJacob Faibussowitsch   PetscCall(DMSNESSetPicard(dm, SNESComputePicard_DMDA, SNESComputePicardJacobian_DMDA, dmdasnes));
5369566063dSJacob Faibussowitsch   PetscCall(DMSNESSetMFFunction(dm, SNESComputeFunction_DMDA, dmdasnes));
537dcad5f1cSBarry Smith   PetscFunctionReturn(0);
538dcad5f1cSBarry Smith }
539