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 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 } 444b33d51aSBarry Smith /*------------------------------------------------------------*/ 45277b19d0SLisandro Dalcin 46d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_Euler(TS ts) 47d71ae5a4SJacob Faibussowitsch { 48277b19d0SLisandro Dalcin TS_Euler *euler = (TS_Euler *)ts->data; 49277b19d0SLisandro Dalcin 50277b19d0SLisandro Dalcin PetscFunctionBegin; 519566063dSJacob Faibussowitsch PetscCall(TSCheckImplicitTerm(ts)); 529566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &euler->update)); 539566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &ts->adapt)); 549566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidatesClear(ts->adapt)); 553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 56277b19d0SLisandro Dalcin } 57277b19d0SLisandro Dalcin 58d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_Euler(TS ts) 59d71ae5a4SJacob Faibussowitsch { 608b1af7b3SBarry Smith TS_Euler *euler = (TS_Euler *)ts->data; 618b1af7b3SBarry Smith 623a40ed3dSBarry Smith PetscFunctionBegin; 639566063dSJacob Faibussowitsch PetscCall(VecDestroy(&euler->update)); 643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 65277b19d0SLisandro Dalcin } 66277b19d0SLisandro Dalcin 67d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_Euler(TS ts) 68d71ae5a4SJacob Faibussowitsch { 69277b19d0SLisandro Dalcin PetscFunctionBegin; 709566063dSJacob Faibussowitsch PetscCall(TSReset_Euler(ts)); 719566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 738b1af7b3SBarry Smith } 748b1af7b3SBarry Smith /*------------------------------------------------------------*/ 758b1af7b3SBarry Smith 76d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetFromOptions_Euler(TS ts, PetscOptionItems *PetscOptionsObject) 77d71ae5a4SJacob Faibussowitsch { 783a40ed3dSBarry Smith PetscFunctionBegin; 793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 804b33d51aSBarry Smith } 814b33d51aSBarry Smith 82d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_Euler(TS ts, PetscViewer viewer) 83d71ae5a4SJacob Faibussowitsch { 843a40ed3dSBarry Smith PetscFunctionBegin; 853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 868b1af7b3SBarry Smith } 878b1af7b3SBarry Smith 88d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_Euler(TS ts, PetscReal t, Vec X) 89d71ae5a4SJacob Faibussowitsch { 903354212dSEmil Constantinescu TS_Euler *euler = (TS_Euler *)ts->data; 913354212dSEmil Constantinescu Vec update = euler->update; 926083293cSBarry Smith PetscReal alpha = (ts->ptime - t) / ts->time_step; 936083293cSBarry Smith 946083293cSBarry Smith PetscFunctionBegin; 959566063dSJacob Faibussowitsch PetscCall(VecWAXPY(X, -ts->time_step, update, ts->vec_sol)); 969566063dSJacob Faibussowitsch PetscCall(VecAXPBY(X, 1.0 - alpha, alpha, ts->vec_sol)); 973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 986083293cSBarry Smith } 996083293cSBarry Smith 100d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSComputeLinearStability_Euler(TS ts, PetscReal xr, PetscReal xi, PetscReal *yr, PetscReal *yi) 101d71ae5a4SJacob Faibussowitsch { 102f9c1d6abSBarry Smith PetscFunctionBegin; 103f9c1d6abSBarry Smith *yr = 1.0 + xr; 104f9c1d6abSBarry Smith *yi = xi; 1053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 106f9c1d6abSBarry Smith } 1078b1af7b3SBarry Smith /* ------------------------------------------------------------ */ 108ebe3b25bSBarry Smith 109ebe3b25bSBarry Smith /*MC 1109596e0b4SJed Brown TSEULER - ODE solver using the explicit forward Euler method 111ebe3b25bSBarry Smith 112d41222bbSBarry Smith Level: beginner 113d41222bbSBarry Smith 114*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSBEULER`, `TSType` 115ebe3b25bSBarry Smith M*/ 116d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_Euler(TS ts) 117d71ae5a4SJacob Faibussowitsch { 1188b1af7b3SBarry Smith TS_Euler *euler; 1198b1af7b3SBarry Smith 1203a40ed3dSBarry Smith PetscFunctionBegin; 1214dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&euler)); 1221566a47fSLisandro Dalcin ts->data = (void *)euler; 1231566a47fSLisandro Dalcin 124000e7ae3SMatthew Knepley ts->ops->setup = TSSetUp_Euler; 125000e7ae3SMatthew Knepley ts->ops->step = TSStep_Euler; 126277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Euler; 127000e7ae3SMatthew Knepley ts->ops->destroy = TSDestroy_Euler; 128000e7ae3SMatthew Knepley ts->ops->setfromoptions = TSSetFromOptions_Euler; 129000e7ae3SMatthew Knepley ts->ops->view = TSView_Euler; 1306083293cSBarry Smith ts->ops->interpolate = TSInterpolate_Euler; 131f9c1d6abSBarry Smith ts->ops->linearstability = TSComputeLinearStability_Euler; 1322ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 133efd4aadfSBarry Smith ts->usessnes = PETSC_FALSE; 1343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1354b33d51aSBarry Smith } 136