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