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 3415091d37SBarry Smith Collective on TS 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 44564e8f4eSLois Curfman McInnes Notes: 45564e8f4eSLois Curfman McInnes The routine to be called here to compute the timestep should be 46564e8f4eSLois Curfman McInnes set by calling TSPseudoSetTimeStep(). 47564e8f4eSLois Curfman McInnes 48db781477SPatrick Sanan .seealso: `TSPseudoTimeStepDefault()`, `TSPseudoSetTimeStep()` 497bf11e45SBarry Smith @*/ 50*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoComputeTimeStep(TS ts, PetscReal *dt) 51*d71ae5a4SJacob 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)); 583a40ed3dSBarry Smith PetscFunctionReturn(0); 597bf11e45SBarry Smith } 607bf11e45SBarry Smith 617bf11e45SBarry Smith /* ------------------------------------------------------------------------------*/ 627bf11e45SBarry Smith /*@C 638d359177SBarry Smith TSPseudoVerifyTimeStepDefault - Default code to verify the quality of the last timestep. 647bf11e45SBarry Smith 6515091d37SBarry Smith Collective on TS 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 82db781477SPatrick Sanan .seealso: `TSPseudoSetVerifyTimeStep()`, `TSPseudoVerifyTimeStep()` 837bf11e45SBarry Smith @*/ 84*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoVerifyTimeStepDefault(TS ts, Vec update, void *dtctx, PetscReal *newdt, PetscBool *flag) 85*d71ae5a4SJacob Faibussowitsch { 863a40ed3dSBarry Smith PetscFunctionBegin; 87a7cc72afSBarry Smith *flag = PETSC_TRUE; 883a40ed3dSBarry Smith PetscFunctionReturn(0); 897bf11e45SBarry Smith } 907bf11e45SBarry Smith 917bf11e45SBarry Smith /*@ 92564e8f4eSLois Curfman McInnes TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable. 937bf11e45SBarry Smith 9415091d37SBarry Smith Collective on TS 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 108564e8f4eSLois Curfman McInnes set by calling TSPseudoSetVerifyTimeStep(). 109564e8f4eSLois Curfman McInnes 110db781477SPatrick Sanan .seealso: `TSPseudoSetVerifyTimeStep()`, `TSPseudoVerifyTimeStepDefault()` 1117bf11e45SBarry Smith @*/ 112*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoVerifyTimeStep(TS ts, Vec update, PetscReal *dt, PetscBool *flag) 113*d71ae5a4SJacob 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)); 1193a40ed3dSBarry Smith PetscFunctionReturn(0); 1207bf11e45SBarry Smith } 1217bf11e45SBarry Smith 1227bf11e45SBarry Smith /* --------------------------------------------------------------------------------*/ 1237bf11e45SBarry Smith 124*d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_Pseudo(TS ts) 125*d71ae5a4SJacob 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)); 156cdbf8f93SLisandro Dalcin PetscFunctionReturn(0); 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)); 1713118ae5eSBarry Smith PetscFunctionReturn(0); 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)); 1763118ae5eSBarry Smith PetscFunctionReturn(0); 1773118ae5eSBarry Smith } 1783a40ed3dSBarry Smith PetscFunctionReturn(0); 1792d3f70b5SBarry Smith } 1802d3f70b5SBarry Smith 1812d3f70b5SBarry Smith /*------------------------------------------------------------*/ 182*d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_Pseudo(TS ts) 183*d71ae5a4SJacob 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)); 1903a40ed3dSBarry Smith PetscFunctionReturn(0); 1912d3f70b5SBarry Smith } 1922d3f70b5SBarry Smith 193*d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_Pseudo(TS ts) 194*d71ae5a4SJacob 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)); 203277b19d0SLisandro Dalcin PetscFunctionReturn(0); 204277b19d0SLisandro Dalcin } 2052d3f70b5SBarry Smith 2062d3f70b5SBarry Smith /*------------------------------------------------------------*/ 2072d3f70b5SBarry Smith 2086f2d6a7bSJed Brown /* 2096f2d6a7bSJed Brown Compute Xdot = (X^{n+1}-X^n)/dt) = 0 2106f2d6a7bSJed Brown */ 211*d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoGetXdot(TS ts, Vec X, Vec *Xdot) 212*d71ae5a4SJacob 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; 2293a40ed3dSBarry Smith PetscFunctionReturn(0); 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 */ 249*d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes, Vec X, Vec Y, TS ts) 250*d71ae5a4SJacob 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)); 2566f2d6a7bSJed Brown PetscFunctionReturn(0); 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 */ 268*d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes, Vec X, Mat AA, Mat BB, TS ts) 269*d71ae5a4SJacob 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)); 2753a40ed3dSBarry Smith PetscFunctionReturn(0); 2762d3f70b5SBarry Smith } 2772d3f70b5SBarry Smith 278*d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_Pseudo(TS ts) 279*d71ae5a4SJacob 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)); 2863a40ed3dSBarry Smith PetscFunctionReturn(0); 2872d3f70b5SBarry Smith } 2882d3f70b5SBarry Smith /*------------------------------------------------------------*/ 2892d3f70b5SBarry Smith 290*d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoMonitorDefault(TS ts, PetscInt step, PetscReal ptime, Vec v, void *dummy) 291*d71ae5a4SJacob 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)); 3043a40ed3dSBarry Smith PetscFunctionReturn(0); 3052d3f70b5SBarry Smith } 3062d3f70b5SBarry Smith 307*d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetFromOptions_Pseudo(TS ts, PetscOptionItems *PetscOptionsObject) 308*d71ae5a4SJacob 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(); 3283a40ed3dSBarry Smith PetscFunctionReturn(0); 3292d3f70b5SBarry Smith } 3302d3f70b5SBarry Smith 331*d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_Pseudo(TS ts, PetscViewer viewer) 332*d71ae5a4SJacob 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 } 3453a40ed3dSBarry Smith PetscFunctionReturn(0); 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 3533f9fe445SBarry Smith Logically Collective on TS 35415091d37SBarry Smith 35582bf6240SBarry Smith Input Parameters: 35615091d37SBarry Smith + ts - timestep context 35782bf6240SBarry Smith . dt - user-defined function to verify timestep 35815091d37SBarry Smith - ctx - [optional] user-defined context for private data 3590298fd71SBarry Smith for the timestep verification routine (may be NULL) 36082bf6240SBarry Smith 36115091d37SBarry Smith Level: advanced 362fee21e36SBarry Smith 36382bf6240SBarry Smith Calling sequence of func: 364a2b725a8SWilliam Gropp $ func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscBool *flag); 36582bf6240SBarry Smith 366a2b725a8SWilliam Gropp + update - latest solution vector 36782bf6240SBarry Smith . ctx - [optional] timestep context 36882bf6240SBarry Smith . newdt - the timestep to use for the next step 369a2b725a8SWilliam Gropp - flag - flag indicating whether the last time step was acceptable 37082bf6240SBarry Smith 37182bf6240SBarry Smith Notes: 37282bf6240SBarry Smith The routine set here will be called by TSPseudoVerifyTimeStep() 37382bf6240SBarry Smith during the timestepping process. 37482bf6240SBarry Smith 375db781477SPatrick Sanan .seealso: `TSPseudoVerifyTimeStepDefault()`, `TSPseudoVerifyTimeStep()` 37682bf6240SBarry Smith @*/ 377*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoSetVerifyTimeStep(TS ts, PetscErrorCode (*dt)(TS, Vec, void *, PetscReal *, PetscBool *), void *ctx) 378*d71ae5a4SJacob Faibussowitsch { 37982bf6240SBarry Smith PetscFunctionBegin; 3800700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 381cac4c232SBarry Smith PetscTryMethod(ts, "TSPseudoSetVerifyTimeStep_C", (TS, PetscErrorCode(*)(TS, Vec, void *, PetscReal *, PetscBool *), void *), (ts, dt, ctx)); 38282bf6240SBarry Smith PetscFunctionReturn(0); 38382bf6240SBarry Smith } 38482bf6240SBarry Smith 38582bf6240SBarry Smith /*@ 38682bf6240SBarry Smith TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 3878d359177SBarry Smith dt when using the TSPseudoTimeStepDefault() routine. 38882bf6240SBarry Smith 3893f9fe445SBarry Smith Logically Collective on TS 390fee21e36SBarry Smith 39115091d37SBarry Smith Input Parameters: 39215091d37SBarry Smith + ts - the timestep context 39315091d37SBarry Smith - inc - the scaling factor >= 1.0 39415091d37SBarry Smith 39582bf6240SBarry Smith Options Database Key: 39667b8a455SSatish Balay . -ts_pseudo_increment <increment> - set pseudo increment 39782bf6240SBarry Smith 39815091d37SBarry Smith Level: advanced 39915091d37SBarry Smith 400db781477SPatrick Sanan .seealso: `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()` 40182bf6240SBarry Smith @*/ 402*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoSetTimeStepIncrement(TS ts, PetscReal inc) 403*d71ae5a4SJacob Faibussowitsch { 40482bf6240SBarry Smith PetscFunctionBegin; 4050700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 406c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts, inc, 2); 407cac4c232SBarry Smith PetscTryMethod(ts, "TSPseudoSetTimeStepIncrement_C", (TS, PetscReal), (ts, inc)); 40882bf6240SBarry Smith PetscFunctionReturn(0); 40982bf6240SBarry Smith } 41082bf6240SBarry Smith 41186534af1SJed Brown /*@ 41286534af1SJed Brown TSPseudoSetMaxTimeStep - Sets the maximum time step 4138d359177SBarry Smith when using the TSPseudoTimeStepDefault() routine. 41486534af1SJed Brown 41586534af1SJed Brown Logically Collective on TS 41686534af1SJed Brown 41786534af1SJed Brown Input Parameters: 41886534af1SJed Brown + ts - the timestep context 41986534af1SJed Brown - maxdt - the maximum time step, use a non-positive value to deactivate 42086534af1SJed Brown 42186534af1SJed Brown Options Database Key: 42267b8a455SSatish Balay . -ts_pseudo_max_dt <increment> - set pseudo max dt 42386534af1SJed Brown 42486534af1SJed Brown Level: advanced 42586534af1SJed Brown 426db781477SPatrick Sanan .seealso: `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()` 42786534af1SJed Brown @*/ 428*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoSetMaxTimeStep(TS ts, PetscReal maxdt) 429*d71ae5a4SJacob Faibussowitsch { 43086534af1SJed Brown PetscFunctionBegin; 43186534af1SJed Brown PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 43286534af1SJed Brown PetscValidLogicalCollectiveReal(ts, maxdt, 2); 433cac4c232SBarry Smith PetscTryMethod(ts, "TSPseudoSetMaxTimeStep_C", (TS, PetscReal), (ts, maxdt)); 43486534af1SJed Brown PetscFunctionReturn(0); 43586534af1SJed Brown } 43686534af1SJed Brown 43782bf6240SBarry Smith /*@ 43882bf6240SBarry Smith TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep 43982bf6240SBarry Smith is computed via the formula 44082bf6240SBarry Smith $ dt = initial_dt*initial_fnorm/current_fnorm 44182bf6240SBarry Smith rather than the default update, 44282bf6240SBarry Smith $ dt = current_dt*previous_fnorm/current_fnorm. 44382bf6240SBarry Smith 4443f9fe445SBarry Smith Logically Collective on TS 44515091d37SBarry Smith 44682bf6240SBarry Smith Input Parameter: 44782bf6240SBarry Smith . ts - the timestep context 44882bf6240SBarry Smith 44982bf6240SBarry Smith Options Database Key: 45067b8a455SSatish Balay . -ts_pseudo_increment_dt_from_initial_dt <true,false> - use the initial dt to determine increment 45182bf6240SBarry Smith 45215091d37SBarry Smith Level: advanced 45315091d37SBarry Smith 454db781477SPatrick Sanan .seealso: `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()` 45582bf6240SBarry Smith @*/ 456*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS ts) 457*d71ae5a4SJacob Faibussowitsch { 45882bf6240SBarry Smith PetscFunctionBegin; 4590700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 460cac4c232SBarry Smith PetscTryMethod(ts, "TSPseudoIncrementDtFromInitialDt_C", (TS), (ts)); 46182bf6240SBarry Smith PetscFunctionReturn(0); 46282bf6240SBarry Smith } 46382bf6240SBarry Smith 464ac226902SBarry Smith /*@C 46582bf6240SBarry Smith TSPseudoSetTimeStep - Sets the user-defined routine to be 46682bf6240SBarry Smith called at each pseudo-timestep to update the timestep. 46782bf6240SBarry Smith 4683f9fe445SBarry Smith Logically Collective on TS 46915091d37SBarry Smith 47082bf6240SBarry Smith Input Parameters: 47115091d37SBarry Smith + ts - timestep context 47282bf6240SBarry Smith . dt - function to compute timestep 47315091d37SBarry Smith - ctx - [optional] user-defined context for private data 4740298fd71SBarry Smith required by the function (may be NULL) 47582bf6240SBarry Smith 47615091d37SBarry Smith Level: intermediate 477fee21e36SBarry Smith 47882bf6240SBarry Smith Calling sequence of func: 479a2b725a8SWilliam Gropp $ func (TS ts,PetscReal *newdt,void *ctx); 48082bf6240SBarry Smith 481a2b725a8SWilliam Gropp + newdt - the newly computed timestep 482a2b725a8SWilliam Gropp - ctx - [optional] timestep context 48382bf6240SBarry Smith 48482bf6240SBarry Smith Notes: 48582bf6240SBarry Smith The routine set here will be called by TSPseudoComputeTimeStep() 48682bf6240SBarry Smith during the timestepping process. 4878d359177SBarry Smith If not set then TSPseudoTimeStepDefault() is automatically used 48882bf6240SBarry Smith 489db781477SPatrick Sanan .seealso: `TSPseudoTimeStepDefault()`, `TSPseudoComputeTimeStep()` 49082bf6240SBarry Smith @*/ 491*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoSetTimeStep(TS ts, PetscErrorCode (*dt)(TS, PetscReal *, void *), void *ctx) 492*d71ae5a4SJacob Faibussowitsch { 49382bf6240SBarry Smith PetscFunctionBegin; 4940700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 495cac4c232SBarry Smith PetscTryMethod(ts, "TSPseudoSetTimeStep_C", (TS, PetscErrorCode(*)(TS, PetscReal *, void *), void *), (ts, dt, ctx)); 49682bf6240SBarry Smith PetscFunctionReturn(0); 49782bf6240SBarry Smith } 49882bf6240SBarry Smith 49982bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 50082bf6240SBarry Smith 501ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(TS, Vec, void *, PetscReal *, PetscBool *); /* force argument to next function to not be extern C*/ 502*d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetVerifyTimeStep_Pseudo(TS ts, FCN1 dt, void *ctx) 503*d71ae5a4SJacob Faibussowitsch { 504be5899b3SLisandro Dalcin TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 50582bf6240SBarry Smith 50682bf6240SBarry Smith PetscFunctionBegin; 50782bf6240SBarry Smith pseudo->verify = dt; 50882bf6240SBarry Smith pseudo->verifyctx = ctx; 50982bf6240SBarry Smith PetscFunctionReturn(0); 51082bf6240SBarry Smith } 51182bf6240SBarry Smith 512*d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetTimeStepIncrement_Pseudo(TS ts, PetscReal inc) 513*d71ae5a4SJacob Faibussowitsch { 5144bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 51582bf6240SBarry Smith 51682bf6240SBarry Smith PetscFunctionBegin; 51782bf6240SBarry Smith pseudo->dt_increment = inc; 51882bf6240SBarry Smith PetscFunctionReturn(0); 51982bf6240SBarry Smith } 52082bf6240SBarry Smith 521*d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetMaxTimeStep_Pseudo(TS ts, PetscReal maxdt) 522*d71ae5a4SJacob Faibussowitsch { 52386534af1SJed Brown TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 52486534af1SJed Brown 52586534af1SJed Brown PetscFunctionBegin; 52686534af1SJed Brown pseudo->dt_max = maxdt; 52786534af1SJed Brown PetscFunctionReturn(0); 52886534af1SJed Brown } 52986534af1SJed Brown 530*d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts) 531*d71ae5a4SJacob Faibussowitsch { 5324bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 53382bf6240SBarry Smith 53482bf6240SBarry Smith PetscFunctionBegin; 5354bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_TRUE; 53682bf6240SBarry Smith PetscFunctionReturn(0); 53782bf6240SBarry Smith } 53882bf6240SBarry Smith 5396849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS, PetscReal *, void *); /* force argument to next function to not be extern C*/ 540*d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetTimeStep_Pseudo(TS ts, FCN2 dt, void *ctx) 541*d71ae5a4SJacob Faibussowitsch { 5424bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 54382bf6240SBarry Smith 54482bf6240SBarry Smith PetscFunctionBegin; 54582bf6240SBarry Smith pseudo->dt = dt; 54682bf6240SBarry Smith pseudo->dtctx = ctx; 54782bf6240SBarry Smith PetscFunctionReturn(0); 54882bf6240SBarry Smith } 54982bf6240SBarry Smith 55082bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 55110e6a065SJed Brown /*MC 55210e6a065SJed Brown TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping 55382bf6240SBarry Smith 55410e6a065SJed Brown This method solves equations of the form 55510e6a065SJed Brown 55610e6a065SJed Brown $ F(X,Xdot) = 0 55710e6a065SJed Brown 55810e6a065SJed Brown for steady state using the iteration 55910e6a065SJed Brown 56010e6a065SJed Brown $ [G'] S = -F(X,0) 56110e6a065SJed Brown $ X += S 56210e6a065SJed Brown 56310e6a065SJed Brown where 56410e6a065SJed Brown 56510e6a065SJed Brown $ G(Y) = F(Y,(Y-X)/dt) 56610e6a065SJed Brown 5676f2d6a7bSJed Brown This is linearly-implicit Euler with the residual always evaluated "at steady 5686f2d6a7bSJed Brown state". See note below. 56910e6a065SJed Brown 57010e6a065SJed Brown Options database keys: 57110e6a065SJed Brown + -ts_pseudo_increment <real> - ratio of increase dt 5723118ae5eSBarry Smith . -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt 5733118ae5eSBarry Smith . -ts_pseudo_fatol <atol> - stop iterating when the function norm is less than atol 5743118ae5eSBarry Smith - -ts_pseudo_frtol <rtol> - stop iterating when the function norm divided by the initial function norm is less than rtol 57510e6a065SJed Brown 57610e6a065SJed Brown Level: beginner 57710e6a065SJed Brown 57810e6a065SJed Brown References: 579606c0280SSatish Balay + * - Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential Algebraic Equations, 2003. 580606c0280SSatish Balay - * - C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998. 58110e6a065SJed Brown 58210e6a065SJed Brown Notes: 5836f2d6a7bSJed Brown The residual computed by this method includes the transient term (Xdot is computed instead of 5846f2d6a7bSJed Brown always being zero), but since the prediction from the last step is always the solution from the 5856f2d6a7bSJed Brown last step, on the first Newton iteration we have 5866f2d6a7bSJed Brown 5876f2d6a7bSJed Brown $ Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0 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 593db781477SPatrick Sanan .seealso: `TSCreate()`, `TS`, `TSSetType()` 59410e6a065SJed Brown 59510e6a065SJed Brown M*/ 596*d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_Pseudo(TS ts) 597*d71ae5a4SJacob Faibussowitsch { 5987bf11e45SBarry Smith TS_Pseudo *pseudo; 599193ac0bcSJed Brown SNES snes; 60019fd82e9SBarry Smith SNESType stype; 6012d3f70b5SBarry Smith 6023a40ed3dSBarry Smith PetscFunctionBegin; 603277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Pseudo; 604000e7ae3SMatthew Knepley ts->ops->destroy = TSDestroy_Pseudo; 605000e7ae3SMatthew Knepley ts->ops->view = TSView_Pseudo; 606000e7ae3SMatthew Knepley ts->ops->setup = TSSetUp_Pseudo; 607000e7ae3SMatthew Knepley ts->ops->step = TSStep_Pseudo; 608000e7ae3SMatthew Knepley ts->ops->setfromoptions = TSSetFromOptions_Pseudo; 6090f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Pseudo; 6100f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Pseudo; 6112ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 612825ab935SBarry Smith ts->usessnes = PETSC_TRUE; 6137bf11e45SBarry Smith 6149566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 6159566063dSJacob Faibussowitsch PetscCall(SNESGetType(snes, &stype)); 6169566063dSJacob Faibussowitsch if (!stype) PetscCall(SNESSetType(snes, SNESKSPONLY)); 6172d3f70b5SBarry Smith 6184dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&pseudo)); 6197bf11e45SBarry Smith ts->data = (void *)pseudo; 6202d3f70b5SBarry Smith 621be5899b3SLisandro Dalcin pseudo->dt = TSPseudoTimeStepDefault; 622be5899b3SLisandro Dalcin pseudo->dtctx = NULL; 62328aa8177SBarry Smith pseudo->dt_increment = 1.1; 6244bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_FALSE; 625193ac0bcSJed Brown pseudo->fnorm = -1; 626be5899b3SLisandro Dalcin pseudo->fnorm_initial = -1; 627be5899b3SLisandro Dalcin pseudo->fnorm_previous = -1; 6283118ae5eSBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 6293118ae5eSBarry Smith pseudo->fatol = 1.e-25; 6303118ae5eSBarry Smith pseudo->frtol = 1.e-5; 6313118ae5eSBarry Smith #else 6323118ae5eSBarry Smith pseudo->fatol = 1.e-50; 6333118ae5eSBarry Smith pseudo->frtol = 1.e-12; 6343118ae5eSBarry Smith #endif 6359566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetVerifyTimeStep_C", TSPseudoSetVerifyTimeStep_Pseudo)); 6369566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStepIncrement_C", TSPseudoSetTimeStepIncrement_Pseudo)); 6379566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetMaxTimeStep_C", TSPseudoSetMaxTimeStep_Pseudo)); 6389566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoIncrementDtFromInitialDt_C", TSPseudoIncrementDtFromInitialDt_Pseudo)); 6399566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStep_C", TSPseudoSetTimeStep_Pseudo)); 6403a40ed3dSBarry Smith PetscFunctionReturn(0); 6412d3f70b5SBarry Smith } 6422d3f70b5SBarry Smith 64382bf6240SBarry Smith /*@C 6448d359177SBarry Smith TSPseudoTimeStepDefault - Default code to compute pseudo-timestepping. 64582bf6240SBarry Smith Use with TSPseudoSetTimeStep(). 64628aa8177SBarry Smith 64715091d37SBarry Smith Collective on TS 64815091d37SBarry Smith 64928aa8177SBarry Smith Input Parameters: 650a2b725a8SWilliam Gropp + ts - the timestep context 651a2b725a8SWilliam Gropp - dtctx - unused timestep context 65228aa8177SBarry Smith 65382bf6240SBarry Smith Output Parameter: 65482bf6240SBarry Smith . newdt - the timestep to use for the next step 65528aa8177SBarry Smith 65615091d37SBarry Smith Level: advanced 65715091d37SBarry Smith 658db781477SPatrick Sanan .seealso: `TSPseudoSetTimeStep()`, `TSPseudoComputeTimeStep()` 65928aa8177SBarry Smith @*/ 660*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoTimeStepDefault(TS ts, PetscReal *newdt, void *dtctx) 661*d71ae5a4SJacob Faibussowitsch { 66282bf6240SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 663be5899b3SLisandro Dalcin PetscReal inc = pseudo->dt_increment; 66428aa8177SBarry Smith 6653a40ed3dSBarry Smith PetscFunctionBegin; 6669566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(pseudo->xdot)); 6679566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts, ts->ptime, ts->vec_sol, pseudo->xdot, pseudo->func, PETSC_FALSE)); 6689566063dSJacob Faibussowitsch PetscCall(VecNorm(pseudo->func, NORM_2, &pseudo->fnorm)); 669be5899b3SLisandro Dalcin if (pseudo->fnorm_initial < 0) { 67082bf6240SBarry Smith /* first time through so compute initial function norm */ 671cdbf8f93SLisandro Dalcin pseudo->fnorm_initial = pseudo->fnorm; 672be5899b3SLisandro Dalcin pseudo->fnorm_previous = pseudo->fnorm; 67382bf6240SBarry Smith } 674bbd56ea5SKarl Rupp if (pseudo->fnorm == 0.0) *newdt = 1.e12 * inc * ts->time_step; 675bbd56ea5SKarl Rupp else if (pseudo->increment_dt_from_initial_dt) *newdt = inc * pseudo->dt_initial * pseudo->fnorm_initial / pseudo->fnorm; 676be5899b3SLisandro Dalcin else *newdt = inc * ts->time_step * pseudo->fnorm_previous / pseudo->fnorm; 67786534af1SJed Brown if (pseudo->dt_max > 0) *newdt = PetscMin(*newdt, pseudo->dt_max); 67882bf6240SBarry Smith pseudo->fnorm_previous = pseudo->fnorm; 6793a40ed3dSBarry Smith PetscFunctionReturn(0); 68028aa8177SBarry Smith } 681