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 func; /* work vector where F(t[i],u[i]) is stored */ 106f2d6a7bSJed Brown Vec xdot; /* work vector for time derivative of state */ 112d3f70b5SBarry Smith 122d3f70b5SBarry Smith /* information used for Pseudo-timestepping */ 132d3f70b5SBarry Smith 146849ba73SBarry Smith PetscErrorCode (*dt)(TS, PetscReal *, void *); /* compute next timestep, and related context */ 152d3f70b5SBarry Smith void *dtctx; 16ace3abfcSBarry Smith PetscErrorCode (*verify)(TS, Vec, void *, PetscReal *, PetscBool *); /* verify previous timestep and related context */ 177bf11e45SBarry Smith void *verifyctx; 182d3f70b5SBarry Smith 19cdbf8f93SLisandro Dalcin PetscReal fnorm_initial, fnorm; /* original and current norm of F(u) */ 2087828ca2SBarry Smith PetscReal fnorm_previous; 2128aa8177SBarry Smith 22cdbf8f93SLisandro Dalcin PetscReal dt_initial; /* initial time-step */ 2387828ca2SBarry Smith PetscReal dt_increment; /* scaling that dt is incremented each time-step */ 2486534af1SJed Brown PetscReal dt_max; /* maximum time step */ 25ace3abfcSBarry Smith PetscBool increment_dt_from_initial_dt; 263118ae5eSBarry Smith PetscReal fatol, frtol; 27eb636060SBarry Smith 28eb636060SBarry Smith PetscObjectState Xstate; /* state of input vector to TSComputeIFunction() with 0 Xdot */ 29940ebdd4SJames Wright 30940ebdd4SJames Wright TSStepStatus status; 317bf11e45SBarry Smith } TS_Pseudo; 322d3f70b5SBarry Smith 330fe2f5c2SJames Wright /*@C 340fe2f5c2SJames Wright TSPseudoComputeFunction - Compute nonlinear residual for pseudo-timestepping 350fe2f5c2SJames Wright 360fe2f5c2SJames Wright This computes the residual for $\dot U = 0$, i.e. $F(U, 0)$ for the IFunction. 370fe2f5c2SJames Wright 380fe2f5c2SJames Wright Collective 390fe2f5c2SJames Wright 400fe2f5c2SJames Wright Input Parameters: 410fe2f5c2SJames Wright + ts - the timestep context 420fe2f5c2SJames Wright - solution - the solution vector 430fe2f5c2SJames Wright 440fe2f5c2SJames Wright Output Parameter: 450fe2f5c2SJames Wright + residual - the nonlinear residual 460fe2f5c2SJames Wright - fnorm - the norm of the nonlinear residual 470fe2f5c2SJames Wright 480fe2f5c2SJames Wright Level: advanced 490fe2f5c2SJames Wright 500fe2f5c2SJames Wright Note: 510fe2f5c2SJames 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. 520fe2f5c2SJames Wright 530fe2f5c2SJames 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. 540fe2f5c2SJames Wright 55f7782ccaSJames Wright To correctly get the residual reuse behavior, `solution` must be the same `Vec` that returned by `TSGetSolution()`. 56f7782ccaSJames Wright 570fe2f5c2SJames Wright .seealso: [](ch_ts), `TSPSEUDO` 580fe2f5c2SJames Wright @*/ 590fe2f5c2SJames Wright PetscErrorCode TSPseudoComputeFunction(TS ts, Vec solution, Vec *residual, PetscReal *fnorm) 600fe2f5c2SJames Wright { 610fe2f5c2SJames Wright TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 620fe2f5c2SJames Wright PetscObjectState Xstate; 630fe2f5c2SJames Wright 640fe2f5c2SJames Wright PetscFunctionBegin; 650fe2f5c2SJames Wright PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 660fe2f5c2SJames Wright PetscValidHeaderSpecific(solution, VEC_CLASSID, 2); 670fe2f5c2SJames Wright if (residual) PetscValidHeaderSpecific(*residual, VEC_CLASSID, 3); 680fe2f5c2SJames Wright if (fnorm) PetscAssertPointer(fnorm, 4); 690fe2f5c2SJames Wright 700fe2f5c2SJames Wright PetscCall(PetscObjectStateGet((PetscObject)solution, &Xstate)); 710fe2f5c2SJames Wright if (Xstate != pseudo->Xstate || pseudo->fnorm < 0) { 720fe2f5c2SJames Wright PetscCall(VecZeroEntries(pseudo->xdot)); 730fe2f5c2SJames Wright PetscCall(TSComputeIFunction(ts, ts->ptime, solution, pseudo->xdot, pseudo->func, PETSC_FALSE)); 740fe2f5c2SJames Wright pseudo->Xstate = Xstate; 750fe2f5c2SJames Wright PetscCall(VecNorm(pseudo->func, NORM_2, &pseudo->fnorm)); 760fe2f5c2SJames Wright } 770fe2f5c2SJames Wright if (residual) *residual = pseudo->func; 780fe2f5c2SJames Wright if (fnorm) *fnorm = pseudo->fnorm; 790fe2f5c2SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 800fe2f5c2SJames Wright } 810fe2f5c2SJames Wright 82d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_Pseudo(TS ts) 83d71ae5a4SJacob Faibussowitsch { 84277b19d0SLisandro Dalcin TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 85940ebdd4SJames Wright PetscInt nits, lits, rejections = 0; 86940ebdd4SJames Wright PetscBool accept; 87*c2cf51afSJames Wright PetscReal next_time_step = ts->time_step, fnorm; 88c558785fSJames Wright TSAdapt adapt; 892d3f70b5SBarry Smith 903a40ed3dSBarry Smith PetscFunctionBegin; 91f7782ccaSJames Wright if (ts->steps == 0) pseudo->dt_initial = ts->time_step; 92940ebdd4SJames Wright 93940ebdd4SJames Wright pseudo->status = TS_STEP_INCOMPLETE; 94940ebdd4SJames Wright while (!ts->reason && pseudo->status != TS_STEP_COMPLETE) { 959566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, ts->ptime + ts->time_step)); 96f7782ccaSJames Wright PetscCall(SNESSolve(ts->snes, NULL, ts->vec_sol)); 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 */ 102f7782ccaSJames Wright PetscCall(TSPostStage(ts, ts->ptime + ts->time_step, 0, &ts->vec_sol)); 103c558785fSJames Wright PetscCall(TSGetAdapt(ts, &adapt)); 104f7782ccaSJames Wright PetscCall(TSAdaptCheckStage(adapt, ts, ts->ptime + ts->time_step, ts->vec_sol, &accept)); 105940ebdd4SJames Wright if (!accept) goto reject_step; 106be5899b3SLisandro Dalcin 107940ebdd4SJames Wright pseudo->status = TS_STEP_PENDING; 108940ebdd4SJames Wright PetscCall(TSAdaptChoose(adapt, ts, ts->time_step, NULL, &next_time_step, &accept)); 109940ebdd4SJames Wright pseudo->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 110940ebdd4SJames Wright if (!accept) { 111940ebdd4SJames Wright ts->time_step = next_time_step; 112940ebdd4SJames Wright goto reject_step; 113940ebdd4SJames Wright } 114be5899b3SLisandro Dalcin ts->ptime += ts->time_step; 115be5899b3SLisandro Dalcin ts->time_step = next_time_step; 116940ebdd4SJames Wright break; 117be5899b3SLisandro Dalcin 118940ebdd4SJames Wright reject_step: 119940ebdd4SJames Wright ts->reject++; 120940ebdd4SJames Wright accept = PETSC_FALSE; 121f7782ccaSJames Wright PetscCall(VecCopy(ts->vec_sol0, ts->vec_sol)); 122940ebdd4SJames Wright if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 123940ebdd4SJames Wright ts->reason = TS_DIVERGED_STEP_REJECTED; 124940ebdd4SJames Wright PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections)); 125940ebdd4SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 126940ebdd4SJames Wright } 127940ebdd4SJames Wright } 128940ebdd4SJames Wright 129940ebdd4SJames Wright // Check solution convergence 130*c2cf51afSJames Wright PetscCall(TSPseudoComputeFunction(ts, ts->vec_sol, NULL, &fnorm)); 131eb636060SBarry Smith 132*c2cf51afSJames Wright if (fnorm < pseudo->fatol) { 1333118ae5eSBarry Smith ts->reason = TS_CONVERGED_PSEUDO_FATOL; 13463a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", converged since fnorm %g < fatol %g\n", ts->steps, (double)pseudo->fnorm, (double)pseudo->frtol)); 1353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1363118ae5eSBarry Smith } 137*c2cf51afSJames Wright if (fnorm / pseudo->fnorm_initial < pseudo->frtol) { 1383118ae5eSBarry Smith ts->reason = TS_CONVERGED_PSEUDO_FRTOL; 13963a3b9bcSJacob 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)); 1403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1413118ae5eSBarry Smith } 1423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1432d3f70b5SBarry Smith } 1442d3f70b5SBarry Smith 145d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_Pseudo(TS ts) 146d71ae5a4SJacob Faibussowitsch { 1477bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 1482d3f70b5SBarry Smith 1493a40ed3dSBarry Smith PetscFunctionBegin; 1509566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pseudo->func)); 1519566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pseudo->xdot)); 1523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1532d3f70b5SBarry Smith } 1542d3f70b5SBarry Smith 155d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_Pseudo(TS ts) 156d71ae5a4SJacob Faibussowitsch { 157277b19d0SLisandro Dalcin PetscFunctionBegin; 1589566063dSJacob Faibussowitsch PetscCall(TSReset_Pseudo(ts)); 1599566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 1609566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetVerifyTimeStep_C", NULL)); 1619566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStepIncrement_C", NULL)); 1629566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetMaxTimeStep_C", NULL)); 1639566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoIncrementDtFromInitialDt_C", NULL)); 1649566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStep_C", NULL)); 1653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 166277b19d0SLisandro Dalcin } 1672d3f70b5SBarry Smith 168d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoGetXdot(TS ts, Vec X, Vec *Xdot) 169d71ae5a4SJacob Faibussowitsch { 1706f2d6a7bSJed Brown TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 1712d3f70b5SBarry Smith 1723a40ed3dSBarry Smith PetscFunctionBegin; 1736f2d6a7bSJed Brown *Xdot = pseudo->xdot; 1743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1752d3f70b5SBarry Smith } 1762d3f70b5SBarry Smith 1776f2d6a7bSJed Brown /* 1786f2d6a7bSJed Brown The transient residual is 1796f2d6a7bSJed Brown 1806f2d6a7bSJed Brown F(U^{n+1},(U^{n+1}-U^n)/dt) = 0 1816f2d6a7bSJed Brown 1826f2d6a7bSJed Brown or for ODE, 1836f2d6a7bSJed Brown 1846f2d6a7bSJed Brown (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0 1856f2d6a7bSJed Brown 1866f2d6a7bSJed Brown This is the function that must be evaluated for transient simulation and for 1876f2d6a7bSJed Brown finite difference Jacobians. On the first Newton step, this algorithm uses 1886f2d6a7bSJed Brown a guess of U^{n+1} = U^n in which case the transient term vanishes and the 1896f2d6a7bSJed Brown residual is actually the steady state residual. Pseudotransient 1906f2d6a7bSJed Brown continuation as described in the literature is a linearly implicit 191eb636060SBarry Smith algorithm, it only takes this one full Newton step with the steady state 1926f2d6a7bSJed Brown residual, and then advances to the next time step. 193eb636060SBarry Smith 194eb636060SBarry Smith This routine avoids a repeated identical call to TSComputeRHSFunction() when that result 195eb636060SBarry Smith is already contained in pseudo->func due to the call to TSComputeIFunction() in TSStep_Pseudo() 196eb636060SBarry Smith 1976f2d6a7bSJed Brown */ 198d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes, Vec X, Vec Y, TS ts) 199d71ae5a4SJacob Faibussowitsch { 200eb636060SBarry Smith TSIFunctionFn *ifunction; 201eb636060SBarry Smith TSRHSFunctionFn *rhsfunction; 202eb636060SBarry Smith void *ctx; 203eb636060SBarry Smith DM dm; 204eb636060SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 20581775c15SStefano Zampini const PetscScalar mdt = 1.0 / ts->time_step; 206eb636060SBarry Smith PetscObjectState Xstate; 207e344a936SJames Wright PetscInt snes_it; 2082d3f70b5SBarry Smith 2093a40ed3dSBarry Smith PetscFunctionBegin; 210eb636060SBarry Smith PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 211eb636060SBarry Smith PetscValidHeaderSpecific(X, VEC_CLASSID, 2); 212eb636060SBarry Smith PetscValidHeaderSpecific(Y, VEC_CLASSID, 3); 213eb636060SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID, 4); 214eb636060SBarry Smith 215eb636060SBarry Smith PetscCall(TSGetDM(ts, &dm)); 216eb636060SBarry Smith PetscCall(DMTSGetIFunction(dm, &ifunction, &ctx)); 217eb636060SBarry Smith PetscCall(DMTSGetRHSFunction(dm, &rhsfunction, NULL)); 218eb636060SBarry Smith PetscCheck(rhsfunction || ifunction, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Must call TSSetRHSFunction() and / or TSSetIFunction()"); 219eb636060SBarry Smith 220eb636060SBarry Smith PetscCall(PetscObjectStateGet((PetscObject)X, &Xstate)); 221e344a936SJames Wright PetscCall(SNESGetIterationNumber(snes, &snes_it)); 222e344a936SJames Wright if (Xstate == pseudo->Xstate && snes_it == 1) { 223e344a936SJames Wright /* reuse the TSComputeIFunction() result performed inside TSStep_Pseudo() */ 224e344a936SJames Wright if (ifunction) PetscCall(VecCopy(pseudo->func, Y)); 225e344a936SJames Wright /* note that pseudo->func contains the negation of TSComputeRHSFunction() */ 226e344a936SJames Wright else PetscCall(VecWAXPY(Y, 1, pseudo->func, pseudo->xdot)); 227e344a936SJames Wright } else { 228f7782ccaSJames Wright PetscCall(VecAXPBYPCZ(pseudo->xdot, -mdt, mdt, 0, ts->vec_sol0, X)); 22981775c15SStefano Zampini PetscCall(TSComputeIFunction(ts, ts->ptime + ts->time_step, X, pseudo->xdot, Y, PETSC_FALSE)); 230eb636060SBarry Smith } 2313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2326f2d6a7bSJed Brown } 2332d3f70b5SBarry Smith 2346f2d6a7bSJed Brown /* 2356f2d6a7bSJed Brown This constructs the Jacobian needed for SNES. For DAE, this is 2366f2d6a7bSJed Brown 2376f2d6a7bSJed Brown dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot 2386f2d6a7bSJed Brown 2396f2d6a7bSJed Brown and for ODE: 2406f2d6a7bSJed Brown 2416f2d6a7bSJed Brown J = I/dt - J_{Frhs} where J_{Frhs} is the given Jacobian of Frhs. 2426f2d6a7bSJed Brown */ 243d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes, Vec X, Mat AA, Mat BB, TS ts) 244d71ae5a4SJacob Faibussowitsch { 2456f2d6a7bSJed Brown Vec Xdot; 2466f2d6a7bSJed Brown 2476f2d6a7bSJed Brown PetscFunctionBegin; 24881775c15SStefano Zampini /* Xdot has already been computed in SNESTSFormFunction_Pseudo (SNES guarantees this) */ 2499566063dSJacob Faibussowitsch PetscCall(TSPseudoGetXdot(ts, X, &Xdot)); 2509566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts, ts->ptime + ts->time_step, X, Xdot, 1. / ts->time_step, AA, BB, PETSC_FALSE)); 2513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2522d3f70b5SBarry Smith } 2532d3f70b5SBarry Smith 254d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_Pseudo(TS ts) 255d71ae5a4SJacob Faibussowitsch { 2567bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 2572d3f70b5SBarry Smith 2583a40ed3dSBarry Smith PetscFunctionBegin; 2599566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &pseudo->func)); 2609566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &pseudo->xdot)); 2613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2622d3f70b5SBarry Smith } 26381775c15SStefano Zampini 264d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoMonitorDefault(TS ts, PetscInt step, PetscReal ptime, Vec v, void *dummy) 265d71ae5a4SJacob Faibussowitsch { 2667bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 267ce94432eSBarry Smith PetscViewer viewer = (PetscViewer)dummy; 2682d3f70b5SBarry Smith 2693a40ed3dSBarry Smith PetscFunctionBegin; 2700fe2f5c2SJames Wright PetscCall(TSPseudoComputeFunction(ts, ts->vec_sol, NULL, NULL)); 2719566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel)); 27263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "TS %" PetscInt_FMT " dt %g time %g fnorm %g\n", step, (double)ts->time_step, (double)ptime, (double)pseudo->fnorm)); 2739566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel)); 2743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2752d3f70b5SBarry Smith } 2762d3f70b5SBarry Smith 277ce78bad3SBarry Smith static PetscErrorCode TSSetFromOptions_Pseudo(TS ts, PetscOptionItems PetscOptionsObject) 278d71ae5a4SJacob Faibussowitsch { 2794bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 280ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE; 281649052a6SBarry Smith PetscViewer viewer; 2822d3f70b5SBarry Smith 2833a40ed3dSBarry Smith PetscFunctionBegin; 284d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Pseudo-timestepping options"); 2859566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_monitor_pseudo", "Monitor convergence", "", flg, &flg, NULL)); 2862d3f70b5SBarry Smith if (flg) { 2879566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ts), "stdout", &viewer)); 28849abdd8aSBarry Smith PetscCall(TSMonitorSet(ts, TSPseudoMonitorDefault, viewer, (PetscCtxDestroyFn *)PetscViewerDestroy)); 28928aa8177SBarry Smith } 290be5899b3SLisandro Dalcin flg = pseudo->increment_dt_from_initial_dt; 2919566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_pseudo_increment_dt_from_initial_dt", "Increase dt as a ratio from original dt", "TSPseudoIncrementDtFromInitialDt", flg, &flg, NULL)); 292be5899b3SLisandro Dalcin pseudo->increment_dt_from_initial_dt = flg; 2939566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_pseudo_increment", "Ratio to increase dt", "TSPseudoSetTimeStepIncrement", pseudo->dt_increment, &pseudo->dt_increment, NULL)); 2949566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_pseudo_max_dt", "Maximum value for dt", "TSPseudoSetMaxTimeStep", pseudo->dt_max, &pseudo->dt_max, NULL)); 2959566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_pseudo_fatol", "Tolerance for norm of function", "", pseudo->fatol, &pseudo->fatol, NULL)); 2969566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_pseudo_frtol", "Relative tolerance for norm of function", "", pseudo->frtol, &pseudo->frtol, NULL)); 297d0609cedSBarry Smith PetscOptionsHeadEnd(); 2983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2992d3f70b5SBarry Smith } 3002d3f70b5SBarry Smith 301d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_Pseudo(TS ts, PetscViewer viewer) 302d71ae5a4SJacob Faibussowitsch { 3033118ae5eSBarry Smith PetscBool isascii; 304d52bd9f3SBarry Smith 3053a40ed3dSBarry Smith PetscFunctionBegin; 3069566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 3073118ae5eSBarry Smith if (isascii) { 3083118ae5eSBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 3099566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " frtol - relative tolerance in function value %g\n", (double)pseudo->frtol)); 3109566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " fatol - absolute tolerance in function value %g\n", (double)pseudo->fatol)); 3119566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " dt_initial - initial timestep %g\n", (double)pseudo->dt_initial)); 3129566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " dt_increment - increase in timestep on successful step %g\n", (double)pseudo->dt_increment)); 3139566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " dt_max - maximum time %g\n", (double)pseudo->dt_max)); 3143118ae5eSBarry Smith } 3153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3162d3f70b5SBarry Smith } 3172d3f70b5SBarry Smith 318ac226902SBarry Smith /*@C 319c558785fSJames Wright TSPseudoTimeStepDefault - Default code to compute pseudo-timestepping. Use with `TSPseudoSetTimeStep()`. 320c558785fSJames Wright 321c558785fSJames Wright Collective, No Fortran Support 322c558785fSJames Wright 323c558785fSJames Wright Input Parameters: 324c558785fSJames Wright + ts - the timestep context 325c558785fSJames Wright - dtctx - unused timestep context 326c558785fSJames Wright 327c558785fSJames Wright Output Parameter: 328c558785fSJames Wright . newdt - the timestep to use for the next step 329c558785fSJames Wright 330c558785fSJames Wright Level: advanced 331c558785fSJames Wright 332c558785fSJames Wright .seealso: [](ch_ts), `TSPseudoSetTimeStep()`, `TSPseudoComputeTimeStep()`, `TSPSEUDO` 333c558785fSJames Wright @*/ 334c558785fSJames Wright PetscErrorCode TSPseudoTimeStepDefault(TS ts, PetscReal *newdt, void *dtctx) 335c558785fSJames Wright { 336c558785fSJames Wright TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 3370fe2f5c2SJames Wright PetscReal inc = pseudo->dt_increment, fnorm; 338c558785fSJames Wright 339c558785fSJames Wright PetscFunctionBegin; 3400fe2f5c2SJames Wright PetscCall(TSPseudoComputeFunction(ts, ts->vec_sol, NULL, &fnorm)); 341c558785fSJames Wright if (pseudo->fnorm_initial < 0) { 342c558785fSJames Wright /* first time through so compute initial function norm */ 3430fe2f5c2SJames Wright pseudo->fnorm_initial = fnorm; 3440fe2f5c2SJames Wright pseudo->fnorm_previous = fnorm; 345c558785fSJames Wright } 3460fe2f5c2SJames Wright if (fnorm == 0.0) *newdt = 1.e12 * inc * ts->time_step; 3470fe2f5c2SJames Wright else if (pseudo->increment_dt_from_initial_dt) *newdt = inc * pseudo->dt_initial * pseudo->fnorm_initial / fnorm; 3480fe2f5c2SJames Wright else *newdt = inc * ts->time_step * pseudo->fnorm_previous / fnorm; 349c558785fSJames Wright if (pseudo->dt_max > 0) *newdt = PetscMin(*newdt, pseudo->dt_max); 3500fe2f5c2SJames Wright pseudo->fnorm_previous = fnorm; 351c558785fSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 352c558785fSJames Wright } 353c558785fSJames Wright 354c558785fSJames 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) 355c558785fSJames Wright { 356c558785fSJames Wright TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 357c558785fSJames Wright 358c558785fSJames Wright PetscFunctionBegin; 359c558785fSJames Wright PetscCall(PetscLogEventBegin(TS_PseudoComputeTimeStep, ts, 0, 0, 0)); 360f6ceb9cfSJames Wright PetscCallBack("TSPSEUDO callback time step", (*pseudo->dt)(ts, next_h, pseudo->dtctx)); 361c558785fSJames Wright PetscCall(PetscLogEventEnd(TS_PseudoComputeTimeStep, ts, 0, 0, 0)); 362c558785fSJames Wright 363c558785fSJames Wright *next_sc = 0; 364c558785fSJames Wright *wlte = -1; /* Weighted local truncation error was not evaluated */ 365c558785fSJames Wright *wltea = -1; /* Weighted absolute local truncation error was not evaluated */ 366c558785fSJames Wright *wlter = -1; /* Weighted relative local truncation error was not evaluated */ 367c558785fSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 368c558785fSJames Wright } 369c558785fSJames Wright 370c558785fSJames Wright static PetscErrorCode TSAdaptCheckStage_TSPseudo(TSAdapt adapt, TS ts, PetscReal t, Vec Y, PetscBool *accept) 371c558785fSJames Wright { 372c558785fSJames Wright TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 373c558785fSJames Wright 374c558785fSJames Wright PetscFunctionBegin; 375c558785fSJames Wright if (pseudo->verify) { 376c558785fSJames Wright PetscReal dt; 377c558785fSJames Wright PetscCall(TSGetTimeStep(ts, &dt)); 378f6ceb9cfSJames Wright PetscCallBack("TSPSEUDO callback verify time step", (*pseudo->verify)(ts, Y, pseudo->verifyctx, &dt, accept)); 379c558785fSJames Wright PetscCall(TSSetTimeStep(ts, dt)); 380c558785fSJames Wright } 381c558785fSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 382c558785fSJames Wright } 383c558785fSJames Wright 384c558785fSJames Wright /*MC 385c558785fSJames Wright TSADAPTTSPSEUDO - TSPseudo adaptive controller for time stepping 386c558785fSJames Wright 387c558785fSJames Wright Level: developer 388c558785fSJames Wright 389c558785fSJames Wright Note: 390c558785fSJames Wright This is only meant to be used with `TSPSEUDO` time integrator. 391c558785fSJames Wright 392c558785fSJames Wright .seealso: [](ch_ts), `TS`, `TSAdapt`, `TSGetAdapt()`, `TSAdaptType`, `TSPSEUDO` 393c558785fSJames Wright M*/ 394c558785fSJames Wright static PetscErrorCode TSAdaptCreate_TSPseudo(TSAdapt adapt) 395c558785fSJames Wright { 396c558785fSJames Wright PetscFunctionBegin; 397c558785fSJames Wright adapt->ops->choose = TSAdaptChoose_TSPseudo; 398c558785fSJames Wright adapt->checkstage = TSAdaptCheckStage_TSPseudo; 399c558785fSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 400c558785fSJames Wright } 401c558785fSJames Wright 402c558785fSJames Wright /*@C 40382bf6240SBarry Smith TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the 40482bf6240SBarry Smith last timestep. 40582bf6240SBarry Smith 406c3339decSBarry Smith Logically Collective 40715091d37SBarry Smith 40882bf6240SBarry Smith Input Parameters: 40915091d37SBarry Smith + ts - timestep context 41082bf6240SBarry Smith . dt - user-defined function to verify timestep 4112fe279fdSBarry Smith - ctx - [optional] user-defined context for private data for the timestep verification routine (may be `NULL`) 41282bf6240SBarry Smith 41320f4b53cSBarry Smith Calling sequence of `func`: 4142fe279fdSBarry Smith + ts - the time-step context 4152fe279fdSBarry Smith . update - latest solution vector 41614d0ab18SJacob Faibussowitsch . ctx - [optional] user-defined timestep context 41782bf6240SBarry Smith . newdt - the timestep to use for the next step 418a2b725a8SWilliam Gropp - flag - flag indicating whether the last time step was acceptable 41982bf6240SBarry Smith 420bcf0153eSBarry Smith Level: advanced 421bcf0153eSBarry Smith 422bcf0153eSBarry Smith Note: 423bcf0153eSBarry Smith The routine set here will be called by `TSPseudoVerifyTimeStep()` 42482bf6240SBarry Smith during the timestepping process. 42582bf6240SBarry Smith 4261cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoVerifyTimeStepDefault()`, `TSPseudoVerifyTimeStep()` 42782bf6240SBarry Smith @*/ 42814d0ab18SJacob Faibussowitsch PetscErrorCode TSPseudoSetVerifyTimeStep(TS ts, PetscErrorCode (*dt)(TS ts, Vec update, void *ctx, PetscReal *newdt, PetscBool *flag), void *ctx) 429d71ae5a4SJacob Faibussowitsch { 43082bf6240SBarry Smith PetscFunctionBegin; 4310700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 432cac4c232SBarry Smith PetscTryMethod(ts, "TSPseudoSetVerifyTimeStep_C", (TS, PetscErrorCode (*)(TS, Vec, void *, PetscReal *, PetscBool *), void *), (ts, dt, ctx)); 4333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 43482bf6240SBarry Smith } 43582bf6240SBarry Smith 43682bf6240SBarry Smith /*@ 43782bf6240SBarry Smith TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 4388d359177SBarry Smith dt when using the TSPseudoTimeStepDefault() routine. 43982bf6240SBarry Smith 440c3339decSBarry Smith Logically Collective 441fee21e36SBarry Smith 44215091d37SBarry Smith Input Parameters: 44315091d37SBarry Smith + ts - the timestep context 44415091d37SBarry Smith - inc - the scaling factor >= 1.0 44515091d37SBarry Smith 44682bf6240SBarry Smith Options Database Key: 44767b8a455SSatish Balay . -ts_pseudo_increment <increment> - set pseudo increment 44882bf6240SBarry Smith 44915091d37SBarry Smith Level: advanced 45015091d37SBarry Smith 4511cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()` 45282bf6240SBarry Smith @*/ 453d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoSetTimeStepIncrement(TS ts, PetscReal inc) 454d71ae5a4SJacob Faibussowitsch { 45582bf6240SBarry Smith PetscFunctionBegin; 4560700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 457c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts, inc, 2); 458cac4c232SBarry Smith PetscTryMethod(ts, "TSPseudoSetTimeStepIncrement_C", (TS, PetscReal), (ts, inc)); 4593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 46082bf6240SBarry Smith } 46182bf6240SBarry Smith 46286534af1SJed Brown /*@ 46386534af1SJed Brown TSPseudoSetMaxTimeStep - Sets the maximum time step 4648d359177SBarry Smith when using the TSPseudoTimeStepDefault() routine. 46586534af1SJed Brown 466c3339decSBarry Smith Logically Collective 46786534af1SJed Brown 46886534af1SJed Brown Input Parameters: 46986534af1SJed Brown + ts - the timestep context 47086534af1SJed Brown - maxdt - the maximum time step, use a non-positive value to deactivate 47186534af1SJed Brown 47286534af1SJed Brown Options Database Key: 47367b8a455SSatish Balay . -ts_pseudo_max_dt <increment> - set pseudo max dt 47486534af1SJed Brown 47586534af1SJed Brown Level: advanced 47686534af1SJed Brown 4771cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()` 47886534af1SJed Brown @*/ 479d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoSetMaxTimeStep(TS ts, PetscReal maxdt) 480d71ae5a4SJacob Faibussowitsch { 48186534af1SJed Brown PetscFunctionBegin; 48286534af1SJed Brown PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 48386534af1SJed Brown PetscValidLogicalCollectiveReal(ts, maxdt, 2); 484cac4c232SBarry Smith PetscTryMethod(ts, "TSPseudoSetMaxTimeStep_C", (TS, PetscReal), (ts, maxdt)); 4853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 48686534af1SJed Brown } 48786534af1SJed Brown 48882bf6240SBarry Smith /*@ 48982bf6240SBarry Smith TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep 490b44f4de4SBarry 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.$ 49182bf6240SBarry Smith 492c3339decSBarry Smith Logically Collective 49315091d37SBarry Smith 49482bf6240SBarry Smith Input Parameter: 49582bf6240SBarry Smith . ts - the timestep context 49682bf6240SBarry Smith 49782bf6240SBarry Smith Options Database Key: 498b44f4de4SBarry Smith . -ts_pseudo_increment_dt_from_initial_dt <true,false> - use the initial $dt$ to determine increment 49982bf6240SBarry Smith 50015091d37SBarry Smith Level: advanced 50115091d37SBarry Smith 5021cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()` 50382bf6240SBarry Smith @*/ 504d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS ts) 505d71ae5a4SJacob Faibussowitsch { 50682bf6240SBarry Smith PetscFunctionBegin; 5070700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 508cac4c232SBarry Smith PetscTryMethod(ts, "TSPseudoIncrementDtFromInitialDt_C", (TS), (ts)); 5093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 51082bf6240SBarry Smith } 51182bf6240SBarry Smith 512ac226902SBarry Smith /*@C 51382bf6240SBarry Smith TSPseudoSetTimeStep - Sets the user-defined routine to be 51482bf6240SBarry Smith called at each pseudo-timestep to update the timestep. 51582bf6240SBarry Smith 516c3339decSBarry Smith Logically Collective 51715091d37SBarry Smith 51882bf6240SBarry Smith Input Parameters: 51915091d37SBarry Smith + ts - timestep context 52082bf6240SBarry Smith . dt - function to compute timestep 5212fe279fdSBarry Smith - ctx - [optional] user-defined context for private data required by the function (may be `NULL`) 52282bf6240SBarry Smith 52320f4b53cSBarry Smith Calling sequence of `dt`: 52414d0ab18SJacob Faibussowitsch + ts - the `TS` context 52514d0ab18SJacob Faibussowitsch . newdt - the newly computed timestep 52614d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined context 52782bf6240SBarry Smith 528bcf0153eSBarry Smith Level: intermediate 52982bf6240SBarry Smith 530bcf0153eSBarry Smith Notes: 531bcf0153eSBarry Smith The routine set here will be called by `TSPseudoComputeTimeStep()` 532bcf0153eSBarry Smith during the timestepping process. 533bcf0153eSBarry Smith 534bcf0153eSBarry Smith If not set then `TSPseudoTimeStepDefault()` is automatically used 535bcf0153eSBarry Smith 5361cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoTimeStepDefault()`, `TSPseudoComputeTimeStep()` 53782bf6240SBarry Smith @*/ 53814d0ab18SJacob Faibussowitsch PetscErrorCode TSPseudoSetTimeStep(TS ts, PetscErrorCode (*dt)(TS ts, PetscReal *newdt, void *ctx), void *ctx) 539d71ae5a4SJacob Faibussowitsch { 54082bf6240SBarry Smith PetscFunctionBegin; 5410700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 542cac4c232SBarry Smith PetscTryMethod(ts, "TSPseudoSetTimeStep_C", (TS, PetscErrorCode (*)(TS, PetscReal *, void *), void *), (ts, dt, ctx)); 5433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 54482bf6240SBarry Smith } 54582bf6240SBarry Smith 546ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(TS, Vec, void *, PetscReal *, PetscBool *); /* force argument to next function to not be extern C*/ 547d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetVerifyTimeStep_Pseudo(TS ts, FCN1 dt, void *ctx) 548d71ae5a4SJacob Faibussowitsch { 549be5899b3SLisandro Dalcin TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 55082bf6240SBarry Smith 55182bf6240SBarry Smith PetscFunctionBegin; 55282bf6240SBarry Smith pseudo->verify = dt; 55382bf6240SBarry Smith pseudo->verifyctx = ctx; 5543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 55582bf6240SBarry Smith } 55682bf6240SBarry Smith 557d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetTimeStepIncrement_Pseudo(TS ts, PetscReal inc) 558d71ae5a4SJacob Faibussowitsch { 5594bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 56082bf6240SBarry Smith 56182bf6240SBarry Smith PetscFunctionBegin; 56282bf6240SBarry Smith pseudo->dt_increment = inc; 5633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 56482bf6240SBarry Smith } 56582bf6240SBarry Smith 566d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetMaxTimeStep_Pseudo(TS ts, PetscReal maxdt) 567d71ae5a4SJacob Faibussowitsch { 56886534af1SJed Brown TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 56986534af1SJed Brown 57086534af1SJed Brown PetscFunctionBegin; 57186534af1SJed Brown pseudo->dt_max = maxdt; 5723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 57386534af1SJed Brown } 57486534af1SJed Brown 575d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts) 576d71ae5a4SJacob Faibussowitsch { 5774bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 57882bf6240SBarry Smith 57982bf6240SBarry Smith PetscFunctionBegin; 5804bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_TRUE; 5813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 58282bf6240SBarry Smith } 58382bf6240SBarry Smith 5846849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS, PetscReal *, void *); /* force argument to next function to not be extern C*/ 585d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetTimeStep_Pseudo(TS ts, FCN2 dt, void *ctx) 586d71ae5a4SJacob Faibussowitsch { 5874bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 58882bf6240SBarry Smith 58982bf6240SBarry Smith PetscFunctionBegin; 59082bf6240SBarry Smith pseudo->dt = dt; 59182bf6240SBarry Smith pseudo->dtctx = ctx; 5923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 59382bf6240SBarry Smith } 59482bf6240SBarry Smith 59510e6a065SJed Brown /*MC 5961d27aa22SBarry Smith TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping {cite}`ckk02` {cite}`kk97` 59782bf6240SBarry Smith 59810e6a065SJed Brown This method solves equations of the form 59910e6a065SJed Brown 6001d27aa22SBarry Smith $$ 6011d27aa22SBarry Smith F(X,Xdot) = 0 6021d27aa22SBarry Smith $$ 60310e6a065SJed Brown 60410e6a065SJed Brown for steady state using the iteration 60510e6a065SJed Brown 6061d27aa22SBarry Smith $$ 60737fdd005SBarry Smith [G'] S = -F(X,0) 60837fdd005SBarry Smith X += S 6091d27aa22SBarry Smith $$ 61010e6a065SJed Brown 61110e6a065SJed Brown where 61210e6a065SJed Brown 6131d27aa22SBarry Smith $$ 6141d27aa22SBarry Smith G(Y) = F(Y,(Y-X)/dt) 6151d27aa22SBarry Smith $$ 61610e6a065SJed Brown 6176f2d6a7bSJed Brown This is linearly-implicit Euler with the residual always evaluated "at steady 6186f2d6a7bSJed Brown state". See note below. 61910e6a065SJed Brown 62093a54799SPierre Jolivet In addition to the modified solve, a dedicated adaptive timestepping scheme is implemented, mimicking the switched evolution relaxation in {cite}`ckk02`. 621fb59f417SJames Wright It determines the next timestep via 622fb59f417SJames Wright 623fb59f417SJames Wright $$ 624fb59f417SJames Wright dt_{n} = r dt_{n-1} \frac{\Vert F(X_{n},0) \Vert}{\Vert F(X_{n-1},0) \Vert} 625fb59f417SJames Wright $$ 626fb59f417SJames Wright 627fb59f417SJames Wright where $r$ is an additional growth factor (set by `-ts_pseudo_increment`). 628fb59f417SJames Wright An alternative formulation is also available that uses the initial timestep and function norm. 629fb59f417SJames Wright 630fb59f417SJames Wright $$ 631fb59f417SJames Wright dt_{n} = r dt_{0} \frac{\Vert F(X_{n},0) \Vert}{\Vert F(X_{0},0) \Vert} 632fb59f417SJames Wright $$ 633fb59f417SJames Wright 634fb59f417SJames Wright This is chosen by setting `-ts_pseudo_increment_dt_from_initial_dt`. 635fb59f417SJames Wright For either method, an upper limit on the timestep can be set by `-ts_pseudo_max_dt`. 636fb59f417SJames Wright 637bcf0153eSBarry Smith Options Database Keys: 63810e6a065SJed Brown + -ts_pseudo_increment <real> - ratio of increase dt 6393118ae5eSBarry Smith . -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt 640fb59f417SJames Wright . -ts_pseudo_max_dt - Maximum dt for adaptive timestepping algorithm 641fb59f417SJames Wright . -ts_pseudo_monitor - Monitor convergence of the function norm 642fb59f417SJames Wright . -ts_pseudo_fatol <atol> - stop iterating when the function norm is less than `atol` 643fb59f417SJames Wright - -ts_pseudo_frtol <rtol> - stop iterating when the function norm divided by the initial function norm is less than `rtol` 64410e6a065SJed Brown 64510e6a065SJed Brown Level: beginner 64610e6a065SJed Brown 64710e6a065SJed Brown Notes: 6486f2d6a7bSJed Brown The residual computed by this method includes the transient term (Xdot is computed instead of 6496f2d6a7bSJed Brown always being zero), but since the prediction from the last step is always the solution from the 6506f2d6a7bSJed Brown last step, on the first Newton iteration we have 6516f2d6a7bSJed Brown 6521d27aa22SBarry Smith $$ 6531d27aa22SBarry Smith Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0 6541d27aa22SBarry Smith $$ 6556f2d6a7bSJed Brown 656eb636060SBarry Smith The Jacobian $ dF/dX + shift*dF/dXdot $ contains a non-zero $ dF/dXdot$ term. In the $ X' = F(X) $ case it 657eb636060SBarry Smith is $ dF/dX + shift*1/dt $. Hence still contains the $ 1/dt $ term so the Jacobian is not simply the 658eb636060SBarry Smith Jacobian of $ F $ and thus this pseudo-transient continuation is not just Newton on $F(x)=0$. 659eb636060SBarry Smith 6606f2d6a7bSJed Brown Therefore, the linear system solved by the first Newton iteration is equivalent to the one 6616f2d6a7bSJed Brown described above and in the papers. If the user chooses to perform multiple Newton iterations, the 6626f2d6a7bSJed Brown algorithm is no longer the one described in the referenced papers. 663fb59f417SJames Wright By default, the `SNESType` is set to `SNESKSPONLY` to match the algorithm from the referenced papers. 664fb59f417SJames Wright Setting the `SNESType` via `-snes_type` will override this default setting. 66510e6a065SJed Brown 6661cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()` 66710e6a065SJed Brown M*/ 668d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_Pseudo(TS ts) 669d71ae5a4SJacob Faibussowitsch { 6707bf11e45SBarry Smith TS_Pseudo *pseudo; 671193ac0bcSJed Brown SNES snes; 67219fd82e9SBarry Smith SNESType stype; 6732d3f70b5SBarry Smith 6743a40ed3dSBarry Smith PetscFunctionBegin; 675277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Pseudo; 676000e7ae3SMatthew Knepley ts->ops->destroy = TSDestroy_Pseudo; 677000e7ae3SMatthew Knepley ts->ops->view = TSView_Pseudo; 678000e7ae3SMatthew Knepley ts->ops->setup = TSSetUp_Pseudo; 679000e7ae3SMatthew Knepley ts->ops->step = TSStep_Pseudo; 680000e7ae3SMatthew Knepley ts->ops->setfromoptions = TSSetFromOptions_Pseudo; 6810f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Pseudo; 6820f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Pseudo; 683c558785fSJames Wright ts->default_adapt_type = TSADAPTTSPSEUDO; 684825ab935SBarry Smith ts->usessnes = PETSC_TRUE; 6857bf11e45SBarry Smith 686c558785fSJames Wright PetscCall(TSAdaptRegister(TSADAPTTSPSEUDO, TSAdaptCreate_TSPseudo)); 687c558785fSJames Wright 6889566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 6899566063dSJacob Faibussowitsch PetscCall(SNESGetType(snes, &stype)); 6909566063dSJacob Faibussowitsch if (!stype) PetscCall(SNESSetType(snes, SNESKSPONLY)); 6912d3f70b5SBarry Smith 6924dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&pseudo)); 6937bf11e45SBarry Smith ts->data = (void *)pseudo; 6942d3f70b5SBarry Smith 695be5899b3SLisandro Dalcin pseudo->dt = TSPseudoTimeStepDefault; 696be5899b3SLisandro Dalcin pseudo->dtctx = NULL; 69728aa8177SBarry Smith pseudo->dt_increment = 1.1; 6984bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_FALSE; 699193ac0bcSJed Brown pseudo->fnorm = -1; 700be5899b3SLisandro Dalcin pseudo->fnorm_initial = -1; 701be5899b3SLisandro Dalcin pseudo->fnorm_previous = -1; 7023118ae5eSBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 7033118ae5eSBarry Smith pseudo->fatol = 1.e-25; 7043118ae5eSBarry Smith pseudo->frtol = 1.e-5; 7053118ae5eSBarry Smith #else 7063118ae5eSBarry Smith pseudo->fatol = 1.e-50; 7073118ae5eSBarry Smith pseudo->frtol = 1.e-12; 7083118ae5eSBarry Smith #endif 7099566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetVerifyTimeStep_C", TSPseudoSetVerifyTimeStep_Pseudo)); 7109566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStepIncrement_C", TSPseudoSetTimeStepIncrement_Pseudo)); 7119566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetMaxTimeStep_C", TSPseudoSetMaxTimeStep_Pseudo)); 7129566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoIncrementDtFromInitialDt_C", TSPseudoIncrementDtFromInitialDt_Pseudo)); 7139566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStep_C", TSPseudoSetTimeStep_Pseudo)); 7143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7152d3f70b5SBarry Smith } 716