xref: /petsc/src/snes/utils/dmdasnes.c (revision e3c51ba28ac7705cb252437d122da1ac93d8fb1e)
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;
146cab3a1bSJed Brown } DM_DA_SNES;
156cab3a1bSJed Brown 
166cab3a1bSJed Brown #undef __FUNCT__
176cab3a1bSJed Brown #define __FUNCT__ "SNESDMDestroy_DMDA"
186cab3a1bSJed Brown static PetscErrorCode SNESDMDestroy_DMDA(SNESDM sdm)
196cab3a1bSJed Brown {
206cab3a1bSJed Brown   PetscErrorCode ierr;
216cab3a1bSJed Brown 
226cab3a1bSJed Brown   PetscFunctionBegin;
236cab3a1bSJed Brown   ierr = PetscFree(sdm->data);CHKERRQ(ierr);
246cab3a1bSJed Brown   PetscFunctionReturn(0);
256cab3a1bSJed Brown }
266cab3a1bSJed Brown 
276cab3a1bSJed Brown #undef __FUNCT__
28e3c0b8a2SPeter Brune #define __FUNCT__ "SNESDMDuplicate_DMDA"
29e3c0b8a2SPeter Brune static PetscErrorCode SNESDMDuplicate_DMDA(SNESDM oldsdm,DM dm)
30e3c0b8a2SPeter Brune {
31e3c0b8a2SPeter Brune   PetscErrorCode ierr;
32e3c0b8a2SPeter Brune   SNESDM         sdm;
33e3c0b8a2SPeter Brune 
34e3c0b8a2SPeter Brune   PetscFunctionBegin;
35e3c0b8a2SPeter Brune   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
36e3c0b8a2SPeter Brune   ierr = PetscNewLog(dm,DM_DA_SNES,&sdm->data);CHKERRQ(ierr);
37e3c0b8a2SPeter Brune   if (oldsdm->data)ierr = PetscMemcpy(sdm->data,oldsdm->data,sizeof(DM_DA_SNES));CHKERRQ(ierr);
38e3c0b8a2SPeter Brune   PetscFunctionReturn(0);
39e3c0b8a2SPeter Brune }
40e3c0b8a2SPeter Brune 
41e3c0b8a2SPeter Brune 
42e3c0b8a2SPeter Brune #undef __FUNCT__
436cab3a1bSJed Brown #define __FUNCT__ "DMDASNESGetContext"
446cab3a1bSJed Brown static PetscErrorCode DMDASNESGetContext(DM dm,SNESDM sdm,DM_DA_SNES **dmdasnes)
456cab3a1bSJed Brown {
466cab3a1bSJed Brown   PetscErrorCode ierr;
476cab3a1bSJed Brown 
486cab3a1bSJed Brown   PetscFunctionBegin;
49b8d0cc33SJed Brown   *dmdasnes = PETSC_NULL;
506cab3a1bSJed Brown   if (!sdm->data) {
516cab3a1bSJed Brown     ierr = PetscNewLog(dm,DM_DA_SNES,&sdm->data);CHKERRQ(ierr);
526cab3a1bSJed Brown     sdm->destroy = SNESDMDestroy_DMDA;
53e3c0b8a2SPeter Brune     sdm->duplicate = SNESDMDuplicate_DMDA;
546cab3a1bSJed Brown   }
556cab3a1bSJed Brown   *dmdasnes = (DM_DA_SNES*)sdm->data;
566cab3a1bSJed Brown   PetscFunctionReturn(0);
576cab3a1bSJed Brown }
586cab3a1bSJed Brown 
596cab3a1bSJed Brown #undef __FUNCT__
606cab3a1bSJed Brown #define __FUNCT__ "SNESComputeFunction_DMDA"
616cab3a1bSJed Brown /*
626cab3a1bSJed Brown   This function should eventually replace:
636cab3a1bSJed Brown     DMDAComputeFunction() and DMDAComputeFunction1()
646cab3a1bSJed Brown  */
656cab3a1bSJed Brown static PetscErrorCode SNESComputeFunction_DMDA(SNES snes,Vec X,Vec F,void *ctx)
666cab3a1bSJed Brown {
676cab3a1bSJed Brown   PetscErrorCode ierr;
686cab3a1bSJed Brown   DM             dm;
696cab3a1bSJed Brown   DM_DA_SNES     *dmdasnes = (DM_DA_SNES*)ctx;
706cab3a1bSJed Brown   DMDALocalInfo  info;
716cab3a1bSJed Brown   Vec            Xloc;
726cab3a1bSJed Brown   void           *x,*f;
736cab3a1bSJed Brown 
746cab3a1bSJed Brown   PetscFunctionBegin;
752961e153SJed Brown   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
762961e153SJed Brown   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
772961e153SJed Brown   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
786cab3a1bSJed Brown   if (!dmdasnes->residuallocal) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Corrupt context");
796cab3a1bSJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
806cab3a1bSJed Brown   ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
816cab3a1bSJed Brown   ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
826cab3a1bSJed Brown   ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
836cab3a1bSJed Brown   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
846cab3a1bSJed Brown   ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
85709ac0d2SJed Brown   switch (dmdasnes->residuallocalimode) {
86709ac0d2SJed Brown   case INSERT_VALUES: {
876cab3a1bSJed Brown     ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr);
886cab3a1bSJed Brown     CHKMEMQ;
896cab3a1bSJed Brown     ierr = (*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx);CHKERRQ(ierr);
906cab3a1bSJed Brown     CHKMEMQ;
916cab3a1bSJed Brown     ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr);
92709ac0d2SJed Brown   } break;
93709ac0d2SJed Brown   case ADD_VALUES: {
94709ac0d2SJed Brown     Vec Floc;
95709ac0d2SJed Brown     ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr);
96709ac0d2SJed Brown     ierr = VecZeroEntries(Floc);CHKERRQ(ierr);
97709ac0d2SJed Brown     ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr);
98709ac0d2SJed Brown     CHKMEMQ;
99709ac0d2SJed Brown     ierr = (*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx);CHKERRQ(ierr);
100709ac0d2SJed Brown     CHKMEMQ;
101709ac0d2SJed Brown     ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr);
102709ac0d2SJed Brown     ierr = VecZeroEntries(F);CHKERRQ(ierr);
103709ac0d2SJed Brown     ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
104709ac0d2SJed Brown     ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
105709ac0d2SJed Brown     ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr);
106709ac0d2SJed Brown   } break;
107709ac0d2SJed Brown   default: SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdasnes->residuallocalimode);
108709ac0d2SJed Brown   }
109709ac0d2SJed Brown   ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
1106cab3a1bSJed Brown   ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
1116cab3a1bSJed Brown   PetscFunctionReturn(0);
1126cab3a1bSJed Brown }
1136cab3a1bSJed Brown 
1146cab3a1bSJed Brown #undef __FUNCT__
1152a4ee8f2SPeter Brune #define __FUNCT__ "SNESComputeObjective_DMDA"
1162a4ee8f2SPeter Brune /*
1172a4ee8f2SPeter Brune   This function should eventually replace:
1182a4ee8f2SPeter Brune     DMDAComputeFunction() and DMDAComputeFunction1()
1192a4ee8f2SPeter Brune  */
1202a4ee8f2SPeter Brune static PetscErrorCode SNESComputeObjective_DMDA(SNES snes,Vec X,PetscReal *ob,void *ctx)
1212a4ee8f2SPeter Brune {
1222a4ee8f2SPeter Brune   PetscErrorCode ierr;
1232a4ee8f2SPeter Brune   DM             dm;
1242a4ee8f2SPeter Brune   DM_DA_SNES     *dmdasnes = (DM_DA_SNES*)ctx;
1252a4ee8f2SPeter Brune   DMDALocalInfo  info;
1262a4ee8f2SPeter Brune   Vec            Xloc;
1272a4ee8f2SPeter Brune   void           *x;
1282a4ee8f2SPeter Brune 
1292a4ee8f2SPeter Brune   PetscFunctionBegin;
1302a4ee8f2SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1312a4ee8f2SPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
1322a4ee8f2SPeter Brune   PetscValidPointer(ob,3);
1332a4ee8f2SPeter Brune   if (!dmdasnes->objectivelocal) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Corrupt context");
1342a4ee8f2SPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
1352a4ee8f2SPeter Brune   ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
1362a4ee8f2SPeter Brune   ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
1372a4ee8f2SPeter Brune   ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
1382a4ee8f2SPeter Brune   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
1392a4ee8f2SPeter Brune   ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
1402a4ee8f2SPeter Brune   CHKMEMQ;
1412a4ee8f2SPeter Brune   ierr = (*dmdasnes->objectivelocal)(&info,x,ob,dmdasnes->objectivelocalctx);CHKERRQ(ierr);
1422a4ee8f2SPeter Brune   CHKMEMQ;
1432a4ee8f2SPeter Brune   ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
1442a4ee8f2SPeter Brune   ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
1452a4ee8f2SPeter Brune   PetscFunctionReturn(0);
1462a4ee8f2SPeter Brune }
1472a4ee8f2SPeter Brune 
1482a4ee8f2SPeter Brune #undef __FUNCT__
1496427ad5bSPeter Brune #define __FUNCT__ "SNESComputeLocalBlockFunction_DMDA"
1506427ad5bSPeter Brune /*
1516427ad5bSPeter Brune   This function should eventually replace:
1526427ad5bSPeter Brune     DMDAComputeFunction() and DMDAComputeFunction1()
1536427ad5bSPeter Brune  */
1546427ad5bSPeter Brune PetscErrorCode SNESComputeLocalBlockFunction_DMDA(SNES snes,Vec X,Vec F,void *dmctx)
1556427ad5bSPeter Brune {
1566427ad5bSPeter Brune   PetscErrorCode ierr;
1576427ad5bSPeter Brune   DM             dm = (DM)dmctx; /* the global DMDA context */
1586427ad5bSPeter Brune   DMDALocalInfo  info;
1596427ad5bSPeter Brune   void           *x,*f;
1606427ad5bSPeter Brune   DM_DA_SNES     *dmdasnes;
1616427ad5bSPeter Brune   SNESDM         sdm;
1626427ad5bSPeter Brune 
1636427ad5bSPeter Brune   PetscFunctionBegin;
1646427ad5bSPeter Brune 
1656427ad5bSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1666427ad5bSPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
1676427ad5bSPeter Brune   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
1686427ad5bSPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,4);
1696427ad5bSPeter Brune 
1706427ad5bSPeter Brune   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
1716427ad5bSPeter Brune   ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr);
172ec31c420SJed Brown   if (!dmdasnes->residuallocal) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Corrupt context");
1736427ad5bSPeter Brune   ierr = DMDAVecGetArray(dm,X,&x);CHKERRQ(ierr);
1746427ad5bSPeter Brune   ierr = DMDAGetLocalBlockInfo(dm,&info);CHKERRQ(ierr);
1756427ad5bSPeter Brune   switch (dmdasnes->residuallocalimode) {
1766427ad5bSPeter Brune   case INSERT_VALUES: {
1776427ad5bSPeter Brune     ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr);
1786427ad5bSPeter Brune     CHKMEMQ;
1796427ad5bSPeter Brune     ierr = (*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx);CHKERRQ(ierr);
1806427ad5bSPeter Brune     CHKMEMQ;
1816427ad5bSPeter Brune     ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr);
1826427ad5bSPeter Brune   } break;
1836427ad5bSPeter Brune   case ADD_VALUES: {
1846427ad5bSPeter Brune     Vec Floc;
1856427ad5bSPeter Brune     ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr);
1866427ad5bSPeter Brune     ierr = VecZeroEntries(Floc);CHKERRQ(ierr);
1876427ad5bSPeter Brune     ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr);
1886427ad5bSPeter Brune     CHKMEMQ;
1896427ad5bSPeter Brune     ierr = (*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx);CHKERRQ(ierr);
1906427ad5bSPeter Brune     CHKMEMQ;
1916427ad5bSPeter Brune     ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr);
1926427ad5bSPeter Brune     ierr = VecZeroEntries(F);CHKERRQ(ierr);
1936427ad5bSPeter Brune     ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
1946427ad5bSPeter Brune     ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr);
1956427ad5bSPeter Brune     ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr);
1966427ad5bSPeter Brune   } break;
1976427ad5bSPeter Brune   default: SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdasnes->residuallocalimode);
1986427ad5bSPeter Brune   }
1996427ad5bSPeter Brune   ierr = DMDAVecRestoreArray(dm,X,&x);CHKERRQ(ierr);
2006427ad5bSPeter Brune   PetscFunctionReturn(0);
2016427ad5bSPeter Brune }
2026427ad5bSPeter Brune 
2036427ad5bSPeter Brune #undef __FUNCT__
2042961e153SJed Brown #define __FUNCT__ "SNESComputeJacobian_DMDA"
2052961e153SJed Brown /*
2062961e153SJed Brown   This function should eventually replace:
2072961e153SJed Brown     DMComputeJacobian() and DMDAComputeJacobian1()
2082961e153SJed Brown  */
2092961e153SJed Brown static PetscErrorCode SNESComputeJacobian_DMDA(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *ctx)
2102961e153SJed Brown {
2112961e153SJed Brown   PetscErrorCode ierr;
2122961e153SJed Brown   DM             dm;
2132961e153SJed Brown   DM_DA_SNES     *dmdasnes = (DM_DA_SNES*)ctx;
2142961e153SJed Brown   DMDALocalInfo  info;
2152961e153SJed Brown   Vec            Xloc;
2162961e153SJed Brown   void           *x;
2172961e153SJed Brown 
2182961e153SJed Brown   PetscFunctionBegin;
2192961e153SJed Brown   if (!dmdasnes->residuallocal) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Corrupt context");
2202961e153SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
2212961e153SJed Brown 
2222961e153SJed Brown   if (dmdasnes->jacobianlocal) {
2232961e153SJed Brown     ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
2242961e153SJed Brown     ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
2252961e153SJed Brown     ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
2262961e153SJed Brown     ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
2272961e153SJed Brown     ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
2282961e153SJed Brown     CHKMEMQ;
2292961e153SJed Brown     ierr = (*dmdasnes->jacobianlocal)(&info,x,*A,*B,mstr,dmdasnes->jacobianlocalctx);CHKERRQ(ierr);
2302961e153SJed Brown     CHKMEMQ;
2312961e153SJed Brown     ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
2322961e153SJed Brown     ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
2332961e153SJed Brown   } else {
2342961e153SJed Brown     MatFDColoring fdcoloring;
2352961e153SJed Brown     ierr = PetscObjectQuery((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject*)&fdcoloring);CHKERRQ(ierr);
2362961e153SJed Brown     if (!fdcoloring) {
2372961e153SJed Brown       ISColoring     coloring;
2382961e153SJed Brown 
2392961e153SJed Brown       ierr = DMCreateColoring(dm,dm->coloringtype,MATAIJ,&coloring);CHKERRQ(ierr);
2402961e153SJed Brown       ierr = MatFDColoringCreate(*B,coloring,&fdcoloring);CHKERRQ(ierr);
2412961e153SJed Brown       ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr);
2422961e153SJed Brown       switch (dm->coloringtype) {
2432961e153SJed Brown       case IS_COLORING_GLOBAL:
2442961e153SJed Brown         ierr = MatFDColoringSetFunction(fdcoloring,(PetscErrorCode(*)(void))SNESComputeFunction_DMDA,dmdasnes);CHKERRQ(ierr);
2452961e153SJed Brown         break;
2462961e153SJed Brown       default: SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"No support for coloring type '%s'",ISColoringTypes[dm->coloringtype]);
2472961e153SJed Brown       }
2482961e153SJed Brown       ierr = PetscObjectSetOptionsPrefix((PetscObject)fdcoloring,((PetscObject)dm)->prefix);CHKERRQ(ierr);
2492961e153SJed Brown       ierr = MatFDColoringSetFromOptions(fdcoloring);CHKERRQ(ierr);
2502961e153SJed Brown       ierr = PetscObjectCompose((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject)fdcoloring);CHKERRQ(ierr);
2512961e153SJed Brown       ierr = PetscObjectDereference((PetscObject)fdcoloring);CHKERRQ(ierr);
2522961e153SJed Brown 
2532961e153SJed Brown       /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call
2542961e153SJed Brown        * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the
2552961e153SJed Brown        * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually
2562961e153SJed Brown        * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be
2572961e153SJed Brown        * taken when the PetscOList for the Vec inside MatFDColoring is destroyed.
2582961e153SJed Brown        */
2592961e153SJed Brown       ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
2602961e153SJed Brown     }
2612961e153SJed Brown     *mstr = SAME_NONZERO_PATTERN;
2622961e153SJed Brown     ierr = MatFDColoringApply(*B,fdcoloring,X,mstr,snes);CHKERRQ(ierr);
2632961e153SJed Brown   }
2642961e153SJed Brown   /* This will be redundant if the user called both, but it's too common to forget. */
2652961e153SJed Brown   if (*A != *B) {
2662961e153SJed Brown     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2672961e153SJed Brown     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2682961e153SJed Brown   }
2692961e153SJed Brown   PetscFunctionReturn(0);
2702961e153SJed Brown }
2712961e153SJed Brown 
2722961e153SJed Brown #undef __FUNCT__
2736427ad5bSPeter Brune #define __FUNCT__ "SNESComputeJacobian_DMDA"
2746427ad5bSPeter Brune /*
2756427ad5bSPeter Brune   This function should eventually replace:
2766427ad5bSPeter Brune     DMComputeJacobian() and DMDAComputeJacobian1()
2776427ad5bSPeter Brune  */
2786427ad5bSPeter Brune PetscErrorCode SNESComputeLocalBlockJacobian_DMDA(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *dmctx)
2796427ad5bSPeter Brune {
2806427ad5bSPeter Brune   PetscErrorCode ierr;
2816427ad5bSPeter Brune   DM             dm = (DM)dmctx; /* the global DMDA context */
2826427ad5bSPeter Brune   DM_DA_SNES     *dmdasnes;
2836427ad5bSPeter Brune   SNESDM         sdm;
2846427ad5bSPeter Brune   DMDALocalInfo  info;
2856427ad5bSPeter Brune   void           *x;
2866427ad5bSPeter Brune 
2876427ad5bSPeter Brune   PetscFunctionBegin;
2886427ad5bSPeter Brune   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
2896427ad5bSPeter Brune   ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr);
290ec31c420SJed Brown   if (!dmdasnes->residuallocal) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Corrupt context");
2916427ad5bSPeter Brune   if (dmdasnes->jacobianlocal) {
2926427ad5bSPeter Brune     ierr = DMDAVecGetArray(dm,X,&x);CHKERRQ(ierr);
2936427ad5bSPeter Brune     CHKMEMQ;
2946427ad5bSPeter Brune     ierr = (*dmdasnes->jacobianlocal)(&info,x,*A,*B,mstr,dmdasnes->jacobianlocalctx);CHKERRQ(ierr);
2956427ad5bSPeter Brune     CHKMEMQ;
2966427ad5bSPeter Brune     ierr = DMDAVecRestoreArray(dm,X,&x);CHKERRQ(ierr);
2976427ad5bSPeter Brune   }
2986427ad5bSPeter Brune   /* This will be redundant if the user called both, but it's too common to forget. */
2996427ad5bSPeter Brune   if (*A != *B) {
3006427ad5bSPeter Brune     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3016427ad5bSPeter Brune     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3026427ad5bSPeter Brune   }
3036427ad5bSPeter Brune   PetscFunctionReturn(0);
3046427ad5bSPeter Brune }
3056427ad5bSPeter Brune 
3066427ad5bSPeter Brune 
3076427ad5bSPeter Brune #undef __FUNCT__
3086cab3a1bSJed Brown #define __FUNCT__ "DMDASNESSetFunctionLocal"
3096cab3a1bSJed Brown /*@C
3106cab3a1bSJed Brown    DMDASNESSetFunctionLocal - set a local residual evaluation function
3116cab3a1bSJed Brown 
3126cab3a1bSJed Brown    Logically Collective
3136cab3a1bSJed Brown 
3146cab3a1bSJed Brown    Input Arguments:
3156cab3a1bSJed Brown +  dm - DM to associate callback with
316*e3c51ba2SJed Brown .  imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part
3176cab3a1bSJed Brown .  func - local residual evaluation
3186cab3a1bSJed Brown -  ctx - optional context for local residual evaluation
3196cab3a1bSJed Brown 
3206cab3a1bSJed Brown    Calling sequence for func:
3216cab3a1bSJed Brown +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
3226cab3a1bSJed Brown .  x - dimensional pointer to state at which to evaluate residual
3236cab3a1bSJed Brown .  f - dimensional pointer to residual, write the residual here
3246cab3a1bSJed Brown -  ctx - optional context passed above
3256cab3a1bSJed Brown 
3266cab3a1bSJed Brown    Level: beginner
3276cab3a1bSJed Brown 
3286cab3a1bSJed Brown .seealso: DMSNESSetFunction(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
3296cab3a1bSJed Brown @*/
330709ac0d2SJed Brown PetscErrorCode DMDASNESSetFunctionLocal(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,void*,void*,void*),void *ctx)
3316cab3a1bSJed Brown {
3326cab3a1bSJed Brown   PetscErrorCode ierr;
3336cab3a1bSJed Brown   SNESDM         sdm;
3346cab3a1bSJed Brown   DM_DA_SNES     *dmdasnes;
3356cab3a1bSJed Brown 
3366cab3a1bSJed Brown   PetscFunctionBegin;
3376cab3a1bSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3386cab3a1bSJed Brown   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
3396cab3a1bSJed Brown   ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr);
340709ac0d2SJed Brown   dmdasnes->residuallocalimode = imode;
3416cab3a1bSJed Brown   dmdasnes->residuallocal = func;
3426cab3a1bSJed Brown   dmdasnes->residuallocalctx = ctx;
3436cab3a1bSJed Brown   ierr = DMSNESSetFunction(dm,SNESComputeFunction_DMDA,dmdasnes);CHKERRQ(ierr);
3446427ad5bSPeter Brune   ierr = DMSNESSetBlockFunction(dm,SNESComputeLocalBlockFunction_DMDA,dm);CHKERRQ(ierr);
3452961e153SJed Brown   if (!sdm->computejacobian) {  /* Call us for the Jacobian too, can be overridden by the user. */
3462961e153SJed Brown     ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes);CHKERRQ(ierr);
3472961e153SJed Brown   }
3482961e153SJed Brown   PetscFunctionReturn(0);
3492961e153SJed Brown }
3502961e153SJed Brown 
3512961e153SJed Brown #undef __FUNCT__
3522961e153SJed Brown #define __FUNCT__ "DMDASNESSetJacobianLocal"
3532961e153SJed Brown /*@C
3542961e153SJed Brown    DMDASNESSetJacobianLocal - set a local residual evaluation function
3552961e153SJed Brown 
3562961e153SJed Brown    Logically Collective
3572961e153SJed Brown 
3582961e153SJed Brown    Input Arguments:
3592961e153SJed Brown +  dm - DM to associate callback with
3602961e153SJed Brown .  func - local residual evaluation
3612961e153SJed Brown -  ctx - optional context for local residual evaluation
3622961e153SJed Brown 
3632961e153SJed Brown    Calling sequence for func:
3642961e153SJed Brown +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
3652961e153SJed Brown .  x - dimensional pointer to state at which to evaluate residual
3662961e153SJed Brown .  f - dimensional pointer to residual, write the residual here
3672961e153SJed Brown -  ctx - optional context passed above
3682961e153SJed Brown 
3692961e153SJed Brown    Level: beginner
3702961e153SJed Brown 
3712961e153SJed Brown .seealso: DMSNESSetJacobian(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
3722961e153SJed Brown @*/
3732961e153SJed Brown PetscErrorCode DMDASNESSetJacobianLocal(DM dm,PetscErrorCode (*func)(DMDALocalInfo*,void*,Mat,Mat,MatStructure*,void*),void *ctx)
3742961e153SJed Brown {
3752961e153SJed Brown   PetscErrorCode ierr;
3762961e153SJed Brown   SNESDM         sdm;
3772961e153SJed Brown   DM_DA_SNES     *dmdasnes;
3782961e153SJed Brown 
3792961e153SJed Brown   PetscFunctionBegin;
3802961e153SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3812961e153SJed Brown   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
3822961e153SJed Brown   ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr);
3832961e153SJed Brown   dmdasnes->jacobianlocal = func;
3842961e153SJed Brown   dmdasnes->jacobianlocalctx = ctx;
3852961e153SJed Brown   ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes);CHKERRQ(ierr);
3866427ad5bSPeter Brune   ierr = DMSNESSetBlockJacobian(dm,SNESComputeLocalBlockJacobian_DMDA,dm);CHKERRQ(ierr);
3876cab3a1bSJed Brown   PetscFunctionReturn(0);
3886cab3a1bSJed Brown }
3892a4ee8f2SPeter Brune 
3902a4ee8f2SPeter Brune 
3912a4ee8f2SPeter Brune #undef __FUNCT__
3922a4ee8f2SPeter Brune #define __FUNCT__ "DMDASNESSetObjectiveLocal"
3932a4ee8f2SPeter Brune /*@C
3942a4ee8f2SPeter Brune    DMDASNESSetObjectiveLocal - set a local residual evaluation function
3952a4ee8f2SPeter Brune 
3962a4ee8f2SPeter Brune    Logically Collective
3972a4ee8f2SPeter Brune 
3982a4ee8f2SPeter Brune    Input Arguments:
3992a4ee8f2SPeter Brune +  dm - DM to associate callback with
4002a4ee8f2SPeter Brune .  func - local objective evaluation
4012a4ee8f2SPeter Brune -  ctx - optional context for local residual evaluation
4022a4ee8f2SPeter Brune 
4032a4ee8f2SPeter Brune    Calling sequence for func:
4042a4ee8f2SPeter Brune +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
4052a4ee8f2SPeter Brune .  x - dimensional pointer to state at which to evaluate residual
4062a4ee8f2SPeter Brune .  ob - eventual objective value
4072a4ee8f2SPeter Brune -  ctx - optional context passed above
4082a4ee8f2SPeter Brune 
4092a4ee8f2SPeter Brune    Level: beginner
4102a4ee8f2SPeter Brune 
4112a4ee8f2SPeter Brune .seealso: DMSNESSetFunction(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
4122a4ee8f2SPeter Brune @*/
4132a4ee8f2SPeter Brune PetscErrorCode DMDASNESSetObjectiveLocal(DM dm,DMDASNESObjective func,void *ctx)
4142a4ee8f2SPeter Brune {
4152a4ee8f2SPeter Brune   PetscErrorCode ierr;
4162a4ee8f2SPeter Brune   SNESDM         sdm;
4172a4ee8f2SPeter Brune   DM_DA_SNES     *dmdasnes;
4182a4ee8f2SPeter Brune 
4192a4ee8f2SPeter Brune   PetscFunctionBegin;
4202a4ee8f2SPeter Brune   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4212a4ee8f2SPeter Brune   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
4222a4ee8f2SPeter Brune   ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr);
4232a4ee8f2SPeter Brune   dmdasnes->objectivelocal = func;
4242a4ee8f2SPeter Brune   dmdasnes->objectivelocalctx = ctx;
4252a4ee8f2SPeter Brune   ierr = DMSNESSetObjective(dm,SNESComputeObjective_DMDA,dmdasnes);CHKERRQ(ierr);
4262a4ee8f2SPeter Brune   PetscFunctionReturn(0);
4272a4ee8f2SPeter Brune }
428