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 304a2ae208SSatish Balay #undef __FUNCT__ 314a2ae208SSatish Balay #define __FUNCT__ "TSPseudoComputeTimeStep" 328d359177SBarry Smith /*@C 337bf11e45SBarry Smith TSPseudoComputeTimeStep - Computes the next timestep for a currently running 34564e8f4eSLois Curfman McInnes pseudo-timestepping process. 352d3f70b5SBarry Smith 3615091d37SBarry Smith Collective on TS 3715091d37SBarry Smith 387bf11e45SBarry Smith Input Parameter: 397bf11e45SBarry Smith . ts - timestep context 407bf11e45SBarry Smith 417bf11e45SBarry Smith Output Parameter: 42fb4a63b6SLois Curfman McInnes . dt - newly computed timestep 43fb4a63b6SLois Curfman McInnes 448d359177SBarry Smith Level: developer 45564e8f4eSLois Curfman McInnes 46564e8f4eSLois Curfman McInnes Notes: 47564e8f4eSLois Curfman McInnes The routine to be called here to compute the timestep should be 48564e8f4eSLois Curfman McInnes set by calling TSPseudoSetTimeStep(). 49564e8f4eSLois Curfman McInnes 50fb4a63b6SLois Curfman McInnes .keywords: timestep, pseudo, compute 51564e8f4eSLois Curfman McInnes 528d359177SBarry Smith .seealso: TSPseudoTimeStepDefault(), TSPseudoSetTimeStep() 537bf11e45SBarry Smith @*/ 547087cfbeSBarry Smith PetscErrorCode TSPseudoComputeTimeStep(TS ts,PetscReal *dt) 557bf11e45SBarry Smith { 567bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 57dfbe8321SBarry Smith PetscErrorCode ierr; 587bf11e45SBarry Smith 593a40ed3dSBarry Smith PetscFunctionBegin; 60d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr); 617bf11e45SBarry Smith ierr = (*pseudo->dt)(ts,dt,pseudo->dtctx);CHKERRQ(ierr); 62d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr); 633a40ed3dSBarry Smith PetscFunctionReturn(0); 647bf11e45SBarry Smith } 657bf11e45SBarry Smith 667bf11e45SBarry Smith 677bf11e45SBarry Smith /* ------------------------------------------------------------------------------*/ 684a2ae208SSatish Balay #undef __FUNCT__ 698d359177SBarry Smith #define __FUNCT__ "TSPseudoVerifyTimeStepDefault" 707bf11e45SBarry Smith /*@C 718d359177SBarry Smith TSPseudoVerifyTimeStepDefault - Default code to verify the quality of the last timestep. 727bf11e45SBarry Smith 7315091d37SBarry Smith Collective on TS 7415091d37SBarry Smith 757bf11e45SBarry Smith Input Parameters: 7615091d37SBarry Smith + ts - the timestep context 777bf11e45SBarry Smith . dtctx - unused timestep context 7815091d37SBarry Smith - update - latest solution vector 797bf11e45SBarry Smith 80564e8f4eSLois Curfman McInnes Output Parameters: 8115091d37SBarry Smith + newdt - the timestep to use for the next step 8215091d37SBarry Smith - flag - flag indicating whether the last time step was acceptable 837bf11e45SBarry Smith 8415091d37SBarry Smith Level: advanced 85fee21e36SBarry Smith 86564e8f4eSLois Curfman McInnes Note: 87564e8f4eSLois Curfman McInnes This routine always returns a flag of 1, indicating an acceptable 88564e8f4eSLois Curfman McInnes timestep. 89564e8f4eSLois Curfman McInnes 90564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, default, verify 91564e8f4eSLois Curfman McInnes 92564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStep() 937bf11e45SBarry Smith @*/ 948d359177SBarry Smith PetscErrorCode TSPseudoVerifyTimeStepDefault(TS ts,Vec update,void *dtctx,PetscReal *newdt,PetscBool *flag) 957bf11e45SBarry Smith { 963a40ed3dSBarry Smith PetscFunctionBegin; 97a7cc72afSBarry Smith *flag = PETSC_TRUE; 983a40ed3dSBarry Smith PetscFunctionReturn(0); 997bf11e45SBarry Smith } 1007bf11e45SBarry Smith 1017bf11e45SBarry Smith 1024a2ae208SSatish Balay #undef __FUNCT__ 1034a2ae208SSatish Balay #define __FUNCT__ "TSPseudoVerifyTimeStep" 1047bf11e45SBarry Smith /*@ 105564e8f4eSLois Curfman McInnes TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable. 1067bf11e45SBarry Smith 10715091d37SBarry Smith Collective on TS 10815091d37SBarry Smith 109fb4a63b6SLois Curfman McInnes Input Parameters: 11015091d37SBarry Smith + ts - timestep context 11115091d37SBarry Smith - update - latest solution vector 1127bf11e45SBarry Smith 113fb4a63b6SLois Curfman McInnes Output Parameters: 11415091d37SBarry Smith + dt - newly computed timestep (if it had to shrink) 11515091d37SBarry Smith - flag - indicates if current timestep was ok 1167bf11e45SBarry Smith 11715091d37SBarry Smith Level: advanced 118fee21e36SBarry Smith 119564e8f4eSLois Curfman McInnes Notes: 120564e8f4eSLois Curfman McInnes The routine to be called here to compute the timestep should be 121564e8f4eSLois Curfman McInnes set by calling TSPseudoSetVerifyTimeStep(). 122564e8f4eSLois Curfman McInnes 123564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, verify 124564e8f4eSLois Curfman McInnes 1258d359177SBarry Smith .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStepDefault() 1267bf11e45SBarry Smith @*/ 1277087cfbeSBarry Smith PetscErrorCode TSPseudoVerifyTimeStep(TS ts,Vec update,PetscReal *dt,PetscBool *flag) 1287bf11e45SBarry Smith { 1297bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 130dfbe8321SBarry Smith PetscErrorCode ierr; 1317bf11e45SBarry Smith 1323a40ed3dSBarry Smith PetscFunctionBegin; 1337bf11e45SBarry Smith 134cb9d8021SPierre Barbier de Reuille *flag = PETSC_TRUE; 135d183316bSPierre Barbier de Reuille ierr = TSFunctionDomainError(ts,ts->ptime,update,flag);CHKERRQ(ierr); 136cb9d8021SPierre Barbier de Reuille if(*flag && pseudo->verify) { 1377bf11e45SBarry Smith ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag);CHKERRQ(ierr); 138cb9d8021SPierre Barbier de Reuille } 1393a40ed3dSBarry Smith PetscFunctionReturn(0); 1407bf11e45SBarry Smith } 1417bf11e45SBarry Smith 1427bf11e45SBarry Smith /* --------------------------------------------------------------------------------*/ 1437bf11e45SBarry Smith 1444a2ae208SSatish Balay #undef __FUNCT__ 1454a2ae208SSatish Balay #define __FUNCT__ "TSStep_Pseudo" 146193ac0bcSJed Brown static PetscErrorCode TSStep_Pseudo(TS ts) 1472d3f70b5SBarry Smith { 148277b19d0SLisandro Dalcin TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 149cdbf8f93SLisandro Dalcin PetscInt its,lits,reject; 150cdbf8f93SLisandro Dalcin PetscBool stepok; 151cdbf8f93SLisandro Dalcin PetscReal next_time_step; 152cdbf8f93SLisandro Dalcin SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING; 153dfbe8321SBarry Smith PetscErrorCode ierr; 1542d3f70b5SBarry Smith 1553a40ed3dSBarry Smith PetscFunctionBegin; 156bbd56ea5SKarl Rupp if (ts->steps == 0) pseudo->dt_initial = ts->time_step; 157193ac0bcSJed Brown ierr = VecCopy(ts->vec_sol,pseudo->update);CHKERRQ(ierr); 158cdbf8f93SLisandro Dalcin next_time_step = ts->time_step; 159cdbf8f93SLisandro Dalcin ierr = TSPseudoComputeTimeStep(ts,&next_time_step);CHKERRQ(ierr); 160cdbf8f93SLisandro Dalcin for (reject=0; reject<ts->max_reject; reject++,ts->reject++) { 161cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 162b8123daeSJed Brown ierr = TSPreStep(ts);CHKERRQ(ierr); 163b8123daeSJed Brown ierr = TSPreStage(ts,ts->ptime+ts->time_step);CHKERRQ(ierr); 1640298fd71SBarry Smith ierr = SNESSolve(ts->snes,NULL,pseudo->update);CHKERRQ(ierr); 165cdbf8f93SLisandro Dalcin ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr); 166b850b91aSLisandro Dalcin ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 1679f954729SBarry Smith ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr); 1689be3e283SDebojyoti Ghosh ierr = TSPostStage(ts,ts->ptime+ts->time_step,0,&(pseudo->update));CHKERRQ(ierr); 1695ef26d82SJed Brown ts->snes_its += its; ts->ksp_its += lits; 170cdbf8f93SLisandro Dalcin ierr = PetscInfo3(ts,"step=%D, nonlinear solve iterations=%D, linear solve iterations=%D\n",ts->steps,its,lits);CHKERRQ(ierr); 171193ac0bcSJed Brown pseudo->fnorm = -1; /* The current norm is no longer valid, monitor must recompute it. */ 172cdbf8f93SLisandro Dalcin ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&next_time_step,&stepok);CHKERRQ(ierr); 173cdbf8f93SLisandro Dalcin if (stepok) break; 174cdbf8f93SLisandro Dalcin } 175f1b97656SJed Brown if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) { 176193ac0bcSJed Brown ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 177193ac0bcSJed 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); 178cdbf8f93SLisandro Dalcin PetscFunctionReturn(0); 1797bf11e45SBarry Smith } 1803118ae5eSBarry Smith if (pseudo->fnorm < 0) { 1813118ae5eSBarry Smith ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr); 1823118ae5eSBarry Smith ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr); 1833118ae5eSBarry Smith ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 1843118ae5eSBarry Smith } 1853118ae5eSBarry Smith if (pseudo->fnorm < pseudo->fatol) { 1863118ae5eSBarry Smith ts->reason = TS_CONVERGED_PSEUDO_FATOL; 1873118ae5eSBarry Smith ierr = PetscInfo3(ts,"step=%D, converged since fnorm %g < fatol %g\n",ts->steps,pseudo->fnorm,pseudo->frtol);CHKERRQ(ierr); 1883118ae5eSBarry Smith PetscFunctionReturn(0); 1893118ae5eSBarry Smith } 1903118ae5eSBarry Smith if (pseudo->fnorm/pseudo->fnorm_initial < pseudo->frtol) { 1913118ae5eSBarry Smith ts->reason = TS_CONVERGED_PSEUDO_FRTOL; 1923118ae5eSBarry Smith 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); 1933118ae5eSBarry Smith PetscFunctionReturn(0); 1943118ae5eSBarry Smith } 195cdbf8f93SLisandro Dalcin if (reject >= ts->max_reject) { 196193ac0bcSJed Brown ts->reason = TS_DIVERGED_STEP_REJECTED; 197cdbf8f93SLisandro Dalcin ierr = PetscInfo2(ts,"step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,reject);CHKERRQ(ierr); 198cdbf8f93SLisandro Dalcin PetscFunctionReturn(0); 199193ac0bcSJed Brown } 200cdbf8f93SLisandro Dalcin ierr = VecCopy(pseudo->update,ts->vec_sol);CHKERRQ(ierr); 201193ac0bcSJed Brown ts->ptime += ts->time_step; 202cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 2032d3f70b5SBarry Smith ts->steps++; 2043a40ed3dSBarry Smith PetscFunctionReturn(0); 2052d3f70b5SBarry Smith } 2062d3f70b5SBarry Smith 2072d3f70b5SBarry Smith /*------------------------------------------------------------*/ 2084a2ae208SSatish Balay #undef __FUNCT__ 209277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Pseudo" 210277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Pseudo(TS ts) 2112d3f70b5SBarry Smith { 2127bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 213dfbe8321SBarry Smith PetscErrorCode ierr; 2142d3f70b5SBarry Smith 2153a40ed3dSBarry Smith PetscFunctionBegin; 2166bf464f9SBarry Smith ierr = VecDestroy(&pseudo->update);CHKERRQ(ierr); 2176bf464f9SBarry Smith ierr = VecDestroy(&pseudo->func);CHKERRQ(ierr); 2186bf464f9SBarry Smith ierr = VecDestroy(&pseudo->xdot);CHKERRQ(ierr); 2193a40ed3dSBarry Smith PetscFunctionReturn(0); 2202d3f70b5SBarry Smith } 2212d3f70b5SBarry Smith 222277b19d0SLisandro Dalcin #undef __FUNCT__ 223277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Pseudo" 224277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Pseudo(TS ts) 225277b19d0SLisandro Dalcin { 226277b19d0SLisandro Dalcin PetscErrorCode ierr; 227277b19d0SLisandro Dalcin 228277b19d0SLisandro Dalcin PetscFunctionBegin; 229277b19d0SLisandro Dalcin ierr = TSReset_Pseudo(ts);CHKERRQ(ierr); 230277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 231bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",NULL);CHKERRQ(ierr); 232bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",NULL);CHKERRQ(ierr); 233bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetMaxTimeStep_C",NULL);CHKERRQ(ierr); 234bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",NULL);CHKERRQ(ierr); 235bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",NULL);CHKERRQ(ierr); 236277b19d0SLisandro Dalcin PetscFunctionReturn(0); 237277b19d0SLisandro Dalcin } 2382d3f70b5SBarry Smith 2392d3f70b5SBarry Smith /*------------------------------------------------------------*/ 2402d3f70b5SBarry Smith 2414a2ae208SSatish Balay #undef __FUNCT__ 2426f2d6a7bSJed Brown #define __FUNCT__ "TSPseudoGetXdot" 2436f2d6a7bSJed Brown /* 2446f2d6a7bSJed Brown Compute Xdot = (X^{n+1}-X^n)/dt) = 0 2456f2d6a7bSJed Brown */ 2466f2d6a7bSJed Brown static PetscErrorCode TSPseudoGetXdot(TS ts,Vec X,Vec *Xdot) 2472d3f70b5SBarry Smith { 2486f2d6a7bSJed Brown TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 249193ac0bcSJed Brown const PetscScalar mdt = 1.0/ts->time_step,*xnp1,*xn; 250193ac0bcSJed Brown PetscScalar *xdot; 251dfbe8321SBarry Smith PetscErrorCode ierr; 252a7cc72afSBarry Smith PetscInt i,n; 2532d3f70b5SBarry Smith 2543a40ed3dSBarry Smith PetscFunctionBegin; 255193ac0bcSJed Brown ierr = VecGetArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr); 256193ac0bcSJed Brown ierr = VecGetArrayRead(X,&xnp1);CHKERRQ(ierr); 2576f2d6a7bSJed Brown ierr = VecGetArray(pseudo->xdot,&xdot);CHKERRQ(ierr); 2586f2d6a7bSJed Brown ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr); 259bbd56ea5SKarl Rupp for (i=0; i<n; i++) xdot[i] = mdt*(xnp1[i] - xn[i]); 260193ac0bcSJed Brown ierr = VecRestoreArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr); 261193ac0bcSJed Brown ierr = VecRestoreArrayRead(X,&xnp1);CHKERRQ(ierr); 2626f2d6a7bSJed Brown ierr = VecRestoreArray(pseudo->xdot,&xdot);CHKERRQ(ierr); 2636f2d6a7bSJed Brown *Xdot = pseudo->xdot; 2643a40ed3dSBarry Smith PetscFunctionReturn(0); 2652d3f70b5SBarry Smith } 2662d3f70b5SBarry Smith 2674a2ae208SSatish Balay #undef __FUNCT__ 2680f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Pseudo" 2696f2d6a7bSJed Brown /* 2706f2d6a7bSJed Brown The transient residual is 2716f2d6a7bSJed Brown 2726f2d6a7bSJed Brown F(U^{n+1},(U^{n+1}-U^n)/dt) = 0 2736f2d6a7bSJed Brown 2746f2d6a7bSJed Brown or for ODE, 2756f2d6a7bSJed Brown 2766f2d6a7bSJed Brown (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0 2776f2d6a7bSJed Brown 2786f2d6a7bSJed Brown This is the function that must be evaluated for transient simulation and for 2796f2d6a7bSJed Brown finite difference Jacobians. On the first Newton step, this algorithm uses 2806f2d6a7bSJed Brown a guess of U^{n+1} = U^n in which case the transient term vanishes and the 2816f2d6a7bSJed Brown residual is actually the steady state residual. Pseudotransient 2826f2d6a7bSJed Brown continuation as described in the literature is a linearly implicit 2836f2d6a7bSJed Brown algorithm, it only takes this one Newton step with the steady state 2846f2d6a7bSJed Brown residual, and then advances to the next time step. 2856f2d6a7bSJed Brown */ 2860f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes,Vec X,Vec Y,TS ts) 2872d3f70b5SBarry Smith { 2886f2d6a7bSJed Brown Vec Xdot; 289dfbe8321SBarry Smith PetscErrorCode ierr; 2902d3f70b5SBarry Smith 2913a40ed3dSBarry Smith PetscFunctionBegin; 2926f2d6a7bSJed Brown ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr); 293193ac0bcSJed Brown ierr = TSComputeIFunction(ts,ts->ptime+ts->time_step,X,Xdot,Y,PETSC_FALSE);CHKERRQ(ierr); 2946f2d6a7bSJed Brown PetscFunctionReturn(0); 2956f2d6a7bSJed Brown } 2962d3f70b5SBarry Smith 2976f2d6a7bSJed Brown #undef __FUNCT__ 2980f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Pseudo" 2996f2d6a7bSJed Brown /* 3006f2d6a7bSJed Brown This constructs the Jacobian needed for SNES. For DAE, this is 3016f2d6a7bSJed Brown 3026f2d6a7bSJed Brown dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot 3036f2d6a7bSJed Brown 3046f2d6a7bSJed Brown and for ODE: 3056f2d6a7bSJed Brown 3066f2d6a7bSJed Brown J = I/dt - J_{Frhs} where J_{Frhs} is the given Jacobian of Frhs. 3076f2d6a7bSJed Brown */ 308d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes,Vec X,Mat AA,Mat BB,TS ts) 3096f2d6a7bSJed Brown { 3106f2d6a7bSJed Brown Vec Xdot; 3116f2d6a7bSJed Brown PetscErrorCode ierr; 3126f2d6a7bSJed Brown 3136f2d6a7bSJed Brown PetscFunctionBegin; 3146f2d6a7bSJed Brown ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr); 315d1e9a80fSBarry Smith ierr = TSComputeIJacobian(ts,ts->ptime+ts->time_step,X,Xdot,1./ts->time_step,AA,BB,PETSC_FALSE);CHKERRQ(ierr); 3163a40ed3dSBarry Smith PetscFunctionReturn(0); 3172d3f70b5SBarry Smith } 3182d3f70b5SBarry Smith 3192d3f70b5SBarry Smith 3204a2ae208SSatish Balay #undef __FUNCT__ 3214a2ae208SSatish Balay #define __FUNCT__ "TSSetUp_Pseudo" 3226849ba73SBarry Smith static PetscErrorCode TSSetUp_Pseudo(TS ts) 3232d3f70b5SBarry Smith { 3247bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 325dfbe8321SBarry Smith PetscErrorCode ierr; 3262d3f70b5SBarry Smith 3273a40ed3dSBarry Smith PetscFunctionBegin; 3287bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr); 3297bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr); 3306f2d6a7bSJed Brown ierr = VecDuplicate(ts->vec_sol,&pseudo->xdot);CHKERRQ(ierr); 3313a40ed3dSBarry Smith PetscFunctionReturn(0); 3322d3f70b5SBarry Smith } 3332d3f70b5SBarry Smith /*------------------------------------------------------------*/ 3342d3f70b5SBarry Smith 3354a2ae208SSatish Balay #undef __FUNCT__ 336a6570f20SBarry Smith #define __FUNCT__ "TSPseudoMonitorDefault" 337649052a6SBarry Smith PetscErrorCode TSPseudoMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy) 3382d3f70b5SBarry Smith { 3397bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 340dfbe8321SBarry Smith PetscErrorCode ierr; 341ce94432eSBarry Smith PetscViewer viewer = (PetscViewer) dummy; 3422d3f70b5SBarry Smith 3433a40ed3dSBarry Smith PetscFunctionBegin; 344193ac0bcSJed Brown if (pseudo->fnorm < 0) { /* The last computed norm is stale, recompute */ 345193ac0bcSJed Brown ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr); 346193ac0bcSJed Brown ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr); 347193ac0bcSJed Brown ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 348193ac0bcSJed Brown } 349649052a6SBarry Smith ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 3507c8652ddSBarry 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); 351649052a6SBarry Smith ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 3523a40ed3dSBarry Smith PetscFunctionReturn(0); 3532d3f70b5SBarry Smith } 3542d3f70b5SBarry Smith 3554a2ae208SSatish Balay #undef __FUNCT__ 3564a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Pseudo" 3574416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Pseudo(PetscOptionItems *PetscOptionsObject,TS ts) 3582d3f70b5SBarry Smith { 3594bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 360dfbe8321SBarry Smith PetscErrorCode ierr; 361ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE; 362649052a6SBarry Smith PetscViewer viewer; 3632d3f70b5SBarry Smith 3643a40ed3dSBarry Smith PetscFunctionBegin; 365e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Pseudo-timestepping options");CHKERRQ(ierr); 3660298fd71SBarry Smith ierr = PetscOptionsBool("-ts_monitor_pseudo","Monitor convergence","TSPseudoMonitorDefault",flg,&flg,NULL);CHKERRQ(ierr); 3672d3f70b5SBarry Smith if (flg) { 368ce94432eSBarry Smith ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ts),"stdout",&viewer);CHKERRQ(ierr); 369649052a6SBarry Smith ierr = TSMonitorSet(ts,TSPseudoMonitorDefault,viewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 37028aa8177SBarry Smith } 37190d69ab7SBarry Smith flg = PETSC_FALSE; 3720298fd71SBarry Smith ierr = PetscOptionsBool("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",flg,&flg,NULL);CHKERRQ(ierr); 373ca4b7087SBarry Smith if (flg) { 374ca4b7087SBarry Smith ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr); 375ca4b7087SBarry Smith } 37694ae4db5SBarry Smith ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,NULL);CHKERRQ(ierr); 37794ae4db5SBarry Smith ierr = PetscOptionsReal("-ts_pseudo_max_dt","Maximum value for dt","TSPseudoSetMaxTimeStep",pseudo->dt_max,&pseudo->dt_max,NULL);CHKERRQ(ierr); 3783118ae5eSBarry Smith ierr = PetscOptionsReal("-ts_pseudo_fatol","Tolerance for norm of function","",pseudo->fatol,&pseudo->fatol,NULL);CHKERRQ(ierr); 3793118ae5eSBarry Smith ierr = PetscOptionsReal("-ts_pseudo_frtol","Relative tolerance for norm of function","",pseudo->frtol,&pseudo->frtol,NULL);CHKERRQ(ierr); 380b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 3813a40ed3dSBarry Smith PetscFunctionReturn(0); 3822d3f70b5SBarry Smith } 3832d3f70b5SBarry Smith 3844a2ae208SSatish Balay #undef __FUNCT__ 3854a2ae208SSatish Balay #define __FUNCT__ "TSView_Pseudo" 3866849ba73SBarry Smith static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer) 3872d3f70b5SBarry Smith { 388d52bd9f3SBarry Smith PetscErrorCode ierr; 3893118ae5eSBarry Smith PetscBool isascii; 390d52bd9f3SBarry Smith 3913a40ed3dSBarry Smith PetscFunctionBegin; 3923118ae5eSBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 3933118ae5eSBarry Smith if (isascii) { 3943118ae5eSBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 3953118ae5eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Parameters for pseudo timestepping\n");CHKERRQ(ierr); 3963118ae5eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," frtol - relative tolerance in function value %g\n",(double)pseudo->frtol);CHKERRQ(ierr); 3973118ae5eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," fatol - absolute tolerance in function value %g\n",(double)pseudo->fatol);CHKERRQ(ierr); 3983118ae5eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," dt_initial - initial timestep %g\n",(double)pseudo->dt_initial);CHKERRQ(ierr); 3993118ae5eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," dt_increment - increase in timestep on successful step %g\n",(double)pseudo->dt_increment);CHKERRQ(ierr); 4003118ae5eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," dt_max - maximum time %g\n",(double)pseudo->dt_max);CHKERRQ(ierr); 4013118ae5eSBarry Smith } 402d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 4033a40ed3dSBarry Smith PetscFunctionReturn(0); 4042d3f70b5SBarry Smith } 4052d3f70b5SBarry Smith 40682bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 4074a2ae208SSatish Balay #undef __FUNCT__ 4084a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep" 409ac226902SBarry Smith /*@C 41082bf6240SBarry Smith TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the 41182bf6240SBarry Smith last timestep. 41282bf6240SBarry Smith 4133f9fe445SBarry Smith Logically Collective on TS 41415091d37SBarry Smith 41582bf6240SBarry Smith Input Parameters: 41615091d37SBarry Smith + ts - timestep context 41782bf6240SBarry Smith . dt - user-defined function to verify timestep 41815091d37SBarry Smith - ctx - [optional] user-defined context for private data 4190298fd71SBarry Smith for the timestep verification routine (may be NULL) 42082bf6240SBarry Smith 42115091d37SBarry Smith Level: advanced 422fee21e36SBarry Smith 42382bf6240SBarry Smith Calling sequence of func: 424ace3abfcSBarry Smith . func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscBool *flag); 42582bf6240SBarry Smith 42682bf6240SBarry Smith . update - latest solution vector 42782bf6240SBarry Smith . ctx - [optional] timestep context 42882bf6240SBarry Smith . newdt - the timestep to use for the next step 42982bf6240SBarry Smith . flag - flag indicating whether the last time step was acceptable 43082bf6240SBarry Smith 43182bf6240SBarry Smith Notes: 43282bf6240SBarry Smith The routine set here will be called by TSPseudoVerifyTimeStep() 43382bf6240SBarry Smith during the timestepping process. 43482bf6240SBarry Smith 43582bf6240SBarry Smith .keywords: timestep, pseudo, set, verify 43682bf6240SBarry Smith 4378d359177SBarry Smith .seealso: TSPseudoVerifyTimeStepDefault(), TSPseudoVerifyTimeStep() 43882bf6240SBarry Smith @*/ 4397087cfbeSBarry Smith PetscErrorCode TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,PetscBool*),void *ctx) 44082bf6240SBarry Smith { 4414ac538c5SBarry Smith PetscErrorCode ierr; 44282bf6240SBarry Smith 44382bf6240SBarry Smith PetscFunctionBegin; 4440700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 4454ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSPseudoSetVerifyTimeStep_C",(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal*,PetscBool*),void*),(ts,dt,ctx));CHKERRQ(ierr); 44682bf6240SBarry Smith PetscFunctionReturn(0); 44782bf6240SBarry Smith } 44882bf6240SBarry Smith 4494a2ae208SSatish Balay #undef __FUNCT__ 4504a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement" 45182bf6240SBarry Smith /*@ 45282bf6240SBarry Smith TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 4538d359177SBarry Smith dt when using the TSPseudoTimeStepDefault() routine. 45482bf6240SBarry Smith 4553f9fe445SBarry Smith Logically Collective on TS 456fee21e36SBarry Smith 45715091d37SBarry Smith Input Parameters: 45815091d37SBarry Smith + ts - the timestep context 45915091d37SBarry Smith - inc - the scaling factor >= 1.0 46015091d37SBarry Smith 46182bf6240SBarry Smith Options Database Key: 462e1bc860dSBarry Smith . -ts_pseudo_increment <increment> 46382bf6240SBarry Smith 46415091d37SBarry Smith Level: advanced 46515091d37SBarry Smith 46682bf6240SBarry Smith .keywords: timestep, pseudo, set, increment 46782bf6240SBarry Smith 4688d359177SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault() 46982bf6240SBarry Smith @*/ 4707087cfbeSBarry Smith PetscErrorCode TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc) 47182bf6240SBarry Smith { 4724ac538c5SBarry Smith PetscErrorCode ierr; 47382bf6240SBarry Smith 47482bf6240SBarry Smith PetscFunctionBegin; 4750700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 476c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts,inc,2); 4774ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSPseudoSetTimeStepIncrement_C",(TS,PetscReal),(ts,inc));CHKERRQ(ierr); 47882bf6240SBarry Smith PetscFunctionReturn(0); 47982bf6240SBarry Smith } 48082bf6240SBarry Smith 4814a2ae208SSatish Balay #undef __FUNCT__ 48286534af1SJed Brown #define __FUNCT__ "TSPseudoSetMaxTimeStep" 48386534af1SJed Brown /*@ 48486534af1SJed Brown TSPseudoSetMaxTimeStep - Sets the maximum time step 4858d359177SBarry Smith when using the TSPseudoTimeStepDefault() routine. 48686534af1SJed Brown 48786534af1SJed Brown Logically Collective on TS 48886534af1SJed Brown 48986534af1SJed Brown Input Parameters: 49086534af1SJed Brown + ts - the timestep context 49186534af1SJed Brown - maxdt - the maximum time step, use a non-positive value to deactivate 49286534af1SJed Brown 49386534af1SJed Brown Options Database Key: 494e1bc860dSBarry Smith . -ts_pseudo_max_dt <increment> 49586534af1SJed Brown 49686534af1SJed Brown Level: advanced 49786534af1SJed Brown 49886534af1SJed Brown .keywords: timestep, pseudo, set 49986534af1SJed Brown 5008d359177SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault() 50186534af1SJed Brown @*/ 50286534af1SJed Brown PetscErrorCode TSPseudoSetMaxTimeStep(TS ts,PetscReal maxdt) 50386534af1SJed Brown { 50486534af1SJed Brown PetscErrorCode ierr; 50586534af1SJed Brown 50686534af1SJed Brown PetscFunctionBegin; 50786534af1SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 50886534af1SJed Brown PetscValidLogicalCollectiveReal(ts,maxdt,2); 50986534af1SJed Brown ierr = PetscTryMethod(ts,"TSPseudoSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));CHKERRQ(ierr); 51086534af1SJed Brown PetscFunctionReturn(0); 51186534af1SJed Brown } 51286534af1SJed Brown 51386534af1SJed Brown #undef __FUNCT__ 5144a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt" 51582bf6240SBarry Smith /*@ 51682bf6240SBarry Smith TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep 51782bf6240SBarry Smith is computed via the formula 51882bf6240SBarry Smith $ dt = initial_dt*initial_fnorm/current_fnorm 51982bf6240SBarry Smith rather than the default update, 52082bf6240SBarry Smith $ dt = current_dt*previous_fnorm/current_fnorm. 52182bf6240SBarry Smith 5223f9fe445SBarry Smith Logically Collective on TS 52315091d37SBarry Smith 52482bf6240SBarry Smith Input Parameter: 52582bf6240SBarry Smith . ts - the timestep context 52682bf6240SBarry Smith 52782bf6240SBarry Smith Options Database Key: 528e1bc860dSBarry Smith . -ts_pseudo_increment_dt_from_initial_dt 52982bf6240SBarry Smith 53015091d37SBarry Smith Level: advanced 53115091d37SBarry Smith 53282bf6240SBarry Smith .keywords: timestep, pseudo, set, increment 53382bf6240SBarry Smith 5348d359177SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault() 53582bf6240SBarry Smith @*/ 5367087cfbeSBarry Smith PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS ts) 53782bf6240SBarry Smith { 5384ac538c5SBarry Smith PetscErrorCode ierr; 53982bf6240SBarry Smith 54082bf6240SBarry Smith PetscFunctionBegin; 5410700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 5424ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSPseudoIncrementDtFromInitialDt_C",(TS),(ts));CHKERRQ(ierr); 54382bf6240SBarry Smith PetscFunctionReturn(0); 54482bf6240SBarry Smith } 54582bf6240SBarry Smith 54682bf6240SBarry Smith 5474a2ae208SSatish Balay #undef __FUNCT__ 5484a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep" 549ac226902SBarry Smith /*@C 55082bf6240SBarry Smith TSPseudoSetTimeStep - Sets the user-defined routine to be 55182bf6240SBarry Smith called at each pseudo-timestep to update the timestep. 55282bf6240SBarry Smith 5533f9fe445SBarry Smith Logically Collective on TS 55415091d37SBarry Smith 55582bf6240SBarry Smith Input Parameters: 55615091d37SBarry Smith + ts - timestep context 55782bf6240SBarry Smith . dt - function to compute timestep 55815091d37SBarry Smith - ctx - [optional] user-defined context for private data 5590298fd71SBarry Smith required by the function (may be NULL) 56082bf6240SBarry Smith 56115091d37SBarry Smith Level: intermediate 562fee21e36SBarry Smith 56382bf6240SBarry Smith Calling sequence of func: 56487828ca2SBarry Smith . func (TS ts,PetscReal *newdt,void *ctx); 56582bf6240SBarry Smith 56682bf6240SBarry Smith . newdt - the newly computed timestep 56782bf6240SBarry Smith . ctx - [optional] timestep context 56882bf6240SBarry Smith 56982bf6240SBarry Smith Notes: 57082bf6240SBarry Smith The routine set here will be called by TSPseudoComputeTimeStep() 57182bf6240SBarry Smith during the timestepping process. 5728d359177SBarry Smith If not set then TSPseudoTimeStepDefault() is automatically used 57382bf6240SBarry Smith 57482bf6240SBarry Smith .keywords: timestep, pseudo, set 57582bf6240SBarry Smith 5768d359177SBarry Smith .seealso: TSPseudoTimeStepDefault(), TSPseudoComputeTimeStep() 57782bf6240SBarry Smith @*/ 5787087cfbeSBarry Smith PetscErrorCode TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void *ctx) 57982bf6240SBarry Smith { 5804ac538c5SBarry Smith PetscErrorCode ierr; 58182bf6240SBarry Smith 58282bf6240SBarry Smith PetscFunctionBegin; 5830700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 5844ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSPseudoSetTimeStep_C",(TS,PetscErrorCode (*)(TS,PetscReal*,void*),void*),(ts,dt,ctx));CHKERRQ(ierr); 58582bf6240SBarry Smith PetscFunctionReturn(0); 58682bf6240SBarry Smith } 58782bf6240SBarry Smith 58882bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 58982bf6240SBarry Smith 590ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscBool*); /* force argument to next function to not be extern C*/ 5914a2ae208SSatish Balay #undef __FUNCT__ 5924a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo" 5937087cfbeSBarry Smith PetscErrorCode TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void *ctx) 59482bf6240SBarry Smith { 59582bf6240SBarry Smith TS_Pseudo *pseudo; 59682bf6240SBarry Smith 59782bf6240SBarry Smith PetscFunctionBegin; 59882bf6240SBarry Smith pseudo = (TS_Pseudo*)ts->data; 59982bf6240SBarry Smith pseudo->verify = dt; 60082bf6240SBarry Smith pseudo->verifyctx = ctx; 60182bf6240SBarry Smith PetscFunctionReturn(0); 60282bf6240SBarry Smith } 60382bf6240SBarry Smith 6044a2ae208SSatish Balay #undef __FUNCT__ 6054a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement_Pseudo" 6067087cfbeSBarry Smith PetscErrorCode TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc) 60782bf6240SBarry Smith { 6084bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 60982bf6240SBarry Smith 61082bf6240SBarry Smith PetscFunctionBegin; 61182bf6240SBarry Smith pseudo->dt_increment = inc; 61282bf6240SBarry Smith PetscFunctionReturn(0); 61382bf6240SBarry Smith } 61482bf6240SBarry Smith 6154a2ae208SSatish Balay #undef __FUNCT__ 61686534af1SJed Brown #define __FUNCT__ "TSPseudoSetMaxTimeStep_Pseudo" 61786534af1SJed Brown PetscErrorCode TSPseudoSetMaxTimeStep_Pseudo(TS ts,PetscReal maxdt) 61886534af1SJed Brown { 61986534af1SJed Brown TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 62086534af1SJed Brown 62186534af1SJed Brown PetscFunctionBegin; 62286534af1SJed Brown pseudo->dt_max = maxdt; 62386534af1SJed Brown PetscFunctionReturn(0); 62486534af1SJed Brown } 62586534af1SJed Brown 62686534af1SJed Brown #undef __FUNCT__ 6274a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt_Pseudo" 6287087cfbeSBarry Smith PetscErrorCode TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts) 62982bf6240SBarry Smith { 6304bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 63182bf6240SBarry Smith 63282bf6240SBarry Smith PetscFunctionBegin; 6334bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_TRUE; 63482bf6240SBarry Smith PetscFunctionReturn(0); 63582bf6240SBarry Smith } 63682bf6240SBarry Smith 6376849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/ 6384a2ae208SSatish Balay #undef __FUNCT__ 6394a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo" 6407087cfbeSBarry Smith PetscErrorCode TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void *ctx) 64182bf6240SBarry Smith { 6424bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 64382bf6240SBarry Smith 64482bf6240SBarry Smith PetscFunctionBegin; 64582bf6240SBarry Smith pseudo->dt = dt; 64682bf6240SBarry Smith pseudo->dtctx = ctx; 64782bf6240SBarry Smith PetscFunctionReturn(0); 64882bf6240SBarry Smith } 64982bf6240SBarry Smith 65082bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 65110e6a065SJed Brown /*MC 65210e6a065SJed Brown TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping 65382bf6240SBarry Smith 65410e6a065SJed Brown This method solves equations of the form 65510e6a065SJed Brown 65610e6a065SJed Brown $ F(X,Xdot) = 0 65710e6a065SJed Brown 65810e6a065SJed Brown for steady state using the iteration 65910e6a065SJed Brown 66010e6a065SJed Brown $ [G'] S = -F(X,0) 66110e6a065SJed Brown $ X += S 66210e6a065SJed Brown 66310e6a065SJed Brown where 66410e6a065SJed Brown 66510e6a065SJed Brown $ G(Y) = F(Y,(Y-X)/dt) 66610e6a065SJed Brown 6676f2d6a7bSJed Brown This is linearly-implicit Euler with the residual always evaluated "at steady 6686f2d6a7bSJed Brown state". See note below. 66910e6a065SJed Brown 67010e6a065SJed Brown Options database keys: 67110e6a065SJed Brown + -ts_pseudo_increment <real> - ratio of increase dt 6723118ae5eSBarry Smith . -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt 6733118ae5eSBarry Smith . -ts_pseudo_fatol <atol> - stop iterating when the function norm is less than atol 6743118ae5eSBarry Smith - -ts_pseudo_frtol <rtol> - stop iterating when the function norm divided by the initial function norm is less than rtol 67510e6a065SJed Brown 67610e6a065SJed Brown Level: beginner 67710e6a065SJed Brown 67810e6a065SJed Brown References: 679*96a0c994SBarry Smith + 1. - Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential Algebraic Equations, 2003. 680*96a0c994SBarry Smith - 2. - C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998. 68110e6a065SJed Brown 68210e6a065SJed Brown Notes: 6836f2d6a7bSJed Brown The residual computed by this method includes the transient term (Xdot is computed instead of 6846f2d6a7bSJed Brown always being zero), but since the prediction from the last step is always the solution from the 6856f2d6a7bSJed Brown last step, on the first Newton iteration we have 6866f2d6a7bSJed Brown 6876f2d6a7bSJed Brown $ Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0 6886f2d6a7bSJed Brown 6896f2d6a7bSJed Brown Therefore, the linear system solved by the first Newton iteration is equivalent to the one 6906f2d6a7bSJed Brown described above and in the papers. If the user chooses to perform multiple Newton iterations, the 6916f2d6a7bSJed Brown algorithm is no longer the one described in the referenced papers. 69210e6a065SJed Brown 69310e6a065SJed Brown .seealso: TSCreate(), TS, TSSetType() 69410e6a065SJed Brown 69510e6a065SJed Brown M*/ 6964a2ae208SSatish Balay #undef __FUNCT__ 6974a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Pseudo" 6988cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Pseudo(TS ts) 6992d3f70b5SBarry Smith { 7007bf11e45SBarry Smith TS_Pseudo *pseudo; 701dfbe8321SBarry Smith PetscErrorCode ierr; 702193ac0bcSJed Brown SNES snes; 70319fd82e9SBarry Smith SNESType stype; 7042d3f70b5SBarry Smith 7053a40ed3dSBarry Smith PetscFunctionBegin; 706277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Pseudo; 707000e7ae3SMatthew Knepley ts->ops->destroy = TSDestroy_Pseudo; 708000e7ae3SMatthew Knepley ts->ops->view = TSView_Pseudo; 7092d3f70b5SBarry Smith 710000e7ae3SMatthew Knepley ts->ops->setup = TSSetUp_Pseudo; 711000e7ae3SMatthew Knepley ts->ops->step = TSStep_Pseudo; 712000e7ae3SMatthew Knepley ts->ops->setfromoptions = TSSetFromOptions_Pseudo; 7130f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Pseudo; 7140f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Pseudo; 7157bf11e45SBarry Smith 716193ac0bcSJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 717193ac0bcSJed Brown ierr = SNESGetType(snes,&stype);CHKERRQ(ierr); 718193ac0bcSJed Brown if (!stype) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);} 7192d3f70b5SBarry Smith 720b00a9115SJed Brown ierr = PetscNewLog(ts,&pseudo);CHKERRQ(ierr); 7217bf11e45SBarry Smith ts->data = (void*)pseudo; 7222d3f70b5SBarry Smith 72328aa8177SBarry Smith pseudo->dt_increment = 1.1; 7244bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_FALSE; 7258d359177SBarry Smith pseudo->dt = TSPseudoTimeStepDefault; 726193ac0bcSJed Brown pseudo->fnorm = -1; 7273118ae5eSBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 7283118ae5eSBarry Smith pseudo->fatol = 1.e-25; 7293118ae5eSBarry Smith pseudo->frtol = 1.e-5; 7303118ae5eSBarry Smith #else 7313118ae5eSBarry Smith pseudo->fatol = 1.e-50; 7323118ae5eSBarry Smith pseudo->frtol = 1.e-12; 7333118ae5eSBarry Smith #endif 734bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr); 735bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr); 736bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetMaxTimeStep_C",TSPseudoSetMaxTimeStep_Pseudo);CHKERRQ(ierr); 737bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr); 738bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr); 7393a40ed3dSBarry Smith PetscFunctionReturn(0); 7402d3f70b5SBarry Smith } 7412d3f70b5SBarry Smith 7424a2ae208SSatish Balay #undef __FUNCT__ 7438d359177SBarry Smith #define __FUNCT__ "TSPseudoTimeStepDefault" 74482bf6240SBarry Smith /*@C 7458d359177SBarry Smith TSPseudoTimeStepDefault - Default code to compute pseudo-timestepping. 74682bf6240SBarry Smith Use with TSPseudoSetTimeStep(). 74728aa8177SBarry Smith 74815091d37SBarry Smith Collective on TS 74915091d37SBarry Smith 75028aa8177SBarry Smith Input Parameters: 75128aa8177SBarry Smith . ts - the timestep context 75282bf6240SBarry Smith . dtctx - unused timestep context 75328aa8177SBarry Smith 75482bf6240SBarry Smith Output Parameter: 75582bf6240SBarry Smith . newdt - the timestep to use for the next step 75628aa8177SBarry Smith 75715091d37SBarry Smith Level: advanced 75815091d37SBarry Smith 75982bf6240SBarry Smith .keywords: timestep, pseudo, default 760564e8f4eSLois Curfman McInnes 76182bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep() 76228aa8177SBarry Smith @*/ 7638d359177SBarry Smith PetscErrorCode TSPseudoTimeStepDefault(TS ts,PetscReal *newdt,void *dtctx) 76428aa8177SBarry Smith { 76582bf6240SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 76687828ca2SBarry Smith PetscReal inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous; 767dfbe8321SBarry Smith PetscErrorCode ierr; 76828aa8177SBarry Smith 7693a40ed3dSBarry Smith PetscFunctionBegin; 770bbd7b040SJed Brown ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr); 771214bc6a2SJed Brown ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr); 77282bf6240SBarry Smith ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 773cdbf8f93SLisandro Dalcin if (pseudo->fnorm_initial == 0.0) { 77482bf6240SBarry Smith /* first time through so compute initial function norm */ 775cdbf8f93SLisandro Dalcin pseudo->fnorm_initial = pseudo->fnorm; 77682bf6240SBarry Smith fnorm_previous = pseudo->fnorm; 77782bf6240SBarry Smith } 778bbd56ea5SKarl Rupp if (pseudo->fnorm == 0.0) *newdt = 1.e12*inc*ts->time_step; 779bbd56ea5SKarl Rupp else if (pseudo->increment_dt_from_initial_dt) *newdt = inc*pseudo->dt_initial*pseudo->fnorm_initial/pseudo->fnorm; 780bbd56ea5SKarl Rupp else *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm; 78186534af1SJed Brown if (pseudo->dt_max > 0) *newdt = PetscMin(*newdt,pseudo->dt_max); 78282bf6240SBarry Smith pseudo->fnorm_previous = pseudo->fnorm; 7833a40ed3dSBarry Smith PetscFunctionReturn(0); 78428aa8177SBarry Smith } 785