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