xref: /petsc/src/snes/utils/dmdasnes.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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 
21942e3340SBarry Smith static PetscErrorCode DMSNESDestroy_DMDA(DMSNES sdm)
226cab3a1bSJed Brown {
236cab3a1bSJed Brown   PetscFunctionBegin;
24*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(sdm->data));
256cab3a1bSJed Brown   PetscFunctionReturn(0);
266cab3a1bSJed Brown }
276cab3a1bSJed Brown 
2822c6f798SBarry Smith static PetscErrorCode DMSNESDuplicate_DMDA(DMSNES oldsdm,DMSNES sdm)
29e3c0b8a2SPeter Brune {
30e3c0b8a2SPeter Brune   PetscFunctionBegin;
31*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscNewLog(sdm,(DMSNES_DA**)&sdm->data));
32dcad5f1cSBarry Smith   if (oldsdm->data) {
33*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMemcpy(sdm->data,oldsdm->data,sizeof(DMSNES_DA)));
34dcad5f1cSBarry Smith   }
35e3c0b8a2SPeter Brune   PetscFunctionReturn(0);
36e3c0b8a2SPeter Brune }
37e3c0b8a2SPeter Brune 
38942e3340SBarry Smith static PetscErrorCode DMDASNESGetContext(DM dm,DMSNES sdm,DMSNES_DA  **dmdasnes)
396cab3a1bSJed Brown {
406cab3a1bSJed Brown   PetscFunctionBegin;
410298fd71SBarry Smith   *dmdasnes = NULL;
426cab3a1bSJed Brown   if (!sdm->data) {
43*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscNewLog(dm,(DMSNES_DA**)&sdm->data));
4422c6f798SBarry Smith     sdm->ops->destroy   = DMSNESDestroy_DMDA;
4522c6f798SBarry Smith     sdm->ops->duplicate = DMSNESDuplicate_DMDA;
466cab3a1bSJed Brown   }
47942e3340SBarry Smith   *dmdasnes = (DMSNES_DA*)sdm->data;
486cab3a1bSJed Brown   PetscFunctionReturn(0);
496cab3a1bSJed Brown }
506cab3a1bSJed Brown 
516cab3a1bSJed Brown static PetscErrorCode SNESComputeFunction_DMDA(SNES snes,Vec X,Vec F,void *ctx)
526cab3a1bSJed Brown {
536cab3a1bSJed Brown   DM             dm;
54942e3340SBarry Smith   DMSNES_DA      *dmdasnes = (DMSNES_DA*)ctx;
556cab3a1bSJed Brown   DMDALocalInfo  info;
566cab3a1bSJed Brown   Vec            Xloc;
576cab3a1bSJed Brown   void           *x,*f;
586cab3a1bSJed Brown 
596cab3a1bSJed Brown   PetscFunctionBegin;
602961e153SJed Brown   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
612961e153SJed Brown   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
622961e153SJed Brown   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
63*5f80ce2aSJacob Faibussowitsch   PetscCheck(dmdasnes->residuallocal,PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context");
64*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetDM(snes,&dm));
65*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(dm,&Xloc));
66*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc));
67*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc));
68*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetLocalInfo(dm,&info));
69*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(dm,Xloc,&x));
70709ac0d2SJed Brown   switch (dmdasnes->residuallocalimode) {
71709ac0d2SJed Brown   case INSERT_VALUES: {
72*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecGetArray(dm,F,&f));
73*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLogEventBegin(SNES_FunctionEval,snes,X,F,0));
746cab3a1bSJed Brown     CHKMEMQ;
75*5f80ce2aSJacob Faibussowitsch     CHKERRQ((*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx));
766cab3a1bSJed Brown     CHKMEMQ;
77*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLogEventEnd(SNES_FunctionEval,snes,X,F,0));
78*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecRestoreArray(dm,F,&f));
79709ac0d2SJed Brown   } break;
80709ac0d2SJed Brown   case ADD_VALUES: {
81709ac0d2SJed Brown     Vec Floc;
82*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetLocalVector(dm,&Floc));
83*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecZeroEntries(Floc));
84*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecGetArray(dm,Floc,&f));
85*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLogEventBegin(SNES_FunctionEval,snes,X,F,0));
86709ac0d2SJed Brown     CHKMEMQ;
87*5f80ce2aSJacob Faibussowitsch     CHKERRQ((*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx));
88709ac0d2SJed Brown     CHKMEMQ;
89*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLogEventEnd(SNES_FunctionEval,snes,X,F,0));
90*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecRestoreArray(dm,Floc,&f));
91*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecZeroEntries(F));
92*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F));
93*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F));
94*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMRestoreLocalVector(dm,&Floc));
95709ac0d2SJed Brown   } break;
9698921bdaSJacob Faibussowitsch   default: SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdasnes->residuallocalimode);
97709ac0d2SJed Brown   }
98*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(dm,Xloc,&x));
99*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(dm,&Xloc));
100422a814eSBarry Smith   if (snes->domainerror) {
101*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetInf(F));
102422a814eSBarry Smith   }
1036cab3a1bSJed Brown   PetscFunctionReturn(0);
1046cab3a1bSJed Brown }
1056cab3a1bSJed Brown 
1062a4ee8f2SPeter Brune static PetscErrorCode SNESComputeObjective_DMDA(SNES snes,Vec X,PetscReal *ob,void *ctx)
1072a4ee8f2SPeter Brune {
1082a4ee8f2SPeter Brune   DM             dm;
109942e3340SBarry Smith   DMSNES_DA      *dmdasnes = (DMSNES_DA*)ctx;
1102a4ee8f2SPeter Brune   DMDALocalInfo  info;
1112a4ee8f2SPeter Brune   Vec            Xloc;
1122a4ee8f2SPeter Brune   void           *x;
1132a4ee8f2SPeter Brune 
1142a4ee8f2SPeter Brune   PetscFunctionBegin;
1152a4ee8f2SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1162a4ee8f2SPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
1172a4ee8f2SPeter Brune   PetscValidPointer(ob,3);
1182c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!dmdasnes->objectivelocal,PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context");
119*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetDM(snes,&dm));
120*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(dm,&Xloc));
121*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc));
122*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc));
123*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetLocalInfo(dm,&info));
124*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(dm,Xloc,&x));
1252a4ee8f2SPeter Brune   CHKMEMQ;
126*5f80ce2aSJacob Faibussowitsch   CHKERRQ((*dmdasnes->objectivelocal)(&info,x,ob,dmdasnes->objectivelocalctx));
1272a4ee8f2SPeter Brune   CHKMEMQ;
128*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(dm,Xloc,&x));
129*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(dm,&Xloc));
1302a4ee8f2SPeter Brune   PetscFunctionReturn(0);
1312a4ee8f2SPeter Brune }
1322a4ee8f2SPeter Brune 
133cc2e6a90SBarry Smith /* Routine is called by example, hence must be labeled PETSC_EXTERN */
134cc2e6a90SBarry Smith PETSC_EXTERN PetscErrorCode SNESComputeJacobian_DMDA(SNES snes,Vec X,Mat A,Mat B,void *ctx)
1352961e153SJed Brown {
1362961e153SJed Brown   DM             dm;
137942e3340SBarry Smith   DMSNES_DA      *dmdasnes = (DMSNES_DA*)ctx;
1382961e153SJed Brown   DMDALocalInfo  info;
1392961e153SJed Brown   Vec            Xloc;
1402961e153SJed Brown   void           *x;
1412961e153SJed Brown 
1422961e153SJed Brown   PetscFunctionBegin;
1432c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!dmdasnes->residuallocal,PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context");
144*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetDM(snes,&dm));
1452961e153SJed Brown 
1462961e153SJed Brown   if (dmdasnes->jacobianlocal) {
147*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetLocalVector(dm,&Xloc));
148*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc));
149*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc));
150*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAGetLocalInfo(dm,&info));
151*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecGetArray(dm,Xloc,&x));
1522961e153SJed Brown     CHKMEMQ;
153*5f80ce2aSJacob Faibussowitsch     CHKERRQ((*dmdasnes->jacobianlocal)(&info,x,A,B,dmdasnes->jacobianlocalctx));
1542961e153SJed Brown     CHKMEMQ;
155*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecRestoreArray(dm,Xloc,&x));
156*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMRestoreLocalVector(dm,&Xloc));
1572961e153SJed Brown   } else {
1582961e153SJed Brown     MatFDColoring fdcoloring;
159*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectQuery((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject*)&fdcoloring));
1602961e153SJed Brown     if (!fdcoloring) {
1612961e153SJed Brown       ISColoring coloring;
1622961e153SJed Brown 
163*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMCreateColoring(dm,dm->coloringtype,&coloring));
164*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatFDColoringCreate(B,coloring,&fdcoloring));
1652961e153SJed Brown       switch (dm->coloringtype) {
1662961e153SJed Brown       case IS_COLORING_GLOBAL:
167*5f80ce2aSJacob Faibussowitsch         CHKERRQ(MatFDColoringSetFunction(fdcoloring,(PetscErrorCode (*)(void))SNESComputeFunction_DMDA,dmdasnes));
1682961e153SJed Brown         break;
16998921bdaSJacob Faibussowitsch       default: SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"No support for coloring type '%s'",ISColoringTypes[dm->coloringtype]);
1702961e153SJed Brown       }
171*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject)fdcoloring,((PetscObject)dm)->prefix));
172*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatFDColoringSetFromOptions(fdcoloring));
173*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatFDColoringSetUp(B,coloring,fdcoloring));
174*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISColoringDestroy(&coloring));
175*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscObjectCompose((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject)fdcoloring));
176*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscObjectDereference((PetscObject)fdcoloring));
1772961e153SJed Brown 
1782961e153SJed Brown       /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call
1792961e153SJed Brown        * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the
1802961e153SJed Brown        * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually
1812961e153SJed Brown        * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be
182140e18c1SBarry Smith        * taken when the PetscObjectList for the Vec inside MatFDColoring is destroyed.
1832961e153SJed Brown        */
184*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscObjectDereference((PetscObject)dm));
1852961e153SJed Brown     }
186*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatFDColoringApply(B,fdcoloring,X,snes));
1872961e153SJed Brown   }
1882961e153SJed Brown   /* This will be redundant if the user called both, but it's too common to forget. */
18994ab13aaSBarry Smith   if (A != B) {
190*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
191*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
1922961e153SJed Brown   }
1932961e153SJed Brown   PetscFunctionReturn(0);
1942961e153SJed Brown }
1952961e153SJed Brown 
1966cab3a1bSJed Brown /*@C
1976cab3a1bSJed Brown    DMDASNESSetFunctionLocal - set a local residual evaluation function
1986cab3a1bSJed Brown 
1996cab3a1bSJed Brown    Logically Collective
2006cab3a1bSJed Brown 
2014165533cSJose E. Roman    Input Parameters:
2026cab3a1bSJed Brown +  dm - DM to associate callback with
203e3c51ba2SJed Brown .  imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part
2046cab3a1bSJed Brown .  func - local residual evaluation
2056cab3a1bSJed Brown -  ctx - optional context for local residual evaluation
2066cab3a1bSJed Brown 
2070038884cSBarry Smith    Calling sequence:
2080038884cSBarry Smith    For PetscErrorCode (*func)(DMDALocalInfo *info,void *x, void *f, void *ctx),
2096cab3a1bSJed Brown +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
2100038884cSBarry Smith .  x - dimensional pointer to state at which to evaluate residual (e.g. PetscScalar *x or **x or ***x)
2110038884cSBarry Smith .  f - dimensional pointer to residual, write the residual here (e.g. PetscScalar *f or **f or ***f)
2126cab3a1bSJed Brown -  ctx - optional context passed above
2136cab3a1bSJed Brown 
2146cab3a1bSJed Brown    Level: beginner
2156cab3a1bSJed Brown 
2160038884cSBarry Smith .seealso: DMDASNESSetJacobianLocal(), DMSNESSetFunction(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
2176cab3a1bSJed Brown @*/
218709ac0d2SJed Brown PetscErrorCode DMDASNESSetFunctionLocal(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,void*,void*,void*),void *ctx)
2196cab3a1bSJed Brown {
220942e3340SBarry Smith   DMSNES         sdm;
221942e3340SBarry Smith   DMSNES_DA      *dmdasnes;
2226cab3a1bSJed Brown 
2236cab3a1bSJed Brown   PetscFunctionBegin;
2246cab3a1bSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
225*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDMSNESWrite(dm,&sdm));
226*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASNESGetContext(dm,sdm,&dmdasnes));
2271aa26658SKarl Rupp 
228709ac0d2SJed Brown   dmdasnes->residuallocalimode = imode;
2296cab3a1bSJed Brown   dmdasnes->residuallocal      = func;
2306cab3a1bSJed Brown   dmdasnes->residuallocalctx   = ctx;
2311aa26658SKarl Rupp 
232*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSNESSetFunction(dm,SNESComputeFunction_DMDA,dmdasnes));
23322c6f798SBarry Smith   if (!sdm->ops->computejacobian) {  /* Call us for the Jacobian too, can be overridden by the user. */
234*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes));
2352961e153SJed Brown   }
2362961e153SJed Brown   PetscFunctionReturn(0);
2372961e153SJed Brown }
2382961e153SJed Brown 
2392961e153SJed Brown /*@C
240e79ed307SJed Brown    DMDASNESSetJacobianLocal - set a local Jacobian evaluation function
2412961e153SJed Brown 
2422961e153SJed Brown    Logically Collective
2432961e153SJed Brown 
2444165533cSJose E. Roman    Input Parameters:
2452961e153SJed Brown +  dm - DM to associate callback with
2466b7eb5ceSMatthew G. Knepley .  func - local Jacobian evaluation
2476b7eb5ceSMatthew G. Knepley -  ctx - optional context for local Jacobian evaluation
2482961e153SJed Brown 
2490038884cSBarry Smith    Calling sequence:
2500038884cSBarry Smith    For PetscErrorCode (*func)(DMDALocalInfo *info,void *x,Mat J,Mat M,void *ctx),
2516b7eb5ceSMatthew G. Knepley +  info - DMDALocalInfo defining the subdomain to evaluate the Jacobian at
2520038884cSBarry Smith .  x - dimensional pointer to state at which to evaluate Jacobian (e.g. PetscScalar *x or **x or ***x)
2536b7eb5ceSMatthew G. Knepley .  J - Mat object for the Jacobian
2546b7eb5ceSMatthew G. Knepley .  M - Mat object for the Jacobian preconditioner matrix
2552961e153SJed Brown -  ctx - optional context passed above
2562961e153SJed Brown 
2572961e153SJed Brown    Level: beginner
2582961e153SJed Brown 
2590038884cSBarry Smith .seealso: DMDASNESSetFunctionLocal(), DMSNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
2602961e153SJed Brown @*/
261d1e9a80fSBarry Smith PetscErrorCode DMDASNESSetJacobianLocal(DM dm,PetscErrorCode (*func)(DMDALocalInfo*,void*,Mat,Mat,void*),void *ctx)
2622961e153SJed Brown {
263942e3340SBarry Smith   DMSNES         sdm;
264942e3340SBarry Smith   DMSNES_DA      *dmdasnes;
2652961e153SJed Brown 
2662961e153SJed Brown   PetscFunctionBegin;
2672961e153SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
268*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDMSNESWrite(dm,&sdm));
269*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASNESGetContext(dm,sdm,&dmdasnes));
2701aa26658SKarl Rupp 
2712961e153SJed Brown   dmdasnes->jacobianlocal    = func;
2722961e153SJed Brown   dmdasnes->jacobianlocalctx = ctx;
2731aa26658SKarl Rupp 
274*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes));
2756cab3a1bSJed Brown   PetscFunctionReturn(0);
2766cab3a1bSJed Brown }
2772a4ee8f2SPeter Brune 
2782a4ee8f2SPeter Brune /*@C
2792a4ee8f2SPeter Brune    DMDASNESSetObjectiveLocal - set a local residual evaluation function
2802a4ee8f2SPeter Brune 
2812a4ee8f2SPeter Brune    Logically Collective
2822a4ee8f2SPeter Brune 
2834165533cSJose E. Roman    Input Parameters:
2842a4ee8f2SPeter Brune +  dm - DM to associate callback with
2852a4ee8f2SPeter Brune .  func - local objective evaluation
2862a4ee8f2SPeter Brune -  ctx - optional context for local residual evaluation
2872a4ee8f2SPeter Brune 
2882a4ee8f2SPeter Brune    Calling sequence for func:
2892a4ee8f2SPeter Brune +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
2902a4ee8f2SPeter Brune .  x - dimensional pointer to state at which to evaluate residual
2912a4ee8f2SPeter Brune .  ob - eventual objective value
2922a4ee8f2SPeter Brune -  ctx - optional context passed above
2932a4ee8f2SPeter Brune 
2942a4ee8f2SPeter Brune    Level: beginner
2952a4ee8f2SPeter Brune 
2960038884cSBarry Smith .seealso: DMSNESSetFunction(), DMDASNESSetJacobianLocal(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
2972a4ee8f2SPeter Brune @*/
2982a4ee8f2SPeter Brune PetscErrorCode DMDASNESSetObjectiveLocal(DM dm,DMDASNESObjective func,void *ctx)
2992a4ee8f2SPeter Brune {
300942e3340SBarry Smith   DMSNES         sdm;
301942e3340SBarry Smith   DMSNES_DA      *dmdasnes;
3022a4ee8f2SPeter Brune 
3032a4ee8f2SPeter Brune   PetscFunctionBegin;
3042a4ee8f2SPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
305*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDMSNESWrite(dm,&sdm));
306*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASNESGetContext(dm,sdm,&dmdasnes));
3071aa26658SKarl Rupp 
3082a4ee8f2SPeter Brune   dmdasnes->objectivelocal    = func;
3092a4ee8f2SPeter Brune   dmdasnes->objectivelocalctx = ctx;
3101aa26658SKarl Rupp 
311*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSNESSetObjective(dm,SNESComputeObjective_DMDA,dmdasnes));
3122a4ee8f2SPeter Brune   PetscFunctionReturn(0);
3132a4ee8f2SPeter Brune }
314dcad5f1cSBarry Smith 
315dcad5f1cSBarry Smith static PetscErrorCode SNESComputePicard_DMDA(SNES snes,Vec X,Vec F,void *ctx)
316dcad5f1cSBarry Smith {
317dcad5f1cSBarry Smith   DM             dm;
318942e3340SBarry Smith   DMSNES_DA      *dmdasnes = (DMSNES_DA*)ctx;
319dcad5f1cSBarry Smith   DMDALocalInfo  info;
320dcad5f1cSBarry Smith   Vec            Xloc;
321dcad5f1cSBarry Smith   void           *x,*f;
322dcad5f1cSBarry Smith 
323dcad5f1cSBarry Smith   PetscFunctionBegin;
324dcad5f1cSBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
325dcad5f1cSBarry Smith   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
326dcad5f1cSBarry Smith   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
3272c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!dmdasnes->rhsplocal,PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context");
328*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetDM(snes,&dm));
329*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(dm,&Xloc));
330*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc));
331*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc));
332*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetLocalInfo(dm,&info));
333*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(dm,Xloc,&x));
334dcad5f1cSBarry Smith   switch (dmdasnes->residuallocalimode) {
335dcad5f1cSBarry Smith   case INSERT_VALUES: {
336*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecGetArray(dm,F,&f));
337dcad5f1cSBarry Smith     CHKMEMQ;
338*5f80ce2aSJacob Faibussowitsch     CHKERRQ((*dmdasnes->rhsplocal)(&info,x,f,dmdasnes->picardlocalctx));
339dcad5f1cSBarry Smith     CHKMEMQ;
340*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecRestoreArray(dm,F,&f));
341dcad5f1cSBarry Smith   } break;
342dcad5f1cSBarry Smith   case ADD_VALUES: {
343dcad5f1cSBarry Smith     Vec Floc;
344*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetLocalVector(dm,&Floc));
345*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecZeroEntries(Floc));
346*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecGetArray(dm,Floc,&f));
347dcad5f1cSBarry Smith     CHKMEMQ;
348*5f80ce2aSJacob Faibussowitsch     CHKERRQ((*dmdasnes->rhsplocal)(&info,x,f,dmdasnes->picardlocalctx));
349dcad5f1cSBarry Smith     CHKMEMQ;
350*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecRestoreArray(dm,Floc,&f));
351*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecZeroEntries(F));
352*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F));
353*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F));
354*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMRestoreLocalVector(dm,&Floc));
355dcad5f1cSBarry Smith   } break;
35698921bdaSJacob Faibussowitsch   default: SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdasnes->residuallocalimode);
357dcad5f1cSBarry Smith   }
358*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(dm,Xloc,&x));
359*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(dm,&Xloc));
360dcad5f1cSBarry Smith   PetscFunctionReturn(0);
361dcad5f1cSBarry Smith }
362dcad5f1cSBarry Smith 
363d1e9a80fSBarry Smith static PetscErrorCode SNESComputePicardJacobian_DMDA(SNES snes,Vec X,Mat A,Mat B,void *ctx)
364dcad5f1cSBarry Smith {
365dcad5f1cSBarry Smith   DM             dm;
366942e3340SBarry Smith   DMSNES_DA      *dmdasnes = (DMSNES_DA*)ctx;
367dcad5f1cSBarry Smith   DMDALocalInfo  info;
368dcad5f1cSBarry Smith   Vec            Xloc;
369dcad5f1cSBarry Smith   void           *x;
370dcad5f1cSBarry Smith 
371dcad5f1cSBarry Smith   PetscFunctionBegin;
3722c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!dmdasnes->jacobianplocal,PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context");
373*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetDM(snes,&dm));
374dcad5f1cSBarry Smith 
375*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(dm,&Xloc));
376*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc));
377*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc));
378*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetLocalInfo(dm,&info));
379*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(dm,Xloc,&x));
380dcad5f1cSBarry Smith   CHKMEMQ;
381*5f80ce2aSJacob Faibussowitsch   CHKERRQ((*dmdasnes->jacobianplocal)(&info,x,A,B,dmdasnes->picardlocalctx));
382dcad5f1cSBarry Smith   CHKMEMQ;
383*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(dm,Xloc,&x));
384*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(dm,&Xloc));
38594ab13aaSBarry Smith   if (A != B) {
386*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
387*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
388dcad5f1cSBarry Smith   }
389dcad5f1cSBarry Smith   PetscFunctionReturn(0);
390dcad5f1cSBarry Smith }
391dcad5f1cSBarry Smith 
392dcad5f1cSBarry Smith /*@C
393dcad5f1cSBarry Smith    DMDASNESSetPicardLocal - set a local right hand side and matrix evaluation function for Picard iteration
394dcad5f1cSBarry Smith 
395dcad5f1cSBarry Smith    Logically Collective
396dcad5f1cSBarry Smith 
3974165533cSJose E. Roman    Input Parameters:
398dcad5f1cSBarry Smith +  dm - DM to associate callback with
399dcad5f1cSBarry Smith .  imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part
400dcad5f1cSBarry Smith .  func - local residual evaluation
401dcad5f1cSBarry Smith -  ctx - optional context for local residual evaluation
402dcad5f1cSBarry Smith 
403dcad5f1cSBarry Smith    Calling sequence for func:
404dcad5f1cSBarry Smith +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
405dcad5f1cSBarry Smith .  x - dimensional pointer to state at which to evaluate residual
406dcad5f1cSBarry Smith .  f - dimensional pointer to residual, write the residual here
407dcad5f1cSBarry Smith -  ctx - optional context passed above
408dcad5f1cSBarry Smith 
40995452b02SPatrick Sanan    Notes:
41095452b02SPatrick Sanan     The user must use
411*5f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESSetFunction(snes,NULL,SNESPicardComputeFunction,&user));
412dcad5f1cSBarry Smith     in their code before calling this routine.
413dcad5f1cSBarry Smith 
414dcad5f1cSBarry Smith    Level: beginner
415dcad5f1cSBarry Smith 
416dcad5f1cSBarry Smith .seealso: DMSNESSetFunction(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
417dcad5f1cSBarry Smith @*/
418dcad5f1cSBarry Smith PetscErrorCode DMDASNESSetPicardLocal(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,void*,void*,void*),
419d1e9a80fSBarry Smith                                       PetscErrorCode (*jac)(DMDALocalInfo*,void*,Mat,Mat,void*),void *ctx)
420dcad5f1cSBarry Smith {
421942e3340SBarry Smith   DMSNES         sdm;
422942e3340SBarry Smith   DMSNES_DA      *dmdasnes;
423dcad5f1cSBarry Smith 
424dcad5f1cSBarry Smith   PetscFunctionBegin;
425dcad5f1cSBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
426*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDMSNESWrite(dm,&sdm));
427*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASNESGetContext(dm,sdm,&dmdasnes));
4281aa26658SKarl Rupp 
429dcad5f1cSBarry Smith   dmdasnes->residuallocalimode = imode;
430dcad5f1cSBarry Smith   dmdasnes->rhsplocal          = func;
431dcad5f1cSBarry Smith   dmdasnes->jacobianplocal     = jac;
432dcad5f1cSBarry Smith   dmdasnes->picardlocalctx     = ctx;
4331aa26658SKarl Rupp 
434*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSNESSetPicard(dm,SNESComputePicard_DMDA,SNESComputePicardJacobian_DMDA,dmdasnes));
435*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSNESSetMFFunction(dm,SNESComputeFunction_DMDA,dmdasnes));
436dcad5f1cSBarry Smith   PetscFunctionReturn(0);
437dcad5f1cSBarry Smith }
438