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); 19ff35dfedSBarry Smith PetscFunctionReturn(0); 20ff35dfedSBarry Smith } 21ff35dfedSBarry Smith 2222c6f798SBarry Smith static PetscErrorCode DMSNESDuplicate_DMLocal(DMSNES oldsdm,DMSNES sdm) 23ff35dfedSBarry Smith { 24ff35dfedSBarry Smith PetscErrorCode ierr; 25ff35dfedSBarry Smith 26ff35dfedSBarry Smith PetscFunctionBegin; 27b00a9115SJed Brown ierr = PetscNewLog(sdm,(DMSNES_Local**)&sdm->data);CHKERRQ(ierr); 28ff35dfedSBarry Smith if (oldsdm->data) { 29942e3340SBarry Smith ierr = PetscMemcpy(sdm->data,oldsdm->data,sizeof(DMSNES_Local));CHKERRQ(ierr); 30ff35dfedSBarry Smith } 31ff35dfedSBarry Smith PetscFunctionReturn(0); 32ff35dfedSBarry Smith } 33ff35dfedSBarry Smith 34942e3340SBarry Smith static PetscErrorCode DMLocalSNESGetContext(DM dm,DMSNES sdm,DMSNES_Local **dmlocalsnes) 35ff35dfedSBarry Smith { 36ff35dfedSBarry Smith PetscErrorCode ierr; 37ff35dfedSBarry Smith 38ff35dfedSBarry Smith PetscFunctionBegin; 390298fd71SBarry Smith *dmlocalsnes = NULL; 40ff35dfedSBarry Smith if (!sdm->data) { 41b00a9115SJed Brown ierr = PetscNewLog(dm,(DMSNES_Local**)&sdm->data);CHKERRQ(ierr); 421aa26658SKarl Rupp 4322c6f798SBarry Smith sdm->ops->destroy = DMSNESDestroy_DMLocal; 4422c6f798SBarry Smith sdm->ops->duplicate = DMSNESDuplicate_DMLocal; 45ff35dfedSBarry Smith } 46942e3340SBarry Smith *dmlocalsnes = (DMSNES_Local*)sdm->data; 47ff35dfedSBarry Smith PetscFunctionReturn(0); 48ff35dfedSBarry Smith } 49ff35dfedSBarry Smith 50b2270762SBarry Smith static PetscErrorCode SNESComputeFunction_DMLocal(SNES snes,Vec X,Vec F,void *ctx) 51ff35dfedSBarry Smith { 52ff35dfedSBarry Smith PetscErrorCode ierr; 53ff35dfedSBarry Smith DM dm; 54942e3340SBarry Smith DMSNES_Local *dmlocalsnes = (DMSNES_Local*)ctx; 55ff35dfedSBarry Smith Vec Xloc,Floc; 56ff35dfedSBarry Smith 57ff35dfedSBarry Smith PetscFunctionBegin; 58ff35dfedSBarry Smith PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 59ff35dfedSBarry Smith PetscValidHeaderSpecific(X,VEC_CLASSID,2); 60ff35dfedSBarry Smith PetscValidHeaderSpecific(F,VEC_CLASSID,3); 61ff35dfedSBarry Smith ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 62ff35dfedSBarry Smith ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 63ff35dfedSBarry Smith ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr); 64ff35dfedSBarry Smith ierr = VecZeroEntries(Xloc);CHKERRQ(ierr); 65bdd6f66aSToby Isaac if (dmlocalsnes->boundarylocal) {ierr = (*dmlocalsnes->boundarylocal)(dm,Xloc,dmlocalsnes->boundarylocalctx);CHKERRQ(ierr);} 66ff35dfedSBarry Smith ierr = VecZeroEntries(Floc);CHKERRQ(ierr); 67ff35dfedSBarry Smith ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 68ff35dfedSBarry Smith ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 69ff35dfedSBarry Smith CHKMEMQ; 70ff35dfedSBarry Smith ierr = (*dmlocalsnes->residuallocal)(dm,Xloc,Floc,dmlocalsnes->residuallocalctx);CHKERRQ(ierr); 71ff35dfedSBarry Smith CHKMEMQ; 72ff35dfedSBarry Smith ierr = VecZeroEntries(F);CHKERRQ(ierr); 73ff35dfedSBarry Smith ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 74ff35dfedSBarry Smith ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 75ff35dfedSBarry Smith ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr); 76ff35dfedSBarry Smith ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 77afa9dec6SMatthew G. Knepley { 78afa9dec6SMatthew G. Knepley char name[PETSC_MAX_PATH_LEN]; 79*f2f0fc17SMatthew G. Knepley char oldname[PETSC_MAX_PATH_LEN]; 80*f2f0fc17SMatthew G. Knepley const char *tmp; 81afa9dec6SMatthew G. Knepley PetscInt it; 82afa9dec6SMatthew G. Knepley 83afa9dec6SMatthew G. Knepley ierr = SNESGetIterationNumber(snes, &it);CHKERRQ(ierr); 84afa9dec6SMatthew G. Knepley ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "Solution, Iterate %d", (int) it);CHKERRQ(ierr); 85*f2f0fc17SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) X, &tmp);CHKERRQ(ierr); 86*f2f0fc17SMatthew G. Knepley ierr = PetscStrncpy(oldname, tmp, PETSC_MAX_PATH_LEN-1);CHKERRQ(ierr); 87afa9dec6SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) X, name);CHKERRQ(ierr); 8885997439SMatthew G. Knepley ierr = VecViewFromOptions(X, (PetscObject) snes, "-dmsnes_solution_vec_view");CHKERRQ(ierr); 89*f2f0fc17SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) X, oldname);CHKERRQ(ierr); 90afa9dec6SMatthew G. Knepley ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "Residual, Iterate %d", (int) it);CHKERRQ(ierr); 91afa9dec6SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) F, name);CHKERRQ(ierr); 9285997439SMatthew G. Knepley ierr = VecViewFromOptions(F, (PetscObject) snes, "-dmsnes_residual_vec_view");CHKERRQ(ierr); 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 { 99ff35dfedSBarry Smith PetscErrorCode ierr; 100ff35dfedSBarry Smith DM dm; 101942e3340SBarry Smith DMSNES_Local *dmlocalsnes = (DMSNES_Local*)ctx; 102ff35dfedSBarry Smith Vec Xloc; 103ff35dfedSBarry Smith 104ff35dfedSBarry Smith PetscFunctionBegin; 105ff35dfedSBarry Smith ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 106ff35dfedSBarry Smith if (dmlocalsnes->jacobianlocal) { 107ff35dfedSBarry Smith ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 108ff35dfedSBarry Smith ierr = VecZeroEntries(Xloc);CHKERRQ(ierr); 109bdd6f66aSToby Isaac if (dmlocalsnes->boundarylocal) {ierr = (*dmlocalsnes->boundarylocal)(dm,Xloc,dmlocalsnes->boundarylocalctx);CHKERRQ(ierr);} 110ff35dfedSBarry Smith ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 111ff35dfedSBarry Smith ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 112ff35dfedSBarry Smith CHKMEMQ; 113d1e9a80fSBarry Smith ierr = (*dmlocalsnes->jacobianlocal)(dm,Xloc,A,B,dmlocalsnes->jacobianlocalctx);CHKERRQ(ierr); 114ff35dfedSBarry Smith CHKMEMQ; 115ff35dfedSBarry Smith ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 116ff35dfedSBarry Smith } else { 117ff35dfedSBarry Smith MatFDColoring fdcoloring; 118ff35dfedSBarry Smith ierr = PetscObjectQuery((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject*)&fdcoloring);CHKERRQ(ierr); 119ff35dfedSBarry Smith if (!fdcoloring) { 120ff35dfedSBarry Smith ISColoring coloring; 121ff35dfedSBarry Smith 122b412c318SBarry Smith ierr = DMCreateColoring(dm,dm->coloringtype,&coloring);CHKERRQ(ierr); 12394ab13aaSBarry Smith ierr = MatFDColoringCreate(B,coloring,&fdcoloring);CHKERRQ(ierr); 124ff35dfedSBarry Smith ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr); 125ff35dfedSBarry Smith switch (dm->coloringtype) { 126ff35dfedSBarry Smith case IS_COLORING_GLOBAL: 127ff35dfedSBarry Smith ierr = MatFDColoringSetFunction(fdcoloring,(PetscErrorCode (*)(void))SNESComputeFunction_DMLocal,dmlocalsnes);CHKERRQ(ierr); 128ff35dfedSBarry Smith break; 129ce94432eSBarry Smith default: SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"No support for coloring type '%s'",ISColoringTypes[dm->coloringtype]); 130ff35dfedSBarry Smith } 131ff35dfedSBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)fdcoloring,((PetscObject)dm)->prefix);CHKERRQ(ierr); 132ff35dfedSBarry Smith ierr = MatFDColoringSetFromOptions(fdcoloring);CHKERRQ(ierr); 13394ab13aaSBarry Smith ierr = MatFDColoringSetUp(B,coloring,fdcoloring);CHKERRQ(ierr); 134ff35dfedSBarry Smith ierr = PetscObjectCompose((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject)fdcoloring);CHKERRQ(ierr); 135ff35dfedSBarry Smith ierr = PetscObjectDereference((PetscObject)fdcoloring);CHKERRQ(ierr); 136ff35dfedSBarry Smith 137ff35dfedSBarry Smith /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call 138ff35dfedSBarry Smith * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the 139ff35dfedSBarry Smith * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually 140ff35dfedSBarry Smith * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be 141140e18c1SBarry Smith * taken when the PetscObjectList for the Vec inside MatFDColoring is destroyed. 142ff35dfedSBarry Smith */ 143ff35dfedSBarry Smith ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr); 144ff35dfedSBarry Smith } 145d1e9a80fSBarry Smith ierr = MatFDColoringApply(B,fdcoloring,X,snes);CHKERRQ(ierr); 146ff35dfedSBarry Smith } 147ff35dfedSBarry Smith /* This will be redundant if the user called both, but it's too common to forget. */ 14894ab13aaSBarry Smith if (A != B) { 14994ab13aaSBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 15094ab13aaSBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 151ff35dfedSBarry Smith } 152ff35dfedSBarry Smith PetscFunctionReturn(0); 153ff35dfedSBarry Smith } 154ff35dfedSBarry Smith 155ff35dfedSBarry Smith /*@C 156ff35dfedSBarry Smith DMSNESSetFunctionLocal - set a local residual evaluation function. This function is called with local vector 157ff35dfedSBarry Smith containing the local vector information PLUS ghost point information. It should compute a result for all local 158ff35dfedSBarry Smith elements and DMSNES will automatically accumulate the overlapping values. 159ff35dfedSBarry Smith 160ff35dfedSBarry Smith Logically Collective 161ff35dfedSBarry Smith 162ff35dfedSBarry Smith Input Arguments: 163ff35dfedSBarry Smith + dm - DM to associate callback with 164ff35dfedSBarry Smith . func - local residual evaluation 165ff35dfedSBarry Smith - ctx - optional context for local residual evaluation 166ff35dfedSBarry Smith 167ff35dfedSBarry Smith Level: beginner 168ff35dfedSBarry Smith 169ff35dfedSBarry Smith .seealso: DMSNESSetFunction(), DMDASNESSetJacobianLocal(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 170ff35dfedSBarry Smith @*/ 171ff35dfedSBarry Smith PetscErrorCode DMSNESSetFunctionLocal(DM dm,PetscErrorCode (*func)(DM,Vec,Vec,void*),void *ctx) 172ff35dfedSBarry Smith { 173ff35dfedSBarry Smith PetscErrorCode ierr; 174942e3340SBarry Smith DMSNES sdm; 175942e3340SBarry Smith DMSNES_Local *dmlocalsnes; 176ff35dfedSBarry Smith 177ff35dfedSBarry Smith PetscFunctionBegin; 178ff35dfedSBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 179942e3340SBarry Smith ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 180ff35dfedSBarry Smith ierr = DMLocalSNESGetContext(dm,sdm,&dmlocalsnes);CHKERRQ(ierr); 1811aa26658SKarl Rupp 182ff35dfedSBarry Smith dmlocalsnes->residuallocal = func; 183ff35dfedSBarry Smith dmlocalsnes->residuallocalctx = ctx; 1841aa26658SKarl Rupp 185ff35dfedSBarry Smith ierr = DMSNESSetFunction(dm,SNESComputeFunction_DMLocal,dmlocalsnes);CHKERRQ(ierr); 18622c6f798SBarry Smith if (!sdm->ops->computejacobian) { /* Call us for the Jacobian too, can be overridden by the user. */ 187ff35dfedSBarry Smith ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMLocal,dmlocalsnes);CHKERRQ(ierr); 188ff35dfedSBarry Smith } 189ff35dfedSBarry Smith PetscFunctionReturn(0); 190ff35dfedSBarry Smith } 191ff35dfedSBarry Smith 192bdd6f66aSToby Isaac /*@C 193bdd6f66aSToby Isaac DMSNESSetBoundaryLocal - set a local boundary value function. This function is called with local vector 194bdd6f66aSToby Isaac containing the local vector information PLUS ghost point information. It should insert values into the local 195bdd6f66aSToby Isaac vector that do not come from the global vector, such as essential boundary condition data. 196bdd6f66aSToby Isaac 197bdd6f66aSToby Isaac Logically Collective 198bdd6f66aSToby Isaac 199bdd6f66aSToby Isaac Input Arguments: 200bdd6f66aSToby Isaac + dm - DM to associate callback with 201bdd6f66aSToby Isaac . func - local boundary value evaluation 202bdd6f66aSToby Isaac - ctx - optional context for local boundary value evaluation 203bdd6f66aSToby Isaac 204bdd6f66aSToby Isaac Level: intermediate 205bdd6f66aSToby Isaac 206bdd6f66aSToby Isaac .seealso: DMSNESSetFunctionLocal(), DMDASNESSetJacobianLocal() 207bdd6f66aSToby Isaac @*/ 208bdd6f66aSToby Isaac PetscErrorCode DMSNESSetBoundaryLocal(DM dm,PetscErrorCode (*func)(DM,Vec,void*),void *ctx) 209bdd6f66aSToby Isaac { 210bdd6f66aSToby Isaac PetscErrorCode ierr; 211bdd6f66aSToby Isaac DMSNES sdm; 212bdd6f66aSToby Isaac DMSNES_Local *dmlocalsnes; 213bdd6f66aSToby Isaac 214bdd6f66aSToby Isaac PetscFunctionBegin; 215bdd6f66aSToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 216bdd6f66aSToby Isaac ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 217bdd6f66aSToby Isaac ierr = DMLocalSNESGetContext(dm,sdm,&dmlocalsnes);CHKERRQ(ierr); 218bdd6f66aSToby Isaac 219bdd6f66aSToby Isaac dmlocalsnes->boundarylocal = func; 220bdd6f66aSToby Isaac dmlocalsnes->boundarylocalctx = ctx; 221bdd6f66aSToby Isaac 222bdd6f66aSToby Isaac PetscFunctionReturn(0); 223bdd6f66aSToby Isaac } 224bdd6f66aSToby Isaac 225ff35dfedSBarry Smith /*@C 226ff35dfedSBarry Smith DMSNESSetJacobianLocal - set a local Jacobian evaluation function 227ff35dfedSBarry Smith 228ff35dfedSBarry Smith Logically Collective 229ff35dfedSBarry Smith 230ff35dfedSBarry Smith Input Arguments: 231ff35dfedSBarry Smith + dm - DM to associate callback with 232ff35dfedSBarry Smith . func - local Jacobian evaluation 233ff35dfedSBarry Smith - ctx - optional context for local Jacobian evaluation 234ff35dfedSBarry Smith 235ff35dfedSBarry Smith Level: beginner 236ff35dfedSBarry Smith 237ff35dfedSBarry Smith .seealso: DMSNESSetJacobian(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 238ff35dfedSBarry Smith @*/ 239d1e9a80fSBarry Smith PetscErrorCode DMSNESSetJacobianLocal(DM dm,PetscErrorCode (*func)(DM,Vec,Mat,Mat,void*),void *ctx) 240ff35dfedSBarry Smith { 241ff35dfedSBarry Smith PetscErrorCode ierr; 242942e3340SBarry Smith DMSNES sdm; 243942e3340SBarry Smith DMSNES_Local *dmlocalsnes; 244ff35dfedSBarry Smith 245ff35dfedSBarry Smith PetscFunctionBegin; 246ff35dfedSBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 247942e3340SBarry Smith ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 248ff35dfedSBarry Smith ierr = DMLocalSNESGetContext(dm,sdm,&dmlocalsnes);CHKERRQ(ierr); 2491aa26658SKarl Rupp 250ff35dfedSBarry Smith dmlocalsnes->jacobianlocal = func; 251ff35dfedSBarry Smith dmlocalsnes->jacobianlocalctx = ctx; 2521aa26658SKarl Rupp 253ff35dfedSBarry Smith ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMLocal,dmlocalsnes);CHKERRQ(ierr); 254ff35dfedSBarry Smith PetscFunctionReturn(0); 255ff35dfedSBarry Smith } 256ff35dfedSBarry Smith 257