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 488d359177SBarry Smith .seealso: TSPseudoTimeStepDefault(), TSPseudoSetTimeStep() 497bf11e45SBarry Smith @*/ 507087cfbeSBarry Smith PetscErrorCode TSPseudoComputeTimeStep(TS ts,PetscReal *dt) 517bf11e45SBarry Smith { 527bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 53dfbe8321SBarry Smith PetscErrorCode ierr; 547bf11e45SBarry Smith 553a40ed3dSBarry Smith PetscFunctionBegin; 56d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr); 577bf11e45SBarry Smith ierr = (*pseudo->dt)(ts,dt,pseudo->dtctx);CHKERRQ(ierr); 58d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr); 593a40ed3dSBarry Smith PetscFunctionReturn(0); 607bf11e45SBarry Smith } 617bf11e45SBarry Smith 627bf11e45SBarry Smith /* ------------------------------------------------------------------------------*/ 637bf11e45SBarry Smith /*@C 648d359177SBarry Smith TSPseudoVerifyTimeStepDefault - Default code to verify the quality of the last timestep. 657bf11e45SBarry Smith 6615091d37SBarry Smith Collective on TS 6715091d37SBarry Smith 687bf11e45SBarry Smith Input Parameters: 6915091d37SBarry Smith + ts - the timestep context 707bf11e45SBarry Smith . dtctx - unused timestep context 7115091d37SBarry Smith - update - latest solution vector 727bf11e45SBarry Smith 73564e8f4eSLois Curfman McInnes Output Parameters: 7415091d37SBarry Smith + newdt - the timestep to use for the next step 7515091d37SBarry Smith - flag - flag indicating whether the last time step was acceptable 767bf11e45SBarry Smith 7715091d37SBarry Smith Level: advanced 78fee21e36SBarry Smith 79564e8f4eSLois Curfman McInnes Note: 80564e8f4eSLois Curfman McInnes This routine always returns a flag of 1, indicating an acceptable 81564e8f4eSLois Curfman McInnes timestep. 82564e8f4eSLois Curfman McInnes 83564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStep() 847bf11e45SBarry Smith @*/ 858d359177SBarry Smith PetscErrorCode TSPseudoVerifyTimeStepDefault(TS ts,Vec update,void *dtctx,PetscReal *newdt,PetscBool *flag) 867bf11e45SBarry Smith { 873a40ed3dSBarry Smith PetscFunctionBegin; 88a7cc72afSBarry Smith *flag = PETSC_TRUE; 893a40ed3dSBarry Smith PetscFunctionReturn(0); 907bf11e45SBarry Smith } 917bf11e45SBarry Smith 927bf11e45SBarry Smith /*@ 93564e8f4eSLois Curfman McInnes TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable. 947bf11e45SBarry Smith 9515091d37SBarry Smith Collective on TS 9615091d37SBarry Smith 97fb4a63b6SLois Curfman McInnes Input Parameters: 9815091d37SBarry Smith + ts - timestep context 9915091d37SBarry Smith - update - latest solution vector 1007bf11e45SBarry Smith 101fb4a63b6SLois Curfman McInnes Output Parameters: 10215091d37SBarry Smith + dt - newly computed timestep (if it had to shrink) 10315091d37SBarry Smith - flag - indicates if current timestep was ok 1047bf11e45SBarry Smith 10515091d37SBarry Smith Level: advanced 106fee21e36SBarry Smith 107564e8f4eSLois Curfman McInnes Notes: 108564e8f4eSLois Curfman McInnes The routine to be called here to compute the timestep should be 109564e8f4eSLois Curfman McInnes set by calling TSPseudoSetVerifyTimeStep(). 110564e8f4eSLois Curfman McInnes 1118d359177SBarry Smith .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStepDefault() 1127bf11e45SBarry Smith @*/ 1137087cfbeSBarry Smith PetscErrorCode TSPseudoVerifyTimeStep(TS ts,Vec update,PetscReal *dt,PetscBool *flag) 1147bf11e45SBarry Smith { 1157bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 116dfbe8321SBarry Smith PetscErrorCode ierr; 1177bf11e45SBarry Smith 1183a40ed3dSBarry Smith PetscFunctionBegin; 119cb9d8021SPierre Barbier de Reuille *flag = PETSC_TRUE; 120be5899b3SLisandro Dalcin if (pseudo->verify) { 1217bf11e45SBarry Smith ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag);CHKERRQ(ierr); 122cb9d8021SPierre Barbier de Reuille } 1233a40ed3dSBarry Smith PetscFunctionReturn(0); 1247bf11e45SBarry Smith } 1257bf11e45SBarry Smith 1267bf11e45SBarry Smith /* --------------------------------------------------------------------------------*/ 1277bf11e45SBarry Smith 128193ac0bcSJed Brown static PetscErrorCode TSStep_Pseudo(TS ts) 1292d3f70b5SBarry Smith { 130277b19d0SLisandro Dalcin TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 131be5899b3SLisandro Dalcin PetscInt nits,lits,reject; 132cdbf8f93SLisandro Dalcin PetscBool stepok; 133be5899b3SLisandro Dalcin PetscReal next_time_step = ts->time_step; 134dfbe8321SBarry Smith PetscErrorCode ierr; 1352d3f70b5SBarry Smith 1363a40ed3dSBarry Smith PetscFunctionBegin; 137bbd56ea5SKarl Rupp if (ts->steps == 0) pseudo->dt_initial = ts->time_step; 138193ac0bcSJed Brown ierr = VecCopy(ts->vec_sol,pseudo->update);CHKERRQ(ierr); 139cdbf8f93SLisandro Dalcin ierr = TSPseudoComputeTimeStep(ts,&next_time_step);CHKERRQ(ierr); 140cdbf8f93SLisandro Dalcin for (reject=0; reject<ts->max_reject; reject++,ts->reject++) { 141cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 142b8123daeSJed Brown ierr = TSPreStage(ts,ts->ptime+ts->time_step);CHKERRQ(ierr); 1430298fd71SBarry Smith ierr = SNESSolve(ts->snes,NULL,pseudo->update);CHKERRQ(ierr); 144be5899b3SLisandro Dalcin ierr = SNESGetIterationNumber(ts->snes,&nits);CHKERRQ(ierr); 145b850b91aSLisandro Dalcin ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 146be5899b3SLisandro Dalcin ts->snes_its += nits; ts->ksp_its += lits; 1479be3e283SDebojyoti Ghosh ierr = TSPostStage(ts,ts->ptime+ts->time_step,0,&(pseudo->update));CHKERRQ(ierr); 148be5899b3SLisandro Dalcin ierr = TSAdaptCheckStage(ts->adapt,ts,ts->ptime+ts->time_step,pseudo->update,&stepok);CHKERRQ(ierr); 149be5899b3SLisandro Dalcin if (!stepok) {next_time_step = ts->time_step; continue;} 150193ac0bcSJed Brown pseudo->fnorm = -1; /* The current norm is no longer valid, monitor must recompute it. */ 151cdbf8f93SLisandro Dalcin ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&next_time_step,&stepok);CHKERRQ(ierr); 152cdbf8f93SLisandro Dalcin if (stepok) break; 153cdbf8f93SLisandro Dalcin } 154be5899b3SLisandro Dalcin if (reject >= ts->max_reject) { 155be5899b3SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 1567d3de750SJacob Faibussowitsch ierr = PetscInfo(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,reject);CHKERRQ(ierr); 157cdbf8f93SLisandro Dalcin PetscFunctionReturn(0); 1587bf11e45SBarry Smith } 159be5899b3SLisandro Dalcin 160be5899b3SLisandro Dalcin ierr = VecCopy(pseudo->update,ts->vec_sol);CHKERRQ(ierr); 161be5899b3SLisandro Dalcin ts->ptime += ts->time_step; 162be5899b3SLisandro Dalcin ts->time_step = next_time_step; 163be5899b3SLisandro Dalcin 1643118ae5eSBarry Smith if (pseudo->fnorm < 0) { 1653118ae5eSBarry Smith ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr); 1663118ae5eSBarry Smith ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr); 1673118ae5eSBarry Smith ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 1683118ae5eSBarry Smith } 1693118ae5eSBarry Smith if (pseudo->fnorm < pseudo->fatol) { 1703118ae5eSBarry Smith ts->reason = TS_CONVERGED_PSEUDO_FATOL; 1717d3de750SJacob Faibussowitsch ierr = PetscInfo(ts,"Step=%D, converged since fnorm %g < fatol %g\n",ts->steps,pseudo->fnorm,pseudo->frtol);CHKERRQ(ierr); 1723118ae5eSBarry Smith PetscFunctionReturn(0); 1733118ae5eSBarry Smith } 1743118ae5eSBarry Smith if (pseudo->fnorm/pseudo->fnorm_initial < pseudo->frtol) { 1753118ae5eSBarry Smith ts->reason = TS_CONVERGED_PSEUDO_FRTOL; 1767d3de750SJacob Faibussowitsch ierr = PetscInfo(ts,"Step=%D, converged since fnorm %g / fnorm_initial %g < frtol %g\n",ts->steps,pseudo->fnorm,pseudo->fnorm_initial,pseudo->fatol);CHKERRQ(ierr); 1773118ae5eSBarry Smith PetscFunctionReturn(0); 1783118ae5eSBarry Smith } 1793a40ed3dSBarry Smith PetscFunctionReturn(0); 1802d3f70b5SBarry Smith } 1812d3f70b5SBarry Smith 1822d3f70b5SBarry Smith /*------------------------------------------------------------*/ 183277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Pseudo(TS ts) 1842d3f70b5SBarry Smith { 1857bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 186dfbe8321SBarry Smith PetscErrorCode ierr; 1872d3f70b5SBarry Smith 1883a40ed3dSBarry Smith PetscFunctionBegin; 1896bf464f9SBarry Smith ierr = VecDestroy(&pseudo->update);CHKERRQ(ierr); 1906bf464f9SBarry Smith ierr = VecDestroy(&pseudo->func);CHKERRQ(ierr); 1916bf464f9SBarry Smith ierr = VecDestroy(&pseudo->xdot);CHKERRQ(ierr); 1923a40ed3dSBarry Smith PetscFunctionReturn(0); 1932d3f70b5SBarry Smith } 1942d3f70b5SBarry Smith 195277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Pseudo(TS ts) 196277b19d0SLisandro Dalcin { 197277b19d0SLisandro Dalcin PetscErrorCode ierr; 198277b19d0SLisandro Dalcin 199277b19d0SLisandro Dalcin PetscFunctionBegin; 200277b19d0SLisandro Dalcin ierr = TSReset_Pseudo(ts);CHKERRQ(ierr); 201277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 202bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",NULL);CHKERRQ(ierr); 203bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",NULL);CHKERRQ(ierr); 204bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetMaxTimeStep_C",NULL);CHKERRQ(ierr); 205bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",NULL);CHKERRQ(ierr); 206bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",NULL);CHKERRQ(ierr); 207277b19d0SLisandro Dalcin PetscFunctionReturn(0); 208277b19d0SLisandro Dalcin } 2092d3f70b5SBarry Smith 2102d3f70b5SBarry Smith /*------------------------------------------------------------*/ 2112d3f70b5SBarry Smith 2126f2d6a7bSJed Brown /* 2136f2d6a7bSJed Brown Compute Xdot = (X^{n+1}-X^n)/dt) = 0 2146f2d6a7bSJed Brown */ 2156f2d6a7bSJed Brown static PetscErrorCode TSPseudoGetXdot(TS ts,Vec X,Vec *Xdot) 2162d3f70b5SBarry Smith { 2176f2d6a7bSJed Brown TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 218193ac0bcSJed Brown const PetscScalar mdt = 1.0/ts->time_step,*xnp1,*xn; 219193ac0bcSJed Brown PetscScalar *xdot; 220dfbe8321SBarry Smith PetscErrorCode ierr; 221a7cc72afSBarry Smith PetscInt i,n; 2222d3f70b5SBarry Smith 2233a40ed3dSBarry Smith PetscFunctionBegin; 224aab5bcd8SJed Brown *Xdot = NULL; 225193ac0bcSJed Brown ierr = VecGetArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr); 226193ac0bcSJed Brown ierr = VecGetArrayRead(X,&xnp1);CHKERRQ(ierr); 2276f2d6a7bSJed Brown ierr = VecGetArray(pseudo->xdot,&xdot);CHKERRQ(ierr); 2286f2d6a7bSJed Brown ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr); 229bbd56ea5SKarl Rupp for (i=0; i<n; i++) xdot[i] = mdt*(xnp1[i] - xn[i]); 230193ac0bcSJed Brown ierr = VecRestoreArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr); 231193ac0bcSJed Brown ierr = VecRestoreArrayRead(X,&xnp1);CHKERRQ(ierr); 2326f2d6a7bSJed Brown ierr = VecRestoreArray(pseudo->xdot,&xdot);CHKERRQ(ierr); 2336f2d6a7bSJed Brown *Xdot = pseudo->xdot; 2343a40ed3dSBarry Smith PetscFunctionReturn(0); 2352d3f70b5SBarry Smith } 2362d3f70b5SBarry Smith 2376f2d6a7bSJed Brown /* 2386f2d6a7bSJed Brown The transient residual is 2396f2d6a7bSJed Brown 2406f2d6a7bSJed Brown F(U^{n+1},(U^{n+1}-U^n)/dt) = 0 2416f2d6a7bSJed Brown 2426f2d6a7bSJed Brown or for ODE, 2436f2d6a7bSJed Brown 2446f2d6a7bSJed Brown (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0 2456f2d6a7bSJed Brown 2466f2d6a7bSJed Brown This is the function that must be evaluated for transient simulation and for 2476f2d6a7bSJed Brown finite difference Jacobians. On the first Newton step, this algorithm uses 2486f2d6a7bSJed Brown a guess of U^{n+1} = U^n in which case the transient term vanishes and the 2496f2d6a7bSJed Brown residual is actually the steady state residual. Pseudotransient 2506f2d6a7bSJed Brown continuation as described in the literature is a linearly implicit 2516f2d6a7bSJed Brown algorithm, it only takes this one Newton step with the steady state 2526f2d6a7bSJed Brown residual, and then advances to the next time step. 2536f2d6a7bSJed Brown */ 2540f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes,Vec X,Vec Y,TS ts) 2552d3f70b5SBarry Smith { 2566f2d6a7bSJed Brown Vec Xdot; 257dfbe8321SBarry Smith PetscErrorCode ierr; 2582d3f70b5SBarry Smith 2593a40ed3dSBarry Smith PetscFunctionBegin; 2606f2d6a7bSJed Brown ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr); 261193ac0bcSJed Brown ierr = TSComputeIFunction(ts,ts->ptime+ts->time_step,X,Xdot,Y,PETSC_FALSE);CHKERRQ(ierr); 2626f2d6a7bSJed Brown PetscFunctionReturn(0); 2636f2d6a7bSJed Brown } 2642d3f70b5SBarry Smith 2656f2d6a7bSJed Brown /* 2666f2d6a7bSJed Brown This constructs the Jacobian needed for SNES. For DAE, this is 2676f2d6a7bSJed Brown 2686f2d6a7bSJed Brown dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot 2696f2d6a7bSJed Brown 2706f2d6a7bSJed Brown and for ODE: 2716f2d6a7bSJed Brown 2726f2d6a7bSJed Brown J = I/dt - J_{Frhs} where J_{Frhs} is the given Jacobian of Frhs. 2736f2d6a7bSJed Brown */ 274d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes,Vec X,Mat AA,Mat BB,TS ts) 2756f2d6a7bSJed Brown { 2766f2d6a7bSJed Brown Vec Xdot; 2776f2d6a7bSJed Brown PetscErrorCode ierr; 2786f2d6a7bSJed Brown 2796f2d6a7bSJed Brown PetscFunctionBegin; 2806f2d6a7bSJed Brown ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr); 281d1e9a80fSBarry Smith ierr = TSComputeIJacobian(ts,ts->ptime+ts->time_step,X,Xdot,1./ts->time_step,AA,BB,PETSC_FALSE);CHKERRQ(ierr); 2823a40ed3dSBarry Smith PetscFunctionReturn(0); 2832d3f70b5SBarry Smith } 2842d3f70b5SBarry Smith 2856849ba73SBarry Smith static PetscErrorCode TSSetUp_Pseudo(TS ts) 2862d3f70b5SBarry Smith { 2877bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 288dfbe8321SBarry Smith PetscErrorCode ierr; 2892d3f70b5SBarry Smith 2903a40ed3dSBarry Smith PetscFunctionBegin; 2917bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr); 2927bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr); 2936f2d6a7bSJed Brown ierr = VecDuplicate(ts->vec_sol,&pseudo->xdot);CHKERRQ(ierr); 2943a40ed3dSBarry Smith PetscFunctionReturn(0); 2952d3f70b5SBarry Smith } 2962d3f70b5SBarry Smith /*------------------------------------------------------------*/ 2972d3f70b5SBarry Smith 298560360afSLisandro Dalcin static PetscErrorCode TSPseudoMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy) 2992d3f70b5SBarry Smith { 3007bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 301dfbe8321SBarry Smith PetscErrorCode ierr; 302ce94432eSBarry Smith PetscViewer viewer = (PetscViewer) dummy; 3032d3f70b5SBarry Smith 3043a40ed3dSBarry Smith PetscFunctionBegin; 305193ac0bcSJed Brown if (pseudo->fnorm < 0) { /* The last computed norm is stale, recompute */ 306193ac0bcSJed Brown ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr); 307193ac0bcSJed Brown ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr); 308193ac0bcSJed Brown ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 309193ac0bcSJed Brown } 310649052a6SBarry Smith ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 3117c8652ddSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"TS %D dt %g time %g fnorm %g\n",step,(double)ts->time_step,(double)ptime,(double)pseudo->fnorm);CHKERRQ(ierr); 312649052a6SBarry Smith ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 3133a40ed3dSBarry Smith PetscFunctionReturn(0); 3142d3f70b5SBarry Smith } 3152d3f70b5SBarry Smith 3164416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Pseudo(PetscOptionItems *PetscOptionsObject,TS ts) 3172d3f70b5SBarry Smith { 3184bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 319dfbe8321SBarry Smith PetscErrorCode ierr; 320ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE; 321649052a6SBarry Smith PetscViewer viewer; 3222d3f70b5SBarry Smith 3233a40ed3dSBarry Smith PetscFunctionBegin; 324e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Pseudo-timestepping options");CHKERRQ(ierr); 325560360afSLisandro Dalcin ierr = PetscOptionsBool("-ts_monitor_pseudo","Monitor convergence","",flg,&flg,NULL);CHKERRQ(ierr); 3262d3f70b5SBarry Smith if (flg) { 327ce94432eSBarry Smith ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ts),"stdout",&viewer);CHKERRQ(ierr); 328649052a6SBarry Smith ierr = TSMonitorSet(ts,TSPseudoMonitorDefault,viewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 32928aa8177SBarry Smith } 330be5899b3SLisandro Dalcin flg = pseudo->increment_dt_from_initial_dt; 3310298fd71SBarry Smith ierr = PetscOptionsBool("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",flg,&flg,NULL);CHKERRQ(ierr); 332be5899b3SLisandro Dalcin pseudo->increment_dt_from_initial_dt = flg; 33394ae4db5SBarry Smith ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,NULL);CHKERRQ(ierr); 33494ae4db5SBarry Smith ierr = PetscOptionsReal("-ts_pseudo_max_dt","Maximum value for dt","TSPseudoSetMaxTimeStep",pseudo->dt_max,&pseudo->dt_max,NULL);CHKERRQ(ierr); 3353118ae5eSBarry Smith ierr = PetscOptionsReal("-ts_pseudo_fatol","Tolerance for norm of function","",pseudo->fatol,&pseudo->fatol,NULL);CHKERRQ(ierr); 3363118ae5eSBarry Smith ierr = PetscOptionsReal("-ts_pseudo_frtol","Relative tolerance for norm of function","",pseudo->frtol,&pseudo->frtol,NULL);CHKERRQ(ierr); 337b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 3383a40ed3dSBarry Smith PetscFunctionReturn(0); 3392d3f70b5SBarry Smith } 3402d3f70b5SBarry Smith 3416849ba73SBarry Smith static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer) 3422d3f70b5SBarry Smith { 343d52bd9f3SBarry Smith PetscErrorCode ierr; 3443118ae5eSBarry Smith PetscBool isascii; 345d52bd9f3SBarry Smith 3463a40ed3dSBarry Smith PetscFunctionBegin; 3473118ae5eSBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 3483118ae5eSBarry Smith if (isascii) { 3493118ae5eSBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 3503118ae5eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," frtol - relative tolerance in function value %g\n",(double)pseudo->frtol);CHKERRQ(ierr); 3513118ae5eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," fatol - absolute tolerance in function value %g\n",(double)pseudo->fatol);CHKERRQ(ierr); 3523118ae5eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," dt_initial - initial timestep %g\n",(double)pseudo->dt_initial);CHKERRQ(ierr); 3533118ae5eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," dt_increment - increase in timestep on successful step %g\n",(double)pseudo->dt_increment);CHKERRQ(ierr); 3543118ae5eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," dt_max - maximum time %g\n",(double)pseudo->dt_max);CHKERRQ(ierr); 3553118ae5eSBarry Smith } 3563a40ed3dSBarry Smith PetscFunctionReturn(0); 3572d3f70b5SBarry Smith } 3582d3f70b5SBarry Smith 35982bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 360ac226902SBarry Smith /*@C 36182bf6240SBarry Smith TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the 36282bf6240SBarry Smith last timestep. 36382bf6240SBarry Smith 3643f9fe445SBarry Smith Logically Collective on TS 36515091d37SBarry Smith 36682bf6240SBarry Smith Input Parameters: 36715091d37SBarry Smith + ts - timestep context 36882bf6240SBarry Smith . dt - user-defined function to verify timestep 36915091d37SBarry Smith - ctx - [optional] user-defined context for private data 3700298fd71SBarry Smith for the timestep verification routine (may be NULL) 37182bf6240SBarry Smith 37215091d37SBarry Smith Level: advanced 373fee21e36SBarry Smith 37482bf6240SBarry Smith Calling sequence of func: 375a2b725a8SWilliam Gropp $ func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscBool *flag); 37682bf6240SBarry Smith 377a2b725a8SWilliam Gropp + update - latest solution vector 37882bf6240SBarry Smith . ctx - [optional] timestep context 37982bf6240SBarry Smith . newdt - the timestep to use for the next step 380a2b725a8SWilliam Gropp - flag - flag indicating whether the last time step was acceptable 38182bf6240SBarry Smith 38282bf6240SBarry Smith Notes: 38382bf6240SBarry Smith The routine set here will be called by TSPseudoVerifyTimeStep() 38482bf6240SBarry Smith during the timestepping process. 38582bf6240SBarry Smith 3868d359177SBarry Smith .seealso: TSPseudoVerifyTimeStepDefault(), TSPseudoVerifyTimeStep() 38782bf6240SBarry Smith @*/ 3887087cfbeSBarry Smith PetscErrorCode TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,PetscBool*),void *ctx) 38982bf6240SBarry Smith { 3904ac538c5SBarry Smith PetscErrorCode ierr; 39182bf6240SBarry Smith 39282bf6240SBarry Smith PetscFunctionBegin; 3930700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3944ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSPseudoSetVerifyTimeStep_C",(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal*,PetscBool*),void*),(ts,dt,ctx));CHKERRQ(ierr); 39582bf6240SBarry Smith PetscFunctionReturn(0); 39682bf6240SBarry Smith } 39782bf6240SBarry Smith 39882bf6240SBarry Smith /*@ 39982bf6240SBarry Smith TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 4008d359177SBarry Smith dt when using the TSPseudoTimeStepDefault() routine. 40182bf6240SBarry Smith 4023f9fe445SBarry Smith Logically Collective on TS 403fee21e36SBarry Smith 40415091d37SBarry Smith Input Parameters: 40515091d37SBarry Smith + ts - the timestep context 40615091d37SBarry Smith - inc - the scaling factor >= 1.0 40715091d37SBarry Smith 40882bf6240SBarry Smith Options Database Key: 409e1bc860dSBarry Smith . -ts_pseudo_increment <increment> 41082bf6240SBarry Smith 41115091d37SBarry Smith Level: advanced 41215091d37SBarry Smith 4138d359177SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault() 41482bf6240SBarry Smith @*/ 4157087cfbeSBarry Smith PetscErrorCode TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc) 41682bf6240SBarry Smith { 4174ac538c5SBarry Smith PetscErrorCode ierr; 41882bf6240SBarry Smith 41982bf6240SBarry Smith PetscFunctionBegin; 4200700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 421c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts,inc,2); 4224ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSPseudoSetTimeStepIncrement_C",(TS,PetscReal),(ts,inc));CHKERRQ(ierr); 42382bf6240SBarry Smith PetscFunctionReturn(0); 42482bf6240SBarry Smith } 42582bf6240SBarry Smith 42686534af1SJed Brown /*@ 42786534af1SJed Brown TSPseudoSetMaxTimeStep - Sets the maximum time step 4288d359177SBarry Smith when using the TSPseudoTimeStepDefault() routine. 42986534af1SJed Brown 43086534af1SJed Brown Logically Collective on TS 43186534af1SJed Brown 43286534af1SJed Brown Input Parameters: 43386534af1SJed Brown + ts - the timestep context 43486534af1SJed Brown - maxdt - the maximum time step, use a non-positive value to deactivate 43586534af1SJed Brown 43686534af1SJed Brown Options Database Key: 437e1bc860dSBarry Smith . -ts_pseudo_max_dt <increment> 43886534af1SJed Brown 43986534af1SJed Brown Level: advanced 44086534af1SJed Brown 4418d359177SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault() 44286534af1SJed Brown @*/ 44386534af1SJed Brown PetscErrorCode TSPseudoSetMaxTimeStep(TS ts,PetscReal maxdt) 44486534af1SJed Brown { 44586534af1SJed Brown PetscErrorCode ierr; 44686534af1SJed Brown 44786534af1SJed Brown PetscFunctionBegin; 44886534af1SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 44986534af1SJed Brown PetscValidLogicalCollectiveReal(ts,maxdt,2); 45086534af1SJed Brown ierr = PetscTryMethod(ts,"TSPseudoSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));CHKERRQ(ierr); 45186534af1SJed Brown PetscFunctionReturn(0); 45286534af1SJed Brown } 45386534af1SJed Brown 45482bf6240SBarry Smith /*@ 45582bf6240SBarry Smith TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep 45682bf6240SBarry Smith is computed via the formula 45782bf6240SBarry Smith $ dt = initial_dt*initial_fnorm/current_fnorm 45882bf6240SBarry Smith rather than the default update, 45982bf6240SBarry Smith $ dt = current_dt*previous_fnorm/current_fnorm. 46082bf6240SBarry Smith 4613f9fe445SBarry Smith Logically Collective on TS 46215091d37SBarry Smith 46382bf6240SBarry Smith Input Parameter: 46482bf6240SBarry Smith . ts - the timestep context 46582bf6240SBarry Smith 46682bf6240SBarry Smith Options Database Key: 467e1bc860dSBarry Smith . -ts_pseudo_increment_dt_from_initial_dt 46882bf6240SBarry Smith 46915091d37SBarry Smith Level: advanced 47015091d37SBarry Smith 4718d359177SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault() 47282bf6240SBarry Smith @*/ 4737087cfbeSBarry Smith PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS ts) 47482bf6240SBarry Smith { 4754ac538c5SBarry Smith PetscErrorCode ierr; 47682bf6240SBarry Smith 47782bf6240SBarry Smith PetscFunctionBegin; 4780700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 4794ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSPseudoIncrementDtFromInitialDt_C",(TS),(ts));CHKERRQ(ierr); 48082bf6240SBarry Smith PetscFunctionReturn(0); 48182bf6240SBarry Smith } 48282bf6240SBarry Smith 483ac226902SBarry Smith /*@C 48482bf6240SBarry Smith TSPseudoSetTimeStep - Sets the user-defined routine to be 48582bf6240SBarry Smith called at each pseudo-timestep to update the timestep. 48682bf6240SBarry Smith 4873f9fe445SBarry Smith Logically Collective on TS 48815091d37SBarry Smith 48982bf6240SBarry Smith Input Parameters: 49015091d37SBarry Smith + ts - timestep context 49182bf6240SBarry Smith . dt - function to compute timestep 49215091d37SBarry Smith - ctx - [optional] user-defined context for private data 4930298fd71SBarry Smith required by the function (may be NULL) 49482bf6240SBarry Smith 49515091d37SBarry Smith Level: intermediate 496fee21e36SBarry Smith 49782bf6240SBarry Smith Calling sequence of func: 498a2b725a8SWilliam Gropp $ func (TS ts,PetscReal *newdt,void *ctx); 49982bf6240SBarry Smith 500a2b725a8SWilliam Gropp + newdt - the newly computed timestep 501a2b725a8SWilliam Gropp - ctx - [optional] timestep context 50282bf6240SBarry Smith 50382bf6240SBarry Smith Notes: 50482bf6240SBarry Smith The routine set here will be called by TSPseudoComputeTimeStep() 50582bf6240SBarry Smith during the timestepping process. 5068d359177SBarry Smith If not set then TSPseudoTimeStepDefault() is automatically used 50782bf6240SBarry Smith 5088d359177SBarry Smith .seealso: TSPseudoTimeStepDefault(), TSPseudoComputeTimeStep() 50982bf6240SBarry Smith @*/ 5107087cfbeSBarry Smith PetscErrorCode TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void *ctx) 51182bf6240SBarry Smith { 5124ac538c5SBarry Smith PetscErrorCode ierr; 51382bf6240SBarry Smith 51482bf6240SBarry Smith PetscFunctionBegin; 5150700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 5164ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSPseudoSetTimeStep_C",(TS,PetscErrorCode (*)(TS,PetscReal*,void*),void*),(ts,dt,ctx));CHKERRQ(ierr); 51782bf6240SBarry Smith PetscFunctionReturn(0); 51882bf6240SBarry Smith } 51982bf6240SBarry Smith 52082bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 52182bf6240SBarry Smith 522ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscBool*); /* force argument to next function to not be extern C*/ 523560360afSLisandro Dalcin static PetscErrorCode TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void *ctx) 52482bf6240SBarry Smith { 525be5899b3SLisandro Dalcin TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 52682bf6240SBarry Smith 52782bf6240SBarry Smith PetscFunctionBegin; 52882bf6240SBarry Smith pseudo->verify = dt; 52982bf6240SBarry Smith pseudo->verifyctx = ctx; 53082bf6240SBarry Smith PetscFunctionReturn(0); 53182bf6240SBarry Smith } 53282bf6240SBarry Smith 533560360afSLisandro Dalcin static PetscErrorCode TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc) 53482bf6240SBarry Smith { 5354bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 53682bf6240SBarry Smith 53782bf6240SBarry Smith PetscFunctionBegin; 53882bf6240SBarry Smith pseudo->dt_increment = inc; 53982bf6240SBarry Smith PetscFunctionReturn(0); 54082bf6240SBarry Smith } 54182bf6240SBarry Smith 542560360afSLisandro Dalcin static PetscErrorCode TSPseudoSetMaxTimeStep_Pseudo(TS ts,PetscReal maxdt) 54386534af1SJed Brown { 54486534af1SJed Brown TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 54586534af1SJed Brown 54686534af1SJed Brown PetscFunctionBegin; 54786534af1SJed Brown pseudo->dt_max = maxdt; 54886534af1SJed Brown PetscFunctionReturn(0); 54986534af1SJed Brown } 55086534af1SJed Brown 551560360afSLisandro Dalcin static PetscErrorCode TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts) 55282bf6240SBarry Smith { 5534bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 55482bf6240SBarry Smith 55582bf6240SBarry Smith PetscFunctionBegin; 5564bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_TRUE; 55782bf6240SBarry Smith PetscFunctionReturn(0); 55882bf6240SBarry Smith } 55982bf6240SBarry Smith 5606849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/ 561560360afSLisandro Dalcin static PetscErrorCode TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void *ctx) 56282bf6240SBarry Smith { 5634bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 56482bf6240SBarry Smith 56582bf6240SBarry Smith PetscFunctionBegin; 56682bf6240SBarry Smith pseudo->dt = dt; 56782bf6240SBarry Smith pseudo->dtctx = ctx; 56882bf6240SBarry Smith PetscFunctionReturn(0); 56982bf6240SBarry Smith } 57082bf6240SBarry Smith 57182bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 57210e6a065SJed Brown /*MC 57310e6a065SJed Brown TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping 57482bf6240SBarry Smith 57510e6a065SJed Brown This method solves equations of the form 57610e6a065SJed Brown 57710e6a065SJed Brown $ F(X,Xdot) = 0 57810e6a065SJed Brown 57910e6a065SJed Brown for steady state using the iteration 58010e6a065SJed Brown 58110e6a065SJed Brown $ [G'] S = -F(X,0) 58210e6a065SJed Brown $ X += S 58310e6a065SJed Brown 58410e6a065SJed Brown where 58510e6a065SJed Brown 58610e6a065SJed Brown $ G(Y) = F(Y,(Y-X)/dt) 58710e6a065SJed Brown 5886f2d6a7bSJed Brown This is linearly-implicit Euler with the residual always evaluated "at steady 5896f2d6a7bSJed Brown state". See note below. 59010e6a065SJed Brown 59110e6a065SJed Brown Options database keys: 59210e6a065SJed Brown + -ts_pseudo_increment <real> - ratio of increase dt 5933118ae5eSBarry Smith . -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt 5943118ae5eSBarry Smith . -ts_pseudo_fatol <atol> - stop iterating when the function norm is less than atol 5953118ae5eSBarry Smith - -ts_pseudo_frtol <rtol> - stop iterating when the function norm divided by the initial function norm is less than rtol 59610e6a065SJed Brown 59710e6a065SJed Brown Level: beginner 59810e6a065SJed Brown 59910e6a065SJed Brown References: 600*606c0280SSatish Balay + * - Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential Algebraic Equations, 2003. 601*606c0280SSatish Balay - * - C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998. 60210e6a065SJed Brown 60310e6a065SJed Brown Notes: 6046f2d6a7bSJed Brown The residual computed by this method includes the transient term (Xdot is computed instead of 6056f2d6a7bSJed Brown always being zero), but since the prediction from the last step is always the solution from the 6066f2d6a7bSJed Brown last step, on the first Newton iteration we have 6076f2d6a7bSJed Brown 6086f2d6a7bSJed Brown $ Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0 6096f2d6a7bSJed Brown 6106f2d6a7bSJed Brown Therefore, the linear system solved by the first Newton iteration is equivalent to the one 6116f2d6a7bSJed Brown described above and in the papers. If the user chooses to perform multiple Newton iterations, the 6126f2d6a7bSJed Brown algorithm is no longer the one described in the referenced papers. 61310e6a065SJed Brown 61410e6a065SJed Brown .seealso: TSCreate(), TS, TSSetType() 61510e6a065SJed Brown 61610e6a065SJed Brown M*/ 6178cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Pseudo(TS ts) 6182d3f70b5SBarry Smith { 6197bf11e45SBarry Smith TS_Pseudo *pseudo; 620dfbe8321SBarry Smith PetscErrorCode ierr; 621193ac0bcSJed Brown SNES snes; 62219fd82e9SBarry Smith SNESType stype; 6232d3f70b5SBarry Smith 6243a40ed3dSBarry Smith PetscFunctionBegin; 625277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Pseudo; 626000e7ae3SMatthew Knepley ts->ops->destroy = TSDestroy_Pseudo; 627000e7ae3SMatthew Knepley ts->ops->view = TSView_Pseudo; 628000e7ae3SMatthew Knepley ts->ops->setup = TSSetUp_Pseudo; 629000e7ae3SMatthew Knepley ts->ops->step = TSStep_Pseudo; 630000e7ae3SMatthew Knepley ts->ops->setfromoptions = TSSetFromOptions_Pseudo; 6310f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Pseudo; 6320f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Pseudo; 6332ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 634825ab935SBarry Smith ts->usessnes = PETSC_TRUE; 6357bf11e45SBarry Smith 636193ac0bcSJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 637193ac0bcSJed Brown ierr = SNESGetType(snes,&stype);CHKERRQ(ierr); 638193ac0bcSJed Brown if (!stype) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);} 6392d3f70b5SBarry Smith 640b00a9115SJed Brown ierr = PetscNewLog(ts,&pseudo);CHKERRQ(ierr); 6417bf11e45SBarry Smith ts->data = (void*)pseudo; 6422d3f70b5SBarry Smith 643be5899b3SLisandro Dalcin pseudo->dt = TSPseudoTimeStepDefault; 644be5899b3SLisandro Dalcin pseudo->dtctx = NULL; 64528aa8177SBarry Smith pseudo->dt_increment = 1.1; 6464bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_FALSE; 647193ac0bcSJed Brown pseudo->fnorm = -1; 648be5899b3SLisandro Dalcin pseudo->fnorm_initial = -1; 649be5899b3SLisandro Dalcin pseudo->fnorm_previous = -1; 6503118ae5eSBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 6513118ae5eSBarry Smith pseudo->fatol = 1.e-25; 6523118ae5eSBarry Smith pseudo->frtol = 1.e-5; 6533118ae5eSBarry Smith #else 6543118ae5eSBarry Smith pseudo->fatol = 1.e-50; 6553118ae5eSBarry Smith pseudo->frtol = 1.e-12; 6563118ae5eSBarry Smith #endif 657bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr); 658bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr); 659bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetMaxTimeStep_C",TSPseudoSetMaxTimeStep_Pseudo);CHKERRQ(ierr); 660bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr); 661bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr); 6623a40ed3dSBarry Smith PetscFunctionReturn(0); 6632d3f70b5SBarry Smith } 6642d3f70b5SBarry Smith 66582bf6240SBarry Smith /*@C 6668d359177SBarry Smith TSPseudoTimeStepDefault - Default code to compute pseudo-timestepping. 66782bf6240SBarry Smith Use with TSPseudoSetTimeStep(). 66828aa8177SBarry Smith 66915091d37SBarry Smith Collective on TS 67015091d37SBarry Smith 67128aa8177SBarry Smith Input Parameters: 672a2b725a8SWilliam Gropp + ts - the timestep context 673a2b725a8SWilliam Gropp - dtctx - unused timestep context 67428aa8177SBarry Smith 67582bf6240SBarry Smith Output Parameter: 67682bf6240SBarry Smith . newdt - the timestep to use for the next step 67728aa8177SBarry Smith 67815091d37SBarry Smith Level: advanced 67915091d37SBarry Smith 68082bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep() 68128aa8177SBarry Smith @*/ 6828d359177SBarry Smith PetscErrorCode TSPseudoTimeStepDefault(TS ts,PetscReal *newdt,void *dtctx) 68328aa8177SBarry Smith { 68482bf6240SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 685be5899b3SLisandro Dalcin PetscReal inc = pseudo->dt_increment; 686dfbe8321SBarry Smith PetscErrorCode ierr; 68728aa8177SBarry Smith 6883a40ed3dSBarry Smith PetscFunctionBegin; 689bbd7b040SJed Brown ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr); 690214bc6a2SJed Brown ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr); 69182bf6240SBarry Smith ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 692be5899b3SLisandro Dalcin if (pseudo->fnorm_initial < 0) { 69382bf6240SBarry Smith /* first time through so compute initial function norm */ 694cdbf8f93SLisandro Dalcin pseudo->fnorm_initial = pseudo->fnorm; 695be5899b3SLisandro Dalcin pseudo->fnorm_previous = pseudo->fnorm; 69682bf6240SBarry Smith } 697bbd56ea5SKarl Rupp if (pseudo->fnorm == 0.0) *newdt = 1.e12*inc*ts->time_step; 698bbd56ea5SKarl Rupp else if (pseudo->increment_dt_from_initial_dt) *newdt = inc*pseudo->dt_initial*pseudo->fnorm_initial/pseudo->fnorm; 699be5899b3SLisandro Dalcin else *newdt = inc*ts->time_step*pseudo->fnorm_previous/pseudo->fnorm; 70086534af1SJed Brown if (pseudo->dt_max > 0) *newdt = PetscMin(*newdt,pseudo->dt_max); 70182bf6240SBarry Smith pseudo->fnorm_previous = pseudo->fnorm; 7023a40ed3dSBarry Smith PetscFunctionReturn(0); 70328aa8177SBarry Smith } 704