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