16cab3a1bSJed Brown #include <petscdmda.h> /*I "petscdmda.h" I*/ 2af0996ceSBarry Smith #include <petsc/private/dmimpl.h> 3af0996ceSBarry Smith #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*); 8d1e9a80fSBarry Smith PetscErrorCode (*jacobianlocal)(DMDALocalInfo*,void*,Mat,Mat,void*); 92a4ee8f2SPeter Brune PetscErrorCode (*objectivelocal)(DMDALocalInfo*,void*,PetscReal*,void*); 105505f7afSJunchao Zhang 115505f7afSJunchao Zhang PetscErrorCode (*residuallocalvec)(DMDALocalInfo*,Vec,Vec,void*); 125505f7afSJunchao Zhang PetscErrorCode (*jacobianlocalvec)(DMDALocalInfo*,Vec,Mat,Mat,void*); 135505f7afSJunchao Zhang PetscErrorCode (*objectivelocalvec)(DMDALocalInfo*,Vec,PetscReal*,void*); 146cab3a1bSJed Brown void *residuallocalctx; 156cab3a1bSJed Brown void *jacobianlocalctx; 162a4ee8f2SPeter Brune void *objectivelocalctx; 17709ac0d2SJed Brown InsertMode residuallocalimode; 18dcad5f1cSBarry Smith 19dcad5f1cSBarry Smith /* For Picard iteration defined locally */ 20dcad5f1cSBarry Smith PetscErrorCode (*rhsplocal)(DMDALocalInfo*,void*,void*,void*); 21d1e9a80fSBarry Smith PetscErrorCode (*jacobianplocal)(DMDALocalInfo*,void*,Mat,Mat,void*); 22dcad5f1cSBarry Smith void *picardlocalctx; 23942e3340SBarry Smith } DMSNES_DA; 246cab3a1bSJed Brown 25942e3340SBarry Smith static PetscErrorCode DMSNESDestroy_DMDA(DMSNES sdm) 266cab3a1bSJed Brown { 276cab3a1bSJed Brown PetscFunctionBegin; 289566063dSJacob Faibussowitsch PetscCall(PetscFree(sdm->data)); 296cab3a1bSJed Brown PetscFunctionReturn(0); 306cab3a1bSJed Brown } 316cab3a1bSJed Brown 3222c6f798SBarry Smith static PetscErrorCode DMSNESDuplicate_DMDA(DMSNES oldsdm,DMSNES sdm) 33e3c0b8a2SPeter Brune { 34e3c0b8a2SPeter Brune PetscFunctionBegin; 359566063dSJacob Faibussowitsch PetscCall(PetscNewLog(sdm,(DMSNES_DA**)&sdm->data)); 361baa6e33SBarry Smith if (oldsdm->data) PetscCall(PetscMemcpy(sdm->data,oldsdm->data,sizeof(DMSNES_DA))); 37e3c0b8a2SPeter Brune PetscFunctionReturn(0); 38e3c0b8a2SPeter Brune } 39e3c0b8a2SPeter Brune 40942e3340SBarry Smith static PetscErrorCode DMDASNESGetContext(DM dm,DMSNES sdm,DMSNES_DA **dmdasnes) 416cab3a1bSJed Brown { 426cab3a1bSJed Brown PetscFunctionBegin; 430298fd71SBarry Smith *dmdasnes = NULL; 446cab3a1bSJed Brown if (!sdm->data) { 459566063dSJacob Faibussowitsch PetscCall(PetscNewLog(dm,(DMSNES_DA**)&sdm->data)); 4622c6f798SBarry Smith sdm->ops->destroy = DMSNESDestroy_DMDA; 4722c6f798SBarry Smith sdm->ops->duplicate = DMSNESDuplicate_DMDA; 486cab3a1bSJed Brown } 49942e3340SBarry Smith *dmdasnes = (DMSNES_DA*)sdm->data; 506cab3a1bSJed Brown PetscFunctionReturn(0); 516cab3a1bSJed Brown } 526cab3a1bSJed Brown 536cab3a1bSJed Brown static PetscErrorCode SNESComputeFunction_DMDA(SNES snes,Vec X,Vec F,void *ctx) 546cab3a1bSJed Brown { 556cab3a1bSJed Brown DM dm; 56942e3340SBarry Smith DMSNES_DA *dmdasnes = (DMSNES_DA*)ctx; 576cab3a1bSJed Brown DMDALocalInfo info; 586cab3a1bSJed Brown Vec Xloc; 596cab3a1bSJed Brown void *x,*f; 606cab3a1bSJed Brown 616cab3a1bSJed Brown PetscFunctionBegin; 622961e153SJed Brown PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 632961e153SJed Brown PetscValidHeaderSpecific(X,VEC_CLASSID,2); 642961e153SJed Brown PetscValidHeaderSpecific(F,VEC_CLASSID,3); 655505f7afSJunchao Zhang PetscCheck(dmdasnes->residuallocal || dmdasnes->residuallocalvec,PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context"); 669566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes,&dm)); 679566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm,&Xloc)); 689566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc)); 699566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc)); 709566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(dm,&info)); 71709ac0d2SJed Brown switch (dmdasnes->residuallocalimode) { 72709ac0d2SJed Brown case INSERT_VALUES: { 739566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(SNES_FunctionEval,snes,X,F,0)); 74*792fecdfSBarry Smith if (dmdasnes->residuallocalvec) PetscCallBack("SNES DMDA local callback function",(*dmdasnes->residuallocalvec)(&info,Xloc,F,dmdasnes->residuallocalctx)); 75ef1023bdSBarry Smith else { 765505f7afSJunchao Zhang PetscCall(DMDAVecGetArray(dm,Xloc,&x)); 775505f7afSJunchao Zhang PetscCall(DMDAVecGetArray(dm,F,&f)); 78*792fecdfSBarry Smith PetscCallBack("SNES DMDA local callback function",(*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx)); 795505f7afSJunchao Zhang PetscCall(DMDAVecRestoreArray(dm,Xloc,&x)); 805505f7afSJunchao Zhang PetscCall(DMDAVecRestoreArray(dm,F,&f)); 815505f7afSJunchao Zhang } 829566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(SNES_FunctionEval,snes,X,F,0)); 83709ac0d2SJed Brown } break; 84709ac0d2SJed Brown case ADD_VALUES: { 85709ac0d2SJed Brown Vec Floc; 869566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm,&Floc)); 879566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Floc)); 889566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(SNES_FunctionEval,snes,X,F,0)); 89*792fecdfSBarry Smith if (dmdasnes->residuallocalvec) PetscCallBack("SNES DMDA local callback function",(*dmdasnes->residuallocalvec)(&info,Xloc,Floc,dmdasnes->residuallocalctx)); 90ef1023bdSBarry Smith else { 915505f7afSJunchao Zhang PetscCall(DMDAVecGetArray(dm,Xloc,&x)); 925505f7afSJunchao Zhang PetscCall(DMDAVecGetArray(dm,Floc,&f)); 93*792fecdfSBarry Smith PetscCallBack("SNES DMDA local callback function",(*dmdasnes->residuallocal)(&info,x,f,dmdasnes->residuallocalctx)); 945505f7afSJunchao Zhang PetscCall(DMDAVecRestoreArray(dm,Xloc,&x)); 955505f7afSJunchao Zhang PetscCall(DMDAVecRestoreArray(dm,Floc,&f)); 965505f7afSJunchao Zhang } 979566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(SNES_FunctionEval,snes,X,F,0)); 989566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 999566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F)); 1009566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F)); 1019566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm,&Floc)); 102709ac0d2SJed Brown } break; 10398921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdasnes->residuallocalimode); 104709ac0d2SJed Brown } 1059566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm,&Xloc)); 1061baa6e33SBarry Smith if (snes->domainerror) PetscCall(VecSetInf(F)); 1076cab3a1bSJed Brown PetscFunctionReturn(0); 1086cab3a1bSJed Brown } 1096cab3a1bSJed Brown 1102a4ee8f2SPeter Brune static PetscErrorCode SNESComputeObjective_DMDA(SNES snes,Vec X,PetscReal *ob,void *ctx) 1112a4ee8f2SPeter Brune { 1122a4ee8f2SPeter Brune DM dm; 113942e3340SBarry Smith DMSNES_DA *dmdasnes = (DMSNES_DA*)ctx; 1142a4ee8f2SPeter Brune DMDALocalInfo info; 1152a4ee8f2SPeter Brune Vec Xloc; 1162a4ee8f2SPeter Brune void *x; 1172a4ee8f2SPeter Brune 1182a4ee8f2SPeter Brune PetscFunctionBegin; 1192a4ee8f2SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1202a4ee8f2SPeter Brune PetscValidHeaderSpecific(X,VEC_CLASSID,2); 121dadcf809SJacob Faibussowitsch PetscValidRealPointer(ob,3); 1225505f7afSJunchao Zhang PetscCheck(dmdasnes->objectivelocal || dmdasnes->objectivelocalvec,PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context"); 1239566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes,&dm)); 1249566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm,&Xloc)); 1259566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc)); 1269566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc)); 1279566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(dm,&info)); 128*792fecdfSBarry Smith if (dmdasnes->objectivelocalvec) PetscCallBack("SNES DMDA local callback objective",(*dmdasnes->objectivelocalvec)(&info,Xloc,ob,dmdasnes->objectivelocalctx)); 129ef1023bdSBarry Smith else { 1309566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(dm,Xloc,&x)); 131*792fecdfSBarry Smith PetscCallBack("SNES DMDA local callback objective",(*dmdasnes->objectivelocal)(&info,x,ob,dmdasnes->objectivelocalctx)); 1329566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dm,Xloc,&x)); 1335505f7afSJunchao Zhang } 1349566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm,&Xloc)); 1352a4ee8f2SPeter Brune PetscFunctionReturn(0); 1362a4ee8f2SPeter Brune } 1372a4ee8f2SPeter Brune 138cc2e6a90SBarry Smith /* Routine is called by example, hence must be labeled PETSC_EXTERN */ 139cc2e6a90SBarry Smith PETSC_EXTERN PetscErrorCode SNESComputeJacobian_DMDA(SNES snes,Vec X,Mat A,Mat B,void *ctx) 1402961e153SJed Brown { 1412961e153SJed Brown DM dm; 142942e3340SBarry Smith DMSNES_DA *dmdasnes = (DMSNES_DA*)ctx; 1432961e153SJed Brown DMDALocalInfo info; 1442961e153SJed Brown Vec Xloc; 1452961e153SJed Brown void *x; 1462961e153SJed Brown 1472961e153SJed Brown PetscFunctionBegin; 1485505f7afSJunchao Zhang PetscCheck(dmdasnes->residuallocal || dmdasnes->residuallocalvec,PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context"); 1499566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes,&dm)); 1502961e153SJed Brown 1515505f7afSJunchao Zhang if (dmdasnes->jacobianlocal || dmdasnes->jacobianlocalctx) { 1529566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm,&Xloc)); 1539566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc)); 1549566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc)); 1559566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(dm,&info)); 156*792fecdfSBarry Smith if (dmdasnes->jacobianlocalvec) PetscCallBack("SNES DMDA local callback Jacobian",(*dmdasnes->jacobianlocalvec)(&info,Xloc,A,B,dmdasnes->jacobianlocalctx)); 157ef1023bdSBarry Smith else { 1589566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(dm,Xloc,&x)); 159*792fecdfSBarry Smith PetscCallBack("SNES DMDA local callback Jacobian",(*dmdasnes->jacobianlocal)(&info,x,A,B,dmdasnes->jacobianlocalctx)); 1609566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dm,Xloc,&x)); 1615505f7afSJunchao Zhang } 1629566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm,&Xloc)); 1632961e153SJed Brown } else { 1642961e153SJed Brown MatFDColoring fdcoloring; 1659566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject*)&fdcoloring)); 1662961e153SJed Brown if (!fdcoloring) { 1672961e153SJed Brown ISColoring coloring; 1682961e153SJed Brown 1699566063dSJacob Faibussowitsch PetscCall(DMCreateColoring(dm,dm->coloringtype,&coloring)); 1709566063dSJacob Faibussowitsch PetscCall(MatFDColoringCreate(B,coloring,&fdcoloring)); 1712961e153SJed Brown switch (dm->coloringtype) { 1722961e153SJed Brown case IS_COLORING_GLOBAL: 1739566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetFunction(fdcoloring,(PetscErrorCode (*)(void))SNESComputeFunction_DMDA,dmdasnes)); 1742961e153SJed Brown break; 17598921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"No support for coloring type '%s'",ISColoringTypes[dm->coloringtype]); 1762961e153SJed Brown } 1779566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fdcoloring,((PetscObject)dm)->prefix)); 1789566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetFromOptions(fdcoloring)); 1799566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetUp(B,coloring,fdcoloring)); 1809566063dSJacob Faibussowitsch PetscCall(ISColoringDestroy(&coloring)); 1819566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject)fdcoloring)); 1829566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)fdcoloring)); 1832961e153SJed Brown 1842961e153SJed Brown /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call 1852961e153SJed Brown * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the 1862961e153SJed Brown * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually 1872961e153SJed Brown * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be 188140e18c1SBarry Smith * taken when the PetscObjectList for the Vec inside MatFDColoring is destroyed. 1892961e153SJed Brown */ 1909566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)dm)); 1912961e153SJed Brown } 1929566063dSJacob Faibussowitsch PetscCall(MatFDColoringApply(B,fdcoloring,X,snes)); 1932961e153SJed Brown } 1942961e153SJed Brown /* This will be redundant if the user called both, but it's too common to forget. */ 19594ab13aaSBarry Smith if (A != B) { 1969566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 1979566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 1982961e153SJed Brown } 1992961e153SJed Brown PetscFunctionReturn(0); 2002961e153SJed Brown } 2012961e153SJed Brown 2026cab3a1bSJed Brown /*@C 2036cab3a1bSJed Brown DMDASNESSetFunctionLocal - set a local residual evaluation function 2046cab3a1bSJed Brown 2056cab3a1bSJed Brown Logically Collective 2066cab3a1bSJed Brown 2074165533cSJose E. Roman Input Parameters: 2086cab3a1bSJed Brown + dm - DM to associate callback with 209e3c51ba2SJed Brown . imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part 2106cab3a1bSJed Brown . func - local residual evaluation 2116cab3a1bSJed Brown - ctx - optional context for local residual evaluation 2126cab3a1bSJed Brown 2130038884cSBarry Smith Calling sequence: 2140038884cSBarry Smith For PetscErrorCode (*func)(DMDALocalInfo *info,void *x, void *f, void *ctx), 2156cab3a1bSJed Brown + info - DMDALocalInfo defining the subdomain to evaluate the residual on 2160038884cSBarry Smith . x - dimensional pointer to state at which to evaluate residual (e.g. PetscScalar *x or **x or ***x) 2170038884cSBarry Smith . f - dimensional pointer to residual, write the residual here (e.g. PetscScalar *f or **f or ***f) 2186cab3a1bSJed Brown - ctx - optional context passed above 2196cab3a1bSJed Brown 2206cab3a1bSJed Brown Level: beginner 2216cab3a1bSJed Brown 222db781477SPatrick Sanan .seealso: `DMDASNESSetJacobianLocal()`, `DMSNESSetFunction()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()` 2236cab3a1bSJed Brown @*/ 224709ac0d2SJed Brown PetscErrorCode DMDASNESSetFunctionLocal(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,void*,void*,void*),void *ctx) 2256cab3a1bSJed Brown { 226942e3340SBarry Smith DMSNES sdm; 227942e3340SBarry Smith DMSNES_DA *dmdasnes; 2286cab3a1bSJed Brown 2296cab3a1bSJed Brown PetscFunctionBegin; 2306cab3a1bSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2319566063dSJacob Faibussowitsch PetscCall(DMGetDMSNESWrite(dm,&sdm)); 2329566063dSJacob Faibussowitsch PetscCall(DMDASNESGetContext(dm,sdm,&dmdasnes)); 2331aa26658SKarl Rupp 234709ac0d2SJed Brown dmdasnes->residuallocalimode = imode; 2356cab3a1bSJed Brown dmdasnes->residuallocal = func; 2366cab3a1bSJed Brown dmdasnes->residuallocalctx = ctx; 2371aa26658SKarl Rupp 2389566063dSJacob Faibussowitsch PetscCall(DMSNESSetFunction(dm,SNESComputeFunction_DMDA,dmdasnes)); 23922c6f798SBarry Smith if (!sdm->ops->computejacobian) { /* Call us for the Jacobian too, can be overridden by the user. */ 2409566063dSJacob Faibussowitsch PetscCall(DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes)); 2412961e153SJed Brown } 2422961e153SJed Brown PetscFunctionReturn(0); 2432961e153SJed Brown } 2442961e153SJed Brown 2452961e153SJed Brown /*@C 2465505f7afSJunchao Zhang DMDASNESSetFunctionLocalVec - set a local residual evaluation function that operates on a local vector 2475505f7afSJunchao Zhang 2485505f7afSJunchao Zhang Logically Collective 2495505f7afSJunchao Zhang 2505505f7afSJunchao Zhang Input Parameters: 2515505f7afSJunchao Zhang + dm - DM to associate callback with 2525505f7afSJunchao Zhang . imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part 2535505f7afSJunchao Zhang . func - local residual evaluation 2545505f7afSJunchao Zhang - ctx - optional context for local residual evaluation 2555505f7afSJunchao Zhang 2565505f7afSJunchao Zhang Calling sequence: 2575505f7afSJunchao Zhang For PetscErrorCode (*func)(DMDALocalInfo *info,Vec x, Vec f, void *ctx), 2585505f7afSJunchao Zhang + info - DMDALocalInfo defining the subdomain to evaluate the residual on 2595505f7afSJunchao Zhang . x - state vector at which to evaluate residual 2605505f7afSJunchao Zhang . f - residual vector 2615505f7afSJunchao Zhang - ctx - optional context passed above 2625505f7afSJunchao Zhang 2635505f7afSJunchao Zhang Level: beginner 2645505f7afSJunchao Zhang 2655505f7afSJunchao Zhang .seealso: `DMDASNESSetFunctionLocal()`, `DMDASNESSetJacobianLocalVec()`, `DMSNESSetFunction()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()` 2665505f7afSJunchao Zhang @*/ 2675505f7afSJunchao Zhang PetscErrorCode DMDASNESSetFunctionLocalVec(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,Vec,Vec,void*),void *ctx) 2685505f7afSJunchao Zhang { 2695505f7afSJunchao Zhang DMSNES sdm; 2705505f7afSJunchao Zhang DMSNES_DA *dmdasnes; 2715505f7afSJunchao Zhang 2725505f7afSJunchao Zhang PetscFunctionBegin; 2735505f7afSJunchao Zhang PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2745505f7afSJunchao Zhang PetscCall(DMGetDMSNESWrite(dm,&sdm)); 2755505f7afSJunchao Zhang PetscCall(DMDASNESGetContext(dm,sdm,&dmdasnes)); 2765505f7afSJunchao Zhang 2775505f7afSJunchao Zhang dmdasnes->residuallocalimode = imode; 2785505f7afSJunchao Zhang dmdasnes->residuallocalvec = func; 2795505f7afSJunchao Zhang dmdasnes->residuallocalctx = ctx; 2805505f7afSJunchao Zhang 2815505f7afSJunchao Zhang PetscCall(DMSNESSetFunction(dm,SNESComputeFunction_DMDA,dmdasnes)); 2825505f7afSJunchao Zhang if (!sdm->ops->computejacobian) { /* Call us for the Jacobian too, can be overridden by the user. */ 2835505f7afSJunchao Zhang PetscCall(DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes)); 2845505f7afSJunchao Zhang } 2855505f7afSJunchao Zhang PetscFunctionReturn(0); 2865505f7afSJunchao Zhang } 2875505f7afSJunchao Zhang 2885505f7afSJunchao Zhang /*@C 289e79ed307SJed Brown DMDASNESSetJacobianLocal - set a local Jacobian evaluation function 2902961e153SJed Brown 2912961e153SJed Brown Logically Collective 2922961e153SJed Brown 2934165533cSJose E. Roman Input Parameters: 2942961e153SJed Brown + dm - DM to associate callback with 2956b7eb5ceSMatthew G. Knepley . func - local Jacobian evaluation 2966b7eb5ceSMatthew G. Knepley - ctx - optional context for local Jacobian evaluation 2972961e153SJed Brown 2980038884cSBarry Smith Calling sequence: 2990038884cSBarry Smith For PetscErrorCode (*func)(DMDALocalInfo *info,void *x,Mat J,Mat M,void *ctx), 3006b7eb5ceSMatthew G. Knepley + info - DMDALocalInfo defining the subdomain to evaluate the Jacobian at 3010038884cSBarry Smith . x - dimensional pointer to state at which to evaluate Jacobian (e.g. PetscScalar *x or **x or ***x) 3026b7eb5ceSMatthew G. Knepley . J - Mat object for the Jacobian 3036b7eb5ceSMatthew G. Knepley . M - Mat object for the Jacobian preconditioner matrix 3042961e153SJed Brown - ctx - optional context passed above 3052961e153SJed Brown 3062961e153SJed Brown Level: beginner 3072961e153SJed Brown 308db781477SPatrick Sanan .seealso: `DMDASNESSetFunctionLocal()`, `DMSNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()` 3092961e153SJed Brown @*/ 310d1e9a80fSBarry Smith PetscErrorCode DMDASNESSetJacobianLocal(DM dm,PetscErrorCode (*func)(DMDALocalInfo*,void*,Mat,Mat,void*),void *ctx) 3112961e153SJed Brown { 312942e3340SBarry Smith DMSNES sdm; 313942e3340SBarry Smith DMSNES_DA *dmdasnes; 3142961e153SJed Brown 3152961e153SJed Brown PetscFunctionBegin; 3162961e153SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3179566063dSJacob Faibussowitsch PetscCall(DMGetDMSNESWrite(dm,&sdm)); 3189566063dSJacob Faibussowitsch PetscCall(DMDASNESGetContext(dm,sdm,&dmdasnes)); 3191aa26658SKarl Rupp 3202961e153SJed Brown dmdasnes->jacobianlocal = func; 3212961e153SJed Brown dmdasnes->jacobianlocalctx = ctx; 3221aa26658SKarl Rupp 3239566063dSJacob Faibussowitsch PetscCall(DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes)); 3246cab3a1bSJed Brown PetscFunctionReturn(0); 3256cab3a1bSJed Brown } 3262a4ee8f2SPeter Brune 3272a4ee8f2SPeter Brune /*@C 3285505f7afSJunchao Zhang DMDASNESSetJacobianLocalVec - set a local Jacobian evaluation function that operates on a local vector 3295505f7afSJunchao Zhang 3305505f7afSJunchao Zhang Logically Collective 3315505f7afSJunchao Zhang 3325505f7afSJunchao Zhang Input Parameters: 3335505f7afSJunchao Zhang + dm - DM to associate callback with 3345505f7afSJunchao Zhang . func - local Jacobian evaluation 3355505f7afSJunchao Zhang - ctx - optional context for local Jacobian evaluation 3365505f7afSJunchao Zhang 3375505f7afSJunchao Zhang Calling sequence: 3385505f7afSJunchao Zhang For PetscErrorCode (*func)(DMDALocalInfo *info,Vec x,Mat J,Mat M,void *ctx), 3395505f7afSJunchao Zhang + info - DMDALocalInfo defining the subdomain to evaluate the Jacobian at 3405505f7afSJunchao Zhang . x - state vector at which to evaluate Jacobian 3415505f7afSJunchao Zhang . J - Mat object for the Jacobian 3425505f7afSJunchao Zhang . M - Mat object for the Jacobian preconditioner matrix 3435505f7afSJunchao Zhang - ctx - optional context passed above 3445505f7afSJunchao Zhang 3455505f7afSJunchao Zhang Level: beginner 3465505f7afSJunchao Zhang 3475505f7afSJunchao Zhang .seealso: `DMDASNESSetJacobianLocal()`, `DMDASNESSetFunctionLocalVec()`, `DMSNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()` 3485505f7afSJunchao Zhang @*/ 3495505f7afSJunchao Zhang PetscErrorCode DMDASNESSetJacobianLocalVec(DM dm,PetscErrorCode (*func)(DMDALocalInfo*,Vec,Mat,Mat,void*),void *ctx) 3505505f7afSJunchao Zhang { 3515505f7afSJunchao Zhang DMSNES sdm; 3525505f7afSJunchao Zhang DMSNES_DA *dmdasnes; 3535505f7afSJunchao Zhang 3545505f7afSJunchao Zhang PetscFunctionBegin; 3555505f7afSJunchao Zhang PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3565505f7afSJunchao Zhang PetscCall(DMGetDMSNESWrite(dm,&sdm)); 3575505f7afSJunchao Zhang PetscCall(DMDASNESGetContext(dm,sdm,&dmdasnes)); 3585505f7afSJunchao Zhang 3595505f7afSJunchao Zhang dmdasnes->jacobianlocalvec = func; 3605505f7afSJunchao Zhang dmdasnes->jacobianlocalctx = ctx; 3615505f7afSJunchao Zhang 3625505f7afSJunchao Zhang PetscCall(DMSNESSetJacobian(dm,SNESComputeJacobian_DMDA,dmdasnes)); 3635505f7afSJunchao Zhang PetscFunctionReturn(0); 3645505f7afSJunchao Zhang } 3655505f7afSJunchao Zhang 3665505f7afSJunchao Zhang /*@C 3672a4ee8f2SPeter Brune DMDASNESSetObjectiveLocal - set a local residual evaluation function 3682a4ee8f2SPeter Brune 3692a4ee8f2SPeter Brune Logically Collective 3702a4ee8f2SPeter Brune 3714165533cSJose E. Roman Input Parameters: 3722a4ee8f2SPeter Brune + dm - DM to associate callback with 3732a4ee8f2SPeter Brune . func - local objective evaluation 3742a4ee8f2SPeter Brune - ctx - optional context for local residual evaluation 3752a4ee8f2SPeter Brune 3762a4ee8f2SPeter Brune Calling sequence for func: 3772a4ee8f2SPeter Brune + info - DMDALocalInfo defining the subdomain to evaluate the residual on 3782a4ee8f2SPeter Brune . x - dimensional pointer to state at which to evaluate residual 3792a4ee8f2SPeter Brune . ob - eventual objective value 3802a4ee8f2SPeter Brune - ctx - optional context passed above 3812a4ee8f2SPeter Brune 3822a4ee8f2SPeter Brune Level: beginner 3832a4ee8f2SPeter Brune 384db781477SPatrick Sanan .seealso: `DMSNESSetFunction()`, `DMDASNESSetJacobianLocal()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()` 3852a4ee8f2SPeter Brune @*/ 3862a4ee8f2SPeter Brune PetscErrorCode DMDASNESSetObjectiveLocal(DM dm,DMDASNESObjective func,void *ctx) 3872a4ee8f2SPeter Brune { 388942e3340SBarry Smith DMSNES sdm; 389942e3340SBarry Smith DMSNES_DA *dmdasnes; 3902a4ee8f2SPeter Brune 3912a4ee8f2SPeter Brune PetscFunctionBegin; 3922a4ee8f2SPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3939566063dSJacob Faibussowitsch PetscCall(DMGetDMSNESWrite(dm,&sdm)); 3949566063dSJacob Faibussowitsch PetscCall(DMDASNESGetContext(dm,sdm,&dmdasnes)); 3951aa26658SKarl Rupp 3962a4ee8f2SPeter Brune dmdasnes->objectivelocal = func; 3972a4ee8f2SPeter Brune dmdasnes->objectivelocalctx = ctx; 3981aa26658SKarl Rupp 3999566063dSJacob Faibussowitsch PetscCall(DMSNESSetObjective(dm,SNESComputeObjective_DMDA,dmdasnes)); 4002a4ee8f2SPeter Brune PetscFunctionReturn(0); 4012a4ee8f2SPeter Brune } 402dcad5f1cSBarry Smith 4035505f7afSJunchao Zhang /*@C 4045505f7afSJunchao Zhang DMDASNESSetObjectiveLocal - set a local residual evaluation function that operates on a local vector 4055505f7afSJunchao Zhang 4065505f7afSJunchao Zhang Logically Collective 4075505f7afSJunchao Zhang 4085505f7afSJunchao Zhang Input Parameters: 4095505f7afSJunchao Zhang + dm - DM to associate callback with 4105505f7afSJunchao Zhang . func - local objective evaluation 4115505f7afSJunchao Zhang - ctx - optional context for local residual evaluation 4125505f7afSJunchao Zhang 4135505f7afSJunchao Zhang Calling sequence 4145505f7afSJunchao Zhang For PetscErrorCode (*func)(DMDALocalInfo *info,Vec x,PetscReal *ob,void *ctx); 4155505f7afSJunchao Zhang + info - DMDALocalInfo defining the subdomain to evaluate the residual on 4165505f7afSJunchao Zhang . x - state vector at which to evaluate residual 4175505f7afSJunchao Zhang . ob - eventual objective value 4185505f7afSJunchao Zhang - ctx - optional context passed above 4195505f7afSJunchao Zhang 4205505f7afSJunchao Zhang Level: beginner 4215505f7afSJunchao Zhang 4225505f7afSJunchao Zhang .seealso: `DMDASNESSetObjectiveLocal()`, `DMSNESSetFunction()`, `DMDASNESSetJacobianLocalVec()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()` 4235505f7afSJunchao Zhang @*/ 4245505f7afSJunchao Zhang PetscErrorCode DMDASNESSetObjectiveLocalVec(DM dm,DMDASNESObjectiveVec func,void *ctx) 4255505f7afSJunchao Zhang { 4265505f7afSJunchao Zhang DMSNES sdm; 4275505f7afSJunchao Zhang DMSNES_DA *dmdasnes; 4285505f7afSJunchao Zhang 4295505f7afSJunchao Zhang PetscFunctionBegin; 4305505f7afSJunchao Zhang PetscValidHeaderSpecific(dm,DM_CLASSID,1); 4315505f7afSJunchao Zhang PetscCall(DMGetDMSNESWrite(dm,&sdm)); 4325505f7afSJunchao Zhang PetscCall(DMDASNESGetContext(dm,sdm,&dmdasnes)); 4335505f7afSJunchao Zhang 4345505f7afSJunchao Zhang dmdasnes->objectivelocalvec = func; 4355505f7afSJunchao Zhang dmdasnes->objectivelocalctx = ctx; 4365505f7afSJunchao Zhang 4375505f7afSJunchao Zhang PetscCall(DMSNESSetObjective(dm,SNESComputeObjective_DMDA,dmdasnes)); 4385505f7afSJunchao Zhang PetscFunctionReturn(0); 4395505f7afSJunchao Zhang } 4405505f7afSJunchao Zhang 441dcad5f1cSBarry Smith static PetscErrorCode SNESComputePicard_DMDA(SNES snes,Vec X,Vec F,void *ctx) 442dcad5f1cSBarry Smith { 443dcad5f1cSBarry Smith DM dm; 444942e3340SBarry Smith DMSNES_DA *dmdasnes = (DMSNES_DA*)ctx; 445dcad5f1cSBarry Smith DMDALocalInfo info; 446dcad5f1cSBarry Smith Vec Xloc; 447dcad5f1cSBarry Smith void *x,*f; 448dcad5f1cSBarry Smith 449dcad5f1cSBarry Smith PetscFunctionBegin; 450dcad5f1cSBarry Smith PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 451dcad5f1cSBarry Smith PetscValidHeaderSpecific(X,VEC_CLASSID,2); 452dcad5f1cSBarry Smith PetscValidHeaderSpecific(F,VEC_CLASSID,3); 45328b400f6SJacob Faibussowitsch PetscCheck(dmdasnes->rhsplocal,PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context"); 4549566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes,&dm)); 4559566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm,&Xloc)); 4569566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc)); 4579566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc)); 4589566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(dm,&info)); 4599566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(dm,Xloc,&x)); 460dcad5f1cSBarry Smith switch (dmdasnes->residuallocalimode) { 461dcad5f1cSBarry Smith case INSERT_VALUES: { 4629566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(dm,F,&f)); 463*792fecdfSBarry Smith PetscCallBack("SNES Picard DMDA local callback function",(*dmdasnes->rhsplocal)(&info,x,f,dmdasnes->picardlocalctx)); 4649566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dm,F,&f)); 465dcad5f1cSBarry Smith } break; 466dcad5f1cSBarry Smith case ADD_VALUES: { 467dcad5f1cSBarry Smith Vec Floc; 4689566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm,&Floc)); 4699566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Floc)); 4709566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(dm,Floc,&f)); 471*792fecdfSBarry Smith PetscCallBack("SNES Picard DMDA local callback function",(*dmdasnes->rhsplocal)(&info,x,f,dmdasnes->picardlocalctx)); 4729566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dm,Floc,&f)); 4739566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 4749566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F)); 4759566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F)); 4769566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm,&Floc)); 477dcad5f1cSBarry Smith } break; 47898921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdasnes->residuallocalimode); 479dcad5f1cSBarry Smith } 4809566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dm,Xloc,&x)); 4819566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm,&Xloc)); 482dcad5f1cSBarry Smith PetscFunctionReturn(0); 483dcad5f1cSBarry Smith } 484dcad5f1cSBarry Smith 485d1e9a80fSBarry Smith static PetscErrorCode SNESComputePicardJacobian_DMDA(SNES snes,Vec X,Mat A,Mat B,void *ctx) 486dcad5f1cSBarry Smith { 487dcad5f1cSBarry Smith DM dm; 488942e3340SBarry Smith DMSNES_DA *dmdasnes = (DMSNES_DA*)ctx; 489dcad5f1cSBarry Smith DMDALocalInfo info; 490dcad5f1cSBarry Smith Vec Xloc; 491dcad5f1cSBarry Smith void *x; 492dcad5f1cSBarry Smith 493dcad5f1cSBarry Smith PetscFunctionBegin; 49428b400f6SJacob Faibussowitsch PetscCheck(dmdasnes->jacobianplocal,PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"Corrupt context"); 4959566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes,&dm)); 496dcad5f1cSBarry Smith 4979566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm,&Xloc)); 4989566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc)); 4999566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc)); 5009566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(dm,&info)); 5019566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(dm,Xloc,&x)); 502*792fecdfSBarry Smith PetscCallBack("SNES Picard DMDA local callback Jacobian",(*dmdasnes->jacobianplocal)(&info,x,A,B,dmdasnes->picardlocalctx)); 5039566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dm,Xloc,&x)); 5049566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm,&Xloc)); 50594ab13aaSBarry Smith if (A != B) { 5069566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 5079566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 508dcad5f1cSBarry Smith } 509dcad5f1cSBarry Smith PetscFunctionReturn(0); 510dcad5f1cSBarry Smith } 511dcad5f1cSBarry Smith 512dcad5f1cSBarry Smith /*@C 513dcad5f1cSBarry Smith DMDASNESSetPicardLocal - set a local right hand side and matrix evaluation function for Picard iteration 514dcad5f1cSBarry Smith 515dcad5f1cSBarry Smith Logically Collective 516dcad5f1cSBarry Smith 5174165533cSJose E. Roman Input Parameters: 518dcad5f1cSBarry Smith + dm - DM to associate callback with 519dcad5f1cSBarry Smith . imode - INSERT_VALUES if local function computes owned part, ADD_VALUES if it contributes to ghosted part 520dcad5f1cSBarry Smith . func - local residual evaluation 521dcad5f1cSBarry Smith - ctx - optional context for local residual evaluation 522dcad5f1cSBarry Smith 523dcad5f1cSBarry Smith Calling sequence for func: 524dcad5f1cSBarry Smith + info - DMDALocalInfo defining the subdomain to evaluate the residual on 525dcad5f1cSBarry Smith . x - dimensional pointer to state at which to evaluate residual 526dcad5f1cSBarry Smith . f - dimensional pointer to residual, write the residual here 527dcad5f1cSBarry Smith - ctx - optional context passed above 528dcad5f1cSBarry Smith 52995452b02SPatrick Sanan Notes: 53095452b02SPatrick Sanan The user must use 5319566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes,NULL,SNESPicardComputeFunction,&user)); 532dcad5f1cSBarry Smith in their code before calling this routine. 533dcad5f1cSBarry Smith 534dcad5f1cSBarry Smith Level: beginner 535dcad5f1cSBarry Smith 536db781477SPatrick Sanan .seealso: `DMSNESSetFunction()`, `DMDASNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()` 537dcad5f1cSBarry Smith @*/ 538dcad5f1cSBarry Smith PetscErrorCode DMDASNESSetPicardLocal(DM dm,InsertMode imode,PetscErrorCode (*func)(DMDALocalInfo*,void*,void*,void*), 539d1e9a80fSBarry Smith PetscErrorCode (*jac)(DMDALocalInfo*,void*,Mat,Mat,void*),void *ctx) 540dcad5f1cSBarry Smith { 541942e3340SBarry Smith DMSNES sdm; 542942e3340SBarry Smith DMSNES_DA *dmdasnes; 543dcad5f1cSBarry Smith 544dcad5f1cSBarry Smith PetscFunctionBegin; 545dcad5f1cSBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 5469566063dSJacob Faibussowitsch PetscCall(DMGetDMSNESWrite(dm,&sdm)); 5479566063dSJacob Faibussowitsch PetscCall(DMDASNESGetContext(dm,sdm,&dmdasnes)); 5481aa26658SKarl Rupp 549dcad5f1cSBarry Smith dmdasnes->residuallocalimode = imode; 550dcad5f1cSBarry Smith dmdasnes->rhsplocal = func; 551dcad5f1cSBarry Smith dmdasnes->jacobianplocal = jac; 552dcad5f1cSBarry Smith dmdasnes->picardlocalctx = ctx; 5531aa26658SKarl Rupp 5549566063dSJacob Faibussowitsch PetscCall(DMSNESSetPicard(dm,SNESComputePicard_DMDA,SNESComputePicardJacobian_DMDA,dmdasnes)); 5559566063dSJacob Faibussowitsch PetscCall(DMSNESSetMFFunction(dm,SNESComputeFunction_DMDA,dmdasnes)); 556dcad5f1cSBarry Smith PetscFunctionReturn(0); 557dcad5f1cSBarry Smith } 558