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 @*/ 507087cfbeSBarry Smith PetscErrorCode TSPseudoComputeTimeStep(TS ts,PetscReal *dt) 517bf11e45SBarry Smith { 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 @*/ 848d359177SBarry Smith PetscErrorCode TSPseudoVerifyTimeStepDefault(TS ts,Vec update,void *dtctx,PetscReal *newdt,PetscBool *flag) 857bf11e45SBarry Smith { 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 @*/ 1127087cfbeSBarry Smith PetscErrorCode TSPseudoVerifyTimeStep(TS ts,Vec update,PetscReal *dt,PetscBool *flag) 1137bf11e45SBarry Smith { 1147bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 1157bf11e45SBarry Smith 1163a40ed3dSBarry Smith PetscFunctionBegin; 117cb9d8021SPierre Barbier de Reuille *flag = PETSC_TRUE; 118*1baa6e33SBarry 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 124193ac0bcSJed Brown static PetscErrorCode TSStep_Pseudo(TS ts) 1252d3f70b5SBarry Smith { 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)); 141be5899b3SLisandro Dalcin ts->snes_its += nits; ts->ksp_its += lits; 1429566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts,ts->ptime+ts->time_step,0,&(pseudo->update))); 1439566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(ts->adapt,ts,ts->ptime+ts->time_step,pseudo->update,&stepok)); 144be5899b3SLisandro Dalcin if (!stepok) {next_time_step = ts->time_step; continue;} 145193ac0bcSJed Brown pseudo->fnorm = -1; /* The current norm is no longer valid, monitor must recompute it. */ 1469566063dSJacob Faibussowitsch PetscCall(TSPseudoVerifyTimeStep(ts,pseudo->update,&next_time_step,&stepok)); 147cdbf8f93SLisandro Dalcin if (stepok) break; 148cdbf8f93SLisandro Dalcin } 149be5899b3SLisandro Dalcin if (reject >= ts->max_reject) { 150be5899b3SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 15163a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts,"Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n",ts->steps,reject)); 152cdbf8f93SLisandro Dalcin PetscFunctionReturn(0); 1537bf11e45SBarry Smith } 154be5899b3SLisandro Dalcin 1559566063dSJacob Faibussowitsch PetscCall(VecCopy(pseudo->update,ts->vec_sol)); 156be5899b3SLisandro Dalcin ts->ptime += ts->time_step; 157be5899b3SLisandro Dalcin ts->time_step = next_time_step; 158be5899b3SLisandro Dalcin 1593118ae5eSBarry Smith if (pseudo->fnorm < 0) { 1609566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(pseudo->xdot)); 1619566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE)); 1629566063dSJacob Faibussowitsch PetscCall(VecNorm(pseudo->func,NORM_2,&pseudo->fnorm)); 1633118ae5eSBarry Smith } 1643118ae5eSBarry Smith if (pseudo->fnorm < pseudo->fatol) { 1653118ae5eSBarry Smith ts->reason = TS_CONVERGED_PSEUDO_FATOL; 16663a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts,"Step=%" PetscInt_FMT ", converged since fnorm %g < fatol %g\n",ts->steps,(double)pseudo->fnorm,(double)pseudo->frtol)); 1673118ae5eSBarry Smith PetscFunctionReturn(0); 1683118ae5eSBarry Smith } 1693118ae5eSBarry Smith if (pseudo->fnorm/pseudo->fnorm_initial < pseudo->frtol) { 1703118ae5eSBarry Smith ts->reason = TS_CONVERGED_PSEUDO_FRTOL; 17163a3b9bcSJacob 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)); 1723118ae5eSBarry Smith PetscFunctionReturn(0); 1733118ae5eSBarry Smith } 1743a40ed3dSBarry Smith PetscFunctionReturn(0); 1752d3f70b5SBarry Smith } 1762d3f70b5SBarry Smith 1772d3f70b5SBarry Smith /*------------------------------------------------------------*/ 178277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Pseudo(TS ts) 1792d3f70b5SBarry Smith { 1807bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 1812d3f70b5SBarry Smith 1823a40ed3dSBarry Smith PetscFunctionBegin; 1839566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pseudo->update)); 1849566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pseudo->func)); 1859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pseudo->xdot)); 1863a40ed3dSBarry Smith PetscFunctionReturn(0); 1872d3f70b5SBarry Smith } 1882d3f70b5SBarry Smith 189277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Pseudo(TS ts) 190277b19d0SLisandro Dalcin { 191277b19d0SLisandro Dalcin PetscFunctionBegin; 1929566063dSJacob Faibussowitsch PetscCall(TSReset_Pseudo(ts)); 1939566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 1949566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",NULL)); 1959566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",NULL)); 1969566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetMaxTimeStep_C",NULL)); 1979566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",NULL)); 1989566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",NULL)); 199277b19d0SLisandro Dalcin PetscFunctionReturn(0); 200277b19d0SLisandro Dalcin } 2012d3f70b5SBarry Smith 2022d3f70b5SBarry Smith /*------------------------------------------------------------*/ 2032d3f70b5SBarry Smith 2046f2d6a7bSJed Brown /* 2056f2d6a7bSJed Brown Compute Xdot = (X^{n+1}-X^n)/dt) = 0 2066f2d6a7bSJed Brown */ 2076f2d6a7bSJed Brown static PetscErrorCode TSPseudoGetXdot(TS ts,Vec X,Vec *Xdot) 2082d3f70b5SBarry Smith { 2096f2d6a7bSJed Brown TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 210193ac0bcSJed Brown const PetscScalar mdt = 1.0/ts->time_step,*xnp1,*xn; 211193ac0bcSJed Brown PetscScalar *xdot; 212a7cc72afSBarry Smith PetscInt i,n; 2132d3f70b5SBarry Smith 2143a40ed3dSBarry Smith PetscFunctionBegin; 215aab5bcd8SJed Brown *Xdot = NULL; 2169566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(ts->vec_sol,&xn)); 2179566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&xnp1)); 2189566063dSJacob Faibussowitsch PetscCall(VecGetArray(pseudo->xdot,&xdot)); 2199566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X,&n)); 220bbd56ea5SKarl Rupp for (i=0; i<n; i++) xdot[i] = mdt*(xnp1[i] - xn[i]); 2219566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(ts->vec_sol,&xn)); 2229566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&xnp1)); 2239566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(pseudo->xdot,&xdot)); 2246f2d6a7bSJed Brown *Xdot = pseudo->xdot; 2253a40ed3dSBarry Smith PetscFunctionReturn(0); 2262d3f70b5SBarry Smith } 2272d3f70b5SBarry Smith 2286f2d6a7bSJed Brown /* 2296f2d6a7bSJed Brown The transient residual is 2306f2d6a7bSJed Brown 2316f2d6a7bSJed Brown F(U^{n+1},(U^{n+1}-U^n)/dt) = 0 2326f2d6a7bSJed Brown 2336f2d6a7bSJed Brown or for ODE, 2346f2d6a7bSJed Brown 2356f2d6a7bSJed Brown (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0 2366f2d6a7bSJed Brown 2376f2d6a7bSJed Brown This is the function that must be evaluated for transient simulation and for 2386f2d6a7bSJed Brown finite difference Jacobians. On the first Newton step, this algorithm uses 2396f2d6a7bSJed Brown a guess of U^{n+1} = U^n in which case the transient term vanishes and the 2406f2d6a7bSJed Brown residual is actually the steady state residual. Pseudotransient 2416f2d6a7bSJed Brown continuation as described in the literature is a linearly implicit 2426f2d6a7bSJed Brown algorithm, it only takes this one Newton step with the steady state 2436f2d6a7bSJed Brown residual, and then advances to the next time step. 2446f2d6a7bSJed Brown */ 2450f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes,Vec X,Vec Y,TS ts) 2462d3f70b5SBarry Smith { 2476f2d6a7bSJed Brown Vec Xdot; 2482d3f70b5SBarry Smith 2493a40ed3dSBarry Smith PetscFunctionBegin; 2509566063dSJacob Faibussowitsch PetscCall(TSPseudoGetXdot(ts,X,&Xdot)); 2519566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts,ts->ptime+ts->time_step,X,Xdot,Y,PETSC_FALSE)); 2526f2d6a7bSJed Brown PetscFunctionReturn(0); 2536f2d6a7bSJed Brown } 2542d3f70b5SBarry Smith 2556f2d6a7bSJed Brown /* 2566f2d6a7bSJed Brown This constructs the Jacobian needed for SNES. For DAE, this is 2576f2d6a7bSJed Brown 2586f2d6a7bSJed Brown dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot 2596f2d6a7bSJed Brown 2606f2d6a7bSJed Brown and for ODE: 2616f2d6a7bSJed Brown 2626f2d6a7bSJed Brown J = I/dt - J_{Frhs} where J_{Frhs} is the given Jacobian of Frhs. 2636f2d6a7bSJed Brown */ 264d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes,Vec X,Mat AA,Mat BB,TS ts) 2656f2d6a7bSJed Brown { 2666f2d6a7bSJed Brown Vec Xdot; 2676f2d6a7bSJed Brown 2686f2d6a7bSJed Brown PetscFunctionBegin; 2699566063dSJacob Faibussowitsch PetscCall(TSPseudoGetXdot(ts,X,&Xdot)); 2709566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts,ts->ptime+ts->time_step,X,Xdot,1./ts->time_step,AA,BB,PETSC_FALSE)); 2713a40ed3dSBarry Smith PetscFunctionReturn(0); 2722d3f70b5SBarry Smith } 2732d3f70b5SBarry Smith 2746849ba73SBarry Smith static PetscErrorCode TSSetUp_Pseudo(TS ts) 2752d3f70b5SBarry Smith { 2767bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 2772d3f70b5SBarry Smith 2783a40ed3dSBarry Smith PetscFunctionBegin; 2799566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&pseudo->update)); 2809566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&pseudo->func)); 2819566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&pseudo->xdot)); 2823a40ed3dSBarry Smith PetscFunctionReturn(0); 2832d3f70b5SBarry Smith } 2842d3f70b5SBarry Smith /*------------------------------------------------------------*/ 2852d3f70b5SBarry Smith 286560360afSLisandro Dalcin static PetscErrorCode TSPseudoMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy) 2872d3f70b5SBarry Smith { 2887bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 289ce94432eSBarry Smith PetscViewer viewer = (PetscViewer) dummy; 2902d3f70b5SBarry Smith 2913a40ed3dSBarry Smith PetscFunctionBegin; 292193ac0bcSJed Brown if (pseudo->fnorm < 0) { /* The last computed norm is stale, recompute */ 2939566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(pseudo->xdot)); 2949566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE)); 2959566063dSJacob Faibussowitsch PetscCall(VecNorm(pseudo->func,NORM_2,&pseudo->fnorm)); 296193ac0bcSJed Brown } 2979566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel)); 29863a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"TS %" PetscInt_FMT " dt %g time %g fnorm %g\n",step,(double)ts->time_step,(double)ptime,(double)pseudo->fnorm)); 2999566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel)); 3003a40ed3dSBarry Smith PetscFunctionReturn(0); 3012d3f70b5SBarry Smith } 3022d3f70b5SBarry Smith 3034416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Pseudo(PetscOptionItems *PetscOptionsObject,TS ts) 3042d3f70b5SBarry Smith { 3054bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 306ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE; 307649052a6SBarry Smith PetscViewer viewer; 3082d3f70b5SBarry Smith 3093a40ed3dSBarry Smith PetscFunctionBegin; 310d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"Pseudo-timestepping options"); 3119566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_monitor_pseudo","Monitor convergence","",flg,&flg,NULL)); 3122d3f70b5SBarry Smith if (flg) { 3139566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ts),"stdout",&viewer)); 3149566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts,TSPseudoMonitorDefault,viewer,(PetscErrorCode (*)(void**))PetscViewerDestroy)); 31528aa8177SBarry Smith } 316be5899b3SLisandro Dalcin flg = pseudo->increment_dt_from_initial_dt; 3179566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",flg,&flg,NULL)); 318be5899b3SLisandro Dalcin pseudo->increment_dt_from_initial_dt = flg; 3199566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,NULL)); 3209566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_pseudo_max_dt","Maximum value for dt","TSPseudoSetMaxTimeStep",pseudo->dt_max,&pseudo->dt_max,NULL)); 3219566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_pseudo_fatol","Tolerance for norm of function","",pseudo->fatol,&pseudo->fatol,NULL)); 3229566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_pseudo_frtol","Relative tolerance for norm of function","",pseudo->frtol,&pseudo->frtol,NULL)); 323d0609cedSBarry Smith PetscOptionsHeadEnd(); 3243a40ed3dSBarry Smith PetscFunctionReturn(0); 3252d3f70b5SBarry Smith } 3262d3f70b5SBarry Smith 3276849ba73SBarry Smith static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer) 3282d3f70b5SBarry Smith { 3293118ae5eSBarry Smith PetscBool isascii; 330d52bd9f3SBarry Smith 3313a40ed3dSBarry Smith PetscFunctionBegin; 3329566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii)); 3333118ae5eSBarry Smith if (isascii) { 3343118ae5eSBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 3359566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," frtol - relative tolerance in function value %g\n",(double)pseudo->frtol)); 3369566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," fatol - absolute tolerance in function value %g\n",(double)pseudo->fatol)); 3379566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," dt_initial - initial timestep %g\n",(double)pseudo->dt_initial)); 3389566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," dt_increment - increase in timestep on successful step %g\n",(double)pseudo->dt_increment)); 3399566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," dt_max - maximum time %g\n",(double)pseudo->dt_max)); 3403118ae5eSBarry Smith } 3413a40ed3dSBarry Smith PetscFunctionReturn(0); 3422d3f70b5SBarry Smith } 3432d3f70b5SBarry Smith 34482bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 345ac226902SBarry Smith /*@C 34682bf6240SBarry Smith TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the 34782bf6240SBarry Smith last timestep. 34882bf6240SBarry Smith 3493f9fe445SBarry Smith Logically Collective on TS 35015091d37SBarry Smith 35182bf6240SBarry Smith Input Parameters: 35215091d37SBarry Smith + ts - timestep context 35382bf6240SBarry Smith . dt - user-defined function to verify timestep 35415091d37SBarry Smith - ctx - [optional] user-defined context for private data 3550298fd71SBarry Smith for the timestep verification routine (may be NULL) 35682bf6240SBarry Smith 35715091d37SBarry Smith Level: advanced 358fee21e36SBarry Smith 35982bf6240SBarry Smith Calling sequence of func: 360a2b725a8SWilliam Gropp $ func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscBool *flag); 36182bf6240SBarry Smith 362a2b725a8SWilliam Gropp + update - latest solution vector 36382bf6240SBarry Smith . ctx - [optional] timestep context 36482bf6240SBarry Smith . newdt - the timestep to use for the next step 365a2b725a8SWilliam Gropp - flag - flag indicating whether the last time step was acceptable 36682bf6240SBarry Smith 36782bf6240SBarry Smith Notes: 36882bf6240SBarry Smith The routine set here will be called by TSPseudoVerifyTimeStep() 36982bf6240SBarry Smith during the timestepping process. 37082bf6240SBarry Smith 371db781477SPatrick Sanan .seealso: `TSPseudoVerifyTimeStepDefault()`, `TSPseudoVerifyTimeStep()` 37282bf6240SBarry Smith @*/ 3737087cfbeSBarry Smith PetscErrorCode TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,PetscBool*),void *ctx) 37482bf6240SBarry Smith { 37582bf6240SBarry Smith PetscFunctionBegin; 3760700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 377cac4c232SBarry Smith PetscTryMethod(ts,"TSPseudoSetVerifyTimeStep_C",(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal*,PetscBool*),void*),(ts,dt,ctx)); 37882bf6240SBarry Smith PetscFunctionReturn(0); 37982bf6240SBarry Smith } 38082bf6240SBarry Smith 38182bf6240SBarry Smith /*@ 38282bf6240SBarry Smith TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 3838d359177SBarry Smith dt when using the TSPseudoTimeStepDefault() routine. 38482bf6240SBarry Smith 3853f9fe445SBarry Smith Logically Collective on TS 386fee21e36SBarry Smith 38715091d37SBarry Smith Input Parameters: 38815091d37SBarry Smith + ts - the timestep context 38915091d37SBarry Smith - inc - the scaling factor >= 1.0 39015091d37SBarry Smith 39182bf6240SBarry Smith Options Database Key: 39267b8a455SSatish Balay . -ts_pseudo_increment <increment> - set pseudo increment 39382bf6240SBarry Smith 39415091d37SBarry Smith Level: advanced 39515091d37SBarry Smith 396db781477SPatrick Sanan .seealso: `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()` 39782bf6240SBarry Smith @*/ 3987087cfbeSBarry Smith PetscErrorCode TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc) 39982bf6240SBarry Smith { 40082bf6240SBarry Smith PetscFunctionBegin; 4010700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 402c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts,inc,2); 403cac4c232SBarry Smith PetscTryMethod(ts,"TSPseudoSetTimeStepIncrement_C",(TS,PetscReal),(ts,inc)); 40482bf6240SBarry Smith PetscFunctionReturn(0); 40582bf6240SBarry Smith } 40682bf6240SBarry Smith 40786534af1SJed Brown /*@ 40886534af1SJed Brown TSPseudoSetMaxTimeStep - Sets the maximum time step 4098d359177SBarry Smith when using the TSPseudoTimeStepDefault() routine. 41086534af1SJed Brown 41186534af1SJed Brown Logically Collective on TS 41286534af1SJed Brown 41386534af1SJed Brown Input Parameters: 41486534af1SJed Brown + ts - the timestep context 41586534af1SJed Brown - maxdt - the maximum time step, use a non-positive value to deactivate 41686534af1SJed Brown 41786534af1SJed Brown Options Database Key: 41867b8a455SSatish Balay . -ts_pseudo_max_dt <increment> - set pseudo max dt 41986534af1SJed Brown 42086534af1SJed Brown Level: advanced 42186534af1SJed Brown 422db781477SPatrick Sanan .seealso: `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()` 42386534af1SJed Brown @*/ 42486534af1SJed Brown PetscErrorCode TSPseudoSetMaxTimeStep(TS ts,PetscReal maxdt) 42586534af1SJed Brown { 42686534af1SJed Brown PetscFunctionBegin; 42786534af1SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 42886534af1SJed Brown PetscValidLogicalCollectiveReal(ts,maxdt,2); 429cac4c232SBarry Smith PetscTryMethod(ts,"TSPseudoSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt)); 43086534af1SJed Brown PetscFunctionReturn(0); 43186534af1SJed Brown } 43286534af1SJed Brown 43382bf6240SBarry Smith /*@ 43482bf6240SBarry Smith TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep 43582bf6240SBarry Smith is computed via the formula 43682bf6240SBarry Smith $ dt = initial_dt*initial_fnorm/current_fnorm 43782bf6240SBarry Smith rather than the default update, 43882bf6240SBarry Smith $ dt = current_dt*previous_fnorm/current_fnorm. 43982bf6240SBarry Smith 4403f9fe445SBarry Smith Logically Collective on TS 44115091d37SBarry Smith 44282bf6240SBarry Smith Input Parameter: 44382bf6240SBarry Smith . ts - the timestep context 44482bf6240SBarry Smith 44582bf6240SBarry Smith Options Database Key: 44667b8a455SSatish Balay . -ts_pseudo_increment_dt_from_initial_dt <true,false> - use the initial dt to determine increment 44782bf6240SBarry Smith 44815091d37SBarry Smith Level: advanced 44915091d37SBarry Smith 450db781477SPatrick Sanan .seealso: `TSPseudoSetTimeStep()`, `TSPseudoTimeStepDefault()` 45182bf6240SBarry Smith @*/ 4527087cfbeSBarry Smith PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS ts) 45382bf6240SBarry Smith { 45482bf6240SBarry Smith PetscFunctionBegin; 4550700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 456cac4c232SBarry Smith PetscTryMethod(ts,"TSPseudoIncrementDtFromInitialDt_C",(TS),(ts)); 45782bf6240SBarry Smith PetscFunctionReturn(0); 45882bf6240SBarry Smith } 45982bf6240SBarry Smith 460ac226902SBarry Smith /*@C 46182bf6240SBarry Smith TSPseudoSetTimeStep - Sets the user-defined routine to be 46282bf6240SBarry Smith called at each pseudo-timestep to update the timestep. 46382bf6240SBarry Smith 4643f9fe445SBarry Smith Logically Collective on TS 46515091d37SBarry Smith 46682bf6240SBarry Smith Input Parameters: 46715091d37SBarry Smith + ts - timestep context 46882bf6240SBarry Smith . dt - function to compute timestep 46915091d37SBarry Smith - ctx - [optional] user-defined context for private data 4700298fd71SBarry Smith required by the function (may be NULL) 47182bf6240SBarry Smith 47215091d37SBarry Smith Level: intermediate 473fee21e36SBarry Smith 47482bf6240SBarry Smith Calling sequence of func: 475a2b725a8SWilliam Gropp $ func (TS ts,PetscReal *newdt,void *ctx); 47682bf6240SBarry Smith 477a2b725a8SWilliam Gropp + newdt - the newly computed timestep 478a2b725a8SWilliam Gropp - ctx - [optional] timestep context 47982bf6240SBarry Smith 48082bf6240SBarry Smith Notes: 48182bf6240SBarry Smith The routine set here will be called by TSPseudoComputeTimeStep() 48282bf6240SBarry Smith during the timestepping process. 4838d359177SBarry Smith If not set then TSPseudoTimeStepDefault() is automatically used 48482bf6240SBarry Smith 485db781477SPatrick Sanan .seealso: `TSPseudoTimeStepDefault()`, `TSPseudoComputeTimeStep()` 48682bf6240SBarry Smith @*/ 4877087cfbeSBarry Smith PetscErrorCode TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void *ctx) 48882bf6240SBarry Smith { 48982bf6240SBarry Smith PetscFunctionBegin; 4900700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 491cac4c232SBarry Smith PetscTryMethod(ts,"TSPseudoSetTimeStep_C",(TS,PetscErrorCode (*)(TS,PetscReal*,void*),void*),(ts,dt,ctx)); 49282bf6240SBarry Smith PetscFunctionReturn(0); 49382bf6240SBarry Smith } 49482bf6240SBarry Smith 49582bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 49682bf6240SBarry Smith 497ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscBool*); /* force argument to next function to not be extern C*/ 498560360afSLisandro Dalcin static PetscErrorCode TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void *ctx) 49982bf6240SBarry Smith { 500be5899b3SLisandro Dalcin TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 50182bf6240SBarry Smith 50282bf6240SBarry Smith PetscFunctionBegin; 50382bf6240SBarry Smith pseudo->verify = dt; 50482bf6240SBarry Smith pseudo->verifyctx = ctx; 50582bf6240SBarry Smith PetscFunctionReturn(0); 50682bf6240SBarry Smith } 50782bf6240SBarry Smith 508560360afSLisandro Dalcin static PetscErrorCode TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc) 50982bf6240SBarry Smith { 5104bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 51182bf6240SBarry Smith 51282bf6240SBarry Smith PetscFunctionBegin; 51382bf6240SBarry Smith pseudo->dt_increment = inc; 51482bf6240SBarry Smith PetscFunctionReturn(0); 51582bf6240SBarry Smith } 51682bf6240SBarry Smith 517560360afSLisandro Dalcin static PetscErrorCode TSPseudoSetMaxTimeStep_Pseudo(TS ts,PetscReal maxdt) 51886534af1SJed Brown { 51986534af1SJed Brown TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 52086534af1SJed Brown 52186534af1SJed Brown PetscFunctionBegin; 52286534af1SJed Brown pseudo->dt_max = maxdt; 52386534af1SJed Brown PetscFunctionReturn(0); 52486534af1SJed Brown } 52586534af1SJed Brown 526560360afSLisandro Dalcin static PetscErrorCode TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts) 52782bf6240SBarry Smith { 5284bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 52982bf6240SBarry Smith 53082bf6240SBarry Smith PetscFunctionBegin; 5314bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_TRUE; 53282bf6240SBarry Smith PetscFunctionReturn(0); 53382bf6240SBarry Smith } 53482bf6240SBarry Smith 5356849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/ 536560360afSLisandro Dalcin static PetscErrorCode TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void *ctx) 53782bf6240SBarry Smith { 5384bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 53982bf6240SBarry Smith 54082bf6240SBarry Smith PetscFunctionBegin; 54182bf6240SBarry Smith pseudo->dt = dt; 54282bf6240SBarry Smith pseudo->dtctx = ctx; 54382bf6240SBarry Smith PetscFunctionReturn(0); 54482bf6240SBarry Smith } 54582bf6240SBarry Smith 54682bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 54710e6a065SJed Brown /*MC 54810e6a065SJed Brown TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping 54982bf6240SBarry Smith 55010e6a065SJed Brown This method solves equations of the form 55110e6a065SJed Brown 55210e6a065SJed Brown $ F(X,Xdot) = 0 55310e6a065SJed Brown 55410e6a065SJed Brown for steady state using the iteration 55510e6a065SJed Brown 55610e6a065SJed Brown $ [G'] S = -F(X,0) 55710e6a065SJed Brown $ X += S 55810e6a065SJed Brown 55910e6a065SJed Brown where 56010e6a065SJed Brown 56110e6a065SJed Brown $ G(Y) = F(Y,(Y-X)/dt) 56210e6a065SJed Brown 5636f2d6a7bSJed Brown This is linearly-implicit Euler with the residual always evaluated "at steady 5646f2d6a7bSJed Brown state". See note below. 56510e6a065SJed Brown 56610e6a065SJed Brown Options database keys: 56710e6a065SJed Brown + -ts_pseudo_increment <real> - ratio of increase dt 5683118ae5eSBarry Smith . -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt 5693118ae5eSBarry Smith . -ts_pseudo_fatol <atol> - stop iterating when the function norm is less than atol 5703118ae5eSBarry Smith - -ts_pseudo_frtol <rtol> - stop iterating when the function norm divided by the initial function norm is less than rtol 57110e6a065SJed Brown 57210e6a065SJed Brown Level: beginner 57310e6a065SJed Brown 57410e6a065SJed Brown References: 575606c0280SSatish Balay + * - Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential Algebraic Equations, 2003. 576606c0280SSatish Balay - * - C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998. 57710e6a065SJed Brown 57810e6a065SJed Brown Notes: 5796f2d6a7bSJed Brown The residual computed by this method includes the transient term (Xdot is computed instead of 5806f2d6a7bSJed Brown always being zero), but since the prediction from the last step is always the solution from the 5816f2d6a7bSJed Brown last step, on the first Newton iteration we have 5826f2d6a7bSJed Brown 5836f2d6a7bSJed Brown $ Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0 5846f2d6a7bSJed Brown 5856f2d6a7bSJed Brown Therefore, the linear system solved by the first Newton iteration is equivalent to the one 5866f2d6a7bSJed Brown described above and in the papers. If the user chooses to perform multiple Newton iterations, the 5876f2d6a7bSJed Brown algorithm is no longer the one described in the referenced papers. 58810e6a065SJed Brown 589db781477SPatrick Sanan .seealso: `TSCreate()`, `TS`, `TSSetType()` 59010e6a065SJed Brown 59110e6a065SJed Brown M*/ 5928cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Pseudo(TS ts) 5932d3f70b5SBarry Smith { 5947bf11e45SBarry Smith TS_Pseudo *pseudo; 595193ac0bcSJed Brown SNES snes; 59619fd82e9SBarry Smith SNESType stype; 5972d3f70b5SBarry Smith 5983a40ed3dSBarry Smith PetscFunctionBegin; 599277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Pseudo; 600000e7ae3SMatthew Knepley ts->ops->destroy = TSDestroy_Pseudo; 601000e7ae3SMatthew Knepley ts->ops->view = TSView_Pseudo; 602000e7ae3SMatthew Knepley ts->ops->setup = TSSetUp_Pseudo; 603000e7ae3SMatthew Knepley ts->ops->step = TSStep_Pseudo; 604000e7ae3SMatthew Knepley ts->ops->setfromoptions = TSSetFromOptions_Pseudo; 6050f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Pseudo; 6060f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Pseudo; 6072ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 608825ab935SBarry Smith ts->usessnes = PETSC_TRUE; 6097bf11e45SBarry Smith 6109566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts,&snes)); 6119566063dSJacob Faibussowitsch PetscCall(SNESGetType(snes,&stype)); 6129566063dSJacob Faibussowitsch if (!stype) PetscCall(SNESSetType(snes,SNESKSPONLY)); 6132d3f70b5SBarry Smith 6149566063dSJacob Faibussowitsch PetscCall(PetscNewLog(ts,&pseudo)); 6157bf11e45SBarry Smith ts->data = (void*)pseudo; 6162d3f70b5SBarry Smith 617be5899b3SLisandro Dalcin pseudo->dt = TSPseudoTimeStepDefault; 618be5899b3SLisandro Dalcin pseudo->dtctx = NULL; 61928aa8177SBarry Smith pseudo->dt_increment = 1.1; 6204bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_FALSE; 621193ac0bcSJed Brown pseudo->fnorm = -1; 622be5899b3SLisandro Dalcin pseudo->fnorm_initial = -1; 623be5899b3SLisandro Dalcin pseudo->fnorm_previous = -1; 6243118ae5eSBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 6253118ae5eSBarry Smith pseudo->fatol = 1.e-25; 6263118ae5eSBarry Smith pseudo->frtol = 1.e-5; 6273118ae5eSBarry Smith #else 6283118ae5eSBarry Smith pseudo->fatol = 1.e-50; 6293118ae5eSBarry Smith pseudo->frtol = 1.e-12; 6303118ae5eSBarry Smith #endif 6319566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",TSPseudoSetVerifyTimeStep_Pseudo)); 6329566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",TSPseudoSetTimeStepIncrement_Pseudo)); 6339566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetMaxTimeStep_C",TSPseudoSetMaxTimeStep_Pseudo)); 6349566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",TSPseudoIncrementDtFromInitialDt_Pseudo)); 6359566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",TSPseudoSetTimeStep_Pseudo)); 6363a40ed3dSBarry Smith PetscFunctionReturn(0); 6372d3f70b5SBarry Smith } 6382d3f70b5SBarry Smith 63982bf6240SBarry Smith /*@C 6408d359177SBarry Smith TSPseudoTimeStepDefault - Default code to compute pseudo-timestepping. 64182bf6240SBarry Smith Use with TSPseudoSetTimeStep(). 64228aa8177SBarry Smith 64315091d37SBarry Smith Collective on TS 64415091d37SBarry Smith 64528aa8177SBarry Smith Input Parameters: 646a2b725a8SWilliam Gropp + ts - the timestep context 647a2b725a8SWilliam Gropp - dtctx - unused timestep context 64828aa8177SBarry Smith 64982bf6240SBarry Smith Output Parameter: 65082bf6240SBarry Smith . newdt - the timestep to use for the next step 65128aa8177SBarry Smith 65215091d37SBarry Smith Level: advanced 65315091d37SBarry Smith 654db781477SPatrick Sanan .seealso: `TSPseudoSetTimeStep()`, `TSPseudoComputeTimeStep()` 65528aa8177SBarry Smith @*/ 6568d359177SBarry Smith PetscErrorCode TSPseudoTimeStepDefault(TS ts,PetscReal *newdt,void *dtctx) 65728aa8177SBarry Smith { 65882bf6240SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 659be5899b3SLisandro Dalcin PetscReal inc = pseudo->dt_increment; 66028aa8177SBarry Smith 6613a40ed3dSBarry Smith PetscFunctionBegin; 6629566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(pseudo->xdot)); 6639566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE)); 6649566063dSJacob Faibussowitsch PetscCall(VecNorm(pseudo->func,NORM_2,&pseudo->fnorm)); 665be5899b3SLisandro Dalcin if (pseudo->fnorm_initial < 0) { 66682bf6240SBarry Smith /* first time through so compute initial function norm */ 667cdbf8f93SLisandro Dalcin pseudo->fnorm_initial = pseudo->fnorm; 668be5899b3SLisandro Dalcin pseudo->fnorm_previous = pseudo->fnorm; 66982bf6240SBarry Smith } 670bbd56ea5SKarl Rupp if (pseudo->fnorm == 0.0) *newdt = 1.e12*inc*ts->time_step; 671bbd56ea5SKarl Rupp else if (pseudo->increment_dt_from_initial_dt) *newdt = inc*pseudo->dt_initial*pseudo->fnorm_initial/pseudo->fnorm; 672be5899b3SLisandro Dalcin else *newdt = inc*ts->time_step*pseudo->fnorm_previous/pseudo->fnorm; 67386534af1SJed Brown if (pseudo->dt_max > 0) *newdt = PetscMin(*newdt,pseudo->dt_max); 67482bf6240SBarry Smith pseudo->fnorm_previous = pseudo->fnorm; 6753a40ed3dSBarry Smith PetscFunctionReturn(0); 67628aa8177SBarry Smith } 677