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