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); 19*5ed5e208SMatthew 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; 28*5ed5e208SMatthew G. Knepley if (sdm->data != oldsdm->data) { 29*5ed5e208SMatthew G. Knepley ierr = PetscFree(sdm->data);CHKERRQ(ierr); 30b00a9115SJed Brown ierr = PetscNewLog(sdm,(DMSNES_Local**)&sdm->data);CHKERRQ(ierr); 31*5ed5e208SMatthew 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 { 54ff35dfedSBarry Smith PetscErrorCode ierr; 55ff35dfedSBarry Smith DM dm; 56942e3340SBarry Smith DMSNES_Local *dmlocalsnes = (DMSNES_Local*)ctx; 57ff35dfedSBarry Smith Vec Xloc,Floc; 58ff35dfedSBarry Smith 59ff35dfedSBarry Smith PetscFunctionBegin; 60ff35dfedSBarry Smith PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 61ff35dfedSBarry Smith PetscValidHeaderSpecific(X,VEC_CLASSID,2); 62ff35dfedSBarry Smith PetscValidHeaderSpecific(F,VEC_CLASSID,3); 63ff35dfedSBarry Smith ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 64ff35dfedSBarry Smith ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 65ff35dfedSBarry Smith ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr); 66ff35dfedSBarry Smith ierr = VecZeroEntries(Xloc);CHKERRQ(ierr); 67bdd6f66aSToby Isaac if (dmlocalsnes->boundarylocal) {ierr = (*dmlocalsnes->boundarylocal)(dm,Xloc,dmlocalsnes->boundarylocalctx);CHKERRQ(ierr);} 68ff35dfedSBarry Smith ierr = VecZeroEntries(Floc);CHKERRQ(ierr); 69ff35dfedSBarry Smith ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 70ff35dfedSBarry Smith ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 71ff35dfedSBarry Smith CHKMEMQ; 72ff35dfedSBarry Smith ierr = (*dmlocalsnes->residuallocal)(dm,Xloc,Floc,dmlocalsnes->residuallocalctx);CHKERRQ(ierr); 73ff35dfedSBarry Smith CHKMEMQ; 74ff35dfedSBarry Smith ierr = VecZeroEntries(F);CHKERRQ(ierr); 75ff35dfedSBarry Smith ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 76ff35dfedSBarry Smith ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 77ff35dfedSBarry Smith ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr); 78ff35dfedSBarry Smith ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 79afa9dec6SMatthew G. Knepley { 80afa9dec6SMatthew G. Knepley char name[PETSC_MAX_PATH_LEN]; 81f2f0fc17SMatthew G. Knepley char oldname[PETSC_MAX_PATH_LEN]; 82f2f0fc17SMatthew G. Knepley const char *tmp; 83afa9dec6SMatthew G. Knepley PetscInt it; 84afa9dec6SMatthew G. Knepley 85afa9dec6SMatthew G. Knepley ierr = SNESGetIterationNumber(snes, &it);CHKERRQ(ierr); 86afa9dec6SMatthew G. Knepley ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "Solution, Iterate %d", (int) it);CHKERRQ(ierr); 87f2f0fc17SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) X, &tmp);CHKERRQ(ierr); 88f2f0fc17SMatthew G. Knepley ierr = PetscStrncpy(oldname, tmp, PETSC_MAX_PATH_LEN-1);CHKERRQ(ierr); 89afa9dec6SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) X, name);CHKERRQ(ierr); 9085997439SMatthew G. Knepley ierr = VecViewFromOptions(X, (PetscObject) snes, "-dmsnes_solution_vec_view");CHKERRQ(ierr); 91f2f0fc17SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) X, oldname);CHKERRQ(ierr); 92afa9dec6SMatthew G. Knepley ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "Residual, Iterate %d", (int) it);CHKERRQ(ierr); 93afa9dec6SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) F, name);CHKERRQ(ierr); 9485997439SMatthew G. Knepley ierr = VecViewFromOptions(F, (PetscObject) snes, "-dmsnes_residual_vec_view");CHKERRQ(ierr); 95afa9dec6SMatthew G. Knepley } 96ff35dfedSBarry Smith PetscFunctionReturn(0); 97ff35dfedSBarry Smith } 98ff35dfedSBarry Smith 99d1e9a80fSBarry Smith static PetscErrorCode SNESComputeJacobian_DMLocal(SNES snes,Vec X,Mat A,Mat B,void *ctx) 100ff35dfedSBarry Smith { 101ff35dfedSBarry Smith PetscErrorCode ierr; 102ff35dfedSBarry Smith DM dm; 103942e3340SBarry Smith DMSNES_Local *dmlocalsnes = (DMSNES_Local*)ctx; 104ff35dfedSBarry Smith Vec Xloc; 105ff35dfedSBarry Smith 106ff35dfedSBarry Smith PetscFunctionBegin; 107ff35dfedSBarry Smith ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 108ff35dfedSBarry Smith if (dmlocalsnes->jacobianlocal) { 109ff35dfedSBarry Smith ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 110ff35dfedSBarry Smith ierr = VecZeroEntries(Xloc);CHKERRQ(ierr); 111bdd6f66aSToby Isaac if (dmlocalsnes->boundarylocal) {ierr = (*dmlocalsnes->boundarylocal)(dm,Xloc,dmlocalsnes->boundarylocalctx);CHKERRQ(ierr);} 112ff35dfedSBarry Smith ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 113ff35dfedSBarry Smith ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 114ff35dfedSBarry Smith CHKMEMQ; 115d1e9a80fSBarry Smith ierr = (*dmlocalsnes->jacobianlocal)(dm,Xloc,A,B,dmlocalsnes->jacobianlocalctx);CHKERRQ(ierr); 116ff35dfedSBarry Smith CHKMEMQ; 117ff35dfedSBarry Smith ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 118ff35dfedSBarry Smith } else { 119ff35dfedSBarry Smith MatFDColoring fdcoloring; 120ff35dfedSBarry Smith ierr = PetscObjectQuery((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject*)&fdcoloring);CHKERRQ(ierr); 121ff35dfedSBarry Smith if (!fdcoloring) { 122ff35dfedSBarry Smith ISColoring coloring; 123ff35dfedSBarry Smith 124b412c318SBarry Smith ierr = DMCreateColoring(dm,dm->coloringtype,&coloring);CHKERRQ(ierr); 12594ab13aaSBarry Smith ierr = MatFDColoringCreate(B,coloring,&fdcoloring);CHKERRQ(ierr); 126ff35dfedSBarry Smith ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr); 127ff35dfedSBarry Smith switch (dm->coloringtype) { 128ff35dfedSBarry Smith case IS_COLORING_GLOBAL: 129ff35dfedSBarry Smith ierr = MatFDColoringSetFunction(fdcoloring,(PetscErrorCode (*)(void))SNESComputeFunction_DMLocal,dmlocalsnes);CHKERRQ(ierr); 130ff35dfedSBarry Smith break; 131ce94432eSBarry Smith default: SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"No support for coloring type '%s'",ISColoringTypes[dm->coloringtype]); 132ff35dfedSBarry Smith } 133ff35dfedSBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)fdcoloring,((PetscObject)dm)->prefix);CHKERRQ(ierr); 134ff35dfedSBarry Smith ierr = MatFDColoringSetFromOptions(fdcoloring);CHKERRQ(ierr); 13594ab13aaSBarry Smith ierr = MatFDColoringSetUp(B,coloring,fdcoloring);CHKERRQ(ierr); 136ff35dfedSBarry Smith ierr = PetscObjectCompose((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject)fdcoloring);CHKERRQ(ierr); 137ff35dfedSBarry Smith ierr = PetscObjectDereference((PetscObject)fdcoloring);CHKERRQ(ierr); 138ff35dfedSBarry Smith 139ff35dfedSBarry Smith /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call 140ff35dfedSBarry Smith * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the 141ff35dfedSBarry Smith * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually 142ff35dfedSBarry Smith * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be 143140e18c1SBarry Smith * taken when the PetscObjectList for the Vec inside MatFDColoring is destroyed. 144ff35dfedSBarry Smith */ 145ff35dfedSBarry Smith ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr); 146ff35dfedSBarry Smith } 147d1e9a80fSBarry Smith ierr = MatFDColoringApply(B,fdcoloring,X,snes);CHKERRQ(ierr); 148ff35dfedSBarry Smith } 149ff35dfedSBarry Smith /* This will be redundant if the user called both, but it's too common to forget. */ 15094ab13aaSBarry Smith if (A != B) { 15194ab13aaSBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 15294ab13aaSBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 153ff35dfedSBarry Smith } 154ff35dfedSBarry Smith PetscFunctionReturn(0); 155ff35dfedSBarry Smith } 156ff35dfedSBarry Smith 157ff35dfedSBarry Smith /*@C 158ff35dfedSBarry Smith DMSNESSetFunctionLocal - set a local residual evaluation function. This function is called with local vector 159ff35dfedSBarry Smith containing the local vector information PLUS ghost point information. It should compute a result for all local 160ff35dfedSBarry Smith elements and DMSNES will automatically accumulate the overlapping values. 161ff35dfedSBarry Smith 162ff35dfedSBarry Smith Logically Collective 163ff35dfedSBarry Smith 164ff35dfedSBarry Smith Input Arguments: 165ff35dfedSBarry Smith + dm - DM to associate callback with 166ff35dfedSBarry Smith . func - local residual evaluation 167ff35dfedSBarry Smith - ctx - optional context for local residual evaluation 168ff35dfedSBarry Smith 169ff35dfedSBarry Smith Level: beginner 170ff35dfedSBarry Smith 171ff35dfedSBarry Smith .seealso: DMSNESSetFunction(), DMDASNESSetJacobianLocal(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 172ff35dfedSBarry Smith @*/ 173ff35dfedSBarry Smith PetscErrorCode DMSNESSetFunctionLocal(DM dm,PetscErrorCode (*func)(DM,Vec,Vec,void*),void *ctx) 174ff35dfedSBarry Smith { 175ff35dfedSBarry Smith PetscErrorCode ierr; 176942e3340SBarry Smith DMSNES sdm; 177942e3340SBarry Smith DMSNES_Local *dmlocalsnes; 178ff35dfedSBarry Smith 179ff35dfedSBarry Smith PetscFunctionBegin; 180ff35dfedSBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 181942e3340SBarry Smith ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 182ff35dfedSBarry Smith ierr = DMLocalSNESGetContext(dm,sdm,&dmlocalsnes);CHKERRQ(ierr); 1831aa26658SKarl Rupp 184ff35dfedSBarry Smith dmlocalsnes->residuallocal = func; 185ff35dfedSBarry Smith dmlocalsnes->residuallocalctx = ctx; 1861aa26658SKarl Rupp 187ff35dfedSBarry Smith ierr = DMSNESSetFunction(dm,SNESComputeFunction_DMLocal,dmlocalsnes);CHKERRQ(ierr); 18822c6f798SBarry Smith if (!sdm->ops->computejacobian) { /* Call us for the Jacobian too, can be overridden by the user. */ 189ff35dfedSBarry Smith ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMLocal,dmlocalsnes);CHKERRQ(ierr); 190ff35dfedSBarry Smith } 191ff35dfedSBarry Smith PetscFunctionReturn(0); 192ff35dfedSBarry Smith } 193ff35dfedSBarry Smith 194bdd6f66aSToby Isaac /*@C 195bdd6f66aSToby Isaac DMSNESSetBoundaryLocal - set a local boundary value function. This function is called with local vector 196bdd6f66aSToby Isaac containing the local vector information PLUS ghost point information. It should insert values into the local 197bdd6f66aSToby Isaac vector that do not come from the global vector, such as essential boundary condition data. 198bdd6f66aSToby Isaac 199bdd6f66aSToby Isaac Logically Collective 200bdd6f66aSToby Isaac 201bdd6f66aSToby Isaac Input Arguments: 202bdd6f66aSToby Isaac + dm - DM to associate callback with 203bdd6f66aSToby Isaac . func - local boundary value evaluation 204bdd6f66aSToby Isaac - ctx - optional context for local boundary value evaluation 205bdd6f66aSToby Isaac 206bdd6f66aSToby Isaac Level: intermediate 207bdd6f66aSToby Isaac 208bdd6f66aSToby Isaac .seealso: DMSNESSetFunctionLocal(), DMDASNESSetJacobianLocal() 209bdd6f66aSToby Isaac @*/ 210bdd6f66aSToby Isaac PetscErrorCode DMSNESSetBoundaryLocal(DM dm,PetscErrorCode (*func)(DM,Vec,void*),void *ctx) 211bdd6f66aSToby Isaac { 212bdd6f66aSToby Isaac PetscErrorCode ierr; 213bdd6f66aSToby Isaac DMSNES sdm; 214bdd6f66aSToby Isaac DMSNES_Local *dmlocalsnes; 215bdd6f66aSToby Isaac 216bdd6f66aSToby Isaac PetscFunctionBegin; 217bdd6f66aSToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 218bdd6f66aSToby Isaac ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 219bdd6f66aSToby Isaac ierr = DMLocalSNESGetContext(dm,sdm,&dmlocalsnes);CHKERRQ(ierr); 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 232ff35dfedSBarry Smith Input Arguments: 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 239ff35dfedSBarry Smith .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 { 243ff35dfedSBarry Smith PetscErrorCode ierr; 244942e3340SBarry Smith DMSNES sdm; 245942e3340SBarry Smith DMSNES_Local *dmlocalsnes; 246ff35dfedSBarry Smith 247ff35dfedSBarry Smith PetscFunctionBegin; 248ff35dfedSBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 249942e3340SBarry Smith ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 250ff35dfedSBarry Smith ierr = DMLocalSNESGetContext(dm,sdm,&dmlocalsnes);CHKERRQ(ierr); 2511aa26658SKarl Rupp 252ff35dfedSBarry Smith dmlocalsnes->jacobianlocal = func; 253ff35dfedSBarry Smith dmlocalsnes->jacobianlocalctx = ctx; 2541aa26658SKarl Rupp 255ff35dfedSBarry Smith ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMLocal,dmlocalsnes);CHKERRQ(ierr); 256ff35dfedSBarry Smith PetscFunctionReturn(0); 257ff35dfedSBarry Smith } 258ff35dfedSBarry Smith 259