1 2 #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 3 #include <petscdm.h> 4 5 static const PetscInt TSEIMEXDefault = 3; 6 7 typedef struct { 8 PetscInt row_ind; /* Return the term T[row_ind][col_ind] */ 9 PetscInt col_ind; /* Return the term T[row_ind][col_ind] */ 10 PetscInt nstages; /* Numbers of stages in current scheme */ 11 PetscInt max_rows; /* Maximum number of rows */ 12 PetscInt *N; /* Harmonic sequence N[max_rows] */ 13 Vec Y; /* States computed during the step, used to complete the step */ 14 Vec Z; /* For shift*(Y-Z) */ 15 Vec *T; /* Working table, size determined by nstages */ 16 Vec YdotRHS; /* f(x) Work vector holding YdotRHS during residual evaluation */ 17 Vec YdotI; /* xdot-g(x) Work vector holding YdotI = G(t,x,xdot) when xdot =0 */ 18 Vec Ydot; /* f(x)+g(x) Work vector */ 19 Vec VecSolPrev; /* Work vector holding the solution from the previous step (used for interpolation) */ 20 PetscReal shift; 21 PetscReal ctime; 22 PetscBool recompute_jacobian; /* Recompute the Jacobian at each stage, default is to freeze the Jacobian at the start of each step */ 23 PetscBool ord_adapt; /* order adapativity */ 24 TSStepStatus status; 25 } TS_EIMEX; 26 27 /* This function is pure */ 28 static PetscInt Map(PetscInt i, PetscInt j, PetscInt s) 29 { 30 return ((2*s-j+1)*j/2+i-j); 31 } 32 33 static PetscErrorCode TSEvaluateStep_EIMEX(TS ts,PetscInt order,Vec X,PetscBool *done) 34 { 35 TS_EIMEX *ext = (TS_EIMEX*)ts->data; 36 const PetscInt ns = ext->nstages; 37 PetscErrorCode ierr; 38 PetscFunctionBegin; 39 ierr = VecCopy(ext->T[Map(ext->row_ind,ext->col_ind,ns)],X);CHKERRQ(ierr); 40 PetscFunctionReturn(0); 41 } 42 43 static PetscErrorCode TSStage_EIMEX(TS ts,PetscInt istage) 44 { 45 TS_EIMEX *ext = (TS_EIMEX*)ts->data; 46 PetscReal h; 47 Vec Y=ext->Y, Z=ext->Z; 48 SNES snes; 49 TSAdapt adapt; 50 PetscInt i,its,lits; 51 PetscBool accept; 52 PetscErrorCode ierr; 53 54 PetscFunctionBegin; 55 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 56 h = ts->time_step/ext->N[istage];/* step size for the istage-th stage */ 57 ext->shift = 1./h; 58 ierr = SNESSetLagJacobian(snes,-2);CHKERRQ(ierr); /* Recompute the Jacobian on this solve, but not again */ 59 ierr = VecCopy(ext->VecSolPrev,Y);CHKERRQ(ierr); /* Take the previous solution as initial step */ 60 61 for (i=0; i<ext->N[istage]; i++) { 62 ext->ctime = ts->ptime + h*i; 63 ierr = VecCopy(Y,Z);CHKERRQ(ierr);/* Save the solution of the previous substep */ 64 ierr = SNESSolve(snes,NULL,Y);CHKERRQ(ierr); 65 ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr); 66 ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr); 67 ts->snes_its += its; ts->ksp_its += lits; 68 ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 69 ierr = TSAdaptCheckStage(adapt,ts,ext->ctime,Y,&accept);CHKERRQ(ierr); 70 } 71 PetscFunctionReturn(0); 72 } 73 74 static PetscErrorCode TSStep_EIMEX(TS ts) 75 { 76 TS_EIMEX *ext = (TS_EIMEX*)ts->data; 77 const PetscInt ns = ext->nstages; 78 Vec *T=ext->T, Y=ext->Y; 79 80 SNES snes; 81 PetscInt i,j; 82 PetscBool accept = PETSC_FALSE; 83 PetscErrorCode ierr; 84 PetscReal alpha,local_error,local_error_a,local_error_r; 85 PetscFunctionBegin; 86 87 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 88 ierr = SNESSetType(snes,"ksponly");CHKERRQ(ierr); 89 ext->status = TS_STEP_INCOMPLETE; 90 91 ierr = VecCopy(ts->vec_sol,ext->VecSolPrev);CHKERRQ(ierr); 92 93 /* Apply n_j steps of the base method to obtain solutions of T(j,1),1<=j<=s */ 94 for (j=0; j<ns; j++) { 95 ierr = TSStage_EIMEX(ts,j);CHKERRQ(ierr); 96 ierr = VecCopy(Y,T[j]);CHKERRQ(ierr); 97 } 98 99 for (i=1;i<ns;i++) { 100 for (j=i;j<ns;j++) { 101 alpha = -(PetscReal)ext->N[j]/ext->N[j-i]; 102 ierr = VecAXPBYPCZ(T[Map(j,i,ns)],alpha,1.0,0,T[Map(j,i-1,ns)],T[Map(j-1,i-1,ns)]);/* T[j][i]=alpha*T[j][i-1]+T[j-1][i-1] */CHKERRQ(ierr); 103 alpha = 1.0/(1.0 + alpha); 104 ierr = VecScale(T[Map(j,i,ns)],alpha);CHKERRQ(ierr); 105 } 106 } 107 108 ierr = TSEvaluateStep(ts,ns,ts->vec_sol,NULL);CHKERRQ(ierr);/*update ts solution */ 109 110 if (ext->ord_adapt && ext->nstages < ext->max_rows) { 111 accept = PETSC_FALSE; 112 while (!accept && ext->nstages < ext->max_rows) { 113 ierr = TSErrorWeightedNorm(ts,ts->vec_sol,T[Map(ext->nstages-1,ext->nstages-2,ext->nstages)],ts->adapt->wnormtype,&local_error,&local_error_a,&local_error_r);CHKERRQ(ierr); 114 accept = (local_error < 1.0)? PETSC_TRUE : PETSC_FALSE; 115 116 if (!accept) {/* add one more stage*/ 117 ierr = TSStage_EIMEX(ts,ext->nstages);CHKERRQ(ierr); 118 ext->nstages++; ext->row_ind++; ext->col_ind++; 119 /*T table need to be recycled*/ 120 ierr = VecDuplicateVecs(ts->vec_sol,(1+ext->nstages)*ext->nstages/2,&ext->T);CHKERRQ(ierr); 121 for (i=0; i<ext->nstages-1; i++) { 122 for (j=0; j<=i; j++) { 123 ierr = VecCopy(T[Map(i,j,ext->nstages-1)],ext->T[Map(i,j,ext->nstages)]);CHKERRQ(ierr); 124 } 125 } 126 ierr = VecDestroyVecs(ext->nstages*(ext->nstages-1)/2,&T);CHKERRQ(ierr); 127 T = ext->T; /*reset the pointer*/ 128 /*recycling finished, store the new solution*/ 129 ierr = VecCopy(Y,T[ext->nstages-1]);CHKERRQ(ierr); 130 /*extrapolation for the newly added stage*/ 131 for (i=1;i<ext->nstages;i++) { 132 alpha = -(PetscReal)ext->N[ext->nstages-1]/ext->N[ext->nstages-1-i]; 133 ierr = VecAXPBYPCZ(T[Map(ext->nstages-1,i,ext->nstages)],alpha,1.0,0,T[Map(ext->nstages-1,i-1,ext->nstages)],T[Map(ext->nstages-1-1,i-1,ext->nstages)]);/*T[ext->nstages-1][i]=alpha*T[ext->nstages-1][i-1]+T[ext->nstages-1-1][i-1]*/CHKERRQ(ierr); 134 alpha = 1.0/(1.0 + alpha); 135 ierr = VecScale(T[Map(ext->nstages-1,i,ext->nstages)],alpha);CHKERRQ(ierr); 136 } 137 /*update ts solution */ 138 ierr = TSEvaluateStep(ts,ext->nstages,ts->vec_sol,NULL);CHKERRQ(ierr); 139 }/*end if !accept*/ 140 }/*end while*/ 141 142 if (ext->nstages == ext->max_rows) { 143 ierr = PetscInfo(ts,"Max number of rows has been used\n");CHKERRQ(ierr); 144 } 145 }/*end if ext->ord_adapt*/ 146 ts->ptime += ts->time_step; 147 ext->status = TS_STEP_COMPLETE; 148 149 if (ext->status != TS_STEP_COMPLETE && !ts->reason) ts->reason = TS_DIVERGED_STEP_REJECTED; 150 PetscFunctionReturn(0); 151 } 152 153 /* cubic Hermit spline */ 154 static PetscErrorCode TSInterpolate_EIMEX(TS ts,PetscReal itime,Vec X) 155 { 156 TS_EIMEX *ext = (TS_EIMEX*)ts->data; 157 PetscReal t,a,b; 158 Vec Y0=ext->VecSolPrev,Y1=ext->Y,Ydot=ext->Ydot,YdotI=ext->YdotI; 159 const PetscReal h = ts->ptime - ts->ptime_prev; 160 PetscErrorCode ierr; 161 PetscFunctionBegin; 162 t = (itime -ts->ptime + h)/h; 163 /* YdotI = -f(x)-g(x) */ 164 165 ierr = VecZeroEntries(Ydot);CHKERRQ(ierr); 166 ierr = TSComputeIFunction(ts,ts->ptime-h,Y0,Ydot,YdotI,PETSC_FALSE);CHKERRQ(ierr); 167 168 a = 2.0*t*t*t - 3.0*t*t + 1.0; 169 b = -(t*t*t - 2.0*t*t + t)*h; 170 ierr = VecAXPBYPCZ(X,a,b,0.0,Y0,YdotI);CHKERRQ(ierr); 171 172 ierr = TSComputeIFunction(ts,ts->ptime,Y1,Ydot,YdotI,PETSC_FALSE);CHKERRQ(ierr); 173 a = -2.0*t*t*t+3.0*t*t; 174 b = -(t*t*t - t*t)*h; 175 ierr = VecAXPBYPCZ(X,a,b,1.0,Y1,YdotI);CHKERRQ(ierr); 176 177 PetscFunctionReturn(0); 178 } 179 180 static PetscErrorCode TSReset_EIMEX(TS ts) 181 { 182 TS_EIMEX *ext = (TS_EIMEX*)ts->data; 183 PetscInt ns; 184 PetscErrorCode ierr; 185 186 PetscFunctionBegin; 187 ns = ext->nstages; 188 ierr = VecDestroyVecs((1+ns)*ns/2,&ext->T);CHKERRQ(ierr); 189 ierr = VecDestroy(&ext->Y);CHKERRQ(ierr); 190 ierr = VecDestroy(&ext->Z);CHKERRQ(ierr); 191 ierr = VecDestroy(&ext->YdotRHS);CHKERRQ(ierr); 192 ierr = VecDestroy(&ext->YdotI);CHKERRQ(ierr); 193 ierr = VecDestroy(&ext->Ydot);CHKERRQ(ierr); 194 ierr = VecDestroy(&ext->VecSolPrev);CHKERRQ(ierr); 195 ierr = PetscFree(ext->N);CHKERRQ(ierr); 196 PetscFunctionReturn(0); 197 } 198 199 static PetscErrorCode TSDestroy_EIMEX(TS ts) 200 { 201 PetscErrorCode ierr; 202 203 PetscFunctionBegin; 204 ierr = TSReset_EIMEX(ts);CHKERRQ(ierr); 205 ierr = PetscFree(ts->data);CHKERRQ(ierr); 206 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSEIMEXSetMaxRows_C",NULL);CHKERRQ(ierr); 207 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSEIMEXSetRowCol_C",NULL);CHKERRQ(ierr); 208 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSEIMEXSetOrdAdapt_C",NULL);CHKERRQ(ierr); 209 PetscFunctionReturn(0); 210 } 211 212 static PetscErrorCode TSEIMEXGetVecs(TS ts,DM dm,Vec *Z,Vec *Ydot,Vec *YdotI, Vec *YdotRHS) 213 { 214 TS_EIMEX *ext = (TS_EIMEX*)ts->data; 215 PetscErrorCode ierr; 216 217 PetscFunctionBegin; 218 if (Z) { 219 if (dm && dm != ts->dm) { 220 ierr = DMGetNamedGlobalVector(dm,"TSEIMEX_Z",Z);CHKERRQ(ierr); 221 } else *Z = ext->Z; 222 } 223 if (Ydot) { 224 if (dm && dm != ts->dm) { 225 ierr = DMGetNamedGlobalVector(dm,"TSEIMEX_Ydot",Ydot);CHKERRQ(ierr); 226 } else *Ydot = ext->Ydot; 227 } 228 if (YdotI) { 229 if (dm && dm != ts->dm) { 230 ierr = DMGetNamedGlobalVector(dm,"TSEIMEX_YdotI",YdotI);CHKERRQ(ierr); 231 } else *YdotI = ext->YdotI; 232 } 233 if (YdotRHS) { 234 if (dm && dm != ts->dm) { 235 ierr = DMGetNamedGlobalVector(dm,"TSEIMEX_YdotRHS",YdotRHS);CHKERRQ(ierr); 236 } else *YdotRHS = ext->YdotRHS; 237 } 238 PetscFunctionReturn(0); 239 } 240 241 static PetscErrorCode TSEIMEXRestoreVecs(TS ts,DM dm,Vec *Z,Vec *Ydot,Vec *YdotI,Vec *YdotRHS) 242 { 243 PetscErrorCode ierr; 244 245 PetscFunctionBegin; 246 if (Z) { 247 if (dm && dm != ts->dm) { 248 ierr = DMRestoreNamedGlobalVector(dm,"TSEIMEX_Z",Z);CHKERRQ(ierr); 249 } 250 } 251 if (Ydot) { 252 if (dm && dm != ts->dm) { 253 ierr = DMRestoreNamedGlobalVector(dm,"TSEIMEX_Ydot",Ydot);CHKERRQ(ierr); 254 } 255 } 256 if (YdotI) { 257 if (dm && dm != ts->dm) { 258 ierr = DMRestoreNamedGlobalVector(dm,"TSEIMEX_YdotI",YdotI);CHKERRQ(ierr); 259 } 260 } 261 if (YdotRHS) { 262 if (dm && dm != ts->dm) { 263 ierr = DMRestoreNamedGlobalVector(dm,"TSEIMEX_YdotRHS",YdotRHS);CHKERRQ(ierr); 264 } 265 } 266 PetscFunctionReturn(0); 267 } 268 269 /* 270 This defines the nonlinear equation that is to be solved with SNES 271 Fn[t0+Theta*dt, U, (U-U0)*shift] = 0 272 In the case of Backward Euler, Fn = (U-U0)/h-g(t1,U)) 273 Since FormIFunction calculates G = ydot - g(t,y), ydot will be set to (U-U0)/h 274 */ 275 static PetscErrorCode SNESTSFormFunction_EIMEX(SNES snes,Vec X,Vec G,TS ts) 276 { 277 TS_EIMEX *ext = (TS_EIMEX*)ts->data; 278 PetscErrorCode ierr; 279 Vec Ydot,Z; 280 DM dm,dmsave; 281 282 PetscFunctionBegin; 283 ierr = VecZeroEntries(G);CHKERRQ(ierr); 284 285 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 286 ierr = TSEIMEXGetVecs(ts,dm,&Z,&Ydot,NULL,NULL);CHKERRQ(ierr); 287 ierr = VecZeroEntries(Ydot);CHKERRQ(ierr); 288 dmsave = ts->dm; 289 ts->dm = dm; 290 ierr = TSComputeIFunction(ts,ext->ctime,X,Ydot,G,PETSC_FALSE);CHKERRQ(ierr); 291 /* PETSC_FALSE indicates non-imex, adding explicit RHS to the implicit I function. */ 292 ierr = VecCopy(G,Ydot);CHKERRQ(ierr); 293 ts->dm = dmsave; 294 ierr = TSEIMEXRestoreVecs(ts,dm,&Z,&Ydot,NULL,NULL);CHKERRQ(ierr); 295 296 PetscFunctionReturn(0); 297 } 298 299 /* 300 This defined the Jacobian matrix for SNES. Jn = (I/h-g'(t,y)) 301 */ 302 static PetscErrorCode SNESTSFormJacobian_EIMEX(SNES snes,Vec X,Mat A,Mat B,TS ts) 303 { 304 TS_EIMEX *ext = (TS_EIMEX*)ts->data; 305 Vec Ydot; 306 PetscErrorCode ierr; 307 DM dm,dmsave; 308 PetscFunctionBegin; 309 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 310 ierr = TSEIMEXGetVecs(ts,dm,NULL,&Ydot,NULL,NULL);CHKERRQ(ierr); 311 /* ierr = VecZeroEntries(Ydot);CHKERRQ(ierr); */ 312 /* ext->Ydot have already been computed in SNESTSFormFunction_EIMEX (SNES guarantees this) */ 313 dmsave = ts->dm; 314 ts->dm = dm; 315 ierr = TSComputeIJacobian(ts,ts->ptime,X,Ydot,ext->shift,A,B,PETSC_TRUE);CHKERRQ(ierr); 316 ts->dm = dmsave; 317 ierr = TSEIMEXRestoreVecs(ts,dm,NULL,&Ydot,NULL,NULL);CHKERRQ(ierr); 318 PetscFunctionReturn(0); 319 } 320 321 static PetscErrorCode DMCoarsenHook_TSEIMEX(DM fine,DM coarse,void *ctx) 322 { 323 324 PetscFunctionBegin; 325 PetscFunctionReturn(0); 326 } 327 328 static PetscErrorCode DMRestrictHook_TSEIMEX(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 329 { 330 TS ts = (TS)ctx; 331 PetscErrorCode ierr; 332 Vec Z,Z_c; 333 334 PetscFunctionBegin; 335 ierr = TSEIMEXGetVecs(ts,fine,&Z,NULL,NULL,NULL);CHKERRQ(ierr); 336 ierr = TSEIMEXGetVecs(ts,coarse,&Z_c,NULL,NULL,NULL);CHKERRQ(ierr); 337 ierr = MatRestrict(restrct,Z,Z_c);CHKERRQ(ierr); 338 ierr = VecPointwiseMult(Z_c,rscale,Z_c);CHKERRQ(ierr); 339 ierr = TSEIMEXRestoreVecs(ts,fine,&Z,NULL,NULL,NULL);CHKERRQ(ierr); 340 ierr = TSEIMEXRestoreVecs(ts,coarse,&Z_c,NULL,NULL,NULL);CHKERRQ(ierr); 341 PetscFunctionReturn(0); 342 } 343 344 static PetscErrorCode TSSetUp_EIMEX(TS ts) 345 { 346 TS_EIMEX *ext = (TS_EIMEX*)ts->data; 347 PetscErrorCode ierr; 348 DM dm; 349 350 PetscFunctionBegin; 351 if (!ext->N) { /* ext->max_rows not set */ 352 ierr = TSEIMEXSetMaxRows(ts,TSEIMEXDefault);CHKERRQ(ierr); 353 } 354 if (-1 == ext->row_ind && -1 == ext->col_ind) { 355 ierr = TSEIMEXSetRowCol(ts,ext->max_rows,ext->max_rows);CHKERRQ(ierr); 356 } else{/* ext->row_ind and col_ind already set */ 357 if (ext->ord_adapt) { 358 ierr = PetscInfo(ts,"Order adaptivity is enabled and TSEIMEXSetRowCol or -ts_eimex_row_col option will take no effect\n");CHKERRQ(ierr); 359 } 360 } 361 362 if (ext->ord_adapt) { 363 ext->nstages = 2; /* Start with the 2-stage scheme */ 364 ierr = TSEIMEXSetRowCol(ts,ext->nstages,ext->nstages);CHKERRQ(ierr); 365 } else{ 366 ext->nstages = ext->max_rows; /* by default nstages is the same as max_rows, this can be changed by setting order adaptivity */ 367 } 368 369 ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr); 370 371 ierr = VecDuplicateVecs(ts->vec_sol,(1+ext->nstages)*ext->nstages/2,&ext->T);CHKERRQ(ierr);/* full T table */ 372 ierr = VecDuplicate(ts->vec_sol,&ext->YdotI);CHKERRQ(ierr); 373 ierr = VecDuplicate(ts->vec_sol,&ext->YdotRHS);CHKERRQ(ierr); 374 ierr = VecDuplicate(ts->vec_sol,&ext->Ydot);CHKERRQ(ierr); 375 ierr = VecDuplicate(ts->vec_sol,&ext->VecSolPrev);CHKERRQ(ierr); 376 ierr = VecDuplicate(ts->vec_sol,&ext->Y);CHKERRQ(ierr); 377 ierr = VecDuplicate(ts->vec_sol,&ext->Z);CHKERRQ(ierr); 378 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 379 if (dm) { 380 ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSEIMEX,DMRestrictHook_TSEIMEX,ts);CHKERRQ(ierr); 381 } 382 PetscFunctionReturn(0); 383 } 384 385 static PetscErrorCode TSSetFromOptions_EIMEX(PetscOptionItems *PetscOptionsObject,TS ts) 386 { 387 TS_EIMEX *ext = (TS_EIMEX*)ts->data; 388 PetscErrorCode ierr; 389 PetscInt tindex[2]; 390 PetscInt np = 2, nrows=TSEIMEXDefault; 391 392 PetscFunctionBegin; 393 tindex[0] = TSEIMEXDefault; 394 tindex[1] = TSEIMEXDefault; 395 ierr = PetscOptionsHead(PetscOptionsObject,"EIMEX ODE solver options");CHKERRQ(ierr); 396 { 397 PetscBool flg; 398 ierr = PetscOptionsInt("-ts_eimex_max_rows","Define the maximum number of rows used","TSEIMEXSetMaxRows",nrows,&nrows,&flg);CHKERRQ(ierr); /* default value 3 */ 399 if (flg) { 400 ierr = TSEIMEXSetMaxRows(ts,nrows);CHKERRQ(ierr); 401 } 402 ierr = PetscOptionsIntArray("-ts_eimex_row_col","Return the specific term in the T table","TSEIMEXSetRowCol",tindex,&np,&flg);CHKERRQ(ierr); 403 if (flg) { 404 ierr = TSEIMEXSetRowCol(ts,tindex[0],tindex[1]);CHKERRQ(ierr); 405 } 406 ierr = PetscOptionsBool("-ts_eimex_order_adapt","Solve the problem with adaptive order","TSEIMEXSetOrdAdapt",ext->ord_adapt,&ext->ord_adapt,NULL);CHKERRQ(ierr); 407 } 408 ierr = PetscOptionsTail();CHKERRQ(ierr); 409 PetscFunctionReturn(0); 410 } 411 412 static PetscErrorCode TSView_EIMEX(TS ts,PetscViewer viewer) 413 { 414 PetscFunctionBegin; 415 PetscFunctionReturn(0); 416 } 417 418 /*@C 419 TSEIMEXSetMaxRows - Set the maximum number of rows for EIMEX schemes 420 421 Logically collective 422 423 Input Parameters: 424 + ts - timestepping context 425 - nrows - maximum number of rows 426 427 Level: intermediate 428 429 .seealso: TSEIMEXSetRowCol(), TSEIMEXSetOrdAdapt(), TSEIMEX 430 @*/ 431 PetscErrorCode TSEIMEXSetMaxRows(TS ts, PetscInt nrows) 432 { 433 PetscErrorCode ierr; 434 PetscFunctionBegin; 435 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 436 ierr = PetscTryMethod(ts,"TSEIMEXSetMaxRows_C",(TS,PetscInt),(ts,nrows));CHKERRQ(ierr); 437 PetscFunctionReturn(0); 438 } 439 440 /*@C 441 TSEIMEXSetRowCol - Set the type index in the T table for the return value 442 443 Logically collective 444 445 Input Parameters: 446 + ts - timestepping context 447 - tindex - index in the T table 448 449 Level: intermediate 450 451 .seealso: TSEIMEXSetMaxRows(), TSEIMEXSetOrdAdapt(), TSEIMEX 452 @*/ 453 PetscErrorCode TSEIMEXSetRowCol(TS ts, PetscInt row, PetscInt col) 454 { 455 PetscErrorCode ierr; 456 PetscFunctionBegin; 457 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 458 ierr = PetscTryMethod(ts,"TSEIMEXSetRowCol_C",(TS,PetscInt, PetscInt),(ts,row,col));CHKERRQ(ierr); 459 PetscFunctionReturn(0); 460 } 461 462 /*@C 463 TSEIMEXSetOrdAdapt - Set the order adaptativity 464 465 Logically collective 466 467 Input Parameters: 468 + ts - timestepping context 469 - tindex - index in the T table 470 471 Level: intermediate 472 473 .seealso: TSEIMEXSetRowCol(), TSEIMEXSetOrdAdapt(), TSEIMEX 474 @*/ 475 PetscErrorCode TSEIMEXSetOrdAdapt(TS ts, PetscBool flg) 476 { 477 PetscErrorCode ierr; 478 PetscFunctionBegin; 479 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 480 ierr = PetscTryMethod(ts,"TSEIMEXSetOrdAdapt_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 481 PetscFunctionReturn(0); 482 } 483 484 static PetscErrorCode TSEIMEXSetMaxRows_EIMEX(TS ts,PetscInt nrows) 485 { 486 TS_EIMEX *ext = (TS_EIMEX*)ts->data; 487 PetscErrorCode ierr; 488 PetscInt i; 489 490 PetscFunctionBegin; 491 if (nrows < 0 || nrows > 100) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Max number of rows (current value %D) should be an integer number between 1 and 100\n",nrows); 492 ierr = PetscFree(ext->N);CHKERRQ(ierr); 493 ext->max_rows = nrows; 494 ierr = PetscMalloc1(nrows,&ext->N);CHKERRQ(ierr); 495 for (i=0;i<nrows;i++) ext->N[i]=i+1; 496 PetscFunctionReturn(0); 497 } 498 499 static PetscErrorCode TSEIMEXSetRowCol_EIMEX(TS ts,PetscInt row,PetscInt col) 500 { 501 TS_EIMEX *ext = (TS_EIMEX*)ts->data; 502 503 PetscFunctionBegin; 504 if (row < 1 || col < 1) SETERRQ2(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"The row or column index (current value %d,%d) should not be less than 1 \n",row,col); 505 if (row > ext->max_rows || col > ext->max_rows) SETERRQ3(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"The row or column index (current value %d,%d) exceeds the maximum number of rows %d\n",row,col,ext->max_rows); 506 if (col > row) SETERRQ2(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"The column index (%d) exceeds the row index (%d)\n",col,row); 507 508 ext->row_ind = row - 1; 509 ext->col_ind = col - 1; /* Array index in C starts from 0 */ 510 PetscFunctionReturn(0); 511 } 512 513 static PetscErrorCode TSEIMEXSetOrdAdapt_EIMEX(TS ts,PetscBool flg) 514 { 515 TS_EIMEX *ext = (TS_EIMEX*)ts->data; 516 PetscFunctionBegin; 517 ext->ord_adapt = flg; 518 PetscFunctionReturn(0); 519 } 520 521 /*MC 522 TSEIMEX - Time stepping with Extrapolated IMEX methods. 523 524 These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly nonlinear such that it 525 is expensive to solve with a fully implicit method. The user should provide the stiff part of the equation using TSSetIFunction() and the 526 non-stiff part with TSSetRHSFunction(). 527 528 Notes: 529 The default is a 3-stage scheme, it can be changed with TSEIMEXSetMaxRows() or -ts_eimex_max_rows 530 531 This method currently only works with ODE, for which the stiff part G(t,X,Xdot) has the form Xdot + Ghat(t,X). 532 533 The general system is written as 534 535 G(t,X,Xdot) = F(t,X) 536 537 where G represents the stiff part and F represents the non-stiff part. The user should provide the stiff part 538 of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction(). 539 This method is designed to be linearly implicit on G and can use an approximate and lagged Jacobian. 540 541 Another common form for the system is 542 543 y'=f(x)+g(x) 544 545 The relationship between F,G and f,g is 546 547 G = y'-g(x), F = f(x) 548 549 References 550 E. Constantinescu and A. Sandu, Extrapolated implicit-explicit time stepping, SIAM Journal on Scientific 551 Computing, 31 (2010), pp. 4452-4477. 552 553 Level: beginner 554 555 .seealso: TSCreate(), TS, TSSetType(), TSEIMEXSetMaxRows(), TSEIMEXSetRowCol(), TSEIMEXSetOrdAdapt() 556 557 M*/ 558 PETSC_EXTERN PetscErrorCode TSCreate_EIMEX(TS ts) 559 { 560 TS_EIMEX *ext; 561 PetscErrorCode ierr; 562 563 PetscFunctionBegin; 564 565 ts->ops->reset = TSReset_EIMEX; 566 ts->ops->destroy = TSDestroy_EIMEX; 567 ts->ops->view = TSView_EIMEX; 568 ts->ops->setup = TSSetUp_EIMEX; 569 ts->ops->step = TSStep_EIMEX; 570 ts->ops->interpolate = TSInterpolate_EIMEX; 571 ts->ops->evaluatestep = TSEvaluateStep_EIMEX; 572 ts->ops->setfromoptions = TSSetFromOptions_EIMEX; 573 ts->ops->snesfunction = SNESTSFormFunction_EIMEX; 574 ts->ops->snesjacobian = SNESTSFormJacobian_EIMEX; 575 ts->default_adapt_type = TSADAPTNONE; 576 577 ts->usessnes = PETSC_TRUE; 578 579 ierr = PetscNewLog(ts,&ext);CHKERRQ(ierr); 580 ts->data = (void*)ext; 581 582 ext->ord_adapt = PETSC_FALSE; /* By default, no order adapativity */ 583 ext->row_ind = -1; 584 ext->col_ind = -1; 585 ext->max_rows = TSEIMEXDefault; 586 ext->nstages = TSEIMEXDefault; 587 588 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSEIMEXSetMaxRows_C", TSEIMEXSetMaxRows_EIMEX);CHKERRQ(ierr); 589 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSEIMEXSetRowCol_C", TSEIMEXSetRowCol_EIMEX);CHKERRQ(ierr); 590 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSEIMEXSetOrdAdapt_C",TSEIMEXSetOrdAdapt_EIMEX);CHKERRQ(ierr); 591 PetscFunctionReturn(0); 592 } 593