1 2 #include <private/tsimpl.h> /*I "petscts.h" I*/ 3 4 /* Logging support */ 5 PetscClassId TS_CLASSID; 6 PetscLogEvent TS_Step, TS_PseudoComputeTimeStep, TS_FunctionEval, TS_JacobianEval; 7 8 #undef __FUNCT__ 9 #define __FUNCT__ "TSSetTypeFromOptions" 10 /* 11 TSSetTypeFromOptions - Sets the type of ts from user options. 12 13 Collective on TS 14 15 Input Parameter: 16 . ts - The ts 17 18 Level: intermediate 19 20 .keywords: TS, set, options, database, type 21 .seealso: TSSetFromOptions(), TSSetType() 22 */ 23 static PetscErrorCode TSSetTypeFromOptions(TS ts) 24 { 25 PetscBool opt; 26 const char *defaultType; 27 char typeName[256]; 28 PetscErrorCode ierr; 29 30 PetscFunctionBegin; 31 if (((PetscObject)ts)->type_name) { 32 defaultType = ((PetscObject)ts)->type_name; 33 } else { 34 defaultType = TSEULER; 35 } 36 37 if (!TSRegisterAllCalled) {ierr = TSRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 38 ierr = PetscOptionsList("-ts_type", "TS method"," TSSetType", TSList, defaultType, typeName, 256, &opt);CHKERRQ(ierr); 39 if (opt) { 40 ierr = TSSetType(ts, typeName);CHKERRQ(ierr); 41 } else { 42 ierr = TSSetType(ts, defaultType);CHKERRQ(ierr); 43 } 44 PetscFunctionReturn(0); 45 } 46 47 #undef __FUNCT__ 48 #define __FUNCT__ "TSSetFromOptions" 49 /*@ 50 TSSetFromOptions - Sets various TS parameters from user options. 51 52 Collective on TS 53 54 Input Parameter: 55 . ts - the TS context obtained from TSCreate() 56 57 Options Database Keys: 58 + -ts_type <type> - TSEULER, TSBEULER, TSSUNDIALS, TSPSEUDO, TSCN, TSRK, TSTHETA, TSGL, TSSSP 59 . -ts_max_steps maxsteps - maximum number of time-steps to take 60 . -ts_max_time time - maximum time to compute to 61 . -ts_dt dt - initial time step 62 . -ts_monitor - print information at each timestep 63 - -ts_monitor_draw - plot information at each timestep 64 65 Level: beginner 66 67 .keywords: TS, timestep, set, options, database 68 69 .seealso: TSGetType() 70 @*/ 71 PetscErrorCode TSSetFromOptions(TS ts) 72 { 73 PetscReal dt; 74 PetscBool opt,flg; 75 PetscErrorCode ierr; 76 PetscViewer monviewer; 77 char monfilename[PETSC_MAX_PATH_LEN]; 78 SNES snes; 79 80 PetscFunctionBegin; 81 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 82 ierr = PetscOptionsBegin(((PetscObject)ts)->comm, ((PetscObject)ts)->prefix, "Time step options", "TS");CHKERRQ(ierr); 83 84 /* Handle generic TS options */ 85 ierr = PetscOptionsInt("-ts_max_steps","Maximum number of time steps","TSSetDuration",ts->max_steps,&ts->max_steps,PETSC_NULL);CHKERRQ(ierr); 86 ierr = PetscOptionsReal("-ts_max_time","Time to run to","TSSetDuration",ts->max_time,&ts->max_time,PETSC_NULL);CHKERRQ(ierr); 87 ierr = PetscOptionsReal("-ts_init_time","Initial time","TSSetInitialTime", ts->ptime, &ts->ptime, PETSC_NULL);CHKERRQ(ierr); 88 ierr = PetscOptionsReal("-ts_dt","Initial time step","TSSetInitialTimeStep",ts->initial_time_step,&dt,&opt);CHKERRQ(ierr); 89 if (opt) {ierr = TSSetInitialTimeStep(ts,ts->ptime,dt);CHKERRQ(ierr);} 90 ierr = PetscOptionsInt("-ts_max_snes_failures","Maximum number of nonlinear solve failures","",ts->max_snes_failures,&ts->max_snes_failures,PETSC_NULL);CHKERRQ(ierr); 91 ierr = PetscOptionsInt("-ts_max_reject","Maximum number of step rejections","",ts->max_reject,&ts->max_reject,PETSC_NULL);CHKERRQ(ierr); 92 ierr = PetscOptionsBool("-ts_error_if_step_failed","Error if no step succeeds","",ts->errorifstepfailed,&ts->errorifstepfailed,PETSC_NULL);CHKERRQ(ierr); 93 94 /* Monitor options */ 95 ierr = PetscOptionsString("-ts_monitor","Monitor timestep size","TSMonitorDefault","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 96 if (flg) { 97 ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,monfilename,&monviewer);CHKERRQ(ierr); 98 ierr = TSMonitorSet(ts,TSMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 99 } 100 opt = PETSC_FALSE; 101 ierr = PetscOptionsBool("-ts_monitor_draw","Monitor timestep size graphically","TSMonitorLG",opt,&opt,PETSC_NULL);CHKERRQ(ierr); 102 if (opt) { 103 ierr = TSMonitorSet(ts,TSMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 104 } 105 opt = PETSC_FALSE; 106 ierr = PetscOptionsBool("-ts_monitor_solution","Monitor solution graphically","TSMonitorSolution",opt,&opt,PETSC_NULL);CHKERRQ(ierr); 107 if (opt) { 108 ierr = TSMonitorSet(ts,TSMonitorSolution,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 109 } 110 111 /* Handle TS type options */ 112 ierr = TSSetTypeFromOptions(ts);CHKERRQ(ierr); 113 114 /* Handle specific TS options */ 115 if (ts->ops->setfromoptions) { 116 ierr = (*ts->ops->setfromoptions)(ts);CHKERRQ(ierr); 117 } 118 119 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 120 ierr = PetscObjectProcessOptionsHandlers((PetscObject)ts);CHKERRQ(ierr); 121 ierr = PetscOptionsEnd();CHKERRQ(ierr); 122 123 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 124 /* Handle subobject options */ 125 if (ts->problem_type == TS_LINEAR) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);} 126 ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 127 PetscFunctionReturn(0); 128 } 129 130 #undef __FUNCT__ 131 #define __FUNCT__ "TSViewFromOptions" 132 /*@ 133 TSViewFromOptions - This function visualizes the ts based upon user options. 134 135 Collective on TS 136 137 Input Parameter: 138 . ts - The ts 139 140 Level: intermediate 141 142 .keywords: TS, view, options, database 143 .seealso: TSSetFromOptions(), TSView() 144 @*/ 145 PetscErrorCode TSViewFromOptions(TS ts,const char title[]) 146 { 147 PetscViewer viewer; 148 PetscDraw draw; 149 PetscBool opt = PETSC_FALSE; 150 char fileName[PETSC_MAX_PATH_LEN]; 151 PetscErrorCode ierr; 152 153 PetscFunctionBegin; 154 ierr = PetscOptionsGetString(((PetscObject)ts)->prefix, "-ts_view", fileName, PETSC_MAX_PATH_LEN, &opt);CHKERRQ(ierr); 155 if (opt && !PetscPreLoadingOn) { 156 ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,fileName,&viewer);CHKERRQ(ierr); 157 ierr = TSView(ts, viewer);CHKERRQ(ierr); 158 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 159 } 160 opt = PETSC_FALSE; 161 ierr = PetscOptionsGetBool(((PetscObject)ts)->prefix, "-ts_view_draw", &opt,PETSC_NULL);CHKERRQ(ierr); 162 if (opt) { 163 ierr = PetscViewerDrawOpen(((PetscObject)ts)->comm, 0, 0, 0, 0, 300, 300, &viewer);CHKERRQ(ierr); 164 ierr = PetscViewerDrawGetDraw(viewer, 0, &draw);CHKERRQ(ierr); 165 if (title) { 166 ierr = PetscDrawSetTitle(draw, (char *)title);CHKERRQ(ierr); 167 } else { 168 ierr = PetscObjectName((PetscObject)ts);CHKERRQ(ierr); 169 ierr = PetscDrawSetTitle(draw, ((PetscObject)ts)->name);CHKERRQ(ierr); 170 } 171 ierr = TSView(ts, viewer);CHKERRQ(ierr); 172 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 173 ierr = PetscDrawPause(draw);CHKERRQ(ierr); 174 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 175 } 176 PetscFunctionReturn(0); 177 } 178 179 #undef __FUNCT__ 180 #define __FUNCT__ "TSComputeRHSJacobian" 181 /*@ 182 TSComputeRHSJacobian - Computes the Jacobian matrix that has been 183 set with TSSetRHSJacobian(). 184 185 Collective on TS and Vec 186 187 Input Parameters: 188 + ts - the TS context 189 . t - current timestep 190 - x - input vector 191 192 Output Parameters: 193 + A - Jacobian matrix 194 . B - optional preconditioning matrix 195 - flag - flag indicating matrix structure 196 197 Notes: 198 Most users should not need to explicitly call this routine, as it 199 is used internally within the nonlinear solvers. 200 201 See KSPSetOperators() for important information about setting the 202 flag parameter. 203 204 Level: developer 205 206 .keywords: SNES, compute, Jacobian, matrix 207 208 .seealso: TSSetRHSJacobian(), KSPSetOperators() 209 @*/ 210 PetscErrorCode TSComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg) 211 { 212 PetscErrorCode ierr; 213 PetscInt Xstate; 214 215 PetscFunctionBegin; 216 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 217 PetscValidHeaderSpecific(X,VEC_CLASSID,3); 218 PetscCheckSameComm(ts,1,X,3); 219 ierr = PetscObjectStateQuery((PetscObject)X,&Xstate);CHKERRQ(ierr); 220 if (ts->rhsjacobian.time == t && (ts->problem_type == TS_LINEAR || (ts->rhsjacobian.X == X && ts->rhsjacobian.Xstate == Xstate))) { 221 *flg = ts->rhsjacobian.mstructure; 222 PetscFunctionReturn(0); 223 } 224 225 if (!ts->userops->rhsjacobian && !ts->userops->ijacobian) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSJacobian() and / or TSSetIJacobian()"); 226 227 if (ts->userops->rhsjacobian) { 228 ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 229 *flg = DIFFERENT_NONZERO_PATTERN; 230 PetscStackPush("TS user Jacobian function"); 231 ierr = (*ts->userops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr); 232 PetscStackPop; 233 ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 234 /* make sure user returned a correct Jacobian and preconditioner */ 235 PetscValidHeaderSpecific(*A,MAT_CLASSID,4); 236 PetscValidHeaderSpecific(*B,MAT_CLASSID,5); 237 } else { 238 ierr = MatZeroEntries(*A);CHKERRQ(ierr); 239 if (*A != *B) {ierr = MatZeroEntries(*B);CHKERRQ(ierr);} 240 *flg = SAME_NONZERO_PATTERN; 241 } 242 ts->rhsjacobian.time = t; 243 ts->rhsjacobian.X = X; 244 ierr = PetscObjectStateQuery((PetscObject)X,&ts->rhsjacobian.Xstate);CHKERRQ(ierr); 245 ts->rhsjacobian.mstructure = *flg; 246 PetscFunctionReturn(0); 247 } 248 249 #undef __FUNCT__ 250 #define __FUNCT__ "TSComputeRHSFunction" 251 /*@ 252 TSComputeRHSFunction - Evaluates the right-hand-side function. 253 254 Collective on TS and Vec 255 256 Input Parameters: 257 + ts - the TS context 258 . t - current time 259 - x - state vector 260 261 Output Parameter: 262 . y - right hand side 263 264 Note: 265 Most users should not need to explicitly call this routine, as it 266 is used internally within the nonlinear solvers. 267 268 Level: developer 269 270 .keywords: TS, compute 271 272 .seealso: TSSetRHSFunction(), TSComputeIFunction() 273 @*/ 274 PetscErrorCode TSComputeRHSFunction(TS ts,PetscReal t,Vec x,Vec y) 275 { 276 PetscErrorCode ierr; 277 278 PetscFunctionBegin; 279 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 280 PetscValidHeaderSpecific(x,VEC_CLASSID,3); 281 PetscValidHeaderSpecific(y,VEC_CLASSID,4); 282 283 if (!ts->userops->rhsfunction && !ts->userops->ifunction) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSFunction() and / or TSSetIFunction()"); 284 285 ierr = PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr); 286 if (ts->userops->rhsfunction) { 287 PetscStackPush("TS user right-hand-side function"); 288 ierr = (*ts->userops->rhsfunction)(ts,t,x,y,ts->funP);CHKERRQ(ierr); 289 PetscStackPop; 290 } else { 291 ierr = VecZeroEntries(y);CHKERRQ(ierr); 292 } 293 294 ierr = PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr); 295 PetscFunctionReturn(0); 296 } 297 298 #undef __FUNCT__ 299 #define __FUNCT__ "TSGetRHSVec_Private" 300 static PetscErrorCode TSGetRHSVec_Private(TS ts,Vec *Frhs) 301 { 302 PetscErrorCode ierr; 303 304 PetscFunctionBegin; 305 if (!ts->Frhs) { 306 ierr = VecDuplicate(ts->vec_sol,&ts->Frhs);CHKERRQ(ierr); 307 } 308 *Frhs = ts->Frhs; 309 PetscFunctionReturn(0); 310 } 311 312 #undef __FUNCT__ 313 #define __FUNCT__ "TSGetRHSMats_Private" 314 static PetscErrorCode TSGetRHSMats_Private(TS ts,Mat *Arhs,Mat *Brhs) 315 { 316 PetscErrorCode ierr; 317 Mat A,B; 318 319 PetscFunctionBegin; 320 ierr = TSGetIJacobian(ts,&A,&B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 321 if (Arhs) { 322 if (!ts->Arhs) { 323 ierr = MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&ts->Arhs);CHKERRQ(ierr); 324 } 325 *Arhs = ts->Arhs; 326 } 327 if (Brhs) { 328 if (!ts->Brhs) { 329 ierr = MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&ts->Brhs);CHKERRQ(ierr); 330 } 331 *Brhs = ts->Brhs; 332 } 333 PetscFunctionReturn(0); 334 } 335 336 #undef __FUNCT__ 337 #define __FUNCT__ "TSComputeIFunction" 338 /*@ 339 TSComputeIFunction - Evaluates the DAE residual written in implicit form F(t,X,Xdot)=0 340 341 Collective on TS and Vec 342 343 Input Parameters: 344 + ts - the TS context 345 . t - current time 346 . X - state vector 347 . Xdot - time derivative of state vector 348 - imex - flag indicates if the method is IMEX so that the RHSFunction should be kept separate 349 350 Output Parameter: 351 . Y - right hand side 352 353 Note: 354 Most users should not need to explicitly call this routine, as it 355 is used internally within the nonlinear solvers. 356 357 If the user did did not write their equations in implicit form, this 358 function recasts them in implicit form. 359 360 Level: developer 361 362 .keywords: TS, compute 363 364 .seealso: TSSetIFunction(), TSComputeRHSFunction() 365 @*/ 366 PetscErrorCode TSComputeIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec Y,PetscBool imex) 367 { 368 PetscErrorCode ierr; 369 370 PetscFunctionBegin; 371 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 372 PetscValidHeaderSpecific(X,VEC_CLASSID,3); 373 PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4); 374 PetscValidHeaderSpecific(Y,VEC_CLASSID,5); 375 376 if (!ts->userops->rhsfunction && !ts->userops->ifunction) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSFunction() and / or TSSetIFunction()"); 377 378 ierr = PetscLogEventBegin(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr); 379 if (ts->userops->ifunction) { 380 PetscStackPush("TS user implicit function"); 381 ierr = (*ts->userops->ifunction)(ts,t,X,Xdot,Y,ts->funP);CHKERRQ(ierr); 382 PetscStackPop; 383 } 384 if (imex) { 385 if (!ts->userops->ifunction) {ierr = VecCopy(Xdot,Y);CHKERRQ(ierr);} 386 } else { 387 if (!ts->userops->ifunction) { 388 ierr = TSComputeRHSFunction(ts,t,X,Y);CHKERRQ(ierr); 389 ierr = VecAYPX(Y,-1,Xdot);CHKERRQ(ierr); 390 } else { 391 Vec Frhs; 392 ierr = TSGetRHSVec_Private(ts,&Frhs);CHKERRQ(ierr); 393 ierr = TSComputeRHSFunction(ts,t,X,Frhs);CHKERRQ(ierr); 394 ierr = VecAXPY(Y,-1,Frhs);CHKERRQ(ierr); 395 } 396 } 397 ierr = PetscLogEventEnd(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr); 398 PetscFunctionReturn(0); 399 } 400 401 #undef __FUNCT__ 402 #define __FUNCT__ "TSComputeIJacobian" 403 /*@ 404 TSComputeIJacobian - Evaluates the Jacobian of the DAE 405 406 Collective on TS and Vec 407 408 Input 409 Input Parameters: 410 + ts - the TS context 411 . t - current timestep 412 . X - state vector 413 . Xdot - time derivative of state vector 414 . shift - shift to apply, see note below 415 - imex - flag indicates if the method is IMEX so that the RHSJacobian should be kept separate 416 417 Output Parameters: 418 + A - Jacobian matrix 419 . B - optional preconditioning matrix 420 - flag - flag indicating matrix structure 421 422 Notes: 423 If F(t,X,Xdot)=0 is the DAE, the required Jacobian is 424 425 dF/dX + shift*dF/dXdot 426 427 Most users should not need to explicitly call this routine, as it 428 is used internally within the nonlinear solvers. 429 430 Level: developer 431 432 .keywords: TS, compute, Jacobian, matrix 433 434 .seealso: TSSetIJacobian() 435 @*/ 436 PetscErrorCode TSComputeIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg,PetscBool imex) 437 { 438 PetscInt Xstate, Xdotstate; 439 PetscErrorCode ierr; 440 441 PetscFunctionBegin; 442 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 443 PetscValidHeaderSpecific(X,VEC_CLASSID,3); 444 PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4); 445 PetscValidPointer(A,6); 446 PetscValidHeaderSpecific(*A,MAT_CLASSID,6); 447 PetscValidPointer(B,7); 448 PetscValidHeaderSpecific(*B,MAT_CLASSID,7); 449 PetscValidPointer(flg,8); 450 ierr = PetscObjectStateQuery((PetscObject)X,&Xstate);CHKERRQ(ierr); 451 ierr = PetscObjectStateQuery((PetscObject)Xdot,&Xdotstate);CHKERRQ(ierr); 452 if (ts->ijacobian.time == t && (ts->problem_type == TS_LINEAR || (ts->ijacobian.X == X && ts->ijacobian.Xstate == Xstate && ts->ijacobian.Xdot == Xdot && ts->ijacobian.Xdotstate == Xdotstate && ts->ijacobian.imex == imex))) { 453 *flg = ts->ijacobian.mstructure; 454 ierr = MatScale(*A, shift / ts->ijacobian.shift);CHKERRQ(ierr); 455 PetscFunctionReturn(0); 456 } 457 458 if (!ts->userops->rhsjacobian && !ts->userops->ijacobian) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSJacobian() and / or TSSetIJacobian()"); 459 460 *flg = SAME_NONZERO_PATTERN; /* In case we're solving a linear problem in which case it wouldn't get initialized below. */ 461 ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 462 if (ts->userops->ijacobian) { 463 PetscStackPush("TS user implicit Jacobian"); 464 ierr = (*ts->userops->ijacobian)(ts,t,X,Xdot,shift,A,B,flg,ts->jacP);CHKERRQ(ierr); 465 PetscStackPop; 466 /* make sure user returned a correct Jacobian and preconditioner */ 467 PetscValidHeaderSpecific(*A,MAT_CLASSID,4); 468 PetscValidHeaderSpecific(*B,MAT_CLASSID,5); 469 } 470 if (imex) { 471 if (!ts->userops->ijacobian) { /* system was written as Xdot = F(t,X) */ 472 ierr = MatZeroEntries(*A);CHKERRQ(ierr); 473 ierr = MatShift(*A,1.0);CHKERRQ(ierr); 474 if (*A != *B) { 475 ierr = MatZeroEntries(*B);CHKERRQ(ierr); 476 ierr = MatShift(*B,shift);CHKERRQ(ierr); 477 } 478 *flg = SAME_PRECONDITIONER; 479 } 480 } else { 481 if (!ts->userops->ijacobian) { 482 ierr = TSComputeRHSJacobian(ts,t,X,A,B,flg);CHKERRQ(ierr); 483 ierr = MatScale(*A,-1);CHKERRQ(ierr); 484 ierr = MatShift(*A,shift);CHKERRQ(ierr); 485 if (*A != *B) { 486 ierr = MatScale(*B,-1);CHKERRQ(ierr); 487 ierr = MatShift(*B,shift);CHKERRQ(ierr); 488 } 489 } else if (ts->userops->rhsjacobian) { 490 Mat Arhs,Brhs; 491 MatStructure axpy,flg2 = DIFFERENT_NONZERO_PATTERN; 492 ierr = TSGetRHSMats_Private(ts,&Arhs,&Brhs);CHKERRQ(ierr); 493 ierr = TSComputeRHSJacobian(ts,t,X,&Arhs,&Brhs,&flg2);CHKERRQ(ierr); 494 axpy = (*flg == flg2) ? SAME_NONZERO_PATTERN : DIFFERENT_NONZERO_PATTERN; 495 ierr = MatAXPY(*A,-1,Arhs,axpy);CHKERRQ(ierr); 496 if (*A != *B) { 497 ierr = MatAXPY(*B,-1,Brhs,axpy);CHKERRQ(ierr); 498 } 499 *flg = PetscMin(*flg,flg2); 500 } 501 } 502 503 ts->ijacobian.time = t; 504 ts->ijacobian.X = X; 505 ts->ijacobian.Xdot = Xdot; 506 ierr = PetscObjectStateQuery((PetscObject)X,&ts->ijacobian.Xstate);CHKERRQ(ierr); 507 ierr = PetscObjectStateQuery((PetscObject)Xdot,&ts->ijacobian.Xdotstate);CHKERRQ(ierr); 508 ts->ijacobian.shift = shift; 509 ts->ijacobian.imex = imex; 510 ts->ijacobian.mstructure = *flg; 511 ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 512 PetscFunctionReturn(0); 513 } 514 515 #undef __FUNCT__ 516 #define __FUNCT__ "TSSetRHSFunction" 517 /*@C 518 TSSetRHSFunction - Sets the routine for evaluating the function, 519 F(t,u), where U_t = F(t,u). 520 521 Logically Collective on TS 522 523 Input Parameters: 524 + ts - the TS context obtained from TSCreate() 525 . r - vector to put the computed right hand side (or PETSC_NULL to have it created) 526 . f - routine for evaluating the right-hand-side function 527 - ctx - [optional] user-defined context for private data for the 528 function evaluation routine (may be PETSC_NULL) 529 530 Calling sequence of func: 531 $ func (TS ts,PetscReal t,Vec u,Vec F,void *ctx); 532 533 + t - current timestep 534 . u - input vector 535 . F - function vector 536 - ctx - [optional] user-defined function context 537 538 Important: 539 The user MUST call either this routine or TSSetMatrices(). 540 541 Level: beginner 542 543 .keywords: TS, timestep, set, right-hand-side, function 544 545 .seealso: TSSetMatrices() 546 @*/ 547 PetscErrorCode TSSetRHSFunction(TS ts,Vec r,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx) 548 { 549 PetscErrorCode ierr; 550 SNES snes; 551 552 PetscFunctionBegin; 553 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 554 if (r) PetscValidHeaderSpecific(r,VEC_CLASSID,2); 555 if (f) ts->userops->rhsfunction = f; 556 if (ctx) ts->funP = ctx; 557 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 558 ierr = SNESSetFunction(snes,r,SNESTSFormFunction,ts);CHKERRQ(ierr); 559 PetscFunctionReturn(0); 560 } 561 562 #undef __FUNCT__ 563 #define __FUNCT__ "TSSetRHSJacobian" 564 /*@C 565 TSSetRHSJacobian - Sets the function to compute the Jacobian of F, 566 where U_t = F(U,t), as well as the location to store the matrix. 567 Use TSSetMatrices() for linear problems. 568 569 Logically Collective on TS 570 571 Input Parameters: 572 + ts - the TS context obtained from TSCreate() 573 . A - Jacobian matrix 574 . B - preconditioner matrix (usually same as A) 575 . f - the Jacobian evaluation routine 576 - ctx - [optional] user-defined context for private data for the 577 Jacobian evaluation routine (may be PETSC_NULL) 578 579 Calling sequence of func: 580 $ func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx); 581 582 + t - current timestep 583 . u - input vector 584 . A - matrix A, where U_t = A(t)u 585 . B - preconditioner matrix, usually the same as A 586 . flag - flag indicating information about the preconditioner matrix 587 structure (same as flag in KSPSetOperators()) 588 - ctx - [optional] user-defined context for matrix evaluation routine 589 590 Notes: 591 See KSPSetOperators() for important information about setting the flag 592 output parameter in the routine func(). Be sure to read this information! 593 594 The routine func() takes Mat * as the matrix arguments rather than Mat. 595 This allows the matrix evaluation routine to replace A and/or B with a 596 completely new matrix structure (not just different matrix elements) 597 when appropriate, for instance, if the nonzero structure is changing 598 throughout the global iterations. 599 600 Level: beginner 601 602 .keywords: TS, timestep, set, right-hand-side, Jacobian 603 604 .seealso: TSDefaultComputeJacobianColor(), 605 SNESDefaultComputeJacobianColor(), TSSetRHSFunction(), TSSetMatrices() 606 607 @*/ 608 PetscErrorCode TSSetRHSJacobian(TS ts,Mat A,Mat B,TSRHSJacobian f,void *ctx) 609 { 610 PetscErrorCode ierr; 611 SNES snes; 612 613 PetscFunctionBegin; 614 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 615 if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2); 616 if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3); 617 if (A) PetscCheckSameComm(ts,1,A,2); 618 if (B) PetscCheckSameComm(ts,1,B,3); 619 620 if (f) ts->userops->rhsjacobian = f; 621 if (ctx) ts->jacP = ctx; 622 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 623 if (!ts->userops->ijacobian) { 624 ierr = SNESSetJacobian(snes,A,B,SNESTSFormJacobian,ts);CHKERRQ(ierr); 625 } 626 if (A) { 627 ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 628 ierr = MatDestroy(&ts->Arhs);CHKERRQ(ierr); 629 ts->Arhs = A; 630 } 631 if (B) { 632 ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 633 ierr = MatDestroy(&ts->Brhs);CHKERRQ(ierr); 634 ts->Brhs = B; 635 } 636 PetscFunctionReturn(0); 637 } 638 639 640 #undef __FUNCT__ 641 #define __FUNCT__ "TSSetIFunction" 642 /*@C 643 TSSetIFunction - Set the function to compute F(t,U,U_t) where F = 0 is the DAE to be solved. 644 645 Logically Collective on TS 646 647 Input Parameters: 648 + ts - the TS context obtained from TSCreate() 649 . r - vector to hold the residual (or PETSC_NULL to have it created internally) 650 . f - the function evaluation routine 651 - ctx - user-defined context for private data for the function evaluation routine (may be PETSC_NULL) 652 653 Calling sequence of f: 654 $ f(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx); 655 656 + t - time at step/stage being solved 657 . u - state vector 658 . u_t - time derivative of state vector 659 . F - function vector 660 - ctx - [optional] user-defined context for matrix evaluation routine 661 662 Important: 663 The user MUST call either this routine, TSSetRHSFunction(), or TSSetMatrices(). This routine must be used when not solving an ODE. 664 665 Level: beginner 666 667 .keywords: TS, timestep, set, DAE, Jacobian 668 669 .seealso: TSSetMatrices(), TSSetRHSFunction(), TSSetIJacobian() 670 @*/ 671 PetscErrorCode TSSetIFunction(TS ts,Vec res,TSIFunction f,void *ctx) 672 { 673 PetscErrorCode ierr; 674 SNES snes; 675 676 PetscFunctionBegin; 677 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 678 if (res) PetscValidHeaderSpecific(res,VEC_CLASSID,2); 679 if (f) ts->userops->ifunction = f; 680 if (ctx) ts->funP = ctx; 681 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 682 ierr = SNESSetFunction(snes,res,SNESTSFormFunction,ts);CHKERRQ(ierr); 683 PetscFunctionReturn(0); 684 } 685 686 #undef __FUNCT__ 687 #define __FUNCT__ "TSGetIFunction" 688 /*@C 689 TSGetIFunction - Returns the vector where the implicit residual is stored and the function/contex to compute it. 690 691 Not Collective 692 693 Input Parameter: 694 . ts - the TS context 695 696 Output Parameter: 697 + r - vector to hold residual (or PETSC_NULL) 698 . func - the function to compute residual (or PETSC_NULL) 699 - ctx - the function context (or PETSC_NULL) 700 701 Level: advanced 702 703 .keywords: TS, nonlinear, get, function 704 705 .seealso: TSSetIFunction(), SNESGetFunction() 706 @*/ 707 PetscErrorCode TSGetIFunction(TS ts,Vec *r,TSIFunction *func,void **ctx) 708 { 709 PetscErrorCode ierr; 710 SNES snes; 711 712 PetscFunctionBegin; 713 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 714 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 715 ierr = SNESGetFunction(snes,r,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 716 if (func) *func = ts->userops->ifunction; 717 if (ctx) *ctx = ts->funP; 718 PetscFunctionReturn(0); 719 } 720 721 #undef __FUNCT__ 722 #define __FUNCT__ "TSGetRHSFunction" 723 /*@C 724 TSGetRHSFunction - Returns the vector where the right hand side is stored and the function/context to compute it. 725 726 Not Collective 727 728 Input Parameter: 729 . ts - the TS context 730 731 Output Parameter: 732 + r - vector to hold computed right hand side (or PETSC_NULL) 733 . func - the function to compute right hand side (or PETSC_NULL) 734 - ctx - the function context (or PETSC_NULL) 735 736 Level: advanced 737 738 .keywords: TS, nonlinear, get, function 739 740 .seealso: TSSetRhsfunction(), SNESGetFunction() 741 @*/ 742 PetscErrorCode TSGetRHSFunction(TS ts,Vec *r,TSRHSFunction *func,void **ctx) 743 { 744 PetscErrorCode ierr; 745 SNES snes; 746 747 PetscFunctionBegin; 748 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 749 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 750 ierr = SNESGetFunction(snes,r,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 751 if (func) *func = ts->userops->rhsfunction; 752 if (ctx) *ctx = ts->funP; 753 PetscFunctionReturn(0); 754 } 755 756 #undef __FUNCT__ 757 #define __FUNCT__ "TSSetIJacobian" 758 /*@C 759 TSSetIJacobian - Set the function to compute the matrix dF/dU + a*dF/dU_t where F(t,U,U_t) is the function 760 you provided with TSSetIFunction(). 761 762 Logically Collective on TS 763 764 Input Parameters: 765 + ts - the TS context obtained from TSCreate() 766 . A - Jacobian matrix 767 . B - preconditioning matrix for A (may be same as A) 768 . f - the Jacobian evaluation routine 769 - ctx - user-defined context for private data for the Jacobian evaluation routine (may be PETSC_NULL) 770 771 Calling sequence of f: 772 $ f(TS ts,PetscReal t,Vec U,Vec U_t,PetscReal a,Mat *A,Mat *B,MatStructure *flag,void *ctx); 773 774 + t - time at step/stage being solved 775 . U - state vector 776 . U_t - time derivative of state vector 777 . a - shift 778 . A - Jacobian of G(U) = F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t 779 . B - preconditioning matrix for A, may be same as A 780 . flag - flag indicating information about the preconditioner matrix 781 structure (same as flag in KSPSetOperators()) 782 - ctx - [optional] user-defined context for matrix evaluation routine 783 784 Notes: 785 The matrices A and B are exactly the matrices that are used by SNES for the nonlinear solve. 786 787 The matrix dF/dU + a*dF/dU_t you provide turns out to be 788 the Jacobian of G(U) = F(t,U,W+a*U) where F(t,U,U_t) = 0 is the DAE to be solved. 789 The time integrator internally approximates U_t by W+a*U where the positive "shift" 790 a and vector W depend on the integration method, step size, and past states. For example with 791 the backward Euler method a = 1/dt and W = -a*U(previous timestep) so 792 W + a*U = a*(U - U(previous timestep)) = (U - U(previous timestep))/dt 793 794 Level: beginner 795 796 .keywords: TS, timestep, DAE, Jacobian 797 798 .seealso: TSSetIFunction(), TSSetRHSJacobian() 799 800 @*/ 801 PetscErrorCode TSSetIJacobian(TS ts,Mat A,Mat B,TSIJacobian f,void *ctx) 802 { 803 PetscErrorCode ierr; 804 SNES snes; 805 806 PetscFunctionBegin; 807 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 808 if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2); 809 if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3); 810 if (A) PetscCheckSameComm(ts,1,A,2); 811 if (B) PetscCheckSameComm(ts,1,B,3); 812 if (f) ts->userops->ijacobian = f; 813 if (ctx) ts->jacP = ctx; 814 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 815 ierr = SNESSetJacobian(snes,A,B,SNESTSFormJacobian,ts);CHKERRQ(ierr); 816 PetscFunctionReturn(0); 817 } 818 819 #undef __FUNCT__ 820 #define __FUNCT__ "TSView" 821 /*@C 822 TSView - Prints the TS data structure. 823 824 Collective on TS 825 826 Input Parameters: 827 + ts - the TS context obtained from TSCreate() 828 - viewer - visualization context 829 830 Options Database Key: 831 . -ts_view - calls TSView() at end of TSStep() 832 833 Notes: 834 The available visualization contexts include 835 + PETSC_VIEWER_STDOUT_SELF - standard output (default) 836 - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 837 output where only the first processor opens 838 the file. All other processors send their 839 data to the first processor to print. 840 841 The user can open an alternative visualization context with 842 PetscViewerASCIIOpen() - output to a specified file. 843 844 Level: beginner 845 846 .keywords: TS, timestep, view 847 848 .seealso: PetscViewerASCIIOpen() 849 @*/ 850 PetscErrorCode TSView(TS ts,PetscViewer viewer) 851 { 852 PetscErrorCode ierr; 853 const TSType type; 854 PetscBool iascii,isstring,isundials; 855 856 PetscFunctionBegin; 857 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 858 if (!viewer) { 859 ierr = PetscViewerASCIIGetStdout(((PetscObject)ts)->comm,&viewer);CHKERRQ(ierr); 860 } 861 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 862 PetscCheckSameComm(ts,1,viewer,2); 863 864 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 865 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 866 if (iascii) { 867 ierr = PetscObjectPrintClassNamePrefixType((PetscObject)ts,viewer,"TS Object");CHKERRQ(ierr); 868 if (ts->ops->view) { 869 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 870 ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr); 871 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 872 } 873 ierr = PetscViewerASCIIPrintf(viewer," maximum steps=%D\n",ts->max_steps);CHKERRQ(ierr); 874 ierr = PetscViewerASCIIPrintf(viewer," maximum time=%G\n",ts->max_time);CHKERRQ(ierr); 875 if (ts->problem_type == TS_NONLINEAR) { 876 ierr = PetscViewerASCIIPrintf(viewer," total number of nonlinear solver iterations=%D\n",ts->nonlinear_its);CHKERRQ(ierr); 877 ierr = PetscViewerASCIIPrintf(viewer," total number of nonlinear solve failures=%D\n",ts->max_snes_failures);CHKERRQ(ierr); 878 } 879 ierr = PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%D\n",ts->linear_its);CHKERRQ(ierr); 880 ierr = PetscViewerASCIIPrintf(viewer," total number of rejected steps=%D\n",ts->reject);CHKERRQ(ierr); 881 } else if (isstring) { 882 ierr = TSGetType(ts,&type);CHKERRQ(ierr); 883 ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr); 884 } 885 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 886 ierr = PetscTypeCompare((PetscObject)ts,TSSUNDIALS,&isundials);CHKERRQ(ierr); 887 if (!isundials && ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);} 888 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 889 PetscFunctionReturn(0); 890 } 891 892 893 #undef __FUNCT__ 894 #define __FUNCT__ "TSSetApplicationContext" 895 /*@ 896 TSSetApplicationContext - Sets an optional user-defined context for 897 the timesteppers. 898 899 Logically Collective on TS 900 901 Input Parameters: 902 + ts - the TS context obtained from TSCreate() 903 - usrP - optional user context 904 905 Level: intermediate 906 907 .keywords: TS, timestep, set, application, context 908 909 .seealso: TSGetApplicationContext() 910 @*/ 911 PetscErrorCode TSSetApplicationContext(TS ts,void *usrP) 912 { 913 PetscFunctionBegin; 914 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 915 ts->user = usrP; 916 PetscFunctionReturn(0); 917 } 918 919 #undef __FUNCT__ 920 #define __FUNCT__ "TSGetApplicationContext" 921 /*@ 922 TSGetApplicationContext - Gets the user-defined context for the 923 timestepper. 924 925 Not Collective 926 927 Input Parameter: 928 . ts - the TS context obtained from TSCreate() 929 930 Output Parameter: 931 . usrP - user context 932 933 Level: intermediate 934 935 .keywords: TS, timestep, get, application, context 936 937 .seealso: TSSetApplicationContext() 938 @*/ 939 PetscErrorCode TSGetApplicationContext(TS ts,void *usrP) 940 { 941 PetscFunctionBegin; 942 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 943 *(void**)usrP = ts->user; 944 PetscFunctionReturn(0); 945 } 946 947 #undef __FUNCT__ 948 #define __FUNCT__ "TSGetTimeStepNumber" 949 /*@ 950 TSGetTimeStepNumber - Gets the current number of timesteps. 951 952 Not Collective 953 954 Input Parameter: 955 . ts - the TS context obtained from TSCreate() 956 957 Output Parameter: 958 . iter - number steps so far 959 960 Level: intermediate 961 962 .keywords: TS, timestep, get, iteration, number 963 @*/ 964 PetscErrorCode TSGetTimeStepNumber(TS ts,PetscInt* iter) 965 { 966 PetscFunctionBegin; 967 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 968 PetscValidIntPointer(iter,2); 969 *iter = ts->steps; 970 PetscFunctionReturn(0); 971 } 972 973 #undef __FUNCT__ 974 #define __FUNCT__ "TSSetInitialTimeStep" 975 /*@ 976 TSSetInitialTimeStep - Sets the initial timestep to be used, 977 as well as the initial time. 978 979 Logically Collective on TS 980 981 Input Parameters: 982 + ts - the TS context obtained from TSCreate() 983 . initial_time - the initial time 984 - time_step - the size of the timestep 985 986 Level: intermediate 987 988 .seealso: TSSetTimeStep(), TSGetTimeStep() 989 990 .keywords: TS, set, initial, timestep 991 @*/ 992 PetscErrorCode TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step) 993 { 994 PetscErrorCode ierr; 995 996 PetscFunctionBegin; 997 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 998 ierr = TSSetTimeStep(ts,time_step);CHKERRQ(ierr); 999 ts->initial_time_step = time_step; 1000 ts->ptime = initial_time; 1001 PetscFunctionReturn(0); 1002 } 1003 1004 #undef __FUNCT__ 1005 #define __FUNCT__ "TSSetTimeStep" 1006 /*@ 1007 TSSetTimeStep - Allows one to reset the timestep at any time, 1008 useful for simple pseudo-timestepping codes. 1009 1010 Logically Collective on TS 1011 1012 Input Parameters: 1013 + ts - the TS context obtained from TSCreate() 1014 - time_step - the size of the timestep 1015 1016 Level: intermediate 1017 1018 .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 1019 1020 .keywords: TS, set, timestep 1021 @*/ 1022 PetscErrorCode TSSetTimeStep(TS ts,PetscReal time_step) 1023 { 1024 PetscFunctionBegin; 1025 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1026 PetscValidLogicalCollectiveReal(ts,time_step,2); 1027 ts->time_step = time_step; 1028 ts->next_time_step = time_step; 1029 PetscFunctionReturn(0); 1030 } 1031 1032 #undef __FUNCT__ 1033 #define __FUNCT__ "TSGetTimeStep" 1034 /*@ 1035 TSGetTimeStep - Gets the current timestep size. 1036 1037 Not Collective 1038 1039 Input Parameter: 1040 . ts - the TS context obtained from TSCreate() 1041 1042 Output Parameter: 1043 . dt - the current timestep size 1044 1045 Level: intermediate 1046 1047 .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 1048 1049 .keywords: TS, get, timestep 1050 @*/ 1051 PetscErrorCode TSGetTimeStep(TS ts,PetscReal* dt) 1052 { 1053 PetscFunctionBegin; 1054 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1055 PetscValidDoublePointer(dt,2); 1056 *dt = ts->time_step; 1057 PetscFunctionReturn(0); 1058 } 1059 1060 #undef __FUNCT__ 1061 #define __FUNCT__ "TSGetSolution" 1062 /*@ 1063 TSGetSolution - Returns the solution at the present timestep. It 1064 is valid to call this routine inside the function that you are evaluating 1065 in order to move to the new timestep. This vector not changed until 1066 the solution at the next timestep has been calculated. 1067 1068 Not Collective, but Vec returned is parallel if TS is parallel 1069 1070 Input Parameter: 1071 . ts - the TS context obtained from TSCreate() 1072 1073 Output Parameter: 1074 . v - the vector containing the solution 1075 1076 Level: intermediate 1077 1078 .seealso: TSGetTimeStep() 1079 1080 .keywords: TS, timestep, get, solution 1081 @*/ 1082 PetscErrorCode TSGetSolution(TS ts,Vec *v) 1083 { 1084 PetscFunctionBegin; 1085 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1086 PetscValidPointer(v,2); 1087 *v = ts->vec_sol; 1088 PetscFunctionReturn(0); 1089 } 1090 1091 /* ----- Routines to initialize and destroy a timestepper ---- */ 1092 #undef __FUNCT__ 1093 #define __FUNCT__ "TSSetProblemType" 1094 /*@ 1095 TSSetProblemType - Sets the type of problem to be solved. 1096 1097 Not collective 1098 1099 Input Parameters: 1100 + ts - The TS 1101 - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms 1102 .vb 1103 U_t = A U 1104 U_t = A(t) U 1105 U_t = F(t,U) 1106 .ve 1107 1108 Level: beginner 1109 1110 .keywords: TS, problem type 1111 .seealso: TSSetUp(), TSProblemType, TS 1112 @*/ 1113 PetscErrorCode TSSetProblemType(TS ts, TSProblemType type) 1114 { 1115 PetscErrorCode ierr; 1116 1117 PetscFunctionBegin; 1118 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1119 ts->problem_type = type; 1120 if (type == TS_LINEAR) { 1121 SNES snes; 1122 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 1123 ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr); 1124 } 1125 PetscFunctionReturn(0); 1126 } 1127 1128 #undef __FUNCT__ 1129 #define __FUNCT__ "TSGetProblemType" 1130 /*@C 1131 TSGetProblemType - Gets the type of problem to be solved. 1132 1133 Not collective 1134 1135 Input Parameter: 1136 . ts - The TS 1137 1138 Output Parameter: 1139 . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms 1140 .vb 1141 M U_t = A U 1142 M(t) U_t = A(t) U 1143 U_t = F(t,U) 1144 .ve 1145 1146 Level: beginner 1147 1148 .keywords: TS, problem type 1149 .seealso: TSSetUp(), TSProblemType, TS 1150 @*/ 1151 PetscErrorCode TSGetProblemType(TS ts, TSProblemType *type) 1152 { 1153 PetscFunctionBegin; 1154 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1155 PetscValidIntPointer(type,2); 1156 *type = ts->problem_type; 1157 PetscFunctionReturn(0); 1158 } 1159 1160 #undef __FUNCT__ 1161 #define __FUNCT__ "TSSetUp" 1162 /*@ 1163 TSSetUp - Sets up the internal data structures for the later use 1164 of a timestepper. 1165 1166 Collective on TS 1167 1168 Input Parameter: 1169 . ts - the TS context obtained from TSCreate() 1170 1171 Notes: 1172 For basic use of the TS solvers the user need not explicitly call 1173 TSSetUp(), since these actions will automatically occur during 1174 the call to TSStep(). However, if one wishes to control this 1175 phase separately, TSSetUp() should be called after TSCreate() 1176 and optional routines of the form TSSetXXX(), but before TSStep(). 1177 1178 Level: advanced 1179 1180 .keywords: TS, timestep, setup 1181 1182 .seealso: TSCreate(), TSStep(), TSDestroy() 1183 @*/ 1184 PetscErrorCode TSSetUp(TS ts) 1185 { 1186 PetscErrorCode ierr; 1187 1188 PetscFunctionBegin; 1189 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1190 if (ts->setupcalled) PetscFunctionReturn(0); 1191 1192 if (!((PetscObject)ts)->type_name) { 1193 ierr = TSSetType(ts,TSEULER);CHKERRQ(ierr); 1194 } 1195 1196 if (!ts->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first"); 1197 1198 if (ts->ops->setup) { 1199 ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr); 1200 } 1201 1202 ts->setupcalled = PETSC_TRUE; 1203 PetscFunctionReturn(0); 1204 } 1205 1206 #undef __FUNCT__ 1207 #define __FUNCT__ "TSReset" 1208 /*@ 1209 TSReset - Resets a TS context and removes any allocated Vecs and Mats. 1210 1211 Collective on TS 1212 1213 Input Parameter: 1214 . ts - the TS context obtained from TSCreate() 1215 1216 Level: beginner 1217 1218 .keywords: TS, timestep, reset 1219 1220 .seealso: TSCreate(), TSSetup(), TSDestroy() 1221 @*/ 1222 PetscErrorCode TSReset(TS ts) 1223 { 1224 PetscErrorCode ierr; 1225 1226 PetscFunctionBegin; 1227 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1228 if (ts->ops->reset) { 1229 ierr = (*ts->ops->reset)(ts);CHKERRQ(ierr); 1230 } 1231 if (ts->snes) {ierr = SNESReset(ts->snes);CHKERRQ(ierr);} 1232 ierr = MatDestroy(&ts->Arhs);CHKERRQ(ierr); 1233 ierr = MatDestroy(&ts->Brhs);CHKERRQ(ierr); 1234 ierr = VecDestroy(&ts->Frhs);CHKERRQ(ierr); 1235 ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr); 1236 if (ts->work) {ierr = VecDestroyVecs(ts->nwork,&ts->work);CHKERRQ(ierr);} 1237 ts->setupcalled = PETSC_FALSE; 1238 PetscFunctionReturn(0); 1239 } 1240 1241 #undef __FUNCT__ 1242 #define __FUNCT__ "TSDestroy" 1243 /*@ 1244 TSDestroy - Destroys the timestepper context that was created 1245 with TSCreate(). 1246 1247 Collective on TS 1248 1249 Input Parameter: 1250 . ts - the TS context obtained from TSCreate() 1251 1252 Level: beginner 1253 1254 .keywords: TS, timestepper, destroy 1255 1256 .seealso: TSCreate(), TSSetUp(), TSSolve() 1257 @*/ 1258 PetscErrorCode TSDestroy(TS *ts) 1259 { 1260 PetscErrorCode ierr; 1261 1262 PetscFunctionBegin; 1263 if (!*ts) PetscFunctionReturn(0); 1264 PetscValidHeaderSpecific((*ts),TS_CLASSID,1); 1265 if (--((PetscObject)(*ts))->refct > 0) {*ts = 0; PetscFunctionReturn(0);} 1266 1267 ierr = TSReset((*ts));CHKERRQ(ierr); 1268 1269 /* if memory was published with AMS then destroy it */ 1270 ierr = PetscObjectDepublish((*ts));CHKERRQ(ierr); 1271 if ((*ts)->ops->destroy) {ierr = (*(*ts)->ops->destroy)((*ts));CHKERRQ(ierr);} 1272 1273 ierr = SNESDestroy(&(*ts)->snes);CHKERRQ(ierr); 1274 ierr = DMDestroy(&(*ts)->dm);CHKERRQ(ierr); 1275 ierr = TSMonitorCancel((*ts));CHKERRQ(ierr); 1276 1277 ierr = PetscFree((*ts)->userops); 1278 1279 ierr = PetscHeaderDestroy(ts);CHKERRQ(ierr); 1280 PetscFunctionReturn(0); 1281 } 1282 1283 #undef __FUNCT__ 1284 #define __FUNCT__ "TSGetSNES" 1285 /*@ 1286 TSGetSNES - Returns the SNES (nonlinear solver) associated with 1287 a TS (timestepper) context. Valid only for nonlinear problems. 1288 1289 Not Collective, but SNES is parallel if TS is parallel 1290 1291 Input Parameter: 1292 . ts - the TS context obtained from TSCreate() 1293 1294 Output Parameter: 1295 . snes - the nonlinear solver context 1296 1297 Notes: 1298 The user can then directly manipulate the SNES context to set various 1299 options, etc. Likewise, the user can then extract and manipulate the 1300 KSP, KSP, and PC contexts as well. 1301 1302 TSGetSNES() does not work for integrators that do not use SNES; in 1303 this case TSGetSNES() returns PETSC_NULL in snes. 1304 1305 Level: beginner 1306 1307 .keywords: timestep, get, SNES 1308 @*/ 1309 PetscErrorCode TSGetSNES(TS ts,SNES *snes) 1310 { 1311 PetscErrorCode ierr; 1312 1313 PetscFunctionBegin; 1314 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1315 PetscValidPointer(snes,2); 1316 if (!ts->snes) { 1317 ierr = SNESCreate(((PetscObject)ts)->comm,&ts->snes);CHKERRQ(ierr); 1318 ierr = PetscLogObjectParent(ts,ts->snes);CHKERRQ(ierr); 1319 ierr = PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);CHKERRQ(ierr); 1320 if (ts->problem_type == TS_LINEAR) { 1321 ierr = SNESSetType(ts->snes,SNESKSPONLY);CHKERRQ(ierr); 1322 } 1323 } 1324 *snes = ts->snes; 1325 PetscFunctionReturn(0); 1326 } 1327 1328 #undef __FUNCT__ 1329 #define __FUNCT__ "TSGetKSP" 1330 /*@ 1331 TSGetKSP - Returns the KSP (linear solver) associated with 1332 a TS (timestepper) context. 1333 1334 Not Collective, but KSP is parallel if TS is parallel 1335 1336 Input Parameter: 1337 . ts - the TS context obtained from TSCreate() 1338 1339 Output Parameter: 1340 . ksp - the nonlinear solver context 1341 1342 Notes: 1343 The user can then directly manipulate the KSP context to set various 1344 options, etc. Likewise, the user can then extract and manipulate the 1345 KSP and PC contexts as well. 1346 1347 TSGetKSP() does not work for integrators that do not use KSP; 1348 in this case TSGetKSP() returns PETSC_NULL in ksp. 1349 1350 Level: beginner 1351 1352 .keywords: timestep, get, KSP 1353 @*/ 1354 PetscErrorCode TSGetKSP(TS ts,KSP *ksp) 1355 { 1356 PetscErrorCode ierr; 1357 SNES snes; 1358 1359 PetscFunctionBegin; 1360 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1361 PetscValidPointer(ksp,2); 1362 if (!((PetscObject)ts)->type_name) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"KSP is not created yet. Call TSSetType() first"); 1363 if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()"); 1364 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 1365 ierr = SNESGetKSP(snes,ksp);CHKERRQ(ierr); 1366 PetscFunctionReturn(0); 1367 } 1368 1369 /* ----------- Routines to set solver parameters ---------- */ 1370 1371 #undef __FUNCT__ 1372 #define __FUNCT__ "TSGetDuration" 1373 /*@ 1374 TSGetDuration - Gets the maximum number of timesteps to use and 1375 maximum time for iteration. 1376 1377 Not Collective 1378 1379 Input Parameters: 1380 + ts - the TS context obtained from TSCreate() 1381 . maxsteps - maximum number of iterations to use, or PETSC_NULL 1382 - maxtime - final time to iterate to, or PETSC_NULL 1383 1384 Level: intermediate 1385 1386 .keywords: TS, timestep, get, maximum, iterations, time 1387 @*/ 1388 PetscErrorCode TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime) 1389 { 1390 PetscFunctionBegin; 1391 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1392 if (maxsteps) { 1393 PetscValidIntPointer(maxsteps,2); 1394 *maxsteps = ts->max_steps; 1395 } 1396 if (maxtime ) { 1397 PetscValidScalarPointer(maxtime,3); 1398 *maxtime = ts->max_time; 1399 } 1400 PetscFunctionReturn(0); 1401 } 1402 1403 #undef __FUNCT__ 1404 #define __FUNCT__ "TSSetDuration" 1405 /*@ 1406 TSSetDuration - Sets the maximum number of timesteps to use and 1407 maximum time for iteration. 1408 1409 Logically Collective on TS 1410 1411 Input Parameters: 1412 + ts - the TS context obtained from TSCreate() 1413 . maxsteps - maximum number of iterations to use 1414 - maxtime - final time to iterate to 1415 1416 Options Database Keys: 1417 . -ts_max_steps <maxsteps> - Sets maxsteps 1418 . -ts_max_time <maxtime> - Sets maxtime 1419 1420 Notes: 1421 The default maximum number of iterations is 5000. Default time is 5.0 1422 1423 Level: intermediate 1424 1425 .keywords: TS, timestep, set, maximum, iterations 1426 @*/ 1427 PetscErrorCode TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime) 1428 { 1429 PetscFunctionBegin; 1430 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1431 PetscValidLogicalCollectiveInt(ts,maxsteps,2); 1432 PetscValidLogicalCollectiveReal(ts,maxtime,2); 1433 if (maxsteps >= 0) ts->max_steps = maxsteps; 1434 if (maxtime != PETSC_DEFAULT) ts->max_time = maxtime; 1435 PetscFunctionReturn(0); 1436 } 1437 1438 #undef __FUNCT__ 1439 #define __FUNCT__ "TSSetSolution" 1440 /*@ 1441 TSSetSolution - Sets the initial solution vector 1442 for use by the TS routines. 1443 1444 Logically Collective on TS and Vec 1445 1446 Input Parameters: 1447 + ts - the TS context obtained from TSCreate() 1448 - x - the solution vector 1449 1450 Level: beginner 1451 1452 .keywords: TS, timestep, set, solution, initial conditions 1453 @*/ 1454 PetscErrorCode TSSetSolution(TS ts,Vec x) 1455 { 1456 PetscErrorCode ierr; 1457 1458 PetscFunctionBegin; 1459 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1460 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 1461 ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr); 1462 ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr); 1463 ts->vec_sol = x; 1464 PetscFunctionReturn(0); 1465 } 1466 1467 #undef __FUNCT__ 1468 #define __FUNCT__ "TSSetPreStep" 1469 /*@C 1470 TSSetPreStep - Sets the general-purpose function 1471 called once at the beginning of each time step. 1472 1473 Logically Collective on TS 1474 1475 Input Parameters: 1476 + ts - The TS context obtained from TSCreate() 1477 - func - The function 1478 1479 Calling sequence of func: 1480 . func (TS ts); 1481 1482 Level: intermediate 1483 1484 .keywords: TS, timestep 1485 @*/ 1486 PetscErrorCode TSSetPreStep(TS ts, PetscErrorCode (*func)(TS)) 1487 { 1488 PetscFunctionBegin; 1489 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1490 ts->ops->prestep = func; 1491 PetscFunctionReturn(0); 1492 } 1493 1494 #undef __FUNCT__ 1495 #define __FUNCT__ "TSPreStep" 1496 /*@C 1497 TSPreStep - Runs the user-defined pre-step function. 1498 1499 Collective on TS 1500 1501 Input Parameters: 1502 . ts - The TS context obtained from TSCreate() 1503 1504 Notes: 1505 TSPreStep() is typically used within time stepping implementations, 1506 so most users would not generally call this routine themselves. 1507 1508 Level: developer 1509 1510 .keywords: TS, timestep 1511 @*/ 1512 PetscErrorCode TSPreStep(TS ts) 1513 { 1514 PetscErrorCode ierr; 1515 1516 PetscFunctionBegin; 1517 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1518 if (ts->ops->prestep) { 1519 PetscStackPush("TS PreStep function"); 1520 ierr = (*ts->ops->prestep)(ts);CHKERRQ(ierr); 1521 PetscStackPop; 1522 } 1523 PetscFunctionReturn(0); 1524 } 1525 1526 #undef __FUNCT__ 1527 #define __FUNCT__ "TSSetPostStep" 1528 /*@C 1529 TSSetPostStep - Sets the general-purpose function 1530 called once at the end of each time step. 1531 1532 Logically Collective on TS 1533 1534 Input Parameters: 1535 + ts - The TS context obtained from TSCreate() 1536 - func - The function 1537 1538 Calling sequence of func: 1539 . func (TS ts); 1540 1541 Level: intermediate 1542 1543 .keywords: TS, timestep 1544 @*/ 1545 PetscErrorCode TSSetPostStep(TS ts, PetscErrorCode (*func)(TS)) 1546 { 1547 PetscFunctionBegin; 1548 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1549 ts->ops->poststep = func; 1550 PetscFunctionReturn(0); 1551 } 1552 1553 #undef __FUNCT__ 1554 #define __FUNCT__ "TSPostStep" 1555 /*@C 1556 TSPostStep - Runs the user-defined post-step function. 1557 1558 Collective on TS 1559 1560 Input Parameters: 1561 . ts - The TS context obtained from TSCreate() 1562 1563 Notes: 1564 TSPostStep() is typically used within time stepping implementations, 1565 so most users would not generally call this routine themselves. 1566 1567 Level: developer 1568 1569 .keywords: TS, timestep 1570 @*/ 1571 PetscErrorCode TSPostStep(TS ts) 1572 { 1573 PetscErrorCode ierr; 1574 1575 PetscFunctionBegin; 1576 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1577 if (ts->ops->poststep) { 1578 PetscStackPush("TS PostStep function"); 1579 ierr = (*ts->ops->poststep)(ts);CHKERRQ(ierr); 1580 PetscStackPop; 1581 } 1582 PetscFunctionReturn(0); 1583 } 1584 1585 /* ------------ Routines to set performance monitoring options ----------- */ 1586 1587 #undef __FUNCT__ 1588 #define __FUNCT__ "TSMonitorSet" 1589 /*@C 1590 TSMonitorSet - Sets an ADDITIONAL function that is to be used at every 1591 timestep to display the iteration's progress. 1592 1593 Logically Collective on TS 1594 1595 Input Parameters: 1596 + ts - the TS context obtained from TSCreate() 1597 . func - monitoring routine 1598 . mctx - [optional] user-defined context for private data for the 1599 monitor routine (use PETSC_NULL if no context is desired) 1600 - monitordestroy - [optional] routine that frees monitor context 1601 (may be PETSC_NULL) 1602 1603 Calling sequence of func: 1604 $ int func(TS ts,PetscInt steps,PetscReal time,Vec x,void *mctx) 1605 1606 + ts - the TS context 1607 . steps - iteration number 1608 . time - current time 1609 . x - current iterate 1610 - mctx - [optional] monitoring context 1611 1612 Notes: 1613 This routine adds an additional monitor to the list of monitors that 1614 already has been loaded. 1615 1616 Fortran notes: Only a single monitor function can be set for each TS object 1617 1618 Level: intermediate 1619 1620 .keywords: TS, timestep, set, monitor 1621 1622 .seealso: TSMonitorDefault(), TSMonitorCancel() 1623 @*/ 1624 PetscErrorCode TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void**)) 1625 { 1626 PetscFunctionBegin; 1627 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1628 if (ts->numbermonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1629 ts->monitor[ts->numbermonitors] = monitor; 1630 ts->mdestroy[ts->numbermonitors] = mdestroy; 1631 ts->monitorcontext[ts->numbermonitors++] = (void*)mctx; 1632 PetscFunctionReturn(0); 1633 } 1634 1635 #undef __FUNCT__ 1636 #define __FUNCT__ "TSMonitorCancel" 1637 /*@C 1638 TSMonitorCancel - Clears all the monitors that have been set on a time-step object. 1639 1640 Logically Collective on TS 1641 1642 Input Parameters: 1643 . ts - the TS context obtained from TSCreate() 1644 1645 Notes: 1646 There is no way to remove a single, specific monitor. 1647 1648 Level: intermediate 1649 1650 .keywords: TS, timestep, set, monitor 1651 1652 .seealso: TSMonitorDefault(), TSMonitorSet() 1653 @*/ 1654 PetscErrorCode TSMonitorCancel(TS ts) 1655 { 1656 PetscErrorCode ierr; 1657 PetscInt i; 1658 1659 PetscFunctionBegin; 1660 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1661 for (i=0; i<ts->numbermonitors; i++) { 1662 if (ts->mdestroy[i]) { 1663 ierr = (*ts->mdestroy[i])(&ts->monitorcontext[i]);CHKERRQ(ierr); 1664 } 1665 } 1666 ts->numbermonitors = 0; 1667 PetscFunctionReturn(0); 1668 } 1669 1670 #undef __FUNCT__ 1671 #define __FUNCT__ "TSMonitorDefault" 1672 /*@ 1673 TSMonitorDefault - Sets the Default monitor 1674 1675 Level: intermediate 1676 1677 .keywords: TS, set, monitor 1678 1679 .seealso: TSMonitorDefault(), TSMonitorSet() 1680 @*/ 1681 PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy) 1682 { 1683 PetscErrorCode ierr; 1684 PetscViewer viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ts)->comm); 1685 1686 PetscFunctionBegin; 1687 ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 1688 ierr = PetscViewerASCIIPrintf(viewer,"%D TS dt %G time %G\n",step,ts->time_step,ptime);CHKERRQ(ierr); 1689 ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 1690 PetscFunctionReturn(0); 1691 } 1692 1693 #undef __FUNCT__ 1694 #define __FUNCT__ "TSStep" 1695 /*@ 1696 TSStep - Steps the requested number of timesteps. 1697 1698 Collective on TS 1699 1700 Input Parameter: 1701 . ts - the TS context obtained from TSCreate() 1702 1703 Output Parameters: 1704 + steps - number of iterations until termination 1705 - ptime - time until termination 1706 1707 Level: beginner 1708 1709 .keywords: TS, timestep, solve 1710 1711 .seealso: TSCreate(), TSSetUp(), TSDestroy() 1712 @*/ 1713 PetscErrorCode TSStep(TS ts) 1714 { 1715 PetscErrorCode ierr; 1716 1717 PetscFunctionBegin; 1718 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1719 1720 ierr = TSSetUp(ts);CHKERRQ(ierr); 1721 1722 ierr = PetscLogEventBegin(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr); 1723 ierr = (*ts->ops->step)(ts);CHKERRQ(ierr); 1724 ierr = PetscLogEventEnd(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr); 1725 PetscFunctionReturn(0); 1726 } 1727 1728 #undef __FUNCT__ 1729 #define __FUNCT__ "TSSolve" 1730 /*@ 1731 TSSolve - Steps the requested number of timesteps. 1732 1733 Collective on TS 1734 1735 Input Parameter: 1736 + ts - the TS context obtained from TSCreate() 1737 - x - the solution vector, or PETSC_NULL if it was set with TSSetSolution() 1738 1739 Level: beginner 1740 1741 .keywords: TS, timestep, solve 1742 1743 .seealso: TSCreate(), TSSetSolution(), TSStep() 1744 @*/ 1745 PetscErrorCode TSSolve(TS ts, Vec x) 1746 { 1747 PetscInt i; 1748 PetscErrorCode ierr; 1749 1750 PetscFunctionBegin; 1751 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1752 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 1753 ierr = TSSetSolution(ts,x); CHKERRQ(ierr); 1754 ierr = TSSetUp(ts); CHKERRQ(ierr); 1755 /* reset time step and iteration counters */ 1756 ts->steps = 0; 1757 ts->linear_its = 0; 1758 ts->nonlinear_its = 0; 1759 ts->reason = TS_CONVERGED_ITERATING; 1760 ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 1761 1762 if (ts->ops->solve) { /* This private interface is transitional and should be removed when all implementations are updated. */ 1763 ierr = (*ts->ops->solve)(ts);CHKERRQ(ierr); 1764 } else { 1765 if (++i >= ts->max_steps) { 1766 ts->reason = TS_CONVERGED_ITS; 1767 } else if (ts->ptime >= ts->max_time) { 1768 ts->reason = TS_CONVERGED_TIME; 1769 } 1770 /* steps the requested number of timesteps. */ 1771 for (i=0; !ts->reason; ) { 1772 ierr = TSPreStep(ts);CHKERRQ(ierr); 1773 ierr = TSStep(ts);CHKERRQ(ierr); 1774 if (ts->reason < 0) { 1775 if (ts->errorifstepfailed) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_NOT_CONVERGED,"TSStep has failed"); 1776 } else if (++i >= ts->max_steps) { 1777 ts->reason = TS_CONVERGED_ITS; 1778 } else if (ts->ptime >= ts->max_time) { 1779 ts->reason = TS_CONVERGED_TIME; 1780 } 1781 ierr = TSPostStep(ts);CHKERRQ(ierr); 1782 ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 1783 } 1784 } 1785 if (!PetscPreLoadingOn) { 1786 ierr = TSViewFromOptions(ts,((PetscObject)ts)->name);CHKERRQ(ierr); 1787 } 1788 PetscFunctionReturn(0); 1789 } 1790 1791 #undef __FUNCT__ 1792 #define __FUNCT__ "TSMonitor" 1793 /* 1794 Runs the user provided monitor routines, if they exists. 1795 */ 1796 PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x) 1797 { 1798 PetscErrorCode ierr; 1799 PetscInt i,n = ts->numbermonitors; 1800 1801 PetscFunctionBegin; 1802 for (i=0; i<n; i++) { 1803 ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr); 1804 } 1805 PetscFunctionReturn(0); 1806 } 1807 1808 /* ------------------------------------------------------------------------*/ 1809 1810 #undef __FUNCT__ 1811 #define __FUNCT__ "TSMonitorLGCreate" 1812 /*@C 1813 TSMonitorLGCreate - Creates a line graph context for use with 1814 TS to monitor convergence of preconditioned residual norms. 1815 1816 Collective on TS 1817 1818 Input Parameters: 1819 + host - the X display to open, or null for the local machine 1820 . label - the title to put in the title bar 1821 . x, y - the screen coordinates of the upper left coordinate of the window 1822 - m, n - the screen width and height in pixels 1823 1824 Output Parameter: 1825 . draw - the drawing context 1826 1827 Options Database Key: 1828 . -ts_monitor_draw - automatically sets line graph monitor 1829 1830 Notes: 1831 Use TSMonitorLGDestroy() to destroy this line graph, not PetscDrawLGDestroy(). 1832 1833 Level: intermediate 1834 1835 .keywords: TS, monitor, line graph, residual, seealso 1836 1837 .seealso: TSMonitorLGDestroy(), TSMonitorSet() 1838 1839 @*/ 1840 PetscErrorCode TSMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1841 { 1842 PetscDraw win; 1843 PetscErrorCode ierr; 1844 1845 PetscFunctionBegin; 1846 ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr); 1847 ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr); 1848 ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr); 1849 ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr); 1850 1851 ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr); 1852 PetscFunctionReturn(0); 1853 } 1854 1855 #undef __FUNCT__ 1856 #define __FUNCT__ "TSMonitorLG" 1857 PetscErrorCode TSMonitorLG(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx) 1858 { 1859 PetscDrawLG lg = (PetscDrawLG) monctx; 1860 PetscReal x,y = ptime; 1861 PetscErrorCode ierr; 1862 1863 PetscFunctionBegin; 1864 if (!monctx) { 1865 MPI_Comm comm; 1866 PetscViewer viewer; 1867 1868 ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr); 1869 viewer = PETSC_VIEWER_DRAW_(comm); 1870 ierr = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr); 1871 } 1872 1873 if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 1874 x = (PetscReal)n; 1875 ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 1876 if (n < 20 || (n % 5)) { 1877 ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 1878 } 1879 PetscFunctionReturn(0); 1880 } 1881 1882 #undef __FUNCT__ 1883 #define __FUNCT__ "TSMonitorLGDestroy" 1884 /*@C 1885 TSMonitorLGDestroy - Destroys a line graph context that was created 1886 with TSMonitorLGCreate(). 1887 1888 Collective on PetscDrawLG 1889 1890 Input Parameter: 1891 . draw - the drawing context 1892 1893 Level: intermediate 1894 1895 .keywords: TS, monitor, line graph, destroy 1896 1897 .seealso: TSMonitorLGCreate(), TSMonitorSet(), TSMonitorLG(); 1898 @*/ 1899 PetscErrorCode TSMonitorLGDestroy(PetscDrawLG *drawlg) 1900 { 1901 PetscDraw draw; 1902 PetscErrorCode ierr; 1903 1904 PetscFunctionBegin; 1905 ierr = PetscDrawLGGetDraw(*drawlg,&draw);CHKERRQ(ierr); 1906 ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr); 1907 ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr); 1908 PetscFunctionReturn(0); 1909 } 1910 1911 #undef __FUNCT__ 1912 #define __FUNCT__ "TSGetTime" 1913 /*@ 1914 TSGetTime - Gets the current time. 1915 1916 Not Collective 1917 1918 Input Parameter: 1919 . ts - the TS context obtained from TSCreate() 1920 1921 Output Parameter: 1922 . t - the current time 1923 1924 Level: beginner 1925 1926 .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 1927 1928 .keywords: TS, get, time 1929 @*/ 1930 PetscErrorCode TSGetTime(TS ts,PetscReal* t) 1931 { 1932 PetscFunctionBegin; 1933 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1934 PetscValidDoublePointer(t,2); 1935 *t = ts->ptime; 1936 PetscFunctionReturn(0); 1937 } 1938 1939 #undef __FUNCT__ 1940 #define __FUNCT__ "TSSetTime" 1941 /*@ 1942 TSSetTime - Allows one to reset the time. 1943 1944 Logically Collective on TS 1945 1946 Input Parameters: 1947 + ts - the TS context obtained from TSCreate() 1948 - time - the time 1949 1950 Level: intermediate 1951 1952 .seealso: TSGetTime(), TSSetDuration() 1953 1954 .keywords: TS, set, time 1955 @*/ 1956 PetscErrorCode TSSetTime(TS ts, PetscReal t) 1957 { 1958 PetscFunctionBegin; 1959 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1960 PetscValidLogicalCollectiveReal(ts,t,2); 1961 ts->ptime = t; 1962 PetscFunctionReturn(0); 1963 } 1964 1965 #undef __FUNCT__ 1966 #define __FUNCT__ "TSSetOptionsPrefix" 1967 /*@C 1968 TSSetOptionsPrefix - Sets the prefix used for searching for all 1969 TS options in the database. 1970 1971 Logically Collective on TS 1972 1973 Input Parameter: 1974 + ts - The TS context 1975 - prefix - The prefix to prepend to all option names 1976 1977 Notes: 1978 A hyphen (-) must NOT be given at the beginning of the prefix name. 1979 The first character of all runtime options is AUTOMATICALLY the 1980 hyphen. 1981 1982 Level: advanced 1983 1984 .keywords: TS, set, options, prefix, database 1985 1986 .seealso: TSSetFromOptions() 1987 1988 @*/ 1989 PetscErrorCode TSSetOptionsPrefix(TS ts,const char prefix[]) 1990 { 1991 PetscErrorCode ierr; 1992 SNES snes; 1993 1994 PetscFunctionBegin; 1995 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1996 ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 1997 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 1998 ierr = SNESSetOptionsPrefix(snes,prefix);CHKERRQ(ierr); 1999 PetscFunctionReturn(0); 2000 } 2001 2002 2003 #undef __FUNCT__ 2004 #define __FUNCT__ "TSAppendOptionsPrefix" 2005 /*@C 2006 TSAppendOptionsPrefix - Appends to the prefix used for searching for all 2007 TS options in the database. 2008 2009 Logically Collective on TS 2010 2011 Input Parameter: 2012 + ts - The TS context 2013 - prefix - The prefix to prepend to all option names 2014 2015 Notes: 2016 A hyphen (-) must NOT be given at the beginning of the prefix name. 2017 The first character of all runtime options is AUTOMATICALLY the 2018 hyphen. 2019 2020 Level: advanced 2021 2022 .keywords: TS, append, options, prefix, database 2023 2024 .seealso: TSGetOptionsPrefix() 2025 2026 @*/ 2027 PetscErrorCode TSAppendOptionsPrefix(TS ts,const char prefix[]) 2028 { 2029 PetscErrorCode ierr; 2030 SNES snes; 2031 2032 PetscFunctionBegin; 2033 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2034 ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 2035 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2036 ierr = SNESAppendOptionsPrefix(snes,prefix);CHKERRQ(ierr); 2037 PetscFunctionReturn(0); 2038 } 2039 2040 #undef __FUNCT__ 2041 #define __FUNCT__ "TSGetOptionsPrefix" 2042 /*@C 2043 TSGetOptionsPrefix - Sets the prefix used for searching for all 2044 TS options in the database. 2045 2046 Not Collective 2047 2048 Input Parameter: 2049 . ts - The TS context 2050 2051 Output Parameter: 2052 . prefix - A pointer to the prefix string used 2053 2054 Notes: On the fortran side, the user should pass in a string 'prifix' of 2055 sufficient length to hold the prefix. 2056 2057 Level: intermediate 2058 2059 .keywords: TS, get, options, prefix, database 2060 2061 .seealso: TSAppendOptionsPrefix() 2062 @*/ 2063 PetscErrorCode TSGetOptionsPrefix(TS ts,const char *prefix[]) 2064 { 2065 PetscErrorCode ierr; 2066 2067 PetscFunctionBegin; 2068 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2069 PetscValidPointer(prefix,2); 2070 ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 2071 PetscFunctionReturn(0); 2072 } 2073 2074 #undef __FUNCT__ 2075 #define __FUNCT__ "TSGetRHSJacobian" 2076 /*@C 2077 TSGetRHSJacobian - Returns the Jacobian J at the present timestep. 2078 2079 Not Collective, but parallel objects are returned if TS is parallel 2080 2081 Input Parameter: 2082 . ts - The TS context obtained from TSCreate() 2083 2084 Output Parameters: 2085 + J - The Jacobian J of F, where U_t = F(U,t) 2086 . M - The preconditioner matrix, usually the same as J 2087 . func - Function to compute the Jacobian of the RHS 2088 - ctx - User-defined context for Jacobian evaluation routine 2089 2090 Notes: You can pass in PETSC_NULL for any return argument you do not need. 2091 2092 Level: intermediate 2093 2094 .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber() 2095 2096 .keywords: TS, timestep, get, matrix, Jacobian 2097 @*/ 2098 PetscErrorCode TSGetRHSJacobian(TS ts,Mat *J,Mat *M,TSRHSJacobian *func,void **ctx) 2099 { 2100 PetscErrorCode ierr; 2101 SNES snes; 2102 2103 PetscFunctionBegin; 2104 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2105 ierr = SNESGetJacobian(snes,J,M,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 2106 if (func) *func = ts->userops->rhsjacobian; 2107 if (ctx) *ctx = ts->jacP; 2108 PetscFunctionReturn(0); 2109 } 2110 2111 #undef __FUNCT__ 2112 #define __FUNCT__ "TSGetIJacobian" 2113 /*@C 2114 TSGetIJacobian - Returns the implicit Jacobian at the present timestep. 2115 2116 Not Collective, but parallel objects are returned if TS is parallel 2117 2118 Input Parameter: 2119 . ts - The TS context obtained from TSCreate() 2120 2121 Output Parameters: 2122 + A - The Jacobian of F(t,U,U_t) 2123 . B - The preconditioner matrix, often the same as A 2124 . f - The function to compute the matrices 2125 - ctx - User-defined context for Jacobian evaluation routine 2126 2127 Notes: You can pass in PETSC_NULL for any return argument you do not need. 2128 2129 Level: advanced 2130 2131 .seealso: TSGetTimeStep(), TSGetRHSJacobian(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber() 2132 2133 .keywords: TS, timestep, get, matrix, Jacobian 2134 @*/ 2135 PetscErrorCode TSGetIJacobian(TS ts,Mat *A,Mat *B,TSIJacobian *f,void **ctx) 2136 { 2137 PetscErrorCode ierr; 2138 SNES snes; 2139 2140 PetscFunctionBegin; 2141 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2142 ierr = SNESGetJacobian(snes,A,B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 2143 if (f) *f = ts->userops->ijacobian; 2144 if (ctx) *ctx = ts->jacP; 2145 PetscFunctionReturn(0); 2146 } 2147 2148 #undef __FUNCT__ 2149 #define __FUNCT__ "TSMonitorSolution" 2150 /*@C 2151 TSMonitorSolution - Monitors progress of the TS solvers by calling 2152 VecView() for the solution at each timestep 2153 2154 Collective on TS 2155 2156 Input Parameters: 2157 + ts - the TS context 2158 . step - current time-step 2159 . ptime - current time 2160 - dummy - either a viewer or PETSC_NULL 2161 2162 Level: intermediate 2163 2164 .keywords: TS, vector, monitor, view 2165 2166 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView() 2167 @*/ 2168 PetscErrorCode TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy) 2169 { 2170 PetscErrorCode ierr; 2171 PetscViewer viewer = (PetscViewer) dummy; 2172 2173 PetscFunctionBegin; 2174 if (!dummy) { 2175 viewer = PETSC_VIEWER_DRAW_(((PetscObject)ts)->comm); 2176 } 2177 ierr = VecView(x,viewer);CHKERRQ(ierr); 2178 PetscFunctionReturn(0); 2179 } 2180 2181 2182 #undef __FUNCT__ 2183 #define __FUNCT__ "TSSetDM" 2184 /*@ 2185 TSSetDM - Sets the DM that may be used by some preconditioners 2186 2187 Logically Collective on TS and DM 2188 2189 Input Parameters: 2190 + ts - the preconditioner context 2191 - dm - the dm 2192 2193 Level: intermediate 2194 2195 2196 .seealso: TSGetDM(), SNESSetDM(), SNESGetDM() 2197 @*/ 2198 PetscErrorCode TSSetDM(TS ts,DM dm) 2199 { 2200 PetscErrorCode ierr; 2201 SNES snes; 2202 2203 PetscFunctionBegin; 2204 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2205 ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr); 2206 ierr = DMDestroy(&ts->dm);CHKERRQ(ierr); 2207 ts->dm = dm; 2208 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2209 ierr = SNESSetDM(snes,dm);CHKERRQ(ierr); 2210 PetscFunctionReturn(0); 2211 } 2212 2213 #undef __FUNCT__ 2214 #define __FUNCT__ "TSGetDM" 2215 /*@ 2216 TSGetDM - Gets the DM that may be used by some preconditioners 2217 2218 Not Collective 2219 2220 Input Parameter: 2221 . ts - the preconditioner context 2222 2223 Output Parameter: 2224 . dm - the dm 2225 2226 Level: intermediate 2227 2228 2229 .seealso: TSSetDM(), SNESSetDM(), SNESGetDM() 2230 @*/ 2231 PetscErrorCode TSGetDM(TS ts,DM *dm) 2232 { 2233 PetscFunctionBegin; 2234 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2235 *dm = ts->dm; 2236 PetscFunctionReturn(0); 2237 } 2238 2239 #undef __FUNCT__ 2240 #define __FUNCT__ "SNESTSFormFunction" 2241 /*@ 2242 SNESTSFormFunction - Function to evaluate nonlinear residual 2243 2244 Logically Collective on SNES 2245 2246 Input Parameter: 2247 + snes - nonlinear solver 2248 . X - the current state at which to evaluate the residual 2249 - ctx - user context, must be a TS 2250 2251 Output Parameter: 2252 . F - the nonlinear residual 2253 2254 Notes: 2255 This function is not normally called by users and is automatically registered with the SNES used by TS. 2256 It is most frequently passed to MatFDColoringSetFunction(). 2257 2258 Level: advanced 2259 2260 .seealso: SNESSetFunction(), MatFDColoringSetFunction() 2261 @*/ 2262 PetscErrorCode SNESTSFormFunction(SNES snes,Vec X,Vec F,void *ctx) 2263 { 2264 TS ts = (TS)ctx; 2265 PetscErrorCode ierr; 2266 2267 PetscFunctionBegin; 2268 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2269 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 2270 PetscValidHeaderSpecific(F,VEC_CLASSID,3); 2271 PetscValidHeaderSpecific(ts,TS_CLASSID,4); 2272 ierr = (ts->ops->snesfunction)(snes,X,F,ts);CHKERRQ(ierr); 2273 PetscFunctionReturn(0); 2274 } 2275 2276 #undef __FUNCT__ 2277 #define __FUNCT__ "SNESTSFormJacobian" 2278 /*@ 2279 SNESTSFormJacobian - Function to evaluate the Jacobian 2280 2281 Collective on SNES 2282 2283 Input Parameter: 2284 + snes - nonlinear solver 2285 . X - the current state at which to evaluate the residual 2286 - ctx - user context, must be a TS 2287 2288 Output Parameter: 2289 + A - the Jacobian 2290 . B - the preconditioning matrix (may be the same as A) 2291 - flag - indicates any structure change in the matrix 2292 2293 Notes: 2294 This function is not normally called by users and is automatically registered with the SNES used by TS. 2295 2296 Level: developer 2297 2298 .seealso: SNESSetJacobian() 2299 @*/ 2300 PetscErrorCode SNESTSFormJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flag,void *ctx) 2301 { 2302 TS ts = (TS)ctx; 2303 PetscErrorCode ierr; 2304 2305 PetscFunctionBegin; 2306 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2307 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 2308 PetscValidPointer(A,3); 2309 PetscValidHeaderSpecific(*A,MAT_CLASSID,3); 2310 PetscValidPointer(B,4); 2311 PetscValidHeaderSpecific(*B,MAT_CLASSID,4); 2312 PetscValidPointer(flag,5); 2313 PetscValidHeaderSpecific(ts,TS_CLASSID,6); 2314 ierr = (ts->ops->snesjacobian)(snes,X,A,B,flag,ts);CHKERRQ(ierr); 2315 PetscFunctionReturn(0); 2316 } 2317 2318 #undef __FUNCT__ 2319 #define __FUNCT__ "TSComputeRHSFunctionLinear" 2320 /*@C 2321 TSComputeRHSFunctionLinear - Evaluate the right hand side via the user-provided Jacobian, for linear problems only 2322 2323 Collective on TS 2324 2325 Input Arguments: 2326 + ts - time stepping context 2327 . t - time at which to evaluate 2328 . X - state at which to evaluate 2329 - ctx - context 2330 2331 Output Arguments: 2332 . F - right hand side 2333 2334 Level: intermediate 2335 2336 Notes: 2337 This function is intended to be passed to TSSetRHSFunction() to evaluate the right hand side for linear problems. 2338 The matrix (and optionally the evaluation context) should be passed to TSSetRHSJacobian(). 2339 2340 .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSJacobianConstant() 2341 @*/ 2342 PetscErrorCode TSComputeRHSFunctionLinear(TS ts,PetscReal t,Vec X,Vec F,void *ctx) 2343 { 2344 PetscErrorCode ierr; 2345 Mat Arhs,Brhs; 2346 MatStructure flg2; 2347 2348 PetscFunctionBegin; 2349 ierr = TSGetRHSMats_Private(ts,&Arhs,&Brhs);CHKERRQ(ierr); 2350 ierr = TSComputeRHSJacobian(ts,t,X,&Arhs,&Brhs,&flg2);CHKERRQ(ierr); 2351 ierr = MatMult(Arhs,X,F);CHKERRQ(ierr); 2352 PetscFunctionReturn(0); 2353 } 2354 2355 #undef __FUNCT__ 2356 #define __FUNCT__ "TSComputeRHSJacobianConstant" 2357 /*@C 2358 TSComputeRHSJacobianConstant - Reuses a Jacobian that is time-independent. 2359 2360 Collective on TS 2361 2362 Input Arguments: 2363 + ts - time stepping context 2364 . t - time at which to evaluate 2365 . X - state at which to evaluate 2366 - ctx - context 2367 2368 Output Arguments: 2369 + A - pointer to operator 2370 . B - pointer to preconditioning matrix 2371 - flg - matrix structure flag 2372 2373 Level: intermediate 2374 2375 Notes: 2376 This function is intended to be passed to TSSetRHSJacobian() to evaluate the Jacobian for linear time-independent problems. 2377 2378 .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSFunctionLinear() 2379 @*/ 2380 PetscErrorCode TSComputeRHSJacobianConstant(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg,void *ctx) 2381 { 2382 2383 PetscFunctionBegin; 2384 *flg = SAME_PRECONDITIONER; 2385 PetscFunctionReturn(0); 2386 } 2387 2388 #undef __FUNCT__ 2389 #define __FUNCT__ "TSComputeIFunctionLinear" 2390 /*@C 2391 TSComputeIFunctionLinear - Evaluate the left hand side via the user-provided Jacobian, for linear problems only 2392 2393 Collective on TS 2394 2395 Input Arguments: 2396 + ts - time stepping context 2397 . t - time at which to evaluate 2398 . X - state at which to evaluate 2399 . Xdot - time derivative of state vector 2400 - ctx - context 2401 2402 Output Arguments: 2403 . F - left hand side 2404 2405 Level: intermediate 2406 2407 Notes: 2408 The assumption here is that the left hand side is of the form A*Xdot (and not A*Xdot + B*X). For other cases, the 2409 user is required to write their own TSComputeIFunction. 2410 This function is intended to be passed to TSSetIFunction() to evaluate the left hand side for linear problems. 2411 The matrix (and optionally the evaluation context) should be passed to TSSetIJacobian(). 2412 2413 .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIJacobianConstant() 2414 @*/ 2415 PetscErrorCode TSComputeIFunctionLinear(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx) 2416 { 2417 PetscErrorCode ierr; 2418 Mat A,B; 2419 MatStructure flg2; 2420 2421 PetscFunctionBegin; 2422 ierr = TSGetIJacobian(ts,&A,&B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 2423 ierr = TSComputeIJacobian(ts,t,X,Xdot,1.0,&A,&B,&flg2,PETSC_TRUE);CHKERRQ(ierr); 2424 ierr = MatMult(A,Xdot,F);CHKERRQ(ierr); 2425 PetscFunctionReturn(0); 2426 } 2427 2428 #undef __FUNCT__ 2429 #define __FUNCT__ "TSComputeIJacobianConstant" 2430 /*@C 2431 TSComputeRHSJacobianConstant - Reuses a Jacobian that is time-independent. 2432 2433 Collective on TS 2434 2435 Input Arguments: 2436 + ts - time stepping context 2437 . t - time at which to evaluate 2438 . X - state at which to evaluate 2439 . Xdot - time derivative of state vector 2440 . shift - shift to apply 2441 - ctx - context 2442 2443 Output Arguments: 2444 + A - pointer to operator 2445 . B - pointer to preconditioning matrix 2446 - flg - matrix structure flag 2447 2448 Level: intermediate 2449 2450 Notes: 2451 This function is intended to be passed to TSSetIJacobian() to evaluate the Jacobian for linear time-independent problems. 2452 2453 .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIFunctionLinear() 2454 @*/ 2455 PetscErrorCode TSComputeIJacobianConstant(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg,void *ctx) 2456 { 2457 2458 PetscFunctionBegin; 2459 *flg = SAME_PRECONDITIONER; 2460 PetscFunctionReturn(0); 2461 } 2462 2463 #if defined(PETSC_HAVE_MATLAB_ENGINE) 2464 #include <mex.h> 2465 2466 typedef struct {char *funcname; mxArray *ctx;} TSMatlabContext; 2467 2468 #undef __FUNCT__ 2469 #define __FUNCT__ "TSComputeFunction_Matlab" 2470 /* 2471 TSComputeFunction_Matlab - Calls the function that has been set with 2472 TSSetFunctionMatlab(). 2473 2474 Collective on TS 2475 2476 Input Parameters: 2477 + snes - the TS context 2478 - x - input vector 2479 2480 Output Parameter: 2481 . y - function vector, as set by TSSetFunction() 2482 2483 Notes: 2484 TSComputeFunction() is typically used within nonlinear solvers 2485 implementations, so most users would not generally call this routine 2486 themselves. 2487 2488 Level: developer 2489 2490 .keywords: TS, nonlinear, compute, function 2491 2492 .seealso: TSSetFunction(), TSGetFunction() 2493 */ 2494 PetscErrorCode TSComputeFunction_Matlab(TS snes,PetscReal time,Vec x,Vec xdot,Vec y, void *ctx) 2495 { 2496 PetscErrorCode ierr; 2497 TSMatlabContext *sctx = (TSMatlabContext *)ctx; 2498 int nlhs = 1,nrhs = 7; 2499 mxArray *plhs[1],*prhs[7]; 2500 long long int lx = 0,lxdot = 0,ly = 0,ls = 0; 2501 2502 PetscFunctionBegin; 2503 PetscValidHeaderSpecific(snes,TS_CLASSID,1); 2504 PetscValidHeaderSpecific(x,VEC_CLASSID,3); 2505 PetscValidHeaderSpecific(xdot,VEC_CLASSID,4); 2506 PetscValidHeaderSpecific(y,VEC_CLASSID,5); 2507 PetscCheckSameComm(snes,1,x,3); 2508 PetscCheckSameComm(snes,1,y,5); 2509 2510 ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 2511 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 2512 ierr = PetscMemcpy(&lxdot,&xdot,sizeof(xdot));CHKERRQ(ierr); 2513 ierr = PetscMemcpy(&ly,&y,sizeof(x));CHKERRQ(ierr); 2514 prhs[0] = mxCreateDoubleScalar((double)ls); 2515 prhs[1] = mxCreateDoubleScalar(time); 2516 prhs[2] = mxCreateDoubleScalar((double)lx); 2517 prhs[3] = mxCreateDoubleScalar((double)lxdot); 2518 prhs[4] = mxCreateDoubleScalar((double)ly); 2519 prhs[5] = mxCreateString(sctx->funcname); 2520 prhs[6] = sctx->ctx; 2521 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeFunctionInternal");CHKERRQ(ierr); 2522 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 2523 mxDestroyArray(prhs[0]); 2524 mxDestroyArray(prhs[1]); 2525 mxDestroyArray(prhs[2]); 2526 mxDestroyArray(prhs[3]); 2527 mxDestroyArray(prhs[4]); 2528 mxDestroyArray(prhs[5]); 2529 mxDestroyArray(plhs[0]); 2530 PetscFunctionReturn(0); 2531 } 2532 2533 2534 #undef __FUNCT__ 2535 #define __FUNCT__ "TSSetFunctionMatlab" 2536 /* 2537 TSSetFunctionMatlab - Sets the function evaluation routine and function 2538 vector for use by the TS routines in solving ODEs 2539 equations from MATLAB. Here the function is a string containing the name of a MATLAB function 2540 2541 Logically Collective on TS 2542 2543 Input Parameters: 2544 + ts - the TS context 2545 - func - function evaluation routine 2546 2547 Calling sequence of func: 2548 $ func (TS ts,PetscReal time,Vec x,Vec xdot,Vec f,void *ctx); 2549 2550 Level: beginner 2551 2552 .keywords: TS, nonlinear, set, function 2553 2554 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction() 2555 */ 2556 PetscErrorCode TSSetFunctionMatlab(TS snes,const char *func,mxArray *ctx) 2557 { 2558 PetscErrorCode ierr; 2559 TSMatlabContext *sctx; 2560 2561 PetscFunctionBegin; 2562 /* currently sctx is memory bleed */ 2563 ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr); 2564 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 2565 /* 2566 This should work, but it doesn't 2567 sctx->ctx = ctx; 2568 mexMakeArrayPersistent(sctx->ctx); 2569 */ 2570 sctx->ctx = mxDuplicateArray(ctx); 2571 ierr = TSSetIFunction(snes,TSComputeFunction_Matlab,sctx);CHKERRQ(ierr); 2572 PetscFunctionReturn(0); 2573 } 2574 2575 #undef __FUNCT__ 2576 #define __FUNCT__ "TSComputeJacobian_Matlab" 2577 /* 2578 TSComputeJacobian_Matlab - Calls the function that has been set with 2579 TSSetJacobianMatlab(). 2580 2581 Collective on TS 2582 2583 Input Parameters: 2584 + snes - the TS context 2585 . x - input vector 2586 . A, B - the matrices 2587 - ctx - user context 2588 2589 Output Parameter: 2590 . flag - structure of the matrix 2591 2592 Level: developer 2593 2594 .keywords: TS, nonlinear, compute, function 2595 2596 .seealso: TSSetFunction(), TSGetFunction() 2597 @*/ 2598 PetscErrorCode TSComputeJacobian_Matlab(TS snes,PetscReal time,Vec x,Vec xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flag, void *ctx) 2599 { 2600 PetscErrorCode ierr; 2601 TSMatlabContext *sctx = (TSMatlabContext *)ctx; 2602 int nlhs = 2,nrhs = 9; 2603 mxArray *plhs[2],*prhs[9]; 2604 long long int lx = 0,lxdot = 0,lA = 0,ls = 0, lB = 0; 2605 2606 PetscFunctionBegin; 2607 PetscValidHeaderSpecific(snes,TS_CLASSID,1); 2608 PetscValidHeaderSpecific(x,VEC_CLASSID,3); 2609 2610 /* call Matlab function in ctx with arguments x and y */ 2611 2612 ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 2613 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 2614 ierr = PetscMemcpy(&lxdot,&xdot,sizeof(x));CHKERRQ(ierr); 2615 ierr = PetscMemcpy(&lA,A,sizeof(x));CHKERRQ(ierr); 2616 ierr = PetscMemcpy(&lB,B,sizeof(x));CHKERRQ(ierr); 2617 prhs[0] = mxCreateDoubleScalar((double)ls); 2618 prhs[1] = mxCreateDoubleScalar((double)time); 2619 prhs[2] = mxCreateDoubleScalar((double)lx); 2620 prhs[3] = mxCreateDoubleScalar((double)lxdot); 2621 prhs[4] = mxCreateDoubleScalar((double)shift); 2622 prhs[5] = mxCreateDoubleScalar((double)lA); 2623 prhs[6] = mxCreateDoubleScalar((double)lB); 2624 prhs[7] = mxCreateString(sctx->funcname); 2625 prhs[8] = sctx->ctx; 2626 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeJacobianInternal");CHKERRQ(ierr); 2627 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 2628 *flag = (MatStructure) mxGetScalar(plhs[1]);CHKERRQ(ierr); 2629 mxDestroyArray(prhs[0]); 2630 mxDestroyArray(prhs[1]); 2631 mxDestroyArray(prhs[2]); 2632 mxDestroyArray(prhs[3]); 2633 mxDestroyArray(prhs[4]); 2634 mxDestroyArray(prhs[5]); 2635 mxDestroyArray(prhs[6]); 2636 mxDestroyArray(prhs[7]); 2637 mxDestroyArray(plhs[0]); 2638 mxDestroyArray(plhs[1]); 2639 PetscFunctionReturn(0); 2640 } 2641 2642 2643 #undef __FUNCT__ 2644 #define __FUNCT__ "TSSetJacobianMatlab" 2645 /* 2646 TSSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices 2647 vector for use by the TS routines in solving ODEs from MATLAB. Here the function is a string containing the name of a MATLAB function 2648 2649 Logically Collective on TS 2650 2651 Input Parameters: 2652 + snes - the TS context 2653 . A,B - Jacobian matrices 2654 . func - function evaluation routine 2655 - ctx - user context 2656 2657 Calling sequence of func: 2658 $ flag = func (TS snes,PetscReal time,Vec x,Vec xdot,Mat A,Mat B,void *ctx); 2659 2660 2661 Level: developer 2662 2663 .keywords: TS, nonlinear, set, function 2664 2665 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction() 2666 */ 2667 PetscErrorCode TSSetJacobianMatlab(TS snes,Mat A,Mat B,const char *func,mxArray *ctx) 2668 { 2669 PetscErrorCode ierr; 2670 TSMatlabContext *sctx; 2671 2672 PetscFunctionBegin; 2673 /* currently sctx is memory bleed */ 2674 ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr); 2675 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 2676 /* 2677 This should work, but it doesn't 2678 sctx->ctx = ctx; 2679 mexMakeArrayPersistent(sctx->ctx); 2680 */ 2681 sctx->ctx = mxDuplicateArray(ctx); 2682 ierr = TSSetIJacobian(snes,A,B,TSComputeJacobian_Matlab,sctx);CHKERRQ(ierr); 2683 PetscFunctionReturn(0); 2684 } 2685 2686 #undef __FUNCT__ 2687 #define __FUNCT__ "TSMonitor_Matlab" 2688 /* 2689 TSMonitor_Matlab - Calls the function that has been set with TSMonitorSetMatlab(). 2690 2691 Collective on TS 2692 2693 .seealso: TSSetFunction(), TSGetFunction() 2694 @*/ 2695 PetscErrorCode TSMonitor_Matlab(TS snes,PetscInt it, PetscReal time,Vec x, void *ctx) 2696 { 2697 PetscErrorCode ierr; 2698 TSMatlabContext *sctx = (TSMatlabContext *)ctx; 2699 int nlhs = 1,nrhs = 6; 2700 mxArray *plhs[1],*prhs[6]; 2701 long long int lx = 0,ls = 0; 2702 2703 PetscFunctionBegin; 2704 PetscValidHeaderSpecific(snes,TS_CLASSID,1); 2705 PetscValidHeaderSpecific(x,VEC_CLASSID,4); 2706 2707 ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 2708 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 2709 prhs[0] = mxCreateDoubleScalar((double)ls); 2710 prhs[1] = mxCreateDoubleScalar((double)it); 2711 prhs[2] = mxCreateDoubleScalar((double)time); 2712 prhs[3] = mxCreateDoubleScalar((double)lx); 2713 prhs[4] = mxCreateString(sctx->funcname); 2714 prhs[5] = sctx->ctx; 2715 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSMonitorInternal");CHKERRQ(ierr); 2716 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 2717 mxDestroyArray(prhs[0]); 2718 mxDestroyArray(prhs[1]); 2719 mxDestroyArray(prhs[2]); 2720 mxDestroyArray(prhs[3]); 2721 mxDestroyArray(prhs[4]); 2722 mxDestroyArray(plhs[0]); 2723 PetscFunctionReturn(0); 2724 } 2725 2726 2727 #undef __FUNCT__ 2728 #define __FUNCT__ "TSMonitorSetMatlab" 2729 /* 2730 TSMonitorSetMatlab - Sets the monitor function from Matlab 2731 2732 Level: developer 2733 2734 .keywords: TS, nonlinear, set, function 2735 2736 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction() 2737 */ 2738 PetscErrorCode TSMonitorSetMatlab(TS snes,const char *func,mxArray *ctx) 2739 { 2740 PetscErrorCode ierr; 2741 TSMatlabContext *sctx; 2742 2743 PetscFunctionBegin; 2744 /* currently sctx is memory bleed */ 2745 ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr); 2746 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 2747 /* 2748 This should work, but it doesn't 2749 sctx->ctx = ctx; 2750 mexMakeArrayPersistent(sctx->ctx); 2751 */ 2752 sctx->ctx = mxDuplicateArray(ctx); 2753 ierr = TSMonitorSet(snes,TSMonitor_Matlab,sctx,PETSC_NULL);CHKERRQ(ierr); 2754 PetscFunctionReturn(0); 2755 } 2756 2757 #endif 2758