163dd3a1aSKris Buschelman #define PETSCTS_DLL 263dd3a1aSKris Buschelman 34b33d51aSBarry Smith /* 4fb4a63b6SLois Curfman McInnes Code for Timestepping with explicit Euler. 54b33d51aSBarry Smith */ 6b9147fbbSdalcinl #include "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; 32ea709b57SSatish Balay PetscScalar dt = ts->time_step; 334b33d51aSBarry Smith 343a40ed3dSBarry Smith PetscFunctionBegin; 358b1af7b3SBarry Smith *steps = -ts->steps; 36c3e30b67SBarry Smith ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 374b33d51aSBarry Smith 388b1af7b3SBarry Smith for (i=0; i<max_steps; i++) { 398b1af7b3SBarry Smith ts->ptime += ts->time_step; 40c3e30b67SBarry Smith ierr = TSComputeRHSFunction(ts,ts->ptime,sol,update);CHKERRQ(ierr); 412dcb1b2aSMatthew Knepley ierr = VecAXPY(sol,dt,update);CHKERRQ(ierr); 428b1af7b3SBarry Smith ts->steps++; 43c3e30b67SBarry Smith ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 448b1af7b3SBarry Smith if (ts->ptime > ts->max_time) break; 458b1af7b3SBarry Smith } 464b33d51aSBarry Smith 478b1af7b3SBarry Smith *steps += ts->steps; 48142b95e3SSatish Balay *ptime = ts->ptime; 493a40ed3dSBarry Smith PetscFunctionReturn(0); 504b33d51aSBarry Smith } 514b33d51aSBarry Smith /*------------------------------------------------------------*/ 524a2ae208SSatish Balay #undef __FUNCT__ 534a2ae208SSatish Balay #define __FUNCT__ "TSDestroy_Euler" 546849ba73SBarry Smith static PetscErrorCode TSDestroy_Euler(TS ts) 554b33d51aSBarry Smith { 568b1af7b3SBarry Smith TS_Euler *euler = (TS_Euler*)ts->data; 57dfbe8321SBarry Smith PetscErrorCode ierr; 588b1af7b3SBarry Smith 593a40ed3dSBarry Smith PetscFunctionBegin; 601713a123SBarry Smith if (euler->update) {ierr = VecDestroy(euler->update);CHKERRQ(ierr);} 61606d414cSSatish Balay ierr = PetscFree(euler);CHKERRQ(ierr); 623a40ed3dSBarry Smith PetscFunctionReturn(0); 638b1af7b3SBarry Smith } 648b1af7b3SBarry Smith /*------------------------------------------------------------*/ 658b1af7b3SBarry Smith 664a2ae208SSatish Balay #undef __FUNCT__ 674a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Euler" 686849ba73SBarry Smith static PetscErrorCode TSSetFromOptions_Euler(TS ts) 698b1af7b3SBarry Smith { 703a40ed3dSBarry Smith PetscFunctionBegin; 713a40ed3dSBarry Smith PetscFunctionReturn(0); 724b33d51aSBarry Smith } 734b33d51aSBarry Smith 744a2ae208SSatish Balay #undef __FUNCT__ 754a2ae208SSatish Balay #define __FUNCT__ "TSView_Euler" 766849ba73SBarry Smith static PetscErrorCode TSView_Euler(TS ts,PetscViewer viewer) 778b1af7b3SBarry Smith { 783a40ed3dSBarry Smith PetscFunctionBegin; 793a40ed3dSBarry Smith PetscFunctionReturn(0); 808b1af7b3SBarry Smith } 818b1af7b3SBarry Smith 828b1af7b3SBarry Smith /* ------------------------------------------------------------ */ 83ebe3b25bSBarry Smith 84ebe3b25bSBarry Smith /*MC 85ebe3b25bSBarry Smith TS_EULER - ODE solver using the explicit forward Euler method 86ebe3b25bSBarry Smith 87d41222bbSBarry Smith Level: beginner 88d41222bbSBarry Smith 89ebe3b25bSBarry Smith .seealso: TSCreate(), TS, TSSetType(), TS_BEULER 90ebe3b25bSBarry Smith 91ebe3b25bSBarry Smith M*/ 92fb2e594dSBarry Smith EXTERN_C_BEGIN 934a2ae208SSatish Balay #undef __FUNCT__ 944a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Euler" 9563dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSCreate_Euler(TS ts) 968b1af7b3SBarry Smith { 978b1af7b3SBarry Smith TS_Euler *euler; 98dfbe8321SBarry Smith PetscErrorCode ierr; 998b1af7b3SBarry Smith 1003a40ed3dSBarry Smith PetscFunctionBegin; 101000e7ae3SMatthew Knepley ts->ops->setup = TSSetUp_Euler; 102000e7ae3SMatthew Knepley ts->ops->step = TSStep_Euler; 103000e7ae3SMatthew Knepley ts->ops->destroy = TSDestroy_Euler; 104000e7ae3SMatthew Knepley ts->ops->setfromoptions = TSSetFromOptions_Euler; 105000e7ae3SMatthew Knepley ts->ops->view = TSView_Euler; 1068b1af7b3SBarry Smith 107*38f2d2fdSLisandro Dalcin ierr = PetscNewLog(ts,TS_Euler,&euler);CHKERRQ(ierr); 1088b1af7b3SBarry Smith ts->data = (void*)euler; 1094b33d51aSBarry Smith 1103a40ed3dSBarry Smith PetscFunctionReturn(0); 1114b33d51aSBarry Smith } 112fb2e594dSBarry Smith EXTERN_C_END 113c3e30b67SBarry Smith 114c3e30b67SBarry Smith 115c3e30b67SBarry Smith 116c3e30b67SBarry Smith 117