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 40a057ae4aSMatthew Knepley ts->ptime += dt; 41c3e30b67SBarry Smith ierr = TSComputeRHSFunction(ts,ts->ptime,sol,update);CHKERRQ(ierr); 422dcb1b2aSMatthew Knepley ierr = VecAXPY(sol,dt,update);CHKERRQ(ierr); 438b1af7b3SBarry Smith ts->steps++; 44c3e30b67SBarry Smith ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 458b1af7b3SBarry Smith if (ts->ptime > ts->max_time) break; 468b1af7b3SBarry Smith } 474b33d51aSBarry Smith 488b1af7b3SBarry Smith *steps += ts->steps; 49142b95e3SSatish Balay *ptime = ts->ptime; 503a40ed3dSBarry Smith PetscFunctionReturn(0); 514b33d51aSBarry Smith } 524b33d51aSBarry Smith /*------------------------------------------------------------*/ 534a2ae208SSatish Balay #undef __FUNCT__ 544a2ae208SSatish Balay #define __FUNCT__ "TSDestroy_Euler" 556849ba73SBarry Smith static PetscErrorCode TSDestroy_Euler(TS ts) 564b33d51aSBarry Smith { 578b1af7b3SBarry Smith TS_Euler *euler = (TS_Euler*)ts->data; 58dfbe8321SBarry Smith PetscErrorCode ierr; 598b1af7b3SBarry Smith 603a40ed3dSBarry Smith PetscFunctionBegin; 611713a123SBarry Smith if (euler->update) {ierr = VecDestroy(euler->update);CHKERRQ(ierr);} 62606d414cSSatish Balay ierr = PetscFree(euler);CHKERRQ(ierr); 633a40ed3dSBarry Smith PetscFunctionReturn(0); 648b1af7b3SBarry Smith } 658b1af7b3SBarry Smith /*------------------------------------------------------------*/ 668b1af7b3SBarry Smith 674a2ae208SSatish Balay #undef __FUNCT__ 684a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Euler" 696849ba73SBarry Smith static PetscErrorCode TSSetFromOptions_Euler(TS ts) 708b1af7b3SBarry Smith { 713a40ed3dSBarry Smith PetscFunctionBegin; 723a40ed3dSBarry Smith PetscFunctionReturn(0); 734b33d51aSBarry Smith } 744b33d51aSBarry Smith 754a2ae208SSatish Balay #undef __FUNCT__ 764a2ae208SSatish Balay #define __FUNCT__ "TSView_Euler" 776849ba73SBarry Smith static PetscErrorCode TSView_Euler(TS ts,PetscViewer viewer) 788b1af7b3SBarry Smith { 793a40ed3dSBarry Smith PetscFunctionBegin; 803a40ed3dSBarry Smith PetscFunctionReturn(0); 818b1af7b3SBarry Smith } 828b1af7b3SBarry Smith 838b1af7b3SBarry Smith /* ------------------------------------------------------------ */ 84ebe3b25bSBarry Smith 85ebe3b25bSBarry Smith /*MC 86*9596e0b4SJed Brown TSEULER - ODE solver using the explicit forward Euler method 87ebe3b25bSBarry Smith 88d41222bbSBarry Smith Level: beginner 89d41222bbSBarry Smith 90*9596e0b4SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSBEULER 91ebe3b25bSBarry Smith 92ebe3b25bSBarry Smith M*/ 93fb2e594dSBarry Smith EXTERN_C_BEGIN 944a2ae208SSatish Balay #undef __FUNCT__ 954a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Euler" 9663dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSCreate_Euler(TS ts) 978b1af7b3SBarry Smith { 988b1af7b3SBarry Smith TS_Euler *euler; 99dfbe8321SBarry Smith PetscErrorCode ierr; 1008b1af7b3SBarry Smith 1013a40ed3dSBarry Smith PetscFunctionBegin; 102000e7ae3SMatthew Knepley ts->ops->setup = TSSetUp_Euler; 103000e7ae3SMatthew Knepley ts->ops->step = TSStep_Euler; 104000e7ae3SMatthew Knepley ts->ops->destroy = TSDestroy_Euler; 105000e7ae3SMatthew Knepley ts->ops->setfromoptions = TSSetFromOptions_Euler; 106000e7ae3SMatthew Knepley ts->ops->view = TSView_Euler; 1078b1af7b3SBarry Smith 10838f2d2fdSLisandro Dalcin ierr = PetscNewLog(ts,TS_Euler,&euler);CHKERRQ(ierr); 1098b1af7b3SBarry Smith ts->data = (void*)euler; 1104b33d51aSBarry Smith 1113a40ed3dSBarry Smith PetscFunctionReturn(0); 1124b33d51aSBarry Smith } 113fb2e594dSBarry Smith EXTERN_C_END 114c3e30b67SBarry Smith 115c3e30b67SBarry Smith 116c3e30b67SBarry Smith 117c3e30b67SBarry Smith 118