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 201 /* cubic Hermit spline */ 202 #undef __FUNCT__ 203 #define __FUNCT__ "TSInterpolate_EIMEX" 204 static PetscErrorCode TSInterpolate_EIMEX(TS ts,PetscReal itime,Vec X) 205 { 206 TS_EIMEX *ext = (TS_EIMEX*)ts->data; 207 PetscReal t,a,b; 208 Vec Y0=ext->VecSolPrev,Y1=ext->Y, 209 Ydot=ext->Ydot,YdotI=ext->YdotI; 210 const PetscReal h = ts->time_step_prev; 211 PetscErrorCode ierr; 212 PetscFunctionBegin; 213 t = (itime -ts->ptime + h)/h; 214 /* YdotI = -f(x)-g(x) */ 215 216 ierr = VecZeroEntries(Ydot);CHKERRQ(ierr); 217 ierr = TSComputeIFunction(ts,ts->ptime-h,Y0,Ydot,YdotI,PETSC_FALSE);CHKERRQ(ierr); 218 219 a = 2.0*t*t*t - 3.0*t*t + 1.0; 220 b = -(t*t*t - 2.0*t*t + t)*h; 221 ierr = VecAXPBYPCZ(X,a,b,0.0,Y0,YdotI);CHKERRQ(ierr); 222 223 ierr = TSComputeIFunction(ts,ts->ptime,Y1,Ydot,YdotI,PETSC_FALSE);CHKERRQ(ierr); 224 a = -2.0*t*t*t+3.0*t*t; 225 b = -(t*t*t - t*t)*h; 226 ierr = VecAXPBYPCZ(X,a,b,1.0,Y1,YdotI);CHKERRQ(ierr); 227 228 PetscFunctionReturn(0); 229 } 230 231 232 #undef __FUNCT__ 233 #define __FUNCT__ "TSReset_EIMEX" 234 static PetscErrorCode TSReset_EIMEX(TS ts) 235 { 236 TS_EIMEX *ext = (TS_EIMEX*)ts->data; 237 PetscInt ns; 238 PetscErrorCode ierr; 239 240 PetscFunctionBegin; 241 ns = ext->nstages; 242 ierr = VecDestroyVecs((1+ns)*ns/2,&ext->T);CHKERRQ(ierr); 243 ierr = VecDestroy(&ext->Y);CHKERRQ(ierr); 244 ierr = VecDestroy(&ext->Z);CHKERRQ(ierr); 245 ierr = VecDestroy(&ext->YdotRHS);CHKERRQ(ierr); 246 ierr = VecDestroy(&ext->YdotI);CHKERRQ(ierr); 247 ierr = VecDestroy(&ext->Ydot);CHKERRQ(ierr); 248 ierr = VecDestroy(&ext->VecSolPrev);CHKERRQ(ierr); 249 ierr = PetscFree(ext->N);CHKERRQ(ierr); 250 PetscFunctionReturn(0); 251 } 252 253 #undef __FUNCT__ 254 #define __FUNCT__ "TSDestroy_EIMEX" 255 static PetscErrorCode TSDestroy_EIMEX(TS ts) 256 { 257 PetscErrorCode ierr; 258 259 PetscFunctionBegin; 260 ierr = TSReset_EIMEX(ts);CHKERRQ(ierr); 261 ierr = PetscFree(ts->data);CHKERRQ(ierr); 262 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSEIMEXSetMaxRows_C",NULL);CHKERRQ(ierr); 263 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSEIMEXSetRowCol_C",NULL);CHKERRQ(ierr); 264 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSEIMEXSetOrdAdapt_C",NULL);CHKERRQ(ierr); 265 266 PetscFunctionReturn(0); 267 } 268 269 270 #undef __FUNCT__ 271 #define __FUNCT__ "TSEIMEXGetVecs" 272 static PetscErrorCode TSEIMEXGetVecs(TS ts,DM dm,Vec *Z,Vec *Ydot,Vec *YdotI, Vec *YdotRHS) 273 { 274 TS_EIMEX *ext = (TS_EIMEX*)ts->data; 275 PetscErrorCode ierr; 276 277 PetscFunctionBegin; 278 if (Z) { 279 if (dm && dm != ts->dm) { 280 ierr = DMGetNamedGlobalVector(dm,"TSEIMEX_Z",Z);CHKERRQ(ierr); 281 } else *Z = ext->Z; 282 } 283 if (Ydot) { 284 if (dm && dm != ts->dm) { 285 ierr = DMGetNamedGlobalVector(dm,"TSEIMEX_Ydot",Ydot);CHKERRQ(ierr); 286 } else *Ydot = ext->Ydot; 287 } 288 if (YdotI) { 289 if (dm && dm != ts->dm) { 290 ierr = DMGetNamedGlobalVector(dm,"TSEIMEX_YdotI",YdotI);CHKERRQ(ierr); 291 } else *YdotI = ext->YdotI; 292 } 293 if (YdotRHS) { 294 if (dm && dm != ts->dm) { 295 ierr = DMGetNamedGlobalVector(dm,"TSEIMEX_YdotRHS",YdotRHS);CHKERRQ(ierr); 296 } else *YdotRHS = ext->YdotRHS; 297 } 298 PetscFunctionReturn(0); 299 } 300 301 302 #undef __FUNCT__ 303 #define __FUNCT__ "TSEIMEXRestoreVecs" 304 static PetscErrorCode TSEIMEXRestoreVecs(TS ts,DM dm,Vec *Z,Vec *Ydot,Vec *YdotI,Vec *YdotRHS) 305 { 306 PetscErrorCode ierr; 307 308 PetscFunctionBegin; 309 if (Z) { 310 if (dm && dm != ts->dm) { 311 ierr = DMRestoreNamedGlobalVector(dm,"TSEIMEX_Z",Z);CHKERRQ(ierr); 312 } 313 } 314 if (Ydot) { 315 if (dm && dm != ts->dm) { 316 ierr = DMRestoreNamedGlobalVector(dm,"TSEIMEX_Ydot",Ydot);CHKERRQ(ierr); 317 } 318 } 319 if (YdotI) { 320 if (dm && dm != ts->dm) { 321 ierr = DMRestoreNamedGlobalVector(dm,"TSEIMEX_YdotI",YdotI);CHKERRQ(ierr); 322 } 323 } 324 if (YdotRHS) { 325 if (dm && dm != ts->dm) { 326 ierr = DMRestoreNamedGlobalVector(dm,"TSEIMEX_YdotRHS",YdotRHS);CHKERRQ(ierr); 327 } 328 } 329 PetscFunctionReturn(0); 330 } 331 332 333 /* 334 This defines the nonlinear equation that is to be solved with SNES 335 Fn[t0+Theta*dt, U, (U-U0)*shift] = 0 336 In the case of Backward Euler, Fn = (U-U0)/h-g(t1,U)) 337 Since FormIFunction calculates G = ydot - g(t,y), ydot will be set to (U-U0)/h 338 */ 339 #undef __FUNCT__ 340 #define __FUNCT__ "SNESTSFormFunction_EIMEX" 341 static PetscErrorCode SNESTSFormFunction_EIMEX(SNES snes,Vec X,Vec G,TS ts) 342 { 343 TS_EIMEX *ext = (TS_EIMEX*)ts->data; 344 PetscErrorCode ierr; 345 Vec Ydot,Z; 346 DM dm,dmsave; 347 348 PetscFunctionBegin; 349 ierr = VecZeroEntries(G);CHKERRQ(ierr); 350 351 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 352 ierr = TSEIMEXGetVecs(ts,dm,&Z,&Ydot,NULL,NULL);CHKERRQ(ierr); 353 ierr = VecZeroEntries(Ydot);CHKERRQ(ierr); 354 dmsave = ts->dm; 355 ts->dm = dm; 356 ierr = TSComputeIFunction(ts,ext->ctime,X,Ydot,G,PETSC_FALSE);CHKERRQ(ierr); 357 /* PETSC_FALSE indicates non-imex, adding explicit RHS to the implicit I function. */ 358 ierr = VecCopy(G,Ydot);CHKERRQ(ierr); 359 ts->dm = dmsave; 360 ierr = TSEIMEXRestoreVecs(ts,dm,&Z,&Ydot,NULL,NULL);CHKERRQ(ierr); 361 362 PetscFunctionReturn(0); 363 } 364 365 /* 366 This defined the Jacobian matrix for SNES. Jn = (I/h-g'(t,y)) 367 */ 368 #undef __FUNCT__ 369 #define __FUNCT__ "SNESTSFormJacobian_EIMEX" 370 static PetscErrorCode SNESTSFormJacobian_EIMEX(SNES snes,Vec X,Mat A,Mat B,TS ts) 371 { 372 TS_EIMEX *ext = (TS_EIMEX*)ts->data; 373 Vec Ydot; 374 PetscErrorCode ierr; 375 DM dm,dmsave; 376 PetscFunctionBegin; 377 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 378 ierr = TSEIMEXGetVecs(ts,dm,NULL,&Ydot,NULL,NULL);CHKERRQ(ierr); 379 /* ierr = VecZeroEntries(Ydot);CHKERRQ(ierr); */ 380 /* ext->Ydot have already been computed in SNESTSFormFunction_EIMEX (SNES guarantees this) */ 381 dmsave = ts->dm; 382 ts->dm = dm; 383 ierr = TSComputeIJacobian(ts,ts->ptime,X,Ydot,ext->shift,A,B,PETSC_TRUE);CHKERRQ(ierr); 384 ts->dm = dmsave; 385 ierr = TSEIMEXRestoreVecs(ts,dm,NULL,&Ydot,NULL,NULL);CHKERRQ(ierr); 386 PetscFunctionReturn(0); 387 } 388 389 #undef __FUNCT__ 390 #define __FUNCT__ "DMCoarsenHook_TSEIMEX" 391 static PetscErrorCode DMCoarsenHook_TSEIMEX(DM fine,DM coarse,void *ctx) 392 { 393 394 PetscFunctionBegin; 395 PetscFunctionReturn(0); 396 } 397 398 #undef __FUNCT__ 399 #define __FUNCT__ "DMRestrictHook_TSEIMEX" 400 static PetscErrorCode DMRestrictHook_TSEIMEX(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 401 { 402 TS ts = (TS)ctx; 403 PetscErrorCode ierr; 404 Vec Z,Z_c; 405 406 PetscFunctionBegin; 407 ierr = TSEIMEXGetVecs(ts,fine,&Z,NULL,NULL,NULL);CHKERRQ(ierr); 408 ierr = TSEIMEXGetVecs(ts,coarse,&Z_c,NULL,NULL,NULL);CHKERRQ(ierr); 409 ierr = MatRestrict(restrct,Z,Z_c);CHKERRQ(ierr); 410 ierr = VecPointwiseMult(Z_c,rscale,Z_c);CHKERRQ(ierr); 411 ierr = TSEIMEXRestoreVecs(ts,fine,&Z,NULL,NULL,NULL);CHKERRQ(ierr); 412 ierr = TSEIMEXRestoreVecs(ts,coarse,&Z_c,NULL,NULL,NULL);CHKERRQ(ierr); 413 PetscFunctionReturn(0); 414 } 415 416 417 #undef __FUNCT__ 418 #define __FUNCT__ "TSSetUp_EIMEX" 419 static PetscErrorCode TSSetUp_EIMEX(TS ts) 420 { 421 TS_EIMEX *ext = (TS_EIMEX*)ts->data; 422 PetscErrorCode ierr; 423 DM dm; 424 PetscFunctionBegin; 425 426 if (!ext->N){ /* ext->max_rows not set */ 427 ierr = TSEIMEXSetMaxRows(ts,TSEIMEXDefault);CHKERRQ(ierr); 428 } 429 if(-1 == ext->row_ind && -1 == ext->col_ind){ 430 ierr = TSEIMEXSetRowCol(ts,ext->max_rows,ext->max_rows);CHKERRQ(ierr); 431 }else{/* ext->row_ind and col_ind already set */ 432 if(ext->ord_adapt){ 433 ierr = PetscInfo(ts,"Order adaptivity is enabled and TSEIMEXSetRowCol or -ts_eimex_row_col option will take no effect\n");CHKERRQ(ierr); 434 } 435 } 436 437 if(ext->ord_adapt){ 438 ext->nstages = 2; /* Start with the 2-stage scheme */ 439 ierr = TSEIMEXSetRowCol(ts,ext->nstages,ext->nstages);CHKERRQ(ierr); 440 }else{ 441 ext->nstages = ext->max_rows; /* by default nstages is the same as max_rows, this can be changed by setting order adaptivity */ 442 } 443 444 ierr = VecDuplicateVecs(ts->vec_sol,(1+ext->nstages)*ext->nstages/2,&ext->T);CHKERRQ(ierr);/* full T table */ 445 ierr = VecDuplicate(ts->vec_sol,&ext->YdotI);CHKERRQ(ierr); 446 ierr = VecDuplicate(ts->vec_sol,&ext->YdotRHS);CHKERRQ(ierr); 447 ierr = VecDuplicate(ts->vec_sol,&ext->Ydot);CHKERRQ(ierr); 448 ierr = VecDuplicate(ts->vec_sol,&ext->VecSolPrev);CHKERRQ(ierr); 449 ierr = VecDuplicate(ts->vec_sol,&ext->Y);CHKERRQ(ierr); 450 ierr = VecDuplicate(ts->vec_sol,&ext->Z);CHKERRQ(ierr); 451 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 452 if (dm) { 453 ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSEIMEX,DMRestrictHook_TSEIMEX,ts);CHKERRQ(ierr); 454 } 455 PetscFunctionReturn(0); 456 } 457 458 #undef __FUNCT__ 459 #define __FUNCT__ "TSSetFromOptions_EIMEX" 460 static PetscErrorCode TSSetFromOptions_EIMEX(TS ts) 461 { 462 TS_EIMEX *ext = (TS_EIMEX*)ts->data; 463 PetscErrorCode ierr; 464 PetscInt tindex[2]; 465 PetscInt np = 2, nrows=TSEIMEXDefault; 466 467 PetscFunctionBegin; 468 tindex[0] = TSEIMEXDefault; 469 tindex[1] = TSEIMEXDefault; 470 ierr = PetscOptionsHead("EIMEX ODE solver options");CHKERRQ(ierr); 471 { 472 PetscBool flg; 473 flg = PETSC_FALSE; 474 ierr = PetscOptionsInt("-ts_eimex_max_rows","Define the maximum number of rows used","TSEIMEXSetMaxRows",nrows,&nrows,&flg);CHKERRQ(ierr); /* default value 3 */ 475 if(flg){ 476 ierr = TSEIMEXSetMaxRows(ts,nrows);CHKERRQ(ierr); 477 } 478 ierr = PetscOptionsIntArray("-ts_eimex_row_col","Return the specific term in the T table","TSEIMEXSetRowCol",tindex,&np,&flg);CHKERRQ(ierr); 479 if(flg){ 480 ierr = TSEIMEXSetRowCol(ts,tindex[0],tindex[1]);CHKERRQ(ierr); 481 } 482 ierr = PetscOptionsBool("-ts_eimex_order_adapt","Solve the problem with adaptive order","TSEIMEXSetOrdAdapt",ext->ord_adapt,&ext->ord_adapt,NULL);CHKERRQ(ierr); 483 ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); 484 } 485 ierr = PetscOptionsTail();CHKERRQ(ierr); 486 PetscFunctionReturn(0); 487 } 488 489 #undef __FUNCT__ 490 #define __FUNCT__ "TSView_EIMEX" 491 static PetscErrorCode TSView_EIMEX(TS ts,PetscViewer viewer) 492 { 493 /* TS_EIMEX *ext = (TS_EIMEX*)ts->data; */ 494 PetscBool iascii; 495 PetscErrorCode ierr; 496 497 PetscFunctionBegin; 498 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 499 if (iascii) { 500 ierr = PetscViewerASCIIPrintf(viewer," EIMEX\n");CHKERRQ(ierr); 501 } 502 ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 503 PetscFunctionReturn(0); 504 } 505 506 507 #undef __FUNCT__ 508 #define __FUNCT__ "TSEIMEXSetMaxRows" 509 /*@C 510 TSEIMEXSetMaxRows - Set the maximum number of rows for EIMEX schemes 511 512 Logically collective 513 514 Input Parameter: 515 + ts - timestepping context 516 - nrows - maximum number of rows 517 518 Level: intermediate 519 520 .seealso: TSEIMEXSetRowCol(), TSEIMEXSetOrdAdapt(), TSEIMEX 521 @*/ 522 PetscErrorCode TSEIMEXSetMaxRows(TS ts, PetscInt nrows) 523 { 524 PetscErrorCode ierr; 525 PetscFunctionBegin; 526 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 527 ierr = PetscTryMethod(ts,"TSEIMEXSetMaxRows_C",(TS,PetscInt),(ts,nrows));CHKERRQ(ierr); 528 PetscFunctionReturn(0); 529 } 530 531 532 #undef __FUNCT__ 533 #define __FUNCT__ "TSEIMEXSetRowCol" 534 /*@C 535 TSEIMEXSetRowCol - Set the type index in the T table for the return value 536 537 Logically collective 538 539 Input Parameter: 540 + ts - timestepping context 541 - tindex - index in the T table 542 543 Level: intermediate 544 545 .seealso: TSEIMEXSetMaxRows(), TSEIMEXSetOrdAdapt(), TSEIMEX 546 @*/ 547 PetscErrorCode TSEIMEXSetRowCol(TS ts, PetscInt row, PetscInt col) 548 { 549 PetscErrorCode ierr; 550 PetscFunctionBegin; 551 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 552 ierr = PetscTryMethod(ts,"TSEIMEXSetRowCol_C",(TS,PetscInt, PetscInt),(ts,row,col));CHKERRQ(ierr); 553 PetscFunctionReturn(0); 554 } 555 556 557 #undef __FUNCT__ 558 #define __FUNCT__ "TSEIMEXSetOrdAdapt" 559 /*@C 560 TSEIMEXSetOrdAdapt - Set the order adaptativity 561 562 Logically collective 563 564 Input Parameter: 565 + ts - timestepping context 566 - tindex - index in the T table 567 568 Level: intermediate 569 570 .seealso: TSEIMEXSetRowCol(), TSEIMEXSetOrdAdapt(), TSEIMEX 571 @*/ 572 PetscErrorCode TSEIMEXSetOrdAdapt(TS ts, PetscBool flg) 573 { 574 PetscErrorCode ierr; 575 PetscFunctionBegin; 576 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 577 ierr = PetscTryMethod(ts,"TSEIMEXSetOrdAdapt_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 578 PetscFunctionReturn(0); 579 } 580 581 582 #undef __FUNCT__ 583 #define __FUNCT__ "TSEIMEXSetMaxRows_EIMEX" 584 static PetscErrorCode TSEIMEXSetMaxRows_EIMEX(TS ts,PetscInt nrows) 585 { 586 TS_EIMEX *ext = (TS_EIMEX*)ts->data; 587 PetscErrorCode ierr; 588 PetscInt i; 589 590 PetscFunctionBegin; 591 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); 592 ierr = PetscFree(ext->N);CHKERRQ(ierr); 593 ext->max_rows = nrows; 594 ierr = PetscMalloc1(nrows,&ext->N);CHKERRQ(ierr); 595 for(i=0;i<nrows;i++) ext->N[i]=i+1; 596 PetscFunctionReturn(0); 597 } 598 599 #undef __FUNCT__ 600 #define __FUNCT__ "TSEIMEXSetRowCol_EIMEX" 601 static PetscErrorCode TSEIMEXSetRowCol_EIMEX(TS ts,PetscInt row,PetscInt col) 602 { 603 TS_EIMEX *ext = (TS_EIMEX*)ts->data; 604 605 PetscFunctionBegin; 606 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); 607 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); 608 if (col > row) SETERRQ2(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"The column index (%d) exceeds the row index (%d)\n",col,row); 609 610 ext->row_ind = row - 1; 611 ext->col_ind = col - 1; /* Array index in C starts from 0 */ 612 PetscFunctionReturn(0); 613 } 614 615 #undef __FUNCT__ 616 #define __FUNCT__ "TSEIMEXSetOrdAdapt_EIMEX" 617 static PetscErrorCode TSEIMEXSetOrdAdapt_EIMEX(TS ts,PetscBool flg) 618 { 619 TS_EIMEX *ext = (TS_EIMEX*)ts->data; 620 PetscFunctionBegin; 621 ext->ord_adapt = flg; 622 PetscFunctionReturn(0); 623 } 624 625 /* ------------------------------------------------------------ */ 626 /*MC 627 TSEIMEX - ODE solver using extrapolated IMEX schemes 628 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(). 629 630 Notes: 631 The default is a 3-stage scheme, it can be changed with TSEIMEXSetMaxRows() or -ts_eimex_max_rows 632 633 This method currently only works with ODE, for which the stiff part G(t,X,Xdot) has the form Xdot + Ghat(t,X). 634 635 Level: beginner 636 637 .seealso: TSCreate(), TS 638 M*/ 639 #undef __FUNCT__ 640 #define __FUNCT__ "TSCreate_EIMEX" 641 PETSC_EXTERN PetscErrorCode TSCreate_EIMEX(TS ts) 642 { 643 TS_EIMEX *ext; 644 PetscErrorCode ierr; 645 646 PetscFunctionBegin; 647 648 ts->ops->reset = TSReset_EIMEX; 649 ts->ops->destroy = TSDestroy_EIMEX; 650 ts->ops->view = TSView_EIMEX; 651 ts->ops->setup = TSSetUp_EIMEX; 652 ts->ops->step = TSStep_EIMEX; 653 ts->ops->interpolate = TSInterpolate_EIMEX; 654 ts->ops->evaluatestep = TSEvaluateStep_EIMEX; 655 ts->ops->setfromoptions = TSSetFromOptions_EIMEX; 656 ts->ops->snesfunction = SNESTSFormFunction_EIMEX; 657 ts->ops->snesjacobian = SNESTSFormJacobian_EIMEX; 658 659 ierr = PetscNewLog(ts,&ext);CHKERRQ(ierr); 660 ts->data = (void*)ext; 661 662 ext->ord_adapt = PETSC_FALSE; /* By default, no order adapativity */ 663 ext->row_ind = -1; 664 ext->col_ind = -1; 665 ext->max_rows = TSEIMEXDefault; 666 ext->nstages = TSEIMEXDefault; 667 668 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSEIMEXSetMaxRows_C", TSEIMEXSetMaxRows_EIMEX);CHKERRQ(ierr); 669 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSEIMEXSetRowCol_C", TSEIMEXSetRowCol_EIMEX);CHKERRQ(ierr); 670 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSEIMEXSetOrdAdapt_C",TSEIMEXSetOrdAdapt_EIMEX);CHKERRQ(ierr); 671 PetscFunctionReturn(0); 672 } 673