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 PetscErrorCode ierr; 16ff35dfedSBarry Smith 17ff35dfedSBarry Smith PetscFunctionBegin; 18ff35dfedSBarry Smith ierr = PetscFree(sdm->data);CHKERRQ(ierr); 195ed5e208SMatthew G. Knepley sdm->data = NULL; 20ff35dfedSBarry Smith PetscFunctionReturn(0); 21ff35dfedSBarry Smith } 22ff35dfedSBarry Smith 2322c6f798SBarry Smith static PetscErrorCode DMSNESDuplicate_DMLocal(DMSNES oldsdm,DMSNES sdm) 24ff35dfedSBarry Smith { 25ff35dfedSBarry Smith PetscErrorCode ierr; 26ff35dfedSBarry Smith 27ff35dfedSBarry Smith PetscFunctionBegin; 285ed5e208SMatthew G. Knepley if (sdm->data != oldsdm->data) { 295ed5e208SMatthew G. Knepley ierr = PetscFree(sdm->data);CHKERRQ(ierr); 30b00a9115SJed Brown ierr = PetscNewLog(sdm,(DMSNES_Local**)&sdm->data);CHKERRQ(ierr); 315ed5e208SMatthew G. Knepley if (oldsdm->data) {ierr = PetscMemcpy(sdm->data,oldsdm->data,sizeof(DMSNES_Local));CHKERRQ(ierr);} 32ff35dfedSBarry Smith } 33ff35dfedSBarry Smith PetscFunctionReturn(0); 34ff35dfedSBarry Smith } 35ff35dfedSBarry Smith 36942e3340SBarry Smith static PetscErrorCode DMLocalSNESGetContext(DM dm,DMSNES sdm,DMSNES_Local **dmlocalsnes) 37ff35dfedSBarry Smith { 38ff35dfedSBarry Smith PetscErrorCode ierr; 39ff35dfedSBarry Smith 40ff35dfedSBarry Smith PetscFunctionBegin; 410298fd71SBarry Smith *dmlocalsnes = NULL; 42ff35dfedSBarry Smith if (!sdm->data) { 43b00a9115SJed Brown ierr = PetscNewLog(dm,(DMSNES_Local**)&sdm->data);CHKERRQ(ierr); 441aa26658SKarl Rupp 4522c6f798SBarry Smith sdm->ops->destroy = DMSNESDestroy_DMLocal; 4622c6f798SBarry Smith sdm->ops->duplicate = DMSNESDuplicate_DMLocal; 47ff35dfedSBarry Smith } 48942e3340SBarry Smith *dmlocalsnes = (DMSNES_Local*)sdm->data; 49ff35dfedSBarry Smith PetscFunctionReturn(0); 50ff35dfedSBarry Smith } 51ff35dfedSBarry Smith 52b2270762SBarry Smith static PetscErrorCode SNESComputeFunction_DMLocal(SNES snes,Vec X,Vec F,void *ctx) 53ff35dfedSBarry Smith { 54942e3340SBarry Smith DMSNES_Local *dmlocalsnes = (DMSNES_Local *) ctx; 557a1e7fa1SMatthew G. Knepley DM dm; 56ff35dfedSBarry Smith Vec Xloc,Floc; 577a1e7fa1SMatthew G. Knepley PetscBool transform; 587a1e7fa1SMatthew G. Knepley PetscErrorCode ierr; 59ff35dfedSBarry Smith 60ff35dfedSBarry Smith PetscFunctionBegin; 61ff35dfedSBarry Smith PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 62ff35dfedSBarry Smith PetscValidHeaderSpecific(X,VEC_CLASSID,2); 63ff35dfedSBarry Smith PetscValidHeaderSpecific(F,VEC_CLASSID,3); 64ff35dfedSBarry Smith ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 65ff35dfedSBarry Smith ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 66ff35dfedSBarry Smith ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr); 67ff35dfedSBarry Smith ierr = VecZeroEntries(Xloc);CHKERRQ(ierr); 68ff35dfedSBarry Smith ierr = VecZeroEntries(Floc);CHKERRQ(ierr); 697a1e7fa1SMatthew G. Knepley /* Non-conforming routines needs boundary values before G2L */ 707a1e7fa1SMatthew G. Knepley if (dmlocalsnes->boundarylocal) {ierr = (*dmlocalsnes->boundarylocal)(dm,Xloc,dmlocalsnes->boundarylocalctx);CHKERRQ(ierr);} 71ff35dfedSBarry Smith ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 72ff35dfedSBarry Smith ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 737a1e7fa1SMatthew G. Knepley /* Need to reset boundary values if we transformed */ 747a1e7fa1SMatthew G. Knepley ierr = DMHasBasisTransform(dm, &transform);CHKERRQ(ierr); 757a1e7fa1SMatthew G. Knepley if (transform && dmlocalsnes->boundarylocal) {ierr = (*dmlocalsnes->boundarylocal)(dm,Xloc,dmlocalsnes->boundarylocalctx);CHKERRQ(ierr);} 76ff35dfedSBarry Smith CHKMEMQ; 77ff35dfedSBarry Smith ierr = (*dmlocalsnes->residuallocal)(dm,Xloc,Floc,dmlocalsnes->residuallocalctx);CHKERRQ(ierr); 78ff35dfedSBarry Smith CHKMEMQ; 79ff35dfedSBarry Smith ierr = VecZeroEntries(F);CHKERRQ(ierr); 80ff35dfedSBarry Smith ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 81ff35dfedSBarry Smith ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 82ff35dfedSBarry Smith ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr); 83ff35dfedSBarry Smith ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 84afa9dec6SMatthew G. Knepley { 85afa9dec6SMatthew G. Knepley char name[PETSC_MAX_PATH_LEN]; 86f2f0fc17SMatthew G. Knepley char oldname[PETSC_MAX_PATH_LEN]; 87f2f0fc17SMatthew G. Knepley const char *tmp; 88afa9dec6SMatthew G. Knepley PetscInt it; 89afa9dec6SMatthew G. Knepley 90afa9dec6SMatthew G. Knepley ierr = SNESGetIterationNumber(snes, &it);CHKERRQ(ierr); 91afa9dec6SMatthew G. Knepley ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "Solution, Iterate %d", (int) it);CHKERRQ(ierr); 92f2f0fc17SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) X, &tmp);CHKERRQ(ierr); 93f2f0fc17SMatthew G. Knepley ierr = PetscStrncpy(oldname, tmp, PETSC_MAX_PATH_LEN-1);CHKERRQ(ierr); 94afa9dec6SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) X, name);CHKERRQ(ierr); 9585997439SMatthew G. Knepley ierr = VecViewFromOptions(X, (PetscObject) snes, "-dmsnes_solution_vec_view");CHKERRQ(ierr); 96f2f0fc17SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) X, oldname);CHKERRQ(ierr); 97afa9dec6SMatthew G. Knepley ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "Residual, Iterate %d", (int) it);CHKERRQ(ierr); 98afa9dec6SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) F, name);CHKERRQ(ierr); 9985997439SMatthew G. Knepley ierr = VecViewFromOptions(F, (PetscObject) snes, "-dmsnes_residual_vec_view");CHKERRQ(ierr); 100afa9dec6SMatthew G. Knepley } 101ff35dfedSBarry Smith PetscFunctionReturn(0); 102ff35dfedSBarry Smith } 103ff35dfedSBarry Smith 104d1e9a80fSBarry Smith static PetscErrorCode SNESComputeJacobian_DMLocal(SNES snes,Vec X,Mat A,Mat B,void *ctx) 105ff35dfedSBarry Smith { 106942e3340SBarry Smith DMSNES_Local *dmlocalsnes = (DMSNES_Local *) ctx; 1077a1e7fa1SMatthew G. Knepley DM dm; 108ff35dfedSBarry Smith Vec Xloc; 1097a1e7fa1SMatthew G. Knepley PetscBool transform; 1107a1e7fa1SMatthew G. Knepley PetscErrorCode ierr; 111ff35dfedSBarry Smith 112ff35dfedSBarry Smith PetscFunctionBegin; 113ff35dfedSBarry Smith ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 114ff35dfedSBarry Smith if (dmlocalsnes->jacobianlocal) { 115ff35dfedSBarry Smith ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 116ff35dfedSBarry Smith ierr = VecZeroEntries(Xloc);CHKERRQ(ierr); 1177a1e7fa1SMatthew G. Knepley /* Non-conforming routines needs boundary values before G2L */ 1187a1e7fa1SMatthew G. Knepley if (dmlocalsnes->boundarylocal) {ierr = (*dmlocalsnes->boundarylocal)(dm,Xloc,dmlocalsnes->boundarylocalctx);CHKERRQ(ierr);} 119ff35dfedSBarry Smith ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 120ff35dfedSBarry Smith ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 1217a1e7fa1SMatthew G. Knepley /* Need to reset boundary values if we transformed */ 1227a1e7fa1SMatthew G. Knepley ierr = DMHasBasisTransform(dm, &transform);CHKERRQ(ierr); 1237a1e7fa1SMatthew G. Knepley if (transform && dmlocalsnes->boundarylocal) {ierr = (*dmlocalsnes->boundarylocal)(dm,Xloc,dmlocalsnes->boundarylocalctx);CHKERRQ(ierr);} 124ff35dfedSBarry Smith CHKMEMQ; 125d1e9a80fSBarry Smith ierr = (*dmlocalsnes->jacobianlocal)(dm,Xloc,A,B,dmlocalsnes->jacobianlocalctx);CHKERRQ(ierr); 126ff35dfedSBarry Smith CHKMEMQ; 127ff35dfedSBarry Smith ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 128ff35dfedSBarry Smith } else { 129ff35dfedSBarry Smith MatFDColoring fdcoloring; 130ff35dfedSBarry Smith ierr = PetscObjectQuery((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject*)&fdcoloring);CHKERRQ(ierr); 131ff35dfedSBarry Smith if (!fdcoloring) { 132ff35dfedSBarry Smith ISColoring coloring; 133ff35dfedSBarry Smith 134b412c318SBarry Smith ierr = DMCreateColoring(dm,dm->coloringtype,&coloring);CHKERRQ(ierr); 13594ab13aaSBarry Smith ierr = MatFDColoringCreate(B,coloring,&fdcoloring);CHKERRQ(ierr); 136ff35dfedSBarry Smith ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr); 137ff35dfedSBarry Smith switch (dm->coloringtype) { 138ff35dfedSBarry Smith case IS_COLORING_GLOBAL: 139ff35dfedSBarry Smith ierr = MatFDColoringSetFunction(fdcoloring,(PetscErrorCode (*)(void))SNESComputeFunction_DMLocal,dmlocalsnes);CHKERRQ(ierr); 140ff35dfedSBarry Smith break; 141ce94432eSBarry Smith default: SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"No support for coloring type '%s'",ISColoringTypes[dm->coloringtype]); 142ff35dfedSBarry Smith } 143ff35dfedSBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)fdcoloring,((PetscObject)dm)->prefix);CHKERRQ(ierr); 144ff35dfedSBarry Smith ierr = MatFDColoringSetFromOptions(fdcoloring);CHKERRQ(ierr); 14594ab13aaSBarry Smith ierr = MatFDColoringSetUp(B,coloring,fdcoloring);CHKERRQ(ierr); 146ff35dfedSBarry Smith ierr = PetscObjectCompose((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject)fdcoloring);CHKERRQ(ierr); 147ff35dfedSBarry Smith ierr = PetscObjectDereference((PetscObject)fdcoloring);CHKERRQ(ierr); 148ff35dfedSBarry Smith 149ff35dfedSBarry Smith /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call 150ff35dfedSBarry Smith * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the 151ff35dfedSBarry Smith * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually 152ff35dfedSBarry Smith * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be 153140e18c1SBarry Smith * taken when the PetscObjectList for the Vec inside MatFDColoring is destroyed. 154ff35dfedSBarry Smith */ 155ff35dfedSBarry Smith ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr); 156ff35dfedSBarry Smith } 157d1e9a80fSBarry Smith ierr = MatFDColoringApply(B,fdcoloring,X,snes);CHKERRQ(ierr); 158ff35dfedSBarry Smith } 159ff35dfedSBarry Smith /* This will be redundant if the user called both, but it's too common to forget. */ 16094ab13aaSBarry Smith if (A != B) { 16194ab13aaSBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 16294ab13aaSBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 163ff35dfedSBarry Smith } 164ff35dfedSBarry Smith PetscFunctionReturn(0); 165ff35dfedSBarry Smith } 166ff35dfedSBarry Smith 167ff35dfedSBarry Smith /*@C 168ff35dfedSBarry Smith DMSNESSetFunctionLocal - set a local residual evaluation function. This function is called with local vector 169ff35dfedSBarry Smith containing the local vector information PLUS ghost point information. It should compute a result for all local 170ff35dfedSBarry Smith elements and DMSNES will automatically accumulate the overlapping values. 171ff35dfedSBarry Smith 172ff35dfedSBarry Smith Logically Collective 173ff35dfedSBarry Smith 174*4165533cSJose E. Roman Input Parameters: 175ff35dfedSBarry Smith + dm - DM to associate callback with 176ff35dfedSBarry Smith . func - local residual evaluation 177ff35dfedSBarry Smith - ctx - optional context for local residual evaluation 178ff35dfedSBarry Smith 179ff35dfedSBarry Smith Level: beginner 180ff35dfedSBarry Smith 181ff35dfedSBarry Smith .seealso: DMSNESSetFunction(), DMDASNESSetJacobianLocal(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 182ff35dfedSBarry Smith @*/ 183ff35dfedSBarry Smith PetscErrorCode DMSNESSetFunctionLocal(DM dm,PetscErrorCode (*func)(DM,Vec,Vec,void*),void *ctx) 184ff35dfedSBarry Smith { 185ff35dfedSBarry Smith PetscErrorCode ierr; 186942e3340SBarry Smith DMSNES sdm; 187942e3340SBarry Smith DMSNES_Local *dmlocalsnes; 188ff35dfedSBarry Smith 189ff35dfedSBarry Smith PetscFunctionBegin; 190ff35dfedSBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 191942e3340SBarry Smith ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 192ff35dfedSBarry Smith ierr = DMLocalSNESGetContext(dm,sdm,&dmlocalsnes);CHKERRQ(ierr); 1931aa26658SKarl Rupp 194ff35dfedSBarry Smith dmlocalsnes->residuallocal = func; 195ff35dfedSBarry Smith dmlocalsnes->residuallocalctx = ctx; 1961aa26658SKarl Rupp 197ff35dfedSBarry Smith ierr = DMSNESSetFunction(dm,SNESComputeFunction_DMLocal,dmlocalsnes);CHKERRQ(ierr); 19822c6f798SBarry Smith if (!sdm->ops->computejacobian) { /* Call us for the Jacobian too, can be overridden by the user. */ 199ff35dfedSBarry Smith ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMLocal,dmlocalsnes);CHKERRQ(ierr); 200ff35dfedSBarry Smith } 201ff35dfedSBarry Smith PetscFunctionReturn(0); 202ff35dfedSBarry Smith } 203ff35dfedSBarry Smith 204bdd6f66aSToby Isaac /*@C 205bdd6f66aSToby Isaac DMSNESSetBoundaryLocal - set a local boundary value function. This function is called with local vector 206bdd6f66aSToby Isaac containing the local vector information PLUS ghost point information. It should insert values into the local 207bdd6f66aSToby Isaac vector that do not come from the global vector, such as essential boundary condition data. 208bdd6f66aSToby Isaac 209bdd6f66aSToby Isaac Logically Collective 210bdd6f66aSToby Isaac 211*4165533cSJose E. Roman Input Parameters: 212bdd6f66aSToby Isaac + dm - DM to associate callback with 213bdd6f66aSToby Isaac . func - local boundary value evaluation 214bdd6f66aSToby Isaac - ctx - optional context for local boundary value evaluation 215bdd6f66aSToby Isaac 216bdd6f66aSToby Isaac Level: intermediate 217bdd6f66aSToby Isaac 218bdd6f66aSToby Isaac .seealso: DMSNESSetFunctionLocal(), DMDASNESSetJacobianLocal() 219bdd6f66aSToby Isaac @*/ 220bdd6f66aSToby Isaac PetscErrorCode DMSNESSetBoundaryLocal(DM dm,PetscErrorCode (*func)(DM,Vec,void*),void *ctx) 221bdd6f66aSToby Isaac { 222bdd6f66aSToby Isaac PetscErrorCode ierr; 223bdd6f66aSToby Isaac DMSNES sdm; 224bdd6f66aSToby Isaac DMSNES_Local *dmlocalsnes; 225bdd6f66aSToby Isaac 226bdd6f66aSToby Isaac PetscFunctionBegin; 227bdd6f66aSToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 228bdd6f66aSToby Isaac ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 229bdd6f66aSToby Isaac ierr = DMLocalSNESGetContext(dm,sdm,&dmlocalsnes);CHKERRQ(ierr); 230bdd6f66aSToby Isaac 231bdd6f66aSToby Isaac dmlocalsnes->boundarylocal = func; 232bdd6f66aSToby Isaac dmlocalsnes->boundarylocalctx = ctx; 233bdd6f66aSToby Isaac 234bdd6f66aSToby Isaac PetscFunctionReturn(0); 235bdd6f66aSToby Isaac } 236bdd6f66aSToby Isaac 237ff35dfedSBarry Smith /*@C 238ff35dfedSBarry Smith DMSNESSetJacobianLocal - set a local Jacobian evaluation function 239ff35dfedSBarry Smith 240ff35dfedSBarry Smith Logically Collective 241ff35dfedSBarry Smith 242*4165533cSJose E. Roman Input Parameters: 243ff35dfedSBarry Smith + dm - DM to associate callback with 244ff35dfedSBarry Smith . func - local Jacobian evaluation 245ff35dfedSBarry Smith - ctx - optional context for local Jacobian evaluation 246ff35dfedSBarry Smith 247ff35dfedSBarry Smith Level: beginner 248ff35dfedSBarry Smith 249ff35dfedSBarry Smith .seealso: DMSNESSetJacobian(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 250ff35dfedSBarry Smith @*/ 251d1e9a80fSBarry Smith PetscErrorCode DMSNESSetJacobianLocal(DM dm,PetscErrorCode (*func)(DM,Vec,Mat,Mat,void*),void *ctx) 252ff35dfedSBarry Smith { 253ff35dfedSBarry Smith PetscErrorCode ierr; 254942e3340SBarry Smith DMSNES sdm; 255942e3340SBarry Smith DMSNES_Local *dmlocalsnes; 256ff35dfedSBarry Smith 257ff35dfedSBarry Smith PetscFunctionBegin; 258ff35dfedSBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 259942e3340SBarry Smith ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 260ff35dfedSBarry Smith ierr = DMLocalSNESGetContext(dm,sdm,&dmlocalsnes);CHKERRQ(ierr); 2611aa26658SKarl Rupp 262ff35dfedSBarry Smith dmlocalsnes->jacobianlocal = func; 263ff35dfedSBarry Smith dmlocalsnes->jacobianlocalctx = ctx; 2641aa26658SKarl Rupp 265ff35dfedSBarry Smith ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMLocal,dmlocalsnes);CHKERRQ(ierr); 26628d58a37SPierre Jolivet PetscFunctionReturn(0); 26728d58a37SPierre Jolivet } 26828d58a37SPierre Jolivet 26928d58a37SPierre Jolivet /*@C 27028d58a37SPierre Jolivet DMSNESGetFunctionLocal - get the local residual evaluation function information set with DMSNESSetFunctionLocal. 27128d58a37SPierre Jolivet 27228d58a37SPierre Jolivet Not Collective 27328d58a37SPierre Jolivet 274*4165533cSJose E. Roman Input Parameter: 275f0fc11ceSJed Brown . dm - DM with the associated callback 27628d58a37SPierre Jolivet 277*4165533cSJose E. Roman Output Parameters: 27828d58a37SPierre Jolivet + func - local residual evaluation 27928d58a37SPierre Jolivet - ctx - context for local residual evaluation 28028d58a37SPierre Jolivet 28128d58a37SPierre Jolivet Level: beginner 28228d58a37SPierre Jolivet 28328d58a37SPierre Jolivet .seealso: DMSNESSetFunction(), DMSNESSetFunctionLocal(), DMDASNESSetJacobianLocal(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 28428d58a37SPierre Jolivet @*/ 28528d58a37SPierre Jolivet PetscErrorCode DMSNESGetFunctionLocal(DM dm,PetscErrorCode (**func)(DM,Vec,Vec,void*),void **ctx) 28628d58a37SPierre Jolivet { 28728d58a37SPierre Jolivet PetscErrorCode ierr; 28828d58a37SPierre Jolivet DMSNES sdm; 28928d58a37SPierre Jolivet DMSNES_Local *dmlocalsnes; 29028d58a37SPierre Jolivet 29128d58a37SPierre Jolivet PetscFunctionBegin; 29228d58a37SPierre Jolivet PetscValidHeaderSpecific(dm,DM_CLASSID,1); 29328d58a37SPierre Jolivet ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr); 29428d58a37SPierre Jolivet ierr = DMLocalSNESGetContext(dm,sdm,&dmlocalsnes);CHKERRQ(ierr); 29528d58a37SPierre Jolivet if (func) *func = dmlocalsnes->residuallocal; 29628d58a37SPierre Jolivet if (ctx) *ctx = dmlocalsnes->residuallocalctx; 29728d58a37SPierre Jolivet PetscFunctionReturn(0); 29828d58a37SPierre Jolivet } 29928d58a37SPierre Jolivet 30028d58a37SPierre Jolivet /*@C 30128d58a37SPierre Jolivet DMSNESGetBoundaryLocal - get the local boundary value function set with DMSNESSetBoundaryLocal. 30228d58a37SPierre Jolivet 30328d58a37SPierre Jolivet Not Collective 30428d58a37SPierre Jolivet 305*4165533cSJose E. Roman Input Parameter: 306f0fc11ceSJed Brown . dm - DM with the associated callback 30728d58a37SPierre Jolivet 308*4165533cSJose E. Roman Output Parameters: 30928d58a37SPierre Jolivet + func - local boundary value evaluation 31028d58a37SPierre Jolivet - ctx - context for local boundary value evaluation 31128d58a37SPierre Jolivet 31228d58a37SPierre Jolivet Level: intermediate 31328d58a37SPierre Jolivet 31428d58a37SPierre Jolivet .seealso: DMSNESSetFunctionLocal(), DMSNESSetBoundaryLocal(), DMDASNESSetJacobianLocal() 31528d58a37SPierre Jolivet @*/ 31628d58a37SPierre Jolivet PetscErrorCode DMSNESGetBoundaryLocal(DM dm,PetscErrorCode (**func)(DM,Vec,void*),void **ctx) 31728d58a37SPierre Jolivet { 31828d58a37SPierre Jolivet PetscErrorCode ierr; 31928d58a37SPierre Jolivet DMSNES sdm; 32028d58a37SPierre Jolivet DMSNES_Local *dmlocalsnes; 32128d58a37SPierre Jolivet 32228d58a37SPierre Jolivet PetscFunctionBegin; 32328d58a37SPierre Jolivet PetscValidHeaderSpecific(dm,DM_CLASSID,1); 32428d58a37SPierre Jolivet ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr); 32528d58a37SPierre Jolivet ierr = DMLocalSNESGetContext(dm,sdm,&dmlocalsnes);CHKERRQ(ierr); 32628d58a37SPierre Jolivet if (func) *func = dmlocalsnes->boundarylocal; 32728d58a37SPierre Jolivet if (ctx) *ctx = dmlocalsnes->boundarylocalctx; 32828d58a37SPierre Jolivet PetscFunctionReturn(0); 32928d58a37SPierre Jolivet } 33028d58a37SPierre Jolivet 33128d58a37SPierre Jolivet /*@C 33228d58a37SPierre Jolivet DMSNESGetJacobianLocal - the local Jacobian evaluation function set with DMSNESSetJacobianLocal. 33328d58a37SPierre Jolivet 33428d58a37SPierre Jolivet Logically Collective 33528d58a37SPierre Jolivet 336*4165533cSJose E. Roman Input Parameter: 337f0fc11ceSJed Brown . dm - DM with the associated callback 33828d58a37SPierre Jolivet 339*4165533cSJose E. Roman Output Parameters: 34028d58a37SPierre Jolivet + func - local Jacobian evaluation 34128d58a37SPierre Jolivet - ctx - context for local Jacobian evaluation 34228d58a37SPierre Jolivet 34328d58a37SPierre Jolivet Level: beginner 34428d58a37SPierre Jolivet 34528d58a37SPierre Jolivet .seealso: DMSNESSetJacobianLocal(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 34628d58a37SPierre Jolivet @*/ 34728d58a37SPierre Jolivet PetscErrorCode DMSNESGetJacobianLocal(DM dm,PetscErrorCode (**func)(DM,Vec,Mat,Mat,void*),void **ctx) 34828d58a37SPierre Jolivet { 34928d58a37SPierre Jolivet PetscErrorCode ierr; 35028d58a37SPierre Jolivet DMSNES sdm; 35128d58a37SPierre Jolivet DMSNES_Local *dmlocalsnes; 35228d58a37SPierre Jolivet 35328d58a37SPierre Jolivet PetscFunctionBegin; 35428d58a37SPierre Jolivet PetscValidHeaderSpecific(dm,DM_CLASSID,1); 35528d58a37SPierre Jolivet ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr); 35628d58a37SPierre Jolivet ierr = DMLocalSNESGetContext(dm,sdm,&dmlocalsnes);CHKERRQ(ierr); 35728d58a37SPierre Jolivet if (func) *func = dmlocalsnes->jacobianlocal; 35828d58a37SPierre Jolivet if (ctx) *ctx = dmlocalsnes->jacobianlocalctx; 359ff35dfedSBarry Smith PetscFunctionReturn(0); 360ff35dfedSBarry Smith } 361