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; 16277b19d0SLisandro Dalcin PetscErrorCode ierr; 174b33d51aSBarry Smith 183a40ed3dSBarry Smith PetscFunctionBegin; 19b8123daeSJed Brown ierr = TSPreStage(ts,ts->ptime);CHKERRQ(ierr); 201566a47fSLisandro Dalcin ierr = TSComputeRHSFunction(ts,ts->ptime,solution,update);CHKERRQ(ierr); 211566a47fSLisandro Dalcin ierr = VecAYPX(update,ts->time_step,solution);CHKERRQ(ierr); 221566a47fSLisandro Dalcin ierr = TSPostStage(ts,ts->ptime,0,&solution);CHKERRQ(ierr); 231566a47fSLisandro Dalcin ierr = TSAdaptCheckStage(ts->adapt,ts,ts->ptime,solution,&stageok);CHKERRQ(ierr); 241566a47fSLisandro Dalcin if(!stageok) {ts->reason = TS_DIVERGED_STEP_REJECTED; PetscFunctionReturn(0);} 251566a47fSLisandro Dalcin ierr = TSFunctionDomainError(ts,ts->ptime+ts->time_step,update,&stageok);CHKERRQ(ierr); 261566a47fSLisandro Dalcin if(!stageok) {ts->reason = TS_DIVERGED_STEP_REJECTED; PetscFunctionReturn(0);} 271566a47fSLisandro Dalcin 281566a47fSLisandro Dalcin ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 291566a47fSLisandro Dalcin if (!accept) {ts->reason = TS_DIVERGED_STEP_REJECTED; PetscFunctionReturn(0);} 301566a47fSLisandro Dalcin ierr = VecCopy(update,solution);CHKERRQ(ierr); 311566a47fSLisandro Dalcin 32186e87acSLisandro Dalcin ts->ptime += ts->time_step; 331566a47fSLisandro Dalcin ts->time_step = next_time_step; 343a40ed3dSBarry Smith PetscFunctionReturn(0); 354b33d51aSBarry Smith } 364b33d51aSBarry Smith /*------------------------------------------------------------*/ 37277b19d0SLisandro Dalcin 38277b19d0SLisandro Dalcin static PetscErrorCode TSSetUp_Euler(TS ts) 39277b19d0SLisandro Dalcin { 40277b19d0SLisandro Dalcin TS_Euler *euler = (TS_Euler*)ts->data; 41277b19d0SLisandro Dalcin PetscErrorCode ierr; 42277b19d0SLisandro Dalcin 43277b19d0SLisandro Dalcin PetscFunctionBegin; 443dd83b38SBarry Smith ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr); 45277b19d0SLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&euler->update);CHKERRQ(ierr); 461566a47fSLisandro Dalcin ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr); 471566a47fSLisandro Dalcin ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr); 48277b19d0SLisandro Dalcin PetscFunctionReturn(0); 49277b19d0SLisandro Dalcin } 50277b19d0SLisandro Dalcin 51277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Euler(TS ts) 524b33d51aSBarry Smith { 538b1af7b3SBarry Smith TS_Euler *euler = (TS_Euler*)ts->data; 54dfbe8321SBarry Smith PetscErrorCode ierr; 558b1af7b3SBarry Smith 563a40ed3dSBarry Smith PetscFunctionBegin; 576bf464f9SBarry Smith ierr = VecDestroy(&euler->update);CHKERRQ(ierr); 58277b19d0SLisandro Dalcin PetscFunctionReturn(0); 59277b19d0SLisandro Dalcin } 60277b19d0SLisandro Dalcin 61277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Euler(TS ts) 62277b19d0SLisandro Dalcin { 63277b19d0SLisandro Dalcin PetscErrorCode ierr; 64277b19d0SLisandro Dalcin 65277b19d0SLisandro Dalcin PetscFunctionBegin; 66277b19d0SLisandro Dalcin ierr = TSReset_Euler(ts);CHKERRQ(ierr); 67277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 683a40ed3dSBarry Smith PetscFunctionReturn(0); 698b1af7b3SBarry Smith } 708b1af7b3SBarry Smith /*------------------------------------------------------------*/ 718b1af7b3SBarry Smith 724416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Euler(PetscOptionItems *PetscOptionsObject,TS ts) 738b1af7b3SBarry Smith { 743a40ed3dSBarry Smith PetscFunctionBegin; 753a40ed3dSBarry Smith PetscFunctionReturn(0); 764b33d51aSBarry Smith } 774b33d51aSBarry Smith 786849ba73SBarry Smith static PetscErrorCode TSView_Euler(TS ts,PetscViewer viewer) 798b1af7b3SBarry Smith { 803a40ed3dSBarry Smith PetscFunctionBegin; 813a40ed3dSBarry Smith PetscFunctionReturn(0); 828b1af7b3SBarry Smith } 838b1af7b3SBarry Smith 846083293cSBarry Smith static PetscErrorCode TSInterpolate_Euler(TS ts,PetscReal t,Vec X) 856083293cSBarry Smith { 863354212dSEmil Constantinescu TS_Euler *euler = (TS_Euler*)ts->data; 873354212dSEmil Constantinescu Vec update = euler->update; 886083293cSBarry Smith PetscReal alpha = (ts->ptime - t)/ts->time_step; 896083293cSBarry Smith PetscErrorCode ierr; 906083293cSBarry Smith 916083293cSBarry Smith PetscFunctionBegin; 923354212dSEmil Constantinescu ierr = VecWAXPY(X,-ts->time_step,update,ts->vec_sol);CHKERRQ(ierr); 933354212dSEmil Constantinescu ierr = VecAXPBY(X,1.0-alpha,alpha,ts->vec_sol);CHKERRQ(ierr); 946083293cSBarry Smith PetscFunctionReturn(0); 956083293cSBarry Smith } 966083293cSBarry Smith 97560360afSLisandro Dalcin static PetscErrorCode TSComputeLinearStability_Euler(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi) 98f9c1d6abSBarry Smith { 99f9c1d6abSBarry Smith PetscFunctionBegin; 100f9c1d6abSBarry Smith *yr = 1.0 + xr; 101f9c1d6abSBarry Smith *yi = xi; 102f9c1d6abSBarry Smith PetscFunctionReturn(0); 103f9c1d6abSBarry Smith } 1048b1af7b3SBarry 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 1119596e0b4SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSBEULER 112ebe3b25bSBarry Smith 113ebe3b25bSBarry Smith M*/ 1148cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Euler(TS ts) 1158b1af7b3SBarry Smith { 1168b1af7b3SBarry Smith TS_Euler *euler; 117dfbe8321SBarry Smith PetscErrorCode ierr; 1188b1af7b3SBarry Smith 1193a40ed3dSBarry Smith PetscFunctionBegin; 1201566a47fSLisandro Dalcin ierr = PetscNewLog(ts,&euler);CHKERRQ(ierr); 1211566a47fSLisandro Dalcin ts->data = (void*)euler; 1221566a47fSLisandro Dalcin 123000e7ae3SMatthew Knepley ts->ops->setup = TSSetUp_Euler; 124000e7ae3SMatthew Knepley ts->ops->step = TSStep_Euler; 125277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Euler; 126000e7ae3SMatthew Knepley ts->ops->destroy = TSDestroy_Euler; 127000e7ae3SMatthew Knepley ts->ops->setfromoptions = TSSetFromOptions_Euler; 128000e7ae3SMatthew Knepley ts->ops->view = TSView_Euler; 1296083293cSBarry Smith ts->ops->interpolate = TSInterpolate_Euler; 130f9c1d6abSBarry Smith ts->ops->linearstability = TSComputeLinearStability_Euler; 1312ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 132*efd4aadfSBarry Smith ts->usessnes = PETSC_FALSE; 1333a40ed3dSBarry Smith PetscFunctionReturn(0); 1344b33d51aSBarry Smith } 135