1af0996ceSBarry Smith #include <petsc/private/dmimpl.h> 2af0996ceSBarry Smith #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/ 3ff35dfedSBarry Smith 4ff35dfedSBarry Smith typedef struct { 5ff35dfedSBarry Smith PetscErrorCode (*residuallocal)(DM,Vec,Vec,void*); 6d1e9a80fSBarry Smith PetscErrorCode (*jacobianlocal)(DM,Vec,Mat,Mat,void*); 7bdd6f66aSToby Isaac PetscErrorCode (*boundarylocal)(DM,Vec,void*); 8ff35dfedSBarry Smith void *residuallocalctx; 9ff35dfedSBarry Smith void *jacobianlocalctx; 10bdd6f66aSToby Isaac void *boundarylocalctx; 11942e3340SBarry Smith } DMSNES_Local; 12ff35dfedSBarry Smith 13942e3340SBarry Smith static PetscErrorCode DMSNESDestroy_DMLocal(DMSNES sdm) 14ff35dfedSBarry Smith { 15ff35dfedSBarry Smith PetscFunctionBegin; 169566063dSJacob Faibussowitsch PetscCall(PetscFree(sdm->data)); 175ed5e208SMatthew G. Knepley sdm->data = NULL; 18ff35dfedSBarry Smith PetscFunctionReturn(0); 19ff35dfedSBarry Smith } 20ff35dfedSBarry Smith 2122c6f798SBarry Smith static PetscErrorCode DMSNESDuplicate_DMLocal(DMSNES oldsdm,DMSNES sdm) 22ff35dfedSBarry Smith { 23ff35dfedSBarry Smith PetscFunctionBegin; 245ed5e208SMatthew G. Knepley if (sdm->data != oldsdm->data) { 259566063dSJacob Faibussowitsch PetscCall(PetscFree(sdm->data)); 269566063dSJacob Faibussowitsch PetscCall(PetscNewLog(sdm,(DMSNES_Local**)&sdm->data)); 279566063dSJacob Faibussowitsch if (oldsdm->data) PetscCall(PetscMemcpy(sdm->data,oldsdm->data,sizeof(DMSNES_Local))); 28ff35dfedSBarry Smith } 29ff35dfedSBarry Smith PetscFunctionReturn(0); 30ff35dfedSBarry Smith } 31ff35dfedSBarry Smith 32942e3340SBarry Smith static PetscErrorCode DMLocalSNESGetContext(DM dm,DMSNES sdm,DMSNES_Local **dmlocalsnes) 33ff35dfedSBarry Smith { 34ff35dfedSBarry Smith PetscFunctionBegin; 350298fd71SBarry Smith *dmlocalsnes = NULL; 36ff35dfedSBarry Smith if (!sdm->data) { 379566063dSJacob Faibussowitsch PetscCall(PetscNewLog(dm,(DMSNES_Local**)&sdm->data)); 381aa26658SKarl Rupp 3922c6f798SBarry Smith sdm->ops->destroy = DMSNESDestroy_DMLocal; 4022c6f798SBarry Smith sdm->ops->duplicate = DMSNESDuplicate_DMLocal; 41ff35dfedSBarry Smith } 42942e3340SBarry Smith *dmlocalsnes = (DMSNES_Local*)sdm->data; 43ff35dfedSBarry Smith PetscFunctionReturn(0); 44ff35dfedSBarry Smith } 45ff35dfedSBarry Smith 46b2270762SBarry Smith static PetscErrorCode SNESComputeFunction_DMLocal(SNES snes,Vec X,Vec F,void *ctx) 47ff35dfedSBarry Smith { 48942e3340SBarry Smith DMSNES_Local *dmlocalsnes = (DMSNES_Local *) ctx; 497a1e7fa1SMatthew G. Knepley DM dm; 50ff35dfedSBarry Smith Vec Xloc,Floc; 517a1e7fa1SMatthew G. Knepley PetscBool transform; 52ff35dfedSBarry Smith 53ff35dfedSBarry Smith PetscFunctionBegin; 54ff35dfedSBarry Smith PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 55ff35dfedSBarry Smith PetscValidHeaderSpecific(X,VEC_CLASSID,2); 56ff35dfedSBarry Smith PetscValidHeaderSpecific(F,VEC_CLASSID,3); 579566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes,&dm)); 589566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm,&Xloc)); 599566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm,&Floc)); 609566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Xloc)); 619566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Floc)); 627a1e7fa1SMatthew G. Knepley /* Non-conforming routines needs boundary values before G2L */ 639566063dSJacob Faibussowitsch if (dmlocalsnes->boundarylocal) PetscCall((*dmlocalsnes->boundarylocal)(dm,Xloc,dmlocalsnes->boundarylocalctx)); 649566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc)); 659566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc)); 667a1e7fa1SMatthew G. Knepley /* Need to reset boundary values if we transformed */ 679566063dSJacob Faibussowitsch PetscCall(DMHasBasisTransform(dm, &transform)); 689566063dSJacob Faibussowitsch if (transform && dmlocalsnes->boundarylocal) PetscCall((*dmlocalsnes->boundarylocal)(dm,Xloc,dmlocalsnes->boundarylocalctx)); 69ff35dfedSBarry Smith CHKMEMQ; 709566063dSJacob Faibussowitsch PetscCall((*dmlocalsnes->residuallocal)(dm,Xloc,Floc,dmlocalsnes->residuallocalctx)); 71ff35dfedSBarry Smith CHKMEMQ; 729566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 739566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F)); 749566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F)); 759566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm,&Floc)); 769566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm,&Xloc)); 77afa9dec6SMatthew G. Knepley { 78afa9dec6SMatthew G. Knepley char name[PETSC_MAX_PATH_LEN]; 79f2f0fc17SMatthew G. Knepley char oldname[PETSC_MAX_PATH_LEN]; 80f2f0fc17SMatthew G. Knepley const char *tmp; 81afa9dec6SMatthew G. Knepley PetscInt it; 82afa9dec6SMatthew G. Knepley 839566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &it)); 849566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "Solution, Iterate %d", (int) it)); 859566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject) X, &tmp)); 869566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(oldname, tmp, PETSC_MAX_PATH_LEN-1)); 879566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) X, name)); 889566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(X, (PetscObject) snes, "-dmsnes_solution_vec_view")); 899566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) X, oldname)); 909566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "Residual, Iterate %d", (int) it)); 919566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) F, name)); 929566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(F, (PetscObject) snes, "-dmsnes_residual_vec_view")); 93afa9dec6SMatthew G. Knepley } 94ff35dfedSBarry Smith PetscFunctionReturn(0); 95ff35dfedSBarry Smith } 96ff35dfedSBarry Smith 97d1e9a80fSBarry Smith static PetscErrorCode SNESComputeJacobian_DMLocal(SNES snes,Vec X,Mat A,Mat B,void *ctx) 98ff35dfedSBarry Smith { 99942e3340SBarry Smith DMSNES_Local *dmlocalsnes = (DMSNES_Local *) ctx; 1007a1e7fa1SMatthew G. Knepley DM dm; 101ff35dfedSBarry Smith Vec Xloc; 1027a1e7fa1SMatthew G. Knepley PetscBool transform; 103ff35dfedSBarry Smith 104ff35dfedSBarry Smith PetscFunctionBegin; 1059566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes,&dm)); 106ff35dfedSBarry Smith if (dmlocalsnes->jacobianlocal) { 1079566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm,&Xloc)); 1089566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Xloc)); 1097a1e7fa1SMatthew G. Knepley /* Non-conforming routines needs boundary values before G2L */ 1109566063dSJacob Faibussowitsch if (dmlocalsnes->boundarylocal) PetscCall((*dmlocalsnes->boundarylocal)(dm,Xloc,dmlocalsnes->boundarylocalctx)); 1119566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc)); 1129566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc)); 1137a1e7fa1SMatthew G. Knepley /* Need to reset boundary values if we transformed */ 1149566063dSJacob Faibussowitsch PetscCall(DMHasBasisTransform(dm, &transform)); 1159566063dSJacob Faibussowitsch if (transform && dmlocalsnes->boundarylocal) PetscCall((*dmlocalsnes->boundarylocal)(dm,Xloc,dmlocalsnes->boundarylocalctx)); 116ff35dfedSBarry Smith CHKMEMQ; 1179566063dSJacob Faibussowitsch PetscCall((*dmlocalsnes->jacobianlocal)(dm,Xloc,A,B,dmlocalsnes->jacobianlocalctx)); 118ff35dfedSBarry Smith CHKMEMQ; 1199566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm,&Xloc)); 120ff35dfedSBarry Smith } else { 121ff35dfedSBarry Smith MatFDColoring fdcoloring; 1229566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject*)&fdcoloring)); 123ff35dfedSBarry Smith if (!fdcoloring) { 124ff35dfedSBarry Smith ISColoring coloring; 125ff35dfedSBarry Smith 1269566063dSJacob Faibussowitsch PetscCall(DMCreateColoring(dm,dm->coloringtype,&coloring)); 1279566063dSJacob Faibussowitsch PetscCall(MatFDColoringCreate(B,coloring,&fdcoloring)); 1289566063dSJacob Faibussowitsch PetscCall(ISColoringDestroy(&coloring)); 129ff35dfedSBarry Smith switch (dm->coloringtype) { 130ff35dfedSBarry Smith case IS_COLORING_GLOBAL: 1319566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetFunction(fdcoloring,(PetscErrorCode (*)(void))SNESComputeFunction_DMLocal,dmlocalsnes)); 132ff35dfedSBarry Smith break; 13398921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"No support for coloring type '%s'",ISColoringTypes[dm->coloringtype]); 134ff35dfedSBarry Smith } 1359566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fdcoloring,((PetscObject)dm)->prefix)); 1369566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetFromOptions(fdcoloring)); 1379566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetUp(B,coloring,fdcoloring)); 1389566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject)fdcoloring)); 1399566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)fdcoloring)); 140ff35dfedSBarry Smith 141ff35dfedSBarry Smith /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call 142ff35dfedSBarry Smith * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the 143ff35dfedSBarry Smith * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually 144ff35dfedSBarry Smith * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be 145140e18c1SBarry Smith * taken when the PetscObjectList for the Vec inside MatFDColoring is destroyed. 146ff35dfedSBarry Smith */ 1479566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)dm)); 148ff35dfedSBarry Smith } 1499566063dSJacob Faibussowitsch PetscCall(MatFDColoringApply(B,fdcoloring,X,snes)); 150ff35dfedSBarry Smith } 151ff35dfedSBarry Smith /* This will be redundant if the user called both, but it's too common to forget. */ 15294ab13aaSBarry Smith if (A != B) { 1539566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 1549566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 155ff35dfedSBarry Smith } 156ff35dfedSBarry Smith PetscFunctionReturn(0); 157ff35dfedSBarry Smith } 158ff35dfedSBarry Smith 159ff35dfedSBarry Smith /*@C 160ff35dfedSBarry Smith DMSNESSetFunctionLocal - set a local residual evaluation function. This function is called with local vector 161ff35dfedSBarry Smith containing the local vector information PLUS ghost point information. It should compute a result for all local 162ff35dfedSBarry Smith elements and DMSNES will automatically accumulate the overlapping values. 163ff35dfedSBarry Smith 164ff35dfedSBarry Smith Logically Collective 165ff35dfedSBarry Smith 1664165533cSJose E. Roman Input Parameters: 167ff35dfedSBarry Smith + dm - DM to associate callback with 168ff35dfedSBarry Smith . func - local residual evaluation 169ff35dfedSBarry Smith - ctx - optional context for local residual evaluation 170ff35dfedSBarry Smith 171ff35dfedSBarry Smith Level: beginner 172ff35dfedSBarry Smith 173*db781477SPatrick Sanan .seealso: `DMSNESSetFunction()`, `DMDASNESSetJacobianLocal()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()` 174ff35dfedSBarry Smith @*/ 175ff35dfedSBarry Smith PetscErrorCode DMSNESSetFunctionLocal(DM dm,PetscErrorCode (*func)(DM,Vec,Vec,void*),void *ctx) 176ff35dfedSBarry Smith { 177942e3340SBarry Smith DMSNES sdm; 178942e3340SBarry Smith DMSNES_Local *dmlocalsnes; 179ff35dfedSBarry Smith 180ff35dfedSBarry Smith PetscFunctionBegin; 181ff35dfedSBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1829566063dSJacob Faibussowitsch PetscCall(DMGetDMSNESWrite(dm,&sdm)); 1839566063dSJacob Faibussowitsch PetscCall(DMLocalSNESGetContext(dm,sdm,&dmlocalsnes)); 1841aa26658SKarl Rupp 185ff35dfedSBarry Smith dmlocalsnes->residuallocal = func; 186ff35dfedSBarry Smith dmlocalsnes->residuallocalctx = ctx; 1871aa26658SKarl Rupp 1889566063dSJacob Faibussowitsch PetscCall(DMSNESSetFunction(dm,SNESComputeFunction_DMLocal,dmlocalsnes)); 18922c6f798SBarry Smith if (!sdm->ops->computejacobian) { /* Call us for the Jacobian too, can be overridden by the user. */ 1909566063dSJacob Faibussowitsch PetscCall(DMSNESSetJacobian(dm,SNESComputeJacobian_DMLocal,dmlocalsnes)); 191ff35dfedSBarry Smith } 192ff35dfedSBarry Smith PetscFunctionReturn(0); 193ff35dfedSBarry Smith } 194ff35dfedSBarry Smith 195bdd6f66aSToby Isaac /*@C 196bdd6f66aSToby Isaac DMSNESSetBoundaryLocal - set a local boundary value function. This function is called with local vector 197bdd6f66aSToby Isaac containing the local vector information PLUS ghost point information. It should insert values into the local 198bdd6f66aSToby Isaac vector that do not come from the global vector, such as essential boundary condition data. 199bdd6f66aSToby Isaac 200bdd6f66aSToby Isaac Logically Collective 201bdd6f66aSToby Isaac 2024165533cSJose E. Roman Input Parameters: 203bdd6f66aSToby Isaac + dm - DM to associate callback with 204bdd6f66aSToby Isaac . func - local boundary value evaluation 205bdd6f66aSToby Isaac - ctx - optional context for local boundary value evaluation 206bdd6f66aSToby Isaac 207bdd6f66aSToby Isaac Level: intermediate 208bdd6f66aSToby Isaac 209*db781477SPatrick Sanan .seealso: `DMSNESSetFunctionLocal()`, `DMDASNESSetJacobianLocal()` 210bdd6f66aSToby Isaac @*/ 211bdd6f66aSToby Isaac PetscErrorCode DMSNESSetBoundaryLocal(DM dm,PetscErrorCode (*func)(DM,Vec,void*),void *ctx) 212bdd6f66aSToby Isaac { 213bdd6f66aSToby Isaac DMSNES sdm; 214bdd6f66aSToby Isaac DMSNES_Local *dmlocalsnes; 215bdd6f66aSToby Isaac 216bdd6f66aSToby Isaac PetscFunctionBegin; 217bdd6f66aSToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2189566063dSJacob Faibussowitsch PetscCall(DMGetDMSNESWrite(dm,&sdm)); 2199566063dSJacob Faibussowitsch PetscCall(DMLocalSNESGetContext(dm,sdm,&dmlocalsnes)); 220bdd6f66aSToby Isaac 221bdd6f66aSToby Isaac dmlocalsnes->boundarylocal = func; 222bdd6f66aSToby Isaac dmlocalsnes->boundarylocalctx = ctx; 223bdd6f66aSToby Isaac 224bdd6f66aSToby Isaac PetscFunctionReturn(0); 225bdd6f66aSToby Isaac } 226bdd6f66aSToby Isaac 227ff35dfedSBarry Smith /*@C 228ff35dfedSBarry Smith DMSNESSetJacobianLocal - set a local Jacobian evaluation function 229ff35dfedSBarry Smith 230ff35dfedSBarry Smith Logically Collective 231ff35dfedSBarry Smith 2324165533cSJose E. Roman Input Parameters: 233ff35dfedSBarry Smith + dm - DM to associate callback with 234ff35dfedSBarry Smith . func - local Jacobian evaluation 235ff35dfedSBarry Smith - ctx - optional context for local Jacobian evaluation 236ff35dfedSBarry Smith 237ff35dfedSBarry Smith Level: beginner 238ff35dfedSBarry Smith 239*db781477SPatrick Sanan .seealso: `DMSNESSetJacobian()`, `DMDASNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()` 240ff35dfedSBarry Smith @*/ 241d1e9a80fSBarry Smith PetscErrorCode DMSNESSetJacobianLocal(DM dm,PetscErrorCode (*func)(DM,Vec,Mat,Mat,void*),void *ctx) 242ff35dfedSBarry Smith { 243942e3340SBarry Smith DMSNES sdm; 244942e3340SBarry Smith DMSNES_Local *dmlocalsnes; 245ff35dfedSBarry Smith 246ff35dfedSBarry Smith PetscFunctionBegin; 247ff35dfedSBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2489566063dSJacob Faibussowitsch PetscCall(DMGetDMSNESWrite(dm,&sdm)); 2499566063dSJacob Faibussowitsch PetscCall(DMLocalSNESGetContext(dm,sdm,&dmlocalsnes)); 2501aa26658SKarl Rupp 251ff35dfedSBarry Smith dmlocalsnes->jacobianlocal = func; 252ff35dfedSBarry Smith dmlocalsnes->jacobianlocalctx = ctx; 2531aa26658SKarl Rupp 2549566063dSJacob Faibussowitsch PetscCall(DMSNESSetJacobian(dm,SNESComputeJacobian_DMLocal,dmlocalsnes)); 25528d58a37SPierre Jolivet PetscFunctionReturn(0); 25628d58a37SPierre Jolivet } 25728d58a37SPierre Jolivet 25828d58a37SPierre Jolivet /*@C 25928d58a37SPierre Jolivet DMSNESGetFunctionLocal - get the local residual evaluation function information set with DMSNESSetFunctionLocal. 26028d58a37SPierre Jolivet 26128d58a37SPierre Jolivet Not Collective 26228d58a37SPierre Jolivet 2634165533cSJose E. Roman Input Parameter: 264f0fc11ceSJed Brown . dm - DM with the associated callback 26528d58a37SPierre Jolivet 2664165533cSJose E. Roman Output Parameters: 26728d58a37SPierre Jolivet + func - local residual evaluation 26828d58a37SPierre Jolivet - ctx - context for local residual evaluation 26928d58a37SPierre Jolivet 27028d58a37SPierre Jolivet Level: beginner 27128d58a37SPierre Jolivet 272*db781477SPatrick Sanan .seealso: `DMSNESSetFunction()`, `DMSNESSetFunctionLocal()`, `DMDASNESSetJacobianLocal()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()` 27328d58a37SPierre Jolivet @*/ 27428d58a37SPierre Jolivet PetscErrorCode DMSNESGetFunctionLocal(DM dm,PetscErrorCode (**func)(DM,Vec,Vec,void*),void **ctx) 27528d58a37SPierre Jolivet { 27628d58a37SPierre Jolivet DMSNES sdm; 27728d58a37SPierre Jolivet DMSNES_Local *dmlocalsnes; 27828d58a37SPierre Jolivet 27928d58a37SPierre Jolivet PetscFunctionBegin; 28028d58a37SPierre Jolivet PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2819566063dSJacob Faibussowitsch PetscCall(DMGetDMSNES(dm,&sdm)); 2829566063dSJacob Faibussowitsch PetscCall(DMLocalSNESGetContext(dm,sdm,&dmlocalsnes)); 28328d58a37SPierre Jolivet if (func) *func = dmlocalsnes->residuallocal; 28428d58a37SPierre Jolivet if (ctx) *ctx = dmlocalsnes->residuallocalctx; 28528d58a37SPierre Jolivet PetscFunctionReturn(0); 28628d58a37SPierre Jolivet } 28728d58a37SPierre Jolivet 28828d58a37SPierre Jolivet /*@C 28928d58a37SPierre Jolivet DMSNESGetBoundaryLocal - get the local boundary value function set with DMSNESSetBoundaryLocal. 29028d58a37SPierre Jolivet 29128d58a37SPierre Jolivet Not Collective 29228d58a37SPierre Jolivet 2934165533cSJose E. Roman Input Parameter: 294f0fc11ceSJed Brown . dm - DM with the associated callback 29528d58a37SPierre Jolivet 2964165533cSJose E. Roman Output Parameters: 29728d58a37SPierre Jolivet + func - local boundary value evaluation 29828d58a37SPierre Jolivet - ctx - context for local boundary value evaluation 29928d58a37SPierre Jolivet 30028d58a37SPierre Jolivet Level: intermediate 30128d58a37SPierre Jolivet 302*db781477SPatrick Sanan .seealso: `DMSNESSetFunctionLocal()`, `DMSNESSetBoundaryLocal()`, `DMDASNESSetJacobianLocal()` 30328d58a37SPierre Jolivet @*/ 30428d58a37SPierre Jolivet PetscErrorCode DMSNESGetBoundaryLocal(DM dm,PetscErrorCode (**func)(DM,Vec,void*),void **ctx) 30528d58a37SPierre Jolivet { 30628d58a37SPierre Jolivet DMSNES sdm; 30728d58a37SPierre Jolivet DMSNES_Local *dmlocalsnes; 30828d58a37SPierre Jolivet 30928d58a37SPierre Jolivet PetscFunctionBegin; 31028d58a37SPierre Jolivet PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3119566063dSJacob Faibussowitsch PetscCall(DMGetDMSNES(dm,&sdm)); 3129566063dSJacob Faibussowitsch PetscCall(DMLocalSNESGetContext(dm,sdm,&dmlocalsnes)); 31328d58a37SPierre Jolivet if (func) *func = dmlocalsnes->boundarylocal; 31428d58a37SPierre Jolivet if (ctx) *ctx = dmlocalsnes->boundarylocalctx; 31528d58a37SPierre Jolivet PetscFunctionReturn(0); 31628d58a37SPierre Jolivet } 31728d58a37SPierre Jolivet 31828d58a37SPierre Jolivet /*@C 31928d58a37SPierre Jolivet DMSNESGetJacobianLocal - the local Jacobian evaluation function set with DMSNESSetJacobianLocal. 32028d58a37SPierre Jolivet 32128d58a37SPierre Jolivet Logically Collective 32228d58a37SPierre Jolivet 3234165533cSJose E. Roman Input Parameter: 324f0fc11ceSJed Brown . dm - DM with the associated callback 32528d58a37SPierre Jolivet 3264165533cSJose E. Roman Output Parameters: 32728d58a37SPierre Jolivet + func - local Jacobian evaluation 32828d58a37SPierre Jolivet - ctx - context for local Jacobian evaluation 32928d58a37SPierre Jolivet 33028d58a37SPierre Jolivet Level: beginner 33128d58a37SPierre Jolivet 332*db781477SPatrick Sanan .seealso: `DMSNESSetJacobianLocal()`, `DMDASNESSetJacobian()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()` 33328d58a37SPierre Jolivet @*/ 33428d58a37SPierre Jolivet PetscErrorCode DMSNESGetJacobianLocal(DM dm,PetscErrorCode (**func)(DM,Vec,Mat,Mat,void*),void **ctx) 33528d58a37SPierre Jolivet { 33628d58a37SPierre Jolivet DMSNES sdm; 33728d58a37SPierre Jolivet DMSNES_Local *dmlocalsnes; 33828d58a37SPierre Jolivet 33928d58a37SPierre Jolivet PetscFunctionBegin; 34028d58a37SPierre Jolivet PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3419566063dSJacob Faibussowitsch PetscCall(DMGetDMSNES(dm,&sdm)); 3429566063dSJacob Faibussowitsch PetscCall(DMLocalSNESGetContext(dm,sdm,&dmlocalsnes)); 34328d58a37SPierre Jolivet if (func) *func = dmlocalsnes->jacobianlocal; 34428d58a37SPierre Jolivet if (ctx) *ctx = dmlocalsnes->jacobianlocalctx; 345ff35dfedSBarry Smith PetscFunctionReturn(0); 346ff35dfedSBarry Smith } 347