14b33d51aSBarry Smith /* 2fb4a63b6SLois Curfman McInnes Code for Timestepping with explicit Euler. 34b33d51aSBarry Smith */ 4af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 54b33d51aSBarry Smith 68b1af7b3SBarry Smith typedef struct { 7277b19d0SLisandro Dalcin Vec update; /* work vector where new solution is formed */ 88b1af7b3SBarry Smith } TS_Euler; 98b1af7b3SBarry Smith 104a2ae208SSatish Balay #undef __FUNCT__ 114a2ae208SSatish Balay #define __FUNCT__ "TSStep_Euler" 12193ac0bcSJed Brown static PetscErrorCode TSStep_Euler(TS ts) 134b33d51aSBarry Smith { 148b1af7b3SBarry Smith TS_Euler *euler = (TS_Euler*)ts->data; 158b1af7b3SBarry Smith Vec sol = ts->vec_sol,update = euler->update; 16277b19d0SLisandro Dalcin PetscErrorCode ierr; 174b33d51aSBarry Smith 183a40ed3dSBarry Smith PetscFunctionBegin; 19b8123daeSJed Brown ierr = TSPreStep(ts);CHKERRQ(ierr); 20b8123daeSJed Brown ierr = TSPreStage(ts,ts->ptime);CHKERRQ(ierr); 21c3e30b67SBarry Smith ierr = TSComputeRHSFunction(ts,ts->ptime,sol,update);CHKERRQ(ierr); 22186e87acSLisandro Dalcin ierr = VecAXPY(sol,ts->time_step,update);CHKERRQ(ierr); 239be3e283SDebojyoti Ghosh ierr = TSPostStage(ts,ts->ptime,0,&sol);CHKERRQ(ierr); 24186e87acSLisandro Dalcin ts->ptime += ts->time_step; 258b1af7b3SBarry Smith ts->steps++; 263a40ed3dSBarry Smith PetscFunctionReturn(0); 274b33d51aSBarry Smith } 284b33d51aSBarry Smith /*------------------------------------------------------------*/ 29277b19d0SLisandro Dalcin 304a2ae208SSatish Balay #undef __FUNCT__ 31277b19d0SLisandro Dalcin #define __FUNCT__ "TSSetUp_Euler" 32277b19d0SLisandro Dalcin static PetscErrorCode TSSetUp_Euler(TS ts) 33277b19d0SLisandro Dalcin { 34277b19d0SLisandro Dalcin TS_Euler *euler = (TS_Euler*)ts->data; 35277b19d0SLisandro Dalcin PetscErrorCode ierr; 36277b19d0SLisandro Dalcin 37277b19d0SLisandro Dalcin PetscFunctionBegin; 38277b19d0SLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&euler->update);CHKERRQ(ierr); 39277b19d0SLisandro Dalcin PetscFunctionReturn(0); 40277b19d0SLisandro Dalcin } 41277b19d0SLisandro Dalcin 42277b19d0SLisandro Dalcin #undef __FUNCT__ 43277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Euler" 44277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Euler(TS ts) 454b33d51aSBarry Smith { 468b1af7b3SBarry Smith TS_Euler *euler = (TS_Euler*)ts->data; 47dfbe8321SBarry Smith PetscErrorCode ierr; 488b1af7b3SBarry Smith 493a40ed3dSBarry Smith PetscFunctionBegin; 506bf464f9SBarry Smith ierr = VecDestroy(&euler->update);CHKERRQ(ierr); 51277b19d0SLisandro Dalcin PetscFunctionReturn(0); 52277b19d0SLisandro Dalcin } 53277b19d0SLisandro Dalcin 54277b19d0SLisandro Dalcin #undef __FUNCT__ 55277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Euler" 56277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Euler(TS ts) 57277b19d0SLisandro Dalcin { 58277b19d0SLisandro Dalcin PetscErrorCode ierr; 59277b19d0SLisandro Dalcin 60277b19d0SLisandro Dalcin PetscFunctionBegin; 61277b19d0SLisandro Dalcin ierr = TSReset_Euler(ts);CHKERRQ(ierr); 62277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 633a40ed3dSBarry Smith PetscFunctionReturn(0); 648b1af7b3SBarry Smith } 658b1af7b3SBarry Smith /*------------------------------------------------------------*/ 668b1af7b3SBarry Smith 674a2ae208SSatish Balay #undef __FUNCT__ 684a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Euler" 698c34d3f5SBarry Smith static PetscErrorCode TSSetFromOptions_Euler(PetscOptions *PetscOptionsObject,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 836083293cSBarry Smith #undef __FUNCT__ 846083293cSBarry Smith #define __FUNCT__ "TSInterpolate_Euler" 856083293cSBarry Smith static PetscErrorCode TSInterpolate_Euler(TS ts,PetscReal t,Vec X) 866083293cSBarry Smith { 87*3354212dSEmil Constantinescu TS_Euler *euler = (TS_Euler*)ts->data; 88*3354212dSEmil Constantinescu Vec update = euler->update; 896083293cSBarry Smith PetscReal alpha = (ts->ptime - t)/ts->time_step; 906083293cSBarry Smith PetscErrorCode ierr; 916083293cSBarry Smith 926083293cSBarry Smith PetscFunctionBegin; 93*3354212dSEmil Constantinescu ierr = VecWAXPY(X,-ts->time_step,update,ts->vec_sol);CHKERRQ(ierr); 94*3354212dSEmil Constantinescu ierr = VecAXPBY(X,1.0-alpha,alpha,ts->vec_sol);CHKERRQ(ierr); 956083293cSBarry Smith PetscFunctionReturn(0); 966083293cSBarry Smith } 976083293cSBarry Smith 98f9c1d6abSBarry Smith #undef __FUNCT__ 99f9c1d6abSBarry Smith #define __FUNCT__ "TSComputeLinearStability_Euler" 100f9c1d6abSBarry Smith PetscErrorCode TSComputeLinearStability_Euler(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi) 101f9c1d6abSBarry Smith { 102f9c1d6abSBarry Smith PetscFunctionBegin; 103f9c1d6abSBarry Smith *yr = 1.0 + xr; 104f9c1d6abSBarry Smith *yi = xi; 105f9c1d6abSBarry Smith PetscFunctionReturn(0); 106f9c1d6abSBarry Smith } 1078b1af7b3SBarry Smith /* ------------------------------------------------------------ */ 108ebe3b25bSBarry Smith 109ebe3b25bSBarry Smith /*MC 1109596e0b4SJed Brown TSEULER - ODE solver using the explicit forward Euler method 111ebe3b25bSBarry Smith 112d41222bbSBarry Smith Level: beginner 113d41222bbSBarry Smith 1149596e0b4SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSBEULER 115ebe3b25bSBarry Smith 116ebe3b25bSBarry Smith M*/ 1174a2ae208SSatish Balay #undef __FUNCT__ 1184a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Euler" 1198cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Euler(TS ts) 1208b1af7b3SBarry Smith { 1218b1af7b3SBarry Smith TS_Euler *euler; 122dfbe8321SBarry Smith PetscErrorCode ierr; 1238b1af7b3SBarry Smith 1243a40ed3dSBarry Smith PetscFunctionBegin; 125000e7ae3SMatthew Knepley ts->ops->setup = TSSetUp_Euler; 126000e7ae3SMatthew Knepley ts->ops->step = TSStep_Euler; 127277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Euler; 128000e7ae3SMatthew Knepley ts->ops->destroy = TSDestroy_Euler; 129000e7ae3SMatthew Knepley ts->ops->setfromoptions = TSSetFromOptions_Euler; 130000e7ae3SMatthew Knepley ts->ops->view = TSView_Euler; 1316083293cSBarry Smith ts->ops->interpolate = TSInterpolate_Euler; 132f9c1d6abSBarry Smith ts->ops->linearstability = TSComputeLinearStability_Euler; 1338b1af7b3SBarry Smith 134b00a9115SJed Brown ierr = PetscNewLog(ts,&euler);CHKERRQ(ierr); 1358b1af7b3SBarry Smith ts->data = (void*)euler; 1363a40ed3dSBarry Smith PetscFunctionReturn(0); 1374b33d51aSBarry Smith } 138