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 PetscFunctionBegin; 578 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 579 if (Arhs){ 580 PetscValidHeaderSpecific(Arhs,MAT_CLASSID,2); 581 PetscCheckSameComm(ts,1,Arhs,2); 582 ts->Arhs = Arhs; 583 ts->ops->rhsmatrix = frhs; 584 } 585 if (Alhs){ 586 PetscValidHeaderSpecific(Alhs,MAT_CLASSID,4); 587 PetscCheckSameComm(ts,1,Alhs,4); 588 ts->Alhs = Alhs; 589 ts->ops->lhsmatrix = flhs; 590 } 591 592 ts->jacP = ctx; 593 ts->matflg = flag; 594 PetscFunctionReturn(0); 595 } 596 597 #undef __FUNCT__ 598 #define __FUNCT__ "TSGetMatrices" 599 /*@C 600 TSGetMatrices - Returns the matrices Arhs and Alhs at the present timestep, 601 where Alhs(t) U_t = Arhs(t) U. 602 603 Not Collective, but parallel objects are returned if TS is parallel 604 605 Input Parameter: 606 . ts - The TS context obtained from TSCreate() 607 608 Output Parameters: 609 + Arhs - The right-hand side matrix 610 . Alhs - The left-hand side matrix 611 - ctx - User-defined context for matrix evaluation routine 612 613 Notes: You can pass in PETSC_NULL for any return argument you do not need. 614 615 Level: intermediate 616 617 .seealso: TSSetMatrices(), TSGetTimeStep(), TSGetTime(), TSGetTimeStepNumber(), TSGetRHSJacobian() 618 619 .keywords: TS, timestep, get, matrix 620 621 @*/ 622 PetscErrorCode TSGetMatrices(TS ts,Mat *Arhs,Mat *Alhs,void **ctx) 623 { 624 PetscFunctionBegin; 625 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 626 if (Arhs) *Arhs = ts->Arhs; 627 if (Alhs) *Alhs = ts->Alhs; 628 if (ctx) *ctx = ts->jacP; 629 PetscFunctionReturn(0); 630 } 631 632 #undef __FUNCT__ 633 #define __FUNCT__ "TSSetRHSJacobian" 634 /*@C 635 TSSetRHSJacobian - Sets the function to compute the Jacobian of F, 636 where U_t = F(U,t), as well as the location to store the matrix. 637 Use TSSetMatrices() for linear problems. 638 639 Logically Collective on TS 640 641 Input Parameters: 642 + ts - the TS context obtained from TSCreate() 643 . A - Jacobian matrix 644 . B - preconditioner matrix (usually same as A) 645 . f - the Jacobian evaluation routine 646 - ctx - [optional] user-defined context for private data for the 647 Jacobian evaluation routine (may be PETSC_NULL) 648 649 Calling sequence of func: 650 $ func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx); 651 652 + t - current timestep 653 . u - input vector 654 . A - matrix A, where U_t = A(t)u 655 . B - preconditioner matrix, usually the same as A 656 . flag - flag indicating information about the preconditioner matrix 657 structure (same as flag in KSPSetOperators()) 658 - ctx - [optional] user-defined context for matrix evaluation routine 659 660 Notes: 661 See KSPSetOperators() for important information about setting the flag 662 output parameter in the routine func(). Be sure to read this information! 663 664 The routine func() takes Mat * as the matrix arguments rather than Mat. 665 This allows the matrix evaluation routine to replace A and/or B with a 666 completely new matrix structure (not just different matrix elements) 667 when appropriate, for instance, if the nonzero structure is changing 668 throughout the global iterations. 669 670 Level: beginner 671 672 .keywords: TS, timestep, set, right-hand-side, Jacobian 673 674 .seealso: TSDefaultComputeJacobianColor(), 675 SNESDefaultComputeJacobianColor(), TSSetRHSFunction(), TSSetMatrices() 676 677 @*/ 678 PetscErrorCode TSSetRHSJacobian(TS ts,Mat A,Mat B,PetscErrorCode (*f)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 679 { 680 PetscFunctionBegin; 681 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 682 PetscValidHeaderSpecific(A,MAT_CLASSID,2); 683 PetscValidHeaderSpecific(B,MAT_CLASSID,3); 684 PetscCheckSameComm(ts,1,A,2); 685 PetscCheckSameComm(ts,1,B,3); 686 if (ts->problem_type != TS_NONLINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Not for linear problems; use TSSetMatrices()"); 687 688 ts->ops->rhsjacobian = f; 689 ts->jacP = ctx; 690 ts->Arhs = A; 691 ts->B = B; 692 PetscFunctionReturn(0); 693 } 694 695 696 #undef __FUNCT__ 697 #define __FUNCT__ "TSSetIFunction" 698 /*@C 699 TSSetIFunction - Set the function to compute F(t,U,U_t) where F = 0 is the DAE to be solved. 700 701 Logically Collective on TS 702 703 Input Parameters: 704 + ts - the TS context obtained from TSCreate() 705 . f - the function evaluation routine 706 - ctx - user-defined context for private data for the function evaluation routine (may be PETSC_NULL) 707 708 Calling sequence of f: 709 $ f(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx); 710 711 + t - time at step/stage being solved 712 . u - state vector 713 . u_t - time derivative of state vector 714 . F - function vector 715 - ctx - [optional] user-defined context for matrix evaluation routine 716 717 Important: 718 The user MUST call either this routine, TSSetRHSFunction(), or TSSetMatrices(). This routine must be used when not solving an ODE. 719 720 Level: beginner 721 722 .keywords: TS, timestep, set, DAE, Jacobian 723 724 .seealso: TSSetMatrices(), TSSetRHSFunction(), TSSetIJacobian() 725 @*/ 726 PetscErrorCode TSSetIFunction(TS ts,TSIFunction f,void *ctx) 727 { 728 729 PetscFunctionBegin; 730 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 731 ts->ops->ifunction = f; 732 ts->funP = ctx; 733 PetscFunctionReturn(0); 734 } 735 736 #undef __FUNCT__ 737 #define __FUNCT__ "TSSetIJacobian" 738 /*@C 739 TSSetIJacobian - Set the function to compute the matrix dF/dU + a*dF/dU_t which is 740 the Jacobian of G(U) = F(t,U,W+a*U) where F(t,U,U_t) = 0 is the DAE to be solved. 741 The time integrator internally approximates U_t by W+a*U where the positive "shift" 742 a and vector W depend on the integration method, step size, and past states. 743 744 Logically Collective on TS 745 746 Input Parameters: 747 + ts - the TS context obtained from TSCreate() 748 . A - Jacobian matrix 749 . B - preconditioning matrix for A (may be same as A) 750 . f - the Jacobian evaluation routine 751 - ctx - user-defined context for private data for the Jacobian evaluation routine (may be PETSC_NULL) 752 753 Calling sequence of f: 754 $ f(TS ts,PetscReal t,Vec U,Vec U_t,PetscReal a,Mat *A,Mat *B,MatStructure *flag,void *ctx); 755 756 + t - time at step/stage being solved 757 . U - state vector 758 . U_t - time derivative of state vector 759 . a - shift 760 . A - Jacobian of G(U) = F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t 761 . B - preconditioning matrix for A, may be same as A 762 . flag - flag indicating information about the preconditioner matrix 763 structure (same as flag in KSPSetOperators()) 764 - ctx - [optional] user-defined context for matrix evaluation routine 765 766 Notes: 767 The matrices A and B are exactly the matrices that are used by SNES for the nonlinear solve. 768 769 Level: beginner 770 771 .keywords: TS, timestep, DAE, Jacobian 772 773 .seealso: TSSetIFunction(), TSSetRHSJacobian() 774 775 @*/ 776 PetscErrorCode TSSetIJacobian(TS ts,Mat A,Mat B,TSIJacobian f,void *ctx) 777 { 778 PetscErrorCode ierr; 779 780 PetscFunctionBegin; 781 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 782 if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2); 783 if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3); 784 if (A) PetscCheckSameComm(ts,1,A,2); 785 if (B) PetscCheckSameComm(ts,1,B,3); 786 if (f) ts->ops->ijacobian = f; 787 if (ctx) ts->jacP = ctx; 788 if (A) { 789 ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 790 if (ts->A) {ierr = MatDestroy(ts->A);CHKERRQ(ierr);} 791 ts->A = A; 792 } 793 #if 0 794 /* The sane and consistent alternative */ 795 if (B) { 796 ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 797 if (ts->B) {ierr = MatDestroy(ts->B);CHKERRQ(ierr);} 798 ts->B = B; 799 } 800 #else 801 /* Don't reference B because TSDestroy() doesn't destroy it. These ownership semantics are awkward and inconsistent. */ 802 if (B) ts->B = B; 803 #endif 804 PetscFunctionReturn(0); 805 } 806 807 #undef __FUNCT__ 808 #define __FUNCT__ "TSView" 809 /*@C 810 TSView - Prints the TS data structure. 811 812 Collective on TS 813 814 Input Parameters: 815 + ts - the TS context obtained from TSCreate() 816 - viewer - visualization context 817 818 Options Database Key: 819 . -ts_view - calls TSView() at end of TSStep() 820 821 Notes: 822 The available visualization contexts include 823 + PETSC_VIEWER_STDOUT_SELF - standard output (default) 824 - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 825 output where only the first processor opens 826 the file. All other processors send their 827 data to the first processor to print. 828 829 The user can open an alternative visualization context with 830 PetscViewerASCIIOpen() - output to a specified file. 831 832 Level: beginner 833 834 .keywords: TS, timestep, view 835 836 .seealso: PetscViewerASCIIOpen() 837 @*/ 838 PetscErrorCode TSView(TS ts,PetscViewer viewer) 839 { 840 PetscErrorCode ierr; 841 const TSType type; 842 PetscBool iascii,isstring; 843 844 PetscFunctionBegin; 845 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 846 if (!viewer) { 847 ierr = PetscViewerASCIIGetStdout(((PetscObject)ts)->comm,&viewer);CHKERRQ(ierr); 848 } 849 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 850 PetscCheckSameComm(ts,1,viewer,2); 851 852 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 853 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 854 if (iascii) { 855 ierr = PetscObjectPrintClassNamePrefixType((PetscObject)ts,viewer,"TS Object");CHKERRQ(ierr); 856 if (ts->ops->view) { 857 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 858 ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr); 859 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 860 } 861 ierr = PetscViewerASCIIPrintf(viewer," maximum steps=%D\n",ts->max_steps);CHKERRQ(ierr); 862 ierr = PetscViewerASCIIPrintf(viewer," maximum time=%G\n",ts->max_time);CHKERRQ(ierr); 863 if (ts->problem_type == TS_NONLINEAR) { 864 ierr = PetscViewerASCIIPrintf(viewer," total number of nonlinear solver iterations=%D\n",ts->nonlinear_its);CHKERRQ(ierr); 865 } 866 ierr = PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%D\n",ts->linear_its);CHKERRQ(ierr); 867 } else if (isstring) { 868 ierr = TSGetType(ts,&type);CHKERRQ(ierr); 869 ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr); 870 } 871 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 872 if (ts->ksp) {ierr = KSPView(ts->ksp,viewer);CHKERRQ(ierr);} 873 if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);} 874 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 875 PetscFunctionReturn(0); 876 } 877 878 879 #undef __FUNCT__ 880 #define __FUNCT__ "TSSetApplicationContext" 881 /*@C 882 TSSetApplicationContext - Sets an optional user-defined context for 883 the timesteppers. 884 885 Logically Collective on TS 886 887 Input Parameters: 888 + ts - the TS context obtained from TSCreate() 889 - usrP - optional user context 890 891 Level: intermediate 892 893 .keywords: TS, timestep, set, application, context 894 895 .seealso: TSGetApplicationContext() 896 @*/ 897 PetscErrorCode TSSetApplicationContext(TS ts,void *usrP) 898 { 899 PetscFunctionBegin; 900 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 901 ts->user = usrP; 902 PetscFunctionReturn(0); 903 } 904 905 #undef __FUNCT__ 906 #define __FUNCT__ "TSGetApplicationContext" 907 /*@C 908 TSGetApplicationContext - Gets the user-defined context for the 909 timestepper. 910 911 Not Collective 912 913 Input Parameter: 914 . ts - the TS context obtained from TSCreate() 915 916 Output Parameter: 917 . usrP - user context 918 919 Level: intermediate 920 921 .keywords: TS, timestep, get, application, context 922 923 .seealso: TSSetApplicationContext() 924 @*/ 925 PetscErrorCode TSGetApplicationContext(TS ts,void **usrP) 926 { 927 PetscFunctionBegin; 928 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 929 *usrP = ts->user; 930 PetscFunctionReturn(0); 931 } 932 933 #undef __FUNCT__ 934 #define __FUNCT__ "TSGetTimeStepNumber" 935 /*@ 936 TSGetTimeStepNumber - Gets the current number of timesteps. 937 938 Not Collective 939 940 Input Parameter: 941 . ts - the TS context obtained from TSCreate() 942 943 Output Parameter: 944 . iter - number steps so far 945 946 Level: intermediate 947 948 .keywords: TS, timestep, get, iteration, number 949 @*/ 950 PetscErrorCode TSGetTimeStepNumber(TS ts,PetscInt* iter) 951 { 952 PetscFunctionBegin; 953 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 954 PetscValidIntPointer(iter,2); 955 *iter = ts->steps; 956 PetscFunctionReturn(0); 957 } 958 959 #undef __FUNCT__ 960 #define __FUNCT__ "TSSetInitialTimeStep" 961 /*@ 962 TSSetInitialTimeStep - Sets the initial timestep to be used, 963 as well as the initial time. 964 965 Logically Collective on TS 966 967 Input Parameters: 968 + ts - the TS context obtained from TSCreate() 969 . initial_time - the initial time 970 - time_step - the size of the timestep 971 972 Level: intermediate 973 974 .seealso: TSSetTimeStep(), TSGetTimeStep() 975 976 .keywords: TS, set, initial, timestep 977 @*/ 978 PetscErrorCode TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step) 979 { 980 PetscFunctionBegin; 981 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 982 ts->time_step = time_step; 983 ts->initial_time_step = time_step; 984 ts->ptime = initial_time; 985 PetscFunctionReturn(0); 986 } 987 988 #undef __FUNCT__ 989 #define __FUNCT__ "TSSetTimeStep" 990 /*@ 991 TSSetTimeStep - Allows one to reset the timestep at any time, 992 useful for simple pseudo-timestepping codes. 993 994 Logically Collective on TS 995 996 Input Parameters: 997 + ts - the TS context obtained from TSCreate() 998 - time_step - the size of the timestep 999 1000 Level: intermediate 1001 1002 .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 1003 1004 .keywords: TS, set, timestep 1005 @*/ 1006 PetscErrorCode TSSetTimeStep(TS ts,PetscReal time_step) 1007 { 1008 PetscFunctionBegin; 1009 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1010 PetscValidLogicalCollectiveReal(ts,time_step,2); 1011 ts->time_step = time_step; 1012 PetscFunctionReturn(0); 1013 } 1014 1015 #undef __FUNCT__ 1016 #define __FUNCT__ "TSGetTimeStep" 1017 /*@ 1018 TSGetTimeStep - Gets the current timestep size. 1019 1020 Not Collective 1021 1022 Input Parameter: 1023 . ts - the TS context obtained from TSCreate() 1024 1025 Output Parameter: 1026 . dt - the current timestep size 1027 1028 Level: intermediate 1029 1030 .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 1031 1032 .keywords: TS, get, timestep 1033 @*/ 1034 PetscErrorCode TSGetTimeStep(TS ts,PetscReal* dt) 1035 { 1036 PetscFunctionBegin; 1037 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1038 PetscValidDoublePointer(dt,2); 1039 *dt = ts->time_step; 1040 PetscFunctionReturn(0); 1041 } 1042 1043 #undef __FUNCT__ 1044 #define __FUNCT__ "TSGetSolution" 1045 /*@ 1046 TSGetSolution - Returns the solution at the present timestep. It 1047 is valid to call this routine inside the function that you are evaluating 1048 in order to move to the new timestep. This vector not changed until 1049 the solution at the next timestep has been calculated. 1050 1051 Not Collective, but Vec returned is parallel if TS is parallel 1052 1053 Input Parameter: 1054 . ts - the TS context obtained from TSCreate() 1055 1056 Output Parameter: 1057 . v - the vector containing the solution 1058 1059 Level: intermediate 1060 1061 .seealso: TSGetTimeStep() 1062 1063 .keywords: TS, timestep, get, solution 1064 @*/ 1065 PetscErrorCode TSGetSolution(TS ts,Vec *v) 1066 { 1067 PetscFunctionBegin; 1068 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1069 PetscValidPointer(v,2); 1070 *v = ts->vec_sol; 1071 PetscFunctionReturn(0); 1072 } 1073 1074 /* ----- Routines to initialize and destroy a timestepper ---- */ 1075 #undef __FUNCT__ 1076 #define __FUNCT__ "TSSetProblemType" 1077 /*@ 1078 TSSetProblemType - Sets the type of problem to be solved. 1079 1080 Not collective 1081 1082 Input Parameters: 1083 + ts - The TS 1084 - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms 1085 .vb 1086 U_t = A U 1087 U_t = A(t) U 1088 U_t = F(t,U) 1089 .ve 1090 1091 Level: beginner 1092 1093 .keywords: TS, problem type 1094 .seealso: TSSetUp(), TSProblemType, TS 1095 @*/ 1096 PetscErrorCode TSSetProblemType(TS ts, TSProblemType type) 1097 { 1098 PetscFunctionBegin; 1099 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1100 ts->problem_type = type; 1101 PetscFunctionReturn(0); 1102 } 1103 1104 #undef __FUNCT__ 1105 #define __FUNCT__ "TSGetProblemType" 1106 /*@C 1107 TSGetProblemType - Gets the type of problem to be solved. 1108 1109 Not collective 1110 1111 Input Parameter: 1112 . ts - The TS 1113 1114 Output Parameter: 1115 . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms 1116 .vb 1117 U_t = A U 1118 U_t = A(t) U 1119 U_t = F(t,U) 1120 .ve 1121 1122 Level: beginner 1123 1124 .keywords: TS, problem type 1125 .seealso: TSSetUp(), TSProblemType, TS 1126 @*/ 1127 PetscErrorCode TSGetProblemType(TS ts, TSProblemType *type) 1128 { 1129 PetscFunctionBegin; 1130 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1131 PetscValidIntPointer(type,2); 1132 *type = ts->problem_type; 1133 PetscFunctionReturn(0); 1134 } 1135 1136 #undef __FUNCT__ 1137 #define __FUNCT__ "TSSetUp" 1138 /*@ 1139 TSSetUp - Sets up the internal data structures for the later use 1140 of a timestepper. 1141 1142 Collective on TS 1143 1144 Input Parameter: 1145 . ts - the TS context obtained from TSCreate() 1146 1147 Notes: 1148 For basic use of the TS solvers the user need not explicitly call 1149 TSSetUp(), since these actions will automatically occur during 1150 the call to TSStep(). However, if one wishes to control this 1151 phase separately, TSSetUp() should be called after TSCreate() 1152 and optional routines of the form TSSetXXX(), but before TSStep(). 1153 1154 Level: advanced 1155 1156 .keywords: TS, timestep, setup 1157 1158 .seealso: TSCreate(), TSStep(), TSDestroy() 1159 @*/ 1160 PetscErrorCode TSSetUp(TS ts) 1161 { 1162 PetscErrorCode ierr; 1163 1164 PetscFunctionBegin; 1165 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1166 if (!ts->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first"); 1167 if (!((PetscObject)ts)->type_name) { 1168 ierr = TSSetType(ts,TSEULER);CHKERRQ(ierr); 1169 } 1170 ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr); 1171 ts->setupcalled = 1; 1172 PetscFunctionReturn(0); 1173 } 1174 1175 #undef __FUNCT__ 1176 #define __FUNCT__ "TSDestroy" 1177 /*@ 1178 TSDestroy - Destroys the timestepper context that was created 1179 with TSCreate(). 1180 1181 Collective on TS 1182 1183 Input Parameter: 1184 . ts - the TS context obtained from TSCreate() 1185 1186 Level: beginner 1187 1188 .keywords: TS, timestepper, destroy 1189 1190 .seealso: TSCreate(), TSSetUp(), TSSolve() 1191 @*/ 1192 PetscErrorCode TSDestroy(TS ts) 1193 { 1194 PetscErrorCode ierr; 1195 1196 PetscFunctionBegin; 1197 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1198 if (--((PetscObject)ts)->refct > 0) PetscFunctionReturn(0); 1199 1200 /* if memory was published with AMS then destroy it */ 1201 ierr = PetscObjectDepublish(ts);CHKERRQ(ierr); 1202 if (ts->ops->destroy) {ierr = (*(ts)->ops->destroy)(ts);CHKERRQ(ierr);} 1203 1204 if (ts->vec_sol) {ierr = VecDestroy(ts->vec_sol);CHKERRQ(ierr);} 1205 if (ts->dm) {ierr = DMDestroy(ts->dm);CHKERRQ(ierr);} 1206 if (ts->A) {ierr = MatDestroy(ts->A);CHKERRQ(ierr);} 1207 if (ts->ksp) {ierr = KSPDestroy(ts->ksp);CHKERRQ(ierr);} 1208 if (ts->snes) {ierr = SNESDestroy(ts->snes);CHKERRQ(ierr);} 1209 ierr = TSMonitorCancel(ts);CHKERRQ(ierr); 1210 1211 ierr = PetscHeaderDestroy(ts);CHKERRQ(ierr); 1212 PetscFunctionReturn(0); 1213 } 1214 1215 #undef __FUNCT__ 1216 #define __FUNCT__ "TSGetSNES" 1217 /*@ 1218 TSGetSNES - Returns the SNES (nonlinear solver) associated with 1219 a TS (timestepper) context. Valid only for nonlinear problems. 1220 1221 Not Collective, but SNES is parallel if TS is parallel 1222 1223 Input Parameter: 1224 . ts - the TS context obtained from TSCreate() 1225 1226 Output Parameter: 1227 . snes - the nonlinear solver context 1228 1229 Notes: 1230 The user can then directly manipulate the SNES context to set various 1231 options, etc. Likewise, the user can then extract and manipulate the 1232 KSP, KSP, and PC contexts as well. 1233 1234 TSGetSNES() does not work for integrators that do not use SNES; in 1235 this case TSGetSNES() returns PETSC_NULL in snes. 1236 1237 Level: beginner 1238 1239 .keywords: timestep, get, SNES 1240 @*/ 1241 PetscErrorCode TSGetSNES(TS ts,SNES *snes) 1242 { 1243 PetscFunctionBegin; 1244 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1245 PetscValidPointer(snes,2); 1246 if (!((PetscObject)ts)->type_name) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"SNES is not created yet. Call TSSetType() first"); 1247 if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Nonlinear only; use TSGetKSP()"); 1248 *snes = ts->snes; 1249 PetscFunctionReturn(0); 1250 } 1251 1252 #undef __FUNCT__ 1253 #define __FUNCT__ "TSGetKSP" 1254 /*@ 1255 TSGetKSP - Returns the KSP (linear solver) associated with 1256 a TS (timestepper) context. 1257 1258 Not Collective, but KSP is parallel if TS is parallel 1259 1260 Input Parameter: 1261 . ts - the TS context obtained from TSCreate() 1262 1263 Output Parameter: 1264 . ksp - the nonlinear solver context 1265 1266 Notes: 1267 The user can then directly manipulate the KSP context to set various 1268 options, etc. Likewise, the user can then extract and manipulate the 1269 KSP and PC contexts as well. 1270 1271 TSGetKSP() does not work for integrators that do not use KSP; 1272 in this case TSGetKSP() returns PETSC_NULL in ksp. 1273 1274 Level: beginner 1275 1276 .keywords: timestep, get, KSP 1277 @*/ 1278 PetscErrorCode TSGetKSP(TS ts,KSP *ksp) 1279 { 1280 PetscFunctionBegin; 1281 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1282 PetscValidPointer(ksp,2); 1283 if (!((PetscObject)ts)->type_name) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"KSP is not created yet. Call TSSetType() first"); 1284 if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()"); 1285 *ksp = ts->ksp; 1286 PetscFunctionReturn(0); 1287 } 1288 1289 /* ----------- Routines to set solver parameters ---------- */ 1290 1291 #undef __FUNCT__ 1292 #define __FUNCT__ "TSGetDuration" 1293 /*@ 1294 TSGetDuration - Gets the maximum number of timesteps to use and 1295 maximum time for iteration. 1296 1297 Not Collective 1298 1299 Input Parameters: 1300 + ts - the TS context obtained from TSCreate() 1301 . maxsteps - maximum number of iterations to use, or PETSC_NULL 1302 - maxtime - final time to iterate to, or PETSC_NULL 1303 1304 Level: intermediate 1305 1306 .keywords: TS, timestep, get, maximum, iterations, time 1307 @*/ 1308 PetscErrorCode TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime) 1309 { 1310 PetscFunctionBegin; 1311 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1312 if (maxsteps) { 1313 PetscValidIntPointer(maxsteps,2); 1314 *maxsteps = ts->max_steps; 1315 } 1316 if (maxtime ) { 1317 PetscValidScalarPointer(maxtime,3); 1318 *maxtime = ts->max_time; 1319 } 1320 PetscFunctionReturn(0); 1321 } 1322 1323 #undef __FUNCT__ 1324 #define __FUNCT__ "TSSetDuration" 1325 /*@ 1326 TSSetDuration - Sets the maximum number of timesteps to use and 1327 maximum time for iteration. 1328 1329 Logically Collective on TS 1330 1331 Input Parameters: 1332 + ts - the TS context obtained from TSCreate() 1333 . maxsteps - maximum number of iterations to use 1334 - maxtime - final time to iterate to 1335 1336 Options Database Keys: 1337 . -ts_max_steps <maxsteps> - Sets maxsteps 1338 . -ts_max_time <maxtime> - Sets maxtime 1339 1340 Notes: 1341 The default maximum number of iterations is 5000. Default time is 5.0 1342 1343 Level: intermediate 1344 1345 .keywords: TS, timestep, set, maximum, iterations 1346 @*/ 1347 PetscErrorCode TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime) 1348 { 1349 PetscFunctionBegin; 1350 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1351 PetscValidLogicalCollectiveInt(ts,maxsteps,2); 1352 PetscValidLogicalCollectiveReal(ts,maxtime,2); 1353 ts->max_steps = maxsteps; 1354 ts->max_time = maxtime; 1355 PetscFunctionReturn(0); 1356 } 1357 1358 #undef __FUNCT__ 1359 #define __FUNCT__ "TSSetSolution" 1360 /*@ 1361 TSSetSolution - Sets the initial solution vector 1362 for use by the TS routines. 1363 1364 Logically Collective on TS and Vec 1365 1366 Input Parameters: 1367 + ts - the TS context obtained from TSCreate() 1368 - x - the solution vector 1369 1370 Level: beginner 1371 1372 .keywords: TS, timestep, set, solution, initial conditions 1373 @*/ 1374 PetscErrorCode TSSetSolution(TS ts,Vec x) 1375 { 1376 PetscErrorCode ierr; 1377 1378 PetscFunctionBegin; 1379 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1380 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 1381 ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr); 1382 if (ts->vec_sol) {ierr = VecDestroy(ts->vec_sol);CHKERRQ(ierr);} 1383 ts->vec_sol = x; 1384 PetscFunctionReturn(0); 1385 } 1386 1387 #undef __FUNCT__ 1388 #define __FUNCT__ "TSSetPreStep" 1389 /*@C 1390 TSSetPreStep - Sets the general-purpose function 1391 called once at the beginning of each time step. 1392 1393 Logically Collective on TS 1394 1395 Input Parameters: 1396 + ts - The TS context obtained from TSCreate() 1397 - func - The function 1398 1399 Calling sequence of func: 1400 . func (TS ts); 1401 1402 Level: intermediate 1403 1404 .keywords: TS, timestep 1405 @*/ 1406 PetscErrorCode TSSetPreStep(TS ts, PetscErrorCode (*func)(TS)) 1407 { 1408 PetscFunctionBegin; 1409 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1410 ts->ops->prestep = func; 1411 PetscFunctionReturn(0); 1412 } 1413 1414 #undef __FUNCT__ 1415 #define __FUNCT__ "TSPreStep" 1416 /*@C 1417 TSPreStep - Runs the user-defined pre-step function. 1418 1419 Collective on TS 1420 1421 Input Parameters: 1422 . ts - The TS context obtained from TSCreate() 1423 1424 Notes: 1425 TSPreStep() is typically used within time stepping implementations, 1426 so most users would not generally call this routine themselves. 1427 1428 Level: developer 1429 1430 .keywords: TS, timestep 1431 @*/ 1432 PetscErrorCode TSPreStep(TS ts) 1433 { 1434 PetscErrorCode ierr; 1435 1436 PetscFunctionBegin; 1437 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1438 if (ts->ops->prestep) { 1439 PetscStackPush("TS PreStep function"); 1440 ierr = (*ts->ops->prestep)(ts);CHKERRQ(ierr); 1441 PetscStackPop; 1442 } 1443 PetscFunctionReturn(0); 1444 } 1445 1446 #undef __FUNCT__ 1447 #define __FUNCT__ "TSSetPostStep" 1448 /*@C 1449 TSSetPostStep - Sets the general-purpose function 1450 called once at the end of each time step. 1451 1452 Logically Collective on TS 1453 1454 Input Parameters: 1455 + ts - The TS context obtained from TSCreate() 1456 - func - The function 1457 1458 Calling sequence of func: 1459 . func (TS ts); 1460 1461 Level: intermediate 1462 1463 .keywords: TS, timestep 1464 @*/ 1465 PetscErrorCode TSSetPostStep(TS ts, PetscErrorCode (*func)(TS)) 1466 { 1467 PetscFunctionBegin; 1468 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1469 ts->ops->poststep = func; 1470 PetscFunctionReturn(0); 1471 } 1472 1473 #undef __FUNCT__ 1474 #define __FUNCT__ "TSPostStep" 1475 /*@C 1476 TSPostStep - Runs the user-defined post-step function. 1477 1478 Collective on TS 1479 1480 Input Parameters: 1481 . ts - The TS context obtained from TSCreate() 1482 1483 Notes: 1484 TSPostStep() is typically used within time stepping implementations, 1485 so most users would not generally call this routine themselves. 1486 1487 Level: developer 1488 1489 .keywords: TS, timestep 1490 @*/ 1491 PetscErrorCode TSPostStep(TS ts) 1492 { 1493 PetscErrorCode ierr; 1494 1495 PetscFunctionBegin; 1496 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1497 if (ts->ops->poststep) { 1498 PetscStackPush("TS PostStep function"); 1499 ierr = (*ts->ops->poststep)(ts);CHKERRQ(ierr); 1500 PetscStackPop; 1501 } 1502 PetscFunctionReturn(0); 1503 } 1504 1505 /* ------------ Routines to set performance monitoring options ----------- */ 1506 1507 #undef __FUNCT__ 1508 #define __FUNCT__ "TSMonitorSet" 1509 /*@C 1510 TSMonitorSet - Sets an ADDITIONAL function that is to be used at every 1511 timestep to display the iteration's progress. 1512 1513 Logically Collective on TS 1514 1515 Input Parameters: 1516 + ts - the TS context obtained from TSCreate() 1517 . func - monitoring routine 1518 . mctx - [optional] user-defined context for private data for the 1519 monitor routine (use PETSC_NULL if no context is desired) 1520 - monitordestroy - [optional] routine that frees monitor context 1521 (may be PETSC_NULL) 1522 1523 Calling sequence of func: 1524 $ int func(TS ts,PetscInt steps,PetscReal time,Vec x,void *mctx) 1525 1526 + ts - the TS context 1527 . steps - iteration number 1528 . time - current time 1529 . x - current iterate 1530 - mctx - [optional] monitoring context 1531 1532 Notes: 1533 This routine adds an additional monitor to the list of monitors that 1534 already has been loaded. 1535 1536 Fortran notes: Only a single monitor function can be set for each TS object 1537 1538 Level: intermediate 1539 1540 .keywords: TS, timestep, set, monitor 1541 1542 .seealso: TSMonitorDefault(), TSMonitorCancel() 1543 @*/ 1544 PetscErrorCode TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void*)) 1545 { 1546 PetscFunctionBegin; 1547 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1548 if (ts->numbermonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1549 ts->monitor[ts->numbermonitors] = monitor; 1550 ts->mdestroy[ts->numbermonitors] = mdestroy; 1551 ts->monitorcontext[ts->numbermonitors++] = (void*)mctx; 1552 PetscFunctionReturn(0); 1553 } 1554 1555 #undef __FUNCT__ 1556 #define __FUNCT__ "TSMonitorCancel" 1557 /*@C 1558 TSMonitorCancel - Clears all the monitors that have been set on a time-step object. 1559 1560 Logically Collective on TS 1561 1562 Input Parameters: 1563 . ts - the TS context obtained from TSCreate() 1564 1565 Notes: 1566 There is no way to remove a single, specific monitor. 1567 1568 Level: intermediate 1569 1570 .keywords: TS, timestep, set, monitor 1571 1572 .seealso: TSMonitorDefault(), TSMonitorSet() 1573 @*/ 1574 PetscErrorCode TSMonitorCancel(TS ts) 1575 { 1576 PetscErrorCode ierr; 1577 PetscInt i; 1578 1579 PetscFunctionBegin; 1580 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1581 for (i=0; i<ts->numbermonitors; i++) { 1582 if (ts->mdestroy[i]) { 1583 ierr = (*ts->mdestroy[i])(ts->monitorcontext[i]);CHKERRQ(ierr); 1584 } 1585 } 1586 ts->numbermonitors = 0; 1587 PetscFunctionReturn(0); 1588 } 1589 1590 #undef __FUNCT__ 1591 #define __FUNCT__ "TSMonitorDefault" 1592 /*@ 1593 TSMonitorDefault - Sets the Default monitor 1594 1595 Level: intermediate 1596 1597 .keywords: TS, set, monitor 1598 1599 .seealso: TSMonitorDefault(), TSMonitorSet() 1600 @*/ 1601 PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *ctx) 1602 { 1603 PetscErrorCode ierr; 1604 PetscViewerASCIIMonitor viewer = (PetscViewerASCIIMonitor)ctx; 1605 1606 PetscFunctionBegin; 1607 if (!ctx) { 1608 ierr = PetscViewerASCIIMonitorCreate(((PetscObject)ts)->comm,"stdout",0,&viewer);CHKERRQ(ierr); 1609 } 1610 ierr = PetscViewerASCIIMonitorPrintf(viewer,"%D TS dt %G time %G\n",step,ts->time_step,ptime);CHKERRQ(ierr); 1611 if (!ctx) { 1612 ierr = PetscViewerASCIIMonitorDestroy(viewer);CHKERRQ(ierr); 1613 } 1614 PetscFunctionReturn(0); 1615 } 1616 1617 #undef __FUNCT__ 1618 #define __FUNCT__ "TSStep" 1619 /*@ 1620 TSStep - Steps the requested number of timesteps. 1621 1622 Collective on TS 1623 1624 Input Parameter: 1625 . ts - the TS context obtained from TSCreate() 1626 1627 Output Parameters: 1628 + steps - number of iterations until termination 1629 - ptime - time until termination 1630 1631 Level: beginner 1632 1633 .keywords: TS, timestep, solve 1634 1635 .seealso: TSCreate(), TSSetUp(), TSDestroy() 1636 @*/ 1637 PetscErrorCode TSStep(TS ts,PetscInt *steps,PetscReal *ptime) 1638 { 1639 PetscErrorCode ierr; 1640 1641 PetscFunctionBegin; 1642 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1643 if (!ts->setupcalled) { 1644 ierr = TSSetUp(ts);CHKERRQ(ierr); 1645 } 1646 1647 ierr = PetscLogEventBegin(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr); 1648 ierr = (*ts->ops->step)(ts, steps, ptime);CHKERRQ(ierr); 1649 ierr = PetscLogEventEnd(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr); 1650 1651 if (!PetscPreLoadingOn) { 1652 ierr = TSViewFromOptions(ts,((PetscObject)ts)->name);CHKERRQ(ierr); 1653 } 1654 PetscFunctionReturn(0); 1655 } 1656 1657 #undef __FUNCT__ 1658 #define __FUNCT__ "TSSolve" 1659 /*@ 1660 TSSolve - Steps the requested number of timesteps. 1661 1662 Collective on TS 1663 1664 Input Parameter: 1665 + ts - the TS context obtained from TSCreate() 1666 - x - the solution vector, or PETSC_NULL if it was set with TSSetSolution() 1667 1668 Level: beginner 1669 1670 .keywords: TS, timestep, solve 1671 1672 .seealso: TSCreate(), TSSetSolution(), TSStep() 1673 @*/ 1674 PetscErrorCode TSSolve(TS ts, Vec x) 1675 { 1676 PetscInt steps; 1677 PetscReal ptime; 1678 PetscErrorCode ierr; 1679 1680 PetscFunctionBegin; 1681 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1682 /* set solution vector if provided */ 1683 if (x) { ierr = TSSetSolution(ts, x); CHKERRQ(ierr); } 1684 /* reset time step and iteration counters */ 1685 ts->steps = 0; ts->linear_its = 0; ts->nonlinear_its = 0; 1686 /* steps the requested number of timesteps. */ 1687 ierr = TSStep(ts, &steps, &ptime);CHKERRQ(ierr); 1688 PetscFunctionReturn(0); 1689 } 1690 1691 #undef __FUNCT__ 1692 #define __FUNCT__ "TSMonitor" 1693 /* 1694 Runs the user provided monitor routines, if they exists. 1695 */ 1696 PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x) 1697 { 1698 PetscErrorCode ierr; 1699 PetscInt i,n = ts->numbermonitors; 1700 1701 PetscFunctionBegin; 1702 for (i=0; i<n; i++) { 1703 ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr); 1704 } 1705 PetscFunctionReturn(0); 1706 } 1707 1708 /* ------------------------------------------------------------------------*/ 1709 1710 #undef __FUNCT__ 1711 #define __FUNCT__ "TSMonitorLGCreate" 1712 /*@C 1713 TSMonitorLGCreate - Creates a line graph context for use with 1714 TS to monitor convergence of preconditioned residual norms. 1715 1716 Collective on TS 1717 1718 Input Parameters: 1719 + host - the X display to open, or null for the local machine 1720 . label - the title to put in the title bar 1721 . x, y - the screen coordinates of the upper left coordinate of the window 1722 - m, n - the screen width and height in pixels 1723 1724 Output Parameter: 1725 . draw - the drawing context 1726 1727 Options Database Key: 1728 . -ts_monitor_draw - automatically sets line graph monitor 1729 1730 Notes: 1731 Use TSMonitorLGDestroy() to destroy this line graph, not PetscDrawLGDestroy(). 1732 1733 Level: intermediate 1734 1735 .keywords: TS, monitor, line graph, residual, seealso 1736 1737 .seealso: TSMonitorLGDestroy(), TSMonitorSet() 1738 1739 @*/ 1740 PetscErrorCode TSMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1741 { 1742 PetscDraw win; 1743 PetscErrorCode ierr; 1744 1745 PetscFunctionBegin; 1746 ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr); 1747 ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr); 1748 ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr); 1749 ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr); 1750 1751 ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr); 1752 PetscFunctionReturn(0); 1753 } 1754 1755 #undef __FUNCT__ 1756 #define __FUNCT__ "TSMonitorLG" 1757 PetscErrorCode TSMonitorLG(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx) 1758 { 1759 PetscDrawLG lg = (PetscDrawLG) monctx; 1760 PetscReal x,y = ptime; 1761 PetscErrorCode ierr; 1762 1763 PetscFunctionBegin; 1764 if (!monctx) { 1765 MPI_Comm comm; 1766 PetscViewer viewer; 1767 1768 ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr); 1769 viewer = PETSC_VIEWER_DRAW_(comm); 1770 ierr = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr); 1771 } 1772 1773 if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 1774 x = (PetscReal)n; 1775 ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 1776 if (n < 20 || (n % 5)) { 1777 ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 1778 } 1779 PetscFunctionReturn(0); 1780 } 1781 1782 #undef __FUNCT__ 1783 #define __FUNCT__ "TSMonitorLGDestroy" 1784 /*@C 1785 TSMonitorLGDestroy - Destroys a line graph context that was created 1786 with TSMonitorLGCreate(). 1787 1788 Collective on PetscDrawLG 1789 1790 Input Parameter: 1791 . draw - the drawing context 1792 1793 Level: intermediate 1794 1795 .keywords: TS, monitor, line graph, destroy 1796 1797 .seealso: TSMonitorLGCreate(), TSMonitorSet(), TSMonitorLG(); 1798 @*/ 1799 PetscErrorCode TSMonitorLGDestroy(PetscDrawLG drawlg) 1800 { 1801 PetscDraw draw; 1802 PetscErrorCode ierr; 1803 1804 PetscFunctionBegin; 1805 ierr = PetscDrawLGGetDraw(drawlg,&draw);CHKERRQ(ierr); 1806 ierr = PetscDrawDestroy(draw);CHKERRQ(ierr); 1807 ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr); 1808 PetscFunctionReturn(0); 1809 } 1810 1811 #undef __FUNCT__ 1812 #define __FUNCT__ "TSGetTime" 1813 /*@ 1814 TSGetTime - Gets the current time. 1815 1816 Not Collective 1817 1818 Input Parameter: 1819 . ts - the TS context obtained from TSCreate() 1820 1821 Output Parameter: 1822 . t - the current time 1823 1824 Level: beginner 1825 1826 .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 1827 1828 .keywords: TS, get, time 1829 @*/ 1830 PetscErrorCode TSGetTime(TS ts,PetscReal* t) 1831 { 1832 PetscFunctionBegin; 1833 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1834 PetscValidDoublePointer(t,2); 1835 *t = ts->ptime; 1836 PetscFunctionReturn(0); 1837 } 1838 1839 #undef __FUNCT__ 1840 #define __FUNCT__ "TSSetTime" 1841 /*@ 1842 TSSetTime - Allows one to reset the time. 1843 1844 Logically Collective on TS 1845 1846 Input Parameters: 1847 + ts - the TS context obtained from TSCreate() 1848 - time - the time 1849 1850 Level: intermediate 1851 1852 .seealso: TSGetTime(), TSSetDuration() 1853 1854 .keywords: TS, set, time 1855 @*/ 1856 PetscErrorCode TSSetTime(TS ts, PetscReal t) 1857 { 1858 PetscFunctionBegin; 1859 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1860 PetscValidLogicalCollectiveReal(ts,t,2); 1861 ts->ptime = t; 1862 PetscFunctionReturn(0); 1863 } 1864 1865 #undef __FUNCT__ 1866 #define __FUNCT__ "TSSetOptionsPrefix" 1867 /*@C 1868 TSSetOptionsPrefix - Sets the prefix used for searching for all 1869 TS options in the database. 1870 1871 Logically Collective on TS 1872 1873 Input Parameter: 1874 + ts - The TS context 1875 - prefix - The prefix to prepend to all option names 1876 1877 Notes: 1878 A hyphen (-) must NOT be given at the beginning of the prefix name. 1879 The first character of all runtime options is AUTOMATICALLY the 1880 hyphen. 1881 1882 Level: advanced 1883 1884 .keywords: TS, set, options, prefix, database 1885 1886 .seealso: TSSetFromOptions() 1887 1888 @*/ 1889 PetscErrorCode TSSetOptionsPrefix(TS ts,const char prefix[]) 1890 { 1891 PetscErrorCode ierr; 1892 1893 PetscFunctionBegin; 1894 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1895 ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 1896 switch(ts->problem_type) { 1897 case TS_NONLINEAR: 1898 if (ts->snes) { 1899 ierr = SNESSetOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr); 1900 } 1901 break; 1902 case TS_LINEAR: 1903 if (ts->ksp) { 1904 ierr = KSPSetOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr); 1905 } 1906 break; 1907 } 1908 PetscFunctionReturn(0); 1909 } 1910 1911 1912 #undef __FUNCT__ 1913 #define __FUNCT__ "TSAppendOptionsPrefix" 1914 /*@C 1915 TSAppendOptionsPrefix - Appends to the prefix used for searching for all 1916 TS options in the database. 1917 1918 Logically Collective on TS 1919 1920 Input Parameter: 1921 + ts - The TS context 1922 - prefix - The prefix to prepend to all option names 1923 1924 Notes: 1925 A hyphen (-) must NOT be given at the beginning of the prefix name. 1926 The first character of all runtime options is AUTOMATICALLY the 1927 hyphen. 1928 1929 Level: advanced 1930 1931 .keywords: TS, append, options, prefix, database 1932 1933 .seealso: TSGetOptionsPrefix() 1934 1935 @*/ 1936 PetscErrorCode TSAppendOptionsPrefix(TS ts,const char prefix[]) 1937 { 1938 PetscErrorCode ierr; 1939 1940 PetscFunctionBegin; 1941 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1942 ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 1943 switch(ts->problem_type) { 1944 case TS_NONLINEAR: 1945 if (ts->snes) { 1946 ierr = SNESAppendOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr); 1947 } 1948 break; 1949 case TS_LINEAR: 1950 if (ts->ksp) { 1951 ierr = KSPAppendOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr); 1952 } 1953 break; 1954 } 1955 PetscFunctionReturn(0); 1956 } 1957 1958 #undef __FUNCT__ 1959 #define __FUNCT__ "TSGetOptionsPrefix" 1960 /*@C 1961 TSGetOptionsPrefix - Sets the prefix used for searching for all 1962 TS options in the database. 1963 1964 Not Collective 1965 1966 Input Parameter: 1967 . ts - The TS context 1968 1969 Output Parameter: 1970 . prefix - A pointer to the prefix string used 1971 1972 Notes: On the fortran side, the user should pass in a string 'prifix' of 1973 sufficient length to hold the prefix. 1974 1975 Level: intermediate 1976 1977 .keywords: TS, get, options, prefix, database 1978 1979 .seealso: TSAppendOptionsPrefix() 1980 @*/ 1981 PetscErrorCode TSGetOptionsPrefix(TS ts,const char *prefix[]) 1982 { 1983 PetscErrorCode ierr; 1984 1985 PetscFunctionBegin; 1986 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1987 PetscValidPointer(prefix,2); 1988 ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 1989 PetscFunctionReturn(0); 1990 } 1991 1992 #undef __FUNCT__ 1993 #define __FUNCT__ "TSGetRHSJacobian" 1994 /*@C 1995 TSGetRHSJacobian - Returns the Jacobian J at the present timestep. 1996 1997 Not Collective, but parallel objects are returned if TS is parallel 1998 1999 Input Parameter: 2000 . ts - The TS context obtained from TSCreate() 2001 2002 Output Parameters: 2003 + J - The Jacobian J of F, where U_t = F(U,t) 2004 . M - The preconditioner matrix, usually the same as J 2005 - ctx - User-defined context for Jacobian evaluation routine 2006 2007 Notes: You can pass in PETSC_NULL for any return argument you do not need. 2008 2009 Level: intermediate 2010 2011 .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber() 2012 2013 .keywords: TS, timestep, get, matrix, Jacobian 2014 @*/ 2015 PetscErrorCode TSGetRHSJacobian(TS ts,Mat *J,Mat *M,void **ctx) 2016 { 2017 PetscFunctionBegin; 2018 if (J) *J = ts->Arhs; 2019 if (M) *M = ts->B; 2020 if (ctx) *ctx = ts->jacP; 2021 PetscFunctionReturn(0); 2022 } 2023 2024 #undef __FUNCT__ 2025 #define __FUNCT__ "TSGetIJacobian" 2026 /*@C 2027 TSGetIJacobian - Returns the implicit Jacobian at the present timestep. 2028 2029 Not Collective, but parallel objects are returned if TS is parallel 2030 2031 Input Parameter: 2032 . ts - The TS context obtained from TSCreate() 2033 2034 Output Parameters: 2035 + A - The Jacobian of F(t,U,U_t) 2036 . B - The preconditioner matrix, often the same as A 2037 . f - The function to compute the matrices 2038 - ctx - User-defined context for Jacobian evaluation routine 2039 2040 Notes: You can pass in PETSC_NULL for any return argument you do not need. 2041 2042 Level: advanced 2043 2044 .seealso: TSGetTimeStep(), TSGetRHSJacobian(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber() 2045 2046 .keywords: TS, timestep, get, matrix, Jacobian 2047 @*/ 2048 PetscErrorCode TSGetIJacobian(TS ts,Mat *A,Mat *B,TSIJacobian *f,void **ctx) 2049 { 2050 PetscFunctionBegin; 2051 if (A) *A = ts->A; 2052 if (B) *B = ts->B; 2053 if (f) *f = ts->ops->ijacobian; 2054 if (ctx) *ctx = ts->jacP; 2055 PetscFunctionReturn(0); 2056 } 2057 2058 #undef __FUNCT__ 2059 #define __FUNCT__ "TSMonitorSolution" 2060 /*@C 2061 TSMonitorSolution - Monitors progress of the TS solvers by calling 2062 VecView() for the solution at each timestep 2063 2064 Collective on TS 2065 2066 Input Parameters: 2067 + ts - the TS context 2068 . step - current time-step 2069 . ptime - current time 2070 - dummy - either a viewer or PETSC_NULL 2071 2072 Level: intermediate 2073 2074 .keywords: TS, vector, monitor, view 2075 2076 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView() 2077 @*/ 2078 PetscErrorCode TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy) 2079 { 2080 PetscErrorCode ierr; 2081 PetscViewer viewer = (PetscViewer) dummy; 2082 2083 PetscFunctionBegin; 2084 if (!dummy) { 2085 viewer = PETSC_VIEWER_DRAW_(((PetscObject)ts)->comm); 2086 } 2087 ierr = VecView(x,viewer);CHKERRQ(ierr); 2088 PetscFunctionReturn(0); 2089 } 2090 2091 2092 #undef __FUNCT__ 2093 #define __FUNCT__ "TSSetDM" 2094 /*@ 2095 TSSetDM - Sets the DM that may be used by some preconditioners 2096 2097 Logically Collective on TS and DM 2098 2099 Input Parameters: 2100 + ts - the preconditioner context 2101 - dm - the dm 2102 2103 Level: intermediate 2104 2105 2106 .seealso: TSGetDM(), SNESSetDM(), SNESGetDM() 2107 @*/ 2108 PetscErrorCode TSSetDM(TS ts,DM dm) 2109 { 2110 PetscErrorCode ierr; 2111 2112 PetscFunctionBegin; 2113 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2114 ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr); 2115 if (ts->dm) {ierr = DMDestroy(ts->dm);CHKERRQ(ierr);} 2116 ts->dm = dm; 2117 if (ts->snes) { 2118 ierr = SNESSetDM(ts->snes,dm);CHKERRQ(ierr); 2119 } 2120 if (ts->ksp) { 2121 ierr = KSPSetDM(ts->ksp,dm);CHKERRQ(ierr); 2122 ierr = KSPSetDMActive(ts->ksp,PETSC_FALSE);CHKERRQ(ierr); 2123 } 2124 PetscFunctionReturn(0); 2125 } 2126 2127 #undef __FUNCT__ 2128 #define __FUNCT__ "TSGetDM" 2129 /*@ 2130 TSGetDM - Gets the DM that may be used by some preconditioners 2131 2132 Not Collective 2133 2134 Input Parameter: 2135 . ts - the preconditioner context 2136 2137 Output Parameter: 2138 . dm - the dm 2139 2140 Level: intermediate 2141 2142 2143 .seealso: TSSetDM(), SNESSetDM(), SNESGetDM() 2144 @*/ 2145 PetscErrorCode TSGetDM(TS ts,DM *dm) 2146 { 2147 PetscFunctionBegin; 2148 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2149 *dm = ts->dm; 2150 PetscFunctionReturn(0); 2151 } 2152 2153 #undef __FUNCT__ 2154 #define __FUNCT__ "SNESTSFormFunction" 2155 /*@ 2156 SNESTSFormFunction - Function to evaluate nonlinear residual 2157 2158 Logically Collective on SNES 2159 2160 Input Parameter: 2161 + snes - nonlinear solver 2162 . X - the current state at which to evaluate the residual 2163 - ctx - user context, must be a TS 2164 2165 Output Parameter: 2166 . F - the nonlinear residual 2167 2168 Notes: 2169 This function is not normally called by users and is automatically registered with the SNES used by TS. 2170 It is most frequently passed to MatFDColoringSetFunction(). 2171 2172 Level: advanced 2173 2174 .seealso: SNESSetFunction(), MatFDColoringSetFunction() 2175 @*/ 2176 PetscErrorCode SNESTSFormFunction(SNES snes,Vec X,Vec F,void *ctx) 2177 { 2178 TS ts = (TS)ctx; 2179 PetscErrorCode ierr; 2180 2181 PetscFunctionBegin; 2182 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2183 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 2184 PetscValidHeaderSpecific(F,VEC_CLASSID,3); 2185 PetscValidHeaderSpecific(ts,TS_CLASSID,4); 2186 ierr = (ts->ops->snesfunction)(snes,X,F,ts);CHKERRQ(ierr); 2187 PetscFunctionReturn(0); 2188 } 2189 2190 #undef __FUNCT__ 2191 #define __FUNCT__ "SNESTSFormJacobian" 2192 /*@ 2193 SNESTSFormJacobian - Function to evaluate the Jacobian 2194 2195 Collective on SNES 2196 2197 Input Parameter: 2198 + snes - nonlinear solver 2199 . X - the current state at which to evaluate the residual 2200 - ctx - user context, must be a TS 2201 2202 Output Parameter: 2203 + A - the Jacobian 2204 . B - the preconditioning matrix (may be the same as A) 2205 - flag - indicates any structure change in the matrix 2206 2207 Notes: 2208 This function is not normally called by users and is automatically registered with the SNES used by TS. 2209 2210 Level: developer 2211 2212 .seealso: SNESSetJacobian() 2213 @*/ 2214 PetscErrorCode SNESTSFormJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flag,void *ctx) 2215 { 2216 TS ts = (TS)ctx; 2217 PetscErrorCode ierr; 2218 2219 PetscFunctionBegin; 2220 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2221 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 2222 PetscValidPointer(A,3); 2223 PetscValidHeaderSpecific(*A,MAT_CLASSID,3); 2224 PetscValidPointer(B,4); 2225 PetscValidHeaderSpecific(*B,MAT_CLASSID,4); 2226 PetscValidPointer(flag,5); 2227 PetscValidHeaderSpecific(ts,TS_CLASSID,6); 2228 ierr = (ts->ops->snesjacobian)(snes,X,A,B,flag,ts);CHKERRQ(ierr); 2229 PetscFunctionReturn(0); 2230 } 2231 2232 #if defined(PETSC_HAVE_MATLAB_ENGINE) 2233 #include <mex.h> 2234 2235 typedef struct {char *funcname; mxArray *ctx;} TSMatlabContext; 2236 2237 #undef __FUNCT__ 2238 #define __FUNCT__ "TSComputeFunction_Matlab" 2239 /* 2240 TSComputeFunction_Matlab - Calls the function that has been set with 2241 TSSetFunctionMatlab(). 2242 2243 Collective on TS 2244 2245 Input Parameters: 2246 + snes - the TS context 2247 - x - input vector 2248 2249 Output Parameter: 2250 . y - function vector, as set by TSSetFunction() 2251 2252 Notes: 2253 TSComputeFunction() is typically used within nonlinear solvers 2254 implementations, so most users would not generally call this routine 2255 themselves. 2256 2257 Level: developer 2258 2259 .keywords: TS, nonlinear, compute, function 2260 2261 .seealso: TSSetFunction(), TSGetFunction() 2262 */ 2263 PetscErrorCode TSComputeFunction_Matlab(TS snes,PetscReal time,Vec x,Vec xdot,Vec y, void *ctx) 2264 { 2265 PetscErrorCode ierr; 2266 TSMatlabContext *sctx = (TSMatlabContext *)ctx; 2267 int nlhs = 1,nrhs = 7; 2268 mxArray *plhs[1],*prhs[7]; 2269 long long int lx = 0,lxdot = 0,ly = 0,ls = 0; 2270 2271 PetscFunctionBegin; 2272 PetscValidHeaderSpecific(snes,TS_CLASSID,1); 2273 PetscValidHeaderSpecific(x,VEC_CLASSID,3); 2274 PetscValidHeaderSpecific(xdot,VEC_CLASSID,4); 2275 PetscValidHeaderSpecific(y,VEC_CLASSID,5); 2276 PetscCheckSameComm(snes,1,x,3); 2277 PetscCheckSameComm(snes,1,y,5); 2278 2279 ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 2280 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 2281 ierr = PetscMemcpy(&lxdot,&xdot,sizeof(xdot));CHKERRQ(ierr); 2282 ierr = PetscMemcpy(&ly,&y,sizeof(x));CHKERRQ(ierr); 2283 prhs[0] = mxCreateDoubleScalar((double)ls); 2284 prhs[1] = mxCreateDoubleScalar(time); 2285 prhs[2] = mxCreateDoubleScalar((double)lx); 2286 prhs[3] = mxCreateDoubleScalar((double)lxdot); 2287 prhs[4] = mxCreateDoubleScalar((double)ly); 2288 prhs[5] = mxCreateString(sctx->funcname); 2289 prhs[6] = sctx->ctx; 2290 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeFunctionInternal");CHKERRQ(ierr); 2291 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 2292 mxDestroyArray(prhs[0]); 2293 mxDestroyArray(prhs[1]); 2294 mxDestroyArray(prhs[2]); 2295 mxDestroyArray(prhs[3]); 2296 mxDestroyArray(prhs[4]); 2297 mxDestroyArray(prhs[5]); 2298 mxDestroyArray(plhs[0]); 2299 PetscFunctionReturn(0); 2300 } 2301 2302 2303 #undef __FUNCT__ 2304 #define __FUNCT__ "TSSetFunctionMatlab" 2305 /* 2306 TSSetFunctionMatlab - Sets the function evaluation routine and function 2307 vector for use by the TS routines in solving ODEs 2308 equations from MATLAB. Here the function is a string containing the name of a MATLAB function 2309 2310 Logically Collective on TS 2311 2312 Input Parameters: 2313 + ts - the TS context 2314 - func - function evaluation routine 2315 2316 Calling sequence of func: 2317 $ func (TS ts,PetscReal time,Vec x,Vec xdot,Vec f,void *ctx); 2318 2319 Level: beginner 2320 2321 .keywords: TS, nonlinear, set, function 2322 2323 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction() 2324 */ 2325 PetscErrorCode TSSetFunctionMatlab(TS snes,const char *func,mxArray *ctx) 2326 { 2327 PetscErrorCode ierr; 2328 TSMatlabContext *sctx; 2329 2330 PetscFunctionBegin; 2331 /* currently sctx is memory bleed */ 2332 ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr); 2333 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 2334 /* 2335 This should work, but it doesn't 2336 sctx->ctx = ctx; 2337 mexMakeArrayPersistent(sctx->ctx); 2338 */ 2339 sctx->ctx = mxDuplicateArray(ctx); 2340 ierr = TSSetIFunction(snes,TSComputeFunction_Matlab,sctx);CHKERRQ(ierr); 2341 PetscFunctionReturn(0); 2342 } 2343 2344 #undef __FUNCT__ 2345 #define __FUNCT__ "TSComputeJacobian_Matlab" 2346 /* 2347 TSComputeJacobian_Matlab - Calls the function that has been set with 2348 TSSetJacobianMatlab(). 2349 2350 Collective on TS 2351 2352 Input Parameters: 2353 + snes - the TS context 2354 . x - input vector 2355 . A, B - the matrices 2356 - ctx - user context 2357 2358 Output Parameter: 2359 . flag - structure of the matrix 2360 2361 Level: developer 2362 2363 .keywords: TS, nonlinear, compute, function 2364 2365 .seealso: TSSetFunction(), TSGetFunction() 2366 @*/ 2367 PetscErrorCode TSComputeJacobian_Matlab(TS snes,PetscReal time,Vec x,Vec xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flag, void *ctx) 2368 { 2369 PetscErrorCode ierr; 2370 TSMatlabContext *sctx = (TSMatlabContext *)ctx; 2371 int nlhs = 2,nrhs = 9; 2372 mxArray *plhs[2],*prhs[9]; 2373 long long int lx = 0,lxdot = 0,lA = 0,ls = 0, lB = 0; 2374 2375 PetscFunctionBegin; 2376 PetscValidHeaderSpecific(snes,TS_CLASSID,1); 2377 PetscValidHeaderSpecific(x,VEC_CLASSID,3); 2378 2379 /* call Matlab function in ctx with arguments x and y */ 2380 2381 ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 2382 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 2383 ierr = PetscMemcpy(&lxdot,&xdot,sizeof(x));CHKERRQ(ierr); 2384 ierr = PetscMemcpy(&lA,A,sizeof(x));CHKERRQ(ierr); 2385 ierr = PetscMemcpy(&lB,B,sizeof(x));CHKERRQ(ierr); 2386 prhs[0] = mxCreateDoubleScalar((double)ls); 2387 prhs[1] = mxCreateDoubleScalar((double)time); 2388 prhs[2] = mxCreateDoubleScalar((double)lx); 2389 prhs[3] = mxCreateDoubleScalar((double)lxdot); 2390 prhs[4] = mxCreateDoubleScalar((double)shift); 2391 prhs[5] = mxCreateDoubleScalar((double)lA); 2392 prhs[6] = mxCreateDoubleScalar((double)lB); 2393 prhs[7] = mxCreateString(sctx->funcname); 2394 prhs[8] = sctx->ctx; 2395 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeJacobianInternal");CHKERRQ(ierr); 2396 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 2397 *flag = (MatStructure) mxGetScalar(plhs[1]);CHKERRQ(ierr); 2398 mxDestroyArray(prhs[0]); 2399 mxDestroyArray(prhs[1]); 2400 mxDestroyArray(prhs[2]); 2401 mxDestroyArray(prhs[3]); 2402 mxDestroyArray(prhs[4]); 2403 mxDestroyArray(prhs[5]); 2404 mxDestroyArray(prhs[6]); 2405 mxDestroyArray(prhs[7]); 2406 mxDestroyArray(plhs[0]); 2407 mxDestroyArray(plhs[1]); 2408 PetscFunctionReturn(0); 2409 } 2410 2411 2412 #undef __FUNCT__ 2413 #define __FUNCT__ "TSSetJacobianMatlab" 2414 /* 2415 TSSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices 2416 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 2417 2418 Logically Collective on TS 2419 2420 Input Parameters: 2421 + snes - the TS context 2422 . A,B - Jacobian matrices 2423 . func - function evaluation routine 2424 - ctx - user context 2425 2426 Calling sequence of func: 2427 $ flag = func (TS snes,PetscReal time,Vec x,Vec xdot,Mat A,Mat B,void *ctx); 2428 2429 2430 Level: developer 2431 2432 .keywords: TS, nonlinear, set, function 2433 2434 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction() 2435 */ 2436 PetscErrorCode TSSetJacobianMatlab(TS snes,Mat A,Mat B,const char *func,mxArray *ctx) 2437 { 2438 PetscErrorCode ierr; 2439 TSMatlabContext *sctx; 2440 2441 PetscFunctionBegin; 2442 /* currently sctx is memory bleed */ 2443 ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr); 2444 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 2445 /* 2446 This should work, but it doesn't 2447 sctx->ctx = ctx; 2448 mexMakeArrayPersistent(sctx->ctx); 2449 */ 2450 sctx->ctx = mxDuplicateArray(ctx); 2451 ierr = TSSetIJacobian(snes,A,B,TSComputeJacobian_Matlab,sctx);CHKERRQ(ierr); 2452 PetscFunctionReturn(0); 2453 } 2454 2455 #undef __FUNCT__ 2456 #define __FUNCT__ "TSMonitor_Matlab" 2457 /* 2458 TSMonitor_Matlab - Calls the function that has been set with TSMonitorSetMatlab(). 2459 2460 Collective on TS 2461 2462 .seealso: TSSetFunction(), TSGetFunction() 2463 @*/ 2464 PetscErrorCode TSMonitor_Matlab(TS snes,PetscInt it, PetscReal time,Vec x, void *ctx) 2465 { 2466 PetscErrorCode ierr; 2467 TSMatlabContext *sctx = (TSMatlabContext *)ctx; 2468 int nlhs = 1,nrhs = 6; 2469 mxArray *plhs[1],*prhs[6]; 2470 long long int lx = 0,ls = 0; 2471 2472 PetscFunctionBegin; 2473 PetscValidHeaderSpecific(snes,TS_CLASSID,1); 2474 PetscValidHeaderSpecific(x,VEC_CLASSID,4); 2475 2476 ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 2477 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 2478 prhs[0] = mxCreateDoubleScalar((double)ls); 2479 prhs[1] = mxCreateDoubleScalar((double)it); 2480 prhs[2] = mxCreateDoubleScalar((double)time); 2481 prhs[3] = mxCreateDoubleScalar((double)lx); 2482 prhs[4] = mxCreateString(sctx->funcname); 2483 prhs[5] = sctx->ctx; 2484 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSMonitorInternal");CHKERRQ(ierr); 2485 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 2486 mxDestroyArray(prhs[0]); 2487 mxDestroyArray(prhs[1]); 2488 mxDestroyArray(prhs[2]); 2489 mxDestroyArray(prhs[3]); 2490 mxDestroyArray(prhs[4]); 2491 mxDestroyArray(plhs[0]); 2492 PetscFunctionReturn(0); 2493 } 2494 2495 2496 #undef __FUNCT__ 2497 #define __FUNCT__ "TSMonitorSetMatlab" 2498 /* 2499 TSMonitorSetMatlab - Sets the monitor function from Matlab 2500 2501 Level: developer 2502 2503 .keywords: TS, nonlinear, set, function 2504 2505 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction() 2506 */ 2507 PetscErrorCode TSMonitorSetMatlab(TS snes,const char *func,mxArray *ctx) 2508 { 2509 PetscErrorCode ierr; 2510 TSMatlabContext *sctx; 2511 2512 PetscFunctionBegin; 2513 /* currently sctx is memory bleed */ 2514 ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr); 2515 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 2516 /* 2517 This should work, but it doesn't 2518 sctx->ctx = ctx; 2519 mexMakeArrayPersistent(sctx->ctx); 2520 */ 2521 sctx->ctx = mxDuplicateArray(ctx); 2522 ierr = TSMonitorSet(snes,TSMonitor_Matlab,sctx,PETSC_NULL);CHKERRQ(ierr); 2523 PetscFunctionReturn(0); 2524 } 2525 2526 #endif 2527