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" 318d359177SBarry Smith /*@C 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 438d359177SBarry Smith Level: developer 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 518d359177SBarry Smith .seealso: TSPseudoTimeStepDefault(), 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__ 688d359177SBarry Smith #define __FUNCT__ "TSPseudoVerifyTimeStepDefault" 697bf11e45SBarry Smith /*@C 708d359177SBarry Smith TSPseudoVerifyTimeStepDefault - 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 @*/ 938d359177SBarry Smith PetscErrorCode TSPseudoVerifyTimeStepDefault(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 1248d359177SBarry Smith .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStepDefault() 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); 1649be3e283SDebojyoti Ghosh ierr = TSPostStage(ts,ts->ptime+ts->time_step,0,&(pseudo->update));CHKERRQ(ierr); 1655ef26d82SJed Brown ts->snes_its += its; ts->ksp_its += lits; 166cdbf8f93SLisandro Dalcin ierr = PetscInfo3(ts,"step=%D, nonlinear solve iterations=%D, linear solve iterations=%D\n",ts->steps,its,lits);CHKERRQ(ierr); 167193ac0bcSJed Brown pseudo->fnorm = -1; /* The current norm is no longer valid, monitor must recompute it. */ 168cdbf8f93SLisandro Dalcin ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&next_time_step,&stepok);CHKERRQ(ierr); 169cdbf8f93SLisandro Dalcin if (stepok) break; 170cdbf8f93SLisandro Dalcin } 171f1b97656SJed Brown if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) { 172193ac0bcSJed Brown ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 173193ac0bcSJed 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); 174cdbf8f93SLisandro Dalcin PetscFunctionReturn(0); 1757bf11e45SBarry Smith } 176cdbf8f93SLisandro Dalcin if (reject >= ts->max_reject) { 177193ac0bcSJed Brown ts->reason = TS_DIVERGED_STEP_REJECTED; 178cdbf8f93SLisandro Dalcin ierr = PetscInfo2(ts,"step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,reject);CHKERRQ(ierr); 179cdbf8f93SLisandro Dalcin PetscFunctionReturn(0); 180193ac0bcSJed Brown } 181cdbf8f93SLisandro Dalcin ierr = VecCopy(pseudo->update,ts->vec_sol);CHKERRQ(ierr); 182193ac0bcSJed Brown ts->ptime += ts->time_step; 183cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 1842d3f70b5SBarry Smith ts->steps++; 1853a40ed3dSBarry Smith PetscFunctionReturn(0); 1862d3f70b5SBarry Smith } 1872d3f70b5SBarry Smith 1882d3f70b5SBarry Smith /*------------------------------------------------------------*/ 1894a2ae208SSatish Balay #undef __FUNCT__ 190277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Pseudo" 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 #undef __FUNCT__ 204277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Pseudo" 205277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Pseudo(TS ts) 206277b19d0SLisandro Dalcin { 207277b19d0SLisandro Dalcin PetscErrorCode ierr; 208277b19d0SLisandro Dalcin 209277b19d0SLisandro Dalcin PetscFunctionBegin; 210277b19d0SLisandro Dalcin ierr = TSReset_Pseudo(ts);CHKERRQ(ierr); 211277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 212bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",NULL);CHKERRQ(ierr); 213bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",NULL);CHKERRQ(ierr); 214bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetMaxTimeStep_C",NULL);CHKERRQ(ierr); 215bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",NULL);CHKERRQ(ierr); 216bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",NULL);CHKERRQ(ierr); 217277b19d0SLisandro Dalcin PetscFunctionReturn(0); 218277b19d0SLisandro Dalcin } 2192d3f70b5SBarry Smith 2202d3f70b5SBarry Smith /*------------------------------------------------------------*/ 2212d3f70b5SBarry Smith 2224a2ae208SSatish Balay #undef __FUNCT__ 2236f2d6a7bSJed Brown #define __FUNCT__ "TSPseudoGetXdot" 2246f2d6a7bSJed Brown /* 2256f2d6a7bSJed Brown Compute Xdot = (X^{n+1}-X^n)/dt) = 0 2266f2d6a7bSJed Brown */ 2276f2d6a7bSJed Brown static PetscErrorCode TSPseudoGetXdot(TS ts,Vec X,Vec *Xdot) 2282d3f70b5SBarry Smith { 2296f2d6a7bSJed Brown TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 230193ac0bcSJed Brown const PetscScalar mdt = 1.0/ts->time_step,*xnp1,*xn; 231193ac0bcSJed Brown PetscScalar *xdot; 232dfbe8321SBarry Smith PetscErrorCode ierr; 233a7cc72afSBarry Smith PetscInt i,n; 2342d3f70b5SBarry Smith 2353a40ed3dSBarry Smith PetscFunctionBegin; 236193ac0bcSJed Brown ierr = VecGetArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr); 237193ac0bcSJed Brown ierr = VecGetArrayRead(X,&xnp1);CHKERRQ(ierr); 2386f2d6a7bSJed Brown ierr = VecGetArray(pseudo->xdot,&xdot);CHKERRQ(ierr); 2396f2d6a7bSJed Brown ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr); 240bbd56ea5SKarl Rupp for (i=0; i<n; i++) xdot[i] = mdt*(xnp1[i] - xn[i]); 241193ac0bcSJed Brown ierr = VecRestoreArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr); 242193ac0bcSJed Brown ierr = VecRestoreArrayRead(X,&xnp1);CHKERRQ(ierr); 2436f2d6a7bSJed Brown ierr = VecRestoreArray(pseudo->xdot,&xdot);CHKERRQ(ierr); 2446f2d6a7bSJed Brown *Xdot = pseudo->xdot; 2453a40ed3dSBarry Smith PetscFunctionReturn(0); 2462d3f70b5SBarry Smith } 2472d3f70b5SBarry Smith 2484a2ae208SSatish Balay #undef __FUNCT__ 2490f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Pseudo" 2506f2d6a7bSJed Brown /* 2516f2d6a7bSJed Brown The transient residual is 2526f2d6a7bSJed Brown 2536f2d6a7bSJed Brown F(U^{n+1},(U^{n+1}-U^n)/dt) = 0 2546f2d6a7bSJed Brown 2556f2d6a7bSJed Brown or for ODE, 2566f2d6a7bSJed Brown 2576f2d6a7bSJed Brown (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0 2586f2d6a7bSJed Brown 2596f2d6a7bSJed Brown This is the function that must be evaluated for transient simulation and for 2606f2d6a7bSJed Brown finite difference Jacobians. On the first Newton step, this algorithm uses 2616f2d6a7bSJed Brown a guess of U^{n+1} = U^n in which case the transient term vanishes and the 2626f2d6a7bSJed Brown residual is actually the steady state residual. Pseudotransient 2636f2d6a7bSJed Brown continuation as described in the literature is a linearly implicit 2646f2d6a7bSJed Brown algorithm, it only takes this one Newton step with the steady state 2656f2d6a7bSJed Brown residual, and then advances to the next time step. 2666f2d6a7bSJed Brown */ 2670f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes,Vec X,Vec Y,TS ts) 2682d3f70b5SBarry Smith { 2696f2d6a7bSJed Brown Vec Xdot; 270dfbe8321SBarry Smith PetscErrorCode ierr; 2712d3f70b5SBarry Smith 2723a40ed3dSBarry Smith PetscFunctionBegin; 2736f2d6a7bSJed Brown ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr); 274193ac0bcSJed Brown ierr = TSComputeIFunction(ts,ts->ptime+ts->time_step,X,Xdot,Y,PETSC_FALSE);CHKERRQ(ierr); 2756f2d6a7bSJed Brown PetscFunctionReturn(0); 2766f2d6a7bSJed Brown } 2772d3f70b5SBarry Smith 2786f2d6a7bSJed Brown #undef __FUNCT__ 2790f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Pseudo" 2806f2d6a7bSJed Brown /* 2816f2d6a7bSJed Brown This constructs the Jacobian needed for SNES. For DAE, this is 2826f2d6a7bSJed Brown 2836f2d6a7bSJed Brown dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot 2846f2d6a7bSJed Brown 2856f2d6a7bSJed Brown and for ODE: 2866f2d6a7bSJed Brown 2876f2d6a7bSJed Brown J = I/dt - J_{Frhs} where J_{Frhs} is the given Jacobian of Frhs. 2886f2d6a7bSJed Brown */ 2890f5c6efeSJed Brown static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes,Vec X,Mat *AA,Mat *BB,MatStructure *str,TS ts) 2906f2d6a7bSJed Brown { 2916f2d6a7bSJed Brown Vec Xdot; 2926f2d6a7bSJed Brown PetscErrorCode ierr; 2936f2d6a7bSJed Brown 2946f2d6a7bSJed Brown PetscFunctionBegin; 2956f2d6a7bSJed Brown ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr); 296193ac0bcSJed Brown ierr = TSComputeIJacobian(ts,ts->ptime+ts->time_step,X,Xdot,1./ts->time_step,AA,BB,str,PETSC_FALSE);CHKERRQ(ierr); 2973a40ed3dSBarry Smith PetscFunctionReturn(0); 2982d3f70b5SBarry Smith } 2992d3f70b5SBarry Smith 3002d3f70b5SBarry Smith 3014a2ae208SSatish Balay #undef __FUNCT__ 3024a2ae208SSatish Balay #define __FUNCT__ "TSSetUp_Pseudo" 3036849ba73SBarry Smith static PetscErrorCode TSSetUp_Pseudo(TS ts) 3042d3f70b5SBarry Smith { 3057bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 306dfbe8321SBarry Smith PetscErrorCode ierr; 3072d3f70b5SBarry Smith 3083a40ed3dSBarry Smith PetscFunctionBegin; 3097bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr); 3107bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr); 3116f2d6a7bSJed Brown ierr = VecDuplicate(ts->vec_sol,&pseudo->xdot);CHKERRQ(ierr); 3123a40ed3dSBarry Smith PetscFunctionReturn(0); 3132d3f70b5SBarry Smith } 3142d3f70b5SBarry Smith /*------------------------------------------------------------*/ 3152d3f70b5SBarry Smith 3164a2ae208SSatish Balay #undef __FUNCT__ 317a6570f20SBarry Smith #define __FUNCT__ "TSPseudoMonitorDefault" 318649052a6SBarry Smith PetscErrorCode TSPseudoMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy) 3192d3f70b5SBarry Smith { 3207bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 321dfbe8321SBarry Smith PetscErrorCode ierr; 322ce94432eSBarry Smith PetscViewer viewer = (PetscViewer) dummy; 3232d3f70b5SBarry Smith 3243a40ed3dSBarry Smith PetscFunctionBegin; 325ce94432eSBarry Smith if (!viewer) { 326ce94432eSBarry Smith ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ts),&viewer);CHKERRQ(ierr); 327ce94432eSBarry Smith } 328193ac0bcSJed Brown if (pseudo->fnorm < 0) { /* The last computed norm is stale, recompute */ 329193ac0bcSJed Brown ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr); 330193ac0bcSJed Brown ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr); 331193ac0bcSJed Brown ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 332193ac0bcSJed Brown } 333649052a6SBarry Smith ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 334*7c8652ddSBarry 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); 335649052a6SBarry Smith ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 3363a40ed3dSBarry Smith PetscFunctionReturn(0); 3372d3f70b5SBarry Smith } 3382d3f70b5SBarry Smith 3394a2ae208SSatish Balay #undef __FUNCT__ 3404a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Pseudo" 3416849ba73SBarry Smith static PetscErrorCode TSSetFromOptions_Pseudo(TS ts) 3422d3f70b5SBarry Smith { 3434bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 344dfbe8321SBarry Smith PetscErrorCode ierr; 345ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE; 346649052a6SBarry Smith PetscViewer viewer; 3472d3f70b5SBarry Smith 3483a40ed3dSBarry Smith PetscFunctionBegin; 349b0a32e0cSBarry Smith ierr = PetscOptionsHead("Pseudo-timestepping options");CHKERRQ(ierr); 3500298fd71SBarry Smith ierr = PetscOptionsBool("-ts_monitor_pseudo","Monitor convergence","TSPseudoMonitorDefault",flg,&flg,NULL);CHKERRQ(ierr); 3512d3f70b5SBarry Smith if (flg) { 352ce94432eSBarry Smith ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ts),"stdout",&viewer);CHKERRQ(ierr); 353649052a6SBarry Smith ierr = TSMonitorSet(ts,TSPseudoMonitorDefault,viewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 35428aa8177SBarry Smith } 35590d69ab7SBarry Smith flg = PETSC_FALSE; 3560298fd71SBarry Smith ierr = PetscOptionsBool("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",flg,&flg,NULL);CHKERRQ(ierr); 357ca4b7087SBarry Smith if (flg) { 358ca4b7087SBarry Smith ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr); 359ca4b7087SBarry Smith } 360419fbf4bSSatish Balay ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,0);CHKERRQ(ierr); 36186534af1SJed Brown ierr = PetscOptionsReal("-ts_pseudo_max_dt","Maximum value for dt","TSPseudoSetMaxTimeStep",pseudo->dt_max,&pseudo->dt_max,0);CHKERRQ(ierr); 362d52bd9f3SBarry Smith 363d52bd9f3SBarry Smith ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); 364b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 3653a40ed3dSBarry Smith PetscFunctionReturn(0); 3662d3f70b5SBarry Smith } 3672d3f70b5SBarry Smith 3684a2ae208SSatish Balay #undef __FUNCT__ 3694a2ae208SSatish Balay #define __FUNCT__ "TSView_Pseudo" 3706849ba73SBarry Smith static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer) 3712d3f70b5SBarry Smith { 372d52bd9f3SBarry Smith PetscErrorCode ierr; 373d52bd9f3SBarry Smith 3743a40ed3dSBarry Smith PetscFunctionBegin; 375d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 3763a40ed3dSBarry Smith PetscFunctionReturn(0); 3772d3f70b5SBarry Smith } 3782d3f70b5SBarry Smith 37982bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 3804a2ae208SSatish Balay #undef __FUNCT__ 3814a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep" 382ac226902SBarry Smith /*@C 38382bf6240SBarry Smith TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the 38482bf6240SBarry Smith last timestep. 38582bf6240SBarry Smith 3863f9fe445SBarry Smith Logically Collective on TS 38715091d37SBarry Smith 38882bf6240SBarry Smith Input Parameters: 38915091d37SBarry Smith + ts - timestep context 39082bf6240SBarry Smith . dt - user-defined function to verify timestep 39115091d37SBarry Smith - ctx - [optional] user-defined context for private data 3920298fd71SBarry Smith for the timestep verification routine (may be NULL) 39382bf6240SBarry Smith 39415091d37SBarry Smith Level: advanced 395fee21e36SBarry Smith 39682bf6240SBarry Smith Calling sequence of func: 397ace3abfcSBarry Smith . func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscBool *flag); 39882bf6240SBarry Smith 39982bf6240SBarry Smith . update - latest solution vector 40082bf6240SBarry Smith . ctx - [optional] timestep context 40182bf6240SBarry Smith . newdt - the timestep to use for the next step 40282bf6240SBarry Smith . flag - flag indicating whether the last time step was acceptable 40382bf6240SBarry Smith 40482bf6240SBarry Smith Notes: 40582bf6240SBarry Smith The routine set here will be called by TSPseudoVerifyTimeStep() 40682bf6240SBarry Smith during the timestepping process. 40782bf6240SBarry Smith 40882bf6240SBarry Smith .keywords: timestep, pseudo, set, verify 40982bf6240SBarry Smith 4108d359177SBarry Smith .seealso: TSPseudoVerifyTimeStepDefault(), TSPseudoVerifyTimeStep() 41182bf6240SBarry Smith @*/ 4127087cfbeSBarry Smith PetscErrorCode TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,PetscBool*),void *ctx) 41382bf6240SBarry Smith { 4144ac538c5SBarry Smith PetscErrorCode ierr; 41582bf6240SBarry Smith 41682bf6240SBarry Smith PetscFunctionBegin; 4170700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 4184ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSPseudoSetVerifyTimeStep_C",(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal*,PetscBool*),void*),(ts,dt,ctx));CHKERRQ(ierr); 41982bf6240SBarry Smith PetscFunctionReturn(0); 42082bf6240SBarry Smith } 42182bf6240SBarry Smith 4224a2ae208SSatish Balay #undef __FUNCT__ 4234a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement" 42482bf6240SBarry Smith /*@ 42582bf6240SBarry Smith TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 4268d359177SBarry Smith dt when using the TSPseudoTimeStepDefault() routine. 42782bf6240SBarry Smith 4283f9fe445SBarry Smith Logically Collective on TS 429fee21e36SBarry Smith 43015091d37SBarry Smith Input Parameters: 43115091d37SBarry Smith + ts - the timestep context 43215091d37SBarry Smith - inc - the scaling factor >= 1.0 43315091d37SBarry Smith 43482bf6240SBarry Smith Options Database Key: 43582bf6240SBarry Smith $ -ts_pseudo_increment <increment> 43682bf6240SBarry Smith 43715091d37SBarry Smith Level: advanced 43815091d37SBarry Smith 43982bf6240SBarry Smith .keywords: timestep, pseudo, set, increment 44082bf6240SBarry Smith 4418d359177SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault() 44282bf6240SBarry Smith @*/ 4437087cfbeSBarry Smith PetscErrorCode TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc) 44482bf6240SBarry Smith { 4454ac538c5SBarry Smith PetscErrorCode ierr; 44682bf6240SBarry Smith 44782bf6240SBarry Smith PetscFunctionBegin; 4480700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 449c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts,inc,2); 4504ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSPseudoSetTimeStepIncrement_C",(TS,PetscReal),(ts,inc));CHKERRQ(ierr); 45182bf6240SBarry Smith PetscFunctionReturn(0); 45282bf6240SBarry Smith } 45382bf6240SBarry Smith 4544a2ae208SSatish Balay #undef __FUNCT__ 45586534af1SJed Brown #define __FUNCT__ "TSPseudoSetMaxTimeStep" 45686534af1SJed Brown /*@ 45786534af1SJed Brown TSPseudoSetMaxTimeStep - Sets the maximum time step 4588d359177SBarry Smith when using the TSPseudoTimeStepDefault() routine. 45986534af1SJed Brown 46086534af1SJed Brown Logically Collective on TS 46186534af1SJed Brown 46286534af1SJed Brown Input Parameters: 46386534af1SJed Brown + ts - the timestep context 46486534af1SJed Brown - maxdt - the maximum time step, use a non-positive value to deactivate 46586534af1SJed Brown 46686534af1SJed Brown Options Database Key: 46786534af1SJed Brown $ -ts_pseudo_max_dt <increment> 46886534af1SJed Brown 46986534af1SJed Brown Level: advanced 47086534af1SJed Brown 47186534af1SJed Brown .keywords: timestep, pseudo, set 47286534af1SJed Brown 4738d359177SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault() 47486534af1SJed Brown @*/ 47586534af1SJed Brown PetscErrorCode TSPseudoSetMaxTimeStep(TS ts,PetscReal maxdt) 47686534af1SJed Brown { 47786534af1SJed Brown PetscErrorCode ierr; 47886534af1SJed Brown 47986534af1SJed Brown PetscFunctionBegin; 48086534af1SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 48186534af1SJed Brown PetscValidLogicalCollectiveReal(ts,maxdt,2); 48286534af1SJed Brown ierr = PetscTryMethod(ts,"TSPseudoSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));CHKERRQ(ierr); 48386534af1SJed Brown PetscFunctionReturn(0); 48486534af1SJed Brown } 48586534af1SJed Brown 48686534af1SJed Brown #undef __FUNCT__ 4874a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt" 48882bf6240SBarry Smith /*@ 48982bf6240SBarry Smith TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep 49082bf6240SBarry Smith is computed via the formula 49182bf6240SBarry Smith $ dt = initial_dt*initial_fnorm/current_fnorm 49282bf6240SBarry Smith rather than the default update, 49382bf6240SBarry Smith $ dt = current_dt*previous_fnorm/current_fnorm. 49482bf6240SBarry Smith 4953f9fe445SBarry Smith Logically Collective on TS 49615091d37SBarry Smith 49782bf6240SBarry Smith Input Parameter: 49882bf6240SBarry Smith . ts - the timestep context 49982bf6240SBarry Smith 50082bf6240SBarry Smith Options Database Key: 50182bf6240SBarry Smith $ -ts_pseudo_increment_dt_from_initial_dt 50282bf6240SBarry Smith 50315091d37SBarry Smith Level: advanced 50415091d37SBarry Smith 50582bf6240SBarry Smith .keywords: timestep, pseudo, set, increment 50682bf6240SBarry Smith 5078d359177SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoTimeStepDefault() 50882bf6240SBarry Smith @*/ 5097087cfbeSBarry Smith PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS ts) 51082bf6240SBarry Smith { 5114ac538c5SBarry Smith PetscErrorCode ierr; 51282bf6240SBarry Smith 51382bf6240SBarry Smith PetscFunctionBegin; 5140700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 5154ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSPseudoIncrementDtFromInitialDt_C",(TS),(ts));CHKERRQ(ierr); 51682bf6240SBarry Smith PetscFunctionReturn(0); 51782bf6240SBarry Smith } 51882bf6240SBarry Smith 51982bf6240SBarry Smith 5204a2ae208SSatish Balay #undef __FUNCT__ 5214a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep" 522ac226902SBarry Smith /*@C 52382bf6240SBarry Smith TSPseudoSetTimeStep - Sets the user-defined routine to be 52482bf6240SBarry Smith called at each pseudo-timestep to update the timestep. 52582bf6240SBarry Smith 5263f9fe445SBarry Smith Logically Collective on TS 52715091d37SBarry Smith 52882bf6240SBarry Smith Input Parameters: 52915091d37SBarry Smith + ts - timestep context 53082bf6240SBarry Smith . dt - function to compute timestep 53115091d37SBarry Smith - ctx - [optional] user-defined context for private data 5320298fd71SBarry Smith required by the function (may be NULL) 53382bf6240SBarry Smith 53415091d37SBarry Smith Level: intermediate 535fee21e36SBarry Smith 53682bf6240SBarry Smith Calling sequence of func: 53787828ca2SBarry Smith . func (TS ts,PetscReal *newdt,void *ctx); 53882bf6240SBarry Smith 53982bf6240SBarry Smith . newdt - the newly computed timestep 54082bf6240SBarry Smith . ctx - [optional] timestep context 54182bf6240SBarry Smith 54282bf6240SBarry Smith Notes: 54382bf6240SBarry Smith The routine set here will be called by TSPseudoComputeTimeStep() 54482bf6240SBarry Smith during the timestepping process. 5458d359177SBarry Smith If not set then TSPseudoTimeStepDefault() is automatically used 54682bf6240SBarry Smith 54782bf6240SBarry Smith .keywords: timestep, pseudo, set 54882bf6240SBarry Smith 5498d359177SBarry Smith .seealso: TSPseudoTimeStepDefault(), TSPseudoComputeTimeStep() 55082bf6240SBarry Smith @*/ 5517087cfbeSBarry Smith PetscErrorCode TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void *ctx) 55282bf6240SBarry Smith { 5534ac538c5SBarry Smith PetscErrorCode ierr; 55482bf6240SBarry Smith 55582bf6240SBarry Smith PetscFunctionBegin; 5560700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 5574ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSPseudoSetTimeStep_C",(TS,PetscErrorCode (*)(TS,PetscReal*,void*),void*),(ts,dt,ctx));CHKERRQ(ierr); 55882bf6240SBarry Smith PetscFunctionReturn(0); 55982bf6240SBarry Smith } 56082bf6240SBarry Smith 56182bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 56282bf6240SBarry Smith 563ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscBool*); /* force argument to next function to not be extern C*/ 5644a2ae208SSatish Balay #undef __FUNCT__ 5654a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo" 5667087cfbeSBarry Smith PetscErrorCode TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void *ctx) 56782bf6240SBarry Smith { 56882bf6240SBarry Smith TS_Pseudo *pseudo; 56982bf6240SBarry Smith 57082bf6240SBarry Smith PetscFunctionBegin; 57182bf6240SBarry Smith pseudo = (TS_Pseudo*)ts->data; 57282bf6240SBarry Smith pseudo->verify = dt; 57382bf6240SBarry Smith pseudo->verifyctx = ctx; 57482bf6240SBarry Smith PetscFunctionReturn(0); 57582bf6240SBarry Smith } 57682bf6240SBarry Smith 5774a2ae208SSatish Balay #undef __FUNCT__ 5784a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement_Pseudo" 5797087cfbeSBarry Smith PetscErrorCode TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc) 58082bf6240SBarry Smith { 5814bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 58282bf6240SBarry Smith 58382bf6240SBarry Smith PetscFunctionBegin; 58482bf6240SBarry Smith pseudo->dt_increment = inc; 58582bf6240SBarry Smith PetscFunctionReturn(0); 58682bf6240SBarry Smith } 58782bf6240SBarry Smith 5884a2ae208SSatish Balay #undef __FUNCT__ 58986534af1SJed Brown #define __FUNCT__ "TSPseudoSetMaxTimeStep_Pseudo" 59086534af1SJed Brown PetscErrorCode TSPseudoSetMaxTimeStep_Pseudo(TS ts,PetscReal maxdt) 59186534af1SJed Brown { 59286534af1SJed Brown TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 59386534af1SJed Brown 59486534af1SJed Brown PetscFunctionBegin; 59586534af1SJed Brown pseudo->dt_max = maxdt; 59686534af1SJed Brown PetscFunctionReturn(0); 59786534af1SJed Brown } 59886534af1SJed Brown 59986534af1SJed Brown #undef __FUNCT__ 6004a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt_Pseudo" 6017087cfbeSBarry Smith PetscErrorCode TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts) 60282bf6240SBarry Smith { 6034bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 60482bf6240SBarry Smith 60582bf6240SBarry Smith PetscFunctionBegin; 6064bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_TRUE; 60782bf6240SBarry Smith PetscFunctionReturn(0); 60882bf6240SBarry Smith } 60982bf6240SBarry Smith 6106849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/ 6114a2ae208SSatish Balay #undef __FUNCT__ 6124a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo" 6137087cfbeSBarry Smith PetscErrorCode TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void *ctx) 61482bf6240SBarry Smith { 6154bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 61682bf6240SBarry Smith 61782bf6240SBarry Smith PetscFunctionBegin; 61882bf6240SBarry Smith pseudo->dt = dt; 61982bf6240SBarry Smith pseudo->dtctx = ctx; 62082bf6240SBarry Smith PetscFunctionReturn(0); 62182bf6240SBarry Smith } 62282bf6240SBarry Smith 62382bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 62410e6a065SJed Brown /*MC 62510e6a065SJed Brown TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping 62682bf6240SBarry Smith 62710e6a065SJed Brown This method solves equations of the form 62810e6a065SJed Brown 62910e6a065SJed Brown $ F(X,Xdot) = 0 63010e6a065SJed Brown 63110e6a065SJed Brown for steady state using the iteration 63210e6a065SJed Brown 63310e6a065SJed Brown $ [G'] S = -F(X,0) 63410e6a065SJed Brown $ X += S 63510e6a065SJed Brown 63610e6a065SJed Brown where 63710e6a065SJed Brown 63810e6a065SJed Brown $ G(Y) = F(Y,(Y-X)/dt) 63910e6a065SJed Brown 6406f2d6a7bSJed Brown This is linearly-implicit Euler with the residual always evaluated "at steady 6416f2d6a7bSJed Brown state". See note below. 64210e6a065SJed Brown 64310e6a065SJed Brown Options database keys: 64410e6a065SJed Brown + -ts_pseudo_increment <real> - ratio of increase dt 64510e6a065SJed Brown - -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt 64610e6a065SJed Brown 64710e6a065SJed Brown Level: beginner 64810e6a065SJed Brown 64910e6a065SJed Brown References: 65010e6a065SJed Brown Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential-Algebraic Equations, 2003. 65110e6a065SJed Brown C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998. 65210e6a065SJed Brown 65310e6a065SJed Brown Notes: 6546f2d6a7bSJed Brown The residual computed by this method includes the transient term (Xdot is computed instead of 6556f2d6a7bSJed Brown always being zero), but since the prediction from the last step is always the solution from the 6566f2d6a7bSJed Brown last step, on the first Newton iteration we have 6576f2d6a7bSJed Brown 6586f2d6a7bSJed Brown $ Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0 6596f2d6a7bSJed Brown 6606f2d6a7bSJed Brown Therefore, the linear system solved by the first Newton iteration is equivalent to the one 6616f2d6a7bSJed Brown described above and in the papers. If the user chooses to perform multiple Newton iterations, the 6626f2d6a7bSJed Brown algorithm is no longer the one described in the referenced papers. 66310e6a065SJed Brown 66410e6a065SJed Brown .seealso: TSCreate(), TS, TSSetType() 66510e6a065SJed Brown 66610e6a065SJed Brown M*/ 6674a2ae208SSatish Balay #undef __FUNCT__ 6684a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Pseudo" 6698cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Pseudo(TS ts) 6702d3f70b5SBarry Smith { 6717bf11e45SBarry Smith TS_Pseudo *pseudo; 672dfbe8321SBarry Smith PetscErrorCode ierr; 673193ac0bcSJed Brown SNES snes; 67419fd82e9SBarry Smith SNESType stype; 6752d3f70b5SBarry Smith 6763a40ed3dSBarry Smith PetscFunctionBegin; 677277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Pseudo; 678000e7ae3SMatthew Knepley ts->ops->destroy = TSDestroy_Pseudo; 679000e7ae3SMatthew Knepley ts->ops->view = TSView_Pseudo; 6802d3f70b5SBarry Smith 681000e7ae3SMatthew Knepley ts->ops->setup = TSSetUp_Pseudo; 682000e7ae3SMatthew Knepley ts->ops->step = TSStep_Pseudo; 683000e7ae3SMatthew Knepley ts->ops->setfromoptions = TSSetFromOptions_Pseudo; 6840f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Pseudo; 6850f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Pseudo; 6867bf11e45SBarry Smith 687193ac0bcSJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 688193ac0bcSJed Brown ierr = SNESGetType(snes,&stype);CHKERRQ(ierr); 689193ac0bcSJed Brown if (!stype) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);} 6902d3f70b5SBarry Smith 691b00a9115SJed Brown ierr = PetscNewLog(ts,&pseudo);CHKERRQ(ierr); 6927bf11e45SBarry Smith ts->data = (void*)pseudo; 6932d3f70b5SBarry Smith 69428aa8177SBarry Smith pseudo->dt_increment = 1.1; 6954bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_FALSE; 6968d359177SBarry Smith pseudo->dt = TSPseudoTimeStepDefault; 697193ac0bcSJed Brown pseudo->fnorm = -1; 69882bf6240SBarry Smith 699bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr); 700bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr); 701bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetMaxTimeStep_C",TSPseudoSetMaxTimeStep_Pseudo);CHKERRQ(ierr); 702bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr); 703bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr); 7043a40ed3dSBarry Smith PetscFunctionReturn(0); 7052d3f70b5SBarry Smith } 7062d3f70b5SBarry Smith 7074a2ae208SSatish Balay #undef __FUNCT__ 7088d359177SBarry Smith #define __FUNCT__ "TSPseudoTimeStepDefault" 70982bf6240SBarry Smith /*@C 7108d359177SBarry Smith TSPseudoTimeStepDefault - Default code to compute pseudo-timestepping. 71182bf6240SBarry Smith Use with TSPseudoSetTimeStep(). 71228aa8177SBarry Smith 71315091d37SBarry Smith Collective on TS 71415091d37SBarry Smith 71528aa8177SBarry Smith Input Parameters: 71628aa8177SBarry Smith . ts - the timestep context 71782bf6240SBarry Smith . dtctx - unused timestep context 71828aa8177SBarry Smith 71982bf6240SBarry Smith Output Parameter: 72082bf6240SBarry Smith . newdt - the timestep to use for the next step 72128aa8177SBarry Smith 72215091d37SBarry Smith Level: advanced 72315091d37SBarry Smith 72482bf6240SBarry Smith .keywords: timestep, pseudo, default 725564e8f4eSLois Curfman McInnes 72682bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep() 72728aa8177SBarry Smith @*/ 7288d359177SBarry Smith PetscErrorCode TSPseudoTimeStepDefault(TS ts,PetscReal *newdt,void *dtctx) 72928aa8177SBarry Smith { 73082bf6240SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 73187828ca2SBarry Smith PetscReal inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous; 732dfbe8321SBarry Smith PetscErrorCode ierr; 73328aa8177SBarry Smith 7343a40ed3dSBarry Smith PetscFunctionBegin; 735bbd7b040SJed Brown ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr); 736214bc6a2SJed Brown ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr); 73782bf6240SBarry Smith ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 738cdbf8f93SLisandro Dalcin if (pseudo->fnorm_initial == 0.0) { 73982bf6240SBarry Smith /* first time through so compute initial function norm */ 740cdbf8f93SLisandro Dalcin pseudo->fnorm_initial = pseudo->fnorm; 74182bf6240SBarry Smith fnorm_previous = pseudo->fnorm; 74282bf6240SBarry Smith } 743bbd56ea5SKarl Rupp if (pseudo->fnorm == 0.0) *newdt = 1.e12*inc*ts->time_step; 744bbd56ea5SKarl Rupp else if (pseudo->increment_dt_from_initial_dt) *newdt = inc*pseudo->dt_initial*pseudo->fnorm_initial/pseudo->fnorm; 745bbd56ea5SKarl Rupp else *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm; 74686534af1SJed Brown if (pseudo->dt_max > 0) *newdt = PetscMin(*newdt,pseudo->dt_max); 74782bf6240SBarry Smith pseudo->fnorm_previous = pseudo->fnorm; 7483a40ed3dSBarry Smith PetscFunctionReturn(0); 74928aa8177SBarry Smith } 750