xref: /petsc/src/snes/utils/dmdasnes.c (revision b8d0cc3360b4ee2189edf50fdf404674eeb79938)
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*);
96cab3a1bSJed Brown   void *residuallocalctx;
106cab3a1bSJed Brown   void *jacobianlocalctx;
116cab3a1bSJed Brown } DM_DA_SNES;
126cab3a1bSJed Brown 
136cab3a1bSJed Brown #undef __FUNCT__
146cab3a1bSJed Brown #define __FUNCT__ "SNESDMDestroy_DMDA"
156cab3a1bSJed Brown static PetscErrorCode SNESDMDestroy_DMDA(SNESDM sdm)
166cab3a1bSJed Brown {
176cab3a1bSJed Brown   PetscErrorCode ierr;
186cab3a1bSJed Brown 
196cab3a1bSJed Brown   PetscFunctionBegin;
206cab3a1bSJed Brown   ierr = PetscFree(sdm->data);CHKERRQ(ierr);
216cab3a1bSJed Brown   PetscFunctionReturn(0);
226cab3a1bSJed Brown }
236cab3a1bSJed Brown 
246cab3a1bSJed Brown #undef __FUNCT__
256cab3a1bSJed Brown #define __FUNCT__ "DMDASNESGetContext"
266cab3a1bSJed Brown static PetscErrorCode DMDASNESGetContext(DM dm,SNESDM sdm,DM_DA_SNES **dmdasnes)
276cab3a1bSJed Brown {
286cab3a1bSJed Brown   PetscErrorCode ierr;
296cab3a1bSJed Brown 
306cab3a1bSJed Brown   PetscFunctionBegin;
31*b8d0cc33SJed Brown   *dmdasnes = PETSC_NULL;
326cab3a1bSJed Brown   if (!sdm->data) {
336cab3a1bSJed Brown     ierr = PetscNewLog(dm,DM_DA_SNES,&sdm->data);CHKERRQ(ierr);
346cab3a1bSJed Brown     sdm->destroy = SNESDMDestroy_DMDA;
356cab3a1bSJed Brown   }
366cab3a1bSJed Brown   *dmdasnes = (DM_DA_SNES*)sdm->data;
376cab3a1bSJed Brown   PetscFunctionReturn(0);
386cab3a1bSJed Brown }
396cab3a1bSJed Brown 
406cab3a1bSJed Brown #undef __FUNCT__
416cab3a1bSJed Brown #define __FUNCT__ "SNESComputeFunction_DMDA"
426cab3a1bSJed Brown /*
436cab3a1bSJed Brown   This function should eventually replace:
446cab3a1bSJed Brown     DMDAComputeFunction() and DMDAComputeFunction1()
456cab3a1bSJed Brown  */
466cab3a1bSJed Brown static PetscErrorCode SNESComputeFunction_DMDA(SNES snes,Vec X,Vec F,void *ctx)
476cab3a1bSJed Brown {
486cab3a1bSJed Brown   PetscErrorCode ierr;
496cab3a1bSJed Brown   DM             dm;
506cab3a1bSJed Brown   DM_DA_SNES     *dmdasnes = (DM_DA_SNES*)ctx;
516cab3a1bSJed Brown   DMDALocalInfo  info;
526cab3a1bSJed Brown   Vec            Xloc;
536cab3a1bSJed Brown   void           *x,*f;
546cab3a1bSJed Brown 
556cab3a1bSJed Brown   PetscFunctionBegin;
562961e153SJed Brown   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
572961e153SJed Brown   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
582961e153SJed Brown   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
596cab3a1bSJed Brown   if (!dmdasnes->residuallocal) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Corrupt context");
606cab3a1bSJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
616cab3a1bSJed Brown   ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
626cab3a1bSJed Brown   ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
636cab3a1bSJed Brown   ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
646cab3a1bSJed Brown   ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
656cab3a1bSJed Brown   ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
666cab3a1bSJed Brown   ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr);
676cab3a1bSJed Brown   CHKMEMQ;
686cab3a1bSJed Brown   ierr = (*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx);CHKERRQ(ierr);
696cab3a1bSJed Brown   CHKMEMQ;
706cab3a1bSJed Brown   ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
716cab3a1bSJed Brown   ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr);
726cab3a1bSJed Brown   ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
736cab3a1bSJed Brown   PetscFunctionReturn(0);
746cab3a1bSJed Brown }
756cab3a1bSJed Brown 
766cab3a1bSJed Brown #undef __FUNCT__
772961e153SJed Brown #define __FUNCT__ "SNESComputeJacobian_DMDA"
782961e153SJed Brown /*
792961e153SJed Brown   This function should eventually replace:
802961e153SJed Brown     DMComputeJacobian() and DMDAComputeJacobian1()
812961e153SJed Brown  */
822961e153SJed Brown static PetscErrorCode SNESComputeJacobian_DMDA(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *ctx)
832961e153SJed Brown {
842961e153SJed Brown   PetscErrorCode ierr;
852961e153SJed Brown   DM             dm;
862961e153SJed Brown   DM_DA_SNES     *dmdasnes = (DM_DA_SNES*)ctx;
872961e153SJed Brown   DMDALocalInfo  info;
882961e153SJed Brown   Vec            Xloc;
892961e153SJed Brown   void           *x;
902961e153SJed Brown 
912961e153SJed Brown   PetscFunctionBegin;
922961e153SJed Brown   if (!dmdasnes->residuallocal) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Corrupt context");
932961e153SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
942961e153SJed Brown 
952961e153SJed Brown   if (dmdasnes->jacobianlocal) {
962961e153SJed Brown     ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr);
972961e153SJed Brown     ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
982961e153SJed Brown     ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr);
992961e153SJed Brown     ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr);
1002961e153SJed Brown     ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr);
1012961e153SJed Brown     CHKMEMQ;
1022961e153SJed Brown     ierr = (*dmdasnes->jacobianlocal)(&info,x,*A,*B,mstr,dmdasnes->jacobianlocalctx);CHKERRQ(ierr);
1032961e153SJed Brown     CHKMEMQ;
1042961e153SJed Brown     ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr);
1052961e153SJed Brown     ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr);
1062961e153SJed Brown   } else {
1072961e153SJed Brown     MatFDColoring fdcoloring;
1082961e153SJed Brown     ierr = PetscObjectQuery((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject*)&fdcoloring);CHKERRQ(ierr);
1092961e153SJed Brown     if (!fdcoloring) {
1102961e153SJed Brown       ISColoring     coloring;
1112961e153SJed Brown 
1122961e153SJed Brown       ierr = DMCreateColoring(dm,dm->coloringtype,MATAIJ,&coloring);CHKERRQ(ierr);
1132961e153SJed Brown       ierr = MatFDColoringCreate(*B,coloring,&fdcoloring);CHKERRQ(ierr);
1142961e153SJed Brown       ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr);
1152961e153SJed Brown       switch (dm->coloringtype) {
1162961e153SJed Brown       case IS_COLORING_GLOBAL:
1172961e153SJed Brown         ierr = MatFDColoringSetFunction(fdcoloring,(PetscErrorCode(*)(void))SNESComputeFunction_DMDA,dmdasnes);CHKERRQ(ierr);
1182961e153SJed Brown         break;
1192961e153SJed Brown       default: SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"No support for coloring type '%s'",ISColoringTypes[dm->coloringtype]);
1202961e153SJed Brown       }
1212961e153SJed Brown       ierr = PetscObjectSetOptionsPrefix((PetscObject)fdcoloring,((PetscObject)dm)->prefix);CHKERRQ(ierr);
1222961e153SJed Brown       ierr = MatFDColoringSetFromOptions(fdcoloring);CHKERRQ(ierr);
1232961e153SJed Brown       ierr = PetscObjectCompose((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject)fdcoloring);CHKERRQ(ierr);
1242961e153SJed Brown       ierr = PetscObjectDereference((PetscObject)fdcoloring);CHKERRQ(ierr);
1252961e153SJed Brown 
1262961e153SJed Brown       /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call
1272961e153SJed Brown        * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the
1282961e153SJed Brown        * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually
1292961e153SJed Brown        * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be
1302961e153SJed Brown        * taken when the PetscOList for the Vec inside MatFDColoring is destroyed.
1312961e153SJed Brown        */
1322961e153SJed Brown       ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
1332961e153SJed Brown     }
1342961e153SJed Brown     *mstr = SAME_NONZERO_PATTERN;
1352961e153SJed Brown     ierr = MatFDColoringApply(*B,fdcoloring,X,mstr,snes);CHKERRQ(ierr);
1362961e153SJed Brown   }
1372961e153SJed Brown   /* This will be redundant if the user called both, but it's too common to forget. */
1382961e153SJed Brown   if (*A != *B) {
1392961e153SJed Brown     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1402961e153SJed Brown     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1412961e153SJed Brown   }
1422961e153SJed Brown   PetscFunctionReturn(0);
1432961e153SJed Brown }
1442961e153SJed Brown 
1452961e153SJed Brown #undef __FUNCT__
1466cab3a1bSJed Brown #define __FUNCT__ "DMDASNESSetFunctionLocal"
1476cab3a1bSJed Brown /*@C
1486cab3a1bSJed Brown    DMDASNESSetFunctionLocal - set a local residual evaluation function
1496cab3a1bSJed Brown 
1506cab3a1bSJed Brown    Logically Collective
1516cab3a1bSJed Brown 
1526cab3a1bSJed Brown    Input Arguments:
1536cab3a1bSJed Brown +  dm - DM to associate callback with
1546cab3a1bSJed Brown .  func - local residual evaluation
1556cab3a1bSJed Brown -  ctx - optional context for local residual evaluation
1566cab3a1bSJed Brown 
1576cab3a1bSJed Brown    Calling sequence for func:
1586cab3a1bSJed Brown +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
1596cab3a1bSJed Brown .  x - dimensional pointer to state at which to evaluate residual
1606cab3a1bSJed Brown .  f - dimensional pointer to residual, write the residual here
1616cab3a1bSJed Brown -  ctx - optional context passed above
1626cab3a1bSJed Brown 
1636cab3a1bSJed Brown    Level: beginner
1646cab3a1bSJed Brown 
1656cab3a1bSJed Brown .seealso: DMSNESSetFunction(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
1666cab3a1bSJed Brown @*/
1676cab3a1bSJed Brown PetscErrorCode DMDASNESSetFunctionLocal(DM dm,PetscErrorCode (*func)(DMDALocalInfo*,void*,void*,void*),void *ctx)
1686cab3a1bSJed Brown {
1696cab3a1bSJed Brown   PetscErrorCode ierr;
1706cab3a1bSJed Brown   SNESDM         sdm;
1716cab3a1bSJed Brown   DM_DA_SNES     *dmdasnes;
1726cab3a1bSJed Brown 
1736cab3a1bSJed Brown   PetscFunctionBegin;
1746cab3a1bSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1756cab3a1bSJed Brown   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
1766cab3a1bSJed Brown   ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr);
1776cab3a1bSJed Brown   dmdasnes->residuallocal = func;
1786cab3a1bSJed Brown   dmdasnes->residuallocalctx = ctx;
1796cab3a1bSJed Brown   ierr = DMSNESSetFunction(dm,SNESComputeFunction_DMDA,dmdasnes);CHKERRQ(ierr);
1802961e153SJed Brown   if (!sdm->computejacobian) {  /* Call us for the Jacobian too, can be overridden by the user. */
1812961e153SJed Brown     ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes);CHKERRQ(ierr);
1822961e153SJed Brown   }
1832961e153SJed Brown   PetscFunctionReturn(0);
1842961e153SJed Brown }
1852961e153SJed Brown 
1862961e153SJed Brown #undef __FUNCT__
1872961e153SJed Brown #define __FUNCT__ "DMDASNESSetJacobianLocal"
1882961e153SJed Brown /*@C
1892961e153SJed Brown    DMDASNESSetJacobianLocal - set a local residual evaluation function
1902961e153SJed Brown 
1912961e153SJed Brown    Logically Collective
1922961e153SJed Brown 
1932961e153SJed Brown    Input Arguments:
1942961e153SJed Brown +  dm - DM to associate callback with
1952961e153SJed Brown .  func - local residual evaluation
1962961e153SJed Brown -  ctx - optional context for local residual evaluation
1972961e153SJed Brown 
1982961e153SJed Brown    Calling sequence for func:
1992961e153SJed Brown +  info - DMDALocalInfo defining the subdomain to evaluate the residual on
2002961e153SJed Brown .  x - dimensional pointer to state at which to evaluate residual
2012961e153SJed Brown .  f - dimensional pointer to residual, write the residual here
2022961e153SJed Brown -  ctx - optional context passed above
2032961e153SJed Brown 
2042961e153SJed Brown    Level: beginner
2052961e153SJed Brown 
2062961e153SJed Brown .seealso: DMSNESSetJacobian(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
2072961e153SJed Brown @*/
2082961e153SJed Brown PetscErrorCode DMDASNESSetJacobianLocal(DM dm,PetscErrorCode (*func)(DMDALocalInfo*,void*,Mat,Mat,MatStructure*,void*),void *ctx)
2092961e153SJed Brown {
2102961e153SJed Brown   PetscErrorCode ierr;
2112961e153SJed Brown   SNESDM         sdm;
2122961e153SJed Brown   DM_DA_SNES     *dmdasnes;
2132961e153SJed Brown 
2142961e153SJed Brown   PetscFunctionBegin;
2152961e153SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2162961e153SJed Brown   ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr);
2172961e153SJed Brown   ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr);
2182961e153SJed Brown   dmdasnes->jacobianlocal = func;
2192961e153SJed Brown   dmdasnes->jacobianlocalctx = ctx;
2202961e153SJed Brown   ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes);CHKERRQ(ierr);
2216cab3a1bSJed Brown   PetscFunctionReturn(0);
2226cab3a1bSJed Brown }
223