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; 26eb636060SBarry Smith 27eb636060SBarry Smith PetscObjectState Xstate; /* state of input vector to TSComputeIFunction() with 0 Xdot */ 287bf11e45SBarry Smith } TS_Pseudo; 292d3f70b5SBarry Smith 302d3f70b5SBarry Smith /* ------------------------------------------------------------------------------*/ 312d3f70b5SBarry Smith 32cc4c1da9SBarry Smith /*@ 337bf11e45SBarry Smith TSPseudoComputeTimeStep - Computes the next timestep for a currently running 34564e8f4eSLois Curfman McInnes pseudo-timestepping process. 352d3f70b5SBarry Smith 36c3339decSBarry Smith Collective 3715091d37SBarry Smith 387bf11e45SBarry Smith Input Parameter: 397bf11e45SBarry Smith . ts - timestep context 407bf11e45SBarry Smith 417bf11e45SBarry Smith Output Parameter: 42fb4a63b6SLois Curfman McInnes . dt - newly computed timestep 43fb4a63b6SLois Curfman McInnes 448d359177SBarry Smith Level: developer 45564e8f4eSLois Curfman McInnes 46bcf0153eSBarry Smith Note: 47564e8f4eSLois Curfman McInnes The routine to be called here to compute the timestep should be 48bcf0153eSBarry Smith set by calling `TSPseudoSetTimeStep()`. 49564e8f4eSLois Curfman McInnes 501cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoTimeStepDefault()`, `TSPseudoSetTimeStep()` 517bf11e45SBarry Smith @*/ 52d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoComputeTimeStep(TS ts, PetscReal *dt) 53d71ae5a4SJacob Faibussowitsch { 547bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 557bf11e45SBarry Smith 563a40ed3dSBarry Smith PetscFunctionBegin; 579566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TS_PseudoComputeTimeStep, ts, 0, 0, 0)); 589566063dSJacob Faibussowitsch PetscCall((*pseudo->dt)(ts, dt, pseudo->dtctx)); 599566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TS_PseudoComputeTimeStep, ts, 0, 0, 0)); 603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 617bf11e45SBarry Smith } 627bf11e45SBarry Smith 637bf11e45SBarry Smith /* ------------------------------------------------------------------------------*/ 647bf11e45SBarry Smith /*@C 658d359177SBarry Smith TSPseudoVerifyTimeStepDefault - Default code to verify the quality of the last timestep. 667bf11e45SBarry Smith 67cc4c1da9SBarry Smith Collective, No Fortran Support 6815091d37SBarry Smith 697bf11e45SBarry Smith Input Parameters: 7015091d37SBarry Smith + ts - the timestep context 717bf11e45SBarry Smith . dtctx - unused timestep context 7215091d37SBarry Smith - update - latest solution vector 737bf11e45SBarry Smith 74564e8f4eSLois Curfman McInnes Output Parameters: 7515091d37SBarry Smith + newdt - the timestep to use for the next step 7615091d37SBarry Smith - flag - flag indicating whether the last time step was acceptable 777bf11e45SBarry Smith 7815091d37SBarry Smith Level: advanced 79fee21e36SBarry Smith 80564e8f4eSLois Curfman McInnes Note: 81564e8f4eSLois Curfman McInnes This routine always returns a flag of 1, indicating an acceptable 82564e8f4eSLois Curfman McInnes timestep. 83564e8f4eSLois Curfman McInnes 841cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetVerifyTimeStep()`, `TSPseudoVerifyTimeStep()` 857bf11e45SBarry Smith @*/ 86d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoVerifyTimeStepDefault(TS ts, Vec update, void *dtctx, PetscReal *newdt, PetscBool *flag) 87d71ae5a4SJacob Faibussowitsch { 883a40ed3dSBarry Smith PetscFunctionBegin; 89a7cc72afSBarry Smith *flag = PETSC_TRUE; 903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 917bf11e45SBarry Smith } 927bf11e45SBarry Smith 937bf11e45SBarry Smith /*@ 94564e8f4eSLois Curfman McInnes TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable. 957bf11e45SBarry Smith 96c3339decSBarry Smith Collective 9715091d37SBarry Smith 98fb4a63b6SLois Curfman McInnes Input Parameters: 9915091d37SBarry Smith + ts - timestep context 10015091d37SBarry Smith - update - latest solution vector 1017bf11e45SBarry Smith 102fb4a63b6SLois Curfman McInnes Output Parameters: 10315091d37SBarry Smith + dt - newly computed timestep (if it had to shrink) 10415091d37SBarry Smith - flag - indicates if current timestep was ok 1057bf11e45SBarry Smith 10615091d37SBarry Smith Level: advanced 107fee21e36SBarry Smith 108564e8f4eSLois Curfman McInnes Notes: 109564e8f4eSLois Curfman McInnes The routine to be called here to compute the timestep should be 110bcf0153eSBarry Smith set by calling `TSPseudoSetVerifyTimeStep()`. 111564e8f4eSLois Curfman McInnes 1121cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetVerifyTimeStep()`, `TSPseudoVerifyTimeStepDefault()` 1137bf11e45SBarry Smith @*/ 114d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoVerifyTimeStep(TS ts, Vec update, PetscReal *dt, PetscBool *flag) 115d71ae5a4SJacob Faibussowitsch { 1167bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 1177bf11e45SBarry Smith 1183a40ed3dSBarry Smith PetscFunctionBegin; 119cb9d8021SPierre Barbier de Reuille *flag = PETSC_TRUE; 1201baa6e33SBarry Smith if (pseudo->verify) PetscCall((*pseudo->verify)(ts, update, pseudo->verifyctx, dt, flag)); 1213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1227bf11e45SBarry Smith } 1237bf11e45SBarry Smith 1247bf11e45SBarry Smith /* --------------------------------------------------------------------------------*/ 1257bf11e45SBarry Smith 126d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_Pseudo(TS ts) 127d71ae5a4SJacob Faibussowitsch { 128277b19d0SLisandro Dalcin TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 129be5899b3SLisandro Dalcin PetscInt nits, lits, reject; 130cdbf8f93SLisandro Dalcin PetscBool stepok; 131be5899b3SLisandro Dalcin PetscReal next_time_step = ts->time_step; 1322d3f70b5SBarry Smith 1333a40ed3dSBarry Smith PetscFunctionBegin; 134eb636060SBarry Smith if (ts->steps == 0) { 135eb636060SBarry Smith pseudo->dt_initial = ts->time_step; 136eb636060SBarry Smith PetscCall(VecCopy(ts->vec_sol, pseudo->update)); /* in all future updates pseudo->update already contains the current time solution */ 137eb636060SBarry Smith } 1389566063dSJacob Faibussowitsch PetscCall(TSPseudoComputeTimeStep(ts, &next_time_step)); 139cdbf8f93SLisandro Dalcin for (reject = 0; reject < ts->max_reject; reject++, ts->reject++) { 140cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 141*81775c15SStefano Zampini if (reject) PetscCall(VecCopy(ts->vec_sol, pseudo->update)); /* need to copy the solution because we got a rejection and pseudo->update contains intermediate computation */ 1429566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, ts->ptime + ts->time_step)); 1439566063dSJacob Faibussowitsch PetscCall(SNESSolve(ts->snes, NULL, pseudo->update)); 1449566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(ts->snes, &nits)); 1459566063dSJacob Faibussowitsch PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits)); 1469371c9d4SSatish Balay ts->snes_its += nits; 1479371c9d4SSatish Balay ts->ksp_its += lits; 148f4f49eeaSPierre Jolivet PetscCall(TSPostStage(ts, ts->ptime + ts->time_step, 0, &pseudo->update)); 1499566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(ts->adapt, ts, ts->ptime + ts->time_step, pseudo->update, &stepok)); 1509371c9d4SSatish Balay if (!stepok) { 1519371c9d4SSatish Balay next_time_step = ts->time_step; 1529371c9d4SSatish Balay continue; 1539371c9d4SSatish Balay } 154eb636060SBarry Smith pseudo->fnorm = -1; /* The current norm is no longer valid */ 1559566063dSJacob Faibussowitsch PetscCall(TSPseudoVerifyTimeStep(ts, pseudo->update, &next_time_step, &stepok)); 156cdbf8f93SLisandro Dalcin if (stepok) break; 157cdbf8f93SLisandro Dalcin } 158be5899b3SLisandro Dalcin if (reject >= ts->max_reject) { 159be5899b3SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 16063a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, reject)); 1613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1627bf11e45SBarry Smith } 163be5899b3SLisandro Dalcin 1649566063dSJacob Faibussowitsch PetscCall(VecCopy(pseudo->update, ts->vec_sol)); 165be5899b3SLisandro Dalcin ts->ptime += ts->time_step; 166be5899b3SLisandro Dalcin ts->time_step = next_time_step; 167be5899b3SLisandro Dalcin 1689566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(pseudo->xdot)); 169eb636060SBarry Smith PetscCall(TSComputeIFunction(ts, ts->ptime, pseudo->update, pseudo->xdot, pseudo->func, PETSC_FALSE)); 170eb636060SBarry Smith PetscCall(PetscObjectStateGet((PetscObject)pseudo->update, &pseudo->Xstate)); 1719566063dSJacob Faibussowitsch PetscCall(VecNorm(pseudo->func, NORM_2, &pseudo->fnorm)); 172eb636060SBarry Smith 1733118ae5eSBarry Smith if (pseudo->fnorm < pseudo->fatol) { 1743118ae5eSBarry Smith ts->reason = TS_CONVERGED_PSEUDO_FATOL; 17563a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", converged since fnorm %g < fatol %g\n", ts->steps, (double)pseudo->fnorm, (double)pseudo->frtol)); 1763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1773118ae5eSBarry Smith } 1783118ae5eSBarry Smith if (pseudo->fnorm / pseudo->fnorm_initial < pseudo->frtol) { 1793118ae5eSBarry Smith ts->reason = TS_CONVERGED_PSEUDO_FRTOL; 18063a3b9bcSJacob 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)); 1813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1823118ae5eSBarry Smith } 1833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1842d3f70b5SBarry Smith } 1852d3f70b5SBarry Smith 1862d3f70b5SBarry Smith /*------------------------------------------------------------*/ 187*81775c15SStefano Zampini 188d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_Pseudo(TS ts) 189d71ae5a4SJacob Faibussowitsch { 1907bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 1912d3f70b5SBarry Smith 1923a40ed3dSBarry Smith PetscFunctionBegin; 1939566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pseudo->update)); 1949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pseudo->func)); 1959566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pseudo->xdot)); 1963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1972d3f70b5SBarry Smith } 1982d3f70b5SBarry Smith 199d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_Pseudo(TS ts) 200d71ae5a4SJacob Faibussowitsch { 201277b19d0SLisandro Dalcin PetscFunctionBegin; 2029566063dSJacob Faibussowitsch PetscCall(TSReset_Pseudo(ts)); 2039566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 2049566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetVerifyTimeStep_C", NULL)); 2059566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStepIncrement_C", NULL)); 2069566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetMaxTimeStep_C", NULL)); 2079566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoIncrementDtFromInitialDt_C", NULL)); 2089566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStep_C", NULL)); 2093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 210277b19d0SLisandro Dalcin } 2112d3f70b5SBarry Smith 2122d3f70b5SBarry Smith /*------------------------------------------------------------*/ 2132d3f70b5SBarry Smith 214d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoGetXdot(TS ts, Vec X, Vec *Xdot) 215d71ae5a4SJacob Faibussowitsch { 2166f2d6a7bSJed Brown TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 2172d3f70b5SBarry Smith 2183a40ed3dSBarry Smith PetscFunctionBegin; 2196f2d6a7bSJed Brown *Xdot = pseudo->xdot; 2203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2212d3f70b5SBarry Smith } 2222d3f70b5SBarry Smith 2236f2d6a7bSJed Brown /* 2246f2d6a7bSJed Brown The transient residual is 2256f2d6a7bSJed Brown 2266f2d6a7bSJed Brown F(U^{n+1},(U^{n+1}-U^n)/dt) = 0 2276f2d6a7bSJed Brown 2286f2d6a7bSJed Brown or for ODE, 2296f2d6a7bSJed Brown 2306f2d6a7bSJed Brown (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0 2316f2d6a7bSJed Brown 2326f2d6a7bSJed Brown This is the function that must be evaluated for transient simulation and for 2336f2d6a7bSJed Brown finite difference Jacobians. On the first Newton step, this algorithm uses 2346f2d6a7bSJed Brown a guess of U^{n+1} = U^n in which case the transient term vanishes and the 2356f2d6a7bSJed Brown residual is actually the steady state residual. Pseudotransient 2366f2d6a7bSJed Brown continuation as described in the literature is a linearly implicit 237eb636060SBarry Smith algorithm, it only takes this one full Newton step with the steady state 2386f2d6a7bSJed Brown residual, and then advances to the next time step. 239eb636060SBarry Smith 240eb636060SBarry Smith This routine avoids a repeated identical call to TSComputeRHSFunction() when that result 241eb636060SBarry Smith is already contained in pseudo->func due to the call to TSComputeIFunction() in TSStep_Pseudo() 242eb636060SBarry Smith 2436f2d6a7bSJed Brown */ 244d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes, Vec X, Vec Y, TS ts) 245d71ae5a4SJacob Faibussowitsch { 246eb636060SBarry Smith TSIFunctionFn *ifunction; 247eb636060SBarry Smith TSRHSFunctionFn *rhsfunction; 248eb636060SBarry Smith void *ctx; 249eb636060SBarry Smith DM dm; 250eb636060SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 251*81775c15SStefano Zampini const PetscScalar mdt = 1.0 / ts->time_step; 252eb636060SBarry Smith PetscBool KSPSNES; 253eb636060SBarry Smith PetscObjectState Xstate; 2542d3f70b5SBarry Smith 2553a40ed3dSBarry Smith PetscFunctionBegin; 256eb636060SBarry Smith PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 257eb636060SBarry Smith PetscValidHeaderSpecific(X, VEC_CLASSID, 2); 258eb636060SBarry Smith PetscValidHeaderSpecific(Y, VEC_CLASSID, 3); 259eb636060SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID, 4); 260eb636060SBarry Smith 261eb636060SBarry Smith PetscCall(TSGetDM(ts, &dm)); 262eb636060SBarry Smith PetscCall(DMTSGetIFunction(dm, &ifunction, &ctx)); 263eb636060SBarry Smith PetscCall(DMTSGetRHSFunction(dm, &rhsfunction, NULL)); 264eb636060SBarry Smith PetscCheck(rhsfunction || ifunction, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Must call TSSetRHSFunction() and / or TSSetIFunction()"); 265eb636060SBarry Smith 266eb636060SBarry Smith PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESKSPONLY, &KSPSNES)); 267eb636060SBarry Smith PetscCall(PetscObjectStateGet((PetscObject)X, &Xstate)); 268eb636060SBarry Smith if (Xstate != pseudo->Xstate || ifunction || !KSPSNES) { 269*81775c15SStefano Zampini PetscCall(VecAXPBYPCZ(pseudo->xdot, -mdt, mdt, 0, ts->vec_sol, X)); 270*81775c15SStefano Zampini PetscCall(TSComputeIFunction(ts, ts->ptime + ts->time_step, X, pseudo->xdot, Y, PETSC_FALSE)); 271eb636060SBarry Smith } else { 272eb636060SBarry Smith /* reuse the TSComputeIFunction() result performed inside TSStep_Pseudo() */ 273eb636060SBarry Smith /* note that pseudo->func contains the negation of TSComputeRHSFunction() */ 274*81775c15SStefano Zampini PetscCall(VecWAXPY(Y, 1, pseudo->func, pseudo->xdot)); 275eb636060SBarry Smith } 2763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2776f2d6a7bSJed Brown } 2782d3f70b5SBarry Smith 2796f2d6a7bSJed Brown /* 2806f2d6a7bSJed Brown This constructs the Jacobian needed for SNES. For DAE, this is 2816f2d6a7bSJed Brown 2826f2d6a7bSJed Brown dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot 2836f2d6a7bSJed Brown 2846f2d6a7bSJed Brown and for ODE: 2856f2d6a7bSJed Brown 2866f2d6a7bSJed Brown J = I/dt - J_{Frhs} where J_{Frhs} is the given Jacobian of Frhs. 2876f2d6a7bSJed Brown */ 288d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes, Vec X, Mat AA, Mat BB, TS ts) 289d71ae5a4SJacob Faibussowitsch { 2906f2d6a7bSJed Brown Vec Xdot; 2916f2d6a7bSJed Brown 2926f2d6a7bSJed Brown PetscFunctionBegin; 293*81775c15SStefano Zampini /* Xdot has already been computed in SNESTSFormFunction_Pseudo (SNES guarantees this) */ 2949566063dSJacob Faibussowitsch PetscCall(TSPseudoGetXdot(ts, X, &Xdot)); 2959566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts, ts->ptime + ts->time_step, X, Xdot, 1. / ts->time_step, AA, BB, PETSC_FALSE)); 2963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2972d3f70b5SBarry Smith } 2982d3f70b5SBarry Smith 299d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_Pseudo(TS ts) 300d71ae5a4SJacob Faibussowitsch { 3017bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 3022d3f70b5SBarry Smith 3033a40ed3dSBarry Smith PetscFunctionBegin; 3049566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &pseudo->update)); 3059566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &pseudo->func)); 3069566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &pseudo->xdot)); 3073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3082d3f70b5SBarry Smith } 309*81775c15SStefano Zampini 3102d3f70b5SBarry Smith /*------------------------------------------------------------*/ 3112d3f70b5SBarry Smith 312d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoMonitorDefault(TS ts, PetscInt step, PetscReal ptime, Vec v, void *dummy) 313d71ae5a4SJacob Faibussowitsch { 3147bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 315ce94432eSBarry Smith PetscViewer viewer = (PetscViewer)dummy; 3162d3f70b5SBarry Smith 3173a40ed3dSBarry Smith PetscFunctionBegin; 318193ac0bcSJed Brown if (pseudo->fnorm < 0) { /* The last computed norm is stale, recompute */ 3199566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(pseudo->xdot)); 3209566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts, ts->ptime, ts->vec_sol, pseudo->xdot, pseudo->func, PETSC_FALSE)); 321eb636060SBarry Smith PetscCall(PetscObjectStateGet((PetscObject)ts->vec_sol, &pseudo->Xstate)); 3229566063dSJacob Faibussowitsch PetscCall(VecNorm(pseudo->func, NORM_2, &pseudo->fnorm)); 323193ac0bcSJed Brown } 3249566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel)); 32563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "TS %" PetscInt_FMT " dt %g time %g fnorm %g\n", step, (double)ts->time_step, (double)ptime, (double)pseudo->fnorm)); 3269566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel)); 3273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3282d3f70b5SBarry Smith } 3292d3f70b5SBarry Smith 330ce78bad3SBarry Smith static PetscErrorCode TSSetFromOptions_Pseudo(TS ts, PetscOptionItems PetscOptionsObject) 331d71ae5a4SJacob Faibussowitsch { 3324bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 333ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE; 334649052a6SBarry Smith PetscViewer viewer; 3352d3f70b5SBarry Smith 3363a40ed3dSBarry Smith PetscFunctionBegin; 337d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Pseudo-timestepping options"); 3389566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_monitor_pseudo", "Monitor convergence", "", flg, &flg, NULL)); 3392d3f70b5SBarry Smith if (flg) { 3409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ts), "stdout", &viewer)); 34149abdd8aSBarry Smith PetscCall(TSMonitorSet(ts, TSPseudoMonitorDefault, viewer, (PetscCtxDestroyFn *)PetscViewerDestroy)); 34228aa8177SBarry Smith } 343be5899b3SLisandro Dalcin flg = pseudo->increment_dt_from_initial_dt; 3449566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_pseudo_increment_dt_from_initial_dt", "Increase dt as a ratio from original dt", "TSPseudoIncrementDtFromInitialDt", flg, &flg, NULL)); 345be5899b3SLisandro Dalcin pseudo->increment_dt_from_initial_dt = flg; 3469566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_pseudo_increment", "Ratio to increase dt", "TSPseudoSetTimeStepIncrement", pseudo->dt_increment, &pseudo->dt_increment, NULL)); 3479566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_pseudo_max_dt", "Maximum value for dt", "TSPseudoSetMaxTimeStep", pseudo->dt_max, &pseudo->dt_max, NULL)); 3489566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_pseudo_fatol", "Tolerance for norm of function", "", pseudo->fatol, &pseudo->fatol, NULL)); 3499566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_pseudo_frtol", "Relative tolerance for norm of function", "", pseudo->frtol, &pseudo->frtol, NULL)); 350d0609cedSBarry Smith PetscOptionsHeadEnd(); 3513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3522d3f70b5SBarry Smith } 3532d3f70b5SBarry Smith 354d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_Pseudo(TS ts, PetscViewer viewer) 355d71ae5a4SJacob Faibussowitsch { 3563118ae5eSBarry Smith PetscBool isascii; 357d52bd9f3SBarry Smith 3583a40ed3dSBarry Smith PetscFunctionBegin; 3599566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 3603118ae5eSBarry Smith if (isascii) { 3613118ae5eSBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 3629566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " frtol - relative tolerance in function value %g\n", (double)pseudo->frtol)); 3639566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " fatol - absolute tolerance in function value %g\n", (double)pseudo->fatol)); 3649566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " dt_initial - initial timestep %g\n", (double)pseudo->dt_initial)); 3659566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " dt_increment - increase in timestep on successful step %g\n", (double)pseudo->dt_increment)); 3669566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " dt_max - maximum time %g\n", (double)pseudo->dt_max)); 3673118ae5eSBarry Smith } 3683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3692d3f70b5SBarry Smith } 3702d3f70b5SBarry Smith 371ac226902SBarry Smith /*@C 37282bf6240SBarry Smith TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the 37382bf6240SBarry Smith last timestep. 37482bf6240SBarry Smith 375c3339decSBarry Smith Logically Collective 37615091d37SBarry Smith 37782bf6240SBarry Smith Input Parameters: 37815091d37SBarry Smith + ts - timestep context 37982bf6240SBarry Smith . dt - user-defined function to verify timestep 3802fe279fdSBarry Smith - ctx - [optional] user-defined context for private data for the timestep verification routine (may be `NULL`) 38182bf6240SBarry Smith 38220f4b53cSBarry Smith Calling sequence of `func`: 3832fe279fdSBarry Smith + ts - the time-step context 3842fe279fdSBarry Smith . update - latest solution vector 38514d0ab18SJacob Faibussowitsch . ctx - [optional] user-defined timestep context 38682bf6240SBarry Smith . newdt - the timestep to use for the next step 387a2b725a8SWilliam Gropp - flag - flag indicating whether the last time step was acceptable 38882bf6240SBarry Smith 389bcf0153eSBarry Smith Level: advanced 390bcf0153eSBarry Smith 391bcf0153eSBarry Smith Note: 392bcf0153eSBarry Smith The routine set here will be called by `TSPseudoVerifyTimeStep()` 39382bf6240SBarry Smith during the timestepping process. 39482bf6240SBarry Smith 3951cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoVerifyTimeStepDefault()`, `TSPseudoVerifyTimeStep()` 39682bf6240SBarry Smith @*/ 39714d0ab18SJacob Faibussowitsch PetscErrorCode TSPseudoSetVerifyTimeStep(TS ts, PetscErrorCode (*dt)(TS ts, Vec update, void *ctx, PetscReal *newdt, PetscBool *flag), void *ctx) 398d71ae5a4SJacob Faibussowitsch { 39982bf6240SBarry Smith PetscFunctionBegin; 4000700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 401cac4c232SBarry Smith PetscTryMethod(ts, "TSPseudoSetVerifyTimeStep_C", (TS, PetscErrorCode (*)(TS, Vec, void *, PetscReal *, PetscBool *), void *), (ts, dt, ctx)); 4023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 40382bf6240SBarry Smith } 40482bf6240SBarry Smith 40582bf6240SBarry Smith /*@ 40682bf6240SBarry Smith TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 4078d359177SBarry Smith dt when using the TSPseudoTimeStepDefault() routine. 40882bf6240SBarry Smith 409c3339decSBarry Smith Logically Collective 410fee21e36SBarry Smith 41115091d37SBarry Smith Input Parameters: 41215091d37SBarry Smith + ts - the timestep context 41315091d37SBarry Smith - inc - the scaling factor >= 1.0 41415091d37SBarry Smith 41582bf6240SBarry Smith Options Database Key: 41667b8a455SSatish Balay . -ts_pseudo_increment <increment> - set pseudo increment 41782bf6240SBarry Smith 41815091d37SBarry Smith Level: advanced 41915091d37SBarry Smith 4201cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()` 42182bf6240SBarry Smith @*/ 422d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoSetTimeStepIncrement(TS ts, PetscReal inc) 423d71ae5a4SJacob Faibussowitsch { 42482bf6240SBarry Smith PetscFunctionBegin; 4250700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 426c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts, inc, 2); 427cac4c232SBarry Smith PetscTryMethod(ts, "TSPseudoSetTimeStepIncrement_C", (TS, PetscReal), (ts, inc)); 4283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 42982bf6240SBarry Smith } 43082bf6240SBarry Smith 43186534af1SJed Brown /*@ 43286534af1SJed Brown TSPseudoSetMaxTimeStep - Sets the maximum time step 4338d359177SBarry Smith when using the TSPseudoTimeStepDefault() routine. 43486534af1SJed Brown 435c3339decSBarry Smith Logically Collective 43686534af1SJed Brown 43786534af1SJed Brown Input Parameters: 43886534af1SJed Brown + ts - the timestep context 43986534af1SJed Brown - maxdt - the maximum time step, use a non-positive value to deactivate 44086534af1SJed Brown 44186534af1SJed Brown Options Database Key: 44267b8a455SSatish Balay . -ts_pseudo_max_dt <increment> - set pseudo max dt 44386534af1SJed Brown 44486534af1SJed Brown Level: advanced 44586534af1SJed Brown 4461cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()` 44786534af1SJed Brown @*/ 448d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoSetMaxTimeStep(TS ts, PetscReal maxdt) 449d71ae5a4SJacob Faibussowitsch { 45086534af1SJed Brown PetscFunctionBegin; 45186534af1SJed Brown PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 45286534af1SJed Brown PetscValidLogicalCollectiveReal(ts, maxdt, 2); 453cac4c232SBarry Smith PetscTryMethod(ts, "TSPseudoSetMaxTimeStep_C", (TS, PetscReal), (ts, maxdt)); 4543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 45586534af1SJed Brown } 45686534af1SJed Brown 45782bf6240SBarry Smith /*@ 45882bf6240SBarry Smith TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep 459b44f4de4SBarry 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.$ 46082bf6240SBarry Smith 461c3339decSBarry Smith Logically Collective 46215091d37SBarry Smith 46382bf6240SBarry Smith Input Parameter: 46482bf6240SBarry Smith . ts - the timestep context 46582bf6240SBarry Smith 46682bf6240SBarry Smith Options Database Key: 467b44f4de4SBarry Smith . -ts_pseudo_increment_dt_from_initial_dt <true,false> - use the initial $dt$ to determine increment 46882bf6240SBarry Smith 46915091d37SBarry Smith Level: advanced 47015091d37SBarry Smith 4711cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()` 47282bf6240SBarry Smith @*/ 473d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS ts) 474d71ae5a4SJacob Faibussowitsch { 47582bf6240SBarry Smith PetscFunctionBegin; 4760700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 477cac4c232SBarry Smith PetscTryMethod(ts, "TSPseudoIncrementDtFromInitialDt_C", (TS), (ts)); 4783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 47982bf6240SBarry Smith } 48082bf6240SBarry Smith 481ac226902SBarry Smith /*@C 48282bf6240SBarry Smith TSPseudoSetTimeStep - Sets the user-defined routine to be 48382bf6240SBarry Smith called at each pseudo-timestep to update the timestep. 48482bf6240SBarry Smith 485c3339decSBarry Smith Logically Collective 48615091d37SBarry Smith 48782bf6240SBarry Smith Input Parameters: 48815091d37SBarry Smith + ts - timestep context 48982bf6240SBarry Smith . dt - function to compute timestep 4902fe279fdSBarry Smith - ctx - [optional] user-defined context for private data required by the function (may be `NULL`) 49182bf6240SBarry Smith 49220f4b53cSBarry Smith Calling sequence of `dt`: 49314d0ab18SJacob Faibussowitsch + ts - the `TS` context 49414d0ab18SJacob Faibussowitsch . newdt - the newly computed timestep 49514d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined context 49682bf6240SBarry Smith 497bcf0153eSBarry Smith Level: intermediate 49882bf6240SBarry Smith 499bcf0153eSBarry Smith Notes: 500bcf0153eSBarry Smith The routine set here will be called by `TSPseudoComputeTimeStep()` 501bcf0153eSBarry Smith during the timestepping process. 502bcf0153eSBarry Smith 503bcf0153eSBarry Smith If not set then `TSPseudoTimeStepDefault()` is automatically used 504bcf0153eSBarry Smith 5051cc06b55SBarry Smith .seealso: [](ch_ts), `TSPSEUDO`, `TSPseudoTimeStepDefault()`, `TSPseudoComputeTimeStep()` 50682bf6240SBarry Smith @*/ 50714d0ab18SJacob Faibussowitsch PetscErrorCode TSPseudoSetTimeStep(TS ts, PetscErrorCode (*dt)(TS ts, PetscReal *newdt, void *ctx), void *ctx) 508d71ae5a4SJacob Faibussowitsch { 50982bf6240SBarry Smith PetscFunctionBegin; 5100700a824SBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 511cac4c232SBarry Smith PetscTryMethod(ts, "TSPseudoSetTimeStep_C", (TS, PetscErrorCode (*)(TS, PetscReal *, void *), void *), (ts, dt, ctx)); 5123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 51382bf6240SBarry Smith } 51482bf6240SBarry Smith 51582bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 51682bf6240SBarry Smith 517ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(TS, Vec, void *, PetscReal *, PetscBool *); /* force argument to next function to not be extern C*/ 518d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetVerifyTimeStep_Pseudo(TS ts, FCN1 dt, void *ctx) 519d71ae5a4SJacob Faibussowitsch { 520be5899b3SLisandro Dalcin TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 52182bf6240SBarry Smith 52282bf6240SBarry Smith PetscFunctionBegin; 52382bf6240SBarry Smith pseudo->verify = dt; 52482bf6240SBarry Smith pseudo->verifyctx = ctx; 5253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 52682bf6240SBarry Smith } 52782bf6240SBarry Smith 528d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetTimeStepIncrement_Pseudo(TS ts, PetscReal inc) 529d71ae5a4SJacob Faibussowitsch { 5304bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 53182bf6240SBarry Smith 53282bf6240SBarry Smith PetscFunctionBegin; 53382bf6240SBarry Smith pseudo->dt_increment = inc; 5343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 53582bf6240SBarry Smith } 53682bf6240SBarry Smith 537d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetMaxTimeStep_Pseudo(TS ts, PetscReal maxdt) 538d71ae5a4SJacob Faibussowitsch { 53986534af1SJed Brown TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 54086534af1SJed Brown 54186534af1SJed Brown PetscFunctionBegin; 54286534af1SJed Brown pseudo->dt_max = maxdt; 5433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 54486534af1SJed Brown } 54586534af1SJed Brown 546d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts) 547d71ae5a4SJacob Faibussowitsch { 5484bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 54982bf6240SBarry Smith 55082bf6240SBarry Smith PetscFunctionBegin; 5514bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_TRUE; 5523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 55382bf6240SBarry Smith } 55482bf6240SBarry Smith 5556849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS, PetscReal *, void *); /* force argument to next function to not be extern C*/ 556d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPseudoSetTimeStep_Pseudo(TS ts, FCN2 dt, void *ctx) 557d71ae5a4SJacob Faibussowitsch { 5584bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 55982bf6240SBarry Smith 56082bf6240SBarry Smith PetscFunctionBegin; 56182bf6240SBarry Smith pseudo->dt = dt; 56282bf6240SBarry Smith pseudo->dtctx = ctx; 5633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 56482bf6240SBarry Smith } 56582bf6240SBarry Smith 56610e6a065SJed Brown /*MC 5671d27aa22SBarry Smith TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping {cite}`ckk02` {cite}`kk97` 56882bf6240SBarry Smith 56910e6a065SJed Brown This method solves equations of the form 57010e6a065SJed Brown 5711d27aa22SBarry Smith $$ 5721d27aa22SBarry Smith F(X,Xdot) = 0 5731d27aa22SBarry Smith $$ 57410e6a065SJed Brown 57510e6a065SJed Brown for steady state using the iteration 57610e6a065SJed Brown 5771d27aa22SBarry Smith $$ 57837fdd005SBarry Smith [G'] S = -F(X,0) 57937fdd005SBarry Smith X += S 5801d27aa22SBarry Smith $$ 58110e6a065SJed Brown 58210e6a065SJed Brown where 58310e6a065SJed Brown 5841d27aa22SBarry Smith $$ 5851d27aa22SBarry Smith G(Y) = F(Y,(Y-X)/dt) 5861d27aa22SBarry Smith $$ 58710e6a065SJed Brown 5886f2d6a7bSJed Brown This is linearly-implicit Euler with the residual always evaluated "at steady 5896f2d6a7bSJed Brown state". See note below. 59010e6a065SJed Brown 59193a54799SPierre Jolivet In addition to the modified solve, a dedicated adaptive timestepping scheme is implemented, mimicking the switched evolution relaxation in {cite}`ckk02`. 592fb59f417SJames Wright It determines the next timestep via 593fb59f417SJames Wright 594fb59f417SJames Wright $$ 595fb59f417SJames Wright dt_{n} = r dt_{n-1} \frac{\Vert F(X_{n},0) \Vert}{\Vert F(X_{n-1},0) \Vert} 596fb59f417SJames Wright $$ 597fb59f417SJames Wright 598fb59f417SJames Wright where $r$ is an additional growth factor (set by `-ts_pseudo_increment`). 599fb59f417SJames Wright An alternative formulation is also available that uses the initial timestep and function norm. 600fb59f417SJames Wright 601fb59f417SJames Wright $$ 602fb59f417SJames Wright dt_{n} = r dt_{0} \frac{\Vert F(X_{n},0) \Vert}{\Vert F(X_{0},0) \Vert} 603fb59f417SJames Wright $$ 604fb59f417SJames Wright 605fb59f417SJames Wright This is chosen by setting `-ts_pseudo_increment_dt_from_initial_dt`. 606fb59f417SJames Wright For either method, an upper limit on the timestep can be set by `-ts_pseudo_max_dt`. 607fb59f417SJames Wright 608bcf0153eSBarry Smith Options Database Keys: 60910e6a065SJed Brown + -ts_pseudo_increment <real> - ratio of increase dt 6103118ae5eSBarry Smith . -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt 611fb59f417SJames Wright . -ts_pseudo_max_dt - Maximum dt for adaptive timestepping algorithm 612fb59f417SJames Wright . -ts_pseudo_monitor - Monitor convergence of the function norm 613fb59f417SJames Wright . -ts_pseudo_fatol <atol> - stop iterating when the function norm is less than `atol` 614fb59f417SJames Wright - -ts_pseudo_frtol <rtol> - stop iterating when the function norm divided by the initial function norm is less than `rtol` 61510e6a065SJed Brown 61610e6a065SJed Brown Level: beginner 61710e6a065SJed Brown 61810e6a065SJed Brown Notes: 6196f2d6a7bSJed Brown The residual computed by this method includes the transient term (Xdot is computed instead of 6206f2d6a7bSJed Brown always being zero), but since the prediction from the last step is always the solution from the 6216f2d6a7bSJed Brown last step, on the first Newton iteration we have 6226f2d6a7bSJed Brown 6231d27aa22SBarry Smith $$ 6241d27aa22SBarry Smith Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0 6251d27aa22SBarry Smith $$ 6266f2d6a7bSJed Brown 627eb636060SBarry Smith The Jacobian $ dF/dX + shift*dF/dXdot $ contains a non-zero $ dF/dXdot$ term. In the $ X' = F(X) $ case it 628eb636060SBarry Smith is $ dF/dX + shift*1/dt $. Hence still contains the $ 1/dt $ term so the Jacobian is not simply the 629eb636060SBarry Smith Jacobian of $ F $ and thus this pseudo-transient continuation is not just Newton on $F(x)=0$. 630eb636060SBarry Smith 6316f2d6a7bSJed Brown Therefore, the linear system solved by the first Newton iteration is equivalent to the one 6326f2d6a7bSJed Brown described above and in the papers. If the user chooses to perform multiple Newton iterations, the 6336f2d6a7bSJed Brown algorithm is no longer the one described in the referenced papers. 634fb59f417SJames Wright By default, the `SNESType` is set to `SNESKSPONLY` to match the algorithm from the referenced papers. 635fb59f417SJames Wright Setting the `SNESType` via `-snes_type` will override this default setting. 63610e6a065SJed Brown 6371cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()` 63810e6a065SJed Brown M*/ 639d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_Pseudo(TS ts) 640d71ae5a4SJacob Faibussowitsch { 6417bf11e45SBarry Smith TS_Pseudo *pseudo; 642193ac0bcSJed Brown SNES snes; 64319fd82e9SBarry Smith SNESType stype; 6442d3f70b5SBarry Smith 6453a40ed3dSBarry Smith PetscFunctionBegin; 646277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Pseudo; 647000e7ae3SMatthew Knepley ts->ops->destroy = TSDestroy_Pseudo; 648000e7ae3SMatthew Knepley ts->ops->view = TSView_Pseudo; 649000e7ae3SMatthew Knepley ts->ops->setup = TSSetUp_Pseudo; 650000e7ae3SMatthew Knepley ts->ops->step = TSStep_Pseudo; 651000e7ae3SMatthew Knepley ts->ops->setfromoptions = TSSetFromOptions_Pseudo; 6520f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Pseudo; 6530f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Pseudo; 6542ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 655825ab935SBarry Smith ts->usessnes = PETSC_TRUE; 6567bf11e45SBarry Smith 6579566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 6589566063dSJacob Faibussowitsch PetscCall(SNESGetType(snes, &stype)); 6599566063dSJacob Faibussowitsch if (!stype) PetscCall(SNESSetType(snes, SNESKSPONLY)); 6602d3f70b5SBarry Smith 6614dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&pseudo)); 6627bf11e45SBarry Smith ts->data = (void *)pseudo; 6632d3f70b5SBarry Smith 664be5899b3SLisandro Dalcin pseudo->dt = TSPseudoTimeStepDefault; 665be5899b3SLisandro Dalcin pseudo->dtctx = NULL; 66628aa8177SBarry Smith pseudo->dt_increment = 1.1; 6674bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_FALSE; 668193ac0bcSJed Brown pseudo->fnorm = -1; 669be5899b3SLisandro Dalcin pseudo->fnorm_initial = -1; 670be5899b3SLisandro Dalcin pseudo->fnorm_previous = -1; 6713118ae5eSBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 6723118ae5eSBarry Smith pseudo->fatol = 1.e-25; 6733118ae5eSBarry Smith pseudo->frtol = 1.e-5; 6743118ae5eSBarry Smith #else 6753118ae5eSBarry Smith pseudo->fatol = 1.e-50; 6763118ae5eSBarry Smith pseudo->frtol = 1.e-12; 6773118ae5eSBarry Smith #endif 6789566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetVerifyTimeStep_C", TSPseudoSetVerifyTimeStep_Pseudo)); 6799566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStepIncrement_C", TSPseudoSetTimeStepIncrement_Pseudo)); 6809566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetMaxTimeStep_C", TSPseudoSetMaxTimeStep_Pseudo)); 6819566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoIncrementDtFromInitialDt_C", TSPseudoIncrementDtFromInitialDt_Pseudo)); 6829566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSPseudoSetTimeStep_C", TSPseudoSetTimeStep_Pseudo)); 6833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6842d3f70b5SBarry Smith } 6852d3f70b5SBarry Smith 68682bf6240SBarry Smith /*@C 687bcf0153eSBarry Smith TSPseudoTimeStepDefault - Default code to compute pseudo-timestepping. Use with `TSPseudoSetTimeStep()`. 68828aa8177SBarry Smith 689cc4c1da9SBarry Smith Collective, No Fortran Support 69015091d37SBarry Smith 69128aa8177SBarry Smith Input Parameters: 692a2b725a8SWilliam Gropp + ts - the timestep context 693a2b725a8SWilliam Gropp - dtctx - unused timestep context 69428aa8177SBarry Smith 69582bf6240SBarry Smith Output Parameter: 69682bf6240SBarry Smith . newdt - the timestep to use for the next step 69728aa8177SBarry Smith 69815091d37SBarry Smith Level: advanced 69915091d37SBarry Smith 7001cc06b55SBarry Smith .seealso: [](ch_ts), `TSPseudoSetTimeStep()`, `TSPseudoComputeTimeStep()`, `TSPSEUDO` 70128aa8177SBarry Smith @*/ 702d71ae5a4SJacob Faibussowitsch PetscErrorCode TSPseudoTimeStepDefault(TS ts, PetscReal *newdt, void *dtctx) 703d71ae5a4SJacob Faibussowitsch { 70482bf6240SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo *)ts->data; 705be5899b3SLisandro Dalcin PetscReal inc = pseudo->dt_increment; 70628aa8177SBarry Smith 7073a40ed3dSBarry Smith PetscFunctionBegin; 70898473560SBarry Smith if (pseudo->fnorm < 0.0) { 7099566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(pseudo->xdot)); 7109566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts, ts->ptime, ts->vec_sol, pseudo->xdot, pseudo->func, PETSC_FALSE)); 711eb636060SBarry Smith PetscCall(PetscObjectStateGet((PetscObject)ts->vec_sol, &pseudo->Xstate)); 7129566063dSJacob Faibussowitsch PetscCall(VecNorm(pseudo->func, NORM_2, &pseudo->fnorm)); 71398473560SBarry Smith } 714be5899b3SLisandro Dalcin if (pseudo->fnorm_initial < 0) { 71582bf6240SBarry Smith /* first time through so compute initial function norm */ 716cdbf8f93SLisandro Dalcin pseudo->fnorm_initial = pseudo->fnorm; 717be5899b3SLisandro Dalcin pseudo->fnorm_previous = pseudo->fnorm; 71882bf6240SBarry Smith } 719bbd56ea5SKarl Rupp if (pseudo->fnorm == 0.0) *newdt = 1.e12 * inc * ts->time_step; 720bbd56ea5SKarl Rupp else if (pseudo->increment_dt_from_initial_dt) *newdt = inc * pseudo->dt_initial * pseudo->fnorm_initial / pseudo->fnorm; 721be5899b3SLisandro Dalcin else *newdt = inc * ts->time_step * pseudo->fnorm_previous / pseudo->fnorm; 72286534af1SJed Brown if (pseudo->dt_max > 0) *newdt = PetscMin(*newdt, pseudo->dt_max); 72382bf6240SBarry Smith pseudo->fnorm_previous = pseudo->fnorm; 7243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 72528aa8177SBarry Smith } 726