xref: /petsc/src/ts/impls/explicit/euler/euler.c (revision d183316b7f3374e75ef313649153ab4029cc8afb)
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;
17cb9d8021SPierre Barbier de Reuille   PetscBool      accept;
184b33d51aSBarry Smith 
193a40ed3dSBarry Smith   PetscFunctionBegin;
20b8123daeSJed Brown   ierr = TSPreStep(ts);CHKERRQ(ierr);
21b8123daeSJed Brown   ierr = TSPreStage(ts,ts->ptime);CHKERRQ(ierr);
22c3e30b67SBarry Smith   ierr = TSComputeRHSFunction(ts,ts->ptime,sol,update);CHKERRQ(ierr);
23186e87acSLisandro Dalcin   ierr = VecAXPY(sol,ts->time_step,update);CHKERRQ(ierr);
24*d183316bSPierre Barbier de Reuille   ierr = TSFunctionDomainError(ts,ts->ptime,sol,&accept);CHKERRQ(ierr);
25cb9d8021SPierre Barbier de Reuille   if(!accept) {
26cb9d8021SPierre Barbier de Reuille     ts->reason = TS_DIVERGED_STEP_REJECTED;
27cb9d8021SPierre Barbier de Reuille     PetscFunctionReturn(0);
28cb9d8021SPierre Barbier de Reuille   }
299be3e283SDebojyoti Ghosh   ierr = TSPostStage(ts,ts->ptime,0,&sol);CHKERRQ(ierr);
30186e87acSLisandro Dalcin   ts->ptime += ts->time_step;
318b1af7b3SBarry Smith   ts->steps++;
323a40ed3dSBarry Smith   PetscFunctionReturn(0);
334b33d51aSBarry Smith }
344b33d51aSBarry Smith /*------------------------------------------------------------*/
35277b19d0SLisandro Dalcin 
364a2ae208SSatish Balay #undef __FUNCT__
37277b19d0SLisandro Dalcin #define __FUNCT__ "TSSetUp_Euler"
38277b19d0SLisandro Dalcin static PetscErrorCode TSSetUp_Euler(TS ts)
39277b19d0SLisandro Dalcin {
40277b19d0SLisandro Dalcin   TS_Euler       *euler = (TS_Euler*)ts->data;
41277b19d0SLisandro Dalcin   PetscErrorCode ierr;
42277b19d0SLisandro Dalcin 
43277b19d0SLisandro Dalcin   PetscFunctionBegin;
44277b19d0SLisandro Dalcin   ierr = VecDuplicate(ts->vec_sol,&euler->update);CHKERRQ(ierr);
45277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
46277b19d0SLisandro Dalcin }
47277b19d0SLisandro Dalcin 
48277b19d0SLisandro Dalcin #undef __FUNCT__
49277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Euler"
50277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Euler(TS ts)
514b33d51aSBarry Smith {
528b1af7b3SBarry Smith   TS_Euler       *euler = (TS_Euler*)ts->data;
53dfbe8321SBarry Smith   PetscErrorCode ierr;
548b1af7b3SBarry Smith 
553a40ed3dSBarry Smith   PetscFunctionBegin;
566bf464f9SBarry Smith   ierr = VecDestroy(&euler->update);CHKERRQ(ierr);
57277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
58277b19d0SLisandro Dalcin }
59277b19d0SLisandro Dalcin 
60277b19d0SLisandro Dalcin #undef __FUNCT__
61277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Euler"
62277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Euler(TS ts)
63277b19d0SLisandro Dalcin {
64277b19d0SLisandro Dalcin   PetscErrorCode ierr;
65277b19d0SLisandro Dalcin 
66277b19d0SLisandro Dalcin   PetscFunctionBegin;
67277b19d0SLisandro Dalcin   ierr = TSReset_Euler(ts);CHKERRQ(ierr);
68277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
693a40ed3dSBarry Smith   PetscFunctionReturn(0);
708b1af7b3SBarry Smith }
718b1af7b3SBarry Smith /*------------------------------------------------------------*/
728b1af7b3SBarry Smith 
734a2ae208SSatish Balay #undef __FUNCT__
744a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Euler"
758c34d3f5SBarry Smith static PetscErrorCode TSSetFromOptions_Euler(PetscOptions *PetscOptionsObject,TS ts)
768b1af7b3SBarry Smith {
773a40ed3dSBarry Smith   PetscFunctionBegin;
783a40ed3dSBarry Smith   PetscFunctionReturn(0);
794b33d51aSBarry Smith }
804b33d51aSBarry Smith 
814a2ae208SSatish Balay #undef __FUNCT__
824a2ae208SSatish Balay #define __FUNCT__ "TSView_Euler"
836849ba73SBarry Smith static PetscErrorCode TSView_Euler(TS ts,PetscViewer viewer)
848b1af7b3SBarry Smith {
853a40ed3dSBarry Smith   PetscFunctionBegin;
863a40ed3dSBarry Smith   PetscFunctionReturn(0);
878b1af7b3SBarry Smith }
888b1af7b3SBarry Smith 
896083293cSBarry Smith #undef __FUNCT__
906083293cSBarry Smith #define __FUNCT__ "TSInterpolate_Euler"
916083293cSBarry Smith static PetscErrorCode TSInterpolate_Euler(TS ts,PetscReal t,Vec X)
926083293cSBarry Smith {
936083293cSBarry Smith   PetscReal      alpha = (ts->ptime - t)/ts->time_step;
946083293cSBarry Smith   PetscErrorCode ierr;
956083293cSBarry Smith 
966083293cSBarry Smith   PetscFunctionBegin;
976083293cSBarry Smith   ierr = VecAXPBY(ts->vec_sol,1.0-alpha,alpha,X);CHKERRQ(ierr);
986083293cSBarry Smith   PetscFunctionReturn(0);
996083293cSBarry Smith }
1006083293cSBarry Smith 
101f9c1d6abSBarry Smith #undef __FUNCT__
102f9c1d6abSBarry Smith #define __FUNCT__ "TSComputeLinearStability_Euler"
103f9c1d6abSBarry Smith PetscErrorCode TSComputeLinearStability_Euler(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
104f9c1d6abSBarry Smith {
105f9c1d6abSBarry Smith   PetscFunctionBegin;
106f9c1d6abSBarry Smith   *yr = 1.0 + xr;
107f9c1d6abSBarry Smith   *yi = xi;
108f9c1d6abSBarry Smith   PetscFunctionReturn(0);
109f9c1d6abSBarry Smith }
1108b1af7b3SBarry Smith /* ------------------------------------------------------------ */
111ebe3b25bSBarry Smith 
112ebe3b25bSBarry Smith /*MC
1139596e0b4SJed Brown       TSEULER - ODE solver using the explicit forward Euler method
114ebe3b25bSBarry Smith 
115d41222bbSBarry Smith   Level: beginner
116d41222bbSBarry Smith 
1179596e0b4SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSBEULER
118ebe3b25bSBarry Smith 
119ebe3b25bSBarry Smith M*/
1204a2ae208SSatish Balay #undef __FUNCT__
1214a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Euler"
1228cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Euler(TS ts)
1238b1af7b3SBarry Smith {
1248b1af7b3SBarry Smith   TS_Euler       *euler;
125dfbe8321SBarry Smith   PetscErrorCode ierr;
1268b1af7b3SBarry Smith 
1273a40ed3dSBarry Smith   PetscFunctionBegin;
128000e7ae3SMatthew Knepley   ts->ops->setup           = TSSetUp_Euler;
129000e7ae3SMatthew Knepley   ts->ops->step            = TSStep_Euler;
130277b19d0SLisandro Dalcin   ts->ops->reset           = TSReset_Euler;
131000e7ae3SMatthew Knepley   ts->ops->destroy         = TSDestroy_Euler;
132000e7ae3SMatthew Knepley   ts->ops->setfromoptions  = TSSetFromOptions_Euler;
133000e7ae3SMatthew Knepley   ts->ops->view            = TSView_Euler;
1346083293cSBarry Smith   ts->ops->interpolate     = TSInterpolate_Euler;
135f9c1d6abSBarry Smith   ts->ops->linearstability = TSComputeLinearStability_Euler;
1368b1af7b3SBarry Smith 
137b00a9115SJed Brown   ierr = PetscNewLog(ts,&euler);CHKERRQ(ierr);
1388b1af7b3SBarry Smith   ts->data = (void*)euler;
1393a40ed3dSBarry Smith   PetscFunctionReturn(0);
1404b33d51aSBarry Smith }
141