1 /* 2 Code for Timestepping with explicit Euler. 3 */ 4 #include <petsc-private/tsimpl.h> /*I "petscts.h" I*/ 5 6 typedef struct { 7 Vec update; /* work vector where new solution is formed */ 8 } TS_Euler; 9 10 #undef __FUNCT__ 11 #define __FUNCT__ "TSStep_Euler" 12 static PetscErrorCode TSStep_Euler(TS ts) 13 { 14 TS_Euler *euler = (TS_Euler*)ts->data; 15 Vec sol = ts->vec_sol,update = euler->update; 16 PetscErrorCode ierr; 17 18 PetscFunctionBegin; 19 ierr = TSPreStep(ts);CHKERRQ(ierr); 20 ierr = TSPreStage(ts,ts->ptime);CHKERRQ(ierr); 21 ierr = TSComputeRHSFunction(ts,ts->ptime,sol,update);CHKERRQ(ierr); 22 ierr = VecAXPY(sol,ts->time_step,update);CHKERRQ(ierr); 23 ts->ptime += ts->time_step; 24 ts->steps++; 25 PetscFunctionReturn(0); 26 } 27 /*------------------------------------------------------------*/ 28 29 #undef __FUNCT__ 30 #define __FUNCT__ "TSSetUp_Euler" 31 static PetscErrorCode TSSetUp_Euler(TS ts) 32 { 33 TS_Euler *euler = (TS_Euler*)ts->data; 34 PetscErrorCode ierr; 35 36 PetscFunctionBegin; 37 ierr = VecDuplicate(ts->vec_sol,&euler->update);CHKERRQ(ierr); 38 PetscFunctionReturn(0); 39 } 40 41 #undef __FUNCT__ 42 #define __FUNCT__ "TSReset_Euler" 43 static PetscErrorCode TSReset_Euler(TS ts) 44 { 45 TS_Euler *euler = (TS_Euler*)ts->data; 46 PetscErrorCode ierr; 47 48 PetscFunctionBegin; 49 ierr = VecDestroy(&euler->update);CHKERRQ(ierr); 50 PetscFunctionReturn(0); 51 } 52 53 #undef __FUNCT__ 54 #define __FUNCT__ "TSDestroy_Euler" 55 static PetscErrorCode TSDestroy_Euler(TS ts) 56 { 57 PetscErrorCode ierr; 58 59 PetscFunctionBegin; 60 ierr = TSReset_Euler(ts);CHKERRQ(ierr); 61 ierr = PetscFree(ts->data);CHKERRQ(ierr); 62 PetscFunctionReturn(0); 63 } 64 /*------------------------------------------------------------*/ 65 66 #undef __FUNCT__ 67 #define __FUNCT__ "TSSetFromOptions_Euler" 68 static PetscErrorCode TSSetFromOptions_Euler(TS ts) 69 { 70 PetscFunctionBegin; 71 PetscFunctionReturn(0); 72 } 73 74 #undef __FUNCT__ 75 #define __FUNCT__ "TSView_Euler" 76 static PetscErrorCode TSView_Euler(TS ts,PetscViewer viewer) 77 { 78 PetscFunctionBegin; 79 PetscFunctionReturn(0); 80 } 81 82 #undef __FUNCT__ 83 #define __FUNCT__ "TSInterpolate_Euler" 84 static PetscErrorCode TSInterpolate_Euler(TS ts,PetscReal t,Vec X) 85 { 86 PetscReal alpha = (ts->ptime - t)/ts->time_step; 87 PetscErrorCode ierr; 88 89 PetscFunctionBegin; 90 ierr = VecAXPBY(ts->vec_sol,1.0-alpha,alpha,X);CHKERRQ(ierr); 91 PetscFunctionReturn(0); 92 } 93 94 #undef __FUNCT__ 95 #define __FUNCT__ "TSComputeLinearStability_Euler" 96 PetscErrorCode TSComputeLinearStability_Euler(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi) 97 { 98 PetscFunctionBegin; 99 *yr = 1.0 + xr; 100 *yi = xi; 101 PetscFunctionReturn(0); 102 } 103 /* ------------------------------------------------------------ */ 104 105 /*MC 106 TSEULER - ODE solver using the explicit forward Euler method 107 108 Level: beginner 109 110 .seealso: TSCreate(), TS, TSSetType(), TSBEULER 111 112 M*/ 113 EXTERN_C_BEGIN 114 #undef __FUNCT__ 115 #define __FUNCT__ "TSCreate_Euler" 116 PetscErrorCode TSCreate_Euler(TS ts) 117 { 118 TS_Euler *euler; 119 PetscErrorCode ierr; 120 121 PetscFunctionBegin; 122 ts->ops->setup = TSSetUp_Euler; 123 ts->ops->step = TSStep_Euler; 124 ts->ops->reset = TSReset_Euler; 125 ts->ops->destroy = TSDestroy_Euler; 126 ts->ops->setfromoptions = TSSetFromOptions_Euler; 127 ts->ops->view = TSView_Euler; 128 ts->ops->interpolate = TSInterpolate_Euler; 129 ts->ops->linearstability = TSComputeLinearStability_Euler; 130 131 ierr = PetscNewLog(ts,TS_Euler,&euler);CHKERRQ(ierr); 132 ts->data = (void*)euler; 133 PetscFunctionReturn(0); 134 } 135 EXTERN_C_END 136