xref: /petsc/src/ts/impls/explicit/euler/euler.c (revision 560360af1e0197128c3ad6271bbfa2e76ad345d4)
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;
368b1af7b3SBarry Smith   ts->steps++;
373a40ed3dSBarry Smith   PetscFunctionReturn(0);
384b33d51aSBarry Smith }
394b33d51aSBarry Smith /*------------------------------------------------------------*/
40277b19d0SLisandro Dalcin 
414a2ae208SSatish Balay #undef __FUNCT__
42277b19d0SLisandro Dalcin #define __FUNCT__ "TSSetUp_Euler"
43277b19d0SLisandro Dalcin static PetscErrorCode TSSetUp_Euler(TS ts)
44277b19d0SLisandro Dalcin {
45277b19d0SLisandro Dalcin   TS_Euler       *euler = (TS_Euler*)ts->data;
46277b19d0SLisandro Dalcin   PetscErrorCode ierr;
471dca01adSEmil Constantinescu   TSRHSFunction  rhsfunction;
481dca01adSEmil Constantinescu   TSIFunction    ifunction;
49277b19d0SLisandro Dalcin 
50277b19d0SLisandro Dalcin   PetscFunctionBegin;
511dca01adSEmil Constantinescu   ierr =  TSGetIFunction(ts,NULL,&ifunction,NULL);CHKERRQ(ierr);
521dca01adSEmil Constantinescu   ierr =  TSGetRHSFunction(ts,NULL,&rhsfunction,NULL);CHKERRQ(ierr);
53a02a3cfeSEmil Constantinescu   if (!rhsfunction || ifunction) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must define RHSFunction() and leave IFunction() undefined in order to use -ts_type euler");
54277b19d0SLisandro Dalcin   ierr = VecDuplicate(ts->vec_sol,&euler->update);CHKERRQ(ierr);
551566a47fSLisandro Dalcin 
561566a47fSLisandro Dalcin   ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr);
571566a47fSLisandro Dalcin   ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr);
58277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
59277b19d0SLisandro Dalcin }
60277b19d0SLisandro Dalcin 
61277b19d0SLisandro Dalcin #undef __FUNCT__
62277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Euler"
63277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Euler(TS ts)
644b33d51aSBarry Smith {
658b1af7b3SBarry Smith   TS_Euler       *euler = (TS_Euler*)ts->data;
66dfbe8321SBarry Smith   PetscErrorCode ierr;
678b1af7b3SBarry Smith 
683a40ed3dSBarry Smith   PetscFunctionBegin;
696bf464f9SBarry Smith   ierr = VecDestroy(&euler->update);CHKERRQ(ierr);
70277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
71277b19d0SLisandro Dalcin }
72277b19d0SLisandro Dalcin 
73277b19d0SLisandro Dalcin #undef __FUNCT__
74277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Euler"
75277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Euler(TS ts)
76277b19d0SLisandro Dalcin {
77277b19d0SLisandro Dalcin   PetscErrorCode ierr;
78277b19d0SLisandro Dalcin 
79277b19d0SLisandro Dalcin   PetscFunctionBegin;
80277b19d0SLisandro Dalcin   ierr = TSReset_Euler(ts);CHKERRQ(ierr);
81277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
823a40ed3dSBarry Smith   PetscFunctionReturn(0);
838b1af7b3SBarry Smith }
848b1af7b3SBarry Smith /*------------------------------------------------------------*/
858b1af7b3SBarry Smith 
864a2ae208SSatish Balay #undef __FUNCT__
874a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Euler"
884416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Euler(PetscOptionItems *PetscOptionsObject,TS ts)
898b1af7b3SBarry Smith {
903a40ed3dSBarry Smith   PetscFunctionBegin;
913a40ed3dSBarry Smith   PetscFunctionReturn(0);
924b33d51aSBarry Smith }
934b33d51aSBarry Smith 
944a2ae208SSatish Balay #undef __FUNCT__
954a2ae208SSatish Balay #define __FUNCT__ "TSView_Euler"
966849ba73SBarry Smith static PetscErrorCode TSView_Euler(TS ts,PetscViewer viewer)
978b1af7b3SBarry Smith {
983a40ed3dSBarry Smith   PetscFunctionBegin;
993a40ed3dSBarry Smith   PetscFunctionReturn(0);
1008b1af7b3SBarry Smith }
1018b1af7b3SBarry Smith 
1026083293cSBarry Smith #undef __FUNCT__
1036083293cSBarry Smith #define __FUNCT__ "TSInterpolate_Euler"
1046083293cSBarry Smith static PetscErrorCode TSInterpolate_Euler(TS ts,PetscReal t,Vec X)
1056083293cSBarry Smith {
1063354212dSEmil Constantinescu   TS_Euler       *euler = (TS_Euler*)ts->data;
1073354212dSEmil Constantinescu   Vec            update = euler->update;
1086083293cSBarry Smith   PetscReal      alpha = (ts->ptime - t)/ts->time_step;
1096083293cSBarry Smith   PetscErrorCode ierr;
1106083293cSBarry Smith 
1116083293cSBarry Smith   PetscFunctionBegin;
1123354212dSEmil Constantinescu   ierr = VecWAXPY(X,-ts->time_step,update,ts->vec_sol);CHKERRQ(ierr);
1133354212dSEmil Constantinescu   ierr = VecAXPBY(X,1.0-alpha,alpha,ts->vec_sol);CHKERRQ(ierr);
1146083293cSBarry Smith   PetscFunctionReturn(0);
1156083293cSBarry Smith }
1166083293cSBarry Smith 
117f9c1d6abSBarry Smith #undef __FUNCT__
118f9c1d6abSBarry Smith #define __FUNCT__ "TSComputeLinearStability_Euler"
119*560360afSLisandro Dalcin static PetscErrorCode TSComputeLinearStability_Euler(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
120f9c1d6abSBarry Smith {
121f9c1d6abSBarry Smith   PetscFunctionBegin;
122f9c1d6abSBarry Smith   *yr = 1.0 + xr;
123f9c1d6abSBarry Smith   *yi = xi;
124f9c1d6abSBarry Smith   PetscFunctionReturn(0);
125f9c1d6abSBarry Smith }
1268b1af7b3SBarry Smith /* ------------------------------------------------------------ */
127ebe3b25bSBarry Smith 
128ebe3b25bSBarry Smith /*MC
1299596e0b4SJed Brown       TSEULER - ODE solver using the explicit forward Euler method
130ebe3b25bSBarry Smith 
131d41222bbSBarry Smith   Level: beginner
132d41222bbSBarry Smith 
1339596e0b4SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSBEULER
134ebe3b25bSBarry Smith 
135ebe3b25bSBarry Smith M*/
1364a2ae208SSatish Balay #undef __FUNCT__
1374a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Euler"
1388cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Euler(TS ts)
1398b1af7b3SBarry Smith {
1408b1af7b3SBarry Smith   TS_Euler       *euler;
141dfbe8321SBarry Smith   PetscErrorCode ierr;
1428b1af7b3SBarry Smith 
1433a40ed3dSBarry Smith   PetscFunctionBegin;
1441566a47fSLisandro Dalcin   ierr = PetscNewLog(ts,&euler);CHKERRQ(ierr);
1451566a47fSLisandro Dalcin   ts->data = (void*)euler;
1461566a47fSLisandro Dalcin 
147000e7ae3SMatthew Knepley   ts->ops->setup           = TSSetUp_Euler;
148000e7ae3SMatthew Knepley   ts->ops->step            = TSStep_Euler;
149277b19d0SLisandro Dalcin   ts->ops->reset           = TSReset_Euler;
150000e7ae3SMatthew Knepley   ts->ops->destroy         = TSDestroy_Euler;
151000e7ae3SMatthew Knepley   ts->ops->setfromoptions  = TSSetFromOptions_Euler;
152000e7ae3SMatthew Knepley   ts->ops->view            = TSView_Euler;
1536083293cSBarry Smith   ts->ops->interpolate     = TSInterpolate_Euler;
154f9c1d6abSBarry Smith   ts->ops->linearstability = TSComputeLinearStability_Euler;
1553a40ed3dSBarry Smith   PetscFunctionReturn(0);
1564b33d51aSBarry Smith }
157