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