xref: /petsc/src/ts/impls/explicit/euler/euler.c (revision 9be3e283b542d23e45ac09f5a477feae8bd90e67)
14b33d51aSBarry Smith /*
2fb4a63b6SLois Curfman McInnes        Code for Timestepping with explicit Euler.
34b33d51aSBarry Smith */
4b45d2f2cSJed Brown #include <petsc-private/tsimpl.h>                /*I   "petscts.h"   I*/
54b33d51aSBarry Smith 
68b1af7b3SBarry Smith typedef struct {
7277b19d0SLisandro Dalcin   Vec update;     /* work vector where new solution is formed  */
88b1af7b3SBarry Smith } TS_Euler;
98b1af7b3SBarry Smith 
104a2ae208SSatish Balay #undef __FUNCT__
114a2ae208SSatish Balay #define __FUNCT__ "TSStep_Euler"
12193ac0bcSJed Brown static PetscErrorCode TSStep_Euler(TS ts)
134b33d51aSBarry Smith {
148b1af7b3SBarry Smith   TS_Euler       *euler = (TS_Euler*)ts->data;
158b1af7b3SBarry Smith   Vec            sol    = ts->vec_sol,update = euler->update;
16277b19d0SLisandro Dalcin   PetscErrorCode ierr;
174b33d51aSBarry Smith 
183a40ed3dSBarry Smith   PetscFunctionBegin;
19b8123daeSJed Brown   ierr = TSPreStep(ts);CHKERRQ(ierr);
20b8123daeSJed Brown   ierr = TSPreStage(ts,ts->ptime);CHKERRQ(ierr);
21c3e30b67SBarry Smith   ierr = TSComputeRHSFunction(ts,ts->ptime,sol,update);CHKERRQ(ierr);
22186e87acSLisandro Dalcin   ierr = VecAXPY(sol,ts->time_step,update);CHKERRQ(ierr);
23*9be3e283SDebojyoti Ghosh   ierr = TSPostStage(ts,ts->ptime,0,&sol);CHKERRQ(ierr);
24186e87acSLisandro Dalcin   ts->ptime += ts->time_step;
258b1af7b3SBarry Smith   ts->steps++;
263a40ed3dSBarry Smith   PetscFunctionReturn(0);
274b33d51aSBarry Smith }
284b33d51aSBarry Smith /*------------------------------------------------------------*/
29277b19d0SLisandro Dalcin 
304a2ae208SSatish Balay #undef __FUNCT__
31277b19d0SLisandro Dalcin #define __FUNCT__ "TSSetUp_Euler"
32277b19d0SLisandro Dalcin static PetscErrorCode TSSetUp_Euler(TS ts)
33277b19d0SLisandro Dalcin {
34277b19d0SLisandro Dalcin   TS_Euler       *euler = (TS_Euler*)ts->data;
35277b19d0SLisandro Dalcin   PetscErrorCode ierr;
36277b19d0SLisandro Dalcin 
37277b19d0SLisandro Dalcin   PetscFunctionBegin;
38277b19d0SLisandro Dalcin   ierr = VecDuplicate(ts->vec_sol,&euler->update);CHKERRQ(ierr);
39277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
40277b19d0SLisandro Dalcin }
41277b19d0SLisandro Dalcin 
42277b19d0SLisandro Dalcin #undef __FUNCT__
43277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Euler"
44277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Euler(TS ts)
454b33d51aSBarry Smith {
468b1af7b3SBarry Smith   TS_Euler       *euler = (TS_Euler*)ts->data;
47dfbe8321SBarry Smith   PetscErrorCode ierr;
488b1af7b3SBarry Smith 
493a40ed3dSBarry Smith   PetscFunctionBegin;
506bf464f9SBarry Smith   ierr = VecDestroy(&euler->update);CHKERRQ(ierr);
51277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
52277b19d0SLisandro Dalcin }
53277b19d0SLisandro Dalcin 
54277b19d0SLisandro Dalcin #undef __FUNCT__
55277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Euler"
56277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Euler(TS ts)
57277b19d0SLisandro Dalcin {
58277b19d0SLisandro Dalcin   PetscErrorCode ierr;
59277b19d0SLisandro Dalcin 
60277b19d0SLisandro Dalcin   PetscFunctionBegin;
61277b19d0SLisandro Dalcin   ierr = TSReset_Euler(ts);CHKERRQ(ierr);
62277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
633a40ed3dSBarry Smith   PetscFunctionReturn(0);
648b1af7b3SBarry Smith }
658b1af7b3SBarry Smith /*------------------------------------------------------------*/
668b1af7b3SBarry Smith 
674a2ae208SSatish Balay #undef __FUNCT__
684a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Euler"
696849ba73SBarry Smith static PetscErrorCode TSSetFromOptions_Euler(TS ts)
708b1af7b3SBarry Smith {
713a40ed3dSBarry Smith   PetscFunctionBegin;
723a40ed3dSBarry Smith   PetscFunctionReturn(0);
734b33d51aSBarry Smith }
744b33d51aSBarry Smith 
754a2ae208SSatish Balay #undef __FUNCT__
764a2ae208SSatish Balay #define __FUNCT__ "TSView_Euler"
776849ba73SBarry Smith static PetscErrorCode TSView_Euler(TS ts,PetscViewer viewer)
788b1af7b3SBarry Smith {
793a40ed3dSBarry Smith   PetscFunctionBegin;
803a40ed3dSBarry Smith   PetscFunctionReturn(0);
818b1af7b3SBarry Smith }
828b1af7b3SBarry Smith 
836083293cSBarry Smith #undef __FUNCT__
846083293cSBarry Smith #define __FUNCT__ "TSInterpolate_Euler"
856083293cSBarry Smith static PetscErrorCode TSInterpolate_Euler(TS ts,PetscReal t,Vec X)
866083293cSBarry Smith {
876083293cSBarry Smith   PetscReal      alpha = (ts->ptime - t)/ts->time_step;
886083293cSBarry Smith   PetscErrorCode ierr;
896083293cSBarry Smith 
906083293cSBarry Smith   PetscFunctionBegin;
916083293cSBarry Smith   ierr = VecAXPBY(ts->vec_sol,1.0-alpha,alpha,X);CHKERRQ(ierr);
926083293cSBarry Smith   PetscFunctionReturn(0);
936083293cSBarry Smith }
946083293cSBarry Smith 
95f9c1d6abSBarry Smith #undef __FUNCT__
96f9c1d6abSBarry Smith #define __FUNCT__ "TSComputeLinearStability_Euler"
97f9c1d6abSBarry Smith PetscErrorCode TSComputeLinearStability_Euler(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
98f9c1d6abSBarry Smith {
99f9c1d6abSBarry Smith   PetscFunctionBegin;
100f9c1d6abSBarry Smith   *yr = 1.0 + xr;
101f9c1d6abSBarry Smith   *yi = xi;
102f9c1d6abSBarry Smith   PetscFunctionReturn(0);
103f9c1d6abSBarry Smith }
1048b1af7b3SBarry Smith /* ------------------------------------------------------------ */
105ebe3b25bSBarry Smith 
106ebe3b25bSBarry Smith /*MC
1079596e0b4SJed Brown       TSEULER - ODE solver using the explicit forward Euler method
108ebe3b25bSBarry Smith 
109d41222bbSBarry Smith   Level: beginner
110d41222bbSBarry Smith 
1119596e0b4SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSBEULER
112ebe3b25bSBarry Smith 
113ebe3b25bSBarry Smith M*/
1144a2ae208SSatish Balay #undef __FUNCT__
1154a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Euler"
1168cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Euler(TS ts)
1178b1af7b3SBarry Smith {
1188b1af7b3SBarry Smith   TS_Euler       *euler;
119dfbe8321SBarry Smith   PetscErrorCode ierr;
1208b1af7b3SBarry Smith 
1213a40ed3dSBarry Smith   PetscFunctionBegin;
122000e7ae3SMatthew Knepley   ts->ops->setup           = TSSetUp_Euler;
123000e7ae3SMatthew Knepley   ts->ops->step            = TSStep_Euler;
124277b19d0SLisandro Dalcin   ts->ops->reset           = TSReset_Euler;
125000e7ae3SMatthew Knepley   ts->ops->destroy         = TSDestroy_Euler;
126000e7ae3SMatthew Knepley   ts->ops->setfromoptions  = TSSetFromOptions_Euler;
127000e7ae3SMatthew Knepley   ts->ops->view            = TSView_Euler;
1286083293cSBarry Smith   ts->ops->interpolate     = TSInterpolate_Euler;
129f9c1d6abSBarry Smith   ts->ops->linearstability = TSComputeLinearStability_Euler;
1308b1af7b3SBarry Smith 
13138f2d2fdSLisandro Dalcin   ierr = PetscNewLog(ts,TS_Euler,&euler);CHKERRQ(ierr);
1328b1af7b3SBarry Smith   ts->data = (void*)euler;
1333a40ed3dSBarry Smith   PetscFunctionReturn(0);
1344b33d51aSBarry Smith }
135