1 #define PETSCTS_DLL 2 3 /* 4 Code for timestepping with implicit Theta method 5 6 Notes: 7 This method can be applied to DAE. 8 9 This method is cast as a 1-stage implicit Runge-Kutta method. 10 11 Theta | Theta 12 ------------- 13 | 1 14 15 To apply a diagonally implicit RK method to DAE, the stage formula 16 17 X_i = x + h sum_j a_ij X'_j 18 19 is interpreted as a formula for X'_i in terms of X_i and known stuff (X'_j, j<i) 20 */ 21 #include "private/tsimpl.h" /*I "petscts.h" I*/ 22 23 typedef struct { 24 Vec X,Xdot; /* Storage for one stage */ 25 PetscTruth extrapolate; 26 PetscReal Theta; 27 PetscReal shift; 28 PetscReal stage_time; 29 } TS_Theta; 30 31 #undef __FUNCT__ 32 #define __FUNCT__ "TSStep_Theta" 33 static PetscErrorCode TSStep_Theta(TS ts,PetscInt *steps,PetscReal *ptime) 34 { 35 TS_Theta *th = (TS_Theta*)ts->data; 36 PetscInt i,max_steps = ts->max_steps,its,lits; 37 PetscErrorCode ierr; 38 39 PetscFunctionBegin; 40 *steps = -ts->steps; 41 *ptime = ts->ptime; 42 43 ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 44 45 for (i=0; i<max_steps; i++) { 46 if (ts->ptime + ts->time_step > ts->max_time) break; 47 ierr = TSPreStep(ts);CHKERRQ(ierr); 48 49 th->stage_time = ts->ptime + th->Theta*ts->time_step; 50 th->shift = 1./(th->Theta*ts->time_step); 51 52 if (th->extrapolate) { 53 ierr = VecWAXPY(th->X,1./th->shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr); 54 } else { 55 ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 56 } 57 ierr = SNESSolve(ts->snes,PETSC_NULL,th->X);CHKERRQ(ierr); 58 ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr); 59 ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 60 ts->nonlinear_its += its; ts->linear_its += lits; 61 62 ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,th->X);CHKERRQ(ierr); 63 ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr); 64 ts->ptime += ts->time_step; 65 ts->steps++; 66 67 ierr = TSPostStep(ts);CHKERRQ(ierr); 68 ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 69 } 70 71 *steps += ts->steps; 72 *ptime = ts->ptime; 73 PetscFunctionReturn(0); 74 } 75 76 /*------------------------------------------------------------*/ 77 #undef __FUNCT__ 78 #define __FUNCT__ "TSDestroy_Theta" 79 static PetscErrorCode TSDestroy_Theta(TS ts) 80 { 81 TS_Theta *th = (TS_Theta*)ts->data; 82 PetscErrorCode ierr; 83 84 PetscFunctionBegin; 85 if (th->X) {ierr = VecDestroy(th->X);CHKERRQ(ierr);} 86 if (th->Xdot) {ierr = VecDestroy(th->Xdot);CHKERRQ(ierr);} 87 ierr = PetscFree(th);CHKERRQ(ierr); 88 PetscFunctionReturn(0); 89 } 90 91 /* 92 This defines the nonlinear equation that is to be solved with SNES 93 G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 94 */ 95 #undef __FUNCT__ 96 #define __FUNCT__ "SNESTSFormFunction_Theta" 97 static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 98 { 99 TS_Theta *th = (TS_Theta*)ts->data; 100 PetscErrorCode ierr; 101 102 PetscFunctionBegin; 103 ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,x);CHKERRQ(ierr); 104 ierr = TSComputeIFunction(ts,th->stage_time,x,th->Xdot,y);CHKERRQ(ierr); 105 PetscFunctionReturn(0); 106 } 107 108 #undef __FUNCT__ 109 #define __FUNCT__ "SNESTSFormJacobian_Theta" 110 static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,TS ts) 111 { 112 TS_Theta *th = (TS_Theta*)ts->data; 113 PetscErrorCode ierr; 114 115 PetscFunctionBegin; 116 /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 117 ierr = TSComputeIJacobian(ts,th->stage_time,x,th->Xdot,th->shift,A,B,str);CHKERRQ(ierr); 118 PetscFunctionReturn(0); 119 } 120 121 122 #undef __FUNCT__ 123 #define __FUNCT__ "TSSetUp_Theta" 124 static PetscErrorCode TSSetUp_Theta(TS ts) 125 { 126 TS_Theta *th = (TS_Theta*)ts->data; 127 PetscErrorCode ierr; 128 Vec res; 129 130 PetscFunctionBegin; 131 if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Only for nonlinear problems"); 132 ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 133 ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 134 ierr = VecDuplicate(ts->vec_sol,&res);CHKERRQ(ierr); 135 ierr = SNESSetFunction(ts->snes,res,SNESTSFormFunction,ts);CHKERRQ(ierr); 136 ierr = VecDestroy(res);CHKERRQ(ierr); 137 /* This is nasty. SNESSetFromOptions() is usually called in TSSetFromOptions(). With -snes_mf_operator, it will 138 replace A and we don't want to mess with that. With -snes_mf, A and B will be replaced as well as the function and 139 context. Note that SNESSetFunction() normally has not been called before SNESSetFromOptions(), so when -snes_mf sets 140 the Jacobian user context to snes->funP, it will actually be NULL. This is not a problem because both snes->funP and 141 snes->jacP should be the TS. */ 142 { 143 Mat A,B; 144 PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*); 145 void *ctx; 146 ierr = SNESGetJacobian(ts->snes,&A,&B,&func,&ctx);CHKERRQ(ierr); 147 ierr = SNESSetJacobian(ts->snes,A?A:ts->A,B?B:ts->B,func?func:&SNESTSFormJacobian,ctx?ctx:ts);CHKERRQ(ierr); 148 } 149 PetscFunctionReturn(0); 150 } 151 /*------------------------------------------------------------*/ 152 153 #undef __FUNCT__ 154 #define __FUNCT__ "TSSetFromOptions_Theta" 155 static PetscErrorCode TSSetFromOptions_Theta(TS ts) 156 { 157 TS_Theta *th = (TS_Theta*)ts->data; 158 PetscErrorCode ierr; 159 160 PetscFunctionBegin; 161 ierr = PetscOptionsHead("Theta ODE solver options");CHKERRQ(ierr); 162 { 163 ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,PETSC_NULL);CHKERRQ(ierr); 164 ierr = PetscOptionsTruth("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,PETSC_NULL);CHKERRQ(ierr); 165 } 166 ierr = PetscOptionsTail();CHKERRQ(ierr); 167 PetscFunctionReturn(0); 168 } 169 170 #undef __FUNCT__ 171 #define __FUNCT__ "TSView_Theta" 172 static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 173 { 174 TS_Theta *th = (TS_Theta*)ts->data; 175 PetscTruth iascii; 176 PetscErrorCode ierr; 177 178 PetscFunctionBegin; 179 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 180 if (iascii) { 181 ierr = PetscViewerASCIIPrintf(viewer," Theta=%G\n",th->Theta);CHKERRQ(ierr); 182 ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate?"yes":"no");CHKERRQ(ierr); 183 } else { 184 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for TS_Theta",((PetscObject)viewer)->type_name); 185 } 186 PetscFunctionReturn(0); 187 } 188 189 EXTERN_C_BEGIN 190 #undef __FUNCT__ 191 #define __FUNCT__ "TSThetaGetTheta_Theta" 192 PetscErrorCode PETSCTS_DLLEXPORT TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 193 { 194 TS_Theta *th = (TS_Theta*)ts->data; 195 196 PetscFunctionBegin; 197 *theta = th->Theta; 198 PetscFunctionReturn(0); 199 } 200 201 #undef __FUNCT__ 202 #define __FUNCT__ "TSThetaSetTheta_Theta" 203 PetscErrorCode PETSCTS_DLLEXPORT TSThetaSetTheta_Theta(TS ts,PetscReal theta) 204 { 205 TS_Theta *th = (TS_Theta*)ts->data; 206 207 PetscFunctionBegin; 208 if (theta <= 0 || 1 < theta) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Theta %G not in range (0,1]",theta); 209 th->Theta = theta; 210 PetscFunctionReturn(0); 211 } 212 EXTERN_C_END 213 214 /* ------------------------------------------------------------ */ 215 /*MC 216 TSTHETA - DAE solver using the implicit Theta method 217 218 Level: beginner 219 220 .seealso: TSCreate(), TS, TSSetType() 221 222 M*/ 223 EXTERN_C_BEGIN 224 #undef __FUNCT__ 225 #define __FUNCT__ "TSCreate_Theta" 226 PetscErrorCode PETSCTS_DLLEXPORT TSCreate_Theta(TS ts) 227 { 228 TS_Theta *th; 229 PetscErrorCode ierr; 230 231 PetscFunctionBegin; 232 ts->ops->destroy = TSDestroy_Theta; 233 ts->ops->view = TSView_Theta; 234 ts->ops->setup = TSSetUp_Theta; 235 ts->ops->step = TSStep_Theta; 236 ts->ops->setfromoptions = TSSetFromOptions_Theta; 237 ts->ops->snesfunction = SNESTSFormFunction_Theta; 238 ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 239 240 ts->problem_type = TS_NONLINEAR; 241 ierr = SNESCreate(((PetscObject)ts)->comm,&ts->snes);CHKERRQ(ierr); 242 ierr = PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);CHKERRQ(ierr); 243 244 ierr = PetscNewLog(ts,TS_Theta,&th);CHKERRQ(ierr); 245 ts->data = (void*)th; 246 247 th->extrapolate = PETSC_FALSE; 248 th->Theta = 0.5; 249 250 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","TSThetaGetTheta_Theta",TSThetaGetTheta_Theta);CHKERRQ(ierr); 251 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","TSThetaSetTheta_Theta",TSThetaSetTheta_Theta);CHKERRQ(ierr); 252 PetscFunctionReturn(0); 253 } 254 EXTERN_C_END 255 256 #undef __FUNCT__ 257 #define __FUNCT__ "TSThetaGetTheta" 258 /*@ 259 TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 260 261 Not Collective 262 263 Input Parameter: 264 . ts - timestepping context 265 266 Output Parameter: 267 . theta - stage abscissa 268 269 Note: 270 Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 271 272 Level: Advanced 273 274 .seealso: TSThetaSetTheta() 275 @*/ 276 PetscErrorCode PETSCTS_DLLEXPORT TSThetaGetTheta(TS ts,PetscReal *theta) 277 { 278 PetscErrorCode ierr,(*f)(TS,PetscReal*); 279 280 PetscFunctionBegin; 281 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 282 PetscValidPointer(theta,2); 283 ierr = PetscObjectQueryFunction((PetscObject)ts,"TSThetaGetTheta_C",(void(**)(void))&f);CHKERRQ(ierr); 284 if (!f) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"TS type %s",((PetscObject)ts)->type_name); 285 ierr = (*f)(ts,theta);CHKERRQ(ierr); 286 PetscFunctionReturn(0); 287 } 288 289 #undef __FUNCT__ 290 #define __FUNCT__ "TSThetaSetTheta" 291 /*@ 292 TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 293 294 Not Collective 295 296 Input Parameter: 297 + ts - timestepping context 298 - theta - stage abscissa 299 300 Options Database: 301 . -ts_theta_theta <theta> 302 303 Level: Intermediate 304 305 .seealso: TSThetaGetTheta() 306 @*/ 307 PetscErrorCode PETSCTS_DLLEXPORT TSThetaSetTheta(TS ts,PetscReal theta) 308 { 309 PetscErrorCode ierr,(*f)(TS,PetscReal); 310 311 PetscFunctionBegin; 312 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 313 ierr = PetscObjectQueryFunction((PetscObject)ts,"TSThetaSetTheta_C",(void(**)(void))&f);CHKERRQ(ierr); 314 if (f) {ierr = (*f)(ts,theta);CHKERRQ(ierr);} 315 PetscFunctionReturn(0); 316 } 317