1 #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 2 #include <petsc/private/tshistoryimpl.h> 3 #include <petscdm.h> 4 5 PetscFunctionList TSTrajectoryList = NULL; 6 PetscBool TSTrajectoryRegisterAllCalled = PETSC_FALSE; 7 PetscClassId TSTRAJECTORY_CLASSID; 8 PetscLogEvent TSTrajectory_Set, TSTrajectory_Get, TSTrajectory_GetVecs, TSTrajectory_SetUp; 9 10 /*@C 11 TSTrajectoryRegister - Adds a way of storing trajectories to the TS package 12 13 Not Collective 14 15 Input Parameters: 16 + name - the name of a new user-defined creation routine 17 - create_func - the creation routine itself 18 19 Notes: 20 TSTrajectoryRegister() may be called multiple times to add several user-defined tses. 21 22 Level: developer 23 24 .seealso: TSTrajectoryRegisterAll() 25 @*/ 26 PetscErrorCode TSTrajectoryRegister(const char sname[],PetscErrorCode (*function)(TSTrajectory,TS)) 27 { 28 PetscFunctionBegin; 29 CHKERRQ(PetscFunctionListAdd(&TSTrajectoryList,sname,function)); 30 PetscFunctionReturn(0); 31 } 32 33 /*@ 34 TSTrajectorySet - Sets a vector of state in the trajectory object 35 36 Collective on TSTrajectory 37 38 Input Parameters: 39 + tj - the trajectory object 40 . ts - the time stepper object (optional) 41 . stepnum - the step number 42 . time - the current time 43 - X - the current solution 44 45 Level: developer 46 47 Notes: Usually one does not call this routine, it is called automatically during TSSolve() 48 49 .seealso: TSTrajectorySetUp(), TSTrajectoryDestroy(), TSTrajectorySetType(), TSTrajectorySetVariableNames(), TSGetTrajectory(), TSTrajectoryGet(), TSTrajectoryGetVecs() 50 @*/ 51 PetscErrorCode TSTrajectorySet(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal time,Vec X) 52 { 53 PetscFunctionBegin; 54 if (!tj) PetscFunctionReturn(0); 55 PetscValidHeaderSpecific(tj,TSTRAJECTORY_CLASSID,1); 56 if (ts) PetscValidHeaderSpecific(ts,TS_CLASSID,2); 57 PetscValidLogicalCollectiveInt(tj,stepnum,3); 58 PetscValidLogicalCollectiveReal(tj,time,4); 59 PetscValidHeaderSpecific(X,VEC_CLASSID,5); 60 PetscCheck(tj->ops->set,PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"TSTrajectory type %s",((PetscObject)tj)->type_name); 61 PetscCheck(tj->setupcalled,PetscObjectComm((PetscObject)tj),PETSC_ERR_ORDER,"TSTrajectorySetUp should be called first"); 62 if (tj->monitor) { 63 CHKERRQ(PetscViewerASCIIPrintf(tj->monitor,"TSTrajectorySet: stepnum %D, time %g (stages %D)\n",stepnum,(double)time,(PetscInt)!tj->solution_only)); 64 } 65 CHKERRQ(PetscLogEventBegin(TSTrajectory_Set,tj,ts,0,0)); 66 CHKERRQ((*tj->ops->set)(tj,ts,stepnum,time,X)); 67 CHKERRQ(PetscLogEventEnd(TSTrajectory_Set,tj,ts,0,0)); 68 if (tj->usehistory) { 69 CHKERRQ(TSHistoryUpdate(tj->tsh,stepnum,time)); 70 } 71 if (tj->lag.caching) tj->lag.Udotcached.time = PETSC_MIN_REAL; 72 PetscFunctionReturn(0); 73 } 74 75 /*@ 76 TSTrajectoryGetNumSteps - Return the number of steps registered in the TSTrajectory via TSTrajectorySet(). 77 78 Not collective. 79 80 Input Parameters: 81 . tj - the trajectory object 82 83 Output Parameter: 84 . steps - the number of steps 85 86 Level: developer 87 88 .seealso: TSTrajectorySet() 89 @*/ 90 PetscErrorCode TSTrajectoryGetNumSteps(TSTrajectory tj, PetscInt *steps) 91 { 92 PetscFunctionBegin; 93 PetscValidHeaderSpecific(tj,TSTRAJECTORY_CLASSID,1); 94 PetscValidIntPointer(steps,2); 95 CHKERRQ(TSHistoryGetNumSteps(tj->tsh,steps)); 96 PetscFunctionReturn(0); 97 } 98 99 /*@ 100 TSTrajectoryGet - Updates the solution vector of a time stepper object by inquiring the TSTrajectory 101 102 Collective on TS 103 104 Input Parameters: 105 + tj - the trajectory object 106 . ts - the time stepper object 107 - stepnum - the step number 108 109 Output Parameter: 110 . time - the time associated with the step number 111 112 Level: developer 113 114 Notes: Usually one does not call this routine, it is called automatically during TSSolve() 115 116 .seealso: TSTrajectorySetUp(), TSTrajectoryDestroy(), TSTrajectorySetType(), TSTrajectorySetVariableNames(), TSGetTrajectory(), TSTrajectorySet(), TSTrajectoryGetVecs(), TSGetSolution() 117 @*/ 118 PetscErrorCode TSTrajectoryGet(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal *time) 119 { 120 PetscFunctionBegin; 121 PetscCheck(tj,PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"TS solver did not save trajectory"); 122 PetscValidHeaderSpecific(tj,TSTRAJECTORY_CLASSID,1); 123 PetscValidHeaderSpecific(ts,TS_CLASSID,2); 124 PetscValidLogicalCollectiveInt(tj,stepnum,3); 125 PetscValidPointer(time,4); 126 PetscCheck(tj->ops->get,PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"TSTrajectory type %s",((PetscObject)tj)->type_name); 127 PetscCheck(tj->setupcalled,PetscObjectComm((PetscObject)tj),PETSC_ERR_ORDER,"TSTrajectorySetUp should be called first"); 128 PetscCheck(stepnum >= 0,PetscObjectComm((PetscObject)tj),PETSC_ERR_PLIB,"Requesting negative step number"); 129 if (tj->monitor) { 130 CHKERRQ(PetscViewerASCIIPrintf(tj->monitor,"TSTrajectoryGet: stepnum %D, stages %D\n",stepnum,(PetscInt)!tj->solution_only)); 131 CHKERRQ(PetscViewerFlush(tj->monitor)); 132 } 133 CHKERRQ(PetscLogEventBegin(TSTrajectory_Get,tj,ts,0,0)); 134 CHKERRQ((*tj->ops->get)(tj,ts,stepnum,time)); 135 CHKERRQ(PetscLogEventEnd(TSTrajectory_Get,tj,ts,0,0)); 136 PetscFunctionReturn(0); 137 } 138 139 /*@ 140 TSTrajectoryGetVecs - Reconstructs the vector of state and its time derivative using information from the TSTrajectory and, possibly, from the TS 141 142 Collective on TS 143 144 Input Parameters: 145 + tj - the trajectory object 146 . ts - the time stepper object (optional) 147 - stepnum - the requested step number 148 149 Input/Output Parameter: 150 151 Output Parameters: 152 + time - On input time for the step if step number is PETSC_DECIDE, on output the time associated with the step number 153 . U - state vector (can be NULL) 154 - Udot - time derivative of state vector (can be NULL) 155 156 Level: developer 157 158 Notes: If the step number is PETSC_DECIDE, the time argument is used to inquire the trajectory. 159 If the requested time does not match any in the trajectory, Lagrangian interpolations are returned. 160 161 .seealso: TSTrajectorySetUp(), TSTrajectoryDestroy(), TSTrajectorySetType(), TSTrajectorySetVariableNames(), TSGetTrajectory(), TSTrajectorySet(), TSTrajectoryGet() 162 @*/ 163 PetscErrorCode TSTrajectoryGetVecs(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal *time,Vec U,Vec Udot) 164 { 165 PetscFunctionBegin; 166 PetscCheck(tj,PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"TS solver did not save trajectory"); 167 PetscValidHeaderSpecific(tj,TSTRAJECTORY_CLASSID,1); 168 if (ts) PetscValidHeaderSpecific(ts,TS_CLASSID,2); 169 PetscValidLogicalCollectiveInt(tj,stepnum,3); 170 PetscValidPointer(time,4); 171 if (U) PetscValidHeaderSpecific(U,VEC_CLASSID,5); 172 if (Udot) PetscValidHeaderSpecific(Udot,VEC_CLASSID,6); 173 if (!U && !Udot) PetscFunctionReturn(0); 174 PetscCheck(tj->setupcalled,PetscObjectComm((PetscObject)tj),PETSC_ERR_ORDER,"TSTrajectorySetUp should be called first"); 175 CHKERRQ(PetscLogEventBegin(TSTrajectory_GetVecs,tj,ts,0,0)); 176 if (tj->monitor) { 177 PetscInt pU,pUdot; 178 pU = U ? 1 : 0; 179 pUdot = Udot ? 1 : 0; 180 CHKERRQ(PetscViewerASCIIPrintf(tj->monitor,"Requested by GetVecs %D %D: stepnum %D, time %g\n",pU,pUdot,stepnum,(double)*time)); 181 CHKERRQ(PetscViewerFlush(tj->monitor)); 182 } 183 if (U && tj->lag.caching) { 184 PetscObjectId id; 185 PetscObjectState state; 186 187 CHKERRQ(PetscObjectStateGet((PetscObject)U,&state)); 188 CHKERRQ(PetscObjectGetId((PetscObject)U,&id)); 189 if (stepnum == PETSC_DECIDE) { 190 if (id == tj->lag.Ucached.id && *time == tj->lag.Ucached.time && state == tj->lag.Ucached.state) U = NULL; 191 } else { 192 if (id == tj->lag.Ucached.id && stepnum == tj->lag.Ucached.step && state == tj->lag.Ucached.state) U = NULL; 193 } 194 if (tj->monitor && !U) { 195 CHKERRQ(PetscViewerASCIIPushTab(tj->monitor)); 196 CHKERRQ(PetscViewerASCIIPrintf(tj->monitor,"State vector cached\n")); 197 CHKERRQ(PetscViewerASCIIPopTab(tj->monitor)); 198 CHKERRQ(PetscViewerFlush(tj->monitor)); 199 } 200 } 201 if (Udot && tj->lag.caching) { 202 PetscObjectId id; 203 PetscObjectState state; 204 205 CHKERRQ(PetscObjectStateGet((PetscObject)Udot,&state)); 206 CHKERRQ(PetscObjectGetId((PetscObject)Udot,&id)); 207 if (stepnum == PETSC_DECIDE) { 208 if (id == tj->lag.Udotcached.id && *time == tj->lag.Udotcached.time && state == tj->lag.Udotcached.state) Udot = NULL; 209 } else { 210 if (id == tj->lag.Udotcached.id && stepnum == tj->lag.Udotcached.step && state == tj->lag.Udotcached.state) Udot = NULL; 211 } 212 if (tj->monitor && !Udot) { 213 CHKERRQ(PetscViewerASCIIPushTab(tj->monitor)); 214 CHKERRQ(PetscViewerASCIIPrintf(tj->monitor,"Derivative vector cached\n")); 215 CHKERRQ(PetscViewerASCIIPopTab(tj->monitor)); 216 CHKERRQ(PetscViewerFlush(tj->monitor)); 217 } 218 } 219 if (!U && !Udot) { 220 CHKERRQ(PetscLogEventEnd(TSTrajectory_GetVecs,tj,ts,0,0)); 221 PetscFunctionReturn(0); 222 } 223 224 if (stepnum == PETSC_DECIDE || Udot) { /* reverse search for requested time in TSHistory */ 225 if (tj->monitor) { 226 CHKERRQ(PetscViewerASCIIPushTab(tj->monitor)); 227 } 228 /* cached states will be updated in the function */ 229 CHKERRQ(TSTrajectoryReconstruct_Private(tj,ts,*time,U,Udot)); 230 if (tj->monitor) { 231 CHKERRQ(PetscViewerASCIIPopTab(tj->monitor)); 232 CHKERRQ(PetscViewerFlush(tj->monitor)); 233 } 234 } else if (U) { /* we were asked to load from stepnum, use TSTrajectoryGet */ 235 TS fakets = ts; 236 Vec U2; 237 238 /* use a fake TS if ts is missing */ 239 if (!ts) { 240 CHKERRQ(PetscObjectQuery((PetscObject)tj,"__fake_ts",(PetscObject*)&fakets)); 241 if (!fakets) { 242 CHKERRQ(TSCreate(PetscObjectComm((PetscObject)tj),&fakets)); 243 CHKERRQ(PetscObjectCompose((PetscObject)tj,"__fake_ts",(PetscObject)fakets)); 244 CHKERRQ(PetscObjectDereference((PetscObject)fakets)); 245 CHKERRQ(VecDuplicate(U,&U2)); 246 CHKERRQ(TSSetSolution(fakets,U2)); 247 CHKERRQ(PetscObjectDereference((PetscObject)U2)); 248 } 249 } 250 CHKERRQ(TSTrajectoryGet(tj,fakets,stepnum,time)); 251 CHKERRQ(TSGetSolution(fakets,&U2)); 252 CHKERRQ(VecCopy(U2,U)); 253 CHKERRQ(PetscObjectStateGet((PetscObject)U,&tj->lag.Ucached.state)); 254 CHKERRQ(PetscObjectGetId((PetscObject)U,&tj->lag.Ucached.id)); 255 tj->lag.Ucached.time = *time; 256 tj->lag.Ucached.step = stepnum; 257 } 258 CHKERRQ(PetscLogEventEnd(TSTrajectory_GetVecs,tj,ts,0,0)); 259 PetscFunctionReturn(0); 260 } 261 262 /*@C 263 TSTrajectoryViewFromOptions - View from Options 264 265 Collective on TSTrajectory 266 267 Input Parameters: 268 + A - the TSTrajectory context 269 . obj - Optional object 270 - name - command line option 271 272 Level: intermediate 273 .seealso: TSTrajectory, TSTrajectoryView, PetscObjectViewFromOptions(), TSTrajectoryCreate() 274 @*/ 275 PetscErrorCode TSTrajectoryViewFromOptions(TSTrajectory A,PetscObject obj,const char name[]) 276 { 277 PetscFunctionBegin; 278 PetscValidHeaderSpecific(A,TSTRAJECTORY_CLASSID,1); 279 CHKERRQ(PetscObjectViewFromOptions((PetscObject)A,obj,name)); 280 PetscFunctionReturn(0); 281 } 282 283 /*@C 284 TSTrajectoryView - Prints information about the trajectory object 285 286 Collective on TSTrajectory 287 288 Input Parameters: 289 + tj - the TSTrajectory context obtained from TSTrajectoryCreate() 290 - viewer - visualization context 291 292 Options Database Key: 293 . -ts_trajectory_view - calls TSTrajectoryView() at end of TSAdjointStep() 294 295 Notes: 296 The available visualization contexts include 297 + PETSC_VIEWER_STDOUT_SELF - standard output (default) 298 - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 299 output where only the first processor opens 300 the file. All other processors send their 301 data to the first processor to print. 302 303 The user can open an alternative visualization context with 304 PetscViewerASCIIOpen() - output to a specified file. 305 306 Level: developer 307 308 .seealso: PetscViewerASCIIOpen() 309 @*/ 310 PetscErrorCode TSTrajectoryView(TSTrajectory tj,PetscViewer viewer) 311 { 312 PetscBool iascii; 313 314 PetscFunctionBegin; 315 PetscValidHeaderSpecific(tj,TSTRAJECTORY_CLASSID,1); 316 if (!viewer) { 317 CHKERRQ(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)tj),&viewer)); 318 } 319 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 320 PetscCheckSameComm(tj,1,viewer,2); 321 322 CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 323 if (iascii) { 324 CHKERRQ(PetscObjectPrintClassNamePrefixType((PetscObject)tj,viewer)); 325 CHKERRQ(PetscViewerASCIIPrintf(viewer," total number of recomputations for adjoint calculation = %D\n",tj->recomps)); 326 CHKERRQ(PetscViewerASCIIPrintf(viewer," disk checkpoint reads = %D\n",tj->diskreads)); 327 CHKERRQ(PetscViewerASCIIPrintf(viewer," disk checkpoint writes = %D\n",tj->diskwrites)); 328 if (tj->ops->view) { 329 CHKERRQ(PetscViewerASCIIPushTab(viewer)); 330 CHKERRQ((*tj->ops->view)(tj,viewer)); 331 CHKERRQ(PetscViewerASCIIPopTab(viewer)); 332 } 333 } 334 PetscFunctionReturn(0); 335 } 336 337 /*@C 338 TSTrajectorySetVariableNames - Sets the name of each component in the solution vector so that it may be saved with the trajectory 339 340 Collective on TSTrajectory 341 342 Input Parameters: 343 + tr - the trajectory context 344 - names - the names of the components, final string must be NULL 345 346 Level: intermediate 347 348 Note: Fortran interface is not possible because of the string array argument 349 350 .seealso: TSTrajectory, TSGetTrajectory() 351 @*/ 352 PetscErrorCode TSTrajectorySetVariableNames(TSTrajectory ctx,const char * const *names) 353 { 354 PetscFunctionBegin; 355 PetscValidHeaderSpecific(ctx,TSTRAJECTORY_CLASSID,1); 356 PetscValidPointer(names,2); 357 CHKERRQ(PetscStrArrayDestroy(&ctx->names)); 358 CHKERRQ(PetscStrArrayallocpy(names,&ctx->names)); 359 PetscFunctionReturn(0); 360 } 361 362 /*@C 363 TSTrajectorySetTransform - Solution vector will be transformed by provided function before being saved to disk 364 365 Collective on TSLGCtx 366 367 Input Parameters: 368 + tj - the TSTrajectory context 369 . transform - the transform function 370 . destroy - function to destroy the optional context 371 - ctx - optional context used by transform function 372 373 Level: intermediate 374 375 .seealso: TSTrajectorySetVariableNames(), TSTrajectory, TSMonitorLGSetTransform() 376 @*/ 377 PetscErrorCode TSTrajectorySetTransform(TSTrajectory tj,PetscErrorCode (*transform)(void*,Vec,Vec*),PetscErrorCode (*destroy)(void*),void *tctx) 378 { 379 PetscFunctionBegin; 380 PetscValidHeaderSpecific(tj,TSTRAJECTORY_CLASSID,1); 381 tj->transform = transform; 382 tj->transformdestroy = destroy; 383 tj->transformctx = tctx; 384 PetscFunctionReturn(0); 385 } 386 387 /*@ 388 TSTrajectoryCreate - This function creates an empty trajectory object used to store the time dependent solution of an ODE/DAE 389 390 Collective 391 392 Input Parameter: 393 . comm - the communicator 394 395 Output Parameter: 396 . tj - the trajectory object 397 398 Level: developer 399 400 Notes: 401 Usually one does not call this routine, it is called automatically when one calls TSSetSaveTrajectory(). 402 403 .seealso: TSTrajectorySetUp(), TSTrajectoryDestroy(), TSTrajectorySetType(), TSTrajectorySetVariableNames(), TSGetTrajectory(), TSTrajectorySetKeepFiles() 404 @*/ 405 PetscErrorCode TSTrajectoryCreate(MPI_Comm comm,TSTrajectory *tj) 406 { 407 TSTrajectory t; 408 409 PetscFunctionBegin; 410 PetscValidPointer(tj,2); 411 *tj = NULL; 412 CHKERRQ(TSInitializePackage()); 413 414 CHKERRQ(PetscHeaderCreate(t,TSTRAJECTORY_CLASSID,"TSTrajectory","Time stepping","TS",comm,TSTrajectoryDestroy,TSTrajectoryView)); 415 t->setupcalled = PETSC_FALSE; 416 CHKERRQ(TSHistoryCreate(comm,&t->tsh)); 417 418 t->lag.order = 1; 419 t->lag.L = NULL; 420 t->lag.T = NULL; 421 t->lag.W = NULL; 422 t->lag.WW = NULL; 423 t->lag.TW = NULL; 424 t->lag.TT = NULL; 425 t->lag.caching = PETSC_TRUE; 426 t->lag.Ucached.id = 0; 427 t->lag.Ucached.state = -1; 428 t->lag.Ucached.time = PETSC_MIN_REAL; 429 t->lag.Ucached.step = PETSC_MAX_INT; 430 t->lag.Udotcached.id = 0; 431 t->lag.Udotcached.state = -1; 432 t->lag.Udotcached.time = PETSC_MIN_REAL; 433 t->lag.Udotcached.step = PETSC_MAX_INT; 434 t->adjoint_solve_mode = PETSC_TRUE; 435 t->solution_only = PETSC_FALSE; 436 t->keepfiles = PETSC_FALSE; 437 t->usehistory = PETSC_TRUE; 438 *tj = t; 439 CHKERRQ(TSTrajectorySetFiletemplate(t,"TS-%06D.bin")); 440 PetscFunctionReturn(0); 441 } 442 443 /*@C 444 TSTrajectorySetType - Sets the storage method to be used as in a trajectory 445 446 Collective on TS 447 448 Input Parameters: 449 + tj - the TSTrajectory context 450 . ts - the TS context 451 - type - a known method 452 453 Options Database Command: 454 . -ts_trajectory_type <type> - Sets the method; use -help for a list of available methods (for instance, basic) 455 456 Level: developer 457 458 .seealso: TS, TSTrajectoryCreate(), TSTrajectorySetFromOptions(), TSTrajectoryDestroy(), TSTrajectoryGetType() 459 460 @*/ 461 PetscErrorCode TSTrajectorySetType(TSTrajectory tj,TS ts,TSTrajectoryType type) 462 { 463 PetscErrorCode (*r)(TSTrajectory,TS); 464 PetscBool match; 465 466 PetscFunctionBegin; 467 PetscValidHeaderSpecific(tj,TSTRAJECTORY_CLASSID,1); 468 CHKERRQ(PetscObjectTypeCompare((PetscObject)tj,type,&match)); 469 if (match) PetscFunctionReturn(0); 470 471 CHKERRQ(PetscFunctionListFind(TSTrajectoryList,type,&r)); 472 PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSTrajectory type: %s",type); 473 if (tj->ops->destroy) { 474 CHKERRQ((*(tj)->ops->destroy)(tj)); 475 476 tj->ops->destroy = NULL; 477 } 478 CHKERRQ(PetscMemzero(tj->ops,sizeof(*tj->ops))); 479 480 CHKERRQ(PetscObjectChangeTypeName((PetscObject)tj,type)); 481 CHKERRQ((*r)(tj,ts)); 482 PetscFunctionReturn(0); 483 } 484 485 /*@C 486 TSTrajectoryGetType - Gets the trajectory type 487 488 Collective on TS 489 490 Input Parameters: 491 + tj - the TSTrajectory context 492 - ts - the TS context 493 494 Output Parameters: 495 . type - a known method 496 497 Level: developer 498 499 .seealso: TS, TSTrajectoryCreate(), TSTrajectorySetFromOptions(), TSTrajectoryDestroy(), TSTrajectorySetType() 500 501 @*/ 502 PetscErrorCode TSTrajectoryGetType(TSTrajectory tj,TS ts,TSTrajectoryType *type) 503 { 504 PetscFunctionBegin; 505 PetscValidHeaderSpecific(tj,TSTRAJECTORY_CLASSID,1); 506 if (type) *type = ((PetscObject)tj)->type_name; 507 PetscFunctionReturn(0); 508 } 509 510 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Basic(TSTrajectory,TS); 511 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Singlefile(TSTrajectory,TS); 512 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Memory(TSTrajectory,TS); 513 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Visualization(TSTrajectory,TS); 514 515 /*@C 516 TSTrajectoryRegisterAll - Registers all of the trajectory storage schecmes in the TS package. 517 518 Not Collective 519 520 Level: developer 521 522 .seealso: TSTrajectoryRegister() 523 @*/ 524 PetscErrorCode TSTrajectoryRegisterAll(void) 525 { 526 PetscFunctionBegin; 527 if (TSTrajectoryRegisterAllCalled) PetscFunctionReturn(0); 528 TSTrajectoryRegisterAllCalled = PETSC_TRUE; 529 530 CHKERRQ(TSTrajectoryRegister(TSTRAJECTORYBASIC,TSTrajectoryCreate_Basic)); 531 CHKERRQ(TSTrajectoryRegister(TSTRAJECTORYSINGLEFILE,TSTrajectoryCreate_Singlefile)); 532 CHKERRQ(TSTrajectoryRegister(TSTRAJECTORYMEMORY,TSTrajectoryCreate_Memory)); 533 CHKERRQ(TSTrajectoryRegister(TSTRAJECTORYVISUALIZATION,TSTrajectoryCreate_Visualization)); 534 PetscFunctionReturn(0); 535 } 536 537 /*@ 538 TSTrajectoryReset - Resets a trajectory context 539 540 Collective on TSTrajectory 541 542 Input Parameter: 543 . tj - the TSTrajectory context obtained from TSTrajectoryCreate() 544 545 Level: developer 546 547 .seealso: TSTrajectoryCreate(), TSTrajectorySetUp() 548 @*/ 549 PetscErrorCode TSTrajectoryReset(TSTrajectory tj) 550 { 551 PetscFunctionBegin; 552 if (!tj) PetscFunctionReturn(0); 553 PetscValidHeaderSpecific(tj,TSTRAJECTORY_CLASSID,1); 554 if (tj->ops->reset) { 555 CHKERRQ((*tj->ops->reset)(tj)); 556 } 557 CHKERRQ(PetscFree(tj->dirfiletemplate)); 558 CHKERRQ(TSHistoryDestroy(&tj->tsh)); 559 CHKERRQ(TSHistoryCreate(PetscObjectComm((PetscObject)tj),&tj->tsh)); 560 tj->setupcalled = PETSC_FALSE; 561 PetscFunctionReturn(0); 562 } 563 564 /*@ 565 TSTrajectoryDestroy - Destroys a trajectory context 566 567 Collective on TSTrajectory 568 569 Input Parameter: 570 . tj - the TSTrajectory context obtained from TSTrajectoryCreate() 571 572 Level: developer 573 574 .seealso: TSTrajectoryCreate(), TSTrajectorySetUp() 575 @*/ 576 PetscErrorCode TSTrajectoryDestroy(TSTrajectory *tj) 577 { 578 PetscFunctionBegin; 579 if (!*tj) PetscFunctionReturn(0); 580 PetscValidHeaderSpecific((*tj),TSTRAJECTORY_CLASSID,1); 581 if (--((PetscObject)(*tj))->refct > 0) {*tj = NULL; PetscFunctionReturn(0);} 582 583 CHKERRQ(TSTrajectoryReset(*tj)); 584 CHKERRQ(TSHistoryDestroy(&(*tj)->tsh)); 585 CHKERRQ(VecDestroyVecs((*tj)->lag.order+1,&(*tj)->lag.W)); 586 CHKERRQ(PetscFree5((*tj)->lag.L,(*tj)->lag.T,(*tj)->lag.WW,(*tj)->lag.TT,(*tj)->lag.TW)); 587 CHKERRQ(VecDestroy(&(*tj)->U)); 588 CHKERRQ(VecDestroy(&(*tj)->Udot)); 589 590 if ((*tj)->transformdestroy) CHKERRQ((*(*tj)->transformdestroy)((*tj)->transformctx)); 591 if ((*tj)->ops->destroy) CHKERRQ((*(*tj)->ops->destroy)((*tj))); 592 if (!((*tj)->keepfiles)) { 593 PetscMPIInt rank; 594 MPI_Comm comm; 595 596 CHKERRQ(PetscObjectGetComm((PetscObject)(*tj),&comm)); 597 CHKERRMPI(MPI_Comm_rank(comm,&rank)); 598 if (rank == 0 && (*tj)->dirname) { /* we own the directory, so we run PetscRMTree on it */ 599 CHKERRQ(PetscRMTree((*tj)->dirname)); 600 } 601 } 602 CHKERRQ(PetscStrArrayDestroy(&(*tj)->names)); 603 CHKERRQ(PetscFree((*tj)->dirname)); 604 CHKERRQ(PetscFree((*tj)->filetemplate)); 605 CHKERRQ(PetscHeaderDestroy(tj)); 606 PetscFunctionReturn(0); 607 } 608 609 /* 610 TSTrajectorySetTypeFromOptions_Private - Sets the type of ts from user options. 611 612 Collective on TSTrajectory 613 614 Input Parameter: 615 + tj - the TSTrajectory context 616 - ts - the TS context 617 618 Options Database Keys: 619 . -ts_trajectory_type <type> - TSTRAJECTORYBASIC, TSTRAJECTORYMEMORY, TSTRAJECTORYSINGLEFILE, TSTRAJECTORYVISUALIZATION 620 621 Level: developer 622 623 .seealso: TSTrajectorySetFromOptions(), TSTrajectorySetType() 624 */ 625 static PetscErrorCode TSTrajectorySetTypeFromOptions_Private(PetscOptionItems *PetscOptionsObject,TSTrajectory tj,TS ts) 626 { 627 PetscBool opt; 628 const char *defaultType; 629 char typeName[256]; 630 631 PetscFunctionBegin; 632 if (((PetscObject)tj)->type_name) defaultType = ((PetscObject)tj)->type_name; 633 else defaultType = TSTRAJECTORYBASIC; 634 635 CHKERRQ(TSTrajectoryRegisterAll()); 636 CHKERRQ(PetscOptionsFList("-ts_trajectory_type","TSTrajectory method","TSTrajectorySetType",TSTrajectoryList,defaultType,typeName,256,&opt)); 637 if (opt) { 638 CHKERRQ(TSTrajectorySetType(tj,ts,typeName)); 639 } else { 640 CHKERRQ(TSTrajectorySetType(tj,ts,defaultType)); 641 } 642 PetscFunctionReturn(0); 643 } 644 645 /*@ 646 TSTrajectorySetUseHistory - Use TSHistory in TSTrajectory 647 648 Collective on TSTrajectory 649 650 Input Parameters: 651 + tj - the TSTrajectory context 652 - flg - PETSC_TRUE to save, PETSC_FALSE to disable 653 654 Options Database Keys: 655 . -ts_trajectory_use_history - have it use TSHistory 656 657 Level: advanced 658 659 .seealso: TSTrajectoryCreate(), TSTrajectoryDestroy(), TSTrajectorySetUp() 660 @*/ 661 PetscErrorCode TSTrajectorySetUseHistory(TSTrajectory tj,PetscBool flg) 662 { 663 PetscFunctionBegin; 664 PetscValidHeaderSpecific(tj,TSTRAJECTORY_CLASSID,1); 665 PetscValidLogicalCollectiveBool(tj,flg,2); 666 tj->usehistory = flg; 667 PetscFunctionReturn(0); 668 } 669 670 /*@ 671 TSTrajectorySetMonitor - Monitor the schedules generated by the checkpointing controller 672 673 Collective on TSTrajectory 674 675 Input Parameters: 676 + tj - the TSTrajectory context 677 - flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable 678 679 Options Database Keys: 680 . -ts_trajectory_monitor - print TSTrajectory information 681 682 Level: developer 683 684 .seealso: TSTrajectoryCreate(), TSTrajectoryDestroy(), TSTrajectorySetUp() 685 @*/ 686 PetscErrorCode TSTrajectorySetMonitor(TSTrajectory tj,PetscBool flg) 687 { 688 PetscFunctionBegin; 689 PetscValidHeaderSpecific(tj,TSTRAJECTORY_CLASSID,1); 690 PetscValidLogicalCollectiveBool(tj,flg,2); 691 if (flg) tj->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)tj)); 692 else tj->monitor = NULL; 693 PetscFunctionReturn(0); 694 } 695 696 /*@ 697 TSTrajectorySetKeepFiles - Keep the files generated by the TSTrajectory 698 699 Collective on TSTrajectory 700 701 Input Parameters: 702 + tj - the TSTrajectory context 703 - flg - PETSC_TRUE to save, PETSC_FALSE to disable 704 705 Options Database Keys: 706 . -ts_trajectory_keep_files - have it keep the files 707 708 Notes: 709 By default the TSTrajectory used for adjoint computations, TSTRAJECTORYBASIC, removes the files it generates at the end of the run. This causes the files to be kept. 710 711 Level: advanced 712 713 .seealso: TSTrajectoryCreate(), TSTrajectoryDestroy(), TSTrajectorySetUp(), TSTrajectorySetMonitor() 714 @*/ 715 PetscErrorCode TSTrajectorySetKeepFiles(TSTrajectory tj,PetscBool flg) 716 { 717 PetscFunctionBegin; 718 PetscValidHeaderSpecific(tj,TSTRAJECTORY_CLASSID,1); 719 PetscValidLogicalCollectiveBool(tj,flg,2); 720 tj->keepfiles = flg; 721 PetscFunctionReturn(0); 722 } 723 724 /*@C 725 TSTrajectorySetDirname - Specify the name of the directory where disk checkpoints are stored. 726 727 Collective on TSTrajectory 728 729 Input Parameters: 730 + tj - the TSTrajectory context 731 - dirname - the directory name 732 733 Options Database Keys: 734 . -ts_trajectory_dirname - set the directory name 735 736 Notes: 737 The final location of the files is determined by dirname/filetemplate where filetemplate was provided by TSTrajectorySetFiletemplate() 738 739 Level: developer 740 741 .seealso: TSTrajectorySetFiletemplate(),TSTrajectorySetUp() 742 @*/ 743 PetscErrorCode TSTrajectorySetDirname(TSTrajectory tj,const char dirname[]) 744 { 745 PetscBool flg; 746 747 PetscFunctionBegin; 748 PetscValidHeaderSpecific(tj,TSTRAJECTORY_CLASSID,1); 749 CHKERRQ(PetscStrcmp(tj->dirname,dirname,&flg)); 750 if (!flg && tj->dirfiletemplate) { 751 SETERRQ(PetscObjectComm((PetscObject)tj),PETSC_ERR_ARG_WRONGSTATE,"Cannot set directoryname after TSTrajectory has been setup"); 752 } 753 CHKERRQ(PetscFree(tj->dirname)); 754 CHKERRQ(PetscStrallocpy(dirname,&tj->dirname)); 755 PetscFunctionReturn(0); 756 } 757 758 /*@C 759 TSTrajectorySetFiletemplate - Specify the name template for the files storing checkpoints. 760 761 Collective on TSTrajectory 762 763 Input Parameters: 764 + tj - the TSTrajectory context 765 - filetemplate - the template 766 767 Options Database Keys: 768 . -ts_trajectory_file_template - set the file name template 769 770 Notes: 771 The name template should be of the form, for example filename-%06D.bin It should not begin with a leading / 772 773 The final location of the files is determined by dirname/filetemplate where dirname was provided by TSTrajectorySetDirname(). The %06D is replaced by the 774 timestep counter 775 776 Level: developer 777 778 .seealso: TSTrajectorySetDirname(),TSTrajectorySetUp() 779 @*/ 780 PetscErrorCode TSTrajectorySetFiletemplate(TSTrajectory tj,const char filetemplate[]) 781 { 782 const char *ptr,*ptr2; 783 784 PetscFunctionBegin; 785 PetscValidHeaderSpecific(tj,TSTRAJECTORY_CLASSID,1); 786 PetscCheck(!tj->dirfiletemplate,PetscObjectComm((PetscObject)tj),PETSC_ERR_ARG_WRONGSTATE,"Cannot set filetemplate after TSTrajectory has been setup"); 787 788 PetscCheck(filetemplate[0],PetscObjectComm((PetscObject)tj),PETSC_ERR_USER,"-ts_trajectory_file_template requires a file name template, e.g. filename-%%06D.bin"); 789 /* Do some cursory validation of the input. */ 790 CHKERRQ(PetscStrstr(filetemplate,"%",(char**)&ptr)); 791 PetscCheck(ptr,PetscObjectComm((PetscObject)tj),PETSC_ERR_USER,"-ts_trajectory_file_template requires a file name template, e.g. filename-%%06D.bin"); 792 for (ptr++; ptr && *ptr; ptr++) { 793 CHKERRQ(PetscStrchr("DdiouxX",*ptr,(char**)&ptr2)); 794 PetscCheck(ptr2 || (*ptr >= '0' && *ptr <= '9'),PetscObjectComm((PetscObject)tj),PETSC_ERR_USER,"Invalid file template argument to -ts_trajectory_file_template, should look like filename-%%06D.bin"); 795 if (ptr2) break; 796 } 797 CHKERRQ(PetscFree(tj->filetemplate)); 798 CHKERRQ(PetscStrallocpy(filetemplate,&tj->filetemplate)); 799 PetscFunctionReturn(0); 800 } 801 802 /*@ 803 TSTrajectorySetFromOptions - Sets various TSTrajectory parameters from user options. 804 805 Collective on TSTrajectory 806 807 Input Parameters: 808 + tj - the TSTrajectory context obtained from TSTrajectoryCreate() 809 - ts - the TS context 810 811 Options Database Keys: 812 + -ts_trajectory_type <type> - TSTRAJECTORYBASIC, TSTRAJECTORYMEMORY, TSTRAJECTORYSINGLEFILE, TSTRAJECTORYVISUALIZATION 813 . -ts_trajectory_keep_files <true,false> - keep the files generated by the code after the program ends. This is true by default for TSTRAJECTORYSINGLEFILE, TSTRAJECTORYVISUALIZATION 814 - -ts_trajectory_monitor - print TSTrajectory information 815 816 Level: developer 817 818 Notes: 819 This is not normally called directly by users 820 821 .seealso: TSSetSaveTrajectory(), TSTrajectorySetUp() 822 @*/ 823 PetscErrorCode TSTrajectorySetFromOptions(TSTrajectory tj,TS ts) 824 { 825 PetscBool set,flg; 826 char dirname[PETSC_MAX_PATH_LEN],filetemplate[PETSC_MAX_PATH_LEN]; 827 PetscErrorCode ierr; 828 829 PetscFunctionBegin; 830 PetscValidHeaderSpecific(tj,TSTRAJECTORY_CLASSID,1); 831 if (ts) PetscValidHeaderSpecific(ts,TS_CLASSID,2); 832 ierr = PetscObjectOptionsBegin((PetscObject)tj);CHKERRQ(ierr); 833 CHKERRQ(TSTrajectorySetTypeFromOptions_Private(PetscOptionsObject,tj,ts)); 834 CHKERRQ(PetscOptionsBool("-ts_trajectory_use_history","Turn on/off usage of TSHistory",NULL,tj->usehistory,&tj->usehistory,NULL)); 835 CHKERRQ(PetscOptionsBool("-ts_trajectory_monitor","Print checkpointing schedules","TSTrajectorySetMonitor",tj->monitor ? PETSC_TRUE:PETSC_FALSE,&flg,&set)); 836 if (set) CHKERRQ(TSTrajectorySetMonitor(tj,flg)); 837 CHKERRQ(PetscOptionsInt("-ts_trajectory_reconstruction_order","Interpolation order for reconstruction",NULL,tj->lag.order,&tj->lag.order,NULL)); 838 CHKERRQ(PetscOptionsBool("-ts_trajectory_reconstruction_caching","Turn on/off caching of TSTrajectoryGetVecs input",NULL,tj->lag.caching,&tj->lag.caching,NULL)); 839 CHKERRQ(PetscOptionsBool("-ts_trajectory_adjointmode","Instruct the trajectory that will be used in a TSAdjointSolve()",NULL,tj->adjoint_solve_mode,&tj->adjoint_solve_mode,NULL)); 840 CHKERRQ(PetscOptionsBool("-ts_trajectory_solution_only","Checkpoint solution only","TSTrajectorySetSolutionOnly",tj->solution_only,&tj->solution_only,NULL)); 841 CHKERRQ(PetscOptionsBool("-ts_trajectory_keep_files","Keep any trajectory files generated during the run","TSTrajectorySetKeepFiles",tj->keepfiles,&flg,&set)); 842 if (set) CHKERRQ(TSTrajectorySetKeepFiles(tj,flg)); 843 844 CHKERRQ(PetscOptionsString("-ts_trajectory_dirname","Directory name for TSTrajectory file","TSTrajectorySetDirname",NULL,dirname,sizeof(dirname)-14,&set)); 845 if (set) { 846 CHKERRQ(TSTrajectorySetDirname(tj,dirname)); 847 } 848 849 CHKERRQ(PetscOptionsString("-ts_trajectory_file_template","Template for TSTrajectory file name, use filename-%06D.bin","TSTrajectorySetFiletemplate",NULL,filetemplate,sizeof(filetemplate),&set)); 850 if (set) { 851 CHKERRQ(TSTrajectorySetFiletemplate(tj,filetemplate)); 852 } 853 854 /* Handle specific TSTrajectory options */ 855 if (tj->ops->setfromoptions) { 856 CHKERRQ((*tj->ops->setfromoptions)(PetscOptionsObject,tj)); 857 } 858 ierr = PetscOptionsEnd();CHKERRQ(ierr); 859 PetscFunctionReturn(0); 860 } 861 862 /*@ 863 TSTrajectorySetUp - Sets up the internal data structures, e.g. stacks, for the later use 864 of a TS trajectory. 865 866 Collective on TS 867 868 Input Parameters: 869 + ts - the TS context obtained from TSCreate() 870 - tj - the TS trajectory context 871 872 Level: developer 873 874 .seealso: TSSetSaveTrajectory(), TSTrajectoryCreate(), TSTrajectoryDestroy() 875 @*/ 876 PetscErrorCode TSTrajectorySetUp(TSTrajectory tj,TS ts) 877 { 878 size_t s1,s2; 879 880 PetscFunctionBegin; 881 if (!tj) PetscFunctionReturn(0); 882 PetscValidHeaderSpecific(tj,TSTRAJECTORY_CLASSID,1); 883 if (ts) PetscValidHeaderSpecific(ts,TS_CLASSID,2); 884 if (tj->setupcalled) PetscFunctionReturn(0); 885 886 CHKERRQ(PetscLogEventBegin(TSTrajectory_SetUp,tj,ts,0,0)); 887 if (!((PetscObject)tj)->type_name) { 888 CHKERRQ(TSTrajectorySetType(tj,ts,TSTRAJECTORYBASIC)); 889 } 890 if (tj->ops->setup) { 891 CHKERRQ((*tj->ops->setup)(tj,ts)); 892 } 893 894 tj->setupcalled = PETSC_TRUE; 895 896 /* Set the counters to zero */ 897 tj->recomps = 0; 898 tj->diskreads = 0; 899 tj->diskwrites = 0; 900 CHKERRQ(PetscStrlen(tj->dirname,&s1)); 901 CHKERRQ(PetscStrlen(tj->filetemplate,&s2)); 902 CHKERRQ(PetscFree(tj->dirfiletemplate)); 903 CHKERRQ(PetscMalloc((s1 + s2 + 10)*sizeof(char),&tj->dirfiletemplate)); 904 CHKERRQ(PetscSNPrintf(tj->dirfiletemplate,s1+s2+10,"%s/%s",tj->dirname,tj->filetemplate)); 905 CHKERRQ(PetscLogEventEnd(TSTrajectory_SetUp,tj,ts,0,0)); 906 PetscFunctionReturn(0); 907 } 908 909 /*@ 910 TSTrajectorySetSolutionOnly - Tells the trajectory to store just the solution, and not any intermediate stage also. 911 912 Collective on TSTrajectory 913 914 Input Parameters: 915 + tj - the TS trajectory context 916 - flg - the boolean flag 917 918 Level: developer 919 920 .seealso: TSSetSaveTrajectory(), TSTrajectoryCreate(), TSTrajectoryDestroy(), TSTrajectoryGetSolutionOnly() 921 @*/ 922 PetscErrorCode TSTrajectorySetSolutionOnly(TSTrajectory tj,PetscBool solution_only) 923 { 924 PetscFunctionBegin; 925 PetscValidHeaderSpecific(tj,TSTRAJECTORY_CLASSID,1); 926 PetscValidLogicalCollectiveBool(tj,solution_only,2); 927 tj->solution_only = solution_only; 928 PetscFunctionReturn(0); 929 } 930 931 /*@ 932 TSTrajectoryGetSolutionOnly - Gets the value set with TSTrajectorySetSolutionOnly. 933 934 Logically collective on TSTrajectory 935 936 Input Parameter: 937 . tj - the TS trajectory context 938 939 Output Parameter: 940 . flg - the boolean flag 941 942 Level: developer 943 944 .seealso: TSSetSaveTrajectory(), TSTrajectoryCreate(), TSTrajectoryDestroy(), TSTrajectorySetSolutionOnly() 945 @*/ 946 PetscErrorCode TSTrajectoryGetSolutionOnly(TSTrajectory tj,PetscBool *solution_only) 947 { 948 PetscFunctionBegin; 949 PetscValidHeaderSpecific(tj,TSTRAJECTORY_CLASSID,1); 950 PetscValidPointer(solution_only,2); 951 *solution_only = tj->solution_only; 952 PetscFunctionReturn(0); 953 } 954 955 /*@ 956 TSTrajectoryGetUpdatedHistoryVecs - Get updated state and time-derivative history vectors. 957 958 Collective on TSTrajectory 959 960 Input Parameters: 961 + tj - the TS trajectory context 962 . ts - the TS solver context 963 - time - the requested time 964 965 Output Parameters: 966 + U - state vector at given time (can be interpolated) 967 - Udot - time-derivative vector at given time (can be interpolated) 968 969 Level: developer 970 971 Notes: The vectors are interpolated if time does not match any time step stored in the TSTrajectory(). Pass NULL to not request a vector. 972 This function differs from TSTrajectoryGetVecs since the vectors obtained cannot be modified, and they need to be returned by 973 calling TSTrajectoryRestoreUpdatedHistoryVecs(). 974 975 .seealso: TSSetSaveTrajectory(), TSTrajectoryCreate(), TSTrajectoryDestroy(), TSTrajectoryRestoreUpdatedHistoryVecs(), TSTrajectoryGetVecs() 976 @*/ 977 PetscErrorCode TSTrajectoryGetUpdatedHistoryVecs(TSTrajectory tj, TS ts, PetscReal time, Vec *U, Vec *Udot) 978 { 979 PetscFunctionBegin; 980 PetscValidHeaderSpecific(tj,TSTRAJECTORY_CLASSID,1); 981 PetscValidHeaderSpecific(ts,TS_CLASSID,2); 982 PetscValidLogicalCollectiveReal(tj,time,3); 983 if (U) PetscValidPointer(U,4); 984 if (Udot) PetscValidPointer(Udot,5); 985 if (U && !tj->U) { 986 DM dm; 987 988 CHKERRQ(TSGetDM(ts,&dm)); 989 CHKERRQ(DMCreateGlobalVector(dm,&tj->U)); 990 } 991 if (Udot && !tj->Udot) { 992 DM dm; 993 994 CHKERRQ(TSGetDM(ts,&dm)); 995 CHKERRQ(DMCreateGlobalVector(dm,&tj->Udot)); 996 } 997 CHKERRQ(TSTrajectoryGetVecs(tj,ts,PETSC_DECIDE,&time,U ? tj->U : NULL,Udot ? tj->Udot : NULL)); 998 if (U) { 999 CHKERRQ(VecLockReadPush(tj->U)); 1000 *U = tj->U; 1001 } 1002 if (Udot) { 1003 CHKERRQ(VecLockReadPush(tj->Udot)); 1004 *Udot = tj->Udot; 1005 } 1006 PetscFunctionReturn(0); 1007 } 1008 1009 /*@ 1010 TSTrajectoryRestoreUpdatedHistoryVecs - Restores updated state and time-derivative history vectors obtained with TSTrajectoryGetUpdatedHistoryVecs(). 1011 1012 Collective on TSTrajectory 1013 1014 Input Parameters: 1015 + tj - the TS trajectory context 1016 . U - state vector at given time (can be interpolated) 1017 - Udot - time-derivative vector at given time (can be interpolated) 1018 1019 Level: developer 1020 1021 .seealso: TSTrajectoryGetUpdatedHistoryVecs() 1022 @*/ 1023 PetscErrorCode TSTrajectoryRestoreUpdatedHistoryVecs(TSTrajectory tj, Vec *U, Vec *Udot) 1024 { 1025 PetscFunctionBegin; 1026 PetscValidHeaderSpecific(tj,TSTRAJECTORY_CLASSID,1); 1027 if (U) PetscValidHeaderSpecific(*U,VEC_CLASSID,2); 1028 if (Udot) PetscValidHeaderSpecific(*Udot,VEC_CLASSID,3); 1029 PetscCheck(!U || *U == tj->U,PetscObjectComm((PetscObject)*U),PETSC_ERR_USER,"U was not obtained from TSTrajectoryGetUpdatedHistoryVecs()"); 1030 PetscCheck(!Udot || *Udot == tj->Udot,PetscObjectComm((PetscObject)*Udot),PETSC_ERR_USER,"Udot was not obtained from TSTrajectoryGetUpdatedHistoryVecs()"); 1031 if (U) { 1032 CHKERRQ(VecLockReadPop(tj->U)); 1033 *U = NULL; 1034 } 1035 if (Udot) { 1036 CHKERRQ(VecLockReadPop(tj->Udot)); 1037 *Udot = NULL; 1038 } 1039 PetscFunctionReturn(0); 1040 } 1041