xref: /petsc/src/ts/impls/explicit/euler/euler.c (revision da9f1d6b25924a16baf1fafcd5e58fa8eaafd3cf)
1 #define PETSCTS_DLL
2 
3 /*
4        Code for Timestepping with explicit Euler.
5 */
6 #include "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   PetscScalar    dt = ts->time_step;
33 
34   PetscFunctionBegin;
35   *steps = -ts->steps;
36   ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
37 
38   for (i=0; i<max_steps; i++) {
39     ts->ptime += ts->time_step;
40     ierr = TSComputeRHSFunction(ts,ts->ptime,sol,update);CHKERRQ(ierr);
41     ierr = VecAXPY(sol,dt,update);CHKERRQ(ierr);
42     ts->steps++;
43     ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
44     if (ts->ptime > ts->max_time) break;
45   }
46 
47   *steps += ts->steps;
48   *ptime  = ts->ptime;
49   PetscFunctionReturn(0);
50 }
51 /*------------------------------------------------------------*/
52 #undef __FUNCT__
53 #define __FUNCT__ "TSDestroy_Euler"
54 static PetscErrorCode TSDestroy_Euler(TS ts)
55 {
56   TS_Euler       *euler = (TS_Euler*)ts->data;
57   PetscErrorCode ierr;
58 
59   PetscFunctionBegin;
60   if (euler->update) {ierr = VecDestroy(euler->update);CHKERRQ(ierr);}
61   ierr = PetscFree(euler);CHKERRQ(ierr);
62   PetscFunctionReturn(0);
63 }
64 /*------------------------------------------------------------*/
65 
66 #undef __FUNCT__
67 #define __FUNCT__ "TSSetFromOptions_Euler"
68 static PetscErrorCode TSSetFromOptions_Euler(TS ts)
69 {
70   PetscFunctionBegin;
71   PetscFunctionReturn(0);
72 }
73 
74 #undef __FUNCT__
75 #define __FUNCT__ "TSView_Euler"
76 static PetscErrorCode TSView_Euler(TS ts,PetscViewer viewer)
77 {
78   PetscFunctionBegin;
79   PetscFunctionReturn(0);
80 }
81 
82 /* ------------------------------------------------------------ */
83 
84 /*MC
85       TS_EULER - ODE solver using the explicit forward Euler method
86 
87   Level: beginner
88 
89 .seealso:  TSCreate(), TS, TSSetType(), TS_BEULER
90 
91 M*/
92 EXTERN_C_BEGIN
93 #undef __FUNCT__
94 #define __FUNCT__ "TSCreate_Euler"
95 PetscErrorCode PETSCTS_DLLEXPORT TSCreate_Euler(TS ts)
96 {
97   TS_Euler       *euler;
98   PetscErrorCode ierr;
99 
100   PetscFunctionBegin;
101   ts->ops->setup           = TSSetUp_Euler;
102   ts->ops->step            = TSStep_Euler;
103   ts->ops->destroy         = TSDestroy_Euler;
104   ts->ops->setfromoptions  = TSSetFromOptions_Euler;
105   ts->ops->view            = TSView_Euler;
106 
107   ierr = PetscNewLog(ts,TS_Euler,&euler);CHKERRQ(ierr);
108   ts->data = (void*)euler;
109 
110   PetscFunctionReturn(0);
111 }
112 EXTERN_C_END
113 
114 
115 
116 
117