xref: /petsc/src/ts/impls/explicit/euler/euler.c (revision f44c6289809fde6fa30cd063773daee08b23f33b)
1 /*
2        Code for Timestepping with explicit Euler.
3 */
4 #include <petsc/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   PetscBool      accept;
18 
19   PetscFunctionBegin;
20   ierr = TSPreStep(ts);CHKERRQ(ierr);
21   ierr = TSPreStage(ts,ts->ptime);CHKERRQ(ierr);
22   ierr = TSComputeRHSFunction(ts,ts->ptime,sol,update);CHKERRQ(ierr);
23   ierr = VecAXPY(sol,ts->time_step,update);CHKERRQ(ierr);
24   ierr = TSFunctionDomainError(ts,ts->ptime,sol,&accept);CHKERRQ(ierr);
25   if(!accept) {
26     ts->reason = TS_DIVERGED_STEP_REJECTED;
27     PetscFunctionReturn(0);
28   }
29   ierr = TSPostStage(ts,ts->ptime,0,&sol);CHKERRQ(ierr);
30   ts->ptime += ts->time_step;
31   ts->steps++;
32   PetscFunctionReturn(0);
33 }
34 /*------------------------------------------------------------*/
35 
36 #undef __FUNCT__
37 #define __FUNCT__ "TSSetUp_Euler"
38 static PetscErrorCode TSSetUp_Euler(TS ts)
39 {
40   TS_Euler       *euler = (TS_Euler*)ts->data;
41   PetscErrorCode ierr;
42   TSRHSFunction  rhsfunction;
43   TSIFunction    ifunction;
44 
45   PetscFunctionBegin;
46   ierr =  TSGetIFunction(ts,NULL,&ifunction,NULL);CHKERRQ(ierr);
47   ierr =  TSGetRHSFunction(ts,NULL,&rhsfunction,NULL);CHKERRQ(ierr);
48   if (!rhsfunction || ifunction) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must define RHSFunction() and leave IFunction() undefined in order to use -ts_type euler");
49   ierr = VecDuplicate(ts->vec_sol,&euler->update);CHKERRQ(ierr);
50   PetscFunctionReturn(0);
51 }
52 
53 #undef __FUNCT__
54 #define __FUNCT__ "TSReset_Euler"
55 static PetscErrorCode TSReset_Euler(TS ts)
56 {
57   TS_Euler       *euler = (TS_Euler*)ts->data;
58   PetscErrorCode ierr;
59 
60   PetscFunctionBegin;
61   ierr = VecDestroy(&euler->update);CHKERRQ(ierr);
62   PetscFunctionReturn(0);
63 }
64 
65 #undef __FUNCT__
66 #define __FUNCT__ "TSDestroy_Euler"
67 static PetscErrorCode TSDestroy_Euler(TS ts)
68 {
69   PetscErrorCode ierr;
70 
71   PetscFunctionBegin;
72   ierr = TSReset_Euler(ts);CHKERRQ(ierr);
73   ierr = PetscFree(ts->data);CHKERRQ(ierr);
74   PetscFunctionReturn(0);
75 }
76 /*------------------------------------------------------------*/
77 
78 #undef __FUNCT__
79 #define __FUNCT__ "TSSetFromOptions_Euler"
80 static PetscErrorCode TSSetFromOptions_Euler(PetscOptionItems *PetscOptionsObject,TS ts)
81 {
82   PetscFunctionBegin;
83   PetscFunctionReturn(0);
84 }
85 
86 #undef __FUNCT__
87 #define __FUNCT__ "TSView_Euler"
88 static PetscErrorCode TSView_Euler(TS ts,PetscViewer viewer)
89 {
90   PetscFunctionBegin;
91   PetscFunctionReturn(0);
92 }
93 
94 #undef __FUNCT__
95 #define __FUNCT__ "TSInterpolate_Euler"
96 static PetscErrorCode TSInterpolate_Euler(TS ts,PetscReal t,Vec X)
97 {
98   PetscReal      alpha = (ts->ptime - t)/ts->time_step;
99   PetscErrorCode ierr;
100 
101   PetscFunctionBegin;
102   ierr = VecAXPBY(ts->vec_sol,1.0-alpha,alpha,X);CHKERRQ(ierr);
103   PetscFunctionReturn(0);
104 }
105 
106 #undef __FUNCT__
107 #define __FUNCT__ "TSComputeLinearStability_Euler"
108 PetscErrorCode TSComputeLinearStability_Euler(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
109 {
110   PetscFunctionBegin;
111   *yr = 1.0 + xr;
112   *yi = xi;
113   PetscFunctionReturn(0);
114 }
115 /* ------------------------------------------------------------ */
116 
117 /*MC
118       TSEULER - ODE solver using the explicit forward Euler method
119 
120   Level: beginner
121 
122 .seealso:  TSCreate(), TS, TSSetType(), TSBEULER
123 
124 M*/
125 #undef __FUNCT__
126 #define __FUNCT__ "TSCreate_Euler"
127 PETSC_EXTERN PetscErrorCode TSCreate_Euler(TS ts)
128 {
129   TS_Euler       *euler;
130   PetscErrorCode ierr;
131 
132   PetscFunctionBegin;
133   ts->ops->setup           = TSSetUp_Euler;
134   ts->ops->step            = TSStep_Euler;
135   ts->ops->reset           = TSReset_Euler;
136   ts->ops->destroy         = TSDestroy_Euler;
137   ts->ops->setfromoptions  = TSSetFromOptions_Euler;
138   ts->ops->view            = TSView_Euler;
139   ts->ops->interpolate     = TSInterpolate_Euler;
140   ts->ops->linearstability = TSComputeLinearStability_Euler;
141 
142   ierr = PetscNewLog(ts,&euler);CHKERRQ(ierr);
143   ts->data = (void*)euler;
144   PetscFunctionReturn(0);
145 }
146