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