Actual source code: euler.c
petsc-3.7.3 2016-08-01
1: /*
2: Code for Timestepping with explicit Euler.
3: */
4: #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/
6: typedef struct {
7: Vec update; /* work vector where new solution is formed */
8: } TS_Euler;
12: static PetscErrorCode TSStep_Euler(TS ts)
13: {
14: TS_Euler *euler = (TS_Euler*)ts->data;
15: Vec solution = ts->vec_sol,update = euler->update;
16: PetscBool stageok,accept = PETSC_TRUE;
17: PetscReal next_time_step = ts->time_step;
21: TSPreStage(ts,ts->ptime);
22: TSComputeRHSFunction(ts,ts->ptime,solution,update);
23: VecAYPX(update,ts->time_step,solution);
24: TSPostStage(ts,ts->ptime,0,&solution);
25: TSAdaptCheckStage(ts->adapt,ts,ts->ptime,solution,&stageok);
26: if(!stageok) {ts->reason = TS_DIVERGED_STEP_REJECTED; return(0);}
27: TSFunctionDomainError(ts,ts->ptime+ts->time_step,update,&stageok);
28: if(!stageok) {ts->reason = TS_DIVERGED_STEP_REJECTED; return(0);}
30: TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);
31: if (!accept) {ts->reason = TS_DIVERGED_STEP_REJECTED; return(0);}
32: VecCopy(update,solution);
34: ts->ptime += ts->time_step;
35: ts->time_step = next_time_step;
36: return(0);
37: }
38: /*------------------------------------------------------------*/
42: static PetscErrorCode TSSetUp_Euler(TS ts)
43: {
44: TS_Euler *euler = (TS_Euler*)ts->data;
46: TSRHSFunction rhsfunction;
47: TSIFunction ifunction;
50: TSGetIFunction(ts,NULL,&ifunction,NULL);
51: TSGetRHSFunction(ts,NULL,&rhsfunction,NULL);
52: if (!rhsfunction || ifunction) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must define RHSFunction() and leave IFunction() undefined in order to use -ts_type euler");
53: VecDuplicate(ts->vec_sol,&euler->update);
55: TSGetAdapt(ts,&ts->adapt);
56: TSAdaptCandidatesClear(ts->adapt);
57: return(0);
58: }
62: static PetscErrorCode TSReset_Euler(TS ts)
63: {
64: TS_Euler *euler = (TS_Euler*)ts->data;
68: VecDestroy(&euler->update);
69: return(0);
70: }
74: static PetscErrorCode TSDestroy_Euler(TS ts)
75: {
79: TSReset_Euler(ts);
80: PetscFree(ts->data);
81: return(0);
82: }
83: /*------------------------------------------------------------*/
87: static PetscErrorCode TSSetFromOptions_Euler(PetscOptionItems *PetscOptionsObject,TS ts)
88: {
90: return(0);
91: }
95: static PetscErrorCode TSView_Euler(TS ts,PetscViewer viewer)
96: {
98: return(0);
99: }
103: static PetscErrorCode TSInterpolate_Euler(TS ts,PetscReal t,Vec X)
104: {
105: TS_Euler *euler = (TS_Euler*)ts->data;
106: Vec update = euler->update;
107: PetscReal alpha = (ts->ptime - t)/ts->time_step;
111: VecWAXPY(X,-ts->time_step,update,ts->vec_sol);
112: VecAXPBY(X,1.0-alpha,alpha,ts->vec_sol);
113: return(0);
114: }
118: static PetscErrorCode TSComputeLinearStability_Euler(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
119: {
121: *yr = 1.0 + xr;
122: *yi = xi;
123: return(0);
124: }
125: /* ------------------------------------------------------------ */
127: /*MC
128: TSEULER - ODE solver using the explicit forward Euler method
130: Level: beginner
132: .seealso: TSCreate(), TS, TSSetType(), TSBEULER
134: M*/
137: PETSC_EXTERN PetscErrorCode TSCreate_Euler(TS ts)
138: {
139: TS_Euler *euler;
143: PetscNewLog(ts,&euler);
144: ts->data = (void*)euler;
146: ts->ops->setup = TSSetUp_Euler;
147: ts->ops->step = TSStep_Euler;
148: ts->ops->reset = TSReset_Euler;
149: ts->ops->destroy = TSDestroy_Euler;
150: ts->ops->setfromoptions = TSSetFromOptions_Euler;
151: ts->ops->view = TSView_Euler;
152: ts->ops->interpolate = TSInterpolate_Euler;
153: ts->ops->linearstability = TSComputeLinearStability_Euler;
154: return(0);
155: }