14b33d51aSBarry Smith /* 2fb4a63b6SLois Curfman McInnes Code for Timestepping with explicit Euler. 34b33d51aSBarry Smith */ 4c6db04a5SJed Brown #include <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" 12a7cc72afSBarry Smith static PetscErrorCode TSStep_Euler(TS ts,PetscInt *steps,PetscReal *ptime) 134b33d51aSBarry Smith { 148b1af7b3SBarry Smith TS_Euler *euler = (TS_Euler*)ts->data; 158b1af7b3SBarry Smith Vec sol = ts->vec_sol,update = euler->update; 16*186e87acSLisandro Dalcin PetscInt i; 17277b19d0SLisandro Dalcin PetscErrorCode ierr; 184b33d51aSBarry Smith 193a40ed3dSBarry Smith PetscFunctionBegin; 208b1af7b3SBarry Smith *steps = -ts->steps; 21*186e87acSLisandro Dalcin *ptime = ts->ptime; 22*186e87acSLisandro Dalcin 23c3e30b67SBarry Smith ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 244b33d51aSBarry Smith 25*186e87acSLisandro Dalcin for (i=0; i<ts->max_steps; i++) { 26*186e87acSLisandro Dalcin if (ts->ptime + ts->time_step > ts->max_time) break; 273f2090d5SJed Brown ierr = TSPreStep(ts);CHKERRQ(ierr); 28*186e87acSLisandro Dalcin 29c3e30b67SBarry Smith ierr = TSComputeRHSFunction(ts,ts->ptime,sol,update);CHKERRQ(ierr); 30*186e87acSLisandro Dalcin 31*186e87acSLisandro Dalcin ierr = VecAXPY(sol,ts->time_step,update);CHKERRQ(ierr); 32*186e87acSLisandro Dalcin ts->ptime += ts->time_step; 338b1af7b3SBarry Smith ts->steps++; 34*186e87acSLisandro Dalcin 353f2090d5SJed Brown ierr = TSPostStep(ts);CHKERRQ(ierr); 36c3e30b67SBarry Smith ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 378b1af7b3SBarry Smith } 384b33d51aSBarry Smith 398b1af7b3SBarry Smith *steps += ts->steps; 40142b95e3SSatish Balay *ptime = ts->ptime; 413a40ed3dSBarry Smith PetscFunctionReturn(0); 424b33d51aSBarry Smith } 434b33d51aSBarry Smith /*------------------------------------------------------------*/ 44277b19d0SLisandro Dalcin 454a2ae208SSatish Balay #undef __FUNCT__ 46277b19d0SLisandro Dalcin #define __FUNCT__ "TSSetUp_Euler" 47277b19d0SLisandro Dalcin static PetscErrorCode TSSetUp_Euler(TS ts) 48277b19d0SLisandro Dalcin { 49277b19d0SLisandro Dalcin TS_Euler *euler = (TS_Euler*)ts->data; 50277b19d0SLisandro Dalcin PetscErrorCode ierr; 51277b19d0SLisandro Dalcin 52277b19d0SLisandro Dalcin PetscFunctionBegin; 53277b19d0SLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&euler->update);CHKERRQ(ierr); 54277b19d0SLisandro Dalcin PetscFunctionReturn(0); 55277b19d0SLisandro Dalcin } 56277b19d0SLisandro Dalcin 57277b19d0SLisandro Dalcin #undef __FUNCT__ 58277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Euler" 59277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Euler(TS ts) 604b33d51aSBarry Smith { 618b1af7b3SBarry Smith TS_Euler *euler = (TS_Euler*)ts->data; 62dfbe8321SBarry Smith PetscErrorCode ierr; 638b1af7b3SBarry Smith 643a40ed3dSBarry Smith PetscFunctionBegin; 651713a123SBarry Smith if (euler->update) {ierr = VecDestroy(euler->update);CHKERRQ(ierr);} 66277b19d0SLisandro Dalcin PetscFunctionReturn(0); 67277b19d0SLisandro Dalcin } 68277b19d0SLisandro Dalcin 69277b19d0SLisandro Dalcin #undef __FUNCT__ 70277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Euler" 71277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Euler(TS ts) 72277b19d0SLisandro Dalcin { 73277b19d0SLisandro Dalcin PetscErrorCode ierr; 74277b19d0SLisandro Dalcin 75277b19d0SLisandro Dalcin PetscFunctionBegin; 76277b19d0SLisandro Dalcin ierr = TSReset_Euler(ts);CHKERRQ(ierr); 77277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 783a40ed3dSBarry Smith PetscFunctionReturn(0); 798b1af7b3SBarry Smith } 808b1af7b3SBarry Smith /*------------------------------------------------------------*/ 818b1af7b3SBarry Smith 824a2ae208SSatish Balay #undef __FUNCT__ 834a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Euler" 846849ba73SBarry Smith static PetscErrorCode TSSetFromOptions_Euler(TS ts) 858b1af7b3SBarry Smith { 863a40ed3dSBarry Smith PetscFunctionBegin; 873a40ed3dSBarry Smith PetscFunctionReturn(0); 884b33d51aSBarry Smith } 894b33d51aSBarry Smith 904a2ae208SSatish Balay #undef __FUNCT__ 914a2ae208SSatish Balay #define __FUNCT__ "TSView_Euler" 926849ba73SBarry Smith static PetscErrorCode TSView_Euler(TS ts,PetscViewer viewer) 938b1af7b3SBarry Smith { 943a40ed3dSBarry Smith PetscFunctionBegin; 953a40ed3dSBarry Smith PetscFunctionReturn(0); 968b1af7b3SBarry Smith } 978b1af7b3SBarry 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*/ 108fb2e594dSBarry Smith EXTERN_C_BEGIN 1094a2ae208SSatish Balay #undef __FUNCT__ 1104a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Euler" 1117087cfbeSBarry Smith PetscErrorCode TSCreate_Euler(TS ts) 1128b1af7b3SBarry Smith { 1138b1af7b3SBarry Smith TS_Euler *euler; 114dfbe8321SBarry Smith PetscErrorCode ierr; 1158b1af7b3SBarry Smith 1163a40ed3dSBarry Smith PetscFunctionBegin; 117000e7ae3SMatthew Knepley ts->ops->setup = TSSetUp_Euler; 118000e7ae3SMatthew Knepley ts->ops->step = TSStep_Euler; 119277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Euler; 120000e7ae3SMatthew Knepley ts->ops->destroy = TSDestroy_Euler; 121000e7ae3SMatthew Knepley ts->ops->setfromoptions = TSSetFromOptions_Euler; 122000e7ae3SMatthew Knepley ts->ops->view = TSView_Euler; 1238b1af7b3SBarry Smith 12438f2d2fdSLisandro Dalcin ierr = PetscNewLog(ts,TS_Euler,&euler);CHKERRQ(ierr); 1258b1af7b3SBarry Smith ts->data = (void*)euler; 1264b33d51aSBarry Smith 1273a40ed3dSBarry Smith PetscFunctionReturn(0); 1284b33d51aSBarry Smith } 129fb2e594dSBarry Smith EXTERN_C_END 130