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