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; 17cb9d8021SPierre Barbier de Reuille PetscBool accept; 184b33d51aSBarry Smith 193a40ed3dSBarry Smith PetscFunctionBegin; 20b8123daeSJed Brown ierr = TSPreStep(ts);CHKERRQ(ierr); 21b8123daeSJed Brown ierr = TSPreStage(ts,ts->ptime);CHKERRQ(ierr); 22c3e30b67SBarry Smith ierr = TSComputeRHSFunction(ts,ts->ptime,sol,update);CHKERRQ(ierr); 23186e87acSLisandro Dalcin ierr = VecAXPY(sol,ts->time_step,update);CHKERRQ(ierr); 24d183316bSPierre Barbier de Reuille ierr = TSFunctionDomainError(ts,ts->ptime,sol,&accept);CHKERRQ(ierr); 25cb9d8021SPierre Barbier de Reuille if(!accept) { 26cb9d8021SPierre Barbier de Reuille ts->reason = TS_DIVERGED_STEP_REJECTED; 27cb9d8021SPierre Barbier de Reuille PetscFunctionReturn(0); 28cb9d8021SPierre Barbier de Reuille } 299be3e283SDebojyoti Ghosh ierr = TSPostStage(ts,ts->ptime,0,&sol);CHKERRQ(ierr); 30186e87acSLisandro Dalcin ts->ptime += ts->time_step; 318b1af7b3SBarry Smith ts->steps++; 323a40ed3dSBarry Smith PetscFunctionReturn(0); 334b33d51aSBarry Smith } 344b33d51aSBarry Smith /*------------------------------------------------------------*/ 35277b19d0SLisandro Dalcin 364a2ae208SSatish Balay #undef __FUNCT__ 37277b19d0SLisandro Dalcin #define __FUNCT__ "TSSetUp_Euler" 38277b19d0SLisandro Dalcin static PetscErrorCode TSSetUp_Euler(TS ts) 39277b19d0SLisandro Dalcin { 40277b19d0SLisandro Dalcin TS_Euler *euler = (TS_Euler*)ts->data; 41277b19d0SLisandro Dalcin PetscErrorCode ierr; 42*1dca01adSEmil Constantinescu TSRHSFunction rhsfunction; 43*1dca01adSEmil Constantinescu TSIFunction ifunction; 44277b19d0SLisandro Dalcin 45277b19d0SLisandro Dalcin PetscFunctionBegin; 46*1dca01adSEmil Constantinescu ierr = TSGetIFunction(ts,NULL,&ifunction,NULL);CHKERRQ(ierr); 47*1dca01adSEmil Constantinescu ierr = TSGetRHSFunction(ts,NULL,&rhsfunction,NULL);CHKERRQ(ierr); 48*1dca01adSEmil Constantinescu if (!rhsfunction || ifunction!=0) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must define RHSFunction() and leave IFunction() undefined in order to use -ts_type euler"); 49277b19d0SLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&euler->update);CHKERRQ(ierr); 50277b19d0SLisandro Dalcin PetscFunctionReturn(0); 51277b19d0SLisandro Dalcin } 52277b19d0SLisandro Dalcin 53277b19d0SLisandro Dalcin #undef __FUNCT__ 54277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Euler" 55277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Euler(TS ts) 564b33d51aSBarry Smith { 578b1af7b3SBarry Smith TS_Euler *euler = (TS_Euler*)ts->data; 58dfbe8321SBarry Smith PetscErrorCode ierr; 598b1af7b3SBarry Smith 603a40ed3dSBarry Smith PetscFunctionBegin; 616bf464f9SBarry Smith ierr = VecDestroy(&euler->update);CHKERRQ(ierr); 62277b19d0SLisandro Dalcin PetscFunctionReturn(0); 63277b19d0SLisandro Dalcin } 64277b19d0SLisandro Dalcin 65277b19d0SLisandro Dalcin #undef __FUNCT__ 66277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Euler" 67277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Euler(TS ts) 68277b19d0SLisandro Dalcin { 69277b19d0SLisandro Dalcin PetscErrorCode ierr; 70277b19d0SLisandro Dalcin 71277b19d0SLisandro Dalcin PetscFunctionBegin; 72277b19d0SLisandro Dalcin ierr = TSReset_Euler(ts);CHKERRQ(ierr); 73277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 743a40ed3dSBarry Smith PetscFunctionReturn(0); 758b1af7b3SBarry Smith } 768b1af7b3SBarry Smith /*------------------------------------------------------------*/ 778b1af7b3SBarry Smith 784a2ae208SSatish Balay #undef __FUNCT__ 794a2ae208SSatish Balay #define __FUNCT__ "TSSetFromOptions_Euler" 808c34d3f5SBarry Smith static PetscErrorCode TSSetFromOptions_Euler(PetscOptions *PetscOptionsObject,TS ts) 818b1af7b3SBarry Smith { 823a40ed3dSBarry Smith PetscFunctionBegin; 833a40ed3dSBarry Smith PetscFunctionReturn(0); 844b33d51aSBarry Smith } 854b33d51aSBarry Smith 864a2ae208SSatish Balay #undef __FUNCT__ 874a2ae208SSatish Balay #define __FUNCT__ "TSView_Euler" 886849ba73SBarry Smith static PetscErrorCode TSView_Euler(TS ts,PetscViewer viewer) 898b1af7b3SBarry Smith { 903a40ed3dSBarry Smith PetscFunctionBegin; 913a40ed3dSBarry Smith PetscFunctionReturn(0); 928b1af7b3SBarry Smith } 938b1af7b3SBarry Smith 946083293cSBarry Smith #undef __FUNCT__ 956083293cSBarry Smith #define __FUNCT__ "TSInterpolate_Euler" 966083293cSBarry Smith static PetscErrorCode TSInterpolate_Euler(TS ts,PetscReal t,Vec X) 976083293cSBarry Smith { 986083293cSBarry Smith PetscReal alpha = (ts->ptime - t)/ts->time_step; 996083293cSBarry Smith PetscErrorCode ierr; 1006083293cSBarry Smith 1016083293cSBarry Smith PetscFunctionBegin; 1026083293cSBarry Smith ierr = VecAXPBY(ts->vec_sol,1.0-alpha,alpha,X);CHKERRQ(ierr); 1036083293cSBarry Smith PetscFunctionReturn(0); 1046083293cSBarry Smith } 1056083293cSBarry Smith 106f9c1d6abSBarry Smith #undef __FUNCT__ 107f9c1d6abSBarry Smith #define __FUNCT__ "TSComputeLinearStability_Euler" 108f9c1d6abSBarry Smith PetscErrorCode TSComputeLinearStability_Euler(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi) 109f9c1d6abSBarry Smith { 110f9c1d6abSBarry Smith PetscFunctionBegin; 111f9c1d6abSBarry Smith *yr = 1.0 + xr; 112f9c1d6abSBarry Smith *yi = xi; 113f9c1d6abSBarry Smith PetscFunctionReturn(0); 114f9c1d6abSBarry Smith } 1158b1af7b3SBarry Smith /* ------------------------------------------------------------ */ 116ebe3b25bSBarry Smith 117ebe3b25bSBarry Smith /*MC 1189596e0b4SJed Brown TSEULER - ODE solver using the explicit forward Euler method 119ebe3b25bSBarry Smith 120d41222bbSBarry Smith Level: beginner 121d41222bbSBarry Smith 1229596e0b4SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSBEULER 123ebe3b25bSBarry Smith 124ebe3b25bSBarry Smith M*/ 1254a2ae208SSatish Balay #undef __FUNCT__ 1264a2ae208SSatish Balay #define __FUNCT__ "TSCreate_Euler" 1278cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Euler(TS ts) 1288b1af7b3SBarry Smith { 1298b1af7b3SBarry Smith TS_Euler *euler; 130dfbe8321SBarry Smith PetscErrorCode ierr; 1318b1af7b3SBarry Smith 1323a40ed3dSBarry Smith PetscFunctionBegin; 133000e7ae3SMatthew Knepley ts->ops->setup = TSSetUp_Euler; 134000e7ae3SMatthew Knepley ts->ops->step = TSStep_Euler; 135277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Euler; 136000e7ae3SMatthew Knepley ts->ops->destroy = TSDestroy_Euler; 137000e7ae3SMatthew Knepley ts->ops->setfromoptions = TSSetFromOptions_Euler; 138000e7ae3SMatthew Knepley ts->ops->view = TSView_Euler; 1396083293cSBarry Smith ts->ops->interpolate = TSInterpolate_Euler; 140f9c1d6abSBarry Smith ts->ops->linearstability = TSComputeLinearStability_Euler; 1418b1af7b3SBarry Smith 142b00a9115SJed Brown ierr = PetscNewLog(ts,&euler);CHKERRQ(ierr); 1438b1af7b3SBarry Smith ts->data = (void*)euler; 1443a40ed3dSBarry Smith PetscFunctionReturn(0); 1454b33d51aSBarry Smith } 146