xref: /petsc/src/ts/impls/explicit/euler/euler.c (revision 193ac0bc5128976c4ec2c1a67f3f9cb026b77f22)
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)
13 {
14   TS_Euler       *euler = (TS_Euler*)ts->data;
15   Vec            sol = ts->vec_sol,update = euler->update;
16   PetscErrorCode ierr;
17 
18   PetscFunctionBegin;
19   ierr = TSComputeRHSFunction(ts,ts->ptime,sol,update);CHKERRQ(ierr);
20   ierr = VecAXPY(sol,ts->time_step,update);CHKERRQ(ierr);
21   ts->ptime += ts->time_step;
22   ts->steps++;
23   PetscFunctionReturn(0);
24 }
25 /*------------------------------------------------------------*/
26 
27 #undef __FUNCT__
28 #define __FUNCT__ "TSSetUp_Euler"
29 static PetscErrorCode TSSetUp_Euler(TS ts)
30 {
31   TS_Euler       *euler = (TS_Euler*)ts->data;
32   PetscErrorCode ierr;
33 
34   PetscFunctionBegin;
35   ierr = VecDuplicate(ts->vec_sol,&euler->update);CHKERRQ(ierr);
36   PetscFunctionReturn(0);
37 }
38 
39 #undef __FUNCT__
40 #define __FUNCT__ "TSReset_Euler"
41 static PetscErrorCode TSReset_Euler(TS ts)
42 {
43   TS_Euler       *euler = (TS_Euler*)ts->data;
44   PetscErrorCode ierr;
45 
46   PetscFunctionBegin;
47   ierr = VecDestroy(&euler->update);CHKERRQ(ierr);
48   PetscFunctionReturn(0);
49 }
50 
51 #undef __FUNCT__
52 #define __FUNCT__ "TSDestroy_Euler"
53 static PetscErrorCode TSDestroy_Euler(TS ts)
54 {
55   PetscErrorCode ierr;
56 
57   PetscFunctionBegin;
58   ierr = TSReset_Euler(ts);CHKERRQ(ierr);
59   ierr = PetscFree(ts->data);CHKERRQ(ierr);
60   PetscFunctionReturn(0);
61 }
62 /*------------------------------------------------------------*/
63 
64 #undef __FUNCT__
65 #define __FUNCT__ "TSSetFromOptions_Euler"
66 static PetscErrorCode TSSetFromOptions_Euler(TS ts)
67 {
68   PetscFunctionBegin;
69   PetscFunctionReturn(0);
70 }
71 
72 #undef __FUNCT__
73 #define __FUNCT__ "TSView_Euler"
74 static PetscErrorCode TSView_Euler(TS ts,PetscViewer viewer)
75 {
76   PetscFunctionBegin;
77   PetscFunctionReturn(0);
78 }
79 
80 /* ------------------------------------------------------------ */
81 
82 /*MC
83       TSEULER - ODE solver using the explicit forward Euler method
84 
85   Level: beginner
86 
87 .seealso:  TSCreate(), TS, TSSetType(), TSBEULER
88 
89 M*/
90 EXTERN_C_BEGIN
91 #undef __FUNCT__
92 #define __FUNCT__ "TSCreate_Euler"
93 PetscErrorCode  TSCreate_Euler(TS ts)
94 {
95   TS_Euler       *euler;
96   PetscErrorCode ierr;
97 
98   PetscFunctionBegin;
99   ts->ops->setup           = TSSetUp_Euler;
100   ts->ops->step            = TSStep_Euler;
101   ts->ops->reset           = TSReset_Euler;
102   ts->ops->destroy         = TSDestroy_Euler;
103   ts->ops->setfromoptions  = TSSetFromOptions_Euler;
104   ts->ops->view            = TSView_Euler;
105 
106   ierr = PetscNewLog(ts,TS_Euler,&euler);CHKERRQ(ierr);
107   ts->data = (void*)euler;
108 
109   PetscFunctionReturn(0);
110 }
111 EXTERN_C_END
112