14b33d51aSBarry Smith /*
2fb4a63b6SLois Curfman McInnes Code for Timestepping with explicit Euler.
34b33d51aSBarry Smith */
4af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/
54b33d51aSBarry Smith
68b1af7b3SBarry Smith typedef struct {
7277b19d0SLisandro Dalcin Vec update; /* work vector where new solution is formed */
88b1af7b3SBarry Smith } TS_Euler;
98b1af7b3SBarry Smith
TSStep_Euler(TS ts)10d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_Euler(TS ts)
11d71ae5a4SJacob Faibussowitsch {
128b1af7b3SBarry Smith TS_Euler *euler = (TS_Euler *)ts->data;
131566a47fSLisandro Dalcin Vec solution = ts->vec_sol, update = euler->update;
141566a47fSLisandro Dalcin PetscBool stageok, accept = PETSC_TRUE;
151566a47fSLisandro Dalcin PetscReal next_time_step = ts->time_step;
164b33d51aSBarry Smith
173a40ed3dSBarry Smith PetscFunctionBegin;
189566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, ts->ptime));
199566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts, ts->ptime, solution, update));
209566063dSJacob Faibussowitsch PetscCall(VecAYPX(update, ts->time_step, solution));
219566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, ts->ptime, 0, &solution));
229566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(ts->adapt, ts, ts->ptime, solution, &stageok));
239371c9d4SSatish Balay if (!stageok) {
249371c9d4SSatish Balay ts->reason = TS_DIVERGED_STEP_REJECTED;
253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
269371c9d4SSatish Balay }
279566063dSJacob Faibussowitsch PetscCall(TSFunctionDomainError(ts, ts->ptime + ts->time_step, update, &stageok));
289371c9d4SSatish Balay if (!stageok) {
299371c9d4SSatish Balay ts->reason = TS_DIVERGED_STEP_REJECTED;
303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
319371c9d4SSatish Balay }
321566a47fSLisandro Dalcin
339566063dSJacob Faibussowitsch PetscCall(TSAdaptChoose(ts->adapt, ts, ts->time_step, NULL, &next_time_step, &accept));
349371c9d4SSatish Balay if (!accept) {
359371c9d4SSatish Balay ts->reason = TS_DIVERGED_STEP_REJECTED;
363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
379371c9d4SSatish Balay }
389566063dSJacob Faibussowitsch PetscCall(VecCopy(update, solution));
391566a47fSLisandro Dalcin
40186e87acSLisandro Dalcin ts->ptime += ts->time_step;
411566a47fSLisandro Dalcin ts->time_step = next_time_step;
423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
434b33d51aSBarry Smith }
44277b19d0SLisandro Dalcin
TSSetUp_Euler(TS ts)45d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_Euler(TS ts)
46d71ae5a4SJacob Faibussowitsch {
47277b19d0SLisandro Dalcin TS_Euler *euler = (TS_Euler *)ts->data;
48277b19d0SLisandro Dalcin
49277b19d0SLisandro Dalcin PetscFunctionBegin;
509566063dSJacob Faibussowitsch PetscCall(TSCheckImplicitTerm(ts));
519566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &euler->update));
529566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &ts->adapt));
539566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidatesClear(ts->adapt));
543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
55277b19d0SLisandro Dalcin }
56277b19d0SLisandro Dalcin
TSReset_Euler(TS ts)57d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_Euler(TS ts)
58d71ae5a4SJacob Faibussowitsch {
598b1af7b3SBarry Smith TS_Euler *euler = (TS_Euler *)ts->data;
608b1af7b3SBarry Smith
613a40ed3dSBarry Smith PetscFunctionBegin;
629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&euler->update));
633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
64277b19d0SLisandro Dalcin }
65277b19d0SLisandro Dalcin
TSDestroy_Euler(TS ts)66d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_Euler(TS ts)
67d71ae5a4SJacob Faibussowitsch {
68277b19d0SLisandro Dalcin PetscFunctionBegin;
699566063dSJacob Faibussowitsch PetscCall(TSReset_Euler(ts));
709566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data));
713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
728b1af7b3SBarry Smith }
738b1af7b3SBarry Smith
TSSetFromOptions_Euler(TS ts,PetscOptionItems PetscOptionsObject)74*ce78bad3SBarry Smith static PetscErrorCode TSSetFromOptions_Euler(TS ts, PetscOptionItems PetscOptionsObject)
75d71ae5a4SJacob Faibussowitsch {
763a40ed3dSBarry Smith PetscFunctionBegin;
773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
784b33d51aSBarry Smith }
794b33d51aSBarry Smith
TSView_Euler(TS ts,PetscViewer viewer)80d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_Euler(TS ts, PetscViewer viewer)
81d71ae5a4SJacob Faibussowitsch {
823a40ed3dSBarry Smith PetscFunctionBegin;
833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
848b1af7b3SBarry Smith }
858b1af7b3SBarry Smith
TSInterpolate_Euler(TS ts,PetscReal t,Vec X)86d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_Euler(TS ts, PetscReal t, Vec X)
87d71ae5a4SJacob Faibussowitsch {
883354212dSEmil Constantinescu TS_Euler *euler = (TS_Euler *)ts->data;
893354212dSEmil Constantinescu Vec update = euler->update;
906083293cSBarry Smith PetscReal alpha = (ts->ptime - t) / ts->time_step;
916083293cSBarry Smith
926083293cSBarry Smith PetscFunctionBegin;
939566063dSJacob Faibussowitsch PetscCall(VecWAXPY(X, -ts->time_step, update, ts->vec_sol));
949566063dSJacob Faibussowitsch PetscCall(VecAXPBY(X, 1.0 - alpha, alpha, ts->vec_sol));
953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
966083293cSBarry Smith }
976083293cSBarry Smith
TSComputeLinearStability_Euler(TS ts,PetscReal xr,PetscReal xi,PetscReal * yr,PetscReal * yi)98d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSComputeLinearStability_Euler(TS ts, PetscReal xr, PetscReal xi, PetscReal *yr, PetscReal *yi)
99d71ae5a4SJacob Faibussowitsch {
100f9c1d6abSBarry Smith PetscFunctionBegin;
101f9c1d6abSBarry Smith *yr = 1.0 + xr;
102f9c1d6abSBarry Smith *yi = xi;
1033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
104f9c1d6abSBarry Smith }
105ebe3b25bSBarry Smith
106ebe3b25bSBarry Smith /*MC
1079596e0b4SJed Brown TSEULER - ODE solver using the explicit forward Euler method
108ebe3b25bSBarry Smith
109d41222bbSBarry Smith Level: beginner
110d41222bbSBarry Smith
1111cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSBEULER`, `TSType`
112ebe3b25bSBarry Smith M*/
TSCreate_Euler(TS ts)113d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_Euler(TS ts)
114d71ae5a4SJacob Faibussowitsch {
1158b1af7b3SBarry Smith TS_Euler *euler;
1168b1af7b3SBarry Smith
1173a40ed3dSBarry Smith PetscFunctionBegin;
1184dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&euler));
1191566a47fSLisandro Dalcin ts->data = (void *)euler;
1201566a47fSLisandro Dalcin
121000e7ae3SMatthew Knepley ts->ops->setup = TSSetUp_Euler;
122000e7ae3SMatthew Knepley ts->ops->step = TSStep_Euler;
123277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Euler;
124000e7ae3SMatthew Knepley ts->ops->destroy = TSDestroy_Euler;
125000e7ae3SMatthew Knepley ts->ops->setfromoptions = TSSetFromOptions_Euler;
126000e7ae3SMatthew Knepley ts->ops->view = TSView_Euler;
1276083293cSBarry Smith ts->ops->interpolate = TSInterpolate_Euler;
128f9c1d6abSBarry Smith ts->ops->linearstability = TSComputeLinearStability_Euler;
1292ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE;
130efd4aadfSBarry Smith ts->usessnes = PETSC_FALSE;
1313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1324b33d51aSBarry Smith }
133