Actual source code: euler.c
1: #define PETSCTS_DLL
3: /*
4: Code for Timestepping with explicit Euler.
5: */
6: #include private/tsimpl.h
8: typedef struct {
9: Vec update; /* work vector where F(t[i],u[i]) is stored */
10: } TS_Euler;
14: static PetscErrorCode TSSetUp_Euler(TS ts)
15: {
16: TS_Euler *euler = (TS_Euler*)ts->data;
20: VecDuplicate(ts->vec_sol,&euler->update);
21: return(0);
22: }
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;
31: PetscInt i,max_steps = ts->max_steps;
32:
34: *steps = -ts->steps;
35: TSMonitor(ts,ts->steps,ts->ptime,sol);
37: for (i=0; i<max_steps; i++) {
38: PetscReal dt = ts->time_step;
40: TSPreStep(ts);
41: ts->ptime += dt;
42: TSComputeRHSFunction(ts,ts->ptime,sol,update);
43: VecAXPY(sol,dt,update);
44: ts->steps++;
45: TSPostStep(ts);
46: TSMonitor(ts,ts->steps,ts->ptime,sol);
47: if (ts->ptime > ts->max_time) break;
48: }
50: *steps += ts->steps;
51: *ptime = ts->ptime;
52: return(0);
53: }
54: /*------------------------------------------------------------*/
57: static PetscErrorCode TSDestroy_Euler(TS ts)
58: {
59: TS_Euler *euler = (TS_Euler*)ts->data;
63: if (euler->update) {VecDestroy(euler->update);}
64: PetscFree(euler);
65: return(0);
66: }
67: /*------------------------------------------------------------*/
71: static PetscErrorCode TSSetFromOptions_Euler(TS ts)
72: {
74: return(0);
75: }
79: static PetscErrorCode TSView_Euler(TS ts,PetscViewer viewer)
80: {
82: return(0);
83: }
85: /* ------------------------------------------------------------ */
87: /*MC
88: TSEULER - ODE solver using the explicit forward Euler method
90: Level: beginner
92: .seealso: TSCreate(), TS, TSSetType(), TSBEULER
94: M*/
98: PetscErrorCode TSCreate_Euler(TS ts)
99: {
100: TS_Euler *euler;
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;
110: PetscNewLog(ts,TS_Euler,&euler);
111: ts->data = (void*)euler;
113: return(0);
114: }