163dd3a1aSKris Buschelman 22d3f70b5SBarry Smith /* 3fb4a63b6SLois Curfman McInnes Code for Timestepping with implicit backwards Euler. 42d3f70b5SBarry Smith */ 5c6db04a5SJed Brown #include <private/tsimpl.h> /*I "petscts.h" I*/ 62d3f70b5SBarry Smith 72d3f70b5SBarry Smith typedef struct { 82d3f70b5SBarry Smith Vec update; /* work vector where new solution is formed */ 92d3f70b5SBarry Smith Vec func; /* work vector where F(t[i],u[i]) is stored */ 106f2d6a7bSJed Brown Vec xdot; /* work vector for time derivative of state */ 112d3f70b5SBarry Smith 122d3f70b5SBarry Smith /* information used for Pseudo-timestepping */ 132d3f70b5SBarry Smith 146849ba73SBarry Smith PetscErrorCode (*dt)(TS,PetscReal*,void*); /* compute next timestep, and related context */ 152d3f70b5SBarry Smith void *dtctx; 16ace3abfcSBarry Smith PetscErrorCode (*verify)(TS,Vec,void*,PetscReal*,PetscBool *); /* verify previous timestep and related context */ 177bf11e45SBarry Smith void *verifyctx; 182d3f70b5SBarry Smith 1987828ca2SBarry Smith PetscReal initial_fnorm,fnorm; /* original and current norm of F(u) */ 2087828ca2SBarry Smith PetscReal fnorm_previous; 2128aa8177SBarry Smith 2287828ca2SBarry Smith PetscReal dt_increment; /* scaling that dt is incremented each time-step */ 23ace3abfcSBarry Smith PetscBool increment_dt_from_initial_dt; 247bf11e45SBarry Smith } TS_Pseudo; 252d3f70b5SBarry Smith 262d3f70b5SBarry Smith /* ------------------------------------------------------------------------------*/ 272d3f70b5SBarry Smith 284a2ae208SSatish Balay #undef __FUNCT__ 294a2ae208SSatish Balay #define __FUNCT__ "TSPseudoComputeTimeStep" 307bf11e45SBarry Smith /*@ 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 4215091d37SBarry Smith Level: advanced 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 48fb4a63b6SLois Curfman McInnes .keywords: timestep, pseudo, compute 49564e8f4eSLois Curfman McInnes 50564e8f4eSLois Curfman McInnes .seealso: TSPseudoDefaultTimeStep(), TSPseudoSetTimeStep() 517bf11e45SBarry Smith @*/ 527087cfbeSBarry Smith PetscErrorCode TSPseudoComputeTimeStep(TS ts,PetscReal *dt) 537bf11e45SBarry Smith { 547bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 55dfbe8321SBarry Smith PetscErrorCode ierr; 567bf11e45SBarry Smith 573a40ed3dSBarry Smith PetscFunctionBegin; 58d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr); 597bf11e45SBarry Smith ierr = (*pseudo->dt)(ts,dt,pseudo->dtctx);CHKERRQ(ierr); 60d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr); 613a40ed3dSBarry Smith PetscFunctionReturn(0); 627bf11e45SBarry Smith } 637bf11e45SBarry Smith 647bf11e45SBarry Smith 657bf11e45SBarry Smith /* ------------------------------------------------------------------------------*/ 664a2ae208SSatish Balay #undef __FUNCT__ 674a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultVerifyTimeStep" 687bf11e45SBarry Smith /*@C 69639f9d9dSBarry Smith TSPseudoDefaultVerifyTimeStep - Default code to verify the quality of the last timestep. 707bf11e45SBarry Smith 7115091d37SBarry Smith Collective on TS 7215091d37SBarry Smith 737bf11e45SBarry Smith Input Parameters: 7415091d37SBarry Smith + ts - the timestep context 757bf11e45SBarry Smith . dtctx - unused timestep context 7615091d37SBarry Smith - update - latest solution vector 777bf11e45SBarry Smith 78564e8f4eSLois Curfman McInnes Output Parameters: 7915091d37SBarry Smith + newdt - the timestep to use for the next step 8015091d37SBarry Smith - flag - flag indicating whether the last time step was acceptable 817bf11e45SBarry Smith 8215091d37SBarry Smith Level: advanced 83fee21e36SBarry Smith 84564e8f4eSLois Curfman McInnes Note: 85564e8f4eSLois Curfman McInnes This routine always returns a flag of 1, indicating an acceptable 86564e8f4eSLois Curfman McInnes timestep. 87564e8f4eSLois Curfman McInnes 88564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, default, verify 89564e8f4eSLois Curfman McInnes 90564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStep() 917bf11e45SBarry Smith @*/ 927087cfbeSBarry Smith PetscErrorCode TSPseudoDefaultVerifyTimeStep(TS ts,Vec update,void *dtctx,PetscReal *newdt,PetscBool *flag) 937bf11e45SBarry Smith { 943a40ed3dSBarry Smith PetscFunctionBegin; 95a7cc72afSBarry Smith *flag = PETSC_TRUE; 963a40ed3dSBarry Smith PetscFunctionReturn(0); 977bf11e45SBarry Smith } 987bf11e45SBarry Smith 997bf11e45SBarry Smith 1004a2ae208SSatish Balay #undef __FUNCT__ 1014a2ae208SSatish Balay #define __FUNCT__ "TSPseudoVerifyTimeStep" 1027bf11e45SBarry Smith /*@ 103564e8f4eSLois Curfman McInnes TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable. 1047bf11e45SBarry Smith 10515091d37SBarry Smith Collective on TS 10615091d37SBarry Smith 107fb4a63b6SLois Curfman McInnes Input Parameters: 10815091d37SBarry Smith + ts - timestep context 10915091d37SBarry Smith - update - latest solution vector 1107bf11e45SBarry Smith 111fb4a63b6SLois Curfman McInnes Output Parameters: 11215091d37SBarry Smith + dt - newly computed timestep (if it had to shrink) 11315091d37SBarry Smith - flag - indicates if current timestep was ok 1147bf11e45SBarry Smith 11515091d37SBarry Smith Level: advanced 116fee21e36SBarry Smith 117564e8f4eSLois Curfman McInnes Notes: 118564e8f4eSLois Curfman McInnes The routine to be called here to compute the timestep should be 119564e8f4eSLois Curfman McInnes set by calling TSPseudoSetVerifyTimeStep(). 120564e8f4eSLois Curfman McInnes 121564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, verify 122564e8f4eSLois Curfman McInnes 123564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoDefaultVerifyTimeStep() 1247bf11e45SBarry Smith @*/ 1257087cfbeSBarry Smith PetscErrorCode TSPseudoVerifyTimeStep(TS ts,Vec update,PetscReal *dt,PetscBool *flag) 1267bf11e45SBarry Smith { 1277bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 128dfbe8321SBarry Smith PetscErrorCode ierr; 1297bf11e45SBarry Smith 1303a40ed3dSBarry Smith PetscFunctionBegin; 131a7cc72afSBarry Smith if (!pseudo->verify) {*flag = PETSC_TRUE; PetscFunctionReturn(0);} 1327bf11e45SBarry Smith 1337bf11e45SBarry Smith ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag);CHKERRQ(ierr); 1347bf11e45SBarry Smith 1353a40ed3dSBarry Smith PetscFunctionReturn(0); 1367bf11e45SBarry Smith } 1377bf11e45SBarry Smith 1387bf11e45SBarry Smith /* --------------------------------------------------------------------------------*/ 1397bf11e45SBarry Smith 1404a2ae208SSatish Balay #undef __FUNCT__ 1414a2ae208SSatish Balay #define __FUNCT__ "TSStep_Pseudo" 142a7cc72afSBarry Smith static PetscErrorCode TSStep_Pseudo(TS ts,PetscInt *steps,PetscReal *ptime) 1432d3f70b5SBarry Smith { 144*277b19d0SLisandro Dalcin TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 145*277b19d0SLisandro Dalcin Vec sol = ts->vec_sol, update = pseudo->update; 146dfbe8321SBarry Smith PetscErrorCode ierr; 147a7cc72afSBarry Smith PetscInt i,max_steps = ts->max_steps,its,lits; 148ace3abfcSBarry Smith PetscBool ok; 14987828ca2SBarry Smith PetscReal current_time_step; 1502d3f70b5SBarry Smith 1513a40ed3dSBarry Smith PetscFunctionBegin; 1522d3f70b5SBarry Smith *steps = -ts->steps; 1532d3f70b5SBarry Smith 154*277b19d0SLisandro Dalcin ierr = VecCopy(sol,update);CHKERRQ(ierr); 1552d3f70b5SBarry Smith for (i=0; i<max_steps && ts->ptime < ts->max_time; i++) { 1567bf11e45SBarry Smith ierr = TSPseudoComputeTimeStep(ts,&ts->time_step);CHKERRQ(ierr); 157e6e5fe25SBarry Smith ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 1587bf11e45SBarry Smith current_time_step = ts->time_step; 1593f2090d5SJed Brown ierr = TSPreStep(ts);CHKERRQ(ierr); 160e6e5fe25SBarry Smith while (PETSC_TRUE) { 1617bf11e45SBarry Smith ts->ptime += current_time_step; 162*277b19d0SLisandro Dalcin ierr = SNESSolve(ts->snes,PETSC_NULL,update);CHKERRQ(ierr); 163b850b91aSLisandro Dalcin ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 1649f954729SBarry Smith ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr); 165c9780b6fSBarry Smith ts->nonlinear_its += its; ts->linear_its += lits; 166*277b19d0SLisandro Dalcin ierr = TSPseudoVerifyTimeStep(ts,update,&ts->time_step,&ok);CHKERRQ(ierr); 1677bf11e45SBarry Smith if (ok) break; 1687bf11e45SBarry Smith ts->ptime -= current_time_step; 1697bf11e45SBarry Smith current_time_step = ts->time_step; 1707bf11e45SBarry Smith } 171*277b19d0SLisandro Dalcin ierr = VecCopy(update,sol);CHKERRQ(ierr); 1722d3f70b5SBarry Smith ts->steps++; 1733f2090d5SJed Brown ierr = TSPostStep(ts);CHKERRQ(ierr); 1742d3f70b5SBarry Smith } 175bbd7b040SJed Brown ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr); 176bbd7b040SJed Brown ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func);CHKERRQ(ierr); 177e6e5fe25SBarry Smith ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 178e6e5fe25SBarry Smith ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 1792d3f70b5SBarry Smith 1802d3f70b5SBarry Smith *steps += ts->steps; 181142b95e3SSatish Balay *ptime = ts->ptime; 1823a40ed3dSBarry Smith PetscFunctionReturn(0); 1832d3f70b5SBarry Smith } 1842d3f70b5SBarry Smith 1852d3f70b5SBarry Smith /*------------------------------------------------------------*/ 1864a2ae208SSatish Balay #undef __FUNCT__ 187*277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Pseudo" 188*277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Pseudo(TS ts) 1892d3f70b5SBarry Smith { 1907bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 191dfbe8321SBarry Smith PetscErrorCode ierr; 1922d3f70b5SBarry Smith 1933a40ed3dSBarry Smith PetscFunctionBegin; 194e4ef1970SSatish Balay if (pseudo->update) {ierr = VecDestroy(pseudo->update);CHKERRQ(ierr);} 1957bf11e45SBarry Smith if (pseudo->func) {ierr = VecDestroy(pseudo->func);CHKERRQ(ierr);} 1966f2d6a7bSJed Brown if (pseudo->xdot) {ierr = VecDestroy(pseudo->xdot);CHKERRQ(ierr);} 1973a40ed3dSBarry Smith PetscFunctionReturn(0); 1982d3f70b5SBarry Smith } 1992d3f70b5SBarry Smith 200*277b19d0SLisandro Dalcin #undef __FUNCT__ 201*277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Pseudo" 202*277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Pseudo(TS ts) 203*277b19d0SLisandro Dalcin { 204*277b19d0SLisandro Dalcin PetscErrorCode ierr; 205*277b19d0SLisandro Dalcin 206*277b19d0SLisandro Dalcin PetscFunctionBegin; 207*277b19d0SLisandro Dalcin ierr = TSReset_Pseudo(ts);CHKERRQ(ierr); 208*277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 209*277b19d0SLisandro Dalcin PetscFunctionReturn(0); 210*277b19d0SLisandro Dalcin } 2112d3f70b5SBarry Smith 2122d3f70b5SBarry Smith /*------------------------------------------------------------*/ 2132d3f70b5SBarry Smith 2144a2ae208SSatish Balay #undef __FUNCT__ 2156f2d6a7bSJed Brown #define __FUNCT__ "TSPseudoGetXdot" 2166f2d6a7bSJed Brown /* 2176f2d6a7bSJed Brown Compute Xdot = (X^{n+1}-X^n)/dt) = 0 2186f2d6a7bSJed Brown */ 2196f2d6a7bSJed Brown static PetscErrorCode TSPseudoGetXdot(TS ts,Vec X,Vec *Xdot) 2202d3f70b5SBarry Smith { 2216f2d6a7bSJed Brown TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 2226f2d6a7bSJed Brown PetscScalar mdt = 1.0/ts->time_step,*xnp1,*xn,*xdot; 223dfbe8321SBarry Smith PetscErrorCode ierr; 224a7cc72afSBarry Smith PetscInt i,n; 2252d3f70b5SBarry Smith 2263a40ed3dSBarry Smith PetscFunctionBegin; 2276f2d6a7bSJed Brown ierr = VecGetArray(ts->vec_sol,&xn);CHKERRQ(ierr); 2286f2d6a7bSJed Brown ierr = VecGetArray(X,&xnp1);CHKERRQ(ierr); 2296f2d6a7bSJed Brown ierr = VecGetArray(pseudo->xdot,&xdot);CHKERRQ(ierr); 2306f2d6a7bSJed Brown ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr); 2312d3f70b5SBarry Smith for (i=0; i<n; i++) { 2326f2d6a7bSJed Brown xdot[i] = mdt*(xnp1[i] - xn[i]); 2332d3f70b5SBarry Smith } 2346f2d6a7bSJed Brown ierr = VecRestoreArray(ts->vec_sol,&xn);CHKERRQ(ierr); 2356f2d6a7bSJed Brown ierr = VecRestoreArray(X,&xnp1);CHKERRQ(ierr); 2366f2d6a7bSJed Brown ierr = VecRestoreArray(pseudo->xdot,&xdot);CHKERRQ(ierr); 2376f2d6a7bSJed Brown *Xdot = pseudo->xdot; 2383a40ed3dSBarry Smith PetscFunctionReturn(0); 2392d3f70b5SBarry Smith } 2402d3f70b5SBarry Smith 2414a2ae208SSatish Balay #undef __FUNCT__ 2420f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Pseudo" 2436f2d6a7bSJed Brown /* 2446f2d6a7bSJed Brown The transient residual is 2456f2d6a7bSJed Brown 2466f2d6a7bSJed Brown F(U^{n+1},(U^{n+1}-U^n)/dt) = 0 2476f2d6a7bSJed Brown 2486f2d6a7bSJed Brown or for ODE, 2496f2d6a7bSJed Brown 2506f2d6a7bSJed Brown (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0 2516f2d6a7bSJed Brown 2526f2d6a7bSJed Brown This is the function that must be evaluated for transient simulation and for 2536f2d6a7bSJed Brown finite difference Jacobians. On the first Newton step, this algorithm uses 2546f2d6a7bSJed Brown a guess of U^{n+1} = U^n in which case the transient term vanishes and the 2556f2d6a7bSJed Brown residual is actually the steady state residual. Pseudotransient 2566f2d6a7bSJed Brown continuation as described in the literature is a linearly implicit 2576f2d6a7bSJed Brown algorithm, it only takes this one Newton step with the steady state 2586f2d6a7bSJed Brown residual, and then advances to the next time step. 2596f2d6a7bSJed Brown */ 2600f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes,Vec X,Vec Y,TS ts) 2612d3f70b5SBarry Smith { 2626f2d6a7bSJed Brown Vec Xdot; 263dfbe8321SBarry Smith PetscErrorCode ierr; 2642d3f70b5SBarry Smith 2653a40ed3dSBarry Smith PetscFunctionBegin; 2666f2d6a7bSJed Brown ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr); 2676f2d6a7bSJed Brown ierr = TSComputeIFunction(ts,ts->ptime,X,Xdot,Y);CHKERRQ(ierr); 2686f2d6a7bSJed Brown PetscFunctionReturn(0); 2696f2d6a7bSJed Brown } 2702d3f70b5SBarry Smith 2716f2d6a7bSJed Brown #undef __FUNCT__ 2720f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Pseudo" 2736f2d6a7bSJed Brown /* 2746f2d6a7bSJed Brown This constructs the Jacobian needed for SNES. For DAE, this is 2756f2d6a7bSJed Brown 2766f2d6a7bSJed Brown dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot 2776f2d6a7bSJed Brown 2786f2d6a7bSJed Brown and for ODE: 2796f2d6a7bSJed Brown 2806f2d6a7bSJed Brown J = I/dt - J_{Frhs} where J_{Frhs} is the given Jacobian of Frhs. 2816f2d6a7bSJed Brown */ 2820f5c6efeSJed Brown static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes,Vec X,Mat *AA,Mat *BB,MatStructure *str,TS ts) 2836f2d6a7bSJed Brown { 2846f2d6a7bSJed Brown Vec Xdot; 2856f2d6a7bSJed Brown PetscErrorCode ierr; 2866f2d6a7bSJed Brown 2876f2d6a7bSJed Brown PetscFunctionBegin; 2886f2d6a7bSJed Brown ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr); 2896f2d6a7bSJed Brown ierr = TSComputeIJacobian(ts,ts->ptime,X,Xdot,1./ts->time_step,AA,BB,str);CHKERRQ(ierr); 2903a40ed3dSBarry Smith PetscFunctionReturn(0); 2912d3f70b5SBarry Smith } 2922d3f70b5SBarry Smith 2932d3f70b5SBarry Smith 2944a2ae208SSatish Balay #undef __FUNCT__ 2954a2ae208SSatish Balay #define __FUNCT__ "TSSetUp_Pseudo" 2966849ba73SBarry Smith static PetscErrorCode TSSetUp_Pseudo(TS ts) 2972d3f70b5SBarry Smith { 2987bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 299dfbe8321SBarry Smith PetscErrorCode ierr; 3002d3f70b5SBarry Smith 3013a40ed3dSBarry Smith PetscFunctionBegin; 3027bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr); 3037bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr); 3046f2d6a7bSJed Brown ierr = VecDuplicate(ts->vec_sol,&pseudo->xdot);CHKERRQ(ierr); 3050f5c6efeSJed Brown ierr = SNESSetFunction(ts->snes,pseudo->func,SNESTSFormFunction,ts);CHKERRQ(ierr); 30644c54cacSJed Brown /* This is nasty. SNESSetFromOptions() is usually called in TSSetFromOptions(). With -snes_mf_operator, it will 30744c54cacSJed 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 30844c54cacSJed Brown context. Note that SNESSetFunction() normally has not been called before SNESSetFromOptions(), so when -snes_mf sets 30944c54cacSJed Brown the Jacobian user context to snes->funP, it will actually be NULL. This is not a problem because both snes->funP and 31044c54cacSJed Brown snes->jacP should be the TS. */ 31144c54cacSJed Brown { 31244c54cacSJed Brown Mat A,B; 31344c54cacSJed Brown PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*); 31444c54cacSJed Brown void *ctx; 31544c54cacSJed Brown ierr = SNESGetJacobian(ts->snes,&A,&B,&func,&ctx);CHKERRQ(ierr); 31644c54cacSJed Brown ierr = SNESSetJacobian(ts->snes,A?A:ts->A,B?B:ts->B,func?func:SNESTSFormJacobian,ctx?ctx:ts);CHKERRQ(ierr); 31744c54cacSJed Brown } 3183a40ed3dSBarry Smith PetscFunctionReturn(0); 3192d3f70b5SBarry Smith } 3202d3f70b5SBarry Smith /*------------------------------------------------------------*/ 3212d3f70b5SBarry Smith 3224a2ae208SSatish Balay #undef __FUNCT__ 323a6570f20SBarry Smith #define __FUNCT__ "TSPseudoMonitorDefault" 324a6570f20SBarry Smith PetscErrorCode TSPseudoMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *ctx) 3252d3f70b5SBarry Smith { 3267bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 327dfbe8321SBarry Smith PetscErrorCode ierr; 328a34d58ebSBarry Smith PetscViewerASCIIMonitor viewer = (PetscViewerASCIIMonitor)ctx; 3292d3f70b5SBarry Smith 3303a40ed3dSBarry Smith PetscFunctionBegin; 331a34d58ebSBarry Smith if (!ctx) { 3327adad957SLisandro Dalcin ierr = PetscViewerASCIIMonitorCreate(((PetscObject)ts)->comm,"stdout",0,&viewer);CHKERRQ(ierr); 333a34d58ebSBarry Smith } 334a34d58ebSBarry Smith ierr = PetscViewerASCIIMonitorPrintf(viewer,"TS %D dt %G time %G fnorm %G\n",step,ts->time_step,ptime,pseudo->fnorm);CHKERRQ(ierr); 335a34d58ebSBarry Smith if (!ctx) { 336a34d58ebSBarry Smith ierr = PetscViewerASCIIMonitorDestroy(viewer);CHKERRQ(ierr); 337a34d58ebSBarry Smith } 3383a40ed3dSBarry Smith PetscFunctionReturn(0); 3392d3f70b5SBarry Smith } 3402d3f70b5SBarry Smith 3414a2ae208SSatish Balay #undef __FUNCT__ 3424a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Pseudo" 3436849ba73SBarry Smith static PetscErrorCode TSSetFromOptions_Pseudo(TS ts) 3442d3f70b5SBarry Smith { 3454bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 346dfbe8321SBarry Smith PetscErrorCode ierr; 347ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE; 348a34d58ebSBarry Smith PetscViewerASCIIMonitor viewer; 3492d3f70b5SBarry Smith 3503a40ed3dSBarry Smith PetscFunctionBegin; 351b0a32e0cSBarry Smith ierr = PetscOptionsHead("Pseudo-timestepping options");CHKERRQ(ierr); 352acfcf0e5SJed Brown ierr = PetscOptionsBool("-ts_monitor","Monitor convergence","TSPseudoMonitorDefault",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 3532d3f70b5SBarry Smith if (flg) { 3547adad957SLisandro Dalcin ierr = PetscViewerASCIIMonitorCreate(((PetscObject)ts)->comm,"stdout",0,&viewer);CHKERRQ(ierr); 355a34d58ebSBarry Smith ierr = TSMonitorSet(ts,TSPseudoMonitorDefault,viewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 35628aa8177SBarry Smith } 35790d69ab7SBarry Smith flg = PETSC_FALSE; 358acfcf0e5SJed 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); 359ca4b7087SBarry Smith if (flg) { 360ca4b7087SBarry Smith ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr); 361ca4b7087SBarry Smith } 362419fbf4bSSatish Balay ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,0);CHKERRQ(ierr); 363b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 3643a40ed3dSBarry Smith PetscFunctionReturn(0); 3652d3f70b5SBarry Smith } 3662d3f70b5SBarry Smith 3674a2ae208SSatish Balay #undef __FUNCT__ 3684a2ae208SSatish Balay #define __FUNCT__ "TSView_Pseudo" 3696849ba73SBarry Smith static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer) 3702d3f70b5SBarry Smith { 3713a40ed3dSBarry Smith PetscFunctionBegin; 3723a40ed3dSBarry Smith PetscFunctionReturn(0); 3732d3f70b5SBarry Smith } 3742d3f70b5SBarry Smith 37582bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 3764a2ae208SSatish Balay #undef __FUNCT__ 3774a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep" 378ac226902SBarry Smith /*@C 37982bf6240SBarry Smith TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the 38082bf6240SBarry Smith last timestep. 38182bf6240SBarry Smith 3823f9fe445SBarry Smith Logically Collective on TS 38315091d37SBarry Smith 38482bf6240SBarry Smith Input Parameters: 38515091d37SBarry Smith + ts - timestep context 38682bf6240SBarry Smith . dt - user-defined function to verify timestep 38715091d37SBarry Smith - ctx - [optional] user-defined context for private data 38882bf6240SBarry Smith for the timestep verification routine (may be PETSC_NULL) 38982bf6240SBarry Smith 39015091d37SBarry Smith Level: advanced 391fee21e36SBarry Smith 39282bf6240SBarry Smith Calling sequence of func: 393ace3abfcSBarry Smith . func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscBool *flag); 39482bf6240SBarry Smith 39582bf6240SBarry Smith . update - latest solution vector 39682bf6240SBarry Smith . ctx - [optional] timestep context 39782bf6240SBarry Smith . newdt - the timestep to use for the next step 39882bf6240SBarry Smith . flag - flag indicating whether the last time step was acceptable 39982bf6240SBarry Smith 40082bf6240SBarry Smith Notes: 40182bf6240SBarry Smith The routine set here will be called by TSPseudoVerifyTimeStep() 40282bf6240SBarry Smith during the timestepping process. 40382bf6240SBarry Smith 40482bf6240SBarry Smith .keywords: timestep, pseudo, set, verify 40582bf6240SBarry Smith 40682bf6240SBarry Smith .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep() 40782bf6240SBarry Smith @*/ 4087087cfbeSBarry Smith PetscErrorCode TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,PetscBool *),void* ctx) 40982bf6240SBarry Smith { 4104ac538c5SBarry Smith PetscErrorCode ierr; 41182bf6240SBarry Smith 41282bf6240SBarry Smith PetscFunctionBegin; 4130700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 4144ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSPseudoSetVerifyTimeStep_C",(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal *,PetscBool *),void *),(ts,dt,ctx));CHKERRQ(ierr); 41582bf6240SBarry Smith PetscFunctionReturn(0); 41682bf6240SBarry Smith } 41782bf6240SBarry Smith 4184a2ae208SSatish Balay #undef __FUNCT__ 4194a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement" 42082bf6240SBarry Smith /*@ 42182bf6240SBarry Smith TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 42282bf6240SBarry Smith dt when using the TSPseudoDefaultTimeStep() routine. 42382bf6240SBarry Smith 4243f9fe445SBarry Smith Logically Collective on TS 425fee21e36SBarry Smith 42615091d37SBarry Smith Input Parameters: 42715091d37SBarry Smith + ts - the timestep context 42815091d37SBarry Smith - inc - the scaling factor >= 1.0 42915091d37SBarry Smith 43082bf6240SBarry Smith Options Database Key: 43182bf6240SBarry Smith $ -ts_pseudo_increment <increment> 43282bf6240SBarry Smith 43315091d37SBarry Smith Level: advanced 43415091d37SBarry Smith 43582bf6240SBarry Smith .keywords: timestep, pseudo, set, increment 43682bf6240SBarry Smith 43782bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 43882bf6240SBarry Smith @*/ 4397087cfbeSBarry Smith PetscErrorCode TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc) 44082bf6240SBarry Smith { 4414ac538c5SBarry Smith PetscErrorCode ierr; 44282bf6240SBarry Smith 44382bf6240SBarry Smith PetscFunctionBegin; 4440700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 445c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts,inc,2); 4464ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSPseudoSetTimeStepIncrement_C",(TS,PetscReal),(ts,inc));CHKERRQ(ierr); 44782bf6240SBarry Smith PetscFunctionReturn(0); 44882bf6240SBarry Smith } 44982bf6240SBarry Smith 4504a2ae208SSatish Balay #undef __FUNCT__ 4514a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt" 45282bf6240SBarry Smith /*@ 45382bf6240SBarry Smith TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep 45482bf6240SBarry Smith is computed via the formula 45582bf6240SBarry Smith $ dt = initial_dt*initial_fnorm/current_fnorm 45682bf6240SBarry Smith rather than the default update, 45782bf6240SBarry Smith $ dt = current_dt*previous_fnorm/current_fnorm. 45882bf6240SBarry Smith 4593f9fe445SBarry Smith Logically Collective on TS 46015091d37SBarry Smith 46182bf6240SBarry Smith Input Parameter: 46282bf6240SBarry Smith . ts - the timestep context 46382bf6240SBarry Smith 46482bf6240SBarry Smith Options Database Key: 46582bf6240SBarry Smith $ -ts_pseudo_increment_dt_from_initial_dt 46682bf6240SBarry Smith 46715091d37SBarry Smith Level: advanced 46815091d37SBarry Smith 46982bf6240SBarry Smith .keywords: timestep, pseudo, set, increment 47082bf6240SBarry Smith 47182bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 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 48382bf6240SBarry Smith 4844a2ae208SSatish Balay #undef __FUNCT__ 4854a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep" 486ac226902SBarry Smith /*@C 48782bf6240SBarry Smith TSPseudoSetTimeStep - Sets the user-defined routine to be 48882bf6240SBarry Smith called at each pseudo-timestep to update the timestep. 48982bf6240SBarry Smith 4903f9fe445SBarry Smith Logically Collective on TS 49115091d37SBarry Smith 49282bf6240SBarry Smith Input Parameters: 49315091d37SBarry Smith + ts - timestep context 49482bf6240SBarry Smith . dt - function to compute timestep 49515091d37SBarry Smith - ctx - [optional] user-defined context for private data 49682bf6240SBarry Smith required by the function (may be PETSC_NULL) 49782bf6240SBarry Smith 49815091d37SBarry Smith Level: intermediate 499fee21e36SBarry Smith 50082bf6240SBarry Smith Calling sequence of func: 50187828ca2SBarry Smith . func (TS ts,PetscReal *newdt,void *ctx); 50282bf6240SBarry Smith 50382bf6240SBarry Smith . newdt - the newly computed timestep 50482bf6240SBarry Smith . ctx - [optional] timestep context 50582bf6240SBarry Smith 50682bf6240SBarry Smith Notes: 50782bf6240SBarry Smith The routine set here will be called by TSPseudoComputeTimeStep() 50882bf6240SBarry Smith during the timestepping process. 50982bf6240SBarry Smith 51082bf6240SBarry Smith .keywords: timestep, pseudo, set 51182bf6240SBarry Smith 51282bf6240SBarry Smith .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep() 51382bf6240SBarry Smith @*/ 5147087cfbeSBarry Smith PetscErrorCode TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void* ctx) 51582bf6240SBarry Smith { 5164ac538c5SBarry Smith PetscErrorCode ierr; 51782bf6240SBarry Smith 51882bf6240SBarry Smith PetscFunctionBegin; 5190700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 5204ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSPseudoSetTimeStep_C",(TS,PetscErrorCode (*)(TS,PetscReal *,void *),void *),(ts,dt,ctx));CHKERRQ(ierr); 52182bf6240SBarry Smith PetscFunctionReturn(0); 52282bf6240SBarry Smith } 52382bf6240SBarry Smith 52482bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 52582bf6240SBarry Smith 526ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscBool *); /* force argument to next function to not be extern C*/ 527fb2e594dSBarry Smith EXTERN_C_BEGIN 5284a2ae208SSatish Balay #undef __FUNCT__ 5294a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo" 5307087cfbeSBarry Smith PetscErrorCode TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void* ctx) 53182bf6240SBarry Smith { 53282bf6240SBarry Smith TS_Pseudo *pseudo; 53382bf6240SBarry Smith 53482bf6240SBarry Smith PetscFunctionBegin; 53582bf6240SBarry Smith pseudo = (TS_Pseudo*)ts->data; 53682bf6240SBarry Smith pseudo->verify = dt; 53782bf6240SBarry Smith pseudo->verifyctx = ctx; 53882bf6240SBarry Smith PetscFunctionReturn(0); 53982bf6240SBarry Smith } 540fb2e594dSBarry Smith EXTERN_C_END 54182bf6240SBarry Smith 542fb2e594dSBarry Smith EXTERN_C_BEGIN 5434a2ae208SSatish Balay #undef __FUNCT__ 5444a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement_Pseudo" 5457087cfbeSBarry Smith PetscErrorCode TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc) 54682bf6240SBarry Smith { 5474bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 54882bf6240SBarry Smith 54982bf6240SBarry Smith PetscFunctionBegin; 55082bf6240SBarry Smith pseudo->dt_increment = inc; 55182bf6240SBarry Smith PetscFunctionReturn(0); 55282bf6240SBarry Smith } 553fb2e594dSBarry Smith EXTERN_C_END 55482bf6240SBarry Smith 555fb2e594dSBarry Smith EXTERN_C_BEGIN 5564a2ae208SSatish Balay #undef __FUNCT__ 5574a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt_Pseudo" 5587087cfbeSBarry Smith PetscErrorCode TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts) 55982bf6240SBarry Smith { 5604bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 56182bf6240SBarry Smith 56282bf6240SBarry Smith PetscFunctionBegin; 5634bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_TRUE; 56482bf6240SBarry Smith PetscFunctionReturn(0); 56582bf6240SBarry Smith } 566fb2e594dSBarry Smith EXTERN_C_END 56782bf6240SBarry Smith 5686849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/ 569fb2e594dSBarry Smith EXTERN_C_BEGIN 5704a2ae208SSatish Balay #undef __FUNCT__ 5714a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo" 5727087cfbeSBarry Smith PetscErrorCode TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void* ctx) 57382bf6240SBarry Smith { 5744bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 57582bf6240SBarry Smith 57682bf6240SBarry Smith PetscFunctionBegin; 57782bf6240SBarry Smith pseudo->dt = dt; 57882bf6240SBarry Smith pseudo->dtctx = ctx; 57982bf6240SBarry Smith PetscFunctionReturn(0); 58082bf6240SBarry Smith } 581fb2e594dSBarry Smith EXTERN_C_END 58282bf6240SBarry Smith 58382bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 58410e6a065SJed Brown /*MC 58510e6a065SJed Brown TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping 58682bf6240SBarry Smith 58710e6a065SJed Brown This method solves equations of the form 58810e6a065SJed Brown 58910e6a065SJed Brown $ F(X,Xdot) = 0 59010e6a065SJed Brown 59110e6a065SJed Brown for steady state using the iteration 59210e6a065SJed Brown 59310e6a065SJed Brown $ [G'] S = -F(X,0) 59410e6a065SJed Brown $ X += S 59510e6a065SJed Brown 59610e6a065SJed Brown where 59710e6a065SJed Brown 59810e6a065SJed Brown $ G(Y) = F(Y,(Y-X)/dt) 59910e6a065SJed Brown 6006f2d6a7bSJed Brown This is linearly-implicit Euler with the residual always evaluated "at steady 6016f2d6a7bSJed Brown state". See note below. 60210e6a065SJed Brown 60310e6a065SJed Brown Options database keys: 60410e6a065SJed Brown + -ts_pseudo_increment <real> - ratio of increase dt 60510e6a065SJed Brown - -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt 60610e6a065SJed Brown 60710e6a065SJed Brown Level: beginner 60810e6a065SJed Brown 60910e6a065SJed Brown References: 61010e6a065SJed Brown Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential-Algebraic Equations, 2003. 61110e6a065SJed Brown C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998. 61210e6a065SJed Brown 61310e6a065SJed Brown Notes: 6146f2d6a7bSJed Brown The residual computed by this method includes the transient term (Xdot is computed instead of 6156f2d6a7bSJed Brown always being zero), but since the prediction from the last step is always the solution from the 6166f2d6a7bSJed Brown last step, on the first Newton iteration we have 6176f2d6a7bSJed Brown 6186f2d6a7bSJed Brown $ Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0 6196f2d6a7bSJed Brown 6206f2d6a7bSJed Brown Therefore, the linear system solved by the first Newton iteration is equivalent to the one 6216f2d6a7bSJed Brown described above and in the papers. If the user chooses to perform multiple Newton iterations, the 6226f2d6a7bSJed Brown algorithm is no longer the one described in the referenced papers. 62310e6a065SJed Brown 62410e6a065SJed Brown .seealso: TSCreate(), TS, TSSetType() 62510e6a065SJed Brown 62610e6a065SJed Brown M*/ 627fb2e594dSBarry Smith EXTERN_C_BEGIN 6284a2ae208SSatish Balay #undef __FUNCT__ 6294a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Pseudo" 6307087cfbeSBarry Smith PetscErrorCode TSCreate_Pseudo(TS ts) 6312d3f70b5SBarry Smith { 6327bf11e45SBarry Smith TS_Pseudo *pseudo; 633dfbe8321SBarry Smith PetscErrorCode ierr; 6342d3f70b5SBarry Smith 6353a40ed3dSBarry Smith PetscFunctionBegin; 636*277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Pseudo; 637000e7ae3SMatthew Knepley ts->ops->destroy = TSDestroy_Pseudo; 638000e7ae3SMatthew Knepley ts->ops->view = TSView_Pseudo; 6392d3f70b5SBarry Smith 64017186662SBarry Smith if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Only for nonlinear problems"); 641000e7ae3SMatthew Knepley ts->ops->setup = TSSetUp_Pseudo; 642000e7ae3SMatthew Knepley ts->ops->step = TSStep_Pseudo; 643000e7ae3SMatthew Knepley ts->ops->setfromoptions = TSSetFromOptions_Pseudo; 6440f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Pseudo; 6450f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Pseudo; 6467bf11e45SBarry Smith 6477bf11e45SBarry Smith /* create the required nonlinear solver context */ 6487adad957SLisandro Dalcin ierr = SNESCreate(((PetscObject)ts)->comm,&ts->snes);CHKERRQ(ierr); 64917efb626SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);CHKERRQ(ierr); 6502d3f70b5SBarry Smith 65138f2d2fdSLisandro Dalcin ierr = PetscNewLog(ts,TS_Pseudo,&pseudo);CHKERRQ(ierr); 6527bf11e45SBarry Smith ts->data = (void*)pseudo; 6532d3f70b5SBarry Smith 65428aa8177SBarry Smith pseudo->dt_increment = 1.1; 6554bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_FALSE; 65628aa8177SBarry Smith pseudo->dt = TSPseudoDefaultTimeStep; 65782bf6240SBarry Smith 658f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C", 659e1311b90SBarry Smith "TSPseudoSetVerifyTimeStep_Pseudo", 6600c97f09dSSatish Balay TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr); 661f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C", 662e1311b90SBarry Smith "TSPseudoSetTimeStepIncrement_Pseudo", 6630c97f09dSSatish Balay TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr); 664f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C", 665e1311b90SBarry Smith "TSPseudoIncrementDtFromInitialDt_Pseudo", 6660c97f09dSSatish Balay TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr); 6670c97f09dSSatish Balay ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C", 6680c97f09dSSatish Balay "TSPseudoSetTimeStep_Pseudo", 6690c97f09dSSatish Balay TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr); 6703a40ed3dSBarry Smith PetscFunctionReturn(0); 6712d3f70b5SBarry Smith } 672fb2e594dSBarry Smith EXTERN_C_END 6732d3f70b5SBarry Smith 6744a2ae208SSatish Balay #undef __FUNCT__ 6754a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultTimeStep" 67682bf6240SBarry Smith /*@C 67782bf6240SBarry Smith TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping. 67882bf6240SBarry Smith Use with TSPseudoSetTimeStep(). 67928aa8177SBarry Smith 68015091d37SBarry Smith Collective on TS 68115091d37SBarry Smith 68228aa8177SBarry Smith Input Parameters: 68328aa8177SBarry Smith . ts - the timestep context 68482bf6240SBarry Smith . dtctx - unused timestep context 68528aa8177SBarry Smith 68682bf6240SBarry Smith Output Parameter: 68782bf6240SBarry Smith . newdt - the timestep to use for the next step 68828aa8177SBarry Smith 68915091d37SBarry Smith Level: advanced 69015091d37SBarry Smith 69182bf6240SBarry Smith .keywords: timestep, pseudo, default 692564e8f4eSLois Curfman McInnes 69382bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep() 69428aa8177SBarry Smith @*/ 6957087cfbeSBarry Smith PetscErrorCode TSPseudoDefaultTimeStep(TS ts,PetscReal* newdt,void* dtctx) 69628aa8177SBarry Smith { 69782bf6240SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 69887828ca2SBarry Smith PetscReal inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous; 699dfbe8321SBarry Smith PetscErrorCode ierr; 70028aa8177SBarry Smith 7013a40ed3dSBarry Smith PetscFunctionBegin; 702bbd7b040SJed Brown ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr); 703bbd7b040SJed Brown ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func);CHKERRQ(ierr); 70482bf6240SBarry Smith ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 70582bf6240SBarry Smith if (pseudo->initial_fnorm == 0.0) { 70682bf6240SBarry Smith /* first time through so compute initial function norm */ 70782bf6240SBarry Smith pseudo->initial_fnorm = pseudo->fnorm; 70882bf6240SBarry Smith fnorm_previous = pseudo->fnorm; 70982bf6240SBarry Smith } 71082bf6240SBarry Smith if (pseudo->fnorm == 0.0) { 71182bf6240SBarry Smith *newdt = 1.e12*inc*ts->time_step; 712004884caSBarry Smith } else if (pseudo->increment_dt_from_initial_dt) { 71382bf6240SBarry Smith *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm; 71482bf6240SBarry Smith } else { 71582bf6240SBarry Smith *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm; 71682bf6240SBarry Smith } 71782bf6240SBarry Smith pseudo->fnorm_previous = pseudo->fnorm; 7183a40ed3dSBarry Smith PetscFunctionReturn(0); 71928aa8177SBarry Smith } 720