xref: /petsc/src/ts/impls/explicit/euler/euler.c (revision 7087cfbefd1a42b179f217f9994fb6cb0d0c1824)
163dd3a1aSKris Buschelman #define PETSCTS_DLL
263dd3a1aSKris Buschelman 
34b33d51aSBarry Smith /*
4fb4a63b6SLois Curfman McInnes        Code for Timestepping with explicit Euler.
54b33d51aSBarry Smith */
67c4f633dSBarry Smith #include "private/tsimpl.h"                /*I   "petscts.h"   I*/
74b33d51aSBarry Smith 
88b1af7b3SBarry Smith typedef struct {
98b1af7b3SBarry Smith   Vec update;     /* work vector where F(t[i],u[i]) is stored */
108b1af7b3SBarry Smith } TS_Euler;
118b1af7b3SBarry Smith 
124a2ae208SSatish Balay #undef __FUNCT__
134a2ae208SSatish Balay #define __FUNCT__ "TSSetUp_Euler"
146849ba73SBarry Smith static PetscErrorCode TSSetUp_Euler(TS ts)
154b33d51aSBarry Smith {
168b1af7b3SBarry Smith   TS_Euler       *euler = (TS_Euler*)ts->data;
17dfbe8321SBarry Smith   PetscErrorCode ierr;
184b33d51aSBarry Smith 
193a40ed3dSBarry Smith   PetscFunctionBegin;
208b1af7b3SBarry Smith   ierr = VecDuplicate(ts->vec_sol,&euler->update);CHKERRQ(ierr);
213a40ed3dSBarry Smith   PetscFunctionReturn(0);
224b33d51aSBarry Smith }
234b33d51aSBarry Smith 
244a2ae208SSatish Balay #undef __FUNCT__
254a2ae208SSatish Balay #define __FUNCT__ "TSStep_Euler"
26a7cc72afSBarry Smith static PetscErrorCode TSStep_Euler(TS ts,PetscInt *steps,PetscReal *ptime)
274b33d51aSBarry Smith {
288b1af7b3SBarry Smith   TS_Euler       *euler = (TS_Euler*)ts->data;
298b1af7b3SBarry Smith   Vec            sol = ts->vec_sol,update = euler->update;
30dfbe8321SBarry Smith   PetscErrorCode ierr;
31a7cc72afSBarry Smith   PetscInt       i,max_steps = ts->max_steps;
324b33d51aSBarry Smith 
333a40ed3dSBarry Smith   PetscFunctionBegin;
348b1af7b3SBarry Smith   *steps = -ts->steps;
35c3e30b67SBarry Smith   ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
364b33d51aSBarry Smith 
378b1af7b3SBarry Smith   for (i=0; i<max_steps; i++) {
38a057ae4aSMatthew Knepley     PetscReal dt = ts->time_step;
39a057ae4aSMatthew Knepley 
403f2090d5SJed Brown     ierr = TSPreStep(ts);CHKERRQ(ierr);
41a057ae4aSMatthew Knepley     ts->ptime += dt;
42c3e30b67SBarry Smith     ierr = TSComputeRHSFunction(ts,ts->ptime,sol,update);CHKERRQ(ierr);
432dcb1b2aSMatthew Knepley     ierr = VecAXPY(sol,dt,update);CHKERRQ(ierr);
448b1af7b3SBarry Smith     ts->steps++;
453f2090d5SJed Brown     ierr = TSPostStep(ts);CHKERRQ(ierr);
46c3e30b67SBarry Smith     ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
478b1af7b3SBarry Smith     if (ts->ptime > ts->max_time) break;
488b1af7b3SBarry Smith   }
494b33d51aSBarry Smith 
508b1af7b3SBarry Smith   *steps += ts->steps;
51142b95e3SSatish Balay   *ptime  = ts->ptime;
523a40ed3dSBarry Smith   PetscFunctionReturn(0);
534b33d51aSBarry Smith }
544b33d51aSBarry Smith /*------------------------------------------------------------*/
554a2ae208SSatish Balay #undef __FUNCT__
564a2ae208SSatish Balay #define __FUNCT__ "TSDestroy_Euler"
576849ba73SBarry Smith static PetscErrorCode TSDestroy_Euler(TS ts)
584b33d51aSBarry Smith {
598b1af7b3SBarry Smith   TS_Euler       *euler = (TS_Euler*)ts->data;
60dfbe8321SBarry Smith   PetscErrorCode ierr;
618b1af7b3SBarry Smith 
623a40ed3dSBarry Smith   PetscFunctionBegin;
631713a123SBarry Smith   if (euler->update) {ierr = VecDestroy(euler->update);CHKERRQ(ierr);}
64606d414cSSatish Balay   ierr = PetscFree(euler);CHKERRQ(ierr);
653a40ed3dSBarry Smith   PetscFunctionReturn(0);
668b1af7b3SBarry Smith }
678b1af7b3SBarry Smith /*------------------------------------------------------------*/
688b1af7b3SBarry Smith 
694a2ae208SSatish Balay #undef __FUNCT__
704a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Euler"
716849ba73SBarry Smith static PetscErrorCode TSSetFromOptions_Euler(TS ts)
728b1af7b3SBarry Smith {
733a40ed3dSBarry Smith   PetscFunctionBegin;
743a40ed3dSBarry Smith   PetscFunctionReturn(0);
754b33d51aSBarry Smith }
764b33d51aSBarry Smith 
774a2ae208SSatish Balay #undef __FUNCT__
784a2ae208SSatish Balay #define __FUNCT__ "TSView_Euler"
796849ba73SBarry Smith static PetscErrorCode TSView_Euler(TS ts,PetscViewer viewer)
808b1af7b3SBarry Smith {
813a40ed3dSBarry Smith   PetscFunctionBegin;
823a40ed3dSBarry Smith   PetscFunctionReturn(0);
838b1af7b3SBarry Smith }
848b1af7b3SBarry Smith 
858b1af7b3SBarry Smith /* ------------------------------------------------------------ */
86ebe3b25bSBarry Smith 
87ebe3b25bSBarry Smith /*MC
889596e0b4SJed Brown       TSEULER - ODE solver using the explicit forward Euler method
89ebe3b25bSBarry Smith 
90d41222bbSBarry Smith   Level: beginner
91d41222bbSBarry Smith 
929596e0b4SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSBEULER
93ebe3b25bSBarry Smith 
94ebe3b25bSBarry Smith M*/
95fb2e594dSBarry Smith EXTERN_C_BEGIN
964a2ae208SSatish Balay #undef __FUNCT__
974a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Euler"
98*7087cfbeSBarry Smith PetscErrorCode  TSCreate_Euler(TS ts)
998b1af7b3SBarry Smith {
1008b1af7b3SBarry Smith   TS_Euler       *euler;
101dfbe8321SBarry Smith   PetscErrorCode ierr;
1028b1af7b3SBarry Smith 
1033a40ed3dSBarry Smith   PetscFunctionBegin;
104000e7ae3SMatthew Knepley   ts->ops->setup           = TSSetUp_Euler;
105000e7ae3SMatthew Knepley   ts->ops->step            = TSStep_Euler;
106000e7ae3SMatthew Knepley   ts->ops->destroy         = TSDestroy_Euler;
107000e7ae3SMatthew Knepley   ts->ops->setfromoptions  = TSSetFromOptions_Euler;
108000e7ae3SMatthew Knepley   ts->ops->view            = TSView_Euler;
1098b1af7b3SBarry Smith 
11038f2d2fdSLisandro Dalcin   ierr = PetscNewLog(ts,TS_Euler,&euler);CHKERRQ(ierr);
1118b1af7b3SBarry Smith   ts->data = (void*)euler;
1124b33d51aSBarry Smith 
1133a40ed3dSBarry Smith   PetscFunctionReturn(0);
1144b33d51aSBarry Smith }
115fb2e594dSBarry Smith EXTERN_C_END
116c3e30b67SBarry Smith 
117c3e30b67SBarry Smith 
118c3e30b67SBarry Smith 
119c3e30b67SBarry Smith 
120