140be0ff1SMatthew G. Knepley /* 240be0ff1SMatthew G. Knepley Code for timestepping with discrete gradient integrators 340be0ff1SMatthew G. Knepley */ 440be0ff1SMatthew G. Knepley #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 540be0ff1SMatthew G. Knepley #include <petscdm.h> 640be0ff1SMatthew G. Knepley 740be0ff1SMatthew G. Knepley PetscBool DGCite = PETSC_FALSE; 840be0ff1SMatthew G. Knepley const char DGCitation[] = "@article{Gonzalez1996,\n" 940be0ff1SMatthew G. Knepley " title = {Time integration and discrete Hamiltonian systems},\n" 1040be0ff1SMatthew G. Knepley " author = {Oscar Gonzalez},\n" 1140be0ff1SMatthew G. Knepley " journal = {Journal of Nonlinear Science},\n" 1240be0ff1SMatthew G. Knepley " volume = {6},\n" 1340be0ff1SMatthew G. Knepley " pages = {449--467},\n" 1440be0ff1SMatthew G. Knepley " doi = {10.1007/978-1-4612-1246-1_10},\n" 1540be0ff1SMatthew G. Knepley " year = {1996}\n}\n"; 1640be0ff1SMatthew G. Knepley 17f940b0e3Sdanofinn const char *DGTypes[] = {"gonzalez", "average", "none", "TSDGType", "DG_", NULL}; 18f940b0e3Sdanofinn 1940be0ff1SMatthew G. Knepley typedef struct { 2040be0ff1SMatthew G. Knepley PetscReal stage_time; 2140be0ff1SMatthew G. Knepley Vec X0, X, Xdot; 22db4653c2SDaniel Finn void *funcCtx; 23f940b0e3Sdanofinn TSDGType discgrad; /* Type of electrostatic model */ 2440be0ff1SMatthew G. Knepley PetscErrorCode (*Sfunc)(TS, PetscReal, Vec, Mat, void *); 2540be0ff1SMatthew G. Knepley PetscErrorCode (*Ffunc)(TS, PetscReal, Vec, PetscScalar *, void *); 2640be0ff1SMatthew G. Knepley PetscErrorCode (*Gfunc)(TS, PetscReal, Vec, Vec, void *); 2740be0ff1SMatthew G. Knepley } TS_DiscGrad; 2840be0ff1SMatthew G. Knepley 29d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDiscGradGetX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot) 30d71ae5a4SJacob Faibussowitsch { 3140be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *)ts->data; 3240be0ff1SMatthew G. Knepley 3340be0ff1SMatthew G. Knepley PetscFunctionBegin; 3440be0ff1SMatthew G. Knepley if (X0) { 359566063dSJacob Faibussowitsch if (dm && dm != ts->dm) PetscCall(DMGetNamedGlobalVector(dm, "TSDiscGrad_X0", X0)); 36ad540459SPierre Jolivet else *X0 = ts->vec_sol; 3740be0ff1SMatthew G. Knepley } 3840be0ff1SMatthew G. Knepley if (Xdot) { 399566063dSJacob Faibussowitsch if (dm && dm != ts->dm) PetscCall(DMGetNamedGlobalVector(dm, "TSDiscGrad_Xdot", Xdot)); 40ad540459SPierre Jolivet else *Xdot = dg->Xdot; 4140be0ff1SMatthew G. Knepley } 423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4340be0ff1SMatthew G. Knepley } 4440be0ff1SMatthew G. Knepley 45d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDiscGradRestoreX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot) 46d71ae5a4SJacob Faibussowitsch { 4740be0ff1SMatthew G. Knepley PetscFunctionBegin; 4840be0ff1SMatthew G. Knepley if (X0) { 499566063dSJacob Faibussowitsch if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSDiscGrad_X0", X0)); 5040be0ff1SMatthew G. Knepley } 5140be0ff1SMatthew G. Knepley if (Xdot) { 529566063dSJacob Faibussowitsch if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSDiscGrad_Xdot", Xdot)); 5340be0ff1SMatthew G. Knepley } 543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5540be0ff1SMatthew G. Knepley } 5640be0ff1SMatthew G. Knepley 57*2a8381b2SBarry Smith static PetscErrorCode DMCoarsenHook_TSDiscGrad(DM fine, DM coarse, PetscCtx ctx) 58d71ae5a4SJacob Faibussowitsch { 5940be0ff1SMatthew G. Knepley PetscFunctionBegin; 603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6140be0ff1SMatthew G. Knepley } 6240be0ff1SMatthew G. Knepley 63*2a8381b2SBarry Smith static PetscErrorCode DMRestrictHook_TSDiscGrad(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, PetscCtx ctx) 64d71ae5a4SJacob Faibussowitsch { 6540be0ff1SMatthew G. Knepley TS ts = (TS)ctx; 6640be0ff1SMatthew G. Knepley Vec X0, Xdot, X0_c, Xdot_c; 6740be0ff1SMatthew G. Knepley 6840be0ff1SMatthew G. Knepley PetscFunctionBegin; 699566063dSJacob Faibussowitsch PetscCall(TSDiscGradGetX0AndXdot(ts, fine, &X0, &Xdot)); 709566063dSJacob Faibussowitsch PetscCall(TSDiscGradGetX0AndXdot(ts, coarse, &X0_c, &Xdot_c)); 719566063dSJacob Faibussowitsch PetscCall(MatRestrict(restrct, X0, X0_c)); 729566063dSJacob Faibussowitsch PetscCall(MatRestrict(restrct, Xdot, Xdot_c)); 739566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(X0_c, rscale, X0_c)); 749566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Xdot_c, rscale, Xdot_c)); 759566063dSJacob Faibussowitsch PetscCall(TSDiscGradRestoreX0AndXdot(ts, fine, &X0, &Xdot)); 769566063dSJacob Faibussowitsch PetscCall(TSDiscGradRestoreX0AndXdot(ts, coarse, &X0_c, &Xdot_c)); 773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7840be0ff1SMatthew G. Knepley } 7940be0ff1SMatthew G. Knepley 80*2a8381b2SBarry Smith static PetscErrorCode DMSubDomainHook_TSDiscGrad(DM dm, DM subdm, PetscCtx ctx) 81d71ae5a4SJacob Faibussowitsch { 8240be0ff1SMatthew G. Knepley PetscFunctionBegin; 833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8440be0ff1SMatthew G. Knepley } 8540be0ff1SMatthew G. Knepley 86*2a8381b2SBarry Smith static PetscErrorCode DMSubDomainRestrictHook_TSDiscGrad(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, PetscCtx ctx) 87d71ae5a4SJacob Faibussowitsch { 8840be0ff1SMatthew G. Knepley TS ts = (TS)ctx; 8940be0ff1SMatthew G. Knepley Vec X0, Xdot, X0_sub, Xdot_sub; 9040be0ff1SMatthew G. Knepley 9140be0ff1SMatthew G. Knepley PetscFunctionBegin; 929566063dSJacob Faibussowitsch PetscCall(TSDiscGradGetX0AndXdot(ts, dm, &X0, &Xdot)); 939566063dSJacob Faibussowitsch PetscCall(TSDiscGradGetX0AndXdot(ts, subdm, &X0_sub, &Xdot_sub)); 9440be0ff1SMatthew G. Knepley 959566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat, X0, X0_sub, INSERT_VALUES, SCATTER_FORWARD)); 969566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(gscat, X0, X0_sub, INSERT_VALUES, SCATTER_FORWARD)); 9740be0ff1SMatthew G. Knepley 989566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat, Xdot, Xdot_sub, INSERT_VALUES, SCATTER_FORWARD)); 999566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(gscat, Xdot, Xdot_sub, INSERT_VALUES, SCATTER_FORWARD)); 10040be0ff1SMatthew G. Knepley 1019566063dSJacob Faibussowitsch PetscCall(TSDiscGradRestoreX0AndXdot(ts, dm, &X0, &Xdot)); 1029566063dSJacob Faibussowitsch PetscCall(TSDiscGradRestoreX0AndXdot(ts, subdm, &X0_sub, &Xdot_sub)); 1033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10440be0ff1SMatthew G. Knepley } 10540be0ff1SMatthew G. Knepley 106d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_DiscGrad(TS ts) 107d71ae5a4SJacob Faibussowitsch { 10840be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *)ts->data; 10940be0ff1SMatthew G. Knepley DM dm; 11040be0ff1SMatthew G. Knepley 11140be0ff1SMatthew G. Knepley PetscFunctionBegin; 1129566063dSJacob Faibussowitsch if (!dg->X) PetscCall(VecDuplicate(ts->vec_sol, &dg->X)); 1139566063dSJacob Faibussowitsch if (!dg->X0) PetscCall(VecDuplicate(ts->vec_sol, &dg->X0)); 1149566063dSJacob Faibussowitsch if (!dg->Xdot) PetscCall(VecDuplicate(ts->vec_sol, &dg->Xdot)); 11540be0ff1SMatthew G. Knepley 1169566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 1179566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_TSDiscGrad, DMRestrictHook_TSDiscGrad, ts)); 1189566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_TSDiscGrad, DMSubDomainRestrictHook_TSDiscGrad, ts)); 1193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12040be0ff1SMatthew G. Knepley } 12140be0ff1SMatthew G. Knepley 122ce78bad3SBarry Smith static PetscErrorCode TSSetFromOptions_DiscGrad(TS ts, PetscOptionItems PetscOptionsObject) 123d71ae5a4SJacob Faibussowitsch { 124db4653c2SDaniel Finn TS_DiscGrad *dg = (TS_DiscGrad *)ts->data; 12540be0ff1SMatthew G. Knepley 12640be0ff1SMatthew G. Knepley PetscFunctionBegin; 127d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Discrete Gradients ODE solver options"); 128d71ae5a4SJacob Faibussowitsch { 129f940b0e3Sdanofinn PetscCall(PetscOptionsEnum("-ts_discgrad_type", "Type of discrete gradient solver", "TSDiscGradSetDGType", DGTypes, (PetscEnum)dg->discgrad, (PetscEnum *)&dg->discgrad, NULL)); 130d71ae5a4SJacob Faibussowitsch } 131d0609cedSBarry Smith PetscOptionsHeadEnd(); 1323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13340be0ff1SMatthew G. Knepley } 13440be0ff1SMatthew G. Knepley 135d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_DiscGrad(TS ts, PetscViewer viewer) 136d71ae5a4SJacob Faibussowitsch { 1379f196a02SMartin Diehl PetscBool isascii; 13840be0ff1SMatthew G. Knepley 13940be0ff1SMatthew G. Knepley PetscFunctionBegin; 1409f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 1419f196a02SMartin Diehl if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer, " Discrete Gradients\n")); 1423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14340be0ff1SMatthew G. Knepley } 14440be0ff1SMatthew G. Knepley 145f940b0e3Sdanofinn static PetscErrorCode TSDiscGradGetType_DiscGrad(TS ts, TSDGType *dgtype) 146d71ae5a4SJacob Faibussowitsch { 147db4653c2SDaniel Finn TS_DiscGrad *dg = (TS_DiscGrad *)ts->data; 148db4653c2SDaniel Finn 149db4653c2SDaniel Finn PetscFunctionBegin; 150f940b0e3Sdanofinn *dgtype = dg->discgrad; 1513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 152db4653c2SDaniel Finn } 153db4653c2SDaniel Finn 154f940b0e3Sdanofinn static PetscErrorCode TSDiscGradSetType_DiscGrad(TS ts, TSDGType dgtype) 155d71ae5a4SJacob Faibussowitsch { 156db4653c2SDaniel Finn TS_DiscGrad *dg = (TS_DiscGrad *)ts->data; 157db4653c2SDaniel Finn 158db4653c2SDaniel Finn PetscFunctionBegin; 159f940b0e3Sdanofinn dg->discgrad = dgtype; 1603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 161db4653c2SDaniel Finn } 162db4653c2SDaniel Finn 163d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_DiscGrad(TS ts) 164d71ae5a4SJacob Faibussowitsch { 16540be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *)ts->data; 16640be0ff1SMatthew G. Knepley 16740be0ff1SMatthew G. Knepley PetscFunctionBegin; 1689566063dSJacob Faibussowitsch PetscCall(VecDestroy(&dg->X)); 1699566063dSJacob Faibussowitsch PetscCall(VecDestroy(&dg->X0)); 1709566063dSJacob Faibussowitsch PetscCall(VecDestroy(&dg->Xdot)); 1713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 17240be0ff1SMatthew G. Knepley } 17340be0ff1SMatthew G. Knepley 174d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_DiscGrad(TS ts) 175d71ae5a4SJacob Faibussowitsch { 17640be0ff1SMatthew G. Knepley DM dm; 17740be0ff1SMatthew G. Knepley 17840be0ff1SMatthew G. Knepley PetscFunctionBegin; 1799566063dSJacob Faibussowitsch PetscCall(TSReset_DiscGrad(ts)); 1809566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 18140be0ff1SMatthew G. Knepley if (dm) { 1829566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookRemove(dm, DMCoarsenHook_TSDiscGrad, DMRestrictHook_TSDiscGrad, ts)); 1839566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookRemove(dm, DMSubDomainHook_TSDiscGrad, DMSubDomainRestrictHook_TSDiscGrad, ts)); 18440be0ff1SMatthew G. Knepley } 1859566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 1869566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSDiscGradGetFormulation_C", NULL)); 1879566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSDiscGradSetFormulation_C", NULL)); 188f940b0e3Sdanofinn PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSDiscGradGetType_C", NULL)); 189f940b0e3Sdanofinn PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSDiscGradSetType_C", NULL)); 1903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 19140be0ff1SMatthew G. Knepley } 19240be0ff1SMatthew G. Knepley 193d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_DiscGrad(TS ts, PetscReal t, Vec X) 194d71ae5a4SJacob Faibussowitsch { 19540be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *)ts->data; 19640be0ff1SMatthew G. Knepley PetscReal dt = t - ts->ptime; 19740be0ff1SMatthew G. Knepley 19840be0ff1SMatthew G. Knepley PetscFunctionBegin; 1999566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, dg->X)); 2009566063dSJacob Faibussowitsch PetscCall(VecWAXPY(X, dt, dg->Xdot, dg->X)); 2013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 20240be0ff1SMatthew G. Knepley } 20340be0ff1SMatthew G. Knepley 204d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDiscGrad_SNESSolve(TS ts, Vec b, Vec x) 205d71ae5a4SJacob Faibussowitsch { 20640be0ff1SMatthew G. Knepley SNES snes; 20740be0ff1SMatthew G. Knepley PetscInt nits, lits; 20840be0ff1SMatthew G. Knepley 20940be0ff1SMatthew G. Knepley PetscFunctionBegin; 2109566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 2119566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, b, x)); 2129566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &nits)); 2139566063dSJacob Faibussowitsch PetscCall(SNESGetLinearSolveIterations(snes, &lits)); 21440be0ff1SMatthew G. Knepley ts->snes_its += nits; 21540be0ff1SMatthew G. Knepley ts->ksp_its += lits; 2163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 21740be0ff1SMatthew G. Knepley } 21840be0ff1SMatthew G. Knepley 219d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_DiscGrad(TS ts) 220d71ae5a4SJacob Faibussowitsch { 22140be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *)ts->data; 22240be0ff1SMatthew G. Knepley TSAdapt adapt; 22340be0ff1SMatthew G. Knepley TSStepStatus status = TS_STEP_INCOMPLETE; 22440be0ff1SMatthew G. Knepley PetscInt rejections = 0; 22540be0ff1SMatthew G. Knepley PetscBool stageok, accept = PETSC_TRUE; 22640be0ff1SMatthew G. Knepley PetscReal next_time_step = ts->time_step; 22740be0ff1SMatthew G. Knepley 22840be0ff1SMatthew G. Knepley PetscFunctionBegin; 2299566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt)); 2309566063dSJacob Faibussowitsch if (!ts->steprollback) PetscCall(VecCopy(ts->vec_sol, dg->X0)); 23140be0ff1SMatthew G. Knepley 23240be0ff1SMatthew G. Knepley while (!ts->reason && status != TS_STEP_COMPLETE) { 23340be0ff1SMatthew G. Knepley PetscReal shift = 1 / (0.5 * ts->time_step); 23440be0ff1SMatthew G. Knepley 23540be0ff1SMatthew G. Knepley dg->stage_time = ts->ptime + 0.5 * ts->time_step; 23640be0ff1SMatthew G. Knepley 2379566063dSJacob Faibussowitsch PetscCall(VecCopy(dg->X0, dg->X)); 2389566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, dg->stage_time)); 2399566063dSJacob Faibussowitsch PetscCall(TSDiscGrad_SNESSolve(ts, NULL, dg->X)); 2409566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, dg->stage_time, 0, &dg->X)); 2419566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(adapt, ts, dg->stage_time, dg->X, &stageok)); 24240be0ff1SMatthew G. Knepley if (!stageok) goto reject_step; 24340be0ff1SMatthew G. Knepley 24440be0ff1SMatthew G. Knepley status = TS_STEP_PENDING; 2459566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(dg->Xdot, -shift, shift, 0, dg->X0, dg->X)); 2469566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vec_sol, ts->time_step, dg->Xdot)); 2479566063dSJacob Faibussowitsch PetscCall(TSAdaptChoose(adapt, ts, ts->time_step, NULL, &next_time_step, &accept)); 24840be0ff1SMatthew G. Knepley status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 24940be0ff1SMatthew G. Knepley if (!accept) { 2509566063dSJacob Faibussowitsch PetscCall(VecCopy(dg->X0, ts->vec_sol)); 25140be0ff1SMatthew G. Knepley ts->time_step = next_time_step; 25240be0ff1SMatthew G. Knepley goto reject_step; 25340be0ff1SMatthew G. Knepley } 25440be0ff1SMatthew G. Knepley ts->ptime += ts->time_step; 25540be0ff1SMatthew G. Knepley ts->time_step = next_time_step; 25640be0ff1SMatthew G. Knepley break; 25740be0ff1SMatthew G. Knepley 25840be0ff1SMatthew G. Knepley reject_step: 2599371c9d4SSatish Balay ts->reject++; 2609371c9d4SSatish Balay accept = PETSC_FALSE; 26140be0ff1SMatthew G. Knepley if (!ts->reason && ts->max_reject >= 0 && ++rejections > ts->max_reject) { 26240be0ff1SMatthew G. Knepley ts->reason = TS_DIVERGED_STEP_REJECTED; 26363a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections)); 26440be0ff1SMatthew G. Knepley } 26540be0ff1SMatthew G. Knepley } 2663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 26740be0ff1SMatthew G. Knepley } 26840be0ff1SMatthew G. Knepley 269d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGetStages_DiscGrad(TS ts, PetscInt *ns, Vec **Y) 270d71ae5a4SJacob Faibussowitsch { 27140be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *)ts->data; 27240be0ff1SMatthew G. Knepley 27340be0ff1SMatthew G. Knepley PetscFunctionBegin; 27440be0ff1SMatthew G. Knepley if (ns) *ns = 1; 275f4f49eeaSPierre Jolivet if (Y) *Y = &dg->X; 2763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 27740be0ff1SMatthew G. Knepley } 27840be0ff1SMatthew G. Knepley 27940be0ff1SMatthew G. Knepley /* 28040be0ff1SMatthew G. Knepley This defines the nonlinear equation that is to be solved with SNES 28140be0ff1SMatthew G. Knepley G(U) = F[t0 + 0.5*dt, U, (U-U0)/dt] = 0 28240be0ff1SMatthew G. Knepley */ 283db4653c2SDaniel Finn 284db4653c2SDaniel Finn /* x = (x+x')/2 */ 285db4653c2SDaniel Finn /* NEED TO CALCULATE x_{n+1} from x and x_{n}*/ 286d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_DiscGrad(SNES snes, Vec x, Vec y, TS ts) 287d71ae5a4SJacob Faibussowitsch { 28840be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *)ts->data; 289db4653c2SDaniel Finn PetscReal norm, shift = 1 / (0.5 * ts->time_step); 290f940b0e3Sdanofinn PetscInt n, dim; 291db4653c2SDaniel Finn Vec X0, Xdot, Xp, Xdiff; 292db4653c2SDaniel Finn Mat S; 293db4653c2SDaniel Finn PetscScalar F = 0, F0 = 0, Gp; 294db4653c2SDaniel Finn Vec G, SgF; 29540be0ff1SMatthew G. Knepley DM dm, dmsave; 29640be0ff1SMatthew G. Knepley 29740be0ff1SMatthew G. Knepley PetscFunctionBegin; 2989566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 299f940b0e3Sdanofinn PetscCall(DMGetDimension(dm, &dim)); 30040be0ff1SMatthew G. Knepley 3019566063dSJacob Faibussowitsch PetscCall(VecDuplicate(y, &Xp)); 3029566063dSJacob Faibussowitsch PetscCall(VecDuplicate(y, &Xdiff)); 3039566063dSJacob Faibussowitsch PetscCall(VecDuplicate(y, &SgF)); 3049566063dSJacob Faibussowitsch PetscCall(VecDuplicate(y, &G)); 305db4653c2SDaniel Finn 306f940b0e3Sdanofinn PetscCall(PetscObjectSetName((PetscObject)x, "x")); 307f940b0e3Sdanofinn PetscCall(VecViewFromOptions(x, NULL, "-x_view")); 308f940b0e3Sdanofinn 3099566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(y, &n)); 3109566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &S)); 311f940b0e3Sdanofinn PetscCall(MatSetSizes(S, n, n, PETSC_DECIDE, PETSC_DECIDE)); 3129566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(S)); 313f940b0e3Sdanofinn PetscInt *S_prealloc_arr; 314f940b0e3Sdanofinn PetscCall(PetscMalloc1(n, &S_prealloc_arr)); 315ac530a7eSPierre Jolivet for (PetscInt i = 0; i < n; ++i) S_prealloc_arr[i] = 2; 316f940b0e3Sdanofinn PetscCall(MatXAIJSetPreallocation(S, 1, S_prealloc_arr, NULL, NULL, NULL)); 3179566063dSJacob Faibussowitsch PetscCall(MatSetUp(S)); 318f940b0e3Sdanofinn PetscCall((*dg->Sfunc)(ts, dg->stage_time, x, S, dg->funcCtx)); 319f940b0e3Sdanofinn PetscCall(PetscFree(S_prealloc_arr)); 320f940b0e3Sdanofinn PetscCall(PetscObjectSetName((PetscObject)S, "S")); 321f940b0e3Sdanofinn PetscCall(MatViewFromOptions(S, NULL, "-S_view")); 3229566063dSJacob Faibussowitsch PetscCall(TSDiscGradGetX0AndXdot(ts, dm, &X0, &Xdot)); 3239566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(Xdot, -shift, shift, 0, X0, x)); /* Xdot = shift (x - X0) */ 324db4653c2SDaniel Finn 325f940b0e3Sdanofinn PetscCall(VecAXPBYPCZ(Xp, -1, 2, 0, X0, x)); /* Xp = 2*x - X0 + (0)*Xp */ 3269566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(Xdiff, -1, 1, 0, X0, Xp)); /* Xdiff = xp - X0 + (0)*Xdiff */ 327db4653c2SDaniel Finn 328f940b0e3Sdanofinn PetscCall(PetscObjectSetName((PetscObject)X0, "X0")); 329f940b0e3Sdanofinn PetscCall(PetscObjectSetName((PetscObject)Xp, "Xp")); 330f940b0e3Sdanofinn PetscCall(VecViewFromOptions(X0, NULL, "-X0_view")); 331f940b0e3Sdanofinn PetscCall(VecViewFromOptions(Xp, NULL, "-Xp_view")); 332f940b0e3Sdanofinn 333f940b0e3Sdanofinn if (dg->discgrad == TS_DG_AVERAGE) { 334f940b0e3Sdanofinn /* Average Value DG: 335f940b0e3Sdanofinn \overline{\nabla} F (x_{n+1},x_{n}) = \int_0^1 \nabla F ((1-\xi)*x_{n+1} + \xi*x_{n}) d \xi */ 336f940b0e3Sdanofinn PetscQuadrature quad; 337f940b0e3Sdanofinn PetscInt Nq; 338f940b0e3Sdanofinn const PetscReal *wq, *xq; 339f940b0e3Sdanofinn Vec Xquad, den; 340f940b0e3Sdanofinn 341f940b0e3Sdanofinn PetscCall(PetscObjectSetName((PetscObject)G, "G")); 342f940b0e3Sdanofinn PetscCall(VecDuplicate(G, &Xquad)); 343f940b0e3Sdanofinn PetscCall(VecDuplicate(G, &den)); 344f940b0e3Sdanofinn PetscCall(VecZeroEntries(G)); 345f940b0e3Sdanofinn 346f940b0e3Sdanofinn /* \overline{\nabla} F = \nabla F ((1-\xi) x_{n} + \xi x_{n+1})*/ 347f940b0e3Sdanofinn PetscCall(PetscDTGaussTensorQuadrature(dim, 1, 2, 0.0, 1.0, &quad)); 348f940b0e3Sdanofinn PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, &xq, &wq)); 349f940b0e3Sdanofinn for (PetscInt q = 0; q < Nq; ++q) { 350f940b0e3Sdanofinn PetscReal xi = xq[q], xim1 = 1 - xq[q]; 351f940b0e3Sdanofinn PetscCall(VecZeroEntries(Xquad)); 352f940b0e3Sdanofinn PetscCall(VecAXPBYPCZ(Xquad, xi, xim1, 1.0, X0, Xp)); 353f940b0e3Sdanofinn PetscCall((*dg->Gfunc)(ts, dg->stage_time, Xquad, den, dg->funcCtx)); 354f940b0e3Sdanofinn PetscCall(VecAXPY(G, wq[q], den)); 355f940b0e3Sdanofinn PetscCall(PetscObjectSetName((PetscObject)den, "den")); 356f940b0e3Sdanofinn PetscCall(VecViewFromOptions(den, NULL, "-den_view")); 357f940b0e3Sdanofinn } 358f940b0e3Sdanofinn PetscCall(VecDestroy(&Xquad)); 359f940b0e3Sdanofinn PetscCall(VecDestroy(&den)); 360f940b0e3Sdanofinn PetscCall(PetscQuadratureDestroy(&quad)); 361f940b0e3Sdanofinn } else if (dg->discgrad == TS_DG_GONZALEZ) { 3629566063dSJacob Faibussowitsch PetscCall((*dg->Ffunc)(ts, dg->stage_time, Xp, &F, dg->funcCtx)); 3639566063dSJacob Faibussowitsch PetscCall((*dg->Ffunc)(ts, dg->stage_time, X0, &F0, dg->funcCtx)); 3649566063dSJacob Faibussowitsch PetscCall((*dg->Gfunc)(ts, dg->stage_time, x, G, dg->funcCtx)); 365db4653c2SDaniel Finn 366db4653c2SDaniel Finn /* Adding Extra Gonzalez Term */ 3679566063dSJacob Faibussowitsch PetscCall(VecDot(Xdiff, G, &Gp)); 3689566063dSJacob Faibussowitsch PetscCall(VecNorm(Xdiff, NORM_2, &norm)); 369db4653c2SDaniel Finn if (norm < PETSC_SQRT_MACHINE_EPSILON) { 370db4653c2SDaniel Finn Gp = 0; 371db4653c2SDaniel Finn } else { 372db4653c2SDaniel Finn /* Gp = (1/|xn+1 - xn|^2) * (F(xn+1) - F(xn) - Gp) */ 373db4653c2SDaniel Finn Gp = (F - F0 - Gp) / PetscSqr(norm); 374db4653c2SDaniel Finn } 3759566063dSJacob Faibussowitsch PetscCall(VecAXPY(G, Gp, Xdiff)); 376f940b0e3Sdanofinn } else if (dg->discgrad == TS_DG_NONE) { 3779566063dSJacob Faibussowitsch PetscCall((*dg->Gfunc)(ts, dg->stage_time, x, G, dg->funcCtx)); 378f940b0e3Sdanofinn } else { 379f940b0e3Sdanofinn SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "DG type not supported."); 380db4653c2SDaniel Finn } 381f940b0e3Sdanofinn PetscCall(MatMult(S, G, SgF)); /* Xdot = S*gradF */ 382f940b0e3Sdanofinn 383f940b0e3Sdanofinn PetscCall(PetscObjectSetName((PetscObject)G, "G")); 384f940b0e3Sdanofinn PetscCall(VecViewFromOptions(G, NULL, "-G_view")); 385f940b0e3Sdanofinn PetscCall(PetscObjectSetName((PetscObject)SgF, "SgF")); 386f940b0e3Sdanofinn PetscCall(VecViewFromOptions(SgF, NULL, "-SgF_view")); 38740be0ff1SMatthew G. Knepley /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 38840be0ff1SMatthew G. Knepley dmsave = ts->dm; 38940be0ff1SMatthew G. Knepley ts->dm = dm; 3909566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(y, 1, -1, 0, Xdot, SgF)); 391f940b0e3Sdanofinn 39240be0ff1SMatthew G. Knepley ts->dm = dmsave; 3939566063dSJacob Faibussowitsch PetscCall(TSDiscGradRestoreX0AndXdot(ts, dm, &X0, &Xdot)); 394db4653c2SDaniel Finn 3959566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Xp)); 3969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Xdiff)); 3979566063dSJacob Faibussowitsch PetscCall(VecDestroy(&SgF)); 3989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&G)); 3999566063dSJacob Faibussowitsch PetscCall(MatDestroy(&S)); 4003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 40140be0ff1SMatthew G. Knepley } 40240be0ff1SMatthew G. Knepley 403d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_DiscGrad(SNES snes, Vec x, Mat A, Mat B, TS ts) 404d71ae5a4SJacob Faibussowitsch { 40540be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *)ts->data; 40640be0ff1SMatthew G. Knepley PetscReal shift = 1 / (0.5 * ts->time_step); 40740be0ff1SMatthew G. Knepley Vec Xdot; 40840be0ff1SMatthew G. Knepley DM dm, dmsave; 40940be0ff1SMatthew G. Knepley 41040be0ff1SMatthew G. Knepley PetscFunctionBegin; 4119566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 41240be0ff1SMatthew G. Knepley /* Xdot has already been computed in SNESTSFormFunction_DiscGrad (SNES guarantees this) */ 4139566063dSJacob Faibussowitsch PetscCall(TSDiscGradGetX0AndXdot(ts, dm, NULL, &Xdot)); 41440be0ff1SMatthew G. Knepley 41540be0ff1SMatthew G. Knepley dmsave = ts->dm; 41640be0ff1SMatthew G. Knepley ts->dm = dm; 4179566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts, dg->stage_time, x, Xdot, shift, A, B, PETSC_FALSE)); 41840be0ff1SMatthew G. Knepley ts->dm = dmsave; 4199566063dSJacob Faibussowitsch PetscCall(TSDiscGradRestoreX0AndXdot(ts, dm, NULL, &Xdot)); 4203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 42140be0ff1SMatthew G. Knepley } 42240be0ff1SMatthew G. Knepley 423*2a8381b2SBarry Smith static PetscErrorCode TSDiscGradGetFormulation_DiscGrad(TS ts, PetscErrorCode (**Sfunc)(TS, PetscReal, Vec, Mat, void *), PetscErrorCode (**Ffunc)(TS, PetscReal, Vec, PetscScalar *, void *), PetscErrorCode (**Gfunc)(TS, PetscReal, Vec, Vec, void *), PetscCtx ctx) 424d71ae5a4SJacob Faibussowitsch { 42540be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *)ts->data; 42640be0ff1SMatthew G. Knepley 42740be0ff1SMatthew G. Knepley PetscFunctionBegin; 42840be0ff1SMatthew G. Knepley *Sfunc = dg->Sfunc; 42940be0ff1SMatthew G. Knepley *Ffunc = dg->Ffunc; 43040be0ff1SMatthew G. Knepley *Gfunc = dg->Gfunc; 4313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 43240be0ff1SMatthew G. Knepley } 43340be0ff1SMatthew G. Knepley 434*2a8381b2SBarry Smith static PetscErrorCode TSDiscGradSetFormulation_DiscGrad(TS ts, PetscErrorCode (*Sfunc)(TS, PetscReal, Vec, Mat, void *), PetscErrorCode (*Ffunc)(TS, PetscReal, Vec, PetscScalar *, void *), PetscErrorCode (*Gfunc)(TS, PetscReal, Vec, Vec, void *), PetscCtx ctx) 435d71ae5a4SJacob Faibussowitsch { 43640be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *)ts->data; 43740be0ff1SMatthew G. Knepley 43840be0ff1SMatthew G. Knepley PetscFunctionBegin; 43940be0ff1SMatthew G. Knepley dg->Sfunc = Sfunc; 44040be0ff1SMatthew G. Knepley dg->Ffunc = Ffunc; 44140be0ff1SMatthew G. Knepley dg->Gfunc = Gfunc; 442db4653c2SDaniel Finn dg->funcCtx = ctx; 4433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 44440be0ff1SMatthew G. Knepley } 44540be0ff1SMatthew G. Knepley 44640be0ff1SMatthew G. Knepley /*MC 44740be0ff1SMatthew G. Knepley TSDISCGRAD - ODE solver using the discrete gradients version of the implicit midpoint method 44840be0ff1SMatthew G. Knepley 449bcf0153eSBarry Smith Level: intermediate 450bcf0153eSBarry Smith 4512ceb51e8SPatrick Sanan Notes: 45214d0ab18SJacob Faibussowitsch This is the implicit midpoint rule, with an optional term that guarantees the discrete 45314d0ab18SJacob Faibussowitsch gradient property. This timestepper applies to systems of the form $u_t = S(u) \nabla F(u)$ 45414d0ab18SJacob Faibussowitsch where $S(u)$ is a linear operator, and $F$ is a functional of $u$. 45540be0ff1SMatthew G. Knepley 4561cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TSSetType()`, `TS`, `TSDISCGRAD`, `TSDiscGradSetFormulation()` 45740be0ff1SMatthew G. Knepley M*/ 458d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_DiscGrad(TS ts) 459d71ae5a4SJacob Faibussowitsch { 46040be0ff1SMatthew G. Knepley TS_DiscGrad *th; 46140be0ff1SMatthew G. Knepley 46240be0ff1SMatthew G. Knepley PetscFunctionBegin; 4639566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(DGCitation, &DGCite)); 46440be0ff1SMatthew G. Knepley ts->ops->reset = TSReset_DiscGrad; 46540be0ff1SMatthew G. Knepley ts->ops->destroy = TSDestroy_DiscGrad; 46640be0ff1SMatthew G. Knepley ts->ops->view = TSView_DiscGrad; 46740be0ff1SMatthew G. Knepley ts->ops->setfromoptions = TSSetFromOptions_DiscGrad; 46840be0ff1SMatthew G. Knepley ts->ops->setup = TSSetUp_DiscGrad; 46940be0ff1SMatthew G. Knepley ts->ops->step = TSStep_DiscGrad; 47040be0ff1SMatthew G. Knepley ts->ops->interpolate = TSInterpolate_DiscGrad; 47140be0ff1SMatthew G. Knepley ts->ops->getstages = TSGetStages_DiscGrad; 47240be0ff1SMatthew G. Knepley ts->ops->snesfunction = SNESTSFormFunction_DiscGrad; 47340be0ff1SMatthew G. Knepley ts->ops->snesjacobian = SNESTSFormJacobian_DiscGrad; 47440be0ff1SMatthew G. Knepley ts->default_adapt_type = TSADAPTNONE; 47540be0ff1SMatthew G. Knepley 47640be0ff1SMatthew G. Knepley ts->usessnes = PETSC_TRUE; 47740be0ff1SMatthew G. Knepley 4784dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&th)); 47940be0ff1SMatthew G. Knepley ts->data = (void *)th; 480db4653c2SDaniel Finn 481f940b0e3Sdanofinn th->discgrad = TS_DG_NONE; 482db4653c2SDaniel Finn 4839566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSDiscGradGetFormulation_C", TSDiscGradGetFormulation_DiscGrad)); 4849566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSDiscGradSetFormulation_C", TSDiscGradSetFormulation_DiscGrad)); 485f940b0e3Sdanofinn PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSDiscGradGetType_C", TSDiscGradGetType_DiscGrad)); 486f940b0e3Sdanofinn PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSDiscGradSetType_C", TSDiscGradSetType_DiscGrad)); 4873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 48840be0ff1SMatthew G. Knepley } 48940be0ff1SMatthew G. Knepley 49040be0ff1SMatthew G. Knepley /*@C 49114d0ab18SJacob Faibussowitsch TSDiscGradGetFormulation - Get the construction method for S, F, and grad F from the 49214d0ab18SJacob Faibussowitsch formulation $u_t = S \nabla F$ for `TSDISCGRAD` 49340be0ff1SMatthew G. Knepley 49440be0ff1SMatthew G. Knepley Not Collective 49540be0ff1SMatthew G. Knepley 49640be0ff1SMatthew G. Knepley Input Parameter: 49740be0ff1SMatthew G. Knepley . ts - timestepping context 49840be0ff1SMatthew G. Knepley 49940be0ff1SMatthew G. Knepley Output Parameters: 50040be0ff1SMatthew G. Knepley + Sfunc - constructor for the S matrix from the formulation 50140be0ff1SMatthew G. Knepley . Ffunc - functional F from the formulation 502f1a722f8SMatthew G. Knepley . Gfunc - constructor for the gradient of F from the formulation 503f1a722f8SMatthew G. Knepley - ctx - the user context 50440be0ff1SMatthew G. Knepley 50520f4b53cSBarry Smith Calling sequence of `Sfunc`: 5068847d985SBarry Smith + ts - the integrator 5078847d985SBarry Smith . time - the current time 5088847d985SBarry Smith . u - the solution 5098847d985SBarry Smith . S - the S-matrix from the formulation 5108847d985SBarry Smith - ctx - the user context 51140be0ff1SMatthew G. Knepley 51220f4b53cSBarry Smith Calling sequence of `Ffunc`: 5138847d985SBarry Smith + ts - the integrator 5148847d985SBarry Smith . time - the current time 5158847d985SBarry Smith . u - the solution 5168847d985SBarry Smith . F - the computed function from the formulation 5178847d985SBarry Smith - ctx - the user context 51840be0ff1SMatthew G. Knepley 51920f4b53cSBarry Smith Calling sequence of `Gfunc`: 5208847d985SBarry Smith + ts - the integrator 5218847d985SBarry Smith . time - the current time 5228847d985SBarry Smith . u - the solution 5238847d985SBarry Smith . G - the gradient of the computed function from the formulation 5248847d985SBarry Smith - ctx - the user context 52540be0ff1SMatthew G. Knepley 52640be0ff1SMatthew G. Knepley Level: intermediate 52740be0ff1SMatthew G. Knepley 5281cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSDISCGRAD`, `TSDiscGradSetFormulation()` 52940be0ff1SMatthew G. Knepley @*/ 530*2a8381b2SBarry Smith PetscErrorCode TSDiscGradGetFormulation(TS ts, PetscErrorCode (**Sfunc)(TS ts, PetscReal time, Vec u, Mat S, PetscCtx ctx), PetscErrorCode (**Ffunc)(TS ts, PetscReal time, Vec u, PetscScalar *F, PetscCtx ctx), PetscErrorCode (**Gfunc)(TS ts, PetscReal time, Vec u, Vec G, PetscCtx ctx), PetscCtx ctx) 531d71ae5a4SJacob Faibussowitsch { 53240be0ff1SMatthew G. Knepley PetscFunctionBegin; 53340be0ff1SMatthew G. Knepley PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 5344f572ea9SToby Isaac PetscAssertPointer(Sfunc, 2); 5354f572ea9SToby Isaac PetscAssertPointer(Ffunc, 3); 5364f572ea9SToby Isaac PetscAssertPointer(Gfunc, 4); 537cac4c232SBarry Smith PetscUseMethod(ts, "TSDiscGradGetFormulation_C", (TS, PetscErrorCode (**Sfunc)(TS, PetscReal, Vec, Mat, void *), PetscErrorCode (**Ffunc)(TS, PetscReal, Vec, PetscScalar *, void *), PetscErrorCode (**Gfunc)(TS, PetscReal, Vec, Vec, void *), void *), (ts, Sfunc, Ffunc, Gfunc, ctx)); 5383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 53940be0ff1SMatthew G. Knepley } 54040be0ff1SMatthew G. Knepley 54140be0ff1SMatthew G. Knepley /*@C 54214d0ab18SJacob Faibussowitsch TSDiscGradSetFormulation - Set the construction method for S, F, and grad F from the 54314d0ab18SJacob Faibussowitsch formulation $u_t = S(u) \nabla F(u)$ for `TSDISCGRAD` 54440be0ff1SMatthew G. Knepley 54540be0ff1SMatthew G. Knepley Not Collective 54640be0ff1SMatthew G. Knepley 54740be0ff1SMatthew G. Knepley Input Parameters: 54840be0ff1SMatthew G. Knepley + ts - timestepping context 54940be0ff1SMatthew G. Knepley . Sfunc - constructor for the S matrix from the formulation 55040be0ff1SMatthew G. Knepley . Ffunc - functional F from the formulation 5512fe279fdSBarry Smith . Gfunc - constructor for the gradient of F from the formulation 5522fe279fdSBarry Smith - ctx - optional context for the functions 55340be0ff1SMatthew G. Knepley 55420f4b53cSBarry Smith Calling sequence of `Sfunc`: 5558847d985SBarry Smith + ts - the integrator 5568847d985SBarry Smith . time - the current time 5578847d985SBarry Smith . u - the solution 5588847d985SBarry Smith . S - the S-matrix from the formulation 5598847d985SBarry Smith - ctx - the user context 56040be0ff1SMatthew G. Knepley 56120f4b53cSBarry Smith Calling sequence of `Ffunc`: 5628847d985SBarry Smith + ts - the integrator 5638847d985SBarry Smith . time - the current time 5648847d985SBarry Smith . u - the solution 5658847d985SBarry Smith . F - the computed function from the formulation 5668847d985SBarry Smith - ctx - the user context 56720f4b53cSBarry Smith 56820f4b53cSBarry Smith Calling sequence of `Gfunc`: 5698847d985SBarry Smith + ts - the integrator 5708847d985SBarry Smith . time - the current time 5718847d985SBarry Smith . u - the solution 5728847d985SBarry Smith . G - the gradient of the computed function from the formulation 5738847d985SBarry Smith - ctx - the user context 57440be0ff1SMatthew G. Knepley 575b43aa488SJacob Faibussowitsch Level: intermediate 57640be0ff1SMatthew G. Knepley 5771cc06b55SBarry Smith .seealso: [](ch_ts), `TSDISCGRAD`, `TSDiscGradGetFormulation()` 57840be0ff1SMatthew G. Knepley @*/ 579*2a8381b2SBarry Smith PetscErrorCode TSDiscGradSetFormulation(TS ts, PetscErrorCode (*Sfunc)(TS ts, PetscReal time, Vec u, Mat S, PetscCtx ctx), PetscErrorCode (*Ffunc)(TS ts, PetscReal time, Vec u, PetscScalar *F, PetscCtx ctx), PetscErrorCode (*Gfunc)(TS ts, PetscReal time, Vec u, Vec G, PetscCtx ctx), PetscCtx ctx) 580d71ae5a4SJacob Faibussowitsch { 58140be0ff1SMatthew G. Knepley PetscFunctionBegin; 58240be0ff1SMatthew G. Knepley PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 58340be0ff1SMatthew G. Knepley PetscValidFunction(Sfunc, 2); 58440be0ff1SMatthew G. Knepley PetscValidFunction(Ffunc, 3); 58540be0ff1SMatthew G. Knepley PetscValidFunction(Gfunc, 4); 586cac4c232SBarry Smith PetscTryMethod(ts, "TSDiscGradSetFormulation_C", (TS, PetscErrorCode (*Sfunc)(TS, PetscReal, Vec, Mat, void *), PetscErrorCode (*Ffunc)(TS, PetscReal, Vec, PetscScalar *, void *), PetscErrorCode (*Gfunc)(TS, PetscReal, Vec, Vec, void *), void *), (ts, Sfunc, Ffunc, Gfunc, ctx)); 5873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 588db4653c2SDaniel Finn } 589db4653c2SDaniel Finn 590db4653c2SDaniel Finn /*@ 591f940b0e3Sdanofinn TSDiscGradGetType - Checks for which discrete gradient to use in formulation for `TSDISCGRAD` 592db4653c2SDaniel Finn 593db4653c2SDaniel Finn Not Collective 594db4653c2SDaniel Finn 595db4653c2SDaniel Finn Input Parameter: 596db4653c2SDaniel Finn . ts - timestepping context 597db4653c2SDaniel Finn 598db4653c2SDaniel Finn Output Parameter: 599f940b0e3Sdanofinn . dgtype - Discrete gradient type <none, gonzalez, average> 600db4653c2SDaniel Finn 601b43aa488SJacob Faibussowitsch Level: advanced 602db4653c2SDaniel Finn 603f940b0e3Sdanofinn .seealso: [](ch_ts), `TSDISCGRAD`, `TSDiscGradSetType()` 604db4653c2SDaniel Finn @*/ 605f940b0e3Sdanofinn PetscErrorCode TSDiscGradGetType(TS ts, TSDGType *dgtype) 606d71ae5a4SJacob Faibussowitsch { 607db4653c2SDaniel Finn PetscFunctionBegin; 608db4653c2SDaniel Finn PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 609f940b0e3Sdanofinn PetscAssertPointer(dgtype, 2); 610f940b0e3Sdanofinn PetscUseMethod(ts, "TSDiscGradGetType_C", (TS, TSDGType *), (ts, dgtype)); 6113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 612db4653c2SDaniel Finn } 613db4653c2SDaniel Finn 614db4653c2SDaniel Finn /*@ 615f940b0e3Sdanofinn TSDiscGradSetType - Sets discrete gradient formulation. 616db4653c2SDaniel Finn 617db4653c2SDaniel Finn Not Collective 618db4653c2SDaniel Finn 619f1a722f8SMatthew G. Knepley Input Parameters: 620db4653c2SDaniel Finn + ts - timestepping context 621f940b0e3Sdanofinn - dgtype - Discrete gradient type <none, gonzalez, average> 622db4653c2SDaniel Finn 623bcf0153eSBarry Smith Options Database Key: 624f940b0e3Sdanofinn . -ts_discgrad_type <type> - flag to choose discrete gradient type 625db4653c2SDaniel Finn 626b43aa488SJacob Faibussowitsch Level: intermediate 627db4653c2SDaniel Finn 62814d0ab18SJacob Faibussowitsch Notes: 629f940b0e3Sdanofinn Without `dgtype` or with type `none`, the discrete gradients timestepper is just implicit midpoint. 63014d0ab18SJacob Faibussowitsch 6311cc06b55SBarry Smith .seealso: [](ch_ts), `TSDISCGRAD` 632db4653c2SDaniel Finn @*/ 633f940b0e3Sdanofinn PetscErrorCode TSDiscGradSetType(TS ts, TSDGType dgtype) 634d71ae5a4SJacob Faibussowitsch { 635db4653c2SDaniel Finn PetscFunctionBegin; 636db4653c2SDaniel Finn PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 637f940b0e3Sdanofinn PetscTryMethod(ts, "TSDiscGradSetType_C", (TS, TSDGType), (ts, dgtype)); 6383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 63940be0ff1SMatthew G. Knepley } 640