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