12d3f70b5SBarry Smith /* 2fb4a63b6SLois Curfman McInnes Code for Timestepping with implicit backwards Euler. 32d3f70b5SBarry Smith */ 4b45d2f2cSJed Brown #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; 257bf11e45SBarry Smith } TS_Pseudo; 262d3f70b5SBarry Smith 272d3f70b5SBarry Smith /* ------------------------------------------------------------------------------*/ 282d3f70b5SBarry Smith 294a2ae208SSatish Balay #undef __FUNCT__ 304a2ae208SSatish Balay #define __FUNCT__ "TSPseudoComputeTimeStep" 317bf11e45SBarry Smith /*@ 327bf11e45SBarry Smith TSPseudoComputeTimeStep - Computes the next timestep for a currently running 33564e8f4eSLois Curfman McInnes pseudo-timestepping process. 342d3f70b5SBarry Smith 3515091d37SBarry Smith Collective on TS 3615091d37SBarry Smith 377bf11e45SBarry Smith Input Parameter: 387bf11e45SBarry Smith . ts - timestep context 397bf11e45SBarry Smith 407bf11e45SBarry Smith Output Parameter: 41fb4a63b6SLois Curfman McInnes . dt - newly computed timestep 42fb4a63b6SLois Curfman McInnes 4315091d37SBarry Smith Level: advanced 44564e8f4eSLois Curfman McInnes 45564e8f4eSLois Curfman McInnes Notes: 46564e8f4eSLois Curfman McInnes The routine to be called here to compute the timestep should be 47564e8f4eSLois Curfman McInnes set by calling TSPseudoSetTimeStep(). 48564e8f4eSLois Curfman McInnes 49fb4a63b6SLois Curfman McInnes .keywords: timestep, pseudo, compute 50564e8f4eSLois Curfman McInnes 51564e8f4eSLois Curfman McInnes .seealso: TSPseudoDefaultTimeStep(), TSPseudoSetTimeStep() 527bf11e45SBarry Smith @*/ 537087cfbeSBarry Smith PetscErrorCode TSPseudoComputeTimeStep(TS ts,PetscReal *dt) 547bf11e45SBarry Smith { 557bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 56dfbe8321SBarry Smith PetscErrorCode ierr; 577bf11e45SBarry Smith 583a40ed3dSBarry Smith PetscFunctionBegin; 59d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr); 607bf11e45SBarry Smith ierr = (*pseudo->dt)(ts,dt,pseudo->dtctx);CHKERRQ(ierr); 61d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr); 623a40ed3dSBarry Smith PetscFunctionReturn(0); 637bf11e45SBarry Smith } 647bf11e45SBarry Smith 657bf11e45SBarry Smith 667bf11e45SBarry Smith /* ------------------------------------------------------------------------------*/ 674a2ae208SSatish Balay #undef __FUNCT__ 684a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultVerifyTimeStep" 697bf11e45SBarry Smith /*@C 70639f9d9dSBarry Smith TSPseudoDefaultVerifyTimeStep - Default code to verify the quality of the last timestep. 717bf11e45SBarry Smith 7215091d37SBarry Smith Collective on TS 7315091d37SBarry Smith 747bf11e45SBarry Smith Input Parameters: 7515091d37SBarry Smith + ts - the timestep context 767bf11e45SBarry Smith . dtctx - unused timestep context 7715091d37SBarry Smith - update - latest solution vector 787bf11e45SBarry Smith 79564e8f4eSLois Curfman McInnes Output Parameters: 8015091d37SBarry Smith + newdt - the timestep to use for the next step 8115091d37SBarry Smith - flag - flag indicating whether the last time step was acceptable 827bf11e45SBarry Smith 8315091d37SBarry Smith Level: advanced 84fee21e36SBarry Smith 85564e8f4eSLois Curfman McInnes Note: 86564e8f4eSLois Curfman McInnes This routine always returns a flag of 1, indicating an acceptable 87564e8f4eSLois Curfman McInnes timestep. 88564e8f4eSLois Curfman McInnes 89564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, default, verify 90564e8f4eSLois Curfman McInnes 91564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStep() 927bf11e45SBarry Smith @*/ 937087cfbeSBarry Smith PetscErrorCode TSPseudoDefaultVerifyTimeStep(TS ts,Vec update,void *dtctx,PetscReal *newdt,PetscBool *flag) 947bf11e45SBarry Smith { 953a40ed3dSBarry Smith PetscFunctionBegin; 96a7cc72afSBarry Smith *flag = PETSC_TRUE; 973a40ed3dSBarry Smith PetscFunctionReturn(0); 987bf11e45SBarry Smith } 997bf11e45SBarry Smith 1007bf11e45SBarry Smith 1014a2ae208SSatish Balay #undef __FUNCT__ 1024a2ae208SSatish Balay #define __FUNCT__ "TSPseudoVerifyTimeStep" 1037bf11e45SBarry Smith /*@ 104564e8f4eSLois Curfman McInnes TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable. 1057bf11e45SBarry Smith 10615091d37SBarry Smith Collective on TS 10715091d37SBarry Smith 108fb4a63b6SLois Curfman McInnes Input Parameters: 10915091d37SBarry Smith + ts - timestep context 11015091d37SBarry Smith - update - latest solution vector 1117bf11e45SBarry Smith 112fb4a63b6SLois Curfman McInnes Output Parameters: 11315091d37SBarry Smith + dt - newly computed timestep (if it had to shrink) 11415091d37SBarry Smith - flag - indicates if current timestep was ok 1157bf11e45SBarry Smith 11615091d37SBarry Smith Level: advanced 117fee21e36SBarry Smith 118564e8f4eSLois Curfman McInnes Notes: 119564e8f4eSLois Curfman McInnes The routine to be called here to compute the timestep should be 120564e8f4eSLois Curfman McInnes set by calling TSPseudoSetVerifyTimeStep(). 121564e8f4eSLois Curfman McInnes 122564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, verify 123564e8f4eSLois Curfman McInnes 124564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoDefaultVerifyTimeStep() 1257bf11e45SBarry Smith @*/ 1267087cfbeSBarry Smith PetscErrorCode TSPseudoVerifyTimeStep(TS ts,Vec update,PetscReal *dt,PetscBool *flag) 1277bf11e45SBarry Smith { 1287bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 129dfbe8321SBarry Smith PetscErrorCode ierr; 1307bf11e45SBarry Smith 1313a40ed3dSBarry Smith PetscFunctionBegin; 132a7cc72afSBarry Smith if (!pseudo->verify) {*flag = PETSC_TRUE; PetscFunctionReturn(0);} 1337bf11e45SBarry Smith 1347bf11e45SBarry Smith ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag);CHKERRQ(ierr); 1353a40ed3dSBarry Smith PetscFunctionReturn(0); 1367bf11e45SBarry Smith } 1377bf11e45SBarry Smith 1387bf11e45SBarry Smith /* --------------------------------------------------------------------------------*/ 1397bf11e45SBarry Smith 1404a2ae208SSatish Balay #undef __FUNCT__ 1414a2ae208SSatish Balay #define __FUNCT__ "TSStep_Pseudo" 142193ac0bcSJed Brown static PetscErrorCode TSStep_Pseudo(TS ts) 1432d3f70b5SBarry Smith { 144277b19d0SLisandro Dalcin TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 145cdbf8f93SLisandro Dalcin PetscInt its,lits,reject; 146cdbf8f93SLisandro Dalcin PetscBool stepok; 147cdbf8f93SLisandro Dalcin PetscReal next_time_step; 148cdbf8f93SLisandro Dalcin SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING; 149dfbe8321SBarry Smith PetscErrorCode ierr; 1502d3f70b5SBarry Smith 1513a40ed3dSBarry Smith PetscFunctionBegin; 152bbd56ea5SKarl Rupp if (ts->steps == 0) pseudo->dt_initial = ts->time_step; 153193ac0bcSJed Brown ierr = VecCopy(ts->vec_sol,pseudo->update);CHKERRQ(ierr); 154cdbf8f93SLisandro Dalcin next_time_step = ts->time_step; 155cdbf8f93SLisandro Dalcin ierr = TSPseudoComputeTimeStep(ts,&next_time_step);CHKERRQ(ierr); 156cdbf8f93SLisandro Dalcin for (reject=0; reject<ts->max_reject; reject++,ts->reject++) { 157cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 158b8123daeSJed Brown ierr = TSPreStep(ts);CHKERRQ(ierr); 159b8123daeSJed Brown ierr = TSPreStage(ts,ts->ptime+ts->time_step);CHKERRQ(ierr); 1600298fd71SBarry Smith ierr = SNESSolve(ts->snes,NULL,pseudo->update);CHKERRQ(ierr); 161cdbf8f93SLisandro Dalcin ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr); 162b850b91aSLisandro Dalcin ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 1639f954729SBarry Smith ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr); 1645ef26d82SJed Brown ts->snes_its += its; ts->ksp_its += lits; 165cdbf8f93SLisandro Dalcin ierr = PetscInfo3(ts,"step=%D, nonlinear solve iterations=%D, linear solve iterations=%D\n",ts->steps,its,lits);CHKERRQ(ierr); 166193ac0bcSJed Brown pseudo->fnorm = -1; /* The current norm is no longer valid, monitor must recompute it. */ 167cdbf8f93SLisandro Dalcin ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&next_time_step,&stepok);CHKERRQ(ierr); 168cdbf8f93SLisandro Dalcin if (stepok) break; 169cdbf8f93SLisandro Dalcin } 170f1b97656SJed Brown if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) { 171193ac0bcSJed Brown ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 172193ac0bcSJed Brown ierr = PetscInfo2(ts,"step=%D, nonlinear solve solve failures %D greater than current TS allowed, stopping solve\n",ts->steps,ts->num_snes_failures);CHKERRQ(ierr); 173cdbf8f93SLisandro Dalcin PetscFunctionReturn(0); 1747bf11e45SBarry Smith } 175cdbf8f93SLisandro Dalcin if (reject >= ts->max_reject) { 176193ac0bcSJed Brown ts->reason = TS_DIVERGED_STEP_REJECTED; 177cdbf8f93SLisandro Dalcin ierr = PetscInfo2(ts,"step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,reject);CHKERRQ(ierr); 178cdbf8f93SLisandro Dalcin PetscFunctionReturn(0); 179193ac0bcSJed Brown } 180cdbf8f93SLisandro Dalcin ierr = VecCopy(pseudo->update,ts->vec_sol);CHKERRQ(ierr); 181193ac0bcSJed Brown ts->ptime += ts->time_step; 182cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 1832d3f70b5SBarry Smith ts->steps++; 1843a40ed3dSBarry Smith PetscFunctionReturn(0); 1852d3f70b5SBarry Smith } 1862d3f70b5SBarry Smith 1872d3f70b5SBarry Smith /*------------------------------------------------------------*/ 1884a2ae208SSatish Balay #undef __FUNCT__ 189277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Pseudo" 190277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Pseudo(TS ts) 1912d3f70b5SBarry Smith { 1927bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 193dfbe8321SBarry Smith PetscErrorCode ierr; 1942d3f70b5SBarry Smith 1953a40ed3dSBarry Smith PetscFunctionBegin; 1966bf464f9SBarry Smith ierr = VecDestroy(&pseudo->update);CHKERRQ(ierr); 1976bf464f9SBarry Smith ierr = VecDestroy(&pseudo->func);CHKERRQ(ierr); 1986bf464f9SBarry Smith ierr = VecDestroy(&pseudo->xdot);CHKERRQ(ierr); 1993a40ed3dSBarry Smith PetscFunctionReturn(0); 2002d3f70b5SBarry Smith } 2012d3f70b5SBarry Smith 202277b19d0SLisandro Dalcin #undef __FUNCT__ 203277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Pseudo" 204277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Pseudo(TS ts) 205277b19d0SLisandro Dalcin { 206277b19d0SLisandro Dalcin PetscErrorCode ierr; 207277b19d0SLisandro Dalcin 208277b19d0SLisandro Dalcin PetscFunctionBegin; 209277b19d0SLisandro Dalcin ierr = TSReset_Pseudo(ts);CHKERRQ(ierr); 210277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 21100de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C","",NULL);CHKERRQ(ierr); 21200de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C","",NULL);CHKERRQ(ierr); 21300de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetMaxTimeStep_C","",NULL);CHKERRQ(ierr); 21400de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C","",NULL);CHKERRQ(ierr); 21500de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStep_C","",NULL);CHKERRQ(ierr); 216277b19d0SLisandro Dalcin PetscFunctionReturn(0); 217277b19d0SLisandro Dalcin } 2182d3f70b5SBarry Smith 2192d3f70b5SBarry Smith /*------------------------------------------------------------*/ 2202d3f70b5SBarry Smith 2214a2ae208SSatish Balay #undef __FUNCT__ 2226f2d6a7bSJed Brown #define __FUNCT__ "TSPseudoGetXdot" 2236f2d6a7bSJed Brown /* 2246f2d6a7bSJed Brown Compute Xdot = (X^{n+1}-X^n)/dt) = 0 2256f2d6a7bSJed Brown */ 2266f2d6a7bSJed Brown static PetscErrorCode TSPseudoGetXdot(TS ts,Vec X,Vec *Xdot) 2272d3f70b5SBarry Smith { 2286f2d6a7bSJed Brown TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 229193ac0bcSJed Brown const PetscScalar mdt = 1.0/ts->time_step,*xnp1,*xn; 230193ac0bcSJed Brown PetscScalar *xdot; 231dfbe8321SBarry Smith PetscErrorCode ierr; 232a7cc72afSBarry Smith PetscInt i,n; 2332d3f70b5SBarry Smith 2343a40ed3dSBarry Smith PetscFunctionBegin; 235193ac0bcSJed Brown ierr = VecGetArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr); 236193ac0bcSJed Brown ierr = VecGetArrayRead(X,&xnp1);CHKERRQ(ierr); 2376f2d6a7bSJed Brown ierr = VecGetArray(pseudo->xdot,&xdot);CHKERRQ(ierr); 2386f2d6a7bSJed Brown ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr); 239bbd56ea5SKarl Rupp for (i=0; i<n; i++) xdot[i] = mdt*(xnp1[i] - xn[i]); 240193ac0bcSJed Brown ierr = VecRestoreArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr); 241193ac0bcSJed Brown ierr = VecRestoreArrayRead(X,&xnp1);CHKERRQ(ierr); 2426f2d6a7bSJed Brown ierr = VecRestoreArray(pseudo->xdot,&xdot);CHKERRQ(ierr); 2436f2d6a7bSJed Brown *Xdot = pseudo->xdot; 2443a40ed3dSBarry Smith PetscFunctionReturn(0); 2452d3f70b5SBarry Smith } 2462d3f70b5SBarry Smith 2474a2ae208SSatish Balay #undef __FUNCT__ 2480f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Pseudo" 2496f2d6a7bSJed Brown /* 2506f2d6a7bSJed Brown The transient residual is 2516f2d6a7bSJed Brown 2526f2d6a7bSJed Brown F(U^{n+1},(U^{n+1}-U^n)/dt) = 0 2536f2d6a7bSJed Brown 2546f2d6a7bSJed Brown or for ODE, 2556f2d6a7bSJed Brown 2566f2d6a7bSJed Brown (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0 2576f2d6a7bSJed Brown 2586f2d6a7bSJed Brown This is the function that must be evaluated for transient simulation and for 2596f2d6a7bSJed Brown finite difference Jacobians. On the first Newton step, this algorithm uses 2606f2d6a7bSJed Brown a guess of U^{n+1} = U^n in which case the transient term vanishes and the 2616f2d6a7bSJed Brown residual is actually the steady state residual. Pseudotransient 2626f2d6a7bSJed Brown continuation as described in the literature is a linearly implicit 2636f2d6a7bSJed Brown algorithm, it only takes this one Newton step with the steady state 2646f2d6a7bSJed Brown residual, and then advances to the next time step. 2656f2d6a7bSJed Brown */ 2660f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes,Vec X,Vec Y,TS ts) 2672d3f70b5SBarry Smith { 2686f2d6a7bSJed Brown Vec Xdot; 269dfbe8321SBarry Smith PetscErrorCode ierr; 2702d3f70b5SBarry Smith 2713a40ed3dSBarry Smith PetscFunctionBegin; 2726f2d6a7bSJed Brown ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr); 273193ac0bcSJed Brown ierr = TSComputeIFunction(ts,ts->ptime+ts->time_step,X,Xdot,Y,PETSC_FALSE);CHKERRQ(ierr); 2746f2d6a7bSJed Brown PetscFunctionReturn(0); 2756f2d6a7bSJed Brown } 2762d3f70b5SBarry Smith 2776f2d6a7bSJed Brown #undef __FUNCT__ 2780f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Pseudo" 2796f2d6a7bSJed Brown /* 2806f2d6a7bSJed Brown This constructs the Jacobian needed for SNES. For DAE, this is 2816f2d6a7bSJed Brown 2826f2d6a7bSJed Brown dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot 2836f2d6a7bSJed Brown 2846f2d6a7bSJed Brown and for ODE: 2856f2d6a7bSJed Brown 2866f2d6a7bSJed Brown J = I/dt - J_{Frhs} where J_{Frhs} is the given Jacobian of Frhs. 2876f2d6a7bSJed Brown */ 2880f5c6efeSJed Brown static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes,Vec X,Mat *AA,Mat *BB,MatStructure *str,TS ts) 2896f2d6a7bSJed Brown { 2906f2d6a7bSJed Brown Vec Xdot; 2916f2d6a7bSJed Brown PetscErrorCode ierr; 2926f2d6a7bSJed Brown 2936f2d6a7bSJed Brown PetscFunctionBegin; 2946f2d6a7bSJed Brown ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr); 295193ac0bcSJed Brown ierr = TSComputeIJacobian(ts,ts->ptime+ts->time_step,X,Xdot,1./ts->time_step,AA,BB,str,PETSC_FALSE);CHKERRQ(ierr); 2963a40ed3dSBarry Smith PetscFunctionReturn(0); 2972d3f70b5SBarry Smith } 2982d3f70b5SBarry Smith 2992d3f70b5SBarry Smith 3004a2ae208SSatish Balay #undef __FUNCT__ 3014a2ae208SSatish Balay #define __FUNCT__ "TSSetUp_Pseudo" 3026849ba73SBarry Smith static PetscErrorCode TSSetUp_Pseudo(TS ts) 3032d3f70b5SBarry Smith { 3047bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 305dfbe8321SBarry Smith PetscErrorCode ierr; 3062d3f70b5SBarry Smith 3073a40ed3dSBarry Smith PetscFunctionBegin; 3087bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr); 3097bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr); 3106f2d6a7bSJed Brown ierr = VecDuplicate(ts->vec_sol,&pseudo->xdot);CHKERRQ(ierr); 3113a40ed3dSBarry Smith PetscFunctionReturn(0); 3122d3f70b5SBarry Smith } 3132d3f70b5SBarry Smith /*------------------------------------------------------------*/ 3142d3f70b5SBarry Smith 3154a2ae208SSatish Balay #undef __FUNCT__ 316a6570f20SBarry Smith #define __FUNCT__ "TSPseudoMonitorDefault" 317649052a6SBarry Smith PetscErrorCode TSPseudoMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy) 3182d3f70b5SBarry Smith { 3197bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 320dfbe8321SBarry Smith PetscErrorCode ierr; 321ce94432eSBarry Smith PetscViewer viewer = (PetscViewer) dummy; 3222d3f70b5SBarry Smith 3233a40ed3dSBarry Smith PetscFunctionBegin; 324ce94432eSBarry Smith if (!viewer) { 325ce94432eSBarry Smith ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ts),&viewer);CHKERRQ(ierr); 326ce94432eSBarry Smith } 327193ac0bcSJed Brown if (pseudo->fnorm < 0) { /* The last computed norm is stale, recompute */ 328193ac0bcSJed Brown ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr); 329193ac0bcSJed Brown ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr); 330193ac0bcSJed Brown ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 331193ac0bcSJed Brown } 332649052a6SBarry Smith ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 333649052a6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"TS %D dt %G time %G fnorm %G\n",step,ts->time_step,ptime,pseudo->fnorm);CHKERRQ(ierr); 334649052a6SBarry Smith ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 3353a40ed3dSBarry Smith PetscFunctionReturn(0); 3362d3f70b5SBarry Smith } 3372d3f70b5SBarry Smith 3384a2ae208SSatish Balay #undef __FUNCT__ 3394a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Pseudo" 3406849ba73SBarry Smith static PetscErrorCode TSSetFromOptions_Pseudo(TS ts) 3412d3f70b5SBarry Smith { 3424bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 343dfbe8321SBarry Smith PetscErrorCode ierr; 344ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE; 345649052a6SBarry Smith PetscViewer viewer; 3462d3f70b5SBarry Smith 3473a40ed3dSBarry Smith PetscFunctionBegin; 348b0a32e0cSBarry Smith ierr = PetscOptionsHead("Pseudo-timestepping options");CHKERRQ(ierr); 3490298fd71SBarry Smith ierr = PetscOptionsBool("-ts_monitor_pseudo","Monitor convergence","TSPseudoMonitorDefault",flg,&flg,NULL);CHKERRQ(ierr); 3502d3f70b5SBarry Smith if (flg) { 351ce94432eSBarry Smith ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ts),"stdout",&viewer);CHKERRQ(ierr); 352649052a6SBarry Smith ierr = TSMonitorSet(ts,TSPseudoMonitorDefault,viewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 35328aa8177SBarry Smith } 35490d69ab7SBarry Smith flg = PETSC_FALSE; 3550298fd71SBarry Smith ierr = PetscOptionsBool("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",flg,&flg,NULL);CHKERRQ(ierr); 356ca4b7087SBarry Smith if (flg) { 357ca4b7087SBarry Smith ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr); 358ca4b7087SBarry Smith } 359419fbf4bSSatish Balay ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,0);CHKERRQ(ierr); 36086534af1SJed Brown ierr = PetscOptionsReal("-ts_pseudo_max_dt","Maximum value for dt","TSPseudoSetMaxTimeStep",pseudo->dt_max,&pseudo->dt_max,0);CHKERRQ(ierr); 361d52bd9f3SBarry Smith 362d52bd9f3SBarry Smith ierr = SNESSetFromOptions(ts->snes);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 { 371d52bd9f3SBarry Smith PetscErrorCode ierr; 372d52bd9f3SBarry Smith 3733a40ed3dSBarry Smith PetscFunctionBegin; 374d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 3753a40ed3dSBarry Smith PetscFunctionReturn(0); 3762d3f70b5SBarry Smith } 3772d3f70b5SBarry Smith 37882bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 3794a2ae208SSatish Balay #undef __FUNCT__ 3804a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep" 381ac226902SBarry Smith /*@C 38282bf6240SBarry Smith TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the 38382bf6240SBarry Smith last timestep. 38482bf6240SBarry Smith 3853f9fe445SBarry Smith Logically Collective on TS 38615091d37SBarry Smith 38782bf6240SBarry Smith Input Parameters: 38815091d37SBarry Smith + ts - timestep context 38982bf6240SBarry Smith . dt - user-defined function to verify timestep 39015091d37SBarry Smith - ctx - [optional] user-defined context for private data 3910298fd71SBarry Smith for the timestep verification routine (may be NULL) 39282bf6240SBarry Smith 39315091d37SBarry Smith Level: advanced 394fee21e36SBarry Smith 39582bf6240SBarry Smith Calling sequence of func: 396ace3abfcSBarry Smith . func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscBool *flag); 39782bf6240SBarry Smith 39882bf6240SBarry Smith . update - latest solution vector 39982bf6240SBarry Smith . ctx - [optional] timestep context 40082bf6240SBarry Smith . newdt - the timestep to use for the next step 40182bf6240SBarry Smith . flag - flag indicating whether the last time step was acceptable 40282bf6240SBarry Smith 40382bf6240SBarry Smith Notes: 40482bf6240SBarry Smith The routine set here will be called by TSPseudoVerifyTimeStep() 40582bf6240SBarry Smith during the timestepping process. 40682bf6240SBarry Smith 40782bf6240SBarry Smith .keywords: timestep, pseudo, set, verify 40882bf6240SBarry Smith 40982bf6240SBarry Smith .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep() 41082bf6240SBarry Smith @*/ 4117087cfbeSBarry Smith PetscErrorCode TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,PetscBool*),void *ctx) 41282bf6240SBarry Smith { 4134ac538c5SBarry Smith PetscErrorCode ierr; 41482bf6240SBarry Smith 41582bf6240SBarry Smith PetscFunctionBegin; 4160700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 4174ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSPseudoSetVerifyTimeStep_C",(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal*,PetscBool*),void*),(ts,dt,ctx));CHKERRQ(ierr); 41882bf6240SBarry Smith PetscFunctionReturn(0); 41982bf6240SBarry Smith } 42082bf6240SBarry Smith 4214a2ae208SSatish Balay #undef __FUNCT__ 4224a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement" 42382bf6240SBarry Smith /*@ 42482bf6240SBarry Smith TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 42582bf6240SBarry Smith dt when using the TSPseudoDefaultTimeStep() routine. 42682bf6240SBarry Smith 4273f9fe445SBarry Smith Logically Collective on TS 428fee21e36SBarry Smith 42915091d37SBarry Smith Input Parameters: 43015091d37SBarry Smith + ts - the timestep context 43115091d37SBarry Smith - inc - the scaling factor >= 1.0 43215091d37SBarry Smith 43382bf6240SBarry Smith Options Database Key: 43482bf6240SBarry Smith $ -ts_pseudo_increment <increment> 43582bf6240SBarry Smith 43615091d37SBarry Smith Level: advanced 43715091d37SBarry Smith 43882bf6240SBarry Smith .keywords: timestep, pseudo, set, increment 43982bf6240SBarry Smith 44082bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 44182bf6240SBarry Smith @*/ 4427087cfbeSBarry Smith PetscErrorCode TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc) 44382bf6240SBarry Smith { 4444ac538c5SBarry Smith PetscErrorCode ierr; 44582bf6240SBarry Smith 44682bf6240SBarry Smith PetscFunctionBegin; 4470700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 448c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts,inc,2); 4494ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSPseudoSetTimeStepIncrement_C",(TS,PetscReal),(ts,inc));CHKERRQ(ierr); 45082bf6240SBarry Smith PetscFunctionReturn(0); 45182bf6240SBarry Smith } 45282bf6240SBarry Smith 4534a2ae208SSatish Balay #undef __FUNCT__ 45486534af1SJed Brown #define __FUNCT__ "TSPseudoSetMaxTimeStep" 45586534af1SJed Brown /*@ 45686534af1SJed Brown TSPseudoSetMaxTimeStep - Sets the maximum time step 45786534af1SJed Brown when using the TSPseudoDefaultTimeStep() routine. 45886534af1SJed Brown 45986534af1SJed Brown Logically Collective on TS 46086534af1SJed Brown 46186534af1SJed Brown Input Parameters: 46286534af1SJed Brown + ts - the timestep context 46386534af1SJed Brown - maxdt - the maximum time step, use a non-positive value to deactivate 46486534af1SJed Brown 46586534af1SJed Brown Options Database Key: 46686534af1SJed Brown $ -ts_pseudo_max_dt <increment> 46786534af1SJed Brown 46886534af1SJed Brown Level: advanced 46986534af1SJed Brown 47086534af1SJed Brown .keywords: timestep, pseudo, set 47186534af1SJed Brown 47286534af1SJed Brown .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 47386534af1SJed Brown @*/ 47486534af1SJed Brown PetscErrorCode TSPseudoSetMaxTimeStep(TS ts,PetscReal maxdt) 47586534af1SJed Brown { 47686534af1SJed Brown PetscErrorCode ierr; 47786534af1SJed Brown 47886534af1SJed Brown PetscFunctionBegin; 47986534af1SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 48086534af1SJed Brown PetscValidLogicalCollectiveReal(ts,maxdt,2); 48186534af1SJed Brown ierr = PetscTryMethod(ts,"TSPseudoSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));CHKERRQ(ierr); 48286534af1SJed Brown PetscFunctionReturn(0); 48386534af1SJed Brown } 48486534af1SJed Brown 48586534af1SJed Brown #undef __FUNCT__ 4864a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt" 48782bf6240SBarry Smith /*@ 48882bf6240SBarry Smith TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep 48982bf6240SBarry Smith is computed via the formula 49082bf6240SBarry Smith $ dt = initial_dt*initial_fnorm/current_fnorm 49182bf6240SBarry Smith rather than the default update, 49282bf6240SBarry Smith $ dt = current_dt*previous_fnorm/current_fnorm. 49382bf6240SBarry Smith 4943f9fe445SBarry Smith Logically Collective on TS 49515091d37SBarry Smith 49682bf6240SBarry Smith Input Parameter: 49782bf6240SBarry Smith . ts - the timestep context 49882bf6240SBarry Smith 49982bf6240SBarry Smith Options Database Key: 50082bf6240SBarry Smith $ -ts_pseudo_increment_dt_from_initial_dt 50182bf6240SBarry Smith 50215091d37SBarry Smith Level: advanced 50315091d37SBarry Smith 50482bf6240SBarry Smith .keywords: timestep, pseudo, set, increment 50582bf6240SBarry Smith 50682bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 50782bf6240SBarry Smith @*/ 5087087cfbeSBarry Smith PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS ts) 50982bf6240SBarry Smith { 5104ac538c5SBarry Smith PetscErrorCode ierr; 51182bf6240SBarry Smith 51282bf6240SBarry Smith PetscFunctionBegin; 5130700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 5144ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSPseudoIncrementDtFromInitialDt_C",(TS),(ts));CHKERRQ(ierr); 51582bf6240SBarry Smith PetscFunctionReturn(0); 51682bf6240SBarry Smith } 51782bf6240SBarry Smith 51882bf6240SBarry Smith 5194a2ae208SSatish Balay #undef __FUNCT__ 5204a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep" 521ac226902SBarry Smith /*@C 52282bf6240SBarry Smith TSPseudoSetTimeStep - Sets the user-defined routine to be 52382bf6240SBarry Smith called at each pseudo-timestep to update the timestep. 52482bf6240SBarry Smith 5253f9fe445SBarry Smith Logically Collective on TS 52615091d37SBarry Smith 52782bf6240SBarry Smith Input Parameters: 52815091d37SBarry Smith + ts - timestep context 52982bf6240SBarry Smith . dt - function to compute timestep 53015091d37SBarry Smith - ctx - [optional] user-defined context for private data 5310298fd71SBarry Smith required by the function (may be NULL) 53282bf6240SBarry Smith 53315091d37SBarry Smith Level: intermediate 534fee21e36SBarry Smith 53582bf6240SBarry Smith Calling sequence of func: 53687828ca2SBarry Smith . func (TS ts,PetscReal *newdt,void *ctx); 53782bf6240SBarry Smith 53882bf6240SBarry Smith . newdt - the newly computed timestep 53982bf6240SBarry Smith . ctx - [optional] timestep context 54082bf6240SBarry Smith 54182bf6240SBarry Smith Notes: 54282bf6240SBarry Smith The routine set here will be called by TSPseudoComputeTimeStep() 54382bf6240SBarry Smith during the timestepping process. 54482bf6240SBarry Smith 54582bf6240SBarry Smith .keywords: timestep, pseudo, set 54682bf6240SBarry Smith 54782bf6240SBarry Smith .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep() 54882bf6240SBarry Smith @*/ 5497087cfbeSBarry Smith PetscErrorCode TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void *ctx) 55082bf6240SBarry Smith { 5514ac538c5SBarry Smith PetscErrorCode ierr; 55282bf6240SBarry Smith 55382bf6240SBarry Smith PetscFunctionBegin; 5540700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 5554ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSPseudoSetTimeStep_C",(TS,PetscErrorCode (*)(TS,PetscReal*,void*),void*),(ts,dt,ctx));CHKERRQ(ierr); 55682bf6240SBarry Smith PetscFunctionReturn(0); 55782bf6240SBarry Smith } 55882bf6240SBarry Smith 55982bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 56082bf6240SBarry Smith 561ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscBool*); /* force argument to next function to not be extern C*/ 5624a2ae208SSatish Balay #undef __FUNCT__ 5634a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo" 5647087cfbeSBarry Smith PetscErrorCode TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void *ctx) 56582bf6240SBarry Smith { 56682bf6240SBarry Smith TS_Pseudo *pseudo; 56782bf6240SBarry Smith 56882bf6240SBarry Smith PetscFunctionBegin; 56982bf6240SBarry Smith pseudo = (TS_Pseudo*)ts->data; 57082bf6240SBarry Smith pseudo->verify = dt; 57182bf6240SBarry Smith pseudo->verifyctx = ctx; 57282bf6240SBarry Smith PetscFunctionReturn(0); 57382bf6240SBarry Smith } 57482bf6240SBarry Smith 5754a2ae208SSatish Balay #undef __FUNCT__ 5764a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement_Pseudo" 5777087cfbeSBarry Smith PetscErrorCode TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc) 57882bf6240SBarry Smith { 5794bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 58082bf6240SBarry Smith 58182bf6240SBarry Smith PetscFunctionBegin; 58282bf6240SBarry Smith pseudo->dt_increment = inc; 58382bf6240SBarry Smith PetscFunctionReturn(0); 58482bf6240SBarry Smith } 58582bf6240SBarry Smith 5864a2ae208SSatish Balay #undef __FUNCT__ 58786534af1SJed Brown #define __FUNCT__ "TSPseudoSetMaxTimeStep_Pseudo" 58886534af1SJed Brown PetscErrorCode TSPseudoSetMaxTimeStep_Pseudo(TS ts,PetscReal maxdt) 58986534af1SJed Brown { 59086534af1SJed Brown TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 59186534af1SJed Brown 59286534af1SJed Brown PetscFunctionBegin; 59386534af1SJed Brown pseudo->dt_max = maxdt; 59486534af1SJed Brown PetscFunctionReturn(0); 59586534af1SJed Brown } 59686534af1SJed Brown 59786534af1SJed Brown #undef __FUNCT__ 5984a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt_Pseudo" 5997087cfbeSBarry Smith PetscErrorCode TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts) 60082bf6240SBarry Smith { 6014bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 60282bf6240SBarry Smith 60382bf6240SBarry Smith PetscFunctionBegin; 6044bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_TRUE; 60582bf6240SBarry Smith PetscFunctionReturn(0); 60682bf6240SBarry Smith } 60782bf6240SBarry Smith 6086849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/ 6094a2ae208SSatish Balay #undef __FUNCT__ 6104a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo" 6117087cfbeSBarry Smith PetscErrorCode TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void *ctx) 61282bf6240SBarry Smith { 6134bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 61482bf6240SBarry Smith 61582bf6240SBarry Smith PetscFunctionBegin; 61682bf6240SBarry Smith pseudo->dt = dt; 61782bf6240SBarry Smith pseudo->dtctx = ctx; 61882bf6240SBarry Smith PetscFunctionReturn(0); 61982bf6240SBarry Smith } 62082bf6240SBarry Smith 62182bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 62210e6a065SJed Brown /*MC 62310e6a065SJed Brown TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping 62482bf6240SBarry Smith 62510e6a065SJed Brown This method solves equations of the form 62610e6a065SJed Brown 62710e6a065SJed Brown $ F(X,Xdot) = 0 62810e6a065SJed Brown 62910e6a065SJed Brown for steady state using the iteration 63010e6a065SJed Brown 63110e6a065SJed Brown $ [G'] S = -F(X,0) 63210e6a065SJed Brown $ X += S 63310e6a065SJed Brown 63410e6a065SJed Brown where 63510e6a065SJed Brown 63610e6a065SJed Brown $ G(Y) = F(Y,(Y-X)/dt) 63710e6a065SJed Brown 6386f2d6a7bSJed Brown This is linearly-implicit Euler with the residual always evaluated "at steady 6396f2d6a7bSJed Brown state". See note below. 64010e6a065SJed Brown 64110e6a065SJed Brown Options database keys: 64210e6a065SJed Brown + -ts_pseudo_increment <real> - ratio of increase dt 64310e6a065SJed Brown - -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt 64410e6a065SJed Brown 64510e6a065SJed Brown Level: beginner 64610e6a065SJed Brown 64710e6a065SJed Brown References: 64810e6a065SJed Brown Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential-Algebraic Equations, 2003. 64910e6a065SJed Brown C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998. 65010e6a065SJed Brown 65110e6a065SJed Brown Notes: 6526f2d6a7bSJed Brown The residual computed by this method includes the transient term (Xdot is computed instead of 6536f2d6a7bSJed Brown always being zero), but since the prediction from the last step is always the solution from the 6546f2d6a7bSJed Brown last step, on the first Newton iteration we have 6556f2d6a7bSJed Brown 6566f2d6a7bSJed Brown $ Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0 6576f2d6a7bSJed Brown 6586f2d6a7bSJed Brown Therefore, the linear system solved by the first Newton iteration is equivalent to the one 6596f2d6a7bSJed Brown described above and in the papers. If the user chooses to perform multiple Newton iterations, the 6606f2d6a7bSJed Brown algorithm is no longer the one described in the referenced papers. 66110e6a065SJed Brown 66210e6a065SJed Brown .seealso: TSCreate(), TS, TSSetType() 66310e6a065SJed Brown 66410e6a065SJed Brown M*/ 6654a2ae208SSatish Balay #undef __FUNCT__ 6664a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Pseudo" 667*8cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Pseudo(TS ts) 6682d3f70b5SBarry Smith { 6697bf11e45SBarry Smith TS_Pseudo *pseudo; 670dfbe8321SBarry Smith PetscErrorCode ierr; 671193ac0bcSJed Brown SNES snes; 67219fd82e9SBarry Smith SNESType stype; 6732d3f70b5SBarry Smith 6743a40ed3dSBarry Smith PetscFunctionBegin; 675277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Pseudo; 676000e7ae3SMatthew Knepley ts->ops->destroy = TSDestroy_Pseudo; 677000e7ae3SMatthew Knepley ts->ops->view = TSView_Pseudo; 6782d3f70b5SBarry Smith 679000e7ae3SMatthew Knepley ts->ops->setup = TSSetUp_Pseudo; 680000e7ae3SMatthew Knepley ts->ops->step = TSStep_Pseudo; 681000e7ae3SMatthew Knepley ts->ops->setfromoptions = TSSetFromOptions_Pseudo; 6820f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Pseudo; 6830f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Pseudo; 6847bf11e45SBarry Smith 685193ac0bcSJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 686193ac0bcSJed Brown ierr = SNESGetType(snes,&stype);CHKERRQ(ierr); 687193ac0bcSJed Brown if (!stype) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);} 6882d3f70b5SBarry Smith 68938f2d2fdSLisandro Dalcin ierr = PetscNewLog(ts,TS_Pseudo,&pseudo);CHKERRQ(ierr); 6907bf11e45SBarry Smith ts->data = (void*)pseudo; 6912d3f70b5SBarry Smith 69228aa8177SBarry Smith pseudo->dt_increment = 1.1; 6934bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_FALSE; 69428aa8177SBarry Smith pseudo->dt = TSPseudoDefaultTimeStep; 695193ac0bcSJed Brown pseudo->fnorm = -1; 69682bf6240SBarry Smith 69700de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C","TSPseudoSetVerifyTimeStep_Pseudo",TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr); 69800de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C","TSPseudoSetTimeStepIncrement_Pseudo",TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr); 69900de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetMaxTimeStep_C","TSPseudoSetMaxTimeStep_Pseudo",TSPseudoSetMaxTimeStep_Pseudo);CHKERRQ(ierr); 70000de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C","TSPseudoIncrementDtFromInitialDt_Pseudo",TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr); 70100de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStep_C","TSPseudoSetTimeStep_Pseudo",TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr); 7023a40ed3dSBarry Smith PetscFunctionReturn(0); 7032d3f70b5SBarry Smith } 7042d3f70b5SBarry Smith 7054a2ae208SSatish Balay #undef __FUNCT__ 7064a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultTimeStep" 70782bf6240SBarry Smith /*@C 70882bf6240SBarry Smith TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping. 70982bf6240SBarry Smith Use with TSPseudoSetTimeStep(). 71028aa8177SBarry Smith 71115091d37SBarry Smith Collective on TS 71215091d37SBarry Smith 71328aa8177SBarry Smith Input Parameters: 71428aa8177SBarry Smith . ts - the timestep context 71582bf6240SBarry Smith . dtctx - unused timestep context 71628aa8177SBarry Smith 71782bf6240SBarry Smith Output Parameter: 71882bf6240SBarry Smith . newdt - the timestep to use for the next step 71928aa8177SBarry Smith 72015091d37SBarry Smith Level: advanced 72115091d37SBarry Smith 72282bf6240SBarry Smith .keywords: timestep, pseudo, default 723564e8f4eSLois Curfman McInnes 72482bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep() 72528aa8177SBarry Smith @*/ 7267087cfbeSBarry Smith PetscErrorCode TSPseudoDefaultTimeStep(TS ts,PetscReal *newdt,void *dtctx) 72728aa8177SBarry Smith { 72882bf6240SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 72987828ca2SBarry Smith PetscReal inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous; 730dfbe8321SBarry Smith PetscErrorCode ierr; 73128aa8177SBarry Smith 7323a40ed3dSBarry Smith PetscFunctionBegin; 733bbd7b040SJed Brown ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr); 734214bc6a2SJed Brown ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr); 73582bf6240SBarry Smith ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 736cdbf8f93SLisandro Dalcin if (pseudo->fnorm_initial == 0.0) { 73782bf6240SBarry Smith /* first time through so compute initial function norm */ 738cdbf8f93SLisandro Dalcin pseudo->fnorm_initial = pseudo->fnorm; 73982bf6240SBarry Smith fnorm_previous = pseudo->fnorm; 74082bf6240SBarry Smith } 741bbd56ea5SKarl Rupp if (pseudo->fnorm == 0.0) *newdt = 1.e12*inc*ts->time_step; 742bbd56ea5SKarl Rupp else if (pseudo->increment_dt_from_initial_dt) *newdt = inc*pseudo->dt_initial*pseudo->fnorm_initial/pseudo->fnorm; 743bbd56ea5SKarl Rupp else *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm; 74486534af1SJed Brown if (pseudo->dt_max > 0) *newdt = PetscMin(*newdt,pseudo->dt_max); 74582bf6240SBarry Smith pseudo->fnorm_previous = pseudo->fnorm; 7463a40ed3dSBarry Smith PetscFunctionReturn(0); 74728aa8177SBarry Smith } 748