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