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 PetscBool opt,flg; 74 PetscErrorCode ierr; 75 PetscViewer monviewer; 76 char monfilename[PETSC_MAX_PATH_LEN]; 77 SNES snes; 78 79 PetscFunctionBegin; 80 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 81 ierr = PetscObjectOptionsBegin((PetscObject)ts);CHKERRQ(ierr); 82 /* Handle TS type options */ 83 ierr = TSSetTypeFromOptions(ts);CHKERRQ(ierr); 84 85 /* Handle generic TS options */ 86 ierr = PetscOptionsInt("-ts_max_steps","Maximum number of time steps","TSSetDuration",ts->max_steps,&ts->max_steps,PETSC_NULL);CHKERRQ(ierr); 87 ierr = PetscOptionsReal("-ts_max_time","Time to run to","TSSetDuration",ts->max_time,&ts->max_time,PETSC_NULL);CHKERRQ(ierr); 88 ierr = PetscOptionsReal("-ts_init_time","Initial time","TSSetTime",ts->ptime,&ts->ptime,PETSC_NULL);CHKERRQ(ierr); 89 ierr = PetscOptionsReal("-ts_dt","Initial time step","TSSetTimeStep",ts->time_step,&ts->time_step,PETSC_NULL);CHKERRQ(ierr); 90 opt = ts->exact_final_time == PETSC_DECIDE ? PETSC_FALSE : (PetscBool)ts->exact_final_time; 91 ierr = PetscOptionsBool("-ts_exact_final_time","Interpolate output to stop exactly at the final time","TSSetExactFinalTime",opt,&opt,&flg);CHKERRQ(ierr); 92 if (flg) {ierr = TSSetExactFinalTime(ts,opt);CHKERRQ(ierr);} 93 ierr = PetscOptionsInt("-ts_max_snes_failures","Maximum number of nonlinear solve failures","",ts->max_snes_failures,&ts->max_snes_failures,PETSC_NULL);CHKERRQ(ierr); 94 ierr = PetscOptionsInt("-ts_max_reject","Maximum number of step rejections","",ts->max_reject,&ts->max_reject,PETSC_NULL);CHKERRQ(ierr); 95 ierr = PetscOptionsBool("-ts_error_if_step_failed","Error if no step succeeds","",ts->errorifstepfailed,&ts->errorifstepfailed,PETSC_NULL);CHKERRQ(ierr); 96 97 /* Monitor options */ 98 ierr = PetscOptionsString("-ts_monitor","Monitor timestep size","TSMonitorDefault","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 99 if (flg) { 100 ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,monfilename,&monviewer);CHKERRQ(ierr); 101 ierr = TSMonitorSet(ts,TSMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 102 } 103 ierr = PetscOptionsString("-ts_monitor_python","Use Python function","TSMonitorSet",0,monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 104 if (flg) {ierr = PetscPythonMonitorSet((PetscObject)ts,monfilename);CHKERRQ(ierr);} 105 106 opt = PETSC_FALSE; 107 ierr = PetscOptionsBool("-ts_monitor_draw","Monitor timestep size graphically","TSMonitorLG",opt,&opt,PETSC_NULL);CHKERRQ(ierr); 108 if (opt) { 109 ierr = TSMonitorSet(ts,TSMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 110 } 111 opt = PETSC_FALSE; 112 ierr = PetscOptionsBool("-ts_monitor_solution","Monitor solution graphically","TSMonitorSolution",opt,&opt,PETSC_NULL);CHKERRQ(ierr); 113 if (opt) { 114 void *ctx; 115 ierr = TSMonitorSolutionCreate(ts,PETSC_NULL,&ctx);CHKERRQ(ierr); 116 ierr = TSMonitorSet(ts,TSMonitorSolution,ctx,TSMonitorSolutionDestroy);CHKERRQ(ierr); 117 } 118 119 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 120 if (ts->problem_type == TS_LINEAR) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);} 121 122 /* Handle specific TS options */ 123 if (ts->ops->setfromoptions) { 124 ierr = (*ts->ops->setfromoptions)(ts);CHKERRQ(ierr); 125 } 126 127 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 128 ierr = PetscObjectProcessOptionsHandlers((PetscObject)ts);CHKERRQ(ierr); 129 ierr = PetscOptionsEnd();CHKERRQ(ierr); 130 PetscFunctionReturn(0); 131 } 132 133 #undef __FUNCT__ 134 #undef __FUNCT__ 135 #define __FUNCT__ "TSComputeRHSJacobian" 136 /*@ 137 TSComputeRHSJacobian - Computes the Jacobian matrix that has been 138 set with TSSetRHSJacobian(). 139 140 Collective on TS and Vec 141 142 Input Parameters: 143 + ts - the TS context 144 . t - current timestep 145 - x - input vector 146 147 Output Parameters: 148 + A - Jacobian matrix 149 . B - optional preconditioning matrix 150 - flag - flag indicating matrix structure 151 152 Notes: 153 Most users should not need to explicitly call this routine, as it 154 is used internally within the nonlinear solvers. 155 156 See KSPSetOperators() for important information about setting the 157 flag parameter. 158 159 Level: developer 160 161 .keywords: SNES, compute, Jacobian, matrix 162 163 .seealso: TSSetRHSJacobian(), KSPSetOperators() 164 @*/ 165 PetscErrorCode TSComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg) 166 { 167 PetscErrorCode ierr; 168 PetscInt Xstate; 169 170 PetscFunctionBegin; 171 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 172 PetscValidHeaderSpecific(X,VEC_CLASSID,3); 173 PetscCheckSameComm(ts,1,X,3); 174 ierr = PetscObjectStateQuery((PetscObject)X,&Xstate);CHKERRQ(ierr); 175 if (ts->rhsjacobian.time == t && (ts->problem_type == TS_LINEAR || (ts->rhsjacobian.X == X && ts->rhsjacobian.Xstate == Xstate))) { 176 *flg = ts->rhsjacobian.mstructure; 177 PetscFunctionReturn(0); 178 } 179 180 if (!ts->userops->rhsjacobian && !ts->userops->ijacobian) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSJacobian() and / or TSSetIJacobian()"); 181 182 if (ts->userops->rhsjacobian) { 183 ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 184 *flg = DIFFERENT_NONZERO_PATTERN; 185 PetscStackPush("TS user Jacobian function"); 186 ierr = (*ts->userops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr); 187 PetscStackPop; 188 ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 189 /* make sure user returned a correct Jacobian and preconditioner */ 190 PetscValidHeaderSpecific(*A,MAT_CLASSID,4); 191 PetscValidHeaderSpecific(*B,MAT_CLASSID,5); 192 } else { 193 ierr = MatZeroEntries(*A);CHKERRQ(ierr); 194 if (*A != *B) {ierr = MatZeroEntries(*B);CHKERRQ(ierr);} 195 *flg = SAME_NONZERO_PATTERN; 196 } 197 ts->rhsjacobian.time = t; 198 ts->rhsjacobian.X = X; 199 ierr = PetscObjectStateQuery((PetscObject)X,&ts->rhsjacobian.Xstate);CHKERRQ(ierr); 200 ts->rhsjacobian.mstructure = *flg; 201 PetscFunctionReturn(0); 202 } 203 204 #undef __FUNCT__ 205 #define __FUNCT__ "TSComputeRHSFunction" 206 /*@ 207 TSComputeRHSFunction - Evaluates the right-hand-side function. 208 209 Collective on TS and Vec 210 211 Input Parameters: 212 + ts - the TS context 213 . t - current time 214 - x - state vector 215 216 Output Parameter: 217 . y - right hand side 218 219 Note: 220 Most users should not need to explicitly call this routine, as it 221 is used internally within the nonlinear solvers. 222 223 Level: developer 224 225 .keywords: TS, compute 226 227 .seealso: TSSetRHSFunction(), TSComputeIFunction() 228 @*/ 229 PetscErrorCode TSComputeRHSFunction(TS ts,PetscReal t,Vec x,Vec y) 230 { 231 PetscErrorCode ierr; 232 233 PetscFunctionBegin; 234 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 235 PetscValidHeaderSpecific(x,VEC_CLASSID,3); 236 PetscValidHeaderSpecific(y,VEC_CLASSID,4); 237 238 if (!ts->userops->rhsfunction && !ts->userops->ifunction) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSFunction() and / or TSSetIFunction()"); 239 240 ierr = PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr); 241 if (ts->userops->rhsfunction) { 242 PetscStackPush("TS user right-hand-side function"); 243 ierr = (*ts->userops->rhsfunction)(ts,t,x,y,ts->funP);CHKERRQ(ierr); 244 PetscStackPop; 245 } else { 246 ierr = VecZeroEntries(y);CHKERRQ(ierr); 247 } 248 249 ierr = PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr); 250 PetscFunctionReturn(0); 251 } 252 253 #undef __FUNCT__ 254 #define __FUNCT__ "TSGetRHSVec_Private" 255 static PetscErrorCode TSGetRHSVec_Private(TS ts,Vec *Frhs) 256 { 257 Vec F; 258 PetscErrorCode ierr; 259 260 PetscFunctionBegin; 261 *Frhs = PETSC_NULL; 262 ierr = TSGetIFunction(ts,&F,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 263 if (!ts->Frhs) { 264 ierr = VecDuplicate(F,&ts->Frhs);CHKERRQ(ierr); 265 } 266 *Frhs = ts->Frhs; 267 PetscFunctionReturn(0); 268 } 269 270 #undef __FUNCT__ 271 #define __FUNCT__ "TSGetRHSMats_Private" 272 static PetscErrorCode TSGetRHSMats_Private(TS ts,Mat *Arhs,Mat *Brhs) 273 { 274 Mat A,B; 275 PetscErrorCode ierr; 276 277 PetscFunctionBegin; 278 ierr = TSGetIJacobian(ts,&A,&B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 279 if (Arhs) { 280 if (!ts->Arhs) { 281 ierr = MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&ts->Arhs);CHKERRQ(ierr); 282 } 283 *Arhs = ts->Arhs; 284 } 285 if (Brhs) { 286 if (!ts->Brhs) { 287 ierr = MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&ts->Brhs);CHKERRQ(ierr); 288 } 289 *Brhs = ts->Brhs; 290 } 291 PetscFunctionReturn(0); 292 } 293 294 #undef __FUNCT__ 295 #define __FUNCT__ "TSComputeIFunction" 296 /*@ 297 TSComputeIFunction - Evaluates the DAE residual written in implicit form F(t,X,Xdot)=0 298 299 Collective on TS and Vec 300 301 Input Parameters: 302 + ts - the TS context 303 . t - current time 304 . X - state vector 305 . Xdot - time derivative of state vector 306 - imex - flag indicates if the method is IMEX so that the RHSFunction should be kept separate 307 308 Output Parameter: 309 . Y - right hand side 310 311 Note: 312 Most users should not need to explicitly call this routine, as it 313 is used internally within the nonlinear solvers. 314 315 If the user did did not write their equations in implicit form, this 316 function recasts them in implicit form. 317 318 Level: developer 319 320 .keywords: TS, compute 321 322 .seealso: TSSetIFunction(), TSComputeRHSFunction() 323 @*/ 324 PetscErrorCode TSComputeIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec Y,PetscBool imex) 325 { 326 PetscErrorCode ierr; 327 328 PetscFunctionBegin; 329 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 330 PetscValidHeaderSpecific(X,VEC_CLASSID,3); 331 PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4); 332 PetscValidHeaderSpecific(Y,VEC_CLASSID,5); 333 334 if (!ts->userops->rhsfunction && !ts->userops->ifunction) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSFunction() and / or TSSetIFunction()"); 335 336 ierr = PetscLogEventBegin(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr); 337 if (ts->userops->ifunction) { 338 PetscStackPush("TS user implicit function"); 339 ierr = (*ts->userops->ifunction)(ts,t,X,Xdot,Y,ts->funP);CHKERRQ(ierr); 340 PetscStackPop; 341 } 342 if (imex) { 343 if (!ts->userops->ifunction) { 344 ierr = VecCopy(Xdot,Y);CHKERRQ(ierr); 345 } 346 } else if (ts->userops->rhsfunction) { 347 if (ts->userops->ifunction) { 348 Vec Frhs; 349 ierr = TSGetRHSVec_Private(ts,&Frhs);CHKERRQ(ierr); 350 ierr = TSComputeRHSFunction(ts,t,X,Frhs);CHKERRQ(ierr); 351 ierr = VecAXPY(Y,-1,Frhs);CHKERRQ(ierr); 352 } else { 353 ierr = TSComputeRHSFunction(ts,t,X,Y);CHKERRQ(ierr); 354 ierr = VecAYPX(Y,-1,Xdot);CHKERRQ(ierr); 355 } 356 } 357 ierr = PetscLogEventEnd(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr); 358 PetscFunctionReturn(0); 359 } 360 361 #undef __FUNCT__ 362 #define __FUNCT__ "TSComputeIJacobian" 363 /*@ 364 TSComputeIJacobian - Evaluates the Jacobian of the DAE 365 366 Collective on TS and Vec 367 368 Input 369 Input Parameters: 370 + ts - the TS context 371 . t - current timestep 372 . X - state vector 373 . Xdot - time derivative of state vector 374 . shift - shift to apply, see note below 375 - imex - flag indicates if the method is IMEX so that the RHSJacobian should be kept separate 376 377 Output Parameters: 378 + A - Jacobian matrix 379 . B - optional preconditioning matrix 380 - flag - flag indicating matrix structure 381 382 Notes: 383 If F(t,X,Xdot)=0 is the DAE, the required Jacobian is 384 385 dF/dX + shift*dF/dXdot 386 387 Most users should not need to explicitly call this routine, as it 388 is used internally within the nonlinear solvers. 389 390 Level: developer 391 392 .keywords: TS, compute, Jacobian, matrix 393 394 .seealso: TSSetIJacobian() 395 @*/ 396 PetscErrorCode TSComputeIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg,PetscBool imex) 397 { 398 PetscInt Xstate, Xdotstate; 399 PetscErrorCode ierr; 400 401 PetscFunctionBegin; 402 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 403 PetscValidHeaderSpecific(X,VEC_CLASSID,3); 404 PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4); 405 PetscValidPointer(A,6); 406 PetscValidHeaderSpecific(*A,MAT_CLASSID,6); 407 PetscValidPointer(B,7); 408 PetscValidHeaderSpecific(*B,MAT_CLASSID,7); 409 PetscValidPointer(flg,8); 410 ierr = PetscObjectStateQuery((PetscObject)X,&Xstate);CHKERRQ(ierr); 411 ierr = PetscObjectStateQuery((PetscObject)Xdot,&Xdotstate);CHKERRQ(ierr); 412 if (ts->ijacobian.time == t && (ts->problem_type == TS_LINEAR || (ts->ijacobian.X == X && ts->ijacobian.Xstate == Xstate && ts->ijacobian.Xdot == Xdot && ts->ijacobian.Xdotstate == Xdotstate && ts->ijacobian.imex == imex))) { 413 *flg = ts->ijacobian.mstructure; 414 ierr = MatScale(*A, shift / ts->ijacobian.shift);CHKERRQ(ierr); 415 PetscFunctionReturn(0); 416 } 417 418 if (!ts->userops->rhsjacobian && !ts->userops->ijacobian) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSJacobian() and / or TSSetIJacobian()"); 419 420 *flg = SAME_NONZERO_PATTERN; /* In case we're solving a linear problem in which case it wouldn't get initialized below. */ 421 ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 422 if (ts->userops->ijacobian) { 423 *flg = DIFFERENT_NONZERO_PATTERN; 424 PetscStackPush("TS user implicit Jacobian"); 425 ierr = (*ts->userops->ijacobian)(ts,t,X,Xdot,shift,A,B,flg,ts->jacP);CHKERRQ(ierr); 426 PetscStackPop; 427 /* make sure user returned a correct Jacobian and preconditioner */ 428 PetscValidHeaderSpecific(*A,MAT_CLASSID,4); 429 PetscValidHeaderSpecific(*B,MAT_CLASSID,5); 430 } 431 if (imex) { 432 if (!ts->userops->ijacobian) { /* system was written as Xdot = F(t,X) */ 433 ierr = MatZeroEntries(*A);CHKERRQ(ierr); 434 ierr = MatShift(*A,shift);CHKERRQ(ierr); 435 if (*A != *B) { 436 ierr = MatZeroEntries(*B);CHKERRQ(ierr); 437 ierr = MatShift(*B,shift);CHKERRQ(ierr); 438 } 439 *flg = SAME_PRECONDITIONER; 440 } 441 } else { 442 if (!ts->userops->ijacobian) { 443 ierr = TSComputeRHSJacobian(ts,t,X,A,B,flg);CHKERRQ(ierr); 444 ierr = MatScale(*A,-1);CHKERRQ(ierr); 445 ierr = MatShift(*A,shift);CHKERRQ(ierr); 446 if (*A != *B) { 447 ierr = MatScale(*B,-1);CHKERRQ(ierr); 448 ierr = MatShift(*B,shift);CHKERRQ(ierr); 449 } 450 } else if (ts->userops->rhsjacobian) { 451 Mat Arhs,Brhs; 452 MatStructure axpy,flg2 = DIFFERENT_NONZERO_PATTERN; 453 ierr = TSGetRHSMats_Private(ts,&Arhs,&Brhs);CHKERRQ(ierr); 454 ierr = TSComputeRHSJacobian(ts,t,X,&Arhs,&Brhs,&flg2);CHKERRQ(ierr); 455 axpy = (*flg == flg2) ? SAME_NONZERO_PATTERN : DIFFERENT_NONZERO_PATTERN; 456 ierr = MatAXPY(*A,-1,Arhs,axpy);CHKERRQ(ierr); 457 if (*A != *B) { 458 ierr = MatAXPY(*B,-1,Brhs,axpy);CHKERRQ(ierr); 459 } 460 *flg = PetscMin(*flg,flg2); 461 } 462 } 463 464 ts->ijacobian.time = t; 465 ts->ijacobian.X = X; 466 ts->ijacobian.Xdot = Xdot; 467 ierr = PetscObjectStateQuery((PetscObject)X,&ts->ijacobian.Xstate);CHKERRQ(ierr); 468 ierr = PetscObjectStateQuery((PetscObject)Xdot,&ts->ijacobian.Xdotstate);CHKERRQ(ierr); 469 ts->ijacobian.shift = shift; 470 ts->ijacobian.imex = imex; 471 ts->ijacobian.mstructure = *flg; 472 ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 473 PetscFunctionReturn(0); 474 } 475 476 #undef __FUNCT__ 477 #define __FUNCT__ "TSSetRHSFunction" 478 /*@C 479 TSSetRHSFunction - Sets the routine for evaluating the function, 480 F(t,u), where U_t = F(t,u). 481 482 Logically Collective on TS 483 484 Input Parameters: 485 + ts - the TS context obtained from TSCreate() 486 . r - vector to put the computed right hand side (or PETSC_NULL to have it created) 487 . f - routine for evaluating the right-hand-side function 488 - ctx - [optional] user-defined context for private data for the 489 function evaluation routine (may be PETSC_NULL) 490 491 Calling sequence of func: 492 $ func (TS ts,PetscReal t,Vec u,Vec F,void *ctx); 493 494 + t - current timestep 495 . u - input vector 496 . F - function vector 497 - ctx - [optional] user-defined function context 498 499 Important: 500 The user MUST call either this routine or TSSetMatrices(). 501 502 Level: beginner 503 504 .keywords: TS, timestep, set, right-hand-side, function 505 506 .seealso: TSSetMatrices() 507 @*/ 508 PetscErrorCode TSSetRHSFunction(TS ts,Vec r,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx) 509 { 510 PetscErrorCode ierr; 511 SNES snes; 512 513 PetscFunctionBegin; 514 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 515 if (r) PetscValidHeaderSpecific(r,VEC_CLASSID,2); 516 if (f) ts->userops->rhsfunction = f; 517 if (ctx) ts->funP = ctx; 518 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 519 ierr = SNESSetFunction(snes,r,SNESTSFormFunction,ts);CHKERRQ(ierr); 520 PetscFunctionReturn(0); 521 } 522 523 #undef __FUNCT__ 524 #define __FUNCT__ "TSSetRHSJacobian" 525 /*@C 526 TSSetRHSJacobian - Sets the function to compute the Jacobian of F, 527 where U_t = F(U,t), as well as the location to store the matrix. 528 Use TSSetMatrices() for linear problems. 529 530 Logically Collective on TS 531 532 Input Parameters: 533 + ts - the TS context obtained from TSCreate() 534 . A - Jacobian matrix 535 . B - preconditioner matrix (usually same as A) 536 . f - the Jacobian evaluation routine 537 - ctx - [optional] user-defined context for private data for the 538 Jacobian evaluation routine (may be PETSC_NULL) 539 540 Calling sequence of func: 541 $ func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx); 542 543 + t - current timestep 544 . u - input vector 545 . A - matrix A, where U_t = A(t)u 546 . B - preconditioner matrix, usually the same as A 547 . flag - flag indicating information about the preconditioner matrix 548 structure (same as flag in KSPSetOperators()) 549 - ctx - [optional] user-defined context for matrix evaluation routine 550 551 Notes: 552 See KSPSetOperators() for important information about setting the flag 553 output parameter in the routine func(). Be sure to read this information! 554 555 The routine func() takes Mat * as the matrix arguments rather than Mat. 556 This allows the matrix evaluation routine to replace A and/or B with a 557 completely new matrix structure (not just different matrix elements) 558 when appropriate, for instance, if the nonzero structure is changing 559 throughout the global iterations. 560 561 Level: beginner 562 563 .keywords: TS, timestep, set, right-hand-side, Jacobian 564 565 .seealso: TSDefaultComputeJacobianColor(), 566 SNESDefaultComputeJacobianColor(), TSSetRHSFunction(), TSSetMatrices() 567 568 @*/ 569 PetscErrorCode TSSetRHSJacobian(TS ts,Mat A,Mat B,TSRHSJacobian f,void *ctx) 570 { 571 PetscErrorCode ierr; 572 SNES snes; 573 574 PetscFunctionBegin; 575 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 576 if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2); 577 if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3); 578 if (A) PetscCheckSameComm(ts,1,A,2); 579 if (B) PetscCheckSameComm(ts,1,B,3); 580 581 if (f) ts->userops->rhsjacobian = f; 582 if (ctx) ts->jacP = ctx; 583 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 584 if (!ts->userops->ijacobian) { 585 ierr = SNESSetJacobian(snes,A,B,SNESTSFormJacobian,ts);CHKERRQ(ierr); 586 } 587 if (A) { 588 ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 589 ierr = MatDestroy(&ts->Arhs);CHKERRQ(ierr); 590 ts->Arhs = A; 591 } 592 if (B) { 593 ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 594 ierr = MatDestroy(&ts->Brhs);CHKERRQ(ierr); 595 ts->Brhs = B; 596 } 597 PetscFunctionReturn(0); 598 } 599 600 601 #undef __FUNCT__ 602 #define __FUNCT__ "TSSetIFunction" 603 /*@C 604 TSSetIFunction - Set the function to compute F(t,U,U_t) where F = 0 is the DAE to be solved. 605 606 Logically Collective on TS 607 608 Input Parameters: 609 + ts - the TS context obtained from TSCreate() 610 . r - vector to hold the residual (or PETSC_NULL to have it created internally) 611 . f - the function evaluation routine 612 - ctx - user-defined context for private data for the function evaluation routine (may be PETSC_NULL) 613 614 Calling sequence of f: 615 $ f(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx); 616 617 + t - time at step/stage being solved 618 . u - state vector 619 . u_t - time derivative of state vector 620 . F - function vector 621 - ctx - [optional] user-defined context for matrix evaluation routine 622 623 Important: 624 The user MUST call either this routine, TSSetRHSFunction(), or TSSetMatrices(). This routine must be used when not solving an ODE. 625 626 Level: beginner 627 628 .keywords: TS, timestep, set, DAE, Jacobian 629 630 .seealso: TSSetMatrices(), TSSetRHSFunction(), TSSetIJacobian() 631 @*/ 632 PetscErrorCode TSSetIFunction(TS ts,Vec res,TSIFunction f,void *ctx) 633 { 634 PetscErrorCode ierr; 635 SNES snes; 636 637 PetscFunctionBegin; 638 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 639 if (res) PetscValidHeaderSpecific(res,VEC_CLASSID,2); 640 if (f) ts->userops->ifunction = f; 641 if (ctx) ts->funP = ctx; 642 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 643 ierr = SNESSetFunction(snes,res,SNESTSFormFunction,ts);CHKERRQ(ierr); 644 PetscFunctionReturn(0); 645 } 646 647 #undef __FUNCT__ 648 #define __FUNCT__ "TSGetIFunction" 649 /*@C 650 TSGetIFunction - Returns the vector where the implicit residual is stored and the function/contex to compute it. 651 652 Not Collective 653 654 Input Parameter: 655 . ts - the TS context 656 657 Output Parameter: 658 + r - vector to hold residual (or PETSC_NULL) 659 . func - the function to compute residual (or PETSC_NULL) 660 - ctx - the function context (or PETSC_NULL) 661 662 Level: advanced 663 664 .keywords: TS, nonlinear, get, function 665 666 .seealso: TSSetIFunction(), SNESGetFunction() 667 @*/ 668 PetscErrorCode TSGetIFunction(TS ts,Vec *r,TSIFunction *func,void **ctx) 669 { 670 PetscErrorCode ierr; 671 SNES snes; 672 673 PetscFunctionBegin; 674 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 675 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 676 ierr = SNESGetFunction(snes,r,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 677 if (func) *func = ts->userops->ifunction; 678 if (ctx) *ctx = ts->funP; 679 PetscFunctionReturn(0); 680 } 681 682 #undef __FUNCT__ 683 #define __FUNCT__ "TSGetRHSFunction" 684 /*@C 685 TSGetRHSFunction - Returns the vector where the right hand side is stored and the function/context to compute it. 686 687 Not Collective 688 689 Input Parameter: 690 . ts - the TS context 691 692 Output Parameter: 693 + r - vector to hold computed right hand side (or PETSC_NULL) 694 . func - the function to compute right hand side (or PETSC_NULL) 695 - ctx - the function context (or PETSC_NULL) 696 697 Level: advanced 698 699 .keywords: TS, nonlinear, get, function 700 701 .seealso: TSSetRhsfunction(), SNESGetFunction() 702 @*/ 703 PetscErrorCode TSGetRHSFunction(TS ts,Vec *r,TSRHSFunction *func,void **ctx) 704 { 705 PetscErrorCode ierr; 706 SNES snes; 707 708 PetscFunctionBegin; 709 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 710 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 711 ierr = SNESGetFunction(snes,r,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 712 if (func) *func = ts->userops->rhsfunction; 713 if (ctx) *ctx = ts->funP; 714 PetscFunctionReturn(0); 715 } 716 717 #undef __FUNCT__ 718 #define __FUNCT__ "TSSetIJacobian" 719 /*@C 720 TSSetIJacobian - Set the function to compute the matrix dF/dU + a*dF/dU_t where F(t,U,U_t) is the function 721 you provided with TSSetIFunction(). 722 723 Logically Collective on TS 724 725 Input Parameters: 726 + ts - the TS context obtained from TSCreate() 727 . A - Jacobian matrix 728 . B - preconditioning matrix for A (may be same as A) 729 . f - the Jacobian evaluation routine 730 - ctx - user-defined context for private data for the Jacobian evaluation routine (may be PETSC_NULL) 731 732 Calling sequence of f: 733 $ f(TS ts,PetscReal t,Vec U,Vec U_t,PetscReal a,Mat *A,Mat *B,MatStructure *flag,void *ctx); 734 735 + t - time at step/stage being solved 736 . U - state vector 737 . U_t - time derivative of state vector 738 . a - shift 739 . A - Jacobian of G(U) = F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t 740 . B - preconditioning matrix for A, may be same as A 741 . flag - flag indicating information about the preconditioner matrix 742 structure (same as flag in KSPSetOperators()) 743 - ctx - [optional] user-defined context for matrix evaluation routine 744 745 Notes: 746 The matrices A and B are exactly the matrices that are used by SNES for the nonlinear solve. 747 748 The matrix dF/dU + a*dF/dU_t you provide turns out to be 749 the Jacobian of G(U) = F(t,U,W+a*U) where F(t,U,U_t) = 0 is the DAE to be solved. 750 The time integrator internally approximates U_t by W+a*U where the positive "shift" 751 a and vector W depend on the integration method, step size, and past states. For example with 752 the backward Euler method a = 1/dt and W = -a*U(previous timestep) so 753 W + a*U = a*(U - U(previous timestep)) = (U - U(previous timestep))/dt 754 755 Level: beginner 756 757 .keywords: TS, timestep, DAE, Jacobian 758 759 .seealso: TSSetIFunction(), TSSetRHSJacobian() 760 761 @*/ 762 PetscErrorCode TSSetIJacobian(TS ts,Mat A,Mat B,TSIJacobian f,void *ctx) 763 { 764 PetscErrorCode ierr; 765 SNES snes; 766 767 PetscFunctionBegin; 768 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 769 if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2); 770 if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3); 771 if (A) PetscCheckSameComm(ts,1,A,2); 772 if (B) PetscCheckSameComm(ts,1,B,3); 773 if (f) ts->userops->ijacobian = f; 774 if (ctx) ts->jacP = ctx; 775 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 776 ierr = SNESSetJacobian(snes,A,B,SNESTSFormJacobian,ts);CHKERRQ(ierr); 777 PetscFunctionReturn(0); 778 } 779 780 #undef __FUNCT__ 781 #define __FUNCT__ "TSView" 782 /*@C 783 TSView - Prints the TS data structure. 784 785 Collective on TS 786 787 Input Parameters: 788 + ts - the TS context obtained from TSCreate() 789 - viewer - visualization context 790 791 Options Database Key: 792 . -ts_view - calls TSView() at end of TSStep() 793 794 Notes: 795 The available visualization contexts include 796 + PETSC_VIEWER_STDOUT_SELF - standard output (default) 797 - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 798 output where only the first processor opens 799 the file. All other processors send their 800 data to the first processor to print. 801 802 The user can open an alternative visualization context with 803 PetscViewerASCIIOpen() - output to a specified file. 804 805 Level: beginner 806 807 .keywords: TS, timestep, view 808 809 .seealso: PetscViewerASCIIOpen() 810 @*/ 811 PetscErrorCode TSView(TS ts,PetscViewer viewer) 812 { 813 PetscErrorCode ierr; 814 const TSType type; 815 PetscBool iascii,isstring,isundials; 816 817 PetscFunctionBegin; 818 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 819 if (!viewer) { 820 ierr = PetscViewerASCIIGetStdout(((PetscObject)ts)->comm,&viewer);CHKERRQ(ierr); 821 } 822 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 823 PetscCheckSameComm(ts,1,viewer,2); 824 825 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 826 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 827 if (iascii) { 828 ierr = PetscObjectPrintClassNamePrefixType((PetscObject)ts,viewer,"TS Object");CHKERRQ(ierr); 829 ierr = PetscViewerASCIIPrintf(viewer," maximum steps=%D\n",ts->max_steps);CHKERRQ(ierr); 830 ierr = PetscViewerASCIIPrintf(viewer," maximum time=%G\n",ts->max_time);CHKERRQ(ierr); 831 if (ts->problem_type == TS_NONLINEAR) { 832 ierr = PetscViewerASCIIPrintf(viewer," total number of nonlinear solver iterations=%D\n",ts->nonlinear_its);CHKERRQ(ierr); 833 ierr = PetscViewerASCIIPrintf(viewer," total number of nonlinear solve failures=%D\n",ts->num_snes_failures);CHKERRQ(ierr); 834 } 835 ierr = PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%D\n",ts->linear_its);CHKERRQ(ierr); 836 ierr = PetscViewerASCIIPrintf(viewer," total number of rejected steps=%D\n",ts->reject);CHKERRQ(ierr); 837 if (ts->ops->view) { 838 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 839 ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr); 840 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 841 } 842 } else if (isstring) { 843 ierr = TSGetType(ts,&type);CHKERRQ(ierr); 844 ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr); 845 } 846 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 847 ierr = PetscTypeCompare((PetscObject)ts,TSSUNDIALS,&isundials);CHKERRQ(ierr); 848 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 849 PetscFunctionReturn(0); 850 } 851 852 853 #undef __FUNCT__ 854 #define __FUNCT__ "TSSetApplicationContext" 855 /*@ 856 TSSetApplicationContext - Sets an optional user-defined context for 857 the timesteppers. 858 859 Logically Collective on TS 860 861 Input Parameters: 862 + ts - the TS context obtained from TSCreate() 863 - usrP - optional user context 864 865 Level: intermediate 866 867 .keywords: TS, timestep, set, application, context 868 869 .seealso: TSGetApplicationContext() 870 @*/ 871 PetscErrorCode TSSetApplicationContext(TS ts,void *usrP) 872 { 873 PetscFunctionBegin; 874 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 875 ts->user = usrP; 876 PetscFunctionReturn(0); 877 } 878 879 #undef __FUNCT__ 880 #define __FUNCT__ "TSGetApplicationContext" 881 /*@ 882 TSGetApplicationContext - Gets the user-defined context for the 883 timestepper. 884 885 Not Collective 886 887 Input Parameter: 888 . ts - the TS context obtained from TSCreate() 889 890 Output Parameter: 891 . usrP - user context 892 893 Level: intermediate 894 895 .keywords: TS, timestep, get, application, context 896 897 .seealso: TSSetApplicationContext() 898 @*/ 899 PetscErrorCode TSGetApplicationContext(TS ts,void *usrP) 900 { 901 PetscFunctionBegin; 902 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 903 *(void**)usrP = ts->user; 904 PetscFunctionReturn(0); 905 } 906 907 #undef __FUNCT__ 908 #define __FUNCT__ "TSGetTimeStepNumber" 909 /*@ 910 TSGetTimeStepNumber - Gets the current number of timesteps. 911 912 Not Collective 913 914 Input Parameter: 915 . ts - the TS context obtained from TSCreate() 916 917 Output Parameter: 918 . iter - number steps so far 919 920 Level: intermediate 921 922 .keywords: TS, timestep, get, iteration, number 923 @*/ 924 PetscErrorCode TSGetTimeStepNumber(TS ts,PetscInt* iter) 925 { 926 PetscFunctionBegin; 927 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 928 PetscValidIntPointer(iter,2); 929 *iter = ts->steps; 930 PetscFunctionReturn(0); 931 } 932 933 #undef __FUNCT__ 934 #define __FUNCT__ "TSSetInitialTimeStep" 935 /*@ 936 TSSetInitialTimeStep - Sets the initial timestep to be used, 937 as well as the initial time. 938 939 Logically Collective on TS 940 941 Input Parameters: 942 + ts - the TS context obtained from TSCreate() 943 . initial_time - the initial time 944 - time_step - the size of the timestep 945 946 Level: intermediate 947 948 .seealso: TSSetTimeStep(), TSGetTimeStep() 949 950 .keywords: TS, set, initial, timestep 951 @*/ 952 PetscErrorCode TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step) 953 { 954 PetscErrorCode ierr; 955 956 PetscFunctionBegin; 957 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 958 ierr = TSSetTimeStep(ts,time_step);CHKERRQ(ierr); 959 ierr = TSSetTime(ts,initial_time);CHKERRQ(ierr); 960 PetscFunctionReturn(0); 961 } 962 963 #undef __FUNCT__ 964 #define __FUNCT__ "TSSetTimeStep" 965 /*@ 966 TSSetTimeStep - Allows one to reset the timestep at any time, 967 useful for simple pseudo-timestepping codes. 968 969 Logically Collective on TS 970 971 Input Parameters: 972 + ts - the TS context obtained from TSCreate() 973 - time_step - the size of the timestep 974 975 Level: intermediate 976 977 .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 978 979 .keywords: TS, set, timestep 980 @*/ 981 PetscErrorCode TSSetTimeStep(TS ts,PetscReal time_step) 982 { 983 PetscFunctionBegin; 984 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 985 PetscValidLogicalCollectiveReal(ts,time_step,2); 986 ts->time_step = time_step; 987 PetscFunctionReturn(0); 988 } 989 990 #undef __FUNCT__ 991 #define __FUNCT__ "TSSetExactFinalTime" 992 /*@ 993 TSSetExactFinalTime - Determines whether to interpolate solution to the 994 exact final time requested by the user or just returns it at the final time 995 it computed. 996 997 Logically Collective on TS 998 999 Input Parameter: 1000 + ts - the time-step context 1001 - ft - PETSC_TRUE if interpolates, else PETSC_FALSE 1002 1003 Level: beginner 1004 1005 .seealso: TSSetDuration() 1006 @*/ 1007 PetscErrorCode TSSetExactFinalTime(TS ts,PetscBool flg) 1008 { 1009 1010 PetscFunctionBegin; 1011 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1012 PetscValidLogicalCollectiveBool(ts,flg,2); 1013 ts->exact_final_time = flg; 1014 PetscFunctionReturn(0); 1015 } 1016 1017 #undef __FUNCT__ 1018 #define __FUNCT__ "TSGetTimeStep" 1019 /*@ 1020 TSGetTimeStep - Gets the current timestep size. 1021 1022 Not Collective 1023 1024 Input Parameter: 1025 . ts - the TS context obtained from TSCreate() 1026 1027 Output Parameter: 1028 . dt - the current timestep size 1029 1030 Level: intermediate 1031 1032 .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 1033 1034 .keywords: TS, get, timestep 1035 @*/ 1036 PetscErrorCode TSGetTimeStep(TS ts,PetscReal* dt) 1037 { 1038 PetscFunctionBegin; 1039 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1040 PetscValidDoublePointer(dt,2); 1041 *dt = ts->time_step; 1042 PetscFunctionReturn(0); 1043 } 1044 1045 #undef __FUNCT__ 1046 #define __FUNCT__ "TSGetSolution" 1047 /*@ 1048 TSGetSolution - Returns the solution at the present timestep. It 1049 is valid to call this routine inside the function that you are evaluating 1050 in order to move to the new timestep. This vector not changed until 1051 the solution at the next timestep has been calculated. 1052 1053 Not Collective, but Vec returned is parallel if TS is parallel 1054 1055 Input Parameter: 1056 . ts - the TS context obtained from TSCreate() 1057 1058 Output Parameter: 1059 . v - the vector containing the solution 1060 1061 Level: intermediate 1062 1063 .seealso: TSGetTimeStep() 1064 1065 .keywords: TS, timestep, get, solution 1066 @*/ 1067 PetscErrorCode TSGetSolution(TS ts,Vec *v) 1068 { 1069 PetscFunctionBegin; 1070 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1071 PetscValidPointer(v,2); 1072 *v = ts->vec_sol; 1073 PetscFunctionReturn(0); 1074 } 1075 1076 /* ----- Routines to initialize and destroy a timestepper ---- */ 1077 #undef __FUNCT__ 1078 #define __FUNCT__ "TSSetProblemType" 1079 /*@ 1080 TSSetProblemType - Sets the type of problem to be solved. 1081 1082 Not collective 1083 1084 Input Parameters: 1085 + ts - The TS 1086 - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms 1087 .vb 1088 U_t = A U 1089 U_t = A(t) U 1090 U_t = F(t,U) 1091 .ve 1092 1093 Level: beginner 1094 1095 .keywords: TS, problem type 1096 .seealso: TSSetUp(), TSProblemType, TS 1097 @*/ 1098 PetscErrorCode TSSetProblemType(TS ts, TSProblemType type) 1099 { 1100 PetscErrorCode ierr; 1101 1102 PetscFunctionBegin; 1103 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1104 ts->problem_type = type; 1105 if (type == TS_LINEAR) { 1106 SNES snes; 1107 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 1108 ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr); 1109 } 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 M U_t = A U 1127 M(t) 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 if (ts->exact_final_time == PETSC_DECIDE) ts->exact_final_time = PETSC_FALSE; 1181 1182 if (!ts->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first"); 1183 1184 if (ts->ops->setup) { 1185 ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr); 1186 } 1187 1188 ts->setupcalled = PETSC_TRUE; 1189 PetscFunctionReturn(0); 1190 } 1191 1192 #undef __FUNCT__ 1193 #define __FUNCT__ "TSReset" 1194 /*@ 1195 TSReset - Resets a TS context and removes any allocated Vecs and Mats. 1196 1197 Collective on TS 1198 1199 Input Parameter: 1200 . ts - the TS context obtained from TSCreate() 1201 1202 Level: beginner 1203 1204 .keywords: TS, timestep, reset 1205 1206 .seealso: TSCreate(), TSSetup(), TSDestroy() 1207 @*/ 1208 PetscErrorCode TSReset(TS ts) 1209 { 1210 PetscErrorCode ierr; 1211 1212 PetscFunctionBegin; 1213 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1214 if (ts->ops->reset) { 1215 ierr = (*ts->ops->reset)(ts);CHKERRQ(ierr); 1216 } 1217 if (ts->snes) {ierr = SNESReset(ts->snes);CHKERRQ(ierr);} 1218 ierr = MatDestroy(&ts->Arhs);CHKERRQ(ierr); 1219 ierr = MatDestroy(&ts->Brhs);CHKERRQ(ierr); 1220 ierr = VecDestroy(&ts->Frhs);CHKERRQ(ierr); 1221 ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr); 1222 ierr = VecDestroyVecs(ts->nwork,&ts->work);CHKERRQ(ierr); 1223 ts->setupcalled = PETSC_FALSE; 1224 PetscFunctionReturn(0); 1225 } 1226 1227 #undef __FUNCT__ 1228 #define __FUNCT__ "TSDestroy" 1229 /*@ 1230 TSDestroy - Destroys the timestepper context that was created 1231 with TSCreate(). 1232 1233 Collective on TS 1234 1235 Input Parameter: 1236 . ts - the TS context obtained from TSCreate() 1237 1238 Level: beginner 1239 1240 .keywords: TS, timestepper, destroy 1241 1242 .seealso: TSCreate(), TSSetUp(), TSSolve() 1243 @*/ 1244 PetscErrorCode TSDestroy(TS *ts) 1245 { 1246 PetscErrorCode ierr; 1247 1248 PetscFunctionBegin; 1249 if (!*ts) PetscFunctionReturn(0); 1250 PetscValidHeaderSpecific((*ts),TS_CLASSID,1); 1251 if (--((PetscObject)(*ts))->refct > 0) {*ts = 0; PetscFunctionReturn(0);} 1252 1253 ierr = TSReset((*ts));CHKERRQ(ierr); 1254 1255 /* if memory was published with AMS then destroy it */ 1256 ierr = PetscObjectDepublish((*ts));CHKERRQ(ierr); 1257 if ((*ts)->ops->destroy) {ierr = (*(*ts)->ops->destroy)((*ts));CHKERRQ(ierr);} 1258 1259 ierr = SNESDestroy(&(*ts)->snes);CHKERRQ(ierr); 1260 ierr = DMDestroy(&(*ts)->dm);CHKERRQ(ierr); 1261 ierr = TSMonitorCancel((*ts));CHKERRQ(ierr); 1262 1263 ierr = PetscFree((*ts)->userops); 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 PetscErrorCode ierr; 1298 1299 PetscFunctionBegin; 1300 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1301 PetscValidPointer(snes,2); 1302 if (!ts->snes) { 1303 ierr = SNESCreate(((PetscObject)ts)->comm,&ts->snes);CHKERRQ(ierr); 1304 ierr = PetscLogObjectParent(ts,ts->snes);CHKERRQ(ierr); 1305 ierr = PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);CHKERRQ(ierr); 1306 if (ts->problem_type == TS_LINEAR) { 1307 ierr = SNESSetType(ts->snes,SNESKSPONLY);CHKERRQ(ierr); 1308 } 1309 } 1310 *snes = ts->snes; 1311 PetscFunctionReturn(0); 1312 } 1313 1314 #undef __FUNCT__ 1315 #define __FUNCT__ "TSGetKSP" 1316 /*@ 1317 TSGetKSP - Returns the KSP (linear solver) associated with 1318 a TS (timestepper) context. 1319 1320 Not Collective, but KSP is parallel if TS is parallel 1321 1322 Input Parameter: 1323 . ts - the TS context obtained from TSCreate() 1324 1325 Output Parameter: 1326 . ksp - the nonlinear solver context 1327 1328 Notes: 1329 The user can then directly manipulate the KSP context to set various 1330 options, etc. Likewise, the user can then extract and manipulate the 1331 KSP and PC contexts as well. 1332 1333 TSGetKSP() does not work for integrators that do not use KSP; 1334 in this case TSGetKSP() returns PETSC_NULL in ksp. 1335 1336 Level: beginner 1337 1338 .keywords: timestep, get, KSP 1339 @*/ 1340 PetscErrorCode TSGetKSP(TS ts,KSP *ksp) 1341 { 1342 PetscErrorCode ierr; 1343 SNES snes; 1344 1345 PetscFunctionBegin; 1346 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1347 PetscValidPointer(ksp,2); 1348 if (!((PetscObject)ts)->type_name) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"KSP is not created yet. Call TSSetType() first"); 1349 if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()"); 1350 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 1351 ierr = SNESGetKSP(snes,ksp);CHKERRQ(ierr); 1352 PetscFunctionReturn(0); 1353 } 1354 1355 /* ----------- Routines to set solver parameters ---------- */ 1356 1357 #undef __FUNCT__ 1358 #define __FUNCT__ "TSGetDuration" 1359 /*@ 1360 TSGetDuration - Gets the maximum number of timesteps to use and 1361 maximum time for iteration. 1362 1363 Not Collective 1364 1365 Input Parameters: 1366 + ts - the TS context obtained from TSCreate() 1367 . maxsteps - maximum number of iterations to use, or PETSC_NULL 1368 - maxtime - final time to iterate to, or PETSC_NULL 1369 1370 Level: intermediate 1371 1372 .keywords: TS, timestep, get, maximum, iterations, time 1373 @*/ 1374 PetscErrorCode TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime) 1375 { 1376 PetscFunctionBegin; 1377 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1378 if (maxsteps) { 1379 PetscValidIntPointer(maxsteps,2); 1380 *maxsteps = ts->max_steps; 1381 } 1382 if (maxtime) { 1383 PetscValidScalarPointer(maxtime,3); 1384 *maxtime = ts->max_time; 1385 } 1386 PetscFunctionReturn(0); 1387 } 1388 1389 #undef __FUNCT__ 1390 #define __FUNCT__ "TSSetDuration" 1391 /*@ 1392 TSSetDuration - Sets the maximum number of timesteps to use and 1393 maximum time for iteration. 1394 1395 Logically Collective on TS 1396 1397 Input Parameters: 1398 + ts - the TS context obtained from TSCreate() 1399 . maxsteps - maximum number of iterations to use 1400 - maxtime - final time to iterate to 1401 1402 Options Database Keys: 1403 . -ts_max_steps <maxsteps> - Sets maxsteps 1404 . -ts_max_time <maxtime> - Sets maxtime 1405 1406 Notes: 1407 The default maximum number of iterations is 5000. Default time is 5.0 1408 1409 Level: intermediate 1410 1411 .keywords: TS, timestep, set, maximum, iterations 1412 1413 .seealso: TSSetExactFinalTime() 1414 @*/ 1415 PetscErrorCode TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime) 1416 { 1417 PetscFunctionBegin; 1418 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1419 PetscValidLogicalCollectiveInt(ts,maxsteps,2); 1420 PetscValidLogicalCollectiveReal(ts,maxtime,2); 1421 if (maxsteps >= 0) ts->max_steps = maxsteps; 1422 if (maxtime != PETSC_DEFAULT) ts->max_time = maxtime; 1423 PetscFunctionReturn(0); 1424 } 1425 1426 #undef __FUNCT__ 1427 #define __FUNCT__ "TSSetSolution" 1428 /*@ 1429 TSSetSolution - Sets the initial solution vector 1430 for use by the TS routines. 1431 1432 Logically Collective on TS and Vec 1433 1434 Input Parameters: 1435 + ts - the TS context obtained from TSCreate() 1436 - x - the solution vector 1437 1438 Level: beginner 1439 1440 .keywords: TS, timestep, set, solution, initial conditions 1441 @*/ 1442 PetscErrorCode TSSetSolution(TS ts,Vec x) 1443 { 1444 PetscErrorCode ierr; 1445 1446 PetscFunctionBegin; 1447 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1448 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 1449 ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr); 1450 ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr); 1451 ts->vec_sol = x; 1452 PetscFunctionReturn(0); 1453 } 1454 1455 #undef __FUNCT__ 1456 #define __FUNCT__ "TSSetPreStep" 1457 /*@C 1458 TSSetPreStep - Sets the general-purpose function 1459 called once at the beginning of each time step. 1460 1461 Logically Collective on TS 1462 1463 Input Parameters: 1464 + ts - The TS context obtained from TSCreate() 1465 - func - The function 1466 1467 Calling sequence of func: 1468 . func (TS ts); 1469 1470 Level: intermediate 1471 1472 .keywords: TS, timestep 1473 @*/ 1474 PetscErrorCode TSSetPreStep(TS ts, PetscErrorCode (*func)(TS)) 1475 { 1476 PetscFunctionBegin; 1477 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1478 ts->ops->prestep = func; 1479 PetscFunctionReturn(0); 1480 } 1481 1482 #undef __FUNCT__ 1483 #define __FUNCT__ "TSPreStep" 1484 /*@ 1485 TSPreStep - Runs the user-defined pre-step function. 1486 1487 Collective on TS 1488 1489 Input Parameters: 1490 . ts - The TS context obtained from TSCreate() 1491 1492 Notes: 1493 TSPreStep() is typically used within time stepping implementations, 1494 so most users would not generally call this routine themselves. 1495 1496 Level: developer 1497 1498 .keywords: TS, timestep 1499 @*/ 1500 PetscErrorCode TSPreStep(TS ts) 1501 { 1502 PetscErrorCode ierr; 1503 1504 PetscFunctionBegin; 1505 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1506 if (ts->ops->prestep) { 1507 PetscStackPush("TS PreStep function"); 1508 ierr = (*ts->ops->prestep)(ts);CHKERRQ(ierr); 1509 PetscStackPop; 1510 } 1511 PetscFunctionReturn(0); 1512 } 1513 1514 #undef __FUNCT__ 1515 #define __FUNCT__ "TSSetPostStep" 1516 /*@C 1517 TSSetPostStep - Sets the general-purpose function 1518 called once at the end of each time step. 1519 1520 Logically Collective on TS 1521 1522 Input Parameters: 1523 + ts - The TS context obtained from TSCreate() 1524 - func - The function 1525 1526 Calling sequence of func: 1527 . func (TS ts); 1528 1529 Level: intermediate 1530 1531 .keywords: TS, timestep 1532 @*/ 1533 PetscErrorCode TSSetPostStep(TS ts, PetscErrorCode (*func)(TS)) 1534 { 1535 PetscFunctionBegin; 1536 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1537 ts->ops->poststep = func; 1538 PetscFunctionReturn(0); 1539 } 1540 1541 #undef __FUNCT__ 1542 #define __FUNCT__ "TSPostStep" 1543 /*@ 1544 TSPostStep - Runs the user-defined post-step function. 1545 1546 Collective on TS 1547 1548 Input Parameters: 1549 . ts - The TS context obtained from TSCreate() 1550 1551 Notes: 1552 TSPostStep() is typically used within time stepping implementations, 1553 so most users would not generally call this routine themselves. 1554 1555 Level: developer 1556 1557 .keywords: TS, timestep 1558 @*/ 1559 PetscErrorCode TSPostStep(TS ts) 1560 { 1561 PetscErrorCode ierr; 1562 1563 PetscFunctionBegin; 1564 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1565 if (ts->ops->poststep) { 1566 PetscStackPush("TS PostStep function"); 1567 ierr = (*ts->ops->poststep)(ts);CHKERRQ(ierr); 1568 PetscStackPop; 1569 } 1570 PetscFunctionReturn(0); 1571 } 1572 1573 /* ------------ Routines to set performance monitoring options ----------- */ 1574 1575 #undef __FUNCT__ 1576 #define __FUNCT__ "TSMonitorSet" 1577 /*@C 1578 TSMonitorSet - Sets an ADDITIONAL function that is to be used at every 1579 timestep to display the iteration's progress. 1580 1581 Logically Collective on TS 1582 1583 Input Parameters: 1584 + ts - the TS context obtained from TSCreate() 1585 . monitor - monitoring routine 1586 . mctx - [optional] user-defined context for private data for the 1587 monitor routine (use PETSC_NULL if no context is desired) 1588 - monitordestroy - [optional] routine that frees monitor context 1589 (may be PETSC_NULL) 1590 1591 Calling sequence of monitor: 1592 $ int monitor(TS ts,PetscInt steps,PetscReal time,Vec x,void *mctx) 1593 1594 + ts - the TS context 1595 . steps - iteration number 1596 . time - current time 1597 . x - current iterate 1598 - mctx - [optional] monitoring context 1599 1600 Notes: 1601 This routine adds an additional monitor to the list of monitors that 1602 already has been loaded. 1603 1604 Fortran notes: Only a single monitor function can be set for each TS object 1605 1606 Level: intermediate 1607 1608 .keywords: TS, timestep, set, monitor 1609 1610 .seealso: TSMonitorDefault(), TSMonitorCancel() 1611 @*/ 1612 PetscErrorCode TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void**)) 1613 { 1614 PetscFunctionBegin; 1615 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1616 if (ts->numbermonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1617 ts->monitor[ts->numbermonitors] = monitor; 1618 ts->mdestroy[ts->numbermonitors] = mdestroy; 1619 ts->monitorcontext[ts->numbermonitors++] = (void*)mctx; 1620 PetscFunctionReturn(0); 1621 } 1622 1623 #undef __FUNCT__ 1624 #define __FUNCT__ "TSMonitorCancel" 1625 /*@C 1626 TSMonitorCancel - Clears all the monitors that have been set on a time-step object. 1627 1628 Logically Collective on TS 1629 1630 Input Parameters: 1631 . ts - the TS context obtained from TSCreate() 1632 1633 Notes: 1634 There is no way to remove a single, specific monitor. 1635 1636 Level: intermediate 1637 1638 .keywords: TS, timestep, set, monitor 1639 1640 .seealso: TSMonitorDefault(), TSMonitorSet() 1641 @*/ 1642 PetscErrorCode TSMonitorCancel(TS ts) 1643 { 1644 PetscErrorCode ierr; 1645 PetscInt i; 1646 1647 PetscFunctionBegin; 1648 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1649 for (i=0; i<ts->numbermonitors; i++) { 1650 if (ts->mdestroy[i]) { 1651 ierr = (*ts->mdestroy[i])(&ts->monitorcontext[i]);CHKERRQ(ierr); 1652 } 1653 } 1654 ts->numbermonitors = 0; 1655 PetscFunctionReturn(0); 1656 } 1657 1658 #undef __FUNCT__ 1659 #define __FUNCT__ "TSMonitorDefault" 1660 /*@ 1661 TSMonitorDefault - Sets the Default monitor 1662 1663 Level: intermediate 1664 1665 .keywords: TS, set, monitor 1666 1667 .seealso: TSMonitorDefault(), TSMonitorSet() 1668 @*/ 1669 PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy) 1670 { 1671 PetscErrorCode ierr; 1672 PetscViewer viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ts)->comm); 1673 1674 PetscFunctionBegin; 1675 ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 1676 ierr = PetscViewerASCIIPrintf(viewer,"%D TS dt %g time %g\n",step,(double)ts->time_step,(double)ptime);CHKERRQ(ierr); 1677 ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 1678 PetscFunctionReturn(0); 1679 } 1680 1681 #undef __FUNCT__ 1682 #define __FUNCT__ "TSSetRetainStages" 1683 /*@ 1684 TSSetRetainStages - Request that all stages in the upcoming step be stored so that interpolation will be available. 1685 1686 Logically Collective on TS 1687 1688 Input Argument: 1689 . ts - time stepping context 1690 1691 Output Argument: 1692 . flg - PETSC_TRUE or PETSC_FALSE 1693 1694 Level: intermediate 1695 1696 .keywords: TS, set 1697 1698 .seealso: TSInterpolate(), TSSetPostStep() 1699 @*/ 1700 PetscErrorCode TSSetRetainStages(TS ts,PetscBool flg) 1701 { 1702 1703 PetscFunctionBegin; 1704 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1705 ts->retain_stages = flg; 1706 PetscFunctionReturn(0); 1707 } 1708 1709 #undef __FUNCT__ 1710 #define __FUNCT__ "TSInterpolate" 1711 /*@ 1712 TSInterpolate - Interpolate the solution computed during the previous step to an arbitrary location in the interval 1713 1714 Collective on TS 1715 1716 Input Argument: 1717 + ts - time stepping context 1718 - t - time to interpolate to 1719 1720 Output Argument: 1721 . X - state at given time 1722 1723 Notes: 1724 The user should call TSSetRetainStages() before taking a step in which interpolation will be requested. 1725 1726 Level: intermediate 1727 1728 Developer Notes: 1729 TSInterpolate() and the storing of previous steps/stages should be generalized to support delay differential equations and continuous adjoints. 1730 1731 .keywords: TS, set 1732 1733 .seealso: TSSetRetainStages(), TSSetPostStep() 1734 @*/ 1735 PetscErrorCode TSInterpolate(TS ts,PetscReal t,Vec X) 1736 { 1737 PetscErrorCode ierr; 1738 1739 PetscFunctionBegin; 1740 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1741 if (t < ts->ptime - ts->time_step_prev || t > ts->ptime) SETERRQ3(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Requested time %G not in last time steps [%G,%G]",t,ts->ptime-ts->time_step_prev,ts->ptime); 1742 if (!ts->ops->interpolate) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"%s does not provide interpolation",((PetscObject)ts)->type_name); 1743 ierr = (*ts->ops->interpolate)(ts,t,X);CHKERRQ(ierr); 1744 PetscFunctionReturn(0); 1745 } 1746 1747 #undef __FUNCT__ 1748 #define __FUNCT__ "TSStep" 1749 /*@ 1750 TSStep - Steps one time step 1751 1752 Collective on TS 1753 1754 Input Parameter: 1755 . ts - the TS context obtained from TSCreate() 1756 1757 Level: intermediate 1758 1759 .keywords: TS, timestep, solve 1760 1761 .seealso: TSCreate(), TSSetUp(), TSDestroy(), TSSolve() 1762 @*/ 1763 PetscErrorCode TSStep(TS ts) 1764 { 1765 PetscReal ptime_prev; 1766 PetscErrorCode ierr; 1767 1768 PetscFunctionBegin; 1769 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1770 ierr = TSSetUp(ts);CHKERRQ(ierr); 1771 1772 ts->reason = TS_CONVERGED_ITERATING; 1773 1774 ptime_prev = ts->ptime; 1775 ierr = PetscLogEventBegin(TS_Step,ts,0,0,0);CHKERRQ(ierr); 1776 ierr = (*ts->ops->step)(ts);CHKERRQ(ierr); 1777 ierr = PetscLogEventEnd(TS_Step,ts,0,0,0);CHKERRQ(ierr); 1778 ts->time_step_prev = ts->ptime - ptime_prev; 1779 1780 if (ts->reason < 0) { 1781 if (ts->errorifstepfailed) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_NOT_CONVERGED,"TSStep has failed"); 1782 } else if (!ts->reason) { 1783 if (ts->steps >= ts->max_steps) 1784 ts->reason = TS_CONVERGED_ITS; 1785 else if (ts->ptime >= ts->max_time) 1786 ts->reason = TS_CONVERGED_TIME; 1787 } 1788 1789 PetscFunctionReturn(0); 1790 } 1791 1792 #undef __FUNCT__ 1793 #define __FUNCT__ "TSSolve" 1794 /*@ 1795 TSSolve - Steps the requested number of timesteps. 1796 1797 Collective on TS 1798 1799 Input Parameter: 1800 + ts - the TS context obtained from TSCreate() 1801 - x - the solution vector 1802 1803 Output Parameter: 1804 . ftime - time of the state vector x upon completion 1805 1806 Level: beginner 1807 1808 Notes: 1809 The final time returned by this function may be different from the time of the internally 1810 held state accessible by TSGetSolution() and TSGetTime() because the method may have 1811 stepped over the final time. 1812 1813 .keywords: TS, timestep, solve 1814 1815 .seealso: TSCreate(), TSSetSolution(), TSStep() 1816 @*/ 1817 PetscErrorCode TSSolve(TS ts,Vec x,PetscReal *ftime) 1818 { 1819 PetscBool flg; 1820 char filename[PETSC_MAX_PATH_LEN]; 1821 PetscViewer viewer; 1822 PetscErrorCode ierr; 1823 1824 PetscFunctionBegin; 1825 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1826 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 1827 if (ts->exact_final_time) { /* Need ts->vec_sol to be distinct so it is not overwritten when we interpolate at the end */ 1828 if (!ts->vec_sol || x == ts->vec_sol) { 1829 Vec y; 1830 ierr = VecDuplicate(x,&y);CHKERRQ(ierr); 1831 ierr = TSSetSolution(ts,y);CHKERRQ(ierr); 1832 ierr = VecDestroy(&y);CHKERRQ(ierr); /* grant ownership */ 1833 } 1834 ierr = VecCopy(x,ts->vec_sol);CHKERRQ(ierr); 1835 } else { 1836 ierr = TSSetSolution(ts,x);CHKERRQ(ierr); 1837 } 1838 ierr = TSSetUp(ts);CHKERRQ(ierr); 1839 /* reset time step and iteration counters */ 1840 ts->steps = 0; 1841 ts->linear_its = 0; 1842 ts->nonlinear_its = 0; 1843 ts->num_snes_failures = 0; 1844 ts->reject = 0; 1845 ts->reason = TS_CONVERGED_ITERATING; 1846 1847 if (ts->ops->solve) { /* This private interface is transitional and should be removed when all implementations are updated. */ 1848 ierr = (*ts->ops->solve)(ts);CHKERRQ(ierr); 1849 ierr = VecCopy(ts->vec_sol,x);CHKERRQ(ierr); 1850 if (ftime) *ftime = ts->ptime; 1851 } else { 1852 /* steps the requested number of timesteps. */ 1853 ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 1854 if (ts->steps >= ts->max_steps) 1855 ts->reason = TS_CONVERGED_ITS; 1856 else if (ts->ptime >= ts->max_time) 1857 ts->reason = TS_CONVERGED_TIME; 1858 while (!ts->reason) { 1859 ierr = TSPreStep(ts);CHKERRQ(ierr); 1860 ierr = TSStep(ts);CHKERRQ(ierr); 1861 ierr = TSPostStep(ts);CHKERRQ(ierr); 1862 ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 1863 } 1864 if (ts->exact_final_time && ts->ptime > ts->max_time) { 1865 ierr = TSInterpolate(ts,ts->max_time,x);CHKERRQ(ierr); 1866 if (ftime) *ftime = ts->max_time; 1867 } else { 1868 ierr = VecCopy(ts->vec_sol,x);CHKERRQ(ierr); 1869 if (ftime) *ftime = ts->ptime; 1870 } 1871 } 1872 ierr = PetscOptionsGetString(((PetscObject)ts)->prefix,"-ts_view",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 1873 if (flg && !PetscPreLoadingOn) { 1874 ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,filename,&viewer);CHKERRQ(ierr); 1875 ierr = TSView(ts,viewer);CHKERRQ(ierr); 1876 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 1877 } 1878 PetscFunctionReturn(0); 1879 } 1880 1881 #undef __FUNCT__ 1882 #define __FUNCT__ "TSMonitor" 1883 /*@ 1884 TSMonitor - Runs all user-provided monitor routines set using TSMonitorSet() 1885 1886 Collective on TS 1887 1888 Input Parameters: 1889 + ts - time stepping context obtained from TSCreate() 1890 . step - step number that has just completed 1891 . ptime - model time of the state 1892 - x - state at the current model time 1893 1894 Notes: 1895 TSMonitor() is typically used within the time stepping implementations. 1896 Users might call this function when using the TSStep() interface instead of TSSolve(). 1897 1898 Level: advanced 1899 1900 .keywords: TS, timestep 1901 @*/ 1902 PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x) 1903 { 1904 PetscErrorCode ierr; 1905 PetscInt i,n = ts->numbermonitors; 1906 1907 PetscFunctionBegin; 1908 for (i=0; i<n; i++) { 1909 ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr); 1910 } 1911 PetscFunctionReturn(0); 1912 } 1913 1914 /* ------------------------------------------------------------------------*/ 1915 1916 #undef __FUNCT__ 1917 #define __FUNCT__ "TSMonitorLGCreate" 1918 /*@C 1919 TSMonitorLGCreate - Creates a line graph context for use with 1920 TS to monitor convergence of preconditioned residual norms. 1921 1922 Collective on TS 1923 1924 Input Parameters: 1925 + host - the X display to open, or null for the local machine 1926 . label - the title to put in the title bar 1927 . x, y - the screen coordinates of the upper left coordinate of the window 1928 - m, n - the screen width and height in pixels 1929 1930 Output Parameter: 1931 . draw - the drawing context 1932 1933 Options Database Key: 1934 . -ts_monitor_draw - automatically sets line graph monitor 1935 1936 Notes: 1937 Use TSMonitorLGDestroy() to destroy this line graph, not PetscDrawLGDestroy(). 1938 1939 Level: intermediate 1940 1941 .keywords: TS, monitor, line graph, residual, seealso 1942 1943 .seealso: TSMonitorLGDestroy(), TSMonitorSet() 1944 1945 @*/ 1946 PetscErrorCode TSMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1947 { 1948 PetscDraw win; 1949 PetscErrorCode ierr; 1950 1951 PetscFunctionBegin; 1952 ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr); 1953 ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr); 1954 ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr); 1955 ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr); 1956 1957 ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr); 1958 PetscFunctionReturn(0); 1959 } 1960 1961 #undef __FUNCT__ 1962 #define __FUNCT__ "TSMonitorLG" 1963 PetscErrorCode TSMonitorLG(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx) 1964 { 1965 PetscDrawLG lg = (PetscDrawLG) monctx; 1966 PetscReal x,y = ptime; 1967 PetscErrorCode ierr; 1968 1969 PetscFunctionBegin; 1970 if (!monctx) { 1971 MPI_Comm comm; 1972 PetscViewer viewer; 1973 1974 ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr); 1975 viewer = PETSC_VIEWER_DRAW_(comm); 1976 ierr = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr); 1977 } 1978 1979 if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 1980 x = (PetscReal)n; 1981 ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 1982 if (n < 20 || (n % 5)) { 1983 ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 1984 } 1985 PetscFunctionReturn(0); 1986 } 1987 1988 #undef __FUNCT__ 1989 #define __FUNCT__ "TSMonitorLGDestroy" 1990 /*@C 1991 TSMonitorLGDestroy - Destroys a line graph context that was created 1992 with TSMonitorLGCreate(). 1993 1994 Collective on PetscDrawLG 1995 1996 Input Parameter: 1997 . draw - the drawing context 1998 1999 Level: intermediate 2000 2001 .keywords: TS, monitor, line graph, destroy 2002 2003 .seealso: TSMonitorLGCreate(), TSMonitorSet(), TSMonitorLG(); 2004 @*/ 2005 PetscErrorCode TSMonitorLGDestroy(PetscDrawLG *drawlg) 2006 { 2007 PetscDraw draw; 2008 PetscErrorCode ierr; 2009 2010 PetscFunctionBegin; 2011 ierr = PetscDrawLGGetDraw(*drawlg,&draw);CHKERRQ(ierr); 2012 ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr); 2013 ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr); 2014 PetscFunctionReturn(0); 2015 } 2016 2017 #undef __FUNCT__ 2018 #define __FUNCT__ "TSGetTime" 2019 /*@ 2020 TSGetTime - Gets the current time. 2021 2022 Not Collective 2023 2024 Input Parameter: 2025 . ts - the TS context obtained from TSCreate() 2026 2027 Output Parameter: 2028 . t - the current time 2029 2030 Level: beginner 2031 2032 .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 2033 2034 .keywords: TS, get, time 2035 @*/ 2036 PetscErrorCode TSGetTime(TS ts,PetscReal* t) 2037 { 2038 PetscFunctionBegin; 2039 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2040 PetscValidDoublePointer(t,2); 2041 *t = ts->ptime; 2042 PetscFunctionReturn(0); 2043 } 2044 2045 #undef __FUNCT__ 2046 #define __FUNCT__ "TSSetTime" 2047 /*@ 2048 TSSetTime - Allows one to reset the time. 2049 2050 Logically Collective on TS 2051 2052 Input Parameters: 2053 + ts - the TS context obtained from TSCreate() 2054 - time - the time 2055 2056 Level: intermediate 2057 2058 .seealso: TSGetTime(), TSSetDuration() 2059 2060 .keywords: TS, set, time 2061 @*/ 2062 PetscErrorCode TSSetTime(TS ts, PetscReal t) 2063 { 2064 PetscFunctionBegin; 2065 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2066 PetscValidLogicalCollectiveReal(ts,t,2); 2067 ts->ptime = t; 2068 PetscFunctionReturn(0); 2069 } 2070 2071 #undef __FUNCT__ 2072 #define __FUNCT__ "TSSetOptionsPrefix" 2073 /*@C 2074 TSSetOptionsPrefix - Sets the prefix used for searching for all 2075 TS options in the database. 2076 2077 Logically Collective on TS 2078 2079 Input Parameter: 2080 + ts - The TS context 2081 - prefix - The prefix to prepend to all option names 2082 2083 Notes: 2084 A hyphen (-) must NOT be given at the beginning of the prefix name. 2085 The first character of all runtime options is AUTOMATICALLY the 2086 hyphen. 2087 2088 Level: advanced 2089 2090 .keywords: TS, set, options, prefix, database 2091 2092 .seealso: TSSetFromOptions() 2093 2094 @*/ 2095 PetscErrorCode TSSetOptionsPrefix(TS ts,const char prefix[]) 2096 { 2097 PetscErrorCode ierr; 2098 SNES snes; 2099 2100 PetscFunctionBegin; 2101 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2102 ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 2103 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2104 ierr = SNESSetOptionsPrefix(snes,prefix);CHKERRQ(ierr); 2105 PetscFunctionReturn(0); 2106 } 2107 2108 2109 #undef __FUNCT__ 2110 #define __FUNCT__ "TSAppendOptionsPrefix" 2111 /*@C 2112 TSAppendOptionsPrefix - Appends to the prefix used for searching for all 2113 TS options in the database. 2114 2115 Logically Collective on TS 2116 2117 Input Parameter: 2118 + ts - The TS context 2119 - prefix - The prefix to prepend to all option names 2120 2121 Notes: 2122 A hyphen (-) must NOT be given at the beginning of the prefix name. 2123 The first character of all runtime options is AUTOMATICALLY the 2124 hyphen. 2125 2126 Level: advanced 2127 2128 .keywords: TS, append, options, prefix, database 2129 2130 .seealso: TSGetOptionsPrefix() 2131 2132 @*/ 2133 PetscErrorCode TSAppendOptionsPrefix(TS ts,const char prefix[]) 2134 { 2135 PetscErrorCode ierr; 2136 SNES snes; 2137 2138 PetscFunctionBegin; 2139 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2140 ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 2141 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2142 ierr = SNESAppendOptionsPrefix(snes,prefix);CHKERRQ(ierr); 2143 PetscFunctionReturn(0); 2144 } 2145 2146 #undef __FUNCT__ 2147 #define __FUNCT__ "TSGetOptionsPrefix" 2148 /*@C 2149 TSGetOptionsPrefix - Sets the prefix used for searching for all 2150 TS options in the database. 2151 2152 Not Collective 2153 2154 Input Parameter: 2155 . ts - The TS context 2156 2157 Output Parameter: 2158 . prefix - A pointer to the prefix string used 2159 2160 Notes: On the fortran side, the user should pass in a string 'prifix' of 2161 sufficient length to hold the prefix. 2162 2163 Level: intermediate 2164 2165 .keywords: TS, get, options, prefix, database 2166 2167 .seealso: TSAppendOptionsPrefix() 2168 @*/ 2169 PetscErrorCode TSGetOptionsPrefix(TS ts,const char *prefix[]) 2170 { 2171 PetscErrorCode ierr; 2172 2173 PetscFunctionBegin; 2174 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2175 PetscValidPointer(prefix,2); 2176 ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 2177 PetscFunctionReturn(0); 2178 } 2179 2180 #undef __FUNCT__ 2181 #define __FUNCT__ "TSGetRHSJacobian" 2182 /*@C 2183 TSGetRHSJacobian - Returns the Jacobian J at the present timestep. 2184 2185 Not Collective, but parallel objects are returned if TS is parallel 2186 2187 Input Parameter: 2188 . ts - The TS context obtained from TSCreate() 2189 2190 Output Parameters: 2191 + J - The Jacobian J of F, where U_t = F(U,t) 2192 . M - The preconditioner matrix, usually the same as J 2193 . func - Function to compute the Jacobian of the RHS 2194 - ctx - User-defined context for Jacobian evaluation routine 2195 2196 Notes: You can pass in PETSC_NULL for any return argument you do not need. 2197 2198 Level: intermediate 2199 2200 .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber() 2201 2202 .keywords: TS, timestep, get, matrix, Jacobian 2203 @*/ 2204 PetscErrorCode TSGetRHSJacobian(TS ts,Mat *J,Mat *M,TSRHSJacobian *func,void **ctx) 2205 { 2206 PetscErrorCode ierr; 2207 SNES snes; 2208 2209 PetscFunctionBegin; 2210 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2211 ierr = SNESGetJacobian(snes,J,M,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 2212 if (func) *func = ts->userops->rhsjacobian; 2213 if (ctx) *ctx = ts->jacP; 2214 PetscFunctionReturn(0); 2215 } 2216 2217 #undef __FUNCT__ 2218 #define __FUNCT__ "TSGetIJacobian" 2219 /*@C 2220 TSGetIJacobian - Returns the implicit Jacobian at the present timestep. 2221 2222 Not Collective, but parallel objects are returned if TS is parallel 2223 2224 Input Parameter: 2225 . ts - The TS context obtained from TSCreate() 2226 2227 Output Parameters: 2228 + A - The Jacobian of F(t,U,U_t) 2229 . B - The preconditioner matrix, often the same as A 2230 . f - The function to compute the matrices 2231 - ctx - User-defined context for Jacobian evaluation routine 2232 2233 Notes: You can pass in PETSC_NULL for any return argument you do not need. 2234 2235 Level: advanced 2236 2237 .seealso: TSGetTimeStep(), TSGetRHSJacobian(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber() 2238 2239 .keywords: TS, timestep, get, matrix, Jacobian 2240 @*/ 2241 PetscErrorCode TSGetIJacobian(TS ts,Mat *A,Mat *B,TSIJacobian *f,void **ctx) 2242 { 2243 PetscErrorCode ierr; 2244 SNES snes; 2245 2246 PetscFunctionBegin; 2247 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2248 ierr = SNESGetJacobian(snes,A,B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 2249 if (f) *f = ts->userops->ijacobian; 2250 if (ctx) *ctx = ts->jacP; 2251 PetscFunctionReturn(0); 2252 } 2253 2254 typedef struct { 2255 PetscViewer viewer; 2256 Vec initialsolution; 2257 PetscBool showinitial; 2258 } TSMonitorSolutionCtx; 2259 2260 #undef __FUNCT__ 2261 #define __FUNCT__ "TSMonitorSolution" 2262 /*@C 2263 TSMonitorSolution - Monitors progress of the TS solvers by calling 2264 VecView() for the solution at each timestep 2265 2266 Collective on TS 2267 2268 Input Parameters: 2269 + ts - the TS context 2270 . step - current time-step 2271 . ptime - current time 2272 - dummy - either a viewer or PETSC_NULL 2273 2274 Level: intermediate 2275 2276 .keywords: TS, vector, monitor, view 2277 2278 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView() 2279 @*/ 2280 PetscErrorCode TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy) 2281 { 2282 PetscErrorCode ierr; 2283 TSMonitorSolutionCtx *ictx = (TSMonitorSolutionCtx*)dummy; 2284 2285 PetscFunctionBegin; 2286 if (!step && ictx->showinitial) { 2287 if (!ictx->initialsolution) { 2288 ierr = VecDuplicate(x,&ictx->initialsolution);CHKERRQ(ierr); 2289 } 2290 ierr = VecCopy(x,ictx->initialsolution);CHKERRQ(ierr); 2291 } 2292 if (ictx->showinitial) { 2293 PetscReal pause; 2294 ierr = PetscViewerDrawGetPause(ictx->viewer,&pause);CHKERRQ(ierr); 2295 ierr = PetscViewerDrawSetPause(ictx->viewer,0.0);CHKERRQ(ierr); 2296 ierr = VecView(ictx->initialsolution,ictx->viewer);CHKERRQ(ierr); 2297 ierr = PetscViewerDrawSetPause(ictx->viewer,pause);CHKERRQ(ierr); 2298 ierr = PetscViewerDrawSetHold(ictx->viewer,PETSC_TRUE);CHKERRQ(ierr); 2299 } 2300 ierr = VecView(x,ictx->viewer);CHKERRQ(ierr); 2301 if (ictx->showinitial) { 2302 ierr = PetscViewerDrawSetHold(ictx->viewer,PETSC_FALSE);CHKERRQ(ierr); 2303 } 2304 PetscFunctionReturn(0); 2305 } 2306 2307 2308 #undef __FUNCT__ 2309 #define __FUNCT__ "TSMonitorSolutionDestroy" 2310 /*@C 2311 TSMonitorSolutionDestroy - Destroys the monitor context for TSMonitorSolution 2312 2313 Collective on TS 2314 2315 Input Parameters: 2316 . ctx - the monitor context 2317 2318 Level: intermediate 2319 2320 .keywords: TS, vector, monitor, view 2321 2322 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorSolution() 2323 @*/ 2324 PetscErrorCode TSMonitorSolutionDestroy(void **ctx) 2325 { 2326 PetscErrorCode ierr; 2327 TSMonitorSolutionCtx *ictx = *(TSMonitorSolutionCtx**)ctx; 2328 2329 PetscFunctionBegin; 2330 ierr = PetscViewerDestroy(&ictx->viewer);CHKERRQ(ierr); 2331 ierr = VecDestroy(&ictx->initialsolution);CHKERRQ(ierr); 2332 ierr = PetscFree(ictx);CHKERRQ(ierr); 2333 PetscFunctionReturn(0); 2334 } 2335 2336 #undef __FUNCT__ 2337 #define __FUNCT__ "TSMonitorSolutionCreate" 2338 /*@C 2339 TSMonitorSolutionCreate - Creates the monitor context for TSMonitorSolution 2340 2341 Collective on TS 2342 2343 Input Parameter: 2344 . ts - time-step context 2345 2346 Output Patameter: 2347 . ctx - the monitor context 2348 2349 Level: intermediate 2350 2351 .keywords: TS, vector, monitor, view 2352 2353 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorSolution() 2354 @*/ 2355 PetscErrorCode TSMonitorSolutionCreate(TS ts,PetscViewer viewer,void **ctx) 2356 { 2357 PetscErrorCode ierr; 2358 TSMonitorSolutionCtx *ictx; 2359 2360 PetscFunctionBegin; 2361 ierr = PetscNew(TSMonitorSolutionCtx,&ictx);CHKERRQ(ierr); 2362 *ctx = (void*)ictx; 2363 if (!viewer) { 2364 viewer = PETSC_VIEWER_DRAW_(((PetscObject)ts)->comm); 2365 } 2366 ierr = PetscObjectReference((PetscObject)viewer);CHKERRQ(ierr); 2367 ictx->viewer = viewer; 2368 ictx->showinitial = PETSC_FALSE; 2369 ierr = PetscOptionsGetBool(((PetscObject)ts)->prefix,"-ts_monitor_solution_initial",&ictx->showinitial,PETSC_NULL);CHKERRQ(ierr); 2370 PetscFunctionReturn(0); 2371 } 2372 2373 #undef __FUNCT__ 2374 #define __FUNCT__ "TSSetDM" 2375 /*@ 2376 TSSetDM - Sets the DM that may be used by some preconditioners 2377 2378 Logically Collective on TS and DM 2379 2380 Input Parameters: 2381 + ts - the preconditioner context 2382 - dm - the dm 2383 2384 Level: intermediate 2385 2386 2387 .seealso: TSGetDM(), SNESSetDM(), SNESGetDM() 2388 @*/ 2389 PetscErrorCode TSSetDM(TS ts,DM dm) 2390 { 2391 PetscErrorCode ierr; 2392 SNES snes; 2393 2394 PetscFunctionBegin; 2395 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2396 ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr); 2397 ierr = DMDestroy(&ts->dm);CHKERRQ(ierr); 2398 ts->dm = dm; 2399 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2400 ierr = SNESSetDM(snes,dm);CHKERRQ(ierr); 2401 PetscFunctionReturn(0); 2402 } 2403 2404 #undef __FUNCT__ 2405 #define __FUNCT__ "TSGetDM" 2406 /*@ 2407 TSGetDM - Gets the DM that may be used by some preconditioners 2408 2409 Not Collective 2410 2411 Input Parameter: 2412 . ts - the preconditioner context 2413 2414 Output Parameter: 2415 . dm - the dm 2416 2417 Level: intermediate 2418 2419 2420 .seealso: TSSetDM(), SNESSetDM(), SNESGetDM() 2421 @*/ 2422 PetscErrorCode TSGetDM(TS ts,DM *dm) 2423 { 2424 PetscFunctionBegin; 2425 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2426 *dm = ts->dm; 2427 PetscFunctionReturn(0); 2428 } 2429 2430 #undef __FUNCT__ 2431 #define __FUNCT__ "SNESTSFormFunction" 2432 /*@ 2433 SNESTSFormFunction - Function to evaluate nonlinear residual 2434 2435 Logically Collective on SNES 2436 2437 Input Parameter: 2438 + snes - nonlinear solver 2439 . X - the current state at which to evaluate the residual 2440 - ctx - user context, must be a TS 2441 2442 Output Parameter: 2443 . F - the nonlinear residual 2444 2445 Notes: 2446 This function is not normally called by users and is automatically registered with the SNES used by TS. 2447 It is most frequently passed to MatFDColoringSetFunction(). 2448 2449 Level: advanced 2450 2451 .seealso: SNESSetFunction(), MatFDColoringSetFunction() 2452 @*/ 2453 PetscErrorCode SNESTSFormFunction(SNES snes,Vec X,Vec F,void *ctx) 2454 { 2455 TS ts = (TS)ctx; 2456 PetscErrorCode ierr; 2457 2458 PetscFunctionBegin; 2459 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2460 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 2461 PetscValidHeaderSpecific(F,VEC_CLASSID,3); 2462 PetscValidHeaderSpecific(ts,TS_CLASSID,4); 2463 ierr = (ts->ops->snesfunction)(snes,X,F,ts);CHKERRQ(ierr); 2464 PetscFunctionReturn(0); 2465 } 2466 2467 #undef __FUNCT__ 2468 #define __FUNCT__ "SNESTSFormJacobian" 2469 /*@ 2470 SNESTSFormJacobian - Function to evaluate the Jacobian 2471 2472 Collective on SNES 2473 2474 Input Parameter: 2475 + snes - nonlinear solver 2476 . X - the current state at which to evaluate the residual 2477 - ctx - user context, must be a TS 2478 2479 Output Parameter: 2480 + A - the Jacobian 2481 . B - the preconditioning matrix (may be the same as A) 2482 - flag - indicates any structure change in the matrix 2483 2484 Notes: 2485 This function is not normally called by users and is automatically registered with the SNES used by TS. 2486 2487 Level: developer 2488 2489 .seealso: SNESSetJacobian() 2490 @*/ 2491 PetscErrorCode SNESTSFormJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flag,void *ctx) 2492 { 2493 TS ts = (TS)ctx; 2494 PetscErrorCode ierr; 2495 2496 PetscFunctionBegin; 2497 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2498 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 2499 PetscValidPointer(A,3); 2500 PetscValidHeaderSpecific(*A,MAT_CLASSID,3); 2501 PetscValidPointer(B,4); 2502 PetscValidHeaderSpecific(*B,MAT_CLASSID,4); 2503 PetscValidPointer(flag,5); 2504 PetscValidHeaderSpecific(ts,TS_CLASSID,6); 2505 ierr = (ts->ops->snesjacobian)(snes,X,A,B,flag,ts);CHKERRQ(ierr); 2506 PetscFunctionReturn(0); 2507 } 2508 2509 #undef __FUNCT__ 2510 #define __FUNCT__ "TSComputeRHSFunctionLinear" 2511 /*@C 2512 TSComputeRHSFunctionLinear - Evaluate the right hand side via the user-provided Jacobian, for linear problems only 2513 2514 Collective on TS 2515 2516 Input Arguments: 2517 + ts - time stepping context 2518 . t - time at which to evaluate 2519 . X - state at which to evaluate 2520 - ctx - context 2521 2522 Output Arguments: 2523 . F - right hand side 2524 2525 Level: intermediate 2526 2527 Notes: 2528 This function is intended to be passed to TSSetRHSFunction() to evaluate the right hand side for linear problems. 2529 The matrix (and optionally the evaluation context) should be passed to TSSetRHSJacobian(). 2530 2531 .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSJacobianConstant() 2532 @*/ 2533 PetscErrorCode TSComputeRHSFunctionLinear(TS ts,PetscReal t,Vec X,Vec F,void *ctx) 2534 { 2535 PetscErrorCode ierr; 2536 Mat Arhs,Brhs; 2537 MatStructure flg2; 2538 2539 PetscFunctionBegin; 2540 ierr = TSGetRHSMats_Private(ts,&Arhs,&Brhs);CHKERRQ(ierr); 2541 ierr = TSComputeRHSJacobian(ts,t,X,&Arhs,&Brhs,&flg2);CHKERRQ(ierr); 2542 ierr = MatMult(Arhs,X,F);CHKERRQ(ierr); 2543 PetscFunctionReturn(0); 2544 } 2545 2546 #undef __FUNCT__ 2547 #define __FUNCT__ "TSComputeRHSJacobianConstant" 2548 /*@C 2549 TSComputeRHSJacobianConstant - Reuses a Jacobian that is time-independent. 2550 2551 Collective on TS 2552 2553 Input Arguments: 2554 + ts - time stepping context 2555 . t - time at which to evaluate 2556 . X - state at which to evaluate 2557 - ctx - context 2558 2559 Output Arguments: 2560 + A - pointer to operator 2561 . B - pointer to preconditioning matrix 2562 - flg - matrix structure flag 2563 2564 Level: intermediate 2565 2566 Notes: 2567 This function is intended to be passed to TSSetRHSJacobian() to evaluate the Jacobian for linear time-independent problems. 2568 2569 .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSFunctionLinear() 2570 @*/ 2571 PetscErrorCode TSComputeRHSJacobianConstant(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg,void *ctx) 2572 { 2573 2574 PetscFunctionBegin; 2575 *flg = SAME_PRECONDITIONER; 2576 PetscFunctionReturn(0); 2577 } 2578 2579 #undef __FUNCT__ 2580 #define __FUNCT__ "TSComputeIFunctionLinear" 2581 /*@C 2582 TSComputeIFunctionLinear - Evaluate the left hand side via the user-provided Jacobian, for linear problems only 2583 2584 Collective on TS 2585 2586 Input Arguments: 2587 + ts - time stepping context 2588 . t - time at which to evaluate 2589 . X - state at which to evaluate 2590 . Xdot - time derivative of state vector 2591 - ctx - context 2592 2593 Output Arguments: 2594 . F - left hand side 2595 2596 Level: intermediate 2597 2598 Notes: 2599 The assumption here is that the left hand side is of the form A*Xdot (and not A*Xdot + B*X). For other cases, the 2600 user is required to write their own TSComputeIFunction. 2601 This function is intended to be passed to TSSetIFunction() to evaluate the left hand side for linear problems. 2602 The matrix (and optionally the evaluation context) should be passed to TSSetIJacobian(). 2603 2604 .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIJacobianConstant() 2605 @*/ 2606 PetscErrorCode TSComputeIFunctionLinear(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx) 2607 { 2608 PetscErrorCode ierr; 2609 Mat A,B; 2610 MatStructure flg2; 2611 2612 PetscFunctionBegin; 2613 ierr = TSGetIJacobian(ts,&A,&B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 2614 ierr = TSComputeIJacobian(ts,t,X,Xdot,1.0,&A,&B,&flg2,PETSC_TRUE);CHKERRQ(ierr); 2615 ierr = MatMult(A,Xdot,F);CHKERRQ(ierr); 2616 PetscFunctionReturn(0); 2617 } 2618 2619 #undef __FUNCT__ 2620 #define __FUNCT__ "TSComputeIJacobianConstant" 2621 /*@C 2622 TSComputeRHSJacobianConstant - Reuses a Jacobian that is time-independent. 2623 2624 Collective on TS 2625 2626 Input Arguments: 2627 + ts - time stepping context 2628 . t - time at which to evaluate 2629 . X - state at which to evaluate 2630 . Xdot - time derivative of state vector 2631 . shift - shift to apply 2632 - ctx - context 2633 2634 Output Arguments: 2635 + A - pointer to operator 2636 . B - pointer to preconditioning matrix 2637 - flg - matrix structure flag 2638 2639 Level: intermediate 2640 2641 Notes: 2642 This function is intended to be passed to TSSetIJacobian() to evaluate the Jacobian for linear time-independent problems. 2643 2644 .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIFunctionLinear() 2645 @*/ 2646 PetscErrorCode TSComputeIJacobianConstant(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg,void *ctx) 2647 { 2648 2649 PetscFunctionBegin; 2650 *flg = SAME_PRECONDITIONER; 2651 PetscFunctionReturn(0); 2652 } 2653 2654 2655 #undef __FUNCT__ 2656 #define __FUNCT__ "TSGetConvergedReason" 2657 /*@ 2658 TSGetConvergedReason - Gets the reason the TS iteration was stopped. 2659 2660 Not Collective 2661 2662 Input Parameter: 2663 . ts - the TS context 2664 2665 Output Parameter: 2666 . reason - negative value indicates diverged, positive value converged, see TSConvergedReason or the 2667 manual pages for the individual convergence tests for complete lists 2668 2669 Level: intermediate 2670 2671 Notes: 2672 Can only be called after the call to TSSolve() is complete. 2673 2674 .keywords: TS, nonlinear, set, convergence, test 2675 2676 .seealso: TSSetConvergenceTest(), TSConvergedReason 2677 @*/ 2678 PetscErrorCode TSGetConvergedReason(TS ts,TSConvergedReason *reason) 2679 { 2680 PetscFunctionBegin; 2681 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2682 PetscValidPointer(reason,2); 2683 *reason = ts->reason; 2684 PetscFunctionReturn(0); 2685 } 2686 2687 2688 #undef __FUNCT__ 2689 #define __FUNCT__ "TSVISetVariableBounds" 2690 /*@ 2691 TSVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu 2692 2693 Input Parameters: 2694 . ts - the TS context. 2695 . xl - lower bound. 2696 . xu - upper bound. 2697 2698 Notes: 2699 If this routine is not called then the lower and upper bounds are set to 2700 SNES_VI_INF and SNES_VI_NINF respectively during SNESSetUp(). 2701 2702 Level: advanced 2703 2704 @*/ 2705 PetscErrorCode TSVISetVariableBounds(TS ts, Vec xl, Vec xu) 2706 { 2707 PetscErrorCode ierr; 2708 SNES snes; 2709 2710 PetscFunctionBegin; 2711 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2712 ierr = SNESVISetVariableBounds(snes,xl,xu);CHKERRQ(ierr); 2713 PetscFunctionReturn(0); 2714 } 2715 2716 #if defined(PETSC_HAVE_MATLAB_ENGINE) 2717 #include <mex.h> 2718 2719 typedef struct {char *funcname; mxArray *ctx;} TSMatlabContext; 2720 2721 #undef __FUNCT__ 2722 #define __FUNCT__ "TSComputeFunction_Matlab" 2723 /* 2724 TSComputeFunction_Matlab - Calls the function that has been set with 2725 TSSetFunctionMatlab(). 2726 2727 Collective on TS 2728 2729 Input Parameters: 2730 + snes - the TS context 2731 - x - input vector 2732 2733 Output Parameter: 2734 . y - function vector, as set by TSSetFunction() 2735 2736 Notes: 2737 TSComputeFunction() is typically used within nonlinear solvers 2738 implementations, so most users would not generally call this routine 2739 themselves. 2740 2741 Level: developer 2742 2743 .keywords: TS, nonlinear, compute, function 2744 2745 .seealso: TSSetFunction(), TSGetFunction() 2746 */ 2747 PetscErrorCode TSComputeFunction_Matlab(TS snes,PetscReal time,Vec x,Vec xdot,Vec y, void *ctx) 2748 { 2749 PetscErrorCode ierr; 2750 TSMatlabContext *sctx = (TSMatlabContext *)ctx; 2751 int nlhs = 1,nrhs = 7; 2752 mxArray *plhs[1],*prhs[7]; 2753 long long int lx = 0,lxdot = 0,ly = 0,ls = 0; 2754 2755 PetscFunctionBegin; 2756 PetscValidHeaderSpecific(snes,TS_CLASSID,1); 2757 PetscValidHeaderSpecific(x,VEC_CLASSID,3); 2758 PetscValidHeaderSpecific(xdot,VEC_CLASSID,4); 2759 PetscValidHeaderSpecific(y,VEC_CLASSID,5); 2760 PetscCheckSameComm(snes,1,x,3); 2761 PetscCheckSameComm(snes,1,y,5); 2762 2763 ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 2764 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 2765 ierr = PetscMemcpy(&lxdot,&xdot,sizeof(xdot));CHKERRQ(ierr); 2766 ierr = PetscMemcpy(&ly,&y,sizeof(x));CHKERRQ(ierr); 2767 prhs[0] = mxCreateDoubleScalar((double)ls); 2768 prhs[1] = mxCreateDoubleScalar(time); 2769 prhs[2] = mxCreateDoubleScalar((double)lx); 2770 prhs[3] = mxCreateDoubleScalar((double)lxdot); 2771 prhs[4] = mxCreateDoubleScalar((double)ly); 2772 prhs[5] = mxCreateString(sctx->funcname); 2773 prhs[6] = sctx->ctx; 2774 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeFunctionInternal");CHKERRQ(ierr); 2775 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 2776 mxDestroyArray(prhs[0]); 2777 mxDestroyArray(prhs[1]); 2778 mxDestroyArray(prhs[2]); 2779 mxDestroyArray(prhs[3]); 2780 mxDestroyArray(prhs[4]); 2781 mxDestroyArray(prhs[5]); 2782 mxDestroyArray(plhs[0]); 2783 PetscFunctionReturn(0); 2784 } 2785 2786 2787 #undef __FUNCT__ 2788 #define __FUNCT__ "TSSetFunctionMatlab" 2789 /* 2790 TSSetFunctionMatlab - Sets the function evaluation routine and function 2791 vector for use by the TS routines in solving ODEs 2792 equations from MATLAB. Here the function is a string containing the name of a MATLAB function 2793 2794 Logically Collective on TS 2795 2796 Input Parameters: 2797 + ts - the TS context 2798 - func - function evaluation routine 2799 2800 Calling sequence of func: 2801 $ func (TS ts,PetscReal time,Vec x,Vec xdot,Vec f,void *ctx); 2802 2803 Level: beginner 2804 2805 .keywords: TS, nonlinear, set, function 2806 2807 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction() 2808 */ 2809 PetscErrorCode TSSetFunctionMatlab(TS ts,const char *func,mxArray *ctx) 2810 { 2811 PetscErrorCode ierr; 2812 TSMatlabContext *sctx; 2813 2814 PetscFunctionBegin; 2815 /* currently sctx is memory bleed */ 2816 ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr); 2817 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 2818 /* 2819 This should work, but it doesn't 2820 sctx->ctx = ctx; 2821 mexMakeArrayPersistent(sctx->ctx); 2822 */ 2823 sctx->ctx = mxDuplicateArray(ctx); 2824 ierr = TSSetIFunction(ts,PETSC_NULL,TSComputeFunction_Matlab,sctx);CHKERRQ(ierr); 2825 PetscFunctionReturn(0); 2826 } 2827 2828 #undef __FUNCT__ 2829 #define __FUNCT__ "TSComputeJacobian_Matlab" 2830 /* 2831 TSComputeJacobian_Matlab - Calls the function that has been set with 2832 TSSetJacobianMatlab(). 2833 2834 Collective on TS 2835 2836 Input Parameters: 2837 + ts - the TS context 2838 . x - input vector 2839 . A, B - the matrices 2840 - ctx - user context 2841 2842 Output Parameter: 2843 . flag - structure of the matrix 2844 2845 Level: developer 2846 2847 .keywords: TS, nonlinear, compute, function 2848 2849 .seealso: TSSetFunction(), TSGetFunction() 2850 @*/ 2851 PetscErrorCode TSComputeJacobian_Matlab(TS ts,PetscReal time,Vec x,Vec xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flag, void *ctx) 2852 { 2853 PetscErrorCode ierr; 2854 TSMatlabContext *sctx = (TSMatlabContext *)ctx; 2855 int nlhs = 2,nrhs = 9; 2856 mxArray *plhs[2],*prhs[9]; 2857 long long int lx = 0,lxdot = 0,lA = 0,ls = 0, lB = 0; 2858 2859 PetscFunctionBegin; 2860 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2861 PetscValidHeaderSpecific(x,VEC_CLASSID,3); 2862 2863 /* call Matlab function in ctx with arguments x and y */ 2864 2865 ierr = PetscMemcpy(&ls,&ts,sizeof(ts));CHKERRQ(ierr); 2866 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 2867 ierr = PetscMemcpy(&lxdot,&xdot,sizeof(x));CHKERRQ(ierr); 2868 ierr = PetscMemcpy(&lA,A,sizeof(x));CHKERRQ(ierr); 2869 ierr = PetscMemcpy(&lB,B,sizeof(x));CHKERRQ(ierr); 2870 prhs[0] = mxCreateDoubleScalar((double)ls); 2871 prhs[1] = mxCreateDoubleScalar((double)time); 2872 prhs[2] = mxCreateDoubleScalar((double)lx); 2873 prhs[3] = mxCreateDoubleScalar((double)lxdot); 2874 prhs[4] = mxCreateDoubleScalar((double)shift); 2875 prhs[5] = mxCreateDoubleScalar((double)lA); 2876 prhs[6] = mxCreateDoubleScalar((double)lB); 2877 prhs[7] = mxCreateString(sctx->funcname); 2878 prhs[8] = sctx->ctx; 2879 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeJacobianInternal");CHKERRQ(ierr); 2880 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 2881 *flag = (MatStructure) mxGetScalar(plhs[1]);CHKERRQ(ierr); 2882 mxDestroyArray(prhs[0]); 2883 mxDestroyArray(prhs[1]); 2884 mxDestroyArray(prhs[2]); 2885 mxDestroyArray(prhs[3]); 2886 mxDestroyArray(prhs[4]); 2887 mxDestroyArray(prhs[5]); 2888 mxDestroyArray(prhs[6]); 2889 mxDestroyArray(prhs[7]); 2890 mxDestroyArray(plhs[0]); 2891 mxDestroyArray(plhs[1]); 2892 PetscFunctionReturn(0); 2893 } 2894 2895 2896 #undef __FUNCT__ 2897 #define __FUNCT__ "TSSetJacobianMatlab" 2898 /* 2899 TSSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices 2900 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 2901 2902 Logically Collective on TS 2903 2904 Input Parameters: 2905 + ts - the TS context 2906 . A,B - Jacobian matrices 2907 . func - function evaluation routine 2908 - ctx - user context 2909 2910 Calling sequence of func: 2911 $ flag = func (TS ts,PetscReal time,Vec x,Vec xdot,Mat A,Mat B,void *ctx); 2912 2913 2914 Level: developer 2915 2916 .keywords: TS, nonlinear, set, function 2917 2918 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction() 2919 */ 2920 PetscErrorCode TSSetJacobianMatlab(TS ts,Mat A,Mat B,const char *func,mxArray *ctx) 2921 { 2922 PetscErrorCode ierr; 2923 TSMatlabContext *sctx; 2924 2925 PetscFunctionBegin; 2926 /* currently sctx is memory bleed */ 2927 ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr); 2928 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 2929 /* 2930 This should work, but it doesn't 2931 sctx->ctx = ctx; 2932 mexMakeArrayPersistent(sctx->ctx); 2933 */ 2934 sctx->ctx = mxDuplicateArray(ctx); 2935 ierr = TSSetIJacobian(ts,A,B,TSComputeJacobian_Matlab,sctx);CHKERRQ(ierr); 2936 PetscFunctionReturn(0); 2937 } 2938 2939 #undef __FUNCT__ 2940 #define __FUNCT__ "TSMonitor_Matlab" 2941 /* 2942 TSMonitor_Matlab - Calls the function that has been set with TSMonitorSetMatlab(). 2943 2944 Collective on TS 2945 2946 .seealso: TSSetFunction(), TSGetFunction() 2947 @*/ 2948 PetscErrorCode TSMonitor_Matlab(TS ts,PetscInt it, PetscReal time,Vec x, void *ctx) 2949 { 2950 PetscErrorCode ierr; 2951 TSMatlabContext *sctx = (TSMatlabContext *)ctx; 2952 int nlhs = 1,nrhs = 6; 2953 mxArray *plhs[1],*prhs[6]; 2954 long long int lx = 0,ls = 0; 2955 2956 PetscFunctionBegin; 2957 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2958 PetscValidHeaderSpecific(x,VEC_CLASSID,4); 2959 2960 ierr = PetscMemcpy(&ls,&ts,sizeof(ts));CHKERRQ(ierr); 2961 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 2962 prhs[0] = mxCreateDoubleScalar((double)ls); 2963 prhs[1] = mxCreateDoubleScalar((double)it); 2964 prhs[2] = mxCreateDoubleScalar((double)time); 2965 prhs[3] = mxCreateDoubleScalar((double)lx); 2966 prhs[4] = mxCreateString(sctx->funcname); 2967 prhs[5] = sctx->ctx; 2968 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSMonitorInternal");CHKERRQ(ierr); 2969 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 2970 mxDestroyArray(prhs[0]); 2971 mxDestroyArray(prhs[1]); 2972 mxDestroyArray(prhs[2]); 2973 mxDestroyArray(prhs[3]); 2974 mxDestroyArray(prhs[4]); 2975 mxDestroyArray(plhs[0]); 2976 PetscFunctionReturn(0); 2977 } 2978 2979 2980 #undef __FUNCT__ 2981 #define __FUNCT__ "TSMonitorSetMatlab" 2982 /* 2983 TSMonitorSetMatlab - Sets the monitor function from Matlab 2984 2985 Level: developer 2986 2987 .keywords: TS, nonlinear, set, function 2988 2989 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction() 2990 */ 2991 PetscErrorCode TSMonitorSetMatlab(TS ts,const char *func,mxArray *ctx) 2992 { 2993 PetscErrorCode ierr; 2994 TSMatlabContext *sctx; 2995 2996 PetscFunctionBegin; 2997 /* currently sctx is memory bleed */ 2998 ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr); 2999 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 3000 /* 3001 This should work, but it doesn't 3002 sctx->ctx = ctx; 3003 mexMakeArrayPersistent(sctx->ctx); 3004 */ 3005 sctx->ctx = mxDuplicateArray(ctx); 3006 ierr = TSMonitorSet(ts,TSMonitor_Matlab,sctx,PETSC_NULL);CHKERRQ(ierr); 3007 PetscFunctionReturn(0); 3008 } 3009 #endif 3010