xref: /petsc/src/ts/impls/explicit/euler/euler.c (revision 3dd83b3843bef2e44404f8ac40307e9aa6509644)
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 
104a2ae208SSatish Balay #undef __FUNCT__
114a2ae208SSatish Balay #define __FUNCT__ "TSStep_Euler"
12193ac0bcSJed Brown static PetscErrorCode TSStep_Euler(TS ts)
134b33d51aSBarry Smith {
148b1af7b3SBarry Smith   TS_Euler       *euler = (TS_Euler*)ts->data;
151566a47fSLisandro Dalcin   Vec            solution = ts->vec_sol,update = euler->update;
161566a47fSLisandro Dalcin   PetscBool      stageok,accept = PETSC_TRUE;
171566a47fSLisandro Dalcin   PetscReal      next_time_step = ts->time_step;
18277b19d0SLisandro Dalcin   PetscErrorCode ierr;
194b33d51aSBarry Smith 
203a40ed3dSBarry Smith   PetscFunctionBegin;
21b8123daeSJed Brown   ierr = TSPreStage(ts,ts->ptime);CHKERRQ(ierr);
221566a47fSLisandro Dalcin   ierr = TSComputeRHSFunction(ts,ts->ptime,solution,update);CHKERRQ(ierr);
231566a47fSLisandro Dalcin   ierr = VecAYPX(update,ts->time_step,solution);CHKERRQ(ierr);
241566a47fSLisandro Dalcin   ierr = TSPostStage(ts,ts->ptime,0,&solution);CHKERRQ(ierr);
251566a47fSLisandro Dalcin   ierr = TSAdaptCheckStage(ts->adapt,ts,ts->ptime,solution,&stageok);CHKERRQ(ierr);
261566a47fSLisandro Dalcin   if(!stageok) {ts->reason = TS_DIVERGED_STEP_REJECTED; PetscFunctionReturn(0);}
271566a47fSLisandro Dalcin   ierr = TSFunctionDomainError(ts,ts->ptime+ts->time_step,update,&stageok);CHKERRQ(ierr);
281566a47fSLisandro Dalcin   if(!stageok) {ts->reason = TS_DIVERGED_STEP_REJECTED; PetscFunctionReturn(0);}
291566a47fSLisandro Dalcin 
301566a47fSLisandro Dalcin   ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
311566a47fSLisandro Dalcin   if (!accept) {ts->reason = TS_DIVERGED_STEP_REJECTED; PetscFunctionReturn(0);}
321566a47fSLisandro Dalcin   ierr = VecCopy(update,solution);CHKERRQ(ierr);
331566a47fSLisandro Dalcin 
34186e87acSLisandro Dalcin   ts->ptime += ts->time_step;
351566a47fSLisandro Dalcin   ts->time_step = next_time_step;
363a40ed3dSBarry Smith   PetscFunctionReturn(0);
374b33d51aSBarry Smith }
384b33d51aSBarry Smith /*------------------------------------------------------------*/
39277b19d0SLisandro Dalcin 
404a2ae208SSatish Balay #undef __FUNCT__
41277b19d0SLisandro Dalcin #define __FUNCT__ "TSSetUp_Euler"
42277b19d0SLisandro Dalcin static PetscErrorCode TSSetUp_Euler(TS ts)
43277b19d0SLisandro Dalcin {
44277b19d0SLisandro Dalcin   TS_Euler       *euler = (TS_Euler*)ts->data;
45277b19d0SLisandro Dalcin   PetscErrorCode ierr;
46277b19d0SLisandro Dalcin 
47277b19d0SLisandro Dalcin   PetscFunctionBegin;
48*3dd83b38SBarry Smith   ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr);
49277b19d0SLisandro Dalcin   ierr = VecDuplicate(ts->vec_sol,&euler->update);CHKERRQ(ierr);
501566a47fSLisandro Dalcin   ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr);
511566a47fSLisandro Dalcin   ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr);
52277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
53277b19d0SLisandro Dalcin }
54277b19d0SLisandro Dalcin 
55277b19d0SLisandro Dalcin #undef __FUNCT__
56277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Euler"
57277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Euler(TS ts)
584b33d51aSBarry Smith {
598b1af7b3SBarry Smith   TS_Euler       *euler = (TS_Euler*)ts->data;
60dfbe8321SBarry Smith   PetscErrorCode ierr;
618b1af7b3SBarry Smith 
623a40ed3dSBarry Smith   PetscFunctionBegin;
636bf464f9SBarry Smith   ierr = VecDestroy(&euler->update);CHKERRQ(ierr);
64277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
65277b19d0SLisandro Dalcin }
66277b19d0SLisandro Dalcin 
67277b19d0SLisandro Dalcin #undef __FUNCT__
68277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Euler"
69277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Euler(TS ts)
70277b19d0SLisandro Dalcin {
71277b19d0SLisandro Dalcin   PetscErrorCode ierr;
72277b19d0SLisandro Dalcin 
73277b19d0SLisandro Dalcin   PetscFunctionBegin;
74277b19d0SLisandro Dalcin   ierr = TSReset_Euler(ts);CHKERRQ(ierr);
75277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
763a40ed3dSBarry Smith   PetscFunctionReturn(0);
778b1af7b3SBarry Smith }
788b1af7b3SBarry Smith /*------------------------------------------------------------*/
798b1af7b3SBarry Smith 
804a2ae208SSatish Balay #undef __FUNCT__
814a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Euler"
824416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Euler(PetscOptionItems *PetscOptionsObject,TS ts)
838b1af7b3SBarry Smith {
843a40ed3dSBarry Smith   PetscFunctionBegin;
853a40ed3dSBarry Smith   PetscFunctionReturn(0);
864b33d51aSBarry Smith }
874b33d51aSBarry Smith 
884a2ae208SSatish Balay #undef __FUNCT__
894a2ae208SSatish Balay #define __FUNCT__ "TSView_Euler"
906849ba73SBarry Smith static PetscErrorCode TSView_Euler(TS ts,PetscViewer viewer)
918b1af7b3SBarry Smith {
923a40ed3dSBarry Smith   PetscFunctionBegin;
933a40ed3dSBarry Smith   PetscFunctionReturn(0);
948b1af7b3SBarry Smith }
958b1af7b3SBarry Smith 
966083293cSBarry Smith #undef __FUNCT__
976083293cSBarry Smith #define __FUNCT__ "TSInterpolate_Euler"
986083293cSBarry Smith static PetscErrorCode TSInterpolate_Euler(TS ts,PetscReal t,Vec X)
996083293cSBarry Smith {
1003354212dSEmil Constantinescu   TS_Euler       *euler = (TS_Euler*)ts->data;
1013354212dSEmil Constantinescu   Vec            update = euler->update;
1026083293cSBarry Smith   PetscReal      alpha = (ts->ptime - t)/ts->time_step;
1036083293cSBarry Smith   PetscErrorCode ierr;
1046083293cSBarry Smith 
1056083293cSBarry Smith   PetscFunctionBegin;
1063354212dSEmil Constantinescu   ierr = VecWAXPY(X,-ts->time_step,update,ts->vec_sol);CHKERRQ(ierr);
1073354212dSEmil Constantinescu   ierr = VecAXPBY(X,1.0-alpha,alpha,ts->vec_sol);CHKERRQ(ierr);
1086083293cSBarry Smith   PetscFunctionReturn(0);
1096083293cSBarry Smith }
1106083293cSBarry Smith 
111f9c1d6abSBarry Smith #undef __FUNCT__
112f9c1d6abSBarry Smith #define __FUNCT__ "TSComputeLinearStability_Euler"
113560360afSLisandro Dalcin static PetscErrorCode TSComputeLinearStability_Euler(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
114f9c1d6abSBarry Smith {
115f9c1d6abSBarry Smith   PetscFunctionBegin;
116f9c1d6abSBarry Smith   *yr = 1.0 + xr;
117f9c1d6abSBarry Smith   *yi = xi;
118f9c1d6abSBarry Smith   PetscFunctionReturn(0);
119f9c1d6abSBarry Smith }
1208b1af7b3SBarry Smith /* ------------------------------------------------------------ */
121ebe3b25bSBarry Smith 
122ebe3b25bSBarry Smith /*MC
1239596e0b4SJed Brown       TSEULER - ODE solver using the explicit forward Euler method
124ebe3b25bSBarry Smith 
125d41222bbSBarry Smith   Level: beginner
126d41222bbSBarry Smith 
1279596e0b4SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSBEULER
128ebe3b25bSBarry Smith 
129ebe3b25bSBarry Smith M*/
1304a2ae208SSatish Balay #undef __FUNCT__
1314a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Euler"
1328cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Euler(TS ts)
1338b1af7b3SBarry Smith {
1348b1af7b3SBarry Smith   TS_Euler       *euler;
135dfbe8321SBarry Smith   PetscErrorCode ierr;
1368b1af7b3SBarry Smith 
1373a40ed3dSBarry Smith   PetscFunctionBegin;
1381566a47fSLisandro Dalcin   ierr = PetscNewLog(ts,&euler);CHKERRQ(ierr);
1391566a47fSLisandro Dalcin   ts->data = (void*)euler;
1401566a47fSLisandro Dalcin 
141000e7ae3SMatthew Knepley   ts->ops->setup           = TSSetUp_Euler;
142000e7ae3SMatthew Knepley   ts->ops->step            = TSStep_Euler;
143277b19d0SLisandro Dalcin   ts->ops->reset           = TSReset_Euler;
144000e7ae3SMatthew Knepley   ts->ops->destroy         = TSDestroy_Euler;
145000e7ae3SMatthew Knepley   ts->ops->setfromoptions  = TSSetFromOptions_Euler;
146000e7ae3SMatthew Knepley   ts->ops->view            = TSView_Euler;
1476083293cSBarry Smith   ts->ops->interpolate     = TSInterpolate_Euler;
148f9c1d6abSBarry Smith   ts->ops->linearstability = TSComputeLinearStability_Euler;
1493a40ed3dSBarry Smith   PetscFunctionReturn(0);
1504b33d51aSBarry Smith }
151