14b33d51aSBarry Smith /* 2fb4a63b6SLois Curfman McInnes Code for Timestepping with explicit Euler. 34b33d51aSBarry Smith */ 4e090d566SSatish Balay #include "src/ts/tsimpl.h" /*I "petscts.h" I*/ 54b33d51aSBarry Smith 68b1af7b3SBarry Smith typedef struct { 78b1af7b3SBarry Smith Vec update; /* work vector where F(t[i],u[i]) is stored */ 88b1af7b3SBarry Smith } TS_Euler; 98b1af7b3SBarry Smith 104a2ae208SSatish Balay #undef __FUNCT__ 114a2ae208SSatish Balay #define __FUNCT__ "TSSetUp_Euler" 12*6849ba73SBarry Smith static PetscErrorCode TSSetUp_Euler(TS ts) 134b33d51aSBarry Smith { 148b1af7b3SBarry Smith TS_Euler *euler = (TS_Euler*)ts->data; 15dfbe8321SBarry Smith PetscErrorCode ierr; 164b33d51aSBarry Smith 173a40ed3dSBarry Smith PetscFunctionBegin; 188b1af7b3SBarry Smith ierr = VecDuplicate(ts->vec_sol,&euler->update);CHKERRQ(ierr); 193a40ed3dSBarry Smith PetscFunctionReturn(0); 204b33d51aSBarry Smith } 214b33d51aSBarry Smith 224a2ae208SSatish Balay #undef __FUNCT__ 234a2ae208SSatish Balay #define __FUNCT__ "TSStep_Euler" 24*6849ba73SBarry Smith static PetscErrorCode TSStep_Euler(TS ts,int *steps,PetscReal *ptime) 254b33d51aSBarry Smith { 268b1af7b3SBarry Smith TS_Euler *euler = (TS_Euler*)ts->data; 278b1af7b3SBarry Smith Vec sol = ts->vec_sol,update = euler->update; 28dfbe8321SBarry Smith PetscErrorCode ierr; 29dfbe8321SBarry Smith int i,max_steps = ts->max_steps; 30ea709b57SSatish Balay PetscScalar dt = ts->time_step; 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++) { 378b1af7b3SBarry Smith ts->ptime += ts->time_step; 38c3e30b67SBarry Smith ierr = TSComputeRHSFunction(ts,ts->ptime,sol,update);CHKERRQ(ierr); 39c3e30b67SBarry Smith ierr = VecAXPY(&dt,update,sol);CHKERRQ(ierr); 408b1af7b3SBarry Smith ts->steps++; 41c3e30b67SBarry Smith ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 428b1af7b3SBarry Smith if (ts->ptime > ts->max_time) break; 438b1af7b3SBarry Smith } 444b33d51aSBarry Smith 458b1af7b3SBarry Smith *steps += ts->steps; 46142b95e3SSatish Balay *ptime = ts->ptime; 473a40ed3dSBarry Smith PetscFunctionReturn(0); 484b33d51aSBarry Smith } 494b33d51aSBarry Smith /*------------------------------------------------------------*/ 504a2ae208SSatish Balay #undef __FUNCT__ 514a2ae208SSatish Balay #define __FUNCT__ "TSDestroy_Euler" 52*6849ba73SBarry Smith static PetscErrorCode TSDestroy_Euler(TS ts) 534b33d51aSBarry Smith { 548b1af7b3SBarry Smith TS_Euler *euler = (TS_Euler*)ts->data; 55dfbe8321SBarry Smith PetscErrorCode ierr; 568b1af7b3SBarry Smith 573a40ed3dSBarry Smith PetscFunctionBegin; 581713a123SBarry Smith if (euler->update) {ierr = VecDestroy(euler->update);CHKERRQ(ierr);} 59606d414cSSatish Balay ierr = PetscFree(euler);CHKERRQ(ierr); 603a40ed3dSBarry Smith PetscFunctionReturn(0); 618b1af7b3SBarry Smith } 628b1af7b3SBarry Smith /*------------------------------------------------------------*/ 638b1af7b3SBarry Smith 644a2ae208SSatish Balay #undef __FUNCT__ 654a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Euler" 66*6849ba73SBarry Smith static PetscErrorCode TSSetFromOptions_Euler(TS ts) 678b1af7b3SBarry Smith { 683a40ed3dSBarry Smith PetscFunctionBegin; 693a40ed3dSBarry Smith PetscFunctionReturn(0); 704b33d51aSBarry Smith } 714b33d51aSBarry Smith 724a2ae208SSatish Balay #undef __FUNCT__ 734a2ae208SSatish Balay #define __FUNCT__ "TSView_Euler" 74*6849ba73SBarry Smith static PetscErrorCode TSView_Euler(TS ts,PetscViewer viewer) 758b1af7b3SBarry Smith { 763a40ed3dSBarry Smith PetscFunctionBegin; 773a40ed3dSBarry Smith PetscFunctionReturn(0); 788b1af7b3SBarry Smith } 798b1af7b3SBarry Smith 808b1af7b3SBarry Smith /* ------------------------------------------------------------ */ 81ebe3b25bSBarry Smith 82ebe3b25bSBarry Smith /*MC 83ebe3b25bSBarry Smith TS_EULER - ODE solver using the explicit forward Euler method 84ebe3b25bSBarry Smith 85d41222bbSBarry Smith Level: beginner 86d41222bbSBarry Smith 87ebe3b25bSBarry Smith .seealso: TSCreate(), TS, TSSetType(), TS_BEULER 88ebe3b25bSBarry Smith 89ebe3b25bSBarry Smith M*/ 90fb2e594dSBarry Smith EXTERN_C_BEGIN 914a2ae208SSatish Balay #undef __FUNCT__ 924a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Euler" 93dfbe8321SBarry Smith PetscErrorCode TSCreate_Euler(TS ts) 948b1af7b3SBarry Smith { 958b1af7b3SBarry Smith TS_Euler *euler; 96dfbe8321SBarry Smith PetscErrorCode ierr; 978b1af7b3SBarry Smith 983a40ed3dSBarry Smith PetscFunctionBegin; 99000e7ae3SMatthew Knepley ts->ops->setup = TSSetUp_Euler; 100000e7ae3SMatthew Knepley ts->ops->step = TSStep_Euler; 101000e7ae3SMatthew Knepley ts->ops->destroy = TSDestroy_Euler; 102000e7ae3SMatthew Knepley ts->ops->setfromoptions = TSSetFromOptions_Euler; 103000e7ae3SMatthew Knepley ts->ops->view = TSView_Euler; 1048b1af7b3SBarry Smith 105b0a32e0cSBarry Smith ierr = PetscNew(TS_Euler,&euler);CHKERRQ(ierr); 106b0a32e0cSBarry Smith PetscLogObjectMemory(ts,sizeof(TS_Euler)); 1078b1af7b3SBarry Smith ts->data = (void*)euler; 1084b33d51aSBarry Smith 1093a40ed3dSBarry Smith PetscFunctionReturn(0); 1104b33d51aSBarry Smith } 111fb2e594dSBarry Smith EXTERN_C_END 112c3e30b67SBarry Smith 113c3e30b67SBarry Smith 114c3e30b67SBarry Smith 115c3e30b67SBarry Smith 116