xref: /petsc/src/ts/impls/explicit/euler/euler.c (revision 6bf464f92cc51e6fd6163850774a6badb2f63b6b)
1 /*
2        Code for Timestepping with explicit Euler.
3 */
4 #include <private/tsimpl.h>                /*I   "petscts.h"   I*/
5 
6 typedef struct {
7   Vec update;     /* work vector where new solution is formed  */
8 } TS_Euler;
9 
10 #undef __FUNCT__
11 #define __FUNCT__ "TSStep_Euler"
12 static PetscErrorCode TSStep_Euler(TS ts,PetscInt *steps,PetscReal *ptime)
13 {
14   TS_Euler       *euler = (TS_Euler*)ts->data;
15   Vec            sol = ts->vec_sol,update = euler->update;
16   PetscInt       i;
17   PetscErrorCode ierr;
18 
19   PetscFunctionBegin;
20   *steps = -ts->steps;
21   *ptime  = ts->ptime;
22 
23   ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
24 
25   for (i=0; i<ts->max_steps; i++) {
26     if (ts->ptime + ts->time_step > ts->max_time) break;
27     ierr = TSPreStep(ts);CHKERRQ(ierr);
28 
29     ierr = TSComputeRHSFunction(ts,ts->ptime,sol,update);CHKERRQ(ierr);
30 
31     ierr = VecAXPY(sol,ts->time_step,update);CHKERRQ(ierr);
32     ts->ptime += ts->time_step;
33     ts->steps++;
34 
35     ierr = TSPostStep(ts);CHKERRQ(ierr);
36     ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
37   }
38 
39   *steps += ts->steps;
40   *ptime  = ts->ptime;
41   PetscFunctionReturn(0);
42 }
43 /*------------------------------------------------------------*/
44 
45 #undef __FUNCT__
46 #define __FUNCT__ "TSSetUp_Euler"
47 static PetscErrorCode TSSetUp_Euler(TS ts)
48 {
49   TS_Euler       *euler = (TS_Euler*)ts->data;
50   PetscErrorCode ierr;
51 
52   PetscFunctionBegin;
53   ierr = VecDuplicate(ts->vec_sol,&euler->update);CHKERRQ(ierr);
54   PetscFunctionReturn(0);
55 }
56 
57 #undef __FUNCT__
58 #define __FUNCT__ "TSReset_Euler"
59 static PetscErrorCode TSReset_Euler(TS ts)
60 {
61   TS_Euler       *euler = (TS_Euler*)ts->data;
62   PetscErrorCode ierr;
63 
64   PetscFunctionBegin;
65   ierr = VecDestroy(&euler->update);CHKERRQ(ierr);
66   PetscFunctionReturn(0);
67 }
68 
69 #undef __FUNCT__
70 #define __FUNCT__ "TSDestroy_Euler"
71 static PetscErrorCode TSDestroy_Euler(TS ts)
72 {
73   PetscErrorCode ierr;
74 
75   PetscFunctionBegin;
76   ierr = TSReset_Euler(ts);CHKERRQ(ierr);
77   ierr = PetscFree(ts->data);CHKERRQ(ierr);
78   PetscFunctionReturn(0);
79 }
80 /*------------------------------------------------------------*/
81 
82 #undef __FUNCT__
83 #define __FUNCT__ "TSSetFromOptions_Euler"
84 static PetscErrorCode TSSetFromOptions_Euler(TS ts)
85 {
86   PetscFunctionBegin;
87   PetscFunctionReturn(0);
88 }
89 
90 #undef __FUNCT__
91 #define __FUNCT__ "TSView_Euler"
92 static PetscErrorCode TSView_Euler(TS ts,PetscViewer viewer)
93 {
94   PetscFunctionBegin;
95   PetscFunctionReturn(0);
96 }
97 
98 /* ------------------------------------------------------------ */
99 
100 /*MC
101       TSEULER - ODE solver using the explicit forward Euler method
102 
103   Level: beginner
104 
105 .seealso:  TSCreate(), TS, TSSetType(), TSBEULER
106 
107 M*/
108 EXTERN_C_BEGIN
109 #undef __FUNCT__
110 #define __FUNCT__ "TSCreate_Euler"
111 PetscErrorCode  TSCreate_Euler(TS ts)
112 {
113   TS_Euler       *euler;
114   PetscErrorCode ierr;
115 
116   PetscFunctionBegin;
117   ts->ops->setup           = TSSetUp_Euler;
118   ts->ops->step            = TSStep_Euler;
119   ts->ops->reset           = TSReset_Euler;
120   ts->ops->destroy         = TSDestroy_Euler;
121   ts->ops->setfromoptions  = TSSetFromOptions_Euler;
122   ts->ops->view            = TSView_Euler;
123 
124   ierr = PetscNewLog(ts,TS_Euler,&euler);CHKERRQ(ierr);
125   ts->data = (void*)euler;
126 
127   PetscFunctionReturn(0);
128 }
129 EXTERN_C_END
130