1af0996ceSBarry Smith #include <petsc/private/dmimpl.h> 2af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 3d433e6cbSMatthew G. Knepley 4d433e6cbSMatthew G. Knepley typedef struct { 5*f2b97b97SMatthew G. Knepley PetscErrorCode (*boundarylocal)(DM,PetscReal,Vec,Vec,void*); 6d433e6cbSMatthew G. Knepley PetscErrorCode (*ifunctionlocal)(DM,PetscReal,Vec,Vec,Vec,void*); 7d433e6cbSMatthew G. Knepley PetscErrorCode (*ijacobianlocal)(DM,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*); 8d433e6cbSMatthew G. Knepley PetscErrorCode (*rhsfunctionlocal)(DM,PetscReal,Vec,Vec,void*); 9*f2b97b97SMatthew G. Knepley void *boundarylocalctx; 10d433e6cbSMatthew G. Knepley void *ifunctionlocalctx; 11d433e6cbSMatthew G. Knepley void *ijacobianlocalctx; 12d433e6cbSMatthew G. Knepley void *rhsfunctionlocalctx; 13d433e6cbSMatthew G. Knepley } DMTS_Local; 14d433e6cbSMatthew G. Knepley 15d433e6cbSMatthew G. Knepley #undef __FUNCT__ 16d433e6cbSMatthew G. Knepley #define __FUNCT__ "DMTSDestroy_DMLocal" 17d433e6cbSMatthew G. Knepley static PetscErrorCode DMTSDestroy_DMLocal(DMTS tdm) 18d433e6cbSMatthew G. Knepley { 19d433e6cbSMatthew G. Knepley PetscErrorCode ierr; 20d433e6cbSMatthew G. Knepley 21d433e6cbSMatthew G. Knepley PetscFunctionBegin; 22d433e6cbSMatthew G. Knepley ierr = PetscFree(tdm->data);CHKERRQ(ierr); 23d433e6cbSMatthew G. Knepley PetscFunctionReturn(0); 24d433e6cbSMatthew G. Knepley } 25d433e6cbSMatthew G. Knepley 26d433e6cbSMatthew G. Knepley #undef __FUNCT__ 27d433e6cbSMatthew G. Knepley #define __FUNCT__ "DMTSDuplicate_DMLocal" 28d433e6cbSMatthew G. Knepley static PetscErrorCode DMTSDuplicate_DMLocal(DMTS oldtdm, DMTS tdm) 29d433e6cbSMatthew G. Knepley { 30d433e6cbSMatthew G. Knepley PetscErrorCode ierr; 31d433e6cbSMatthew G. Knepley 32d433e6cbSMatthew G. Knepley PetscFunctionBegin; 33d433e6cbSMatthew G. Knepley ierr = PetscNewLog(tdm, (DMTS_Local **) &tdm->data);CHKERRQ(ierr); 34d433e6cbSMatthew G. Knepley if (oldtdm->data) {ierr = PetscMemcpy(tdm->data, oldtdm->data, sizeof(DMTS_Local));CHKERRQ(ierr);} 35d433e6cbSMatthew G. Knepley PetscFunctionReturn(0); 36d433e6cbSMatthew G. Knepley } 37d433e6cbSMatthew G. Knepley 38d433e6cbSMatthew G. Knepley #undef __FUNCT__ 39d433e6cbSMatthew G. Knepley #define __FUNCT__ "DMLocalTSGetContext" 40d433e6cbSMatthew G. Knepley static PetscErrorCode DMLocalTSGetContext(DM dm, DMTS tdm, DMTS_Local **dmlocalts) 41d433e6cbSMatthew G. Knepley { 42d433e6cbSMatthew G. Knepley PetscErrorCode ierr; 43d433e6cbSMatthew G. Knepley 44d433e6cbSMatthew G. Knepley PetscFunctionBegin; 45d433e6cbSMatthew G. Knepley *dmlocalts = NULL; 46d433e6cbSMatthew G. Knepley if (!tdm->data) { 47d433e6cbSMatthew G. Knepley ierr = PetscNewLog(dm, (DMTS_Local **) &tdm->data);CHKERRQ(ierr); 48d433e6cbSMatthew G. Knepley 49d433e6cbSMatthew G. Knepley tdm->ops->destroy = DMTSDestroy_DMLocal; 50d433e6cbSMatthew G. Knepley tdm->ops->duplicate = DMTSDuplicate_DMLocal; 51d433e6cbSMatthew G. Knepley } 52d433e6cbSMatthew G. Knepley *dmlocalts = (DMTS_Local *) tdm->data; 53d433e6cbSMatthew G. Knepley PetscFunctionReturn(0); 54d433e6cbSMatthew G. Knepley } 55d433e6cbSMatthew G. Knepley 56d433e6cbSMatthew G. Knepley #undef __FUNCT__ 57d433e6cbSMatthew G. Knepley #define __FUNCT__ "TSComputeIFunction_DMLocal" 58d433e6cbSMatthew G. Knepley static PetscErrorCode TSComputeIFunction_DMLocal(TS ts, PetscReal time, Vec X, Vec X_t, Vec F, void *ctx) 59d433e6cbSMatthew G. Knepley { 60d433e6cbSMatthew G. Knepley DM dm; 61d433e6cbSMatthew G. Knepley Vec locX, locX_t, locF; 62d433e6cbSMatthew G. Knepley DMTS_Local *dmlocalts = (DMTS_Local *) ctx; 63d433e6cbSMatthew G. Knepley PetscErrorCode ierr; 64d433e6cbSMatthew G. Knepley 65d433e6cbSMatthew G. Knepley PetscFunctionBegin; 66d433e6cbSMatthew G. Knepley PetscValidHeaderSpecific(ts,TS_CLASSID,1); 67d433e6cbSMatthew G. Knepley PetscValidHeaderSpecific(X,VEC_CLASSID,3); 68d433e6cbSMatthew G. Knepley PetscValidHeaderSpecific(X_t,VEC_CLASSID,4); 69d433e6cbSMatthew G. Knepley PetscValidHeaderSpecific(F,VEC_CLASSID,5); 70d433e6cbSMatthew G. Knepley ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 71d433e6cbSMatthew G. Knepley ierr = DMGetLocalVector(dm, &locX);CHKERRQ(ierr); 72d433e6cbSMatthew G. Knepley ierr = DMGetLocalVector(dm, &locX_t);CHKERRQ(ierr); 73d433e6cbSMatthew G. Knepley ierr = DMGetLocalVector(dm, &locF);CHKERRQ(ierr); 74d433e6cbSMatthew G. Knepley ierr = VecZeroEntries(locX);CHKERRQ(ierr); 75d433e6cbSMatthew G. Knepley ierr = VecZeroEntries(locX_t);CHKERRQ(ierr); 76*f2b97b97SMatthew G. Knepley if (dmlocalts->boundarylocal) {ierr = (*dmlocalts->boundarylocal)(dm, time, locX, locX_t,dmlocalts->boundarylocalctx);CHKERRQ(ierr);} 77d433e6cbSMatthew G. Knepley ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX);CHKERRQ(ierr); 78d433e6cbSMatthew G. Knepley ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX);CHKERRQ(ierr); 79d433e6cbSMatthew G. Knepley ierr = DMGlobalToLocalBegin(dm, X_t, INSERT_VALUES, locX_t);CHKERRQ(ierr); 80d433e6cbSMatthew G. Knepley ierr = DMGlobalToLocalEnd(dm, X_t, INSERT_VALUES, locX_t);CHKERRQ(ierr); 81d433e6cbSMatthew G. Knepley ierr = VecZeroEntries(locF);CHKERRQ(ierr); 82d433e6cbSMatthew G. Knepley CHKMEMQ; 83d433e6cbSMatthew G. Knepley ierr = (*dmlocalts->ifunctionlocal)(dm, time, locX, locX_t, locF, dmlocalts->ifunctionlocalctx);CHKERRQ(ierr); 84d433e6cbSMatthew G. Knepley CHKMEMQ; 85d433e6cbSMatthew G. Knepley ierr = VecZeroEntries(F);CHKERRQ(ierr); 86d433e6cbSMatthew G. Knepley ierr = DMLocalToGlobalBegin(dm, locF, ADD_VALUES, F);CHKERRQ(ierr); 87d433e6cbSMatthew G. Knepley ierr = DMLocalToGlobalEnd(dm, locF, ADD_VALUES, F);CHKERRQ(ierr); 88d433e6cbSMatthew G. Knepley ierr = DMRestoreLocalVector(dm, &locX);CHKERRQ(ierr); 89d433e6cbSMatthew G. Knepley ierr = DMRestoreLocalVector(dm, &locX_t);CHKERRQ(ierr); 90d433e6cbSMatthew G. Knepley ierr = DMRestoreLocalVector(dm, &locF);CHKERRQ(ierr); 91d433e6cbSMatthew G. Knepley PetscFunctionReturn(0); 92d433e6cbSMatthew G. Knepley } 93d433e6cbSMatthew G. Knepley 94d433e6cbSMatthew G. Knepley #undef __FUNCT__ 95d433e6cbSMatthew G. Knepley #define __FUNCT__ "TSComputeRHSFunction_DMLocal" 96d433e6cbSMatthew G. Knepley static PetscErrorCode TSComputeRHSFunction_DMLocal(TS ts, PetscReal time, Vec X, Vec F, void *ctx) 97d433e6cbSMatthew G. Knepley { 98d433e6cbSMatthew G. Knepley DM dm; 99d433e6cbSMatthew G. Knepley Vec locX; 100d433e6cbSMatthew G. Knepley DMTS_Local *dmlocalts = (DMTS_Local *) ctx; 101d433e6cbSMatthew G. Knepley PetscErrorCode ierr; 102d433e6cbSMatthew G. Knepley 103d433e6cbSMatthew G. Knepley PetscFunctionBegin; 104d433e6cbSMatthew G. Knepley PetscValidHeaderSpecific(ts,TS_CLASSID,1); 105d433e6cbSMatthew G. Knepley PetscValidHeaderSpecific(X,VEC_CLASSID,3); 106d433e6cbSMatthew G. Knepley PetscValidHeaderSpecific(F,VEC_CLASSID,5); 107d433e6cbSMatthew G. Knepley ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 108d433e6cbSMatthew G. Knepley ierr = DMGetLocalVector(dm, &locX);CHKERRQ(ierr); 109d433e6cbSMatthew G. Knepley ierr = VecZeroEntries(locX);CHKERRQ(ierr); 110*f2b97b97SMatthew G. Knepley if (dmlocalts->boundarylocal) {ierr = (*dmlocalts->boundarylocal)(dm,time,locX,NULL,dmlocalts->boundarylocalctx);CHKERRQ(ierr);} 111d433e6cbSMatthew G. Knepley ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX);CHKERRQ(ierr); 112d433e6cbSMatthew G. Knepley ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX);CHKERRQ(ierr); 113d433e6cbSMatthew G. Knepley ierr = VecZeroEntries(F);CHKERRQ(ierr); 114d433e6cbSMatthew G. Knepley CHKMEMQ; 115d433e6cbSMatthew G. Knepley ierr = (*dmlocalts->rhsfunctionlocal)(dm, time, locX, F, dmlocalts->rhsfunctionlocalctx);CHKERRQ(ierr); 116d433e6cbSMatthew G. Knepley CHKMEMQ; 117d433e6cbSMatthew G. Knepley ierr = DMRestoreLocalVector(dm, &locX);CHKERRQ(ierr); 118d433e6cbSMatthew G. Knepley PetscFunctionReturn(0); 119d433e6cbSMatthew G. Knepley } 120d433e6cbSMatthew G. Knepley 121d433e6cbSMatthew G. Knepley #undef __FUNCT__ 122d433e6cbSMatthew G. Knepley #define __FUNCT__ "TSComputeIJacobian_DMLocal" 123d433e6cbSMatthew G. Knepley static PetscErrorCode TSComputeIJacobian_DMLocal(TS ts, PetscReal time, Vec X, Vec X_t, PetscReal a, Mat A, Mat B, void *ctx) 124d433e6cbSMatthew G. Knepley { 125d433e6cbSMatthew G. Knepley DM dm; 126d433e6cbSMatthew G. Knepley Vec locX, locX_t; 127d433e6cbSMatthew G. Knepley DMTS_Local *dmlocalts = (DMTS_Local *) ctx; 128d433e6cbSMatthew G. Knepley PetscErrorCode ierr; 129d433e6cbSMatthew G. Knepley 130d433e6cbSMatthew G. Knepley PetscFunctionBegin; 131d433e6cbSMatthew G. Knepley ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 132d433e6cbSMatthew G. Knepley if (dmlocalts->ijacobianlocal) { 133d433e6cbSMatthew G. Knepley ierr = DMGetLocalVector(dm, &locX);CHKERRQ(ierr); 134d433e6cbSMatthew G. Knepley ierr = DMGetLocalVector(dm, &locX_t);CHKERRQ(ierr); 135d433e6cbSMatthew G. Knepley ierr = VecZeroEntries(locX);CHKERRQ(ierr); 136*f2b97b97SMatthew G. Knepley if (dmlocalts->boundarylocal) {ierr = (*dmlocalts->boundarylocal)(dm,time,locX,locX_t,dmlocalts->boundarylocalctx);CHKERRQ(ierr);} 137d433e6cbSMatthew G. Knepley ierr = VecZeroEntries(locX_t);CHKERRQ(ierr); 138d433e6cbSMatthew G. Knepley ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX);CHKERRQ(ierr); 139d433e6cbSMatthew G. Knepley ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX);CHKERRQ(ierr); 140d433e6cbSMatthew G. Knepley ierr = DMGlobalToLocalBegin(dm, X_t, INSERT_VALUES, locX_t);CHKERRQ(ierr); 141d433e6cbSMatthew G. Knepley ierr = DMGlobalToLocalEnd(dm, X_t, INSERT_VALUES, locX_t);CHKERRQ(ierr); 142d433e6cbSMatthew G. Knepley CHKMEMQ; 143d433e6cbSMatthew G. Knepley ierr = (*dmlocalts->ijacobianlocal)(dm, time, locX, locX_t, a, A, B, dmlocalts->ijacobianlocalctx);CHKERRQ(ierr); 144d433e6cbSMatthew G. Knepley CHKMEMQ; 145d433e6cbSMatthew G. Knepley ierr = DMRestoreLocalVector(dm, &locX);CHKERRQ(ierr); 146d433e6cbSMatthew G. Knepley ierr = DMRestoreLocalVector(dm, &locX_t);CHKERRQ(ierr); 147d433e6cbSMatthew G. Knepley } else { 148d433e6cbSMatthew G. Knepley MatFDColoring fdcoloring; 149d433e6cbSMatthew G. Knepley ierr = PetscObjectQuery((PetscObject) dm, "DMDASNES_FDCOLORING", (PetscObject *) &fdcoloring);CHKERRQ(ierr); 150d433e6cbSMatthew G. Knepley if (!fdcoloring) { 151d433e6cbSMatthew G. Knepley ISColoring coloring; 152d433e6cbSMatthew G. Knepley 153d433e6cbSMatthew G. Knepley ierr = DMCreateColoring(dm, dm->coloringtype, &coloring);CHKERRQ(ierr); 154d433e6cbSMatthew G. Knepley ierr = MatFDColoringCreate(B, coloring, &fdcoloring);CHKERRQ(ierr); 155d433e6cbSMatthew G. Knepley ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr); 156d433e6cbSMatthew G. Knepley switch (dm->coloringtype) { 157d433e6cbSMatthew G. Knepley case IS_COLORING_GLOBAL: 158d433e6cbSMatthew G. Knepley ierr = MatFDColoringSetFunction(fdcoloring, (PetscErrorCode (*)(void)) TSComputeIFunction_DMLocal, dmlocalts);CHKERRQ(ierr); 159d433e6cbSMatthew G. Knepley break; 160d433e6cbSMatthew G. Knepley default: SETERRQ1(PetscObjectComm((PetscObject) ts), PETSC_ERR_SUP, "No support for coloring type '%s'", ISColoringTypes[dm->coloringtype]); 161d433e6cbSMatthew G. Knepley } 162d433e6cbSMatthew G. Knepley ierr = PetscObjectSetOptionsPrefix((PetscObject) fdcoloring, ((PetscObject) dm)->prefix);CHKERRQ(ierr); 163d433e6cbSMatthew G. Knepley ierr = MatFDColoringSetFromOptions(fdcoloring);CHKERRQ(ierr); 164d433e6cbSMatthew G. Knepley ierr = MatFDColoringSetUp(B, coloring, fdcoloring);CHKERRQ(ierr); 165d433e6cbSMatthew G. Knepley ierr = PetscObjectCompose((PetscObject) dm, "DMDASNES_FDCOLORING", (PetscObject) fdcoloring);CHKERRQ(ierr); 166d433e6cbSMatthew G. Knepley ierr = PetscObjectDereference((PetscObject) fdcoloring);CHKERRQ(ierr); 167d433e6cbSMatthew G. Knepley 168d433e6cbSMatthew G. Knepley /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call 169d433e6cbSMatthew G. Knepley * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the 170d433e6cbSMatthew G. Knepley * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually 171d433e6cbSMatthew G. Knepley * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be 172d433e6cbSMatthew G. Knepley * taken when the PetscObjectList for the Vec inside MatFDColoring is destroyed. 173d433e6cbSMatthew G. Knepley */ 174d433e6cbSMatthew G. Knepley ierr = PetscObjectDereference((PetscObject) dm);CHKERRQ(ierr); 175d433e6cbSMatthew G. Knepley } 176d433e6cbSMatthew G. Knepley ierr = MatFDColoringApply(B, fdcoloring, X, ts);CHKERRQ(ierr); 177d433e6cbSMatthew G. Knepley } 178d433e6cbSMatthew G. Knepley /* This will be redundant if the user called both, but it's too common to forget. */ 179d433e6cbSMatthew G. Knepley if (A != B) { 180d433e6cbSMatthew G. Knepley ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 181d433e6cbSMatthew G. Knepley ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 182d433e6cbSMatthew G. Knepley } 183d433e6cbSMatthew G. Knepley PetscFunctionReturn(0); 184d433e6cbSMatthew G. Knepley } 185d433e6cbSMatthew G. Knepley 186d433e6cbSMatthew G. Knepley #undef __FUNCT__ 1877eb36faeSMatthew G. Knepley #define __FUNCT__ "DMTSSetBoundaryLocal" 1887eb36faeSMatthew G. Knepley /*@C 189*f2b97b97SMatthew G. Knepley DMTSSetBoundaryLocal - set the function for essential boundary data for a local implicit function evaluation. 190a40652d4SToby Isaac It should set the essential boundary data for the local portion of the solution X, as well its time derivative X_t (if it is not NULL). 191a40652d4SToby Isaac Vectors are initialized to zero before this function, so it is only needed for non homogeneous data. 192a40652d4SToby Isaac 193a40652d4SToby Isaac Note that this function is somewhat optional: boundary data could potentially be inserted by a function passed to 194a40652d4SToby Isaac DMTSSetIFunctionLocal(). The use case for this function is for discretizations with constraints (see 195a40652d4SToby Isaac DMGetDefaultConstraints()): this function inserts boundary values before constraint interpolation. 1967eb36faeSMatthew G. Knepley 1977eb36faeSMatthew G. Knepley Logically Collective 1987eb36faeSMatthew G. Knepley 1997eb36faeSMatthew G. Knepley Input Arguments: 2007eb36faeSMatthew G. Knepley + dm - DM to associate callback with 201a40652d4SToby Isaac . func - local function evaluation 202a40652d4SToby Isaac - ctx - context for function evaluation 2037eb36faeSMatthew G. Knepley 2047eb36faeSMatthew G. Knepley Level: intermediate 2057eb36faeSMatthew G. Knepley 206a40652d4SToby Isaac .seealso: DMTSSetIFunction(), DMTSSetIJacobianLocal() 2077eb36faeSMatthew G. Knepley @*/ 208*f2b97b97SMatthew G. Knepley PetscErrorCode DMTSSetBoundaryLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, void *), void *ctx) 2097eb36faeSMatthew G. Knepley { 2107eb36faeSMatthew G. Knepley DMTS tdm; 2117eb36faeSMatthew G. Knepley DMTS_Local *dmlocalts; 2127eb36faeSMatthew G. Knepley PetscErrorCode ierr; 2137eb36faeSMatthew G. Knepley 2147eb36faeSMatthew G. Knepley PetscFunctionBegin; 2157eb36faeSMatthew G. Knepley PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2167eb36faeSMatthew G. Knepley ierr = DMGetDMTSWrite(dm, &tdm);CHKERRQ(ierr); 2177eb36faeSMatthew G. Knepley ierr = DMLocalTSGetContext(dm, tdm, &dmlocalts);CHKERRQ(ierr); 2187eb36faeSMatthew G. Knepley 2197eb36faeSMatthew G. Knepley dmlocalts->boundarylocal = func; 2207eb36faeSMatthew G. Knepley dmlocalts->boundarylocalctx = ctx; 221a40652d4SToby Isaac 2227eb36faeSMatthew G. Knepley PetscFunctionReturn(0); 2237eb36faeSMatthew G. Knepley } 2247eb36faeSMatthew G. Knepley 2257eb36faeSMatthew G. Knepley #undef __FUNCT__ 226d433e6cbSMatthew G. Knepley #define __FUNCT__ "DMTSSetIFunctionLocal" 227d433e6cbSMatthew G. Knepley /*@C 228d433e6cbSMatthew G. Knepley DMTSSetIFunctionLocal - set a local implicit function evaluation function. This function is called with local vector 229d433e6cbSMatthew G. Knepley containing the local vector information PLUS ghost point information. It should compute a result for all local 230d433e6cbSMatthew G. Knepley elements and DMTS will automatically accumulate the overlapping values. 231d433e6cbSMatthew G. Knepley 232d433e6cbSMatthew G. Knepley Logically Collective 233d433e6cbSMatthew G. Knepley 234d433e6cbSMatthew G. Knepley Input Arguments: 235d433e6cbSMatthew G. Knepley + dm - DM to associate callback with 236d433e6cbSMatthew G. Knepley . func - local function evaluation 237d433e6cbSMatthew G. Knepley - ctx - context for function evaluation 238d433e6cbSMatthew G. Knepley 239d433e6cbSMatthew G. Knepley Level: beginner 240d433e6cbSMatthew G. Knepley 241d433e6cbSMatthew G. Knepley .seealso: DMTSSetIFunction(), DMTSSetIJacobianLocal() 242d433e6cbSMatthew G. Knepley @*/ 243d433e6cbSMatthew G. Knepley PetscErrorCode DMTSSetIFunctionLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, Vec, void *), void *ctx) 244d433e6cbSMatthew G. Knepley { 245d433e6cbSMatthew G. Knepley DMTS tdm; 246d433e6cbSMatthew G. Knepley DMTS_Local *dmlocalts; 247d433e6cbSMatthew G. Knepley PetscErrorCode ierr; 248d433e6cbSMatthew G. Knepley 249d433e6cbSMatthew G. Knepley PetscFunctionBegin; 250d433e6cbSMatthew G. Knepley PetscValidHeaderSpecific(dm,DM_CLASSID,1); 251d433e6cbSMatthew G. Knepley ierr = DMGetDMTSWrite(dm, &tdm);CHKERRQ(ierr); 252d433e6cbSMatthew G. Knepley ierr = DMLocalTSGetContext(dm, tdm, &dmlocalts);CHKERRQ(ierr); 253d433e6cbSMatthew G. Knepley 254d433e6cbSMatthew G. Knepley dmlocalts->ifunctionlocal = func; 255d433e6cbSMatthew G. Knepley dmlocalts->ifunctionlocalctx = ctx; 256d433e6cbSMatthew G. Knepley 257d433e6cbSMatthew G. Knepley ierr = DMTSSetIFunction(dm, TSComputeIFunction_DMLocal, dmlocalts);CHKERRQ(ierr); 258d433e6cbSMatthew G. Knepley if (!tdm->ops->ijacobian) { /* Call us for the Jacobian too, can be overridden by the user. */ 259d433e6cbSMatthew G. Knepley ierr = DMTSSetIJacobian(dm, TSComputeIJacobian_DMLocal, dmlocalts);CHKERRQ(ierr); 260d433e6cbSMatthew G. Knepley } 261d433e6cbSMatthew G. Knepley PetscFunctionReturn(0); 262d433e6cbSMatthew G. Knepley } 263d433e6cbSMatthew G. Knepley 264d433e6cbSMatthew G. Knepley #undef __FUNCT__ 265d433e6cbSMatthew G. Knepley #define __FUNCT__ "DMTSSetIJacobianLocal" 266d433e6cbSMatthew G. Knepley /*@C 267d433e6cbSMatthew G. Knepley DMTSSetIJacobianLocal - set a local Jacobian evaluation function 268d433e6cbSMatthew G. Knepley 269d433e6cbSMatthew G. Knepley Logically Collective 270d433e6cbSMatthew G. Knepley 271d433e6cbSMatthew G. Knepley Input Arguments: 272d433e6cbSMatthew G. Knepley + dm - DM to associate callback with 273d433e6cbSMatthew G. Knepley . func - local Jacobian evaluation 274d433e6cbSMatthew G. Knepley - ctx - optional context for local Jacobian evaluation 275d433e6cbSMatthew G. Knepley 276d433e6cbSMatthew G. Knepley Level: beginner 277d433e6cbSMatthew G. Knepley 278d433e6cbSMatthew G. Knepley .seealso: DMTSSetIFunctionLocal(), DMTSSetIJacobian(), DMTSSetIFunction() 279d433e6cbSMatthew G. Knepley @*/ 280d433e6cbSMatthew G. Knepley PetscErrorCode DMTSSetIJacobianLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void *ctx) 281d433e6cbSMatthew G. Knepley { 282d433e6cbSMatthew G. Knepley DMTS tdm; 283d433e6cbSMatthew G. Knepley DMTS_Local *dmlocalts; 284d433e6cbSMatthew G. Knepley PetscErrorCode ierr; 285d433e6cbSMatthew G. Knepley 286d433e6cbSMatthew G. Knepley PetscFunctionBegin; 287d433e6cbSMatthew G. Knepley PetscValidHeaderSpecific(dm,DM_CLASSID,1); 288d433e6cbSMatthew G. Knepley ierr = DMGetDMTSWrite(dm, &tdm);CHKERRQ(ierr); 289d433e6cbSMatthew G. Knepley ierr = DMLocalTSGetContext(dm, tdm, &dmlocalts);CHKERRQ(ierr); 290d433e6cbSMatthew G. Knepley 291d433e6cbSMatthew G. Knepley dmlocalts->ijacobianlocal = func; 292d433e6cbSMatthew G. Knepley dmlocalts->ijacobianlocalctx = ctx; 293d433e6cbSMatthew G. Knepley 294d433e6cbSMatthew G. Knepley ierr = DMTSSetIJacobian(dm, TSComputeIJacobian_DMLocal, dmlocalts);CHKERRQ(ierr); 295d433e6cbSMatthew G. Knepley PetscFunctionReturn(0); 296d433e6cbSMatthew G. Knepley } 297d433e6cbSMatthew G. Knepley 298d433e6cbSMatthew G. Knepley #undef __FUNCT__ 299d433e6cbSMatthew G. Knepley #define __FUNCT__ "DMTSSetRHSFunctionLocal" 300d433e6cbSMatthew G. Knepley /*@C 301d433e6cbSMatthew G. Knepley DMTSSetRHSFunctionLocal - set a local rhs function evaluation function. This function is called with local vector 302d433e6cbSMatthew G. Knepley containing the local vector information PLUS ghost point information. It should compute a result for all local 303d433e6cbSMatthew G. Knepley elements and DMTS will automatically accumulate the overlapping values. 304d433e6cbSMatthew G. Knepley 305d433e6cbSMatthew G. Knepley Logically Collective 306d433e6cbSMatthew G. Knepley 307d433e6cbSMatthew G. Knepley Input Arguments: 308d433e6cbSMatthew G. Knepley + dm - DM to associate callback with 309d433e6cbSMatthew G. Knepley . func - local function evaluation 310d433e6cbSMatthew G. Knepley - ctx - context for function evaluation 311d433e6cbSMatthew G. Knepley 312d433e6cbSMatthew G. Knepley Level: beginner 313d433e6cbSMatthew G. Knepley 314d433e6cbSMatthew G. Knepley .seealso: DMTSSetRHSFunction(), DMTSSetIFunction(), DMTSSetIJacobianLocal() 315d433e6cbSMatthew G. Knepley @*/ 316d433e6cbSMatthew G. Knepley PetscErrorCode DMTSSetRHSFunctionLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, void *), void *ctx) 317d433e6cbSMatthew G. Knepley { 318d433e6cbSMatthew G. Knepley DMTS tdm; 319d433e6cbSMatthew G. Knepley DMTS_Local *dmlocalts; 320d433e6cbSMatthew G. Knepley PetscErrorCode ierr; 321d433e6cbSMatthew G. Knepley 322d433e6cbSMatthew G. Knepley PetscFunctionBegin; 323d433e6cbSMatthew G. Knepley PetscValidHeaderSpecific(dm,DM_CLASSID,1); 324d433e6cbSMatthew G. Knepley ierr = DMGetDMTSWrite(dm, &tdm);CHKERRQ(ierr); 325d433e6cbSMatthew G. Knepley ierr = DMLocalTSGetContext(dm, tdm, &dmlocalts);CHKERRQ(ierr); 326d433e6cbSMatthew G. Knepley 327d433e6cbSMatthew G. Knepley dmlocalts->rhsfunctionlocal = func; 328d433e6cbSMatthew G. Knepley dmlocalts->rhsfunctionlocalctx = ctx; 329d433e6cbSMatthew G. Knepley 330d433e6cbSMatthew G. Knepley ierr = DMTSSetRHSFunction(dm, TSComputeRHSFunction_DMLocal, dmlocalts);CHKERRQ(ierr); 331d433e6cbSMatthew G. Knepley PetscFunctionReturn(0); 332d433e6cbSMatthew G. Knepley } 333