Actual source code: euler.c
petsc-3.11.4 2019-09-28
1: /*
2: Code for Timestepping with explicit Euler.
3: */
4: #include <petsc/private/tsimpl.h>
6: typedef struct {
7: Vec update; /* work vector where new solution is formed */
8: } TS_Euler;
10: static PetscErrorCode TSStep_Euler(TS ts)
11: {
12: TS_Euler *euler = (TS_Euler*)ts->data;
13: Vec solution = ts->vec_sol,update = euler->update;
14: PetscBool stageok,accept = PETSC_TRUE;
15: PetscReal next_time_step = ts->time_step;
19: TSPreStage(ts,ts->ptime);
20: TSComputeRHSFunction(ts,ts->ptime,solution,update);
21: VecAYPX(update,ts->time_step,solution);
22: TSPostStage(ts,ts->ptime,0,&solution);
23: TSAdaptCheckStage(ts->adapt,ts,ts->ptime,solution,&stageok);
24: if(!stageok) {ts->reason = TS_DIVERGED_STEP_REJECTED; return(0);}
25: TSFunctionDomainError(ts,ts->ptime+ts->time_step,update,&stageok);
26: if(!stageok) {ts->reason = TS_DIVERGED_STEP_REJECTED; return(0);}
28: TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);
29: if (!accept) {ts->reason = TS_DIVERGED_STEP_REJECTED; return(0);}
30: VecCopy(update,solution);
32: ts->ptime += ts->time_step;
33: ts->time_step = next_time_step;
34: return(0);
35: }
36: /*------------------------------------------------------------*/
38: static PetscErrorCode TSSetUp_Euler(TS ts)
39: {
40: TS_Euler *euler = (TS_Euler*)ts->data;
44: TSCheckImplicitTerm(ts);
45: VecDuplicate(ts->vec_sol,&euler->update);
46: TSGetAdapt(ts,&ts->adapt);
47: TSAdaptCandidatesClear(ts->adapt);
48: return(0);
49: }
51: static PetscErrorCode TSReset_Euler(TS ts)
52: {
53: TS_Euler *euler = (TS_Euler*)ts->data;
57: VecDestroy(&euler->update);
58: return(0);
59: }
61: static PetscErrorCode TSDestroy_Euler(TS ts)
62: {
66: TSReset_Euler(ts);
67: PetscFree(ts->data);
68: return(0);
69: }
70: /*------------------------------------------------------------*/
72: static PetscErrorCode TSSetFromOptions_Euler(PetscOptionItems *PetscOptionsObject,TS ts)
73: {
75: return(0);
76: }
78: static PetscErrorCode TSView_Euler(TS ts,PetscViewer viewer)
79: {
81: return(0);
82: }
84: static PetscErrorCode TSInterpolate_Euler(TS ts,PetscReal t,Vec X)
85: {
86: TS_Euler *euler = (TS_Euler*)ts->data;
87: Vec update = euler->update;
88: PetscReal alpha = (ts->ptime - t)/ts->time_step;
92: VecWAXPY(X,-ts->time_step,update,ts->vec_sol);
93: VecAXPBY(X,1.0-alpha,alpha,ts->vec_sol);
94: return(0);
95: }
97: static PetscErrorCode TSComputeLinearStability_Euler(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
98: {
100: *yr = 1.0 + xr;
101: *yi = xi;
102: return(0);
103: }
104: /* ------------------------------------------------------------ */
106: /*MC
107: TSEULER - ODE solver using the explicit forward Euler method
109: Level: beginner
111: .seealso: TSCreate(), TS, TSSetType(), TSBEULER
113: M*/
114: PETSC_EXTERN PetscErrorCode TSCreate_Euler(TS ts)
115: {
116: TS_Euler *euler;
120: PetscNewLog(ts,&euler);
121: ts->data = (void*)euler;
123: ts->ops->setup = TSSetUp_Euler;
124: ts->ops->step = TSStep_Euler;
125: ts->ops->reset = TSReset_Euler;
126: ts->ops->destroy = TSDestroy_Euler;
127: ts->ops->setfromoptions = TSSetFromOptions_Euler;
128: ts->ops->view = TSView_Euler;
129: ts->ops->interpolate = TSInterpolate_Euler;
130: ts->ops->linearstability = TSComputeLinearStability_Euler;
131: ts->default_adapt_type = TSADAPTNONE;
132: ts->usessnes = PETSC_FALSE;
133: return(0);
134: }