xref: /petsc/src/ts/impls/explicit/euler/euler.c (revision 193ac0bc5128976c4ec2c1a67f3f9cb026b77f22)
14b33d51aSBarry Smith /*
2fb4a63b6SLois Curfman McInnes        Code for Timestepping with explicit Euler.
34b33d51aSBarry Smith */
4c6db04a5SJed Brown #include <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"
12*193ac0bcSJed 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;
19c3e30b67SBarry Smith   ierr = TSComputeRHSFunction(ts,ts->ptime,sol,update);CHKERRQ(ierr);
20186e87acSLisandro Dalcin   ierr = VecAXPY(sol,ts->time_step,update);CHKERRQ(ierr);
21186e87acSLisandro Dalcin   ts->ptime += ts->time_step;
228b1af7b3SBarry Smith   ts->steps++;
233a40ed3dSBarry Smith   PetscFunctionReturn(0);
244b33d51aSBarry Smith }
254b33d51aSBarry Smith /*------------------------------------------------------------*/
26277b19d0SLisandro Dalcin 
274a2ae208SSatish Balay #undef __FUNCT__
28277b19d0SLisandro Dalcin #define __FUNCT__ "TSSetUp_Euler"
29277b19d0SLisandro Dalcin static PetscErrorCode TSSetUp_Euler(TS ts)
30277b19d0SLisandro Dalcin {
31277b19d0SLisandro Dalcin   TS_Euler       *euler = (TS_Euler*)ts->data;
32277b19d0SLisandro Dalcin   PetscErrorCode ierr;
33277b19d0SLisandro Dalcin 
34277b19d0SLisandro Dalcin   PetscFunctionBegin;
35277b19d0SLisandro Dalcin   ierr = VecDuplicate(ts->vec_sol,&euler->update);CHKERRQ(ierr);
36277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
37277b19d0SLisandro Dalcin }
38277b19d0SLisandro Dalcin 
39277b19d0SLisandro Dalcin #undef __FUNCT__
40277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Euler"
41277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Euler(TS ts)
424b33d51aSBarry Smith {
438b1af7b3SBarry Smith   TS_Euler       *euler = (TS_Euler*)ts->data;
44dfbe8321SBarry Smith   PetscErrorCode ierr;
458b1af7b3SBarry Smith 
463a40ed3dSBarry Smith   PetscFunctionBegin;
476bf464f9SBarry Smith   ierr = VecDestroy(&euler->update);CHKERRQ(ierr);
48277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
49277b19d0SLisandro Dalcin }
50277b19d0SLisandro Dalcin 
51277b19d0SLisandro Dalcin #undef __FUNCT__
52277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Euler"
53277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Euler(TS ts)
54277b19d0SLisandro Dalcin {
55277b19d0SLisandro Dalcin   PetscErrorCode ierr;
56277b19d0SLisandro Dalcin 
57277b19d0SLisandro Dalcin   PetscFunctionBegin;
58277b19d0SLisandro Dalcin   ierr = TSReset_Euler(ts);CHKERRQ(ierr);
59277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
603a40ed3dSBarry Smith   PetscFunctionReturn(0);
618b1af7b3SBarry Smith }
628b1af7b3SBarry Smith /*------------------------------------------------------------*/
638b1af7b3SBarry Smith 
644a2ae208SSatish Balay #undef __FUNCT__
654a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Euler"
666849ba73SBarry Smith static PetscErrorCode TSSetFromOptions_Euler(TS ts)
678b1af7b3SBarry Smith {
683a40ed3dSBarry Smith   PetscFunctionBegin;
693a40ed3dSBarry Smith   PetscFunctionReturn(0);
704b33d51aSBarry Smith }
714b33d51aSBarry Smith 
724a2ae208SSatish Balay #undef __FUNCT__
734a2ae208SSatish Balay #define __FUNCT__ "TSView_Euler"
746849ba73SBarry Smith static PetscErrorCode TSView_Euler(TS ts,PetscViewer viewer)
758b1af7b3SBarry Smith {
763a40ed3dSBarry Smith   PetscFunctionBegin;
773a40ed3dSBarry Smith   PetscFunctionReturn(0);
788b1af7b3SBarry Smith }
798b1af7b3SBarry Smith 
808b1af7b3SBarry Smith /* ------------------------------------------------------------ */
81ebe3b25bSBarry Smith 
82ebe3b25bSBarry Smith /*MC
839596e0b4SJed Brown       TSEULER - ODE solver using the explicit forward Euler method
84ebe3b25bSBarry Smith 
85d41222bbSBarry Smith   Level: beginner
86d41222bbSBarry Smith 
879596e0b4SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSBEULER
88ebe3b25bSBarry Smith 
89ebe3b25bSBarry Smith M*/
90fb2e594dSBarry Smith EXTERN_C_BEGIN
914a2ae208SSatish Balay #undef __FUNCT__
924a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Euler"
937087cfbeSBarry Smith PetscErrorCode  TSCreate_Euler(TS ts)
948b1af7b3SBarry Smith {
958b1af7b3SBarry Smith   TS_Euler       *euler;
96dfbe8321SBarry Smith   PetscErrorCode ierr;
978b1af7b3SBarry Smith 
983a40ed3dSBarry Smith   PetscFunctionBegin;
99000e7ae3SMatthew Knepley   ts->ops->setup           = TSSetUp_Euler;
100000e7ae3SMatthew Knepley   ts->ops->step            = TSStep_Euler;
101277b19d0SLisandro Dalcin   ts->ops->reset           = TSReset_Euler;
102000e7ae3SMatthew Knepley   ts->ops->destroy         = TSDestroy_Euler;
103000e7ae3SMatthew Knepley   ts->ops->setfromoptions  = TSSetFromOptions_Euler;
104000e7ae3SMatthew Knepley   ts->ops->view            = TSView_Euler;
1058b1af7b3SBarry Smith 
10638f2d2fdSLisandro Dalcin   ierr = PetscNewLog(ts,TS_Euler,&euler);CHKERRQ(ierr);
1078b1af7b3SBarry Smith   ts->data = (void*)euler;
1084b33d51aSBarry Smith 
1093a40ed3dSBarry Smith   PetscFunctionReturn(0);
1104b33d51aSBarry Smith }
111fb2e594dSBarry Smith EXTERN_C_END
112