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 10193ac0bcSJed Brown static PetscErrorCode TSStep_Euler(TS ts) 114b33d51aSBarry Smith { 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; 18*9566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts,ts->ptime)); 19*9566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts,ts->ptime,solution,update)); 20*9566063dSJacob Faibussowitsch PetscCall(VecAYPX(update,ts->time_step,solution)); 21*9566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts,ts->ptime,0,&solution)); 22*9566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(ts->adapt,ts,ts->ptime,solution,&stageok)); 231566a47fSLisandro Dalcin if (!stageok) {ts->reason = TS_DIVERGED_STEP_REJECTED; PetscFunctionReturn(0);} 24*9566063dSJacob Faibussowitsch PetscCall(TSFunctionDomainError(ts,ts->ptime+ts->time_step,update,&stageok)); 251566a47fSLisandro Dalcin if (!stageok) {ts->reason = TS_DIVERGED_STEP_REJECTED; PetscFunctionReturn(0);} 261566a47fSLisandro Dalcin 27*9566063dSJacob Faibussowitsch PetscCall(TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept)); 281566a47fSLisandro Dalcin if (!accept) {ts->reason = TS_DIVERGED_STEP_REJECTED; PetscFunctionReturn(0);} 29*9566063dSJacob Faibussowitsch PetscCall(VecCopy(update,solution)); 301566a47fSLisandro Dalcin 31186e87acSLisandro Dalcin ts->ptime += ts->time_step; 321566a47fSLisandro Dalcin ts->time_step = next_time_step; 333a40ed3dSBarry Smith PetscFunctionReturn(0); 344b33d51aSBarry Smith } 354b33d51aSBarry Smith /*------------------------------------------------------------*/ 36277b19d0SLisandro Dalcin 37277b19d0SLisandro Dalcin static PetscErrorCode TSSetUp_Euler(TS ts) 38277b19d0SLisandro Dalcin { 39277b19d0SLisandro Dalcin TS_Euler *euler = (TS_Euler*)ts->data; 40277b19d0SLisandro Dalcin 41277b19d0SLisandro Dalcin PetscFunctionBegin; 42*9566063dSJacob Faibussowitsch PetscCall(TSCheckImplicitTerm(ts)); 43*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&euler->update)); 44*9566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts,&ts->adapt)); 45*9566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidatesClear(ts->adapt)); 46277b19d0SLisandro Dalcin PetscFunctionReturn(0); 47277b19d0SLisandro Dalcin } 48277b19d0SLisandro Dalcin 49277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Euler(TS ts) 504b33d51aSBarry Smith { 518b1af7b3SBarry Smith TS_Euler *euler = (TS_Euler*)ts->data; 528b1af7b3SBarry Smith 533a40ed3dSBarry Smith PetscFunctionBegin; 54*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&euler->update)); 55277b19d0SLisandro Dalcin PetscFunctionReturn(0); 56277b19d0SLisandro Dalcin } 57277b19d0SLisandro Dalcin 58277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Euler(TS ts) 59277b19d0SLisandro Dalcin { 60277b19d0SLisandro Dalcin PetscFunctionBegin; 61*9566063dSJacob Faibussowitsch PetscCall(TSReset_Euler(ts)); 62*9566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 633a40ed3dSBarry Smith PetscFunctionReturn(0); 648b1af7b3SBarry Smith } 658b1af7b3SBarry Smith /*------------------------------------------------------------*/ 668b1af7b3SBarry Smith 674416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Euler(PetscOptionItems *PetscOptionsObject,TS ts) 688b1af7b3SBarry Smith { 693a40ed3dSBarry Smith PetscFunctionBegin; 703a40ed3dSBarry Smith PetscFunctionReturn(0); 714b33d51aSBarry Smith } 724b33d51aSBarry Smith 736849ba73SBarry Smith static PetscErrorCode TSView_Euler(TS ts,PetscViewer viewer) 748b1af7b3SBarry Smith { 753a40ed3dSBarry Smith PetscFunctionBegin; 763a40ed3dSBarry Smith PetscFunctionReturn(0); 778b1af7b3SBarry Smith } 788b1af7b3SBarry Smith 796083293cSBarry Smith static PetscErrorCode TSInterpolate_Euler(TS ts,PetscReal t,Vec X) 806083293cSBarry Smith { 813354212dSEmil Constantinescu TS_Euler *euler = (TS_Euler*)ts->data; 823354212dSEmil Constantinescu Vec update = euler->update; 836083293cSBarry Smith PetscReal alpha = (ts->ptime - t)/ts->time_step; 846083293cSBarry Smith 856083293cSBarry Smith PetscFunctionBegin; 86*9566063dSJacob Faibussowitsch PetscCall(VecWAXPY(X,-ts->time_step,update,ts->vec_sol)); 87*9566063dSJacob Faibussowitsch PetscCall(VecAXPBY(X,1.0-alpha,alpha,ts->vec_sol)); 886083293cSBarry Smith PetscFunctionReturn(0); 896083293cSBarry Smith } 906083293cSBarry Smith 91560360afSLisandro Dalcin static PetscErrorCode TSComputeLinearStability_Euler(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi) 92f9c1d6abSBarry Smith { 93f9c1d6abSBarry Smith PetscFunctionBegin; 94f9c1d6abSBarry Smith *yr = 1.0 + xr; 95f9c1d6abSBarry Smith *yi = xi; 96f9c1d6abSBarry Smith PetscFunctionReturn(0); 97f9c1d6abSBarry Smith } 988b1af7b3SBarry Smith /* ------------------------------------------------------------ */ 99ebe3b25bSBarry Smith 100ebe3b25bSBarry Smith /*MC 1019596e0b4SJed Brown TSEULER - ODE solver using the explicit forward Euler method 102ebe3b25bSBarry Smith 103d41222bbSBarry Smith Level: beginner 104d41222bbSBarry Smith 1059596e0b4SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSBEULER 106ebe3b25bSBarry Smith 107ebe3b25bSBarry Smith M*/ 1088cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Euler(TS ts) 1098b1af7b3SBarry Smith { 1108b1af7b3SBarry Smith TS_Euler *euler; 1118b1af7b3SBarry Smith 1123a40ed3dSBarry Smith PetscFunctionBegin; 113*9566063dSJacob Faibussowitsch PetscCall(PetscNewLog(ts,&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