1ff35dfedSBarry Smith #include <petsc-private/dmimpl.h> 2ff35dfedSBarry Smith #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/ 3ff35dfedSBarry Smith 4ff35dfedSBarry Smith typedef struct { 5ff35dfedSBarry Smith PetscErrorCode (*residuallocal)(DM,Vec,Vec,void*); 6ff35dfedSBarry Smith PetscErrorCode (*jacobianlocal)(DM,Vec,Mat,Mat,MatStructure*,void*); 7ff35dfedSBarry Smith void *residuallocalctx; 8ff35dfedSBarry Smith void *jacobianlocalctx; 9942e3340SBarry Smith } DMSNES_Local; 10ff35dfedSBarry Smith 11ff35dfedSBarry Smith #undef __FUNCT__ 12942e3340SBarry Smith #define __FUNCT__ "DMSNESDestroy_DMLocal" 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 22ff35dfedSBarry Smith #undef __FUNCT__ 23942e3340SBarry Smith #define __FUNCT__ "DMSNESDuplicate_DMLocal" 2422c6f798SBarry Smith static PetscErrorCode DMSNESDuplicate_DMLocal(DMSNES oldsdm,DMSNES sdm) 25ff35dfedSBarry Smith { 26ff35dfedSBarry Smith PetscErrorCode ierr; 27ff35dfedSBarry Smith 28ff35dfedSBarry Smith PetscFunctionBegin; 2922c6f798SBarry Smith ierr = PetscNewLog(sdm,DMSNES_Local,&sdm->data);CHKERRQ(ierr); 30ff35dfedSBarry Smith if (oldsdm->data) { 31942e3340SBarry Smith ierr = PetscMemcpy(sdm->data,oldsdm->data,sizeof(DMSNES_Local));CHKERRQ(ierr); 32ff35dfedSBarry Smith } 33ff35dfedSBarry Smith PetscFunctionReturn(0); 34ff35dfedSBarry Smith } 35ff35dfedSBarry Smith 36ff35dfedSBarry Smith #undef __FUNCT__ 37ff35dfedSBarry Smith #define __FUNCT__ "DMLocalSNESGetContext" 38942e3340SBarry Smith static PetscErrorCode DMLocalSNESGetContext(DM dm,DMSNES sdm,DMSNES_Local **dmlocalsnes) 39ff35dfedSBarry Smith { 40ff35dfedSBarry Smith PetscErrorCode ierr; 41ff35dfedSBarry Smith 42ff35dfedSBarry Smith PetscFunctionBegin; 430298fd71SBarry Smith *dmlocalsnes = NULL; 44ff35dfedSBarry Smith if (!sdm->data) { 45942e3340SBarry Smith ierr = PetscNewLog(dm,DMSNES_Local,&sdm->data);CHKERRQ(ierr); 461aa26658SKarl Rupp 4722c6f798SBarry Smith sdm->ops->destroy = DMSNESDestroy_DMLocal; 4822c6f798SBarry Smith sdm->ops->duplicate = DMSNESDuplicate_DMLocal; 49ff35dfedSBarry Smith } 50942e3340SBarry Smith *dmlocalsnes = (DMSNES_Local*)sdm->data; 51ff35dfedSBarry Smith PetscFunctionReturn(0); 52ff35dfedSBarry Smith } 53ff35dfedSBarry Smith 54ff35dfedSBarry Smith #undef __FUNCT__ 55ff35dfedSBarry Smith #define __FUNCT__ "SNESComputeFunction_DMLocal" 56b2270762SBarry Smith static PetscErrorCode SNESComputeFunction_DMLocal(SNES snes,Vec X,Vec F,void *ctx) 57ff35dfedSBarry Smith { 58ff35dfedSBarry Smith PetscErrorCode ierr; 59ff35dfedSBarry Smith DM dm; 60942e3340SBarry Smith DMSNES_Local *dmlocalsnes = (DMSNES_Local*)ctx; 61ff35dfedSBarry Smith Vec Xloc,Floc; 62ff35dfedSBarry Smith 63ff35dfedSBarry Smith PetscFunctionBegin; 64ff35dfedSBarry Smith PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 65ff35dfedSBarry Smith PetscValidHeaderSpecific(X,VEC_CLASSID,2); 66ff35dfedSBarry Smith PetscValidHeaderSpecific(F,VEC_CLASSID,3); 67ff35dfedSBarry Smith ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 68ff35dfedSBarry Smith ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 69ff35dfedSBarry Smith ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr); 70ff35dfedSBarry Smith ierr = VecZeroEntries(Xloc);CHKERRQ(ierr); 71ff35dfedSBarry Smith ierr = VecZeroEntries(Floc);CHKERRQ(ierr); 72ff35dfedSBarry Smith ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 73ff35dfedSBarry Smith ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 74ff35dfedSBarry Smith CHKMEMQ; 75ff35dfedSBarry Smith ierr = (*dmlocalsnes->residuallocal)(dm,Xloc,Floc,dmlocalsnes->residuallocalctx);CHKERRQ(ierr); 76ff35dfedSBarry Smith CHKMEMQ; 77ff35dfedSBarry Smith ierr = VecZeroEntries(F);CHKERRQ(ierr); 78ff35dfedSBarry Smith ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 79ff35dfedSBarry Smith ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 80ff35dfedSBarry Smith ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr); 81ff35dfedSBarry Smith ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 82ff35dfedSBarry Smith PetscFunctionReturn(0); 83ff35dfedSBarry Smith } 84ff35dfedSBarry Smith 85ff35dfedSBarry Smith #undef __FUNCT__ 86ff35dfedSBarry Smith #define __FUNCT__ "SNESComputeJacobian_DMLocal" 87ff35dfedSBarry Smith static PetscErrorCode SNESComputeJacobian_DMLocal(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *ctx) 88ff35dfedSBarry Smith { 89ff35dfedSBarry Smith PetscErrorCode ierr; 90ff35dfedSBarry Smith DM dm; 91942e3340SBarry Smith DMSNES_Local *dmlocalsnes = (DMSNES_Local*)ctx; 92ff35dfedSBarry Smith Vec Xloc; 93ff35dfedSBarry Smith 94ff35dfedSBarry Smith PetscFunctionBegin; 95ff35dfedSBarry Smith ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 96ff35dfedSBarry Smith if (dmlocalsnes->jacobianlocal) { 97ff35dfedSBarry Smith ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 98ff35dfedSBarry Smith ierr = VecZeroEntries(Xloc);CHKERRQ(ierr); 99ff35dfedSBarry Smith ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 100ff35dfedSBarry Smith ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 101ff35dfedSBarry Smith CHKMEMQ; 102ff35dfedSBarry Smith ierr = (*dmlocalsnes->jacobianlocal)(dm,Xloc,*A,*B,mstr,dmlocalsnes->jacobianlocalctx);CHKERRQ(ierr); 103ff35dfedSBarry Smith CHKMEMQ; 104ff35dfedSBarry Smith ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 105ff35dfedSBarry Smith } else { 106ff35dfedSBarry Smith MatFDColoring fdcoloring; 107ff35dfedSBarry Smith ierr = PetscObjectQuery((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject*)&fdcoloring);CHKERRQ(ierr); 108ff35dfedSBarry Smith if (!fdcoloring) { 109ff35dfedSBarry Smith ISColoring coloring; 110ff35dfedSBarry Smith 111ff35dfedSBarry Smith ierr = DMCreateColoring(dm,dm->coloringtype,dm->mattype,&coloring);CHKERRQ(ierr); 112ff35dfedSBarry Smith ierr = MatFDColoringCreate(*B,coloring,&fdcoloring);CHKERRQ(ierr); 113ff35dfedSBarry Smith ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr); 114ff35dfedSBarry Smith switch (dm->coloringtype) { 115ff35dfedSBarry Smith case IS_COLORING_GLOBAL: 116ff35dfedSBarry Smith ierr = MatFDColoringSetFunction(fdcoloring,(PetscErrorCode (*)(void))SNESComputeFunction_DMLocal,dmlocalsnes);CHKERRQ(ierr); 117ff35dfedSBarry Smith break; 118*ce94432eSBarry Smith default: SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"No support for coloring type '%s'",ISColoringTypes[dm->coloringtype]); 119ff35dfedSBarry Smith } 120ff35dfedSBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)fdcoloring,((PetscObject)dm)->prefix);CHKERRQ(ierr); 121ff35dfedSBarry Smith ierr = MatFDColoringSetFromOptions(fdcoloring);CHKERRQ(ierr); 122ff35dfedSBarry Smith ierr = PetscObjectCompose((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject)fdcoloring);CHKERRQ(ierr); 123ff35dfedSBarry Smith ierr = PetscObjectDereference((PetscObject)fdcoloring);CHKERRQ(ierr); 124ff35dfedSBarry Smith 125ff35dfedSBarry Smith /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call 126ff35dfedSBarry Smith * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the 127ff35dfedSBarry Smith * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually 128ff35dfedSBarry Smith * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be 129140e18c1SBarry Smith * taken when the PetscObjectList for the Vec inside MatFDColoring is destroyed. 130ff35dfedSBarry Smith */ 131ff35dfedSBarry Smith ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr); 132ff35dfedSBarry Smith } 133ff35dfedSBarry Smith *mstr = SAME_NONZERO_PATTERN; 134ff35dfedSBarry Smith ierr = MatFDColoringApply(*B,fdcoloring,X,mstr,snes);CHKERRQ(ierr); 135ff35dfedSBarry Smith } 136ff35dfedSBarry Smith /* This will be redundant if the user called both, but it's too common to forget. */ 137ff35dfedSBarry Smith if (*A != *B) { 138ff35dfedSBarry Smith ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 139ff35dfedSBarry Smith ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 140ff35dfedSBarry Smith } 141ff35dfedSBarry Smith PetscFunctionReturn(0); 142ff35dfedSBarry Smith } 143ff35dfedSBarry Smith 144ff35dfedSBarry Smith #undef __FUNCT__ 145ff35dfedSBarry Smith #define __FUNCT__ "DMSNESSetFunctionLocal" 146ff35dfedSBarry Smith /*@C 147ff35dfedSBarry Smith DMSNESSetFunctionLocal - set a local residual evaluation function. This function is called with local vector 148ff35dfedSBarry Smith containing the local vector information PLUS ghost point information. It should compute a result for all local 149ff35dfedSBarry Smith elements and DMSNES will automatically accumulate the overlapping values. 150ff35dfedSBarry Smith 151ff35dfedSBarry Smith Logically Collective 152ff35dfedSBarry Smith 153ff35dfedSBarry Smith Input Arguments: 154ff35dfedSBarry Smith + dm - DM to associate callback with 155ff35dfedSBarry Smith . func - local residual evaluation 156ff35dfedSBarry Smith - ctx - optional context for local residual evaluation 157ff35dfedSBarry Smith 158ff35dfedSBarry Smith Level: beginner 159ff35dfedSBarry Smith 160ff35dfedSBarry Smith .seealso: DMSNESSetFunction(), DMDASNESSetJacobianLocal(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 161ff35dfedSBarry Smith @*/ 162ff35dfedSBarry Smith PetscErrorCode DMSNESSetFunctionLocal(DM dm,PetscErrorCode (*func)(DM,Vec,Vec,void*),void *ctx) 163ff35dfedSBarry Smith { 164ff35dfedSBarry Smith PetscErrorCode ierr; 165942e3340SBarry Smith DMSNES sdm; 166942e3340SBarry Smith DMSNES_Local *dmlocalsnes; 167ff35dfedSBarry Smith 168ff35dfedSBarry Smith PetscFunctionBegin; 169ff35dfedSBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 170942e3340SBarry Smith ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 171ff35dfedSBarry Smith ierr = DMLocalSNESGetContext(dm,sdm,&dmlocalsnes);CHKERRQ(ierr); 1721aa26658SKarl Rupp 173ff35dfedSBarry Smith dmlocalsnes->residuallocal = func; 174ff35dfedSBarry Smith dmlocalsnes->residuallocalctx = ctx; 1751aa26658SKarl Rupp 176ff35dfedSBarry Smith ierr = DMSNESSetFunction(dm,SNESComputeFunction_DMLocal,dmlocalsnes);CHKERRQ(ierr); 17722c6f798SBarry Smith if (!sdm->ops->computejacobian) { /* Call us for the Jacobian too, can be overridden by the user. */ 178ff35dfedSBarry Smith ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMLocal,dmlocalsnes);CHKERRQ(ierr); 179ff35dfedSBarry Smith } 180ff35dfedSBarry Smith PetscFunctionReturn(0); 181ff35dfedSBarry Smith } 182ff35dfedSBarry Smith 183ff35dfedSBarry Smith #undef __FUNCT__ 184ff35dfedSBarry Smith #define __FUNCT__ "DMSNESSetJacobianLocal" 185ff35dfedSBarry Smith /*@C 186ff35dfedSBarry Smith DMSNESSetJacobianLocal - set a local Jacobian evaluation function 187ff35dfedSBarry Smith 188ff35dfedSBarry Smith Logically Collective 189ff35dfedSBarry Smith 190ff35dfedSBarry Smith Input Arguments: 191ff35dfedSBarry Smith + dm - DM to associate callback with 192ff35dfedSBarry Smith . func - local Jacobian evaluation 193ff35dfedSBarry Smith - ctx - optional context for local Jacobian evaluation 194ff35dfedSBarry Smith 195ff35dfedSBarry Smith Level: beginner 196ff35dfedSBarry Smith 197ff35dfedSBarry Smith .seealso: DMSNESSetJacobian(), DMDASNESSetJacobian(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 198ff35dfedSBarry Smith @*/ 199ff35dfedSBarry Smith PetscErrorCode DMSNESSetJacobianLocal(DM dm,PetscErrorCode (*func)(DM,Vec,Mat,Mat,MatStructure*,void*),void *ctx) 200ff35dfedSBarry Smith { 201ff35dfedSBarry Smith PetscErrorCode ierr; 202942e3340SBarry Smith DMSNES sdm; 203942e3340SBarry Smith DMSNES_Local *dmlocalsnes; 204ff35dfedSBarry Smith 205ff35dfedSBarry Smith PetscFunctionBegin; 206ff35dfedSBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 207942e3340SBarry Smith ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 208ff35dfedSBarry Smith ierr = DMLocalSNESGetContext(dm,sdm,&dmlocalsnes);CHKERRQ(ierr); 2091aa26658SKarl Rupp 210ff35dfedSBarry Smith dmlocalsnes->jacobianlocal = func; 211ff35dfedSBarry Smith dmlocalsnes->jacobianlocalctx = ctx; 2121aa26658SKarl Rupp 213ff35dfedSBarry Smith ierr = DMSNESSetJacobian(dm,SNESComputeJacobian_DMLocal,dmlocalsnes);CHKERRQ(ierr); 214ff35dfedSBarry Smith PetscFunctionReturn(0); 215ff35dfedSBarry Smith } 216ff35dfedSBarry Smith 217