12d3f70b5SBarry Smith /* 2fb4a63b6SLois Curfman McInnes Code for Timestepping with implicit backwards Euler. 32d3f70b5SBarry Smith */ 4c6db04a5SJed Brown #include <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 18*cdbf8f93SLisandro Dalcin PetscReal fnorm_initial,fnorm; /* original and current norm of F(u) */ 1987828ca2SBarry Smith PetscReal fnorm_previous; 2028aa8177SBarry Smith 21*cdbf8f93SLisandro Dalcin PetscReal dt_initial; /* initial time-step */ 2287828ca2SBarry Smith PetscReal dt_increment; /* scaling that dt is incremented each time-step */ 23ace3abfcSBarry Smith PetscBool increment_dt_from_initial_dt; 247bf11e45SBarry Smith } TS_Pseudo; 252d3f70b5SBarry Smith 262d3f70b5SBarry Smith /* ------------------------------------------------------------------------------*/ 272d3f70b5SBarry Smith 284a2ae208SSatish Balay #undef __FUNCT__ 294a2ae208SSatish Balay #define __FUNCT__ "TSPseudoComputeTimeStep" 307bf11e45SBarry Smith /*@ 317bf11e45SBarry Smith TSPseudoComputeTimeStep - Computes the next timestep for a currently running 32564e8f4eSLois Curfman McInnes pseudo-timestepping process. 332d3f70b5SBarry Smith 3415091d37SBarry Smith Collective on TS 3515091d37SBarry Smith 367bf11e45SBarry Smith Input Parameter: 377bf11e45SBarry Smith . ts - timestep context 387bf11e45SBarry Smith 397bf11e45SBarry Smith Output Parameter: 40fb4a63b6SLois Curfman McInnes . dt - newly computed timestep 41fb4a63b6SLois Curfman McInnes 4215091d37SBarry Smith Level: advanced 43564e8f4eSLois Curfman McInnes 44564e8f4eSLois Curfman McInnes Notes: 45564e8f4eSLois Curfman McInnes The routine to be called here to compute the timestep should be 46564e8f4eSLois Curfman McInnes set by calling TSPseudoSetTimeStep(). 47564e8f4eSLois Curfman McInnes 48fb4a63b6SLois Curfman McInnes .keywords: timestep, pseudo, compute 49564e8f4eSLois Curfman McInnes 50564e8f4eSLois Curfman McInnes .seealso: TSPseudoDefaultTimeStep(), TSPseudoSetTimeStep() 517bf11e45SBarry Smith @*/ 527087cfbeSBarry Smith PetscErrorCode TSPseudoComputeTimeStep(TS ts,PetscReal *dt) 537bf11e45SBarry Smith { 547bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 55dfbe8321SBarry Smith PetscErrorCode ierr; 567bf11e45SBarry Smith 573a40ed3dSBarry Smith PetscFunctionBegin; 58d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr); 597bf11e45SBarry Smith ierr = (*pseudo->dt)(ts,dt,pseudo->dtctx);CHKERRQ(ierr); 60d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr); 613a40ed3dSBarry Smith PetscFunctionReturn(0); 627bf11e45SBarry Smith } 637bf11e45SBarry Smith 647bf11e45SBarry Smith 657bf11e45SBarry Smith /* ------------------------------------------------------------------------------*/ 664a2ae208SSatish Balay #undef __FUNCT__ 674a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultVerifyTimeStep" 687bf11e45SBarry Smith /*@C 69639f9d9dSBarry Smith TSPseudoDefaultVerifyTimeStep - Default code to verify the quality of the last timestep. 707bf11e45SBarry Smith 7115091d37SBarry Smith Collective on TS 7215091d37SBarry Smith 737bf11e45SBarry Smith Input Parameters: 7415091d37SBarry Smith + ts - the timestep context 757bf11e45SBarry Smith . dtctx - unused timestep context 7615091d37SBarry Smith - update - latest solution vector 777bf11e45SBarry Smith 78564e8f4eSLois Curfman McInnes Output Parameters: 7915091d37SBarry Smith + newdt - the timestep to use for the next step 8015091d37SBarry Smith - flag - flag indicating whether the last time step was acceptable 817bf11e45SBarry Smith 8215091d37SBarry Smith Level: advanced 83fee21e36SBarry Smith 84564e8f4eSLois Curfman McInnes Note: 85564e8f4eSLois Curfman McInnes This routine always returns a flag of 1, indicating an acceptable 86564e8f4eSLois Curfman McInnes timestep. 87564e8f4eSLois Curfman McInnes 88564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, default, verify 89564e8f4eSLois Curfman McInnes 90564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStep() 917bf11e45SBarry Smith @*/ 927087cfbeSBarry Smith PetscErrorCode TSPseudoDefaultVerifyTimeStep(TS ts,Vec update,void *dtctx,PetscReal *newdt,PetscBool *flag) 937bf11e45SBarry Smith { 943a40ed3dSBarry Smith PetscFunctionBegin; 95a7cc72afSBarry Smith *flag = PETSC_TRUE; 963a40ed3dSBarry Smith PetscFunctionReturn(0); 977bf11e45SBarry Smith } 987bf11e45SBarry Smith 997bf11e45SBarry Smith 1004a2ae208SSatish Balay #undef __FUNCT__ 1014a2ae208SSatish Balay #define __FUNCT__ "TSPseudoVerifyTimeStep" 1027bf11e45SBarry Smith /*@ 103564e8f4eSLois Curfman McInnes TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable. 1047bf11e45SBarry Smith 10515091d37SBarry Smith Collective on TS 10615091d37SBarry Smith 107fb4a63b6SLois Curfman McInnes Input Parameters: 10815091d37SBarry Smith + ts - timestep context 10915091d37SBarry Smith - update - latest solution vector 1107bf11e45SBarry Smith 111fb4a63b6SLois Curfman McInnes Output Parameters: 11215091d37SBarry Smith + dt - newly computed timestep (if it had to shrink) 11315091d37SBarry Smith - flag - indicates if current timestep was ok 1147bf11e45SBarry Smith 11515091d37SBarry Smith Level: advanced 116fee21e36SBarry Smith 117564e8f4eSLois Curfman McInnes Notes: 118564e8f4eSLois Curfman McInnes The routine to be called here to compute the timestep should be 119564e8f4eSLois Curfman McInnes set by calling TSPseudoSetVerifyTimeStep(). 120564e8f4eSLois Curfman McInnes 121564e8f4eSLois Curfman McInnes .keywords: timestep, pseudo, verify 122564e8f4eSLois Curfman McInnes 123564e8f4eSLois Curfman McInnes .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoDefaultVerifyTimeStep() 1247bf11e45SBarry Smith @*/ 1257087cfbeSBarry Smith PetscErrorCode TSPseudoVerifyTimeStep(TS ts,Vec update,PetscReal *dt,PetscBool *flag) 1267bf11e45SBarry Smith { 1277bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 128dfbe8321SBarry Smith PetscErrorCode ierr; 1297bf11e45SBarry Smith 1303a40ed3dSBarry Smith PetscFunctionBegin; 131a7cc72afSBarry Smith if (!pseudo->verify) {*flag = PETSC_TRUE; PetscFunctionReturn(0);} 1327bf11e45SBarry Smith 1337bf11e45SBarry Smith ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag);CHKERRQ(ierr); 1347bf11e45SBarry Smith 1353a40ed3dSBarry Smith PetscFunctionReturn(0); 1367bf11e45SBarry Smith } 1377bf11e45SBarry Smith 1387bf11e45SBarry Smith /* --------------------------------------------------------------------------------*/ 1397bf11e45SBarry Smith 1404a2ae208SSatish Balay #undef __FUNCT__ 1414a2ae208SSatish Balay #define __FUNCT__ "TSStep_Pseudo" 142193ac0bcSJed Brown static PetscErrorCode TSStep_Pseudo(TS ts) 1432d3f70b5SBarry Smith { 144277b19d0SLisandro Dalcin TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 145*cdbf8f93SLisandro Dalcin PetscInt its,lits,reject; 146*cdbf8f93SLisandro Dalcin PetscBool stepok; 147*cdbf8f93SLisandro Dalcin PetscReal next_time_step; 148*cdbf8f93SLisandro Dalcin SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING; 149dfbe8321SBarry Smith PetscErrorCode ierr; 1502d3f70b5SBarry Smith 1513a40ed3dSBarry Smith PetscFunctionBegin; 152*cdbf8f93SLisandro Dalcin if (ts->steps == 0) { 153*cdbf8f93SLisandro Dalcin pseudo->dt_initial = ts->time_step; 154*cdbf8f93SLisandro Dalcin } 155193ac0bcSJed Brown ierr = VecCopy(ts->vec_sol,pseudo->update);CHKERRQ(ierr); 156*cdbf8f93SLisandro Dalcin next_time_step = ts->time_step; 157*cdbf8f93SLisandro Dalcin ierr = TSPseudoComputeTimeStep(ts,&next_time_step);CHKERRQ(ierr); 158*cdbf8f93SLisandro Dalcin for (reject=0; reject<ts->max_reject; reject++,ts->reject++) { 159*cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 160193ac0bcSJed Brown ierr = SNESSolve(ts->snes,PETSC_NULL,pseudo->update);CHKERRQ(ierr); 161*cdbf8f93SLisandro 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); 164c9780b6fSBarry Smith ts->nonlinear_its += its; ts->linear_its += lits; 165*cdbf8f93SLisandro 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. */ 167*cdbf8f93SLisandro Dalcin ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&next_time_step,&stepok);CHKERRQ(ierr); 168*cdbf8f93SLisandro Dalcin if (stepok) break; 169*cdbf8f93SLisandro Dalcin } 170*cdbf8f93SLisandro Dalcin if (snesreason < 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); 173*cdbf8f93SLisandro Dalcin PetscFunctionReturn(0); 1747bf11e45SBarry Smith } 175*cdbf8f93SLisandro Dalcin if (reject >= ts->max_reject) { 176193ac0bcSJed Brown ts->reason = TS_DIVERGED_STEP_REJECTED; 177*cdbf8f93SLisandro Dalcin ierr = PetscInfo2(ts,"step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,reject);CHKERRQ(ierr); 178*cdbf8f93SLisandro Dalcin PetscFunctionReturn(0); 179193ac0bcSJed Brown } 180*cdbf8f93SLisandro Dalcin ierr = VecCopy(pseudo->update,ts->vec_sol);CHKERRQ(ierr); 181193ac0bcSJed Brown ts->ptime += ts->time_step; 182*cdbf8f93SLisandro Dalcin ts->time_step_prev = ts->time_step; 183*cdbf8f93SLisandro 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); 212335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C","",PETSC_NULL);CHKERRQ(ierr); 213335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C","",PETSC_NULL);CHKERRQ(ierr); 214335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C","",PETSC_NULL);CHKERRQ(ierr); 215335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C","",PETSC_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); 2392d3f70b5SBarry Smith for (i=0; i<n; i++) { 2406f2d6a7bSJed Brown xdot[i] = mdt*(xnp1[i] - xn[i]); 2412d3f70b5SBarry Smith } 242193ac0bcSJed Brown ierr = VecRestoreArrayRead(ts->vec_sol,&xn);CHKERRQ(ierr); 243193ac0bcSJed Brown ierr = VecRestoreArrayRead(X,&xnp1);CHKERRQ(ierr); 2446f2d6a7bSJed Brown ierr = VecRestoreArray(pseudo->xdot,&xdot);CHKERRQ(ierr); 2456f2d6a7bSJed Brown *Xdot = pseudo->xdot; 2463a40ed3dSBarry Smith PetscFunctionReturn(0); 2472d3f70b5SBarry Smith } 2482d3f70b5SBarry Smith 2494a2ae208SSatish Balay #undef __FUNCT__ 2500f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Pseudo" 2516f2d6a7bSJed Brown /* 2526f2d6a7bSJed Brown The transient residual is 2536f2d6a7bSJed Brown 2546f2d6a7bSJed Brown F(U^{n+1},(U^{n+1}-U^n)/dt) = 0 2556f2d6a7bSJed Brown 2566f2d6a7bSJed Brown or for ODE, 2576f2d6a7bSJed Brown 2586f2d6a7bSJed Brown (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0 2596f2d6a7bSJed Brown 2606f2d6a7bSJed Brown This is the function that must be evaluated for transient simulation and for 2616f2d6a7bSJed Brown finite difference Jacobians. On the first Newton step, this algorithm uses 2626f2d6a7bSJed Brown a guess of U^{n+1} = U^n in which case the transient term vanishes and the 2636f2d6a7bSJed Brown residual is actually the steady state residual. Pseudotransient 2646f2d6a7bSJed Brown continuation as described in the literature is a linearly implicit 2656f2d6a7bSJed Brown algorithm, it only takes this one Newton step with the steady state 2666f2d6a7bSJed Brown residual, and then advances to the next time step. 2676f2d6a7bSJed Brown */ 2680f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes,Vec X,Vec Y,TS ts) 2692d3f70b5SBarry Smith { 2706f2d6a7bSJed Brown Vec Xdot; 271dfbe8321SBarry Smith PetscErrorCode ierr; 2722d3f70b5SBarry Smith 2733a40ed3dSBarry Smith PetscFunctionBegin; 2746f2d6a7bSJed Brown ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr); 275193ac0bcSJed Brown ierr = TSComputeIFunction(ts,ts->ptime+ts->time_step,X,Xdot,Y,PETSC_FALSE);CHKERRQ(ierr); 2766f2d6a7bSJed Brown PetscFunctionReturn(0); 2776f2d6a7bSJed Brown } 2782d3f70b5SBarry Smith 2796f2d6a7bSJed Brown #undef __FUNCT__ 2800f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Pseudo" 2816f2d6a7bSJed Brown /* 2826f2d6a7bSJed Brown This constructs the Jacobian needed for SNES. For DAE, this is 2836f2d6a7bSJed Brown 2846f2d6a7bSJed Brown dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot 2856f2d6a7bSJed Brown 2866f2d6a7bSJed Brown and for ODE: 2876f2d6a7bSJed Brown 2886f2d6a7bSJed Brown J = I/dt - J_{Frhs} where J_{Frhs} is the given Jacobian of Frhs. 2896f2d6a7bSJed Brown */ 2900f5c6efeSJed Brown static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes,Vec X,Mat *AA,Mat *BB,MatStructure *str,TS ts) 2916f2d6a7bSJed Brown { 2926f2d6a7bSJed Brown Vec Xdot; 2936f2d6a7bSJed Brown PetscErrorCode ierr; 2946f2d6a7bSJed Brown 2956f2d6a7bSJed Brown PetscFunctionBegin; 2966f2d6a7bSJed Brown ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr); 297193ac0bcSJed Brown ierr = TSComputeIJacobian(ts,ts->ptime+ts->time_step,X,Xdot,1./ts->time_step,AA,BB,str,PETSC_FALSE);CHKERRQ(ierr); 2983a40ed3dSBarry Smith PetscFunctionReturn(0); 2992d3f70b5SBarry Smith } 3002d3f70b5SBarry Smith 3012d3f70b5SBarry Smith 3024a2ae208SSatish Balay #undef __FUNCT__ 3034a2ae208SSatish Balay #define __FUNCT__ "TSSetUp_Pseudo" 3046849ba73SBarry Smith static PetscErrorCode TSSetUp_Pseudo(TS ts) 3052d3f70b5SBarry Smith { 3067bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 307dfbe8321SBarry Smith PetscErrorCode ierr; 3082d3f70b5SBarry Smith 3093a40ed3dSBarry Smith PetscFunctionBegin; 3107bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr); 3117bf11e45SBarry Smith ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr); 3126f2d6a7bSJed Brown ierr = VecDuplicate(ts->vec_sol,&pseudo->xdot);CHKERRQ(ierr); 3133a40ed3dSBarry Smith PetscFunctionReturn(0); 3142d3f70b5SBarry Smith } 3152d3f70b5SBarry Smith /*------------------------------------------------------------*/ 3162d3f70b5SBarry Smith 3174a2ae208SSatish Balay #undef __FUNCT__ 318a6570f20SBarry Smith #define __FUNCT__ "TSPseudoMonitorDefault" 319649052a6SBarry Smith PetscErrorCode TSPseudoMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy) 3202d3f70b5SBarry Smith { 3217bf11e45SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 322dfbe8321SBarry Smith PetscErrorCode ierr; 323649052a6SBarry Smith PetscViewer viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ts)->comm); 3242d3f70b5SBarry Smith 3253a40ed3dSBarry Smith PetscFunctionBegin; 326193ac0bcSJed Brown if (pseudo->fnorm < 0) { /* The last computed norm is stale, recompute */ 327193ac0bcSJed Brown ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr); 328193ac0bcSJed Brown ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr); 329193ac0bcSJed Brown ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 330193ac0bcSJed Brown } 331649052a6SBarry Smith ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 332649052a6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"TS %D dt %G time %G fnorm %G\n",step,ts->time_step,ptime,pseudo->fnorm);CHKERRQ(ierr); 333649052a6SBarry Smith ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 3343a40ed3dSBarry Smith PetscFunctionReturn(0); 3352d3f70b5SBarry Smith } 3362d3f70b5SBarry Smith 3374a2ae208SSatish Balay #undef __FUNCT__ 3384a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Pseudo" 3396849ba73SBarry Smith static PetscErrorCode TSSetFromOptions_Pseudo(TS ts) 3402d3f70b5SBarry Smith { 3414bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 342dfbe8321SBarry Smith PetscErrorCode ierr; 343ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE; 344649052a6SBarry Smith PetscViewer viewer; 3452d3f70b5SBarry Smith 3463a40ed3dSBarry Smith PetscFunctionBegin; 347b0a32e0cSBarry Smith ierr = PetscOptionsHead("Pseudo-timestepping options");CHKERRQ(ierr); 348193ac0bcSJed Brown ierr = PetscOptionsBool("-ts_monitor_pseudo","Monitor convergence","TSPseudoMonitorDefault",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 3492d3f70b5SBarry Smith if (flg) { 350649052a6SBarry Smith ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,"stdout",&viewer);CHKERRQ(ierr); 351649052a6SBarry Smith ierr = TSMonitorSet(ts,TSPseudoMonitorDefault,viewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 35228aa8177SBarry Smith } 35390d69ab7SBarry Smith flg = PETSC_FALSE; 354acfcf0e5SJed Brown ierr = PetscOptionsBool("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 355ca4b7087SBarry Smith if (flg) { 356ca4b7087SBarry Smith ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr); 357ca4b7087SBarry Smith } 358419fbf4bSSatish Balay ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,0);CHKERRQ(ierr); 359d52bd9f3SBarry Smith 360d52bd9f3SBarry Smith ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); 361b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 3623a40ed3dSBarry Smith PetscFunctionReturn(0); 3632d3f70b5SBarry Smith } 3642d3f70b5SBarry Smith 3654a2ae208SSatish Balay #undef __FUNCT__ 3664a2ae208SSatish Balay #define __FUNCT__ "TSView_Pseudo" 3676849ba73SBarry Smith static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer) 3682d3f70b5SBarry Smith { 369d52bd9f3SBarry Smith PetscErrorCode ierr; 370d52bd9f3SBarry Smith 3713a40ed3dSBarry Smith PetscFunctionBegin; 372d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 3733a40ed3dSBarry Smith PetscFunctionReturn(0); 3742d3f70b5SBarry Smith } 3752d3f70b5SBarry Smith 37682bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 3774a2ae208SSatish Balay #undef __FUNCT__ 3784a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep" 379ac226902SBarry Smith /*@C 38082bf6240SBarry Smith TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the 38182bf6240SBarry Smith last timestep. 38282bf6240SBarry Smith 3833f9fe445SBarry Smith Logically Collective on TS 38415091d37SBarry Smith 38582bf6240SBarry Smith Input Parameters: 38615091d37SBarry Smith + ts - timestep context 38782bf6240SBarry Smith . dt - user-defined function to verify timestep 38815091d37SBarry Smith - ctx - [optional] user-defined context for private data 38982bf6240SBarry Smith for the timestep verification routine (may be PETSC_NULL) 39082bf6240SBarry Smith 39115091d37SBarry Smith Level: advanced 392fee21e36SBarry Smith 39382bf6240SBarry Smith Calling sequence of func: 394ace3abfcSBarry Smith . func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscBool *flag); 39582bf6240SBarry Smith 39682bf6240SBarry Smith . update - latest solution vector 39782bf6240SBarry Smith . ctx - [optional] timestep context 39882bf6240SBarry Smith . newdt - the timestep to use for the next step 39982bf6240SBarry Smith . flag - flag indicating whether the last time step was acceptable 40082bf6240SBarry Smith 40182bf6240SBarry Smith Notes: 40282bf6240SBarry Smith The routine set here will be called by TSPseudoVerifyTimeStep() 40382bf6240SBarry Smith during the timestepping process. 40482bf6240SBarry Smith 40582bf6240SBarry Smith .keywords: timestep, pseudo, set, verify 40682bf6240SBarry Smith 40782bf6240SBarry Smith .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep() 40882bf6240SBarry Smith @*/ 4097087cfbeSBarry Smith PetscErrorCode TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,PetscBool *),void* ctx) 41082bf6240SBarry Smith { 4114ac538c5SBarry Smith PetscErrorCode ierr; 41282bf6240SBarry Smith 41382bf6240SBarry Smith PetscFunctionBegin; 4140700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 4154ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSPseudoSetVerifyTimeStep_C",(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal *,PetscBool *),void *),(ts,dt,ctx));CHKERRQ(ierr); 41682bf6240SBarry Smith PetscFunctionReturn(0); 41782bf6240SBarry Smith } 41882bf6240SBarry Smith 4194a2ae208SSatish Balay #undef __FUNCT__ 4204a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement" 42182bf6240SBarry Smith /*@ 42282bf6240SBarry Smith TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 42382bf6240SBarry Smith dt when using the TSPseudoDefaultTimeStep() routine. 42482bf6240SBarry Smith 4253f9fe445SBarry Smith Logically Collective on TS 426fee21e36SBarry Smith 42715091d37SBarry Smith Input Parameters: 42815091d37SBarry Smith + ts - the timestep context 42915091d37SBarry Smith - inc - the scaling factor >= 1.0 43015091d37SBarry Smith 43182bf6240SBarry Smith Options Database Key: 43282bf6240SBarry Smith $ -ts_pseudo_increment <increment> 43382bf6240SBarry Smith 43415091d37SBarry Smith Level: advanced 43515091d37SBarry Smith 43682bf6240SBarry Smith .keywords: timestep, pseudo, set, increment 43782bf6240SBarry Smith 43882bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 43982bf6240SBarry Smith @*/ 4407087cfbeSBarry Smith PetscErrorCode TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc) 44182bf6240SBarry Smith { 4424ac538c5SBarry Smith PetscErrorCode ierr; 44382bf6240SBarry Smith 44482bf6240SBarry Smith PetscFunctionBegin; 4450700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 446c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts,inc,2); 4474ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSPseudoSetTimeStepIncrement_C",(TS,PetscReal),(ts,inc));CHKERRQ(ierr); 44882bf6240SBarry Smith PetscFunctionReturn(0); 44982bf6240SBarry Smith } 45082bf6240SBarry Smith 4514a2ae208SSatish Balay #undef __FUNCT__ 4524a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt" 45382bf6240SBarry Smith /*@ 45482bf6240SBarry Smith TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep 45582bf6240SBarry Smith is computed via the formula 45682bf6240SBarry Smith $ dt = initial_dt*initial_fnorm/current_fnorm 45782bf6240SBarry Smith rather than the default update, 45882bf6240SBarry Smith $ dt = current_dt*previous_fnorm/current_fnorm. 45982bf6240SBarry Smith 4603f9fe445SBarry Smith Logically Collective on TS 46115091d37SBarry Smith 46282bf6240SBarry Smith Input Parameter: 46382bf6240SBarry Smith . ts - the timestep context 46482bf6240SBarry Smith 46582bf6240SBarry Smith Options Database Key: 46682bf6240SBarry Smith $ -ts_pseudo_increment_dt_from_initial_dt 46782bf6240SBarry Smith 46815091d37SBarry Smith Level: advanced 46915091d37SBarry Smith 47082bf6240SBarry Smith .keywords: timestep, pseudo, set, increment 47182bf6240SBarry Smith 47282bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 47382bf6240SBarry Smith @*/ 4747087cfbeSBarry Smith PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS ts) 47582bf6240SBarry Smith { 4764ac538c5SBarry Smith PetscErrorCode ierr; 47782bf6240SBarry Smith 47882bf6240SBarry Smith PetscFunctionBegin; 4790700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 4804ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSPseudoIncrementDtFromInitialDt_C",(TS),(ts));CHKERRQ(ierr); 48182bf6240SBarry Smith PetscFunctionReturn(0); 48282bf6240SBarry Smith } 48382bf6240SBarry Smith 48482bf6240SBarry Smith 4854a2ae208SSatish Balay #undef __FUNCT__ 4864a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep" 487ac226902SBarry Smith /*@C 48882bf6240SBarry Smith TSPseudoSetTimeStep - Sets the user-defined routine to be 48982bf6240SBarry Smith called at each pseudo-timestep to update the timestep. 49082bf6240SBarry Smith 4913f9fe445SBarry Smith Logically Collective on TS 49215091d37SBarry Smith 49382bf6240SBarry Smith Input Parameters: 49415091d37SBarry Smith + ts - timestep context 49582bf6240SBarry Smith . dt - function to compute timestep 49615091d37SBarry Smith - ctx - [optional] user-defined context for private data 49782bf6240SBarry Smith required by the function (may be PETSC_NULL) 49882bf6240SBarry Smith 49915091d37SBarry Smith Level: intermediate 500fee21e36SBarry Smith 50182bf6240SBarry Smith Calling sequence of func: 50287828ca2SBarry Smith . func (TS ts,PetscReal *newdt,void *ctx); 50382bf6240SBarry Smith 50482bf6240SBarry Smith . newdt - the newly computed timestep 50582bf6240SBarry Smith . ctx - [optional] timestep context 50682bf6240SBarry Smith 50782bf6240SBarry Smith Notes: 50882bf6240SBarry Smith The routine set here will be called by TSPseudoComputeTimeStep() 50982bf6240SBarry Smith during the timestepping process. 51082bf6240SBarry Smith 51182bf6240SBarry Smith .keywords: timestep, pseudo, set 51282bf6240SBarry Smith 51382bf6240SBarry Smith .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep() 51482bf6240SBarry Smith @*/ 5157087cfbeSBarry Smith PetscErrorCode TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void* ctx) 51682bf6240SBarry Smith { 5174ac538c5SBarry Smith PetscErrorCode ierr; 51882bf6240SBarry Smith 51982bf6240SBarry Smith PetscFunctionBegin; 5200700a824SBarry Smith PetscValidHeaderSpecific(ts,TS_CLASSID,1); 5214ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSPseudoSetTimeStep_C",(TS,PetscErrorCode (*)(TS,PetscReal *,void *),void *),(ts,dt,ctx));CHKERRQ(ierr); 52282bf6240SBarry Smith PetscFunctionReturn(0); 52382bf6240SBarry Smith } 52482bf6240SBarry Smith 52582bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 52682bf6240SBarry Smith 527ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscBool *); /* force argument to next function to not be extern C*/ 528fb2e594dSBarry Smith EXTERN_C_BEGIN 5294a2ae208SSatish Balay #undef __FUNCT__ 5304a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo" 5317087cfbeSBarry Smith PetscErrorCode TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void* ctx) 53282bf6240SBarry Smith { 53382bf6240SBarry Smith TS_Pseudo *pseudo; 53482bf6240SBarry Smith 53582bf6240SBarry Smith PetscFunctionBegin; 53682bf6240SBarry Smith pseudo = (TS_Pseudo*)ts->data; 53782bf6240SBarry Smith pseudo->verify = dt; 53882bf6240SBarry Smith pseudo->verifyctx = ctx; 53982bf6240SBarry Smith PetscFunctionReturn(0); 54082bf6240SBarry Smith } 541fb2e594dSBarry Smith EXTERN_C_END 54282bf6240SBarry Smith 543fb2e594dSBarry Smith EXTERN_C_BEGIN 5444a2ae208SSatish Balay #undef __FUNCT__ 5454a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStepIncrement_Pseudo" 5467087cfbeSBarry Smith PetscErrorCode TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc) 54782bf6240SBarry Smith { 5484bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 54982bf6240SBarry Smith 55082bf6240SBarry Smith PetscFunctionBegin; 55182bf6240SBarry Smith pseudo->dt_increment = inc; 55282bf6240SBarry Smith PetscFunctionReturn(0); 55382bf6240SBarry Smith } 554fb2e594dSBarry Smith EXTERN_C_END 55582bf6240SBarry Smith 556fb2e594dSBarry Smith EXTERN_C_BEGIN 5574a2ae208SSatish Balay #undef __FUNCT__ 5584a2ae208SSatish Balay #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt_Pseudo" 5597087cfbeSBarry Smith PetscErrorCode TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts) 56082bf6240SBarry Smith { 5614bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 56282bf6240SBarry Smith 56382bf6240SBarry Smith PetscFunctionBegin; 5644bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_TRUE; 56582bf6240SBarry Smith PetscFunctionReturn(0); 56682bf6240SBarry Smith } 567fb2e594dSBarry Smith EXTERN_C_END 56882bf6240SBarry Smith 5696849ba73SBarry Smith typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/ 570fb2e594dSBarry Smith EXTERN_C_BEGIN 5714a2ae208SSatish Balay #undef __FUNCT__ 5724a2ae208SSatish Balay #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo" 5737087cfbeSBarry Smith PetscErrorCode TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void* ctx) 57482bf6240SBarry Smith { 5754bbc92c1SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 57682bf6240SBarry Smith 57782bf6240SBarry Smith PetscFunctionBegin; 57882bf6240SBarry Smith pseudo->dt = dt; 57982bf6240SBarry Smith pseudo->dtctx = ctx; 58082bf6240SBarry Smith PetscFunctionReturn(0); 58182bf6240SBarry Smith } 582fb2e594dSBarry Smith EXTERN_C_END 58382bf6240SBarry Smith 58482bf6240SBarry Smith /* ----------------------------------------------------------------------------- */ 58510e6a065SJed Brown /*MC 58610e6a065SJed Brown TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping 58782bf6240SBarry Smith 58810e6a065SJed Brown This method solves equations of the form 58910e6a065SJed Brown 59010e6a065SJed Brown $ F(X,Xdot) = 0 59110e6a065SJed Brown 59210e6a065SJed Brown for steady state using the iteration 59310e6a065SJed Brown 59410e6a065SJed Brown $ [G'] S = -F(X,0) 59510e6a065SJed Brown $ X += S 59610e6a065SJed Brown 59710e6a065SJed Brown where 59810e6a065SJed Brown 59910e6a065SJed Brown $ G(Y) = F(Y,(Y-X)/dt) 60010e6a065SJed Brown 6016f2d6a7bSJed Brown This is linearly-implicit Euler with the residual always evaluated "at steady 6026f2d6a7bSJed Brown state". See note below. 60310e6a065SJed Brown 60410e6a065SJed Brown Options database keys: 60510e6a065SJed Brown + -ts_pseudo_increment <real> - ratio of increase dt 60610e6a065SJed Brown - -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt 60710e6a065SJed Brown 60810e6a065SJed Brown Level: beginner 60910e6a065SJed Brown 61010e6a065SJed Brown References: 61110e6a065SJed Brown Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential-Algebraic Equations, 2003. 61210e6a065SJed Brown C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998. 61310e6a065SJed Brown 61410e6a065SJed Brown Notes: 6156f2d6a7bSJed Brown The residual computed by this method includes the transient term (Xdot is computed instead of 6166f2d6a7bSJed Brown always being zero), but since the prediction from the last step is always the solution from the 6176f2d6a7bSJed Brown last step, on the first Newton iteration we have 6186f2d6a7bSJed Brown 6196f2d6a7bSJed Brown $ Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0 6206f2d6a7bSJed Brown 6216f2d6a7bSJed Brown Therefore, the linear system solved by the first Newton iteration is equivalent to the one 6226f2d6a7bSJed Brown described above and in the papers. If the user chooses to perform multiple Newton iterations, the 6236f2d6a7bSJed Brown algorithm is no longer the one described in the referenced papers. 62410e6a065SJed Brown 62510e6a065SJed Brown .seealso: TSCreate(), TS, TSSetType() 62610e6a065SJed Brown 62710e6a065SJed Brown M*/ 628fb2e594dSBarry Smith EXTERN_C_BEGIN 6294a2ae208SSatish Balay #undef __FUNCT__ 6304a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Pseudo" 6317087cfbeSBarry Smith PetscErrorCode TSCreate_Pseudo(TS ts) 6322d3f70b5SBarry Smith { 6337bf11e45SBarry Smith TS_Pseudo *pseudo; 634dfbe8321SBarry Smith PetscErrorCode ierr; 635193ac0bcSJed Brown SNES snes; 636193ac0bcSJed Brown const SNESType stype; 6372d3f70b5SBarry Smith 6383a40ed3dSBarry Smith PetscFunctionBegin; 639277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Pseudo; 640000e7ae3SMatthew Knepley ts->ops->destroy = TSDestroy_Pseudo; 641000e7ae3SMatthew Knepley ts->ops->view = TSView_Pseudo; 6422d3f70b5SBarry Smith 643000e7ae3SMatthew Knepley ts->ops->setup = TSSetUp_Pseudo; 644000e7ae3SMatthew Knepley ts->ops->step = TSStep_Pseudo; 645000e7ae3SMatthew Knepley ts->ops->setfromoptions = TSSetFromOptions_Pseudo; 6460f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Pseudo; 6470f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Pseudo; 6487bf11e45SBarry Smith 649193ac0bcSJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 650193ac0bcSJed Brown ierr = SNESGetType(snes,&stype);CHKERRQ(ierr); 651193ac0bcSJed Brown if (!stype) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);} 6522d3f70b5SBarry Smith 65338f2d2fdSLisandro Dalcin ierr = PetscNewLog(ts,TS_Pseudo,&pseudo);CHKERRQ(ierr); 6547bf11e45SBarry Smith ts->data = (void*)pseudo; 6552d3f70b5SBarry Smith 65628aa8177SBarry Smith pseudo->dt_increment = 1.1; 6574bbc92c1SBarry Smith pseudo->increment_dt_from_initial_dt = PETSC_FALSE; 65828aa8177SBarry Smith pseudo->dt = TSPseudoDefaultTimeStep; 659193ac0bcSJed Brown pseudo->fnorm = -1; 66082bf6240SBarry Smith 661f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C", 662e1311b90SBarry Smith "TSPseudoSetVerifyTimeStep_Pseudo", 6630c97f09dSSatish Balay TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr); 664f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C", 665e1311b90SBarry Smith "TSPseudoSetTimeStepIncrement_Pseudo", 6660c97f09dSSatish Balay TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr); 667f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C", 668e1311b90SBarry Smith "TSPseudoIncrementDtFromInitialDt_Pseudo", 6690c97f09dSSatish Balay TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr); 6700c97f09dSSatish Balay ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C", 6710c97f09dSSatish Balay "TSPseudoSetTimeStep_Pseudo", 6720c97f09dSSatish Balay TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr); 6733a40ed3dSBarry Smith PetscFunctionReturn(0); 6742d3f70b5SBarry Smith } 675fb2e594dSBarry Smith EXTERN_C_END 6762d3f70b5SBarry Smith 6774a2ae208SSatish Balay #undef __FUNCT__ 6784a2ae208SSatish Balay #define __FUNCT__ "TSPseudoDefaultTimeStep" 67982bf6240SBarry Smith /*@C 68082bf6240SBarry Smith TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping. 68182bf6240SBarry Smith Use with TSPseudoSetTimeStep(). 68228aa8177SBarry Smith 68315091d37SBarry Smith Collective on TS 68415091d37SBarry Smith 68528aa8177SBarry Smith Input Parameters: 68628aa8177SBarry Smith . ts - the timestep context 68782bf6240SBarry Smith . dtctx - unused timestep context 68828aa8177SBarry Smith 68982bf6240SBarry Smith Output Parameter: 69082bf6240SBarry Smith . newdt - the timestep to use for the next step 69128aa8177SBarry Smith 69215091d37SBarry Smith Level: advanced 69315091d37SBarry Smith 69482bf6240SBarry Smith .keywords: timestep, pseudo, default 695564e8f4eSLois Curfman McInnes 69682bf6240SBarry Smith .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep() 69728aa8177SBarry Smith @*/ 6987087cfbeSBarry Smith PetscErrorCode TSPseudoDefaultTimeStep(TS ts,PetscReal* newdt,void* dtctx) 69928aa8177SBarry Smith { 70082bf6240SBarry Smith TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 70187828ca2SBarry Smith PetscReal inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous; 702dfbe8321SBarry Smith PetscErrorCode ierr; 70328aa8177SBarry Smith 7043a40ed3dSBarry Smith PetscFunctionBegin; 705bbd7b040SJed Brown ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr); 706214bc6a2SJed Brown ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func,PETSC_FALSE);CHKERRQ(ierr); 70782bf6240SBarry Smith ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 708*cdbf8f93SLisandro Dalcin if (pseudo->fnorm_initial == 0.0) { 70982bf6240SBarry Smith /* first time through so compute initial function norm */ 710*cdbf8f93SLisandro Dalcin pseudo->fnorm_initial = pseudo->fnorm; 71182bf6240SBarry Smith fnorm_previous = pseudo->fnorm; 71282bf6240SBarry Smith } 71382bf6240SBarry Smith if (pseudo->fnorm == 0.0) { 71482bf6240SBarry Smith *newdt = 1.e12*inc*ts->time_step; 715004884caSBarry Smith } else if (pseudo->increment_dt_from_initial_dt) { 716*cdbf8f93SLisandro Dalcin *newdt = inc*pseudo->dt_initial*pseudo->fnorm_initial/pseudo->fnorm; 71782bf6240SBarry Smith } else { 71882bf6240SBarry Smith *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm; 71982bf6240SBarry Smith } 72082bf6240SBarry Smith pseudo->fnorm_previous = pseudo->fnorm; 7213a40ed3dSBarry Smith PetscFunctionReturn(0); 72228aa8177SBarry Smith } 723