12d3f70b5SBarry Smith /* 2fb4a63b6SLois Curfman McInnes Code for Timestepping with implicit backwards Euler. 32d3f70b5SBarry Smith */ 4af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 52d3f70b5SBarry Smith 6c558785fSJames Wright #define TSADAPTTSPSEUDO "tspseudo" 7c558785fSJames Wright 82d3f70b5SBarry Smith typedef struct { 92d3f70b5SBarry Smith Vec update; /* work vector where new solution is formed */ 102d3f70b5SBarry Smith Vec func; /* work vector where F(t[i],u[i]) is stored */ 116f2d6a7bSJed Brown Vec xdot; /* work vector for time derivative of state */ 122d3f70b5SBarry Smith 132d3f70b5SBarry Smith /* information used for Pseudo-timestepping */ 142d3f70b5SBarry Smith 156849ba73SBarry Smith PetscErrorCode (*dt)(TS, PetscReal *, void *); /* compute next timestep, and related context */ 162d3f70b5SBarry Smith void *dtctx; 17ace3abfcSBarry Smith PetscErrorCode (*verify)(TS, Vec, void *, PetscReal *, PetscBool *); /* verify previous timestep and related context */ 187bf11e45SBarry Smith void *verifyctx; 192d3f70b5SBarry Smith 20cdbf8f93SLisandro Dalcin PetscReal fnorm_initial, fnorm; /* original and current norm of F(u) */ 2187828ca2SBarry Smith PetscReal fnorm_previous; 2228aa8177SBarry Smith 23cdbf8f93SLisandro Dalcin PetscReal dt_initial; /* initial time-step */ 2487828ca2SBarry Smith PetscReal dt_increment; /* scaling that dt is incremented each time-step */ 2586534af1SJed Brown PetscReal dt_max; /* maximum time step */ 26ace3abfcSBarry Smith PetscBool increment_dt_from_initial_dt; 273118ae5eSBarry Smith PetscReal fatol, frtol; 28eb636060SBarry Smith 29eb636060SBarry Smith PetscObjectState Xstate; /* state of input vector to TSComputeIFunction() with 0 Xdot */ 307bf11e45SBarry Smith } TS_Pseudo; 312d3f70b5SBarry Smith 32*0fe2f5c2SJames Wright /*@C 33*0fe2f5c2SJames Wright TSPseudoComputeFunction - Compute nonlinear residual for pseudo-timestepping 34*0fe2f5c2SJames Wright 35*0fe2f5c2SJames Wright This computes the residual for $\dot U = 0$, i.e. $F(U, 0)$ for the IFunction. 36*0fe2f5c2SJames Wright 37*0fe2f5c2SJames Wright Collective 38*0fe2f5c2SJames Wright 39*0fe2f5c2SJames Wright Input Parameters: 40*0fe2f5c2SJames Wright + ts - the timestep context 41*0fe2f5c2SJames Wright - solution - the solution vector 42*0fe2f5c2SJames Wright 43*0fe2f5c2SJames Wright Output Parameter: 44*0fe2f5c2SJames Wright + residual - the nonlinear residual 45*0fe2f5c2SJames Wright - fnorm - the norm of the nonlinear residual 46*0fe2f5c2SJames Wright 47*0fe2f5c2SJames Wright Level: advanced 48*0fe2f5c2SJames Wright 49*0fe2f5c2SJames Wright Note: 50*0fe2f5c2SJames Wright `TSPSEUDO` records the nonlinear residual and the `solution` vector used to generate it. If given the same `solution` vector (as determined by the vectors `PetscObjectState`), this function will return those recorded values. 51*0fe2f5c2SJames Wright 52*0fe2f5c2SJames Wright This would be used in a custom adaptive timestepping implementation that needs access to the residual, but reuses the calculation done by `TSPSEUDO` by default. 53*0fe2f5c2SJames Wright 54*0fe2f5c2SJames Wright .seealso: [](ch_ts), `TSPSEUDO` 55*0fe2f5c2SJames Wright @*/ 56*0fe2f5c2SJames Wright PetscErrorCode TSPseudoComputeFunction(TS ts, Vec solution, Vec *residual, PetscReal *fnorm) 57*0fe2f5c2SJames Wright { 58*0fe2f5c2SJames Wright TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 59*0fe2f5c2SJames Wright PetscObjectState Xstate; 60*0fe2f5c2SJames Wright 61*0fe2f5c2SJames Wright PetscFunctionBegin; 62*0fe2f5c2SJames Wright PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 63*0fe2f5c2SJames Wright PetscValidHeaderSpecific(solution, VEC_CLASSID, 2); 64*0fe2f5c2SJames Wright if (residual) PetscValidHeaderSpecific(*residual, VEC_CLASSID, 3); 65*0fe2f5c2SJames Wright if (fnorm) PetscAssertPointer(fnorm, 4); 66*0fe2f5c2SJames Wright 67*0fe2f5c2SJames Wright PetscCall(PetscObjectStateGet((PetscObject)solution, &Xstate)); 68*0fe2f5c2SJames Wright if (Xstate != pseudo->Xstate || pseudo->fnorm < 0) { 69*0fe2f5c2SJames Wright PetscCall(VecZeroEntries(pseudo->xdot)); 70*0fe2f5c2SJames Wright PetscCall(TSComputeIFunction(ts, ts->ptime, solution, pseudo->xdot, pseudo->func, PETSC_FALSE)); 71*0fe2f5c2SJames Wright pseudo->Xstate = Xstate; 72*0fe2f5c2SJames Wright PetscCall(VecNorm(pseudo->func, NORM_2, &pseudo->fnorm)); 73*0fe2f5c2SJames Wright } 74*0fe2f5c2SJames Wright if (residual) *residual = pseudo->func; 75*0fe2f5c2SJames Wright if (fnorm) *fnorm = pseudo->fnorm; 76*0fe2f5c2SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 77*0fe2f5c2SJames Wright } 78*0fe2f5c2SJames Wright 79d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_Pseudo(TS ts) 80d71ae5a4SJacob Faibussowitsch { 81277b19d0SLisandro Dalcin TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 82be5899b3SLisandro Dalcin PetscInt nits, lits, reject; 83cdbf8f93SLisandro Dalcin PetscBool stepok; 84be5899b3SLisandro Dalcin PetscReal next_time_step = ts->time_step; 85c558785fSJames Wright TSAdapt adapt; 862d3f70b5SBarry Smith 873a40ed3dSBarry Smith PetscFunctionBegin; 88eb636060SBarry Smith if (ts->steps == 0) { 89eb636060SBarry Smith pseudo->dt_initial = ts->time_step; 90eb636060SBarry Smith PetscCall(VecCopy(ts->vec_sol, pseudo->update)); /* in all future updates pseudo->update already contains the current time solution */ 91eb636060SBarry Smith } 92cdbf8f93SLisandro Dalcin for (reject = 0; reject < ts->max_reject; reject++, ts->reject++) { 93cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 94c558785fSJames Wright if (reject) PetscCall(VecCopy(ts->vec_sol0, pseudo->update)); /* need to copy the solution because we got a rejection and pseudo->update contains intermediate computation */ 959566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, ts->ptime + ts->time_step)); 969566063dSJacob Faibussowitsch PetscCall(SNESSolve(ts->snes, NULL, pseudo->update)); 979566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(ts->snes, &nits)); 989566063dSJacob Faibussowitsch PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits)); 999371c9d4SSatish Balay ts->snes_its += nits; 1009371c9d4SSatish Balay ts->ksp_its += lits; 101c558785fSJames Wright pseudo->fnorm = -1; /* The current norm is no longer valid */ 102f4f49eeaSPierre Jolivet PetscCall(TSPostStage(ts, ts->ptime + ts->time_step, 0, &pseudo->update)); 103c558785fSJames Wright PetscCall(TSGetAdapt(ts, &adapt)); 104c558785fSJames Wright PetscCall(TSAdaptCheckStage(adapt, ts, ts->ptime + ts->time_step, pseudo->update, &stepok)); 1059371c9d4SSatish Balay if (!stepok) { 1069371c9d4SSatish Balay next_time_step = ts->time_step; 1079371c9d4SSatish Balay continue; 1089371c9d4SSatish Balay } 109c558785fSJames Wright PetscCall(VecCopy(pseudo->update, ts->vec_sol)); 110c558785fSJames Wright PetscCall(TSAdaptChoose(adapt, ts, ts->time_step, NULL, &next_time_step, &stepok)); 111cdbf8f93SLisandro Dalcin if (stepok) break; 112cdbf8f93SLisandro Dalcin } 113be5899b3SLisandro Dalcin if (reject >= ts->max_reject) { 114be5899b3SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 11563a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, reject)); 1163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1177bf11e45SBarry Smith } 118be5899b3SLisandro Dalcin 119be5899b3SLisandro Dalcin ts->ptime += ts->time_step; 120be5899b3SLisandro Dalcin ts->time_step = next_time_step; 121be5899b3SLisandro Dalcin 122*0fe2f5c2SJames Wright PetscCall(TSPseudoComputeFunction(ts, ts->vec_sol, NULL, NULL)); 123eb636060SBarry Smith 1243118ae5eSBarry Smith if (pseudo->fnorm < pseudo->fatol) { 1253118ae5eSBarry Smith ts->reason = TS_CONVERGED_PSEUDO_FATOL; 12663a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", converged since fnorm %g < fatol %g\n", ts->steps, (double)pseudo->fnorm, (double)pseudo->frtol)); 1273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1283118ae5eSBarry Smith } 1293118ae5eSBarry Smith if (pseudo->fnorm / pseudo->fnorm_initial < pseudo->frtol) { 1303118ae5eSBarry Smith ts->reason = TS_CONVERGED_PSEUDO_FRTOL; 13163a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", converged since fnorm %g / fnorm_initial %g < frtol %g\n", ts->steps, (double)pseudo->fnorm, (double)pseudo->fnorm_initial, (double)pseudo->fatol)); 1323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1333118ae5eSBarry Smith } 1343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1352d3f70b5SBarry Smith } 1362d3f70b5SBarry Smith 137d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_Pseudo(TS ts) 138d71ae5a4SJacob Faibussowitsch { 1397bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 1402d3f70b5SBarry Smith 1413a40ed3dSBarry Smith PetscFunctionBegin; 1429566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pseudo->update)); 1439566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pseudo->func)); 1449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pseudo->xdot)); 1453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1462d3f70b5SBarry Smith } 1472d3f70b5SBarry Smith 148d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_Pseudo(TS ts) 149d71ae5a4SJacob Faibussowitsch { 150277b19d0SLisandro Dalcin PetscFunctionBegin; 1519566063dSJacob Faibussowitsch PetscCall(TSReset_Pseudo(ts)); 1529566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 1539566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetVerifyTimeStep_C", NULL)); 1549566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStepIncrement_C", NULL)); 1559566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetMaxTimeStep_C", NULL)); 1569566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoIncrementDtFromInitialDt_C", NULL)); 1579566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStep_C", NULL)); 1583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 159277b19d0SLisandro Dalcin } 1602d3f70b5SBarry Smith 161d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoGetXdot(TS ts, Vec X, Vec *Xdot) 162d71ae5a4SJacob Faibussowitsch { 1636f2d6a7bSJed Brown TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 1642d3f70b5SBarry Smith 1653a40ed3dSBarry Smith PetscFunctionBegin; 1666f2d6a7bSJed Brown *Xdot = pseudo->xdot; 1673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1682d3f70b5SBarry Smith } 1692d3f70b5SBarry Smith 1706f2d6a7bSJed Brown /* 1716f2d6a7bSJed Brown The transient residual is 1726f2d6a7bSJed Brown 1736f2d6a7bSJed Brown F(U^{n+1},(U^{n+1}-U^n)/dt) = 0 1746f2d6a7bSJed Brown 1756f2d6a7bSJed Brown or for ODE, 1766f2d6a7bSJed Brown 1776f2d6a7bSJed Brown (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0 1786f2d6a7bSJed Brown 1796f2d6a7bSJed Brown This is the function that must be evaluated for transient simulation and for 1806f2d6a7bSJed Brown finite difference Jacobians. On the first Newton step, this algorithm uses 1816f2d6a7bSJed Brown a guess of U^{n+1} = U^n in which case the transient term vanishes and the 1826f2d6a7bSJed Brown residual is actually the steady state residual. Pseudotransient 1836f2d6a7bSJed Brown continuation as described in the literature is a linearly implicit 184eb636060SBarry Smith algorithm, it only takes this one full Newton step with the steady state 1856f2d6a7bSJed Brown residual, and then advances to the next time step. 186eb636060SBarry Smith 187eb636060SBarry Smith This routine avoids a repeated identical call to TSComputeRHSFunction() when that result 188eb636060SBarry Smith is already contained in pseudo->func due to the call to TSComputeIFunction() in TSStep_Pseudo() 189eb636060SBarry Smith 1906f2d6a7bSJed Brown */ 191d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes, Vec X, Vec Y, TS ts) 192d71ae5a4SJacob Faibussowitsch { 193eb636060SBarry Smith TSIFunctionFn *ifunction; 194eb636060SBarry Smith TSRHSFunctionFn *rhsfunction; 195eb636060SBarry Smith void *ctx; 196eb636060SBarry Smith DM dm; 197eb636060SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 19881775c15SStefano Zampini const PetscScalar mdt = 1.0 / ts->time_step; 199eb636060SBarry Smith PetscBool KSPSNES; 200eb636060SBarry Smith PetscObjectState Xstate; 2012d3f70b5SBarry Smith 2023a40ed3dSBarry Smith PetscFunctionBegin; 203eb636060SBarry Smith PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 204eb636060SBarry Smith PetscValidHeaderSpecific(X, VEC_CLASSID, 2); 205eb636060SBarry Smith PetscValidHeaderSpecific(Y, VEC_CLASSID, 3); 206eb636060SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID, 4); 207eb636060SBarry Smith 208eb636060SBarry Smith PetscCall(TSGetDM(ts, &dm)); 209eb636060SBarry Smith PetscCall(DMTSGetIFunction(dm, &ifunction, &ctx)); 210eb636060SBarry Smith PetscCall(DMTSGetRHSFunction(dm, &rhsfunction, NULL)); 211eb636060SBarry Smith PetscCheck(rhsfunction || ifunction, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Must call TSSetRHSFunction() and / or TSSetIFunction()"); 212eb636060SBarry Smith 213eb636060SBarry Smith PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESKSPONLY, &KSPSNES)); 214eb636060SBarry Smith PetscCall(PetscObjectStateGet((PetscObject)X, &Xstate)); 215eb636060SBarry Smith if (Xstate != pseudo->Xstate || ifunction || !KSPSNES) { 216c558785fSJames Wright // X == ts->vec_sol for first SNES iteration, so pseudo->xdot == 0 21781775c15SStefano Zampini PetscCall(VecAXPBYPCZ(pseudo->xdot, -mdt, mdt, 0, ts->vec_sol, X)); 21881775c15SStefano Zampini PetscCall(TSComputeIFunction(ts, ts->ptime + ts->time_step, X, pseudo->xdot, Y, PETSC_FALSE)); 219eb636060SBarry Smith } else { 220eb636060SBarry Smith /* reuse the TSComputeIFunction() result performed inside TSStep_Pseudo() */ 221eb636060SBarry Smith /* note that pseudo->func contains the negation of TSComputeRHSFunction() */ 22281775c15SStefano Zampini PetscCall(VecWAXPY(Y, 1, pseudo->func, pseudo->xdot)); 223eb636060SBarry Smith } 2243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2256f2d6a7bSJed Brown } 2262d3f70b5SBarry Smith 2276f2d6a7bSJed Brown /* 2286f2d6a7bSJed Brown This constructs the Jacobian needed for SNES. For DAE, this is 2296f2d6a7bSJed Brown 2306f2d6a7bSJed Brown dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot 2316f2d6a7bSJed Brown 2326f2d6a7bSJed Brown and for ODE: 2336f2d6a7bSJed Brown 2346f2d6a7bSJed Brown J = I/dt - J_{Frhs} where J_{Frhs} is the given Jacobian of Frhs. 2356f2d6a7bSJed Brown */ 236d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes, Vec X, Mat AA, Mat BB, TS ts) 237d71ae5a4SJacob Faibussowitsch { 2386f2d6a7bSJed Brown Vec Xdot; 2396f2d6a7bSJed Brown 2406f2d6a7bSJed Brown PetscFunctionBegin; 24181775c15SStefano Zampini /* Xdot has already been computed in SNESTSFormFunction_Pseudo (SNES guarantees this) */ 2429566063dSJacob Faibussowitsch PetscCall(TSPseudoGetXdot(ts, X, &Xdot)); 2439566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts, ts->ptime + ts->time_step, X, Xdot, 1. / ts->time_step, AA, BB, PETSC_FALSE)); 2443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2452d3f70b5SBarry Smith } 2462d3f70b5SBarry Smith 247d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_Pseudo(TS ts) 248d71ae5a4SJacob Faibussowitsch { 2497bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 2502d3f70b5SBarry Smith 2513a40ed3dSBarry Smith PetscFunctionBegin; 2529566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &pseudo->update)); 2539566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &pseudo->func)); 2549566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &pseudo->xdot)); 2553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2562d3f70b5SBarry Smith } 25781775c15SStefano Zampini 258d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoMonitorDefault(TS ts, PetscInt step, PetscReal ptime, Vec v, void *dummy) 259d71ae5a4SJacob Faibussowitsch { 2607bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 261ce94432eSBarry Smith PetscViewer viewer = (PetscViewer)dummy; 2622d3f70b5SBarry Smith 2633a40ed3dSBarry Smith PetscFunctionBegin; 264*0fe2f5c2SJames Wright PetscCall(TSPseudoComputeFunction(ts, ts->vec_sol, NULL, NULL)); 2659566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel)); 26663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "TS %" PetscInt_FMT " dt %g time %g fnorm %g\n", step, (double)ts->time_step, (double)ptime, (double)pseudo->fnorm)); 2679566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel)); 2683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2692d3f70b5SBarry Smith } 2702d3f70b5SBarry Smith 271ce78bad3SBarry Smith static PetscErrorCode TSSetFromOptions_Pseudo(TS ts, PetscOptionItems PetscOptionsObject) 272d71ae5a4SJacob Faibussowitsch { 2734bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 274ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE; 275649052a6SBarry Smith PetscViewer viewer; 2762d3f70b5SBarry Smith 2773a40ed3dSBarry Smith PetscFunctionBegin; 278d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Pseudo-timestepping options"); 2799566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_monitor_pseudo", "Monitor convergence", "", flg, &flg, NULL)); 2802d3f70b5SBarry Smith if (flg) { 2819566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ts), "stdout", &viewer)); 28249abdd8aSBarry Smith PetscCall(TSMonitorSet(ts, TSPseudoMonitorDefault, viewer, (PetscCtxDestroyFn *)PetscViewerDestroy)); 28328aa8177SBarry Smith } 284be5899b3SLisandro Dalcin flg = pseudo->increment_dt_from_initial_dt; 2859566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_pseudo_increment_dt_from_initial_dt", "Increase dt as a ratio from original dt", "TSPseudoIncrementDtFromInitialDt", flg, &flg, NULL)); 286be5899b3SLisandro Dalcin pseudo->increment_dt_from_initial_dt = flg; 2879566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_pseudo_increment", "Ratio to increase dt", "TSPseudoSetTimeStepIncrement", pseudo->dt_increment, &pseudo->dt_increment, NULL)); 2889566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_pseudo_max_dt", "Maximum value for dt", "TSPseudoSetMaxTimeStep", pseudo->dt_max, &pseudo->dt_max, NULL)); 2899566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_pseudo_fatol", "Tolerance for norm of function", "", pseudo->fatol, &pseudo->fatol, NULL)); 2909566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_pseudo_frtol", "Relative tolerance for norm of function", "", pseudo->frtol, &pseudo->frtol, NULL)); 291d0609cedSBarry Smith PetscOptionsHeadEnd(); 2923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2932d3f70b5SBarry Smith } 2942d3f70b5SBarry Smith 295d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_Pseudo(TS ts, PetscViewer viewer) 296d71ae5a4SJacob Faibussowitsch { 2973118ae5eSBarry Smith PetscBool isascii; 298d52bd9f3SBarry Smith 2993a40ed3dSBarry Smith PetscFunctionBegin; 3009566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 3013118ae5eSBarry Smith if (isascii) { 3023118ae5eSBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 3039566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " frtol - relative tolerance in function value %g\n", (double)pseudo->frtol)); 3049566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " fatol - absolute tolerance in function value %g\n", (double)pseudo->fatol)); 3059566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " dt_initial - initial timestep %g\n", (double)pseudo->dt_initial)); 3069566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " dt_increment - increase in timestep on successful step %g\n", (double)pseudo->dt_increment)); 3079566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " dt_max - maximum time %g\n", (double)pseudo->dt_max)); 3083118ae5eSBarry Smith } 3093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3102d3f70b5SBarry Smith } 3112d3f70b5SBarry Smith 312ac226902SBarry Smith /*@C 313c558785fSJames Wright TSPseudoTimeStepDefault - Default code to compute pseudo-timestepping. Use with `TSPseudoSetTimeStep()`. 314c558785fSJames Wright 315c558785fSJames Wright Collective, No Fortran Support 316c558785fSJames Wright 317c558785fSJames Wright Input Parameters: 318c558785fSJames Wright + ts - the timestep context 319c558785fSJames Wright - dtctx - unused timestep context 320c558785fSJames Wright 321c558785fSJames Wright Output Parameter: 322c558785fSJames Wright . newdt - the timestep to use for the next step 323c558785fSJames Wright 324c558785fSJames Wright Level: advanced 325c558785fSJames Wright 326c558785fSJames Wright .seealso: [](ch_ts), `TSPseudoSetTimeStep()`, `TSPseudoComputeTimeStep()`, `TSPSEUDO` 327c558785fSJames Wright @*/ 328c558785fSJames Wright PetscErrorCode TSPseudoTimeStepDefault(TS ts, PetscReal *newdt, void *dtctx) 329c558785fSJames Wright { 330c558785fSJames Wright TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 331*0fe2f5c2SJames Wright PetscReal inc = pseudo->dt_increment, fnorm; 332c558785fSJames Wright 333c558785fSJames Wright PetscFunctionBegin; 334*0fe2f5c2SJames Wright PetscCall(TSPseudoComputeFunction(ts, ts->vec_sol, NULL, &fnorm)); 335c558785fSJames Wright if (pseudo->fnorm_initial < 0) { 336c558785fSJames Wright /* first time through so compute initial function norm */ 337*0fe2f5c2SJames Wright pseudo->fnorm_initial = fnorm; 338*0fe2f5c2SJames Wright pseudo->fnorm_previous = fnorm; 339c558785fSJames Wright } 340*0fe2f5c2SJames Wright if (fnorm == 0.0) *newdt = 1.e12 * inc * ts->time_step; 341*0fe2f5c2SJames Wright else if (pseudo->increment_dt_from_initial_dt) *newdt = inc * pseudo->dt_initial * pseudo->fnorm_initial / fnorm; 342*0fe2f5c2SJames Wright else *newdt = inc * ts->time_step * pseudo->fnorm_previous / fnorm; 343c558785fSJames Wright if (pseudo->dt_max > 0) *newdt = PetscMin(*newdt, pseudo->dt_max); 344*0fe2f5c2SJames Wright pseudo->fnorm_previous = fnorm; 345c558785fSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 346c558785fSJames Wright } 347c558785fSJames Wright 348c558785fSJames Wright static PetscErrorCode TSAdaptChoose_TSPseudo(TSAdapt adapt, TS ts, PetscReal h, PetscInt *next_sc, PetscReal *next_h, PetscBool *accept, PetscReal *wlte, PetscReal *wltea, PetscReal *wlter) 349c558785fSJames Wright { 350c558785fSJames Wright TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 351c558785fSJames Wright 352c558785fSJames Wright PetscFunctionBegin; 353c558785fSJames Wright PetscCall(PetscLogEventBegin(TS_PseudoComputeTimeStep, ts, 0, 0, 0)); 354c558785fSJames Wright PetscCall((*pseudo->dt)(ts, next_h, pseudo->dtctx)); 355c558785fSJames Wright PetscCall(PetscLogEventEnd(TS_PseudoComputeTimeStep, ts, 0, 0, 0)); 356c558785fSJames Wright 357c558785fSJames Wright *next_sc = 0; 358c558785fSJames Wright *wlte = -1; /* Weighted local truncation error was not evaluated */ 359c558785fSJames Wright *wltea = -1; /* Weighted absolute local truncation error was not evaluated */ 360c558785fSJames Wright *wlter = -1; /* Weighted relative local truncation error was not evaluated */ 361c558785fSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 362c558785fSJames Wright } 363c558785fSJames Wright 364c558785fSJames Wright static PetscErrorCode TSAdaptCheckStage_TSPseudo(TSAdapt adapt, TS ts, PetscReal t, Vec Y, PetscBool *accept) 365c558785fSJames Wright { 366c558785fSJames Wright TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 367c558785fSJames Wright 368c558785fSJames Wright PetscFunctionBegin; 369c558785fSJames Wright if (pseudo->verify) { 370c558785fSJames Wright PetscReal dt; 371c558785fSJames Wright PetscCall(TSGetTimeStep(ts, &dt)); 372c558785fSJames Wright PetscCall((*pseudo->verify)(ts, Y, pseudo->verifyctx, &dt, accept)); 373c558785fSJames Wright PetscCall(TSSetTimeStep(ts, dt)); 374c558785fSJames Wright } 375c558785fSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 376c558785fSJames Wright } 377c558785fSJames Wright 378c558785fSJames Wright /*MC 379c558785fSJames Wright TSADAPTTSPSEUDO - TSPseudo adaptive controller for time stepping 380c558785fSJames Wright 381c558785fSJames Wright Level: developer 382c558785fSJames Wright 383c558785fSJames Wright Note: 384c558785fSJames Wright This is only meant to be used with `TSPSEUDO` time integrator. 385c558785fSJames Wright 386c558785fSJames Wright .seealso: [](ch_ts), `TS`, `TSAdapt`, `TSGetAdapt()`, `TSAdaptType`, `TSPSEUDO` 387c558785fSJames Wright M*/ 388c558785fSJames Wright static PetscErrorCode TSAdaptCreate_TSPseudo(TSAdapt adapt) 389c558785fSJames Wright { 390c558785fSJames Wright PetscFunctionBegin; 391c558785fSJames Wright adapt->ops->choose = TSAdaptChoose_TSPseudo; 392c558785fSJames Wright adapt->checkstage = TSAdaptCheckStage_TSPseudo; 393c558785fSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 394c558785fSJames Wright } 395c558785fSJames Wright 396c558785fSJames Wright /*@ 397c558785fSJames Wright TSPseudoComputeTimeStep - Computes the next timestep for a currently running 398c558785fSJames Wright pseudo-timestepping process. 399c558785fSJames Wright 400c558785fSJames Wright Collective 401c558785fSJames Wright 402c558785fSJames Wright Input Parameter: 403c558785fSJames Wright . ts - timestep context 404c558785fSJames Wright 405c558785fSJames Wright Output Parameter: 406c558785fSJames Wright . dt - newly computed timestep 407c558785fSJames Wright 408c558785fSJames Wright Level: developer 409c558785fSJames Wright 410c558785fSJames Wright Note: 411c558785fSJames Wright The routine to be called here to compute the timestep should be 412c558785fSJames Wright set by calling `TSPseudoSetTimeStep()`. 413c558785fSJames Wright 414c558785fSJames Wright .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoTimeStepDefault()`, `TSPseudoSetTimeStep()` 415c558785fSJames Wright @*/ 416c558785fSJames Wright PetscErrorCode TSPseudoComputeTimeStep(TS ts, PetscReal *dt) 417c558785fSJames Wright { 418c558785fSJames Wright TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 419c558785fSJames Wright 420c558785fSJames Wright PetscFunctionBegin; 421c558785fSJames Wright PetscCall(PetscLogEventBegin(TS_PseudoComputeTimeStep, ts, 0, 0, 0)); 422c558785fSJames Wright PetscCall((*pseudo->dt)(ts, dt, pseudo->dtctx)); 423c558785fSJames Wright PetscCall(PetscLogEventEnd(TS_PseudoComputeTimeStep, ts, 0, 0, 0)); 424c558785fSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 425c558785fSJames Wright } 426c558785fSJames Wright 427c558785fSJames Wright /*@C 428c558785fSJames Wright TSPseudoVerifyTimeStepDefault - Default code to verify the quality of the last timestep. 429c558785fSJames Wright 430c558785fSJames Wright Collective, No Fortran Support 431c558785fSJames Wright 432c558785fSJames Wright Input Parameters: 433c558785fSJames Wright + ts - the timestep context 434c558785fSJames Wright . dtctx - unused timestep context 435c558785fSJames Wright - update - latest solution vector 436c558785fSJames Wright 437c558785fSJames Wright Output Parameters: 438c558785fSJames Wright + newdt - the timestep to use for the next step 439c558785fSJames Wright - flag - flag indicating whether the last time step was acceptable 440c558785fSJames Wright 441c558785fSJames Wright Level: advanced 442c558785fSJames Wright 443c558785fSJames Wright Note: 444c558785fSJames Wright This routine always returns a flag of 1, indicating an acceptable 445c558785fSJames Wright timestep. 446c558785fSJames Wright 447c558785fSJames Wright .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetVerifyTimeStep()`, `TSPseudoVerifyTimeStep()` 448c558785fSJames Wright @*/ 449c558785fSJames Wright PetscErrorCode TSPseudoVerifyTimeStepDefault(TS ts, Vec update, void *dtctx, PetscReal *newdt, PetscBool *flag) 450c558785fSJames Wright { 451c558785fSJames Wright PetscFunctionBegin; 452c558785fSJames Wright // NOTE: This function was never used 453c558785fSJames Wright *flag = PETSC_TRUE; 454c558785fSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 455c558785fSJames Wright } 456c558785fSJames Wright 457c558785fSJames Wright /*@ 458c558785fSJames Wright TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable. 459c558785fSJames Wright 460c558785fSJames Wright Collective 461c558785fSJames Wright 462c558785fSJames Wright Input Parameters: 463c558785fSJames Wright + ts - timestep context 464c558785fSJames Wright - update - latest solution vector 465c558785fSJames Wright 466c558785fSJames Wright Output Parameters: 467c558785fSJames Wright + dt - newly computed timestep (if it had to shrink) 468c558785fSJames Wright - flag - indicates if current timestep was ok 469c558785fSJames Wright 470c558785fSJames Wright Level: advanced 471c558785fSJames Wright 472c558785fSJames Wright Notes: 473c558785fSJames Wright The routine to be called here to compute the timestep should be 474c558785fSJames Wright set by calling `TSPseudoSetVerifyTimeStep()`. 475c558785fSJames Wright 476c558785fSJames Wright .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetVerifyTimeStep()`, `TSPseudoVerifyTimeStepDefault()` 477c558785fSJames Wright @*/ 478c558785fSJames Wright PetscErrorCode TSPseudoVerifyTimeStep(TS ts, Vec update, PetscReal *dt, PetscBool *flag) 479c558785fSJames Wright { 480c558785fSJames Wright TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 481c558785fSJames Wright 482c558785fSJames Wright PetscFunctionBegin; 483c558785fSJames Wright // NOTE: This function is never used 484c558785fSJames Wright *flag = PETSC_TRUE; 485c558785fSJames Wright if (pseudo->verify) PetscCall((*pseudo->verify)(ts, update, pseudo->verifyctx, dt, flag)); 486c558785fSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 487c558785fSJames Wright } 488c558785fSJames Wright 489c558785fSJames Wright /*@C 49082bf6240SBarry Smith TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the 49182bf6240SBarry Smith last timestep. 49282bf6240SBarry Smith 493c3339decSBarry Smith Logically Collective 49415091d37SBarry Smith 49582bf6240SBarry Smith Input Parameters: 49615091d37SBarry Smith + ts - timestep context 49782bf6240SBarry Smith . dt - user-defined function to verify timestep 4982fe279fdSBarry Smith - ctx - [optional] user-defined context for private data for the timestep verification routine (may be `NULL`) 49982bf6240SBarry Smith 50020f4b53cSBarry Smith Calling sequence of `func`: 5012fe279fdSBarry Smith + ts - the time-step context 5022fe279fdSBarry Smith . update - latest solution vector 50314d0ab18SJacob Faibussowitsch . ctx - [optional] user-defined timestep context 50482bf6240SBarry Smith . newdt - the timestep to use for the next step 505a2b725a8SWilliam Gropp - flag - flag indicating whether the last time step was acceptable 50682bf6240SBarry Smith 507bcf0153eSBarry Smith Level: advanced 508bcf0153eSBarry Smith 509bcf0153eSBarry Smith Note: 510bcf0153eSBarry Smith The routine set here will be called by `TSPseudoVerifyTimeStep()` 51182bf6240SBarry Smith during the timestepping process. 51282bf6240SBarry Smith 5131cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoVerifyTimeStepDefault()`, `TSPseudoVerifyTimeStep()` 51482bf6240SBarry Smith @*/ 51514d0ab18SJacob Faibussowitsch PetscErrorCode TSPseudoSetVerifyTimeStep(TS ts, PetscErrorCode (*dt)(TS ts, Vec update, void *ctx, PetscReal *newdt, PetscBool *flag), void *ctx) 516d71ae5a4SJacob Faibussowitsch { 51782bf6240SBarry Smith PetscFunctionBegin; 5180700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 519cac4c232SBarry Smith PetscTryMethod(ts, "TSPseudoSetVerifyTimeStep_C", (TS, PetscErrorCode (*)(TS, Vec, void *, PetscReal *, PetscBool *), void *), (ts, dt, ctx)); 5203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 52182bf6240SBarry Smith } 52282bf6240SBarry Smith 52382bf6240SBarry Smith /*@ 52482bf6240SBarry Smith TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 5258d359177SBarry Smith dt when using the TSPseudoTimeStepDefault() routine. 52682bf6240SBarry Smith 527c3339decSBarry Smith Logically Collective 528fee21e36SBarry Smith 52915091d37SBarry Smith Input Parameters: 53015091d37SBarry Smith + ts - the timestep context 53115091d37SBarry Smith - inc - the scaling factor >= 1.0 53215091d37SBarry Smith 53382bf6240SBarry Smith Options Database Key: 53467b8a455SSatish Balay . -ts_pseudo_increment <increment> - set pseudo increment 53582bf6240SBarry Smith 53615091d37SBarry Smith Level: advanced 53715091d37SBarry Smith 5381cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()` 53982bf6240SBarry Smith @*/ 540d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoSetTimeStepIncrement(TS ts, PetscReal inc) 541d71ae5a4SJacob Faibussowitsch { 54282bf6240SBarry Smith PetscFunctionBegin; 5430700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 544c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts, inc, 2); 545cac4c232SBarry Smith PetscTryMethod(ts, "TSPseudoSetTimeStepIncrement_C", (TS, PetscReal), (ts, inc)); 5463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 54782bf6240SBarry Smith } 54882bf6240SBarry Smith 54986534af1SJed Brown /*@ 55086534af1SJed Brown TSPseudoSetMaxTimeStep - Sets the maximum time step 5518d359177SBarry Smith when using the TSPseudoTimeStepDefault() routine. 55286534af1SJed Brown 553c3339decSBarry Smith Logically Collective 55486534af1SJed Brown 55586534af1SJed Brown Input Parameters: 55686534af1SJed Brown + ts - the timestep context 55786534af1SJed Brown - maxdt - the maximum time step, use a non-positive value to deactivate 55886534af1SJed Brown 55986534af1SJed Brown Options Database Key: 56067b8a455SSatish Balay . -ts_pseudo_max_dt <increment> - set pseudo max dt 56186534af1SJed Brown 56286534af1SJed Brown Level: advanced 56386534af1SJed Brown 5641cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()` 56586534af1SJed Brown @*/ 566d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoSetMaxTimeStep(TS ts, PetscReal maxdt) 567d71ae5a4SJacob Faibussowitsch { 56886534af1SJed Brown PetscFunctionBegin; 56986534af1SJed Brown PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 57086534af1SJed Brown PetscValidLogicalCollectiveReal(ts, maxdt, 2); 571cac4c232SBarry Smith PetscTryMethod(ts, "TSPseudoSetMaxTimeStep_C", (TS, PetscReal), (ts, maxdt)); 5723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 57386534af1SJed Brown } 57486534af1SJed Brown 57582bf6240SBarry Smith /*@ 57682bf6240SBarry Smith TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep 577b44f4de4SBarry Smith is computed via the formula $ dt = initial\_dt*initial\_fnorm/current\_fnorm $ rather than the default update, $ dt = current\_dt*previous\_fnorm/current\_fnorm.$ 57882bf6240SBarry Smith 579c3339decSBarry Smith Logically Collective 58015091d37SBarry Smith 58182bf6240SBarry Smith Input Parameter: 58282bf6240SBarry Smith . ts - the timestep context 58382bf6240SBarry Smith 58482bf6240SBarry Smith Options Database Key: 585b44f4de4SBarry Smith . -ts_pseudo_increment_dt_from_initial_dt <true,false> - use the initial $dt$ to determine increment 58682bf6240SBarry Smith 58715091d37SBarry Smith Level: advanced 58815091d37SBarry Smith 5891cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()` 59082bf6240SBarry Smith @*/ 591d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS ts) 592d71ae5a4SJacob Faibussowitsch { 59382bf6240SBarry Smith PetscFunctionBegin; 5940700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 595cac4c232SBarry Smith PetscTryMethod(ts, "TSPseudoIncrementDtFromInitialDt_C", (TS), (ts)); 5963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 59782bf6240SBarry Smith } 59882bf6240SBarry Smith 599ac226902SBarry Smith /*@C 60082bf6240SBarry Smith TSPseudoSetTimeStep - Sets the user-defined routine to be 60182bf6240SBarry Smith called at each pseudo-timestep to update the timestep. 60282bf6240SBarry Smith 603c3339decSBarry Smith Logically Collective 60415091d37SBarry Smith 60582bf6240SBarry Smith Input Parameters: 60615091d37SBarry Smith + ts - timestep context 60782bf6240SBarry Smith . dt - function to compute timestep 6082fe279fdSBarry Smith - ctx - [optional] user-defined context for private data required by the function (may be `NULL`) 60982bf6240SBarry Smith 61020f4b53cSBarry Smith Calling sequence of `dt`: 61114d0ab18SJacob Faibussowitsch + ts - the `TS` context 61214d0ab18SJacob Faibussowitsch . newdt - the newly computed timestep 61314d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined context 61482bf6240SBarry Smith 615bcf0153eSBarry Smith Level: intermediate 61682bf6240SBarry Smith 617bcf0153eSBarry Smith Notes: 618bcf0153eSBarry Smith The routine set here will be called by `TSPseudoComputeTimeStep()` 619bcf0153eSBarry Smith during the timestepping process. 620bcf0153eSBarry Smith 621bcf0153eSBarry Smith If not set then `TSPseudoTimeStepDefault()` is automatically used 622bcf0153eSBarry Smith 6231cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoTimeStepDefault()`, `TSPseudoComputeTimeStep()` 62482bf6240SBarry Smith @*/ 62514d0ab18SJacob Faibussowitsch PetscErrorCode TSPseudoSetTimeStep(TS ts, PetscErrorCode (*dt)(TS ts, PetscReal *newdt, void *ctx), void *ctx) 626d71ae5a4SJacob Faibussowitsch { 62782bf6240SBarry Smith PetscFunctionBegin; 6280700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 629cac4c232SBarry Smith PetscTryMethod(ts, "TSPseudoSetTimeStep_C", (TS, PetscErrorCode (*)(TS, PetscReal *, void *), void *), (ts, dt, ctx)); 6303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 63182bf6240SBarry Smith } 63282bf6240SBarry Smith 633ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(TS, Vec, void *, PetscReal *, PetscBool *); /* force argument to next function to not be extern C*/ 634d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetVerifyTimeStep_Pseudo(TS ts, FCN1 dt, void *ctx) 635d71ae5a4SJacob Faibussowitsch { 636be5899b3SLisandro Dalcin TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 63782bf6240SBarry Smith 63882bf6240SBarry Smith PetscFunctionBegin; 63982bf6240SBarry Smith pseudo->verify = dt; 64082bf6240SBarry Smith pseudo->verifyctx = ctx; 6413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 64282bf6240SBarry Smith } 64382bf6240SBarry Smith 644d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetTimeStepIncrement_Pseudo(TS ts, PetscReal inc) 645d71ae5a4SJacob Faibussowitsch { 6464bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 64782bf6240SBarry Smith 64882bf6240SBarry Smith PetscFunctionBegin; 64982bf6240SBarry Smith pseudo->dt_increment = inc; 6503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 65182bf6240SBarry Smith } 65282bf6240SBarry Smith 653d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetMaxTimeStep_Pseudo(TS ts, PetscReal maxdt) 654d71ae5a4SJacob Faibussowitsch { 65586534af1SJed Brown TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 65686534af1SJed Brown 65786534af1SJed Brown PetscFunctionBegin; 65886534af1SJed Brown pseudo->dt_max = maxdt; 6593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 66086534af1SJed Brown } 66186534af1SJed Brown 662d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts) 663d71ae5a4SJacob Faibussowitsch { 6644bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 66582bf6240SBarry Smith 66682bf6240SBarry Smith PetscFunctionBegin; 6674bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_TRUE; 6683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 66982bf6240SBarry Smith } 67082bf6240SBarry Smith 6716849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS, PetscReal *, void *); /* force argument to next function to not be extern C*/ 672d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetTimeStep_Pseudo(TS ts, FCN2 dt, void *ctx) 673d71ae5a4SJacob Faibussowitsch { 6744bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 67582bf6240SBarry Smith 67682bf6240SBarry Smith PetscFunctionBegin; 67782bf6240SBarry Smith pseudo->dt = dt; 67882bf6240SBarry Smith pseudo->dtctx = ctx; 6793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 68082bf6240SBarry Smith } 68182bf6240SBarry Smith 68210e6a065SJed Brown /*MC 6831d27aa22SBarry Smith TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping {cite}`ckk02` {cite}`kk97` 68482bf6240SBarry Smith 68510e6a065SJed Brown This method solves equations of the form 68610e6a065SJed Brown 6871d27aa22SBarry Smith $$ 6881d27aa22SBarry Smith F(X,Xdot) = 0 6891d27aa22SBarry Smith $$ 69010e6a065SJed Brown 69110e6a065SJed Brown for steady state using the iteration 69210e6a065SJed Brown 6931d27aa22SBarry Smith $$ 69437fdd005SBarry Smith [G'] S = -F(X,0) 69537fdd005SBarry Smith X += S 6961d27aa22SBarry Smith $$ 69710e6a065SJed Brown 69810e6a065SJed Brown where 69910e6a065SJed Brown 7001d27aa22SBarry Smith $$ 7011d27aa22SBarry Smith G(Y) = F(Y,(Y-X)/dt) 7021d27aa22SBarry Smith $$ 70310e6a065SJed Brown 7046f2d6a7bSJed Brown This is linearly-implicit Euler with the residual always evaluated "at steady 7056f2d6a7bSJed Brown state". See note below. 70610e6a065SJed Brown 70793a54799SPierre Jolivet In addition to the modified solve, a dedicated adaptive timestepping scheme is implemented, mimicking the switched evolution relaxation in {cite}`ckk02`. 708fb59f417SJames Wright It determines the next timestep via 709fb59f417SJames Wright 710fb59f417SJames Wright $$ 711fb59f417SJames Wright dt_{n} = r dt_{n-1} \frac{\Vert F(X_{n},0) \Vert}{\Vert F(X_{n-1},0) \Vert} 712fb59f417SJames Wright $$ 713fb59f417SJames Wright 714fb59f417SJames Wright where $r$ is an additional growth factor (set by `-ts_pseudo_increment`). 715fb59f417SJames Wright An alternative formulation is also available that uses the initial timestep and function norm. 716fb59f417SJames Wright 717fb59f417SJames Wright $$ 718fb59f417SJames Wright dt_{n} = r dt_{0} \frac{\Vert F(X_{n},0) \Vert}{\Vert F(X_{0},0) \Vert} 719fb59f417SJames Wright $$ 720fb59f417SJames Wright 721fb59f417SJames Wright This is chosen by setting `-ts_pseudo_increment_dt_from_initial_dt`. 722fb59f417SJames Wright For either method, an upper limit on the timestep can be set by `-ts_pseudo_max_dt`. 723fb59f417SJames Wright 724bcf0153eSBarry Smith Options Database Keys: 72510e6a065SJed Brown + -ts_pseudo_increment <real> - ratio of increase dt 7263118ae5eSBarry Smith . -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt 727fb59f417SJames Wright . -ts_pseudo_max_dt - Maximum dt for adaptive timestepping algorithm 728fb59f417SJames Wright . -ts_pseudo_monitor - Monitor convergence of the function norm 729fb59f417SJames Wright . -ts_pseudo_fatol <atol> - stop iterating when the function norm is less than `atol` 730fb59f417SJames Wright - -ts_pseudo_frtol <rtol> - stop iterating when the function norm divided by the initial function norm is less than `rtol` 73110e6a065SJed Brown 73210e6a065SJed Brown Level: beginner 73310e6a065SJed Brown 73410e6a065SJed Brown Notes: 7356f2d6a7bSJed Brown The residual computed by this method includes the transient term (Xdot is computed instead of 7366f2d6a7bSJed Brown always being zero), but since the prediction from the last step is always the solution from the 7376f2d6a7bSJed Brown last step, on the first Newton iteration we have 7386f2d6a7bSJed Brown 7391d27aa22SBarry Smith $$ 7401d27aa22SBarry Smith Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0 7411d27aa22SBarry Smith $$ 7426f2d6a7bSJed Brown 743eb636060SBarry Smith The Jacobian $ dF/dX + shift*dF/dXdot $ contains a non-zero $ dF/dXdot$ term. In the $ X' = F(X) $ case it 744eb636060SBarry Smith is $ dF/dX + shift*1/dt $. Hence still contains the $ 1/dt $ term so the Jacobian is not simply the 745eb636060SBarry Smith Jacobian of $ F $ and thus this pseudo-transient continuation is not just Newton on $F(x)=0$. 746eb636060SBarry Smith 7476f2d6a7bSJed Brown Therefore, the linear system solved by the first Newton iteration is equivalent to the one 7486f2d6a7bSJed Brown described above and in the papers. If the user chooses to perform multiple Newton iterations, the 7496f2d6a7bSJed Brown algorithm is no longer the one described in the referenced papers. 750fb59f417SJames Wright By default, the `SNESType` is set to `SNESKSPONLY` to match the algorithm from the referenced papers. 751fb59f417SJames Wright Setting the `SNESType` via `-snes_type` will override this default setting. 75210e6a065SJed Brown 7531cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()` 75410e6a065SJed Brown M*/ 755d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_Pseudo(TS ts) 756d71ae5a4SJacob Faibussowitsch { 7577bf11e45SBarry Smith TS_Pseudo *pseudo; 758193ac0bcSJed Brown SNES snes; 75919fd82e9SBarry Smith SNESType stype; 7602d3f70b5SBarry Smith 7613a40ed3dSBarry Smith PetscFunctionBegin; 762277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Pseudo; 763000e7ae3SMatthew Knepley ts->ops->destroy = TSDestroy_Pseudo; 764000e7ae3SMatthew Knepley ts->ops->view = TSView_Pseudo; 765000e7ae3SMatthew Knepley ts->ops->setup = TSSetUp_Pseudo; 766000e7ae3SMatthew Knepley ts->ops->step = TSStep_Pseudo; 767000e7ae3SMatthew Knepley ts->ops->setfromoptions = TSSetFromOptions_Pseudo; 7680f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Pseudo; 7690f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Pseudo; 770c558785fSJames Wright ts->default_adapt_type = TSADAPTTSPSEUDO; 771825ab935SBarry Smith ts->usessnes = PETSC_TRUE; 7727bf11e45SBarry Smith 773c558785fSJames Wright PetscCall(TSAdaptRegister(TSADAPTTSPSEUDO, TSAdaptCreate_TSPseudo)); 774c558785fSJames Wright 7759566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 7769566063dSJacob Faibussowitsch PetscCall(SNESGetType(snes, &stype)); 7779566063dSJacob Faibussowitsch if (!stype) PetscCall(SNESSetType(snes, SNESKSPONLY)); 7782d3f70b5SBarry Smith 7794dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&pseudo)); 7807bf11e45SBarry Smith ts->data = (void *)pseudo; 7812d3f70b5SBarry Smith 782be5899b3SLisandro Dalcin pseudo->dt = TSPseudoTimeStepDefault; 783be5899b3SLisandro Dalcin pseudo->dtctx = NULL; 78428aa8177SBarry Smith pseudo->dt_increment = 1.1; 7854bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_FALSE; 786193ac0bcSJed Brown pseudo->fnorm = -1; 787be5899b3SLisandro Dalcin pseudo->fnorm_initial = -1; 788be5899b3SLisandro Dalcin pseudo->fnorm_previous = -1; 7893118ae5eSBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 7903118ae5eSBarry Smith pseudo->fatol = 1.e-25; 7913118ae5eSBarry Smith pseudo->frtol = 1.e-5; 7923118ae5eSBarry Smith #else 7933118ae5eSBarry Smith pseudo->fatol = 1.e-50; 7943118ae5eSBarry Smith pseudo->frtol = 1.e-12; 7953118ae5eSBarry Smith #endif 7969566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetVerifyTimeStep_C", TSPseudoSetVerifyTimeStep_Pseudo)); 7979566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStepIncrement_C", TSPseudoSetTimeStepIncrement_Pseudo)); 7989566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetMaxTimeStep_C", TSPseudoSetMaxTimeStep_Pseudo)); 7999566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoIncrementDtFromInitialDt_C", TSPseudoIncrementDtFromInitialDt_Pseudo)); 8009566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStep_C", TSPseudoSetTimeStep_Pseudo)); 8013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8022d3f70b5SBarry Smith } 803