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 solution = ts->vec_sol,update = euler->update; 16 PetscBool stageok,accept = PETSC_TRUE; 17 PetscReal next_time_step = ts->time_step; 18 PetscErrorCode ierr; 19 20 PetscFunctionBegin; 21 ierr = TSPreStage(ts,ts->ptime);CHKERRQ(ierr); 22 ierr = TSComputeRHSFunction(ts,ts->ptime,solution,update);CHKERRQ(ierr); 23 ierr = VecAYPX(update,ts->time_step,solution);CHKERRQ(ierr); 24 ierr = TSPostStage(ts,ts->ptime,0,&solution);CHKERRQ(ierr); 25 ierr = TSAdaptCheckStage(ts->adapt,ts,ts->ptime,solution,&stageok);CHKERRQ(ierr); 26 if(!stageok) {ts->reason = TS_DIVERGED_STEP_REJECTED; PetscFunctionReturn(0);} 27 ierr = TSFunctionDomainError(ts,ts->ptime+ts->time_step,update,&stageok);CHKERRQ(ierr); 28 if(!stageok) {ts->reason = TS_DIVERGED_STEP_REJECTED; PetscFunctionReturn(0);} 29 30 ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 31 if (!accept) {ts->reason = TS_DIVERGED_STEP_REJECTED; PetscFunctionReturn(0);} 32 ierr = VecCopy(update,solution);CHKERRQ(ierr); 33 34 ts->ptime += ts->time_step; 35 ts->time_step = next_time_step; 36 PetscFunctionReturn(0); 37 } 38 /*------------------------------------------------------------*/ 39 40 #undef __FUNCT__ 41 #define __FUNCT__ "TSSetUp_Euler" 42 static PetscErrorCode TSSetUp_Euler(TS ts) 43 { 44 TS_Euler *euler = (TS_Euler*)ts->data; 45 PetscErrorCode ierr; 46 TSRHSFunction rhsfunction; 47 TSIFunction ifunction; 48 49 PetscFunctionBegin; 50 ierr = TSGetIFunction(ts,NULL,&ifunction,NULL);CHKERRQ(ierr); 51 ierr = TSGetRHSFunction(ts,NULL,&rhsfunction,NULL);CHKERRQ(ierr); 52 if (!rhsfunction || ifunction) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must define RHSFunction() and leave IFunction() undefined in order to use -ts_type euler"); 53 ierr = VecDuplicate(ts->vec_sol,&euler->update);CHKERRQ(ierr); 54 55 ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr); 56 ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr); 57 PetscFunctionReturn(0); 58 } 59 60 #undef __FUNCT__ 61 #define __FUNCT__ "TSReset_Euler" 62 static PetscErrorCode TSReset_Euler(TS ts) 63 { 64 TS_Euler *euler = (TS_Euler*)ts->data; 65 PetscErrorCode ierr; 66 67 PetscFunctionBegin; 68 ierr = VecDestroy(&euler->update);CHKERRQ(ierr); 69 PetscFunctionReturn(0); 70 } 71 72 #undef __FUNCT__ 73 #define __FUNCT__ "TSDestroy_Euler" 74 static PetscErrorCode TSDestroy_Euler(TS ts) 75 { 76 PetscErrorCode ierr; 77 78 PetscFunctionBegin; 79 ierr = TSReset_Euler(ts);CHKERRQ(ierr); 80 ierr = PetscFree(ts->data);CHKERRQ(ierr); 81 PetscFunctionReturn(0); 82 } 83 /*------------------------------------------------------------*/ 84 85 #undef __FUNCT__ 86 #define __FUNCT__ "TSSetFromOptions_Euler" 87 static PetscErrorCode TSSetFromOptions_Euler(PetscOptionItems *PetscOptionsObject,TS ts) 88 { 89 PetscFunctionBegin; 90 PetscFunctionReturn(0); 91 } 92 93 #undef __FUNCT__ 94 #define __FUNCT__ "TSView_Euler" 95 static PetscErrorCode TSView_Euler(TS ts,PetscViewer viewer) 96 { 97 PetscFunctionBegin; 98 PetscFunctionReturn(0); 99 } 100 101 #undef __FUNCT__ 102 #define __FUNCT__ "TSInterpolate_Euler" 103 static PetscErrorCode TSInterpolate_Euler(TS ts,PetscReal t,Vec X) 104 { 105 TS_Euler *euler = (TS_Euler*)ts->data; 106 Vec update = euler->update; 107 PetscReal alpha = (ts->ptime - t)/ts->time_step; 108 PetscErrorCode ierr; 109 110 PetscFunctionBegin; 111 ierr = VecWAXPY(X,-ts->time_step,update,ts->vec_sol);CHKERRQ(ierr); 112 ierr = VecAXPBY(X,1.0-alpha,alpha,ts->vec_sol);CHKERRQ(ierr); 113 PetscFunctionReturn(0); 114 } 115 116 #undef __FUNCT__ 117 #define __FUNCT__ "TSComputeLinearStability_Euler" 118 static PetscErrorCode TSComputeLinearStability_Euler(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi) 119 { 120 PetscFunctionBegin; 121 *yr = 1.0 + xr; 122 *yi = xi; 123 PetscFunctionReturn(0); 124 } 125 /* ------------------------------------------------------------ */ 126 127 /*MC 128 TSEULER - ODE solver using the explicit forward Euler method 129 130 Level: beginner 131 132 .seealso: TSCreate(), TS, TSSetType(), TSBEULER 133 134 M*/ 135 #undef __FUNCT__ 136 #define __FUNCT__ "TSCreate_Euler" 137 PETSC_EXTERN PetscErrorCode TSCreate_Euler(TS ts) 138 { 139 TS_Euler *euler; 140 PetscErrorCode ierr; 141 142 PetscFunctionBegin; 143 ierr = PetscNewLog(ts,&euler);CHKERRQ(ierr); 144 ts->data = (void*)euler; 145 146 ts->ops->setup = TSSetUp_Euler; 147 ts->ops->step = TSStep_Euler; 148 ts->ops->reset = TSReset_Euler; 149 ts->ops->destroy = TSDestroy_Euler; 150 ts->ops->setfromoptions = TSSetFromOptions_Euler; 151 ts->ops->view = TSView_Euler; 152 ts->ops->interpolate = TSInterpolate_Euler; 153 ts->ops->linearstability = TSComputeLinearStability_Euler; 154 PetscFunctionReturn(0); 155 } 156