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; 11709ac0d2SJed 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__ 26*e3c0b8a2SPeter Brune #define __FUNCT__ "SNESDMDuplicate_DMDA" 27*e3c0b8a2SPeter Brune static PetscErrorCode SNESDMDuplicate_DMDA(SNESDM oldsdm,DM dm) 28*e3c0b8a2SPeter Brune { 29*e3c0b8a2SPeter Brune PetscErrorCode ierr; 30*e3c0b8a2SPeter Brune SNESDM sdm; 31*e3c0b8a2SPeter Brune 32*e3c0b8a2SPeter Brune PetscFunctionBegin; 33*e3c0b8a2SPeter Brune ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 34*e3c0b8a2SPeter Brune ierr = PetscNewLog(dm,DM_DA_SNES,&sdm->data);CHKERRQ(ierr); 35*e3c0b8a2SPeter Brune if (oldsdm->data)ierr = PetscMemcpy(sdm->data,oldsdm->data,sizeof(DM_DA_SNES));CHKERRQ(ierr); 36*e3c0b8a2SPeter Brune PetscFunctionReturn(0); 37*e3c0b8a2SPeter Brune } 38*e3c0b8a2SPeter Brune 39*e3c0b8a2SPeter Brune 40*e3c0b8a2SPeter Brune #undef __FUNCT__ 416cab3a1bSJed Brown #define __FUNCT__ "DMDASNESGetContext" 426cab3a1bSJed Brown static PetscErrorCode DMDASNESGetContext(DM dm,SNESDM sdm,DM_DA_SNES **dmdasnes) 436cab3a1bSJed Brown { 446cab3a1bSJed Brown PetscErrorCode ierr; 456cab3a1bSJed Brown 466cab3a1bSJed Brown PetscFunctionBegin; 47b8d0cc33SJed Brown *dmdasnes = PETSC_NULL; 486cab3a1bSJed Brown if (!sdm->data) { 496cab3a1bSJed Brown ierr = PetscNewLog(dm,DM_DA_SNES,&sdm->data);CHKERRQ(ierr); 506cab3a1bSJed Brown sdm->destroy = SNESDMDestroy_DMDA; 51*e3c0b8a2SPeter Brune sdm->duplicate = SNESDMDuplicate_DMDA; 526cab3a1bSJed Brown } 536cab3a1bSJed Brown *dmdasnes = (DM_DA_SNES*)sdm->data; 546cab3a1bSJed Brown PetscFunctionReturn(0); 556cab3a1bSJed Brown } 566cab3a1bSJed Brown 576cab3a1bSJed Brown #undef __FUNCT__ 586cab3a1bSJed Brown #define __FUNCT__ "SNESComputeFunction_DMDA" 596cab3a1bSJed Brown /* 606cab3a1bSJed Brown This function should eventually replace: 616cab3a1bSJed Brown DMDAComputeFunction() and DMDAComputeFunction1() 626cab3a1bSJed Brown */ 636cab3a1bSJed Brown static PetscErrorCode SNESComputeFunction_DMDA(SNES snes,Vec X,Vec F,void *ctx) 646cab3a1bSJed Brown { 656cab3a1bSJed Brown PetscErrorCode ierr; 666cab3a1bSJed Brown DM dm; 676cab3a1bSJed Brown DM_DA_SNES *dmdasnes = (DM_DA_SNES*)ctx; 686cab3a1bSJed Brown DMDALocalInfo info; 696cab3a1bSJed Brown Vec Xloc; 706cab3a1bSJed Brown void *x,*f; 716cab3a1bSJed Brown 726cab3a1bSJed Brown PetscFunctionBegin; 732961e153SJed Brown PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 742961e153SJed Brown PetscValidHeaderSpecific(X,VEC_CLASSID,2); 752961e153SJed Brown PetscValidHeaderSpecific(F,VEC_CLASSID,3); 766cab3a1bSJed Brown if (!dmdasnes->residuallocal) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Corrupt context"); 776cab3a1bSJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 786cab3a1bSJed Brown ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 796cab3a1bSJed Brown ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 806cab3a1bSJed Brown ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 816cab3a1bSJed Brown ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 826cab3a1bSJed Brown ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr); 83709ac0d2SJed Brown switch (dmdasnes->residuallocalimode) { 84709ac0d2SJed Brown case INSERT_VALUES: { 856cab3a1bSJed Brown ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr); 866cab3a1bSJed Brown CHKMEMQ; 876cab3a1bSJed Brown ierr = (*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx);CHKERRQ(ierr); 886cab3a1bSJed Brown CHKMEMQ; 896cab3a1bSJed Brown ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr); 90709ac0d2SJed Brown } break; 91709ac0d2SJed Brown case ADD_VALUES: { 92709ac0d2SJed Brown Vec Floc; 93709ac0d2SJed Brown ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr); 94709ac0d2SJed Brown ierr = VecZeroEntries(Floc);CHKERRQ(ierr); 95709ac0d2SJed Brown ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr); 96709ac0d2SJed Brown CHKMEMQ; 97709ac0d2SJed Brown ierr = (*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx);CHKERRQ(ierr); 98709ac0d2SJed Brown CHKMEMQ; 99709ac0d2SJed Brown ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr); 100709ac0d2SJed Brown ierr = VecZeroEntries(F);CHKERRQ(ierr); 101709ac0d2SJed Brown ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 102709ac0d2SJed Brown ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 103709ac0d2SJed Brown ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr); 104709ac0d2SJed Brown } break; 105709ac0d2SJed Brown default: SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdasnes->residuallocalimode); 106709ac0d2SJed Brown } 107709ac0d2SJed Brown ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr); 1086cab3a1bSJed Brown ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 1096cab3a1bSJed Brown PetscFunctionReturn(0); 1106cab3a1bSJed Brown } 1116cab3a1bSJed Brown 1126cab3a1bSJed Brown #undef __FUNCT__ 1136427ad5bSPeter Brune #define __FUNCT__ "SNESComputeLocalBlockFunction_DMDA" 1146427ad5bSPeter Brune /* 1156427ad5bSPeter Brune This function should eventually replace: 1166427ad5bSPeter Brune DMDAComputeFunction() and DMDAComputeFunction1() 1176427ad5bSPeter Brune */ 1186427ad5bSPeter Brune PetscErrorCode SNESComputeLocalBlockFunction_DMDA(SNES snes,Vec X,Vec F,void *dmctx) 1196427ad5bSPeter Brune { 1206427ad5bSPeter Brune PetscErrorCode ierr; 1216427ad5bSPeter Brune DM dm = (DM)dmctx; /* the global DMDA context */ 1226427ad5bSPeter Brune DMDALocalInfo info; 1236427ad5bSPeter Brune void *x,*f; 1246427ad5bSPeter Brune DM_DA_SNES *dmdasnes; 1256427ad5bSPeter Brune SNESDM sdm; 1266427ad5bSPeter Brune 1276427ad5bSPeter Brune PetscFunctionBegin; 1286427ad5bSPeter Brune 1296427ad5bSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1306427ad5bSPeter Brune PetscValidHeaderSpecific(X,VEC_CLASSID,2); 1316427ad5bSPeter Brune PetscValidHeaderSpecific(F,VEC_CLASSID,3); 1326427ad5bSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,4); 1336427ad5bSPeter Brune 1346427ad5bSPeter Brune if (!dmdasnes->residuallocal) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Corrupt context"); 1356427ad5bSPeter Brune ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 1366427ad5bSPeter Brune ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr); 1376427ad5bSPeter Brune ierr = DMDAVecGetArray(dm,X,&x);CHKERRQ(ierr); 1386427ad5bSPeter Brune ierr = DMDAGetLocalBlockInfo(dm,&info);CHKERRQ(ierr); 1396427ad5bSPeter Brune switch (dmdasnes->residuallocalimode) { 1406427ad5bSPeter Brune case INSERT_VALUES: { 1416427ad5bSPeter Brune ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr); 1426427ad5bSPeter Brune CHKMEMQ; 1436427ad5bSPeter Brune ierr = (*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx);CHKERRQ(ierr); 1446427ad5bSPeter Brune CHKMEMQ; 1456427ad5bSPeter Brune ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr); 1466427ad5bSPeter Brune } break; 1476427ad5bSPeter Brune case ADD_VALUES: { 1486427ad5bSPeter Brune Vec Floc; 1496427ad5bSPeter Brune ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr); 1506427ad5bSPeter Brune ierr = VecZeroEntries(Floc);CHKERRQ(ierr); 1516427ad5bSPeter Brune ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr); 1526427ad5bSPeter Brune CHKMEMQ; 1536427ad5bSPeter Brune ierr = (*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx);CHKERRQ(ierr); 1546427ad5bSPeter Brune CHKMEMQ; 1556427ad5bSPeter Brune ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr); 1566427ad5bSPeter Brune ierr = VecZeroEntries(F);CHKERRQ(ierr); 1576427ad5bSPeter Brune ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 1586427ad5bSPeter Brune ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 1596427ad5bSPeter Brune ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr); 1606427ad5bSPeter Brune } break; 1616427ad5bSPeter Brune default: SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdasnes->residuallocalimode); 1626427ad5bSPeter Brune } 1636427ad5bSPeter Brune ierr = DMDAVecRestoreArray(dm,X,&x);CHKERRQ(ierr); 1646427ad5bSPeter Brune PetscFunctionReturn(0); 1656427ad5bSPeter Brune } 1666427ad5bSPeter Brune 1676427ad5bSPeter Brune #undef __FUNCT__ 1682961e153SJed Brown #define __FUNCT__ "SNESComputeJacobian_DMDA" 1692961e153SJed Brown /* 1702961e153SJed Brown This function should eventually replace: 1712961e153SJed Brown DMComputeJacobian() and DMDAComputeJacobian1() 1722961e153SJed Brown */ 1732961e153SJed Brown static PetscErrorCode SNESComputeJacobian_DMDA(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *ctx) 1742961e153SJed Brown { 1752961e153SJed Brown PetscErrorCode ierr; 1762961e153SJed Brown DM dm; 1772961e153SJed Brown DM_DA_SNES *dmdasnes = (DM_DA_SNES*)ctx; 1782961e153SJed Brown DMDALocalInfo info; 1792961e153SJed Brown Vec Xloc; 1802961e153SJed Brown void *x; 1812961e153SJed Brown 1822961e153SJed Brown PetscFunctionBegin; 1832961e153SJed Brown if (!dmdasnes->residuallocal) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Corrupt context"); 1842961e153SJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 1852961e153SJed Brown 1862961e153SJed Brown if (dmdasnes->jacobianlocal) { 1872961e153SJed Brown ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 1882961e153SJed Brown ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 1892961e153SJed Brown ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 1902961e153SJed Brown ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 1912961e153SJed Brown ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr); 1922961e153SJed Brown CHKMEMQ; 1932961e153SJed Brown ierr = (*dmdasnes->jacobianlocal)(&info,x,*A,*B,mstr,dmdasnes->jacobianlocalctx);CHKERRQ(ierr); 1942961e153SJed Brown CHKMEMQ; 1952961e153SJed Brown ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr); 1962961e153SJed Brown ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 1972961e153SJed Brown } else { 1982961e153SJed Brown MatFDColoring fdcoloring; 1992961e153SJed Brown ierr = PetscObjectQuery((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject*)&fdcoloring);CHKERRQ(ierr); 2002961e153SJed Brown if (!fdcoloring) { 2012961e153SJed Brown ISColoring coloring; 2022961e153SJed Brown 2032961e153SJed Brown ierr = DMCreateColoring(dm,dm->coloringtype,MATAIJ,&coloring);CHKERRQ(ierr); 2042961e153SJed Brown ierr = MatFDColoringCreate(*B,coloring,&fdcoloring);CHKERRQ(ierr); 2052961e153SJed Brown ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr); 2062961e153SJed Brown switch (dm->coloringtype) { 2072961e153SJed Brown case IS_COLORING_GLOBAL: 2082961e153SJed Brown ierr = MatFDColoringSetFunction(fdcoloring,(PetscErrorCode(*)(void))SNESComputeFunction_DMDA,dmdasnes);CHKERRQ(ierr); 2092961e153SJed Brown break; 2102961e153SJed Brown default: SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"No support for coloring type '%s'",ISColoringTypes[dm->coloringtype]); 2112961e153SJed Brown } 2122961e153SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)fdcoloring,((PetscObject)dm)->prefix);CHKERRQ(ierr); 2132961e153SJed Brown ierr = MatFDColoringSetFromOptions(fdcoloring);CHKERRQ(ierr); 2142961e153SJed Brown ierr = PetscObjectCompose((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject)fdcoloring);CHKERRQ(ierr); 2152961e153SJed Brown ierr = PetscObjectDereference((PetscObject)fdcoloring);CHKERRQ(ierr); 2162961e153SJed Brown 2172961e153SJed Brown /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call 2182961e153SJed Brown * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the 2192961e153SJed Brown * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually 2202961e153SJed Brown * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be 2212961e153SJed Brown * taken when the PetscOList for the Vec inside MatFDColoring is destroyed. 2222961e153SJed Brown */ 2232961e153SJed Brown ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr); 2242961e153SJed Brown } 2252961e153SJed Brown *mstr = SAME_NONZERO_PATTERN; 2262961e153SJed Brown ierr = MatFDColoringApply(*B,fdcoloring,X,mstr,snes);CHKERRQ(ierr); 2272961e153SJed Brown } 2282961e153SJed Brown /* This will be redundant if the user called both, but it's too common to forget. */ 2292961e153SJed Brown if (*A != *B) { 2302961e153SJed Brown ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2312961e153SJed Brown ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2322961e153SJed Brown } 2332961e153SJed Brown PetscFunctionReturn(0); 2342961e153SJed Brown } 2352961e153SJed Brown 2362961e153SJed Brown #undef __FUNCT__ 2376427ad5bSPeter Brune #define __FUNCT__ "SNESComputeJacobian_DMDA" 2386427ad5bSPeter Brune /* 2396427ad5bSPeter Brune This function should eventually replace: 2406427ad5bSPeter Brune DMComputeJacobian() and DMDAComputeJacobian1() 2416427ad5bSPeter Brune */ 2426427ad5bSPeter Brune PetscErrorCode SNESComputeLocalBlockJacobian_DMDA(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *dmctx) 2436427ad5bSPeter Brune { 2446427ad5bSPeter Brune PetscErrorCode ierr; 2456427ad5bSPeter Brune DM dm = (DM)dmctx; /* the global DMDA context */ 2466427ad5bSPeter Brune DM_DA_SNES *dmdasnes; 2476427ad5bSPeter Brune SNESDM sdm; 2486427ad5bSPeter Brune DMDALocalInfo info; 2496427ad5bSPeter Brune void *x; 2506427ad5bSPeter Brune 2516427ad5bSPeter Brune PetscFunctionBegin; 2526427ad5bSPeter Brune if (!dmdasnes->residuallocal) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"Corrupt context"); 2536427ad5bSPeter Brune ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 2546427ad5bSPeter Brune ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr); 2556427ad5bSPeter Brune if (dmdasnes->jacobianlocal) { 2566427ad5bSPeter Brune ierr = DMDAVecGetArray(dm,X,&x);CHKERRQ(ierr); 2576427ad5bSPeter Brune CHKMEMQ; 2586427ad5bSPeter Brune ierr = (*dmdasnes->jacobianlocal)(&info,x,*A,*B,mstr,dmdasnes->jacobianlocalctx);CHKERRQ(ierr); 2596427ad5bSPeter Brune CHKMEMQ; 2606427ad5bSPeter Brune ierr = DMDAVecRestoreArray(dm,X,&x);CHKERRQ(ierr); 2616427ad5bSPeter Brune } 2626427ad5bSPeter Brune /* This will be redundant if the user called both, but it's too common to forget. */ 2636427ad5bSPeter Brune if (*A != *B) { 2646427ad5bSPeter Brune ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2656427ad5bSPeter Brune ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2666427ad5bSPeter Brune } 2676427ad5bSPeter Brune PetscFunctionReturn(0); 2686427ad5bSPeter Brune } 2696427ad5bSPeter Brune 2706427ad5bSPeter Brune 2716427ad5bSPeter Brune #undef __FUNCT__ 2726cab3a1bSJed Brown #define __FUNCT__ "DMDASNESSetFunctionLocal" 2736cab3a1bSJed Brown /*@C 2746cab3a1bSJed Brown DMDASNESSetFunctionLocal - set a local residual evaluation function 2756cab3a1bSJed Brown 2766cab3a1bSJed Brown Logically Collective 2776cab3a1bSJed Brown 2786cab3a1bSJed Brown Input Arguments: 2796cab3a1bSJed Brown + dm - DM to associate callback with 2806cab3a1bSJed Brown . func - local residual evaluation 2816cab3a1bSJed Brown - ctx - optional context for local residual evaluation 2826cab3a1bSJed Brown 2836cab3a1bSJed Brown Calling sequence for func: 2846cab3a1bSJed Brown + info - DMDALocalInfo defining the subdomain to evaluate the residual on 285709ac0d2SJed Brown . imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part 2866cab3a1bSJed Brown . x - dimensional pointer to state at which to evaluate residual 2876cab3a1bSJed Brown . f - dimensional pointer to residual, write the residual here 2886cab3a1bSJed Brown - ctx - optional context passed above 2896cab3a1bSJed Brown 2906cab3a1bSJed Brown Level: beginner 2916cab3a1bSJed Brown 2926cab3a1bSJed Brown .seealso: DMSNESSetFunction(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 2936cab3a1bSJed Brown @*/ 294709ac0d2SJed Brown PetscErrorCode DMDASNESSetFunctionLocal(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,void*,void*,void*),void *ctx) 2956cab3a1bSJed Brown { 2966cab3a1bSJed Brown PetscErrorCode ierr; 2976cab3a1bSJed Brown SNESDM sdm; 2986cab3a1bSJed Brown DM_DA_SNES *dmdasnes; 2996cab3a1bSJed Brown 3006cab3a1bSJed Brown PetscFunctionBegin; 3016cab3a1bSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3026cab3a1bSJed Brown ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 3036cab3a1bSJed Brown ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr); 304709ac0d2SJed Brown dmdasnes->residuallocalimode = imode; 3056cab3a1bSJed Brown dmdasnes->residuallocal = func; 3066cab3a1bSJed Brown dmdasnes->residuallocalctx = ctx; 3076cab3a1bSJed Brown ierr = DMSNESSetFunction(dm,SNESComputeFunction_DMDA,dmdasnes);CHKERRQ(ierr); 3086427ad5bSPeter Brune ierr = DMSNESSetBlockFunction(dm,SNESComputeLocalBlockFunction_DMDA,dm);CHKERRQ(ierr); 3092961e153SJed Brown if (!sdm->computejacobian) { /* Call us for the Jacobian too, can be overridden by the user. */ 3102961e153SJed Brown ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes);CHKERRQ(ierr); 3112961e153SJed Brown } 3122961e153SJed Brown PetscFunctionReturn(0); 3132961e153SJed Brown } 3142961e153SJed Brown 3152961e153SJed Brown #undef __FUNCT__ 3162961e153SJed Brown #define __FUNCT__ "DMDASNESSetJacobianLocal" 3172961e153SJed Brown /*@C 3182961e153SJed Brown DMDASNESSetJacobianLocal - set a local residual evaluation function 3192961e153SJed Brown 3202961e153SJed Brown Logically Collective 3212961e153SJed Brown 3222961e153SJed Brown Input Arguments: 3232961e153SJed Brown + dm - DM to associate callback with 3242961e153SJed Brown . func - local residual evaluation 3252961e153SJed Brown - ctx - optional context for local residual evaluation 3262961e153SJed Brown 3272961e153SJed Brown Calling sequence for func: 3282961e153SJed Brown + info - DMDALocalInfo defining the subdomain to evaluate the residual on 3292961e153SJed Brown . x - dimensional pointer to state at which to evaluate residual 3302961e153SJed Brown . f - dimensional pointer to residual, write the residual here 3312961e153SJed Brown - ctx - optional context passed above 3322961e153SJed Brown 3332961e153SJed Brown Level: beginner 3342961e153SJed Brown 3352961e153SJed Brown .seealso: DMSNESSetJacobian(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 3362961e153SJed Brown @*/ 3372961e153SJed Brown PetscErrorCode DMDASNESSetJacobianLocal(DM dm,PetscErrorCode (*func)(DMDALocalInfo*,void*,Mat,Mat,MatStructure*,void*),void *ctx) 3382961e153SJed Brown { 3392961e153SJed Brown PetscErrorCode ierr; 3402961e153SJed Brown SNESDM sdm; 3412961e153SJed Brown DM_DA_SNES *dmdasnes; 3422961e153SJed Brown 3432961e153SJed Brown PetscFunctionBegin; 3442961e153SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3452961e153SJed Brown ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 3462961e153SJed Brown ierr = DMDASNESGetContext(dm,sdm,&dmdasnes);CHKERRQ(ierr); 3472961e153SJed Brown dmdasnes->jacobianlocal = func; 3482961e153SJed Brown dmdasnes->jacobianlocalctx = ctx; 3492961e153SJed Brown ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes);CHKERRQ(ierr); 3506427ad5bSPeter Brune ierr = DMSNESSetBlockJacobian(dm,SNESComputeLocalBlockJacobian_DMDA,dm);CHKERRQ(ierr); 3516cab3a1bSJed Brown PetscFunctionReturn(0); 3526cab3a1bSJed Brown } 353