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*); 7*bdd6f66aSToby Isaac PetscErrorCode (*boundarylocal)(DM,Vec,void*); 8ff35dfedSBarry Smith void *residuallocalctx; 9ff35dfedSBarry Smith void *jacobianlocalctx; 10*bdd6f66aSToby Isaac void *boundarylocalctx; 11942e3340SBarry Smith } DMSNES_Local; 12ff35dfedSBarry Smith 13ff35dfedSBarry Smith #undef __FUNCT__ 14942e3340SBarry Smith #define __FUNCT__ "DMSNESDestroy_DMLocal" 15942e3340SBarry Smith static PetscErrorCode DMSNESDestroy_DMLocal(DMSNES sdm) 16ff35dfedSBarry Smith { 17ff35dfedSBarry Smith PetscErrorCode ierr; 18ff35dfedSBarry Smith 19ff35dfedSBarry Smith PetscFunctionBegin; 20ff35dfedSBarry Smith ierr = PetscFree(sdm->data);CHKERRQ(ierr); 21ff35dfedSBarry Smith PetscFunctionReturn(0); 22ff35dfedSBarry Smith } 23ff35dfedSBarry Smith 24ff35dfedSBarry Smith #undef __FUNCT__ 25942e3340SBarry Smith #define __FUNCT__ "DMSNESDuplicate_DMLocal" 2622c6f798SBarry Smith static PetscErrorCode DMSNESDuplicate_DMLocal(DMSNES oldsdm,DMSNES sdm) 27ff35dfedSBarry Smith { 28ff35dfedSBarry Smith PetscErrorCode ierr; 29ff35dfedSBarry Smith 30ff35dfedSBarry Smith PetscFunctionBegin; 31b00a9115SJed Brown ierr = PetscNewLog(sdm,(DMSNES_Local**)&sdm->data);CHKERRQ(ierr); 32ff35dfedSBarry Smith if (oldsdm->data) { 33942e3340SBarry Smith ierr = PetscMemcpy(sdm->data,oldsdm->data,sizeof(DMSNES_Local));CHKERRQ(ierr); 34ff35dfedSBarry Smith } 35ff35dfedSBarry Smith PetscFunctionReturn(0); 36ff35dfedSBarry Smith } 37ff35dfedSBarry Smith 38ff35dfedSBarry Smith #undef __FUNCT__ 39ff35dfedSBarry Smith #define __FUNCT__ "DMLocalSNESGetContext" 40942e3340SBarry Smith static PetscErrorCode DMLocalSNESGetContext(DM dm,DMSNES sdm,DMSNES_Local **dmlocalsnes) 41ff35dfedSBarry Smith { 42ff35dfedSBarry Smith PetscErrorCode ierr; 43ff35dfedSBarry Smith 44ff35dfedSBarry Smith PetscFunctionBegin; 450298fd71SBarry Smith *dmlocalsnes = NULL; 46ff35dfedSBarry Smith if (!sdm->data) { 47b00a9115SJed Brown ierr = PetscNewLog(dm,(DMSNES_Local**)&sdm->data);CHKERRQ(ierr); 481aa26658SKarl Rupp 4922c6f798SBarry Smith sdm->ops->destroy = DMSNESDestroy_DMLocal; 5022c6f798SBarry Smith sdm->ops->duplicate = DMSNESDuplicate_DMLocal; 51ff35dfedSBarry Smith } 52942e3340SBarry Smith *dmlocalsnes = (DMSNES_Local*)sdm->data; 53ff35dfedSBarry Smith PetscFunctionReturn(0); 54ff35dfedSBarry Smith } 55ff35dfedSBarry Smith 56ff35dfedSBarry Smith #undef __FUNCT__ 57ff35dfedSBarry Smith #define __FUNCT__ "SNESComputeFunction_DMLocal" 58b2270762SBarry Smith static PetscErrorCode SNESComputeFunction_DMLocal(SNES snes,Vec X,Vec F,void *ctx) 59ff35dfedSBarry Smith { 60ff35dfedSBarry Smith PetscErrorCode ierr; 61ff35dfedSBarry Smith DM dm; 62942e3340SBarry Smith DMSNES_Local *dmlocalsnes = (DMSNES_Local*)ctx; 63ff35dfedSBarry Smith Vec Xloc,Floc; 64ff35dfedSBarry Smith 65ff35dfedSBarry Smith PetscFunctionBegin; 66ff35dfedSBarry Smith PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 67ff35dfedSBarry Smith PetscValidHeaderSpecific(X,VEC_CLASSID,2); 68ff35dfedSBarry Smith PetscValidHeaderSpecific(F,VEC_CLASSID,3); 69ff35dfedSBarry Smith ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 70ff35dfedSBarry Smith ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 71ff35dfedSBarry Smith ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr); 72ff35dfedSBarry Smith ierr = VecZeroEntries(Xloc);CHKERRQ(ierr); 73*bdd6f66aSToby Isaac if (dmlocalsnes->boundarylocal) {ierr = (*dmlocalsnes->boundarylocal)(dm,Xloc,dmlocalsnes->boundarylocalctx);CHKERRQ(ierr);} 74ff35dfedSBarry Smith ierr = VecZeroEntries(Floc);CHKERRQ(ierr); 75ff35dfedSBarry Smith ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 76ff35dfedSBarry Smith ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 77ff35dfedSBarry Smith CHKMEMQ; 78ff35dfedSBarry Smith ierr = (*dmlocalsnes->residuallocal)(dm,Xloc,Floc,dmlocalsnes->residuallocalctx);CHKERRQ(ierr); 79ff35dfedSBarry Smith CHKMEMQ; 80ff35dfedSBarry Smith ierr = VecZeroEntries(F);CHKERRQ(ierr); 81ff35dfedSBarry Smith ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 82ff35dfedSBarry Smith ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 83ff35dfedSBarry Smith ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr); 84ff35dfedSBarry Smith ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 85ff35dfedSBarry Smith PetscFunctionReturn(0); 86ff35dfedSBarry Smith } 87ff35dfedSBarry Smith 88ff35dfedSBarry Smith #undef __FUNCT__ 89ff35dfedSBarry Smith #define __FUNCT__ "SNESComputeJacobian_DMLocal" 90d1e9a80fSBarry Smith static PetscErrorCode SNESComputeJacobian_DMLocal(SNES snes,Vec X,Mat A,Mat B,void *ctx) 91ff35dfedSBarry Smith { 92ff35dfedSBarry Smith PetscErrorCode ierr; 93ff35dfedSBarry Smith DM dm; 94942e3340SBarry Smith DMSNES_Local *dmlocalsnes = (DMSNES_Local*)ctx; 95ff35dfedSBarry Smith Vec Xloc; 96ff35dfedSBarry Smith 97ff35dfedSBarry Smith PetscFunctionBegin; 98ff35dfedSBarry Smith ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 99ff35dfedSBarry Smith if (dmlocalsnes->jacobianlocal) { 100ff35dfedSBarry Smith ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 101ff35dfedSBarry Smith ierr = VecZeroEntries(Xloc);CHKERRQ(ierr); 102*bdd6f66aSToby Isaac if (dmlocalsnes->boundarylocal) {ierr = (*dmlocalsnes->boundarylocal)(dm,Xloc,dmlocalsnes->boundarylocalctx);CHKERRQ(ierr);} 103ff35dfedSBarry Smith ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 104ff35dfedSBarry Smith ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 105ff35dfedSBarry Smith CHKMEMQ; 106d1e9a80fSBarry Smith ierr = (*dmlocalsnes->jacobianlocal)(dm,Xloc,A,B,dmlocalsnes->jacobianlocalctx);CHKERRQ(ierr); 107ff35dfedSBarry Smith CHKMEMQ; 108ff35dfedSBarry Smith ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 109ff35dfedSBarry Smith } else { 110ff35dfedSBarry Smith MatFDColoring fdcoloring; 111ff35dfedSBarry Smith ierr = PetscObjectQuery((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject*)&fdcoloring);CHKERRQ(ierr); 112ff35dfedSBarry Smith if (!fdcoloring) { 113ff35dfedSBarry Smith ISColoring coloring; 114ff35dfedSBarry Smith 115b412c318SBarry Smith ierr = DMCreateColoring(dm,dm->coloringtype,&coloring);CHKERRQ(ierr); 11694ab13aaSBarry Smith ierr = MatFDColoringCreate(B,coloring,&fdcoloring);CHKERRQ(ierr); 117ff35dfedSBarry Smith ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr); 118ff35dfedSBarry Smith switch (dm->coloringtype) { 119ff35dfedSBarry Smith case IS_COLORING_GLOBAL: 120ff35dfedSBarry Smith ierr = MatFDColoringSetFunction(fdcoloring,(PetscErrorCode (*)(void))SNESComputeFunction_DMLocal,dmlocalsnes);CHKERRQ(ierr); 121ff35dfedSBarry Smith break; 122ce94432eSBarry Smith default: SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"No support for coloring type '%s'",ISColoringTypes[dm->coloringtype]); 123ff35dfedSBarry Smith } 124ff35dfedSBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)fdcoloring,((PetscObject)dm)->prefix);CHKERRQ(ierr); 125ff35dfedSBarry Smith ierr = MatFDColoringSetFromOptions(fdcoloring);CHKERRQ(ierr); 12694ab13aaSBarry Smith ierr = MatFDColoringSetUp(B,coloring,fdcoloring);CHKERRQ(ierr); 127ff35dfedSBarry Smith ierr = PetscObjectCompose((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject)fdcoloring);CHKERRQ(ierr); 128ff35dfedSBarry Smith ierr = PetscObjectDereference((PetscObject)fdcoloring);CHKERRQ(ierr); 129ff35dfedSBarry Smith 130ff35dfedSBarry Smith /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call 131ff35dfedSBarry Smith * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the 132ff35dfedSBarry Smith * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually 133ff35dfedSBarry Smith * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be 134140e18c1SBarry Smith * taken when the PetscObjectList for the Vec inside MatFDColoring is destroyed. 135ff35dfedSBarry Smith */ 136ff35dfedSBarry Smith ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr); 137ff35dfedSBarry Smith } 138d1e9a80fSBarry Smith ierr = MatFDColoringApply(B,fdcoloring,X,snes);CHKERRQ(ierr); 139ff35dfedSBarry Smith } 140ff35dfedSBarry Smith /* This will be redundant if the user called both, but it's too common to forget. */ 14194ab13aaSBarry Smith if (A != B) { 14294ab13aaSBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14394ab13aaSBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 144ff35dfedSBarry Smith } 145ff35dfedSBarry Smith PetscFunctionReturn(0); 146ff35dfedSBarry Smith } 147ff35dfedSBarry Smith 148ff35dfedSBarry Smith #undef __FUNCT__ 149ff35dfedSBarry Smith #define __FUNCT__ "DMSNESSetFunctionLocal" 150ff35dfedSBarry Smith /*@C 151ff35dfedSBarry Smith DMSNESSetFunctionLocal - set a local residual evaluation function. This function is called with local vector 152ff35dfedSBarry Smith containing the local vector information PLUS ghost point information. It should compute a result for all local 153ff35dfedSBarry Smith elements and DMSNES will automatically accumulate the overlapping values. 154ff35dfedSBarry Smith 155ff35dfedSBarry Smith Logically Collective 156ff35dfedSBarry Smith 157ff35dfedSBarry Smith Input Arguments: 158ff35dfedSBarry Smith + dm - DM to associate callback with 159ff35dfedSBarry Smith . func - local residual evaluation 160ff35dfedSBarry Smith - ctx - optional context for local residual evaluation 161ff35dfedSBarry Smith 162ff35dfedSBarry Smith Level: beginner 163ff35dfedSBarry Smith 164ff35dfedSBarry Smith .seealso: DMSNESSetFunction(), DMDASNESSetJacobianLocal(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 165ff35dfedSBarry Smith @*/ 166ff35dfedSBarry Smith PetscErrorCode DMSNESSetFunctionLocal(DM dm,PetscErrorCode (*func)(DM,Vec,Vec,void*),void *ctx) 167ff35dfedSBarry Smith { 168ff35dfedSBarry Smith PetscErrorCode ierr; 169942e3340SBarry Smith DMSNES sdm; 170942e3340SBarry Smith DMSNES_Local *dmlocalsnes; 171ff35dfedSBarry Smith 172ff35dfedSBarry Smith PetscFunctionBegin; 173ff35dfedSBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 174942e3340SBarry Smith ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 175ff35dfedSBarry Smith ierr = DMLocalSNESGetContext(dm,sdm,&dmlocalsnes);CHKERRQ(ierr); 1761aa26658SKarl Rupp 177ff35dfedSBarry Smith dmlocalsnes->residuallocal = func; 178ff35dfedSBarry Smith dmlocalsnes->residuallocalctx = ctx; 1791aa26658SKarl Rupp 180ff35dfedSBarry Smith ierr = DMSNESSetFunction(dm,SNESComputeFunction_DMLocal,dmlocalsnes);CHKERRQ(ierr); 18122c6f798SBarry Smith if (!sdm->ops->computejacobian) { /* Call us for the Jacobian too, can be overridden by the user. */ 182ff35dfedSBarry Smith ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMLocal,dmlocalsnes);CHKERRQ(ierr); 183ff35dfedSBarry Smith } 184ff35dfedSBarry Smith PetscFunctionReturn(0); 185ff35dfedSBarry Smith } 186ff35dfedSBarry Smith 187ff35dfedSBarry Smith #undef __FUNCT__ 188*bdd6f66aSToby Isaac #define __FUNCT__ "DMSNESSetBoundaryLocal" 189*bdd6f66aSToby Isaac /*@C 190*bdd6f66aSToby Isaac DMSNESSetBoundaryLocal - set a local boundary value function. This function is called with local vector 191*bdd6f66aSToby Isaac containing the local vector information PLUS ghost point information. It should insert values into the local 192*bdd6f66aSToby Isaac vector that do not come from the global vector, such as essential boundary condition data. 193*bdd6f66aSToby Isaac 194*bdd6f66aSToby Isaac Logically Collective 195*bdd6f66aSToby Isaac 196*bdd6f66aSToby Isaac Input Arguments: 197*bdd6f66aSToby Isaac + dm - DM to associate callback with 198*bdd6f66aSToby Isaac . func - local boundary value evaluation 199*bdd6f66aSToby Isaac - ctx - optional context for local boundary value evaluation 200*bdd6f66aSToby Isaac 201*bdd6f66aSToby Isaac Level: intermediate 202*bdd6f66aSToby Isaac 203*bdd6f66aSToby Isaac .seealso: DMSNESSetFunctionLocal(), DMDASNESSetJacobianLocal() 204*bdd6f66aSToby Isaac @*/ 205*bdd6f66aSToby Isaac PetscErrorCode DMSNESSetBoundaryLocal(DM dm,PetscErrorCode (*func)(DM,Vec,void*),void *ctx) 206*bdd6f66aSToby Isaac { 207*bdd6f66aSToby Isaac PetscErrorCode ierr; 208*bdd6f66aSToby Isaac DMSNES sdm; 209*bdd6f66aSToby Isaac DMSNES_Local *dmlocalsnes; 210*bdd6f66aSToby Isaac 211*bdd6f66aSToby Isaac PetscFunctionBegin; 212*bdd6f66aSToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 213*bdd6f66aSToby Isaac ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 214*bdd6f66aSToby Isaac ierr = DMLocalSNESGetContext(dm,sdm,&dmlocalsnes);CHKERRQ(ierr); 215*bdd6f66aSToby Isaac 216*bdd6f66aSToby Isaac dmlocalsnes->boundarylocal = func; 217*bdd6f66aSToby Isaac dmlocalsnes->boundarylocalctx = ctx; 218*bdd6f66aSToby Isaac 219*bdd6f66aSToby Isaac PetscFunctionReturn(0); 220*bdd6f66aSToby Isaac } 221*bdd6f66aSToby Isaac 222*bdd6f66aSToby Isaac #undef __FUNCT__ 223ff35dfedSBarry Smith #define __FUNCT__ "DMSNESSetJacobianLocal" 224ff35dfedSBarry Smith /*@C 225ff35dfedSBarry Smith DMSNESSetJacobianLocal - set a local Jacobian evaluation function 226ff35dfedSBarry Smith 227ff35dfedSBarry Smith Logically Collective 228ff35dfedSBarry Smith 229ff35dfedSBarry Smith Input Arguments: 230ff35dfedSBarry Smith + dm - DM to associate callback with 231ff35dfedSBarry Smith . func - local Jacobian evaluation 232ff35dfedSBarry Smith - ctx - optional context for local Jacobian evaluation 233ff35dfedSBarry Smith 234ff35dfedSBarry Smith Level: beginner 235ff35dfedSBarry Smith 236ff35dfedSBarry Smith .seealso: DMSNESSetJacobian(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 237ff35dfedSBarry Smith @*/ 238d1e9a80fSBarry Smith PetscErrorCode DMSNESSetJacobianLocal(DM dm,PetscErrorCode (*func)(DM,Vec,Mat,Mat,void*),void *ctx) 239ff35dfedSBarry Smith { 240ff35dfedSBarry Smith PetscErrorCode ierr; 241942e3340SBarry Smith DMSNES sdm; 242942e3340SBarry Smith DMSNES_Local *dmlocalsnes; 243ff35dfedSBarry Smith 244ff35dfedSBarry Smith PetscFunctionBegin; 245ff35dfedSBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 246942e3340SBarry Smith ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 247ff35dfedSBarry Smith ierr = DMLocalSNESGetContext(dm,sdm,&dmlocalsnes);CHKERRQ(ierr); 2481aa26658SKarl Rupp 249ff35dfedSBarry Smith dmlocalsnes->jacobianlocal = func; 250ff35dfedSBarry Smith dmlocalsnes->jacobianlocalctx = ctx; 2511aa26658SKarl Rupp 252ff35dfedSBarry Smith ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMLocal,dmlocalsnes);CHKERRQ(ierr); 253ff35dfedSBarry Smith PetscFunctionReturn(0); 254ff35dfedSBarry Smith } 255ff35dfedSBarry Smith 256