12d3f70b5SBarry Smith /* 2fb4a63b6SLois Curfman McInnes Code for Timestepping with implicit backwards Euler. 32d3f70b5SBarry Smith */ 4c6db04a5SJed Brown #include <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 1887828ca2SBarry Smith PetscReal initial_fnorm,fnorm; /* original and current norm of F(u) */ 1987828ca2SBarry Smith PetscReal fnorm_previous; 2028aa8177SBarry Smith 2187828ca2SBarry Smith PetscReal dt_increment; /* scaling that dt is incremented each time-step */ 22ace3abfcSBarry Smith PetscBool increment_dt_from_initial_dt; 237bf11e45SBarry Smith } TS_Pseudo; 242d3f70b5SBarry Smith 252d3f70b5SBarry Smith /* ------------------------------------------------------------------------------*/ 262d3f70b5SBarry Smith 274a2ae208SSatish Balay #undef __FUNCT__ 284a2ae208SSatish Balay #define __FUNCT__ "TSPseudoComputeTimeStep" 297bf11e45SBarry Smith /*@ 307bf11e45SBarry Smith TSPseudoComputeTimeStep - Computes the next timestep for a currently running 31564e8f4eSLois Curfman McInnes pseudo-timestepping process. 322d3f70b5SBarry Smith 3315091d37SBarry Smith Collective on TS 3415091d37SBarry Smith 357bf11e45SBarry Smith Input Parameter: 367bf11e45SBarry Smith . ts - timestep context 377bf11e45SBarry Smith 387bf11e45SBarry Smith Output Parameter: 39fb4a63b6SLois Curfman McInnes . dt - newly computed timestep 40fb4a63b6SLois Curfman McInnes 4115091d37SBarry Smith Level: advanced 42564e8f4eSLois Curfman McInnes 43564e8f4eSLois Curfman McInnes Notes: 44564e8f4eSLois Curfman McInnes The routine to be called here to compute the timestep should be 45564e8f4eSLois Curfman McInnes set by calling TSPseudoSetTimeStep(). 46564e8f4eSLois Curfman McInnes 47fb4a63b6SLois Curfman McInnes .keywords: timestep, pseudo, compute 48564e8f4eSLois Curfman McInnes 49564e8f4eSLois Curfman McInnes .seealso: TSPseudoDefaultTimeStep(), TSPseudoSetTimeStep() 507bf11e45SBarry Smith @*/ 517087cfbeSBarry Smith PetscErrorCode TSPseudoComputeTimeStep(TS ts,PetscReal *dt) 527bf11e45SBarry Smith { 537bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 54dfbe8321SBarry Smith PetscErrorCode ierr; 557bf11e45SBarry Smith 563a40ed3dSBarry Smith PetscFunctionBegin; 57d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr); 587bf11e45SBarry Smith ierr = (*pseudo->dt)(ts,dt,pseudo->dtctx);CHKERRQ(ierr); 59d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr); 603a40ed3dSBarry Smith PetscFunctionReturn(0); 617bf11e45SBarry Smith } 627bf11e45SBarry Smith 637bf11e45SBarry Smith 647bf11e45SBarry Smith /* ------------------------------------------------------------------------------*/ 654a2ae208SSatish Balay #undef __FUNCT__ 664a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultVerifyTimeStep" 677bf11e45SBarry Smith /*@C 68639f9d9dSBarry Smith TSPseudoDefaultVerifyTimeStep - Default code to verify the quality of the last timestep. 697bf11e45SBarry Smith 7015091d37SBarry Smith Collective on TS 7115091d37SBarry Smith 727bf11e45SBarry Smith Input Parameters: 7315091d37SBarry Smith + ts - the timestep context 747bf11e45SBarry Smith . dtctx - unused timestep context 7515091d37SBarry Smith - update - latest solution vector 767bf11e45SBarry Smith 77564e8f4eSLois Curfman McInnes Output Parameters: 7815091d37SBarry Smith + newdt - the timestep to use for the next step 7915091d37SBarry Smith - flag - flag indicating whether the last time step was acceptable 807bf11e45SBarry Smith 8115091d37SBarry Smith Level: advanced 82fee21e36SBarry Smith 83564e8f4eSLois Curfman McInnes Note: 84564e8f4eSLois Curfman McInnes This routine always returns a flag of 1, indicating an acceptable 85564e8f4eSLois Curfman McInnes timestep. 86564e8f4eSLois Curfman McInnes 87564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, default, verify 88564e8f4eSLois Curfman McInnes 89564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStep() 907bf11e45SBarry Smith @*/ 917087cfbeSBarry Smith PetscErrorCode TSPseudoDefaultVerifyTimeStep(TS ts,Vec update,void *dtctx,PetscReal *newdt,PetscBool *flag) 927bf11e45SBarry Smith { 933a40ed3dSBarry Smith PetscFunctionBegin; 94a7cc72afSBarry Smith *flag = PETSC_TRUE; 953a40ed3dSBarry Smith PetscFunctionReturn(0); 967bf11e45SBarry Smith } 977bf11e45SBarry Smith 987bf11e45SBarry Smith 994a2ae208SSatish Balay #undef __FUNCT__ 1004a2ae208SSatish Balay #define __FUNCT__ "TSPseudoVerifyTimeStep" 1017bf11e45SBarry Smith /*@ 102564e8f4eSLois Curfman McInnes TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable. 1037bf11e45SBarry Smith 10415091d37SBarry Smith Collective on TS 10515091d37SBarry Smith 106fb4a63b6SLois Curfman McInnes Input Parameters: 10715091d37SBarry Smith + ts - timestep context 10815091d37SBarry Smith - update - latest solution vector 1097bf11e45SBarry Smith 110fb4a63b6SLois Curfman McInnes Output Parameters: 11115091d37SBarry Smith + dt - newly computed timestep (if it had to shrink) 11215091d37SBarry Smith - flag - indicates if current timestep was ok 1137bf11e45SBarry Smith 11415091d37SBarry Smith Level: advanced 115fee21e36SBarry Smith 116564e8f4eSLois Curfman McInnes Notes: 117564e8f4eSLois Curfman McInnes The routine to be called here to compute the timestep should be 118564e8f4eSLois Curfman McInnes set by calling TSPseudoSetVerifyTimeStep(). 119564e8f4eSLois Curfman McInnes 120564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, verify 121564e8f4eSLois Curfman McInnes 122564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoDefaultVerifyTimeStep() 1237bf11e45SBarry Smith @*/ 1247087cfbeSBarry Smith PetscErrorCode TSPseudoVerifyTimeStep(TS ts,Vec update,PetscReal *dt,PetscBool *flag) 1257bf11e45SBarry Smith { 1267bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 127dfbe8321SBarry Smith PetscErrorCode ierr; 1287bf11e45SBarry Smith 1293a40ed3dSBarry Smith PetscFunctionBegin; 130a7cc72afSBarry Smith if (!pseudo->verify) {*flag = PETSC_TRUE; PetscFunctionReturn(0);} 1317bf11e45SBarry Smith 1327bf11e45SBarry Smith ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag);CHKERRQ(ierr); 1337bf11e45SBarry Smith 1343a40ed3dSBarry Smith PetscFunctionReturn(0); 1357bf11e45SBarry Smith } 1367bf11e45SBarry Smith 1377bf11e45SBarry Smith /* --------------------------------------------------------------------------------*/ 1387bf11e45SBarry Smith 1394a2ae208SSatish Balay #undef __FUNCT__ 1404a2ae208SSatish Balay #define __FUNCT__ "TSStep_Pseudo" 141a7cc72afSBarry Smith static PetscErrorCode TSStep_Pseudo(TS ts,PetscInt *steps,PetscReal *ptime) 1422d3f70b5SBarry Smith { 143277b19d0SLisandro Dalcin TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 144277b19d0SLisandro Dalcin Vec sol = ts->vec_sol, update = pseudo->update; 145dfbe8321SBarry Smith PetscErrorCode ierr; 146*186e87acSLisandro Dalcin PetscInt i,its,lits; 147ace3abfcSBarry Smith PetscBool ok; 14887828ca2SBarry Smith PetscReal current_time_step; 1492d3f70b5SBarry Smith 1503a40ed3dSBarry Smith PetscFunctionBegin; 1512d3f70b5SBarry Smith *steps = -ts->steps; 152*186e87acSLisandro Dalcin *ptime = ts->ptime; 1532d3f70b5SBarry Smith 154277b19d0SLisandro Dalcin ierr = VecCopy(sol,update);CHKERRQ(ierr); 155*186e87acSLisandro Dalcin for (i=0; i<ts->max_steps && ts->ptime < ts->max_time; i++) { 156*186e87acSLisandro Dalcin 1577bf11e45SBarry Smith ierr = TSPseudoComputeTimeStep(ts,&ts->time_step);CHKERRQ(ierr); 158e6e5fe25SBarry Smith ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 1597bf11e45SBarry Smith current_time_step = ts->time_step; 1603f2090d5SJed Brown ierr = TSPreStep(ts);CHKERRQ(ierr); 161e6e5fe25SBarry Smith while (PETSC_TRUE) { 1627bf11e45SBarry Smith ts->ptime += current_time_step; 163277b19d0SLisandro Dalcin ierr = SNESSolve(ts->snes,PETSC_NULL,update);CHKERRQ(ierr); 164b850b91aSLisandro Dalcin ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 1659f954729SBarry Smith ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr); 166c9780b6fSBarry Smith ts->nonlinear_its += its; ts->linear_its += lits; 167277b19d0SLisandro Dalcin ierr = TSPseudoVerifyTimeStep(ts,update,&ts->time_step,&ok);CHKERRQ(ierr); 1687bf11e45SBarry Smith if (ok) break; 1697bf11e45SBarry Smith ts->ptime -= current_time_step; 1707bf11e45SBarry Smith current_time_step = ts->time_step; 1717bf11e45SBarry Smith } 172277b19d0SLisandro Dalcin ierr = VecCopy(update,sol);CHKERRQ(ierr); 1732d3f70b5SBarry Smith ts->steps++; 1743f2090d5SJed Brown ierr = TSPostStep(ts);CHKERRQ(ierr); 1752d3f70b5SBarry Smith } 176bbd7b040SJed Brown ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr); 177bbd7b040SJed Brown ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func);CHKERRQ(ierr); 178e6e5fe25SBarry Smith ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 179e6e5fe25SBarry Smith ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 1802d3f70b5SBarry Smith 1812d3f70b5SBarry Smith *steps += ts->steps; 182142b95e3SSatish Balay *ptime = ts->ptime; 1833a40ed3dSBarry Smith PetscFunctionReturn(0); 1842d3f70b5SBarry Smith } 1852d3f70b5SBarry Smith 1862d3f70b5SBarry Smith /*------------------------------------------------------------*/ 1874a2ae208SSatish Balay #undef __FUNCT__ 188277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Pseudo" 189277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Pseudo(TS ts) 1902d3f70b5SBarry Smith { 1917bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 192dfbe8321SBarry Smith PetscErrorCode ierr; 1932d3f70b5SBarry Smith 1943a40ed3dSBarry Smith PetscFunctionBegin; 195e4ef1970SSatish Balay if (pseudo->update) {ierr = VecDestroy(pseudo->update);CHKERRQ(ierr);} 1967bf11e45SBarry Smith if (pseudo->func) {ierr = VecDestroy(pseudo->func);CHKERRQ(ierr);} 1976f2d6a7bSJed Brown if (pseudo->xdot) {ierr = VecDestroy(pseudo->xdot);CHKERRQ(ierr);} 1983a40ed3dSBarry Smith PetscFunctionReturn(0); 1992d3f70b5SBarry Smith } 2002d3f70b5SBarry Smith 201277b19d0SLisandro Dalcin #undef __FUNCT__ 202277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Pseudo" 203277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Pseudo(TS ts) 204277b19d0SLisandro Dalcin { 205277b19d0SLisandro Dalcin PetscErrorCode ierr; 206277b19d0SLisandro Dalcin 207277b19d0SLisandro Dalcin PetscFunctionBegin; 208277b19d0SLisandro Dalcin ierr = TSReset_Pseudo(ts);CHKERRQ(ierr); 209277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 210277b19d0SLisandro Dalcin PetscFunctionReturn(0); 211277b19d0SLisandro Dalcin } 2122d3f70b5SBarry Smith 2132d3f70b5SBarry Smith /*------------------------------------------------------------*/ 2142d3f70b5SBarry Smith 2154a2ae208SSatish Balay #undef __FUNCT__ 2166f2d6a7bSJed Brown #define __FUNCT__ "TSPseudoGetXdot" 2176f2d6a7bSJed Brown /* 2186f2d6a7bSJed Brown Compute Xdot = (X^{n+1}-X^n)/dt) = 0 2196f2d6a7bSJed Brown */ 2206f2d6a7bSJed Brown static PetscErrorCode TSPseudoGetXdot(TS ts,Vec X,Vec *Xdot) 2212d3f70b5SBarry Smith { 2226f2d6a7bSJed Brown TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 2236f2d6a7bSJed Brown PetscScalar mdt = 1.0/ts->time_step,*xnp1,*xn,*xdot; 224dfbe8321SBarry Smith PetscErrorCode ierr; 225a7cc72afSBarry Smith PetscInt i,n; 2262d3f70b5SBarry Smith 2273a40ed3dSBarry Smith PetscFunctionBegin; 2286f2d6a7bSJed Brown ierr = VecGetArray(ts->vec_sol,&xn);CHKERRQ(ierr); 2296f2d6a7bSJed Brown ierr = VecGetArray(X,&xnp1);CHKERRQ(ierr); 2306f2d6a7bSJed Brown ierr = VecGetArray(pseudo->xdot,&xdot);CHKERRQ(ierr); 2316f2d6a7bSJed Brown ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr); 2322d3f70b5SBarry Smith for (i=0; i<n; i++) { 2336f2d6a7bSJed Brown xdot[i] = mdt*(xnp1[i] - xn[i]); 2342d3f70b5SBarry Smith } 2356f2d6a7bSJed Brown ierr = VecRestoreArray(ts->vec_sol,&xn);CHKERRQ(ierr); 2366f2d6a7bSJed Brown ierr = VecRestoreArray(X,&xnp1);CHKERRQ(ierr); 2376f2d6a7bSJed Brown ierr = VecRestoreArray(pseudo->xdot,&xdot);CHKERRQ(ierr); 2386f2d6a7bSJed Brown *Xdot = pseudo->xdot; 2393a40ed3dSBarry Smith PetscFunctionReturn(0); 2402d3f70b5SBarry Smith } 2412d3f70b5SBarry Smith 2424a2ae208SSatish Balay #undef __FUNCT__ 2430f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Pseudo" 2446f2d6a7bSJed Brown /* 2456f2d6a7bSJed Brown The transient residual is 2466f2d6a7bSJed Brown 2476f2d6a7bSJed Brown F(U^{n+1},(U^{n+1}-U^n)/dt) = 0 2486f2d6a7bSJed Brown 2496f2d6a7bSJed Brown or for ODE, 2506f2d6a7bSJed Brown 2516f2d6a7bSJed Brown (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0 2526f2d6a7bSJed Brown 2536f2d6a7bSJed Brown This is the function that must be evaluated for transient simulation and for 2546f2d6a7bSJed Brown finite difference Jacobians. On the first Newton step, this algorithm uses 2556f2d6a7bSJed Brown a guess of U^{n+1} = U^n in which case the transient term vanishes and the 2566f2d6a7bSJed Brown residual is actually the steady state residual. Pseudotransient 2576f2d6a7bSJed Brown continuation as described in the literature is a linearly implicit 2586f2d6a7bSJed Brown algorithm, it only takes this one Newton step with the steady state 2596f2d6a7bSJed Brown residual, and then advances to the next time step. 2606f2d6a7bSJed Brown */ 2610f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes,Vec X,Vec Y,TS ts) 2622d3f70b5SBarry Smith { 2636f2d6a7bSJed Brown Vec Xdot; 264dfbe8321SBarry Smith PetscErrorCode ierr; 2652d3f70b5SBarry Smith 2663a40ed3dSBarry Smith PetscFunctionBegin; 2676f2d6a7bSJed Brown ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr); 2686f2d6a7bSJed Brown ierr = TSComputeIFunction(ts,ts->ptime,X,Xdot,Y);CHKERRQ(ierr); 2696f2d6a7bSJed Brown PetscFunctionReturn(0); 2706f2d6a7bSJed Brown } 2712d3f70b5SBarry Smith 2726f2d6a7bSJed Brown #undef __FUNCT__ 2730f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Pseudo" 2746f2d6a7bSJed Brown /* 2756f2d6a7bSJed Brown This constructs the Jacobian needed for SNES. For DAE, this is 2766f2d6a7bSJed Brown 2776f2d6a7bSJed Brown dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot 2786f2d6a7bSJed Brown 2796f2d6a7bSJed Brown and for ODE: 2806f2d6a7bSJed Brown 2816f2d6a7bSJed Brown J = I/dt - J_{Frhs} where J_{Frhs} is the given Jacobian of Frhs. 2826f2d6a7bSJed Brown */ 2830f5c6efeSJed Brown static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes,Vec X,Mat *AA,Mat *BB,MatStructure *str,TS ts) 2846f2d6a7bSJed Brown { 2856f2d6a7bSJed Brown Vec Xdot; 2866f2d6a7bSJed Brown PetscErrorCode ierr; 2876f2d6a7bSJed Brown 2886f2d6a7bSJed Brown PetscFunctionBegin; 2896f2d6a7bSJed Brown ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr); 2906f2d6a7bSJed Brown ierr = TSComputeIJacobian(ts,ts->ptime,X,Xdot,1./ts->time_step,AA,BB,str);CHKERRQ(ierr); 2913a40ed3dSBarry Smith PetscFunctionReturn(0); 2922d3f70b5SBarry Smith } 2932d3f70b5SBarry Smith 2942d3f70b5SBarry Smith 2954a2ae208SSatish Balay #undef __FUNCT__ 2964a2ae208SSatish Balay #define __FUNCT__ "TSSetUp_Pseudo" 2976849ba73SBarry Smith static PetscErrorCode TSSetUp_Pseudo(TS ts) 2982d3f70b5SBarry Smith { 2997bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 300dfbe8321SBarry Smith PetscErrorCode ierr; 3012d3f70b5SBarry Smith 3023a40ed3dSBarry Smith PetscFunctionBegin; 3037bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr); 3047bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr); 3056f2d6a7bSJed Brown ierr = VecDuplicate(ts->vec_sol,&pseudo->xdot);CHKERRQ(ierr); 3060f5c6efeSJed Brown ierr = SNESSetFunction(ts->snes,pseudo->func,SNESTSFormFunction,ts);CHKERRQ(ierr); 30744c54cacSJed Brown /* This is nasty. SNESSetFromOptions() is usually called in TSSetFromOptions(). With -snes_mf_operator, it will 30844c54cacSJed Brown replace A and we don't want to mess with that. With -snes_mf, A and B will be replaced as well as the function and 30944c54cacSJed Brown context. Note that SNESSetFunction() normally has not been called before SNESSetFromOptions(), so when -snes_mf sets 31044c54cacSJed Brown the Jacobian user context to snes->funP, it will actually be NULL. This is not a problem because both snes->funP and 31144c54cacSJed Brown snes->jacP should be the TS. */ 31244c54cacSJed Brown { 31344c54cacSJed Brown Mat A,B; 31444c54cacSJed Brown PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*); 31544c54cacSJed Brown void *ctx; 31644c54cacSJed Brown ierr = SNESGetJacobian(ts->snes,&A,&B,&func,&ctx);CHKERRQ(ierr); 31744c54cacSJed Brown ierr = SNESSetJacobian(ts->snes,A?A:ts->A,B?B:ts->B,func?func:SNESTSFormJacobian,ctx?ctx:ts);CHKERRQ(ierr); 31844c54cacSJed Brown } 3193a40ed3dSBarry Smith PetscFunctionReturn(0); 3202d3f70b5SBarry Smith } 3212d3f70b5SBarry Smith /*------------------------------------------------------------*/ 3222d3f70b5SBarry Smith 3234a2ae208SSatish Balay #undef __FUNCT__ 324a6570f20SBarry Smith #define __FUNCT__ "TSPseudoMonitorDefault" 325a6570f20SBarry Smith PetscErrorCode TSPseudoMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *ctx) 3262d3f70b5SBarry Smith { 3277bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 328dfbe8321SBarry Smith PetscErrorCode ierr; 329a34d58ebSBarry Smith PetscViewerASCIIMonitor viewer = (PetscViewerASCIIMonitor)ctx; 3302d3f70b5SBarry Smith 3313a40ed3dSBarry Smith PetscFunctionBegin; 332a34d58ebSBarry Smith if (!ctx) { 3337adad957SLisandro Dalcin ierr = PetscViewerASCIIMonitorCreate(((PetscObject)ts)->comm,"stdout",0,&viewer);CHKERRQ(ierr); 334a34d58ebSBarry Smith } 335a34d58ebSBarry Smith ierr = PetscViewerASCIIMonitorPrintf(viewer,"TS %D dt %G time %G fnorm %G\n",step,ts->time_step,ptime,pseudo->fnorm);CHKERRQ(ierr); 336a34d58ebSBarry Smith if (!ctx) { 337a34d58ebSBarry Smith ierr = PetscViewerASCIIMonitorDestroy(viewer);CHKERRQ(ierr); 338a34d58ebSBarry Smith } 3393a40ed3dSBarry Smith PetscFunctionReturn(0); 3402d3f70b5SBarry Smith } 3412d3f70b5SBarry Smith 3424a2ae208SSatish Balay #undef __FUNCT__ 3434a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Pseudo" 3446849ba73SBarry Smith static PetscErrorCode TSSetFromOptions_Pseudo(TS ts) 3452d3f70b5SBarry Smith { 3464bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 347dfbe8321SBarry Smith PetscErrorCode ierr; 348ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE; 349a34d58ebSBarry Smith PetscViewerASCIIMonitor viewer; 3502d3f70b5SBarry Smith 3513a40ed3dSBarry Smith PetscFunctionBegin; 352b0a32e0cSBarry Smith ierr = PetscOptionsHead("Pseudo-timestepping options");CHKERRQ(ierr); 353acfcf0e5SJed Brown ierr = PetscOptionsBool("-ts_monitor","Monitor convergence","TSPseudoMonitorDefault",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 3542d3f70b5SBarry Smith if (flg) { 3557adad957SLisandro Dalcin ierr = PetscViewerASCIIMonitorCreate(((PetscObject)ts)->comm,"stdout",0,&viewer);CHKERRQ(ierr); 356a34d58ebSBarry Smith ierr = TSMonitorSet(ts,TSPseudoMonitorDefault,viewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 35728aa8177SBarry Smith } 35890d69ab7SBarry Smith flg = PETSC_FALSE; 359acfcf0e5SJed Brown ierr = PetscOptionsBool("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 360ca4b7087SBarry Smith if (flg) { 361ca4b7087SBarry Smith ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr); 362ca4b7087SBarry Smith } 363419fbf4bSSatish Balay ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,0);CHKERRQ(ierr); 364b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 3653a40ed3dSBarry Smith PetscFunctionReturn(0); 3662d3f70b5SBarry Smith } 3672d3f70b5SBarry Smith 3684a2ae208SSatish Balay #undef __FUNCT__ 3694a2ae208SSatish Balay #define __FUNCT__ "TSView_Pseudo" 3706849ba73SBarry Smith static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer) 3712d3f70b5SBarry Smith { 3723a40ed3dSBarry Smith PetscFunctionBegin; 3733a40ed3dSBarry Smith PetscFunctionReturn(0); 3742d3f70b5SBarry Smith } 3752d3f70b5SBarry Smith 37682bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 3774a2ae208SSatish Balay #undef __FUNCT__ 3784a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep" 379ac226902SBarry Smith /*@C 38082bf6240SBarry Smith TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the 38182bf6240SBarry Smith last timestep. 38282bf6240SBarry Smith 3833f9fe445SBarry Smith Logically Collective on TS 38415091d37SBarry Smith 38582bf6240SBarry Smith Input Parameters: 38615091d37SBarry Smith + ts - timestep context 38782bf6240SBarry Smith . dt - user-defined function to verify timestep 38815091d37SBarry Smith - ctx - [optional] user-defined context for private data 38982bf6240SBarry Smith for the timestep verification routine (may be PETSC_NULL) 39082bf6240SBarry Smith 39115091d37SBarry Smith Level: advanced 392fee21e36SBarry Smith 39382bf6240SBarry Smith Calling sequence of func: 394ace3abfcSBarry Smith . func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscBool *flag); 39582bf6240SBarry Smith 39682bf6240SBarry Smith . update - latest solution vector 39782bf6240SBarry Smith . ctx - [optional] timestep context 39882bf6240SBarry Smith . newdt - the timestep to use for the next step 39982bf6240SBarry Smith . flag - flag indicating whether the last time step was acceptable 40082bf6240SBarry Smith 40182bf6240SBarry Smith Notes: 40282bf6240SBarry Smith The routine set here will be called by TSPseudoVerifyTimeStep() 40382bf6240SBarry Smith during the timestepping process. 40482bf6240SBarry Smith 40582bf6240SBarry Smith .keywords: timestep, pseudo, set, verify 40682bf6240SBarry Smith 40782bf6240SBarry Smith .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep() 40882bf6240SBarry Smith @*/ 4097087cfbeSBarry Smith PetscErrorCode TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,PetscBool *),void* ctx) 41082bf6240SBarry Smith { 4114ac538c5SBarry Smith PetscErrorCode ierr; 41282bf6240SBarry Smith 41382bf6240SBarry Smith PetscFunctionBegin; 4140700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 4154ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSPseudoSetVerifyTimeStep_C",(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal *,PetscBool *),void *),(ts,dt,ctx));CHKERRQ(ierr); 41682bf6240SBarry Smith PetscFunctionReturn(0); 41782bf6240SBarry Smith } 41882bf6240SBarry Smith 4194a2ae208SSatish Balay #undef __FUNCT__ 4204a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement" 42182bf6240SBarry Smith /*@ 42282bf6240SBarry Smith TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 42382bf6240SBarry Smith dt when using the TSPseudoDefaultTimeStep() routine. 42482bf6240SBarry Smith 4253f9fe445SBarry Smith Logically Collective on TS 426fee21e36SBarry Smith 42715091d37SBarry Smith Input Parameters: 42815091d37SBarry Smith + ts - the timestep context 42915091d37SBarry Smith - inc - the scaling factor >= 1.0 43015091d37SBarry Smith 43182bf6240SBarry Smith Options Database Key: 43282bf6240SBarry Smith $ -ts_pseudo_increment <increment> 43382bf6240SBarry Smith 43415091d37SBarry Smith Level: advanced 43515091d37SBarry Smith 43682bf6240SBarry Smith .keywords: timestep, pseudo, set, increment 43782bf6240SBarry Smith 43882bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 43982bf6240SBarry Smith @*/ 4407087cfbeSBarry Smith PetscErrorCode TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc) 44182bf6240SBarry Smith { 4424ac538c5SBarry Smith PetscErrorCode ierr; 44382bf6240SBarry Smith 44482bf6240SBarry Smith PetscFunctionBegin; 4450700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 446c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts,inc,2); 4474ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSPseudoSetTimeStepIncrement_C",(TS,PetscReal),(ts,inc));CHKERRQ(ierr); 44882bf6240SBarry Smith PetscFunctionReturn(0); 44982bf6240SBarry Smith } 45082bf6240SBarry Smith 4514a2ae208SSatish Balay #undef __FUNCT__ 4524a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt" 45382bf6240SBarry Smith /*@ 45482bf6240SBarry Smith TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep 45582bf6240SBarry Smith is computed via the formula 45682bf6240SBarry Smith $ dt = initial_dt*initial_fnorm/current_fnorm 45782bf6240SBarry Smith rather than the default update, 45882bf6240SBarry Smith $ dt = current_dt*previous_fnorm/current_fnorm. 45982bf6240SBarry Smith 4603f9fe445SBarry Smith Logically Collective on TS 46115091d37SBarry Smith 46282bf6240SBarry Smith Input Parameter: 46382bf6240SBarry Smith . ts - the timestep context 46482bf6240SBarry Smith 46582bf6240SBarry Smith Options Database Key: 46682bf6240SBarry Smith $ -ts_pseudo_increment_dt_from_initial_dt 46782bf6240SBarry Smith 46815091d37SBarry Smith Level: advanced 46915091d37SBarry Smith 47082bf6240SBarry Smith .keywords: timestep, pseudo, set, increment 47182bf6240SBarry Smith 47282bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 47382bf6240SBarry Smith @*/ 4747087cfbeSBarry Smith PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS ts) 47582bf6240SBarry Smith { 4764ac538c5SBarry Smith PetscErrorCode ierr; 47782bf6240SBarry Smith 47882bf6240SBarry Smith PetscFunctionBegin; 4790700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 4804ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSPseudoIncrementDtFromInitialDt_C",(TS),(ts));CHKERRQ(ierr); 48182bf6240SBarry Smith PetscFunctionReturn(0); 48282bf6240SBarry Smith } 48382bf6240SBarry Smith 48482bf6240SBarry Smith 4854a2ae208SSatish Balay #undef __FUNCT__ 4864a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep" 487ac226902SBarry Smith /*@C 48882bf6240SBarry Smith TSPseudoSetTimeStep - Sets the user-defined routine to be 48982bf6240SBarry Smith called at each pseudo-timestep to update the timestep. 49082bf6240SBarry Smith 4913f9fe445SBarry Smith Logically Collective on TS 49215091d37SBarry Smith 49382bf6240SBarry Smith Input Parameters: 49415091d37SBarry Smith + ts - timestep context 49582bf6240SBarry Smith . dt - function to compute timestep 49615091d37SBarry Smith - ctx - [optional] user-defined context for private data 49782bf6240SBarry Smith required by the function (may be PETSC_NULL) 49882bf6240SBarry Smith 49915091d37SBarry Smith Level: intermediate 500fee21e36SBarry Smith 50182bf6240SBarry Smith Calling sequence of func: 50287828ca2SBarry Smith . func (TS ts,PetscReal *newdt,void *ctx); 50382bf6240SBarry Smith 50482bf6240SBarry Smith . newdt - the newly computed timestep 50582bf6240SBarry Smith . ctx - [optional] timestep context 50682bf6240SBarry Smith 50782bf6240SBarry Smith Notes: 50882bf6240SBarry Smith The routine set here will be called by TSPseudoComputeTimeStep() 50982bf6240SBarry Smith during the timestepping process. 51082bf6240SBarry Smith 51182bf6240SBarry Smith .keywords: timestep, pseudo, set 51282bf6240SBarry Smith 51382bf6240SBarry Smith .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep() 51482bf6240SBarry Smith @*/ 5157087cfbeSBarry Smith PetscErrorCode TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void* ctx) 51682bf6240SBarry Smith { 5174ac538c5SBarry Smith PetscErrorCode ierr; 51882bf6240SBarry Smith 51982bf6240SBarry Smith PetscFunctionBegin; 5200700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 5214ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSPseudoSetTimeStep_C",(TS,PetscErrorCode (*)(TS,PetscReal *,void *),void *),(ts,dt,ctx));CHKERRQ(ierr); 52282bf6240SBarry Smith PetscFunctionReturn(0); 52382bf6240SBarry Smith } 52482bf6240SBarry Smith 52582bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 52682bf6240SBarry Smith 527ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscBool *); /* force argument to next function to not be extern C*/ 528fb2e594dSBarry Smith EXTERN_C_BEGIN 5294a2ae208SSatish Balay #undef __FUNCT__ 5304a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo" 5317087cfbeSBarry Smith PetscErrorCode TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void* ctx) 53282bf6240SBarry Smith { 53382bf6240SBarry Smith TS_Pseudo *pseudo; 53482bf6240SBarry Smith 53582bf6240SBarry Smith PetscFunctionBegin; 53682bf6240SBarry Smith pseudo = (TS_Pseudo*)ts->data; 53782bf6240SBarry Smith pseudo->verify = dt; 53882bf6240SBarry Smith pseudo->verifyctx = ctx; 53982bf6240SBarry Smith PetscFunctionReturn(0); 54082bf6240SBarry Smith } 541fb2e594dSBarry Smith EXTERN_C_END 54282bf6240SBarry Smith 543fb2e594dSBarry Smith EXTERN_C_BEGIN 5444a2ae208SSatish Balay #undef __FUNCT__ 5454a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement_Pseudo" 5467087cfbeSBarry Smith PetscErrorCode TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc) 54782bf6240SBarry Smith { 5484bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 54982bf6240SBarry Smith 55082bf6240SBarry Smith PetscFunctionBegin; 55182bf6240SBarry Smith pseudo->dt_increment = inc; 55282bf6240SBarry Smith PetscFunctionReturn(0); 55382bf6240SBarry Smith } 554fb2e594dSBarry Smith EXTERN_C_END 55582bf6240SBarry Smith 556fb2e594dSBarry Smith EXTERN_C_BEGIN 5574a2ae208SSatish Balay #undef __FUNCT__ 5584a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt_Pseudo" 5597087cfbeSBarry Smith PetscErrorCode TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts) 56082bf6240SBarry Smith { 5614bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 56282bf6240SBarry Smith 56382bf6240SBarry Smith PetscFunctionBegin; 5644bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_TRUE; 56582bf6240SBarry Smith PetscFunctionReturn(0); 56682bf6240SBarry Smith } 567fb2e594dSBarry Smith EXTERN_C_END 56882bf6240SBarry Smith 5696849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/ 570fb2e594dSBarry Smith EXTERN_C_BEGIN 5714a2ae208SSatish Balay #undef __FUNCT__ 5724a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo" 5737087cfbeSBarry Smith PetscErrorCode TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void* ctx) 57482bf6240SBarry Smith { 5754bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 57682bf6240SBarry Smith 57782bf6240SBarry Smith PetscFunctionBegin; 57882bf6240SBarry Smith pseudo->dt = dt; 57982bf6240SBarry Smith pseudo->dtctx = ctx; 58082bf6240SBarry Smith PetscFunctionReturn(0); 58182bf6240SBarry Smith } 582fb2e594dSBarry Smith EXTERN_C_END 58382bf6240SBarry Smith 58482bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 58510e6a065SJed Brown /*MC 58610e6a065SJed Brown TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping 58782bf6240SBarry Smith 58810e6a065SJed Brown This method solves equations of the form 58910e6a065SJed Brown 59010e6a065SJed Brown $ F(X,Xdot) = 0 59110e6a065SJed Brown 59210e6a065SJed Brown for steady state using the iteration 59310e6a065SJed Brown 59410e6a065SJed Brown $ [G'] S = -F(X,0) 59510e6a065SJed Brown $ X += S 59610e6a065SJed Brown 59710e6a065SJed Brown where 59810e6a065SJed Brown 59910e6a065SJed Brown $ G(Y) = F(Y,(Y-X)/dt) 60010e6a065SJed Brown 6016f2d6a7bSJed Brown This is linearly-implicit Euler with the residual always evaluated "at steady 6026f2d6a7bSJed Brown state". See note below. 60310e6a065SJed Brown 60410e6a065SJed Brown Options database keys: 60510e6a065SJed Brown + -ts_pseudo_increment <real> - ratio of increase dt 60610e6a065SJed Brown - -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt 60710e6a065SJed Brown 60810e6a065SJed Brown Level: beginner 60910e6a065SJed Brown 61010e6a065SJed Brown References: 61110e6a065SJed Brown Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential-Algebraic Equations, 2003. 61210e6a065SJed Brown C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998. 61310e6a065SJed Brown 61410e6a065SJed Brown Notes: 6156f2d6a7bSJed Brown The residual computed by this method includes the transient term (Xdot is computed instead of 6166f2d6a7bSJed Brown always being zero), but since the prediction from the last step is always the solution from the 6176f2d6a7bSJed Brown last step, on the first Newton iteration we have 6186f2d6a7bSJed Brown 6196f2d6a7bSJed Brown $ Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0 6206f2d6a7bSJed Brown 6216f2d6a7bSJed Brown Therefore, the linear system solved by the first Newton iteration is equivalent to the one 6226f2d6a7bSJed Brown described above and in the papers. If the user chooses to perform multiple Newton iterations, the 6236f2d6a7bSJed Brown algorithm is no longer the one described in the referenced papers. 62410e6a065SJed Brown 62510e6a065SJed Brown .seealso: TSCreate(), TS, TSSetType() 62610e6a065SJed Brown 62710e6a065SJed Brown M*/ 628fb2e594dSBarry Smith EXTERN_C_BEGIN 6294a2ae208SSatish Balay #undef __FUNCT__ 6304a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Pseudo" 6317087cfbeSBarry Smith PetscErrorCode TSCreate_Pseudo(TS ts) 6322d3f70b5SBarry Smith { 6337bf11e45SBarry Smith TS_Pseudo *pseudo; 634dfbe8321SBarry Smith PetscErrorCode ierr; 6352d3f70b5SBarry Smith 6363a40ed3dSBarry Smith PetscFunctionBegin; 637277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Pseudo; 638000e7ae3SMatthew Knepley ts->ops->destroy = TSDestroy_Pseudo; 639000e7ae3SMatthew Knepley ts->ops->view = TSView_Pseudo; 6402d3f70b5SBarry Smith 64117186662SBarry Smith if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Only for nonlinear problems"); 642000e7ae3SMatthew Knepley ts->ops->setup = TSSetUp_Pseudo; 643000e7ae3SMatthew Knepley ts->ops->step = TSStep_Pseudo; 644000e7ae3SMatthew Knepley ts->ops->setfromoptions = TSSetFromOptions_Pseudo; 6450f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Pseudo; 6460f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Pseudo; 6477bf11e45SBarry Smith 6487bf11e45SBarry Smith /* create the required nonlinear solver context */ 6497adad957SLisandro Dalcin ierr = SNESCreate(((PetscObject)ts)->comm,&ts->snes);CHKERRQ(ierr); 65017efb626SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);CHKERRQ(ierr); 6512d3f70b5SBarry Smith 65238f2d2fdSLisandro Dalcin ierr = PetscNewLog(ts,TS_Pseudo,&pseudo);CHKERRQ(ierr); 6537bf11e45SBarry Smith ts->data = (void*)pseudo; 6542d3f70b5SBarry Smith 65528aa8177SBarry Smith pseudo->dt_increment = 1.1; 6564bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_FALSE; 65728aa8177SBarry Smith pseudo->dt = TSPseudoDefaultTimeStep; 65882bf6240SBarry Smith 659f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C", 660e1311b90SBarry Smith "TSPseudoSetVerifyTimeStep_Pseudo", 6610c97f09dSSatish Balay TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr); 662f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C", 663e1311b90SBarry Smith "TSPseudoSetTimeStepIncrement_Pseudo", 6640c97f09dSSatish Balay TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr); 665f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C", 666e1311b90SBarry Smith "TSPseudoIncrementDtFromInitialDt_Pseudo", 6670c97f09dSSatish Balay TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr); 6680c97f09dSSatish Balay ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C", 6690c97f09dSSatish Balay "TSPseudoSetTimeStep_Pseudo", 6700c97f09dSSatish Balay TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr); 6713a40ed3dSBarry Smith PetscFunctionReturn(0); 6722d3f70b5SBarry Smith } 673fb2e594dSBarry Smith EXTERN_C_END 6742d3f70b5SBarry Smith 6754a2ae208SSatish Balay #undef __FUNCT__ 6764a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultTimeStep" 67782bf6240SBarry Smith /*@C 67882bf6240SBarry Smith TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping. 67982bf6240SBarry Smith Use with TSPseudoSetTimeStep(). 68028aa8177SBarry Smith 68115091d37SBarry Smith Collective on TS 68215091d37SBarry Smith 68328aa8177SBarry Smith Input Parameters: 68428aa8177SBarry Smith . ts - the timestep context 68582bf6240SBarry Smith . dtctx - unused timestep context 68628aa8177SBarry Smith 68782bf6240SBarry Smith Output Parameter: 68882bf6240SBarry Smith . newdt - the timestep to use for the next step 68928aa8177SBarry Smith 69015091d37SBarry Smith Level: advanced 69115091d37SBarry Smith 69282bf6240SBarry Smith .keywords: timestep, pseudo, default 693564e8f4eSLois Curfman McInnes 69482bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep() 69528aa8177SBarry Smith @*/ 6967087cfbeSBarry Smith PetscErrorCode TSPseudoDefaultTimeStep(TS ts,PetscReal* newdt,void* dtctx) 69728aa8177SBarry Smith { 69882bf6240SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 69987828ca2SBarry Smith PetscReal inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous; 700dfbe8321SBarry Smith PetscErrorCode ierr; 70128aa8177SBarry Smith 7023a40ed3dSBarry Smith PetscFunctionBegin; 703bbd7b040SJed Brown ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr); 704bbd7b040SJed Brown ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func);CHKERRQ(ierr); 70582bf6240SBarry Smith ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 70682bf6240SBarry Smith if (pseudo->initial_fnorm == 0.0) { 70782bf6240SBarry Smith /* first time through so compute initial function norm */ 70882bf6240SBarry Smith pseudo->initial_fnorm = pseudo->fnorm; 70982bf6240SBarry Smith fnorm_previous = pseudo->fnorm; 71082bf6240SBarry Smith } 71182bf6240SBarry Smith if (pseudo->fnorm == 0.0) { 71282bf6240SBarry Smith *newdt = 1.e12*inc*ts->time_step; 713004884caSBarry Smith } else if (pseudo->increment_dt_from_initial_dt) { 71482bf6240SBarry Smith *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm; 71582bf6240SBarry Smith } else { 71682bf6240SBarry Smith *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm; 71782bf6240SBarry Smith } 71882bf6240SBarry Smith pseudo->fnorm_previous = pseudo->fnorm; 7193a40ed3dSBarry Smith PetscFunctionReturn(0); 72028aa8177SBarry Smith } 721