Actual source code: euler.c

  1: /*
  2:        Code for Timestepping with explicit Euler.
  3: */
  4: #include <private/tsimpl.h>                /*I   "petscts.h"   I*/

  6: typedef struct {
  7:   Vec update;     /* work vector where new solution is formed  */
  8: } TS_Euler;

 12: static PetscErrorCode TSStep_Euler(TS ts)
 13: {
 14:   TS_Euler       *euler = (TS_Euler*)ts->data;
 15:   Vec            sol = ts->vec_sol,update = euler->update;

 19:   TSComputeRHSFunction(ts,ts->ptime,sol,update);
 20:   VecAXPY(sol,ts->time_step,update);
 21:   ts->ptime += ts->time_step;
 22:   ts->steps++;
 23:   return(0);
 24: }
 25: /*------------------------------------------------------------*/

 29: static PetscErrorCode TSSetUp_Euler(TS ts)
 30: {
 31:   TS_Euler       *euler = (TS_Euler*)ts->data;

 35:   VecDuplicate(ts->vec_sol,&euler->update);
 36:   return(0);
 37: }

 41: static PetscErrorCode TSReset_Euler(TS ts)
 42: {
 43:   TS_Euler       *euler = (TS_Euler*)ts->data;

 47:   VecDestroy(&euler->update);
 48:   return(0);
 49: }

 53: static PetscErrorCode TSDestroy_Euler(TS ts)
 54: {

 58:   TSReset_Euler(ts);
 59:   PetscFree(ts->data);
 60:   return(0);
 61: }
 62: /*------------------------------------------------------------*/

 66: static PetscErrorCode TSSetFromOptions_Euler(TS ts)
 67: {
 69:   return(0);
 70: }

 74: static PetscErrorCode TSView_Euler(TS ts,PetscViewer viewer)
 75: {
 77:   return(0);
 78: }

 82: static PetscErrorCode TSInterpolate_Euler(TS ts,PetscReal t,Vec X)
 83: {
 84:   PetscReal      alpha = (ts->ptime - t)/ts->time_step;

 88:   VecAXPBY(ts->vec_sol,1.0-alpha,alpha,X);
 89:   return(0);
 90: }

 92: /* ------------------------------------------------------------ */

 94: /*MC
 95:       TSEULER - ODE solver using the explicit forward Euler method

 97:   Level: beginner

 99: .seealso:  TSCreate(), TS, TSSetType(), TSBEULER

101: M*/
105: PetscErrorCode  TSCreate_Euler(TS ts)
106: {
107:   TS_Euler       *euler;

111:   ts->ops->setup           = TSSetUp_Euler;
112:   ts->ops->step            = TSStep_Euler;
113:   ts->ops->reset           = TSReset_Euler;
114:   ts->ops->destroy         = TSDestroy_Euler;
115:   ts->ops->setfromoptions  = TSSetFromOptions_Euler;
116:   ts->ops->view            = TSView_Euler;
117:   ts->ops->interpolate     = TSInterpolate_Euler;

119:   PetscNewLog(ts,TS_Euler,&euler);
120:   ts->data = (void*)euler;

122:   return(0);
123: }