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 30cc4c1da9SBarry Smith /*@ 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 65cc4c1da9SBarry Smith Collective, No Fortran Support 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; 143f4f49eeaSPierre Jolivet 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 307ce78bad3SBarry Smith 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)); 31849abdd8aSBarry Smith PetscCall(TSMonitorSet(ts, TSPseudoMonitorDefault, viewer, (PetscCtxDestroyFn *)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 348ac226902SBarry Smith /*@C 34982bf6240SBarry Smith TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the 35082bf6240SBarry Smith last timestep. 35182bf6240SBarry Smith 352c3339decSBarry Smith Logically Collective 35315091d37SBarry Smith 35482bf6240SBarry Smith Input Parameters: 35515091d37SBarry Smith + ts - timestep context 35682bf6240SBarry Smith . dt - user-defined function to verify timestep 3572fe279fdSBarry Smith - ctx - [optional] user-defined context for private data for the timestep verification routine (may be `NULL`) 35882bf6240SBarry Smith 35920f4b53cSBarry Smith Calling sequence of `func`: 3602fe279fdSBarry Smith + ts - the time-step context 3612fe279fdSBarry Smith . update - latest solution vector 36214d0ab18SJacob Faibussowitsch . ctx - [optional] user-defined timestep context 36382bf6240SBarry Smith . newdt - the timestep to use for the next step 364a2b725a8SWilliam Gropp - flag - flag indicating whether the last time step was acceptable 36582bf6240SBarry Smith 366bcf0153eSBarry Smith Level: advanced 367bcf0153eSBarry Smith 368bcf0153eSBarry Smith Note: 369bcf0153eSBarry Smith The routine set here will be called by `TSPseudoVerifyTimeStep()` 37082bf6240SBarry Smith during the timestepping process. 37182bf6240SBarry Smith 3721cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoVerifyTimeStepDefault()`, `TSPseudoVerifyTimeStep()` 37382bf6240SBarry Smith @*/ 37414d0ab18SJacob Faibussowitsch PetscErrorCode TSPseudoSetVerifyTimeStep(TS ts, PetscErrorCode (*dt)(TS ts, Vec update, void *ctx, PetscReal *newdt, PetscBool *flag), void *ctx) 375d71ae5a4SJacob Faibussowitsch { 37682bf6240SBarry Smith PetscFunctionBegin; 3770700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 378cac4c232SBarry Smith PetscTryMethod(ts, "TSPseudoSetVerifyTimeStep_C", (TS, PetscErrorCode (*)(TS, Vec, void *, PetscReal *, PetscBool *), void *), (ts, dt, ctx)); 3793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 38082bf6240SBarry Smith } 38182bf6240SBarry Smith 38282bf6240SBarry Smith /*@ 38382bf6240SBarry Smith TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 3848d359177SBarry Smith dt when using the TSPseudoTimeStepDefault() routine. 38582bf6240SBarry Smith 386c3339decSBarry Smith Logically Collective 387fee21e36SBarry Smith 38815091d37SBarry Smith Input Parameters: 38915091d37SBarry Smith + ts - the timestep context 39015091d37SBarry Smith - inc - the scaling factor >= 1.0 39115091d37SBarry Smith 39282bf6240SBarry Smith Options Database Key: 39367b8a455SSatish Balay . -ts_pseudo_increment <increment> - set pseudo increment 39482bf6240SBarry Smith 39515091d37SBarry Smith Level: advanced 39615091d37SBarry Smith 3971cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()` 39882bf6240SBarry Smith @*/ 399d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoSetTimeStepIncrement(TS ts, PetscReal inc) 400d71ae5a4SJacob Faibussowitsch { 40182bf6240SBarry Smith PetscFunctionBegin; 4020700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 403c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts, inc, 2); 404cac4c232SBarry Smith PetscTryMethod(ts, "TSPseudoSetTimeStepIncrement_C", (TS, PetscReal), (ts, inc)); 4053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 40682bf6240SBarry Smith } 40782bf6240SBarry Smith 40886534af1SJed Brown /*@ 40986534af1SJed Brown TSPseudoSetMaxTimeStep - Sets the maximum time step 4108d359177SBarry Smith when using the TSPseudoTimeStepDefault() routine. 41186534af1SJed Brown 412c3339decSBarry Smith Logically Collective 41386534af1SJed Brown 41486534af1SJed Brown Input Parameters: 41586534af1SJed Brown + ts - the timestep context 41686534af1SJed Brown - maxdt - the maximum time step, use a non-positive value to deactivate 41786534af1SJed Brown 41886534af1SJed Brown Options Database Key: 41967b8a455SSatish Balay . -ts_pseudo_max_dt <increment> - set pseudo max dt 42086534af1SJed Brown 42186534af1SJed Brown Level: advanced 42286534af1SJed Brown 4231cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()` 42486534af1SJed Brown @*/ 425d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoSetMaxTimeStep(TS ts, PetscReal maxdt) 426d71ae5a4SJacob Faibussowitsch { 42786534af1SJed Brown PetscFunctionBegin; 42886534af1SJed Brown PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 42986534af1SJed Brown PetscValidLogicalCollectiveReal(ts, maxdt, 2); 430cac4c232SBarry Smith PetscTryMethod(ts, "TSPseudoSetMaxTimeStep_C", (TS, PetscReal), (ts, maxdt)); 4313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 43286534af1SJed Brown } 43386534af1SJed Brown 43482bf6240SBarry Smith /*@ 43582bf6240SBarry Smith TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep 436b44f4de4SBarry 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.$ 43782bf6240SBarry Smith 438c3339decSBarry Smith Logically Collective 43915091d37SBarry Smith 44082bf6240SBarry Smith Input Parameter: 44182bf6240SBarry Smith . ts - the timestep context 44282bf6240SBarry Smith 44382bf6240SBarry Smith Options Database Key: 444b44f4de4SBarry Smith . -ts_pseudo_increment_dt_from_initial_dt <true,false> - use the initial $dt$ to determine increment 44582bf6240SBarry Smith 44615091d37SBarry Smith Level: advanced 44715091d37SBarry Smith 4481cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()` 44982bf6240SBarry Smith @*/ 450d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS ts) 451d71ae5a4SJacob Faibussowitsch { 45282bf6240SBarry Smith PetscFunctionBegin; 4530700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 454cac4c232SBarry Smith PetscTryMethod(ts, "TSPseudoIncrementDtFromInitialDt_C", (TS), (ts)); 4553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 45682bf6240SBarry Smith } 45782bf6240SBarry Smith 458ac226902SBarry Smith /*@C 45982bf6240SBarry Smith TSPseudoSetTimeStep - Sets the user-defined routine to be 46082bf6240SBarry Smith called at each pseudo-timestep to update the timestep. 46182bf6240SBarry Smith 462c3339decSBarry Smith Logically Collective 46315091d37SBarry Smith 46482bf6240SBarry Smith Input Parameters: 46515091d37SBarry Smith + ts - timestep context 46682bf6240SBarry Smith . dt - function to compute timestep 4672fe279fdSBarry Smith - ctx - [optional] user-defined context for private data required by the function (may be `NULL`) 46882bf6240SBarry Smith 46920f4b53cSBarry Smith Calling sequence of `dt`: 47014d0ab18SJacob Faibussowitsch + ts - the `TS` context 47114d0ab18SJacob Faibussowitsch . newdt - the newly computed timestep 47214d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined context 47382bf6240SBarry Smith 474bcf0153eSBarry Smith Level: intermediate 47582bf6240SBarry Smith 476bcf0153eSBarry Smith Notes: 477bcf0153eSBarry Smith The routine set here will be called by `TSPseudoComputeTimeStep()` 478bcf0153eSBarry Smith during the timestepping process. 479bcf0153eSBarry Smith 480bcf0153eSBarry Smith If not set then `TSPseudoTimeStepDefault()` is automatically used 481bcf0153eSBarry Smith 4821cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoTimeStepDefault()`, `TSPseudoComputeTimeStep()` 48382bf6240SBarry Smith @*/ 48414d0ab18SJacob Faibussowitsch PetscErrorCode TSPseudoSetTimeStep(TS ts, PetscErrorCode (*dt)(TS ts, PetscReal *newdt, void *ctx), void *ctx) 485d71ae5a4SJacob Faibussowitsch { 48682bf6240SBarry Smith PetscFunctionBegin; 4870700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 488cac4c232SBarry Smith PetscTryMethod(ts, "TSPseudoSetTimeStep_C", (TS, PetscErrorCode (*)(TS, PetscReal *, void *), void *), (ts, dt, ctx)); 4893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 49082bf6240SBarry Smith } 49182bf6240SBarry Smith 49282bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 49382bf6240SBarry Smith 494ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(TS, Vec, void *, PetscReal *, PetscBool *); /* force argument to next function to not be extern C*/ 495d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetVerifyTimeStep_Pseudo(TS ts, FCN1 dt, void *ctx) 496d71ae5a4SJacob Faibussowitsch { 497be5899b3SLisandro Dalcin TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 49882bf6240SBarry Smith 49982bf6240SBarry Smith PetscFunctionBegin; 50082bf6240SBarry Smith pseudo->verify = dt; 50182bf6240SBarry Smith pseudo->verifyctx = ctx; 5023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 50382bf6240SBarry Smith } 50482bf6240SBarry Smith 505d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetTimeStepIncrement_Pseudo(TS ts, PetscReal inc) 506d71ae5a4SJacob Faibussowitsch { 5074bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 50882bf6240SBarry Smith 50982bf6240SBarry Smith PetscFunctionBegin; 51082bf6240SBarry Smith pseudo->dt_increment = inc; 5113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 51282bf6240SBarry Smith } 51382bf6240SBarry Smith 514d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetMaxTimeStep_Pseudo(TS ts, PetscReal maxdt) 515d71ae5a4SJacob Faibussowitsch { 51686534af1SJed Brown TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 51786534af1SJed Brown 51886534af1SJed Brown PetscFunctionBegin; 51986534af1SJed Brown pseudo->dt_max = maxdt; 5203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 52186534af1SJed Brown } 52286534af1SJed Brown 523d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts) 524d71ae5a4SJacob Faibussowitsch { 5254bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 52682bf6240SBarry Smith 52782bf6240SBarry Smith PetscFunctionBegin; 5284bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_TRUE; 5293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 53082bf6240SBarry Smith } 53182bf6240SBarry Smith 5326849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS, PetscReal *, void *); /* force argument to next function to not be extern C*/ 533d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetTimeStep_Pseudo(TS ts, FCN2 dt, void *ctx) 534d71ae5a4SJacob Faibussowitsch { 5354bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 53682bf6240SBarry Smith 53782bf6240SBarry Smith PetscFunctionBegin; 53882bf6240SBarry Smith pseudo->dt = dt; 53982bf6240SBarry Smith pseudo->dtctx = ctx; 5403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 54182bf6240SBarry Smith } 54282bf6240SBarry Smith 54310e6a065SJed Brown /*MC 5441d27aa22SBarry Smith TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping {cite}`ckk02` {cite}`kk97` 54582bf6240SBarry Smith 54610e6a065SJed Brown This method solves equations of the form 54710e6a065SJed Brown 5481d27aa22SBarry Smith $$ 5491d27aa22SBarry Smith F(X,Xdot) = 0 5501d27aa22SBarry Smith $$ 55110e6a065SJed Brown 55210e6a065SJed Brown for steady state using the iteration 55310e6a065SJed Brown 5541d27aa22SBarry Smith $$ 55537fdd005SBarry Smith [G'] S = -F(X,0) 55637fdd005SBarry Smith X += S 5571d27aa22SBarry Smith $$ 55810e6a065SJed Brown 55910e6a065SJed Brown where 56010e6a065SJed Brown 5611d27aa22SBarry Smith $$ 5621d27aa22SBarry Smith G(Y) = F(Y,(Y-X)/dt) 5631d27aa22SBarry Smith $$ 56410e6a065SJed Brown 5656f2d6a7bSJed Brown This is linearly-implicit Euler with the residual always evaluated "at steady 5666f2d6a7bSJed Brown state". See note below. 56710e6a065SJed Brown 568*fb59f417SJames Wright In addition to the modified solve, a dedicated adaptive timestepping scheme is implemented, mimicing the switched evolution relaxation in {cite}`ckk02`. 569*fb59f417SJames Wright It determines the next timestep via 570*fb59f417SJames Wright 571*fb59f417SJames Wright $$ 572*fb59f417SJames Wright dt_{n} = r dt_{n-1} \frac{\Vert F(X_{n},0) \Vert}{\Vert F(X_{n-1},0) \Vert} 573*fb59f417SJames Wright $$ 574*fb59f417SJames Wright 575*fb59f417SJames Wright where $r$ is an additional growth factor (set by `-ts_pseudo_increment`). 576*fb59f417SJames Wright An alternative formulation is also available that uses the initial timestep and function norm. 577*fb59f417SJames Wright 578*fb59f417SJames Wright $$ 579*fb59f417SJames Wright dt_{n} = r dt_{0} \frac{\Vert F(X_{n},0) \Vert}{\Vert F(X_{0},0) \Vert} 580*fb59f417SJames Wright $$ 581*fb59f417SJames Wright 582*fb59f417SJames Wright This is chosen by setting `-ts_pseudo_increment_dt_from_initial_dt`. 583*fb59f417SJames Wright For either method, an upper limit on the timestep can be set by `-ts_pseudo_max_dt`. 584*fb59f417SJames Wright 585bcf0153eSBarry Smith Options Database Keys: 58610e6a065SJed Brown + -ts_pseudo_increment <real> - ratio of increase dt 5873118ae5eSBarry Smith . -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt 588*fb59f417SJames Wright . -ts_pseudo_max_dt - Maximum dt for adaptive timestepping algorithm 589*fb59f417SJames Wright . -ts_pseudo_monitor - Monitor convergence of the function norm 590*fb59f417SJames Wright . -ts_pseudo_fatol <atol> - stop iterating when the function norm is less than `atol` 591*fb59f417SJames Wright - -ts_pseudo_frtol <rtol> - stop iterating when the function norm divided by the initial function norm is less than `rtol` 59210e6a065SJed Brown 59310e6a065SJed Brown Level: beginner 59410e6a065SJed Brown 59510e6a065SJed Brown Notes: 5966f2d6a7bSJed Brown The residual computed by this method includes the transient term (Xdot is computed instead of 5976f2d6a7bSJed Brown always being zero), but since the prediction from the last step is always the solution from the 5986f2d6a7bSJed Brown last step, on the first Newton iteration we have 5996f2d6a7bSJed Brown 6001d27aa22SBarry Smith $$ 6011d27aa22SBarry Smith Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0 6021d27aa22SBarry Smith $$ 6036f2d6a7bSJed Brown 6046f2d6a7bSJed Brown Therefore, the linear system solved by the first Newton iteration is equivalent to the one 6056f2d6a7bSJed Brown described above and in the papers. If the user chooses to perform multiple Newton iterations, the 6066f2d6a7bSJed Brown algorithm is no longer the one described in the referenced papers. 607*fb59f417SJames Wright By default, the `SNESType` is set to `SNESKSPONLY` to match the algorithm from the referenced papers. 608*fb59f417SJames Wright Setting the `SNESType` via `-snes_type` will override this default setting. 60910e6a065SJed Brown 6101cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()` 61110e6a065SJed Brown M*/ 612d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_Pseudo(TS ts) 613d71ae5a4SJacob Faibussowitsch { 6147bf11e45SBarry Smith TS_Pseudo *pseudo; 615193ac0bcSJed Brown SNES snes; 61619fd82e9SBarry Smith SNESType stype; 6172d3f70b5SBarry Smith 6183a40ed3dSBarry Smith PetscFunctionBegin; 619277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Pseudo; 620000e7ae3SMatthew Knepley ts->ops->destroy = TSDestroy_Pseudo; 621000e7ae3SMatthew Knepley ts->ops->view = TSView_Pseudo; 622000e7ae3SMatthew Knepley ts->ops->setup = TSSetUp_Pseudo; 623000e7ae3SMatthew Knepley ts->ops->step = TSStep_Pseudo; 624000e7ae3SMatthew Knepley ts->ops->setfromoptions = TSSetFromOptions_Pseudo; 6250f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Pseudo; 6260f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Pseudo; 6272ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 628825ab935SBarry Smith ts->usessnes = PETSC_TRUE; 6297bf11e45SBarry Smith 6309566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 6319566063dSJacob Faibussowitsch PetscCall(SNESGetType(snes, &stype)); 6329566063dSJacob Faibussowitsch if (!stype) PetscCall(SNESSetType(snes, SNESKSPONLY)); 6332d3f70b5SBarry Smith 6344dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&pseudo)); 6357bf11e45SBarry Smith ts->data = (void *)pseudo; 6362d3f70b5SBarry Smith 637be5899b3SLisandro Dalcin pseudo->dt = TSPseudoTimeStepDefault; 638be5899b3SLisandro Dalcin pseudo->dtctx = NULL; 63928aa8177SBarry Smith pseudo->dt_increment = 1.1; 6404bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_FALSE; 641193ac0bcSJed Brown pseudo->fnorm = -1; 642be5899b3SLisandro Dalcin pseudo->fnorm_initial = -1; 643be5899b3SLisandro Dalcin pseudo->fnorm_previous = -1; 6443118ae5eSBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 6453118ae5eSBarry Smith pseudo->fatol = 1.e-25; 6463118ae5eSBarry Smith pseudo->frtol = 1.e-5; 6473118ae5eSBarry Smith #else 6483118ae5eSBarry Smith pseudo->fatol = 1.e-50; 6493118ae5eSBarry Smith pseudo->frtol = 1.e-12; 6503118ae5eSBarry Smith #endif 6519566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetVerifyTimeStep_C", TSPseudoSetVerifyTimeStep_Pseudo)); 6529566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStepIncrement_C", TSPseudoSetTimeStepIncrement_Pseudo)); 6539566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetMaxTimeStep_C", TSPseudoSetMaxTimeStep_Pseudo)); 6549566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoIncrementDtFromInitialDt_C", TSPseudoIncrementDtFromInitialDt_Pseudo)); 6559566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStep_C", TSPseudoSetTimeStep_Pseudo)); 6563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6572d3f70b5SBarry Smith } 6582d3f70b5SBarry Smith 65982bf6240SBarry Smith /*@C 660bcf0153eSBarry Smith TSPseudoTimeStepDefault - Default code to compute pseudo-timestepping. Use with `TSPseudoSetTimeStep()`. 66128aa8177SBarry Smith 662cc4c1da9SBarry Smith Collective, No Fortran Support 66315091d37SBarry Smith 66428aa8177SBarry Smith Input Parameters: 665a2b725a8SWilliam Gropp + ts - the timestep context 666a2b725a8SWilliam Gropp - dtctx - unused timestep context 66728aa8177SBarry Smith 66882bf6240SBarry Smith Output Parameter: 66982bf6240SBarry Smith . newdt - the timestep to use for the next step 67028aa8177SBarry Smith 67115091d37SBarry Smith Level: advanced 67215091d37SBarry Smith 6731cc06b55SBarry Smith .seealso: [](ch_ts), `TSPseudoSetTimeStep()`, `TSPseudoComputeTimeStep()`, `TSPSEUDO` 67428aa8177SBarry Smith @*/ 675d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoTimeStepDefault(TS ts, PetscReal *newdt, void *dtctx) 676d71ae5a4SJacob Faibussowitsch { 67782bf6240SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 678be5899b3SLisandro Dalcin PetscReal inc = pseudo->dt_increment; 67928aa8177SBarry Smith 6803a40ed3dSBarry Smith PetscFunctionBegin; 6819566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(pseudo->xdot)); 6829566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts, ts->ptime, ts->vec_sol, pseudo->xdot, pseudo->func, PETSC_FALSE)); 6839566063dSJacob Faibussowitsch PetscCall(VecNorm(pseudo->func, NORM_2, &pseudo->fnorm)); 684be5899b3SLisandro Dalcin if (pseudo->fnorm_initial < 0) { 68582bf6240SBarry Smith /* first time through so compute initial function norm */ 686cdbf8f93SLisandro Dalcin pseudo->fnorm_initial = pseudo->fnorm; 687be5899b3SLisandro Dalcin pseudo->fnorm_previous = pseudo->fnorm; 68882bf6240SBarry Smith } 689bbd56ea5SKarl Rupp if (pseudo->fnorm == 0.0) *newdt = 1.e12 * inc * ts->time_step; 690bbd56ea5SKarl Rupp else if (pseudo->increment_dt_from_initial_dt) *newdt = inc * pseudo->dt_initial * pseudo->fnorm_initial / pseudo->fnorm; 691be5899b3SLisandro Dalcin else *newdt = inc * ts->time_step * pseudo->fnorm_previous / pseudo->fnorm; 69286534af1SJed Brown if (pseudo->dt_max > 0) *newdt = PetscMin(*newdt, pseudo->dt_max); 69382bf6240SBarry Smith pseudo->fnorm_previous = pseudo->fnorm; 6943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 69528aa8177SBarry Smith } 696