163dd3a1aSKris Buschelman 24b33d51aSBarry Smith /* 3fb4a63b6SLois Curfman McInnes Code for Timestepping with explicit Euler. 44b33d51aSBarry Smith */ 5*c6db04a5SJed Brown #include <private/tsimpl.h> /*I "petscts.h" I*/ 64b33d51aSBarry Smith 78b1af7b3SBarry Smith typedef struct { 88b1af7b3SBarry Smith Vec update; /* work vector where F(t[i],u[i]) is stored */ 98b1af7b3SBarry Smith } TS_Euler; 108b1af7b3SBarry Smith 114a2ae208SSatish Balay #undef __FUNCT__ 124a2ae208SSatish Balay #define __FUNCT__ "TSSetUp_Euler" 136849ba73SBarry Smith static PetscErrorCode TSSetUp_Euler(TS ts) 144b33d51aSBarry Smith { 158b1af7b3SBarry Smith TS_Euler *euler = (TS_Euler*)ts->data; 16dfbe8321SBarry Smith PetscErrorCode ierr; 174b33d51aSBarry Smith 183a40ed3dSBarry Smith PetscFunctionBegin; 198b1af7b3SBarry Smith ierr = VecDuplicate(ts->vec_sol,&euler->update);CHKERRQ(ierr); 203a40ed3dSBarry Smith PetscFunctionReturn(0); 214b33d51aSBarry Smith } 224b33d51aSBarry Smith 234a2ae208SSatish Balay #undef __FUNCT__ 244a2ae208SSatish Balay #define __FUNCT__ "TSStep_Euler" 25a7cc72afSBarry Smith static PetscErrorCode TSStep_Euler(TS ts,PetscInt *steps,PetscReal *ptime) 264b33d51aSBarry Smith { 278b1af7b3SBarry Smith TS_Euler *euler = (TS_Euler*)ts->data; 288b1af7b3SBarry Smith Vec sol = ts->vec_sol,update = euler->update; 29dfbe8321SBarry Smith PetscErrorCode ierr; 30a7cc72afSBarry Smith PetscInt i,max_steps = ts->max_steps; 314b33d51aSBarry Smith 323a40ed3dSBarry Smith PetscFunctionBegin; 338b1af7b3SBarry Smith *steps = -ts->steps; 34c3e30b67SBarry Smith ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 354b33d51aSBarry Smith 368b1af7b3SBarry Smith for (i=0; i<max_steps; i++) { 37a057ae4aSMatthew Knepley PetscReal dt = ts->time_step; 38a057ae4aSMatthew Knepley 393f2090d5SJed Brown ierr = TSPreStep(ts);CHKERRQ(ierr); 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++; 443f2090d5SJed Brown ierr = TSPostStep(ts);CHKERRQ(ierr); 45c3e30b67SBarry Smith ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 468b1af7b3SBarry Smith if (ts->ptime > ts->max_time) break; 478b1af7b3SBarry Smith } 484b33d51aSBarry Smith 498b1af7b3SBarry Smith *steps += ts->steps; 50142b95e3SSatish Balay *ptime = ts->ptime; 513a40ed3dSBarry Smith PetscFunctionReturn(0); 524b33d51aSBarry Smith } 534b33d51aSBarry Smith /*------------------------------------------------------------*/ 544a2ae208SSatish Balay #undef __FUNCT__ 554a2ae208SSatish Balay #define __FUNCT__ "TSDestroy_Euler" 566849ba73SBarry Smith static PetscErrorCode TSDestroy_Euler(TS ts) 574b33d51aSBarry Smith { 588b1af7b3SBarry Smith TS_Euler *euler = (TS_Euler*)ts->data; 59dfbe8321SBarry Smith PetscErrorCode ierr; 608b1af7b3SBarry Smith 613a40ed3dSBarry Smith PetscFunctionBegin; 621713a123SBarry Smith if (euler->update) {ierr = VecDestroy(euler->update);CHKERRQ(ierr);} 63606d414cSSatish Balay ierr = PetscFree(euler);CHKERRQ(ierr); 643a40ed3dSBarry Smith PetscFunctionReturn(0); 658b1af7b3SBarry Smith } 668b1af7b3SBarry Smith /*------------------------------------------------------------*/ 678b1af7b3SBarry Smith 684a2ae208SSatish Balay #undef __FUNCT__ 694a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Euler" 706849ba73SBarry Smith static PetscErrorCode TSSetFromOptions_Euler(TS ts) 718b1af7b3SBarry Smith { 723a40ed3dSBarry Smith PetscFunctionBegin; 733a40ed3dSBarry Smith PetscFunctionReturn(0); 744b33d51aSBarry Smith } 754b33d51aSBarry Smith 764a2ae208SSatish Balay #undef __FUNCT__ 774a2ae208SSatish Balay #define __FUNCT__ "TSView_Euler" 786849ba73SBarry Smith static PetscErrorCode TSView_Euler(TS ts,PetscViewer viewer) 798b1af7b3SBarry Smith { 803a40ed3dSBarry Smith PetscFunctionBegin; 813a40ed3dSBarry Smith PetscFunctionReturn(0); 828b1af7b3SBarry Smith } 838b1af7b3SBarry Smith 848b1af7b3SBarry Smith /* ------------------------------------------------------------ */ 85ebe3b25bSBarry Smith 86ebe3b25bSBarry Smith /*MC 879596e0b4SJed Brown TSEULER - ODE solver using the explicit forward Euler method 88ebe3b25bSBarry Smith 89d41222bbSBarry Smith Level: beginner 90d41222bbSBarry Smith 919596e0b4SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSBEULER 92ebe3b25bSBarry Smith 93ebe3b25bSBarry Smith M*/ 94fb2e594dSBarry Smith EXTERN_C_BEGIN 954a2ae208SSatish Balay #undef __FUNCT__ 964a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Euler" 977087cfbeSBarry Smith PetscErrorCode TSCreate_Euler(TS ts) 988b1af7b3SBarry Smith { 998b1af7b3SBarry Smith TS_Euler *euler; 100dfbe8321SBarry Smith PetscErrorCode ierr; 1018b1af7b3SBarry Smith 1023a40ed3dSBarry Smith PetscFunctionBegin; 103000e7ae3SMatthew Knepley ts->ops->setup = TSSetUp_Euler; 104000e7ae3SMatthew Knepley ts->ops->step = TSStep_Euler; 105000e7ae3SMatthew Knepley ts->ops->destroy = TSDestroy_Euler; 106000e7ae3SMatthew Knepley ts->ops->setfromoptions = TSSetFromOptions_Euler; 107000e7ae3SMatthew Knepley ts->ops->view = TSView_Euler; 1088b1af7b3SBarry Smith 10938f2d2fdSLisandro Dalcin ierr = PetscNewLog(ts,TS_Euler,&euler);CHKERRQ(ierr); 1108b1af7b3SBarry Smith ts->data = (void*)euler; 1114b33d51aSBarry Smith 1123a40ed3dSBarry Smith PetscFunctionReturn(0); 1134b33d51aSBarry Smith } 114fb2e594dSBarry Smith EXTERN_C_END 115c3e30b67SBarry Smith 116c3e30b67SBarry Smith 117c3e30b67SBarry Smith 118c3e30b67SBarry Smith 119