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 48fb4a63b6SLois Curfman McInnes .keywords: timestep, pseudo, compute 49564e8f4eSLois Curfman McInnes 508d359177SBarry Smith .seealso: TSPseudoTimeStepDefault(), 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 /* ------------------------------------------------------------------------------*/ 667bf11e45SBarry Smith /*@C 678d359177SBarry Smith TSPseudoVerifyTimeStepDefault - Default code to verify the quality of the last timestep. 687bf11e45SBarry Smith 6915091d37SBarry Smith Collective on TS 7015091d37SBarry Smith 717bf11e45SBarry Smith Input Parameters: 7215091d37SBarry Smith + ts - the timestep context 737bf11e45SBarry Smith . dtctx - unused timestep context 7415091d37SBarry Smith - update - latest solution vector 757bf11e45SBarry Smith 76564e8f4eSLois Curfman McInnes Output Parameters: 7715091d37SBarry Smith + newdt - the timestep to use for the next step 7815091d37SBarry Smith - flag - flag indicating whether the last time step was acceptable 797bf11e45SBarry Smith 8015091d37SBarry Smith Level: advanced 81fee21e36SBarry Smith 82564e8f4eSLois Curfman McInnes Note: 83564e8f4eSLois Curfman McInnes This routine always returns a flag of 1, indicating an acceptable 84564e8f4eSLois Curfman McInnes timestep. 85564e8f4eSLois Curfman McInnes 86564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, default, verify 87564e8f4eSLois Curfman McInnes 88564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStep() 897bf11e45SBarry Smith @*/ 908d359177SBarry Smith PetscErrorCode TSPseudoVerifyTimeStepDefault(TS ts,Vec update,void *dtctx,PetscReal *newdt,PetscBool *flag) 917bf11e45SBarry Smith { 923a40ed3dSBarry Smith PetscFunctionBegin; 93a7cc72afSBarry Smith *flag = PETSC_TRUE; 943a40ed3dSBarry Smith PetscFunctionReturn(0); 957bf11e45SBarry Smith } 967bf11e45SBarry Smith 977bf11e45SBarry Smith 987bf11e45SBarry Smith /*@ 99564e8f4eSLois Curfman McInnes TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable. 1007bf11e45SBarry Smith 10115091d37SBarry Smith Collective on TS 10215091d37SBarry Smith 103fb4a63b6SLois Curfman McInnes Input Parameters: 10415091d37SBarry Smith + ts - timestep context 10515091d37SBarry Smith - update - latest solution vector 1067bf11e45SBarry Smith 107fb4a63b6SLois Curfman McInnes Output Parameters: 10815091d37SBarry Smith + dt - newly computed timestep (if it had to shrink) 10915091d37SBarry Smith - flag - indicates if current timestep was ok 1107bf11e45SBarry Smith 11115091d37SBarry Smith Level: advanced 112fee21e36SBarry Smith 113564e8f4eSLois Curfman McInnes Notes: 114564e8f4eSLois Curfman McInnes The routine to be called here to compute the timestep should be 115564e8f4eSLois Curfman McInnes set by calling TSPseudoSetVerifyTimeStep(). 116564e8f4eSLois Curfman McInnes 117564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, verify 118564e8f4eSLois Curfman McInnes 1198d359177SBarry Smith .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStepDefault() 1207bf11e45SBarry Smith @*/ 1217087cfbeSBarry Smith PetscErrorCode TSPseudoVerifyTimeStep(TS ts,Vec update,PetscReal *dt,PetscBool *flag) 1227bf11e45SBarry Smith { 1237bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 124dfbe8321SBarry Smith PetscErrorCode ierr; 1257bf11e45SBarry Smith 1263a40ed3dSBarry Smith PetscFunctionBegin; 127cb9d8021SPierre Barbier de Reuille *flag = PETSC_TRUE; 128be5899b3SLisandro Dalcin if(pseudo->verify) { 1297bf11e45SBarry Smith ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag);CHKERRQ(ierr); 130cb9d8021SPierre Barbier de Reuille } 1313a40ed3dSBarry Smith PetscFunctionReturn(0); 1327bf11e45SBarry Smith } 1337bf11e45SBarry Smith 1347bf11e45SBarry Smith /* --------------------------------------------------------------------------------*/ 1357bf11e45SBarry Smith 136193ac0bcSJed Brown static PetscErrorCode TSStep_Pseudo(TS ts) 1372d3f70b5SBarry Smith { 138277b19d0SLisandro Dalcin TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 139be5899b3SLisandro Dalcin PetscInt nits,lits,reject; 140cdbf8f93SLisandro Dalcin PetscBool stepok; 141be5899b3SLisandro Dalcin PetscReal next_time_step = ts->time_step; 142dfbe8321SBarry Smith PetscErrorCode ierr; 1432d3f70b5SBarry Smith 1443a40ed3dSBarry Smith PetscFunctionBegin; 145bbd56ea5SKarl Rupp if (ts->steps == 0) pseudo->dt_initial = ts->time_step; 146193ac0bcSJed Brown ierr = VecCopy(ts->vec_sol,pseudo->update);CHKERRQ(ierr); 147cdbf8f93SLisandro Dalcin ierr = TSPseudoComputeTimeStep(ts,&next_time_step);CHKERRQ(ierr); 148cdbf8f93SLisandro Dalcin for (reject=0; reject<ts->max_reject; reject++,ts->reject++) { 149cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 150b8123daeSJed Brown ierr = TSPreStage(ts,ts->ptime+ts->time_step);CHKERRQ(ierr); 1510298fd71SBarry Smith ierr = SNESSolve(ts->snes,NULL,pseudo->update);CHKERRQ(ierr); 152be5899b3SLisandro Dalcin ierr = SNESGetIterationNumber(ts->snes,&nits);CHKERRQ(ierr); 153b850b91aSLisandro Dalcin ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 154be5899b3SLisandro Dalcin ts->snes_its += nits; ts->ksp_its += lits; 1559be3e283SDebojyoti Ghosh ierr = TSPostStage(ts,ts->ptime+ts->time_step,0,&(pseudo->update));CHKERRQ(ierr); 156be5899b3SLisandro Dalcin ierr = TSAdaptCheckStage(ts->adapt,ts,ts->ptime+ts->time_step,pseudo->update,&stepok);CHKERRQ(ierr); 157be5899b3SLisandro Dalcin if (!stepok) {next_time_step = ts->time_step; continue;} 158193ac0bcSJed Brown pseudo->fnorm = -1; /* The current norm is no longer valid, monitor must recompute it. */ 159cdbf8f93SLisandro Dalcin ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&next_time_step,&stepok);CHKERRQ(ierr); 160cdbf8f93SLisandro Dalcin if (stepok) break; 161cdbf8f93SLisandro Dalcin } 162be5899b3SLisandro Dalcin if (reject >= ts->max_reject) { 163be5899b3SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 164be5899b3SLisandro Dalcin ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,reject);CHKERRQ(ierr); 165cdbf8f93SLisandro Dalcin PetscFunctionReturn(0); 1667bf11e45SBarry Smith } 167be5899b3SLisandro Dalcin 168be5899b3SLisandro Dalcin ierr = VecCopy(pseudo->update,ts->vec_sol);CHKERRQ(ierr); 169be5899b3SLisandro Dalcin ts->ptime += ts->time_step; 170be5899b3SLisandro Dalcin ts->time_step = next_time_step; 171be5899b3SLisandro Dalcin 1723118ae5eSBarry Smith if (pseudo->fnorm < 0) { 1733118ae5eSBarry Smith ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr); 1743118ae5eSBarry Smith ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr); 1753118ae5eSBarry Smith ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 1763118ae5eSBarry Smith } 1773118ae5eSBarry Smith if (pseudo->fnorm < pseudo->fatol) { 1783118ae5eSBarry Smith ts->reason = TS_CONVERGED_PSEUDO_FATOL; 179be5899b3SLisandro Dalcin ierr = PetscInfo3(ts,"Step=%D, converged since fnorm %g < fatol %g\n",ts->steps,pseudo->fnorm,pseudo->frtol);CHKERRQ(ierr); 1803118ae5eSBarry Smith PetscFunctionReturn(0); 1813118ae5eSBarry Smith } 1823118ae5eSBarry Smith if (pseudo->fnorm/pseudo->fnorm_initial < pseudo->frtol) { 1833118ae5eSBarry Smith ts->reason = TS_CONVERGED_PSEUDO_FRTOL; 184be5899b3SLisandro Dalcin ierr = PetscInfo4(ts,"Step=%D, converged since fnorm %g / fnorm_initial %g < frtol %g\n",ts->steps,pseudo->fnorm,pseudo->fnorm_initial,pseudo->fatol);CHKERRQ(ierr); 1853118ae5eSBarry Smith PetscFunctionReturn(0); 1863118ae5eSBarry Smith } 1873a40ed3dSBarry Smith PetscFunctionReturn(0); 1882d3f70b5SBarry Smith } 1892d3f70b5SBarry Smith 1902d3f70b5SBarry Smith /*------------------------------------------------------------*/ 191277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Pseudo(TS ts) 1922d3f70b5SBarry Smith { 1937bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 194dfbe8321SBarry Smith PetscErrorCode ierr; 1952d3f70b5SBarry Smith 1963a40ed3dSBarry Smith PetscFunctionBegin; 1976bf464f9SBarry Smith ierr = VecDestroy(&pseudo->update);CHKERRQ(ierr); 1986bf464f9SBarry Smith ierr = VecDestroy(&pseudo->func);CHKERRQ(ierr); 1996bf464f9SBarry Smith ierr = VecDestroy(&pseudo->xdot);CHKERRQ(ierr); 2003a40ed3dSBarry Smith PetscFunctionReturn(0); 2012d3f70b5SBarry Smith } 2022d3f70b5SBarry Smith 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); 210bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",NULL);CHKERRQ(ierr); 211bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",NULL);CHKERRQ(ierr); 212bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetMaxTimeStep_C",NULL);CHKERRQ(ierr); 213bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",NULL);CHKERRQ(ierr); 214bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",NULL);CHKERRQ(ierr); 215277b19d0SLisandro Dalcin PetscFunctionReturn(0); 216277b19d0SLisandro Dalcin } 2172d3f70b5SBarry Smith 2182d3f70b5SBarry Smith /*------------------------------------------------------------*/ 2192d3f70b5SBarry Smith 2206f2d6a7bSJed Brown /* 2216f2d6a7bSJed Brown Compute Xdot = (X^{n+1}-X^n)/dt) = 0 2226f2d6a7bSJed Brown */ 2236f2d6a7bSJed Brown static PetscErrorCode TSPseudoGetXdot(TS ts,Vec X,Vec *Xdot) 2242d3f70b5SBarry Smith { 2256f2d6a7bSJed Brown TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 226193ac0bcSJed Brown const PetscScalar mdt = 1.0/ts->time_step,*xnp1,*xn; 227193ac0bcSJed Brown PetscScalar *xdot; 228dfbe8321SBarry Smith PetscErrorCode ierr; 229a7cc72afSBarry Smith PetscInt i,n; 2302d3f70b5SBarry Smith 2313a40ed3dSBarry Smith PetscFunctionBegin; 232*aab5bcd8SJed Brown *Xdot = NULL; 233193ac0bcSJed Brown ierr = VecGetArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr); 234193ac0bcSJed Brown ierr = VecGetArrayRead(X,&xnp1);CHKERRQ(ierr); 2356f2d6a7bSJed Brown ierr = VecGetArray(pseudo->xdot,&xdot);CHKERRQ(ierr); 2366f2d6a7bSJed Brown ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr); 237bbd56ea5SKarl Rupp for (i=0; i<n; i++) xdot[i] = mdt*(xnp1[i] - xn[i]); 238193ac0bcSJed Brown ierr = VecRestoreArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr); 239193ac0bcSJed Brown ierr = VecRestoreArrayRead(X,&xnp1);CHKERRQ(ierr); 2406f2d6a7bSJed Brown ierr = VecRestoreArray(pseudo->xdot,&xdot);CHKERRQ(ierr); 2416f2d6a7bSJed Brown *Xdot = pseudo->xdot; 2423a40ed3dSBarry Smith PetscFunctionReturn(0); 2432d3f70b5SBarry Smith } 2442d3f70b5SBarry Smith 2456f2d6a7bSJed Brown /* 2466f2d6a7bSJed Brown The transient residual is 2476f2d6a7bSJed Brown 2486f2d6a7bSJed Brown F(U^{n+1},(U^{n+1}-U^n)/dt) = 0 2496f2d6a7bSJed Brown 2506f2d6a7bSJed Brown or for ODE, 2516f2d6a7bSJed Brown 2526f2d6a7bSJed Brown (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0 2536f2d6a7bSJed Brown 2546f2d6a7bSJed Brown This is the function that must be evaluated for transient simulation and for 2556f2d6a7bSJed Brown finite difference Jacobians. On the first Newton step, this algorithm uses 2566f2d6a7bSJed Brown a guess of U^{n+1} = U^n in which case the transient term vanishes and the 2576f2d6a7bSJed Brown residual is actually the steady state residual. Pseudotransient 2586f2d6a7bSJed Brown continuation as described in the literature is a linearly implicit 2596f2d6a7bSJed Brown algorithm, it only takes this one Newton step with the steady state 2606f2d6a7bSJed Brown residual, and then advances to the next time step. 2616f2d6a7bSJed Brown */ 2620f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes,Vec X,Vec Y,TS ts) 2632d3f70b5SBarry Smith { 2646f2d6a7bSJed Brown Vec Xdot; 265dfbe8321SBarry Smith PetscErrorCode ierr; 2662d3f70b5SBarry Smith 2673a40ed3dSBarry Smith PetscFunctionBegin; 2686f2d6a7bSJed Brown ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr); 269193ac0bcSJed Brown ierr = TSComputeIFunction(ts,ts->ptime+ts->time_step,X,Xdot,Y,PETSC_FALSE);CHKERRQ(ierr); 2706f2d6a7bSJed Brown PetscFunctionReturn(0); 2716f2d6a7bSJed Brown } 2722d3f70b5SBarry Smith 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 */ 282d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes,Vec X,Mat AA,Mat BB,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); 289d1e9a80fSBarry Smith ierr = TSComputeIJacobian(ts,ts->ptime+ts->time_step,X,Xdot,1./ts->time_step,AA,BB,PETSC_FALSE);CHKERRQ(ierr); 2903a40ed3dSBarry Smith PetscFunctionReturn(0); 2912d3f70b5SBarry Smith } 2922d3f70b5SBarry Smith 2932d3f70b5SBarry Smith 2946849ba73SBarry Smith static PetscErrorCode TSSetUp_Pseudo(TS ts) 2952d3f70b5SBarry Smith { 2967bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 297dfbe8321SBarry Smith PetscErrorCode ierr; 2982d3f70b5SBarry Smith 2993a40ed3dSBarry Smith PetscFunctionBegin; 3007bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr); 3017bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr); 3026f2d6a7bSJed Brown ierr = VecDuplicate(ts->vec_sol,&pseudo->xdot);CHKERRQ(ierr); 3033a40ed3dSBarry Smith PetscFunctionReturn(0); 3042d3f70b5SBarry Smith } 3052d3f70b5SBarry Smith /*------------------------------------------------------------*/ 3062d3f70b5SBarry Smith 307560360afSLisandro Dalcin static PetscErrorCode TSPseudoMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy) 3082d3f70b5SBarry Smith { 3097bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 310dfbe8321SBarry Smith PetscErrorCode ierr; 311ce94432eSBarry Smith PetscViewer viewer = (PetscViewer) dummy; 3122d3f70b5SBarry Smith 3133a40ed3dSBarry Smith PetscFunctionBegin; 314193ac0bcSJed Brown if (pseudo->fnorm < 0) { /* The last computed norm is stale, recompute */ 315193ac0bcSJed Brown ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr); 316193ac0bcSJed Brown ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr); 317193ac0bcSJed Brown ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 318193ac0bcSJed Brown } 319649052a6SBarry Smith ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 3207c8652ddSBarry 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); 321649052a6SBarry Smith ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 3223a40ed3dSBarry Smith PetscFunctionReturn(0); 3232d3f70b5SBarry Smith } 3242d3f70b5SBarry Smith 3254416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Pseudo(PetscOptionItems *PetscOptionsObject,TS ts) 3262d3f70b5SBarry Smith { 3274bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 328dfbe8321SBarry Smith PetscErrorCode ierr; 329ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE; 330649052a6SBarry Smith PetscViewer viewer; 3312d3f70b5SBarry Smith 3323a40ed3dSBarry Smith PetscFunctionBegin; 333e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Pseudo-timestepping options");CHKERRQ(ierr); 334560360afSLisandro Dalcin ierr = PetscOptionsBool("-ts_monitor_pseudo","Monitor convergence","",flg,&flg,NULL);CHKERRQ(ierr); 3352d3f70b5SBarry Smith if (flg) { 336ce94432eSBarry Smith ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ts),"stdout",&viewer);CHKERRQ(ierr); 337649052a6SBarry Smith ierr = TSMonitorSet(ts,TSPseudoMonitorDefault,viewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 33828aa8177SBarry Smith } 339be5899b3SLisandro Dalcin flg = pseudo->increment_dt_from_initial_dt; 3400298fd71SBarry Smith ierr = PetscOptionsBool("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",flg,&flg,NULL);CHKERRQ(ierr); 341be5899b3SLisandro Dalcin pseudo->increment_dt_from_initial_dt = flg; 34294ae4db5SBarry Smith ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,NULL);CHKERRQ(ierr); 34394ae4db5SBarry Smith ierr = PetscOptionsReal("-ts_pseudo_max_dt","Maximum value for dt","TSPseudoSetMaxTimeStep",pseudo->dt_max,&pseudo->dt_max,NULL);CHKERRQ(ierr); 3443118ae5eSBarry Smith ierr = PetscOptionsReal("-ts_pseudo_fatol","Tolerance for norm of function","",pseudo->fatol,&pseudo->fatol,NULL);CHKERRQ(ierr); 3453118ae5eSBarry Smith ierr = PetscOptionsReal("-ts_pseudo_frtol","Relative tolerance for norm of function","",pseudo->frtol,&pseudo->frtol,NULL);CHKERRQ(ierr); 346b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 3473a40ed3dSBarry Smith PetscFunctionReturn(0); 3482d3f70b5SBarry Smith } 3492d3f70b5SBarry Smith 3506849ba73SBarry Smith static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer) 3512d3f70b5SBarry Smith { 352d52bd9f3SBarry Smith PetscErrorCode ierr; 3533118ae5eSBarry Smith PetscBool isascii; 354d52bd9f3SBarry Smith 3553a40ed3dSBarry Smith PetscFunctionBegin; 3563118ae5eSBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 3573118ae5eSBarry Smith if (isascii) { 3583118ae5eSBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 3593118ae5eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," frtol - relative tolerance in function value %g\n",(double)pseudo->frtol);CHKERRQ(ierr); 3603118ae5eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," fatol - absolute tolerance in function value %g\n",(double)pseudo->fatol);CHKERRQ(ierr); 3613118ae5eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," dt_initial - initial timestep %g\n",(double)pseudo->dt_initial);CHKERRQ(ierr); 3623118ae5eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," dt_increment - increase in timestep on successful step %g\n",(double)pseudo->dt_increment);CHKERRQ(ierr); 3633118ae5eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," dt_max - maximum time %g\n",(double)pseudo->dt_max);CHKERRQ(ierr); 3643118ae5eSBarry Smith } 3653a40ed3dSBarry Smith PetscFunctionReturn(0); 3662d3f70b5SBarry Smith } 3672d3f70b5SBarry Smith 36882bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 369ac226902SBarry Smith /*@C 37082bf6240SBarry Smith TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the 37182bf6240SBarry Smith last timestep. 37282bf6240SBarry Smith 3733f9fe445SBarry Smith Logically Collective on TS 37415091d37SBarry Smith 37582bf6240SBarry Smith Input Parameters: 37615091d37SBarry Smith + ts - timestep context 37782bf6240SBarry Smith . dt - user-defined function to verify timestep 37815091d37SBarry Smith - ctx - [optional] user-defined context for private data 3790298fd71SBarry Smith for the timestep verification routine (may be NULL) 38082bf6240SBarry Smith 38115091d37SBarry Smith Level: advanced 382fee21e36SBarry Smith 38382bf6240SBarry Smith Calling sequence of func: 384ace3abfcSBarry Smith . func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscBool *flag); 38582bf6240SBarry Smith 38682bf6240SBarry Smith . update - latest solution vector 38782bf6240SBarry Smith . ctx - [optional] timestep context 38882bf6240SBarry Smith . newdt - the timestep to use for the next step 38982bf6240SBarry Smith . flag - flag indicating whether the last time step was acceptable 39082bf6240SBarry Smith 39182bf6240SBarry Smith Notes: 39282bf6240SBarry Smith The routine set here will be called by TSPseudoVerifyTimeStep() 39382bf6240SBarry Smith during the timestepping process. 39482bf6240SBarry Smith 39582bf6240SBarry Smith .keywords: timestep, pseudo, set, verify 39682bf6240SBarry Smith 3978d359177SBarry Smith .seealso: TSPseudoVerifyTimeStepDefault(), TSPseudoVerifyTimeStep() 39882bf6240SBarry Smith @*/ 3997087cfbeSBarry Smith PetscErrorCode TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,PetscBool*),void *ctx) 40082bf6240SBarry Smith { 4014ac538c5SBarry Smith PetscErrorCode ierr; 40282bf6240SBarry Smith 40382bf6240SBarry Smith PetscFunctionBegin; 4040700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 4054ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSPseudoSetVerifyTimeStep_C",(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal*,PetscBool*),void*),(ts,dt,ctx));CHKERRQ(ierr); 40682bf6240SBarry Smith PetscFunctionReturn(0); 40782bf6240SBarry Smith } 40882bf6240SBarry Smith 40982bf6240SBarry Smith /*@ 41082bf6240SBarry Smith TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 4118d359177SBarry Smith dt when using the TSPseudoTimeStepDefault() routine. 41282bf6240SBarry Smith 4133f9fe445SBarry Smith Logically Collective on TS 414fee21e36SBarry Smith 41515091d37SBarry Smith Input Parameters: 41615091d37SBarry Smith + ts - the timestep context 41715091d37SBarry Smith - inc - the scaling factor >= 1.0 41815091d37SBarry Smith 41982bf6240SBarry Smith Options Database Key: 420e1bc860dSBarry Smith . -ts_pseudo_increment <increment> 42182bf6240SBarry Smith 42215091d37SBarry Smith Level: advanced 42315091d37SBarry Smith 42482bf6240SBarry Smith .keywords: timestep, pseudo, set, increment 42582bf6240SBarry Smith 4268d359177SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault() 42782bf6240SBarry Smith @*/ 4287087cfbeSBarry Smith PetscErrorCode TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc) 42982bf6240SBarry Smith { 4304ac538c5SBarry Smith PetscErrorCode ierr; 43182bf6240SBarry Smith 43282bf6240SBarry Smith PetscFunctionBegin; 4330700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 434c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts,inc,2); 4354ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSPseudoSetTimeStepIncrement_C",(TS,PetscReal),(ts,inc));CHKERRQ(ierr); 43682bf6240SBarry Smith PetscFunctionReturn(0); 43782bf6240SBarry Smith } 43882bf6240SBarry Smith 43986534af1SJed Brown /*@ 44086534af1SJed Brown TSPseudoSetMaxTimeStep - Sets the maximum time step 4418d359177SBarry Smith when using the TSPseudoTimeStepDefault() routine. 44286534af1SJed Brown 44386534af1SJed Brown Logically Collective on TS 44486534af1SJed Brown 44586534af1SJed Brown Input Parameters: 44686534af1SJed Brown + ts - the timestep context 44786534af1SJed Brown - maxdt - the maximum time step, use a non-positive value to deactivate 44886534af1SJed Brown 44986534af1SJed Brown Options Database Key: 450e1bc860dSBarry Smith . -ts_pseudo_max_dt <increment> 45186534af1SJed Brown 45286534af1SJed Brown Level: advanced 45386534af1SJed Brown 45486534af1SJed Brown .keywords: timestep, pseudo, set 45586534af1SJed Brown 4568d359177SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault() 45786534af1SJed Brown @*/ 45886534af1SJed Brown PetscErrorCode TSPseudoSetMaxTimeStep(TS ts,PetscReal maxdt) 45986534af1SJed Brown { 46086534af1SJed Brown PetscErrorCode ierr; 46186534af1SJed Brown 46286534af1SJed Brown PetscFunctionBegin; 46386534af1SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 46486534af1SJed Brown PetscValidLogicalCollectiveReal(ts,maxdt,2); 46586534af1SJed Brown ierr = PetscTryMethod(ts,"TSPseudoSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));CHKERRQ(ierr); 46686534af1SJed Brown PetscFunctionReturn(0); 46786534af1SJed Brown } 46886534af1SJed Brown 46982bf6240SBarry Smith /*@ 47082bf6240SBarry Smith TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep 47182bf6240SBarry Smith is computed via the formula 47282bf6240SBarry Smith $ dt = initial_dt*initial_fnorm/current_fnorm 47382bf6240SBarry Smith rather than the default update, 47482bf6240SBarry Smith $ dt = current_dt*previous_fnorm/current_fnorm. 47582bf6240SBarry Smith 4763f9fe445SBarry Smith Logically Collective on TS 47715091d37SBarry Smith 47882bf6240SBarry Smith Input Parameter: 47982bf6240SBarry Smith . ts - the timestep context 48082bf6240SBarry Smith 48182bf6240SBarry Smith Options Database Key: 482e1bc860dSBarry Smith . -ts_pseudo_increment_dt_from_initial_dt 48382bf6240SBarry Smith 48415091d37SBarry Smith Level: advanced 48515091d37SBarry Smith 48682bf6240SBarry Smith .keywords: timestep, pseudo, set, increment 48782bf6240SBarry Smith 4888d359177SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault() 48982bf6240SBarry Smith @*/ 4907087cfbeSBarry Smith PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS ts) 49182bf6240SBarry Smith { 4924ac538c5SBarry Smith PetscErrorCode ierr; 49382bf6240SBarry Smith 49482bf6240SBarry Smith PetscFunctionBegin; 4950700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 4964ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSPseudoIncrementDtFromInitialDt_C",(TS),(ts));CHKERRQ(ierr); 49782bf6240SBarry Smith PetscFunctionReturn(0); 49882bf6240SBarry Smith } 49982bf6240SBarry Smith 50082bf6240SBarry Smith 501ac226902SBarry Smith /*@C 50282bf6240SBarry Smith TSPseudoSetTimeStep - Sets the user-defined routine to be 50382bf6240SBarry Smith called at each pseudo-timestep to update the timestep. 50482bf6240SBarry Smith 5053f9fe445SBarry Smith Logically Collective on TS 50615091d37SBarry Smith 50782bf6240SBarry Smith Input Parameters: 50815091d37SBarry Smith + ts - timestep context 50982bf6240SBarry Smith . dt - function to compute timestep 51015091d37SBarry Smith - ctx - [optional] user-defined context for private data 5110298fd71SBarry Smith required by the function (may be NULL) 51282bf6240SBarry Smith 51315091d37SBarry Smith Level: intermediate 514fee21e36SBarry Smith 51582bf6240SBarry Smith Calling sequence of func: 51687828ca2SBarry Smith . func (TS ts,PetscReal *newdt,void *ctx); 51782bf6240SBarry Smith 51882bf6240SBarry Smith . newdt - the newly computed timestep 51982bf6240SBarry Smith . ctx - [optional] timestep context 52082bf6240SBarry Smith 52182bf6240SBarry Smith Notes: 52282bf6240SBarry Smith The routine set here will be called by TSPseudoComputeTimeStep() 52382bf6240SBarry Smith during the timestepping process. 5248d359177SBarry Smith If not set then TSPseudoTimeStepDefault() is automatically used 52582bf6240SBarry Smith 52682bf6240SBarry Smith .keywords: timestep, pseudo, set 52782bf6240SBarry Smith 5288d359177SBarry Smith .seealso: TSPseudoTimeStepDefault(), TSPseudoComputeTimeStep() 52982bf6240SBarry Smith @*/ 5307087cfbeSBarry Smith PetscErrorCode TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void *ctx) 53182bf6240SBarry Smith { 5324ac538c5SBarry Smith PetscErrorCode ierr; 53382bf6240SBarry Smith 53482bf6240SBarry Smith PetscFunctionBegin; 5350700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 5364ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSPseudoSetTimeStep_C",(TS,PetscErrorCode (*)(TS,PetscReal*,void*),void*),(ts,dt,ctx));CHKERRQ(ierr); 53782bf6240SBarry Smith PetscFunctionReturn(0); 53882bf6240SBarry Smith } 53982bf6240SBarry Smith 54082bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 54182bf6240SBarry Smith 542ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscBool*); /* force argument to next function to not be extern C*/ 543560360afSLisandro Dalcin static PetscErrorCode TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void *ctx) 54482bf6240SBarry Smith { 545be5899b3SLisandro Dalcin TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 54682bf6240SBarry Smith 54782bf6240SBarry Smith PetscFunctionBegin; 54882bf6240SBarry Smith pseudo->verify = dt; 54982bf6240SBarry Smith pseudo->verifyctx = ctx; 55082bf6240SBarry Smith PetscFunctionReturn(0); 55182bf6240SBarry Smith } 55282bf6240SBarry Smith 553560360afSLisandro Dalcin static PetscErrorCode TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc) 55482bf6240SBarry Smith { 5554bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 55682bf6240SBarry Smith 55782bf6240SBarry Smith PetscFunctionBegin; 55882bf6240SBarry Smith pseudo->dt_increment = inc; 55982bf6240SBarry Smith PetscFunctionReturn(0); 56082bf6240SBarry Smith } 56182bf6240SBarry Smith 562560360afSLisandro Dalcin static PetscErrorCode TSPseudoSetMaxTimeStep_Pseudo(TS ts,PetscReal maxdt) 56386534af1SJed Brown { 56486534af1SJed Brown TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 56586534af1SJed Brown 56686534af1SJed Brown PetscFunctionBegin; 56786534af1SJed Brown pseudo->dt_max = maxdt; 56886534af1SJed Brown PetscFunctionReturn(0); 56986534af1SJed Brown } 57086534af1SJed Brown 571560360afSLisandro Dalcin static PetscErrorCode TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts) 57282bf6240SBarry Smith { 5734bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 57482bf6240SBarry Smith 57582bf6240SBarry Smith PetscFunctionBegin; 5764bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_TRUE; 57782bf6240SBarry Smith PetscFunctionReturn(0); 57882bf6240SBarry Smith } 57982bf6240SBarry Smith 5806849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/ 581560360afSLisandro Dalcin static PetscErrorCode TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void *ctx) 58282bf6240SBarry Smith { 5834bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 58482bf6240SBarry Smith 58582bf6240SBarry Smith PetscFunctionBegin; 58682bf6240SBarry Smith pseudo->dt = dt; 58782bf6240SBarry Smith pseudo->dtctx = ctx; 58882bf6240SBarry Smith PetscFunctionReturn(0); 58982bf6240SBarry Smith } 59082bf6240SBarry Smith 59182bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 59210e6a065SJed Brown /*MC 59310e6a065SJed Brown TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping 59482bf6240SBarry Smith 59510e6a065SJed Brown This method solves equations of the form 59610e6a065SJed Brown 59710e6a065SJed Brown $ F(X,Xdot) = 0 59810e6a065SJed Brown 59910e6a065SJed Brown for steady state using the iteration 60010e6a065SJed Brown 60110e6a065SJed Brown $ [G'] S = -F(X,0) 60210e6a065SJed Brown $ X += S 60310e6a065SJed Brown 60410e6a065SJed Brown where 60510e6a065SJed Brown 60610e6a065SJed Brown $ G(Y) = F(Y,(Y-X)/dt) 60710e6a065SJed Brown 6086f2d6a7bSJed Brown This is linearly-implicit Euler with the residual always evaluated "at steady 6096f2d6a7bSJed Brown state". See note below. 61010e6a065SJed Brown 61110e6a065SJed Brown Options database keys: 61210e6a065SJed Brown + -ts_pseudo_increment <real> - ratio of increase dt 6133118ae5eSBarry Smith . -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt 6143118ae5eSBarry Smith . -ts_pseudo_fatol <atol> - stop iterating when the function norm is less than atol 6153118ae5eSBarry Smith - -ts_pseudo_frtol <rtol> - stop iterating when the function norm divided by the initial function norm is less than rtol 61610e6a065SJed Brown 61710e6a065SJed Brown Level: beginner 61810e6a065SJed Brown 61910e6a065SJed Brown References: 62096a0c994SBarry Smith + 1. - Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential Algebraic Equations, 2003. 62196a0c994SBarry Smith - 2. - C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998. 62210e6a065SJed Brown 62310e6a065SJed Brown Notes: 6246f2d6a7bSJed Brown The residual computed by this method includes the transient term (Xdot is computed instead of 6256f2d6a7bSJed Brown always being zero), but since the prediction from the last step is always the solution from the 6266f2d6a7bSJed Brown last step, on the first Newton iteration we have 6276f2d6a7bSJed Brown 6286f2d6a7bSJed Brown $ Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0 6296f2d6a7bSJed Brown 6306f2d6a7bSJed Brown Therefore, the linear system solved by the first Newton iteration is equivalent to the one 6316f2d6a7bSJed Brown described above and in the papers. If the user chooses to perform multiple Newton iterations, the 6326f2d6a7bSJed Brown algorithm is no longer the one described in the referenced papers. 63310e6a065SJed Brown 63410e6a065SJed Brown .seealso: TSCreate(), TS, TSSetType() 63510e6a065SJed Brown 63610e6a065SJed Brown M*/ 6378cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Pseudo(TS ts) 6382d3f70b5SBarry Smith { 6397bf11e45SBarry Smith TS_Pseudo *pseudo; 640dfbe8321SBarry Smith PetscErrorCode ierr; 641193ac0bcSJed Brown SNES snes; 64219fd82e9SBarry Smith SNESType stype; 6432d3f70b5SBarry Smith 6443a40ed3dSBarry Smith PetscFunctionBegin; 645277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Pseudo; 646000e7ae3SMatthew Knepley ts->ops->destroy = TSDestroy_Pseudo; 647000e7ae3SMatthew Knepley ts->ops->view = TSView_Pseudo; 648000e7ae3SMatthew Knepley ts->ops->setup = TSSetUp_Pseudo; 649000e7ae3SMatthew Knepley ts->ops->step = TSStep_Pseudo; 650000e7ae3SMatthew Knepley ts->ops->setfromoptions = TSSetFromOptions_Pseudo; 6510f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Pseudo; 6520f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Pseudo; 6532ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 6547bf11e45SBarry Smith 655193ac0bcSJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 656193ac0bcSJed Brown ierr = SNESGetType(snes,&stype);CHKERRQ(ierr); 657193ac0bcSJed Brown if (!stype) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);} 6582d3f70b5SBarry Smith 659b00a9115SJed Brown ierr = PetscNewLog(ts,&pseudo);CHKERRQ(ierr); 6607bf11e45SBarry Smith ts->data = (void*)pseudo; 6612d3f70b5SBarry Smith 662be5899b3SLisandro Dalcin pseudo->dt = TSPseudoTimeStepDefault; 663be5899b3SLisandro Dalcin pseudo->dtctx = NULL; 66428aa8177SBarry Smith pseudo->dt_increment = 1.1; 6654bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_FALSE; 666193ac0bcSJed Brown pseudo->fnorm = -1; 667be5899b3SLisandro Dalcin pseudo->fnorm_initial = -1; 668be5899b3SLisandro Dalcin pseudo->fnorm_previous = -1; 6693118ae5eSBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 6703118ae5eSBarry Smith pseudo->fatol = 1.e-25; 6713118ae5eSBarry Smith pseudo->frtol = 1.e-5; 6723118ae5eSBarry Smith #else 6733118ae5eSBarry Smith pseudo->fatol = 1.e-50; 6743118ae5eSBarry Smith pseudo->frtol = 1.e-12; 6753118ae5eSBarry Smith #endif 676bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr); 677bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr); 678bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetMaxTimeStep_C",TSPseudoSetMaxTimeStep_Pseudo);CHKERRQ(ierr); 679bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr); 680bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr); 6813a40ed3dSBarry Smith PetscFunctionReturn(0); 6822d3f70b5SBarry Smith } 6832d3f70b5SBarry Smith 68482bf6240SBarry Smith /*@C 6858d359177SBarry Smith TSPseudoTimeStepDefault - Default code to compute pseudo-timestepping. 68682bf6240SBarry Smith Use with TSPseudoSetTimeStep(). 68728aa8177SBarry Smith 68815091d37SBarry Smith Collective on TS 68915091d37SBarry Smith 69028aa8177SBarry Smith Input Parameters: 69128aa8177SBarry Smith . ts - the timestep context 69282bf6240SBarry Smith . dtctx - unused timestep context 69328aa8177SBarry Smith 69482bf6240SBarry Smith Output Parameter: 69582bf6240SBarry Smith . newdt - the timestep to use for the next step 69628aa8177SBarry Smith 69715091d37SBarry Smith Level: advanced 69815091d37SBarry Smith 69982bf6240SBarry Smith .keywords: timestep, pseudo, default 700564e8f4eSLois Curfman McInnes 70182bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep() 70228aa8177SBarry Smith @*/ 7038d359177SBarry Smith PetscErrorCode TSPseudoTimeStepDefault(TS ts,PetscReal *newdt,void *dtctx) 70428aa8177SBarry Smith { 70582bf6240SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 706be5899b3SLisandro Dalcin PetscReal inc = pseudo->dt_increment; 707dfbe8321SBarry Smith PetscErrorCode ierr; 70828aa8177SBarry Smith 7093a40ed3dSBarry Smith PetscFunctionBegin; 710bbd7b040SJed Brown ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr); 711214bc6a2SJed Brown ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr); 71282bf6240SBarry Smith ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 713be5899b3SLisandro Dalcin if (pseudo->fnorm_initial < 0) { 71482bf6240SBarry Smith /* first time through so compute initial function norm */ 715cdbf8f93SLisandro Dalcin pseudo->fnorm_initial = pseudo->fnorm; 716be5899b3SLisandro Dalcin pseudo->fnorm_previous = pseudo->fnorm; 71782bf6240SBarry Smith } 718bbd56ea5SKarl Rupp if (pseudo->fnorm == 0.0) *newdt = 1.e12*inc*ts->time_step; 719bbd56ea5SKarl Rupp else if (pseudo->increment_dt_from_initial_dt) *newdt = inc*pseudo->dt_initial*pseudo->fnorm_initial/pseudo->fnorm; 720be5899b3SLisandro Dalcin else *newdt = inc*ts->time_step*pseudo->fnorm_previous/pseudo->fnorm; 72186534af1SJed Brown if (pseudo->dt_max > 0) *newdt = PetscMin(*newdt,pseudo->dt_max); 72282bf6240SBarry Smith pseudo->fnorm_previous = pseudo->fnorm; 7233a40ed3dSBarry Smith PetscFunctionReturn(0); 72428aa8177SBarry Smith } 725