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