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