xref: /petsc/src/snes/utils/dmdasnes.c (revision 6427ad5be45d8802ffd76d36c8ddb5eadd18bb26)
16cab3a1bSJed Brown #include <petscdmda.h>          /*I "petscdmda.h" I*/
2b45d2f2cSJed Brown #include <petsc-private/dmimpl.h>
3b45d2f2cSJed Brown #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*);
86cab3a1bSJed Brown   PetscErrorCode (*jacobianlocal)(DMDALocalInfo*,void*,Mat,Mat,MatStructure*,void*);
96cab3a1bSJed Brown   void *residuallocalctx;
106cab3a1bSJed Brown   void *jacobianlocalctx;
11709ac0d2SJed Brown   InsertMode residuallocalimode;
126cab3a1bSJed Brown } DM_DA_SNES;
136cab3a1bSJed Brown 
146cab3a1bSJed Brown #undef __FUNCT__
156cab3a1bSJed Brown #define __FUNCT__ "SNESDMDestroy_DMDA"
166cab3a1bSJed Brown static PetscErrorCode SNESDMDestroy_DMDA(SNESDM sdm)
176cab3a1bSJed Brown {
186cab3a1bSJed Brown   PetscErrorCode ierr;
196cab3a1bSJed Brown 
206cab3a1bSJed Brown   PetscFunctionBegin;
216cab3a1bSJed Brown   ierr = PetscFree(sdm->data);CHKERRQ(ierr);
226cab3a1bSJed Brown   PetscFunctionReturn(0);
236cab3a1bSJed Brown }
246cab3a1bSJed Brown 
256cab3a1bSJed Brown #undef __FUNCT__
266cab3a1bSJed Brown #define __FUNCT__ "DMDASNESGetContext"
276cab3a1bSJed Brown static PetscErrorCode DMDASNESGetContext(DM dm,SNESDM sdm,DM_DA_SNES **dmdasnes)
286cab3a1bSJed Brown {
296cab3a1bSJed Brown   PetscErrorCode ierr;
306cab3a1bSJed Brown 
316cab3a1bSJed Brown   PetscFunctionBegin;
32b8d0cc33SJed Brown   *dmdasnes = PETSC_NULL;
336cab3a1bSJed Brown   if (!sdm->data) {
346cab3a1bSJed Brown     ierr = PetscNewLog(dm,DM_DA_SNES,&sdm->data);CHKERRQ(ierr);
356cab3a1bSJed Brown     sdm->destroy = SNESDMDestroy_DMDA;
366cab3a1bSJed Brown   }
376cab3a1bSJed Brown   *dmdasnes = (DM_DA_SNES*)sdm->data;
386cab3a1bSJed Brown   PetscFunctionReturn(0);
396cab3a1bSJed Brown }
406cab3a1bSJed Brown 
416cab3a1bSJed Brown #undef __FUNCT__
426cab3a1bSJed Brown #define __FUNCT__ "SNESComputeFunction_DMDA"
436cab3a1bSJed Brown /*
446cab3a1bSJed Brown   This function should eventually replace:
456cab3a1bSJed Brown     DMDAComputeFunction() and DMDAComputeFunction1()
466cab3a1bSJed Brown  */
476cab3a1bSJed Brown static PetscErrorCode SNESComputeFunction_DMDA(SNES snes,Vec X,Vec F,void *ctx)
486cab3a1bSJed Brown {
496cab3a1bSJed Brown   PetscErrorCode ierr;
506cab3a1bSJed Brown   DM             dm;
516cab3a1bSJed Brown   DM_DA_SNES     *dmdasnes = (DM_DA_SNES*)ctx;
526cab3a1bSJed Brown   DMDALocalInfo  info;
536cab3a1bSJed Brown   Vec            Xloc;
546cab3a1bSJed Brown   void           *x,*f;
556cab3a1bSJed Brown 
566cab3a1bSJed Brown   PetscFunctionBegin;
572961e153SJed Brown   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
582961e153SJed Brown   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
592961e153SJed Brown   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
606cab3a1bSJed Brown   if (!dmdasnes->residuallocal) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Corrupt context");
616cab3a1bSJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
626cab3a1bSJed Brown   ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
636cab3a1bSJed Brown   ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
646cab3a1bSJed Brown   ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
656cab3a1bSJed Brown   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
666cab3a1bSJed Brown   ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
67709ac0d2SJed Brown   switch (dmdasnes->residuallocalimode) {
68709ac0d2SJed Brown   case INSERT_VALUES: {
696cab3a1bSJed Brown     ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr);
706cab3a1bSJed Brown     CHKMEMQ;
716cab3a1bSJed Brown     ierr = (*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx);CHKERRQ(ierr);
726cab3a1bSJed Brown     CHKMEMQ;
736cab3a1bSJed Brown     ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr);
74709ac0d2SJed Brown   } break;
75709ac0d2SJed Brown   case ADD_VALUES: {
76709ac0d2SJed Brown     Vec Floc;
77709ac0d2SJed Brown     ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr);
78709ac0d2SJed Brown     ierr = VecZeroEntries(Floc);CHKERRQ(ierr);
79709ac0d2SJed Brown     ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr);
80709ac0d2SJed Brown     CHKMEMQ;
81709ac0d2SJed Brown     ierr = (*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx);CHKERRQ(ierr);
82709ac0d2SJed Brown     CHKMEMQ;
83709ac0d2SJed Brown     ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr);
84709ac0d2SJed Brown     ierr = VecZeroEntries(F);CHKERRQ(ierr);
85709ac0d2SJed Brown     ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
86709ac0d2SJed Brown     ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
87709ac0d2SJed Brown     ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr);
88709ac0d2SJed Brown   } break;
89709ac0d2SJed Brown   default: SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdasnes->residuallocalimode);
90709ac0d2SJed Brown   }
91709ac0d2SJed Brown   ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
926cab3a1bSJed Brown   ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
936cab3a1bSJed Brown   PetscFunctionReturn(0);
946cab3a1bSJed Brown }
956cab3a1bSJed Brown 
966cab3a1bSJed Brown #undef __FUNCT__
97*6427ad5bSPeter Brune #define __FUNCT__ "SNESComputeLocalBlockFunction_DMDA"
98*6427ad5bSPeter Brune /*
99*6427ad5bSPeter Brune   This function should eventually replace:
100*6427ad5bSPeter Brune     DMDAComputeFunction() and DMDAComputeFunction1()
101*6427ad5bSPeter Brune  */
102*6427ad5bSPeter Brune PetscErrorCode SNESComputeLocalBlockFunction_DMDA(SNES snes,Vec X,Vec F,void *dmctx)
103*6427ad5bSPeter Brune {
104*6427ad5bSPeter Brune   PetscErrorCode ierr;
105*6427ad5bSPeter Brune   DM             dm = (DM)dmctx; /* the global DMDA context */
106*6427ad5bSPeter Brune   DMDALocalInfo  info;
107*6427ad5bSPeter Brune   void           *x,*f;
108*6427ad5bSPeter Brune   DM_DA_SNES     *dmdasnes;
109*6427ad5bSPeter Brune   SNESDM         sdm;
110*6427ad5bSPeter Brune 
111*6427ad5bSPeter Brune   PetscFunctionBegin;
112*6427ad5bSPeter Brune 
113*6427ad5bSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
114*6427ad5bSPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
115*6427ad5bSPeter Brune   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
116*6427ad5bSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,4);
117*6427ad5bSPeter Brune 
118*6427ad5bSPeter Brune   if (!dmdasnes->residuallocal) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Corrupt context");
119*6427ad5bSPeter Brune   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
120*6427ad5bSPeter Brune   ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr);
121*6427ad5bSPeter Brune   ierr = DMDAVecGetArray(dm,X,&x);CHKERRQ(ierr);
122*6427ad5bSPeter Brune   ierr = DMDAGetLocalBlockInfo(dm,&info);CHKERRQ(ierr);
123*6427ad5bSPeter Brune   switch (dmdasnes->residuallocalimode) {
124*6427ad5bSPeter Brune   case INSERT_VALUES: {
125*6427ad5bSPeter Brune     ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr);
126*6427ad5bSPeter Brune     CHKMEMQ;
127*6427ad5bSPeter Brune     ierr = (*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx);CHKERRQ(ierr);
128*6427ad5bSPeter Brune     CHKMEMQ;
129*6427ad5bSPeter Brune     ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr);
130*6427ad5bSPeter Brune   } break;
131*6427ad5bSPeter Brune   case ADD_VALUES: {
132*6427ad5bSPeter Brune     Vec Floc;
133*6427ad5bSPeter Brune     ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr);
134*6427ad5bSPeter Brune     ierr = VecZeroEntries(Floc);CHKERRQ(ierr);
135*6427ad5bSPeter Brune     ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr);
136*6427ad5bSPeter Brune     CHKMEMQ;
137*6427ad5bSPeter Brune     ierr = (*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx);CHKERRQ(ierr);
138*6427ad5bSPeter Brune     CHKMEMQ;
139*6427ad5bSPeter Brune     ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr);
140*6427ad5bSPeter Brune     ierr = VecZeroEntries(F);CHKERRQ(ierr);
141*6427ad5bSPeter Brune     ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
142*6427ad5bSPeter Brune     ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
143*6427ad5bSPeter Brune     ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr);
144*6427ad5bSPeter Brune   } break;
145*6427ad5bSPeter Brune   default: SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdasnes->residuallocalimode);
146*6427ad5bSPeter Brune   }
147*6427ad5bSPeter Brune   ierr = DMDAVecRestoreArray(dm,X,&x);CHKERRQ(ierr);
148*6427ad5bSPeter Brune   PetscFunctionReturn(0);
149*6427ad5bSPeter Brune }
150*6427ad5bSPeter Brune 
151*6427ad5bSPeter Brune #undef __FUNCT__
1522961e153SJed Brown #define __FUNCT__ "SNESComputeJacobian_DMDA"
1532961e153SJed Brown /*
1542961e153SJed Brown   This function should eventually replace:
1552961e153SJed Brown     DMComputeJacobian() and DMDAComputeJacobian1()
1562961e153SJed Brown  */
1572961e153SJed Brown static PetscErrorCode SNESComputeJacobian_DMDA(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *ctx)
1582961e153SJed Brown {
1592961e153SJed Brown   PetscErrorCode ierr;
1602961e153SJed Brown   DM             dm;
1612961e153SJed Brown   DM_DA_SNES     *dmdasnes = (DM_DA_SNES*)ctx;
1622961e153SJed Brown   DMDALocalInfo  info;
1632961e153SJed Brown   Vec            Xloc;
1642961e153SJed Brown   void           *x;
1652961e153SJed Brown 
1662961e153SJed Brown   PetscFunctionBegin;
1672961e153SJed Brown   if (!dmdasnes->residuallocal) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Corrupt context");
1682961e153SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
1692961e153SJed Brown 
1702961e153SJed Brown   if (dmdasnes->jacobianlocal) {
1712961e153SJed Brown     ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
1722961e153SJed Brown     ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
1732961e153SJed Brown     ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
1742961e153SJed Brown     ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
1752961e153SJed Brown     ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
1762961e153SJed Brown     CHKMEMQ;
1772961e153SJed Brown     ierr = (*dmdasnes->jacobianlocal)(&info,x,*A,*B,mstr,dmdasnes->jacobianlocalctx);CHKERRQ(ierr);
1782961e153SJed Brown     CHKMEMQ;
1792961e153SJed Brown     ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
1802961e153SJed Brown     ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
1812961e153SJed Brown   } else {
1822961e153SJed Brown     MatFDColoring fdcoloring;
1832961e153SJed Brown     ierr = PetscObjectQuery((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject*)&fdcoloring);CHKERRQ(ierr);
1842961e153SJed Brown     if (!fdcoloring) {
1852961e153SJed Brown       ISColoring     coloring;
1862961e153SJed Brown 
1872961e153SJed Brown       ierr = DMCreateColoring(dm,dm->coloringtype,MATAIJ,&coloring);CHKERRQ(ierr);
1882961e153SJed Brown       ierr = MatFDColoringCreate(*B,coloring,&fdcoloring);CHKERRQ(ierr);
1892961e153SJed Brown       ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr);
1902961e153SJed Brown       switch (dm->coloringtype) {
1912961e153SJed Brown       case IS_COLORING_GLOBAL:
1922961e153SJed Brown         ierr = MatFDColoringSetFunction(fdcoloring,(PetscErrorCode(*)(void))SNESComputeFunction_DMDA,dmdasnes);CHKERRQ(ierr);
1932961e153SJed Brown         break;
1942961e153SJed Brown       default: SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"No support for coloring type '%s'",ISColoringTypes[dm->coloringtype]);
1952961e153SJed Brown       }
1962961e153SJed Brown       ierr = PetscObjectSetOptionsPrefix((PetscObject)fdcoloring,((PetscObject)dm)->prefix);CHKERRQ(ierr);
1972961e153SJed Brown       ierr = MatFDColoringSetFromOptions(fdcoloring);CHKERRQ(ierr);
1982961e153SJed Brown       ierr = PetscObjectCompose((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject)fdcoloring);CHKERRQ(ierr);
1992961e153SJed Brown       ierr = PetscObjectDereference((PetscObject)fdcoloring);CHKERRQ(ierr);
2002961e153SJed Brown 
2012961e153SJed Brown       /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call
2022961e153SJed Brown        * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the
2032961e153SJed Brown        * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually
2042961e153SJed Brown        * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be
2052961e153SJed Brown        * taken when the PetscOList for the Vec inside MatFDColoring is destroyed.
2062961e153SJed Brown        */
2072961e153SJed Brown       ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
2082961e153SJed Brown     }
2092961e153SJed Brown     *mstr = SAME_NONZERO_PATTERN;
2102961e153SJed Brown     ierr = MatFDColoringApply(*B,fdcoloring,X,mstr,snes);CHKERRQ(ierr);
2112961e153SJed Brown   }
2122961e153SJed Brown   /* This will be redundant if the user called both, but it's too common to forget. */
2132961e153SJed Brown   if (*A != *B) {
2142961e153SJed Brown     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2152961e153SJed Brown     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2162961e153SJed Brown   }
2172961e153SJed Brown   PetscFunctionReturn(0);
2182961e153SJed Brown }
2192961e153SJed Brown 
2202961e153SJed Brown #undef __FUNCT__
221*6427ad5bSPeter Brune #define __FUNCT__ "SNESComputeJacobian_DMDA"
222*6427ad5bSPeter Brune /*
223*6427ad5bSPeter Brune   This function should eventually replace:
224*6427ad5bSPeter Brune     DMComputeJacobian() and DMDAComputeJacobian1()
225*6427ad5bSPeter Brune  */
226*6427ad5bSPeter Brune PetscErrorCode SNESComputeLocalBlockJacobian_DMDA(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *dmctx)
227*6427ad5bSPeter Brune {
228*6427ad5bSPeter Brune   PetscErrorCode ierr;
229*6427ad5bSPeter Brune   DM             dm = (DM)dmctx; /* the global DMDA context */
230*6427ad5bSPeter Brune   DM_DA_SNES     *dmdasnes;
231*6427ad5bSPeter Brune   SNESDM         sdm;
232*6427ad5bSPeter Brune   DMDALocalInfo  info;
233*6427ad5bSPeter Brune   void           *x;
234*6427ad5bSPeter Brune 
235*6427ad5bSPeter Brune   PetscFunctionBegin;
236*6427ad5bSPeter Brune   if (!dmdasnes->residuallocal) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Corrupt context");
237*6427ad5bSPeter Brune   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
238*6427ad5bSPeter Brune   ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr);
239*6427ad5bSPeter Brune   if (dmdasnes->jacobianlocal) {
240*6427ad5bSPeter Brune     ierr = DMDAVecGetArray(dm,X,&x);CHKERRQ(ierr);
241*6427ad5bSPeter Brune     CHKMEMQ;
242*6427ad5bSPeter Brune     ierr = (*dmdasnes->jacobianlocal)(&info,x,*A,*B,mstr,dmdasnes->jacobianlocalctx);CHKERRQ(ierr);
243*6427ad5bSPeter Brune     CHKMEMQ;
244*6427ad5bSPeter Brune     ierr = DMDAVecRestoreArray(dm,X,&x);CHKERRQ(ierr);
245*6427ad5bSPeter Brune   }
246*6427ad5bSPeter Brune   /* This will be redundant if the user called both, but it's too common to forget. */
247*6427ad5bSPeter Brune   if (*A != *B) {
248*6427ad5bSPeter Brune     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
249*6427ad5bSPeter Brune     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
250*6427ad5bSPeter Brune   }
251*6427ad5bSPeter Brune   PetscFunctionReturn(0);
252*6427ad5bSPeter Brune }
253*6427ad5bSPeter Brune 
254*6427ad5bSPeter Brune 
255*6427ad5bSPeter Brune #undef __FUNCT__
2566cab3a1bSJed Brown #define __FUNCT__ "DMDASNESSetFunctionLocal"
2576cab3a1bSJed Brown /*@C
2586cab3a1bSJed Brown    DMDASNESSetFunctionLocal - set a local residual evaluation function
2596cab3a1bSJed Brown 
2606cab3a1bSJed Brown    Logically Collective
2616cab3a1bSJed Brown 
2626cab3a1bSJed Brown    Input Arguments:
2636cab3a1bSJed Brown +  dm - DM to associate callback with
2646cab3a1bSJed Brown .  func - local residual evaluation
2656cab3a1bSJed Brown -  ctx - optional context for local residual evaluation
2666cab3a1bSJed Brown 
2676cab3a1bSJed Brown    Calling sequence for func:
2686cab3a1bSJed Brown +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
269709ac0d2SJed Brown .  imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part
2706cab3a1bSJed Brown .  x - dimensional pointer to state at which to evaluate residual
2716cab3a1bSJed Brown .  f - dimensional pointer to residual, write the residual here
2726cab3a1bSJed Brown -  ctx - optional context passed above
2736cab3a1bSJed Brown 
2746cab3a1bSJed Brown    Level: beginner
2756cab3a1bSJed Brown 
2766cab3a1bSJed Brown .seealso: DMSNESSetFunction(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
2776cab3a1bSJed Brown @*/
278709ac0d2SJed Brown PetscErrorCode DMDASNESSetFunctionLocal(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,void*,void*,void*),void *ctx)
2796cab3a1bSJed Brown {
2806cab3a1bSJed Brown   PetscErrorCode ierr;
2816cab3a1bSJed Brown   SNESDM         sdm;
2826cab3a1bSJed Brown   DM_DA_SNES     *dmdasnes;
2836cab3a1bSJed Brown 
2846cab3a1bSJed Brown   PetscFunctionBegin;
2856cab3a1bSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2866cab3a1bSJed Brown   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
2876cab3a1bSJed Brown   ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr);
288709ac0d2SJed Brown   dmdasnes->residuallocalimode = imode;
2896cab3a1bSJed Brown   dmdasnes->residuallocal = func;
2906cab3a1bSJed Brown   dmdasnes->residuallocalctx = ctx;
2916cab3a1bSJed Brown   ierr = DMSNESSetFunction(dm,SNESComputeFunction_DMDA,dmdasnes);CHKERRQ(ierr);
292*6427ad5bSPeter Brune   ierr = DMSNESSetBlockFunction(dm,SNESComputeLocalBlockFunction_DMDA,dm);CHKERRQ(ierr);
2932961e153SJed Brown   if (!sdm->computejacobian) {  /* Call us for the Jacobian too, can be overridden by the user. */
2942961e153SJed Brown     ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes);CHKERRQ(ierr);
2952961e153SJed Brown   }
2962961e153SJed Brown   PetscFunctionReturn(0);
2972961e153SJed Brown }
2982961e153SJed Brown 
2992961e153SJed Brown #undef __FUNCT__
3002961e153SJed Brown #define __FUNCT__ "DMDASNESSetJacobianLocal"
3012961e153SJed Brown /*@C
3022961e153SJed Brown    DMDASNESSetJacobianLocal - set a local residual evaluation function
3032961e153SJed Brown 
3042961e153SJed Brown    Logically Collective
3052961e153SJed Brown 
3062961e153SJed Brown    Input Arguments:
3072961e153SJed Brown +  dm - DM to associate callback with
3082961e153SJed Brown .  func - local residual evaluation
3092961e153SJed Brown -  ctx - optional context for local residual evaluation
3102961e153SJed Brown 
3112961e153SJed Brown    Calling sequence for func:
3122961e153SJed Brown +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
3132961e153SJed Brown .  x - dimensional pointer to state at which to evaluate residual
3142961e153SJed Brown .  f - dimensional pointer to residual, write the residual here
3152961e153SJed Brown -  ctx - optional context passed above
3162961e153SJed Brown 
3172961e153SJed Brown    Level: beginner
3182961e153SJed Brown 
3192961e153SJed Brown .seealso: DMSNESSetJacobian(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
3202961e153SJed Brown @*/
3212961e153SJed Brown PetscErrorCode DMDASNESSetJacobianLocal(DM dm,PetscErrorCode (*func)(DMDALocalInfo*,void*,Mat,Mat,MatStructure*,void*),void *ctx)
3222961e153SJed Brown {
3232961e153SJed Brown   PetscErrorCode ierr;
3242961e153SJed Brown   SNESDM         sdm;
3252961e153SJed Brown   DM_DA_SNES     *dmdasnes;
3262961e153SJed Brown 
3272961e153SJed Brown   PetscFunctionBegin;
3282961e153SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3292961e153SJed Brown   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
3302961e153SJed Brown   ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr);
3312961e153SJed Brown   dmdasnes->jacobianlocal = func;
3322961e153SJed Brown   dmdasnes->jacobianlocalctx = ctx;
3332961e153SJed Brown   ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes);CHKERRQ(ierr);
334*6427ad5bSPeter Brune   ierr = DMSNESSetBlockJacobian(dm,SNESComputeLocalBlockJacobian_DMDA,dm);CHKERRQ(ierr);
3356cab3a1bSJed Brown   PetscFunctionReturn(0);
3366cab3a1bSJed Brown }
337