xref: /petsc/src/snes/utils/dmdasnes.c (revision 98921bda46e76d7aaed9e0138c5ff9d0ce93f355)
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   PetscErrorCode ierr;
246cab3a1bSJed Brown 
256cab3a1bSJed Brown   PetscFunctionBegin;
266cab3a1bSJed Brown   ierr = PetscFree(sdm->data);CHKERRQ(ierr);
276cab3a1bSJed Brown   PetscFunctionReturn(0);
286cab3a1bSJed Brown }
296cab3a1bSJed Brown 
3022c6f798SBarry Smith static PetscErrorCode DMSNESDuplicate_DMDA(DMSNES oldsdm,DMSNES sdm)
31e3c0b8a2SPeter Brune {
32e3c0b8a2SPeter Brune   PetscErrorCode ierr;
33e3c0b8a2SPeter Brune 
34e3c0b8a2SPeter Brune   PetscFunctionBegin;
35b00a9115SJed Brown   ierr = PetscNewLog(sdm,(DMSNES_DA**)&sdm->data);CHKERRQ(ierr);
36dcad5f1cSBarry Smith   if (oldsdm->data) {
37942e3340SBarry Smith     ierr = PetscMemcpy(sdm->data,oldsdm->data,sizeof(DMSNES_DA));CHKERRQ(ierr);
38dcad5f1cSBarry Smith   }
39e3c0b8a2SPeter Brune   PetscFunctionReturn(0);
40e3c0b8a2SPeter Brune }
41e3c0b8a2SPeter Brune 
42942e3340SBarry Smith static PetscErrorCode DMDASNESGetContext(DM dm,DMSNES sdm,DMSNES_DA  **dmdasnes)
436cab3a1bSJed Brown {
446cab3a1bSJed Brown   PetscErrorCode ierr;
456cab3a1bSJed Brown 
466cab3a1bSJed Brown   PetscFunctionBegin;
470298fd71SBarry Smith   *dmdasnes = NULL;
486cab3a1bSJed Brown   if (!sdm->data) {
49b00a9115SJed Brown     ierr = PetscNewLog(dm,(DMSNES_DA**)&sdm->data);CHKERRQ(ierr);
5022c6f798SBarry Smith     sdm->ops->destroy   = DMSNESDestroy_DMDA;
5122c6f798SBarry Smith     sdm->ops->duplicate = DMSNESDuplicate_DMDA;
526cab3a1bSJed Brown   }
53942e3340SBarry Smith   *dmdasnes = (DMSNES_DA*)sdm->data;
546cab3a1bSJed Brown   PetscFunctionReturn(0);
556cab3a1bSJed Brown }
566cab3a1bSJed Brown 
576cab3a1bSJed Brown static PetscErrorCode SNESComputeFunction_DMDA(SNES snes,Vec X,Vec F,void *ctx)
586cab3a1bSJed Brown {
596cab3a1bSJed Brown   PetscErrorCode ierr;
606cab3a1bSJed Brown   DM             dm;
61942e3340SBarry Smith   DMSNES_DA      *dmdasnes = (DMSNES_DA*)ctx;
626cab3a1bSJed Brown   DMDALocalInfo  info;
636cab3a1bSJed Brown   Vec            Xloc;
646cab3a1bSJed Brown   void           *x,*f;
656cab3a1bSJed Brown 
666cab3a1bSJed Brown   PetscFunctionBegin;
672961e153SJed Brown   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
682961e153SJed Brown   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
692961e153SJed Brown   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
70ce94432eSBarry Smith   if (!dmdasnes->residuallocal) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context");
716cab3a1bSJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
726cab3a1bSJed Brown   ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
736cab3a1bSJed Brown   ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
746cab3a1bSJed Brown   ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
756cab3a1bSJed Brown   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
766cab3a1bSJed Brown   ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
77709ac0d2SJed Brown   switch (dmdasnes->residuallocalimode) {
78709ac0d2SJed Brown   case INSERT_VALUES: {
796cab3a1bSJed Brown     ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr);
8062cd874fSBarry Smith     ierr = PetscLogEventBegin(SNES_FunctionEval,snes,X,F,0);CHKERRQ(ierr);
816cab3a1bSJed Brown     CHKMEMQ;
826cab3a1bSJed Brown     ierr = (*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx);CHKERRQ(ierr);
836cab3a1bSJed Brown     CHKMEMQ;
8462cd874fSBarry Smith     ierr = PetscLogEventEnd(SNES_FunctionEval,snes,X,F,0);CHKERRQ(ierr);
856cab3a1bSJed Brown     ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr);
86709ac0d2SJed Brown   } break;
87709ac0d2SJed Brown   case ADD_VALUES: {
88709ac0d2SJed Brown     Vec Floc;
89709ac0d2SJed Brown     ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr);
90709ac0d2SJed Brown     ierr = VecZeroEntries(Floc);CHKERRQ(ierr);
91709ac0d2SJed Brown     ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr);
9262cd874fSBarry Smith     ierr = PetscLogEventBegin(SNES_FunctionEval,snes,X,F,0);CHKERRQ(ierr);
93709ac0d2SJed Brown     CHKMEMQ;
94709ac0d2SJed Brown     ierr = (*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx);CHKERRQ(ierr);
95709ac0d2SJed Brown     CHKMEMQ;
9662cd874fSBarry Smith     ierr = PetscLogEventEnd(SNES_FunctionEval,snes,X,F,0);CHKERRQ(ierr);
97709ac0d2SJed Brown     ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr);
98709ac0d2SJed Brown     ierr = VecZeroEntries(F);CHKERRQ(ierr);
99709ac0d2SJed Brown     ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
100709ac0d2SJed Brown     ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
101709ac0d2SJed Brown     ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr);
102709ac0d2SJed Brown   } break;
103*98921bdaSJacob Faibussowitsch   default: SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdasnes->residuallocalimode);
104709ac0d2SJed Brown   }
105709ac0d2SJed Brown   ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
1066cab3a1bSJed Brown   ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
107422a814eSBarry Smith   if (snes->domainerror) {
108422a814eSBarry Smith     ierr = VecSetInf(F);CHKERRQ(ierr);
109422a814eSBarry Smith   }
1106cab3a1bSJed Brown   PetscFunctionReturn(0);
1116cab3a1bSJed Brown }
1126cab3a1bSJed Brown 
1132a4ee8f2SPeter Brune static PetscErrorCode SNESComputeObjective_DMDA(SNES snes,Vec X,PetscReal *ob,void *ctx)
1142a4ee8f2SPeter Brune {
1152a4ee8f2SPeter Brune   PetscErrorCode ierr;
1162a4ee8f2SPeter Brune   DM             dm;
117942e3340SBarry Smith   DMSNES_DA      *dmdasnes = (DMSNES_DA*)ctx;
1182a4ee8f2SPeter Brune   DMDALocalInfo  info;
1192a4ee8f2SPeter Brune   Vec            Xloc;
1202a4ee8f2SPeter Brune   void           *x;
1212a4ee8f2SPeter Brune 
1222a4ee8f2SPeter Brune   PetscFunctionBegin;
1232a4ee8f2SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1242a4ee8f2SPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
1252a4ee8f2SPeter Brune   PetscValidPointer(ob,3);
126ce94432eSBarry Smith   if (!dmdasnes->objectivelocal) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context");
1272a4ee8f2SPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
1282a4ee8f2SPeter Brune   ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
1292a4ee8f2SPeter Brune   ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
1302a4ee8f2SPeter Brune   ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
1312a4ee8f2SPeter Brune   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
1322a4ee8f2SPeter Brune   ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
1332a4ee8f2SPeter Brune   CHKMEMQ;
1342a4ee8f2SPeter Brune   ierr = (*dmdasnes->objectivelocal)(&info,x,ob,dmdasnes->objectivelocalctx);CHKERRQ(ierr);
1352a4ee8f2SPeter Brune   CHKMEMQ;
1362a4ee8f2SPeter Brune   ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
1372a4ee8f2SPeter Brune   ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
1382a4ee8f2SPeter Brune   PetscFunctionReturn(0);
1392a4ee8f2SPeter Brune }
1402a4ee8f2SPeter Brune 
141cc2e6a90SBarry Smith /* Routine is called by example, hence must be labeled PETSC_EXTERN */
142cc2e6a90SBarry Smith PETSC_EXTERN PetscErrorCode SNESComputeJacobian_DMDA(SNES snes,Vec X,Mat A,Mat B,void *ctx)
1432961e153SJed Brown {
1442961e153SJed Brown   PetscErrorCode ierr;
1452961e153SJed Brown   DM             dm;
146942e3340SBarry Smith   DMSNES_DA      *dmdasnes = (DMSNES_DA*)ctx;
1472961e153SJed Brown   DMDALocalInfo  info;
1482961e153SJed Brown   Vec            Xloc;
1492961e153SJed Brown   void           *x;
1502961e153SJed Brown 
1512961e153SJed Brown   PetscFunctionBegin;
152ce94432eSBarry Smith   if (!dmdasnes->residuallocal) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context");
1532961e153SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
1542961e153SJed Brown 
1552961e153SJed Brown   if (dmdasnes->jacobianlocal) {
1562961e153SJed Brown     ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
1572961e153SJed Brown     ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
1582961e153SJed Brown     ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
1592961e153SJed Brown     ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
1602961e153SJed Brown     ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
1612961e153SJed Brown     CHKMEMQ;
162d1e9a80fSBarry Smith     ierr = (*dmdasnes->jacobianlocal)(&info,x,A,B,dmdasnes->jacobianlocalctx);CHKERRQ(ierr);
1632961e153SJed Brown     CHKMEMQ;
1642961e153SJed Brown     ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
1652961e153SJed Brown     ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
1662961e153SJed Brown   } else {
1672961e153SJed Brown     MatFDColoring fdcoloring;
1682961e153SJed Brown     ierr = PetscObjectQuery((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject*)&fdcoloring);CHKERRQ(ierr);
1692961e153SJed Brown     if (!fdcoloring) {
1702961e153SJed Brown       ISColoring coloring;
1712961e153SJed Brown 
172b412c318SBarry Smith       ierr = DMCreateColoring(dm,dm->coloringtype,&coloring);CHKERRQ(ierr);
17394ab13aaSBarry Smith       ierr = MatFDColoringCreate(B,coloring,&fdcoloring);CHKERRQ(ierr);
1742961e153SJed Brown       switch (dm->coloringtype) {
1752961e153SJed Brown       case IS_COLORING_GLOBAL:
1762961e153SJed Brown         ierr = MatFDColoringSetFunction(fdcoloring,(PetscErrorCode (*)(void))SNESComputeFunction_DMDA,dmdasnes);CHKERRQ(ierr);
1772961e153SJed Brown         break;
178*98921bdaSJacob Faibussowitsch       default: SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"No support for coloring type '%s'",ISColoringTypes[dm->coloringtype]);
1792961e153SJed Brown       }
1802961e153SJed Brown       ierr = PetscObjectSetOptionsPrefix((PetscObject)fdcoloring,((PetscObject)dm)->prefix);CHKERRQ(ierr);
1812961e153SJed Brown       ierr = MatFDColoringSetFromOptions(fdcoloring);CHKERRQ(ierr);
18294ab13aaSBarry Smith       ierr = MatFDColoringSetUp(B,coloring,fdcoloring);CHKERRQ(ierr);
183f86b9fbaSHong Zhang       ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr);
1842961e153SJed Brown       ierr = PetscObjectCompose((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject)fdcoloring);CHKERRQ(ierr);
1852961e153SJed Brown       ierr = PetscObjectDereference((PetscObject)fdcoloring);CHKERRQ(ierr);
1862961e153SJed Brown 
1872961e153SJed Brown       /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call
1882961e153SJed Brown        * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the
1892961e153SJed Brown        * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually
1902961e153SJed Brown        * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be
191140e18c1SBarry Smith        * taken when the PetscObjectList for the Vec inside MatFDColoring is destroyed.
1922961e153SJed Brown        */
1932961e153SJed Brown       ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
1942961e153SJed Brown     }
195d1e9a80fSBarry Smith     ierr = MatFDColoringApply(B,fdcoloring,X,snes);CHKERRQ(ierr);
1962961e153SJed Brown   }
1972961e153SJed Brown   /* This will be redundant if the user called both, but it's too common to forget. */
19894ab13aaSBarry Smith   if (A != B) {
19994ab13aaSBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
20094ab13aaSBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2012961e153SJed Brown   }
2022961e153SJed Brown   PetscFunctionReturn(0);
2032961e153SJed Brown }
2042961e153SJed Brown 
2056cab3a1bSJed Brown /*@C
2066cab3a1bSJed Brown    DMDASNESSetFunctionLocal - set a local residual evaluation function
2076cab3a1bSJed Brown 
2086cab3a1bSJed Brown    Logically Collective
2096cab3a1bSJed Brown 
2104165533cSJose E. Roman    Input Parameters:
2116cab3a1bSJed Brown +  dm - DM to associate callback with
212e3c51ba2SJed Brown .  imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part
2136cab3a1bSJed Brown .  func - local residual evaluation
2146cab3a1bSJed Brown -  ctx - optional context for local residual evaluation
2156cab3a1bSJed Brown 
2160038884cSBarry Smith    Calling sequence:
2170038884cSBarry Smith    For PetscErrorCode (*func)(DMDALocalInfo *info,void *x, void *f, void *ctx),
2186cab3a1bSJed Brown +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
2190038884cSBarry Smith .  x - dimensional pointer to state at which to evaluate residual (e.g. PetscScalar *x or **x or ***x)
2200038884cSBarry Smith .  f - dimensional pointer to residual, write the residual here (e.g. PetscScalar *f or **f or ***f)
2216cab3a1bSJed Brown -  ctx - optional context passed above
2226cab3a1bSJed Brown 
2236cab3a1bSJed Brown    Level: beginner
2246cab3a1bSJed Brown 
2250038884cSBarry Smith .seealso: DMDASNESSetJacobianLocal(), DMSNESSetFunction(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
2266cab3a1bSJed Brown @*/
227709ac0d2SJed Brown PetscErrorCode DMDASNESSetFunctionLocal(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,void*,void*,void*),void *ctx)
2286cab3a1bSJed Brown {
2296cab3a1bSJed Brown   PetscErrorCode ierr;
230942e3340SBarry Smith   DMSNES         sdm;
231942e3340SBarry Smith   DMSNES_DA      *dmdasnes;
2326cab3a1bSJed Brown 
2336cab3a1bSJed Brown   PetscFunctionBegin;
2346cab3a1bSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
235942e3340SBarry Smith   ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr);
2366cab3a1bSJed Brown   ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr);
2371aa26658SKarl Rupp 
238709ac0d2SJed Brown   dmdasnes->residuallocalimode = imode;
2396cab3a1bSJed Brown   dmdasnes->residuallocal      = func;
2406cab3a1bSJed Brown   dmdasnes->residuallocalctx   = ctx;
2411aa26658SKarl Rupp 
2426cab3a1bSJed Brown   ierr = DMSNESSetFunction(dm,SNESComputeFunction_DMDA,dmdasnes);CHKERRQ(ierr);
24322c6f798SBarry Smith   if (!sdm->ops->computejacobian) {  /* Call us for the Jacobian too, can be overridden by the user. */
2442961e153SJed Brown     ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes);CHKERRQ(ierr);
2452961e153SJed Brown   }
2462961e153SJed Brown   PetscFunctionReturn(0);
2472961e153SJed Brown }
2482961e153SJed Brown 
2492961e153SJed Brown /*@C
250e79ed307SJed Brown    DMDASNESSetJacobianLocal - set a local Jacobian evaluation function
2512961e153SJed Brown 
2522961e153SJed Brown    Logically Collective
2532961e153SJed Brown 
2544165533cSJose E. Roman    Input Parameters:
2552961e153SJed Brown +  dm - DM to associate callback with
2566b7eb5ceSMatthew G. Knepley .  func - local Jacobian evaluation
2576b7eb5ceSMatthew G. Knepley -  ctx - optional context for local Jacobian evaluation
2582961e153SJed Brown 
2590038884cSBarry Smith    Calling sequence:
2600038884cSBarry Smith    For PetscErrorCode (*func)(DMDALocalInfo *info,void *x,Mat J,Mat M,void *ctx),
2616b7eb5ceSMatthew G. Knepley +  info - DMDALocalInfo defining the subdomain to evaluate the Jacobian at
2620038884cSBarry Smith .  x - dimensional pointer to state at which to evaluate Jacobian (e.g. PetscScalar *x or **x or ***x)
2636b7eb5ceSMatthew G. Knepley .  J - Mat object for the Jacobian
2646b7eb5ceSMatthew G. Knepley .  M - Mat object for the Jacobian preconditioner matrix
2652961e153SJed Brown -  ctx - optional context passed above
2662961e153SJed Brown 
2672961e153SJed Brown    Level: beginner
2682961e153SJed Brown 
2690038884cSBarry Smith .seealso: DMDASNESSetFunctionLocal(), DMSNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
2702961e153SJed Brown @*/
271d1e9a80fSBarry Smith PetscErrorCode DMDASNESSetJacobianLocal(DM dm,PetscErrorCode (*func)(DMDALocalInfo*,void*,Mat,Mat,void*),void *ctx)
2722961e153SJed Brown {
2732961e153SJed Brown   PetscErrorCode ierr;
274942e3340SBarry Smith   DMSNES         sdm;
275942e3340SBarry Smith   DMSNES_DA      *dmdasnes;
2762961e153SJed Brown 
2772961e153SJed Brown   PetscFunctionBegin;
2782961e153SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
279942e3340SBarry Smith   ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr);
2802961e153SJed Brown   ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr);
2811aa26658SKarl Rupp 
2822961e153SJed Brown   dmdasnes->jacobianlocal    = func;
2832961e153SJed Brown   dmdasnes->jacobianlocalctx = ctx;
2841aa26658SKarl Rupp 
2852961e153SJed Brown   ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes);CHKERRQ(ierr);
2866cab3a1bSJed Brown   PetscFunctionReturn(0);
2876cab3a1bSJed Brown }
2882a4ee8f2SPeter Brune 
2892a4ee8f2SPeter Brune /*@C
2902a4ee8f2SPeter Brune    DMDASNESSetObjectiveLocal - set a local residual evaluation function
2912a4ee8f2SPeter Brune 
2922a4ee8f2SPeter Brune    Logically Collective
2932a4ee8f2SPeter Brune 
2944165533cSJose E. Roman    Input Parameters:
2952a4ee8f2SPeter Brune +  dm - DM to associate callback with
2962a4ee8f2SPeter Brune .  func - local objective evaluation
2972a4ee8f2SPeter Brune -  ctx - optional context for local residual evaluation
2982a4ee8f2SPeter Brune 
2992a4ee8f2SPeter Brune    Calling sequence for func:
3002a4ee8f2SPeter Brune +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
3012a4ee8f2SPeter Brune .  x - dimensional pointer to state at which to evaluate residual
3022a4ee8f2SPeter Brune .  ob - eventual objective value
3032a4ee8f2SPeter Brune -  ctx - optional context passed above
3042a4ee8f2SPeter Brune 
3052a4ee8f2SPeter Brune    Level: beginner
3062a4ee8f2SPeter Brune 
3070038884cSBarry Smith .seealso: DMSNESSetFunction(), DMDASNESSetJacobianLocal(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
3082a4ee8f2SPeter Brune @*/
3092a4ee8f2SPeter Brune PetscErrorCode DMDASNESSetObjectiveLocal(DM dm,DMDASNESObjective func,void *ctx)
3102a4ee8f2SPeter Brune {
3112a4ee8f2SPeter Brune   PetscErrorCode ierr;
312942e3340SBarry Smith   DMSNES         sdm;
313942e3340SBarry Smith   DMSNES_DA      *dmdasnes;
3142a4ee8f2SPeter Brune 
3152a4ee8f2SPeter Brune   PetscFunctionBegin;
3162a4ee8f2SPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
317942e3340SBarry Smith   ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr);
3182a4ee8f2SPeter Brune   ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr);
3191aa26658SKarl Rupp 
3202a4ee8f2SPeter Brune   dmdasnes->objectivelocal    = func;
3212a4ee8f2SPeter Brune   dmdasnes->objectivelocalctx = ctx;
3221aa26658SKarl Rupp 
3232a4ee8f2SPeter Brune   ierr = DMSNESSetObjective(dm,SNESComputeObjective_DMDA,dmdasnes);CHKERRQ(ierr);
3242a4ee8f2SPeter Brune   PetscFunctionReturn(0);
3252a4ee8f2SPeter Brune }
326dcad5f1cSBarry Smith 
327dcad5f1cSBarry Smith static PetscErrorCode SNESComputePicard_DMDA(SNES snes,Vec X,Vec F,void *ctx)
328dcad5f1cSBarry Smith {
329dcad5f1cSBarry Smith   PetscErrorCode ierr;
330dcad5f1cSBarry Smith   DM             dm;
331942e3340SBarry Smith   DMSNES_DA      *dmdasnes = (DMSNES_DA*)ctx;
332dcad5f1cSBarry Smith   DMDALocalInfo  info;
333dcad5f1cSBarry Smith   Vec            Xloc;
334dcad5f1cSBarry Smith   void           *x,*f;
335dcad5f1cSBarry Smith 
336dcad5f1cSBarry Smith   PetscFunctionBegin;
337dcad5f1cSBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
338dcad5f1cSBarry Smith   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
339dcad5f1cSBarry Smith   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
340ce94432eSBarry Smith   if (!dmdasnes->rhsplocal) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context");
341dcad5f1cSBarry Smith   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
342dcad5f1cSBarry Smith   ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
343dcad5f1cSBarry Smith   ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
344dcad5f1cSBarry Smith   ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
345dcad5f1cSBarry Smith   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
346dcad5f1cSBarry Smith   ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
347dcad5f1cSBarry Smith   switch (dmdasnes->residuallocalimode) {
348dcad5f1cSBarry Smith   case INSERT_VALUES: {
349dcad5f1cSBarry Smith     ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr);
350dcad5f1cSBarry Smith     CHKMEMQ;
351dcad5f1cSBarry Smith     ierr = (*dmdasnes->rhsplocal)(&info,x,f,dmdasnes->picardlocalctx);CHKERRQ(ierr);
352dcad5f1cSBarry Smith     CHKMEMQ;
353dcad5f1cSBarry Smith     ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr);
354dcad5f1cSBarry Smith   } break;
355dcad5f1cSBarry Smith   case ADD_VALUES: {
356dcad5f1cSBarry Smith     Vec Floc;
357dcad5f1cSBarry Smith     ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr);
358dcad5f1cSBarry Smith     ierr = VecZeroEntries(Floc);CHKERRQ(ierr);
359dcad5f1cSBarry Smith     ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr);
360dcad5f1cSBarry Smith     CHKMEMQ;
361dcad5f1cSBarry Smith     ierr = (*dmdasnes->rhsplocal)(&info,x,f,dmdasnes->picardlocalctx);CHKERRQ(ierr);
362dcad5f1cSBarry Smith     CHKMEMQ;
363dcad5f1cSBarry Smith     ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr);
364dcad5f1cSBarry Smith     ierr = VecZeroEntries(F);CHKERRQ(ierr);
365dcad5f1cSBarry Smith     ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
366dcad5f1cSBarry Smith     ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
367dcad5f1cSBarry Smith     ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr);
368dcad5f1cSBarry Smith   } break;
369*98921bdaSJacob Faibussowitsch   default: SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdasnes->residuallocalimode);
370dcad5f1cSBarry Smith   }
371dcad5f1cSBarry Smith   ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
372dcad5f1cSBarry Smith   ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
373dcad5f1cSBarry Smith   PetscFunctionReturn(0);
374dcad5f1cSBarry Smith }
375dcad5f1cSBarry Smith 
376d1e9a80fSBarry Smith static PetscErrorCode SNESComputePicardJacobian_DMDA(SNES snes,Vec X,Mat A,Mat B,void *ctx)
377dcad5f1cSBarry Smith {
378dcad5f1cSBarry Smith   PetscErrorCode ierr;
379dcad5f1cSBarry Smith   DM             dm;
380942e3340SBarry Smith   DMSNES_DA      *dmdasnes = (DMSNES_DA*)ctx;
381dcad5f1cSBarry Smith   DMDALocalInfo  info;
382dcad5f1cSBarry Smith   Vec            Xloc;
383dcad5f1cSBarry Smith   void           *x;
384dcad5f1cSBarry Smith 
385dcad5f1cSBarry Smith   PetscFunctionBegin;
386ce94432eSBarry Smith   if (!dmdasnes->jacobianplocal) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context");
387dcad5f1cSBarry Smith   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
388dcad5f1cSBarry Smith 
389dcad5f1cSBarry Smith   ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
390dcad5f1cSBarry Smith   ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
391dcad5f1cSBarry Smith   ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
392dcad5f1cSBarry Smith   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
393dcad5f1cSBarry Smith   ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
394dcad5f1cSBarry Smith   CHKMEMQ;
395d1e9a80fSBarry Smith   ierr = (*dmdasnes->jacobianplocal)(&info,x,A,B,dmdasnes->picardlocalctx);CHKERRQ(ierr);
396dcad5f1cSBarry Smith   CHKMEMQ;
397dcad5f1cSBarry Smith   ierr  = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
398dcad5f1cSBarry Smith   ierr  = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
39994ab13aaSBarry Smith   if (A != B) {
40094ab13aaSBarry Smith     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
40194ab13aaSBarry Smith     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
402dcad5f1cSBarry Smith   }
403dcad5f1cSBarry Smith   PetscFunctionReturn(0);
404dcad5f1cSBarry Smith }
405dcad5f1cSBarry Smith 
406dcad5f1cSBarry Smith /*@C
407dcad5f1cSBarry Smith    DMDASNESSetPicardLocal - set a local right hand side and matrix evaluation function for Picard iteration
408dcad5f1cSBarry Smith 
409dcad5f1cSBarry Smith    Logically Collective
410dcad5f1cSBarry Smith 
4114165533cSJose E. Roman    Input Parameters:
412dcad5f1cSBarry Smith +  dm - DM to associate callback with
413dcad5f1cSBarry Smith .  imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part
414dcad5f1cSBarry Smith .  func - local residual evaluation
415dcad5f1cSBarry Smith -  ctx - optional context for local residual evaluation
416dcad5f1cSBarry Smith 
417dcad5f1cSBarry Smith    Calling sequence for func:
418dcad5f1cSBarry Smith +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
419dcad5f1cSBarry Smith .  x - dimensional pointer to state at which to evaluate residual
420dcad5f1cSBarry Smith .  f - dimensional pointer to residual, write the residual here
421dcad5f1cSBarry Smith -  ctx - optional context passed above
422dcad5f1cSBarry Smith 
42395452b02SPatrick Sanan    Notes:
42495452b02SPatrick Sanan     The user must use
4250298fd71SBarry Smith     ierr = SNESSetFunction(snes,NULL,SNESPicardComputeFunction,&user);CHKERRQ(ierr);
426dcad5f1cSBarry Smith     in their code before calling this routine.
427dcad5f1cSBarry Smith 
428dcad5f1cSBarry Smith    Level: beginner
429dcad5f1cSBarry Smith 
430dcad5f1cSBarry Smith .seealso: DMSNESSetFunction(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
431dcad5f1cSBarry Smith @*/
432dcad5f1cSBarry Smith PetscErrorCode DMDASNESSetPicardLocal(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,void*,void*,void*),
433d1e9a80fSBarry Smith                                       PetscErrorCode (*jac)(DMDALocalInfo*,void*,Mat,Mat,void*),void *ctx)
434dcad5f1cSBarry Smith {
435dcad5f1cSBarry Smith   PetscErrorCode ierr;
436942e3340SBarry Smith   DMSNES         sdm;
437942e3340SBarry Smith   DMSNES_DA      *dmdasnes;
438dcad5f1cSBarry Smith 
439dcad5f1cSBarry Smith   PetscFunctionBegin;
440dcad5f1cSBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
441942e3340SBarry Smith   ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr);
442dcad5f1cSBarry Smith   ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr);
4431aa26658SKarl Rupp 
444dcad5f1cSBarry Smith   dmdasnes->residuallocalimode = imode;
445dcad5f1cSBarry Smith   dmdasnes->rhsplocal          = func;
446dcad5f1cSBarry Smith   dmdasnes->jacobianplocal     = jac;
447dcad5f1cSBarry Smith   dmdasnes->picardlocalctx     = ctx;
4481aa26658SKarl Rupp 
449dcad5f1cSBarry Smith   ierr = DMSNESSetPicard(dm,SNESComputePicard_DMDA,SNESComputePicardJacobian_DMDA,dmdasnes);CHKERRQ(ierr);
450bbc1464cSBarry Smith   ierr = DMSNESSetMFFunction(dm,SNESComputeFunction_DMDA,dmdasnes);CHKERRQ(ierr);
451dcad5f1cSBarry Smith   PetscFunctionReturn(0);
452dcad5f1cSBarry Smith }
453