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 2740be0ff1SMatthew G. Knepley static PetscErrorCode TSDiscGradGetX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot) 2840be0ff1SMatthew G. Knepley { 2940be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *) ts->data; 3040be0ff1SMatthew G. Knepley 3140be0ff1SMatthew G. Knepley PetscFunctionBegin; 3240be0ff1SMatthew G. Knepley if (X0) { 339566063dSJacob Faibussowitsch if (dm && dm != ts->dm) PetscCall(DMGetNamedGlobalVector(dm, "TSDiscGrad_X0", X0)); 3440be0ff1SMatthew G. Knepley else {*X0 = ts->vec_sol;} 3540be0ff1SMatthew G. Knepley } 3640be0ff1SMatthew G. Knepley if (Xdot) { 379566063dSJacob Faibussowitsch if (dm && dm != ts->dm) PetscCall(DMGetNamedGlobalVector(dm, "TSDiscGrad_Xdot", Xdot)); 3840be0ff1SMatthew G. Knepley else {*Xdot = dg->Xdot;} 3940be0ff1SMatthew G. Knepley } 4040be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 4140be0ff1SMatthew G. Knepley } 4240be0ff1SMatthew G. Knepley 4340be0ff1SMatthew G. Knepley static PetscErrorCode TSDiscGradRestoreX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot) 4440be0ff1SMatthew G. Knepley { 4540be0ff1SMatthew G. Knepley PetscFunctionBegin; 4640be0ff1SMatthew G. Knepley if (X0) { 479566063dSJacob Faibussowitsch if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSDiscGrad_X0", X0)); 4840be0ff1SMatthew G. Knepley } 4940be0ff1SMatthew G. Knepley if (Xdot) { 509566063dSJacob Faibussowitsch if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSDiscGrad_Xdot", Xdot)); 5140be0ff1SMatthew G. Knepley } 5240be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 5340be0ff1SMatthew G. Knepley } 5440be0ff1SMatthew G. Knepley 5540be0ff1SMatthew G. Knepley static PetscErrorCode DMCoarsenHook_TSDiscGrad(DM fine, DM coarse, void *ctx) 5640be0ff1SMatthew G. Knepley { 5740be0ff1SMatthew G. Knepley PetscFunctionBegin; 5840be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 5940be0ff1SMatthew G. Knepley } 6040be0ff1SMatthew G. Knepley 6140be0ff1SMatthew G. Knepley static PetscErrorCode DMRestrictHook_TSDiscGrad(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx) 6240be0ff1SMatthew G. Knepley { 6340be0ff1SMatthew G. Knepley TS ts = (TS) ctx; 6440be0ff1SMatthew G. Knepley Vec X0, Xdot, X0_c, Xdot_c; 6540be0ff1SMatthew G. Knepley 6640be0ff1SMatthew G. Knepley PetscFunctionBegin; 679566063dSJacob Faibussowitsch PetscCall(TSDiscGradGetX0AndXdot(ts, fine, &X0, &Xdot)); 689566063dSJacob Faibussowitsch PetscCall(TSDiscGradGetX0AndXdot(ts, coarse, &X0_c, &Xdot_c)); 699566063dSJacob Faibussowitsch PetscCall(MatRestrict(restrct, X0, X0_c)); 709566063dSJacob Faibussowitsch PetscCall(MatRestrict(restrct, Xdot, Xdot_c)); 719566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(X0_c, rscale, X0_c)); 729566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Xdot_c, rscale, Xdot_c)); 739566063dSJacob Faibussowitsch PetscCall(TSDiscGradRestoreX0AndXdot(ts, fine, &X0, &Xdot)); 749566063dSJacob Faibussowitsch PetscCall(TSDiscGradRestoreX0AndXdot(ts, coarse, &X0_c, &Xdot_c)); 7540be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 7640be0ff1SMatthew G. Knepley } 7740be0ff1SMatthew G. Knepley 7840be0ff1SMatthew G. Knepley static PetscErrorCode DMSubDomainHook_TSDiscGrad(DM dm, DM subdm, void *ctx) 7940be0ff1SMatthew G. Knepley { 8040be0ff1SMatthew G. Knepley PetscFunctionBegin; 8140be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 8240be0ff1SMatthew G. Knepley } 8340be0ff1SMatthew G. Knepley 8440be0ff1SMatthew G. Knepley static PetscErrorCode DMSubDomainRestrictHook_TSDiscGrad(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx) 8540be0ff1SMatthew G. Knepley { 8640be0ff1SMatthew G. Knepley TS ts = (TS) ctx; 8740be0ff1SMatthew G. Knepley Vec X0, Xdot, X0_sub, Xdot_sub; 8840be0ff1SMatthew G. Knepley 8940be0ff1SMatthew G. Knepley PetscFunctionBegin; 909566063dSJacob Faibussowitsch PetscCall(TSDiscGradGetX0AndXdot(ts, dm, &X0, &Xdot)); 919566063dSJacob Faibussowitsch PetscCall(TSDiscGradGetX0AndXdot(ts, subdm, &X0_sub, &Xdot_sub)); 9240be0ff1SMatthew G. Knepley 939566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat, X0, X0_sub, INSERT_VALUES, SCATTER_FORWARD)); 949566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(gscat, X0, X0_sub, INSERT_VALUES, SCATTER_FORWARD)); 9540be0ff1SMatthew G. Knepley 969566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat, Xdot, Xdot_sub, INSERT_VALUES, SCATTER_FORWARD)); 979566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(gscat, Xdot, Xdot_sub, INSERT_VALUES, SCATTER_FORWARD)); 9840be0ff1SMatthew G. Knepley 999566063dSJacob Faibussowitsch PetscCall(TSDiscGradRestoreX0AndXdot(ts, dm, &X0, &Xdot)); 1009566063dSJacob Faibussowitsch PetscCall(TSDiscGradRestoreX0AndXdot(ts, subdm, &X0_sub, &Xdot_sub)); 10140be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 10240be0ff1SMatthew G. Knepley } 10340be0ff1SMatthew G. Knepley 10440be0ff1SMatthew G. Knepley static PetscErrorCode TSSetUp_DiscGrad(TS ts) 10540be0ff1SMatthew G. Knepley { 10640be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *) ts->data; 10740be0ff1SMatthew G. Knepley DM dm; 10840be0ff1SMatthew G. Knepley 10940be0ff1SMatthew G. Knepley PetscFunctionBegin; 1109566063dSJacob Faibussowitsch if (!dg->X) PetscCall(VecDuplicate(ts->vec_sol, &dg->X)); 1119566063dSJacob Faibussowitsch if (!dg->X0) PetscCall(VecDuplicate(ts->vec_sol, &dg->X0)); 1129566063dSJacob Faibussowitsch if (!dg->Xdot) PetscCall(VecDuplicate(ts->vec_sol, &dg->Xdot)); 11340be0ff1SMatthew G. Knepley 1149566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 1159566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_TSDiscGrad, DMRestrictHook_TSDiscGrad, ts)); 1169566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_TSDiscGrad, DMSubDomainRestrictHook_TSDiscGrad, ts)); 11740be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 11840be0ff1SMatthew G. Knepley } 11940be0ff1SMatthew G. Knepley 12040be0ff1SMatthew G. Knepley static PetscErrorCode TSSetFromOptions_DiscGrad(PetscOptionItems *PetscOptionsObject, TS ts) 12140be0ff1SMatthew G. Knepley { 122db4653c2SDaniel Finn TS_DiscGrad *dg = (TS_DiscGrad *) ts->data; 12340be0ff1SMatthew G. Knepley 12440be0ff1SMatthew G. Knepley PetscFunctionBegin; 125d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Discrete Gradients ODE solver options"); 12640be0ff1SMatthew G. Knepley { 1279566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_discgrad_gonzalez","Use Gonzalez term in discrete gradients formulation","TSDiscGradUseGonzalez",dg->gonzalez,&dg->gonzalez,NULL)); 12840be0ff1SMatthew G. Knepley } 129d0609cedSBarry Smith PetscOptionsHeadEnd(); 13040be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 13140be0ff1SMatthew G. Knepley } 13240be0ff1SMatthew G. Knepley 13340be0ff1SMatthew G. Knepley static PetscErrorCode TSView_DiscGrad(TS ts,PetscViewer viewer) 13440be0ff1SMatthew G. Knepley { 13540be0ff1SMatthew G. Knepley PetscBool iascii; 13640be0ff1SMatthew G. Knepley 13740be0ff1SMatthew G. Knepley PetscFunctionBegin; 1389566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii)); 13940be0ff1SMatthew G. Knepley if (iascii) { 1409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Discrete Gradients\n")); 14140be0ff1SMatthew G. Knepley } 14240be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 14340be0ff1SMatthew G. Knepley } 14440be0ff1SMatthew G. Knepley 145db4653c2SDaniel Finn static PetscErrorCode TSDiscGradIsGonzalez_DiscGrad(TS ts,PetscBool *gonzalez) 146db4653c2SDaniel Finn { 147db4653c2SDaniel Finn TS_DiscGrad *dg = (TS_DiscGrad*)ts->data; 148db4653c2SDaniel Finn 149db4653c2SDaniel Finn PetscFunctionBegin; 150db4653c2SDaniel Finn *gonzalez = dg->gonzalez; 151db4653c2SDaniel Finn PetscFunctionReturn(0); 152db4653c2SDaniel Finn } 153db4653c2SDaniel Finn 154db4653c2SDaniel Finn static PetscErrorCode TSDiscGradUseGonzalez_DiscGrad(TS ts,PetscBool flg) 155db4653c2SDaniel Finn { 156db4653c2SDaniel Finn TS_DiscGrad *dg = (TS_DiscGrad*)ts->data; 157db4653c2SDaniel Finn 158db4653c2SDaniel Finn PetscFunctionBegin; 159db4653c2SDaniel Finn dg->gonzalez = flg; 160db4653c2SDaniel Finn PetscFunctionReturn(0); 161db4653c2SDaniel Finn } 162db4653c2SDaniel Finn 16340be0ff1SMatthew G. Knepley static PetscErrorCode TSReset_DiscGrad(TS ts) 16440be0ff1SMatthew G. Knepley { 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)); 17140be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 17240be0ff1SMatthew G. Knepley } 17340be0ff1SMatthew G. Knepley 17440be0ff1SMatthew G. Knepley static PetscErrorCode TSDestroy_DiscGrad(TS ts) 17540be0ff1SMatthew G. Knepley { 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)); 188*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSDiscGradIsGonzalez_C",NULL)); 189*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSDiscGradUseGonzalez_C",NULL)); 19040be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 19140be0ff1SMatthew G. Knepley } 19240be0ff1SMatthew G. Knepley 19340be0ff1SMatthew G. Knepley static PetscErrorCode TSInterpolate_DiscGrad(TS ts, PetscReal t, Vec X) 19440be0ff1SMatthew G. Knepley { 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)); 20140be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 20240be0ff1SMatthew G. Knepley } 20340be0ff1SMatthew G. Knepley 20440be0ff1SMatthew G. Knepley static PetscErrorCode TSDiscGrad_SNESSolve(TS ts, Vec b, Vec x) 20540be0ff1SMatthew G. Knepley { 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; 21640be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 21740be0ff1SMatthew G. Knepley } 21840be0ff1SMatthew G. Knepley 21940be0ff1SMatthew G. Knepley static PetscErrorCode TSStep_DiscGrad(TS ts) 22040be0ff1SMatthew G. Knepley { 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: 25940be0ff1SMatthew G. Knepley ts->reject++; accept = PETSC_FALSE; 26040be0ff1SMatthew G. Knepley if (!ts->reason && ts->max_reject >= 0 && ++rejections > ts->max_reject) { 26140be0ff1SMatthew G. Knepley ts->reason = TS_DIVERGED_STEP_REJECTED; 26263a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections)); 26340be0ff1SMatthew G. Knepley } 26440be0ff1SMatthew G. Knepley } 26540be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 26640be0ff1SMatthew G. Knepley } 26740be0ff1SMatthew G. Knepley 26840be0ff1SMatthew G. Knepley static PetscErrorCode TSGetStages_DiscGrad(TS ts, PetscInt *ns, Vec **Y) 26940be0ff1SMatthew G. Knepley { 27040be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *) ts->data; 27140be0ff1SMatthew G. Knepley 27240be0ff1SMatthew G. Knepley PetscFunctionBegin; 27340be0ff1SMatthew G. Knepley if (ns) *ns = 1; 27440be0ff1SMatthew G. Knepley if (Y) *Y = &(dg->X); 27540be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 27640be0ff1SMatthew G. Knepley } 27740be0ff1SMatthew G. Knepley 27840be0ff1SMatthew G. Knepley /* 27940be0ff1SMatthew G. Knepley This defines the nonlinear equation that is to be solved with SNES 28040be0ff1SMatthew G. Knepley G(U) = F[t0 + 0.5*dt, U, (U-U0)/dt] = 0 28140be0ff1SMatthew G. Knepley */ 282db4653c2SDaniel Finn 283db4653c2SDaniel Finn /* x = (x+x')/2 */ 284db4653c2SDaniel Finn /* NEED TO CALCULATE x_{n+1} from x and x_{n}*/ 28540be0ff1SMatthew G. Knepley static PetscErrorCode SNESTSFormFunction_DiscGrad(SNES snes, Vec x, Vec y, TS ts) 28640be0ff1SMatthew G. Knepley { 287db4653c2SDaniel Finn 28840be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *) ts->data; 289db4653c2SDaniel Finn PetscReal norm, shift = 1/(0.5*ts->time_step); 290db4653c2SDaniel Finn PetscInt n; 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)); 29940be0ff1SMatthew G. Knepley 3009566063dSJacob Faibussowitsch PetscCall(VecDuplicate(y, &Xp)); 3019566063dSJacob Faibussowitsch PetscCall(VecDuplicate(y, &Xdiff)); 3029566063dSJacob Faibussowitsch PetscCall(VecDuplicate(y, &SgF)); 3039566063dSJacob Faibussowitsch PetscCall(VecDuplicate(y, &G)); 304db4653c2SDaniel Finn 3059566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(y, &n)); 3069566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &S)); 3079566063dSJacob Faibussowitsch PetscCall(MatSetSizes(S, PETSC_DECIDE, PETSC_DECIDE, n, n)); 3089566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(S)); 3099566063dSJacob Faibussowitsch PetscCall(MatSetUp(S)); 3109566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(S,MAT_FINAL_ASSEMBLY)); 3119566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(S,MAT_FINAL_ASSEMBLY)); 312db4653c2SDaniel Finn 3139566063dSJacob Faibussowitsch PetscCall(TSDiscGradGetX0AndXdot(ts, dm, &X0, &Xdot)); 3149566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(Xdot, -shift, shift, 0, X0, x)); /* Xdot = shift (x - X0) */ 315db4653c2SDaniel Finn 3169566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(Xp, -1, 2, 0, X0, x)); /* Xp = 2*x - X0 + (0)*Xmid */ 3179566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(Xdiff, -1, 1, 0, X0, Xp)); /* Xdiff = xp - X0 + (0)*Xdiff */ 318db4653c2SDaniel Finn 319db4653c2SDaniel Finn if (dg->gonzalez) { 3209566063dSJacob Faibussowitsch PetscCall((*dg->Sfunc)(ts, dg->stage_time, x , S, dg->funcCtx)); 3219566063dSJacob Faibussowitsch PetscCall((*dg->Ffunc)(ts, dg->stage_time, Xp, &F, dg->funcCtx)); 3229566063dSJacob Faibussowitsch PetscCall((*dg->Ffunc)(ts, dg->stage_time, X0, &F0, dg->funcCtx)); 3239566063dSJacob Faibussowitsch PetscCall((*dg->Gfunc)(ts, dg->stage_time, x , G, dg->funcCtx)); 324db4653c2SDaniel Finn 325db4653c2SDaniel Finn /* Adding Extra Gonzalez Term */ 3269566063dSJacob Faibussowitsch PetscCall(VecDot(Xdiff, G, &Gp)); 3279566063dSJacob Faibussowitsch PetscCall(VecNorm(Xdiff, NORM_2, &norm)); 328db4653c2SDaniel Finn if (norm < PETSC_SQRT_MACHINE_EPSILON) { 329db4653c2SDaniel Finn Gp = 0; 330db4653c2SDaniel Finn } else { 331db4653c2SDaniel Finn /* Gp = (1/|xn+1 - xn|^2) * (F(xn+1) - F(xn) - Gp) */ 332db4653c2SDaniel Finn Gp = (F - F0 - Gp)/PetscSqr(norm); 333db4653c2SDaniel Finn } 3349566063dSJacob Faibussowitsch PetscCall(VecAXPY(G, Gp, Xdiff)); 3359566063dSJacob Faibussowitsch PetscCall(MatMult(S, G , SgF)); /* S*gradF */ 336db4653c2SDaniel Finn 337db4653c2SDaniel Finn } else { 3389566063dSJacob Faibussowitsch PetscCall((*dg->Sfunc)(ts, dg->stage_time, x, S, dg->funcCtx)); 3399566063dSJacob Faibussowitsch PetscCall((*dg->Gfunc)(ts, dg->stage_time, x, G, dg->funcCtx)); 340db4653c2SDaniel Finn 3419566063dSJacob Faibussowitsch PetscCall(MatMult(S, G , SgF)); /* Xdot = S*gradF */ 342db4653c2SDaniel Finn } 34340be0ff1SMatthew G. Knepley /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 34440be0ff1SMatthew G. Knepley dmsave = ts->dm; 34540be0ff1SMatthew G. Knepley ts->dm = dm; 3469566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(y, 1, -1, 0, Xdot, SgF)); 34740be0ff1SMatthew G. Knepley ts->dm = dmsave; 3489566063dSJacob Faibussowitsch PetscCall(TSDiscGradRestoreX0AndXdot(ts, dm, &X0, &Xdot)); 349db4653c2SDaniel Finn 3509566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Xp)); 3519566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Xdiff)); 3529566063dSJacob Faibussowitsch PetscCall(VecDestroy(&SgF)); 3539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&G)); 3549566063dSJacob Faibussowitsch PetscCall(MatDestroy(&S)); 355db4653c2SDaniel Finn 35640be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 35740be0ff1SMatthew G. Knepley } 35840be0ff1SMatthew G. Knepley 35940be0ff1SMatthew G. Knepley static PetscErrorCode SNESTSFormJacobian_DiscGrad(SNES snes, Vec x, Mat A, Mat B, TS ts) 36040be0ff1SMatthew G. Knepley { 36140be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *) ts->data; 36240be0ff1SMatthew G. Knepley PetscReal shift = 1/(0.5*ts->time_step); 36340be0ff1SMatthew G. Knepley Vec Xdot; 36440be0ff1SMatthew G. Knepley DM dm,dmsave; 36540be0ff1SMatthew G. Knepley 36640be0ff1SMatthew G. Knepley PetscFunctionBegin; 3679566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 36840be0ff1SMatthew G. Knepley /* Xdot has already been computed in SNESTSFormFunction_DiscGrad (SNES guarantees this) */ 3699566063dSJacob Faibussowitsch PetscCall(TSDiscGradGetX0AndXdot(ts, dm, NULL, &Xdot)); 37040be0ff1SMatthew G. Knepley 37140be0ff1SMatthew G. Knepley dmsave = ts->dm; 37240be0ff1SMatthew G. Knepley ts->dm = dm; 3739566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts, dg->stage_time, x, Xdot, shift, A, B, PETSC_FALSE)); 37440be0ff1SMatthew G. Knepley ts->dm = dmsave; 3759566063dSJacob Faibussowitsch PetscCall(TSDiscGradRestoreX0AndXdot(ts, dm, NULL, &Xdot)); 37640be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 37740be0ff1SMatthew G. Knepley } 37840be0ff1SMatthew G. Knepley 379db4653c2SDaniel Finn 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) 38040be0ff1SMatthew G. Knepley { 38140be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *) ts->data; 38240be0ff1SMatthew G. Knepley 38340be0ff1SMatthew G. Knepley PetscFunctionBegin; 38440be0ff1SMatthew G. Knepley *Sfunc = dg->Sfunc; 38540be0ff1SMatthew G. Knepley *Ffunc = dg->Ffunc; 38640be0ff1SMatthew G. Knepley *Gfunc = dg->Gfunc; 38740be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 38840be0ff1SMatthew G. Knepley } 38940be0ff1SMatthew G. Knepley 390db4653c2SDaniel Finn 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) 39140be0ff1SMatthew G. Knepley { 39240be0ff1SMatthew G. Knepley TS_DiscGrad *dg = (TS_DiscGrad *) ts->data; 39340be0ff1SMatthew G. Knepley 39440be0ff1SMatthew G. Knepley PetscFunctionBegin; 39540be0ff1SMatthew G. Knepley dg->Sfunc = Sfunc; 39640be0ff1SMatthew G. Knepley dg->Ffunc = Ffunc; 39740be0ff1SMatthew G. Knepley dg->Gfunc = Gfunc; 398db4653c2SDaniel Finn dg->funcCtx = ctx; 39940be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 40040be0ff1SMatthew G. Knepley } 40140be0ff1SMatthew G. Knepley 40240be0ff1SMatthew G. Knepley /*MC 40340be0ff1SMatthew G. Knepley TSDISCGRAD - ODE solver using the discrete gradients version of the implicit midpoint method 40440be0ff1SMatthew G. Knepley 4052ceb51e8SPatrick Sanan Notes: 4062ceb51e8SPatrick Sanan This is the implicit midpoint rule, with an optional term that guarantees the discrete gradient property. This 40740be0ff1SMatthew G. Knepley timestepper applies to systems of the form 40840be0ff1SMatthew G. Knepley $ u_t = S(u) grad F(u) 40940be0ff1SMatthew G. Knepley where S(u) is a linear operator, and F is a functional of u. 41040be0ff1SMatthew G. Knepley 411f1a722f8SMatthew G. Knepley Level: intermediate 412f1a722f8SMatthew G. Knepley 413db781477SPatrick Sanan .seealso: `TSCreate()`, `TSSetType()`, `TS`, `TSDISCGRAD`, `TSDiscGradSetFormulation()` 41440be0ff1SMatthew G. Knepley M*/ 41540be0ff1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSCreate_DiscGrad(TS ts) 41640be0ff1SMatthew G. Knepley { 41740be0ff1SMatthew G. Knepley TS_DiscGrad *th; 41840be0ff1SMatthew G. Knepley 41940be0ff1SMatthew G. Knepley PetscFunctionBegin; 4209566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(DGCitation, &DGCite)); 42140be0ff1SMatthew G. Knepley ts->ops->reset = TSReset_DiscGrad; 42240be0ff1SMatthew G. Knepley ts->ops->destroy = TSDestroy_DiscGrad; 42340be0ff1SMatthew G. Knepley ts->ops->view = TSView_DiscGrad; 42440be0ff1SMatthew G. Knepley ts->ops->setfromoptions = TSSetFromOptions_DiscGrad; 42540be0ff1SMatthew G. Knepley ts->ops->setup = TSSetUp_DiscGrad; 42640be0ff1SMatthew G. Knepley ts->ops->step = TSStep_DiscGrad; 42740be0ff1SMatthew G. Knepley ts->ops->interpolate = TSInterpolate_DiscGrad; 42840be0ff1SMatthew G. Knepley ts->ops->getstages = TSGetStages_DiscGrad; 42940be0ff1SMatthew G. Knepley ts->ops->snesfunction = SNESTSFormFunction_DiscGrad; 43040be0ff1SMatthew G. Knepley ts->ops->snesjacobian = SNESTSFormJacobian_DiscGrad; 43140be0ff1SMatthew G. Knepley ts->default_adapt_type = TSADAPTNONE; 43240be0ff1SMatthew G. Knepley 43340be0ff1SMatthew G. Knepley ts->usessnes = PETSC_TRUE; 43440be0ff1SMatthew G. Knepley 4359566063dSJacob Faibussowitsch PetscCall(PetscNewLog(ts,&th)); 43640be0ff1SMatthew G. Knepley ts->data = (void*)th; 437db4653c2SDaniel Finn 438db4653c2SDaniel Finn th->gonzalez = PETSC_FALSE; 439db4653c2SDaniel Finn 4409566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSDiscGradGetFormulation_C",TSDiscGradGetFormulation_DiscGrad)); 4419566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSDiscGradSetFormulation_C",TSDiscGradSetFormulation_DiscGrad)); 4429566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSDiscGradIsGonzalez_C",TSDiscGradIsGonzalez_DiscGrad)); 4439566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSDiscGradUseGonzalez_C",TSDiscGradUseGonzalez_DiscGrad)); 44440be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 44540be0ff1SMatthew G. Knepley } 44640be0ff1SMatthew G. Knepley 44740be0ff1SMatthew G. Knepley /*@C 44840be0ff1SMatthew G. Knepley TSDiscGradGetFormulation - Get the construction method for S, F, and grad F from the formulation u_t = S grad F 44940be0ff1SMatthew G. Knepley 45040be0ff1SMatthew G. Knepley Not Collective 45140be0ff1SMatthew G. Knepley 45240be0ff1SMatthew G. Knepley Input Parameter: 45340be0ff1SMatthew G. Knepley . ts - timestepping context 45440be0ff1SMatthew G. Knepley 45540be0ff1SMatthew G. Knepley Output Parameters: 45640be0ff1SMatthew G. Knepley + Sfunc - constructor for the S matrix from the formulation 45740be0ff1SMatthew G. Knepley . Ffunc - functional F from the formulation 458f1a722f8SMatthew G. Knepley . Gfunc - constructor for the gradient of F from the formulation 459f1a722f8SMatthew G. Knepley - ctx - the user context 46040be0ff1SMatthew G. Knepley 46140be0ff1SMatthew G. Knepley Calling sequence of Sfunc: 46240be0ff1SMatthew G. Knepley $ PetscErrorCode func(TS ts, PetscReal time, Vec u, Mat S, void *) 46340be0ff1SMatthew G. Knepley 46440be0ff1SMatthew G. Knepley Calling sequence of Ffunc: 46540be0ff1SMatthew G. Knepley $ PetscErrorCode func(TS ts, PetscReal time, Vec u, PetscScalar *F, void *) 46640be0ff1SMatthew G. Knepley 46740be0ff1SMatthew G. Knepley Calling sequence of Gfunc: 46840be0ff1SMatthew G. Knepley $ PetscErrorCode func(TS ts, PetscReal time, Vec u, Vec G, void *) 46940be0ff1SMatthew G. Knepley 47040be0ff1SMatthew G. Knepley Level: intermediate 47140be0ff1SMatthew G. Knepley 472db781477SPatrick Sanan .seealso: `TSDiscGradSetFormulation()` 47340be0ff1SMatthew G. Knepley @*/ 474db4653c2SDaniel Finn 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) 47540be0ff1SMatthew G. Knepley { 47640be0ff1SMatthew G. Knepley PetscFunctionBegin; 47740be0ff1SMatthew G. Knepley PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 47840be0ff1SMatthew G. Knepley PetscValidPointer(Sfunc, 2); 47940be0ff1SMatthew G. Knepley PetscValidPointer(Ffunc, 3); 48040be0ff1SMatthew G. Knepley PetscValidPointer(Gfunc, 4); 481cac4c232SBarry 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)); 48240be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 48340be0ff1SMatthew G. Knepley } 48440be0ff1SMatthew G. Knepley 48540be0ff1SMatthew G. Knepley /*@C 48640be0ff1SMatthew G. Knepley TSDiscGradSetFormulation - Set the construction method for S, F, and grad F from the formulation u_t = S(u) grad F(u) 48740be0ff1SMatthew G. Knepley 48840be0ff1SMatthew G. Knepley Not Collective 48940be0ff1SMatthew G. Knepley 49040be0ff1SMatthew G. Knepley Input Parameters: 49140be0ff1SMatthew G. Knepley + ts - timestepping context 49240be0ff1SMatthew G. Knepley . Sfunc - constructor for the S matrix from the formulation 49340be0ff1SMatthew G. Knepley . Ffunc - functional F from the formulation 49440be0ff1SMatthew G. Knepley - Gfunc - constructor for the gradient of F from the formulation 49540be0ff1SMatthew G. Knepley Calling sequence of Sfunc: 49640be0ff1SMatthew G. Knepley $ PetscErrorCode func(TS ts, PetscReal time, Vec u, Mat S, void *) 49740be0ff1SMatthew G. Knepley 49840be0ff1SMatthew G. Knepley Calling sequence of Ffunc: 49940be0ff1SMatthew G. Knepley $ PetscErrorCode func(TS ts, PetscReal time, Vec u, PetscScalar *F, void *) 50040be0ff1SMatthew G. Knepley 50140be0ff1SMatthew G. Knepley Calling sequence of Gfunc: 50240be0ff1SMatthew G. Knepley $ PetscErrorCode func(TS ts, PetscReal time, Vec u, Vec G, void *) 50340be0ff1SMatthew G. Knepley 50440be0ff1SMatthew G. Knepley Level: Intermediate 50540be0ff1SMatthew G. Knepley 506db781477SPatrick Sanan .seealso: `TSDiscGradGetFormulation()` 50740be0ff1SMatthew G. Knepley @*/ 50840be0ff1SMatthew G. Knepley 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) 50940be0ff1SMatthew G. Knepley { 51040be0ff1SMatthew G. Knepley PetscFunctionBegin; 51140be0ff1SMatthew G. Knepley PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 51240be0ff1SMatthew G. Knepley PetscValidFunction(Sfunc, 2); 51340be0ff1SMatthew G. Knepley PetscValidFunction(Ffunc, 3); 51440be0ff1SMatthew G. Knepley PetscValidFunction(Gfunc, 4); 515cac4c232SBarry 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)); 516db4653c2SDaniel Finn PetscFunctionReturn(0); 517db4653c2SDaniel Finn } 518db4653c2SDaniel Finn 519db4653c2SDaniel Finn /*@ 520db4653c2SDaniel Finn TSDiscGradIsGonzalez - Checks flag for whether to use additional conservative terms in discrete gradient formulation. 521db4653c2SDaniel Finn 522db4653c2SDaniel Finn Not Collective 523db4653c2SDaniel Finn 524db4653c2SDaniel Finn Input Parameter: 525db4653c2SDaniel Finn . ts - timestepping context 526db4653c2SDaniel Finn 527db4653c2SDaniel Finn Output Parameter: 528db4653c2SDaniel Finn . gonzalez - PETSC_TRUE when using the Gonzalez term 529db4653c2SDaniel Finn 530db4653c2SDaniel Finn Level: Advanced 531db4653c2SDaniel Finn 532db781477SPatrick Sanan .seealso: `TSDiscGradUseGonzalez()`, `TSDISCGRAD` 533db4653c2SDaniel Finn @*/ 534db4653c2SDaniel Finn PetscErrorCode TSDiscGradIsGonzalez(TS ts,PetscBool *gonzalez) 535db4653c2SDaniel Finn { 536db4653c2SDaniel Finn PetscFunctionBegin; 537db4653c2SDaniel Finn PetscValidHeaderSpecific(ts,TS_CLASSID,1); 538dadcf809SJacob Faibussowitsch PetscValidBoolPointer(gonzalez,2); 539cac4c232SBarry Smith PetscUseMethod(ts,"TSDiscGradIsGonzalez_C",(TS,PetscBool*),(ts,gonzalez)); 540db4653c2SDaniel Finn PetscFunctionReturn(0); 541db4653c2SDaniel Finn } 542db4653c2SDaniel Finn 543db4653c2SDaniel Finn /*@ 544db4653c2SDaniel Finn TSDiscGradUseGonzalez - Sets discrete gradient formulation with or without additional conservative terms. Without flag, the discrete gradients timestepper is just backwards euler 545db4653c2SDaniel Finn 546db4653c2SDaniel Finn Not Collective 547db4653c2SDaniel Finn 548f1a722f8SMatthew G. Knepley Input Parameters: 549db4653c2SDaniel Finn + ts - timestepping context 550db4653c2SDaniel Finn - flg - PETSC_TRUE to use the Gonzalez term 551db4653c2SDaniel Finn 552db4653c2SDaniel Finn Options Database: 55367b8a455SSatish Balay . -ts_discgrad_gonzalez <flg> - use the Gonzalez term for the discrete gradient formulation 554db4653c2SDaniel Finn 555db4653c2SDaniel Finn Level: Intermediate 556db4653c2SDaniel Finn 557db781477SPatrick Sanan .seealso: `TSDISCGRAD` 558db4653c2SDaniel Finn @*/ 559db4653c2SDaniel Finn PetscErrorCode TSDiscGradUseGonzalez(TS ts,PetscBool flg) 560db4653c2SDaniel Finn { 561db4653c2SDaniel Finn PetscFunctionBegin; 562db4653c2SDaniel Finn PetscValidHeaderSpecific(ts,TS_CLASSID,1); 563cac4c232SBarry Smith PetscTryMethod(ts,"TSDiscGradUseGonzalez_C",(TS,PetscBool),(ts,flg)); 56440be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 56540be0ff1SMatthew G. Knepley } 566