xref: /petsc/src/ts/impls/explicit/euler/euler.c (revision b8123dae7f872ca50579e8005f6d2cd0e4534f13)
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;
19*b8123daeSJed Brown   ierr = TSPreStep(ts);CHKERRQ(ierr);
20*b8123daeSJed 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);
23186e87acSLisandro Dalcin   ts->ptime += ts->time_step;
248b1af7b3SBarry Smith   ts->steps++;
253a40ed3dSBarry Smith   PetscFunctionReturn(0);
264b33d51aSBarry Smith }
274b33d51aSBarry Smith /*------------------------------------------------------------*/
28277b19d0SLisandro Dalcin 
294a2ae208SSatish Balay #undef __FUNCT__
30277b19d0SLisandro Dalcin #define __FUNCT__ "TSSetUp_Euler"
31277b19d0SLisandro Dalcin static PetscErrorCode TSSetUp_Euler(TS ts)
32277b19d0SLisandro Dalcin {
33277b19d0SLisandro Dalcin   TS_Euler       *euler = (TS_Euler*)ts->data;
34277b19d0SLisandro Dalcin   PetscErrorCode ierr;
35277b19d0SLisandro Dalcin 
36277b19d0SLisandro Dalcin   PetscFunctionBegin;
37277b19d0SLisandro Dalcin   ierr = VecDuplicate(ts->vec_sol,&euler->update);CHKERRQ(ierr);
38277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
39277b19d0SLisandro Dalcin }
40277b19d0SLisandro Dalcin 
41277b19d0SLisandro Dalcin #undef __FUNCT__
42277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Euler"
43277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Euler(TS ts)
444b33d51aSBarry Smith {
458b1af7b3SBarry Smith   TS_Euler       *euler = (TS_Euler*)ts->data;
46dfbe8321SBarry Smith   PetscErrorCode ierr;
478b1af7b3SBarry Smith 
483a40ed3dSBarry Smith   PetscFunctionBegin;
496bf464f9SBarry Smith   ierr = VecDestroy(&euler->update);CHKERRQ(ierr);
50277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
51277b19d0SLisandro Dalcin }
52277b19d0SLisandro Dalcin 
53277b19d0SLisandro Dalcin #undef __FUNCT__
54277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Euler"
55277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Euler(TS ts)
56277b19d0SLisandro Dalcin {
57277b19d0SLisandro Dalcin   PetscErrorCode ierr;
58277b19d0SLisandro Dalcin 
59277b19d0SLisandro Dalcin   PetscFunctionBegin;
60277b19d0SLisandro Dalcin   ierr = TSReset_Euler(ts);CHKERRQ(ierr);
61277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
623a40ed3dSBarry Smith   PetscFunctionReturn(0);
638b1af7b3SBarry Smith }
648b1af7b3SBarry Smith /*------------------------------------------------------------*/
658b1af7b3SBarry Smith 
664a2ae208SSatish Balay #undef __FUNCT__
674a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Euler"
686849ba73SBarry Smith static PetscErrorCode TSSetFromOptions_Euler(TS ts)
698b1af7b3SBarry Smith {
703a40ed3dSBarry Smith   PetscFunctionBegin;
713a40ed3dSBarry Smith   PetscFunctionReturn(0);
724b33d51aSBarry Smith }
734b33d51aSBarry Smith 
744a2ae208SSatish Balay #undef __FUNCT__
754a2ae208SSatish Balay #define __FUNCT__ "TSView_Euler"
766849ba73SBarry Smith static PetscErrorCode TSView_Euler(TS ts,PetscViewer viewer)
778b1af7b3SBarry Smith {
783a40ed3dSBarry Smith   PetscFunctionBegin;
793a40ed3dSBarry Smith   PetscFunctionReturn(0);
808b1af7b3SBarry Smith }
818b1af7b3SBarry Smith 
826083293cSBarry Smith #undef __FUNCT__
836083293cSBarry Smith #define __FUNCT__ "TSInterpolate_Euler"
846083293cSBarry Smith static PetscErrorCode TSInterpolate_Euler(TS ts,PetscReal t,Vec X)
856083293cSBarry Smith {
866083293cSBarry Smith   PetscReal      alpha = (ts->ptime - t)/ts->time_step;
876083293cSBarry Smith   PetscErrorCode ierr;
886083293cSBarry Smith 
896083293cSBarry Smith   PetscFunctionBegin;
906083293cSBarry Smith   ierr = VecAXPBY(ts->vec_sol,1.0-alpha,alpha,X);CHKERRQ(ierr);
916083293cSBarry Smith   PetscFunctionReturn(0);
926083293cSBarry Smith }
936083293cSBarry Smith 
948b1af7b3SBarry Smith /* ------------------------------------------------------------ */
95ebe3b25bSBarry Smith 
96ebe3b25bSBarry Smith /*MC
979596e0b4SJed Brown       TSEULER - ODE solver using the explicit forward Euler method
98ebe3b25bSBarry Smith 
99d41222bbSBarry Smith   Level: beginner
100d41222bbSBarry Smith 
1019596e0b4SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSBEULER
102ebe3b25bSBarry Smith 
103ebe3b25bSBarry Smith M*/
104fb2e594dSBarry Smith EXTERN_C_BEGIN
1054a2ae208SSatish Balay #undef __FUNCT__
1064a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Euler"
1077087cfbeSBarry Smith PetscErrorCode  TSCreate_Euler(TS ts)
1088b1af7b3SBarry Smith {
1098b1af7b3SBarry Smith   TS_Euler       *euler;
110dfbe8321SBarry Smith   PetscErrorCode ierr;
1118b1af7b3SBarry Smith 
1123a40ed3dSBarry Smith   PetscFunctionBegin;
113000e7ae3SMatthew Knepley   ts->ops->setup           = TSSetUp_Euler;
114000e7ae3SMatthew Knepley   ts->ops->step            = TSStep_Euler;
115277b19d0SLisandro Dalcin   ts->ops->reset           = TSReset_Euler;
116000e7ae3SMatthew Knepley   ts->ops->destroy         = TSDestroy_Euler;
117000e7ae3SMatthew Knepley   ts->ops->setfromoptions  = TSSetFromOptions_Euler;
118000e7ae3SMatthew Knepley   ts->ops->view            = TSView_Euler;
1196083293cSBarry Smith   ts->ops->interpolate     = TSInterpolate_Euler;
1208b1af7b3SBarry Smith 
12138f2d2fdSLisandro Dalcin   ierr = PetscNewLog(ts,TS_Euler,&euler);CHKERRQ(ierr);
1228b1af7b3SBarry Smith   ts->data = (void*)euler;
1234b33d51aSBarry Smith 
1243a40ed3dSBarry Smith   PetscFunctionReturn(0);
1254b33d51aSBarry Smith }
126fb2e594dSBarry Smith EXTERN_C_END
127