xref: /petsc/src/snes/utils/dmdasnes.c (revision 62cd874fff1a4a8f364c92545a705179b653d996)
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*);
106cab3a1bSJed Brown   void       *residuallocalctx;
116cab3a1bSJed Brown   void       *jacobianlocalctx;
122a4ee8f2SPeter Brune   void       *objectivelocalctx;
13709ac0d2SJed Brown   InsertMode residuallocalimode;
14dcad5f1cSBarry Smith 
15dcad5f1cSBarry Smith   /*   For Picard iteration defined locally */
16dcad5f1cSBarry Smith   PetscErrorCode (*rhsplocal)(DMDALocalInfo*,void*,void*,void*);
17d1e9a80fSBarry Smith   PetscErrorCode (*jacobianplocal)(DMDALocalInfo*,void*,Mat,Mat,void*);
18dcad5f1cSBarry Smith   void *picardlocalctx;
19942e3340SBarry Smith } DMSNES_DA;
206cab3a1bSJed Brown 
216cab3a1bSJed Brown #undef __FUNCT__
22942e3340SBarry Smith #define __FUNCT__ "DMSNESDestroy_DMDA"
23942e3340SBarry Smith static PetscErrorCode DMSNESDestroy_DMDA(DMSNES sdm)
246cab3a1bSJed Brown {
256cab3a1bSJed Brown   PetscErrorCode ierr;
266cab3a1bSJed Brown 
276cab3a1bSJed Brown   PetscFunctionBegin;
286cab3a1bSJed Brown   ierr = PetscFree(sdm->data);CHKERRQ(ierr);
296cab3a1bSJed Brown   PetscFunctionReturn(0);
306cab3a1bSJed Brown }
316cab3a1bSJed Brown 
326cab3a1bSJed Brown #undef __FUNCT__
33942e3340SBarry Smith #define __FUNCT__ "DMSNESDuplicate_DMDA"
3422c6f798SBarry Smith static PetscErrorCode DMSNESDuplicate_DMDA(DMSNES oldsdm,DMSNES sdm)
35e3c0b8a2SPeter Brune {
36e3c0b8a2SPeter Brune   PetscErrorCode ierr;
37e3c0b8a2SPeter Brune 
38e3c0b8a2SPeter Brune   PetscFunctionBegin;
39b00a9115SJed Brown   ierr = PetscNewLog(sdm,(DMSNES_DA**)&sdm->data);CHKERRQ(ierr);
40dcad5f1cSBarry Smith   if (oldsdm->data) {
41942e3340SBarry Smith     ierr = PetscMemcpy(sdm->data,oldsdm->data,sizeof(DMSNES_DA));CHKERRQ(ierr);
42dcad5f1cSBarry Smith   }
43e3c0b8a2SPeter Brune   PetscFunctionReturn(0);
44e3c0b8a2SPeter Brune }
45e3c0b8a2SPeter Brune 
46e3c0b8a2SPeter Brune 
47e3c0b8a2SPeter Brune #undef __FUNCT__
486cab3a1bSJed Brown #define __FUNCT__ "DMDASNESGetContext"
49942e3340SBarry Smith static PetscErrorCode DMDASNESGetContext(DM dm,DMSNES sdm,DMSNES_DA  **dmdasnes)
506cab3a1bSJed Brown {
516cab3a1bSJed Brown   PetscErrorCode ierr;
526cab3a1bSJed Brown 
536cab3a1bSJed Brown   PetscFunctionBegin;
540298fd71SBarry Smith   *dmdasnes = NULL;
556cab3a1bSJed Brown   if (!sdm->data) {
56b00a9115SJed Brown     ierr = PetscNewLog(dm,(DMSNES_DA**)&sdm->data);CHKERRQ(ierr);
5722c6f798SBarry Smith     sdm->ops->destroy   = DMSNESDestroy_DMDA;
5822c6f798SBarry Smith     sdm->ops->duplicate = DMSNESDuplicate_DMDA;
596cab3a1bSJed Brown   }
60942e3340SBarry Smith   *dmdasnes = (DMSNES_DA*)sdm->data;
616cab3a1bSJed Brown   PetscFunctionReturn(0);
626cab3a1bSJed Brown }
636cab3a1bSJed Brown 
646cab3a1bSJed Brown #undef __FUNCT__
656cab3a1bSJed Brown #define __FUNCT__ "SNESComputeFunction_DMDA"
666cab3a1bSJed Brown static PetscErrorCode SNESComputeFunction_DMDA(SNES snes,Vec X,Vec F,void *ctx)
676cab3a1bSJed Brown {
686cab3a1bSJed Brown   PetscErrorCode ierr;
696cab3a1bSJed Brown   DM             dm;
70942e3340SBarry Smith   DMSNES_DA      *dmdasnes = (DMSNES_DA*)ctx;
716cab3a1bSJed Brown   DMDALocalInfo  info;
726cab3a1bSJed Brown   Vec            Xloc;
736cab3a1bSJed Brown   void           *x,*f;
746cab3a1bSJed Brown 
756cab3a1bSJed Brown   PetscFunctionBegin;
762961e153SJed Brown   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
772961e153SJed Brown   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
782961e153SJed Brown   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
79ce94432eSBarry Smith   if (!dmdasnes->residuallocal) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context");
806cab3a1bSJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
816cab3a1bSJed Brown   ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
826cab3a1bSJed Brown   ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
836cab3a1bSJed Brown   ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
846cab3a1bSJed Brown   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
856cab3a1bSJed Brown   ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
86709ac0d2SJed Brown   switch (dmdasnes->residuallocalimode) {
87709ac0d2SJed Brown   case INSERT_VALUES: {
886cab3a1bSJed Brown     ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr);
89*62cd874fSBarry Smith     ierr = PetscLogEventBegin(SNES_FunctionEval,snes,X,F,0);CHKERRQ(ierr);
906cab3a1bSJed Brown     CHKMEMQ;
916cab3a1bSJed Brown     ierr = (*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx);CHKERRQ(ierr);
926cab3a1bSJed Brown     CHKMEMQ;
93*62cd874fSBarry Smith     ierr = PetscLogEventEnd(SNES_FunctionEval,snes,X,F,0);CHKERRQ(ierr);
946cab3a1bSJed Brown     ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr);
95709ac0d2SJed Brown   } break;
96709ac0d2SJed Brown   case ADD_VALUES: {
97709ac0d2SJed Brown     Vec Floc;
98709ac0d2SJed Brown     ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr);
99709ac0d2SJed Brown     ierr = VecZeroEntries(Floc);CHKERRQ(ierr);
100709ac0d2SJed Brown     ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr);
101*62cd874fSBarry Smith     ierr = PetscLogEventBegin(SNES_FunctionEval,snes,X,F,0);CHKERRQ(ierr);
102709ac0d2SJed Brown     CHKMEMQ;
103709ac0d2SJed Brown     ierr = (*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx);CHKERRQ(ierr);
104709ac0d2SJed Brown     CHKMEMQ;
105*62cd874fSBarry Smith     ierr = PetscLogEventEnd(SNES_FunctionEval,snes,X,F,0);CHKERRQ(ierr);
106709ac0d2SJed Brown     ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr);
107709ac0d2SJed Brown     ierr = VecZeroEntries(F);CHKERRQ(ierr);
108709ac0d2SJed Brown     ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
109709ac0d2SJed Brown     ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
110709ac0d2SJed Brown     ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr);
111709ac0d2SJed Brown   } break;
112ce94432eSBarry Smith   default: SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdasnes->residuallocalimode);
113709ac0d2SJed Brown   }
114709ac0d2SJed Brown   ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
1156cab3a1bSJed Brown   ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
116422a814eSBarry Smith   if (snes->domainerror) {
117422a814eSBarry Smith     ierr = VecSetInf(F);CHKERRQ(ierr);
118422a814eSBarry Smith   }
1196cab3a1bSJed Brown   PetscFunctionReturn(0);
1206cab3a1bSJed Brown }
1216cab3a1bSJed Brown 
1226cab3a1bSJed Brown #undef __FUNCT__
1232a4ee8f2SPeter Brune #define __FUNCT__ "SNESComputeObjective_DMDA"
1242a4ee8f2SPeter Brune static PetscErrorCode SNESComputeObjective_DMDA(SNES snes,Vec X,PetscReal *ob,void *ctx)
1252a4ee8f2SPeter Brune {
1262a4ee8f2SPeter Brune   PetscErrorCode ierr;
1272a4ee8f2SPeter Brune   DM             dm;
128942e3340SBarry Smith   DMSNES_DA      *dmdasnes = (DMSNES_DA*)ctx;
1292a4ee8f2SPeter Brune   DMDALocalInfo  info;
1302a4ee8f2SPeter Brune   Vec            Xloc;
1312a4ee8f2SPeter Brune   void           *x;
1322a4ee8f2SPeter Brune 
1332a4ee8f2SPeter Brune   PetscFunctionBegin;
1342a4ee8f2SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1352a4ee8f2SPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
1362a4ee8f2SPeter Brune   PetscValidPointer(ob,3);
137ce94432eSBarry Smith   if (!dmdasnes->objectivelocal) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context");
1382a4ee8f2SPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
1392a4ee8f2SPeter Brune   ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
1402a4ee8f2SPeter Brune   ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
1412a4ee8f2SPeter Brune   ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
1422a4ee8f2SPeter Brune   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
1432a4ee8f2SPeter Brune   ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
1442a4ee8f2SPeter Brune   CHKMEMQ;
1452a4ee8f2SPeter Brune   ierr = (*dmdasnes->objectivelocal)(&info,x,ob,dmdasnes->objectivelocalctx);CHKERRQ(ierr);
1462a4ee8f2SPeter Brune   CHKMEMQ;
1472a4ee8f2SPeter Brune   ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
1482a4ee8f2SPeter Brune   ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
1492a4ee8f2SPeter Brune   PetscFunctionReturn(0);
1502a4ee8f2SPeter Brune }
1512a4ee8f2SPeter Brune 
1526427ad5bSPeter Brune 
1536427ad5bSPeter Brune #undef __FUNCT__
1542961e153SJed Brown #define __FUNCT__ "SNESComputeJacobian_DMDA"
155422a814eSBarry Smith PetscErrorCode SNESComputeJacobian_DMDA(SNES snes,Vec X,Mat A,Mat B,void *ctx)
1562961e153SJed Brown {
1572961e153SJed Brown   PetscErrorCode ierr;
1582961e153SJed Brown   DM             dm;
159942e3340SBarry Smith   DMSNES_DA      *dmdasnes = (DMSNES_DA*)ctx;
1602961e153SJed Brown   DMDALocalInfo  info;
1612961e153SJed Brown   Vec            Xloc;
1622961e153SJed Brown   void           *x;
1632961e153SJed Brown 
1642961e153SJed Brown   PetscFunctionBegin;
165ce94432eSBarry Smith   if (!dmdasnes->residuallocal) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context");
1662961e153SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
1672961e153SJed Brown 
1682961e153SJed Brown   if (dmdasnes->jacobianlocal) {
1692961e153SJed Brown     ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
1702961e153SJed Brown     ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
1712961e153SJed Brown     ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
1722961e153SJed Brown     ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
1732961e153SJed Brown     ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
1742961e153SJed Brown     CHKMEMQ;
175d1e9a80fSBarry Smith     ierr = (*dmdasnes->jacobianlocal)(&info,x,A,B,dmdasnes->jacobianlocalctx);CHKERRQ(ierr);
1762961e153SJed Brown     CHKMEMQ;
1772961e153SJed Brown     ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
1782961e153SJed Brown     ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
1792961e153SJed Brown   } else {
1802961e153SJed Brown     MatFDColoring fdcoloring;
1812961e153SJed Brown     ierr = PetscObjectQuery((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject*)&fdcoloring);CHKERRQ(ierr);
1822961e153SJed Brown     if (!fdcoloring) {
1832961e153SJed Brown       ISColoring coloring;
1842961e153SJed Brown 
185b412c318SBarry Smith       ierr = DMCreateColoring(dm,dm->coloringtype,&coloring);CHKERRQ(ierr);
18694ab13aaSBarry Smith       ierr = MatFDColoringCreate(B,coloring,&fdcoloring);CHKERRQ(ierr);
1872961e153SJed Brown       switch (dm->coloringtype) {
1882961e153SJed Brown       case IS_COLORING_GLOBAL:
1892961e153SJed Brown         ierr = MatFDColoringSetFunction(fdcoloring,(PetscErrorCode (*)(void))SNESComputeFunction_DMDA,dmdasnes);CHKERRQ(ierr);
1902961e153SJed Brown         break;
191ce94432eSBarry Smith       default: SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"No support for coloring type '%s'",ISColoringTypes[dm->coloringtype]);
1922961e153SJed Brown       }
1932961e153SJed Brown       ierr = PetscObjectSetOptionsPrefix((PetscObject)fdcoloring,((PetscObject)dm)->prefix);CHKERRQ(ierr);
1942961e153SJed Brown       ierr = MatFDColoringSetFromOptions(fdcoloring);CHKERRQ(ierr);
19594ab13aaSBarry Smith       ierr = MatFDColoringSetUp(B,coloring,fdcoloring);CHKERRQ(ierr);
196f86b9fbaSHong Zhang       ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr);
1972961e153SJed Brown       ierr = PetscObjectCompose((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject)fdcoloring);CHKERRQ(ierr);
1982961e153SJed Brown       ierr = PetscObjectDereference((PetscObject)fdcoloring);CHKERRQ(ierr);
1992961e153SJed Brown 
2002961e153SJed Brown       /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call
2012961e153SJed Brown        * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the
2022961e153SJed Brown        * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually
2032961e153SJed Brown        * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be
204140e18c1SBarry Smith        * taken when the PetscObjectList for the Vec inside MatFDColoring is destroyed.
2052961e153SJed Brown        */
2062961e153SJed Brown       ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
2072961e153SJed Brown     }
208d1e9a80fSBarry Smith     ierr  = MatFDColoringApply(B,fdcoloring,X,snes);CHKERRQ(ierr);
2092961e153SJed Brown   }
2102961e153SJed Brown   /* This will be redundant if the user called both, but it's too common to forget. */
21194ab13aaSBarry Smith   if (A != B) {
21294ab13aaSBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
21394ab13aaSBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2142961e153SJed Brown   }
2152961e153SJed Brown   PetscFunctionReturn(0);
2162961e153SJed Brown }
2172961e153SJed Brown 
2182961e153SJed Brown #undef __FUNCT__
2196cab3a1bSJed Brown #define __FUNCT__ "DMDASNESSetFunctionLocal"
2206cab3a1bSJed Brown /*@C
2216cab3a1bSJed Brown    DMDASNESSetFunctionLocal - set a local residual evaluation function
2226cab3a1bSJed Brown 
2236cab3a1bSJed Brown    Logically Collective
2246cab3a1bSJed Brown 
2256cab3a1bSJed Brown    Input Arguments:
2266cab3a1bSJed Brown +  dm - DM to associate callback with
227e3c51ba2SJed Brown .  imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part
2286cab3a1bSJed Brown .  func - local residual evaluation
2296cab3a1bSJed Brown -  ctx - optional context for local residual evaluation
2306cab3a1bSJed Brown 
2310038884cSBarry Smith    Calling sequence:
2320038884cSBarry Smith    For PetscErrorCode (*func)(DMDALocalInfo *info,void *x, void *f, void *ctx),
2336cab3a1bSJed Brown +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
2340038884cSBarry Smith .  x - dimensional pointer to state at which to evaluate residual (e.g. PetscScalar *x or **x or ***x)
2350038884cSBarry Smith .  f - dimensional pointer to residual, write the residual here (e.g. PetscScalar *f or **f or ***f)
2366cab3a1bSJed Brown -  ctx - optional context passed above
2376cab3a1bSJed Brown 
2386cab3a1bSJed Brown    Level: beginner
2396cab3a1bSJed Brown 
2400038884cSBarry Smith .seealso: DMDASNESSetJacobianLocal(), DMSNESSetFunction(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
2416cab3a1bSJed Brown @*/
242709ac0d2SJed Brown PetscErrorCode DMDASNESSetFunctionLocal(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,void*,void*,void*),void *ctx)
2436cab3a1bSJed Brown {
2446cab3a1bSJed Brown   PetscErrorCode ierr;
245942e3340SBarry Smith   DMSNES         sdm;
246942e3340SBarry Smith   DMSNES_DA      *dmdasnes;
2476cab3a1bSJed Brown 
2486cab3a1bSJed Brown   PetscFunctionBegin;
2496cab3a1bSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
250942e3340SBarry Smith   ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr);
2516cab3a1bSJed Brown   ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr);
2521aa26658SKarl Rupp 
253709ac0d2SJed Brown   dmdasnes->residuallocalimode = imode;
2546cab3a1bSJed Brown   dmdasnes->residuallocal      = func;
2556cab3a1bSJed Brown   dmdasnes->residuallocalctx   = ctx;
2561aa26658SKarl Rupp 
2576cab3a1bSJed Brown   ierr = DMSNESSetFunction(dm,SNESComputeFunction_DMDA,dmdasnes);CHKERRQ(ierr);
25822c6f798SBarry Smith   if (!sdm->ops->computejacobian) {  /* Call us for the Jacobian too, can be overridden by the user. */
2592961e153SJed Brown     ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes);CHKERRQ(ierr);
2602961e153SJed Brown   }
2612961e153SJed Brown   PetscFunctionReturn(0);
2622961e153SJed Brown }
2632961e153SJed Brown 
2642961e153SJed Brown #undef __FUNCT__
2652961e153SJed Brown #define __FUNCT__ "DMDASNESSetJacobianLocal"
2662961e153SJed Brown /*@C
267e79ed307SJed Brown    DMDASNESSetJacobianLocal - set a local Jacobian evaluation function
2682961e153SJed Brown 
2692961e153SJed Brown    Logically Collective
2702961e153SJed Brown 
2712961e153SJed Brown    Input Arguments:
2722961e153SJed Brown +  dm - DM to associate callback with
2736b7eb5ceSMatthew G. Knepley .  func - local Jacobian evaluation
2746b7eb5ceSMatthew G. Knepley -  ctx - optional context for local Jacobian evaluation
2752961e153SJed Brown 
2760038884cSBarry Smith    Calling sequence:
2770038884cSBarry Smith    For PetscErrorCode (*func)(DMDALocalInfo *info,void *x,Mat J,Mat M,void *ctx),
2786b7eb5ceSMatthew G. Knepley +  info - DMDALocalInfo defining the subdomain to evaluate the Jacobian at
2790038884cSBarry Smith .  x - dimensional pointer to state at which to evaluate Jacobian (e.g. PetscScalar *x or **x or ***x)
2806b7eb5ceSMatthew G. Knepley .  J - Mat object for the Jacobian
2816b7eb5ceSMatthew G. Knepley .  M - Mat object for the Jacobian preconditioner matrix
2822961e153SJed Brown -  ctx - optional context passed above
2832961e153SJed Brown 
2842961e153SJed Brown    Level: beginner
2852961e153SJed Brown 
2860038884cSBarry Smith .seealso: DMDASNESSetFunctionLocal(), DMSNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
2872961e153SJed Brown @*/
288d1e9a80fSBarry Smith PetscErrorCode DMDASNESSetJacobianLocal(DM dm,PetscErrorCode (*func)(DMDALocalInfo*,void*,Mat,Mat,void*),void *ctx)
2892961e153SJed Brown {
2902961e153SJed Brown   PetscErrorCode ierr;
291942e3340SBarry Smith   DMSNES         sdm;
292942e3340SBarry Smith   DMSNES_DA      *dmdasnes;
2932961e153SJed Brown 
2942961e153SJed Brown   PetscFunctionBegin;
2952961e153SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
296942e3340SBarry Smith   ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr);
2972961e153SJed Brown   ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr);
2981aa26658SKarl Rupp 
2992961e153SJed Brown   dmdasnes->jacobianlocal    = func;
3002961e153SJed Brown   dmdasnes->jacobianlocalctx = ctx;
3011aa26658SKarl Rupp 
3022961e153SJed Brown   ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes);CHKERRQ(ierr);
3036cab3a1bSJed Brown   PetscFunctionReturn(0);
3046cab3a1bSJed Brown }
3052a4ee8f2SPeter Brune 
3062a4ee8f2SPeter Brune 
3072a4ee8f2SPeter Brune #undef __FUNCT__
3082a4ee8f2SPeter Brune #define __FUNCT__ "DMDASNESSetObjectiveLocal"
3092a4ee8f2SPeter Brune /*@C
3102a4ee8f2SPeter Brune    DMDASNESSetObjectiveLocal - set a local residual evaluation function
3112a4ee8f2SPeter Brune 
3122a4ee8f2SPeter Brune    Logically Collective
3132a4ee8f2SPeter Brune 
3142a4ee8f2SPeter Brune    Input Arguments:
3152a4ee8f2SPeter Brune +  dm - DM to associate callback with
3162a4ee8f2SPeter Brune .  func - local objective evaluation
3172a4ee8f2SPeter Brune -  ctx - optional context for local residual evaluation
3182a4ee8f2SPeter Brune 
3192a4ee8f2SPeter Brune    Calling sequence for func:
3202a4ee8f2SPeter Brune +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
3212a4ee8f2SPeter Brune .  x - dimensional pointer to state at which to evaluate residual
3222a4ee8f2SPeter Brune .  ob - eventual objective value
3232a4ee8f2SPeter Brune -  ctx - optional context passed above
3242a4ee8f2SPeter Brune 
3252a4ee8f2SPeter Brune    Level: beginner
3262a4ee8f2SPeter Brune 
3270038884cSBarry Smith .seealso: DMSNESSetFunction(), DMDASNESSetJacobianLocal(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
3282a4ee8f2SPeter Brune @*/
3292a4ee8f2SPeter Brune PetscErrorCode DMDASNESSetObjectiveLocal(DM dm,DMDASNESObjective func,void *ctx)
3302a4ee8f2SPeter Brune {
3312a4ee8f2SPeter Brune   PetscErrorCode ierr;
332942e3340SBarry Smith   DMSNES         sdm;
333942e3340SBarry Smith   DMSNES_DA      *dmdasnes;
3342a4ee8f2SPeter Brune 
3352a4ee8f2SPeter Brune   PetscFunctionBegin;
3362a4ee8f2SPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
337942e3340SBarry Smith   ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr);
3382a4ee8f2SPeter Brune   ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr);
3391aa26658SKarl Rupp 
3402a4ee8f2SPeter Brune   dmdasnes->objectivelocal    = func;
3412a4ee8f2SPeter Brune   dmdasnes->objectivelocalctx = ctx;
3421aa26658SKarl Rupp 
3432a4ee8f2SPeter Brune   ierr = DMSNESSetObjective(dm,SNESComputeObjective_DMDA,dmdasnes);CHKERRQ(ierr);
3442a4ee8f2SPeter Brune   PetscFunctionReturn(0);
3452a4ee8f2SPeter Brune }
346dcad5f1cSBarry Smith 
347dcad5f1cSBarry Smith #undef __FUNCT__
348dcad5f1cSBarry Smith #define __FUNCT__ "SNESComputePicard_DMDA"
349dcad5f1cSBarry Smith static PetscErrorCode SNESComputePicard_DMDA(SNES snes,Vec X,Vec F,void *ctx)
350dcad5f1cSBarry Smith {
351dcad5f1cSBarry Smith   PetscErrorCode ierr;
352dcad5f1cSBarry Smith   DM             dm;
353942e3340SBarry Smith   DMSNES_DA      *dmdasnes = (DMSNES_DA*)ctx;
354dcad5f1cSBarry Smith   DMDALocalInfo  info;
355dcad5f1cSBarry Smith   Vec            Xloc;
356dcad5f1cSBarry Smith   void           *x,*f;
357dcad5f1cSBarry Smith 
358dcad5f1cSBarry Smith   PetscFunctionBegin;
359dcad5f1cSBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
360dcad5f1cSBarry Smith   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
361dcad5f1cSBarry Smith   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
362ce94432eSBarry Smith   if (!dmdasnes->rhsplocal) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context");
363dcad5f1cSBarry Smith   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
364dcad5f1cSBarry Smith   ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
365dcad5f1cSBarry Smith   ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
366dcad5f1cSBarry Smith   ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
367dcad5f1cSBarry Smith   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
368dcad5f1cSBarry Smith   ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
369dcad5f1cSBarry Smith   switch (dmdasnes->residuallocalimode) {
370dcad5f1cSBarry Smith   case INSERT_VALUES: {
371dcad5f1cSBarry Smith     ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr);
372dcad5f1cSBarry Smith     CHKMEMQ;
373dcad5f1cSBarry Smith     ierr = (*dmdasnes->rhsplocal)(&info,x,f,dmdasnes->picardlocalctx);CHKERRQ(ierr);
374dcad5f1cSBarry Smith     CHKMEMQ;
375dcad5f1cSBarry Smith     ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr);
376dcad5f1cSBarry Smith   } break;
377dcad5f1cSBarry Smith   case ADD_VALUES: {
378dcad5f1cSBarry Smith     Vec Floc;
379dcad5f1cSBarry Smith     ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr);
380dcad5f1cSBarry Smith     ierr = VecZeroEntries(Floc);CHKERRQ(ierr);
381dcad5f1cSBarry Smith     ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr);
382dcad5f1cSBarry Smith     CHKMEMQ;
383dcad5f1cSBarry Smith     ierr = (*dmdasnes->rhsplocal)(&info,x,f,dmdasnes->picardlocalctx);CHKERRQ(ierr);
384dcad5f1cSBarry Smith     CHKMEMQ;
385dcad5f1cSBarry Smith     ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr);
386dcad5f1cSBarry Smith     ierr = VecZeroEntries(F);CHKERRQ(ierr);
387dcad5f1cSBarry Smith     ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
388dcad5f1cSBarry Smith     ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
389dcad5f1cSBarry Smith     ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr);
390dcad5f1cSBarry Smith   } break;
391ce94432eSBarry Smith   default: SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdasnes->residuallocalimode);
392dcad5f1cSBarry Smith   }
393dcad5f1cSBarry Smith   ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
394dcad5f1cSBarry Smith   ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
395dcad5f1cSBarry Smith   PetscFunctionReturn(0);
396dcad5f1cSBarry Smith }
397dcad5f1cSBarry Smith 
398dcad5f1cSBarry Smith #undef __FUNCT__
399dcad5f1cSBarry Smith #define __FUNCT__ "SNESComputePicardJacobian_DMDA"
400d1e9a80fSBarry Smith static PetscErrorCode SNESComputePicardJacobian_DMDA(SNES snes,Vec X,Mat A,Mat B,void *ctx)
401dcad5f1cSBarry Smith {
402dcad5f1cSBarry Smith   PetscErrorCode ierr;
403dcad5f1cSBarry Smith   DM             dm;
404942e3340SBarry Smith   DMSNES_DA      *dmdasnes = (DMSNES_DA*)ctx;
405dcad5f1cSBarry Smith   DMDALocalInfo  info;
406dcad5f1cSBarry Smith   Vec            Xloc;
407dcad5f1cSBarry Smith   void           *x;
408dcad5f1cSBarry Smith 
409dcad5f1cSBarry Smith   PetscFunctionBegin;
410ce94432eSBarry Smith   if (!dmdasnes->jacobianplocal) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context");
411dcad5f1cSBarry Smith   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
412dcad5f1cSBarry Smith 
413dcad5f1cSBarry Smith   ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
414dcad5f1cSBarry Smith   ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
415dcad5f1cSBarry Smith   ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
416dcad5f1cSBarry Smith   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
417dcad5f1cSBarry Smith   ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
418dcad5f1cSBarry Smith   CHKMEMQ;
419d1e9a80fSBarry Smith   ierr = (*dmdasnes->jacobianplocal)(&info,x,A,B,dmdasnes->picardlocalctx);CHKERRQ(ierr);
420dcad5f1cSBarry Smith   CHKMEMQ;
421dcad5f1cSBarry Smith   ierr  = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
422dcad5f1cSBarry Smith   ierr  = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
42394ab13aaSBarry Smith   if (A != B) {
42494ab13aaSBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
42594ab13aaSBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
426dcad5f1cSBarry Smith   }
427dcad5f1cSBarry Smith   PetscFunctionReturn(0);
428dcad5f1cSBarry Smith }
429dcad5f1cSBarry Smith 
430dcad5f1cSBarry Smith #undef __FUNCT__
431dcad5f1cSBarry Smith #define __FUNCT__ "DMDASNESSetPicardLocal"
432dcad5f1cSBarry Smith /*@C
433dcad5f1cSBarry Smith    DMDASNESSetPicardLocal - set a local right hand side and matrix evaluation function for Picard iteration
434dcad5f1cSBarry Smith 
435dcad5f1cSBarry Smith    Logically Collective
436dcad5f1cSBarry Smith 
437dcad5f1cSBarry Smith    Input Arguments:
438dcad5f1cSBarry Smith +  dm - DM to associate callback with
439dcad5f1cSBarry Smith .  imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part
440dcad5f1cSBarry Smith .  func - local residual evaluation
441dcad5f1cSBarry Smith -  ctx - optional context for local residual evaluation
442dcad5f1cSBarry Smith 
443dcad5f1cSBarry Smith    Calling sequence for func:
444dcad5f1cSBarry Smith +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
445dcad5f1cSBarry Smith .  x - dimensional pointer to state at which to evaluate residual
446dcad5f1cSBarry Smith .  f - dimensional pointer to residual, write the residual here
447dcad5f1cSBarry Smith -  ctx - optional context passed above
448dcad5f1cSBarry Smith 
449dcad5f1cSBarry Smith    Notes:  The user must use
450dcad5f1cSBarry Smith     extern PetscErrorCode  SNESPicardComputeFunction(SNES,Vec,Vec,void*);
45194ab13aaSBarry Smith     extern PetscErrorCode  SNESPicardComputeJacobian(SNES,Vec,Mat,Mat,MatStructure*,void*);
4520298fd71SBarry Smith     ierr = SNESSetFunction(snes,NULL,SNESPicardComputeFunction,&user);CHKERRQ(ierr);
453dcad5f1cSBarry Smith     in their code before calling this routine.
454dcad5f1cSBarry Smith 
455dcad5f1cSBarry Smith 
456dcad5f1cSBarry Smith    Level: beginner
457dcad5f1cSBarry Smith 
458dcad5f1cSBarry Smith .seealso: DMSNESSetFunction(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
459dcad5f1cSBarry Smith @*/
460dcad5f1cSBarry Smith PetscErrorCode DMDASNESSetPicardLocal(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,void*,void*,void*),
461d1e9a80fSBarry Smith                                       PetscErrorCode (*jac)(DMDALocalInfo*,void*,Mat,Mat,void*),void *ctx)
462dcad5f1cSBarry Smith {
463dcad5f1cSBarry Smith   PetscErrorCode ierr;
464942e3340SBarry Smith   DMSNES         sdm;
465942e3340SBarry Smith   DMSNES_DA      *dmdasnes;
466dcad5f1cSBarry Smith 
467dcad5f1cSBarry Smith   PetscFunctionBegin;
468dcad5f1cSBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
469942e3340SBarry Smith   ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr);
470dcad5f1cSBarry Smith   ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr);
4711aa26658SKarl Rupp 
472dcad5f1cSBarry Smith   dmdasnes->residuallocalimode = imode;
473dcad5f1cSBarry Smith   dmdasnes->rhsplocal          = func;
474dcad5f1cSBarry Smith   dmdasnes->jacobianplocal     = jac;
475dcad5f1cSBarry Smith   dmdasnes->picardlocalctx     = ctx;
4761aa26658SKarl Rupp 
477dcad5f1cSBarry Smith   ierr = DMSNESSetPicard(dm,SNESComputePicard_DMDA,SNESComputePicardJacobian_DMDA,dmdasnes);CHKERRQ(ierr);
478dcad5f1cSBarry Smith   PetscFunctionReturn(0);
479dcad5f1cSBarry Smith }
480