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