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