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 109371c9d4SSatish Balay static PetscErrorCode TSStep_Euler(TS ts) { 118b1af7b3SBarry Smith TS_Euler *euler = (TS_Euler *)ts->data; 121566a47fSLisandro Dalcin Vec solution = ts->vec_sol, update = euler->update; 131566a47fSLisandro Dalcin PetscBool stageok, accept = PETSC_TRUE; 141566a47fSLisandro Dalcin PetscReal next_time_step = ts->time_step; 154b33d51aSBarry Smith 163a40ed3dSBarry Smith PetscFunctionBegin; 179566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, ts->ptime)); 189566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts, ts->ptime, solution, update)); 199566063dSJacob Faibussowitsch PetscCall(VecAYPX(update, ts->time_step, solution)); 209566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, ts->ptime, 0, &solution)); 219566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(ts->adapt, ts, ts->ptime, solution, &stageok)); 229371c9d4SSatish Balay if (!stageok) { 239371c9d4SSatish Balay ts->reason = TS_DIVERGED_STEP_REJECTED; 249371c9d4SSatish Balay PetscFunctionReturn(0); 259371c9d4SSatish Balay } 269566063dSJacob Faibussowitsch PetscCall(TSFunctionDomainError(ts, ts->ptime + ts->time_step, update, &stageok)); 279371c9d4SSatish Balay if (!stageok) { 289371c9d4SSatish Balay ts->reason = TS_DIVERGED_STEP_REJECTED; 299371c9d4SSatish Balay PetscFunctionReturn(0); 309371c9d4SSatish Balay } 311566a47fSLisandro Dalcin 329566063dSJacob Faibussowitsch PetscCall(TSAdaptChoose(ts->adapt, ts, ts->time_step, NULL, &next_time_step, &accept)); 339371c9d4SSatish Balay if (!accept) { 349371c9d4SSatish Balay ts->reason = TS_DIVERGED_STEP_REJECTED; 359371c9d4SSatish Balay PetscFunctionReturn(0); 369371c9d4SSatish Balay } 379566063dSJacob Faibussowitsch PetscCall(VecCopy(update, solution)); 381566a47fSLisandro Dalcin 39186e87acSLisandro Dalcin ts->ptime += ts->time_step; 401566a47fSLisandro Dalcin ts->time_step = next_time_step; 413a40ed3dSBarry Smith PetscFunctionReturn(0); 424b33d51aSBarry Smith } 434b33d51aSBarry Smith /*------------------------------------------------------------*/ 44277b19d0SLisandro Dalcin 459371c9d4SSatish Balay static PetscErrorCode TSSetUp_Euler(TS ts) { 46277b19d0SLisandro Dalcin TS_Euler *euler = (TS_Euler *)ts->data; 47277b19d0SLisandro Dalcin 48277b19d0SLisandro Dalcin PetscFunctionBegin; 499566063dSJacob Faibussowitsch PetscCall(TSCheckImplicitTerm(ts)); 509566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &euler->update)); 519566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &ts->adapt)); 529566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidatesClear(ts->adapt)); 53277b19d0SLisandro Dalcin PetscFunctionReturn(0); 54277b19d0SLisandro Dalcin } 55277b19d0SLisandro Dalcin 569371c9d4SSatish Balay static PetscErrorCode TSReset_Euler(TS ts) { 578b1af7b3SBarry Smith TS_Euler *euler = (TS_Euler *)ts->data; 588b1af7b3SBarry Smith 593a40ed3dSBarry Smith PetscFunctionBegin; 609566063dSJacob Faibussowitsch PetscCall(VecDestroy(&euler->update)); 61277b19d0SLisandro Dalcin PetscFunctionReturn(0); 62277b19d0SLisandro Dalcin } 63277b19d0SLisandro Dalcin 649371c9d4SSatish Balay static PetscErrorCode TSDestroy_Euler(TS ts) { 65277b19d0SLisandro Dalcin PetscFunctionBegin; 669566063dSJacob Faibussowitsch PetscCall(TSReset_Euler(ts)); 679566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 683a40ed3dSBarry Smith PetscFunctionReturn(0); 698b1af7b3SBarry Smith } 708b1af7b3SBarry Smith /*------------------------------------------------------------*/ 718b1af7b3SBarry Smith 729371c9d4SSatish Balay static PetscErrorCode TSSetFromOptions_Euler(TS ts, PetscOptionItems *PetscOptionsObject) { 733a40ed3dSBarry Smith PetscFunctionBegin; 743a40ed3dSBarry Smith PetscFunctionReturn(0); 754b33d51aSBarry Smith } 764b33d51aSBarry Smith 779371c9d4SSatish Balay static PetscErrorCode TSView_Euler(TS ts, PetscViewer viewer) { 783a40ed3dSBarry Smith PetscFunctionBegin; 793a40ed3dSBarry Smith PetscFunctionReturn(0); 808b1af7b3SBarry Smith } 818b1af7b3SBarry Smith 829371c9d4SSatish Balay static PetscErrorCode TSInterpolate_Euler(TS ts, PetscReal t, Vec X) { 833354212dSEmil Constantinescu TS_Euler *euler = (TS_Euler *)ts->data; 843354212dSEmil Constantinescu Vec update = euler->update; 856083293cSBarry Smith PetscReal alpha = (ts->ptime - t) / ts->time_step; 866083293cSBarry Smith 876083293cSBarry Smith PetscFunctionBegin; 889566063dSJacob Faibussowitsch PetscCall(VecWAXPY(X, -ts->time_step, update, ts->vec_sol)); 899566063dSJacob Faibussowitsch PetscCall(VecAXPBY(X, 1.0 - alpha, alpha, ts->vec_sol)); 906083293cSBarry Smith PetscFunctionReturn(0); 916083293cSBarry Smith } 926083293cSBarry Smith 939371c9d4SSatish Balay static PetscErrorCode TSComputeLinearStability_Euler(TS ts, PetscReal xr, PetscReal xi, PetscReal *yr, PetscReal *yi) { 94f9c1d6abSBarry Smith PetscFunctionBegin; 95f9c1d6abSBarry Smith *yr = 1.0 + xr; 96f9c1d6abSBarry Smith *yi = xi; 97f9c1d6abSBarry Smith PetscFunctionReturn(0); 98f9c1d6abSBarry Smith } 998b1af7b3SBarry Smith /* ------------------------------------------------------------ */ 100ebe3b25bSBarry Smith 101ebe3b25bSBarry Smith /*MC 1029596e0b4SJed Brown TSEULER - ODE solver using the explicit forward Euler method 103ebe3b25bSBarry Smith 104d41222bbSBarry Smith Level: beginner 105d41222bbSBarry Smith 106db781477SPatrick Sanan .seealso: `TSCreate()`, `TS`, `TSSetType()`, `TSBEULER` 107ebe3b25bSBarry Smith 108ebe3b25bSBarry Smith M*/ 1099371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TSCreate_Euler(TS ts) { 1108b1af7b3SBarry Smith TS_Euler *euler; 1118b1af7b3SBarry Smith 1123a40ed3dSBarry Smith PetscFunctionBegin; 113*4dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&euler)); 1141566a47fSLisandro Dalcin ts->data = (void *)euler; 1151566a47fSLisandro Dalcin 116000e7ae3SMatthew Knepley ts->ops->setup = TSSetUp_Euler; 117000e7ae3SMatthew Knepley ts->ops->step = TSStep_Euler; 118277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Euler; 119000e7ae3SMatthew Knepley ts->ops->destroy = TSDestroy_Euler; 120000e7ae3SMatthew Knepley ts->ops->setfromoptions = TSSetFromOptions_Euler; 121000e7ae3SMatthew Knepley ts->ops->view = TSView_Euler; 1226083293cSBarry Smith ts->ops->interpolate = TSInterpolate_Euler; 123f9c1d6abSBarry Smith ts->ops->linearstability = TSComputeLinearStability_Euler; 1242ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 125efd4aadfSBarry Smith ts->usessnes = PETSC_FALSE; 1263a40ed3dSBarry Smith PetscFunctionReturn(0); 1274b33d51aSBarry Smith } 128