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