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_always; 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 1203 if (ts->dm) {ierr = DMDestroy(ts->dm);CHKERRQ(ierr);} 1204 if (ts->A) {ierr = MatDestroy(ts->A);CHKERRQ(ierr);} 1205 if (ts->ksp) {ierr = KSPDestroy(ts->ksp);CHKERRQ(ierr);} 1206 if (ts->snes) {ierr = SNESDestroy(ts->snes);CHKERRQ(ierr);} 1207 if (ts->ops->destroy) {ierr = (*(ts)->ops->destroy)(ts);CHKERRQ(ierr);} 1208 ierr = TSMonitorCancel(ts);CHKERRQ(ierr); 1209 ierr = PetscHeaderDestroy(ts);CHKERRQ(ierr); 1210 PetscFunctionReturn(0); 1211 } 1212 1213 #undef __FUNCT__ 1214 #define __FUNCT__ "TSGetSNES" 1215 /*@ 1216 TSGetSNES - Returns the SNES (nonlinear solver) associated with 1217 a TS (timestepper) context. Valid only for nonlinear problems. 1218 1219 Not Collective, but SNES is parallel if TS is parallel 1220 1221 Input Parameter: 1222 . ts - the TS context obtained from TSCreate() 1223 1224 Output Parameter: 1225 . snes - the nonlinear solver context 1226 1227 Notes: 1228 The user can then directly manipulate the SNES context to set various 1229 options, etc. Likewise, the user can then extract and manipulate the 1230 KSP, KSP, and PC contexts as well. 1231 1232 TSGetSNES() does not work for integrators that do not use SNES; in 1233 this case TSGetSNES() returns PETSC_NULL in snes. 1234 1235 Level: beginner 1236 1237 .keywords: timestep, get, SNES 1238 @*/ 1239 PetscErrorCode TSGetSNES(TS ts,SNES *snes) 1240 { 1241 PetscFunctionBegin; 1242 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1243 PetscValidPointer(snes,2); 1244 if (!((PetscObject)ts)->type_name) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"SNES is not created yet. Call TSSetType() first"); 1245 if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Nonlinear only; use TSGetKSP()"); 1246 *snes = ts->snes; 1247 PetscFunctionReturn(0); 1248 } 1249 1250 #undef __FUNCT__ 1251 #define __FUNCT__ "TSGetKSP" 1252 /*@ 1253 TSGetKSP - Returns the KSP (linear solver) associated with 1254 a TS (timestepper) context. 1255 1256 Not Collective, but KSP is parallel if TS is parallel 1257 1258 Input Parameter: 1259 . ts - the TS context obtained from TSCreate() 1260 1261 Output Parameter: 1262 . ksp - the nonlinear solver context 1263 1264 Notes: 1265 The user can then directly manipulate the KSP context to set various 1266 options, etc. Likewise, the user can then extract and manipulate the 1267 KSP and PC contexts as well. 1268 1269 TSGetKSP() does not work for integrators that do not use KSP; 1270 in this case TSGetKSP() returns PETSC_NULL in ksp. 1271 1272 Level: beginner 1273 1274 .keywords: timestep, get, KSP 1275 @*/ 1276 PetscErrorCode TSGetKSP(TS ts,KSP *ksp) 1277 { 1278 PetscFunctionBegin; 1279 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1280 PetscValidPointer(ksp,2); 1281 if (!((PetscObject)ts)->type_name) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"KSP is not created yet. Call TSSetType() first"); 1282 if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()"); 1283 *ksp = ts->ksp; 1284 PetscFunctionReturn(0); 1285 } 1286 1287 /* ----------- Routines to set solver parameters ---------- */ 1288 1289 #undef __FUNCT__ 1290 #define __FUNCT__ "TSGetDuration" 1291 /*@ 1292 TSGetDuration - Gets the maximum number of timesteps to use and 1293 maximum time for iteration. 1294 1295 Not Collective 1296 1297 Input Parameters: 1298 + ts - the TS context obtained from TSCreate() 1299 . maxsteps - maximum number of iterations to use, or PETSC_NULL 1300 - maxtime - final time to iterate to, or PETSC_NULL 1301 1302 Level: intermediate 1303 1304 .keywords: TS, timestep, get, maximum, iterations, time 1305 @*/ 1306 PetscErrorCode TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime) 1307 { 1308 PetscFunctionBegin; 1309 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1310 if (maxsteps) { 1311 PetscValidIntPointer(maxsteps,2); 1312 *maxsteps = ts->max_steps; 1313 } 1314 if (maxtime ) { 1315 PetscValidScalarPointer(maxtime,3); 1316 *maxtime = ts->max_time; 1317 } 1318 PetscFunctionReturn(0); 1319 } 1320 1321 #undef __FUNCT__ 1322 #define __FUNCT__ "TSSetDuration" 1323 /*@ 1324 TSSetDuration - Sets the maximum number of timesteps to use and 1325 maximum time for iteration. 1326 1327 Logically Collective on TS 1328 1329 Input Parameters: 1330 + ts - the TS context obtained from TSCreate() 1331 . maxsteps - maximum number of iterations to use 1332 - maxtime - final time to iterate to 1333 1334 Options Database Keys: 1335 . -ts_max_steps <maxsteps> - Sets maxsteps 1336 . -ts_max_time <maxtime> - Sets maxtime 1337 1338 Notes: 1339 The default maximum number of iterations is 5000. Default time is 5.0 1340 1341 Level: intermediate 1342 1343 .keywords: TS, timestep, set, maximum, iterations 1344 @*/ 1345 PetscErrorCode TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime) 1346 { 1347 PetscFunctionBegin; 1348 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1349 PetscValidLogicalCollectiveInt(ts,maxsteps,2); 1350 PetscValidLogicalCollectiveReal(ts,maxtime,2); 1351 ts->max_steps = maxsteps; 1352 ts->max_time = maxtime; 1353 PetscFunctionReturn(0); 1354 } 1355 1356 #undef __FUNCT__ 1357 #define __FUNCT__ "TSSetSolution" 1358 /*@ 1359 TSSetSolution - Sets the initial solution vector 1360 for use by the TS routines. 1361 1362 Logically Collective on TS and Vec 1363 1364 Input Parameters: 1365 + ts - the TS context obtained from TSCreate() 1366 - x - the solution vector 1367 1368 Level: beginner 1369 1370 .keywords: TS, timestep, set, solution, initial conditions 1371 @*/ 1372 PetscErrorCode TSSetSolution(TS ts,Vec x) 1373 { 1374 PetscFunctionBegin; 1375 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1376 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 1377 ts->vec_sol = ts->vec_sol_always = x; 1378 PetscFunctionReturn(0); 1379 } 1380 1381 #undef __FUNCT__ 1382 #define __FUNCT__ "TSSetPreStep" 1383 /*@C 1384 TSSetPreStep - Sets the general-purpose function 1385 called once at the beginning of each time step. 1386 1387 Logically Collective on TS 1388 1389 Input Parameters: 1390 + ts - The TS context obtained from TSCreate() 1391 - func - The function 1392 1393 Calling sequence of func: 1394 . func (TS ts); 1395 1396 Level: intermediate 1397 1398 .keywords: TS, timestep 1399 @*/ 1400 PetscErrorCode TSSetPreStep(TS ts, PetscErrorCode (*func)(TS)) 1401 { 1402 PetscFunctionBegin; 1403 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1404 ts->ops->prestep = func; 1405 PetscFunctionReturn(0); 1406 } 1407 1408 #undef __FUNCT__ 1409 #define __FUNCT__ "TSPreStep" 1410 /*@C 1411 TSPreStep - Runs the user-defined pre-step function. 1412 1413 Collective on TS 1414 1415 Input Parameters: 1416 . ts - The TS context obtained from TSCreate() 1417 1418 Notes: 1419 TSPreStep() is typically used within time stepping implementations, 1420 so most users would not generally call this routine themselves. 1421 1422 Level: developer 1423 1424 .keywords: TS, timestep 1425 @*/ 1426 PetscErrorCode TSPreStep(TS ts) 1427 { 1428 PetscErrorCode ierr; 1429 1430 PetscFunctionBegin; 1431 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1432 if (ts->ops->prestep) { 1433 PetscStackPush("TS PreStep function"); 1434 ierr = (*ts->ops->prestep)(ts);CHKERRQ(ierr); 1435 PetscStackPop; 1436 } 1437 PetscFunctionReturn(0); 1438 } 1439 1440 #undef __FUNCT__ 1441 #define __FUNCT__ "TSDefaultPreStep" 1442 /*@ 1443 TSDefaultPreStep - The default pre-stepping function which does nothing. 1444 1445 Collective on TS 1446 1447 Input Parameters: 1448 . ts - The TS context obtained from TSCreate() 1449 1450 Level: developer 1451 1452 .keywords: TS, timestep 1453 @*/ 1454 PetscErrorCode TSDefaultPreStep(TS ts) 1455 { 1456 PetscFunctionBegin; 1457 PetscFunctionReturn(0); 1458 } 1459 1460 #undef __FUNCT__ 1461 #define __FUNCT__ "TSSetPostStep" 1462 /*@C 1463 TSSetPostStep - Sets the general-purpose function 1464 called once at the end of each time step. 1465 1466 Logically Collective on TS 1467 1468 Input Parameters: 1469 + ts - The TS context obtained from TSCreate() 1470 - func - The function 1471 1472 Calling sequence of func: 1473 . func (TS ts); 1474 1475 Level: intermediate 1476 1477 .keywords: TS, timestep 1478 @*/ 1479 PetscErrorCode TSSetPostStep(TS ts, PetscErrorCode (*func)(TS)) 1480 { 1481 PetscFunctionBegin; 1482 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1483 ts->ops->poststep = func; 1484 PetscFunctionReturn(0); 1485 } 1486 1487 #undef __FUNCT__ 1488 #define __FUNCT__ "TSPostStep" 1489 /*@C 1490 TSPostStep - Runs the user-defined post-step function. 1491 1492 Collective on TS 1493 1494 Input Parameters: 1495 . ts - The TS context obtained from TSCreate() 1496 1497 Notes: 1498 TSPostStep() is typically used within time stepping implementations, 1499 so most users would not generally call this routine themselves. 1500 1501 Level: developer 1502 1503 .keywords: TS, timestep 1504 @*/ 1505 PetscErrorCode TSPostStep(TS ts) 1506 { 1507 PetscErrorCode ierr; 1508 1509 PetscFunctionBegin; 1510 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1511 if (ts->ops->poststep) { 1512 PetscStackPush("TS PostStep function"); 1513 ierr = (*ts->ops->poststep)(ts);CHKERRQ(ierr); 1514 PetscStackPop; 1515 } 1516 PetscFunctionReturn(0); 1517 } 1518 1519 #undef __FUNCT__ 1520 #define __FUNCT__ "TSDefaultPostStep" 1521 /*@ 1522 TSDefaultPostStep - The default post-stepping function which does nothing. 1523 1524 Collective on TS 1525 1526 Input Parameters: 1527 . ts - The TS context obtained from TSCreate() 1528 1529 Level: developer 1530 1531 .keywords: TS, timestep 1532 @*/ 1533 PetscErrorCode TSDefaultPostStep(TS ts) 1534 { 1535 PetscFunctionBegin; 1536 PetscFunctionReturn(0); 1537 } 1538 1539 /* ------------ Routines to set performance monitoring options ----------- */ 1540 1541 #undef __FUNCT__ 1542 #define __FUNCT__ "TSMonitorSet" 1543 /*@C 1544 TSMonitorSet - Sets an ADDITIONAL function that is to be used at every 1545 timestep to display the iteration's progress. 1546 1547 Logically Collective on TS 1548 1549 Input Parameters: 1550 + ts - the TS context obtained from TSCreate() 1551 . func - monitoring routine 1552 . mctx - [optional] user-defined context for private data for the 1553 monitor routine (use PETSC_NULL if no context is desired) 1554 - monitordestroy - [optional] routine that frees monitor context 1555 (may be PETSC_NULL) 1556 1557 Calling sequence of func: 1558 $ int func(TS ts,PetscInt steps,PetscReal time,Vec x,void *mctx) 1559 1560 + ts - the TS context 1561 . steps - iteration number 1562 . time - current time 1563 . x - current iterate 1564 - mctx - [optional] monitoring context 1565 1566 Notes: 1567 This routine adds an additional monitor to the list of monitors that 1568 already has been loaded. 1569 1570 Fortran notes: Only a single monitor function can be set for each TS object 1571 1572 Level: intermediate 1573 1574 .keywords: TS, timestep, set, monitor 1575 1576 .seealso: TSMonitorDefault(), TSMonitorCancel() 1577 @*/ 1578 PetscErrorCode TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void*)) 1579 { 1580 PetscFunctionBegin; 1581 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1582 if (ts->numbermonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1583 ts->monitor[ts->numbermonitors] = monitor; 1584 ts->mdestroy[ts->numbermonitors] = mdestroy; 1585 ts->monitorcontext[ts->numbermonitors++] = (void*)mctx; 1586 PetscFunctionReturn(0); 1587 } 1588 1589 #undef __FUNCT__ 1590 #define __FUNCT__ "TSMonitorCancel" 1591 /*@C 1592 TSMonitorCancel - Clears all the monitors that have been set on a time-step object. 1593 1594 Logically Collective on TS 1595 1596 Input Parameters: 1597 . ts - the TS context obtained from TSCreate() 1598 1599 Notes: 1600 There is no way to remove a single, specific monitor. 1601 1602 Level: intermediate 1603 1604 .keywords: TS, timestep, set, monitor 1605 1606 .seealso: TSMonitorDefault(), TSMonitorSet() 1607 @*/ 1608 PetscErrorCode TSMonitorCancel(TS ts) 1609 { 1610 PetscErrorCode ierr; 1611 PetscInt i; 1612 1613 PetscFunctionBegin; 1614 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1615 for (i=0; i<ts->numbermonitors; i++) { 1616 if (ts->mdestroy[i]) { 1617 ierr = (*ts->mdestroy[i])(ts->monitorcontext[i]);CHKERRQ(ierr); 1618 } 1619 } 1620 ts->numbermonitors = 0; 1621 PetscFunctionReturn(0); 1622 } 1623 1624 #undef __FUNCT__ 1625 #define __FUNCT__ "TSMonitorDefault" 1626 /*@ 1627 TSMonitorDefault - Sets the Default monitor 1628 1629 Level: intermediate 1630 1631 .keywords: TS, set, monitor 1632 1633 .seealso: TSMonitorDefault(), TSMonitorSet() 1634 @*/ 1635 PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *ctx) 1636 { 1637 PetscErrorCode ierr; 1638 PetscViewerASCIIMonitor viewer = (PetscViewerASCIIMonitor)ctx; 1639 1640 PetscFunctionBegin; 1641 if (!ctx) { 1642 ierr = PetscViewerASCIIMonitorCreate(((PetscObject)ts)->comm,"stdout",0,&viewer);CHKERRQ(ierr); 1643 } 1644 ierr = PetscViewerASCIIMonitorPrintf(viewer,"%D TS dt %G time %G\n",step,ts->time_step,ptime);CHKERRQ(ierr); 1645 if (!ctx) { 1646 ierr = PetscViewerASCIIMonitorDestroy(viewer);CHKERRQ(ierr); 1647 } 1648 PetscFunctionReturn(0); 1649 } 1650 1651 #undef __FUNCT__ 1652 #define __FUNCT__ "TSStep" 1653 /*@ 1654 TSStep - Steps the requested number of timesteps. 1655 1656 Collective on TS 1657 1658 Input Parameter: 1659 . ts - the TS context obtained from TSCreate() 1660 1661 Output Parameters: 1662 + steps - number of iterations until termination 1663 - ptime - time until termination 1664 1665 Level: beginner 1666 1667 .keywords: TS, timestep, solve 1668 1669 .seealso: TSCreate(), TSSetUp(), TSDestroy() 1670 @*/ 1671 PetscErrorCode TSStep(TS ts,PetscInt *steps,PetscReal *ptime) 1672 { 1673 PetscErrorCode ierr; 1674 1675 PetscFunctionBegin; 1676 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1677 if (!ts->setupcalled) { 1678 ierr = TSSetUp(ts);CHKERRQ(ierr); 1679 } 1680 1681 ierr = PetscLogEventBegin(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr); 1682 ierr = (*ts->ops->step)(ts, steps, ptime);CHKERRQ(ierr); 1683 ierr = PetscLogEventEnd(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr); 1684 1685 if (!PetscPreLoadingOn) { 1686 ierr = TSViewFromOptions(ts,((PetscObject)ts)->name);CHKERRQ(ierr); 1687 } 1688 PetscFunctionReturn(0); 1689 } 1690 1691 #undef __FUNCT__ 1692 #define __FUNCT__ "TSSolve" 1693 /*@ 1694 TSSolve - Steps the requested number of timesteps. 1695 1696 Collective on TS 1697 1698 Input Parameter: 1699 + ts - the TS context obtained from TSCreate() 1700 - x - the solution vector, or PETSC_NULL if it was set with TSSetSolution() 1701 1702 Level: beginner 1703 1704 .keywords: TS, timestep, solve 1705 1706 .seealso: TSCreate(), TSSetSolution(), TSStep() 1707 @*/ 1708 PetscErrorCode TSSolve(TS ts, Vec x) 1709 { 1710 PetscInt steps; 1711 PetscReal ptime; 1712 PetscErrorCode ierr; 1713 1714 PetscFunctionBegin; 1715 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1716 /* set solution vector if provided */ 1717 if (x) { ierr = TSSetSolution(ts, x); CHKERRQ(ierr); } 1718 /* reset time step and iteration counters */ 1719 ts->steps = 0; ts->linear_its = 0; ts->nonlinear_its = 0; 1720 /* steps the requested number of timesteps. */ 1721 ierr = TSStep(ts, &steps, &ptime);CHKERRQ(ierr); 1722 PetscFunctionReturn(0); 1723 } 1724 1725 #undef __FUNCT__ 1726 #define __FUNCT__ "TSMonitor" 1727 /* 1728 Runs the user provided monitor routines, if they exists. 1729 */ 1730 PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x) 1731 { 1732 PetscErrorCode ierr; 1733 PetscInt i,n = ts->numbermonitors; 1734 1735 PetscFunctionBegin; 1736 for (i=0; i<n; i++) { 1737 ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr); 1738 } 1739 PetscFunctionReturn(0); 1740 } 1741 1742 /* ------------------------------------------------------------------------*/ 1743 1744 #undef __FUNCT__ 1745 #define __FUNCT__ "TSMonitorLGCreate" 1746 /*@C 1747 TSMonitorLGCreate - Creates a line graph context for use with 1748 TS to monitor convergence of preconditioned residual norms. 1749 1750 Collective on TS 1751 1752 Input Parameters: 1753 + host - the X display to open, or null for the local machine 1754 . label - the title to put in the title bar 1755 . x, y - the screen coordinates of the upper left coordinate of the window 1756 - m, n - the screen width and height in pixels 1757 1758 Output Parameter: 1759 . draw - the drawing context 1760 1761 Options Database Key: 1762 . -ts_monitor_draw - automatically sets line graph monitor 1763 1764 Notes: 1765 Use TSMonitorLGDestroy() to destroy this line graph, not PetscDrawLGDestroy(). 1766 1767 Level: intermediate 1768 1769 .keywords: TS, monitor, line graph, residual, seealso 1770 1771 .seealso: TSMonitorLGDestroy(), TSMonitorSet() 1772 1773 @*/ 1774 PetscErrorCode TSMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1775 { 1776 PetscDraw win; 1777 PetscErrorCode ierr; 1778 1779 PetscFunctionBegin; 1780 ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr); 1781 ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr); 1782 ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr); 1783 ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr); 1784 1785 ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr); 1786 PetscFunctionReturn(0); 1787 } 1788 1789 #undef __FUNCT__ 1790 #define __FUNCT__ "TSMonitorLG" 1791 PetscErrorCode TSMonitorLG(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx) 1792 { 1793 PetscDrawLG lg = (PetscDrawLG) monctx; 1794 PetscReal x,y = ptime; 1795 PetscErrorCode ierr; 1796 1797 PetscFunctionBegin; 1798 if (!monctx) { 1799 MPI_Comm comm; 1800 PetscViewer viewer; 1801 1802 ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr); 1803 viewer = PETSC_VIEWER_DRAW_(comm); 1804 ierr = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr); 1805 } 1806 1807 if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 1808 x = (PetscReal)n; 1809 ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 1810 if (n < 20 || (n % 5)) { 1811 ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 1812 } 1813 PetscFunctionReturn(0); 1814 } 1815 1816 #undef __FUNCT__ 1817 #define __FUNCT__ "TSMonitorLGDestroy" 1818 /*@C 1819 TSMonitorLGDestroy - Destroys a line graph context that was created 1820 with TSMonitorLGCreate(). 1821 1822 Collective on PetscDrawLG 1823 1824 Input Parameter: 1825 . draw - the drawing context 1826 1827 Level: intermediate 1828 1829 .keywords: TS, monitor, line graph, destroy 1830 1831 .seealso: TSMonitorLGCreate(), TSMonitorSet(), TSMonitorLG(); 1832 @*/ 1833 PetscErrorCode TSMonitorLGDestroy(PetscDrawLG drawlg) 1834 { 1835 PetscDraw draw; 1836 PetscErrorCode ierr; 1837 1838 PetscFunctionBegin; 1839 ierr = PetscDrawLGGetDraw(drawlg,&draw);CHKERRQ(ierr); 1840 ierr = PetscDrawDestroy(draw);CHKERRQ(ierr); 1841 ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr); 1842 PetscFunctionReturn(0); 1843 } 1844 1845 #undef __FUNCT__ 1846 #define __FUNCT__ "TSGetTime" 1847 /*@ 1848 TSGetTime - Gets the current time. 1849 1850 Not Collective 1851 1852 Input Parameter: 1853 . ts - the TS context obtained from TSCreate() 1854 1855 Output Parameter: 1856 . t - the current time 1857 1858 Level: beginner 1859 1860 .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 1861 1862 .keywords: TS, get, time 1863 @*/ 1864 PetscErrorCode TSGetTime(TS ts,PetscReal* t) 1865 { 1866 PetscFunctionBegin; 1867 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1868 PetscValidDoublePointer(t,2); 1869 *t = ts->ptime; 1870 PetscFunctionReturn(0); 1871 } 1872 1873 #undef __FUNCT__ 1874 #define __FUNCT__ "TSSetTime" 1875 /*@ 1876 TSSetTime - Allows one to reset the time. 1877 1878 Logically Collective on TS 1879 1880 Input Parameters: 1881 + ts - the TS context obtained from TSCreate() 1882 - time - the time 1883 1884 Level: intermediate 1885 1886 .seealso: TSGetTime(), TSSetDuration() 1887 1888 .keywords: TS, set, time 1889 @*/ 1890 PetscErrorCode TSSetTime(TS ts, PetscReal t) 1891 { 1892 PetscFunctionBegin; 1893 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1894 PetscValidLogicalCollectiveReal(ts,t,2); 1895 ts->ptime = t; 1896 PetscFunctionReturn(0); 1897 } 1898 1899 #undef __FUNCT__ 1900 #define __FUNCT__ "TSSetOptionsPrefix" 1901 /*@C 1902 TSSetOptionsPrefix - Sets the prefix used for searching for all 1903 TS options in the database. 1904 1905 Logically Collective on TS 1906 1907 Input Parameter: 1908 + ts - The TS context 1909 - prefix - The prefix to prepend to all option names 1910 1911 Notes: 1912 A hyphen (-) must NOT be given at the beginning of the prefix name. 1913 The first character of all runtime options is AUTOMATICALLY the 1914 hyphen. 1915 1916 Level: advanced 1917 1918 .keywords: TS, set, options, prefix, database 1919 1920 .seealso: TSSetFromOptions() 1921 1922 @*/ 1923 PetscErrorCode TSSetOptionsPrefix(TS ts,const char prefix[]) 1924 { 1925 PetscErrorCode ierr; 1926 1927 PetscFunctionBegin; 1928 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1929 ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 1930 switch(ts->problem_type) { 1931 case TS_NONLINEAR: 1932 if (ts->snes) { 1933 ierr = SNESSetOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr); 1934 } 1935 break; 1936 case TS_LINEAR: 1937 if (ts->ksp) { 1938 ierr = KSPSetOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr); 1939 } 1940 break; 1941 } 1942 PetscFunctionReturn(0); 1943 } 1944 1945 1946 #undef __FUNCT__ 1947 #define __FUNCT__ "TSAppendOptionsPrefix" 1948 /*@C 1949 TSAppendOptionsPrefix - Appends to the prefix used for searching for all 1950 TS options in the database. 1951 1952 Logically Collective on TS 1953 1954 Input Parameter: 1955 + ts - The TS context 1956 - prefix - The prefix to prepend to all option names 1957 1958 Notes: 1959 A hyphen (-) must NOT be given at the beginning of the prefix name. 1960 The first character of all runtime options is AUTOMATICALLY the 1961 hyphen. 1962 1963 Level: advanced 1964 1965 .keywords: TS, append, options, prefix, database 1966 1967 .seealso: TSGetOptionsPrefix() 1968 1969 @*/ 1970 PetscErrorCode TSAppendOptionsPrefix(TS ts,const char prefix[]) 1971 { 1972 PetscErrorCode ierr; 1973 1974 PetscFunctionBegin; 1975 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1976 ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 1977 switch(ts->problem_type) { 1978 case TS_NONLINEAR: 1979 if (ts->snes) { 1980 ierr = SNESAppendOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr); 1981 } 1982 break; 1983 case TS_LINEAR: 1984 if (ts->ksp) { 1985 ierr = KSPAppendOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr); 1986 } 1987 break; 1988 } 1989 PetscFunctionReturn(0); 1990 } 1991 1992 #undef __FUNCT__ 1993 #define __FUNCT__ "TSGetOptionsPrefix" 1994 /*@C 1995 TSGetOptionsPrefix - Sets the prefix used for searching for all 1996 TS options in the database. 1997 1998 Not Collective 1999 2000 Input Parameter: 2001 . ts - The TS context 2002 2003 Output Parameter: 2004 . prefix - A pointer to the prefix string used 2005 2006 Notes: On the fortran side, the user should pass in a string 'prifix' of 2007 sufficient length to hold the prefix. 2008 2009 Level: intermediate 2010 2011 .keywords: TS, get, options, prefix, database 2012 2013 .seealso: TSAppendOptionsPrefix() 2014 @*/ 2015 PetscErrorCode TSGetOptionsPrefix(TS ts,const char *prefix[]) 2016 { 2017 PetscErrorCode ierr; 2018 2019 PetscFunctionBegin; 2020 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2021 PetscValidPointer(prefix,2); 2022 ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 2023 PetscFunctionReturn(0); 2024 } 2025 2026 #undef __FUNCT__ 2027 #define __FUNCT__ "TSGetRHSJacobian" 2028 /*@C 2029 TSGetRHSJacobian - Returns the Jacobian J at the present timestep. 2030 2031 Not Collective, but parallel objects are returned if TS is parallel 2032 2033 Input Parameter: 2034 . ts - The TS context obtained from TSCreate() 2035 2036 Output Parameters: 2037 + J - The Jacobian J of F, where U_t = F(U,t) 2038 . M - The preconditioner matrix, usually the same as J 2039 - ctx - User-defined context for Jacobian evaluation routine 2040 2041 Notes: You can pass in PETSC_NULL for any return argument you do not need. 2042 2043 Level: intermediate 2044 2045 .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber() 2046 2047 .keywords: TS, timestep, get, matrix, Jacobian 2048 @*/ 2049 PetscErrorCode TSGetRHSJacobian(TS ts,Mat *J,Mat *M,void **ctx) 2050 { 2051 PetscFunctionBegin; 2052 if (J) *J = ts->Arhs; 2053 if (M) *M = ts->B; 2054 if (ctx) *ctx = ts->jacP; 2055 PetscFunctionReturn(0); 2056 } 2057 2058 #undef __FUNCT__ 2059 #define __FUNCT__ "TSGetIJacobian" 2060 /*@C 2061 TSGetIJacobian - Returns the implicit Jacobian at the present timestep. 2062 2063 Not Collective, but parallel objects are returned if TS is parallel 2064 2065 Input Parameter: 2066 . ts - The TS context obtained from TSCreate() 2067 2068 Output Parameters: 2069 + A - The Jacobian of F(t,U,U_t) 2070 . B - The preconditioner matrix, often the same as A 2071 . f - The function to compute the matrices 2072 - ctx - User-defined context for Jacobian evaluation routine 2073 2074 Notes: You can pass in PETSC_NULL for any return argument you do not need. 2075 2076 Level: advanced 2077 2078 .seealso: TSGetTimeStep(), TSGetRHSJacobian(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber() 2079 2080 .keywords: TS, timestep, get, matrix, Jacobian 2081 @*/ 2082 PetscErrorCode TSGetIJacobian(TS ts,Mat *A,Mat *B,TSIJacobian *f,void **ctx) 2083 { 2084 PetscFunctionBegin; 2085 if (A) *A = ts->A; 2086 if (B) *B = ts->B; 2087 if (f) *f = ts->ops->ijacobian; 2088 if (ctx) *ctx = ts->jacP; 2089 PetscFunctionReturn(0); 2090 } 2091 2092 #undef __FUNCT__ 2093 #define __FUNCT__ "TSMonitorSolution" 2094 /*@C 2095 TSMonitorSolution - Monitors progress of the TS solvers by calling 2096 VecView() for the solution at each timestep 2097 2098 Collective on TS 2099 2100 Input Parameters: 2101 + ts - the TS context 2102 . step - current time-step 2103 . ptime - current time 2104 - dummy - either a viewer or PETSC_NULL 2105 2106 Level: intermediate 2107 2108 .keywords: TS, vector, monitor, view 2109 2110 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView() 2111 @*/ 2112 PetscErrorCode TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy) 2113 { 2114 PetscErrorCode ierr; 2115 PetscViewer viewer = (PetscViewer) dummy; 2116 2117 PetscFunctionBegin; 2118 if (!dummy) { 2119 viewer = PETSC_VIEWER_DRAW_(((PetscObject)ts)->comm); 2120 } 2121 ierr = VecView(x,viewer);CHKERRQ(ierr); 2122 PetscFunctionReturn(0); 2123 } 2124 2125 2126 #undef __FUNCT__ 2127 #define __FUNCT__ "TSSetDM" 2128 /*@ 2129 TSSetDM - Sets the DM that may be used by some preconditioners 2130 2131 Logically Collective on TS and DM 2132 2133 Input Parameters: 2134 + ts - the preconditioner context 2135 - dm - the dm 2136 2137 Level: intermediate 2138 2139 2140 .seealso: TSGetDM(), SNESSetDM(), SNESGetDM() 2141 @*/ 2142 PetscErrorCode TSSetDM(TS ts,DM dm) 2143 { 2144 PetscErrorCode ierr; 2145 2146 PetscFunctionBegin; 2147 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2148 ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr); 2149 if (ts->dm) {ierr = DMDestroy(ts->dm);CHKERRQ(ierr);} 2150 ts->dm = dm; 2151 if (ts->snes) { 2152 ierr = SNESSetDM(ts->snes,dm);CHKERRQ(ierr); 2153 } 2154 if (ts->ksp) { 2155 ierr = KSPSetDM(ts->ksp,dm);CHKERRQ(ierr); 2156 ierr = KSPSetDMActive(ts->ksp,PETSC_FALSE);CHKERRQ(ierr); 2157 } 2158 PetscFunctionReturn(0); 2159 } 2160 2161 #undef __FUNCT__ 2162 #define __FUNCT__ "TSGetDM" 2163 /*@ 2164 TSGetDM - Gets the DM that may be used by some preconditioners 2165 2166 Not Collective 2167 2168 Input Parameter: 2169 . ts - the preconditioner context 2170 2171 Output Parameter: 2172 . dm - the dm 2173 2174 Level: intermediate 2175 2176 2177 .seealso: TSSetDM(), SNESSetDM(), SNESGetDM() 2178 @*/ 2179 PetscErrorCode TSGetDM(TS ts,DM *dm) 2180 { 2181 PetscFunctionBegin; 2182 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2183 *dm = ts->dm; 2184 PetscFunctionReturn(0); 2185 } 2186 2187 #undef __FUNCT__ 2188 #define __FUNCT__ "SNESTSFormFunction" 2189 /*@ 2190 SNESTSFormFunction - Function to evaluate nonlinear residual 2191 2192 Logically Collective on SNES 2193 2194 Input Parameter: 2195 + snes - nonlinear solver 2196 . X - the current state at which to evaluate the residual 2197 - ctx - user context, must be a TS 2198 2199 Output Parameter: 2200 . F - the nonlinear residual 2201 2202 Notes: 2203 This function is not normally called by users and is automatically registered with the SNES used by TS. 2204 It is most frequently passed to MatFDColoringSetFunction(). 2205 2206 Level: advanced 2207 2208 .seealso: SNESSetFunction(), MatFDColoringSetFunction() 2209 @*/ 2210 PetscErrorCode SNESTSFormFunction(SNES snes,Vec X,Vec F,void *ctx) 2211 { 2212 TS ts = (TS)ctx; 2213 PetscErrorCode ierr; 2214 2215 PetscFunctionBegin; 2216 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2217 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 2218 PetscValidHeaderSpecific(F,VEC_CLASSID,3); 2219 PetscValidHeaderSpecific(ts,TS_CLASSID,4); 2220 ierr = (ts->ops->snesfunction)(snes,X,F,ts);CHKERRQ(ierr); 2221 PetscFunctionReturn(0); 2222 } 2223 2224 #undef __FUNCT__ 2225 #define __FUNCT__ "SNESTSFormJacobian" 2226 /*@ 2227 SNESTSFormJacobian - Function to evaluate the Jacobian 2228 2229 Collective on SNES 2230 2231 Input Parameter: 2232 + snes - nonlinear solver 2233 . X - the current state at which to evaluate the residual 2234 - ctx - user context, must be a TS 2235 2236 Output Parameter: 2237 + A - the Jacobian 2238 . B - the preconditioning matrix (may be the same as A) 2239 - flag - indicates any structure change in the matrix 2240 2241 Notes: 2242 This function is not normally called by users and is automatically registered with the SNES used by TS. 2243 2244 Level: developer 2245 2246 .seealso: SNESSetJacobian() 2247 @*/ 2248 PetscErrorCode SNESTSFormJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flag,void *ctx) 2249 { 2250 TS ts = (TS)ctx; 2251 PetscErrorCode ierr; 2252 2253 PetscFunctionBegin; 2254 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2255 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 2256 PetscValidPointer(A,3); 2257 PetscValidHeaderSpecific(*A,MAT_CLASSID,3); 2258 PetscValidPointer(B,4); 2259 PetscValidHeaderSpecific(*B,MAT_CLASSID,4); 2260 PetscValidPointer(flag,5); 2261 PetscValidHeaderSpecific(ts,TS_CLASSID,6); 2262 ierr = (ts->ops->snesjacobian)(snes,X,A,B,flag,ts);CHKERRQ(ierr); 2263 PetscFunctionReturn(0); 2264 } 2265 2266 #if defined(PETSC_HAVE_MATLAB_ENGINE) 2267 #include "mex.h" 2268 2269 typedef struct {char *funcname; mxArray *ctx;} TSMatlabContext; 2270 2271 #undef __FUNCT__ 2272 #define __FUNCT__ "TSComputeFunction_Matlab" 2273 /* 2274 TSComputeFunction_Matlab - Calls the function that has been set with 2275 TSSetFunctionMatlab(). 2276 2277 Collective on TS 2278 2279 Input Parameters: 2280 + snes - the TS context 2281 - x - input vector 2282 2283 Output Parameter: 2284 . y - function vector, as set by TSSetFunction() 2285 2286 Notes: 2287 TSComputeFunction() is typically used within nonlinear solvers 2288 implementations, so most users would not generally call this routine 2289 themselves. 2290 2291 Level: developer 2292 2293 .keywords: TS, nonlinear, compute, function 2294 2295 .seealso: TSSetFunction(), TSGetFunction() 2296 */ 2297 PetscErrorCode TSComputeFunction_Matlab(TS snes,PetscReal time,Vec x,Vec xdot,Vec y, void *ctx) 2298 { 2299 PetscErrorCode ierr; 2300 TSMatlabContext *sctx = (TSMatlabContext *)ctx; 2301 int nlhs = 1,nrhs = 7; 2302 mxArray *plhs[1],*prhs[7]; 2303 long long int lx = 0,lxdot = 0,ly = 0,ls = 0; 2304 2305 PetscFunctionBegin; 2306 PetscValidHeaderSpecific(snes,TS_CLASSID,1); 2307 PetscValidHeaderSpecific(x,VEC_CLASSID,3); 2308 PetscValidHeaderSpecific(xdot,VEC_CLASSID,4); 2309 PetscValidHeaderSpecific(y,VEC_CLASSID,5); 2310 PetscCheckSameComm(snes,1,x,3); 2311 PetscCheckSameComm(snes,1,y,5); 2312 2313 ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 2314 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 2315 ierr = PetscMemcpy(&lxdot,&xdot,sizeof(xdot));CHKERRQ(ierr); 2316 ierr = PetscMemcpy(&ly,&y,sizeof(x));CHKERRQ(ierr); 2317 prhs[0] = mxCreateDoubleScalar((double)ls); 2318 prhs[1] = mxCreateDoubleScalar(time); 2319 prhs[2] = mxCreateDoubleScalar((double)lx); 2320 prhs[3] = mxCreateDoubleScalar((double)lxdot); 2321 prhs[4] = mxCreateDoubleScalar((double)ly); 2322 prhs[5] = mxCreateString(sctx->funcname); 2323 prhs[6] = sctx->ctx; 2324 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeFunctionInternal");CHKERRQ(ierr); 2325 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 2326 mxDestroyArray(prhs[0]); 2327 mxDestroyArray(prhs[1]); 2328 mxDestroyArray(prhs[2]); 2329 mxDestroyArray(prhs[3]); 2330 mxDestroyArray(prhs[4]); 2331 mxDestroyArray(prhs[5]); 2332 mxDestroyArray(plhs[0]); 2333 PetscFunctionReturn(0); 2334 } 2335 2336 2337 #undef __FUNCT__ 2338 #define __FUNCT__ "TSSetFunctionMatlab" 2339 /* 2340 TSSetFunctionMatlab - Sets the function evaluation routine and function 2341 vector for use by the TS routines in solving ODEs 2342 equations from MATLAB. Here the function is a string containing the name of a MATLAB function 2343 2344 Logically Collective on TS 2345 2346 Input Parameters: 2347 + ts - the TS context 2348 - func - function evaluation routine 2349 2350 Calling sequence of func: 2351 $ func (TS ts,PetscReal time,Vec x,Vec xdot,Vec f,void *ctx); 2352 2353 Level: beginner 2354 2355 .keywords: TS, nonlinear, set, function 2356 2357 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction() 2358 */ 2359 PetscErrorCode TSSetFunctionMatlab(TS snes,const char *func,mxArray *ctx) 2360 { 2361 PetscErrorCode ierr; 2362 TSMatlabContext *sctx; 2363 2364 PetscFunctionBegin; 2365 /* currently sctx is memory bleed */ 2366 ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr); 2367 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 2368 /* 2369 This should work, but it doesn't 2370 sctx->ctx = ctx; 2371 mexMakeArrayPersistent(sctx->ctx); 2372 */ 2373 sctx->ctx = mxDuplicateArray(ctx); 2374 ierr = TSSetIFunction(snes,TSComputeFunction_Matlab,sctx);CHKERRQ(ierr); 2375 PetscFunctionReturn(0); 2376 } 2377 2378 #undef __FUNCT__ 2379 #define __FUNCT__ "TSComputeJacobian_Matlab" 2380 /* 2381 TSComputeJacobian_Matlab - Calls the function that has been set with 2382 TSSetJacobianMatlab(). 2383 2384 Collective on TS 2385 2386 Input Parameters: 2387 + snes - the TS context 2388 . x - input vector 2389 . A, B - the matrices 2390 - ctx - user context 2391 2392 Output Parameter: 2393 . flag - structure of the matrix 2394 2395 Level: developer 2396 2397 .keywords: TS, nonlinear, compute, function 2398 2399 .seealso: TSSetFunction(), TSGetFunction() 2400 @*/ 2401 PetscErrorCode TSComputeJacobian_Matlab(TS snes,PetscReal time,Vec x,Vec xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flag, void *ctx) 2402 { 2403 PetscErrorCode ierr; 2404 TSMatlabContext *sctx = (TSMatlabContext *)ctx; 2405 int nlhs = 2,nrhs = 9; 2406 mxArray *plhs[2],*prhs[9]; 2407 long long int lx = 0,lxdot = 0,lA = 0,ls = 0, lB = 0; 2408 2409 PetscFunctionBegin; 2410 PetscValidHeaderSpecific(snes,TS_CLASSID,1); 2411 PetscValidHeaderSpecific(x,VEC_CLASSID,3); 2412 2413 /* call Matlab function in ctx with arguments x and y */ 2414 2415 ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 2416 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 2417 ierr = PetscMemcpy(&lxdot,&xdot,sizeof(x));CHKERRQ(ierr); 2418 ierr = PetscMemcpy(&lA,A,sizeof(x));CHKERRQ(ierr); 2419 ierr = PetscMemcpy(&lB,B,sizeof(x));CHKERRQ(ierr); 2420 prhs[0] = mxCreateDoubleScalar((double)ls); 2421 prhs[1] = mxCreateDoubleScalar((double)time); 2422 prhs[2] = mxCreateDoubleScalar((double)lx); 2423 prhs[3] = mxCreateDoubleScalar((double)lxdot); 2424 prhs[4] = mxCreateDoubleScalar((double)shift); 2425 prhs[5] = mxCreateDoubleScalar((double)lA); 2426 prhs[6] = mxCreateDoubleScalar((double)lB); 2427 prhs[7] = mxCreateString(sctx->funcname); 2428 prhs[8] = sctx->ctx; 2429 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeJacobianInternal");CHKERRQ(ierr); 2430 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 2431 *flag = (MatStructure) mxGetScalar(plhs[1]);CHKERRQ(ierr); 2432 mxDestroyArray(prhs[0]); 2433 mxDestroyArray(prhs[1]); 2434 mxDestroyArray(prhs[2]); 2435 mxDestroyArray(prhs[3]); 2436 mxDestroyArray(prhs[4]); 2437 mxDestroyArray(prhs[5]); 2438 mxDestroyArray(prhs[6]); 2439 mxDestroyArray(prhs[7]); 2440 mxDestroyArray(plhs[0]); 2441 mxDestroyArray(plhs[1]); 2442 PetscFunctionReturn(0); 2443 } 2444 2445 2446 #undef __FUNCT__ 2447 #define __FUNCT__ "TSSetJacobianMatlab" 2448 /* 2449 TSSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices 2450 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 2451 2452 Logically Collective on TS 2453 2454 Input Parameters: 2455 + snes - the TS context 2456 . A,B - Jacobian matrices 2457 . func - function evaluation routine 2458 - ctx - user context 2459 2460 Calling sequence of func: 2461 $ flag = func (TS snes,PetscReal time,Vec x,Vec xdot,Mat A,Mat B,void *ctx); 2462 2463 2464 Level: developer 2465 2466 .keywords: TS, nonlinear, set, function 2467 2468 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction() 2469 */ 2470 PetscErrorCode TSSetJacobianMatlab(TS snes,Mat A,Mat B,const char *func,mxArray *ctx) 2471 { 2472 PetscErrorCode ierr; 2473 TSMatlabContext *sctx; 2474 2475 PetscFunctionBegin; 2476 /* currently sctx is memory bleed */ 2477 ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr); 2478 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 2479 /* 2480 This should work, but it doesn't 2481 sctx->ctx = ctx; 2482 mexMakeArrayPersistent(sctx->ctx); 2483 */ 2484 sctx->ctx = mxDuplicateArray(ctx); 2485 ierr = TSSetIJacobian(snes,A,B,TSComputeJacobian_Matlab,sctx);CHKERRQ(ierr); 2486 PetscFunctionReturn(0); 2487 } 2488 2489 #undef __FUNCT__ 2490 #define __FUNCT__ "TSMonitor_Matlab" 2491 /* 2492 TSMonitor_Matlab - Calls the function that has been set with TSMonitorSetMatlab(). 2493 2494 Collective on TS 2495 2496 .seealso: TSSetFunction(), TSGetFunction() 2497 @*/ 2498 PetscErrorCode TSMonitor_Matlab(TS snes,PetscInt it, PetscReal time,Vec x, void *ctx) 2499 { 2500 PetscErrorCode ierr; 2501 TSMatlabContext *sctx = (TSMatlabContext *)ctx; 2502 int nlhs = 1,nrhs = 6; 2503 mxArray *plhs[1],*prhs[6]; 2504 long long int lx = 0,ls = 0; 2505 2506 PetscFunctionBegin; 2507 PetscValidHeaderSpecific(snes,TS_CLASSID,1); 2508 PetscValidHeaderSpecific(x,VEC_CLASSID,4); 2509 2510 ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 2511 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 2512 prhs[0] = mxCreateDoubleScalar((double)ls); 2513 prhs[1] = mxCreateDoubleScalar((double)it); 2514 prhs[2] = mxCreateDoubleScalar((double)time); 2515 prhs[3] = mxCreateDoubleScalar((double)lx); 2516 prhs[4] = mxCreateString(sctx->funcname); 2517 prhs[5] = sctx->ctx; 2518 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSMonitorInternal");CHKERRQ(ierr); 2519 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 2520 mxDestroyArray(prhs[0]); 2521 mxDestroyArray(prhs[1]); 2522 mxDestroyArray(prhs[2]); 2523 mxDestroyArray(prhs[3]); 2524 mxDestroyArray(prhs[4]); 2525 mxDestroyArray(plhs[0]); 2526 PetscFunctionReturn(0); 2527 } 2528 2529 2530 #undef __FUNCT__ 2531 #define __FUNCT__ "TSMonitorSetMatlab" 2532 /* 2533 TSMonitorSetMatlab - Sets the monitor function from Matlab 2534 2535 Level: developer 2536 2537 .keywords: TS, nonlinear, set, function 2538 2539 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction() 2540 */ 2541 PetscErrorCode TSMonitorSetMatlab(TS snes,const char *func,mxArray *ctx) 2542 { 2543 PetscErrorCode ierr; 2544 TSMatlabContext *sctx; 2545 2546 PetscFunctionBegin; 2547 /* currently sctx is memory bleed */ 2548 ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr); 2549 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 2550 /* 2551 This should work, but it doesn't 2552 sctx->ctx = ctx; 2553 mexMakeArrayPersistent(sctx->ctx); 2554 */ 2555 sctx->ctx = mxDuplicateArray(ctx); 2556 ierr = TSMonitorSet(snes,TSMonitor_Matlab,sctx,PETSC_NULL);CHKERRQ(ierr); 2557 PetscFunctionReturn(0); 2558 } 2559 2560 #endif 2561