1 2 #include <petsc-private/tsimpl.h> /*I "petscts.h" I*/ 3 #include <petscdmshell.h> 4 5 /* Logging support */ 6 PetscClassId TS_CLASSID; 7 PetscLogEvent TS_Step, TS_PseudoComputeTimeStep, TS_FunctionEval, TS_JacobianEval; 8 9 #undef __FUNCT__ 10 #define __FUNCT__ "TSSetTypeFromOptions" 11 /* 12 TSSetTypeFromOptions - Sets the type of ts from user options. 13 14 Collective on TS 15 16 Input Parameter: 17 . ts - The ts 18 19 Level: intermediate 20 21 .keywords: TS, set, options, database, type 22 .seealso: TSSetFromOptions(), TSSetType() 23 */ 24 static PetscErrorCode TSSetTypeFromOptions(TS ts) 25 { 26 PetscBool opt; 27 const char *defaultType; 28 char typeName[256]; 29 PetscErrorCode ierr; 30 31 PetscFunctionBegin; 32 if (((PetscObject)ts)->type_name) { 33 defaultType = ((PetscObject)ts)->type_name; 34 } else { 35 defaultType = TSEULER; 36 } 37 38 if (!TSRegisterAllCalled) {ierr = TSRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 39 ierr = PetscOptionsList("-ts_type", "TS method"," TSSetType", TSList, defaultType, typeName, 256, &opt);CHKERRQ(ierr); 40 if (opt) { 41 ierr = TSSetType(ts, typeName);CHKERRQ(ierr); 42 } else { 43 ierr = TSSetType(ts, defaultType);CHKERRQ(ierr); 44 } 45 PetscFunctionReturn(0); 46 } 47 48 #undef __FUNCT__ 49 #define __FUNCT__ "TSSetFromOptions" 50 /*@ 51 TSSetFromOptions - Sets various TS parameters from user options. 52 53 Collective on TS 54 55 Input Parameter: 56 . ts - the TS context obtained from TSCreate() 57 58 Options Database Keys: 59 + -ts_type <type> - TSEULER, TSBEULER, TSSUNDIALS, TSPSEUDO, TSCN, TSRK, TSTHETA, TSGL, TSSSP 60 . -ts_max_steps maxsteps - maximum number of time-steps to take 61 . -ts_final_time time - maximum time to compute to 62 . -ts_dt dt - initial time step 63 . -ts_monitor - print information at each timestep 64 - -ts_monitor_draw - plot information at each timestep 65 66 Level: beginner 67 68 .keywords: TS, timestep, set, options, database 69 70 .seealso: TSGetType() 71 @*/ 72 PetscErrorCode TSSetFromOptions(TS ts) 73 { 74 PetscBool opt,flg; 75 PetscErrorCode ierr; 76 PetscViewer monviewer; 77 char monfilename[PETSC_MAX_PATH_LEN]; 78 SNES snes; 79 TSAdapt adapt; 80 PetscReal time_step; 81 82 PetscFunctionBegin; 83 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 84 ierr = PetscObjectOptionsBegin((PetscObject)ts);CHKERRQ(ierr); 85 /* Handle TS type options */ 86 ierr = TSSetTypeFromOptions(ts);CHKERRQ(ierr); 87 88 /* Handle generic TS options */ 89 ierr = PetscOptionsInt("-ts_max_steps","Maximum number of time steps","TSSetDuration",ts->max_steps,&ts->max_steps,PETSC_NULL);CHKERRQ(ierr); 90 ierr = PetscOptionsReal("-ts_final_time","Time to run to","TSSetDuration",ts->max_time,&ts->max_time,PETSC_NULL);CHKERRQ(ierr); 91 ierr = PetscOptionsReal("-ts_init_time","Initial time","TSSetTime",ts->ptime,&ts->ptime,PETSC_NULL);CHKERRQ(ierr); 92 ierr = PetscOptionsReal("-ts_dt","Initial time step","TSSetTimeStep",ts->time_step,&time_step,&flg);CHKERRQ(ierr); 93 if (flg) { 94 ierr = TSSetTimeStep(ts,time_step);CHKERRQ(ierr); 95 } 96 opt = ts->exact_final_time == PETSC_DECIDE ? PETSC_FALSE : (PetscBool)ts->exact_final_time; 97 ierr = PetscOptionsBool("-ts_exact_final_time","Interpolate output to stop exactly at the final time","TSSetExactFinalTime",opt,&opt,&flg);CHKERRQ(ierr); 98 if (flg) {ierr = TSSetExactFinalTime(ts,opt);CHKERRQ(ierr);} 99 ierr = PetscOptionsInt("-ts_max_snes_failures","Maximum number of nonlinear solve failures","TSSetMaxSNESFailures",ts->max_snes_failures,&ts->max_snes_failures,PETSC_NULL);CHKERRQ(ierr); 100 ierr = PetscOptionsInt("-ts_max_reject","Maximum number of step rejections before step fails","TSSetMaxStepRejections",ts->max_reject,&ts->max_reject,PETSC_NULL);CHKERRQ(ierr); 101 ierr = PetscOptionsBool("-ts_error_if_step_fails","Error if no step succeeds","TSSetErrorIfStepFails",ts->errorifstepfailed,&ts->errorifstepfailed,PETSC_NULL);CHKERRQ(ierr); 102 ierr = PetscOptionsReal("-ts_rtol","Relative tolerance for local truncation error","TSSetTolerances",ts->rtol,&ts->rtol,PETSC_NULL);CHKERRQ(ierr); 103 ierr = PetscOptionsReal("-ts_atol","Absolute tolerance for local truncation error","TSSetTolerances",ts->atol,&ts->atol,PETSC_NULL);CHKERRQ(ierr); 104 105 /* Monitor options */ 106 ierr = PetscOptionsString("-ts_monitor","Monitor timestep size","TSMonitorDefault","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 107 if (flg) { 108 ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,monfilename,&monviewer);CHKERRQ(ierr); 109 ierr = TSMonitorSet(ts,TSMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 110 } 111 ierr = PetscOptionsString("-ts_monitor_python","Use Python function","TSMonitorSet",0,monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 112 if (flg) {ierr = PetscPythonMonitorSet((PetscObject)ts,monfilename);CHKERRQ(ierr);} 113 114 opt = PETSC_FALSE; 115 ierr = PetscOptionsBool("-ts_monitor_draw","Monitor timestep size graphically","TSMonitorLG",opt,&opt,PETSC_NULL);CHKERRQ(ierr); 116 if (opt) { 117 ierr = TSMonitorSet(ts,TSMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 118 } 119 opt = PETSC_FALSE; 120 ierr = PetscOptionsBool("-ts_monitor_solution","Monitor solution graphically","TSMonitorSolution",opt,&opt,PETSC_NULL);CHKERRQ(ierr); 121 if (opt) { 122 void *ctx; 123 ierr = TSMonitorSolutionCreate(ts,PETSC_NULL,&ctx);CHKERRQ(ierr); 124 ierr = TSMonitorSet(ts,TSMonitorSolution,ctx,TSMonitorSolutionDestroy);CHKERRQ(ierr); 125 } 126 opt = PETSC_FALSE; 127 ierr = PetscOptionsString("-ts_monitor_solution_binary","Save each solution to a binary file","TSMonitorSolutionBinary",0,monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 128 if (flg) { 129 PetscViewer ctx; 130 if (monfilename[0]) { 131 ierr = PetscViewerBinaryOpen(((PetscObject)ts)->comm,monfilename,FILE_MODE_WRITE,&ctx);CHKERRQ(ierr); 132 } else { 133 ctx = PETSC_VIEWER_BINARY_(((PetscObject)ts)->comm); 134 } 135 ierr = TSMonitorSet(ts,TSMonitorSolutionBinary,ctx,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 136 } 137 opt = PETSC_FALSE; 138 ierr = PetscOptionsString("-ts_monitor_solution_vtk","Save each time step to a binary file, use filename-%%03D.vts","TSMonitorSolutionVTK",0,monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 139 if (flg) { 140 const char *ptr,*ptr2; 141 char *filetemplate; 142 if (!monfilename[0]) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"-ts_monitor_solution_vtk requires a file template, e.g. filename-%%03D.vts"); 143 /* Do some cursory validation of the input. */ 144 ierr = PetscStrstr(monfilename,"%",(char**)&ptr);CHKERRQ(ierr); 145 if (!ptr) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"-ts_monitor_solution_vtk requires a file template, e.g. filename-%%03D.vts"); 146 for (ptr++ ; ptr && *ptr; ptr++) { 147 ierr = PetscStrchr("DdiouxX",*ptr,(char**)&ptr2);CHKERRQ(ierr); 148 if (!ptr2 && (*ptr < '0' || '9' < *ptr)) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Invalid file template argument to -ts_monitor_solution_vtk, should look like filename-%%03D.vts"); 149 if (ptr2) break; 150 } 151 ierr = PetscStrallocpy(monfilename,&filetemplate);CHKERRQ(ierr); 152 ierr = TSMonitorSet(ts,TSMonitorSolutionVTK,filetemplate,(PetscErrorCode (*)(void**))TSMonitorSolutionVTKDestroy);CHKERRQ(ierr); 153 } 154 155 ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 156 ierr = TSAdaptSetFromOptions(adapt);CHKERRQ(ierr); 157 158 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 159 if (ts->problem_type == TS_LINEAR) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);} 160 161 /* Handle specific TS options */ 162 if (ts->ops->setfromoptions) { 163 ierr = (*ts->ops->setfromoptions)(ts);CHKERRQ(ierr); 164 } 165 166 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 167 ierr = PetscObjectProcessOptionsHandlers((PetscObject)ts);CHKERRQ(ierr); 168 ierr = PetscOptionsEnd();CHKERRQ(ierr); 169 PetscFunctionReturn(0); 170 } 171 172 #undef __FUNCT__ 173 #undef __FUNCT__ 174 #define __FUNCT__ "TSComputeRHSJacobian" 175 /*@ 176 TSComputeRHSJacobian - Computes the Jacobian matrix that has been 177 set with TSSetRHSJacobian(). 178 179 Collective on TS and Vec 180 181 Input Parameters: 182 + ts - the TS context 183 . t - current timestep 184 - x - input vector 185 186 Output Parameters: 187 + A - Jacobian matrix 188 . B - optional preconditioning matrix 189 - flag - flag indicating matrix structure 190 191 Notes: 192 Most users should not need to explicitly call this routine, as it 193 is used internally within the nonlinear solvers. 194 195 See KSPSetOperators() for important information about setting the 196 flag parameter. 197 198 Level: developer 199 200 .keywords: SNES, compute, Jacobian, matrix 201 202 .seealso: TSSetRHSJacobian(), KSPSetOperators() 203 @*/ 204 PetscErrorCode TSComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg) 205 { 206 PetscErrorCode ierr; 207 PetscInt Xstate; 208 DM dm; 209 TSDM tsdm; 210 TSRHSJacobian rhsjacobianfunc; 211 void *ctx; 212 TSIJacobian ijacobianfunc; 213 214 PetscFunctionBegin; 215 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 216 PetscValidHeaderSpecific(X,VEC_CLASSID,3); 217 PetscCheckSameComm(ts,1,X,3); 218 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 219 ierr = DMTSGetContext(dm,&tsdm);CHKERRQ(ierr); 220 ierr = DMTSGetRHSJacobian(dm,&rhsjacobianfunc,&ctx);CHKERRQ(ierr); 221 ierr = DMTSGetIJacobian(dm,&ijacobianfunc,PETSC_NULL);CHKERRQ(ierr); 222 ierr = PetscObjectStateQuery((PetscObject)X,&Xstate);CHKERRQ(ierr); 223 if (ts->rhsjacobian.time == t && (ts->problem_type == TS_LINEAR || (ts->rhsjacobian.X == X && ts->rhsjacobian.Xstate == Xstate))) { 224 *flg = ts->rhsjacobian.mstructure; 225 PetscFunctionReturn(0); 226 } 227 228 if (!rhsjacobianfunc && !ijacobianfunc) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSJacobian() and / or TSSetIJacobian()"); 229 230 if (rhsjacobianfunc) { 231 ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 232 *flg = DIFFERENT_NONZERO_PATTERN; 233 PetscStackPush("TS user Jacobian function"); 234 ierr = (*rhsjacobianfunc)(ts,t,X,A,B,flg,ctx);CHKERRQ(ierr); 235 PetscStackPop; 236 ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 237 /* make sure user returned a correct Jacobian and preconditioner */ 238 PetscValidHeaderSpecific(*A,MAT_CLASSID,4); 239 PetscValidHeaderSpecific(*B,MAT_CLASSID,5); 240 } else { 241 ierr = MatZeroEntries(*A);CHKERRQ(ierr); 242 if (*A != *B) {ierr = MatZeroEntries(*B);CHKERRQ(ierr);} 243 *flg = SAME_NONZERO_PATTERN; 244 } 245 ts->rhsjacobian.time = t; 246 ts->rhsjacobian.X = X; 247 ierr = PetscObjectStateQuery((PetscObject)X,&ts->rhsjacobian.Xstate);CHKERRQ(ierr); 248 ts->rhsjacobian.mstructure = *flg; 249 PetscFunctionReturn(0); 250 } 251 252 #undef __FUNCT__ 253 #define __FUNCT__ "TSComputeRHSFunction" 254 /*@ 255 TSComputeRHSFunction - Evaluates the right-hand-side function. 256 257 Collective on TS and Vec 258 259 Input Parameters: 260 + ts - the TS context 261 . t - current time 262 - x - state vector 263 264 Output Parameter: 265 . y - right hand side 266 267 Note: 268 Most users should not need to explicitly call this routine, as it 269 is used internally within the nonlinear solvers. 270 271 Level: developer 272 273 .keywords: TS, compute 274 275 .seealso: TSSetRHSFunction(), TSComputeIFunction() 276 @*/ 277 PetscErrorCode TSComputeRHSFunction(TS ts,PetscReal t,Vec x,Vec y) 278 { 279 PetscErrorCode ierr; 280 281 TSRHSFunction rhsfunction; 282 TSIFunction ifunction; 283 void *ctx; 284 DM dm; 285 286 PetscFunctionBegin; 287 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 288 PetscValidHeaderSpecific(x,VEC_CLASSID,3); 289 PetscValidHeaderSpecific(y,VEC_CLASSID,4); 290 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 291 ierr = DMTSGetRHSFunction(dm,&rhsfunction,&ctx);CHKERRQ(ierr); 292 ierr = DMTSGetIFunction(dm,&ifunction,PETSC_NULL);CHKERRQ(ierr); 293 294 if (!rhsfunction && !ifunction) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSFunction() and / or TSSetIFunction()"); 295 296 ierr = PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr); 297 if (rhsfunction) { 298 PetscStackPush("TS user right-hand-side function"); 299 ierr = (*rhsfunction)(ts,t,x,y,ctx);CHKERRQ(ierr); 300 PetscStackPop; 301 } else { 302 ierr = VecZeroEntries(y);CHKERRQ(ierr); 303 } 304 305 ierr = PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr); 306 PetscFunctionReturn(0); 307 } 308 309 #undef __FUNCT__ 310 #define __FUNCT__ "TSGetRHSVec_Private" 311 static PetscErrorCode TSGetRHSVec_Private(TS ts,Vec *Frhs) 312 { 313 Vec F; 314 PetscErrorCode ierr; 315 316 PetscFunctionBegin; 317 *Frhs = PETSC_NULL; 318 ierr = TSGetIFunction(ts,&F,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 319 if (!ts->Frhs) { 320 ierr = VecDuplicate(F,&ts->Frhs);CHKERRQ(ierr); 321 } 322 *Frhs = ts->Frhs; 323 PetscFunctionReturn(0); 324 } 325 326 #undef __FUNCT__ 327 #define __FUNCT__ "TSGetRHSMats_Private" 328 static PetscErrorCode TSGetRHSMats_Private(TS ts,Mat *Arhs,Mat *Brhs) 329 { 330 Mat A,B; 331 PetscErrorCode ierr; 332 333 PetscFunctionBegin; 334 ierr = TSGetIJacobian(ts,&A,&B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 335 if (Arhs) { 336 if (!ts->Arhs) { 337 ierr = MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&ts->Arhs);CHKERRQ(ierr); 338 } 339 *Arhs = ts->Arhs; 340 } 341 if (Brhs) { 342 if (!ts->Brhs) { 343 ierr = MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&ts->Brhs);CHKERRQ(ierr); 344 } 345 *Brhs = ts->Brhs; 346 } 347 PetscFunctionReturn(0); 348 } 349 350 #undef __FUNCT__ 351 #define __FUNCT__ "TSComputeIFunction" 352 /*@ 353 TSComputeIFunction - Evaluates the DAE residual written in implicit form F(t,X,Xdot)=0 354 355 Collective on TS and Vec 356 357 Input Parameters: 358 + ts - the TS context 359 . t - current time 360 . X - state vector 361 . Xdot - time derivative of state vector 362 - imex - flag indicates if the method is IMEX so that the RHSFunction should be kept separate 363 364 Output Parameter: 365 . Y - right hand side 366 367 Note: 368 Most users should not need to explicitly call this routine, as it 369 is used internally within the nonlinear solvers. 370 371 If the user did did not write their equations in implicit form, this 372 function recasts them in implicit form. 373 374 Level: developer 375 376 .keywords: TS, compute 377 378 .seealso: TSSetIFunction(), TSComputeRHSFunction() 379 @*/ 380 PetscErrorCode TSComputeIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec Y,PetscBool imex) 381 { 382 PetscErrorCode ierr; 383 TSIFunction ifunction; 384 TSRHSFunction rhsfunction; 385 void *ctx; 386 DM dm; 387 388 PetscFunctionBegin; 389 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 390 PetscValidHeaderSpecific(X,VEC_CLASSID,3); 391 PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4); 392 PetscValidHeaderSpecific(Y,VEC_CLASSID,5); 393 394 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 395 ierr = DMTSGetIFunction(dm,&ifunction,&ctx);CHKERRQ(ierr); 396 ierr = DMTSGetRHSFunction(dm,&rhsfunction,PETSC_NULL);CHKERRQ(ierr); 397 398 if (!rhsfunction && !ifunction) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSFunction() and / or TSSetIFunction()"); 399 400 ierr = PetscLogEventBegin(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr); 401 if (ifunction) { 402 PetscStackPush("TS user implicit function"); 403 ierr = (*ifunction)(ts,t,X,Xdot,Y,ctx);CHKERRQ(ierr); 404 PetscStackPop; 405 } 406 if (imex) { 407 if (!ifunction) { 408 ierr = VecCopy(Xdot,Y);CHKERRQ(ierr); 409 } 410 } else if (rhsfunction) { 411 if (ifunction) { 412 Vec Frhs; 413 ierr = TSGetRHSVec_Private(ts,&Frhs);CHKERRQ(ierr); 414 ierr = TSComputeRHSFunction(ts,t,X,Frhs);CHKERRQ(ierr); 415 ierr = VecAXPY(Y,-1,Frhs);CHKERRQ(ierr); 416 } else { 417 ierr = TSComputeRHSFunction(ts,t,X,Y);CHKERRQ(ierr); 418 ierr = VecAYPX(Y,-1,Xdot);CHKERRQ(ierr); 419 } 420 } 421 ierr = PetscLogEventEnd(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr); 422 PetscFunctionReturn(0); 423 } 424 425 #undef __FUNCT__ 426 #define __FUNCT__ "TSComputeIJacobian" 427 /*@ 428 TSComputeIJacobian - Evaluates the Jacobian of the DAE 429 430 Collective on TS and Vec 431 432 Input 433 Input Parameters: 434 + ts - the TS context 435 . t - current timestep 436 . X - state vector 437 . Xdot - time derivative of state vector 438 . shift - shift to apply, see note below 439 - imex - flag indicates if the method is IMEX so that the RHSJacobian should be kept separate 440 441 Output Parameters: 442 + A - Jacobian matrix 443 . B - optional preconditioning matrix 444 - flag - flag indicating matrix structure 445 446 Notes: 447 If F(t,X,Xdot)=0 is the DAE, the required Jacobian is 448 449 dF/dX + shift*dF/dXdot 450 451 Most users should not need to explicitly call this routine, as it 452 is used internally within the nonlinear solvers. 453 454 Level: developer 455 456 .keywords: TS, compute, Jacobian, matrix 457 458 .seealso: TSSetIJacobian() 459 @*/ 460 PetscErrorCode TSComputeIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg,PetscBool imex) 461 { 462 PetscInt Xstate, Xdotstate; 463 PetscErrorCode ierr; 464 TSIJacobian ijacobian; 465 TSRHSJacobian rhsjacobian; 466 DM dm; 467 void *ctx; 468 469 PetscFunctionBegin; 470 471 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 472 PetscValidHeaderSpecific(X,VEC_CLASSID,3); 473 PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4); 474 PetscValidPointer(A,6); 475 PetscValidHeaderSpecific(*A,MAT_CLASSID,6); 476 PetscValidPointer(B,7); 477 PetscValidHeaderSpecific(*B,MAT_CLASSID,7); 478 PetscValidPointer(flg,8); 479 480 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 481 ierr = DMTSGetIJacobian(dm,&ijacobian,&ctx);CHKERRQ(ierr); 482 ierr = DMTSGetRHSJacobian(dm,&rhsjacobian,PETSC_NULL);CHKERRQ(ierr); 483 484 ierr = PetscObjectStateQuery((PetscObject)X,&Xstate);CHKERRQ(ierr); 485 ierr = PetscObjectStateQuery((PetscObject)Xdot,&Xdotstate);CHKERRQ(ierr); 486 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))) { 487 *flg = ts->ijacobian.mstructure; 488 ierr = MatScale(*A, shift / ts->ijacobian.shift);CHKERRQ(ierr); 489 PetscFunctionReturn(0); 490 } 491 492 if (!rhsjacobian && !ijacobian) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSJacobian() and / or TSSetIJacobian()"); 493 494 *flg = SAME_NONZERO_PATTERN; /* In case we're solving a linear problem in which case it wouldn't get initialized below. */ 495 ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 496 if (ijacobian) { 497 *flg = DIFFERENT_NONZERO_PATTERN; 498 PetscStackPush("TS user implicit Jacobian"); 499 ierr = (*ijacobian)(ts,t,X,Xdot,shift,A,B,flg,ctx);CHKERRQ(ierr); 500 PetscStackPop; 501 /* make sure user returned a correct Jacobian and preconditioner */ 502 PetscValidHeaderSpecific(*A,MAT_CLASSID,4); 503 PetscValidHeaderSpecific(*B,MAT_CLASSID,5); 504 } 505 if (imex) { 506 if (!ijacobian) { /* system was written as Xdot = F(t,X) */ 507 ierr = MatZeroEntries(*A);CHKERRQ(ierr); 508 ierr = MatShift(*A,shift);CHKERRQ(ierr); 509 if (*A != *B) { 510 ierr = MatZeroEntries(*B);CHKERRQ(ierr); 511 ierr = MatShift(*B,shift);CHKERRQ(ierr); 512 } 513 *flg = SAME_PRECONDITIONER; 514 } 515 } else { 516 if (!ijacobian) { 517 ierr = TSComputeRHSJacobian(ts,t,X,A,B,flg);CHKERRQ(ierr); 518 ierr = MatScale(*A,-1);CHKERRQ(ierr); 519 ierr = MatShift(*A,shift);CHKERRQ(ierr); 520 if (*A != *B) { 521 ierr = MatScale(*B,-1);CHKERRQ(ierr); 522 ierr = MatShift(*B,shift);CHKERRQ(ierr); 523 } 524 } else if (rhsjacobian) { 525 Mat Arhs,Brhs; 526 MatStructure axpy,flg2 = DIFFERENT_NONZERO_PATTERN; 527 ierr = TSGetRHSMats_Private(ts,&Arhs,&Brhs);CHKERRQ(ierr); 528 ierr = TSComputeRHSJacobian(ts,t,X,&Arhs,&Brhs,&flg2);CHKERRQ(ierr); 529 axpy = (*flg == flg2) ? SAME_NONZERO_PATTERN : DIFFERENT_NONZERO_PATTERN; 530 ierr = MatAXPY(*A,-1,Arhs,axpy);CHKERRQ(ierr); 531 if (*A != *B) { 532 ierr = MatAXPY(*B,-1,Brhs,axpy);CHKERRQ(ierr); 533 } 534 *flg = PetscMin(*flg,flg2); 535 } 536 } 537 538 ts->ijacobian.time = t; 539 ts->ijacobian.X = X; 540 ts->ijacobian.Xdot = Xdot; 541 ierr = PetscObjectStateQuery((PetscObject)X,&ts->ijacobian.Xstate);CHKERRQ(ierr); 542 ierr = PetscObjectStateQuery((PetscObject)Xdot,&ts->ijacobian.Xdotstate);CHKERRQ(ierr); 543 ts->ijacobian.shift = shift; 544 ts->ijacobian.imex = imex; 545 ts->ijacobian.mstructure = *flg; 546 ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 547 PetscFunctionReturn(0); 548 } 549 550 #undef __FUNCT__ 551 #define __FUNCT__ "TSSetRHSFunction" 552 /*@C 553 TSSetRHSFunction - Sets the routine for evaluating the function, 554 F(t,u), where U_t = F(t,u). 555 556 Logically Collective on TS 557 558 Input Parameters: 559 + ts - the TS context obtained from TSCreate() 560 . r - vector to put the computed right hand side (or PETSC_NULL to have it created) 561 . f - routine for evaluating the right-hand-side function 562 - ctx - [optional] user-defined context for private data for the 563 function evaluation routine (may be PETSC_NULL) 564 565 Calling sequence of func: 566 $ func (TS ts,PetscReal t,Vec u,Vec F,void *ctx); 567 568 + t - current timestep 569 . u - input vector 570 . F - function vector 571 - ctx - [optional] user-defined function context 572 573 Level: beginner 574 575 .keywords: TS, timestep, set, right-hand-side, function 576 577 .seealso: TSSetRHSJacobian(), TSSetIJacobian() 578 @*/ 579 PetscErrorCode TSSetRHSFunction(TS ts,Vec r,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx) 580 { 581 PetscErrorCode ierr; 582 SNES snes; 583 Vec ralloc = PETSC_NULL; 584 DM dm; 585 586 PetscFunctionBegin; 587 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 588 if (r) PetscValidHeaderSpecific(r,VEC_CLASSID,2); 589 590 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 591 ierr = DMTSSetRHSFunction(dm,f,ctx);CHKERRQ(ierr); 592 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 593 if (!r && !ts->dm && ts->vec_sol) { 594 ierr = VecDuplicate(ts->vec_sol,&ralloc);CHKERRQ(ierr); 595 r = ralloc; 596 } 597 ierr = SNESSetFunction(snes,r,SNESTSFormFunction,ts);CHKERRQ(ierr); 598 ierr = VecDestroy(&ralloc);CHKERRQ(ierr); 599 PetscFunctionReturn(0); 600 } 601 602 #undef __FUNCT__ 603 #define __FUNCT__ "TSSetRHSJacobian" 604 /*@C 605 TSSetRHSJacobian - Sets the function to compute the Jacobian of F, 606 where U_t = F(U,t), as well as the location to store the matrix. 607 608 Logically Collective on TS 609 610 Input Parameters: 611 + ts - the TS context obtained from TSCreate() 612 . A - Jacobian matrix 613 . B - preconditioner matrix (usually same as A) 614 . f - the Jacobian evaluation routine 615 - ctx - [optional] user-defined context for private data for the 616 Jacobian evaluation routine (may be PETSC_NULL) 617 618 Calling sequence of func: 619 $ func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx); 620 621 + t - current timestep 622 . u - input vector 623 . A - matrix A, where U_t = A(t)u 624 . B - preconditioner matrix, usually the same as A 625 . flag - flag indicating information about the preconditioner matrix 626 structure (same as flag in KSPSetOperators()) 627 - ctx - [optional] user-defined context for matrix evaluation routine 628 629 Notes: 630 See KSPSetOperators() for important information about setting the flag 631 output parameter in the routine func(). Be sure to read this information! 632 633 The routine func() takes Mat * as the matrix arguments rather than Mat. 634 This allows the matrix evaluation routine to replace A and/or B with a 635 completely new matrix structure (not just different matrix elements) 636 when appropriate, for instance, if the nonzero structure is changing 637 throughout the global iterations. 638 639 Level: beginner 640 641 .keywords: TS, timestep, set, right-hand-side, Jacobian 642 643 .seealso: SNESDefaultComputeJacobianColor(), TSSetRHSFunction() 644 645 @*/ 646 PetscErrorCode TSSetRHSJacobian(TS ts,Mat A,Mat B,TSRHSJacobian f,void *ctx) 647 { 648 PetscErrorCode ierr; 649 SNES snes; 650 DM dm; 651 TSIJacobian ijacobian; 652 653 PetscFunctionBegin; 654 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 655 if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2); 656 if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3); 657 if (A) PetscCheckSameComm(ts,1,A,2); 658 if (B) PetscCheckSameComm(ts,1,B,3); 659 660 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 661 ierr = DMTSSetRHSJacobian(dm,f,ctx);CHKERRQ(ierr); 662 ierr = DMTSGetIJacobian(dm,&ijacobian,PETSC_NULL);CHKERRQ(ierr); 663 664 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 665 if (!ijacobian) { 666 ierr = SNESSetJacobian(snes,A,B,SNESTSFormJacobian,ts);CHKERRQ(ierr); 667 } 668 if (A) { 669 ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 670 ierr = MatDestroy(&ts->Arhs);CHKERRQ(ierr); 671 ts->Arhs = A; 672 } 673 if (B) { 674 ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 675 ierr = MatDestroy(&ts->Brhs);CHKERRQ(ierr); 676 ts->Brhs = B; 677 } 678 PetscFunctionReturn(0); 679 } 680 681 682 #undef __FUNCT__ 683 #define __FUNCT__ "TSSetIFunction" 684 /*@C 685 TSSetIFunction - Set the function to compute F(t,U,U_t) where F = 0 is the DAE to be solved. 686 687 Logically Collective on TS 688 689 Input Parameters: 690 + ts - the TS context obtained from TSCreate() 691 . r - vector to hold the residual (or PETSC_NULL to have it created internally) 692 . f - the function evaluation routine 693 - ctx - user-defined context for private data for the function evaluation routine (may be PETSC_NULL) 694 695 Calling sequence of f: 696 $ f(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx); 697 698 + t - time at step/stage being solved 699 . u - state vector 700 . u_t - time derivative of state vector 701 . F - function vector 702 - ctx - [optional] user-defined context for matrix evaluation routine 703 704 Important: 705 The user MUST call either this routine, TSSetRHSFunction(). This routine must be used when not solving an ODE, for example a DAE. 706 707 Level: beginner 708 709 .keywords: TS, timestep, set, DAE, Jacobian 710 711 .seealso: TSSetRHSJacobian(), TSSetRHSFunction(), TSSetIJacobian() 712 @*/ 713 PetscErrorCode TSSetIFunction(TS ts,Vec res,TSIFunction f,void *ctx) 714 { 715 PetscErrorCode ierr; 716 SNES snes; 717 Vec resalloc = PETSC_NULL; 718 DM dm; 719 720 PetscFunctionBegin; 721 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 722 if (res) PetscValidHeaderSpecific(res,VEC_CLASSID,2); 723 724 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 725 ierr = DMTSSetIFunction(dm,f,ctx);CHKERRQ(ierr); 726 727 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 728 if (!res && !ts->dm && ts->vec_sol) { 729 ierr = VecDuplicate(ts->vec_sol,&resalloc);CHKERRQ(ierr); 730 res = resalloc; 731 } 732 ierr = SNESSetFunction(snes,res,SNESTSFormFunction,ts);CHKERRQ(ierr); 733 ierr = VecDestroy(&resalloc);CHKERRQ(ierr); 734 735 PetscFunctionReturn(0); 736 } 737 738 #undef __FUNCT__ 739 #define __FUNCT__ "TSGetIFunction" 740 /*@C 741 TSGetIFunction - Returns the vector where the implicit residual is stored and the function/contex to compute it. 742 743 Not Collective 744 745 Input Parameter: 746 . ts - the TS context 747 748 Output Parameter: 749 + r - vector to hold residual (or PETSC_NULL) 750 . func - the function to compute residual (or PETSC_NULL) 751 - ctx - the function context (or PETSC_NULL) 752 753 Level: advanced 754 755 .keywords: TS, nonlinear, get, function 756 757 .seealso: TSSetIFunction(), SNESGetFunction() 758 @*/ 759 PetscErrorCode TSGetIFunction(TS ts,Vec *r,TSIFunction *func,void **ctx) 760 { 761 PetscErrorCode ierr; 762 SNES snes; 763 DM dm; 764 765 PetscFunctionBegin; 766 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 767 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 768 ierr = SNESGetFunction(snes,r,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 769 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 770 ierr = DMTSGetIFunction(dm,func,ctx);CHKERRQ(ierr); 771 PetscFunctionReturn(0); 772 } 773 774 #undef __FUNCT__ 775 #define __FUNCT__ "TSGetRHSFunction" 776 /*@C 777 TSGetRHSFunction - Returns the vector where the right hand side is stored and the function/context to compute it. 778 779 Not Collective 780 781 Input Parameter: 782 . ts - the TS context 783 784 Output Parameter: 785 + r - vector to hold computed right hand side (or PETSC_NULL) 786 . func - the function to compute right hand side (or PETSC_NULL) 787 - ctx - the function context (or PETSC_NULL) 788 789 Level: advanced 790 791 .keywords: TS, nonlinear, get, function 792 793 .seealso: TSSetRhsfunction(), SNESGetFunction() 794 @*/ 795 PetscErrorCode TSGetRHSFunction(TS ts,Vec *r,TSRHSFunction *func,void **ctx) 796 { 797 PetscErrorCode ierr; 798 SNES snes; 799 DM dm; 800 801 PetscFunctionBegin; 802 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 803 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 804 ierr = SNESGetFunction(snes,r,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 805 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 806 ierr = DMTSGetRHSFunction(dm,func,ctx);CHKERRQ(ierr); 807 PetscFunctionReturn(0); 808 } 809 810 #undef __FUNCT__ 811 #define __FUNCT__ "TSSetIJacobian" 812 /*@C 813 TSSetIJacobian - Set the function to compute the matrix dF/dU + a*dF/dU_t where F(t,U,U_t) is the function 814 you provided with TSSetIFunction(). 815 816 Logically Collective on TS 817 818 Input Parameters: 819 + ts - the TS context obtained from TSCreate() 820 . A - Jacobian matrix 821 . B - preconditioning matrix for A (may be same as A) 822 . f - the Jacobian evaluation routine 823 - ctx - user-defined context for private data for the Jacobian evaluation routine (may be PETSC_NULL) 824 825 Calling sequence of f: 826 $ f(TS ts,PetscReal t,Vec U,Vec U_t,PetscReal a,Mat *A,Mat *B,MatStructure *flag,void *ctx); 827 828 + t - time at step/stage being solved 829 . U - state vector 830 . U_t - time derivative of state vector 831 . a - shift 832 . A - Jacobian of G(U) = F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t 833 . B - preconditioning matrix for A, may be same as A 834 . flag - flag indicating information about the preconditioner matrix 835 structure (same as flag in KSPSetOperators()) 836 - ctx - [optional] user-defined context for matrix evaluation routine 837 838 Notes: 839 The matrices A and B are exactly the matrices that are used by SNES for the nonlinear solve. 840 841 The matrix dF/dU + a*dF/dU_t you provide turns out to be 842 the Jacobian of G(U) = F(t,U,W+a*U) where F(t,U,U_t) = 0 is the DAE to be solved. 843 The time integrator internally approximates U_t by W+a*U where the positive "shift" 844 a and vector W depend on the integration method, step size, and past states. For example with 845 the backward Euler method a = 1/dt and W = -a*U(previous timestep) so 846 W + a*U = a*(U - U(previous timestep)) = (U - U(previous timestep))/dt 847 848 Level: beginner 849 850 .keywords: TS, timestep, DAE, Jacobian 851 852 .seealso: TSSetIFunction(), TSSetRHSJacobian(), SNESDefaultComputeJacobianColor(), SNESDefaultComputeJacobian() 853 854 @*/ 855 PetscErrorCode TSSetIJacobian(TS ts,Mat A,Mat B,TSIJacobian f,void *ctx) 856 { 857 PetscErrorCode ierr; 858 SNES snes; 859 DM dm; 860 861 PetscFunctionBegin; 862 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 863 if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2); 864 if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3); 865 if (A) PetscCheckSameComm(ts,1,A,2); 866 if (B) PetscCheckSameComm(ts,1,B,3); 867 868 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 869 ierr = DMTSSetIJacobian(dm,f,ctx);CHKERRQ(ierr); 870 871 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 872 ierr = SNESSetJacobian(snes,A,B,SNESTSFormJacobian,ts);CHKERRQ(ierr); 873 PetscFunctionReturn(0); 874 } 875 876 #undef __FUNCT__ 877 #define __FUNCT__ "TSView" 878 /*@C 879 TSView - Prints the TS data structure. 880 881 Collective on TS 882 883 Input Parameters: 884 + ts - the TS context obtained from TSCreate() 885 - viewer - visualization context 886 887 Options Database Key: 888 . -ts_view - calls TSView() at end of TSStep() 889 890 Notes: 891 The available visualization contexts include 892 + PETSC_VIEWER_STDOUT_SELF - standard output (default) 893 - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 894 output where only the first processor opens 895 the file. All other processors send their 896 data to the first processor to print. 897 898 The user can open an alternative visualization context with 899 PetscViewerASCIIOpen() - output to a specified file. 900 901 Level: beginner 902 903 .keywords: TS, timestep, view 904 905 .seealso: PetscViewerASCIIOpen() 906 @*/ 907 PetscErrorCode TSView(TS ts,PetscViewer viewer) 908 { 909 PetscErrorCode ierr; 910 const TSType type; 911 PetscBool iascii,isstring,isundials; 912 913 PetscFunctionBegin; 914 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 915 if (!viewer) { 916 ierr = PetscViewerASCIIGetStdout(((PetscObject)ts)->comm,&viewer);CHKERRQ(ierr); 917 } 918 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 919 PetscCheckSameComm(ts,1,viewer,2); 920 921 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 922 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 923 if (iascii) { 924 ierr = PetscObjectPrintClassNamePrefixType((PetscObject)ts,viewer,"TS Object");CHKERRQ(ierr); 925 ierr = PetscViewerASCIIPrintf(viewer," maximum steps=%D\n",ts->max_steps);CHKERRQ(ierr); 926 ierr = PetscViewerASCIIPrintf(viewer," maximum time=%G\n",ts->max_time);CHKERRQ(ierr); 927 if (ts->problem_type == TS_NONLINEAR) { 928 ierr = PetscViewerASCIIPrintf(viewer," total number of nonlinear solver iterations=%D\n",ts->snes_its);CHKERRQ(ierr); 929 ierr = PetscViewerASCIIPrintf(viewer," total number of nonlinear solve failures=%D\n",ts->num_snes_failures);CHKERRQ(ierr); 930 } 931 ierr = PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%D\n",ts->ksp_its);CHKERRQ(ierr); 932 ierr = PetscViewerASCIIPrintf(viewer," total number of rejected steps=%D\n",ts->reject);CHKERRQ(ierr); 933 if (ts->ops->view) { 934 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 935 ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr); 936 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 937 } 938 } else if (isstring) { 939 ierr = TSGetType(ts,&type);CHKERRQ(ierr); 940 ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr); 941 } 942 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 943 ierr = PetscObjectTypeCompare((PetscObject)ts,TSSUNDIALS,&isundials);CHKERRQ(ierr); 944 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 945 PetscFunctionReturn(0); 946 } 947 948 949 #undef __FUNCT__ 950 #define __FUNCT__ "TSSetApplicationContext" 951 /*@ 952 TSSetApplicationContext - Sets an optional user-defined context for 953 the timesteppers. 954 955 Logically Collective on TS 956 957 Input Parameters: 958 + ts - the TS context obtained from TSCreate() 959 - usrP - optional user context 960 961 Level: intermediate 962 963 .keywords: TS, timestep, set, application, context 964 965 .seealso: TSGetApplicationContext() 966 @*/ 967 PetscErrorCode TSSetApplicationContext(TS ts,void *usrP) 968 { 969 PetscFunctionBegin; 970 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 971 ts->user = usrP; 972 PetscFunctionReturn(0); 973 } 974 975 #undef __FUNCT__ 976 #define __FUNCT__ "TSGetApplicationContext" 977 /*@ 978 TSGetApplicationContext - Gets the user-defined context for the 979 timestepper. 980 981 Not Collective 982 983 Input Parameter: 984 . ts - the TS context obtained from TSCreate() 985 986 Output Parameter: 987 . usrP - user context 988 989 Level: intermediate 990 991 .keywords: TS, timestep, get, application, context 992 993 .seealso: TSSetApplicationContext() 994 @*/ 995 PetscErrorCode TSGetApplicationContext(TS ts,void *usrP) 996 { 997 PetscFunctionBegin; 998 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 999 *(void**)usrP = ts->user; 1000 PetscFunctionReturn(0); 1001 } 1002 1003 #undef __FUNCT__ 1004 #define __FUNCT__ "TSGetTimeStepNumber" 1005 /*@ 1006 TSGetTimeStepNumber - Gets the number of time steps completed. 1007 1008 Not Collective 1009 1010 Input Parameter: 1011 . ts - the TS context obtained from TSCreate() 1012 1013 Output Parameter: 1014 . iter - number of steps completed so far 1015 1016 Level: intermediate 1017 1018 .keywords: TS, timestep, get, iteration, number 1019 .seealso: TSGetTime(), TSGetTimeStep(), TSSetPreStep(), TSSetPreStage(), TSSetPostStep() 1020 @*/ 1021 PetscErrorCode TSGetTimeStepNumber(TS ts,PetscInt* iter) 1022 { 1023 PetscFunctionBegin; 1024 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1025 PetscValidIntPointer(iter,2); 1026 *iter = ts->steps; 1027 PetscFunctionReturn(0); 1028 } 1029 1030 #undef __FUNCT__ 1031 #define __FUNCT__ "TSSetInitialTimeStep" 1032 /*@ 1033 TSSetInitialTimeStep - Sets the initial timestep to be used, 1034 as well as the initial time. 1035 1036 Logically Collective on TS 1037 1038 Input Parameters: 1039 + ts - the TS context obtained from TSCreate() 1040 . initial_time - the initial time 1041 - time_step - the size of the timestep 1042 1043 Level: intermediate 1044 1045 .seealso: TSSetTimeStep(), TSGetTimeStep() 1046 1047 .keywords: TS, set, initial, timestep 1048 @*/ 1049 PetscErrorCode TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step) 1050 { 1051 PetscErrorCode ierr; 1052 1053 PetscFunctionBegin; 1054 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1055 ierr = TSSetTimeStep(ts,time_step);CHKERRQ(ierr); 1056 ierr = TSSetTime(ts,initial_time);CHKERRQ(ierr); 1057 PetscFunctionReturn(0); 1058 } 1059 1060 #undef __FUNCT__ 1061 #define __FUNCT__ "TSSetTimeStep" 1062 /*@ 1063 TSSetTimeStep - Allows one to reset the timestep at any time, 1064 useful for simple pseudo-timestepping codes. 1065 1066 Logically Collective on TS 1067 1068 Input Parameters: 1069 + ts - the TS context obtained from TSCreate() 1070 - time_step - the size of the timestep 1071 1072 Level: intermediate 1073 1074 .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 1075 1076 .keywords: TS, set, timestep 1077 @*/ 1078 PetscErrorCode TSSetTimeStep(TS ts,PetscReal time_step) 1079 { 1080 PetscFunctionBegin; 1081 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1082 PetscValidLogicalCollectiveReal(ts,time_step,2); 1083 ts->time_step = time_step; 1084 ts->time_step_orig = time_step; 1085 PetscFunctionReturn(0); 1086 } 1087 1088 #undef __FUNCT__ 1089 #define __FUNCT__ "TSSetExactFinalTime" 1090 /*@ 1091 TSSetExactFinalTime - Determines whether to interpolate solution to the 1092 exact final time requested by the user or just returns it at the final time 1093 it computed. 1094 1095 Logically Collective on TS 1096 1097 Input Parameter: 1098 + ts - the time-step context 1099 - ft - PETSC_TRUE if interpolates, else PETSC_FALSE 1100 1101 Level: beginner 1102 1103 .seealso: TSSetDuration() 1104 @*/ 1105 PetscErrorCode TSSetExactFinalTime(TS ts,PetscBool flg) 1106 { 1107 1108 PetscFunctionBegin; 1109 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1110 PetscValidLogicalCollectiveBool(ts,flg,2); 1111 ts->exact_final_time = flg; 1112 PetscFunctionReturn(0); 1113 } 1114 1115 #undef __FUNCT__ 1116 #define __FUNCT__ "TSGetTimeStep" 1117 /*@ 1118 TSGetTimeStep - Gets the current timestep size. 1119 1120 Not Collective 1121 1122 Input Parameter: 1123 . ts - the TS context obtained from TSCreate() 1124 1125 Output Parameter: 1126 . dt - the current timestep size 1127 1128 Level: intermediate 1129 1130 .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 1131 1132 .keywords: TS, get, timestep 1133 @*/ 1134 PetscErrorCode TSGetTimeStep(TS ts,PetscReal* dt) 1135 { 1136 PetscFunctionBegin; 1137 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1138 PetscValidRealPointer(dt,2); 1139 *dt = ts->time_step; 1140 PetscFunctionReturn(0); 1141 } 1142 1143 #undef __FUNCT__ 1144 #define __FUNCT__ "TSGetSolution" 1145 /*@ 1146 TSGetSolution - Returns the solution at the present timestep. It 1147 is valid to call this routine inside the function that you are evaluating 1148 in order to move to the new timestep. This vector not changed until 1149 the solution at the next timestep has been calculated. 1150 1151 Not Collective, but Vec returned is parallel if TS is parallel 1152 1153 Input Parameter: 1154 . ts - the TS context obtained from TSCreate() 1155 1156 Output Parameter: 1157 . v - the vector containing the solution 1158 1159 Level: intermediate 1160 1161 .seealso: TSGetTimeStep() 1162 1163 .keywords: TS, timestep, get, solution 1164 @*/ 1165 PetscErrorCode TSGetSolution(TS ts,Vec *v) 1166 { 1167 PetscFunctionBegin; 1168 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1169 PetscValidPointer(v,2); 1170 *v = ts->vec_sol; 1171 PetscFunctionReturn(0); 1172 } 1173 1174 /* ----- Routines to initialize and destroy a timestepper ---- */ 1175 #undef __FUNCT__ 1176 #define __FUNCT__ "TSSetProblemType" 1177 /*@ 1178 TSSetProblemType - Sets the type of problem to be solved. 1179 1180 Not collective 1181 1182 Input Parameters: 1183 + ts - The TS 1184 - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms 1185 .vb 1186 U_t = A U 1187 U_t = A(t) U 1188 U_t = F(t,U) 1189 .ve 1190 1191 Level: beginner 1192 1193 .keywords: TS, problem type 1194 .seealso: TSSetUp(), TSProblemType, TS 1195 @*/ 1196 PetscErrorCode TSSetProblemType(TS ts, TSProblemType type) 1197 { 1198 PetscErrorCode ierr; 1199 1200 PetscFunctionBegin; 1201 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1202 ts->problem_type = type; 1203 if (type == TS_LINEAR) { 1204 SNES snes; 1205 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 1206 ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr); 1207 } 1208 PetscFunctionReturn(0); 1209 } 1210 1211 #undef __FUNCT__ 1212 #define __FUNCT__ "TSGetProblemType" 1213 /*@C 1214 TSGetProblemType - Gets the type of problem to be solved. 1215 1216 Not collective 1217 1218 Input Parameter: 1219 . ts - The TS 1220 1221 Output Parameter: 1222 . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms 1223 .vb 1224 M U_t = A U 1225 M(t) U_t = A(t) U 1226 U_t = F(t,U) 1227 .ve 1228 1229 Level: beginner 1230 1231 .keywords: TS, problem type 1232 .seealso: TSSetUp(), TSProblemType, TS 1233 @*/ 1234 PetscErrorCode TSGetProblemType(TS ts, TSProblemType *type) 1235 { 1236 PetscFunctionBegin; 1237 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1238 PetscValidIntPointer(type,2); 1239 *type = ts->problem_type; 1240 PetscFunctionReturn(0); 1241 } 1242 1243 #undef __FUNCT__ 1244 #define __FUNCT__ "TSSetUp" 1245 /*@ 1246 TSSetUp - Sets up the internal data structures for the later use 1247 of a timestepper. 1248 1249 Collective on TS 1250 1251 Input Parameter: 1252 . ts - the TS context obtained from TSCreate() 1253 1254 Notes: 1255 For basic use of the TS solvers the user need not explicitly call 1256 TSSetUp(), since these actions will automatically occur during 1257 the call to TSStep(). However, if one wishes to control this 1258 phase separately, TSSetUp() should be called after TSCreate() 1259 and optional routines of the form TSSetXXX(), but before TSStep(). 1260 1261 Level: advanced 1262 1263 .keywords: TS, timestep, setup 1264 1265 .seealso: TSCreate(), TSStep(), TSDestroy() 1266 @*/ 1267 PetscErrorCode TSSetUp(TS ts) 1268 { 1269 PetscErrorCode ierr; 1270 DM dm; 1271 PetscErrorCode (*func)(SNES,Vec,Vec,void*); 1272 PetscErrorCode (*jac)(SNES,Vec,Mat*,Mat*,MatStructure*,void*); 1273 TSIJacobian ijac; 1274 TSRHSJacobian rhsjac; 1275 1276 PetscFunctionBegin; 1277 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1278 if (ts->setupcalled) PetscFunctionReturn(0); 1279 1280 if (!((PetscObject)ts)->type_name) { 1281 ierr = TSSetType(ts,TSEULER);CHKERRQ(ierr); 1282 } 1283 if (ts->exact_final_time == PETSC_DECIDE) ts->exact_final_time = PETSC_FALSE; 1284 1285 if (!ts->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first"); 1286 1287 ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr); 1288 1289 if (ts->ops->setup) { 1290 ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr); 1291 } 1292 1293 /* in the case where we've set a DMTSFunction or what have you, we need the default SNESFunction 1294 to be set right but can't do it elsewhere due to the overreliance on ctx=ts. 1295 */ 1296 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 1297 ierr = DMSNESGetFunction(dm,&func,PETSC_NULL);CHKERRQ(ierr); 1298 if (!func) { 1299 ierr =DMSNESSetFunction(dm,SNESTSFormFunction,ts);CHKERRQ(ierr); 1300 } 1301 /* if the SNES doesn't have a jacobian set and the TS has an ijacobian or rhsjacobian set, set the SNES to use it. 1302 Otherwise, the SNES will use coloring internally to form the Jacobian. 1303 */ 1304 ierr = DMSNESGetJacobian(dm,&jac,PETSC_NULL);CHKERRQ(ierr); 1305 ierr = DMTSGetIJacobian(dm,&ijac,PETSC_NULL);CHKERRQ(ierr); 1306 ierr = DMTSGetRHSJacobian(dm,&rhsjac,PETSC_NULL);CHKERRQ(ierr); 1307 if (!jac && (ijac || rhsjac)) { 1308 ierr = DMSNESSetJacobian(dm,SNESTSFormJacobian,ts);CHKERRQ(ierr); 1309 } 1310 ts->setupcalled = PETSC_TRUE; 1311 PetscFunctionReturn(0); 1312 } 1313 1314 #undef __FUNCT__ 1315 #define __FUNCT__ "TSReset" 1316 /*@ 1317 TSReset - Resets a TS context and removes any allocated Vecs and Mats. 1318 1319 Collective on TS 1320 1321 Input Parameter: 1322 . ts - the TS context obtained from TSCreate() 1323 1324 Level: beginner 1325 1326 .keywords: TS, timestep, reset 1327 1328 .seealso: TSCreate(), TSSetup(), TSDestroy() 1329 @*/ 1330 PetscErrorCode TSReset(TS ts) 1331 { 1332 PetscErrorCode ierr; 1333 1334 PetscFunctionBegin; 1335 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1336 if (ts->ops->reset) { 1337 ierr = (*ts->ops->reset)(ts);CHKERRQ(ierr); 1338 } 1339 if (ts->snes) {ierr = SNESReset(ts->snes);CHKERRQ(ierr);} 1340 ierr = MatDestroy(&ts->Arhs);CHKERRQ(ierr); 1341 ierr = MatDestroy(&ts->Brhs);CHKERRQ(ierr); 1342 ierr = VecDestroy(&ts->Frhs);CHKERRQ(ierr); 1343 ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr); 1344 ierr = VecDestroy(&ts->vatol);CHKERRQ(ierr); 1345 ierr = VecDestroy(&ts->vrtol);CHKERRQ(ierr); 1346 ierr = VecDestroyVecs(ts->nwork,&ts->work);CHKERRQ(ierr); 1347 ts->setupcalled = PETSC_FALSE; 1348 PetscFunctionReturn(0); 1349 } 1350 1351 #undef __FUNCT__ 1352 #define __FUNCT__ "TSDestroy" 1353 /*@ 1354 TSDestroy - Destroys the timestepper context that was created 1355 with TSCreate(). 1356 1357 Collective on TS 1358 1359 Input Parameter: 1360 . ts - the TS context obtained from TSCreate() 1361 1362 Level: beginner 1363 1364 .keywords: TS, timestepper, destroy 1365 1366 .seealso: TSCreate(), TSSetUp(), TSSolve() 1367 @*/ 1368 PetscErrorCode TSDestroy(TS *ts) 1369 { 1370 PetscErrorCode ierr; 1371 1372 PetscFunctionBegin; 1373 if (!*ts) PetscFunctionReturn(0); 1374 PetscValidHeaderSpecific((*ts),TS_CLASSID,1); 1375 if (--((PetscObject)(*ts))->refct > 0) {*ts = 0; PetscFunctionReturn(0);} 1376 1377 ierr = TSReset((*ts));CHKERRQ(ierr); 1378 1379 /* if memory was published with AMS then destroy it */ 1380 ierr = PetscObjectDepublish((*ts));CHKERRQ(ierr); 1381 if ((*ts)->ops->destroy) {ierr = (*(*ts)->ops->destroy)((*ts));CHKERRQ(ierr);} 1382 1383 ierr = TSAdaptDestroy(&(*ts)->adapt);CHKERRQ(ierr); 1384 ierr = SNESDestroy(&(*ts)->snes);CHKERRQ(ierr); 1385 ierr = DMDestroy(&(*ts)->dm);CHKERRQ(ierr); 1386 ierr = TSMonitorCancel((*ts));CHKERRQ(ierr); 1387 1388 ierr = PetscHeaderDestroy(ts);CHKERRQ(ierr); 1389 PetscFunctionReturn(0); 1390 } 1391 1392 #undef __FUNCT__ 1393 #define __FUNCT__ "TSGetSNES" 1394 /*@ 1395 TSGetSNES - Returns the SNES (nonlinear solver) associated with 1396 a TS (timestepper) context. Valid only for nonlinear problems. 1397 1398 Not Collective, but SNES is parallel if TS is parallel 1399 1400 Input Parameter: 1401 . ts - the TS context obtained from TSCreate() 1402 1403 Output Parameter: 1404 . snes - the nonlinear solver context 1405 1406 Notes: 1407 The user can then directly manipulate the SNES context to set various 1408 options, etc. Likewise, the user can then extract and manipulate the 1409 KSP, KSP, and PC contexts as well. 1410 1411 TSGetSNES() does not work for integrators that do not use SNES; in 1412 this case TSGetSNES() returns PETSC_NULL in snes. 1413 1414 Level: beginner 1415 1416 .keywords: timestep, get, SNES 1417 @*/ 1418 PetscErrorCode TSGetSNES(TS ts,SNES *snes) 1419 { 1420 PetscErrorCode ierr; 1421 1422 PetscFunctionBegin; 1423 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1424 PetscValidPointer(snes,2); 1425 if (!ts->snes) { 1426 ierr = SNESCreate(((PetscObject)ts)->comm,&ts->snes);CHKERRQ(ierr); 1427 ierr = SNESSetFunction(ts->snes,PETSC_NULL,SNESTSFormFunction,ts);CHKERRQ(ierr); 1428 ierr = PetscLogObjectParent(ts,ts->snes);CHKERRQ(ierr); 1429 ierr = PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);CHKERRQ(ierr); 1430 if (ts->dm) {ierr = SNESSetDM(ts->snes,ts->dm);CHKERRQ(ierr);} 1431 if (ts->problem_type == TS_LINEAR) { 1432 ierr = SNESSetType(ts->snes,SNESKSPONLY);CHKERRQ(ierr); 1433 } 1434 } 1435 *snes = ts->snes; 1436 PetscFunctionReturn(0); 1437 } 1438 1439 #undef __FUNCT__ 1440 #define __FUNCT__ "TSGetKSP" 1441 /*@ 1442 TSGetKSP - Returns the KSP (linear solver) associated with 1443 a TS (timestepper) context. 1444 1445 Not Collective, but KSP is parallel if TS is parallel 1446 1447 Input Parameter: 1448 . ts - the TS context obtained from TSCreate() 1449 1450 Output Parameter: 1451 . ksp - the nonlinear solver context 1452 1453 Notes: 1454 The user can then directly manipulate the KSP context to set various 1455 options, etc. Likewise, the user can then extract and manipulate the 1456 KSP and PC contexts as well. 1457 1458 TSGetKSP() does not work for integrators that do not use KSP; 1459 in this case TSGetKSP() returns PETSC_NULL in ksp. 1460 1461 Level: beginner 1462 1463 .keywords: timestep, get, KSP 1464 @*/ 1465 PetscErrorCode TSGetKSP(TS ts,KSP *ksp) 1466 { 1467 PetscErrorCode ierr; 1468 SNES snes; 1469 1470 PetscFunctionBegin; 1471 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1472 PetscValidPointer(ksp,2); 1473 if (!((PetscObject)ts)->type_name) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"KSP is not created yet. Call TSSetType() first"); 1474 if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()"); 1475 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 1476 ierr = SNESGetKSP(snes,ksp);CHKERRQ(ierr); 1477 PetscFunctionReturn(0); 1478 } 1479 1480 /* ----------- Routines to set solver parameters ---------- */ 1481 1482 #undef __FUNCT__ 1483 #define __FUNCT__ "TSGetDuration" 1484 /*@ 1485 TSGetDuration - Gets the maximum number of timesteps to use and 1486 maximum time for iteration. 1487 1488 Not Collective 1489 1490 Input Parameters: 1491 + ts - the TS context obtained from TSCreate() 1492 . maxsteps - maximum number of iterations to use, or PETSC_NULL 1493 - maxtime - final time to iterate to, or PETSC_NULL 1494 1495 Level: intermediate 1496 1497 .keywords: TS, timestep, get, maximum, iterations, time 1498 @*/ 1499 PetscErrorCode TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime) 1500 { 1501 PetscFunctionBegin; 1502 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1503 if (maxsteps) { 1504 PetscValidIntPointer(maxsteps,2); 1505 *maxsteps = ts->max_steps; 1506 } 1507 if (maxtime) { 1508 PetscValidScalarPointer(maxtime,3); 1509 *maxtime = ts->max_time; 1510 } 1511 PetscFunctionReturn(0); 1512 } 1513 1514 #undef __FUNCT__ 1515 #define __FUNCT__ "TSSetDuration" 1516 /*@ 1517 TSSetDuration - Sets the maximum number of timesteps to use and 1518 maximum time for iteration. 1519 1520 Logically Collective on TS 1521 1522 Input Parameters: 1523 + ts - the TS context obtained from TSCreate() 1524 . maxsteps - maximum number of iterations to use 1525 - maxtime - final time to iterate to 1526 1527 Options Database Keys: 1528 . -ts_max_steps <maxsteps> - Sets maxsteps 1529 . -ts_final_time <maxtime> - Sets maxtime 1530 1531 Notes: 1532 The default maximum number of iterations is 5000. Default time is 5.0 1533 1534 Level: intermediate 1535 1536 .keywords: TS, timestep, set, maximum, iterations 1537 1538 .seealso: TSSetExactFinalTime() 1539 @*/ 1540 PetscErrorCode TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime) 1541 { 1542 PetscFunctionBegin; 1543 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1544 PetscValidLogicalCollectiveInt(ts,maxsteps,2); 1545 PetscValidLogicalCollectiveReal(ts,maxtime,2); 1546 if (maxsteps >= 0) ts->max_steps = maxsteps; 1547 if (maxtime != PETSC_DEFAULT) ts->max_time = maxtime; 1548 PetscFunctionReturn(0); 1549 } 1550 1551 #undef __FUNCT__ 1552 #define __FUNCT__ "TSSetSolution" 1553 /*@ 1554 TSSetSolution - Sets the initial solution vector 1555 for use by the TS routines. 1556 1557 Logically Collective on TS and Vec 1558 1559 Input Parameters: 1560 + ts - the TS context obtained from TSCreate() 1561 - x - the solution vector 1562 1563 Level: beginner 1564 1565 .keywords: TS, timestep, set, solution, initial conditions 1566 @*/ 1567 PetscErrorCode TSSetSolution(TS ts,Vec x) 1568 { 1569 PetscErrorCode ierr; 1570 DM dm; 1571 1572 PetscFunctionBegin; 1573 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1574 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 1575 ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr); 1576 ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr); 1577 ts->vec_sol = x; 1578 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 1579 ierr = DMShellSetGlobalVector(dm,x);CHKERRQ(ierr); 1580 PetscFunctionReturn(0); 1581 } 1582 1583 #undef __FUNCT__ 1584 #define __FUNCT__ "TSSetPreStep" 1585 /*@C 1586 TSSetPreStep - Sets the general-purpose function 1587 called once at the beginning of each time step. 1588 1589 Logically Collective on TS 1590 1591 Input Parameters: 1592 + ts - The TS context obtained from TSCreate() 1593 - func - The function 1594 1595 Calling sequence of func: 1596 . func (TS ts); 1597 1598 Level: intermediate 1599 1600 Note: 1601 If a step is rejected, TSStep() will call this routine again before each attempt. 1602 The last completed time step number can be queried using TSGetTimeStepNumber(), the 1603 size of the step being attempted can be obtained using TSGetTimeStep(). 1604 1605 .keywords: TS, timestep 1606 .seealso: TSSetPreStage(), TSSetPostStep(), TSStep() 1607 @*/ 1608 PetscErrorCode TSSetPreStep(TS ts, PetscErrorCode (*func)(TS)) 1609 { 1610 PetscFunctionBegin; 1611 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1612 ts->ops->prestep = func; 1613 PetscFunctionReturn(0); 1614 } 1615 1616 #undef __FUNCT__ 1617 #define __FUNCT__ "TSPreStep" 1618 /*@ 1619 TSPreStep - Runs the user-defined pre-step function. 1620 1621 Collective on TS 1622 1623 Input Parameters: 1624 . ts - The TS context obtained from TSCreate() 1625 1626 Notes: 1627 TSPreStep() is typically used within time stepping implementations, 1628 so most users would not generally call this routine themselves. 1629 1630 Level: developer 1631 1632 .keywords: TS, timestep 1633 .seealso: TSSetPreStep(), TSPreStage(), TSPostStep() 1634 @*/ 1635 PetscErrorCode TSPreStep(TS ts) 1636 { 1637 PetscErrorCode ierr; 1638 1639 PetscFunctionBegin; 1640 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1641 if (ts->ops->prestep) { 1642 PetscStackPush("TS PreStep function"); 1643 ierr = (*ts->ops->prestep)(ts);CHKERRQ(ierr); 1644 PetscStackPop; 1645 } 1646 PetscFunctionReturn(0); 1647 } 1648 1649 #undef __FUNCT__ 1650 #define __FUNCT__ "TSSetPreStage" 1651 /*@C 1652 TSSetPreStage - Sets the general-purpose function 1653 called once at the beginning of each stage. 1654 1655 Logically Collective on TS 1656 1657 Input Parameters: 1658 + ts - The TS context obtained from TSCreate() 1659 - func - The function 1660 1661 Calling sequence of func: 1662 . PetscErrorCode func(TS ts, PetscReal stagetime); 1663 1664 Level: intermediate 1665 1666 Note: 1667 There may be several stages per time step. If the solve for a given stage fails, the step may be rejected and retried. 1668 The time step number being computed can be queried using TSGetTimeStepNumber() and the total size of the step being 1669 attempted can be obtained using TSGetTimeStep(). The time at the start of the step is available via TSGetTime(). 1670 1671 .keywords: TS, timestep 1672 .seealso: TSSetPreStep(), TSSetPostStep(), TSGetApplicationContext() 1673 @*/ 1674 PetscErrorCode TSSetPreStage(TS ts, PetscErrorCode (*func)(TS,PetscReal)) 1675 { 1676 PetscFunctionBegin; 1677 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1678 ts->ops->prestage = func; 1679 PetscFunctionReturn(0); 1680 } 1681 1682 #undef __FUNCT__ 1683 #define __FUNCT__ "TSPreStage" 1684 /*@ 1685 TSPreStage - Runs the user-defined pre-stage function set using TSSetPreStage() 1686 1687 Collective on TS 1688 1689 Input Parameters: 1690 . ts - The TS context obtained from TSCreate() 1691 1692 Notes: 1693 TSPreStage() is typically used within time stepping implementations, 1694 most users would not generally call this routine themselves. 1695 1696 Level: developer 1697 1698 .keywords: TS, timestep 1699 .seealso: TSSetPreStep(), TSPreStep(), TSPostStep() 1700 @*/ 1701 PetscErrorCode TSPreStage(TS ts, PetscReal stagetime) 1702 { 1703 PetscErrorCode ierr; 1704 1705 PetscFunctionBegin; 1706 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1707 if (ts->ops->prestage) { 1708 PetscStackPush("TS PreStage function"); 1709 ierr = (*ts->ops->prestage)(ts,stagetime);CHKERRQ(ierr); 1710 PetscStackPop; 1711 } 1712 PetscFunctionReturn(0); 1713 } 1714 1715 #undef __FUNCT__ 1716 #define __FUNCT__ "TSSetPostStep" 1717 /*@C 1718 TSSetPostStep - Sets the general-purpose function 1719 called once at the end of each time step. 1720 1721 Logically Collective on TS 1722 1723 Input Parameters: 1724 + ts - The TS context obtained from TSCreate() 1725 - func - The function 1726 1727 Calling sequence of func: 1728 $ func (TS ts); 1729 1730 Level: intermediate 1731 1732 .keywords: TS, timestep 1733 .seealso: TSSetPreStep(), TSSetPreStage(), TSGetTimeStep(), TSGetTimeStepNumber(), TSGetTime() 1734 @*/ 1735 PetscErrorCode TSSetPostStep(TS ts, PetscErrorCode (*func)(TS)) 1736 { 1737 PetscFunctionBegin; 1738 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1739 ts->ops->poststep = func; 1740 PetscFunctionReturn(0); 1741 } 1742 1743 #undef __FUNCT__ 1744 #define __FUNCT__ "TSPostStep" 1745 /*@ 1746 TSPostStep - Runs the user-defined post-step function. 1747 1748 Collective on TS 1749 1750 Input Parameters: 1751 . ts - The TS context obtained from TSCreate() 1752 1753 Notes: 1754 TSPostStep() is typically used within time stepping implementations, 1755 so most users would not generally call this routine themselves. 1756 1757 Level: developer 1758 1759 .keywords: TS, timestep 1760 @*/ 1761 PetscErrorCode TSPostStep(TS ts) 1762 { 1763 PetscErrorCode ierr; 1764 1765 PetscFunctionBegin; 1766 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1767 if (ts->ops->poststep) { 1768 PetscStackPush("TS PostStep function"); 1769 ierr = (*ts->ops->poststep)(ts);CHKERRQ(ierr); 1770 PetscStackPop; 1771 } 1772 PetscFunctionReturn(0); 1773 } 1774 1775 /* ------------ Routines to set performance monitoring options ----------- */ 1776 1777 #undef __FUNCT__ 1778 #define __FUNCT__ "TSMonitorSet" 1779 /*@C 1780 TSMonitorSet - Sets an ADDITIONAL function that is to be used at every 1781 timestep to display the iteration's progress. 1782 1783 Logically Collective on TS 1784 1785 Input Parameters: 1786 + ts - the TS context obtained from TSCreate() 1787 . monitor - monitoring routine 1788 . mctx - [optional] user-defined context for private data for the 1789 monitor routine (use PETSC_NULL if no context is desired) 1790 - monitordestroy - [optional] routine that frees monitor context 1791 (may be PETSC_NULL) 1792 1793 Calling sequence of monitor: 1794 $ int monitor(TS ts,PetscInt steps,PetscReal time,Vec x,void *mctx) 1795 1796 + ts - the TS context 1797 . steps - iteration number 1798 . time - current time 1799 . x - current iterate 1800 - mctx - [optional] monitoring context 1801 1802 Notes: 1803 This routine adds an additional monitor to the list of monitors that 1804 already has been loaded. 1805 1806 Fortran notes: Only a single monitor function can be set for each TS object 1807 1808 Level: intermediate 1809 1810 .keywords: TS, timestep, set, monitor 1811 1812 .seealso: TSMonitorDefault(), TSMonitorCancel() 1813 @*/ 1814 PetscErrorCode TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void**)) 1815 { 1816 PetscFunctionBegin; 1817 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1818 if (ts->numbermonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1819 ts->monitor[ts->numbermonitors] = monitor; 1820 ts->mdestroy[ts->numbermonitors] = mdestroy; 1821 ts->monitorcontext[ts->numbermonitors++] = (void*)mctx; 1822 PetscFunctionReturn(0); 1823 } 1824 1825 #undef __FUNCT__ 1826 #define __FUNCT__ "TSMonitorCancel" 1827 /*@C 1828 TSMonitorCancel - Clears all the monitors that have been set on a time-step object. 1829 1830 Logically Collective on TS 1831 1832 Input Parameters: 1833 . ts - the TS context obtained from TSCreate() 1834 1835 Notes: 1836 There is no way to remove a single, specific monitor. 1837 1838 Level: intermediate 1839 1840 .keywords: TS, timestep, set, monitor 1841 1842 .seealso: TSMonitorDefault(), TSMonitorSet() 1843 @*/ 1844 PetscErrorCode TSMonitorCancel(TS ts) 1845 { 1846 PetscErrorCode ierr; 1847 PetscInt i; 1848 1849 PetscFunctionBegin; 1850 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1851 for (i=0; i<ts->numbermonitors; i++) { 1852 if (ts->mdestroy[i]) { 1853 ierr = (*ts->mdestroy[i])(&ts->monitorcontext[i]);CHKERRQ(ierr); 1854 } 1855 } 1856 ts->numbermonitors = 0; 1857 PetscFunctionReturn(0); 1858 } 1859 1860 #undef __FUNCT__ 1861 #define __FUNCT__ "TSMonitorDefault" 1862 /*@ 1863 TSMonitorDefault - Sets the Default monitor 1864 1865 Level: intermediate 1866 1867 .keywords: TS, set, monitor 1868 1869 .seealso: TSMonitorDefault(), TSMonitorSet() 1870 @*/ 1871 PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy) 1872 { 1873 PetscErrorCode ierr; 1874 PetscViewer viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ts)->comm); 1875 1876 PetscFunctionBegin; 1877 ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 1878 ierr = PetscViewerASCIIPrintf(viewer,"%D TS dt %g time %g\n",step,(double)ts->time_step,(double)ptime);CHKERRQ(ierr); 1879 ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 1880 PetscFunctionReturn(0); 1881 } 1882 1883 #undef __FUNCT__ 1884 #define __FUNCT__ "TSSetRetainStages" 1885 /*@ 1886 TSSetRetainStages - Request that all stages in the upcoming step be stored so that interpolation will be available. 1887 1888 Logically Collective on TS 1889 1890 Input Argument: 1891 . ts - time stepping context 1892 1893 Output Argument: 1894 . flg - PETSC_TRUE or PETSC_FALSE 1895 1896 Level: intermediate 1897 1898 .keywords: TS, set 1899 1900 .seealso: TSInterpolate(), TSSetPostStep() 1901 @*/ 1902 PetscErrorCode TSSetRetainStages(TS ts,PetscBool flg) 1903 { 1904 1905 PetscFunctionBegin; 1906 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1907 ts->retain_stages = flg; 1908 PetscFunctionReturn(0); 1909 } 1910 1911 #undef __FUNCT__ 1912 #define __FUNCT__ "TSInterpolate" 1913 /*@ 1914 TSInterpolate - Interpolate the solution computed during the previous step to an arbitrary location in the interval 1915 1916 Collective on TS 1917 1918 Input Argument: 1919 + ts - time stepping context 1920 - t - time to interpolate to 1921 1922 Output Argument: 1923 . X - state at given time 1924 1925 Notes: 1926 The user should call TSSetRetainStages() before taking a step in which interpolation will be requested. 1927 1928 Level: intermediate 1929 1930 Developer Notes: 1931 TSInterpolate() and the storing of previous steps/stages should be generalized to support delay differential equations and continuous adjoints. 1932 1933 .keywords: TS, set 1934 1935 .seealso: TSSetRetainStages(), TSSetPostStep() 1936 @*/ 1937 PetscErrorCode TSInterpolate(TS ts,PetscReal t,Vec X) 1938 { 1939 PetscErrorCode ierr; 1940 1941 PetscFunctionBegin; 1942 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1943 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); 1944 if (!ts->ops->interpolate) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"%s does not provide interpolation",((PetscObject)ts)->type_name); 1945 ierr = (*ts->ops->interpolate)(ts,t,X);CHKERRQ(ierr); 1946 PetscFunctionReturn(0); 1947 } 1948 1949 #undef __FUNCT__ 1950 #define __FUNCT__ "TSStep" 1951 /*@ 1952 TSStep - Steps one time step 1953 1954 Collective on TS 1955 1956 Input Parameter: 1957 . ts - the TS context obtained from TSCreate() 1958 1959 Level: intermediate 1960 1961 Notes: 1962 The hook set using TSSetPreStep() is called before each attempt to take the step. In general, the time step size may 1963 be changed due to adaptive error controller or solve failures. Note that steps may contain multiple stages. 1964 1965 .keywords: TS, timestep, solve 1966 1967 .seealso: TSCreate(), TSSetUp(), TSDestroy(), TSSolve(), TSSetPreStep(), TSSetPreStage() 1968 @*/ 1969 PetscErrorCode TSStep(TS ts) 1970 { 1971 PetscReal ptime_prev; 1972 PetscErrorCode ierr; 1973 1974 PetscFunctionBegin; 1975 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1976 ierr = TSSetUp(ts);CHKERRQ(ierr); 1977 1978 ts->reason = TS_CONVERGED_ITERATING; 1979 1980 ptime_prev = ts->ptime; 1981 ierr = PetscLogEventBegin(TS_Step,ts,0,0,0);CHKERRQ(ierr); 1982 ierr = (*ts->ops->step)(ts);CHKERRQ(ierr); 1983 ierr = PetscLogEventEnd(TS_Step,ts,0,0,0);CHKERRQ(ierr); 1984 ts->time_step_prev = ts->ptime - ptime_prev; 1985 1986 if (ts->reason < 0) { 1987 if (ts->errorifstepfailed) { 1988 if (ts->reason == TS_DIVERGED_NONLINEAR_SOLVE) { 1989 SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s, increase -ts_max_snes_failures or make negative to attempt recovery",TSConvergedReasons[ts->reason]); 1990 } else SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s",TSConvergedReasons[ts->reason]); 1991 } 1992 } else if (!ts->reason) { 1993 if (ts->steps >= ts->max_steps) 1994 ts->reason = TS_CONVERGED_ITS; 1995 else if (ts->ptime >= ts->max_time) 1996 ts->reason = TS_CONVERGED_TIME; 1997 } 1998 1999 PetscFunctionReturn(0); 2000 } 2001 2002 #undef __FUNCT__ 2003 #define __FUNCT__ "TSEvaluateStep" 2004 /*@ 2005 TSEvaluateStep - Evaluate the solution at the end of a time step with a given order of accuracy. 2006 2007 Collective on TS 2008 2009 Input Arguments: 2010 + ts - time stepping context 2011 . order - desired order of accuracy 2012 - done - whether the step was evaluated at this order (pass PETSC_NULL to generate an error if not available) 2013 2014 Output Arguments: 2015 . X - state at the end of the current step 2016 2017 Level: advanced 2018 2019 Notes: 2020 This function cannot be called until all stages have been evaluated. 2021 It is normally called by adaptive controllers before a step has been accepted and may also be called by the user after TSStep() has returned. 2022 2023 .seealso: TSStep(), TSAdapt 2024 @*/ 2025 PetscErrorCode TSEvaluateStep(TS ts,PetscInt order,Vec X,PetscBool *done) 2026 { 2027 PetscErrorCode ierr; 2028 2029 PetscFunctionBegin; 2030 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2031 PetscValidType(ts,1); 2032 PetscValidHeaderSpecific(X,VEC_CLASSID,3); 2033 if (!ts->ops->evaluatestep) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSEvaluateStep not implemented for type '%s'",((PetscObject)ts)->type_name); 2034 ierr = (*ts->ops->evaluatestep)(ts,order,X,done);CHKERRQ(ierr); 2035 PetscFunctionReturn(0); 2036 } 2037 2038 #undef __FUNCT__ 2039 #define __FUNCT__ "TSSolve" 2040 /*@ 2041 TSSolve - Steps the requested number of timesteps. 2042 2043 Collective on TS 2044 2045 Input Parameter: 2046 + ts - the TS context obtained from TSCreate() 2047 - x - the solution vector 2048 2049 Output Parameter: 2050 . ftime - time of the state vector x upon completion 2051 2052 Level: beginner 2053 2054 Notes: 2055 The final time returned by this function may be different from the time of the internally 2056 held state accessible by TSGetSolution() and TSGetTime() because the method may have 2057 stepped over the final time. 2058 2059 .keywords: TS, timestep, solve 2060 2061 .seealso: TSCreate(), TSSetSolution(), TSStep() 2062 @*/ 2063 PetscErrorCode TSSolve(TS ts,Vec x,PetscReal *ftime) 2064 { 2065 PetscBool flg; 2066 char filename[PETSC_MAX_PATH_LEN]; 2067 PetscViewer viewer; 2068 PetscErrorCode ierr; 2069 2070 PetscFunctionBegin; 2071 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2072 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 2073 if (ts->exact_final_time) { /* Need ts->vec_sol to be distinct so it is not overwritten when we interpolate at the end */ 2074 if (!ts->vec_sol || x == ts->vec_sol) { 2075 Vec y; 2076 ierr = VecDuplicate(x,&y);CHKERRQ(ierr); 2077 ierr = TSSetSolution(ts,y);CHKERRQ(ierr); 2078 ierr = VecDestroy(&y);CHKERRQ(ierr); /* grant ownership */ 2079 } 2080 ierr = VecCopy(x,ts->vec_sol);CHKERRQ(ierr); 2081 } else { 2082 ierr = TSSetSolution(ts,x);CHKERRQ(ierr); 2083 } 2084 ierr = TSSetUp(ts);CHKERRQ(ierr); 2085 /* reset time step and iteration counters */ 2086 ts->steps = 0; 2087 ts->ksp_its = 0; 2088 ts->snes_its = 0; 2089 ts->num_snes_failures = 0; 2090 ts->reject = 0; 2091 ts->reason = TS_CONVERGED_ITERATING; 2092 2093 if (ts->ops->solve) { /* This private interface is transitional and should be removed when all implementations are updated. */ 2094 ierr = (*ts->ops->solve)(ts);CHKERRQ(ierr); 2095 ierr = VecCopy(ts->vec_sol,x);CHKERRQ(ierr); 2096 if (ftime) *ftime = ts->ptime; 2097 } else { 2098 /* steps the requested number of timesteps. */ 2099 ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 2100 if (ts->steps >= ts->max_steps) 2101 ts->reason = TS_CONVERGED_ITS; 2102 else if (ts->ptime >= ts->max_time) 2103 ts->reason = TS_CONVERGED_TIME; 2104 while (!ts->reason) { 2105 ierr = TSStep(ts);CHKERRQ(ierr); 2106 ierr = TSPostStep(ts);CHKERRQ(ierr); 2107 ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 2108 } 2109 if (ts->exact_final_time && ts->ptime > ts->max_time) { 2110 ierr = TSInterpolate(ts,ts->max_time,x);CHKERRQ(ierr); 2111 if (ftime) *ftime = ts->max_time; 2112 } else { 2113 ierr = VecCopy(ts->vec_sol,x);CHKERRQ(ierr); 2114 if (ftime) *ftime = ts->ptime; 2115 } 2116 } 2117 ierr = PetscOptionsGetString(((PetscObject)ts)->prefix,"-ts_view",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 2118 if (flg && !PetscPreLoadingOn) { 2119 ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,filename,&viewer);CHKERRQ(ierr); 2120 ierr = TSView(ts,viewer);CHKERRQ(ierr); 2121 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 2122 } 2123 PetscFunctionReturn(0); 2124 } 2125 2126 #undef __FUNCT__ 2127 #define __FUNCT__ "TSMonitor" 2128 /*@ 2129 TSMonitor - Runs all user-provided monitor routines set using TSMonitorSet() 2130 2131 Collective on TS 2132 2133 Input Parameters: 2134 + ts - time stepping context obtained from TSCreate() 2135 . step - step number that has just completed 2136 . ptime - model time of the state 2137 - x - state at the current model time 2138 2139 Notes: 2140 TSMonitor() is typically used within the time stepping implementations. 2141 Users might call this function when using the TSStep() interface instead of TSSolve(). 2142 2143 Level: advanced 2144 2145 .keywords: TS, timestep 2146 @*/ 2147 PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x) 2148 { 2149 PetscErrorCode ierr; 2150 PetscInt i,n = ts->numbermonitors; 2151 2152 PetscFunctionBegin; 2153 for (i=0; i<n; i++) { 2154 ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr); 2155 } 2156 PetscFunctionReturn(0); 2157 } 2158 2159 /* ------------------------------------------------------------------------*/ 2160 2161 #undef __FUNCT__ 2162 #define __FUNCT__ "TSMonitorLGCreate" 2163 /*@C 2164 TSMonitorLGCreate - Creates a line graph context for use with 2165 TS to monitor convergence of preconditioned residual norms. 2166 2167 Collective on TS 2168 2169 Input Parameters: 2170 + host - the X display to open, or null for the local machine 2171 . label - the title to put in the title bar 2172 . x, y - the screen coordinates of the upper left coordinate of the window 2173 - m, n - the screen width and height in pixels 2174 2175 Output Parameter: 2176 . draw - the drawing context 2177 2178 Options Database Key: 2179 . -ts_monitor_draw - automatically sets line graph monitor 2180 2181 Notes: 2182 Use TSMonitorLGDestroy() to destroy this line graph, not PetscDrawLGDestroy(). 2183 2184 Level: intermediate 2185 2186 .keywords: TS, monitor, line graph, residual, seealso 2187 2188 .seealso: TSMonitorLGDestroy(), TSMonitorSet() 2189 2190 @*/ 2191 PetscErrorCode TSMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 2192 { 2193 PetscDraw win; 2194 PetscErrorCode ierr; 2195 2196 PetscFunctionBegin; 2197 ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr); 2198 ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr); 2199 ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr); 2200 ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr); 2201 2202 ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr); 2203 PetscFunctionReturn(0); 2204 } 2205 2206 #undef __FUNCT__ 2207 #define __FUNCT__ "TSMonitorLG" 2208 PetscErrorCode TSMonitorLG(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx) 2209 { 2210 PetscDrawLG lg = (PetscDrawLG) monctx; 2211 PetscReal x,y = ptime; 2212 PetscErrorCode ierr; 2213 2214 PetscFunctionBegin; 2215 if (!monctx) { 2216 MPI_Comm comm; 2217 PetscViewer viewer; 2218 2219 ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr); 2220 viewer = PETSC_VIEWER_DRAW_(comm); 2221 ierr = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr); 2222 } 2223 2224 if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 2225 x = (PetscReal)n; 2226 ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 2227 if (n < 20 || (n % 5)) { 2228 ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 2229 } 2230 PetscFunctionReturn(0); 2231 } 2232 2233 #undef __FUNCT__ 2234 #define __FUNCT__ "TSMonitorLGDestroy" 2235 /*@C 2236 TSMonitorLGDestroy - Destroys a line graph context that was created 2237 with TSMonitorLGCreate(). 2238 2239 Collective on PetscDrawLG 2240 2241 Input Parameter: 2242 . draw - the drawing context 2243 2244 Level: intermediate 2245 2246 .keywords: TS, monitor, line graph, destroy 2247 2248 .seealso: TSMonitorLGCreate(), TSMonitorSet(), TSMonitorLG(); 2249 @*/ 2250 PetscErrorCode TSMonitorLGDestroy(PetscDrawLG *drawlg) 2251 { 2252 PetscDraw draw; 2253 PetscErrorCode ierr; 2254 2255 PetscFunctionBegin; 2256 ierr = PetscDrawLGGetDraw(*drawlg,&draw);CHKERRQ(ierr); 2257 ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr); 2258 ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr); 2259 PetscFunctionReturn(0); 2260 } 2261 2262 #undef __FUNCT__ 2263 #define __FUNCT__ "TSGetTime" 2264 /*@ 2265 TSGetTime - Gets the time of the most recently completed step. 2266 2267 Not Collective 2268 2269 Input Parameter: 2270 . ts - the TS context obtained from TSCreate() 2271 2272 Output Parameter: 2273 . t - the current time 2274 2275 Level: beginner 2276 2277 Note: 2278 When called during time step evaluation (e.g. during residual evaluation or via hooks set using TSSetPreStep(), 2279 TSSetPreStage(), or TSSetPostStep()), the time is the time at the start of the step being evaluated. 2280 2281 .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 2282 2283 .keywords: TS, get, time 2284 @*/ 2285 PetscErrorCode TSGetTime(TS ts,PetscReal* t) 2286 { 2287 PetscFunctionBegin; 2288 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2289 PetscValidRealPointer(t,2); 2290 *t = ts->ptime; 2291 PetscFunctionReturn(0); 2292 } 2293 2294 #undef __FUNCT__ 2295 #define __FUNCT__ "TSSetTime" 2296 /*@ 2297 TSSetTime - Allows one to reset the time. 2298 2299 Logically Collective on TS 2300 2301 Input Parameters: 2302 + ts - the TS context obtained from TSCreate() 2303 - time - the time 2304 2305 Level: intermediate 2306 2307 .seealso: TSGetTime(), TSSetDuration() 2308 2309 .keywords: TS, set, time 2310 @*/ 2311 PetscErrorCode TSSetTime(TS ts, PetscReal t) 2312 { 2313 PetscFunctionBegin; 2314 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2315 PetscValidLogicalCollectiveReal(ts,t,2); 2316 ts->ptime = t; 2317 PetscFunctionReturn(0); 2318 } 2319 2320 #undef __FUNCT__ 2321 #define __FUNCT__ "TSSetOptionsPrefix" 2322 /*@C 2323 TSSetOptionsPrefix - Sets the prefix used for searching for all 2324 TS options in the database. 2325 2326 Logically Collective on TS 2327 2328 Input Parameter: 2329 + ts - The TS context 2330 - prefix - The prefix to prepend to all option names 2331 2332 Notes: 2333 A hyphen (-) must NOT be given at the beginning of the prefix name. 2334 The first character of all runtime options is AUTOMATICALLY the 2335 hyphen. 2336 2337 Level: advanced 2338 2339 .keywords: TS, set, options, prefix, database 2340 2341 .seealso: TSSetFromOptions() 2342 2343 @*/ 2344 PetscErrorCode TSSetOptionsPrefix(TS ts,const char prefix[]) 2345 { 2346 PetscErrorCode ierr; 2347 SNES snes; 2348 2349 PetscFunctionBegin; 2350 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2351 ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 2352 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2353 ierr = SNESSetOptionsPrefix(snes,prefix);CHKERRQ(ierr); 2354 PetscFunctionReturn(0); 2355 } 2356 2357 2358 #undef __FUNCT__ 2359 #define __FUNCT__ "TSAppendOptionsPrefix" 2360 /*@C 2361 TSAppendOptionsPrefix - Appends to the prefix used for searching for all 2362 TS options in the database. 2363 2364 Logically Collective on TS 2365 2366 Input Parameter: 2367 + ts - The TS context 2368 - prefix - The prefix to prepend to all option names 2369 2370 Notes: 2371 A hyphen (-) must NOT be given at the beginning of the prefix name. 2372 The first character of all runtime options is AUTOMATICALLY the 2373 hyphen. 2374 2375 Level: advanced 2376 2377 .keywords: TS, append, options, prefix, database 2378 2379 .seealso: TSGetOptionsPrefix() 2380 2381 @*/ 2382 PetscErrorCode TSAppendOptionsPrefix(TS ts,const char prefix[]) 2383 { 2384 PetscErrorCode ierr; 2385 SNES snes; 2386 2387 PetscFunctionBegin; 2388 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2389 ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 2390 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2391 ierr = SNESAppendOptionsPrefix(snes,prefix);CHKERRQ(ierr); 2392 PetscFunctionReturn(0); 2393 } 2394 2395 #undef __FUNCT__ 2396 #define __FUNCT__ "TSGetOptionsPrefix" 2397 /*@C 2398 TSGetOptionsPrefix - Sets the prefix used for searching for all 2399 TS options in the database. 2400 2401 Not Collective 2402 2403 Input Parameter: 2404 . ts - The TS context 2405 2406 Output Parameter: 2407 . prefix - A pointer to the prefix string used 2408 2409 Notes: On the fortran side, the user should pass in a string 'prifix' of 2410 sufficient length to hold the prefix. 2411 2412 Level: intermediate 2413 2414 .keywords: TS, get, options, prefix, database 2415 2416 .seealso: TSAppendOptionsPrefix() 2417 @*/ 2418 PetscErrorCode TSGetOptionsPrefix(TS ts,const char *prefix[]) 2419 { 2420 PetscErrorCode ierr; 2421 2422 PetscFunctionBegin; 2423 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2424 PetscValidPointer(prefix,2); 2425 ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 2426 PetscFunctionReturn(0); 2427 } 2428 2429 #undef __FUNCT__ 2430 #define __FUNCT__ "TSGetRHSJacobian" 2431 /*@C 2432 TSGetRHSJacobian - Returns the Jacobian J at the present timestep. 2433 2434 Not Collective, but parallel objects are returned if TS is parallel 2435 2436 Input Parameter: 2437 . ts - The TS context obtained from TSCreate() 2438 2439 Output Parameters: 2440 + J - The Jacobian J of F, where U_t = F(U,t) 2441 . M - The preconditioner matrix, usually the same as J 2442 . func - Function to compute the Jacobian of the RHS 2443 - ctx - User-defined context for Jacobian evaluation routine 2444 2445 Notes: You can pass in PETSC_NULL for any return argument you do not need. 2446 2447 Level: intermediate 2448 2449 .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber() 2450 2451 .keywords: TS, timestep, get, matrix, Jacobian 2452 @*/ 2453 PetscErrorCode TSGetRHSJacobian(TS ts,Mat *J,Mat *M,TSRHSJacobian *func,void **ctx) 2454 { 2455 PetscErrorCode ierr; 2456 SNES snes; 2457 DM dm; 2458 2459 PetscFunctionBegin; 2460 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2461 ierr = SNESGetJacobian(snes,J,M,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 2462 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 2463 ierr = DMTSGetRHSJacobian(dm,func,ctx);CHKERRQ(ierr); 2464 PetscFunctionReturn(0); 2465 } 2466 2467 #undef __FUNCT__ 2468 #define __FUNCT__ "TSGetIJacobian" 2469 /*@C 2470 TSGetIJacobian - Returns the implicit Jacobian at the present timestep. 2471 2472 Not Collective, but parallel objects are returned if TS is parallel 2473 2474 Input Parameter: 2475 . ts - The TS context obtained from TSCreate() 2476 2477 Output Parameters: 2478 + A - The Jacobian of F(t,U,U_t) 2479 . B - The preconditioner matrix, often the same as A 2480 . f - The function to compute the matrices 2481 - ctx - User-defined context for Jacobian evaluation routine 2482 2483 Notes: You can pass in PETSC_NULL for any return argument you do not need. 2484 2485 Level: advanced 2486 2487 .seealso: TSGetTimeStep(), TSGetRHSJacobian(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber() 2488 2489 .keywords: TS, timestep, get, matrix, Jacobian 2490 @*/ 2491 PetscErrorCode TSGetIJacobian(TS ts,Mat *A,Mat *B,TSIJacobian *f,void **ctx) 2492 { 2493 PetscErrorCode ierr; 2494 SNES snes; 2495 DM dm; 2496 PetscFunctionBegin; 2497 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2498 ierr = SNESGetJacobian(snes,A,B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 2499 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 2500 ierr = DMTSGetIJacobian(dm,f,ctx);CHKERRQ(ierr); 2501 PetscFunctionReturn(0); 2502 } 2503 2504 typedef struct { 2505 PetscViewer viewer; 2506 Vec initialsolution; 2507 PetscBool showinitial; 2508 } TSMonitorSolutionCtx; 2509 2510 #undef __FUNCT__ 2511 #define __FUNCT__ "TSMonitorSolution" 2512 /*@C 2513 TSMonitorSolution - Monitors progress of the TS solvers by calling 2514 VecView() for the solution at each timestep 2515 2516 Collective on TS 2517 2518 Input Parameters: 2519 + ts - the TS context 2520 . step - current time-step 2521 . ptime - current time 2522 - dummy - either a viewer or PETSC_NULL 2523 2524 Level: intermediate 2525 2526 .keywords: TS, vector, monitor, view 2527 2528 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView() 2529 @*/ 2530 PetscErrorCode TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy) 2531 { 2532 PetscErrorCode ierr; 2533 TSMonitorSolutionCtx *ictx = (TSMonitorSolutionCtx*)dummy; 2534 2535 PetscFunctionBegin; 2536 if (!step && ictx->showinitial) { 2537 if (!ictx->initialsolution) { 2538 ierr = VecDuplicate(x,&ictx->initialsolution);CHKERRQ(ierr); 2539 } 2540 ierr = VecCopy(x,ictx->initialsolution);CHKERRQ(ierr); 2541 } 2542 if (ictx->showinitial) { 2543 PetscReal pause; 2544 ierr = PetscViewerDrawGetPause(ictx->viewer,&pause);CHKERRQ(ierr); 2545 ierr = PetscViewerDrawSetPause(ictx->viewer,0.0);CHKERRQ(ierr); 2546 ierr = VecView(ictx->initialsolution,ictx->viewer);CHKERRQ(ierr); 2547 ierr = PetscViewerDrawSetPause(ictx->viewer,pause);CHKERRQ(ierr); 2548 ierr = PetscViewerDrawSetHold(ictx->viewer,PETSC_TRUE);CHKERRQ(ierr); 2549 } 2550 ierr = VecView(x,ictx->viewer);CHKERRQ(ierr); 2551 if (ictx->showinitial) { 2552 ierr = PetscViewerDrawSetHold(ictx->viewer,PETSC_FALSE);CHKERRQ(ierr); 2553 } 2554 PetscFunctionReturn(0); 2555 } 2556 2557 2558 #undef __FUNCT__ 2559 #define __FUNCT__ "TSMonitorSolutionDestroy" 2560 /*@C 2561 TSMonitorSolutionDestroy - Destroys the monitor context for TSMonitorSolution 2562 2563 Collective on TS 2564 2565 Input Parameters: 2566 . ctx - the monitor context 2567 2568 Level: intermediate 2569 2570 .keywords: TS, vector, monitor, view 2571 2572 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorSolution() 2573 @*/ 2574 PetscErrorCode TSMonitorSolutionDestroy(void **ctx) 2575 { 2576 PetscErrorCode ierr; 2577 TSMonitorSolutionCtx *ictx = *(TSMonitorSolutionCtx**)ctx; 2578 2579 PetscFunctionBegin; 2580 ierr = PetscViewerDestroy(&ictx->viewer);CHKERRQ(ierr); 2581 ierr = VecDestroy(&ictx->initialsolution);CHKERRQ(ierr); 2582 ierr = PetscFree(ictx);CHKERRQ(ierr); 2583 PetscFunctionReturn(0); 2584 } 2585 2586 #undef __FUNCT__ 2587 #define __FUNCT__ "TSMonitorSolutionCreate" 2588 /*@C 2589 TSMonitorSolutionCreate - Creates the monitor context for TSMonitorSolution 2590 2591 Collective on TS 2592 2593 Input Parameter: 2594 . ts - time-step context 2595 2596 Output Patameter: 2597 . ctx - the monitor context 2598 2599 Level: intermediate 2600 2601 .keywords: TS, vector, monitor, view 2602 2603 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorSolution() 2604 @*/ 2605 PetscErrorCode TSMonitorSolutionCreate(TS ts,PetscViewer viewer,void **ctx) 2606 { 2607 PetscErrorCode ierr; 2608 TSMonitorSolutionCtx *ictx; 2609 2610 PetscFunctionBegin; 2611 ierr = PetscNew(TSMonitorSolutionCtx,&ictx);CHKERRQ(ierr); 2612 *ctx = (void*)ictx; 2613 if (!viewer) { 2614 viewer = PETSC_VIEWER_DRAW_(((PetscObject)ts)->comm); 2615 } 2616 ierr = PetscObjectReference((PetscObject)viewer);CHKERRQ(ierr); 2617 ictx->viewer = viewer; 2618 ictx->showinitial = PETSC_FALSE; 2619 ierr = PetscOptionsGetBool(((PetscObject)ts)->prefix,"-ts_monitor_solution_initial",&ictx->showinitial,PETSC_NULL);CHKERRQ(ierr); 2620 PetscFunctionReturn(0); 2621 } 2622 2623 #undef __FUNCT__ 2624 #define __FUNCT__ "TSSetDM" 2625 /*@ 2626 TSSetDM - Sets the DM that may be used by some preconditioners 2627 2628 Logically Collective on TS and DM 2629 2630 Input Parameters: 2631 + ts - the preconditioner context 2632 - dm - the dm 2633 2634 Level: intermediate 2635 2636 2637 .seealso: TSGetDM(), SNESSetDM(), SNESGetDM() 2638 @*/ 2639 PetscErrorCode TSSetDM(TS ts,DM dm) 2640 { 2641 PetscErrorCode ierr; 2642 SNES snes; 2643 TSDM tsdm; 2644 2645 PetscFunctionBegin; 2646 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2647 ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr); 2648 if (ts->dm) { /* Move the TSDM context over to the new DM unless the new DM already has one */ 2649 PetscContainer oldcontainer,container; 2650 ierr = PetscObjectQuery((PetscObject)ts->dm,"TSDM",(PetscObject*)&oldcontainer);CHKERRQ(ierr); 2651 ierr = PetscObjectQuery((PetscObject)dm,"TSDM",(PetscObject*)&container);CHKERRQ(ierr); 2652 if (oldcontainer && !container) { 2653 ierr = DMTSCopyContext(ts->dm,dm);CHKERRQ(ierr); 2654 ierr = DMTSGetContext(ts->dm,&tsdm);CHKERRQ(ierr); 2655 if (tsdm->originaldm == ts->dm) { /* Grant write privileges to the replacement DM */ 2656 tsdm->originaldm = dm; 2657 } 2658 } 2659 ierr = DMDestroy(&ts->dm);CHKERRQ(ierr); 2660 } 2661 ts->dm = dm; 2662 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2663 ierr = SNESSetDM(snes,dm);CHKERRQ(ierr); 2664 PetscFunctionReturn(0); 2665 } 2666 2667 #undef __FUNCT__ 2668 #define __FUNCT__ "TSGetDM" 2669 /*@ 2670 TSGetDM - Gets the DM that may be used by some preconditioners 2671 2672 Not Collective 2673 2674 Input Parameter: 2675 . ts - the preconditioner context 2676 2677 Output Parameter: 2678 . dm - the dm 2679 2680 Level: intermediate 2681 2682 2683 .seealso: TSSetDM(), SNESSetDM(), SNESGetDM() 2684 @*/ 2685 PetscErrorCode TSGetDM(TS ts,DM *dm) 2686 { 2687 PetscErrorCode ierr; 2688 2689 PetscFunctionBegin; 2690 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2691 if (!ts->dm) { 2692 ierr = DMShellCreate(((PetscObject)ts)->comm,&ts->dm);CHKERRQ(ierr); 2693 if (ts->snes) {ierr = SNESSetDM(ts->snes,ts->dm);CHKERRQ(ierr);} 2694 } 2695 *dm = ts->dm; 2696 PetscFunctionReturn(0); 2697 } 2698 2699 #undef __FUNCT__ 2700 #define __FUNCT__ "SNESTSFormFunction" 2701 /*@ 2702 SNESTSFormFunction - Function to evaluate nonlinear residual 2703 2704 Logically Collective on SNES 2705 2706 Input Parameter: 2707 + snes - nonlinear solver 2708 . X - the current state at which to evaluate the residual 2709 - ctx - user context, must be a TS 2710 2711 Output Parameter: 2712 . F - the nonlinear residual 2713 2714 Notes: 2715 This function is not normally called by users and is automatically registered with the SNES used by TS. 2716 It is most frequently passed to MatFDColoringSetFunction(). 2717 2718 Level: advanced 2719 2720 .seealso: SNESSetFunction(), MatFDColoringSetFunction() 2721 @*/ 2722 PetscErrorCode SNESTSFormFunction(SNES snes,Vec X,Vec F,void *ctx) 2723 { 2724 TS ts = (TS)ctx; 2725 PetscErrorCode ierr; 2726 2727 PetscFunctionBegin; 2728 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2729 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 2730 PetscValidHeaderSpecific(F,VEC_CLASSID,3); 2731 PetscValidHeaderSpecific(ts,TS_CLASSID,4); 2732 ierr = (ts->ops->snesfunction)(snes,X,F,ts);CHKERRQ(ierr); 2733 PetscFunctionReturn(0); 2734 } 2735 2736 #undef __FUNCT__ 2737 #define __FUNCT__ "SNESTSFormJacobian" 2738 /*@ 2739 SNESTSFormJacobian - Function to evaluate the Jacobian 2740 2741 Collective on SNES 2742 2743 Input Parameter: 2744 + snes - nonlinear solver 2745 . X - the current state at which to evaluate the residual 2746 - ctx - user context, must be a TS 2747 2748 Output Parameter: 2749 + A - the Jacobian 2750 . B - the preconditioning matrix (may be the same as A) 2751 - flag - indicates any structure change in the matrix 2752 2753 Notes: 2754 This function is not normally called by users and is automatically registered with the SNES used by TS. 2755 2756 Level: developer 2757 2758 .seealso: SNESSetJacobian() 2759 @*/ 2760 PetscErrorCode SNESTSFormJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flag,void *ctx) 2761 { 2762 TS ts = (TS)ctx; 2763 PetscErrorCode ierr; 2764 2765 PetscFunctionBegin; 2766 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2767 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 2768 PetscValidPointer(A,3); 2769 PetscValidHeaderSpecific(*A,MAT_CLASSID,3); 2770 PetscValidPointer(B,4); 2771 PetscValidHeaderSpecific(*B,MAT_CLASSID,4); 2772 PetscValidPointer(flag,5); 2773 PetscValidHeaderSpecific(ts,TS_CLASSID,6); 2774 ierr = (ts->ops->snesjacobian)(snes,X,A,B,flag,ts);CHKERRQ(ierr); 2775 PetscFunctionReturn(0); 2776 } 2777 2778 #undef __FUNCT__ 2779 #define __FUNCT__ "TSComputeRHSFunctionLinear" 2780 /*@C 2781 TSComputeRHSFunctionLinear - Evaluate the right hand side via the user-provided Jacobian, for linear problems only 2782 2783 Collective on TS 2784 2785 Input Arguments: 2786 + ts - time stepping context 2787 . t - time at which to evaluate 2788 . X - state at which to evaluate 2789 - ctx - context 2790 2791 Output Arguments: 2792 . F - right hand side 2793 2794 Level: intermediate 2795 2796 Notes: 2797 This function is intended to be passed to TSSetRHSFunction() to evaluate the right hand side for linear problems. 2798 The matrix (and optionally the evaluation context) should be passed to TSSetRHSJacobian(). 2799 2800 .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSJacobianConstant() 2801 @*/ 2802 PetscErrorCode TSComputeRHSFunctionLinear(TS ts,PetscReal t,Vec X,Vec F,void *ctx) 2803 { 2804 PetscErrorCode ierr; 2805 Mat Arhs,Brhs; 2806 MatStructure flg2; 2807 2808 PetscFunctionBegin; 2809 ierr = TSGetRHSMats_Private(ts,&Arhs,&Brhs);CHKERRQ(ierr); 2810 ierr = TSComputeRHSJacobian(ts,t,X,&Arhs,&Brhs,&flg2);CHKERRQ(ierr); 2811 ierr = MatMult(Arhs,X,F);CHKERRQ(ierr); 2812 PetscFunctionReturn(0); 2813 } 2814 2815 #undef __FUNCT__ 2816 #define __FUNCT__ "TSComputeRHSJacobianConstant" 2817 /*@C 2818 TSComputeRHSJacobianConstant - Reuses a Jacobian that is time-independent. 2819 2820 Collective on TS 2821 2822 Input Arguments: 2823 + ts - time stepping context 2824 . t - time at which to evaluate 2825 . X - state at which to evaluate 2826 - ctx - context 2827 2828 Output Arguments: 2829 + A - pointer to operator 2830 . B - pointer to preconditioning matrix 2831 - flg - matrix structure flag 2832 2833 Level: intermediate 2834 2835 Notes: 2836 This function is intended to be passed to TSSetRHSJacobian() to evaluate the Jacobian for linear time-independent problems. 2837 2838 .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSFunctionLinear() 2839 @*/ 2840 PetscErrorCode TSComputeRHSJacobianConstant(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg,void *ctx) 2841 { 2842 2843 PetscFunctionBegin; 2844 *flg = SAME_PRECONDITIONER; 2845 PetscFunctionReturn(0); 2846 } 2847 2848 #undef __FUNCT__ 2849 #define __FUNCT__ "TSComputeIFunctionLinear" 2850 /*@C 2851 TSComputeIFunctionLinear - Evaluate the left hand side via the user-provided Jacobian, for linear problems only 2852 2853 Collective on TS 2854 2855 Input Arguments: 2856 + ts - time stepping context 2857 . t - time at which to evaluate 2858 . X - state at which to evaluate 2859 . Xdot - time derivative of state vector 2860 - ctx - context 2861 2862 Output Arguments: 2863 . F - left hand side 2864 2865 Level: intermediate 2866 2867 Notes: 2868 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 2869 user is required to write their own TSComputeIFunction. 2870 This function is intended to be passed to TSSetIFunction() to evaluate the left hand side for linear problems. 2871 The matrix (and optionally the evaluation context) should be passed to TSSetIJacobian(). 2872 2873 .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIJacobianConstant() 2874 @*/ 2875 PetscErrorCode TSComputeIFunctionLinear(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx) 2876 { 2877 PetscErrorCode ierr; 2878 Mat A,B; 2879 MatStructure flg2; 2880 2881 PetscFunctionBegin; 2882 ierr = TSGetIJacobian(ts,&A,&B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 2883 ierr = TSComputeIJacobian(ts,t,X,Xdot,1.0,&A,&B,&flg2,PETSC_TRUE);CHKERRQ(ierr); 2884 ierr = MatMult(A,Xdot,F);CHKERRQ(ierr); 2885 PetscFunctionReturn(0); 2886 } 2887 2888 #undef __FUNCT__ 2889 #define __FUNCT__ "TSComputeIJacobianConstant" 2890 /*@C 2891 TSComputeIJacobianConstant - Reuses a Jacobian that is time-independent. 2892 2893 Collective on TS 2894 2895 Input Arguments: 2896 + ts - time stepping context 2897 . t - time at which to evaluate 2898 . X - state at which to evaluate 2899 . Xdot - time derivative of state vector 2900 . shift - shift to apply 2901 - ctx - context 2902 2903 Output Arguments: 2904 + A - pointer to operator 2905 . B - pointer to preconditioning matrix 2906 - flg - matrix structure flag 2907 2908 Level: intermediate 2909 2910 Notes: 2911 This function is intended to be passed to TSSetIJacobian() to evaluate the Jacobian for linear time-independent problems. 2912 2913 .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIFunctionLinear() 2914 @*/ 2915 PetscErrorCode TSComputeIJacobianConstant(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg,void *ctx) 2916 { 2917 2918 PetscFunctionBegin; 2919 *flg = SAME_PRECONDITIONER; 2920 PetscFunctionReturn(0); 2921 } 2922 2923 2924 #undef __FUNCT__ 2925 #define __FUNCT__ "TSGetConvergedReason" 2926 /*@ 2927 TSGetConvergedReason - Gets the reason the TS iteration was stopped. 2928 2929 Not Collective 2930 2931 Input Parameter: 2932 . ts - the TS context 2933 2934 Output Parameter: 2935 . reason - negative value indicates diverged, positive value converged, see TSConvergedReason or the 2936 manual pages for the individual convergence tests for complete lists 2937 2938 Level: intermediate 2939 2940 Notes: 2941 Can only be called after the call to TSSolve() is complete. 2942 2943 .keywords: TS, nonlinear, set, convergence, test 2944 2945 .seealso: TSSetConvergenceTest(), TSConvergedReason 2946 @*/ 2947 PetscErrorCode TSGetConvergedReason(TS ts,TSConvergedReason *reason) 2948 { 2949 PetscFunctionBegin; 2950 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2951 PetscValidPointer(reason,2); 2952 *reason = ts->reason; 2953 PetscFunctionReturn(0); 2954 } 2955 2956 #undef __FUNCT__ 2957 #define __FUNCT__ "TSGetSNESIterations" 2958 /*@ 2959 TSGetSNESIterations - Gets the total number of nonlinear iterations 2960 used by the time integrator. 2961 2962 Not Collective 2963 2964 Input Parameter: 2965 . ts - TS context 2966 2967 Output Parameter: 2968 . nits - number of nonlinear iterations 2969 2970 Notes: 2971 This counter is reset to zero for each successive call to TSSolve(). 2972 2973 Level: intermediate 2974 2975 .keywords: TS, get, number, nonlinear, iterations 2976 2977 .seealso: TSGetKSPIterations() 2978 @*/ 2979 PetscErrorCode TSGetSNESIterations(TS ts,PetscInt *nits) 2980 { 2981 PetscFunctionBegin; 2982 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2983 PetscValidIntPointer(nits,2); 2984 *nits = ts->snes_its; 2985 PetscFunctionReturn(0); 2986 } 2987 2988 #undef __FUNCT__ 2989 #define __FUNCT__ "TSGetKSPIterations" 2990 /*@ 2991 TSGetKSPIterations - Gets the total number of linear iterations 2992 used by the time integrator. 2993 2994 Not Collective 2995 2996 Input Parameter: 2997 . ts - TS context 2998 2999 Output Parameter: 3000 . lits - number of linear iterations 3001 3002 Notes: 3003 This counter is reset to zero for each successive call to TSSolve(). 3004 3005 Level: intermediate 3006 3007 .keywords: TS, get, number, linear, iterations 3008 3009 .seealso: TSGetSNESIterations(), SNESGetKSPIterations() 3010 @*/ 3011 PetscErrorCode TSGetKSPIterations(TS ts,PetscInt *lits) 3012 { 3013 PetscFunctionBegin; 3014 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3015 PetscValidIntPointer(lits,2); 3016 *lits = ts->ksp_its; 3017 PetscFunctionReturn(0); 3018 } 3019 3020 #undef __FUNCT__ 3021 #define __FUNCT__ "TSGetStepRejections" 3022 /*@ 3023 TSGetStepRejections - Gets the total number of rejected steps. 3024 3025 Not Collective 3026 3027 Input Parameter: 3028 . ts - TS context 3029 3030 Output Parameter: 3031 . rejects - number of steps rejected 3032 3033 Notes: 3034 This counter is reset to zero for each successive call to TSSolve(). 3035 3036 Level: intermediate 3037 3038 .keywords: TS, get, number 3039 3040 .seealso: TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetSNESFailures(), TSSetMaxSNESFailures(), TSSetErrorIfStepFails() 3041 @*/ 3042 PetscErrorCode TSGetStepRejections(TS ts,PetscInt *rejects) 3043 { 3044 PetscFunctionBegin; 3045 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3046 PetscValidIntPointer(rejects,2); 3047 *rejects = ts->reject; 3048 PetscFunctionReturn(0); 3049 } 3050 3051 #undef __FUNCT__ 3052 #define __FUNCT__ "TSGetSNESFailures" 3053 /*@ 3054 TSGetSNESFailures - Gets the total number of failed SNES solves 3055 3056 Not Collective 3057 3058 Input Parameter: 3059 . ts - TS context 3060 3061 Output Parameter: 3062 . fails - number of failed nonlinear solves 3063 3064 Notes: 3065 This counter is reset to zero for each successive call to TSSolve(). 3066 3067 Level: intermediate 3068 3069 .keywords: TS, get, number 3070 3071 .seealso: TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSSetMaxSNESFailures() 3072 @*/ 3073 PetscErrorCode TSGetSNESFailures(TS ts,PetscInt *fails) 3074 { 3075 PetscFunctionBegin; 3076 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3077 PetscValidIntPointer(fails,2); 3078 *fails = ts->num_snes_failures; 3079 PetscFunctionReturn(0); 3080 } 3081 3082 #undef __FUNCT__ 3083 #define __FUNCT__ "TSSetMaxStepRejections" 3084 /*@ 3085 TSSetMaxStepRejections - Sets the maximum number of step rejections before a step fails 3086 3087 Not Collective 3088 3089 Input Parameter: 3090 + ts - TS context 3091 - rejects - maximum number of rejected steps, pass -1 for unlimited 3092 3093 Notes: 3094 The counter is reset to zero for each step 3095 3096 Options Database Key: 3097 . -ts_max_reject - Maximum number of step rejections before a step fails 3098 3099 Level: intermediate 3100 3101 .keywords: TS, set, maximum, number 3102 3103 .seealso: TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxSNESFailures(), TSGetStepRejections(), TSGetSNESFailures(), TSSetErrorIfStepFails(), TSGetConvergedReason() 3104 @*/ 3105 PetscErrorCode TSSetMaxStepRejections(TS ts,PetscInt rejects) 3106 { 3107 PetscFunctionBegin; 3108 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3109 ts->max_reject = rejects; 3110 PetscFunctionReturn(0); 3111 } 3112 3113 #undef __FUNCT__ 3114 #define __FUNCT__ "TSSetMaxSNESFailures" 3115 /*@ 3116 TSSetMaxSNESFailures - Sets the maximum number of failed SNES solves 3117 3118 Not Collective 3119 3120 Input Parameter: 3121 + ts - TS context 3122 - fails - maximum number of failed nonlinear solves, pass -1 for unlimited 3123 3124 Notes: 3125 The counter is reset to zero for each successive call to TSSolve(). 3126 3127 Options Database Key: 3128 . -ts_max_snes_failures - Maximum number of nonlinear solve failures 3129 3130 Level: intermediate 3131 3132 .keywords: TS, set, maximum, number 3133 3134 .seealso: TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSGetSNESFailures(), SNESGetConvergedReason(), TSGetConvergedReason() 3135 @*/ 3136 PetscErrorCode TSSetMaxSNESFailures(TS ts,PetscInt fails) 3137 { 3138 PetscFunctionBegin; 3139 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3140 ts->max_snes_failures = fails; 3141 PetscFunctionReturn(0); 3142 } 3143 3144 #undef __FUNCT__ 3145 #define __FUNCT__ "TSSetErrorIfStepFails()" 3146 /*@ 3147 TSSetErrorIfStepFails - Error if no step succeeds 3148 3149 Not Collective 3150 3151 Input Parameter: 3152 + ts - TS context 3153 - err - PETSC_TRUE to error if no step succeeds, PETSC_FALSE to return without failure 3154 3155 Options Database Key: 3156 . -ts_error_if_step_fails - Error if no step succeeds 3157 3158 Level: intermediate 3159 3160 .keywords: TS, set, error 3161 3162 .seealso: TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSGetSNESFailures(), TSSetErrorIfStepFails(), TSGetConvergedReason() 3163 @*/ 3164 PetscErrorCode TSSetErrorIfStepFails(TS ts,PetscBool err) 3165 { 3166 PetscFunctionBegin; 3167 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3168 ts->errorifstepfailed = err; 3169 PetscFunctionReturn(0); 3170 } 3171 3172 #undef __FUNCT__ 3173 #define __FUNCT__ "TSMonitorSolutionBinary" 3174 /*@C 3175 TSMonitorSolutionBinary - Monitors progress of the TS solvers by VecView() for the solution at each timestep. Normally the viewer is a binary file 3176 3177 Collective on TS 3178 3179 Input Parameters: 3180 + ts - the TS context 3181 . step - current time-step 3182 . ptime - current time 3183 . x - current state 3184 - viewer - binary viewer 3185 3186 Level: intermediate 3187 3188 .keywords: TS, vector, monitor, view 3189 3190 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView() 3191 @*/ 3192 PetscErrorCode TSMonitorSolutionBinary(TS ts,PetscInt step,PetscReal ptime,Vec x,void *viewer) 3193 { 3194 PetscErrorCode ierr; 3195 PetscViewer v = (PetscViewer)viewer; 3196 3197 PetscFunctionBegin; 3198 ierr = VecView(x,v);CHKERRQ(ierr); 3199 PetscFunctionReturn(0); 3200 } 3201 3202 #undef __FUNCT__ 3203 #define __FUNCT__ "TSMonitorSolutionVTK" 3204 /*@C 3205 TSMonitorSolutionVTK - Monitors progress of the TS solvers by VecView() for the solution at each timestep. 3206 3207 Collective on TS 3208 3209 Input Parameters: 3210 + ts - the TS context 3211 . step - current time-step 3212 . ptime - current time 3213 . x - current state 3214 - filenametemplate - string containing a format specifier for the integer time step (e.g. %03D) 3215 3216 Level: intermediate 3217 3218 Notes: 3219 The VTK format does not allow writing multiple time steps in the same file, therefore a different file will be written for each time step. 3220 These are named according to the file name template. 3221 3222 This function is normally passed as an argument to TSMonitorSet() along with TSMonitorSolutionVTKDestroy(). 3223 3224 .keywords: TS, vector, monitor, view 3225 3226 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView() 3227 @*/ 3228 PetscErrorCode TSMonitorSolutionVTK(TS ts,PetscInt step,PetscReal ptime,Vec x,void *filenametemplate) 3229 { 3230 PetscErrorCode ierr; 3231 char filename[PETSC_MAX_PATH_LEN]; 3232 PetscViewer viewer; 3233 3234 PetscFunctionBegin; 3235 ierr = PetscSNPrintf(filename,sizeof filename,(const char*)filenametemplate,step);CHKERRQ(ierr); 3236 ierr = PetscViewerVTKOpen(((PetscObject)ts)->comm,filename,FILE_MODE_WRITE,&viewer);CHKERRQ(ierr); 3237 ierr = VecView(x,viewer);CHKERRQ(ierr); 3238 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 3239 PetscFunctionReturn(0); 3240 } 3241 3242 #undef __FUNCT__ 3243 #define __FUNCT__ "TSMonitorSolutionVTKDestroy" 3244 /*@C 3245 TSMonitorSolutionVTKDestroy - Destroy context for monitoring 3246 3247 Collective on TS 3248 3249 Input Parameters: 3250 . filenametemplate - string containing a format specifier for the integer time step (e.g. %03D) 3251 3252 Level: intermediate 3253 3254 Note: 3255 This function is normally passed to TSMonitorSet() along with TSMonitorSolutionVTK(). 3256 3257 .keywords: TS, vector, monitor, view 3258 3259 .seealso: TSMonitorSet(), TSMonitorSolutionVTK() 3260 @*/ 3261 PetscErrorCode TSMonitorSolutionVTKDestroy(void *filenametemplate) 3262 { 3263 PetscErrorCode ierr; 3264 3265 PetscFunctionBegin; 3266 ierr = PetscFree(*(char**)filenametemplate);CHKERRQ(ierr); 3267 PetscFunctionReturn(0); 3268 } 3269 3270 #undef __FUNCT__ 3271 #define __FUNCT__ "TSGetAdapt" 3272 /*@ 3273 TSGetAdapt - Get the adaptive controller context for the current method 3274 3275 Collective on TS if controller has not been created yet 3276 3277 Input Arguments: 3278 . ts - time stepping context 3279 3280 Output Arguments: 3281 . adapt - adaptive controller 3282 3283 Level: intermediate 3284 3285 .seealso: TSAdapt, TSAdaptSetType(), TSAdaptChoose() 3286 @*/ 3287 PetscErrorCode TSGetAdapt(TS ts,TSAdapt *adapt) 3288 { 3289 PetscErrorCode ierr; 3290 3291 PetscFunctionBegin; 3292 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3293 PetscValidPointer(adapt,2); 3294 if (!ts->adapt) { 3295 ierr = TSAdaptCreate(((PetscObject)ts)->comm,&ts->adapt);CHKERRQ(ierr); 3296 ierr = PetscLogObjectParent(ts,ts->adapt);CHKERRQ(ierr); 3297 ierr = PetscObjectIncrementTabLevel((PetscObject)ts->adapt,(PetscObject)ts,1);CHKERRQ(ierr); 3298 } 3299 *adapt = ts->adapt; 3300 PetscFunctionReturn(0); 3301 } 3302 3303 #undef __FUNCT__ 3304 #define __FUNCT__ "TSSetTolerances" 3305 /*@ 3306 TSSetTolerances - Set tolerances for local truncation error when using adaptive controller 3307 3308 Logically Collective 3309 3310 Input Arguments: 3311 + ts - time integration context 3312 . atol - scalar absolute tolerances, PETSC_DECIDE to leave current value 3313 . vatol - vector of absolute tolerances or PETSC_NULL, used in preference to atol if present 3314 . rtol - scalar relative tolerances, PETSC_DECIDE to leave current value 3315 - vrtol - vector of relative tolerances or PETSC_NULL, used in preference to atol if present 3316 3317 Level: beginner 3318 3319 .seealso: TS, TSAdapt, TSVecNormWRMS(), TSGetTolerances() 3320 @*/ 3321 PetscErrorCode TSSetTolerances(TS ts,PetscReal atol,Vec vatol,PetscReal rtol,Vec vrtol) 3322 { 3323 PetscErrorCode ierr; 3324 3325 PetscFunctionBegin; 3326 if (atol != PETSC_DECIDE && atol != PETSC_DEFAULT) ts->atol = atol; 3327 if (vatol) { 3328 ierr = PetscObjectReference((PetscObject)vatol);CHKERRQ(ierr); 3329 ierr = VecDestroy(&ts->vatol);CHKERRQ(ierr); 3330 ts->vatol = vatol; 3331 } 3332 if (rtol != PETSC_DECIDE && rtol != PETSC_DEFAULT) ts->rtol = rtol; 3333 if (vrtol) { 3334 ierr = PetscObjectReference((PetscObject)vrtol);CHKERRQ(ierr); 3335 ierr = VecDestroy(&ts->vrtol);CHKERRQ(ierr); 3336 ts->vrtol = vrtol; 3337 } 3338 PetscFunctionReturn(0); 3339 } 3340 3341 #undef __FUNCT__ 3342 #define __FUNCT__ "TSGetTolerances" 3343 /*@ 3344 TSGetTolerances - Get tolerances for local truncation error when using adaptive controller 3345 3346 Logically Collective 3347 3348 Input Arguments: 3349 . ts - time integration context 3350 3351 Output Arguments: 3352 + atol - scalar absolute tolerances, PETSC_NULL to ignore 3353 . vatol - vector of absolute tolerances, PETSC_NULL to ignore 3354 . rtol - scalar relative tolerances, PETSC_NULL to ignore 3355 - vrtol - vector of relative tolerances, PETSC_NULL to ignore 3356 3357 Level: beginner 3358 3359 .seealso: TS, TSAdapt, TSVecNormWRMS(), TSSetTolerances() 3360 @*/ 3361 PetscErrorCode TSGetTolerances(TS ts,PetscReal *atol,Vec *vatol,PetscReal *rtol,Vec *vrtol) 3362 { 3363 3364 PetscFunctionBegin; 3365 if (atol) *atol = ts->atol; 3366 if (vatol) *vatol = ts->vatol; 3367 if (rtol) *rtol = ts->rtol; 3368 if (vrtol) *vrtol = ts->vrtol; 3369 PetscFunctionReturn(0); 3370 } 3371 3372 #undef __FUNCT__ 3373 #define __FUNCT__ "TSErrorNormWRMS" 3374 /*@ 3375 TSErrorNormWRMS - compute a weighted norm of the difference between a vector and the current state 3376 3377 Collective on TS 3378 3379 Input Arguments: 3380 + ts - time stepping context 3381 - Y - state vector to be compared to ts->vec_sol 3382 3383 Output Arguments: 3384 . norm - weighted norm, a value of 1.0 is considered small 3385 3386 Level: developer 3387 3388 .seealso: TSSetTolerances() 3389 @*/ 3390 PetscErrorCode TSErrorNormWRMS(TS ts,Vec Y,PetscReal *norm) 3391 { 3392 PetscErrorCode ierr; 3393 PetscInt i,n,N; 3394 const PetscScalar *x,*y; 3395 Vec X; 3396 PetscReal sum,gsum; 3397 3398 PetscFunctionBegin; 3399 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3400 PetscValidHeaderSpecific(Y,VEC_CLASSID,2); 3401 PetscValidPointer(norm,3); 3402 X = ts->vec_sol; 3403 PetscCheckSameTypeAndComm(X,1,Y,2); 3404 if (X == Y) SETERRQ(((PetscObject)X)->comm,PETSC_ERR_ARG_IDN,"Y cannot be the TS solution vector"); 3405 3406 ierr = VecGetSize(X,&N);CHKERRQ(ierr); 3407 ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr); 3408 ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 3409 ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 3410 sum = 0.; 3411 if (ts->vatol && ts->vrtol) { 3412 const PetscScalar *atol,*rtol; 3413 ierr = VecGetArrayRead(ts->vatol,&atol);CHKERRQ(ierr); 3414 ierr = VecGetArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr); 3415 for (i=0; i<n; i++) { 3416 PetscReal tol = PetscRealPart(atol[i]) + PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i])); 3417 sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol); 3418 } 3419 ierr = VecRestoreArrayRead(ts->vatol,&atol);CHKERRQ(ierr); 3420 ierr = VecRestoreArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr); 3421 } else if (ts->vatol) { /* vector atol, scalar rtol */ 3422 const PetscScalar *atol; 3423 ierr = VecGetArrayRead(ts->vatol,&atol);CHKERRQ(ierr); 3424 for (i=0; i<n; i++) { 3425 PetscReal tol = PetscRealPart(atol[i]) + ts->rtol * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i])); 3426 sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol); 3427 } 3428 ierr = VecRestoreArrayRead(ts->vatol,&atol);CHKERRQ(ierr); 3429 } else if (ts->vrtol) { /* scalar atol, vector rtol */ 3430 const PetscScalar *rtol; 3431 ierr = VecGetArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr); 3432 for (i=0; i<n; i++) { 3433 PetscReal tol = ts->atol + PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i])); 3434 sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol); 3435 } 3436 ierr = VecRestoreArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr); 3437 } else { /* scalar atol, scalar rtol */ 3438 for (i=0; i<n; i++) { 3439 PetscReal tol = ts->atol + ts->rtol * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i])); 3440 sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol); 3441 } 3442 } 3443 ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 3444 ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 3445 3446 ierr = MPI_Allreduce(&sum,&gsum,1,MPIU_REAL,MPIU_SUM,((PetscObject)ts)->comm);CHKERRQ(ierr); 3447 *norm = PetscSqrtReal(gsum / N); 3448 if (PetscIsInfOrNanScalar(*norm)) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 3449 PetscFunctionReturn(0); 3450 } 3451 3452 #undef __FUNCT__ 3453 #define __FUNCT__ "TSSetCFLTimeLocal" 3454 /*@ 3455 TSSetCFLTimeLocal - Set the local CFL constraint relative to forward Euler 3456 3457 Logically Collective on TS 3458 3459 Input Arguments: 3460 + ts - time stepping context 3461 - cfltime - maximum stable time step if using forward Euler (value can be different on each process) 3462 3463 Note: 3464 After calling this function, the global CFL time can be obtained by calling TSGetCFLTime() 3465 3466 Level: intermediate 3467 3468 .seealso: TSGetCFLTime(), TSADAPTCFL 3469 @*/ 3470 PetscErrorCode TSSetCFLTimeLocal(TS ts,PetscReal cfltime) 3471 { 3472 3473 PetscFunctionBegin; 3474 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3475 ts->cfltime_local = cfltime; 3476 ts->cfltime = -1.; 3477 PetscFunctionReturn(0); 3478 } 3479 3480 #undef __FUNCT__ 3481 #define __FUNCT__ "TSGetCFLTime" 3482 /*@ 3483 TSGetCFLTime - Get the maximum stable time step according to CFL criteria applied to forward Euler 3484 3485 Collective on TS 3486 3487 Input Arguments: 3488 . ts - time stepping context 3489 3490 Output Arguments: 3491 . cfltime - maximum stable time step for forward Euler 3492 3493 Level: advanced 3494 3495 .seealso: TSSetCFLTimeLocal() 3496 @*/ 3497 PetscErrorCode TSGetCFLTime(TS ts,PetscReal *cfltime) 3498 { 3499 PetscErrorCode ierr; 3500 3501 PetscFunctionBegin; 3502 if (ts->cfltime < 0) { 3503 ierr = MPI_Allreduce(&ts->cfltime_local,&ts->cfltime,1,MPIU_REAL,MPIU_MIN,((PetscObject)ts)->comm);CHKERRQ(ierr); 3504 } 3505 *cfltime = ts->cfltime; 3506 PetscFunctionReturn(0); 3507 } 3508 3509 #undef __FUNCT__ 3510 #define __FUNCT__ "TSVISetVariableBounds" 3511 /*@ 3512 TSVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu 3513 3514 Input Parameters: 3515 . ts - the TS context. 3516 . xl - lower bound. 3517 . xu - upper bound. 3518 3519 Notes: 3520 If this routine is not called then the lower and upper bounds are set to 3521 SNES_VI_NINF and SNES_VI_INF respectively during SNESSetUp(). 3522 3523 Level: advanced 3524 3525 @*/ 3526 PetscErrorCode TSVISetVariableBounds(TS ts, Vec xl, Vec xu) 3527 { 3528 PetscErrorCode ierr; 3529 SNES snes; 3530 3531 PetscFunctionBegin; 3532 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 3533 ierr = SNESVISetVariableBounds(snes,xl,xu);CHKERRQ(ierr); 3534 PetscFunctionReturn(0); 3535 } 3536 3537 #if defined(PETSC_HAVE_MATLAB_ENGINE) 3538 #include <mex.h> 3539 3540 typedef struct {char *funcname; mxArray *ctx;} TSMatlabContext; 3541 3542 #undef __FUNCT__ 3543 #define __FUNCT__ "TSComputeFunction_Matlab" 3544 /* 3545 TSComputeFunction_Matlab - Calls the function that has been set with 3546 TSSetFunctionMatlab(). 3547 3548 Collective on TS 3549 3550 Input Parameters: 3551 + snes - the TS context 3552 - x - input vector 3553 3554 Output Parameter: 3555 . y - function vector, as set by TSSetFunction() 3556 3557 Notes: 3558 TSComputeFunction() is typically used within nonlinear solvers 3559 implementations, so most users would not generally call this routine 3560 themselves. 3561 3562 Level: developer 3563 3564 .keywords: TS, nonlinear, compute, function 3565 3566 .seealso: TSSetFunction(), TSGetFunction() 3567 */ 3568 PetscErrorCode TSComputeFunction_Matlab(TS snes,PetscReal time,Vec x,Vec xdot,Vec y, void *ctx) 3569 { 3570 PetscErrorCode ierr; 3571 TSMatlabContext *sctx = (TSMatlabContext *)ctx; 3572 int nlhs = 1,nrhs = 7; 3573 mxArray *plhs[1],*prhs[7]; 3574 long long int lx = 0,lxdot = 0,ly = 0,ls = 0; 3575 3576 PetscFunctionBegin; 3577 PetscValidHeaderSpecific(snes,TS_CLASSID,1); 3578 PetscValidHeaderSpecific(x,VEC_CLASSID,3); 3579 PetscValidHeaderSpecific(xdot,VEC_CLASSID,4); 3580 PetscValidHeaderSpecific(y,VEC_CLASSID,5); 3581 PetscCheckSameComm(snes,1,x,3); 3582 PetscCheckSameComm(snes,1,y,5); 3583 3584 ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 3585 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 3586 ierr = PetscMemcpy(&lxdot,&xdot,sizeof(xdot));CHKERRQ(ierr); 3587 ierr = PetscMemcpy(&ly,&y,sizeof(x));CHKERRQ(ierr); 3588 prhs[0] = mxCreateDoubleScalar((double)ls); 3589 prhs[1] = mxCreateDoubleScalar(time); 3590 prhs[2] = mxCreateDoubleScalar((double)lx); 3591 prhs[3] = mxCreateDoubleScalar((double)lxdot); 3592 prhs[4] = mxCreateDoubleScalar((double)ly); 3593 prhs[5] = mxCreateString(sctx->funcname); 3594 prhs[6] = sctx->ctx; 3595 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeFunctionInternal");CHKERRQ(ierr); 3596 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 3597 mxDestroyArray(prhs[0]); 3598 mxDestroyArray(prhs[1]); 3599 mxDestroyArray(prhs[2]); 3600 mxDestroyArray(prhs[3]); 3601 mxDestroyArray(prhs[4]); 3602 mxDestroyArray(prhs[5]); 3603 mxDestroyArray(plhs[0]); 3604 PetscFunctionReturn(0); 3605 } 3606 3607 3608 #undef __FUNCT__ 3609 #define __FUNCT__ "TSSetFunctionMatlab" 3610 /* 3611 TSSetFunctionMatlab - Sets the function evaluation routine and function 3612 vector for use by the TS routines in solving ODEs 3613 equations from MATLAB. Here the function is a string containing the name of a MATLAB function 3614 3615 Logically Collective on TS 3616 3617 Input Parameters: 3618 + ts - the TS context 3619 - func - function evaluation routine 3620 3621 Calling sequence of func: 3622 $ func (TS ts,PetscReal time,Vec x,Vec xdot,Vec f,void *ctx); 3623 3624 Level: beginner 3625 3626 .keywords: TS, nonlinear, set, function 3627 3628 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction() 3629 */ 3630 PetscErrorCode TSSetFunctionMatlab(TS ts,const char *func,mxArray *ctx) 3631 { 3632 PetscErrorCode ierr; 3633 TSMatlabContext *sctx; 3634 3635 PetscFunctionBegin; 3636 /* currently sctx is memory bleed */ 3637 ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr); 3638 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 3639 /* 3640 This should work, but it doesn't 3641 sctx->ctx = ctx; 3642 mexMakeArrayPersistent(sctx->ctx); 3643 */ 3644 sctx->ctx = mxDuplicateArray(ctx); 3645 ierr = TSSetIFunction(ts,PETSC_NULL,TSComputeFunction_Matlab,sctx);CHKERRQ(ierr); 3646 PetscFunctionReturn(0); 3647 } 3648 3649 #undef __FUNCT__ 3650 #define __FUNCT__ "TSComputeJacobian_Matlab" 3651 /* 3652 TSComputeJacobian_Matlab - Calls the function that has been set with 3653 TSSetJacobianMatlab(). 3654 3655 Collective on TS 3656 3657 Input Parameters: 3658 + ts - the TS context 3659 . x - input vector 3660 . A, B - the matrices 3661 - ctx - user context 3662 3663 Output Parameter: 3664 . flag - structure of the matrix 3665 3666 Level: developer 3667 3668 .keywords: TS, nonlinear, compute, function 3669 3670 .seealso: TSSetFunction(), TSGetFunction() 3671 @*/ 3672 PetscErrorCode TSComputeJacobian_Matlab(TS ts,PetscReal time,Vec x,Vec xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flag, void *ctx) 3673 { 3674 PetscErrorCode ierr; 3675 TSMatlabContext *sctx = (TSMatlabContext *)ctx; 3676 int nlhs = 2,nrhs = 9; 3677 mxArray *plhs[2],*prhs[9]; 3678 long long int lx = 0,lxdot = 0,lA = 0,ls = 0, lB = 0; 3679 3680 PetscFunctionBegin; 3681 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3682 PetscValidHeaderSpecific(x,VEC_CLASSID,3); 3683 3684 /* call Matlab function in ctx with arguments x and y */ 3685 3686 ierr = PetscMemcpy(&ls,&ts,sizeof(ts));CHKERRQ(ierr); 3687 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 3688 ierr = PetscMemcpy(&lxdot,&xdot,sizeof(x));CHKERRQ(ierr); 3689 ierr = PetscMemcpy(&lA,A,sizeof(x));CHKERRQ(ierr); 3690 ierr = PetscMemcpy(&lB,B,sizeof(x));CHKERRQ(ierr); 3691 prhs[0] = mxCreateDoubleScalar((double)ls); 3692 prhs[1] = mxCreateDoubleScalar((double)time); 3693 prhs[2] = mxCreateDoubleScalar((double)lx); 3694 prhs[3] = mxCreateDoubleScalar((double)lxdot); 3695 prhs[4] = mxCreateDoubleScalar((double)shift); 3696 prhs[5] = mxCreateDoubleScalar((double)lA); 3697 prhs[6] = mxCreateDoubleScalar((double)lB); 3698 prhs[7] = mxCreateString(sctx->funcname); 3699 prhs[8] = sctx->ctx; 3700 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeJacobianInternal");CHKERRQ(ierr); 3701 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 3702 *flag = (MatStructure) mxGetScalar(plhs[1]);CHKERRQ(ierr); 3703 mxDestroyArray(prhs[0]); 3704 mxDestroyArray(prhs[1]); 3705 mxDestroyArray(prhs[2]); 3706 mxDestroyArray(prhs[3]); 3707 mxDestroyArray(prhs[4]); 3708 mxDestroyArray(prhs[5]); 3709 mxDestroyArray(prhs[6]); 3710 mxDestroyArray(prhs[7]); 3711 mxDestroyArray(plhs[0]); 3712 mxDestroyArray(plhs[1]); 3713 PetscFunctionReturn(0); 3714 } 3715 3716 3717 #undef __FUNCT__ 3718 #define __FUNCT__ "TSSetJacobianMatlab" 3719 /* 3720 TSSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices 3721 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 3722 3723 Logically Collective on TS 3724 3725 Input Parameters: 3726 + ts - the TS context 3727 . A,B - Jacobian matrices 3728 . func - function evaluation routine 3729 - ctx - user context 3730 3731 Calling sequence of func: 3732 $ flag = func (TS ts,PetscReal time,Vec x,Vec xdot,Mat A,Mat B,void *ctx); 3733 3734 3735 Level: developer 3736 3737 .keywords: TS, nonlinear, set, function 3738 3739 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction() 3740 */ 3741 PetscErrorCode TSSetJacobianMatlab(TS ts,Mat A,Mat B,const char *func,mxArray *ctx) 3742 { 3743 PetscErrorCode ierr; 3744 TSMatlabContext *sctx; 3745 3746 PetscFunctionBegin; 3747 /* currently sctx is memory bleed */ 3748 ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr); 3749 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 3750 /* 3751 This should work, but it doesn't 3752 sctx->ctx = ctx; 3753 mexMakeArrayPersistent(sctx->ctx); 3754 */ 3755 sctx->ctx = mxDuplicateArray(ctx); 3756 ierr = TSSetIJacobian(ts,A,B,TSComputeJacobian_Matlab,sctx);CHKERRQ(ierr); 3757 PetscFunctionReturn(0); 3758 } 3759 3760 #undef __FUNCT__ 3761 #define __FUNCT__ "TSMonitor_Matlab" 3762 /* 3763 TSMonitor_Matlab - Calls the function that has been set with TSMonitorSetMatlab(). 3764 3765 Collective on TS 3766 3767 .seealso: TSSetFunction(), TSGetFunction() 3768 @*/ 3769 PetscErrorCode TSMonitor_Matlab(TS ts,PetscInt it, PetscReal time,Vec x, void *ctx) 3770 { 3771 PetscErrorCode ierr; 3772 TSMatlabContext *sctx = (TSMatlabContext *)ctx; 3773 int nlhs = 1,nrhs = 6; 3774 mxArray *plhs[1],*prhs[6]; 3775 long long int lx = 0,ls = 0; 3776 3777 PetscFunctionBegin; 3778 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3779 PetscValidHeaderSpecific(x,VEC_CLASSID,4); 3780 3781 ierr = PetscMemcpy(&ls,&ts,sizeof(ts));CHKERRQ(ierr); 3782 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 3783 prhs[0] = mxCreateDoubleScalar((double)ls); 3784 prhs[1] = mxCreateDoubleScalar((double)it); 3785 prhs[2] = mxCreateDoubleScalar((double)time); 3786 prhs[3] = mxCreateDoubleScalar((double)lx); 3787 prhs[4] = mxCreateString(sctx->funcname); 3788 prhs[5] = sctx->ctx; 3789 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSMonitorInternal");CHKERRQ(ierr); 3790 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 3791 mxDestroyArray(prhs[0]); 3792 mxDestroyArray(prhs[1]); 3793 mxDestroyArray(prhs[2]); 3794 mxDestroyArray(prhs[3]); 3795 mxDestroyArray(prhs[4]); 3796 mxDestroyArray(plhs[0]); 3797 PetscFunctionReturn(0); 3798 } 3799 3800 3801 #undef __FUNCT__ 3802 #define __FUNCT__ "TSMonitorSetMatlab" 3803 /* 3804 TSMonitorSetMatlab - Sets the monitor function from Matlab 3805 3806 Level: developer 3807 3808 .keywords: TS, nonlinear, set, function 3809 3810 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction() 3811 */ 3812 PetscErrorCode TSMonitorSetMatlab(TS ts,const char *func,mxArray *ctx) 3813 { 3814 PetscErrorCode ierr; 3815 TSMatlabContext *sctx; 3816 3817 PetscFunctionBegin; 3818 /* currently sctx is memory bleed */ 3819 ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr); 3820 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 3821 /* 3822 This should work, but it doesn't 3823 sctx->ctx = ctx; 3824 mexMakeArrayPersistent(sctx->ctx); 3825 */ 3826 sctx->ctx = mxDuplicateArray(ctx); 3827 ierr = TSMonitorSet(ts,TSMonitor_Matlab,sctx,PETSC_NULL);CHKERRQ(ierr); 3828 PetscFunctionReturn(0); 3829 } 3830 #endif 3831