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