xref: /petsc/src/ts/impls/explicit/euler/euler.c (revision 277b19d07ec08e548e5816b82c213278d45c3326)
163dd3a1aSKris Buschelman 
24b33d51aSBarry Smith /*
3fb4a63b6SLois Curfman McInnes        Code for Timestepping with explicit Euler.
44b33d51aSBarry Smith */
5c6db04a5SJed Brown #include <private/tsimpl.h>                /*I   "petscts.h"   I*/
64b33d51aSBarry Smith 
78b1af7b3SBarry Smith typedef struct {
8*277b19d0SLisandro Dalcin   Vec update;     /* work vector where new solution is formed  */
98b1af7b3SBarry Smith } TS_Euler;
108b1af7b3SBarry Smith 
114a2ae208SSatish Balay #undef __FUNCT__
124a2ae208SSatish Balay #define __FUNCT__ "TSStep_Euler"
13a7cc72afSBarry Smith static PetscErrorCode TSStep_Euler(TS ts,PetscInt *steps,PetscReal *ptime)
144b33d51aSBarry Smith {
158b1af7b3SBarry Smith   TS_Euler       *euler = (TS_Euler*)ts->data;
168b1af7b3SBarry Smith   Vec            sol = ts->vec_sol,update = euler->update;
17a7cc72afSBarry Smith   PetscInt       i,max_steps = ts->max_steps;
18*277b19d0SLisandro Dalcin   PetscErrorCode ierr;
194b33d51aSBarry Smith 
203a40ed3dSBarry Smith   PetscFunctionBegin;
218b1af7b3SBarry Smith   *steps = -ts->steps;
22c3e30b67SBarry Smith   ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
234b33d51aSBarry Smith 
248b1af7b3SBarry Smith   for (i=0; i<max_steps; i++) {
25a057ae4aSMatthew Knepley     PetscReal dt = ts->time_step;
26a057ae4aSMatthew Knepley 
273f2090d5SJed Brown     ierr = TSPreStep(ts);CHKERRQ(ierr);
28a057ae4aSMatthew Knepley     ts->ptime += dt;
29c3e30b67SBarry Smith     ierr = TSComputeRHSFunction(ts,ts->ptime,sol,update);CHKERRQ(ierr);
302dcb1b2aSMatthew Knepley     ierr = VecAXPY(sol,dt,update);CHKERRQ(ierr);
318b1af7b3SBarry Smith     ts->steps++;
323f2090d5SJed Brown     ierr = TSPostStep(ts);CHKERRQ(ierr);
33c3e30b67SBarry Smith     ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
348b1af7b3SBarry Smith     if (ts->ptime > ts->max_time) break;
358b1af7b3SBarry Smith   }
364b33d51aSBarry Smith 
378b1af7b3SBarry Smith   *steps += ts->steps;
38142b95e3SSatish Balay   *ptime  = ts->ptime;
393a40ed3dSBarry Smith   PetscFunctionReturn(0);
404b33d51aSBarry Smith }
414b33d51aSBarry Smith /*------------------------------------------------------------*/
42*277b19d0SLisandro Dalcin 
434a2ae208SSatish Balay #undef __FUNCT__
44*277b19d0SLisandro Dalcin #define __FUNCT__ "TSSetUp_Euler"
45*277b19d0SLisandro Dalcin static PetscErrorCode TSSetUp_Euler(TS ts)
46*277b19d0SLisandro Dalcin {
47*277b19d0SLisandro Dalcin   TS_Euler       *euler = (TS_Euler*)ts->data;
48*277b19d0SLisandro Dalcin   PetscErrorCode ierr;
49*277b19d0SLisandro Dalcin 
50*277b19d0SLisandro Dalcin   PetscFunctionBegin;
51*277b19d0SLisandro Dalcin   ierr = VecDuplicate(ts->vec_sol,&euler->update);CHKERRQ(ierr);
52*277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
53*277b19d0SLisandro Dalcin }
54*277b19d0SLisandro Dalcin 
55*277b19d0SLisandro Dalcin #undef __FUNCT__
56*277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Euler"
57*277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Euler(TS ts)
584b33d51aSBarry Smith {
598b1af7b3SBarry Smith   TS_Euler       *euler = (TS_Euler*)ts->data;
60dfbe8321SBarry Smith   PetscErrorCode ierr;
618b1af7b3SBarry Smith 
623a40ed3dSBarry Smith   PetscFunctionBegin;
631713a123SBarry Smith   if (euler->update) {ierr = VecDestroy(euler->update);CHKERRQ(ierr);}
64*277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
65*277b19d0SLisandro Dalcin }
66*277b19d0SLisandro Dalcin 
67*277b19d0SLisandro Dalcin #undef __FUNCT__
68*277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Euler"
69*277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Euler(TS ts)
70*277b19d0SLisandro Dalcin {
71*277b19d0SLisandro Dalcin   PetscErrorCode ierr;
72*277b19d0SLisandro Dalcin 
73*277b19d0SLisandro Dalcin   PetscFunctionBegin;
74*277b19d0SLisandro Dalcin   ierr = TSReset_Euler(ts);CHKERRQ(ierr);
75*277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
763a40ed3dSBarry Smith   PetscFunctionReturn(0);
778b1af7b3SBarry Smith }
788b1af7b3SBarry Smith /*------------------------------------------------------------*/
798b1af7b3SBarry Smith 
804a2ae208SSatish Balay #undef __FUNCT__
814a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Euler"
826849ba73SBarry Smith static PetscErrorCode TSSetFromOptions_Euler(TS ts)
838b1af7b3SBarry Smith {
843a40ed3dSBarry Smith   PetscFunctionBegin;
853a40ed3dSBarry Smith   PetscFunctionReturn(0);
864b33d51aSBarry Smith }
874b33d51aSBarry Smith 
884a2ae208SSatish Balay #undef __FUNCT__
894a2ae208SSatish Balay #define __FUNCT__ "TSView_Euler"
906849ba73SBarry Smith static PetscErrorCode TSView_Euler(TS ts,PetscViewer viewer)
918b1af7b3SBarry Smith {
923a40ed3dSBarry Smith   PetscFunctionBegin;
933a40ed3dSBarry Smith   PetscFunctionReturn(0);
948b1af7b3SBarry Smith }
958b1af7b3SBarry Smith 
968b1af7b3SBarry Smith /* ------------------------------------------------------------ */
97ebe3b25bSBarry Smith 
98ebe3b25bSBarry Smith /*MC
999596e0b4SJed Brown       TSEULER - ODE solver using the explicit forward Euler method
100ebe3b25bSBarry Smith 
101d41222bbSBarry Smith   Level: beginner
102d41222bbSBarry Smith 
1039596e0b4SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSBEULER
104ebe3b25bSBarry Smith 
105ebe3b25bSBarry Smith M*/
106fb2e594dSBarry Smith EXTERN_C_BEGIN
1074a2ae208SSatish Balay #undef __FUNCT__
1084a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Euler"
1097087cfbeSBarry Smith PetscErrorCode  TSCreate_Euler(TS ts)
1108b1af7b3SBarry Smith {
1118b1af7b3SBarry Smith   TS_Euler       *euler;
112dfbe8321SBarry Smith   PetscErrorCode ierr;
1138b1af7b3SBarry Smith 
1143a40ed3dSBarry Smith   PetscFunctionBegin;
115000e7ae3SMatthew Knepley   ts->ops->setup           = TSSetUp_Euler;
116000e7ae3SMatthew Knepley   ts->ops->step            = TSStep_Euler;
117*277b19d0SLisandro Dalcin   ts->ops->reset           = TSReset_Euler;
118000e7ae3SMatthew Knepley   ts->ops->destroy         = TSDestroy_Euler;
119000e7ae3SMatthew Knepley   ts->ops->setfromoptions  = TSSetFromOptions_Euler;
120000e7ae3SMatthew Knepley   ts->ops->view            = TSView_Euler;
1218b1af7b3SBarry Smith 
12238f2d2fdSLisandro Dalcin   ierr = PetscNewLog(ts,TS_Euler,&euler);CHKERRQ(ierr);
1238b1af7b3SBarry Smith   ts->data = (void*)euler;
1244b33d51aSBarry Smith 
1253a40ed3dSBarry Smith   PetscFunctionReturn(0);
1264b33d51aSBarry Smith }
127fb2e594dSBarry Smith EXTERN_C_END
128c3e30b67SBarry Smith 
129c3e30b67SBarry Smith 
130c3e30b67SBarry Smith 
131c3e30b67SBarry Smith 
132