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 62d3f70b5SBarry Smith typedef struct { 72d3f70b5SBarry Smith Vec update; /* work vector where new solution is formed */ 82d3f70b5SBarry Smith Vec func; /* work vector where F(t[i],u[i]) is stored */ 96f2d6a7bSJed Brown Vec xdot; /* work vector for time derivative of state */ 102d3f70b5SBarry Smith 112d3f70b5SBarry Smith /* information used for Pseudo-timestepping */ 122d3f70b5SBarry Smith 136849ba73SBarry Smith PetscErrorCode (*dt)(TS, PetscReal *, void *); /* compute next timestep, and related context */ 142d3f70b5SBarry Smith void *dtctx; 15ace3abfcSBarry Smith PetscErrorCode (*verify)(TS, Vec, void *, PetscReal *, PetscBool *); /* verify previous timestep and related context */ 167bf11e45SBarry Smith void *verifyctx; 172d3f70b5SBarry Smith 18cdbf8f93SLisandro Dalcin PetscReal fnorm_initial, fnorm; /* original and current norm of F(u) */ 1987828ca2SBarry Smith PetscReal fnorm_previous; 2028aa8177SBarry Smith 21cdbf8f93SLisandro Dalcin PetscReal dt_initial; /* initial time-step */ 2287828ca2SBarry Smith PetscReal dt_increment; /* scaling that dt is incremented each time-step */ 2386534af1SJed Brown PetscReal dt_max; /* maximum time step */ 24ace3abfcSBarry Smith PetscBool increment_dt_from_initial_dt; 253118ae5eSBarry Smith PetscReal fatol, frtol; 267bf11e45SBarry Smith } TS_Pseudo; 272d3f70b5SBarry Smith 282d3f70b5SBarry Smith /* ------------------------------------------------------------------------------*/ 292d3f70b5SBarry Smith 308d359177SBarry Smith /*@C 317bf11e45SBarry Smith TSPseudoComputeTimeStep - Computes the next timestep for a currently running 32564e8f4eSLois Curfman McInnes pseudo-timestepping process. 332d3f70b5SBarry Smith 34c3339decSBarry Smith Collective 3515091d37SBarry Smith 367bf11e45SBarry Smith Input Parameter: 377bf11e45SBarry Smith . ts - timestep context 387bf11e45SBarry Smith 397bf11e45SBarry Smith Output Parameter: 40fb4a63b6SLois Curfman McInnes . dt - newly computed timestep 41fb4a63b6SLois Curfman McInnes 428d359177SBarry Smith Level: developer 43564e8f4eSLois Curfman McInnes 44bcf0153eSBarry Smith Note: 45564e8f4eSLois Curfman McInnes The routine to be called here to compute the timestep should be 46bcf0153eSBarry Smith set by calling `TSPseudoSetTimeStep()`. 47564e8f4eSLois Curfman McInnes 481cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoTimeStepDefault()`, `TSPseudoSetTimeStep()` 497bf11e45SBarry Smith @*/ 50d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoComputeTimeStep(TS ts, PetscReal *dt) 51d71ae5a4SJacob Faibussowitsch { 527bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 537bf11e45SBarry Smith 543a40ed3dSBarry Smith PetscFunctionBegin; 559566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TS_PseudoComputeTimeStep, ts, 0, 0, 0)); 569566063dSJacob Faibussowitsch PetscCall((*pseudo->dt)(ts, dt, pseudo->dtctx)); 579566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TS_PseudoComputeTimeStep, ts, 0, 0, 0)); 583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 597bf11e45SBarry Smith } 607bf11e45SBarry Smith 617bf11e45SBarry Smith /* ------------------------------------------------------------------------------*/ 627bf11e45SBarry Smith /*@C 638d359177SBarry Smith TSPseudoVerifyTimeStepDefault - Default code to verify the quality of the last timestep. 647bf11e45SBarry Smith 65c3339decSBarry Smith Collective 6615091d37SBarry Smith 677bf11e45SBarry Smith Input Parameters: 6815091d37SBarry Smith + ts - the timestep context 697bf11e45SBarry Smith . dtctx - unused timestep context 7015091d37SBarry Smith - update - latest solution vector 717bf11e45SBarry Smith 72564e8f4eSLois Curfman McInnes Output Parameters: 7315091d37SBarry Smith + newdt - the timestep to use for the next step 7415091d37SBarry Smith - flag - flag indicating whether the last time step was acceptable 757bf11e45SBarry Smith 7615091d37SBarry Smith Level: advanced 77fee21e36SBarry Smith 78564e8f4eSLois Curfman McInnes Note: 79564e8f4eSLois Curfman McInnes This routine always returns a flag of 1, indicating an acceptable 80564e8f4eSLois Curfman McInnes timestep. 81564e8f4eSLois Curfman McInnes 821cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetVerifyTimeStep()`, `TSPseudoVerifyTimeStep()` 837bf11e45SBarry Smith @*/ 84d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoVerifyTimeStepDefault(TS ts, Vec update, void *dtctx, PetscReal *newdt, PetscBool *flag) 85d71ae5a4SJacob Faibussowitsch { 863a40ed3dSBarry Smith PetscFunctionBegin; 87a7cc72afSBarry Smith *flag = PETSC_TRUE; 883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 897bf11e45SBarry Smith } 907bf11e45SBarry Smith 917bf11e45SBarry Smith /*@ 92564e8f4eSLois Curfman McInnes TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable. 937bf11e45SBarry Smith 94c3339decSBarry Smith Collective 9515091d37SBarry Smith 96fb4a63b6SLois Curfman McInnes Input Parameters: 9715091d37SBarry Smith + ts - timestep context 9815091d37SBarry Smith - update - latest solution vector 997bf11e45SBarry Smith 100fb4a63b6SLois Curfman McInnes Output Parameters: 10115091d37SBarry Smith + dt - newly computed timestep (if it had to shrink) 10215091d37SBarry Smith - flag - indicates if current timestep was ok 1037bf11e45SBarry Smith 10415091d37SBarry Smith Level: advanced 105fee21e36SBarry Smith 106564e8f4eSLois Curfman McInnes Notes: 107564e8f4eSLois Curfman McInnes The routine to be called here to compute the timestep should be 108bcf0153eSBarry Smith set by calling `TSPseudoSetVerifyTimeStep()`. 109564e8f4eSLois Curfman McInnes 1101cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetVerifyTimeStep()`, `TSPseudoVerifyTimeStepDefault()` 1117bf11e45SBarry Smith @*/ 112d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoVerifyTimeStep(TS ts, Vec update, PetscReal *dt, PetscBool *flag) 113d71ae5a4SJacob Faibussowitsch { 1147bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 1157bf11e45SBarry Smith 1163a40ed3dSBarry Smith PetscFunctionBegin; 117cb9d8021SPierre Barbier de Reuille *flag = PETSC_TRUE; 1181baa6e33SBarry Smith if (pseudo->verify) PetscCall((*pseudo->verify)(ts, update, pseudo->verifyctx, dt, flag)); 1193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1207bf11e45SBarry Smith } 1217bf11e45SBarry Smith 1227bf11e45SBarry Smith /* --------------------------------------------------------------------------------*/ 1237bf11e45SBarry Smith 124d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_Pseudo(TS ts) 125d71ae5a4SJacob Faibussowitsch { 126277b19d0SLisandro Dalcin TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 127be5899b3SLisandro Dalcin PetscInt nits, lits, reject; 128cdbf8f93SLisandro Dalcin PetscBool stepok; 129be5899b3SLisandro Dalcin PetscReal next_time_step = ts->time_step; 1302d3f70b5SBarry Smith 1313a40ed3dSBarry Smith PetscFunctionBegin; 132bbd56ea5SKarl Rupp if (ts->steps == 0) pseudo->dt_initial = ts->time_step; 1339566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, pseudo->update)); 1349566063dSJacob Faibussowitsch PetscCall(TSPseudoComputeTimeStep(ts, &next_time_step)); 135cdbf8f93SLisandro Dalcin for (reject = 0; reject < ts->max_reject; reject++, ts->reject++) { 136cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 1379566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, ts->ptime + ts->time_step)); 1389566063dSJacob Faibussowitsch PetscCall(SNESSolve(ts->snes, NULL, pseudo->update)); 1399566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(ts->snes, &nits)); 1409566063dSJacob Faibussowitsch PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits)); 1419371c9d4SSatish Balay ts->snes_its += nits; 1429371c9d4SSatish Balay ts->ksp_its += lits; 1439566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, ts->ptime + ts->time_step, 0, &(pseudo->update))); 1449566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(ts->adapt, ts, ts->ptime + ts->time_step, pseudo->update, &stepok)); 1459371c9d4SSatish Balay if (!stepok) { 1469371c9d4SSatish Balay next_time_step = ts->time_step; 1479371c9d4SSatish Balay continue; 1489371c9d4SSatish Balay } 149193ac0bcSJed Brown pseudo->fnorm = -1; /* The current norm is no longer valid, monitor must recompute it. */ 1509566063dSJacob Faibussowitsch PetscCall(TSPseudoVerifyTimeStep(ts, pseudo->update, &next_time_step, &stepok)); 151cdbf8f93SLisandro Dalcin if (stepok) break; 152cdbf8f93SLisandro Dalcin } 153be5899b3SLisandro Dalcin if (reject >= ts->max_reject) { 154be5899b3SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 15563a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, reject)); 1563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1577bf11e45SBarry Smith } 158be5899b3SLisandro Dalcin 1599566063dSJacob Faibussowitsch PetscCall(VecCopy(pseudo->update, ts->vec_sol)); 160be5899b3SLisandro Dalcin ts->ptime += ts->time_step; 161be5899b3SLisandro Dalcin ts->time_step = next_time_step; 162be5899b3SLisandro Dalcin 1633118ae5eSBarry Smith if (pseudo->fnorm < 0) { 1649566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(pseudo->xdot)); 1659566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts, ts->ptime, ts->vec_sol, pseudo->xdot, pseudo->func, PETSC_FALSE)); 1669566063dSJacob Faibussowitsch PetscCall(VecNorm(pseudo->func, NORM_2, &pseudo->fnorm)); 1673118ae5eSBarry Smith } 1683118ae5eSBarry Smith if (pseudo->fnorm < pseudo->fatol) { 1693118ae5eSBarry Smith ts->reason = TS_CONVERGED_PSEUDO_FATOL; 17063a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", converged since fnorm %g < fatol %g\n", ts->steps, (double)pseudo->fnorm, (double)pseudo->frtol)); 1713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1723118ae5eSBarry Smith } 1733118ae5eSBarry Smith if (pseudo->fnorm / pseudo->fnorm_initial < pseudo->frtol) { 1743118ae5eSBarry Smith ts->reason = TS_CONVERGED_PSEUDO_FRTOL; 17563a3b9bcSJacob 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)); 1763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1773118ae5eSBarry Smith } 1783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1792d3f70b5SBarry Smith } 1802d3f70b5SBarry Smith 1812d3f70b5SBarry Smith /*------------------------------------------------------------*/ 182d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_Pseudo(TS ts) 183d71ae5a4SJacob Faibussowitsch { 1847bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 1852d3f70b5SBarry Smith 1863a40ed3dSBarry Smith PetscFunctionBegin; 1879566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pseudo->update)); 1889566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pseudo->func)); 1899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pseudo->xdot)); 1903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1912d3f70b5SBarry Smith } 1922d3f70b5SBarry Smith 193d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_Pseudo(TS ts) 194d71ae5a4SJacob Faibussowitsch { 195277b19d0SLisandro Dalcin PetscFunctionBegin; 1969566063dSJacob Faibussowitsch PetscCall(TSReset_Pseudo(ts)); 1979566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 1989566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetVerifyTimeStep_C", NULL)); 1999566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStepIncrement_C", NULL)); 2009566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetMaxTimeStep_C", NULL)); 2019566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoIncrementDtFromInitialDt_C", NULL)); 2029566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStep_C", NULL)); 2033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 204277b19d0SLisandro Dalcin } 2052d3f70b5SBarry Smith 2062d3f70b5SBarry Smith /*------------------------------------------------------------*/ 2072d3f70b5SBarry Smith 2086f2d6a7bSJed Brown /* 2096f2d6a7bSJed Brown Compute Xdot = (X^{n+1}-X^n)/dt) = 0 2106f2d6a7bSJed Brown */ 211d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoGetXdot(TS ts, Vec X, Vec *Xdot) 212d71ae5a4SJacob Faibussowitsch { 2136f2d6a7bSJed Brown TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 214193ac0bcSJed Brown const PetscScalar mdt = 1.0 / ts->time_step, *xnp1, *xn; 215193ac0bcSJed Brown PetscScalar *xdot; 216a7cc72afSBarry Smith PetscInt i, n; 2172d3f70b5SBarry Smith 2183a40ed3dSBarry Smith PetscFunctionBegin; 219aab5bcd8SJed Brown *Xdot = NULL; 2209566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(ts->vec_sol, &xn)); 2219566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &xnp1)); 2229566063dSJacob Faibussowitsch PetscCall(VecGetArray(pseudo->xdot, &xdot)); 2239566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &n)); 224bbd56ea5SKarl Rupp for (i = 0; i < n; i++) xdot[i] = mdt * (xnp1[i] - xn[i]); 2259566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(ts->vec_sol, &xn)); 2269566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &xnp1)); 2279566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(pseudo->xdot, &xdot)); 2286f2d6a7bSJed Brown *Xdot = pseudo->xdot; 2293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2302d3f70b5SBarry Smith } 2312d3f70b5SBarry Smith 2326f2d6a7bSJed Brown /* 2336f2d6a7bSJed Brown The transient residual is 2346f2d6a7bSJed Brown 2356f2d6a7bSJed Brown F(U^{n+1},(U^{n+1}-U^n)/dt) = 0 2366f2d6a7bSJed Brown 2376f2d6a7bSJed Brown or for ODE, 2386f2d6a7bSJed Brown 2396f2d6a7bSJed Brown (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0 2406f2d6a7bSJed Brown 2416f2d6a7bSJed Brown This is the function that must be evaluated for transient simulation and for 2426f2d6a7bSJed Brown finite difference Jacobians. On the first Newton step, this algorithm uses 2436f2d6a7bSJed Brown a guess of U^{n+1} = U^n in which case the transient term vanishes and the 2446f2d6a7bSJed Brown residual is actually the steady state residual. Pseudotransient 2456f2d6a7bSJed Brown continuation as described in the literature is a linearly implicit 2466f2d6a7bSJed Brown algorithm, it only takes this one Newton step with the steady state 2476f2d6a7bSJed Brown residual, and then advances to the next time step. 2486f2d6a7bSJed Brown */ 249d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes, Vec X, Vec Y, TS ts) 250d71ae5a4SJacob Faibussowitsch { 2516f2d6a7bSJed Brown Vec Xdot; 2522d3f70b5SBarry Smith 2533a40ed3dSBarry Smith PetscFunctionBegin; 2549566063dSJacob Faibussowitsch PetscCall(TSPseudoGetXdot(ts, X, &Xdot)); 2559566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts, ts->ptime + ts->time_step, X, Xdot, Y, PETSC_FALSE)); 2563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2576f2d6a7bSJed Brown } 2582d3f70b5SBarry Smith 2596f2d6a7bSJed Brown /* 2606f2d6a7bSJed Brown This constructs the Jacobian needed for SNES. For DAE, this is 2616f2d6a7bSJed Brown 2626f2d6a7bSJed Brown dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot 2636f2d6a7bSJed Brown 2646f2d6a7bSJed Brown and for ODE: 2656f2d6a7bSJed Brown 2666f2d6a7bSJed Brown J = I/dt - J_{Frhs} where J_{Frhs} is the given Jacobian of Frhs. 2676f2d6a7bSJed Brown */ 268d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes, Vec X, Mat AA, Mat BB, TS ts) 269d71ae5a4SJacob Faibussowitsch { 2706f2d6a7bSJed Brown Vec Xdot; 2716f2d6a7bSJed Brown 2726f2d6a7bSJed Brown PetscFunctionBegin; 2739566063dSJacob Faibussowitsch PetscCall(TSPseudoGetXdot(ts, X, &Xdot)); 2749566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts, ts->ptime + ts->time_step, X, Xdot, 1. / ts->time_step, AA, BB, PETSC_FALSE)); 2753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2762d3f70b5SBarry Smith } 2772d3f70b5SBarry Smith 278d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_Pseudo(TS ts) 279d71ae5a4SJacob Faibussowitsch { 2807bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 2812d3f70b5SBarry Smith 2823a40ed3dSBarry Smith PetscFunctionBegin; 2839566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &pseudo->update)); 2849566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &pseudo->func)); 2859566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &pseudo->xdot)); 2863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2872d3f70b5SBarry Smith } 2882d3f70b5SBarry Smith /*------------------------------------------------------------*/ 2892d3f70b5SBarry Smith 290d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoMonitorDefault(TS ts, PetscInt step, PetscReal ptime, Vec v, void *dummy) 291d71ae5a4SJacob Faibussowitsch { 2927bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 293ce94432eSBarry Smith PetscViewer viewer = (PetscViewer)dummy; 2942d3f70b5SBarry Smith 2953a40ed3dSBarry Smith PetscFunctionBegin; 296193ac0bcSJed Brown if (pseudo->fnorm < 0) { /* The last computed norm is stale, recompute */ 2979566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(pseudo->xdot)); 2989566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts, ts->ptime, ts->vec_sol, pseudo->xdot, pseudo->func, PETSC_FALSE)); 2999566063dSJacob Faibussowitsch PetscCall(VecNorm(pseudo->func, NORM_2, &pseudo->fnorm)); 300193ac0bcSJed Brown } 3019566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel)); 30263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "TS %" PetscInt_FMT " dt %g time %g fnorm %g\n", step, (double)ts->time_step, (double)ptime, (double)pseudo->fnorm)); 3039566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel)); 3043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3052d3f70b5SBarry Smith } 3062d3f70b5SBarry Smith 307d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetFromOptions_Pseudo(TS ts, PetscOptionItems *PetscOptionsObject) 308d71ae5a4SJacob Faibussowitsch { 3094bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 310ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE; 311649052a6SBarry Smith PetscViewer viewer; 3122d3f70b5SBarry Smith 3133a40ed3dSBarry Smith PetscFunctionBegin; 314d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Pseudo-timestepping options"); 3159566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_monitor_pseudo", "Monitor convergence", "", flg, &flg, NULL)); 3162d3f70b5SBarry Smith if (flg) { 3179566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ts), "stdout", &viewer)); 3189566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts, TSPseudoMonitorDefault, viewer, (PetscErrorCode(*)(void **))PetscViewerDestroy)); 31928aa8177SBarry Smith } 320be5899b3SLisandro Dalcin flg = pseudo->increment_dt_from_initial_dt; 3219566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_pseudo_increment_dt_from_initial_dt", "Increase dt as a ratio from original dt", "TSPseudoIncrementDtFromInitialDt", flg, &flg, NULL)); 322be5899b3SLisandro Dalcin pseudo->increment_dt_from_initial_dt = flg; 3239566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_pseudo_increment", "Ratio to increase dt", "TSPseudoSetTimeStepIncrement", pseudo->dt_increment, &pseudo->dt_increment, NULL)); 3249566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_pseudo_max_dt", "Maximum value for dt", "TSPseudoSetMaxTimeStep", pseudo->dt_max, &pseudo->dt_max, NULL)); 3259566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_pseudo_fatol", "Tolerance for norm of function", "", pseudo->fatol, &pseudo->fatol, NULL)); 3269566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_pseudo_frtol", "Relative tolerance for norm of function", "", pseudo->frtol, &pseudo->frtol, NULL)); 327d0609cedSBarry Smith PetscOptionsHeadEnd(); 3283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3292d3f70b5SBarry Smith } 3302d3f70b5SBarry Smith 331d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_Pseudo(TS ts, PetscViewer viewer) 332d71ae5a4SJacob Faibussowitsch { 3333118ae5eSBarry Smith PetscBool isascii; 334d52bd9f3SBarry Smith 3353a40ed3dSBarry Smith PetscFunctionBegin; 3369566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 3373118ae5eSBarry Smith if (isascii) { 3383118ae5eSBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 3399566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " frtol - relative tolerance in function value %g\n", (double)pseudo->frtol)); 3409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " fatol - absolute tolerance in function value %g\n", (double)pseudo->fatol)); 3419566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " dt_initial - initial timestep %g\n", (double)pseudo->dt_initial)); 3429566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " dt_increment - increase in timestep on successful step %g\n", (double)pseudo->dt_increment)); 3439566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " dt_max - maximum time %g\n", (double)pseudo->dt_max)); 3443118ae5eSBarry Smith } 3453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3462d3f70b5SBarry Smith } 3472d3f70b5SBarry Smith 34882bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 349ac226902SBarry Smith /*@C 35082bf6240SBarry Smith TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the 35182bf6240SBarry Smith last timestep. 35282bf6240SBarry Smith 353c3339decSBarry Smith Logically Collective 35415091d37SBarry Smith 35582bf6240SBarry Smith Input Parameters: 35615091d37SBarry Smith + ts - timestep context 35782bf6240SBarry Smith . dt - user-defined function to verify timestep 3582fe279fdSBarry Smith - ctx - [optional] user-defined context for private data for the timestep verification routine (may be `NULL`) 35982bf6240SBarry Smith 36020f4b53cSBarry Smith Calling sequence of `func`: 3612fe279fdSBarry Smith + ts - the time-step context 3622fe279fdSBarry Smith . update - latest solution vector 36314d0ab18SJacob Faibussowitsch . ctx - [optional] user-defined timestep context 36482bf6240SBarry Smith . newdt - the timestep to use for the next step 365a2b725a8SWilliam Gropp - flag - flag indicating whether the last time step was acceptable 36682bf6240SBarry Smith 367bcf0153eSBarry Smith Level: advanced 368bcf0153eSBarry Smith 369bcf0153eSBarry Smith Note: 370bcf0153eSBarry Smith The routine set here will be called by `TSPseudoVerifyTimeStep()` 37182bf6240SBarry Smith during the timestepping process. 37282bf6240SBarry Smith 3731cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoVerifyTimeStepDefault()`, `TSPseudoVerifyTimeStep()` 37482bf6240SBarry Smith @*/ 37514d0ab18SJacob Faibussowitsch PetscErrorCode TSPseudoSetVerifyTimeStep(TS ts, PetscErrorCode (*dt)(TS ts, Vec update, void *ctx, PetscReal *newdt, PetscBool *flag), void *ctx) 376d71ae5a4SJacob Faibussowitsch { 37782bf6240SBarry Smith PetscFunctionBegin; 3780700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 379cac4c232SBarry Smith PetscTryMethod(ts, "TSPseudoSetVerifyTimeStep_C", (TS, PetscErrorCode(*)(TS, Vec, void *, PetscReal *, PetscBool *), void *), (ts, dt, ctx)); 3803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 38182bf6240SBarry Smith } 38282bf6240SBarry Smith 38382bf6240SBarry Smith /*@ 38482bf6240SBarry Smith TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 3858d359177SBarry Smith dt when using the TSPseudoTimeStepDefault() routine. 38682bf6240SBarry Smith 387c3339decSBarry Smith Logically Collective 388fee21e36SBarry Smith 38915091d37SBarry Smith Input Parameters: 39015091d37SBarry Smith + ts - the timestep context 39115091d37SBarry Smith - inc - the scaling factor >= 1.0 39215091d37SBarry Smith 39382bf6240SBarry Smith Options Database Key: 39467b8a455SSatish Balay . -ts_pseudo_increment <increment> - set pseudo increment 39582bf6240SBarry Smith 39615091d37SBarry Smith Level: advanced 39715091d37SBarry Smith 3981cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()` 39982bf6240SBarry Smith @*/ 400d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoSetTimeStepIncrement(TS ts, PetscReal inc) 401d71ae5a4SJacob Faibussowitsch { 40282bf6240SBarry Smith PetscFunctionBegin; 4030700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 404c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts, inc, 2); 405cac4c232SBarry Smith PetscTryMethod(ts, "TSPseudoSetTimeStepIncrement_C", (TS, PetscReal), (ts, inc)); 4063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 40782bf6240SBarry Smith } 40882bf6240SBarry Smith 40986534af1SJed Brown /*@ 41086534af1SJed Brown TSPseudoSetMaxTimeStep - Sets the maximum time step 4118d359177SBarry Smith when using the TSPseudoTimeStepDefault() routine. 41286534af1SJed Brown 413c3339decSBarry Smith Logically Collective 41486534af1SJed Brown 41586534af1SJed Brown Input Parameters: 41686534af1SJed Brown + ts - the timestep context 41786534af1SJed Brown - maxdt - the maximum time step, use a non-positive value to deactivate 41886534af1SJed Brown 41986534af1SJed Brown Options Database Key: 42067b8a455SSatish Balay . -ts_pseudo_max_dt <increment> - set pseudo max dt 42186534af1SJed Brown 42286534af1SJed Brown Level: advanced 42386534af1SJed Brown 4241cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()` 42586534af1SJed Brown @*/ 426d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoSetMaxTimeStep(TS ts, PetscReal maxdt) 427d71ae5a4SJacob Faibussowitsch { 42886534af1SJed Brown PetscFunctionBegin; 42986534af1SJed Brown PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 43086534af1SJed Brown PetscValidLogicalCollectiveReal(ts, maxdt, 2); 431cac4c232SBarry Smith PetscTryMethod(ts, "TSPseudoSetMaxTimeStep_C", (TS, PetscReal), (ts, maxdt)); 4323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 43386534af1SJed Brown } 43486534af1SJed Brown 43582bf6240SBarry Smith /*@ 43682bf6240SBarry Smith TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep 43782bf6240SBarry Smith is computed via the formula 43882bf6240SBarry Smith $ dt = initial_dt*initial_fnorm/current_fnorm 43982bf6240SBarry Smith rather than the default update, 44082bf6240SBarry Smith $ dt = current_dt*previous_fnorm/current_fnorm. 44182bf6240SBarry Smith 442c3339decSBarry Smith Logically Collective 44315091d37SBarry Smith 44482bf6240SBarry Smith Input Parameter: 44582bf6240SBarry Smith . ts - the timestep context 44682bf6240SBarry Smith 44782bf6240SBarry Smith Options Database Key: 44867b8a455SSatish Balay . -ts_pseudo_increment_dt_from_initial_dt <true,false> - use the initial dt to determine increment 44982bf6240SBarry Smith 45015091d37SBarry Smith Level: advanced 45115091d37SBarry Smith 4521cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()` 45382bf6240SBarry Smith @*/ 454d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS ts) 455d71ae5a4SJacob Faibussowitsch { 45682bf6240SBarry Smith PetscFunctionBegin; 4570700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 458cac4c232SBarry Smith PetscTryMethod(ts, "TSPseudoIncrementDtFromInitialDt_C", (TS), (ts)); 4593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 46082bf6240SBarry Smith } 46182bf6240SBarry Smith 462ac226902SBarry Smith /*@C 46382bf6240SBarry Smith TSPseudoSetTimeStep - Sets the user-defined routine to be 46482bf6240SBarry Smith called at each pseudo-timestep to update the timestep. 46582bf6240SBarry Smith 466c3339decSBarry Smith Logically Collective 46715091d37SBarry Smith 46882bf6240SBarry Smith Input Parameters: 46915091d37SBarry Smith + ts - timestep context 47082bf6240SBarry Smith . dt - function to compute timestep 4712fe279fdSBarry Smith - ctx - [optional] user-defined context for private data required by the function (may be `NULL`) 47282bf6240SBarry Smith 47320f4b53cSBarry Smith Calling sequence of `dt`: 47414d0ab18SJacob Faibussowitsch + ts - the `TS` context 47514d0ab18SJacob Faibussowitsch . newdt - the newly computed timestep 47614d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined context 47782bf6240SBarry Smith 478bcf0153eSBarry Smith Level: intermediate 47982bf6240SBarry Smith 480bcf0153eSBarry Smith Notes: 481bcf0153eSBarry Smith The routine set here will be called by `TSPseudoComputeTimeStep()` 482bcf0153eSBarry Smith during the timestepping process. 483bcf0153eSBarry Smith 484bcf0153eSBarry Smith If not set then `TSPseudoTimeStepDefault()` is automatically used 485bcf0153eSBarry Smith 4861cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoTimeStepDefault()`, `TSPseudoComputeTimeStep()` 48782bf6240SBarry Smith @*/ 48814d0ab18SJacob Faibussowitsch PetscErrorCode TSPseudoSetTimeStep(TS ts, PetscErrorCode (*dt)(TS ts, PetscReal *newdt, void *ctx), void *ctx) 489d71ae5a4SJacob Faibussowitsch { 49082bf6240SBarry Smith PetscFunctionBegin; 4910700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 492cac4c232SBarry Smith PetscTryMethod(ts, "TSPseudoSetTimeStep_C", (TS, PetscErrorCode(*)(TS, PetscReal *, void *), void *), (ts, dt, ctx)); 4933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 49482bf6240SBarry Smith } 49582bf6240SBarry Smith 49682bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 49782bf6240SBarry Smith 498ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(TS, Vec, void *, PetscReal *, PetscBool *); /* force argument to next function to not be extern C*/ 499d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetVerifyTimeStep_Pseudo(TS ts, FCN1 dt, void *ctx) 500d71ae5a4SJacob Faibussowitsch { 501be5899b3SLisandro Dalcin TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 50282bf6240SBarry Smith 50382bf6240SBarry Smith PetscFunctionBegin; 50482bf6240SBarry Smith pseudo->verify = dt; 50582bf6240SBarry Smith pseudo->verifyctx = ctx; 5063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 50782bf6240SBarry Smith } 50882bf6240SBarry Smith 509d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetTimeStepIncrement_Pseudo(TS ts, PetscReal inc) 510d71ae5a4SJacob Faibussowitsch { 5114bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 51282bf6240SBarry Smith 51382bf6240SBarry Smith PetscFunctionBegin; 51482bf6240SBarry Smith pseudo->dt_increment = inc; 5153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 51682bf6240SBarry Smith } 51782bf6240SBarry Smith 518d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetMaxTimeStep_Pseudo(TS ts, PetscReal maxdt) 519d71ae5a4SJacob Faibussowitsch { 52086534af1SJed Brown TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 52186534af1SJed Brown 52286534af1SJed Brown PetscFunctionBegin; 52386534af1SJed Brown pseudo->dt_max = maxdt; 5243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 52586534af1SJed Brown } 52686534af1SJed Brown 527d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts) 528d71ae5a4SJacob Faibussowitsch { 5294bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 53082bf6240SBarry Smith 53182bf6240SBarry Smith PetscFunctionBegin; 5324bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_TRUE; 5333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 53482bf6240SBarry Smith } 53582bf6240SBarry Smith 5366849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS, PetscReal *, void *); /* force argument to next function to not be extern C*/ 537d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetTimeStep_Pseudo(TS ts, FCN2 dt, void *ctx) 538d71ae5a4SJacob Faibussowitsch { 5394bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 54082bf6240SBarry Smith 54182bf6240SBarry Smith PetscFunctionBegin; 54282bf6240SBarry Smith pseudo->dt = dt; 54382bf6240SBarry Smith pseudo->dtctx = ctx; 5443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 54582bf6240SBarry Smith } 54682bf6240SBarry Smith 54710e6a065SJed Brown /*MC 548*1d27aa22SBarry Smith TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping {cite}`ckk02` {cite}`kk97` 54982bf6240SBarry Smith 55010e6a065SJed Brown This method solves equations of the form 55110e6a065SJed Brown 552*1d27aa22SBarry Smith $$ 553*1d27aa22SBarry Smith F(X,Xdot) = 0 554*1d27aa22SBarry Smith $$ 55510e6a065SJed Brown 55610e6a065SJed Brown for steady state using the iteration 55710e6a065SJed Brown 558*1d27aa22SBarry Smith $$ 55937fdd005SBarry Smith [G'] S = -F(X,0) 56037fdd005SBarry Smith X += S 561*1d27aa22SBarry Smith $$ 56210e6a065SJed Brown 56310e6a065SJed Brown where 56410e6a065SJed Brown 565*1d27aa22SBarry Smith $$ 566*1d27aa22SBarry Smith G(Y) = F(Y,(Y-X)/dt) 567*1d27aa22SBarry Smith $$ 56810e6a065SJed Brown 5696f2d6a7bSJed Brown This is linearly-implicit Euler with the residual always evaluated "at steady 5706f2d6a7bSJed Brown state". See note below. 57110e6a065SJed Brown 572bcf0153eSBarry Smith Options Database Keys: 57310e6a065SJed Brown + -ts_pseudo_increment <real> - ratio of increase dt 5743118ae5eSBarry Smith . -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt 5753118ae5eSBarry Smith . -ts_pseudo_fatol <atol> - stop iterating when the function norm is less than atol 5763118ae5eSBarry Smith - -ts_pseudo_frtol <rtol> - stop iterating when the function norm divided by the initial function norm is less than rtol 57710e6a065SJed Brown 57810e6a065SJed Brown Level: beginner 57910e6a065SJed Brown 58010e6a065SJed Brown Notes: 5816f2d6a7bSJed Brown The residual computed by this method includes the transient term (Xdot is computed instead of 5826f2d6a7bSJed Brown always being zero), but since the prediction from the last step is always the solution from the 5836f2d6a7bSJed Brown last step, on the first Newton iteration we have 5846f2d6a7bSJed Brown 585*1d27aa22SBarry Smith $$ 586*1d27aa22SBarry Smith Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0 587*1d27aa22SBarry Smith $$ 5886f2d6a7bSJed Brown 5896f2d6a7bSJed Brown Therefore, the linear system solved by the first Newton iteration is equivalent to the one 5906f2d6a7bSJed Brown described above and in the papers. If the user chooses to perform multiple Newton iterations, the 5916f2d6a7bSJed Brown algorithm is no longer the one described in the referenced papers. 59210e6a065SJed Brown 5931cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()` 59410e6a065SJed Brown M*/ 595d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_Pseudo(TS ts) 596d71ae5a4SJacob Faibussowitsch { 5977bf11e45SBarry Smith TS_Pseudo *pseudo; 598193ac0bcSJed Brown SNES snes; 59919fd82e9SBarry Smith SNESType stype; 6002d3f70b5SBarry Smith 6013a40ed3dSBarry Smith PetscFunctionBegin; 602277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Pseudo; 603000e7ae3SMatthew Knepley ts->ops->destroy = TSDestroy_Pseudo; 604000e7ae3SMatthew Knepley ts->ops->view = TSView_Pseudo; 605000e7ae3SMatthew Knepley ts->ops->setup = TSSetUp_Pseudo; 606000e7ae3SMatthew Knepley ts->ops->step = TSStep_Pseudo; 607000e7ae3SMatthew Knepley ts->ops->setfromoptions = TSSetFromOptions_Pseudo; 6080f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Pseudo; 6090f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Pseudo; 6102ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 611825ab935SBarry Smith ts->usessnes = PETSC_TRUE; 6127bf11e45SBarry Smith 6139566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 6149566063dSJacob Faibussowitsch PetscCall(SNESGetType(snes, &stype)); 6159566063dSJacob Faibussowitsch if (!stype) PetscCall(SNESSetType(snes, SNESKSPONLY)); 6162d3f70b5SBarry Smith 6174dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&pseudo)); 6187bf11e45SBarry Smith ts->data = (void *)pseudo; 6192d3f70b5SBarry Smith 620be5899b3SLisandro Dalcin pseudo->dt = TSPseudoTimeStepDefault; 621be5899b3SLisandro Dalcin pseudo->dtctx = NULL; 62228aa8177SBarry Smith pseudo->dt_increment = 1.1; 6234bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_FALSE; 624193ac0bcSJed Brown pseudo->fnorm = -1; 625be5899b3SLisandro Dalcin pseudo->fnorm_initial = -1; 626be5899b3SLisandro Dalcin pseudo->fnorm_previous = -1; 6273118ae5eSBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 6283118ae5eSBarry Smith pseudo->fatol = 1.e-25; 6293118ae5eSBarry Smith pseudo->frtol = 1.e-5; 6303118ae5eSBarry Smith #else 6313118ae5eSBarry Smith pseudo->fatol = 1.e-50; 6323118ae5eSBarry Smith pseudo->frtol = 1.e-12; 6333118ae5eSBarry Smith #endif 6349566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetVerifyTimeStep_C", TSPseudoSetVerifyTimeStep_Pseudo)); 6359566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStepIncrement_C", TSPseudoSetTimeStepIncrement_Pseudo)); 6369566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetMaxTimeStep_C", TSPseudoSetMaxTimeStep_Pseudo)); 6379566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoIncrementDtFromInitialDt_C", TSPseudoIncrementDtFromInitialDt_Pseudo)); 6389566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStep_C", TSPseudoSetTimeStep_Pseudo)); 6393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6402d3f70b5SBarry Smith } 6412d3f70b5SBarry Smith 64282bf6240SBarry Smith /*@C 643bcf0153eSBarry Smith TSPseudoTimeStepDefault - Default code to compute pseudo-timestepping. Use with `TSPseudoSetTimeStep()`. 64428aa8177SBarry Smith 645c3339decSBarry Smith Collective 64615091d37SBarry Smith 64728aa8177SBarry Smith Input Parameters: 648a2b725a8SWilliam Gropp + ts - the timestep context 649a2b725a8SWilliam Gropp - dtctx - unused timestep context 65028aa8177SBarry Smith 65182bf6240SBarry Smith Output Parameter: 65282bf6240SBarry Smith . newdt - the timestep to use for the next step 65328aa8177SBarry Smith 65415091d37SBarry Smith Level: advanced 65515091d37SBarry Smith 6561cc06b55SBarry Smith .seealso: [](ch_ts), `TSPseudoSetTimeStep()`, `TSPseudoComputeTimeStep()`, `TSPSEUDO` 65728aa8177SBarry Smith @*/ 658d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoTimeStepDefault(TS ts, PetscReal *newdt, void *dtctx) 659d71ae5a4SJacob Faibussowitsch { 66082bf6240SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 661be5899b3SLisandro Dalcin PetscReal inc = pseudo->dt_increment; 66228aa8177SBarry Smith 6633a40ed3dSBarry Smith PetscFunctionBegin; 6649566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(pseudo->xdot)); 6659566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts, ts->ptime, ts->vec_sol, pseudo->xdot, pseudo->func, PETSC_FALSE)); 6669566063dSJacob Faibussowitsch PetscCall(VecNorm(pseudo->func, NORM_2, &pseudo->fnorm)); 667be5899b3SLisandro Dalcin if (pseudo->fnorm_initial < 0) { 66882bf6240SBarry Smith /* first time through so compute initial function norm */ 669cdbf8f93SLisandro Dalcin pseudo->fnorm_initial = pseudo->fnorm; 670be5899b3SLisandro Dalcin pseudo->fnorm_previous = pseudo->fnorm; 67182bf6240SBarry Smith } 672bbd56ea5SKarl Rupp if (pseudo->fnorm == 0.0) *newdt = 1.e12 * inc * ts->time_step; 673bbd56ea5SKarl Rupp else if (pseudo->increment_dt_from_initial_dt) *newdt = inc * pseudo->dt_initial * pseudo->fnorm_initial / pseudo->fnorm; 674be5899b3SLisandro Dalcin else *newdt = inc * ts->time_step * pseudo->fnorm_previous / pseudo->fnorm; 67586534af1SJed Brown if (pseudo->dt_max > 0) *newdt = PetscMin(*newdt, pseudo->dt_max); 67682bf6240SBarry Smith pseudo->fnorm_previous = pseudo->fnorm; 6773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 67828aa8177SBarry Smith } 679