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