xref: /petsc/src/snes/utils/dmdasnes.c (revision 942e334039a1a47ff909c032afdf0581bd880e9f)
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*);
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*);
17dcad5f1cSBarry Smith   PetscErrorCode (*jacobianplocal)(DMDALocalInfo*,void*,Mat,Mat,MatStructure*,void*);
18dcad5f1cSBarry Smith   void *picardlocalctx;
19*942e3340SBarry Smith } DMSNES_DA;
206cab3a1bSJed Brown 
216cab3a1bSJed Brown #undef __FUNCT__
22*942e3340SBarry Smith #define __FUNCT__ "DMSNESDestroy_DMDA"
23*942e3340SBarry 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__
33*942e3340SBarry Smith #define __FUNCT__ "DMSNESDuplicate_DMDA"
34*942e3340SBarry Smith static PetscErrorCode DMSNESDuplicate_DMDA(DMSNES oldsdm,DM dm)
35e3c0b8a2SPeter Brune {
36e3c0b8a2SPeter Brune   PetscErrorCode ierr;
37*942e3340SBarry Smith   DMSNES         sdm;
38e3c0b8a2SPeter Brune 
39e3c0b8a2SPeter Brune   PetscFunctionBegin;
40*942e3340SBarry Smith   ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr);
41*942e3340SBarry Smith   ierr = PetscNewLog(dm,DMSNES_DA ,&sdm->data);CHKERRQ(ierr);
42dcad5f1cSBarry Smith   if (oldsdm->data) {
43*942e3340SBarry Smith     ierr = PetscMemcpy(sdm->data,oldsdm->data,sizeof(DMSNES_DA ));CHKERRQ(ierr);
44dcad5f1cSBarry Smith   }
45e3c0b8a2SPeter Brune   PetscFunctionReturn(0);
46e3c0b8a2SPeter Brune }
47e3c0b8a2SPeter Brune 
48e3c0b8a2SPeter Brune 
49e3c0b8a2SPeter Brune #undef __FUNCT__
506cab3a1bSJed Brown #define __FUNCT__ "DMDASNESGetContext"
51*942e3340SBarry Smith static PetscErrorCode DMDASNESGetContext(DM dm,DMSNES sdm,DMSNES_DA  **dmdasnes)
526cab3a1bSJed Brown {
536cab3a1bSJed Brown   PetscErrorCode ierr;
546cab3a1bSJed Brown 
556cab3a1bSJed Brown   PetscFunctionBegin;
56b8d0cc33SJed Brown   *dmdasnes = PETSC_NULL;
576cab3a1bSJed Brown   if (!sdm->data) {
58*942e3340SBarry Smith     ierr = PetscNewLog(dm,DMSNES_DA ,&sdm->data);CHKERRQ(ierr);
59*942e3340SBarry Smith     sdm->destroy   = DMSNESDestroy_DMDA;
60*942e3340SBarry Smith     sdm->duplicate = DMSNESDuplicate_DMDA;
616cab3a1bSJed Brown   }
62*942e3340SBarry Smith   *dmdasnes = (DMSNES_DA *)sdm->data;
636cab3a1bSJed Brown   PetscFunctionReturn(0);
646cab3a1bSJed Brown }
656cab3a1bSJed Brown 
666cab3a1bSJed Brown #undef __FUNCT__
676cab3a1bSJed Brown #define __FUNCT__ "SNESComputeFunction_DMDA"
686cab3a1bSJed Brown /*
696cab3a1bSJed Brown   This function should eventually replace:
706cab3a1bSJed Brown     DMDAComputeFunction() and DMDAComputeFunction1()
716cab3a1bSJed Brown  */
726cab3a1bSJed Brown static PetscErrorCode SNESComputeFunction_DMDA(SNES snes,Vec X,Vec F,void *ctx)
736cab3a1bSJed Brown {
746cab3a1bSJed Brown   PetscErrorCode ierr;
756cab3a1bSJed Brown   DM             dm;
76*942e3340SBarry Smith   DMSNES_DA      *dmdasnes = (DMSNES_DA *)ctx;
776cab3a1bSJed Brown   DMDALocalInfo  info;
786cab3a1bSJed Brown   Vec            Xloc;
796cab3a1bSJed Brown   void           *x,*f;
806cab3a1bSJed Brown 
816cab3a1bSJed Brown   PetscFunctionBegin;
822961e153SJed Brown   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
832961e153SJed Brown   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
842961e153SJed Brown   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
856cab3a1bSJed Brown   if (!dmdasnes->residuallocal) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Corrupt context");
866cab3a1bSJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
876cab3a1bSJed Brown   ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
886cab3a1bSJed Brown   ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
896cab3a1bSJed Brown   ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
906cab3a1bSJed Brown   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
916cab3a1bSJed Brown   ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
92709ac0d2SJed Brown   switch (dmdasnes->residuallocalimode) {
93709ac0d2SJed Brown   case INSERT_VALUES: {
946cab3a1bSJed Brown     ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr);
956cab3a1bSJed Brown     CHKMEMQ;
966cab3a1bSJed Brown     ierr = (*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx);CHKERRQ(ierr);
976cab3a1bSJed Brown     CHKMEMQ;
986cab3a1bSJed Brown     ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr);
99709ac0d2SJed Brown   } break;
100709ac0d2SJed Brown   case ADD_VALUES: {
101709ac0d2SJed Brown     Vec Floc;
102709ac0d2SJed Brown     ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr);
103709ac0d2SJed Brown     ierr = VecZeroEntries(Floc);CHKERRQ(ierr);
104709ac0d2SJed Brown     ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr);
105709ac0d2SJed Brown     CHKMEMQ;
106709ac0d2SJed Brown     ierr = (*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx);CHKERRQ(ierr);
107709ac0d2SJed Brown     CHKMEMQ;
108709ac0d2SJed Brown     ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr);
109709ac0d2SJed Brown     ierr = VecZeroEntries(F);CHKERRQ(ierr);
110709ac0d2SJed Brown     ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
111709ac0d2SJed Brown     ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
112709ac0d2SJed Brown     ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr);
113709ac0d2SJed Brown   } break;
114709ac0d2SJed Brown   default: SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdasnes->residuallocalimode);
115709ac0d2SJed Brown   }
116709ac0d2SJed Brown   ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
1176cab3a1bSJed Brown   ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
1186cab3a1bSJed Brown   PetscFunctionReturn(0);
1196cab3a1bSJed Brown }
1206cab3a1bSJed Brown 
1216cab3a1bSJed Brown #undef __FUNCT__
1222a4ee8f2SPeter Brune #define __FUNCT__ "SNESComputeObjective_DMDA"
1232a4ee8f2SPeter Brune /*
1242a4ee8f2SPeter Brune   This function should eventually replace:
1252a4ee8f2SPeter Brune     DMDAComputeFunction() and DMDAComputeFunction1()
1262a4ee8f2SPeter Brune  */
1272a4ee8f2SPeter Brune static PetscErrorCode SNESComputeObjective_DMDA(SNES snes,Vec X,PetscReal *ob,void *ctx)
1282a4ee8f2SPeter Brune {
1292a4ee8f2SPeter Brune   PetscErrorCode ierr;
1302a4ee8f2SPeter Brune   DM             dm;
131*942e3340SBarry Smith   DMSNES_DA      *dmdasnes = (DMSNES_DA *)ctx;
1322a4ee8f2SPeter Brune   DMDALocalInfo  info;
1332a4ee8f2SPeter Brune   Vec            Xloc;
1342a4ee8f2SPeter Brune   void           *x;
1352a4ee8f2SPeter Brune 
1362a4ee8f2SPeter Brune   PetscFunctionBegin;
1372a4ee8f2SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1382a4ee8f2SPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
1392a4ee8f2SPeter Brune   PetscValidPointer(ob,3);
1402a4ee8f2SPeter Brune   if (!dmdasnes->objectivelocal) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Corrupt context");
1412a4ee8f2SPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
1422a4ee8f2SPeter Brune   ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
1432a4ee8f2SPeter Brune   ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
1442a4ee8f2SPeter Brune   ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
1452a4ee8f2SPeter Brune   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
1462a4ee8f2SPeter Brune   ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
1472a4ee8f2SPeter Brune   CHKMEMQ;
1482a4ee8f2SPeter Brune   ierr = (*dmdasnes->objectivelocal)(&info,x,ob,dmdasnes->objectivelocalctx);CHKERRQ(ierr);
1492a4ee8f2SPeter Brune   CHKMEMQ;
1502a4ee8f2SPeter Brune   ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
1512a4ee8f2SPeter Brune   ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
1522a4ee8f2SPeter Brune   PetscFunctionReturn(0);
1532a4ee8f2SPeter Brune }
1542a4ee8f2SPeter Brune 
1552a4ee8f2SPeter Brune #undef __FUNCT__
1566427ad5bSPeter Brune #define __FUNCT__ "SNESComputeLocalBlockFunction_DMDA"
1576427ad5bSPeter Brune /*
1586427ad5bSPeter Brune   This function should eventually replace:
1596427ad5bSPeter Brune     DMDAComputeFunction() and DMDAComputeFunction1()
1606427ad5bSPeter Brune  */
1616427ad5bSPeter Brune PetscErrorCode SNESComputeLocalBlockFunction_DMDA(SNES snes,Vec X,Vec F,void *dmctx)
1626427ad5bSPeter Brune {
1636427ad5bSPeter Brune   PetscErrorCode ierr;
1646427ad5bSPeter Brune   DM             dm = (DM)dmctx; /* the global DMDA context */
1656427ad5bSPeter Brune   DMDALocalInfo  info;
1666427ad5bSPeter Brune   void           *x,*f;
167*942e3340SBarry Smith   DMSNES_DA      *dmdasnes;
168*942e3340SBarry Smith   DMSNES         sdm;
1696427ad5bSPeter Brune 
1706427ad5bSPeter Brune   PetscFunctionBegin;
1716427ad5bSPeter Brune 
1726427ad5bSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1736427ad5bSPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
1746427ad5bSPeter Brune   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
1756427ad5bSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,4);
1766427ad5bSPeter Brune 
177*942e3340SBarry Smith   ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr);
1786427ad5bSPeter Brune   ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr);
179ec31c420SJed Brown   if (!dmdasnes->residuallocal) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Corrupt context");
1806427ad5bSPeter Brune   ierr = DMDAVecGetArray(dm,X,&x);CHKERRQ(ierr);
1816427ad5bSPeter Brune   ierr = DMDAGetLocalBlockInfo(dm,&info);CHKERRQ(ierr);
1826427ad5bSPeter Brune   switch (dmdasnes->residuallocalimode) {
1836427ad5bSPeter Brune   case INSERT_VALUES: {
1846427ad5bSPeter Brune     ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr);
1856427ad5bSPeter Brune     CHKMEMQ;
1866427ad5bSPeter Brune     ierr = (*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx);CHKERRQ(ierr);
1876427ad5bSPeter Brune     CHKMEMQ;
1886427ad5bSPeter Brune     ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr);
1896427ad5bSPeter Brune   } break;
1906427ad5bSPeter Brune   case ADD_VALUES: {
1916427ad5bSPeter Brune     Vec Floc;
1926427ad5bSPeter Brune     ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr);
1936427ad5bSPeter Brune     ierr = VecZeroEntries(Floc);CHKERRQ(ierr);
1946427ad5bSPeter Brune     ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr);
1956427ad5bSPeter Brune     CHKMEMQ;
1966427ad5bSPeter Brune     ierr = (*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx);CHKERRQ(ierr);
1976427ad5bSPeter Brune     CHKMEMQ;
1986427ad5bSPeter Brune     ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr);
1996427ad5bSPeter Brune     ierr = VecZeroEntries(F);CHKERRQ(ierr);
2006427ad5bSPeter Brune     ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
2016427ad5bSPeter Brune     ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
2026427ad5bSPeter Brune     ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr);
2036427ad5bSPeter Brune   } break;
2046427ad5bSPeter Brune   default: SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdasnes->residuallocalimode);
2056427ad5bSPeter Brune   }
2066427ad5bSPeter Brune   ierr = DMDAVecRestoreArray(dm,X,&x);CHKERRQ(ierr);
2076427ad5bSPeter Brune   PetscFunctionReturn(0);
2086427ad5bSPeter Brune }
2096427ad5bSPeter Brune 
2106427ad5bSPeter Brune #undef __FUNCT__
2112961e153SJed Brown #define __FUNCT__ "SNESComputeJacobian_DMDA"
2122961e153SJed Brown static PetscErrorCode SNESComputeJacobian_DMDA(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *ctx)
2132961e153SJed Brown {
2142961e153SJed Brown   PetscErrorCode ierr;
2152961e153SJed Brown   DM             dm;
216*942e3340SBarry Smith   DMSNES_DA      *dmdasnes = (DMSNES_DA *)ctx;
2172961e153SJed Brown   DMDALocalInfo  info;
2182961e153SJed Brown   Vec            Xloc;
2192961e153SJed Brown   void           *x;
2202961e153SJed Brown 
2212961e153SJed Brown   PetscFunctionBegin;
2222961e153SJed Brown   if (!dmdasnes->residuallocal) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Corrupt context");
2232961e153SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
2242961e153SJed Brown 
2252961e153SJed Brown   if (dmdasnes->jacobianlocal) {
2262961e153SJed Brown     ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
2272961e153SJed Brown     ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
2282961e153SJed Brown     ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
2292961e153SJed Brown     ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
2302961e153SJed Brown     ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
2312961e153SJed Brown     CHKMEMQ;
2322961e153SJed Brown     ierr = (*dmdasnes->jacobianlocal)(&info,x,*A,*B,mstr,dmdasnes->jacobianlocalctx);CHKERRQ(ierr);
2332961e153SJed Brown     CHKMEMQ;
2342961e153SJed Brown     ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
2352961e153SJed Brown     ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
2362961e153SJed Brown   } else {
2372961e153SJed Brown     MatFDColoring fdcoloring;
2382961e153SJed Brown     ierr = PetscObjectQuery((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject*)&fdcoloring);CHKERRQ(ierr);
2392961e153SJed Brown     if (!fdcoloring) {
2402961e153SJed Brown       ISColoring     coloring;
2412961e153SJed Brown 
242c97a5455SBarry Smith       ierr = DMCreateColoring(dm,dm->coloringtype,dm->mattype,&coloring);CHKERRQ(ierr);
2432961e153SJed Brown       ierr = MatFDColoringCreate(*B,coloring,&fdcoloring);CHKERRQ(ierr);
2442961e153SJed Brown       ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr);
2452961e153SJed Brown       switch (dm->coloringtype) {
2462961e153SJed Brown       case IS_COLORING_GLOBAL:
2472961e153SJed Brown         ierr = MatFDColoringSetFunction(fdcoloring,(PetscErrorCode(*)(void))SNESComputeFunction_DMDA,dmdasnes);CHKERRQ(ierr);
2482961e153SJed Brown         break;
2492961e153SJed Brown       default: SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"No support for coloring type '%s'",ISColoringTypes[dm->coloringtype]);
2502961e153SJed Brown       }
2512961e153SJed Brown       ierr = PetscObjectSetOptionsPrefix((PetscObject)fdcoloring,((PetscObject)dm)->prefix);CHKERRQ(ierr);
2522961e153SJed Brown       ierr = MatFDColoringSetFromOptions(fdcoloring);CHKERRQ(ierr);
2532961e153SJed Brown       ierr = PetscObjectCompose((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject)fdcoloring);CHKERRQ(ierr);
2542961e153SJed Brown       ierr = PetscObjectDereference((PetscObject)fdcoloring);CHKERRQ(ierr);
2552961e153SJed Brown 
2562961e153SJed Brown       /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call
2572961e153SJed Brown        * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the
2582961e153SJed Brown        * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually
2592961e153SJed Brown        * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be
2602961e153SJed Brown        * taken when the PetscOList for the Vec inside MatFDColoring is destroyed.
2612961e153SJed Brown        */
2622961e153SJed Brown       ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
2632961e153SJed Brown     }
2642961e153SJed Brown     *mstr = SAME_NONZERO_PATTERN;
2652961e153SJed Brown     ierr = MatFDColoringApply(*B,fdcoloring,X,mstr,snes);CHKERRQ(ierr);
2662961e153SJed Brown   }
2672961e153SJed Brown   /* This will be redundant if the user called both, but it's too common to forget. */
2682961e153SJed Brown   if (*A != *B) {
2692961e153SJed Brown     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2702961e153SJed Brown     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2712961e153SJed Brown   }
2722961e153SJed Brown   PetscFunctionReturn(0);
2732961e153SJed Brown }
2742961e153SJed Brown 
2752961e153SJed Brown #undef __FUNCT__
276dcad5f1cSBarry Smith #define __FUNCT__ "SNESComputeLocalBlockJacobian_DMDA"
2776427ad5bSPeter Brune PetscErrorCode SNESComputeLocalBlockJacobian_DMDA(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *dmctx)
2786427ad5bSPeter Brune {
2796427ad5bSPeter Brune   PetscErrorCode ierr;
2806427ad5bSPeter Brune   DM             dm = (DM)dmctx; /* the global DMDA context */
281*942e3340SBarry Smith   DMSNES_DA      *dmdasnes;
282*942e3340SBarry Smith   DMSNES         sdm;
2836427ad5bSPeter Brune   DMDALocalInfo  info;
2846427ad5bSPeter Brune   void           *x;
2856427ad5bSPeter Brune 
2866427ad5bSPeter Brune   PetscFunctionBegin;
287*942e3340SBarry Smith   ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr);
2886427ad5bSPeter Brune   ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr);
289ec31c420SJed Brown   if (!dmdasnes->residuallocal) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Corrupt context");
2906427ad5bSPeter Brune   if (dmdasnes->jacobianlocal) {
2916427ad5bSPeter Brune     ierr = DMDAVecGetArray(dm,X,&x);CHKERRQ(ierr);
2926427ad5bSPeter Brune     CHKMEMQ;
2936427ad5bSPeter Brune     ierr = (*dmdasnes->jacobianlocal)(&info,x,*A,*B,mstr,dmdasnes->jacobianlocalctx);CHKERRQ(ierr);
2946427ad5bSPeter Brune     CHKMEMQ;
2956427ad5bSPeter Brune     ierr = DMDAVecRestoreArray(dm,X,&x);CHKERRQ(ierr);
2966427ad5bSPeter Brune   }
2976427ad5bSPeter Brune   /* This will be redundant if the user called both, but it's too common to forget. */
2986427ad5bSPeter Brune   if (*A != *B) {
2996427ad5bSPeter Brune     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3006427ad5bSPeter Brune     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3016427ad5bSPeter Brune   }
3026427ad5bSPeter Brune   PetscFunctionReturn(0);
3036427ad5bSPeter Brune }
3046427ad5bSPeter Brune 
3056427ad5bSPeter Brune #undef __FUNCT__
3066cab3a1bSJed Brown #define __FUNCT__ "DMDASNESSetFunctionLocal"
3076cab3a1bSJed Brown /*@C
3086cab3a1bSJed Brown    DMDASNESSetFunctionLocal - set a local residual evaluation function
3096cab3a1bSJed Brown 
3106cab3a1bSJed Brown    Logically Collective
3116cab3a1bSJed Brown 
3126cab3a1bSJed Brown    Input Arguments:
3136cab3a1bSJed Brown +  dm - DM to associate callback with
314e3c51ba2SJed Brown .  imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part
3156cab3a1bSJed Brown .  func - local residual evaluation
3166cab3a1bSJed Brown -  ctx - optional context for local residual evaluation
3176cab3a1bSJed Brown 
3186cab3a1bSJed Brown    Calling sequence for func:
3196cab3a1bSJed Brown +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
3206cab3a1bSJed Brown .  x - dimensional pointer to state at which to evaluate residual
3216cab3a1bSJed Brown .  f - dimensional pointer to residual, write the residual here
3226cab3a1bSJed Brown -  ctx - optional context passed above
3236cab3a1bSJed Brown 
3246cab3a1bSJed Brown    Level: beginner
3256cab3a1bSJed Brown 
3266cab3a1bSJed Brown .seealso: DMSNESSetFunction(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
3276cab3a1bSJed Brown @*/
328709ac0d2SJed Brown PetscErrorCode DMDASNESSetFunctionLocal(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,void*,void*,void*),void *ctx)
3296cab3a1bSJed Brown {
3306cab3a1bSJed Brown   PetscErrorCode ierr;
331*942e3340SBarry Smith   DMSNES         sdm;
332*942e3340SBarry Smith   DMSNES_DA      *dmdasnes;
3336cab3a1bSJed Brown 
3346cab3a1bSJed Brown   PetscFunctionBegin;
3356cab3a1bSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
336*942e3340SBarry Smith   ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr);
3376cab3a1bSJed Brown   ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr);
338709ac0d2SJed Brown   dmdasnes->residuallocalimode = imode;
3396cab3a1bSJed Brown   dmdasnes->residuallocal = func;
3406cab3a1bSJed Brown   dmdasnes->residuallocalctx = ctx;
3416cab3a1bSJed Brown   ierr = DMSNESSetFunction(dm,SNESComputeFunction_DMDA,dmdasnes);CHKERRQ(ierr);
3426427ad5bSPeter Brune   ierr = DMSNESSetBlockFunction(dm,SNESComputeLocalBlockFunction_DMDA,dm);CHKERRQ(ierr);
3432961e153SJed Brown   if (!sdm->computejacobian) {  /* Call us for the Jacobian too, can be overridden by the user. */
3442961e153SJed Brown     ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes);CHKERRQ(ierr);
3452961e153SJed Brown   }
3462961e153SJed Brown   PetscFunctionReturn(0);
3472961e153SJed Brown }
3482961e153SJed Brown 
3492961e153SJed Brown #undef __FUNCT__
3502961e153SJed Brown #define __FUNCT__ "DMDASNESSetJacobianLocal"
3512961e153SJed Brown /*@C
3522961e153SJed Brown    DMDASNESSetJacobianLocal - set a local residual evaluation function
3532961e153SJed Brown 
3542961e153SJed Brown    Logically Collective
3552961e153SJed Brown 
3562961e153SJed Brown    Input Arguments:
3572961e153SJed Brown +  dm - DM to associate callback with
3582961e153SJed Brown .  func - local residual evaluation
3592961e153SJed Brown -  ctx - optional context for local residual evaluation
3602961e153SJed Brown 
3612961e153SJed Brown    Calling sequence for func:
3622961e153SJed Brown +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
3632961e153SJed Brown .  x - dimensional pointer to state at which to evaluate residual
3642961e153SJed Brown .  f - dimensional pointer to residual, write the residual here
3652961e153SJed Brown -  ctx - optional context passed above
3662961e153SJed Brown 
3672961e153SJed Brown    Level: beginner
3682961e153SJed Brown 
3692961e153SJed Brown .seealso: DMSNESSetJacobian(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
3702961e153SJed Brown @*/
3712961e153SJed Brown PetscErrorCode DMDASNESSetJacobianLocal(DM dm,PetscErrorCode (*func)(DMDALocalInfo*,void*,Mat,Mat,MatStructure*,void*),void *ctx)
3722961e153SJed Brown {
3732961e153SJed Brown   PetscErrorCode ierr;
374*942e3340SBarry Smith   DMSNES         sdm;
375*942e3340SBarry Smith   DMSNES_DA      *dmdasnes;
3762961e153SJed Brown 
3772961e153SJed Brown   PetscFunctionBegin;
3782961e153SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
379*942e3340SBarry Smith   ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr);
3802961e153SJed Brown   ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr);
3812961e153SJed Brown   dmdasnes->jacobianlocal = func;
3822961e153SJed Brown   dmdasnes->jacobianlocalctx = ctx;
3832961e153SJed Brown   ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes);CHKERRQ(ierr);
3846427ad5bSPeter Brune   ierr = DMSNESSetBlockJacobian(dm,SNESComputeLocalBlockJacobian_DMDA,dm);CHKERRQ(ierr);
3856cab3a1bSJed Brown   PetscFunctionReturn(0);
3866cab3a1bSJed Brown }
3872a4ee8f2SPeter Brune 
3882a4ee8f2SPeter Brune 
3892a4ee8f2SPeter Brune #undef __FUNCT__
3902a4ee8f2SPeter Brune #define __FUNCT__ "DMDASNESSetObjectiveLocal"
3912a4ee8f2SPeter Brune /*@C
3922a4ee8f2SPeter Brune    DMDASNESSetObjectiveLocal - set a local residual evaluation function
3932a4ee8f2SPeter Brune 
3942a4ee8f2SPeter Brune    Logically Collective
3952a4ee8f2SPeter Brune 
3962a4ee8f2SPeter Brune    Input Arguments:
3972a4ee8f2SPeter Brune +  dm - DM to associate callback with
3982a4ee8f2SPeter Brune .  func - local objective evaluation
3992a4ee8f2SPeter Brune -  ctx - optional context for local residual evaluation
4002a4ee8f2SPeter Brune 
4012a4ee8f2SPeter Brune    Calling sequence for func:
4022a4ee8f2SPeter Brune +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
4032a4ee8f2SPeter Brune .  x - dimensional pointer to state at which to evaluate residual
4042a4ee8f2SPeter Brune .  ob - eventual objective value
4052a4ee8f2SPeter Brune -  ctx - optional context passed above
4062a4ee8f2SPeter Brune 
4072a4ee8f2SPeter Brune    Level: beginner
4082a4ee8f2SPeter Brune 
4092a4ee8f2SPeter Brune .seealso: DMSNESSetFunction(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
4102a4ee8f2SPeter Brune @*/
4112a4ee8f2SPeter Brune PetscErrorCode DMDASNESSetObjectiveLocal(DM dm,DMDASNESObjective func,void *ctx)
4122a4ee8f2SPeter Brune {
4132a4ee8f2SPeter Brune   PetscErrorCode ierr;
414*942e3340SBarry Smith   DMSNES         sdm;
415*942e3340SBarry Smith   DMSNES_DA      *dmdasnes;
4162a4ee8f2SPeter Brune 
4172a4ee8f2SPeter Brune   PetscFunctionBegin;
4182a4ee8f2SPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
419*942e3340SBarry Smith   ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr);
4202a4ee8f2SPeter Brune   ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr);
4212a4ee8f2SPeter Brune   dmdasnes->objectivelocal = func;
4222a4ee8f2SPeter Brune   dmdasnes->objectivelocalctx = ctx;
4232a4ee8f2SPeter Brune   ierr = DMSNESSetObjective(dm,SNESComputeObjective_DMDA,dmdasnes);CHKERRQ(ierr);
4242a4ee8f2SPeter Brune   PetscFunctionReturn(0);
4252a4ee8f2SPeter Brune }
426dcad5f1cSBarry Smith 
427dcad5f1cSBarry Smith #undef __FUNCT__
428dcad5f1cSBarry Smith #define __FUNCT__ "SNESComputePicard_DMDA"
429dcad5f1cSBarry Smith static PetscErrorCode SNESComputePicard_DMDA(SNES snes,Vec X,Vec F,void *ctx)
430dcad5f1cSBarry Smith {
431dcad5f1cSBarry Smith   PetscErrorCode ierr;
432dcad5f1cSBarry Smith   DM             dm;
433*942e3340SBarry Smith   DMSNES_DA      *dmdasnes = (DMSNES_DA *)ctx;
434dcad5f1cSBarry Smith   DMDALocalInfo  info;
435dcad5f1cSBarry Smith   Vec            Xloc;
436dcad5f1cSBarry Smith   void           *x,*f;
437dcad5f1cSBarry Smith 
438dcad5f1cSBarry Smith   PetscFunctionBegin;
439dcad5f1cSBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
440dcad5f1cSBarry Smith   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
441dcad5f1cSBarry Smith   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
442dcad5f1cSBarry Smith   if (!dmdasnes->rhsplocal) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Corrupt context");
443dcad5f1cSBarry Smith   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
444dcad5f1cSBarry Smith   ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
445dcad5f1cSBarry Smith   ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
446dcad5f1cSBarry Smith   ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
447dcad5f1cSBarry Smith   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
448dcad5f1cSBarry Smith   ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
449dcad5f1cSBarry Smith   switch (dmdasnes->residuallocalimode) {
450dcad5f1cSBarry Smith   case INSERT_VALUES: {
451dcad5f1cSBarry Smith     ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr);
452dcad5f1cSBarry Smith     CHKMEMQ;
453dcad5f1cSBarry Smith     ierr = (*dmdasnes->rhsplocal)(&info,x,f,dmdasnes->picardlocalctx);CHKERRQ(ierr);
454dcad5f1cSBarry Smith     CHKMEMQ;
455dcad5f1cSBarry Smith     ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr);
456dcad5f1cSBarry Smith   } break;
457dcad5f1cSBarry Smith   case ADD_VALUES: {
458dcad5f1cSBarry Smith     Vec Floc;
459dcad5f1cSBarry Smith     ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr);
460dcad5f1cSBarry Smith     ierr = VecZeroEntries(Floc);CHKERRQ(ierr);
461dcad5f1cSBarry Smith     ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr);
462dcad5f1cSBarry Smith     CHKMEMQ;
463dcad5f1cSBarry Smith     ierr = (*dmdasnes->rhsplocal)(&info,x,f,dmdasnes->picardlocalctx);CHKERRQ(ierr);
464dcad5f1cSBarry Smith     CHKMEMQ;
465dcad5f1cSBarry Smith     ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr);
466dcad5f1cSBarry Smith     ierr = VecZeroEntries(F);CHKERRQ(ierr);
467dcad5f1cSBarry Smith     ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
468dcad5f1cSBarry Smith     ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
469dcad5f1cSBarry Smith     ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr);
470dcad5f1cSBarry Smith   } break;
471dcad5f1cSBarry Smith   default: SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdasnes->residuallocalimode);
472dcad5f1cSBarry Smith   }
473dcad5f1cSBarry Smith   ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
474dcad5f1cSBarry Smith   ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
475dcad5f1cSBarry Smith   PetscFunctionReturn(0);
476dcad5f1cSBarry Smith }
477dcad5f1cSBarry Smith 
478dcad5f1cSBarry Smith #undef __FUNCT__
479dcad5f1cSBarry Smith #define __FUNCT__ "SNESComputePicardJacobian_DMDA"
480dcad5f1cSBarry Smith static PetscErrorCode SNESComputePicardJacobian_DMDA(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *ctx)
481dcad5f1cSBarry Smith {
482dcad5f1cSBarry Smith   PetscErrorCode ierr;
483dcad5f1cSBarry Smith   DM             dm;
484*942e3340SBarry Smith   DMSNES_DA      *dmdasnes = (DMSNES_DA *)ctx;
485dcad5f1cSBarry Smith   DMDALocalInfo  info;
486dcad5f1cSBarry Smith   Vec            Xloc;
487dcad5f1cSBarry Smith   void           *x;
488dcad5f1cSBarry Smith 
489dcad5f1cSBarry Smith   PetscFunctionBegin;
490dcad5f1cSBarry Smith   if (!dmdasnes->jacobianplocal) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Corrupt context");
491dcad5f1cSBarry Smith   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
492dcad5f1cSBarry Smith 
493dcad5f1cSBarry Smith   ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
494dcad5f1cSBarry Smith   ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
495dcad5f1cSBarry Smith   ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
496dcad5f1cSBarry Smith   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
497dcad5f1cSBarry Smith   ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
498dcad5f1cSBarry Smith   CHKMEMQ;
499dcad5f1cSBarry Smith   ierr = (*dmdasnes->jacobianplocal)(&info,x,*A,*B,mstr,dmdasnes->picardlocalctx);CHKERRQ(ierr);
500dcad5f1cSBarry Smith   CHKMEMQ;
501dcad5f1cSBarry Smith   ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
502dcad5f1cSBarry Smith   ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
503dcad5f1cSBarry Smith   *mstr = SAME_NONZERO_PATTERN;
504dcad5f1cSBarry Smith   if (*A != *B) {
505dcad5f1cSBarry Smith     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
506dcad5f1cSBarry Smith     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
507dcad5f1cSBarry Smith   }
508dcad5f1cSBarry Smith   PetscFunctionReturn(0);
509dcad5f1cSBarry Smith }
510dcad5f1cSBarry Smith 
511dcad5f1cSBarry Smith #undef __FUNCT__
512dcad5f1cSBarry Smith #define __FUNCT__ "DMDASNESSetPicardLocal"
513dcad5f1cSBarry Smith /*@C
514dcad5f1cSBarry Smith    DMDASNESSetPicardLocal - set a local right hand side and matrix evaluation function for Picard iteration
515dcad5f1cSBarry Smith 
516dcad5f1cSBarry Smith    Logically Collective
517dcad5f1cSBarry Smith 
518dcad5f1cSBarry Smith    Input Arguments:
519dcad5f1cSBarry Smith +  dm - DM to associate callback with
520dcad5f1cSBarry Smith .  imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part
521dcad5f1cSBarry Smith .  func - local residual evaluation
522dcad5f1cSBarry Smith -  ctx - optional context for local residual evaluation
523dcad5f1cSBarry Smith 
524dcad5f1cSBarry Smith    Calling sequence for func:
525dcad5f1cSBarry Smith +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
526dcad5f1cSBarry Smith .  x - dimensional pointer to state at which to evaluate residual
527dcad5f1cSBarry Smith .  f - dimensional pointer to residual, write the residual here
528dcad5f1cSBarry Smith -  ctx - optional context passed above
529dcad5f1cSBarry Smith 
530dcad5f1cSBarry Smith    Notes:  The user must use
531dcad5f1cSBarry Smith     extern PetscErrorCode  SNESPicardComputeFunction(SNES,Vec,Vec,void*);
532dcad5f1cSBarry Smith     extern PetscErrorCode  SNESPicardComputeJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
533dcad5f1cSBarry Smith     ierr = SNESSetFunction(snes,PETSC_NULL,SNESPicardComputeFunction,&user);CHKERRQ(ierr);
534dcad5f1cSBarry Smith     in their code before calling this routine.
535dcad5f1cSBarry Smith 
536dcad5f1cSBarry Smith 
537dcad5f1cSBarry Smith    Level: beginner
538dcad5f1cSBarry Smith 
539dcad5f1cSBarry Smith .seealso: DMSNESSetFunction(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
540dcad5f1cSBarry Smith @*/
541dcad5f1cSBarry Smith PetscErrorCode DMDASNESSetPicardLocal(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,void*,void*,void*),
542dcad5f1cSBarry Smith                                         PetscErrorCode (*jac)(DMDALocalInfo*,void*,Mat,Mat,MatStructure*,void*),void *ctx)
543dcad5f1cSBarry Smith {
544dcad5f1cSBarry Smith   PetscErrorCode ierr;
545*942e3340SBarry Smith   DMSNES         sdm;
546*942e3340SBarry Smith   DMSNES_DA      *dmdasnes;
547dcad5f1cSBarry Smith 
548dcad5f1cSBarry Smith   PetscFunctionBegin;
549dcad5f1cSBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
550*942e3340SBarry Smith   ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr);
551dcad5f1cSBarry Smith   ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr);
552dcad5f1cSBarry Smith   dmdasnes->residuallocalimode = imode;
553dcad5f1cSBarry Smith   dmdasnes->rhsplocal          = func;
554dcad5f1cSBarry Smith   dmdasnes->jacobianplocal     = jac;
555dcad5f1cSBarry Smith   dmdasnes->picardlocalctx     = ctx;
556dcad5f1cSBarry Smith   ierr = DMSNESSetPicard(dm,SNESComputePicard_DMDA,SNESComputePicardJacobian_DMDA,dmdasnes);CHKERRQ(ierr);
557dcad5f1cSBarry Smith   PetscFunctionReturn(0);
558dcad5f1cSBarry Smith }
559