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