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; 11*709ac0d2SJed Brown InsertMode residuallocalimode; 126cab3a1bSJed Brown } DM_DA_SNES; 136cab3a1bSJed Brown 146cab3a1bSJed Brown #undef __FUNCT__ 156cab3a1bSJed Brown #define __FUNCT__ "SNESDMDestroy_DMDA" 166cab3a1bSJed Brown static PetscErrorCode SNESDMDestroy_DMDA(SNESDM sdm) 176cab3a1bSJed Brown { 186cab3a1bSJed Brown PetscErrorCode ierr; 196cab3a1bSJed Brown 206cab3a1bSJed Brown PetscFunctionBegin; 216cab3a1bSJed Brown ierr = PetscFree(sdm->data);CHKERRQ(ierr); 226cab3a1bSJed Brown PetscFunctionReturn(0); 236cab3a1bSJed Brown } 246cab3a1bSJed Brown 256cab3a1bSJed Brown #undef __FUNCT__ 266cab3a1bSJed Brown #define __FUNCT__ "DMDASNESGetContext" 276cab3a1bSJed Brown static PetscErrorCode DMDASNESGetContext(DM dm,SNESDM sdm,DM_DA_SNES **dmdasnes) 286cab3a1bSJed Brown { 296cab3a1bSJed Brown PetscErrorCode ierr; 306cab3a1bSJed Brown 316cab3a1bSJed Brown PetscFunctionBegin; 32b8d0cc33SJed Brown *dmdasnes = PETSC_NULL; 336cab3a1bSJed Brown if (!sdm->data) { 346cab3a1bSJed Brown ierr = PetscNewLog(dm,DM_DA_SNES,&sdm->data);CHKERRQ(ierr); 356cab3a1bSJed Brown sdm->destroy = SNESDMDestroy_DMDA; 366cab3a1bSJed Brown } 376cab3a1bSJed Brown *dmdasnes = (DM_DA_SNES*)sdm->data; 386cab3a1bSJed Brown PetscFunctionReturn(0); 396cab3a1bSJed Brown } 406cab3a1bSJed Brown 416cab3a1bSJed Brown #undef __FUNCT__ 426cab3a1bSJed Brown #define __FUNCT__ "SNESComputeFunction_DMDA" 436cab3a1bSJed Brown /* 446cab3a1bSJed Brown This function should eventually replace: 456cab3a1bSJed Brown DMDAComputeFunction() and DMDAComputeFunction1() 466cab3a1bSJed Brown */ 476cab3a1bSJed Brown static PetscErrorCode SNESComputeFunction_DMDA(SNES snes,Vec X,Vec F,void *ctx) 486cab3a1bSJed Brown { 496cab3a1bSJed Brown PetscErrorCode ierr; 506cab3a1bSJed Brown DM dm; 516cab3a1bSJed Brown DM_DA_SNES *dmdasnes = (DM_DA_SNES*)ctx; 526cab3a1bSJed Brown DMDALocalInfo info; 536cab3a1bSJed Brown Vec Xloc; 546cab3a1bSJed Brown void *x,*f; 556cab3a1bSJed Brown 566cab3a1bSJed Brown PetscFunctionBegin; 572961e153SJed Brown PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 582961e153SJed Brown PetscValidHeaderSpecific(X,VEC_CLASSID,2); 592961e153SJed Brown PetscValidHeaderSpecific(F,VEC_CLASSID,3); 606cab3a1bSJed Brown if (!dmdasnes->residuallocal) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Corrupt context"); 616cab3a1bSJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 626cab3a1bSJed Brown ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 636cab3a1bSJed Brown ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 646cab3a1bSJed Brown ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 656cab3a1bSJed Brown ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 666cab3a1bSJed Brown ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr); 67*709ac0d2SJed Brown switch (dmdasnes->residuallocalimode) { 68*709ac0d2SJed Brown case INSERT_VALUES: { 696cab3a1bSJed Brown ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr); 706cab3a1bSJed Brown CHKMEMQ; 716cab3a1bSJed Brown ierr = (*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx);CHKERRQ(ierr); 726cab3a1bSJed Brown CHKMEMQ; 736cab3a1bSJed Brown ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr); 74*709ac0d2SJed Brown } break; 75*709ac0d2SJed Brown case ADD_VALUES: { 76*709ac0d2SJed Brown Vec Floc; 77*709ac0d2SJed Brown ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr); 78*709ac0d2SJed Brown ierr = VecZeroEntries(Floc);CHKERRQ(ierr); 79*709ac0d2SJed Brown ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr); 80*709ac0d2SJed Brown CHKMEMQ; 81*709ac0d2SJed Brown ierr = (*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx);CHKERRQ(ierr); 82*709ac0d2SJed Brown CHKMEMQ; 83*709ac0d2SJed Brown ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr); 84*709ac0d2SJed Brown ierr = VecZeroEntries(F);CHKERRQ(ierr); 85*709ac0d2SJed Brown ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 86*709ac0d2SJed Brown ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 87*709ac0d2SJed Brown ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr); 88*709ac0d2SJed Brown } break; 89*709ac0d2SJed Brown default: SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdasnes->residuallocalimode); 90*709ac0d2SJed Brown } 91*709ac0d2SJed Brown ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr); 926cab3a1bSJed Brown ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 936cab3a1bSJed Brown PetscFunctionReturn(0); 946cab3a1bSJed Brown } 956cab3a1bSJed Brown 966cab3a1bSJed Brown #undef __FUNCT__ 972961e153SJed Brown #define __FUNCT__ "SNESComputeJacobian_DMDA" 982961e153SJed Brown /* 992961e153SJed Brown This function should eventually replace: 1002961e153SJed Brown DMComputeJacobian() and DMDAComputeJacobian1() 1012961e153SJed Brown */ 1022961e153SJed Brown static PetscErrorCode SNESComputeJacobian_DMDA(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *ctx) 1032961e153SJed Brown { 1042961e153SJed Brown PetscErrorCode ierr; 1052961e153SJed Brown DM dm; 1062961e153SJed Brown DM_DA_SNES *dmdasnes = (DM_DA_SNES*)ctx; 1072961e153SJed Brown DMDALocalInfo info; 1082961e153SJed Brown Vec Xloc; 1092961e153SJed Brown void *x; 1102961e153SJed Brown 1112961e153SJed Brown PetscFunctionBegin; 1122961e153SJed Brown if (!dmdasnes->residuallocal) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Corrupt context"); 1132961e153SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 1142961e153SJed Brown 1152961e153SJed Brown if (dmdasnes->jacobianlocal) { 1162961e153SJed Brown ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 1172961e153SJed Brown ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 1182961e153SJed Brown ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 1192961e153SJed Brown ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 1202961e153SJed Brown ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr); 1212961e153SJed Brown CHKMEMQ; 1222961e153SJed Brown ierr = (*dmdasnes->jacobianlocal)(&info,x,*A,*B,mstr,dmdasnes->jacobianlocalctx);CHKERRQ(ierr); 1232961e153SJed Brown CHKMEMQ; 1242961e153SJed Brown ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr); 1252961e153SJed Brown ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 1262961e153SJed Brown } else { 1272961e153SJed Brown MatFDColoring fdcoloring; 1282961e153SJed Brown ierr = PetscObjectQuery((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject*)&fdcoloring);CHKERRQ(ierr); 1292961e153SJed Brown if (!fdcoloring) { 1302961e153SJed Brown ISColoring coloring; 1312961e153SJed Brown 1322961e153SJed Brown ierr = DMCreateColoring(dm,dm->coloringtype,MATAIJ,&coloring);CHKERRQ(ierr); 1332961e153SJed Brown ierr = MatFDColoringCreate(*B,coloring,&fdcoloring);CHKERRQ(ierr); 1342961e153SJed Brown ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr); 1352961e153SJed Brown switch (dm->coloringtype) { 1362961e153SJed Brown case IS_COLORING_GLOBAL: 1372961e153SJed Brown ierr = MatFDColoringSetFunction(fdcoloring,(PetscErrorCode(*)(void))SNESComputeFunction_DMDA,dmdasnes);CHKERRQ(ierr); 1382961e153SJed Brown break; 1392961e153SJed Brown default: SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"No support for coloring type '%s'",ISColoringTypes[dm->coloringtype]); 1402961e153SJed Brown } 1412961e153SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)fdcoloring,((PetscObject)dm)->prefix);CHKERRQ(ierr); 1422961e153SJed Brown ierr = MatFDColoringSetFromOptions(fdcoloring);CHKERRQ(ierr); 1432961e153SJed Brown ierr = PetscObjectCompose((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject)fdcoloring);CHKERRQ(ierr); 1442961e153SJed Brown ierr = PetscObjectDereference((PetscObject)fdcoloring);CHKERRQ(ierr); 1452961e153SJed Brown 1462961e153SJed Brown /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call 1472961e153SJed Brown * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the 1482961e153SJed Brown * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually 1492961e153SJed Brown * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be 1502961e153SJed Brown * taken when the PetscOList for the Vec inside MatFDColoring is destroyed. 1512961e153SJed Brown */ 1522961e153SJed Brown ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr); 1532961e153SJed Brown } 1542961e153SJed Brown *mstr = SAME_NONZERO_PATTERN; 1552961e153SJed Brown ierr = MatFDColoringApply(*B,fdcoloring,X,mstr,snes);CHKERRQ(ierr); 1562961e153SJed Brown } 1572961e153SJed Brown /* This will be redundant if the user called both, but it's too common to forget. */ 1582961e153SJed Brown if (*A != *B) { 1592961e153SJed Brown ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1602961e153SJed Brown ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1612961e153SJed Brown } 1622961e153SJed Brown PetscFunctionReturn(0); 1632961e153SJed Brown } 1642961e153SJed Brown 1652961e153SJed Brown #undef __FUNCT__ 1666cab3a1bSJed Brown #define __FUNCT__ "DMDASNESSetFunctionLocal" 1676cab3a1bSJed Brown /*@C 1686cab3a1bSJed Brown DMDASNESSetFunctionLocal - set a local residual evaluation function 1696cab3a1bSJed Brown 1706cab3a1bSJed Brown Logically Collective 1716cab3a1bSJed Brown 1726cab3a1bSJed Brown Input Arguments: 1736cab3a1bSJed Brown + dm - DM to associate callback with 1746cab3a1bSJed Brown . func - local residual evaluation 1756cab3a1bSJed Brown - ctx - optional context for local residual evaluation 1766cab3a1bSJed Brown 1776cab3a1bSJed Brown Calling sequence for func: 1786cab3a1bSJed Brown + info - DMDALocalInfo defining the subdomain to evaluate the residual on 179*709ac0d2SJed Brown . imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part 1806cab3a1bSJed Brown . x - dimensional pointer to state at which to evaluate residual 1816cab3a1bSJed Brown . f - dimensional pointer to residual, write the residual here 1826cab3a1bSJed Brown - ctx - optional context passed above 1836cab3a1bSJed Brown 1846cab3a1bSJed Brown Level: beginner 1856cab3a1bSJed Brown 1866cab3a1bSJed Brown .seealso: DMSNESSetFunction(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 1876cab3a1bSJed Brown @*/ 188*709ac0d2SJed Brown PetscErrorCode DMDASNESSetFunctionLocal(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,void*,void*,void*),void *ctx) 1896cab3a1bSJed Brown { 1906cab3a1bSJed Brown PetscErrorCode ierr; 1916cab3a1bSJed Brown SNESDM sdm; 1926cab3a1bSJed Brown DM_DA_SNES *dmdasnes; 1936cab3a1bSJed Brown 1946cab3a1bSJed Brown PetscFunctionBegin; 1956cab3a1bSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1966cab3a1bSJed Brown ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 1976cab3a1bSJed Brown ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr); 198*709ac0d2SJed Brown dmdasnes->residuallocalimode = imode; 1996cab3a1bSJed Brown dmdasnes->residuallocal = func; 2006cab3a1bSJed Brown dmdasnes->residuallocalctx = ctx; 2016cab3a1bSJed Brown ierr = DMSNESSetFunction(dm,SNESComputeFunction_DMDA,dmdasnes);CHKERRQ(ierr); 2022961e153SJed Brown if (!sdm->computejacobian) { /* Call us for the Jacobian too, can be overridden by the user. */ 2032961e153SJed Brown ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes);CHKERRQ(ierr); 2042961e153SJed Brown } 2052961e153SJed Brown PetscFunctionReturn(0); 2062961e153SJed Brown } 2072961e153SJed Brown 2082961e153SJed Brown #undef __FUNCT__ 2092961e153SJed Brown #define __FUNCT__ "DMDASNESSetJacobianLocal" 2102961e153SJed Brown /*@C 2112961e153SJed Brown DMDASNESSetJacobianLocal - set a local residual evaluation function 2122961e153SJed Brown 2132961e153SJed Brown Logically Collective 2142961e153SJed Brown 2152961e153SJed Brown Input Arguments: 2162961e153SJed Brown + dm - DM to associate callback with 2172961e153SJed Brown . func - local residual evaluation 2182961e153SJed Brown - ctx - optional context for local residual evaluation 2192961e153SJed Brown 2202961e153SJed Brown Calling sequence for func: 2212961e153SJed Brown + info - DMDALocalInfo defining the subdomain to evaluate the residual on 2222961e153SJed Brown . x - dimensional pointer to state at which to evaluate residual 2232961e153SJed Brown . f - dimensional pointer to residual, write the residual here 2242961e153SJed Brown - ctx - optional context passed above 2252961e153SJed Brown 2262961e153SJed Brown Level: beginner 2272961e153SJed Brown 2282961e153SJed Brown .seealso: DMSNESSetJacobian(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 2292961e153SJed Brown @*/ 2302961e153SJed Brown PetscErrorCode DMDASNESSetJacobianLocal(DM dm,PetscErrorCode (*func)(DMDALocalInfo*,void*,Mat,Mat,MatStructure*,void*),void *ctx) 2312961e153SJed Brown { 2322961e153SJed Brown PetscErrorCode ierr; 2332961e153SJed Brown SNESDM sdm; 2342961e153SJed Brown DM_DA_SNES *dmdasnes; 2352961e153SJed Brown 2362961e153SJed Brown PetscFunctionBegin; 2372961e153SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2382961e153SJed Brown ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 2392961e153SJed Brown ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr); 2402961e153SJed Brown dmdasnes->jacobianlocal = func; 2412961e153SJed Brown dmdasnes->jacobianlocalctx = ctx; 2422961e153SJed Brown ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes);CHKERRQ(ierr); 2436cab3a1bSJed Brown PetscFunctionReturn(0); 2446cab3a1bSJed Brown } 245