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