173f4d377SMatthew Knepley /*$Id: euler.c,v 1.30 2001/08/07 03:04:22 balay Exp $*/ 24b33d51aSBarry Smith /* 3fb4a63b6SLois Curfman McInnes Code for Timestepping with explicit Euler. 44b33d51aSBarry Smith */ 5e090d566SSatish Balay #include "src/ts/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" 138b1af7b3SBarry Smith static int TSSetUp_Euler(TS ts) 144b33d51aSBarry Smith { 158b1af7b3SBarry Smith TS_Euler *euler = (TS_Euler*)ts->data; 164b33d51aSBarry Smith int 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" 2587828ca2SBarry Smith static int TSStep_Euler(TS ts,int *steps,PetscReal *ptime) 264b33d51aSBarry Smith { 278b1af7b3SBarry Smith TS_Euler *euler = (TS_Euler*)ts->data; 288b1af7b3SBarry Smith Vec sol = ts->vec_sol,update = euler->update; 298b1af7b3SBarry Smith int ierr,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" 52e1311b90SBarry Smith static int TSDestroy_Euler(TS ts) 534b33d51aSBarry Smith { 548b1af7b3SBarry Smith TS_Euler *euler = (TS_Euler*)ts->data; 55606d414cSSatish Balay int 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" 668b1af7b3SBarry Smith static int 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" 74b0a32e0cSBarry Smith static int TSView_Euler(TS ts,PetscViewer viewer) 758b1af7b3SBarry Smith { 763a40ed3dSBarry Smith PetscFunctionBegin; 773a40ed3dSBarry Smith PetscFunctionReturn(0); 788b1af7b3SBarry Smith } 798b1af7b3SBarry Smith 808b1af7b3SBarry Smith /* ------------------------------------------------------------ */ 81fb2e594dSBarry Smith EXTERN_C_BEGIN 824a2ae208SSatish Balay #undef __FUNCT__ 834a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Euler" 848b1af7b3SBarry Smith int TSCreate_Euler(TS ts) 858b1af7b3SBarry Smith { 868b1af7b3SBarry Smith TS_Euler *euler; 87b0a32e0cSBarry Smith int ierr; 888b1af7b3SBarry Smith 893a40ed3dSBarry Smith PetscFunctionBegin; 90*000e7ae3SMatthew Knepley ts->ops->setup = TSSetUp_Euler; 91*000e7ae3SMatthew Knepley ts->ops->step = TSStep_Euler; 92*000e7ae3SMatthew Knepley ts->ops->destroy = TSDestroy_Euler; 93*000e7ae3SMatthew Knepley ts->ops->setfromoptions = TSSetFromOptions_Euler; 94*000e7ae3SMatthew Knepley ts->ops->view = TSView_Euler; 958b1af7b3SBarry Smith 96b0a32e0cSBarry Smith ierr = PetscNew(TS_Euler,&euler);CHKERRQ(ierr); 97b0a32e0cSBarry Smith PetscLogObjectMemory(ts,sizeof(TS_Euler)); 988b1af7b3SBarry Smith ts->data = (void*)euler; 994b33d51aSBarry Smith 1003a40ed3dSBarry Smith PetscFunctionReturn(0); 1014b33d51aSBarry Smith } 102fb2e594dSBarry Smith EXTERN_C_END 103c3e30b67SBarry Smith 104c3e30b67SBarry Smith 105c3e30b67SBarry Smith 106c3e30b67SBarry Smith 107