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 1740be0ff1SMatthew G. Knepley typedef struct { 1840be0ff1SMatthew G. Knepley PetscReal stage_time; 1940be0ff1SMatthew G. Knepley Vec X0, X, Xdot; 20db4653c2SDaniel Finn void *funcCtx; 21db4653c2SDaniel Finn PetscBool gonzalez; 2240be0ff1SMatthew G. Knepley PetscErrorCode (*Sfunc)(TS, PetscReal, Vec, Mat, void *); 2340be0ff1SMatthew G. Knepley PetscErrorCode (*Ffunc)(TS, PetscReal, Vec, PetscScalar *, void *); 2440be0ff1SMatthew G. Knepley PetscErrorCode (*Gfunc)(TS, PetscReal, Vec, Vec, void *); 2540be0ff1SMatthew G. Knepley } TS_DiscGrad; 2640be0ff1SMatthew G. Knepley 279371c9d4SSatish Balay static PetscErrorCode TSDiscGradGetX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot) { 2840be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *)ts->data; 2940be0ff1SMatthew G. Knepley 3040be0ff1SMatthew G. Knepley PetscFunctionBegin; 3140be0ff1SMatthew G. Knepley if (X0) { 329566063dSJacob Faibussowitsch if (dm && dm != ts->dm) PetscCall(DMGetNamedGlobalVector(dm, "TSDiscGrad_X0", X0)); 3340be0ff1SMatthew G. Knepley else { *X0 = ts->vec_sol; } 3440be0ff1SMatthew G. Knepley } 3540be0ff1SMatthew G. Knepley if (Xdot) { 369566063dSJacob Faibussowitsch if (dm && dm != ts->dm) PetscCall(DMGetNamedGlobalVector(dm, "TSDiscGrad_Xdot", Xdot)); 3740be0ff1SMatthew G. Knepley else { *Xdot = dg->Xdot; } 3840be0ff1SMatthew G. Knepley } 3940be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 4040be0ff1SMatthew G. Knepley } 4140be0ff1SMatthew G. Knepley 429371c9d4SSatish Balay static PetscErrorCode TSDiscGradRestoreX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot) { 4340be0ff1SMatthew G. Knepley PetscFunctionBegin; 4440be0ff1SMatthew G. Knepley if (X0) { 459566063dSJacob Faibussowitsch if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSDiscGrad_X0", X0)); 4640be0ff1SMatthew G. Knepley } 4740be0ff1SMatthew G. Knepley if (Xdot) { 489566063dSJacob Faibussowitsch if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSDiscGrad_Xdot", Xdot)); 4940be0ff1SMatthew G. Knepley } 5040be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 5140be0ff1SMatthew G. Knepley } 5240be0ff1SMatthew G. Knepley 539371c9d4SSatish Balay static PetscErrorCode DMCoarsenHook_TSDiscGrad(DM fine, DM coarse, void *ctx) { 5440be0ff1SMatthew G. Knepley PetscFunctionBegin; 5540be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 5640be0ff1SMatthew G. Knepley } 5740be0ff1SMatthew G. Knepley 589371c9d4SSatish Balay static PetscErrorCode DMRestrictHook_TSDiscGrad(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx) { 5940be0ff1SMatthew G. Knepley TS ts = (TS)ctx; 6040be0ff1SMatthew G. Knepley Vec X0, Xdot, X0_c, Xdot_c; 6140be0ff1SMatthew G. Knepley 6240be0ff1SMatthew G. Knepley PetscFunctionBegin; 639566063dSJacob Faibussowitsch PetscCall(TSDiscGradGetX0AndXdot(ts, fine, &X0, &Xdot)); 649566063dSJacob Faibussowitsch PetscCall(TSDiscGradGetX0AndXdot(ts, coarse, &X0_c, &Xdot_c)); 659566063dSJacob Faibussowitsch PetscCall(MatRestrict(restrct, X0, X0_c)); 669566063dSJacob Faibussowitsch PetscCall(MatRestrict(restrct, Xdot, Xdot_c)); 679566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(X0_c, rscale, X0_c)); 689566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Xdot_c, rscale, Xdot_c)); 699566063dSJacob Faibussowitsch PetscCall(TSDiscGradRestoreX0AndXdot(ts, fine, &X0, &Xdot)); 709566063dSJacob Faibussowitsch PetscCall(TSDiscGradRestoreX0AndXdot(ts, coarse, &X0_c, &Xdot_c)); 7140be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 7240be0ff1SMatthew G. Knepley } 7340be0ff1SMatthew G. Knepley 749371c9d4SSatish Balay static PetscErrorCode DMSubDomainHook_TSDiscGrad(DM dm, DM subdm, void *ctx) { 7540be0ff1SMatthew G. Knepley PetscFunctionBegin; 7640be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 7740be0ff1SMatthew G. Knepley } 7840be0ff1SMatthew G. Knepley 799371c9d4SSatish Balay static PetscErrorCode DMSubDomainRestrictHook_TSDiscGrad(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx) { 8040be0ff1SMatthew G. Knepley TS ts = (TS)ctx; 8140be0ff1SMatthew G. Knepley Vec X0, Xdot, X0_sub, Xdot_sub; 8240be0ff1SMatthew G. Knepley 8340be0ff1SMatthew G. Knepley PetscFunctionBegin; 849566063dSJacob Faibussowitsch PetscCall(TSDiscGradGetX0AndXdot(ts, dm, &X0, &Xdot)); 859566063dSJacob Faibussowitsch PetscCall(TSDiscGradGetX0AndXdot(ts, subdm, &X0_sub, &Xdot_sub)); 8640be0ff1SMatthew G. Knepley 879566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat, X0, X0_sub, INSERT_VALUES, SCATTER_FORWARD)); 889566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(gscat, X0, X0_sub, INSERT_VALUES, SCATTER_FORWARD)); 8940be0ff1SMatthew G. Knepley 909566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat, Xdot, Xdot_sub, INSERT_VALUES, SCATTER_FORWARD)); 919566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(gscat, Xdot, Xdot_sub, INSERT_VALUES, SCATTER_FORWARD)); 9240be0ff1SMatthew G. Knepley 939566063dSJacob Faibussowitsch PetscCall(TSDiscGradRestoreX0AndXdot(ts, dm, &X0, &Xdot)); 949566063dSJacob Faibussowitsch PetscCall(TSDiscGradRestoreX0AndXdot(ts, subdm, &X0_sub, &Xdot_sub)); 9540be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 9640be0ff1SMatthew G. Knepley } 9740be0ff1SMatthew G. Knepley 989371c9d4SSatish Balay static PetscErrorCode TSSetUp_DiscGrad(TS ts) { 9940be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *)ts->data; 10040be0ff1SMatthew G. Knepley DM dm; 10140be0ff1SMatthew G. Knepley 10240be0ff1SMatthew G. Knepley PetscFunctionBegin; 1039566063dSJacob Faibussowitsch if (!dg->X) PetscCall(VecDuplicate(ts->vec_sol, &dg->X)); 1049566063dSJacob Faibussowitsch if (!dg->X0) PetscCall(VecDuplicate(ts->vec_sol, &dg->X0)); 1059566063dSJacob Faibussowitsch if (!dg->Xdot) PetscCall(VecDuplicate(ts->vec_sol, &dg->Xdot)); 10640be0ff1SMatthew G. Knepley 1079566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 1089566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_TSDiscGrad, DMRestrictHook_TSDiscGrad, ts)); 1099566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_TSDiscGrad, DMSubDomainRestrictHook_TSDiscGrad, ts)); 11040be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 11140be0ff1SMatthew G. Knepley } 11240be0ff1SMatthew G. Knepley 1139371c9d4SSatish Balay static PetscErrorCode TSSetFromOptions_DiscGrad(TS ts, PetscOptionItems *PetscOptionsObject) { 114db4653c2SDaniel Finn TS_DiscGrad *dg = (TS_DiscGrad *)ts->data; 11540be0ff1SMatthew G. Knepley 11640be0ff1SMatthew G. Knepley PetscFunctionBegin; 117d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Discrete Gradients ODE solver options"); 1189371c9d4SSatish Balay { PetscCall(PetscOptionsBool("-ts_discgrad_gonzalez", "Use Gonzalez term in discrete gradients formulation", "TSDiscGradUseGonzalez", dg->gonzalez, &dg->gonzalez, NULL)); } 119d0609cedSBarry Smith PetscOptionsHeadEnd(); 12040be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 12140be0ff1SMatthew G. Knepley } 12240be0ff1SMatthew G. Knepley 1239371c9d4SSatish Balay static PetscErrorCode TSView_DiscGrad(TS ts, PetscViewer viewer) { 12440be0ff1SMatthew G. Knepley PetscBool iascii; 12540be0ff1SMatthew G. Knepley 12640be0ff1SMatthew G. Knepley PetscFunctionBegin; 1279566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 128*48a46eb9SPierre Jolivet if (iascii) PetscCall(PetscViewerASCIIPrintf(viewer, " Discrete Gradients\n")); 12940be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 13040be0ff1SMatthew G. Knepley } 13140be0ff1SMatthew G. Knepley 1329371c9d4SSatish Balay static PetscErrorCode TSDiscGradIsGonzalez_DiscGrad(TS ts, PetscBool *gonzalez) { 133db4653c2SDaniel Finn TS_DiscGrad *dg = (TS_DiscGrad *)ts->data; 134db4653c2SDaniel Finn 135db4653c2SDaniel Finn PetscFunctionBegin; 136db4653c2SDaniel Finn *gonzalez = dg->gonzalez; 137db4653c2SDaniel Finn PetscFunctionReturn(0); 138db4653c2SDaniel Finn } 139db4653c2SDaniel Finn 1409371c9d4SSatish Balay static PetscErrorCode TSDiscGradUseGonzalez_DiscGrad(TS ts, PetscBool flg) { 141db4653c2SDaniel Finn TS_DiscGrad *dg = (TS_DiscGrad *)ts->data; 142db4653c2SDaniel Finn 143db4653c2SDaniel Finn PetscFunctionBegin; 144db4653c2SDaniel Finn dg->gonzalez = flg; 145db4653c2SDaniel Finn PetscFunctionReturn(0); 146db4653c2SDaniel Finn } 147db4653c2SDaniel Finn 1489371c9d4SSatish Balay static PetscErrorCode TSReset_DiscGrad(TS ts) { 14940be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *)ts->data; 15040be0ff1SMatthew G. Knepley 15140be0ff1SMatthew G. Knepley PetscFunctionBegin; 1529566063dSJacob Faibussowitsch PetscCall(VecDestroy(&dg->X)); 1539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&dg->X0)); 1549566063dSJacob Faibussowitsch PetscCall(VecDestroy(&dg->Xdot)); 15540be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 15640be0ff1SMatthew G. Knepley } 15740be0ff1SMatthew G. Knepley 1589371c9d4SSatish Balay static PetscErrorCode TSDestroy_DiscGrad(TS ts) { 15940be0ff1SMatthew G. Knepley DM dm; 16040be0ff1SMatthew G. Knepley 16140be0ff1SMatthew G. Knepley PetscFunctionBegin; 1629566063dSJacob Faibussowitsch PetscCall(TSReset_DiscGrad(ts)); 1639566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 16440be0ff1SMatthew G. Knepley if (dm) { 1659566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookRemove(dm, DMCoarsenHook_TSDiscGrad, DMRestrictHook_TSDiscGrad, ts)); 1669566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookRemove(dm, DMSubDomainHook_TSDiscGrad, DMSubDomainRestrictHook_TSDiscGrad, ts)); 16740be0ff1SMatthew G. Knepley } 1689566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 1699566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSDiscGradGetFormulation_C", NULL)); 1709566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSDiscGradSetFormulation_C", NULL)); 1712e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSDiscGradIsGonzalez_C", NULL)); 1722e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSDiscGradUseGonzalez_C", NULL)); 17340be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 17440be0ff1SMatthew G. Knepley } 17540be0ff1SMatthew G. Knepley 1769371c9d4SSatish Balay static PetscErrorCode TSInterpolate_DiscGrad(TS ts, PetscReal t, Vec X) { 17740be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *)ts->data; 17840be0ff1SMatthew G. Knepley PetscReal dt = t - ts->ptime; 17940be0ff1SMatthew G. Knepley 18040be0ff1SMatthew G. Knepley PetscFunctionBegin; 1819566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, dg->X)); 1829566063dSJacob Faibussowitsch PetscCall(VecWAXPY(X, dt, dg->Xdot, dg->X)); 18340be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 18440be0ff1SMatthew G. Knepley } 18540be0ff1SMatthew G. Knepley 1869371c9d4SSatish Balay static PetscErrorCode TSDiscGrad_SNESSolve(TS ts, Vec b, Vec x) { 18740be0ff1SMatthew G. Knepley SNES snes; 18840be0ff1SMatthew G. Knepley PetscInt nits, lits; 18940be0ff1SMatthew G. Knepley 19040be0ff1SMatthew G. Knepley PetscFunctionBegin; 1919566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 1929566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, b, x)); 1939566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &nits)); 1949566063dSJacob Faibussowitsch PetscCall(SNESGetLinearSolveIterations(snes, &lits)); 19540be0ff1SMatthew G. Knepley ts->snes_its += nits; 19640be0ff1SMatthew G. Knepley ts->ksp_its += lits; 19740be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 19840be0ff1SMatthew G. Knepley } 19940be0ff1SMatthew G. Knepley 2009371c9d4SSatish Balay static PetscErrorCode TSStep_DiscGrad(TS ts) { 20140be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *)ts->data; 20240be0ff1SMatthew G. Knepley TSAdapt adapt; 20340be0ff1SMatthew G. Knepley TSStepStatus status = TS_STEP_INCOMPLETE; 20440be0ff1SMatthew G. Knepley PetscInt rejections = 0; 20540be0ff1SMatthew G. Knepley PetscBool stageok, accept = PETSC_TRUE; 20640be0ff1SMatthew G. Knepley PetscReal next_time_step = ts->time_step; 20740be0ff1SMatthew G. Knepley 20840be0ff1SMatthew G. Knepley PetscFunctionBegin; 2099566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt)); 2109566063dSJacob Faibussowitsch if (!ts->steprollback) PetscCall(VecCopy(ts->vec_sol, dg->X0)); 21140be0ff1SMatthew G. Knepley 21240be0ff1SMatthew G. Knepley while (!ts->reason && status != TS_STEP_COMPLETE) { 21340be0ff1SMatthew G. Knepley PetscReal shift = 1 / (0.5 * ts->time_step); 21440be0ff1SMatthew G. Knepley 21540be0ff1SMatthew G. Knepley dg->stage_time = ts->ptime + 0.5 * ts->time_step; 21640be0ff1SMatthew G. Knepley 2179566063dSJacob Faibussowitsch PetscCall(VecCopy(dg->X0, dg->X)); 2189566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, dg->stage_time)); 2199566063dSJacob Faibussowitsch PetscCall(TSDiscGrad_SNESSolve(ts, NULL, dg->X)); 2209566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, dg->stage_time, 0, &dg->X)); 2219566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(adapt, ts, dg->stage_time, dg->X, &stageok)); 22240be0ff1SMatthew G. Knepley if (!stageok) goto reject_step; 22340be0ff1SMatthew G. Knepley 22440be0ff1SMatthew G. Knepley status = TS_STEP_PENDING; 2259566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(dg->Xdot, -shift, shift, 0, dg->X0, dg->X)); 2269566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vec_sol, ts->time_step, dg->Xdot)); 2279566063dSJacob Faibussowitsch PetscCall(TSAdaptChoose(adapt, ts, ts->time_step, NULL, &next_time_step, &accept)); 22840be0ff1SMatthew G. Knepley status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 22940be0ff1SMatthew G. Knepley if (!accept) { 2309566063dSJacob Faibussowitsch PetscCall(VecCopy(dg->X0, ts->vec_sol)); 23140be0ff1SMatthew G. Knepley ts->time_step = next_time_step; 23240be0ff1SMatthew G. Knepley goto reject_step; 23340be0ff1SMatthew G. Knepley } 23440be0ff1SMatthew G. Knepley ts->ptime += ts->time_step; 23540be0ff1SMatthew G. Knepley ts->time_step = next_time_step; 23640be0ff1SMatthew G. Knepley break; 23740be0ff1SMatthew G. Knepley 23840be0ff1SMatthew G. Knepley reject_step: 2399371c9d4SSatish Balay ts->reject++; 2409371c9d4SSatish Balay accept = PETSC_FALSE; 24140be0ff1SMatthew G. Knepley if (!ts->reason && ts->max_reject >= 0 && ++rejections > ts->max_reject) { 24240be0ff1SMatthew G. Knepley ts->reason = TS_DIVERGED_STEP_REJECTED; 24363a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections)); 24440be0ff1SMatthew G. Knepley } 24540be0ff1SMatthew G. Knepley } 24640be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 24740be0ff1SMatthew G. Knepley } 24840be0ff1SMatthew G. Knepley 2499371c9d4SSatish Balay static PetscErrorCode TSGetStages_DiscGrad(TS ts, PetscInt *ns, Vec **Y) { 25040be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *)ts->data; 25140be0ff1SMatthew G. Knepley 25240be0ff1SMatthew G. Knepley PetscFunctionBegin; 25340be0ff1SMatthew G. Knepley if (ns) *ns = 1; 25440be0ff1SMatthew G. Knepley if (Y) *Y = &(dg->X); 25540be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 25640be0ff1SMatthew G. Knepley } 25740be0ff1SMatthew G. Knepley 25840be0ff1SMatthew G. Knepley /* 25940be0ff1SMatthew G. Knepley This defines the nonlinear equation that is to be solved with SNES 26040be0ff1SMatthew G. Knepley G(U) = F[t0 + 0.5*dt, U, (U-U0)/dt] = 0 26140be0ff1SMatthew G. Knepley */ 262db4653c2SDaniel Finn 263db4653c2SDaniel Finn /* x = (x+x')/2 */ 264db4653c2SDaniel Finn /* NEED TO CALCULATE x_{n+1} from x and x_{n}*/ 2659371c9d4SSatish Balay static PetscErrorCode SNESTSFormFunction_DiscGrad(SNES snes, Vec x, Vec y, TS ts) { 26640be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *)ts->data; 267db4653c2SDaniel Finn PetscReal norm, shift = 1 / (0.5 * ts->time_step); 268db4653c2SDaniel Finn PetscInt n; 269db4653c2SDaniel Finn Vec X0, Xdot, Xp, Xdiff; 270db4653c2SDaniel Finn Mat S; 271db4653c2SDaniel Finn PetscScalar F = 0, F0 = 0, Gp; 272db4653c2SDaniel Finn Vec G, SgF; 27340be0ff1SMatthew G. Knepley DM dm, dmsave; 27440be0ff1SMatthew G. Knepley 27540be0ff1SMatthew G. Knepley PetscFunctionBegin; 2769566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 27740be0ff1SMatthew G. Knepley 2789566063dSJacob Faibussowitsch PetscCall(VecDuplicate(y, &Xp)); 2799566063dSJacob Faibussowitsch PetscCall(VecDuplicate(y, &Xdiff)); 2809566063dSJacob Faibussowitsch PetscCall(VecDuplicate(y, &SgF)); 2819566063dSJacob Faibussowitsch PetscCall(VecDuplicate(y, &G)); 282db4653c2SDaniel Finn 2839566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(y, &n)); 2849566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &S)); 2859566063dSJacob Faibussowitsch PetscCall(MatSetSizes(S, PETSC_DECIDE, PETSC_DECIDE, n, n)); 2869566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(S)); 2879566063dSJacob Faibussowitsch PetscCall(MatSetUp(S)); 2889566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(S, MAT_FINAL_ASSEMBLY)); 2899566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(S, MAT_FINAL_ASSEMBLY)); 290db4653c2SDaniel Finn 2919566063dSJacob Faibussowitsch PetscCall(TSDiscGradGetX0AndXdot(ts, dm, &X0, &Xdot)); 2929566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(Xdot, -shift, shift, 0, X0, x)); /* Xdot = shift (x - X0) */ 293db4653c2SDaniel Finn 2949566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(Xp, -1, 2, 0, X0, x)); /* Xp = 2*x - X0 + (0)*Xmid */ 2959566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(Xdiff, -1, 1, 0, X0, Xp)); /* Xdiff = xp - X0 + (0)*Xdiff */ 296db4653c2SDaniel Finn 297db4653c2SDaniel Finn if (dg->gonzalez) { 2989566063dSJacob Faibussowitsch PetscCall((*dg->Sfunc)(ts, dg->stage_time, x, S, dg->funcCtx)); 2999566063dSJacob Faibussowitsch PetscCall((*dg->Ffunc)(ts, dg->stage_time, Xp, &F, dg->funcCtx)); 3009566063dSJacob Faibussowitsch PetscCall((*dg->Ffunc)(ts, dg->stage_time, X0, &F0, dg->funcCtx)); 3019566063dSJacob Faibussowitsch PetscCall((*dg->Gfunc)(ts, dg->stage_time, x, G, dg->funcCtx)); 302db4653c2SDaniel Finn 303db4653c2SDaniel Finn /* Adding Extra Gonzalez Term */ 3049566063dSJacob Faibussowitsch PetscCall(VecDot(Xdiff, G, &Gp)); 3059566063dSJacob Faibussowitsch PetscCall(VecNorm(Xdiff, NORM_2, &norm)); 306db4653c2SDaniel Finn if (norm < PETSC_SQRT_MACHINE_EPSILON) { 307db4653c2SDaniel Finn Gp = 0; 308db4653c2SDaniel Finn } else { 309db4653c2SDaniel Finn /* Gp = (1/|xn+1 - xn|^2) * (F(xn+1) - F(xn) - Gp) */ 310db4653c2SDaniel Finn Gp = (F - F0 - Gp) / PetscSqr(norm); 311db4653c2SDaniel Finn } 3129566063dSJacob Faibussowitsch PetscCall(VecAXPY(G, Gp, Xdiff)); 3139566063dSJacob Faibussowitsch PetscCall(MatMult(S, G, SgF)); /* S*gradF */ 314db4653c2SDaniel Finn 315db4653c2SDaniel Finn } else { 3169566063dSJacob Faibussowitsch PetscCall((*dg->Sfunc)(ts, dg->stage_time, x, S, dg->funcCtx)); 3179566063dSJacob Faibussowitsch PetscCall((*dg->Gfunc)(ts, dg->stage_time, x, G, dg->funcCtx)); 318db4653c2SDaniel Finn 3199566063dSJacob Faibussowitsch PetscCall(MatMult(S, G, SgF)); /* Xdot = S*gradF */ 320db4653c2SDaniel Finn } 32140be0ff1SMatthew G. Knepley /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 32240be0ff1SMatthew G. Knepley dmsave = ts->dm; 32340be0ff1SMatthew G. Knepley ts->dm = dm; 3249566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(y, 1, -1, 0, Xdot, SgF)); 32540be0ff1SMatthew G. Knepley ts->dm = dmsave; 3269566063dSJacob Faibussowitsch PetscCall(TSDiscGradRestoreX0AndXdot(ts, dm, &X0, &Xdot)); 327db4653c2SDaniel Finn 3289566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Xp)); 3299566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Xdiff)); 3309566063dSJacob Faibussowitsch PetscCall(VecDestroy(&SgF)); 3319566063dSJacob Faibussowitsch PetscCall(VecDestroy(&G)); 3329566063dSJacob Faibussowitsch PetscCall(MatDestroy(&S)); 333db4653c2SDaniel Finn 33440be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 33540be0ff1SMatthew G. Knepley } 33640be0ff1SMatthew G. Knepley 3379371c9d4SSatish Balay static PetscErrorCode SNESTSFormJacobian_DiscGrad(SNES snes, Vec x, Mat A, Mat B, TS ts) { 33840be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *)ts->data; 33940be0ff1SMatthew G. Knepley PetscReal shift = 1 / (0.5 * ts->time_step); 34040be0ff1SMatthew G. Knepley Vec Xdot; 34140be0ff1SMatthew G. Knepley DM dm, dmsave; 34240be0ff1SMatthew G. Knepley 34340be0ff1SMatthew G. Knepley PetscFunctionBegin; 3449566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 34540be0ff1SMatthew G. Knepley /* Xdot has already been computed in SNESTSFormFunction_DiscGrad (SNES guarantees this) */ 3469566063dSJacob Faibussowitsch PetscCall(TSDiscGradGetX0AndXdot(ts, dm, NULL, &Xdot)); 34740be0ff1SMatthew G. Knepley 34840be0ff1SMatthew G. Knepley dmsave = ts->dm; 34940be0ff1SMatthew G. Knepley ts->dm = dm; 3509566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts, dg->stage_time, x, Xdot, shift, A, B, PETSC_FALSE)); 35140be0ff1SMatthew G. Knepley ts->dm = dmsave; 3529566063dSJacob Faibussowitsch PetscCall(TSDiscGradRestoreX0AndXdot(ts, dm, NULL, &Xdot)); 35340be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 35440be0ff1SMatthew G. Knepley } 35540be0ff1SMatthew G. Knepley 3569371c9d4SSatish Balay 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 *), void *ctx) { 35740be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *)ts->data; 35840be0ff1SMatthew G. Knepley 35940be0ff1SMatthew G. Knepley PetscFunctionBegin; 36040be0ff1SMatthew G. Knepley *Sfunc = dg->Sfunc; 36140be0ff1SMatthew G. Knepley *Ffunc = dg->Ffunc; 36240be0ff1SMatthew G. Knepley *Gfunc = dg->Gfunc; 36340be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 36440be0ff1SMatthew G. Knepley } 36540be0ff1SMatthew G. Knepley 3669371c9d4SSatish Balay 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 *), void *ctx) { 36740be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *)ts->data; 36840be0ff1SMatthew G. Knepley 36940be0ff1SMatthew G. Knepley PetscFunctionBegin; 37040be0ff1SMatthew G. Knepley dg->Sfunc = Sfunc; 37140be0ff1SMatthew G. Knepley dg->Ffunc = Ffunc; 37240be0ff1SMatthew G. Knepley dg->Gfunc = Gfunc; 373db4653c2SDaniel Finn dg->funcCtx = ctx; 37440be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 37540be0ff1SMatthew G. Knepley } 37640be0ff1SMatthew G. Knepley 37740be0ff1SMatthew G. Knepley /*MC 37840be0ff1SMatthew G. Knepley TSDISCGRAD - ODE solver using the discrete gradients version of the implicit midpoint method 37940be0ff1SMatthew G. Knepley 3802ceb51e8SPatrick Sanan Notes: 3812ceb51e8SPatrick Sanan This is the implicit midpoint rule, with an optional term that guarantees the discrete gradient property. This 38240be0ff1SMatthew G. Knepley timestepper applies to systems of the form 38340be0ff1SMatthew G. Knepley $ u_t = S(u) grad F(u) 38440be0ff1SMatthew G. Knepley where S(u) is a linear operator, and F is a functional of u. 38540be0ff1SMatthew G. Knepley 386f1a722f8SMatthew G. Knepley Level: intermediate 387f1a722f8SMatthew G. Knepley 388db781477SPatrick Sanan .seealso: `TSCreate()`, `TSSetType()`, `TS`, `TSDISCGRAD`, `TSDiscGradSetFormulation()` 38940be0ff1SMatthew G. Knepley M*/ 3909371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TSCreate_DiscGrad(TS ts) { 39140be0ff1SMatthew G. Knepley TS_DiscGrad *th; 39240be0ff1SMatthew G. Knepley 39340be0ff1SMatthew G. Knepley PetscFunctionBegin; 3949566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(DGCitation, &DGCite)); 39540be0ff1SMatthew G. Knepley ts->ops->reset = TSReset_DiscGrad; 39640be0ff1SMatthew G. Knepley ts->ops->destroy = TSDestroy_DiscGrad; 39740be0ff1SMatthew G. Knepley ts->ops->view = TSView_DiscGrad; 39840be0ff1SMatthew G. Knepley ts->ops->setfromoptions = TSSetFromOptions_DiscGrad; 39940be0ff1SMatthew G. Knepley ts->ops->setup = TSSetUp_DiscGrad; 40040be0ff1SMatthew G. Knepley ts->ops->step = TSStep_DiscGrad; 40140be0ff1SMatthew G. Knepley ts->ops->interpolate = TSInterpolate_DiscGrad; 40240be0ff1SMatthew G. Knepley ts->ops->getstages = TSGetStages_DiscGrad; 40340be0ff1SMatthew G. Knepley ts->ops->snesfunction = SNESTSFormFunction_DiscGrad; 40440be0ff1SMatthew G. Knepley ts->ops->snesjacobian = SNESTSFormJacobian_DiscGrad; 40540be0ff1SMatthew G. Knepley ts->default_adapt_type = TSADAPTNONE; 40640be0ff1SMatthew G. Knepley 40740be0ff1SMatthew G. Knepley ts->usessnes = PETSC_TRUE; 40840be0ff1SMatthew G. Knepley 4099566063dSJacob Faibussowitsch PetscCall(PetscNewLog(ts, &th)); 41040be0ff1SMatthew G. Knepley ts->data = (void *)th; 411db4653c2SDaniel Finn 412db4653c2SDaniel Finn th->gonzalez = PETSC_FALSE; 413db4653c2SDaniel Finn 4149566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSDiscGradGetFormulation_C", TSDiscGradGetFormulation_DiscGrad)); 4159566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSDiscGradSetFormulation_C", TSDiscGradSetFormulation_DiscGrad)); 4169566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSDiscGradIsGonzalez_C", TSDiscGradIsGonzalez_DiscGrad)); 4179566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSDiscGradUseGonzalez_C", TSDiscGradUseGonzalez_DiscGrad)); 41840be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 41940be0ff1SMatthew G. Knepley } 42040be0ff1SMatthew G. Knepley 42140be0ff1SMatthew G. Knepley /*@C 42240be0ff1SMatthew G. Knepley TSDiscGradGetFormulation - Get the construction method for S, F, and grad F from the formulation u_t = S grad F 42340be0ff1SMatthew G. Knepley 42440be0ff1SMatthew G. Knepley Not Collective 42540be0ff1SMatthew G. Knepley 42640be0ff1SMatthew G. Knepley Input Parameter: 42740be0ff1SMatthew G. Knepley . ts - timestepping context 42840be0ff1SMatthew G. Knepley 42940be0ff1SMatthew G. Knepley Output Parameters: 43040be0ff1SMatthew G. Knepley + Sfunc - constructor for the S matrix from the formulation 43140be0ff1SMatthew G. Knepley . Ffunc - functional F from the formulation 432f1a722f8SMatthew G. Knepley . Gfunc - constructor for the gradient of F from the formulation 433f1a722f8SMatthew G. Knepley - ctx - the user context 43440be0ff1SMatthew G. Knepley 43540be0ff1SMatthew G. Knepley Calling sequence of Sfunc: 43640be0ff1SMatthew G. Knepley $ PetscErrorCode func(TS ts, PetscReal time, Vec u, Mat S, void *) 43740be0ff1SMatthew G. Knepley 43840be0ff1SMatthew G. Knepley Calling sequence of Ffunc: 43940be0ff1SMatthew G. Knepley $ PetscErrorCode func(TS ts, PetscReal time, Vec u, PetscScalar *F, void *) 44040be0ff1SMatthew G. Knepley 44140be0ff1SMatthew G. Knepley Calling sequence of Gfunc: 44240be0ff1SMatthew G. Knepley $ PetscErrorCode func(TS ts, PetscReal time, Vec u, Vec G, void *) 44340be0ff1SMatthew G. Knepley 44440be0ff1SMatthew G. Knepley Level: intermediate 44540be0ff1SMatthew G. Knepley 446db781477SPatrick Sanan .seealso: `TSDiscGradSetFormulation()` 44740be0ff1SMatthew G. Knepley @*/ 4489371c9d4SSatish Balay PetscErrorCode TSDiscGradGetFormulation(TS ts, PetscErrorCode (**Sfunc)(TS, PetscReal, Vec, Mat, void *), PetscErrorCode (**Ffunc)(TS, PetscReal, Vec, PetscScalar *, void *), PetscErrorCode (**Gfunc)(TS, PetscReal, Vec, Vec, void *), void *ctx) { 44940be0ff1SMatthew G. Knepley PetscFunctionBegin; 45040be0ff1SMatthew G. Knepley PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 45140be0ff1SMatthew G. Knepley PetscValidPointer(Sfunc, 2); 45240be0ff1SMatthew G. Knepley PetscValidPointer(Ffunc, 3); 45340be0ff1SMatthew G. Knepley PetscValidPointer(Gfunc, 4); 454cac4c232SBarry 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)); 45540be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 45640be0ff1SMatthew G. Knepley } 45740be0ff1SMatthew G. Knepley 45840be0ff1SMatthew G. Knepley /*@C 45940be0ff1SMatthew G. Knepley TSDiscGradSetFormulation - Set the construction method for S, F, and grad F from the formulation u_t = S(u) grad F(u) 46040be0ff1SMatthew G. Knepley 46140be0ff1SMatthew G. Knepley Not Collective 46240be0ff1SMatthew G. Knepley 46340be0ff1SMatthew G. Knepley Input Parameters: 46440be0ff1SMatthew G. Knepley + ts - timestepping context 46540be0ff1SMatthew G. Knepley . Sfunc - constructor for the S matrix from the formulation 46640be0ff1SMatthew G. Knepley . Ffunc - functional F from the formulation 46740be0ff1SMatthew G. Knepley - Gfunc - constructor for the gradient of F from the formulation 46840be0ff1SMatthew G. Knepley Calling sequence of Sfunc: 46940be0ff1SMatthew G. Knepley $ PetscErrorCode func(TS ts, PetscReal time, Vec u, Mat S, void *) 47040be0ff1SMatthew G. Knepley 47140be0ff1SMatthew G. Knepley Calling sequence of Ffunc: 47240be0ff1SMatthew G. Knepley $ PetscErrorCode func(TS ts, PetscReal time, Vec u, PetscScalar *F, void *) 47340be0ff1SMatthew G. Knepley 47440be0ff1SMatthew G. Knepley Calling sequence of Gfunc: 47540be0ff1SMatthew G. Knepley $ PetscErrorCode func(TS ts, PetscReal time, Vec u, Vec G, void *) 47640be0ff1SMatthew G. Knepley 47740be0ff1SMatthew G. Knepley Level: Intermediate 47840be0ff1SMatthew G. Knepley 479db781477SPatrick Sanan .seealso: `TSDiscGradGetFormulation()` 48040be0ff1SMatthew G. Knepley @*/ 4819371c9d4SSatish Balay PetscErrorCode TSDiscGradSetFormulation(TS ts, PetscErrorCode (*Sfunc)(TS, PetscReal, Vec, Mat, void *), PetscErrorCode (*Ffunc)(TS, PetscReal, Vec, PetscScalar *, void *), PetscErrorCode (*Gfunc)(TS, PetscReal, Vec, Vec, void *), void *ctx) { 48240be0ff1SMatthew G. Knepley PetscFunctionBegin; 48340be0ff1SMatthew G. Knepley PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 48440be0ff1SMatthew G. Knepley PetscValidFunction(Sfunc, 2); 48540be0ff1SMatthew G. Knepley PetscValidFunction(Ffunc, 3); 48640be0ff1SMatthew G. Knepley PetscValidFunction(Gfunc, 4); 487cac4c232SBarry 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)); 488db4653c2SDaniel Finn PetscFunctionReturn(0); 489db4653c2SDaniel Finn } 490db4653c2SDaniel Finn 491db4653c2SDaniel Finn /*@ 492db4653c2SDaniel Finn TSDiscGradIsGonzalez - Checks flag for whether to use additional conservative terms in discrete gradient formulation. 493db4653c2SDaniel Finn 494db4653c2SDaniel Finn Not Collective 495db4653c2SDaniel Finn 496db4653c2SDaniel Finn Input Parameter: 497db4653c2SDaniel Finn . ts - timestepping context 498db4653c2SDaniel Finn 499db4653c2SDaniel Finn Output Parameter: 500db4653c2SDaniel Finn . gonzalez - PETSC_TRUE when using the Gonzalez term 501db4653c2SDaniel Finn 502db4653c2SDaniel Finn Level: Advanced 503db4653c2SDaniel Finn 504db781477SPatrick Sanan .seealso: `TSDiscGradUseGonzalez()`, `TSDISCGRAD` 505db4653c2SDaniel Finn @*/ 5069371c9d4SSatish Balay PetscErrorCode TSDiscGradIsGonzalez(TS ts, PetscBool *gonzalez) { 507db4653c2SDaniel Finn PetscFunctionBegin; 508db4653c2SDaniel Finn PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 509dadcf809SJacob Faibussowitsch PetscValidBoolPointer(gonzalez, 2); 510cac4c232SBarry Smith PetscUseMethod(ts, "TSDiscGradIsGonzalez_C", (TS, PetscBool *), (ts, gonzalez)); 511db4653c2SDaniel Finn PetscFunctionReturn(0); 512db4653c2SDaniel Finn } 513db4653c2SDaniel Finn 514db4653c2SDaniel Finn /*@ 515db4653c2SDaniel Finn TSDiscGradUseGonzalez - Sets discrete gradient formulation with or without additional conservative terms. Without flag, the discrete gradients timestepper is just backwards euler 516db4653c2SDaniel Finn 517db4653c2SDaniel Finn Not Collective 518db4653c2SDaniel Finn 519f1a722f8SMatthew G. Knepley Input Parameters: 520db4653c2SDaniel Finn + ts - timestepping context 521db4653c2SDaniel Finn - flg - PETSC_TRUE to use the Gonzalez term 522db4653c2SDaniel Finn 523db4653c2SDaniel Finn Options Database: 52467b8a455SSatish Balay . -ts_discgrad_gonzalez <flg> - use the Gonzalez term for the discrete gradient formulation 525db4653c2SDaniel Finn 526db4653c2SDaniel Finn Level: Intermediate 527db4653c2SDaniel Finn 528db781477SPatrick Sanan .seealso: `TSDISCGRAD` 529db4653c2SDaniel Finn @*/ 5309371c9d4SSatish Balay PetscErrorCode TSDiscGradUseGonzalez(TS ts, PetscBool flg) { 531db4653c2SDaniel Finn PetscFunctionBegin; 532db4653c2SDaniel Finn PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 533cac4c232SBarry Smith PetscTryMethod(ts, "TSDiscGradUseGonzalez_C", (TS, PetscBool), (ts, flg)); 53440be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 53540be0ff1SMatthew G. Knepley } 536