Actual source code: ts.c
petsc-3.12.5 2020-03-29
1: #include <petsc/private/tsimpl.h>
2: #include <petscdmshell.h>
3: #include <petscdmda.h>
4: #include <petscviewer.h>
5: #include <petscdraw.h>
7: #define SkipSmallValue(a,b,tol) if(PetscAbsScalar(a)< tol || PetscAbsScalar(b)< tol) continue;
9: /* Logging support */
10: PetscClassId TS_CLASSID, DMTS_CLASSID;
11: PetscLogEvent TS_Step, TS_PseudoComputeTimeStep, TS_FunctionEval, TS_JacobianEval;
13: const char *const TSExactFinalTimeOptions[] = {"UNSPECIFIED","STEPOVER","INTERPOLATE","MATCHSTEP","TSExactFinalTimeOption","TS_EXACTFINALTIME_",0};
16: /*@C
17: TSMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user
19: Collective on TS
21: Input Parameters:
22: + ts - TS object you wish to monitor
23: . name - the monitor type one is seeking
24: . help - message indicating what monitoring is done
25: . manual - manual page for the monitor
26: . monitor - the monitor function
27: - monitorsetup - a function that is called once ONLY if the user selected this monitor that may set additional features of the TS or PetscViewer objects
29: Level: developer
31: .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
32: PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
33: PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
34: PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
35: PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
36: PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
37: PetscOptionsFList(), PetscOptionsEList()
38: @*/
39: PetscErrorCode TSMonitorSetFromOptions(TS ts,const char name[],const char help[], const char manual[],PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*),PetscErrorCode (*monitorsetup)(TS,PetscViewerAndFormat*))
40: {
41: PetscErrorCode ierr;
42: PetscViewer viewer;
43: PetscViewerFormat format;
44: PetscBool flg;
47: PetscOptionsGetViewer(PetscObjectComm((PetscObject)ts),((PetscObject) ts)->options,((PetscObject)ts)->prefix,name,&viewer,&format,&flg);
48: if (flg) {
49: PetscViewerAndFormat *vf;
50: PetscViewerAndFormatCreate(viewer,format,&vf);
51: PetscObjectDereference((PetscObject)viewer);
52: if (monitorsetup) {
53: (*monitorsetup)(ts,vf);
54: }
55: TSMonitorSet(ts,(PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);
56: }
57: return(0);
58: }
60: static PetscErrorCode TSAdaptSetDefaultType(TSAdapt adapt,TSAdaptType default_type)
61: {
67: if (!((PetscObject)adapt)->type_name) {
68: TSAdaptSetType(adapt,default_type);
69: }
70: return(0);
71: }
73: /*@
74: TSSetFromOptions - Sets various TS parameters from user options.
76: Collective on TS
78: Input Parameter:
79: . ts - the TS context obtained from TSCreate()
81: Options Database Keys:
82: + -ts_type <type> - TSEULER, TSBEULER, TSSUNDIALS, TSPSEUDO, TSCN, TSRK, TSTHETA, TSALPHA, TSGLLE, TSSSP, TSGLEE, TSBSYMP
83: . -ts_save_trajectory - checkpoint the solution at each time-step
84: . -ts_max_time <time> - maximum time to compute to
85: . -ts_max_steps <steps> - maximum number of time-steps to take
86: . -ts_init_time <time> - initial time to start computation
87: . -ts_final_time <time> - final time to compute to (deprecated: use -ts_max_time)
88: . -ts_dt <dt> - initial time step
89: . -ts_exact_final_time <stepover,interpolate,matchstep> - whether to stop at the exact given final time and how to compute the solution at that ti,e
90: . -ts_max_snes_failures <maxfailures> - Maximum number of nonlinear solve failures allowed
91: . -ts_max_reject <maxrejects> - Maximum number of step rejections before step fails
92: . -ts_error_if_step_fails <true,false> - Error if no step succeeds
93: . -ts_rtol <rtol> - relative tolerance for local truncation error
94: . -ts_atol <atol> Absolute tolerance for local truncation error
95: . -ts_rhs_jacobian_test_mult -mat_shell_test_mult_view - test the Jacobian at each iteration against finite difference with RHS function
96: . -ts_rhs_jacobian_test_mult_transpose -mat_shell_test_mult_transpose_view - test the Jacobian at each iteration against finite difference with RHS function
97: . -ts_adjoint_solve <yes,no> After solving the ODE/DAE solve the adjoint problem (requires -ts_save_trajectory)
98: . -ts_fd_color - Use finite differences with coloring to compute IJacobian
99: . -ts_monitor - print information at each timestep
100: . -ts_monitor_lg_solution - Monitor solution graphically
101: . -ts_monitor_lg_error - Monitor error graphically
102: . -ts_monitor_error - Monitors norm of error
103: . -ts_monitor_lg_timestep - Monitor timestep size graphically
104: . -ts_monitor_lg_timestep_log - Monitor log timestep size graphically
105: . -ts_monitor_lg_snes_iterations - Monitor number nonlinear iterations for each timestep graphically
106: . -ts_monitor_lg_ksp_iterations - Monitor number nonlinear iterations for each timestep graphically
107: . -ts_monitor_sp_eig - Monitor eigenvalues of linearized operator graphically
108: . -ts_monitor_draw_solution - Monitor solution graphically
109: . -ts_monitor_draw_solution_phase <xleft,yleft,xright,yright> - Monitor solution graphically with phase diagram, requires problem with exactly 2 degrees of freedom
110: . -ts_monitor_draw_error - Monitor error graphically, requires use to have provided TSSetSolutionFunction()
111: . -ts_monitor_solution [ascii binary draw][:filename][:viewerformat] - monitors the solution at each timestep
112: . -ts_monitor_solution_vtk <filename.vts,filename.vtu> - Save each time step to a binary file, use filename-%%03D.vts (filename-%%03D.vtu)
113: - -ts_monitor_envelope - determine maximum and minimum value of each component of the solution over the solution time
115: Developer Note:
116: We should unify all the -ts_monitor options in the way that -xxx_view has been unified
118: Level: beginner
120: .seealso: TSGetType()
121: @*/
122: PetscErrorCode TSSetFromOptions(TS ts)
123: {
124: PetscBool opt,flg,tflg;
125: PetscErrorCode ierr;
126: char monfilename[PETSC_MAX_PATH_LEN];
127: PetscReal time_step;
128: TSExactFinalTimeOption eftopt;
129: char dir[16];
130: TSIFunction ifun;
131: const char *defaultType;
132: char typeName[256];
137: TSRegisterAll();
138: TSGetIFunction(ts,NULL,&ifun,NULL);
140: PetscObjectOptionsBegin((PetscObject)ts);
141: if (((PetscObject)ts)->type_name) defaultType = ((PetscObject)ts)->type_name;
142: else defaultType = ifun ? TSBEULER : TSEULER;
143: PetscOptionsFList("-ts_type","TS method","TSSetType",TSList,defaultType,typeName,256,&opt);
144: if (opt) {
145: TSSetType(ts,typeName);
146: } else {
147: TSSetType(ts,defaultType);
148: }
150: /* Handle generic TS options */
151: PetscOptionsDeprecated("-ts_final_time","-ts_max_time","3.10",NULL);
152: PetscOptionsReal("-ts_max_time","Maximum time to run to","TSSetMaxTime",ts->max_time,&ts->max_time,NULL);
153: PetscOptionsInt("-ts_max_steps","Maximum number of time steps","TSSetMaxSteps",ts->max_steps,&ts->max_steps,NULL);
154: PetscOptionsReal("-ts_init_time","Initial time","TSSetTime",ts->ptime,&ts->ptime,NULL);
155: PetscOptionsReal("-ts_dt","Initial time step","TSSetTimeStep",ts->time_step,&time_step,&flg);
156: if (flg) {TSSetTimeStep(ts,time_step);}
157: PetscOptionsEnum("-ts_exact_final_time","Option for handling of final time step","TSSetExactFinalTime",TSExactFinalTimeOptions,(PetscEnum)ts->exact_final_time,(PetscEnum*)&eftopt,&flg);
158: if (flg) {TSSetExactFinalTime(ts,eftopt);}
159: PetscOptionsInt("-ts_max_snes_failures","Maximum number of nonlinear solve failures","TSSetMaxSNESFailures",ts->max_snes_failures,&ts->max_snes_failures,NULL);
160: PetscOptionsInt("-ts_max_reject","Maximum number of step rejections before step fails","TSSetMaxStepRejections",ts->max_reject,&ts->max_reject,NULL);
161: PetscOptionsBool("-ts_error_if_step_fails","Error if no step succeeds","TSSetErrorIfStepFails",ts->errorifstepfailed,&ts->errorifstepfailed,NULL);
162: PetscOptionsReal("-ts_rtol","Relative tolerance for local truncation error","TSSetTolerances",ts->rtol,&ts->rtol,NULL);
163: PetscOptionsReal("-ts_atol","Absolute tolerance for local truncation error","TSSetTolerances",ts->atol,&ts->atol,NULL);
165: PetscOptionsBool("-ts_rhs_jacobian_test_mult","Test the RHS Jacobian for consistency with RHS at each solve ","None",ts->testjacobian,&ts->testjacobian,NULL);
166: PetscOptionsBool("-ts_rhs_jacobian_test_mult_transpose","Test the RHS Jacobian transpose for consistency with RHS at each solve ","None",ts->testjacobiantranspose,&ts->testjacobiantranspose,NULL);
167: PetscOptionsBool("-ts_use_splitrhsfunction","Use the split RHS function for multirate solvers ","TSSetUseSplitRHSFunction",ts->use_splitrhsfunction,&ts->use_splitrhsfunction,NULL);
168: #if defined(PETSC_HAVE_SAWS)
169: {
170: PetscBool set;
171: flg = PETSC_FALSE;
172: PetscOptionsBool("-ts_saws_block","Block for SAWs memory snooper at end of TSSolve","PetscObjectSAWsBlock",((PetscObject)ts)->amspublishblock,&flg,&set);
173: if (set) {
174: PetscObjectSAWsSetBlock((PetscObject)ts,flg);
175: }
176: }
177: #endif
179: /* Monitor options */
180: TSMonitorSetFromOptions(ts,"-ts_monitor","Monitor time and timestep size","TSMonitorDefault",TSMonitorDefault,NULL);
181: TSMonitorSetFromOptions(ts,"-ts_monitor_extreme","Monitor extreme values of the solution","TSMonitorExtreme",TSMonitorExtreme,NULL);
182: TSMonitorSetFromOptions(ts,"-ts_monitor_solution","View the solution at each timestep","TSMonitorSolution",TSMonitorSolution,NULL);
184: PetscOptionsString("-ts_monitor_python","Use Python function","TSMonitorSet",0,monfilename,PETSC_MAX_PATH_LEN,&flg);
185: if (flg) {PetscPythonMonitorSet((PetscObject)ts,monfilename);}
187: PetscOptionsName("-ts_monitor_lg_solution","Monitor solution graphically","TSMonitorLGSolution",&opt);
188: if (opt) {
189: TSMonitorLGCtx ctx;
190: PetscInt howoften = 1;
192: PetscOptionsInt("-ts_monitor_lg_solution","Monitor solution graphically","TSMonitorLGSolution",howoften,&howoften,NULL);
193: TSMonitorLGCtxCreate(PETSC_COMM_SELF,0,0,PETSC_DECIDE,PETSC_DECIDE,400,300,howoften,&ctx);
194: TSMonitorSet(ts,TSMonitorLGSolution,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);
195: }
197: PetscOptionsName("-ts_monitor_lg_error","Monitor error graphically","TSMonitorLGError",&opt);
198: if (opt) {
199: TSMonitorLGCtx ctx;
200: PetscInt howoften = 1;
202: PetscOptionsInt("-ts_monitor_lg_error","Monitor error graphically","TSMonitorLGError",howoften,&howoften,NULL);
203: TSMonitorLGCtxCreate(PETSC_COMM_SELF,0,0,PETSC_DECIDE,PETSC_DECIDE,400,300,howoften,&ctx);
204: TSMonitorSet(ts,TSMonitorLGError,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);
205: }
206: TSMonitorSetFromOptions(ts,"-ts_monitor_error","View the error at each timestep","TSMonitorError",TSMonitorError,NULL);
208: PetscOptionsName("-ts_monitor_lg_timestep","Monitor timestep size graphically","TSMonitorLGTimeStep",&opt);
209: if (opt) {
210: TSMonitorLGCtx ctx;
211: PetscInt howoften = 1;
213: PetscOptionsInt("-ts_monitor_lg_timestep","Monitor timestep size graphically","TSMonitorLGTimeStep",howoften,&howoften,NULL);
214: TSMonitorLGCtxCreate(PetscObjectComm((PetscObject)ts),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,howoften,&ctx);
215: TSMonitorSet(ts,TSMonitorLGTimeStep,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);
216: }
217: PetscOptionsName("-ts_monitor_lg_timestep_log","Monitor log timestep size graphically","TSMonitorLGTimeStep",&opt);
218: if (opt) {
219: TSMonitorLGCtx ctx;
220: PetscInt howoften = 1;
222: PetscOptionsInt("-ts_monitor_lg_timestep_log","Monitor log timestep size graphically","TSMonitorLGTimeStep",howoften,&howoften,NULL);
223: TSMonitorLGCtxCreate(PetscObjectComm((PetscObject)ts),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,howoften,&ctx);
224: TSMonitorSet(ts,TSMonitorLGTimeStep,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);
225: ctx->semilogy = PETSC_TRUE;
226: }
228: PetscOptionsName("-ts_monitor_lg_snes_iterations","Monitor number nonlinear iterations for each timestep graphically","TSMonitorLGSNESIterations",&opt);
229: if (opt) {
230: TSMonitorLGCtx ctx;
231: PetscInt howoften = 1;
233: PetscOptionsInt("-ts_monitor_lg_snes_iterations","Monitor number nonlinear iterations for each timestep graphically","TSMonitorLGSNESIterations",howoften,&howoften,NULL);
234: TSMonitorLGCtxCreate(PetscObjectComm((PetscObject)ts),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,howoften,&ctx);
235: TSMonitorSet(ts,TSMonitorLGSNESIterations,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);
236: }
237: PetscOptionsName("-ts_monitor_lg_ksp_iterations","Monitor number nonlinear iterations for each timestep graphically","TSMonitorLGKSPIterations",&opt);
238: if (opt) {
239: TSMonitorLGCtx ctx;
240: PetscInt howoften = 1;
242: PetscOptionsInt("-ts_monitor_lg_ksp_iterations","Monitor number nonlinear iterations for each timestep graphically","TSMonitorLGKSPIterations",howoften,&howoften,NULL);
243: TSMonitorLGCtxCreate(PetscObjectComm((PetscObject)ts),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,howoften,&ctx);
244: TSMonitorSet(ts,TSMonitorLGKSPIterations,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);
245: }
246: PetscOptionsName("-ts_monitor_sp_eig","Monitor eigenvalues of linearized operator graphically","TSMonitorSPEig",&opt);
247: if (opt) {
248: TSMonitorSPEigCtx ctx;
249: PetscInt howoften = 1;
251: PetscOptionsInt("-ts_monitor_sp_eig","Monitor eigenvalues of linearized operator graphically","TSMonitorSPEig",howoften,&howoften,NULL);
252: TSMonitorSPEigCtxCreate(PETSC_COMM_SELF,0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx);
253: TSMonitorSet(ts,TSMonitorSPEig,ctx,(PetscErrorCode (*)(void**))TSMonitorSPEigCtxDestroy);
254: }
255: PetscOptionsName("-ts_monitor_sp_swarm","Display particle phase from the DMSwarm","TSMonitorSPSwarm",&opt);
256: if (opt) {
257: TSMonitorSPCtx ctx;
258: PetscInt howoften = 1;
259: PetscOptionsInt("-ts_monitor_sp_swarm","Display particles phase from the DMSwarm","TSMonitorSPSwarm",howoften,&howoften,NULL);
260: TSMonitorSPCtxCreate(PETSC_COMM_SELF, NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 300, 300, howoften, &ctx);
261: TSMonitorSet(ts, TSMonitorSPSwarmSolution, ctx, (PetscErrorCode (*)(void**))TSMonitorSPCtxDestroy);
262: }
263: opt = PETSC_FALSE;
264: PetscOptionsName("-ts_monitor_draw_solution","Monitor solution graphically","TSMonitorDrawSolution",&opt);
265: if (opt) {
266: TSMonitorDrawCtx ctx;
267: PetscInt howoften = 1;
269: PetscOptionsInt("-ts_monitor_draw_solution","Monitor solution graphically","TSMonitorDrawSolution",howoften,&howoften,NULL);
270: TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),0,"Computed Solution",PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx);
271: TSMonitorSet(ts,TSMonitorDrawSolution,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);
272: }
273: opt = PETSC_FALSE;
274: PetscOptionsName("-ts_monitor_draw_solution_phase","Monitor solution graphically","TSMonitorDrawSolutionPhase",&opt);
275: if (opt) {
276: TSMonitorDrawCtx ctx;
277: PetscReal bounds[4];
278: PetscInt n = 4;
279: PetscDraw draw;
280: PetscDrawAxis axis;
282: PetscOptionsRealArray("-ts_monitor_draw_solution_phase","Monitor solution graphically","TSMonitorDrawSolutionPhase",bounds,&n,NULL);
283: if (n != 4) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONG,"Must provide bounding box of phase field");
284: TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,1,&ctx);
285: PetscViewerDrawGetDraw(ctx->viewer,0,&draw);
286: PetscViewerDrawGetDrawAxis(ctx->viewer,0,&axis);
287: PetscDrawAxisSetLimits(axis,bounds[0],bounds[2],bounds[1],bounds[3]);
288: PetscDrawAxisSetLabels(axis,"Phase Diagram","Variable 1","Variable 2");
289: TSMonitorSet(ts,TSMonitorDrawSolutionPhase,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);
290: }
291: opt = PETSC_FALSE;
292: PetscOptionsName("-ts_monitor_draw_error","Monitor error graphically","TSMonitorDrawError",&opt);
293: if (opt) {
294: TSMonitorDrawCtx ctx;
295: PetscInt howoften = 1;
297: PetscOptionsInt("-ts_monitor_draw_error","Monitor error graphically","TSMonitorDrawError",howoften,&howoften,NULL);
298: TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),0,"Error",PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx);
299: TSMonitorSet(ts,TSMonitorDrawError,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);
300: }
301: opt = PETSC_FALSE;
302: PetscOptionsName("-ts_monitor_draw_solution_function","Monitor solution provided by TSMonitorSetSolutionFunction() graphically","TSMonitorDrawSolutionFunction",&opt);
303: if (opt) {
304: TSMonitorDrawCtx ctx;
305: PetscInt howoften = 1;
307: PetscOptionsInt("-ts_monitor_draw_solution_function","Monitor solution provided by TSMonitorSetSolutionFunction() graphically","TSMonitorDrawSolutionFunction",howoften,&howoften,NULL);
308: TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),0,"Solution provided by user function",PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx);
309: TSMonitorSet(ts,TSMonitorDrawSolutionFunction,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);
310: }
312: opt = PETSC_FALSE;
313: 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);
314: if (flg) {
315: const char *ptr,*ptr2;
316: char *filetemplate;
317: if (!monfilename[0]) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"-ts_monitor_solution_vtk requires a file template, e.g. filename-%%03D.vts");
318: /* Do some cursory validation of the input. */
319: PetscStrstr(monfilename,"%",(char**)&ptr);
320: if (!ptr) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"-ts_monitor_solution_vtk requires a file template, e.g. filename-%%03D.vts");
321: for (ptr++; ptr && *ptr; ptr++) {
322: PetscStrchr("DdiouxX",*ptr,(char**)&ptr2);
323: if (!ptr2 && (*ptr < '0' || '9' < *ptr)) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Invalid file template argument to -ts_monitor_solution_vtk, should look like filename-%%03D.vts");
324: if (ptr2) break;
325: }
326: PetscStrallocpy(monfilename,&filetemplate);
327: TSMonitorSet(ts,TSMonitorSolutionVTK,filetemplate,(PetscErrorCode (*)(void**))TSMonitorSolutionVTKDestroy);
328: }
330: PetscOptionsString("-ts_monitor_dmda_ray","Display a ray of the solution","None","y=0",dir,16,&flg);
331: if (flg) {
332: TSMonitorDMDARayCtx *rayctx;
333: int ray = 0;
334: DMDADirection ddir;
335: DM da;
336: PetscMPIInt rank;
338: if (dir[1] != '=') SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONG,"Unknown ray %s",dir);
339: if (dir[0] == 'x') ddir = DMDA_X;
340: else if (dir[0] == 'y') ddir = DMDA_Y;
341: else SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONG,"Unknown ray %s",dir);
342: sscanf(dir+2,"%d",&ray);
344: PetscInfo2(((PetscObject)ts),"Displaying DMDA ray %c = %d\n",dir[0],ray);
345: PetscNew(&rayctx);
346: TSGetDM(ts,&da);
347: DMDAGetRay(da,ddir,ray,&rayctx->ray,&rayctx->scatter);
348: MPI_Comm_rank(PetscObjectComm((PetscObject)ts),&rank);
349: if (!rank) {
350: PetscViewerDrawOpen(PETSC_COMM_SELF,0,0,0,0,600,300,&rayctx->viewer);
351: }
352: rayctx->lgctx = NULL;
353: TSMonitorSet(ts,TSMonitorDMDARay,rayctx,TSMonitorDMDARayDestroy);
354: }
355: PetscOptionsString("-ts_monitor_lg_dmda_ray","Display a ray of the solution","None","x=0",dir,16,&flg);
356: if (flg) {
357: TSMonitorDMDARayCtx *rayctx;
358: int ray = 0;
359: DMDADirection ddir;
360: DM da;
361: PetscInt howoften = 1;
363: if (dir[1] != '=') SETERRQ1(PetscObjectComm((PetscObject) ts), PETSC_ERR_ARG_WRONG, "Malformed ray %s", dir);
364: if (dir[0] == 'x') ddir = DMDA_X;
365: else if (dir[0] == 'y') ddir = DMDA_Y;
366: else SETERRQ1(PetscObjectComm((PetscObject) ts), PETSC_ERR_ARG_WRONG, "Unknown ray direction %s", dir);
367: sscanf(dir+2, "%d", &ray);
369: PetscInfo2(((PetscObject) ts),"Displaying LG DMDA ray %c = %d\n", dir[0], ray);
370: PetscNew(&rayctx);
371: TSGetDM(ts, &da);
372: DMDAGetRay(da, ddir, ray, &rayctx->ray, &rayctx->scatter);
373: TSMonitorLGCtxCreate(PETSC_COMM_SELF,0,0,PETSC_DECIDE,PETSC_DECIDE,600,400,howoften,&rayctx->lgctx);
374: TSMonitorSet(ts, TSMonitorLGDMDARay, rayctx, TSMonitorDMDARayDestroy);
375: }
377: PetscOptionsName("-ts_monitor_envelope","Monitor maximum and minimum value of each component of the solution","TSMonitorEnvelope",&opt);
378: if (opt) {
379: TSMonitorEnvelopeCtx ctx;
381: TSMonitorEnvelopeCtxCreate(ts,&ctx);
382: TSMonitorSet(ts,TSMonitorEnvelope,ctx,(PetscErrorCode (*)(void**))TSMonitorEnvelopeCtxDestroy);
383: }
385: flg = PETSC_FALSE;
386: PetscOptionsBool("-ts_fd_color", "Use finite differences with coloring to compute IJacobian", "TSComputeJacobianDefaultColor", flg, &flg, NULL);
387: if (flg) {
388: DM dm;
389: DMTS tdm;
391: TSGetDM(ts, &dm);
392: DMGetDMTS(dm, &tdm);
393: tdm->ijacobianctx = NULL;
394: TSSetIJacobian(ts, NULL, NULL, TSComputeIJacobianDefaultColor, 0);
395: PetscInfo(ts, "Setting default finite difference coloring Jacobian matrix\n");
396: }
398: /* Handle specific TS options */
399: if (ts->ops->setfromoptions) {
400: (*ts->ops->setfromoptions)(PetscOptionsObject,ts);
401: }
403: /* Handle TSAdapt options */
404: TSGetAdapt(ts,&ts->adapt);
405: TSAdaptSetDefaultType(ts->adapt,ts->default_adapt_type);
406: TSAdaptSetFromOptions(PetscOptionsObject,ts->adapt);
408: /* TS trajectory must be set after TS, since it may use some TS options above */
409: tflg = ts->trajectory ? PETSC_TRUE : PETSC_FALSE;
410: PetscOptionsBool("-ts_save_trajectory","Save the solution at each timestep","TSSetSaveTrajectory",tflg,&tflg,NULL);
411: if (tflg) {
412: TSSetSaveTrajectory(ts);
413: }
415: TSAdjointSetFromOptions(PetscOptionsObject,ts);
417: /* process any options handlers added with PetscObjectAddOptionsHandler() */
418: PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)ts);
419: PetscOptionsEnd();
421: if (ts->trajectory) {
422: TSTrajectorySetFromOptions(ts->trajectory,ts);
423: }
425: /* why do we have to do this here and not during TSSetUp? */
426: TSGetSNES(ts,&ts->snes);
427: if (ts->problem_type == TS_LINEAR) {
428: PetscObjectTypeCompareAny((PetscObject)ts->snes,&flg,SNESKSPONLY,SNESKSPTRANSPOSEONLY,"");
429: if (!flg) { SNESSetType(ts->snes,SNESKSPONLY); }
430: }
431: SNESSetFromOptions(ts->snes);
432: return(0);
433: }
435: /*@
436: TSGetTrajectory - Gets the trajectory from a TS if it exists
438: Collective on TS
440: Input Parameters:
441: . ts - the TS context obtained from TSCreate()
443: Output Parameters;
444: . tr - the TSTrajectory object, if it exists
446: Note: This routine should be called after all TS options have been set
448: Level: advanced
450: .seealso: TSGetTrajectory(), TSAdjointSolve(), TSTrajectory, TSTrajectoryCreate()
452: @*/
453: PetscErrorCode TSGetTrajectory(TS ts,TSTrajectory *tr)
454: {
457: *tr = ts->trajectory;
458: return(0);
459: }
461: /*@
462: TSSetSaveTrajectory - Causes the TS to save its solutions as it iterates forward in time in a TSTrajectory object
464: Collective on TS
466: Input Parameters:
467: . ts - the TS context obtained from TSCreate()
469: Options Database:
470: + -ts_save_trajectory - saves the trajectory to a file
471: - -ts_trajectory_type type
473: Note: This routine should be called after all TS options have been set
475: The TSTRAJECTORYVISUALIZATION files can be loaded into Python with $PETSC_DIR/lib/petsc/bin/PetscBinaryIOTrajectory.py and
476: MATLAB with $PETSC_DIR/share/petsc/matlab/PetscReadBinaryTrajectory.m
478: Level: intermediate
480: .seealso: TSGetTrajectory(), TSAdjointSolve()
482: @*/
483: PetscErrorCode TSSetSaveTrajectory(TS ts)
484: {
489: if (!ts->trajectory) {
490: TSTrajectoryCreate(PetscObjectComm((PetscObject)ts),&ts->trajectory);
491: }
492: return(0);
493: }
495: /*@
496: TSResetTrajectory - Destroys and recreates the internal TSTrajectory object
498: Collective on TS
500: Input Parameters:
501: . ts - the TS context obtained from TSCreate()
503: Level: intermediate
505: .seealso: TSGetTrajectory(), TSAdjointSolve()
507: @*/
508: PetscErrorCode TSResetTrajectory(TS ts)
509: {
514: if (ts->trajectory) {
515: TSTrajectoryDestroy(&ts->trajectory);
516: TSTrajectoryCreate(PetscObjectComm((PetscObject)ts),&ts->trajectory);
517: }
518: return(0);
519: }
521: /*@
522: TSComputeRHSJacobian - Computes the Jacobian matrix that has been
523: set with TSSetRHSJacobian().
525: Collective on TS
527: Input Parameters:
528: + ts - the TS context
529: . t - current timestep
530: - U - input vector
532: Output Parameters:
533: + A - Jacobian matrix
534: . B - optional preconditioning matrix
535: - flag - flag indicating matrix structure
537: Notes:
538: Most users should not need to explicitly call this routine, as it
539: is used internally within the nonlinear solvers.
541: See KSPSetOperators() for important information about setting the
542: flag parameter.
544: Level: developer
546: .seealso: TSSetRHSJacobian(), KSPSetOperators()
547: @*/
548: PetscErrorCode TSComputeRHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B)
549: {
550: PetscErrorCode ierr;
551: PetscObjectState Ustate;
552: PetscObjectId Uid;
553: DM dm;
554: DMTS tsdm;
555: TSRHSJacobian rhsjacobianfunc;
556: void *ctx;
557: TSIJacobian ijacobianfunc;
558: TSRHSFunction rhsfunction;
564: TSGetDM(ts,&dm);
565: DMGetDMTS(dm,&tsdm);
566: DMTSGetRHSJacobian(dm,&rhsjacobianfunc,&ctx);
567: DMTSGetIJacobian(dm,&ijacobianfunc,NULL);
568: DMTSGetRHSFunction(dm,&rhsfunction,&ctx);
569: PetscObjectStateGet((PetscObject)U,&Ustate);
570: PetscObjectGetId((PetscObject)U,&Uid);
572: if (ts->rhsjacobian.time == t && (ts->problem_type == TS_LINEAR || (ts->rhsjacobian.Xid == Uid && ts->rhsjacobian.Xstate == Ustate)) && (rhsfunction != TSComputeRHSFunctionLinear)) {
573: /* restore back RHS Jacobian matrices if they have been shifted/scaled */
574: if (A == ts->Arhs) {
575: if (ts->rhsjacobian.shift != 0) {
576: MatShift(A,-ts->rhsjacobian.shift);
577: }
578: if (ts->rhsjacobian.scale != 1.) {
579: MatScale(A,1./ts->rhsjacobian.scale);
580: }
581: }
582: if (B && B == ts->Brhs && A != B) {
583: if (ts->rhsjacobian.shift != 0) {
584: MatShift(B,-ts->rhsjacobian.shift);
585: }
586: if (ts->rhsjacobian.scale != 1.) {
587: MatScale(B,1./ts->rhsjacobian.scale);
588: }
589: }
590: ts->rhsjacobian.shift = 0;
591: ts->rhsjacobian.scale = 1.;
592: return(0);
593: }
595: if (!rhsjacobianfunc && !ijacobianfunc) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must call TSSetRHSJacobian() and / or TSSetIJacobian()");
597: if (ts->rhsjacobian.reuse) {
598: if (A == ts->Arhs) {
599: /* MatScale has a short path for this case.
600: However, this code path is taken the first time TSComputeRHSJacobian is called
601: and the matrices have not assembled yet */
602: if (ts->rhsjacobian.shift != 0) {
603: MatShift(A,-ts->rhsjacobian.shift);
604: }
605: if (ts->rhsjacobian.scale != 1.) {
606: MatScale(A,1./ts->rhsjacobian.scale);
607: }
608: }
609: if (B && B == ts->Brhs && A != B) {
610: if (ts->rhsjacobian.shift != 0) {
611: MatShift(B,-ts->rhsjacobian.shift);
612: }
613: if (ts->rhsjacobian.scale != 1.) {
614: MatScale(B,1./ts->rhsjacobian.scale);
615: }
616: }
617: }
619: if (rhsjacobianfunc) {
620: PetscBool missing;
621: PetscLogEventBegin(TS_JacobianEval,ts,U,A,B);
622: PetscStackPush("TS user Jacobian function");
623: (*rhsjacobianfunc)(ts,t,U,A,B,ctx);
624: PetscStackPop;
625: PetscLogEventEnd(TS_JacobianEval,ts,U,A,B);
626: if (A) {
627: MatMissingDiagonal(A,&missing,NULL);
628: if (missing) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Amat passed to TSSetRHSJacobian() must have all diagonal entries set, if they are zero you must still set them with a zero value");
629: }
630: if (B && B != A) {
631: MatMissingDiagonal(B,&missing,NULL);
632: if (missing) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Bmat passed to TSSetRHSJacobian() must have all diagonal entries set, if they are zero you must still set them with a zero value");
633: }
634: } else {
635: MatZeroEntries(A);
636: if (B && A != B) {MatZeroEntries(B);}
637: }
638: ts->rhsjacobian.time = t;
639: ts->rhsjacobian.shift = 0;
640: ts->rhsjacobian.scale = 1.;
641: PetscObjectGetId((PetscObject)U,&ts->rhsjacobian.Xid);
642: PetscObjectStateGet((PetscObject)U,&ts->rhsjacobian.Xstate);
643: return(0);
644: }
646: /*@
647: TSComputeRHSFunction - Evaluates the right-hand-side function.
649: Collective on TS
651: Input Parameters:
652: + ts - the TS context
653: . t - current time
654: - U - state vector
656: Output Parameter:
657: . y - right hand side
659: Note:
660: Most users should not need to explicitly call this routine, as it
661: is used internally within the nonlinear solvers.
663: Level: developer
665: .seealso: TSSetRHSFunction(), TSComputeIFunction()
666: @*/
667: PetscErrorCode TSComputeRHSFunction(TS ts,PetscReal t,Vec U,Vec y)
668: {
670: TSRHSFunction rhsfunction;
671: TSIFunction ifunction;
672: void *ctx;
673: DM dm;
679: TSGetDM(ts,&dm);
680: DMTSGetRHSFunction(dm,&rhsfunction,&ctx);
681: DMTSGetIFunction(dm,&ifunction,NULL);
683: if (!rhsfunction && !ifunction) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must call TSSetRHSFunction() and / or TSSetIFunction()");
685: PetscLogEventBegin(TS_FunctionEval,ts,U,y,0);
686: if (rhsfunction) {
687: PetscStackPush("TS user right-hand-side function");
688: (*rhsfunction)(ts,t,U,y,ctx);
689: PetscStackPop;
690: } else {
691: VecZeroEntries(y);
692: }
694: PetscLogEventEnd(TS_FunctionEval,ts,U,y,0);
695: return(0);
696: }
698: /*@
699: TSComputeSolutionFunction - Evaluates the solution function.
701: Collective on TS
703: Input Parameters:
704: + ts - the TS context
705: - t - current time
707: Output Parameter:
708: . U - the solution
710: Note:
711: Most users should not need to explicitly call this routine, as it
712: is used internally within the nonlinear solvers.
714: Level: developer
716: .seealso: TSSetSolutionFunction(), TSSetRHSFunction(), TSComputeIFunction()
717: @*/
718: PetscErrorCode TSComputeSolutionFunction(TS ts,PetscReal t,Vec U)
719: {
720: PetscErrorCode ierr;
721: TSSolutionFunction solutionfunction;
722: void *ctx;
723: DM dm;
728: TSGetDM(ts,&dm);
729: DMTSGetSolutionFunction(dm,&solutionfunction,&ctx);
731: if (solutionfunction) {
732: PetscStackPush("TS user solution function");
733: (*solutionfunction)(ts,t,U,ctx);
734: PetscStackPop;
735: }
736: return(0);
737: }
738: /*@
739: TSComputeForcingFunction - Evaluates the forcing function.
741: Collective on TS
743: Input Parameters:
744: + ts - the TS context
745: - t - current time
747: Output Parameter:
748: . U - the function value
750: Note:
751: Most users should not need to explicitly call this routine, as it
752: is used internally within the nonlinear solvers.
754: Level: developer
756: .seealso: TSSetSolutionFunction(), TSSetRHSFunction(), TSComputeIFunction()
757: @*/
758: PetscErrorCode TSComputeForcingFunction(TS ts,PetscReal t,Vec U)
759: {
760: PetscErrorCode ierr, (*forcing)(TS,PetscReal,Vec,void*);
761: void *ctx;
762: DM dm;
767: TSGetDM(ts,&dm);
768: DMTSGetForcingFunction(dm,&forcing,&ctx);
770: if (forcing) {
771: PetscStackPush("TS user forcing function");
772: (*forcing)(ts,t,U,ctx);
773: PetscStackPop;
774: }
775: return(0);
776: }
778: static PetscErrorCode TSGetRHSVec_Private(TS ts,Vec *Frhs)
779: {
780: Vec F;
784: *Frhs = NULL;
785: TSGetIFunction(ts,&F,NULL,NULL);
786: if (!ts->Frhs) {
787: VecDuplicate(F,&ts->Frhs);
788: }
789: *Frhs = ts->Frhs;
790: return(0);
791: }
793: PetscErrorCode TSGetRHSMats_Private(TS ts,Mat *Arhs,Mat *Brhs)
794: {
795: Mat A,B;
797: TSIJacobian ijacobian;
800: if (Arhs) *Arhs = NULL;
801: if (Brhs) *Brhs = NULL;
802: TSGetIJacobian(ts,&A,&B,&ijacobian,NULL);
803: if (Arhs) {
804: if (!ts->Arhs) {
805: if (ijacobian) {
806: MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&ts->Arhs);
807: } else {
808: ts->Arhs = A;
809: PetscObjectReference((PetscObject)A);
810: }
811: } else {
812: PetscBool flg;
813: SNESGetUseMatrixFree(ts->snes,NULL,&flg);
814: /* Handle case where user provided only RHSJacobian and used -snes_mf_operator */
815: if (flg && !ijacobian && ts->Arhs == ts->Brhs){
816: PetscObjectDereference((PetscObject)ts->Arhs);
817: ts->Arhs = A;
818: PetscObjectReference((PetscObject)A);
819: }
820: }
821: *Arhs = ts->Arhs;
822: }
823: if (Brhs) {
824: if (!ts->Brhs) {
825: if (A != B) {
826: if (ijacobian) {
827: MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&ts->Brhs);
828: } else {
829: ts->Brhs = B;
830: PetscObjectReference((PetscObject)B);
831: }
832: } else {
833: PetscObjectReference((PetscObject)ts->Arhs);
834: ts->Brhs = ts->Arhs;
835: }
836: }
837: *Brhs = ts->Brhs;
838: }
839: return(0);
840: }
842: /*@
843: TSComputeIFunction - Evaluates the DAE residual written in implicit form F(t,U,Udot)=0
845: Collective on TS
847: Input Parameters:
848: + ts - the TS context
849: . t - current time
850: . U - state vector
851: . Udot - time derivative of state vector
852: - imex - flag indicates if the method is IMEX so that the RHSFunction should be kept separate
854: Output Parameter:
855: . Y - right hand side
857: Note:
858: Most users should not need to explicitly call this routine, as it
859: is used internally within the nonlinear solvers.
861: If the user did did not write their equations in implicit form, this
862: function recasts them in implicit form.
864: Level: developer
866: .seealso: TSSetIFunction(), TSComputeRHSFunction()
867: @*/
868: PetscErrorCode TSComputeIFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec Y,PetscBool imex)
869: {
871: TSIFunction ifunction;
872: TSRHSFunction rhsfunction;
873: void *ctx;
874: DM dm;
882: TSGetDM(ts,&dm);
883: DMTSGetIFunction(dm,&ifunction,&ctx);
884: DMTSGetRHSFunction(dm,&rhsfunction,NULL);
886: if (!rhsfunction && !ifunction) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must call TSSetRHSFunction() and / or TSSetIFunction()");
888: PetscLogEventBegin(TS_FunctionEval,ts,U,Udot,Y);
889: if (ifunction) {
890: PetscStackPush("TS user implicit function");
891: (*ifunction)(ts,t,U,Udot,Y,ctx);
892: PetscStackPop;
893: }
894: if (imex) {
895: if (!ifunction) {
896: VecCopy(Udot,Y);
897: }
898: } else if (rhsfunction) {
899: if (ifunction) {
900: Vec Frhs;
901: TSGetRHSVec_Private(ts,&Frhs);
902: TSComputeRHSFunction(ts,t,U,Frhs);
903: VecAXPY(Y,-1,Frhs);
904: } else {
905: TSComputeRHSFunction(ts,t,U,Y);
906: VecAYPX(Y,-1,Udot);
907: }
908: }
909: PetscLogEventEnd(TS_FunctionEval,ts,U,Udot,Y);
910: return(0);
911: }
913: /*@
914: TSComputeIJacobian - Evaluates the Jacobian of the DAE
916: Collective on TS
918: Input
919: Input Parameters:
920: + ts - the TS context
921: . t - current timestep
922: . U - state vector
923: . Udot - time derivative of state vector
924: . shift - shift to apply, see note below
925: - imex - flag indicates if the method is IMEX so that the RHSJacobian should be kept separate
927: Output Parameters:
928: + A - Jacobian matrix
929: - B - matrix from which the preconditioner is constructed; often the same as A
931: Notes:
932: If F(t,U,Udot)=0 is the DAE, the required Jacobian is
934: dF/dU + shift*dF/dUdot
936: Most users should not need to explicitly call this routine, as it
937: is used internally within the nonlinear solvers.
939: Level: developer
941: .seealso: TSSetIJacobian()
942: @*/
943: PetscErrorCode TSComputeIJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal shift,Mat A,Mat B,PetscBool imex)
944: {
946: TSIJacobian ijacobian;
947: TSRHSJacobian rhsjacobian;
948: DM dm;
949: void *ctx;
960: TSGetDM(ts,&dm);
961: DMTSGetIJacobian(dm,&ijacobian,&ctx);
962: DMTSGetRHSJacobian(dm,&rhsjacobian,NULL);
964: if (!rhsjacobian && !ijacobian) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must call TSSetRHSJacobian() and / or TSSetIJacobian()");
966: PetscLogEventBegin(TS_JacobianEval,ts,U,A,B);
967: if (ijacobian) {
968: PetscBool missing;
969: PetscStackPush("TS user implicit Jacobian");
970: (*ijacobian)(ts,t,U,Udot,shift,A,B,ctx);
971: PetscStackPop;
972: MatMissingDiagonal(A,&missing,NULL);
973: if (missing) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Amat passed to TSSetIJacobian() must have all diagonal entries set, if they are zero you must still set them with a zero value");
974: if (B != A) {
975: MatMissingDiagonal(B,&missing,NULL);
976: if (missing) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Bmat passed to TSSetIJacobian() must have all diagonal entries set, if they are zero you must still set them with a zero value");
977: }
978: }
979: if (imex) {
980: if (!ijacobian) { /* system was written as Udot = G(t,U) */
981: PetscBool assembled;
982: if (rhsjacobian) {
983: Mat Arhs = NULL;
984: TSGetRHSMats_Private(ts,&Arhs,NULL);
985: if (A == Arhs) {
986: if (rhsjacobian == TSComputeRHSJacobianConstant) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Unsupported operation! cannot use TSComputeRHSJacobianConstant");
987: ts->rhsjacobian.time = PETSC_MIN_REAL;
988: }
989: }
990: MatZeroEntries(A);
991: MatAssembled(A,&assembled);
992: if (!assembled) {
993: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
994: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
995: }
996: MatShift(A,shift);
997: if (A != B) {
998: MatZeroEntries(B);
999: MatAssembled(B,&assembled);
1000: if (!assembled) {
1001: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
1002: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
1003: }
1004: MatShift(B,shift);
1005: }
1006: }
1007: } else {
1008: Mat Arhs = NULL,Brhs = NULL;
1009: if (rhsjacobian) {
1010: TSGetRHSMats_Private(ts,&Arhs,&Brhs);
1011: TSComputeRHSJacobian(ts,t,U,Arhs,Brhs);
1012: }
1013: if (Arhs == A) { /* No IJacobian, so we only have the RHS matrix */
1014: PetscBool flg;
1015: ts->rhsjacobian.scale = -1;
1016: ts->rhsjacobian.shift = shift;
1017: SNESGetUseMatrixFree(ts->snes,NULL,&flg);
1018: /* since -snes_mf_operator uses the full SNES function it does not need to be shifted or scaled here */
1019: if (!flg) {
1020: MatScale(A,-1);
1021: MatShift(A,shift);
1022: }
1023: if (A != B) {
1024: MatScale(B,-1);
1025: MatShift(B,shift);
1026: }
1027: } else if (Arhs) { /* Both IJacobian and RHSJacobian */
1028: MatStructure axpy = DIFFERENT_NONZERO_PATTERN;
1029: if (!ijacobian) { /* No IJacobian provided, but we have a separate RHS matrix */
1030: MatZeroEntries(A);
1031: MatShift(A,shift);
1032: if (A != B) {
1033: MatZeroEntries(B);
1034: MatShift(B,shift);
1035: }
1036: }
1037: MatAXPY(A,-1,Arhs,axpy);
1038: if (A != B) {
1039: MatAXPY(B,-1,Brhs,axpy);
1040: }
1041: }
1042: }
1043: PetscLogEventEnd(TS_JacobianEval,ts,U,A,B);
1044: return(0);
1045: }
1047: /*@C
1048: TSSetRHSFunction - Sets the routine for evaluating the function,
1049: where U_t = G(t,u).
1051: Logically Collective on TS
1053: Input Parameters:
1054: + ts - the TS context obtained from TSCreate()
1055: . r - vector to put the computed right hand side (or NULL to have it created)
1056: . f - routine for evaluating the right-hand-side function
1057: - ctx - [optional] user-defined context for private data for the
1058: function evaluation routine (may be NULL)
1060: Calling sequence of func:
1061: $ PetscErrorCode func (TS ts,PetscReal t,Vec u,Vec F,void *ctx);
1063: + t - current timestep
1064: . u - input vector
1065: . F - function vector
1066: - ctx - [optional] user-defined function context
1068: Level: beginner
1070: Notes:
1071: You must call this function or TSSetIFunction() to define your ODE. You cannot use this function when solving a DAE.
1073: .seealso: TSSetRHSJacobian(), TSSetIJacobian(), TSSetIFunction()
1074: @*/
1075: PetscErrorCode TSSetRHSFunction(TS ts,Vec r,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx)
1076: {
1078: SNES snes;
1079: Vec ralloc = NULL;
1080: DM dm;
1086: TSGetDM(ts,&dm);
1087: DMTSSetRHSFunction(dm,f,ctx);
1088: TSGetSNES(ts,&snes);
1089: if (!r && !ts->dm && ts->vec_sol) {
1090: VecDuplicate(ts->vec_sol,&ralloc);
1091: r = ralloc;
1092: }
1093: SNESSetFunction(snes,r,SNESTSFormFunction,ts);
1094: VecDestroy(&ralloc);
1095: return(0);
1096: }
1098: /*@C
1099: TSSetSolutionFunction - Provide a function that computes the solution of the ODE or DAE
1101: Logically Collective on TS
1103: Input Parameters:
1104: + ts - the TS context obtained from TSCreate()
1105: . f - routine for evaluating the solution
1106: - ctx - [optional] user-defined context for private data for the
1107: function evaluation routine (may be NULL)
1109: Calling sequence of func:
1110: $ PetscErrorCode func (TS ts,PetscReal t,Vec u,void *ctx);
1112: + t - current timestep
1113: . u - output vector
1114: - ctx - [optional] user-defined function context
1116: Options Database:
1117: + -ts_monitor_lg_error - create a graphical monitor of error history, requires user to have provided TSSetSolutionFunction()
1118: - -ts_monitor_draw_error - Monitor error graphically, requires user to have provided TSSetSolutionFunction()
1120: Notes:
1121: This routine is used for testing accuracy of time integration schemes when you already know the solution.
1122: If analytic solutions are not known for your system, consider using the Method of Manufactured Solutions to
1123: create closed-form solutions with non-physical forcing terms.
1125: For low-dimensional problems solved in serial, such as small discrete systems, TSMonitorLGError() can be used to monitor the error history.
1127: Level: beginner
1129: .seealso: TSSetRHSJacobian(), TSSetIJacobian(), TSComputeSolutionFunction(), TSSetForcingFunction(), TSSetSolution(), TSGetSolution(), TSMonitorLGError(), TSMonitorDrawError()
1130: @*/
1131: PetscErrorCode TSSetSolutionFunction(TS ts,PetscErrorCode (*f)(TS,PetscReal,Vec,void*),void *ctx)
1132: {
1134: DM dm;
1138: TSGetDM(ts,&dm);
1139: DMTSSetSolutionFunction(dm,f,ctx);
1140: return(0);
1141: }
1143: /*@C
1144: TSSetForcingFunction - Provide a function that computes a forcing term for a ODE or PDE
1146: Logically Collective on TS
1148: Input Parameters:
1149: + ts - the TS context obtained from TSCreate()
1150: . func - routine for evaluating the forcing function
1151: - ctx - [optional] user-defined context for private data for the
1152: function evaluation routine (may be NULL)
1154: Calling sequence of func:
1155: $ PetscErrorCode func (TS ts,PetscReal t,Vec f,void *ctx);
1157: + t - current timestep
1158: . f - output vector
1159: - ctx - [optional] user-defined function context
1161: Notes:
1162: This routine is useful for testing accuracy of time integration schemes when using the Method of Manufactured Solutions to
1163: create closed-form solutions with a non-physical forcing term. It allows you to use the Method of Manufactored Solution without directly editing the
1164: definition of the problem you are solving and hence possibly introducing bugs.
1166: This replaces the ODE F(u,u_t,t) = 0 the TS is solving with F(u,u_t,t) - func(t) = 0
1168: This forcing function does not depend on the solution to the equations, it can only depend on spatial location, time, and possibly parameters, the
1169: parameters can be passed in the ctx variable.
1171: For low-dimensional problems solved in serial, such as small discrete systems, TSMonitorLGError() can be used to monitor the error history.
1173: Level: beginner
1175: .seealso: TSSetRHSJacobian(), TSSetIJacobian(), TSComputeSolutionFunction(), TSSetSolutionFunction()
1176: @*/
1177: PetscErrorCode TSSetForcingFunction(TS ts,TSForcingFunction func,void *ctx)
1178: {
1180: DM dm;
1184: TSGetDM(ts,&dm);
1185: DMTSSetForcingFunction(dm,func,ctx);
1186: return(0);
1187: }
1189: /*@C
1190: TSSetRHSJacobian - Sets the function to compute the Jacobian of G,
1191: where U_t = G(U,t), as well as the location to store the matrix.
1193: Logically Collective on TS
1195: Input Parameters:
1196: + ts - the TS context obtained from TSCreate()
1197: . Amat - (approximate) Jacobian matrix
1198: . Pmat - matrix from which preconditioner is to be constructed (usually the same as Amat)
1199: . f - the Jacobian evaluation routine
1200: - ctx - [optional] user-defined context for private data for the
1201: Jacobian evaluation routine (may be NULL)
1203: Calling sequence of f:
1204: $ PetscErrorCode func (TS ts,PetscReal t,Vec u,Mat A,Mat B,void *ctx);
1206: + t - current timestep
1207: . u - input vector
1208: . Amat - (approximate) Jacobian matrix
1209: . Pmat - matrix from which preconditioner is to be constructed (usually the same as Amat)
1210: - ctx - [optional] user-defined context for matrix evaluation routine
1212: Notes:
1213: You must set all the diagonal entries of the matrices, if they are zero you must still set them with a zero value
1215: The TS solver may modify the nonzero structure and the entries of the matrices Amat and Pmat between the calls to f()
1216: You should not assume the values are the same in the next call to f() as you set them in the previous call.
1218: Level: beginner
1220: .seealso: SNESComputeJacobianDefaultColor(), TSSetRHSFunction(), TSRHSJacobianSetReuse(), TSSetIJacobian()
1222: @*/
1223: PetscErrorCode TSSetRHSJacobian(TS ts,Mat Amat,Mat Pmat,TSRHSJacobian f,void *ctx)
1224: {
1226: SNES snes;
1227: DM dm;
1228: TSIJacobian ijacobian;
1237: TSGetDM(ts,&dm);
1238: DMTSSetRHSJacobian(dm,f,ctx);
1239: if (f == TSComputeRHSJacobianConstant) {
1240: /* Handle this case automatically for the user; otherwise user should call themselves. */
1241: TSRHSJacobianSetReuse(ts,PETSC_TRUE);
1242: }
1243: DMTSGetIJacobian(dm,&ijacobian,NULL);
1244: TSGetSNES(ts,&snes);
1245: if (!ijacobian) {
1246: SNESSetJacobian(snes,Amat,Pmat,SNESTSFormJacobian,ts);
1247: }
1248: if (Amat) {
1249: PetscObjectReference((PetscObject)Amat);
1250: MatDestroy(&ts->Arhs);
1251: ts->Arhs = Amat;
1252: }
1253: if (Pmat) {
1254: PetscObjectReference((PetscObject)Pmat);
1255: MatDestroy(&ts->Brhs);
1256: ts->Brhs = Pmat;
1257: }
1258: return(0);
1259: }
1261: /*@C
1262: TSSetIFunction - Set the function to compute F(t,U,U_t) where F() = 0 is the DAE to be solved.
1264: Logically Collective on TS
1266: Input Parameters:
1267: + ts - the TS context obtained from TSCreate()
1268: . r - vector to hold the residual (or NULL to have it created internally)
1269: . f - the function evaluation routine
1270: - ctx - user-defined context for private data for the function evaluation routine (may be NULL)
1272: Calling sequence of f:
1273: $ PetscErrorCode f(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx);
1275: + t - time at step/stage being solved
1276: . u - state vector
1277: . u_t - time derivative of state vector
1278: . F - function vector
1279: - ctx - [optional] user-defined context for matrix evaluation routine
1281: Important:
1282: The user MUST call either this routine or TSSetRHSFunction() to define the ODE. When solving DAEs you must use this function.
1284: Level: beginner
1286: .seealso: TSSetRHSJacobian(), TSSetRHSFunction(), TSSetIJacobian()
1287: @*/
1288: PetscErrorCode TSSetIFunction(TS ts,Vec r,TSIFunction f,void *ctx)
1289: {
1291: SNES snes;
1292: Vec ralloc = NULL;
1293: DM dm;
1299: TSGetDM(ts,&dm);
1300: DMTSSetIFunction(dm,f,ctx);
1302: TSGetSNES(ts,&snes);
1303: if (!r && !ts->dm && ts->vec_sol) {
1304: VecDuplicate(ts->vec_sol,&ralloc);
1305: r = ralloc;
1306: }
1307: SNESSetFunction(snes,r,SNESTSFormFunction,ts);
1308: VecDestroy(&ralloc);
1309: return(0);
1310: }
1312: /*@C
1313: TSGetIFunction - Returns the vector where the implicit residual is stored and the function/contex to compute it.
1315: Not Collective
1317: Input Parameter:
1318: . ts - the TS context
1320: Output Parameter:
1321: + r - vector to hold residual (or NULL)
1322: . func - the function to compute residual (or NULL)
1323: - ctx - the function context (or NULL)
1325: Level: advanced
1327: .seealso: TSSetIFunction(), SNESGetFunction()
1328: @*/
1329: PetscErrorCode TSGetIFunction(TS ts,Vec *r,TSIFunction *func,void **ctx)
1330: {
1332: SNES snes;
1333: DM dm;
1337: TSGetSNES(ts,&snes);
1338: SNESGetFunction(snes,r,NULL,NULL);
1339: TSGetDM(ts,&dm);
1340: DMTSGetIFunction(dm,func,ctx);
1341: return(0);
1342: }
1344: /*@C
1345: TSGetRHSFunction - Returns the vector where the right hand side is stored and the function/context to compute it.
1347: Not Collective
1349: Input Parameter:
1350: . ts - the TS context
1352: Output Parameter:
1353: + r - vector to hold computed right hand side (or NULL)
1354: . func - the function to compute right hand side (or NULL)
1355: - ctx - the function context (or NULL)
1357: Level: advanced
1359: .seealso: TSSetRHSFunction(), SNESGetFunction()
1360: @*/
1361: PetscErrorCode TSGetRHSFunction(TS ts,Vec *r,TSRHSFunction *func,void **ctx)
1362: {
1364: SNES snes;
1365: DM dm;
1369: TSGetSNES(ts,&snes);
1370: SNESGetFunction(snes,r,NULL,NULL);
1371: TSGetDM(ts,&dm);
1372: DMTSGetRHSFunction(dm,func,ctx);
1373: return(0);
1374: }
1376: /*@C
1377: TSSetIJacobian - Set the function to compute the matrix dF/dU + a*dF/dU_t where F(t,U,U_t) is the function
1378: provided with TSSetIFunction().
1380: Logically Collective on TS
1382: Input Parameters:
1383: + ts - the TS context obtained from TSCreate()
1384: . Amat - (approximate) Jacobian matrix
1385: . Pmat - matrix used to compute preconditioner (usually the same as Amat)
1386: . f - the Jacobian evaluation routine
1387: - ctx - user-defined context for private data for the Jacobian evaluation routine (may be NULL)
1389: Calling sequence of f:
1390: $ PetscErrorCode f(TS ts,PetscReal t,Vec U,Vec U_t,PetscReal a,Mat Amat,Mat Pmat,void *ctx);
1392: + t - time at step/stage being solved
1393: . U - state vector
1394: . U_t - time derivative of state vector
1395: . a - shift
1396: . Amat - (approximate) Jacobian of F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t
1397: . Pmat - matrix used for constructing preconditioner, usually the same as Amat
1398: - ctx - [optional] user-defined context for matrix evaluation routine
1400: Notes:
1401: The matrices Amat and Pmat are exactly the matrices that are used by SNES for the nonlinear solve.
1403: If you know the operator Amat has a null space you can use MatSetNullSpace() and MatSetTransposeNullSpace() to supply the null
1404: space to Amat and the KSP solvers will automatically use that null space as needed during the solution process.
1406: The matrix dF/dU + a*dF/dU_t you provide turns out to be
1407: the Jacobian of F(t,U,W+a*U) where F(t,U,U_t) = 0 is the DAE to be solved.
1408: The time integrator internally approximates U_t by W+a*U where the positive "shift"
1409: a and vector W depend on the integration method, step size, and past states. For example with
1410: the backward Euler method a = 1/dt and W = -a*U(previous timestep) so
1411: W + a*U = a*(U - U(previous timestep)) = (U - U(previous timestep))/dt
1413: You must set all the diagonal entries of the matrices, if they are zero you must still set them with a zero value
1415: The TS solver may modify the nonzero structure and the entries of the matrices Amat and Pmat between the calls to f()
1416: You should not assume the values are the same in the next call to f() as you set them in the previous call.
1418: Level: beginner
1420: .seealso: TSSetIFunction(), TSSetRHSJacobian(), SNESComputeJacobianDefaultColor(), SNESComputeJacobianDefault(), TSSetRHSFunction()
1422: @*/
1423: PetscErrorCode TSSetIJacobian(TS ts,Mat Amat,Mat Pmat,TSIJacobian f,void *ctx)
1424: {
1426: SNES snes;
1427: DM dm;
1436: TSGetDM(ts,&dm);
1437: DMTSSetIJacobian(dm,f,ctx);
1439: TSGetSNES(ts,&snes);
1440: SNESSetJacobian(snes,Amat,Pmat,SNESTSFormJacobian,ts);
1441: return(0);
1442: }
1444: /*@
1445: TSRHSJacobianSetReuse - restore RHS Jacobian before re-evaluating. Without this flag, TS will change the sign and
1446: shift the RHS Jacobian for a finite-time-step implicit solve, in which case the user function will need to recompute
1447: the entire Jacobian. The reuse flag must be set if the evaluation function will assume that the matrix entries have
1448: not been changed by the TS.
1450: Logically Collective
1452: Input Arguments:
1453: + ts - TS context obtained from TSCreate()
1454: - reuse - PETSC_TRUE if the RHS Jacobian
1456: Level: intermediate
1458: .seealso: TSSetRHSJacobian(), TSComputeRHSJacobianConstant()
1459: @*/
1460: PetscErrorCode TSRHSJacobianSetReuse(TS ts,PetscBool reuse)
1461: {
1463: ts->rhsjacobian.reuse = reuse;
1464: return(0);
1465: }
1467: /*@C
1468: TSSetI2Function - Set the function to compute F(t,U,U_t,U_tt) where F = 0 is the DAE to be solved.
1470: Logically Collective on TS
1472: Input Parameters:
1473: + ts - the TS context obtained from TSCreate()
1474: . F - vector to hold the residual (or NULL to have it created internally)
1475: . fun - the function evaluation routine
1476: - ctx - user-defined context for private data for the function evaluation routine (may be NULL)
1478: Calling sequence of fun:
1479: $ PetscErrorCode fun(TS ts,PetscReal t,Vec U,Vec U_t,Vec U_tt,Vec F,ctx);
1481: + t - time at step/stage being solved
1482: . U - state vector
1483: . U_t - time derivative of state vector
1484: . U_tt - second time derivative of state vector
1485: . F - function vector
1486: - ctx - [optional] user-defined context for matrix evaluation routine (may be NULL)
1488: Level: beginner
1490: .seealso: TSSetI2Jacobian()
1491: @*/
1492: PetscErrorCode TSSetI2Function(TS ts,Vec F,TSI2Function fun,void *ctx)
1493: {
1494: DM dm;
1500: TSSetIFunction(ts,F,NULL,NULL);
1501: TSGetDM(ts,&dm);
1502: DMTSSetI2Function(dm,fun,ctx);
1503: return(0);
1504: }
1506: /*@C
1507: TSGetI2Function - Returns the vector where the implicit residual is stored and the function/contex to compute it.
1509: Not Collective
1511: Input Parameter:
1512: . ts - the TS context
1514: Output Parameter:
1515: + r - vector to hold residual (or NULL)
1516: . fun - the function to compute residual (or NULL)
1517: - ctx - the function context (or NULL)
1519: Level: advanced
1521: .seealso: TSSetI2Function(), SNESGetFunction()
1522: @*/
1523: PetscErrorCode TSGetI2Function(TS ts,Vec *r,TSI2Function *fun,void **ctx)
1524: {
1526: SNES snes;
1527: DM dm;
1531: TSGetSNES(ts,&snes);
1532: SNESGetFunction(snes,r,NULL,NULL);
1533: TSGetDM(ts,&dm);
1534: DMTSGetI2Function(dm,fun,ctx);
1535: return(0);
1536: }
1538: /*@C
1539: TSSetI2Jacobian - Set the function to compute the matrix dF/dU + v*dF/dU_t + a*dF/dU_tt
1540: where F(t,U,U_t,U_tt) is the function you provided with TSSetI2Function().
1542: Logically Collective on TS
1544: Input Parameters:
1545: + ts - the TS context obtained from TSCreate()
1546: . J - Jacobian matrix
1547: . P - preconditioning matrix for J (may be same as J)
1548: . jac - the Jacobian evaluation routine
1549: - ctx - user-defined context for private data for the Jacobian evaluation routine (may be NULL)
1551: Calling sequence of jac:
1552: $ PetscErrorCode jac(TS ts,PetscReal t,Vec U,Vec U_t,Vec U_tt,PetscReal v,PetscReal a,Mat J,Mat P,void *ctx);
1554: + t - time at step/stage being solved
1555: . U - state vector
1556: . U_t - time derivative of state vector
1557: . U_tt - second time derivative of state vector
1558: . v - shift for U_t
1559: . a - shift for U_tt
1560: . J - Jacobian of G(U) = F(t,U,W+v*U,W'+a*U), equivalent to dF/dU + v*dF/dU_t + a*dF/dU_tt
1561: . P - preconditioning matrix for J, may be same as J
1562: - ctx - [optional] user-defined context for matrix evaluation routine
1564: Notes:
1565: The matrices J and P are exactly the matrices that are used by SNES for the nonlinear solve.
1567: The matrix dF/dU + v*dF/dU_t + a*dF/dU_tt you provide turns out to be
1568: the Jacobian of G(U) = F(t,U,W+v*U,W'+a*U) where F(t,U,U_t,U_tt) = 0 is the DAE to be solved.
1569: The time integrator internally approximates U_t by W+v*U and U_tt by W'+a*U where the positive "shift"
1570: parameters 'v' and 'a' and vectors W, W' depend on the integration method, step size, and past states.
1572: Level: beginner
1574: .seealso: TSSetI2Function()
1575: @*/
1576: PetscErrorCode TSSetI2Jacobian(TS ts,Mat J,Mat P,TSI2Jacobian jac,void *ctx)
1577: {
1578: DM dm;
1585: TSSetIJacobian(ts,J,P,NULL,NULL);
1586: TSGetDM(ts,&dm);
1587: DMTSSetI2Jacobian(dm,jac,ctx);
1588: return(0);
1589: }
1591: /*@C
1592: TSGetI2Jacobian - Returns the implicit Jacobian at the present timestep.
1594: Not Collective, but parallel objects are returned if TS is parallel
1596: Input Parameter:
1597: . ts - The TS context obtained from TSCreate()
1599: Output Parameters:
1600: + J - The (approximate) Jacobian of F(t,U,U_t,U_tt)
1601: . P - The matrix from which the preconditioner is constructed, often the same as J
1602: . jac - The function to compute the Jacobian matrices
1603: - ctx - User-defined context for Jacobian evaluation routine
1605: Notes:
1606: You can pass in NULL for any return argument you do not need.
1608: Level: advanced
1610: .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetStepNumber()
1612: @*/
1613: PetscErrorCode TSGetI2Jacobian(TS ts,Mat *J,Mat *P,TSI2Jacobian *jac,void **ctx)
1614: {
1616: SNES snes;
1617: DM dm;
1620: TSGetSNES(ts,&snes);
1621: SNESSetUpMatrices(snes);
1622: SNESGetJacobian(snes,J,P,NULL,NULL);
1623: TSGetDM(ts,&dm);
1624: DMTSGetI2Jacobian(dm,jac,ctx);
1625: return(0);
1626: }
1628: /*@
1629: TSComputeI2Function - Evaluates the DAE residual written in implicit form F(t,U,U_t,U_tt) = 0
1631: Collective on TS
1633: Input Parameters:
1634: + ts - the TS context
1635: . t - current time
1636: . U - state vector
1637: . V - time derivative of state vector (U_t)
1638: - A - second time derivative of state vector (U_tt)
1640: Output Parameter:
1641: . F - the residual vector
1643: Note:
1644: Most users should not need to explicitly call this routine, as it
1645: is used internally within the nonlinear solvers.
1647: Level: developer
1649: .seealso: TSSetI2Function()
1650: @*/
1651: PetscErrorCode TSComputeI2Function(TS ts,PetscReal t,Vec U,Vec V,Vec A,Vec F)
1652: {
1653: DM dm;
1654: TSI2Function I2Function;
1655: void *ctx;
1656: TSRHSFunction rhsfunction;
1666: TSGetDM(ts,&dm);
1667: DMTSGetI2Function(dm,&I2Function,&ctx);
1668: DMTSGetRHSFunction(dm,&rhsfunction,NULL);
1670: if (!I2Function) {
1671: TSComputeIFunction(ts,t,U,A,F,PETSC_FALSE);
1672: return(0);
1673: }
1675: PetscLogEventBegin(TS_FunctionEval,ts,U,V,F);
1677: PetscStackPush("TS user implicit function");
1678: I2Function(ts,t,U,V,A,F,ctx);
1679: PetscStackPop;
1681: if (rhsfunction) {
1682: Vec Frhs;
1683: TSGetRHSVec_Private(ts,&Frhs);
1684: TSComputeRHSFunction(ts,t,U,Frhs);
1685: VecAXPY(F,-1,Frhs);
1686: }
1688: PetscLogEventEnd(TS_FunctionEval,ts,U,V,F);
1689: return(0);
1690: }
1692: /*@
1693: TSComputeI2Jacobian - Evaluates the Jacobian of the DAE
1695: Collective on TS
1697: Input Parameters:
1698: + ts - the TS context
1699: . t - current timestep
1700: . U - state vector
1701: . V - time derivative of state vector
1702: . A - second time derivative of state vector
1703: . shiftV - shift to apply, see note below
1704: - shiftA - shift to apply, see note below
1706: Output Parameters:
1707: + J - Jacobian matrix
1708: - P - optional preconditioning matrix
1710: Notes:
1711: If F(t,U,V,A)=0 is the DAE, the required Jacobian is
1713: dF/dU + shiftV*dF/dV + shiftA*dF/dA
1715: Most users should not need to explicitly call this routine, as it
1716: is used internally within the nonlinear solvers.
1718: Level: developer
1720: .seealso: TSSetI2Jacobian()
1721: @*/
1722: PetscErrorCode TSComputeI2Jacobian(TS ts,PetscReal t,Vec U,Vec V,Vec A,PetscReal shiftV,PetscReal shiftA,Mat J,Mat P)
1723: {
1724: DM dm;
1725: TSI2Jacobian I2Jacobian;
1726: void *ctx;
1727: TSRHSJacobian rhsjacobian;
1738: TSGetDM(ts,&dm);
1739: DMTSGetI2Jacobian(dm,&I2Jacobian,&ctx);
1740: DMTSGetRHSJacobian(dm,&rhsjacobian,NULL);
1742: if (!I2Jacobian) {
1743: TSComputeIJacobian(ts,t,U,A,shiftA,J,P,PETSC_FALSE);
1744: return(0);
1745: }
1747: PetscLogEventBegin(TS_JacobianEval,ts,U,J,P);
1749: PetscStackPush("TS user implicit Jacobian");
1750: I2Jacobian(ts,t,U,V,A,shiftV,shiftA,J,P,ctx);
1751: PetscStackPop;
1753: if (rhsjacobian) {
1754: Mat Jrhs,Prhs; MatStructure axpy = DIFFERENT_NONZERO_PATTERN;
1755: TSGetRHSMats_Private(ts,&Jrhs,&Prhs);
1756: TSComputeRHSJacobian(ts,t,U,Jrhs,Prhs);
1757: MatAXPY(J,-1,Jrhs,axpy);
1758: if (P != J) {MatAXPY(P,-1,Prhs,axpy);}
1759: }
1761: PetscLogEventEnd(TS_JacobianEval,ts,U,J,P);
1762: return(0);
1763: }
1765: /*@
1766: TS2SetSolution - Sets the initial solution and time derivative vectors
1767: for use by the TS routines handling second order equations.
1769: Logically Collective on TS
1771: Input Parameters:
1772: + ts - the TS context obtained from TSCreate()
1773: . u - the solution vector
1774: - v - the time derivative vector
1776: Level: beginner
1778: @*/
1779: PetscErrorCode TS2SetSolution(TS ts,Vec u,Vec v)
1780: {
1787: TSSetSolution(ts,u);
1788: PetscObjectReference((PetscObject)v);
1789: VecDestroy(&ts->vec_dot);
1790: ts->vec_dot = v;
1791: return(0);
1792: }
1794: /*@
1795: TS2GetSolution - Returns the solution and time derivative at the present timestep
1796: for second order equations. It is valid to call this routine inside the function
1797: that you are evaluating in order to move to the new timestep. This vector not
1798: changed until the solution at the next timestep has been calculated.
1800: Not Collective, but Vec returned is parallel if TS is parallel
1802: Input Parameter:
1803: . ts - the TS context obtained from TSCreate()
1805: Output Parameter:
1806: + u - the vector containing the solution
1807: - v - the vector containing the time derivative
1809: Level: intermediate
1811: .seealso: TS2SetSolution(), TSGetTimeStep(), TSGetTime()
1813: @*/
1814: PetscErrorCode TS2GetSolution(TS ts,Vec *u,Vec *v)
1815: {
1820: if (u) *u = ts->vec_sol;
1821: if (v) *v = ts->vec_dot;
1822: return(0);
1823: }
1825: /*@C
1826: TSLoad - Loads a KSP that has been stored in binary with KSPView().
1828: Collective on PetscViewer
1830: Input Parameters:
1831: + newdm - the newly loaded TS, this needs to have been created with TSCreate() or
1832: some related function before a call to TSLoad().
1833: - viewer - binary file viewer, obtained from PetscViewerBinaryOpen()
1835: Level: intermediate
1837: Notes:
1838: The type is determined by the data in the file, any type set into the TS before this call is ignored.
1840: Notes for advanced users:
1841: Most users should not need to know the details of the binary storage
1842: format, since TSLoad() and TSView() completely hide these details.
1843: But for anyone who's interested, the standard binary matrix storage
1844: format is
1845: .vb
1846: has not yet been determined
1847: .ve
1849: .seealso: PetscViewerBinaryOpen(), TSView(), MatLoad(), VecLoad()
1850: @*/
1851: PetscErrorCode TSLoad(TS ts, PetscViewer viewer)
1852: {
1854: PetscBool isbinary;
1855: PetscInt classid;
1856: char type[256];
1857: DMTS sdm;
1858: DM dm;
1863: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
1864: if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");
1866: PetscViewerBinaryRead(viewer,&classid,1,NULL,PETSC_INT);
1867: if (classid != TS_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONG,"Not TS next in file");
1868: PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);
1869: TSSetType(ts, type);
1870: if (ts->ops->load) {
1871: (*ts->ops->load)(ts,viewer);
1872: }
1873: DMCreate(PetscObjectComm((PetscObject)ts),&dm);
1874: DMLoad(dm,viewer);
1875: TSSetDM(ts,dm);
1876: DMCreateGlobalVector(ts->dm,&ts->vec_sol);
1877: VecLoad(ts->vec_sol,viewer);
1878: DMGetDMTS(ts->dm,&sdm);
1879: DMTSLoad(sdm,viewer);
1880: return(0);
1881: }
1883: #include <petscdraw.h>
1884: #if defined(PETSC_HAVE_SAWS)
1885: #include <petscviewersaws.h>
1886: #endif
1887: /*@C
1888: TSView - Prints the TS data structure.
1890: Collective on TS
1892: Input Parameters:
1893: + ts - the TS context obtained from TSCreate()
1894: - viewer - visualization context
1896: Options Database Key:
1897: . -ts_view - calls TSView() at end of TSStep()
1899: Notes:
1900: The available visualization contexts include
1901: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
1902: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
1903: output where only the first processor opens
1904: the file. All other processors send their
1905: data to the first processor to print.
1907: The user can open an alternative visualization context with
1908: PetscViewerASCIIOpen() - output to a specified file.
1910: Level: beginner
1912: .seealso: PetscViewerASCIIOpen()
1913: @*/
1914: PetscErrorCode TSView(TS ts,PetscViewer viewer)
1915: {
1917: TSType type;
1918: PetscBool iascii,isstring,isundials,isbinary,isdraw;
1919: DMTS sdm;
1920: #if defined(PETSC_HAVE_SAWS)
1921: PetscBool issaws;
1922: #endif
1926: if (!viewer) {
1927: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ts),&viewer);
1928: }
1932: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1933: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
1934: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
1935: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
1936: #if defined(PETSC_HAVE_SAWS)
1937: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSAWS,&issaws);
1938: #endif
1939: if (iascii) {
1940: PetscObjectPrintClassNamePrefixType((PetscObject)ts,viewer);
1941: if (ts->ops->view) {
1942: PetscViewerASCIIPushTab(viewer);
1943: (*ts->ops->view)(ts,viewer);
1944: PetscViewerASCIIPopTab(viewer);
1945: }
1946: if (ts->max_steps < PETSC_MAX_INT) {
1947: PetscViewerASCIIPrintf(viewer," maximum steps=%D\n",ts->max_steps);
1948: }
1949: if (ts->max_time < PETSC_MAX_REAL) {
1950: PetscViewerASCIIPrintf(viewer," maximum time=%g\n",(double)ts->max_time);
1951: }
1952: if (ts->usessnes) {
1953: PetscBool lin;
1954: if (ts->problem_type == TS_NONLINEAR) {
1955: PetscViewerASCIIPrintf(viewer," total number of nonlinear solver iterations=%D\n",ts->snes_its);
1956: }
1957: PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%D\n",ts->ksp_its);
1958: PetscObjectTypeCompareAny((PetscObject)ts->snes,&lin,SNESKSPONLY,SNESKSPTRANSPOSEONLY,"");
1959: PetscViewerASCIIPrintf(viewer," total number of %slinear solve failures=%D\n",lin ? "" : "non",ts->num_snes_failures);
1960: }
1961: PetscViewerASCIIPrintf(viewer," total number of rejected steps=%D\n",ts->reject);
1962: if (ts->vrtol) {
1963: PetscViewerASCIIPrintf(viewer," using vector of relative error tolerances, ");
1964: } else {
1965: PetscViewerASCIIPrintf(viewer," using relative error tolerance of %g, ",(double)ts->rtol);
1966: }
1967: if (ts->vatol) {
1968: PetscViewerASCIIPrintf(viewer," using vector of absolute error tolerances\n");
1969: } else {
1970: PetscViewerASCIIPrintf(viewer," using absolute error tolerance of %g\n",(double)ts->atol);
1971: }
1972: PetscViewerASCIIPushTab(viewer);
1973: TSAdaptView(ts->adapt,viewer);
1974: PetscViewerASCIIPopTab(viewer);
1975: } else if (isstring) {
1976: TSGetType(ts,&type);
1977: PetscViewerStringSPrintf(viewer," TSType: %-7.7s",type);
1978: if (ts->ops->view) {(*ts->ops->view)(ts,viewer);}
1979: } else if (isbinary) {
1980: PetscInt classid = TS_FILE_CLASSID;
1981: MPI_Comm comm;
1982: PetscMPIInt rank;
1983: char type[256];
1985: PetscObjectGetComm((PetscObject)ts,&comm);
1986: MPI_Comm_rank(comm,&rank);
1987: if (!rank) {
1988: PetscViewerBinaryWrite(viewer,&classid,1,PETSC_INT,PETSC_FALSE);
1989: PetscStrncpy(type,((PetscObject)ts)->type_name,256);
1990: PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);
1991: }
1992: if (ts->ops->view) {
1993: (*ts->ops->view)(ts,viewer);
1994: }
1995: if (ts->adapt) {TSAdaptView(ts->adapt,viewer);}
1996: DMView(ts->dm,viewer);
1997: VecView(ts->vec_sol,viewer);
1998: DMGetDMTS(ts->dm,&sdm);
1999: DMTSView(sdm,viewer);
2000: } else if (isdraw) {
2001: PetscDraw draw;
2002: char str[36];
2003: PetscReal x,y,bottom,h;
2005: PetscViewerDrawGetDraw(viewer,0,&draw);
2006: PetscDrawGetCurrentPoint(draw,&x,&y);
2007: PetscStrcpy(str,"TS: ");
2008: PetscStrcat(str,((PetscObject)ts)->type_name);
2009: PetscDrawStringBoxed(draw,x,y,PETSC_DRAW_BLACK,PETSC_DRAW_BLACK,str,NULL,&h);
2010: bottom = y - h;
2011: PetscDrawPushCurrentPoint(draw,x,bottom);
2012: if (ts->ops->view) {
2013: (*ts->ops->view)(ts,viewer);
2014: }
2015: if (ts->adapt) {TSAdaptView(ts->adapt,viewer);}
2016: if (ts->snes) {SNESView(ts->snes,viewer);}
2017: PetscDrawPopCurrentPoint(draw);
2018: #if defined(PETSC_HAVE_SAWS)
2019: } else if (issaws) {
2020: PetscMPIInt rank;
2021: const char *name;
2023: PetscObjectGetName((PetscObject)ts,&name);
2024: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
2025: if (!((PetscObject)ts)->amsmem && !rank) {
2026: char dir[1024];
2028: PetscObjectViewSAWs((PetscObject)ts,viewer);
2029: PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/time_step",name);
2030: PetscStackCallSAWs(SAWs_Register,(dir,&ts->steps,1,SAWs_READ,SAWs_INT));
2031: PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/time",name);
2032: PetscStackCallSAWs(SAWs_Register,(dir,&ts->ptime,1,SAWs_READ,SAWs_DOUBLE));
2033: }
2034: if (ts->ops->view) {
2035: (*ts->ops->view)(ts,viewer);
2036: }
2037: #endif
2038: }
2039: if (ts->snes && ts->usessnes) {
2040: PetscViewerASCIIPushTab(viewer);
2041: SNESView(ts->snes,viewer);
2042: PetscViewerASCIIPopTab(viewer);
2043: }
2044: DMGetDMTS(ts->dm,&sdm);
2045: DMTSView(sdm,viewer);
2047: PetscViewerASCIIPushTab(viewer);
2048: PetscObjectTypeCompare((PetscObject)ts,TSSUNDIALS,&isundials);
2049: PetscViewerASCIIPopTab(viewer);
2050: return(0);
2051: }
2053: /*@
2054: TSSetApplicationContext - Sets an optional user-defined context for
2055: the timesteppers.
2057: Logically Collective on TS
2059: Input Parameters:
2060: + ts - the TS context obtained from TSCreate()
2061: - usrP - optional user context
2063: Fortran Notes:
2064: To use this from Fortran you must write a Fortran interface definition for this
2065: function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
2067: Level: intermediate
2069: .seealso: TSGetApplicationContext()
2070: @*/
2071: PetscErrorCode TSSetApplicationContext(TS ts,void *usrP)
2072: {
2075: ts->user = usrP;
2076: return(0);
2077: }
2079: /*@
2080: TSGetApplicationContext - Gets the user-defined context for the
2081: timestepper.
2083: Not Collective
2085: Input Parameter:
2086: . ts - the TS context obtained from TSCreate()
2088: Output Parameter:
2089: . usrP - user context
2091: Fortran Notes:
2092: To use this from Fortran you must write a Fortran interface definition for this
2093: function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
2095: Level: intermediate
2097: .seealso: TSSetApplicationContext()
2098: @*/
2099: PetscErrorCode TSGetApplicationContext(TS ts,void *usrP)
2100: {
2103: *(void**)usrP = ts->user;
2104: return(0);
2105: }
2107: /*@
2108: TSGetStepNumber - Gets the number of steps completed.
2110: Not Collective
2112: Input Parameter:
2113: . ts - the TS context obtained from TSCreate()
2115: Output Parameter:
2116: . steps - number of steps completed so far
2118: Level: intermediate
2120: .seealso: TSGetTime(), TSGetTimeStep(), TSSetPreStep(), TSSetPreStage(), TSSetPostStage(), TSSetPostStep()
2121: @*/
2122: PetscErrorCode TSGetStepNumber(TS ts,PetscInt *steps)
2123: {
2127: *steps = ts->steps;
2128: return(0);
2129: }
2131: /*@
2132: TSSetStepNumber - Sets the number of steps completed.
2134: Logically Collective on TS
2136: Input Parameters:
2137: + ts - the TS context
2138: - steps - number of steps completed so far
2140: Notes:
2141: For most uses of the TS solvers the user need not explicitly call
2142: TSSetStepNumber(), as the step counter is appropriately updated in
2143: TSSolve()/TSStep()/TSRollBack(). Power users may call this routine to
2144: reinitialize timestepping by setting the step counter to zero (and time
2145: to the initial time) to solve a similar problem with different initial
2146: conditions or parameters. Other possible use case is to continue
2147: timestepping from a previously interrupted run in such a way that TS
2148: monitors will be called with a initial nonzero step counter.
2150: Level: advanced
2152: .seealso: TSGetStepNumber(), TSSetTime(), TSSetTimeStep(), TSSetSolution()
2153: @*/
2154: PetscErrorCode TSSetStepNumber(TS ts,PetscInt steps)
2155: {
2159: if (steps < 0) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Step number must be non-negative");
2160: ts->steps = steps;
2161: return(0);
2162: }
2164: /*@
2165: TSSetTimeStep - Allows one to reset the timestep at any time,
2166: useful for simple pseudo-timestepping codes.
2168: Logically Collective on TS
2170: Input Parameters:
2171: + ts - the TS context obtained from TSCreate()
2172: - time_step - the size of the timestep
2174: Level: intermediate
2176: .seealso: TSGetTimeStep(), TSSetTime()
2178: @*/
2179: PetscErrorCode TSSetTimeStep(TS ts,PetscReal time_step)
2180: {
2184: ts->time_step = time_step;
2185: return(0);
2186: }
2188: /*@
2189: TSSetExactFinalTime - Determines whether to adapt the final time step to
2190: match the exact final time, interpolate solution to the exact final time,
2191: or just return at the final time TS computed.
2193: Logically Collective on TS
2195: Input Parameter:
2196: + ts - the time-step context
2197: - eftopt - exact final time option
2199: $ TS_EXACTFINALTIME_STEPOVER - Don't do anything if final time is exceeded
2200: $ TS_EXACTFINALTIME_INTERPOLATE - Interpolate back to final time
2201: $ TS_EXACTFINALTIME_MATCHSTEP - Adapt final time step to match the final time
2203: Options Database:
2204: . -ts_exact_final_time <stepover,interpolate,matchstep> - select the final step at runtime
2206: Warning: If you use the option TS_EXACTFINALTIME_STEPOVER the solution may be at a very different time
2207: then the final time you selected.
2209: Level: beginner
2211: .seealso: TSExactFinalTimeOption, TSGetExactFinalTime()
2212: @*/
2213: PetscErrorCode TSSetExactFinalTime(TS ts,TSExactFinalTimeOption eftopt)
2214: {
2218: ts->exact_final_time = eftopt;
2219: return(0);
2220: }
2222: /*@
2223: TSGetExactFinalTime - Gets the exact final time option.
2225: Not Collective
2227: Input Parameter:
2228: . ts - the TS context
2230: Output Parameter:
2231: . eftopt - exact final time option
2233: Level: beginner
2235: .seealso: TSExactFinalTimeOption, TSSetExactFinalTime()
2236: @*/
2237: PetscErrorCode TSGetExactFinalTime(TS ts,TSExactFinalTimeOption *eftopt)
2238: {
2242: *eftopt = ts->exact_final_time;
2243: return(0);
2244: }
2246: /*@
2247: TSGetTimeStep - Gets the current timestep size.
2249: Not Collective
2251: Input Parameter:
2252: . ts - the TS context obtained from TSCreate()
2254: Output Parameter:
2255: . dt - the current timestep size
2257: Level: intermediate
2259: .seealso: TSSetTimeStep(), TSGetTime()
2261: @*/
2262: PetscErrorCode TSGetTimeStep(TS ts,PetscReal *dt)
2263: {
2267: *dt = ts->time_step;
2268: return(0);
2269: }
2271: /*@
2272: TSGetSolution - Returns the solution at the present timestep. It
2273: is valid to call this routine inside the function that you are evaluating
2274: in order to move to the new timestep. This vector not changed until
2275: the solution at the next timestep has been calculated.
2277: Not Collective, but Vec returned is parallel if TS is parallel
2279: Input Parameter:
2280: . ts - the TS context obtained from TSCreate()
2282: Output Parameter:
2283: . v - the vector containing the solution
2285: Note: If you used TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP); this does not return the solution at the requested
2286: final time. It returns the solution at the next timestep.
2288: Level: intermediate
2290: .seealso: TSGetTimeStep(), TSGetTime(), TSGetSolveTime(), TSGetSolutionComponents(), TSSetSolutionFunction()
2292: @*/
2293: PetscErrorCode TSGetSolution(TS ts,Vec *v)
2294: {
2298: *v = ts->vec_sol;
2299: return(0);
2300: }
2302: /*@
2303: TSGetSolutionComponents - Returns any solution components at the present
2304: timestep, if available for the time integration method being used.
2305: Solution components are quantities that share the same size and
2306: structure as the solution vector.
2308: Not Collective, but Vec returned is parallel if TS is parallel
2310: Parameters :
2311: + ts - the TS context obtained from TSCreate() (input parameter).
2312: . n - If v is PETSC_NULL, then the number of solution components is
2313: returned through n, else the n-th solution component is
2314: returned in v.
2315: - v - the vector containing the n-th solution component
2316: (may be PETSC_NULL to use this function to find out
2317: the number of solutions components).
2319: Level: advanced
2321: .seealso: TSGetSolution()
2323: @*/
2324: PetscErrorCode TSGetSolutionComponents(TS ts,PetscInt *n,Vec *v)
2325: {
2330: if (!ts->ops->getsolutioncomponents) *n = 0;
2331: else {
2332: (*ts->ops->getsolutioncomponents)(ts,n,v);
2333: }
2334: return(0);
2335: }
2337: /*@
2338: TSGetAuxSolution - Returns an auxiliary solution at the present
2339: timestep, if available for the time integration method being used.
2341: Not Collective, but Vec returned is parallel if TS is parallel
2343: Parameters :
2344: + ts - the TS context obtained from TSCreate() (input parameter).
2345: - v - the vector containing the auxiliary solution
2347: Level: intermediate
2349: .seealso: TSGetSolution()
2351: @*/
2352: PetscErrorCode TSGetAuxSolution(TS ts,Vec *v)
2353: {
2358: if (ts->ops->getauxsolution) {
2359: (*ts->ops->getauxsolution)(ts,v);
2360: } else {
2361: VecZeroEntries(*v);
2362: }
2363: return(0);
2364: }
2366: /*@
2367: TSGetTimeError - Returns the estimated error vector, if the chosen
2368: TSType has an error estimation functionality.
2370: Not Collective, but Vec returned is parallel if TS is parallel
2372: Note: MUST call after TSSetUp()
2374: Parameters :
2375: + ts - the TS context obtained from TSCreate() (input parameter).
2376: . n - current estimate (n=0) or previous one (n=-1)
2377: - v - the vector containing the error (same size as the solution).
2379: Level: intermediate
2381: .seealso: TSGetSolution(), TSSetTimeError()
2383: @*/
2384: PetscErrorCode TSGetTimeError(TS ts,PetscInt n,Vec *v)
2385: {
2390: if (ts->ops->gettimeerror) {
2391: (*ts->ops->gettimeerror)(ts,n,v);
2392: } else {
2393: VecZeroEntries(*v);
2394: }
2395: return(0);
2396: }
2398: /*@
2399: TSSetTimeError - Sets the estimated error vector, if the chosen
2400: TSType has an error estimation functionality. This can be used
2401: to restart such a time integrator with a given error vector.
2403: Not Collective, but Vec returned is parallel if TS is parallel
2405: Parameters :
2406: + ts - the TS context obtained from TSCreate() (input parameter).
2407: - v - the vector containing the error (same size as the solution).
2409: Level: intermediate
2411: .seealso: TSSetSolution(), TSGetTimeError)
2413: @*/
2414: PetscErrorCode TSSetTimeError(TS ts,Vec v)
2415: {
2420: if (!ts->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetUp() first");
2421: if (ts->ops->settimeerror) {
2422: (*ts->ops->settimeerror)(ts,v);
2423: }
2424: return(0);
2425: }
2427: /* ----- Routines to initialize and destroy a timestepper ---- */
2428: /*@
2429: TSSetProblemType - Sets the type of problem to be solved.
2431: Not collective
2433: Input Parameters:
2434: + ts - The TS
2435: - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
2436: .vb
2437: U_t - A U = 0 (linear)
2438: U_t - A(t) U = 0 (linear)
2439: F(t,U,U_t) = 0 (nonlinear)
2440: .ve
2442: Level: beginner
2444: .seealso: TSSetUp(), TSProblemType, TS
2445: @*/
2446: PetscErrorCode TSSetProblemType(TS ts, TSProblemType type)
2447: {
2452: ts->problem_type = type;
2453: if (type == TS_LINEAR) {
2454: SNES snes;
2455: TSGetSNES(ts,&snes);
2456: SNESSetType(snes,SNESKSPONLY);
2457: }
2458: return(0);
2459: }
2461: /*@C
2462: TSGetProblemType - Gets the type of problem to be solved.
2464: Not collective
2466: Input Parameter:
2467: . ts - The TS
2469: Output Parameter:
2470: . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
2471: .vb
2472: M U_t = A U
2473: M(t) U_t = A(t) U
2474: F(t,U,U_t)
2475: .ve
2477: Level: beginner
2479: .seealso: TSSetUp(), TSProblemType, TS
2480: @*/
2481: PetscErrorCode TSGetProblemType(TS ts, TSProblemType *type)
2482: {
2486: *type = ts->problem_type;
2487: return(0);
2488: }
2490: /*@
2491: TSSetUp - Sets up the internal data structures for the later use
2492: of a timestepper.
2494: Collective on TS
2496: Input Parameter:
2497: . ts - the TS context obtained from TSCreate()
2499: Notes:
2500: For basic use of the TS solvers the user need not explicitly call
2501: TSSetUp(), since these actions will automatically occur during
2502: the call to TSStep() or TSSolve(). However, if one wishes to control this
2503: phase separately, TSSetUp() should be called after TSCreate()
2504: and optional routines of the form TSSetXXX(), but before TSStep() and TSSolve().
2506: Level: advanced
2508: .seealso: TSCreate(), TSStep(), TSDestroy(), TSSolve()
2509: @*/
2510: PetscErrorCode TSSetUp(TS ts)
2511: {
2513: DM dm;
2514: PetscErrorCode (*func)(SNES,Vec,Vec,void*);
2515: PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*);
2516: TSIFunction ifun;
2517: TSIJacobian ijac;
2518: TSI2Jacobian i2jac;
2519: TSRHSJacobian rhsjac;
2520: PetscBool isnone;
2524: if (ts->setupcalled) return(0);
2526: if (!((PetscObject)ts)->type_name) {
2527: TSGetIFunction(ts,NULL,&ifun,NULL);
2528: TSSetType(ts,ifun ? TSBEULER : TSEULER);
2529: }
2531: if (!ts->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first");
2533: if (!ts->Jacp && ts->Jacprhs) { /* IJacobianP shares the same matrix with RHSJacobianP if only RHSJacobianP is provided */
2534: PetscObjectReference((PetscObject)ts->Jacprhs);
2535: ts->Jacp = ts->Jacprhs;
2536: }
2538: if (ts->quadraturets) {
2539: TSSetUp(ts->quadraturets);
2540: VecDestroy(&ts->vec_costintegrand);
2541: VecDuplicate(ts->quadraturets->vec_sol,&ts->vec_costintegrand);
2542: }
2544: TSGetRHSJacobian(ts,NULL,NULL,&rhsjac,NULL);
2545: if (ts->rhsjacobian.reuse && rhsjac == TSComputeRHSJacobianConstant) {
2546: Mat Amat,Pmat;
2547: SNES snes;
2548: TSGetSNES(ts,&snes);
2549: SNESGetJacobian(snes,&Amat,&Pmat,NULL,NULL);
2550: /* Matching matrices implies that an IJacobian is NOT set, because if it had been set, the IJacobian's matrix would
2551: * have displaced the RHS matrix */
2552: if (Amat && Amat == ts->Arhs) {
2553: /* we need to copy the values of the matrix because for the constant Jacobian case the user will never set the numerical values in this new location */
2554: MatDuplicate(ts->Arhs,MAT_COPY_VALUES,&Amat);
2555: SNESSetJacobian(snes,Amat,NULL,NULL,NULL);
2556: MatDestroy(&Amat);
2557: }
2558: if (Pmat && Pmat == ts->Brhs) {
2559: MatDuplicate(ts->Brhs,MAT_COPY_VALUES,&Pmat);
2560: SNESSetJacobian(snes,NULL,Pmat,NULL,NULL);
2561: MatDestroy(&Pmat);
2562: }
2563: }
2565: TSGetAdapt(ts,&ts->adapt);
2566: TSAdaptSetDefaultType(ts->adapt,ts->default_adapt_type);
2568: if (ts->ops->setup) {
2569: (*ts->ops->setup)(ts);
2570: }
2572: /* Attempt to check/preset a default value for the exact final time option */
2573: PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&isnone);
2574: if (!isnone && ts->exact_final_time == TS_EXACTFINALTIME_UNSPECIFIED)
2575: ts->exact_final_time = TS_EXACTFINALTIME_MATCHSTEP;
2577: /* In the case where we've set a DMTSFunction or what have you, we need the default SNESFunction
2578: to be set right but can't do it elsewhere due to the overreliance on ctx=ts.
2579: */
2580: TSGetDM(ts,&dm);
2581: DMSNESGetFunction(dm,&func,NULL);
2582: if (!func) {
2583: DMSNESSetFunction(dm,SNESTSFormFunction,ts);
2584: }
2585: /* If the SNES doesn't have a jacobian set and the TS has an ijacobian or rhsjacobian set, set the SNES to use it.
2586: Otherwise, the SNES will use coloring internally to form the Jacobian.
2587: */
2588: DMSNESGetJacobian(dm,&jac,NULL);
2589: DMTSGetIJacobian(dm,&ijac,NULL);
2590: DMTSGetI2Jacobian(dm,&i2jac,NULL);
2591: DMTSGetRHSJacobian(dm,&rhsjac,NULL);
2592: if (!jac && (ijac || i2jac || rhsjac)) {
2593: DMSNESSetJacobian(dm,SNESTSFormJacobian,ts);
2594: }
2596: /* if time integration scheme has a starting method, call it */
2597: if (ts->ops->startingmethod) {
2598: (*ts->ops->startingmethod)(ts);
2599: }
2601: ts->setupcalled = PETSC_TRUE;
2602: return(0);
2603: }
2605: /*@
2606: TSReset - Resets a TS context and removes any allocated Vecs and Mats.
2608: Collective on TS
2610: Input Parameter:
2611: . ts - the TS context obtained from TSCreate()
2613: Level: beginner
2615: .seealso: TSCreate(), TSSetup(), TSDestroy()
2616: @*/
2617: PetscErrorCode TSReset(TS ts)
2618: {
2619: TS_RHSSplitLink ilink = ts->tsrhssplit,next;
2620: PetscErrorCode ierr;
2625: if (ts->ops->reset) {
2626: (*ts->ops->reset)(ts);
2627: }
2628: if (ts->snes) {SNESReset(ts->snes);}
2629: if (ts->adapt) {TSAdaptReset(ts->adapt);}
2631: MatDestroy(&ts->Arhs);
2632: MatDestroy(&ts->Brhs);
2633: VecDestroy(&ts->Frhs);
2634: VecDestroy(&ts->vec_sol);
2635: VecDestroy(&ts->vec_dot);
2636: VecDestroy(&ts->vatol);
2637: VecDestroy(&ts->vrtol);
2638: VecDestroyVecs(ts->nwork,&ts->work);
2640: MatDestroy(&ts->Jacprhs);
2641: MatDestroy(&ts->Jacp);
2642: if (ts->forward_solve) {
2643: TSForwardReset(ts);
2644: }
2645: if (ts->quadraturets) {
2646: TSReset(ts->quadraturets);
2647: VecDestroy(&ts->vec_costintegrand);
2648: }
2649: while (ilink) {
2650: next = ilink->next;
2651: TSDestroy(&ilink->ts);
2652: PetscFree(ilink->splitname);
2653: ISDestroy(&ilink->is);
2654: PetscFree(ilink);
2655: ilink = next;
2656: }
2657: ts->num_rhs_splits = 0;
2658: ts->setupcalled = PETSC_FALSE;
2659: return(0);
2660: }
2662: /*@
2663: TSDestroy - Destroys the timestepper context that was created
2664: with TSCreate().
2666: Collective on TS
2668: Input Parameter:
2669: . ts - the TS context obtained from TSCreate()
2671: Level: beginner
2673: .seealso: TSCreate(), TSSetUp(), TSSolve()
2674: @*/
2675: PetscErrorCode TSDestroy(TS *ts)
2676: {
2680: if (!*ts) return(0);
2682: if (--((PetscObject)(*ts))->refct > 0) {*ts = 0; return(0);}
2684: TSReset(*ts);
2685: TSAdjointReset(*ts);
2686: if ((*ts)->forward_solve) {
2687: TSForwardReset(*ts);
2688: }
2689: /* if memory was published with SAWs then destroy it */
2690: PetscObjectSAWsViewOff((PetscObject)*ts);
2691: if ((*ts)->ops->destroy) {(*(*ts)->ops->destroy)((*ts));}
2693: TSTrajectoryDestroy(&(*ts)->trajectory);
2695: TSAdaptDestroy(&(*ts)->adapt);
2696: TSEventDestroy(&(*ts)->event);
2698: SNESDestroy(&(*ts)->snes);
2699: DMDestroy(&(*ts)->dm);
2700: TSMonitorCancel((*ts));
2701: TSAdjointMonitorCancel((*ts));
2703: TSDestroy(&(*ts)->quadraturets);
2704: PetscHeaderDestroy(ts);
2705: return(0);
2706: }
2708: /*@
2709: TSGetSNES - Returns the SNES (nonlinear solver) associated with
2710: a TS (timestepper) context. Valid only for nonlinear problems.
2712: Not Collective, but SNES is parallel if TS is parallel
2714: Input Parameter:
2715: . ts - the TS context obtained from TSCreate()
2717: Output Parameter:
2718: . snes - the nonlinear solver context
2720: Notes:
2721: The user can then directly manipulate the SNES context to set various
2722: options, etc. Likewise, the user can then extract and manipulate the
2723: KSP, KSP, and PC contexts as well.
2725: TSGetSNES() does not work for integrators that do not use SNES; in
2726: this case TSGetSNES() returns NULL in snes.
2728: Level: beginner
2730: @*/
2731: PetscErrorCode TSGetSNES(TS ts,SNES *snes)
2732: {
2738: if (!ts->snes) {
2739: SNESCreate(PetscObjectComm((PetscObject)ts),&ts->snes);
2740: PetscObjectSetOptions((PetscObject)ts->snes,((PetscObject)ts)->options);
2741: SNESSetFunction(ts->snes,NULL,SNESTSFormFunction,ts);
2742: PetscLogObjectParent((PetscObject)ts,(PetscObject)ts->snes);
2743: PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);
2744: if (ts->dm) {SNESSetDM(ts->snes,ts->dm);}
2745: if (ts->problem_type == TS_LINEAR) {
2746: SNESSetType(ts->snes,SNESKSPONLY);
2747: }
2748: }
2749: *snes = ts->snes;
2750: return(0);
2751: }
2753: /*@
2754: TSSetSNES - Set the SNES (nonlinear solver) to be used by the timestepping context
2756: Collective
2758: Input Parameter:
2759: + ts - the TS context obtained from TSCreate()
2760: - snes - the nonlinear solver context
2762: Notes:
2763: Most users should have the TS created by calling TSGetSNES()
2765: Level: developer
2767: @*/
2768: PetscErrorCode TSSetSNES(TS ts,SNES snes)
2769: {
2771: PetscErrorCode (*func)(SNES,Vec,Mat,Mat,void*);
2776: PetscObjectReference((PetscObject)snes);
2777: SNESDestroy(&ts->snes);
2779: ts->snes = snes;
2781: SNESSetFunction(ts->snes,NULL,SNESTSFormFunction,ts);
2782: SNESGetJacobian(ts->snes,NULL,NULL,&func,NULL);
2783: if (func == SNESTSFormJacobian) {
2784: SNESSetJacobian(ts->snes,NULL,NULL,SNESTSFormJacobian,ts);
2785: }
2786: return(0);
2787: }
2789: /*@
2790: TSGetKSP - Returns the KSP (linear solver) associated with
2791: a TS (timestepper) context.
2793: Not Collective, but KSP is parallel if TS is parallel
2795: Input Parameter:
2796: . ts - the TS context obtained from TSCreate()
2798: Output Parameter:
2799: . ksp - the nonlinear solver context
2801: Notes:
2802: The user can then directly manipulate the KSP context to set various
2803: options, etc. Likewise, the user can then extract and manipulate the
2804: KSP and PC contexts as well.
2806: TSGetKSP() does not work for integrators that do not use KSP;
2807: in this case TSGetKSP() returns NULL in ksp.
2809: Level: beginner
2811: @*/
2812: PetscErrorCode TSGetKSP(TS ts,KSP *ksp)
2813: {
2815: SNES snes;
2820: if (!((PetscObject)ts)->type_name) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"KSP is not created yet. Call TSSetType() first");
2821: if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()");
2822: TSGetSNES(ts,&snes);
2823: SNESGetKSP(snes,ksp);
2824: return(0);
2825: }
2827: /* ----------- Routines to set solver parameters ---------- */
2829: /*@
2830: TSSetMaxSteps - Sets the maximum number of steps to use.
2832: Logically Collective on TS
2834: Input Parameters:
2835: + ts - the TS context obtained from TSCreate()
2836: - maxsteps - maximum number of steps to use
2838: Options Database Keys:
2839: . -ts_max_steps <maxsteps> - Sets maxsteps
2841: Notes:
2842: The default maximum number of steps is 5000
2844: Level: intermediate
2846: .seealso: TSGetMaxSteps(), TSSetMaxTime(), TSSetExactFinalTime()
2847: @*/
2848: PetscErrorCode TSSetMaxSteps(TS ts,PetscInt maxsteps)
2849: {
2853: if (maxsteps < 0 ) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of steps must be non-negative");
2854: ts->max_steps = maxsteps;
2855: return(0);
2856: }
2858: /*@
2859: TSGetMaxSteps - Gets the maximum number of steps to use.
2861: Not Collective
2863: Input Parameters:
2864: . ts - the TS context obtained from TSCreate()
2866: Output Parameter:
2867: . maxsteps - maximum number of steps to use
2869: Level: advanced
2871: .seealso: TSSetMaxSteps(), TSGetMaxTime(), TSSetMaxTime()
2872: @*/
2873: PetscErrorCode TSGetMaxSteps(TS ts,PetscInt *maxsteps)
2874: {
2878: *maxsteps = ts->max_steps;
2879: return(0);
2880: }
2882: /*@
2883: TSSetMaxTime - Sets the maximum (or final) time for timestepping.
2885: Logically Collective on TS
2887: Input Parameters:
2888: + ts - the TS context obtained from TSCreate()
2889: - maxtime - final time to step to
2891: Options Database Keys:
2892: . -ts_max_time <maxtime> - Sets maxtime
2894: Notes:
2895: The default maximum time is 5.0
2897: Level: intermediate
2899: .seealso: TSGetMaxTime(), TSSetMaxSteps(), TSSetExactFinalTime()
2900: @*/
2901: PetscErrorCode TSSetMaxTime(TS ts,PetscReal maxtime)
2902: {
2906: ts->max_time = maxtime;
2907: return(0);
2908: }
2910: /*@
2911: TSGetMaxTime - Gets the maximum (or final) time for timestepping.
2913: Not Collective
2915: Input Parameters:
2916: . ts - the TS context obtained from TSCreate()
2918: Output Parameter:
2919: . maxtime - final time to step to
2921: Level: advanced
2923: .seealso: TSSetMaxTime(), TSGetMaxSteps(), TSSetMaxSteps()
2924: @*/
2925: PetscErrorCode TSGetMaxTime(TS ts,PetscReal *maxtime)
2926: {
2930: *maxtime = ts->max_time;
2931: return(0);
2932: }
2934: /*@
2935: TSSetInitialTimeStep - Deprecated, use TSSetTime() and TSSetTimeStep().
2937: Level: deprecated
2939: @*/
2940: PetscErrorCode TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step)
2941: {
2945: TSSetTime(ts,initial_time);
2946: TSSetTimeStep(ts,time_step);
2947: return(0);
2948: }
2950: /*@
2951: TSGetDuration - Deprecated, use TSGetMaxSteps() and TSGetMaxTime().
2953: Level: deprecated
2955: @*/
2956: PetscErrorCode TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime)
2957: {
2960: if (maxsteps) {
2962: *maxsteps = ts->max_steps;
2963: }
2964: if (maxtime) {
2966: *maxtime = ts->max_time;
2967: }
2968: return(0);
2969: }
2971: /*@
2972: TSSetDuration - Deprecated, use TSSetMaxSteps() and TSSetMaxTime().
2974: Level: deprecated
2976: @*/
2977: PetscErrorCode TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime)
2978: {
2983: if (maxsteps >= 0) ts->max_steps = maxsteps;
2984: if (maxtime != PETSC_DEFAULT) ts->max_time = maxtime;
2985: return(0);
2986: }
2988: /*@
2989: TSGetTimeStepNumber - Deprecated, use TSGetStepNumber().
2991: Level: deprecated
2993: @*/
2994: PetscErrorCode TSGetTimeStepNumber(TS ts,PetscInt *steps) { return TSGetStepNumber(ts,steps); }
2996: /*@
2997: TSGetTotalSteps - Deprecated, use TSGetStepNumber().
2999: Level: deprecated
3001: @*/
3002: PetscErrorCode TSGetTotalSteps(TS ts,PetscInt *steps) { return TSGetStepNumber(ts,steps); }
3004: /*@
3005: TSSetSolution - Sets the initial solution vector
3006: for use by the TS routines.
3008: Logically Collective on TS
3010: Input Parameters:
3011: + ts - the TS context obtained from TSCreate()
3012: - u - the solution vector
3014: Level: beginner
3016: .seealso: TSSetSolutionFunction(), TSGetSolution(), TSCreate()
3017: @*/
3018: PetscErrorCode TSSetSolution(TS ts,Vec u)
3019: {
3021: DM dm;
3026: PetscObjectReference((PetscObject)u);
3027: VecDestroy(&ts->vec_sol);
3028: ts->vec_sol = u;
3030: TSGetDM(ts,&dm);
3031: DMShellSetGlobalVector(dm,u);
3032: return(0);
3033: }
3035: /*@C
3036: TSSetPreStep - Sets the general-purpose function
3037: called once at the beginning of each time step.
3039: Logically Collective on TS
3041: Input Parameters:
3042: + ts - The TS context obtained from TSCreate()
3043: - func - The function
3045: Calling sequence of func:
3046: . PetscErrorCode func (TS ts);
3048: Level: intermediate
3050: .seealso: TSSetPreStage(), TSSetPostStage(), TSSetPostStep(), TSStep(), TSRestartStep()
3051: @*/
3052: PetscErrorCode TSSetPreStep(TS ts, PetscErrorCode (*func)(TS))
3053: {
3056: ts->prestep = func;
3057: return(0);
3058: }
3060: /*@
3061: TSPreStep - Runs the user-defined pre-step function.
3063: Collective on TS
3065: Input Parameters:
3066: . ts - The TS context obtained from TSCreate()
3068: Notes:
3069: TSPreStep() is typically used within time stepping implementations,
3070: so most users would not generally call this routine themselves.
3072: Level: developer
3074: .seealso: TSSetPreStep(), TSPreStage(), TSPostStage(), TSPostStep()
3075: @*/
3076: PetscErrorCode TSPreStep(TS ts)
3077: {
3082: if (ts->prestep) {
3083: Vec U;
3084: PetscObjectState sprev,spost;
3086: TSGetSolution(ts,&U);
3087: PetscObjectStateGet((PetscObject)U,&sprev);
3088: PetscStackCallStandard((*ts->prestep),(ts));
3089: PetscObjectStateGet((PetscObject)U,&spost);
3090: if (sprev != spost) {TSRestartStep(ts);}
3091: }
3092: return(0);
3093: }
3095: /*@C
3096: TSSetPreStage - Sets the general-purpose function
3097: called once at the beginning of each stage.
3099: Logically Collective on TS
3101: Input Parameters:
3102: + ts - The TS context obtained from TSCreate()
3103: - func - The function
3105: Calling sequence of func:
3106: . PetscErrorCode func(TS ts, PetscReal stagetime);
3108: Level: intermediate
3110: Note:
3111: There may be several stages per time step. If the solve for a given stage fails, the step may be rejected and retried.
3112: The time step number being computed can be queried using TSGetStepNumber() and the total size of the step being
3113: attempted can be obtained using TSGetTimeStep(). The time at the start of the step is available via TSGetTime().
3115: .seealso: TSSetPostStage(), TSSetPreStep(), TSSetPostStep(), TSGetApplicationContext()
3116: @*/
3117: PetscErrorCode TSSetPreStage(TS ts, PetscErrorCode (*func)(TS,PetscReal))
3118: {
3121: ts->prestage = func;
3122: return(0);
3123: }
3125: /*@C
3126: TSSetPostStage - Sets the general-purpose function
3127: called once at the end of each stage.
3129: Logically Collective on TS
3131: Input Parameters:
3132: + ts - The TS context obtained from TSCreate()
3133: - func - The function
3135: Calling sequence of func:
3136: . PetscErrorCode func(TS ts, PetscReal stagetime, PetscInt stageindex, Vec* Y);
3138: Level: intermediate
3140: Note:
3141: There may be several stages per time step. If the solve for a given stage fails, the step may be rejected and retried.
3142: The time step number being computed can be queried using TSGetStepNumber() and the total size of the step being
3143: attempted can be obtained using TSGetTimeStep(). The time at the start of the step is available via TSGetTime().
3145: .seealso: TSSetPreStage(), TSSetPreStep(), TSSetPostStep(), TSGetApplicationContext()
3146: @*/
3147: PetscErrorCode TSSetPostStage(TS ts, PetscErrorCode (*func)(TS,PetscReal,PetscInt,Vec*))
3148: {
3151: ts->poststage = func;
3152: return(0);
3153: }
3155: /*@C
3156: TSSetPostEvaluate - Sets the general-purpose function
3157: called once at the end of each step evaluation.
3159: Logically Collective on TS
3161: Input Parameters:
3162: + ts - The TS context obtained from TSCreate()
3163: - func - The function
3165: Calling sequence of func:
3166: . PetscErrorCode func(TS ts);
3168: Level: intermediate
3170: Note:
3171: Semantically, TSSetPostEvaluate() differs from TSSetPostStep() since the function it sets is called before event-handling
3172: thus guaranteeing the same solution (computed by the time-stepper) will be passed to it. On the other hand, TSPostStep()
3173: may be passed a different solution, possibly changed by the event handler. TSPostEvaluate() is called after the next step
3174: solution is evaluated allowing to modify it, if need be. The solution can be obtained with TSGetSolution(), the time step
3175: with TSGetTimeStep(), and the time at the start of the step is available via TSGetTime()
3177: .seealso: TSSetPreStage(), TSSetPreStep(), TSSetPostStep(), TSGetApplicationContext()
3178: @*/
3179: PetscErrorCode TSSetPostEvaluate(TS ts, PetscErrorCode (*func)(TS))
3180: {
3183: ts->postevaluate = func;
3184: return(0);
3185: }
3187: /*@
3188: TSPreStage - Runs the user-defined pre-stage function set using TSSetPreStage()
3190: Collective on TS
3192: Input Parameters:
3193: . ts - The TS context obtained from TSCreate()
3194: stagetime - The absolute time of the current stage
3196: Notes:
3197: TSPreStage() is typically used within time stepping implementations,
3198: most users would not generally call this routine themselves.
3200: Level: developer
3202: .seealso: TSPostStage(), TSSetPreStep(), TSPreStep(), TSPostStep()
3203: @*/
3204: PetscErrorCode TSPreStage(TS ts, PetscReal stagetime)
3205: {
3208: if (ts->prestage) {
3209: PetscStackCallStandard((*ts->prestage),(ts,stagetime));
3210: }
3211: return(0);
3212: }
3214: /*@
3215: TSPostStage - Runs the user-defined post-stage function set using TSSetPostStage()
3217: Collective on TS
3219: Input Parameters:
3220: . ts - The TS context obtained from TSCreate()
3221: stagetime - The absolute time of the current stage
3222: stageindex - Stage number
3223: Y - Array of vectors (of size = total number
3224: of stages) with the stage solutions
3226: Notes:
3227: TSPostStage() is typically used within time stepping implementations,
3228: most users would not generally call this routine themselves.
3230: Level: developer
3232: .seealso: TSPreStage(), TSSetPreStep(), TSPreStep(), TSPostStep()
3233: @*/
3234: PetscErrorCode TSPostStage(TS ts, PetscReal stagetime, PetscInt stageindex, Vec *Y)
3235: {
3238: if (ts->poststage) {
3239: PetscStackCallStandard((*ts->poststage),(ts,stagetime,stageindex,Y));
3240: }
3241: return(0);
3242: }
3244: /*@
3245: TSPostEvaluate - Runs the user-defined post-evaluate function set using TSSetPostEvaluate()
3247: Collective on TS
3249: Input Parameters:
3250: . ts - The TS context obtained from TSCreate()
3252: Notes:
3253: TSPostEvaluate() is typically used within time stepping implementations,
3254: most users would not generally call this routine themselves.
3256: Level: developer
3258: .seealso: TSSetPostEvaluate(), TSSetPreStep(), TSPreStep(), TSPostStep()
3259: @*/
3260: PetscErrorCode TSPostEvaluate(TS ts)
3261: {
3266: if (ts->postevaluate) {
3267: Vec U;
3268: PetscObjectState sprev,spost;
3270: TSGetSolution(ts,&U);
3271: PetscObjectStateGet((PetscObject)U,&sprev);
3272: PetscStackCallStandard((*ts->postevaluate),(ts));
3273: PetscObjectStateGet((PetscObject)U,&spost);
3274: if (sprev != spost) {TSRestartStep(ts);}
3275: }
3276: return(0);
3277: }
3279: /*@C
3280: TSSetPostStep - Sets the general-purpose function
3281: called once at the end of each time step.
3283: Logically Collective on TS
3285: Input Parameters:
3286: + ts - The TS context obtained from TSCreate()
3287: - func - The function
3289: Calling sequence of func:
3290: $ func (TS ts);
3292: Notes:
3293: The function set by TSSetPostStep() is called after each successful step. The solution vector X
3294: obtained by TSGetSolution() may be different than that computed at the step end if the event handler
3295: locates an event and TSPostEvent() modifies it. Use TSSetPostEvaluate() if an unmodified solution is needed instead.
3297: Level: intermediate
3299: .seealso: TSSetPreStep(), TSSetPreStage(), TSSetPostEvaluate(), TSGetTimeStep(), TSGetStepNumber(), TSGetTime(), TSRestartStep()
3300: @*/
3301: PetscErrorCode TSSetPostStep(TS ts, PetscErrorCode (*func)(TS))
3302: {
3305: ts->poststep = func;
3306: return(0);
3307: }
3309: /*@
3310: TSPostStep - Runs the user-defined post-step function.
3312: Collective on TS
3314: Input Parameters:
3315: . ts - The TS context obtained from TSCreate()
3317: Notes:
3318: TSPostStep() is typically used within time stepping implementations,
3319: so most users would not generally call this routine themselves.
3321: Level: developer
3323: @*/
3324: PetscErrorCode TSPostStep(TS ts)
3325: {
3330: if (ts->poststep) {
3331: Vec U;
3332: PetscObjectState sprev,spost;
3334: TSGetSolution(ts,&U);
3335: PetscObjectStateGet((PetscObject)U,&sprev);
3336: PetscStackCallStandard((*ts->poststep),(ts));
3337: PetscObjectStateGet((PetscObject)U,&spost);
3338: if (sprev != spost) {TSRestartStep(ts);}
3339: }
3340: return(0);
3341: }
3343: /* ------------ Routines to set performance monitoring options ----------- */
3345: /*@C
3346: TSMonitorSet - Sets an ADDITIONAL function that is to be used at every
3347: timestep to display the iteration's progress.
3349: Logically Collective on TS
3351: Input Parameters:
3352: + ts - the TS context obtained from TSCreate()
3353: . monitor - monitoring routine
3354: . mctx - [optional] user-defined context for private data for the
3355: monitor routine (use NULL if no context is desired)
3356: - monitordestroy - [optional] routine that frees monitor context
3357: (may be NULL)
3359: Calling sequence of monitor:
3360: $ PetscErrorCode monitor(TS ts,PetscInt steps,PetscReal time,Vec u,void *mctx)
3362: + ts - the TS context
3363: . steps - iteration number (after the final time step the monitor routine may be called with a step of -1, this indicates the solution has been interpolated to this time)
3364: . time - current time
3365: . u - current iterate
3366: - mctx - [optional] monitoring context
3368: Notes:
3369: This routine adds an additional monitor to the list of monitors that
3370: already has been loaded.
3372: Fortran Notes:
3373: Only a single monitor function can be set for each TS object
3375: Level: intermediate
3377: .seealso: TSMonitorDefault(), TSMonitorCancel()
3378: @*/
3379: PetscErrorCode TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void**))
3380: {
3382: PetscInt i;
3383: PetscBool identical;
3387: for (i=0; i<ts->numbermonitors;i++) {
3388: PetscMonitorCompare((PetscErrorCode (*)(void))monitor,mctx,mdestroy,(PetscErrorCode (*)(void))ts->monitor[i],ts->monitorcontext[i],ts->monitordestroy[i],&identical);
3389: if (identical) return(0);
3390: }
3391: if (ts->numbermonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
3392: ts->monitor[ts->numbermonitors] = monitor;
3393: ts->monitordestroy[ts->numbermonitors] = mdestroy;
3394: ts->monitorcontext[ts->numbermonitors++] = (void*)mctx;
3395: return(0);
3396: }
3398: /*@C
3399: TSMonitorCancel - Clears all the monitors that have been set on a time-step object.
3401: Logically Collective on TS
3403: Input Parameters:
3404: . ts - the TS context obtained from TSCreate()
3406: Notes:
3407: There is no way to remove a single, specific monitor.
3409: Level: intermediate
3411: .seealso: TSMonitorDefault(), TSMonitorSet()
3412: @*/
3413: PetscErrorCode TSMonitorCancel(TS ts)
3414: {
3416: PetscInt i;
3420: for (i=0; i<ts->numbermonitors; i++) {
3421: if (ts->monitordestroy[i]) {
3422: (*ts->monitordestroy[i])(&ts->monitorcontext[i]);
3423: }
3424: }
3425: ts->numbermonitors = 0;
3426: return(0);
3427: }
3429: /*@C
3430: TSMonitorDefault - The Default monitor, prints the timestep and time for each step
3432: Level: intermediate
3434: .seealso: TSMonitorSet()
3435: @*/
3436: PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscViewerAndFormat *vf)
3437: {
3439: PetscViewer viewer = vf->viewer;
3440: PetscBool iascii,ibinary;
3444: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
3445: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&ibinary);
3446: PetscViewerPushFormat(viewer,vf->format);
3447: if (iascii) {
3448: PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);
3449: if (step == -1){ /* this indicates it is an interpolated solution */
3450: PetscViewerASCIIPrintf(viewer,"Interpolated solution at time %g between steps %D and %D\n",(double)ptime,ts->steps-1,ts->steps);
3451: } else {
3452: PetscViewerASCIIPrintf(viewer,"%D TS dt %g time %g%s",step,(double)ts->time_step,(double)ptime,ts->steprollback ? " (r)\n" : "\n");
3453: }
3454: PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);
3455: } else if (ibinary) {
3456: PetscMPIInt rank;
3457: MPI_Comm_rank(PetscObjectComm((PetscObject)viewer),&rank);
3458: if (!rank) {
3459: PetscBool skipHeader;
3460: PetscInt classid = REAL_FILE_CLASSID;
3462: PetscViewerBinaryGetSkipHeader(viewer,&skipHeader);
3463: if (!skipHeader) {
3464: PetscViewerBinaryWrite(viewer,&classid,1,PETSC_INT,PETSC_FALSE);
3465: }
3466: PetscRealView(1,&ptime,viewer);
3467: } else {
3468: PetscRealView(0,&ptime,viewer);
3469: }
3470: }
3471: PetscViewerPopFormat(viewer);
3472: return(0);
3473: }
3475: /*@C
3476: TSMonitorExtreme - Prints the extreme values of the solution at each timestep
3478: Level: intermediate
3480: .seealso: TSMonitorSet()
3481: @*/
3482: PetscErrorCode TSMonitorExtreme(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscViewerAndFormat *vf)
3483: {
3485: PetscViewer viewer = vf->viewer;
3486: PetscBool iascii;
3487: PetscReal max,min;
3489:
3492: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
3493: PetscViewerPushFormat(viewer,vf->format);
3494: if (iascii) {
3495: VecMax(v,NULL,&max);
3496: VecMin(v,NULL,&min);
3497: PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);
3498: PetscViewerASCIIPrintf(viewer,"%D TS dt %g time %g%s max %g min %g\n",step,(double)ts->time_step,(double)ptime,ts->steprollback ? " (r)" : "",(double)max,(double)min);
3499: PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);
3500: }
3501: PetscViewerPopFormat(viewer);
3502: return(0);
3503: }
3505: /*@
3506: TSInterpolate - Interpolate the solution computed during the previous step to an arbitrary location in the interval
3508: Collective on TS
3510: Input Argument:
3511: + ts - time stepping context
3512: - t - time to interpolate to
3514: Output Argument:
3515: . U - state at given time
3517: Level: intermediate
3519: Developer Notes:
3520: TSInterpolate() and the storing of previous steps/stages should be generalized to support delay differential equations and continuous adjoints.
3522: .seealso: TSSetExactFinalTime(), TSSolve()
3523: @*/
3524: PetscErrorCode TSInterpolate(TS ts,PetscReal t,Vec U)
3525: {
3531: if (t < ts->ptime_prev || t > ts->ptime) SETERRQ3(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Requested time %g not in last time steps [%g,%g]",t,(double)ts->ptime_prev,(double)ts->ptime);
3532: if (!ts->ops->interpolate) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"%s does not provide interpolation",((PetscObject)ts)->type_name);
3533: (*ts->ops->interpolate)(ts,t,U);
3534: return(0);
3535: }
3537: /*@
3538: TSStep - Steps one time step
3540: Collective on TS
3542: Input Parameter:
3543: . ts - the TS context obtained from TSCreate()
3545: Level: developer
3547: Notes:
3548: The public interface for the ODE/DAE solvers is TSSolve(), you should almost for sure be using that routine and not this routine.
3550: The hook set using TSSetPreStep() is called before each attempt to take the step. In general, the time step size may
3551: be changed due to adaptive error controller or solve failures. Note that steps may contain multiple stages.
3553: This may over-step the final time provided in TSSetMaxTime() depending on the time-step used. TSSolve() interpolates to exactly the
3554: time provided in TSSetMaxTime(). One can use TSInterpolate() to determine an interpolated solution within the final timestep.
3556: .seealso: TSCreate(), TSSetUp(), TSDestroy(), TSSolve(), TSSetPreStep(), TSSetPreStage(), TSSetPostStage(), TSInterpolate()
3557: @*/
3558: PetscErrorCode TSStep(TS ts)
3559: {
3560: PetscErrorCode ierr;
3561: static PetscBool cite = PETSC_FALSE;
3562: PetscReal ptime;
3566: PetscCitationsRegister("@techreport{tspaper,\n"
3567: " title = {{PETSc/TS}: A Modern Scalable {DAE/ODE} Solver Library},\n"
3568: " author = {Shrirang Abhyankar and Jed Brown and Emil Constantinescu and Debojyoti Ghosh and Barry F. Smith},\n"
3569: " type = {Preprint},\n"
3570: " number = {ANL/MCS-P5061-0114},\n"
3571: " institution = {Argonne National Laboratory},\n"
3572: " year = {2014}\n}\n",&cite);
3574: TSSetUp(ts);
3575: TSTrajectorySetUp(ts->trajectory,ts);
3577: if (ts->max_time >= PETSC_MAX_REAL && ts->max_steps == PETSC_MAX_INT) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"You must call TSSetMaxTime() or TSSetMaxSteps(), or use -ts_max_time <time> or -ts_max_steps <steps>");
3578: if (ts->exact_final_time == TS_EXACTFINALTIME_UNSPECIFIED) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"You must call TSSetExactFinalTime() or use -ts_exact_final_time <stepover,interpolate,matchstep> before calling TSStep()");
3579: if (ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP && !ts->adapt) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Since TS is not adaptive you cannot use TS_EXACTFINALTIME_MATCHSTEP, suggest TS_EXACTFINALTIME_INTERPOLATE");
3581: if (!ts->steps) ts->ptime_prev = ts->ptime;
3582: ptime = ts->ptime; ts->ptime_prev_rollback = ts->ptime_prev;
3583: ts->reason = TS_CONVERGED_ITERATING;
3584: if (!ts->ops->step) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSStep not implemented for type '%s'",((PetscObject)ts)->type_name);
3585: PetscLogEventBegin(TS_Step,ts,0,0,0);
3586: (*ts->ops->step)(ts);
3587: PetscLogEventEnd(TS_Step,ts,0,0,0);
3588: ts->ptime_prev = ptime;
3589: ts->steps++;
3590: ts->steprollback = PETSC_FALSE;
3591: ts->steprestart = PETSC_FALSE;
3593: if (ts->reason < 0) {
3594: if (ts->errorifstepfailed) {
3595: if (ts->reason == TS_DIVERGED_NONLINEAR_SOLVE) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s, increase -ts_max_snes_failures or make negative to attempt recovery",TSConvergedReasons[ts->reason]);
3596: else SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s",TSConvergedReasons[ts->reason]);
3597: }
3598: } else if (!ts->reason) {
3599: if (ts->steps >= ts->max_steps) ts->reason = TS_CONVERGED_ITS;
3600: else if (ts->ptime >= ts->max_time) ts->reason = TS_CONVERGED_TIME;
3601: }
3602: return(0);
3603: }
3605: /*@
3606: TSEvaluateWLTE - Evaluate the weighted local truncation error norm
3607: at the end of a time step with a given order of accuracy.
3609: Collective on TS
3611: Input Arguments:
3612: + ts - time stepping context
3613: . wnormtype - norm type, either NORM_2 or NORM_INFINITY
3614: - order - optional, desired order for the error evaluation or PETSC_DECIDE
3616: Output Arguments:
3617: + order - optional, the actual order of the error evaluation
3618: - wlte - the weighted local truncation error norm
3620: Level: advanced
3622: Notes:
3623: If the timestepper cannot evaluate the error in a particular step
3624: (eg. in the first step or restart steps after event handling),
3625: this routine returns wlte=-1.0 .
3627: .seealso: TSStep(), TSAdapt, TSErrorWeightedNorm()
3628: @*/
3629: PetscErrorCode TSEvaluateWLTE(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte)
3630: {
3640: if (wnormtype != NORM_2 && wnormtype != NORM_INFINITY) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"No support for norm type %s",NormTypes[wnormtype]);
3641: if (!ts->ops->evaluatewlte) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSEvaluateWLTE not implemented for type '%s'",((PetscObject)ts)->type_name);
3642: (*ts->ops->evaluatewlte)(ts,wnormtype,order,wlte);
3643: return(0);
3644: }
3646: /*@
3647: TSEvaluateStep - Evaluate the solution at the end of a time step with a given order of accuracy.
3649: Collective on TS
3651: Input Arguments:
3652: + ts - time stepping context
3653: . order - desired order of accuracy
3654: - done - whether the step was evaluated at this order (pass NULL to generate an error if not available)
3656: Output Arguments:
3657: . U - state at the end of the current step
3659: Level: advanced
3661: Notes:
3662: This function cannot be called until all stages have been evaluated.
3663: 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.
3665: .seealso: TSStep(), TSAdapt
3666: @*/
3667: PetscErrorCode TSEvaluateStep(TS ts,PetscInt order,Vec U,PetscBool *done)
3668: {
3675: if (!ts->ops->evaluatestep) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSEvaluateStep not implemented for type '%s'",((PetscObject)ts)->type_name);
3676: (*ts->ops->evaluatestep)(ts,order,U,done);
3677: return(0);
3678: }
3680: /*@
3681: TSSolve - Steps the requested number of timesteps.
3683: Collective on TS
3685: Input Parameter:
3686: + ts - the TS context obtained from TSCreate()
3687: - u - the solution vector (can be null if TSSetSolution() was used and TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP) was not used,
3688: otherwise must contain the initial conditions and will contain the solution at the final requested time
3690: Level: beginner
3692: Notes:
3693: The final time returned by this function may be different from the time of the internally
3694: held state accessible by TSGetSolution() and TSGetTime() because the method may have
3695: stepped over the final time.
3697: .seealso: TSCreate(), TSSetSolution(), TSStep(), TSGetTime(), TSGetSolveTime()
3698: @*/
3699: PetscErrorCode TSSolve(TS ts,Vec u)
3700: {
3701: Vec solution;
3702: PetscErrorCode ierr;
3708: if (ts->exact_final_time == TS_EXACTFINALTIME_INTERPOLATE && u) { /* Need ts->vec_sol to be distinct so it is not overwritten when we interpolate at the end */
3709: if (!ts->vec_sol || u == ts->vec_sol) {
3710: VecDuplicate(u,&solution);
3711: TSSetSolution(ts,solution);
3712: VecDestroy(&solution); /* grant ownership */
3713: }
3714: VecCopy(u,ts->vec_sol);
3715: if (ts->forward_solve) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Sensitivity analysis does not support the mode TS_EXACTFINALTIME_INTERPOLATE");
3716: } else if (u) {
3717: TSSetSolution(ts,u);
3718: }
3719: TSSetUp(ts);
3720: TSTrajectorySetUp(ts->trajectory,ts);
3722: if (ts->max_time >= PETSC_MAX_REAL && ts->max_steps == PETSC_MAX_INT) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"You must call TSSetMaxTime() or TSSetMaxSteps(), or use -ts_max_time <time> or -ts_max_steps <steps>");
3723: if (ts->exact_final_time == TS_EXACTFINALTIME_UNSPECIFIED) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"You must call TSSetExactFinalTime() or use -ts_exact_final_time <stepover,interpolate,matchstep> before calling TSSolve()");
3724: if (ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP && !ts->adapt) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Since TS is not adaptive you cannot use TS_EXACTFINALTIME_MATCHSTEP, suggest TS_EXACTFINALTIME_INTERPOLATE");
3726: if (ts->forward_solve) {
3727: TSForwardSetUp(ts);
3728: }
3730: /* reset number of steps only when the step is not restarted. ARKIMEX
3731: restarts the step after an event. Resetting these counters in such case causes
3732: TSTrajectory to incorrectly save the output files
3733: */
3734: /* reset time step and iteration counters */
3735: if (!ts->steps) {
3736: ts->ksp_its = 0;
3737: ts->snes_its = 0;
3738: ts->num_snes_failures = 0;
3739: ts->reject = 0;
3740: ts->steprestart = PETSC_TRUE;
3741: ts->steprollback = PETSC_FALSE;
3742: ts->rhsjacobian.time = PETSC_MIN_REAL;
3743: }
3744: if (ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP && ts->ptime < ts->max_time && ts->ptime + ts->time_step > ts->max_time) ts->time_step = ts->max_time - ts->ptime;
3745: ts->reason = TS_CONVERGED_ITERATING;
3747: TSViewFromOptions(ts,NULL,"-ts_view_pre");
3749: if (ts->ops->solve) { /* This private interface is transitional and should be removed when all implementations are updated. */
3750: (*ts->ops->solve)(ts);
3751: if (u) {VecCopy(ts->vec_sol,u);}
3752: ts->solvetime = ts->ptime;
3753: solution = ts->vec_sol;
3754: } else { /* Step the requested number of timesteps. */
3755: if (ts->steps >= ts->max_steps) ts->reason = TS_CONVERGED_ITS;
3756: else if (ts->ptime >= ts->max_time) ts->reason = TS_CONVERGED_TIME;
3758: if (!ts->steps) {
3759: TSTrajectorySet(ts->trajectory,ts,ts->steps,ts->ptime,ts->vec_sol);
3760: TSEventInitialize(ts->event,ts,ts->ptime,ts->vec_sol);
3761: }
3763: while (!ts->reason) {
3764: TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);
3765: if (!ts->steprollback) {
3766: TSPreStep(ts);
3767: }
3768: TSStep(ts);
3769: if (ts->testjacobian) {
3770: TSRHSJacobianTest(ts,NULL);
3771: }
3772: if (ts->testjacobiantranspose) {
3773: TSRHSJacobianTestTranspose(ts,NULL);
3774: }
3775: if (ts->quadraturets && ts->costintegralfwd) { /* Must evaluate the cost integral before event is handled. The cost integral value can also be rolled back. */
3776: TSForwardCostIntegral(ts);
3777: }
3778: if (ts->forward_solve) { /* compute forward sensitivities before event handling because postevent() may change RHS and jump conditions may have to be applied */
3779: TSForwardStep(ts);
3780: }
3781: TSPostEvaluate(ts);
3782: TSEventHandler(ts); /* The right-hand side may be changed due to event. Be careful with Any computation using the RHS information after this point. */
3783: if (ts->steprollback) {
3784: TSPostEvaluate(ts);
3785: }
3786: if (!ts->steprollback) {
3787: TSTrajectorySet(ts->trajectory,ts,ts->steps,ts->ptime,ts->vec_sol);
3788: TSPostStep(ts);
3789: }
3790: }
3791: TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);
3793: if (ts->exact_final_time == TS_EXACTFINALTIME_INTERPOLATE && ts->ptime > ts->max_time) {
3794: TSInterpolate(ts,ts->max_time,u);
3795: ts->solvetime = ts->max_time;
3796: solution = u;
3797: TSMonitor(ts,-1,ts->solvetime,solution);
3798: } else {
3799: if (u) {VecCopy(ts->vec_sol,u);}
3800: ts->solvetime = ts->ptime;
3801: solution = ts->vec_sol;
3802: }
3803: }
3805: TSViewFromOptions(ts,NULL,"-ts_view");
3806: VecViewFromOptions(solution,NULL,"-ts_view_solution");
3807: PetscObjectSAWsBlock((PetscObject)ts);
3808: if (ts->adjoint_solve) {
3809: TSAdjointSolve(ts);
3810: }
3811: return(0);
3812: }
3814: /*@C
3815: TSMonitor - Runs all user-provided monitor routines set using TSMonitorSet()
3817: Collective on TS
3819: Input Parameters:
3820: + ts - time stepping context obtained from TSCreate()
3821: . step - step number that has just completed
3822: . ptime - model time of the state
3823: - u - state at the current model time
3825: Notes:
3826: TSMonitor() is typically used automatically within the time stepping implementations.
3827: Users would almost never call this routine directly.
3829: A step of -1 indicates that the monitor is being called on a solution obtained by interpolating from computed solutions
3831: Level: developer
3833: @*/
3834: PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec u)
3835: {
3836: DM dm;
3837: PetscInt i,n = ts->numbermonitors;
3844: TSGetDM(ts,&dm);
3845: DMSetOutputSequenceNumber(dm,step,ptime);
3847: VecLockReadPush(u);
3848: for (i=0; i<n; i++) {
3849: (*ts->monitor[i])(ts,step,ptime,u,ts->monitorcontext[i]);
3850: }
3851: VecLockReadPop(u);
3852: return(0);
3853: }
3855: /* ------------------------------------------------------------------------*/
3856: /*@C
3857: TSMonitorLGCtxCreate - Creates a TSMonitorLGCtx context for use with
3858: TS to monitor the solution process graphically in various ways
3860: Collective on TS
3862: Input Parameters:
3863: + host - the X display to open, or null for the local machine
3864: . label - the title to put in the title bar
3865: . x, y - the screen coordinates of the upper left coordinate of the window
3866: . m, n - the screen width and height in pixels
3867: - howoften - if positive then determines the frequency of the plotting, if -1 then only at the final time
3869: Output Parameter:
3870: . ctx - the context
3872: Options Database Key:
3873: + -ts_monitor_lg_timestep - automatically sets line graph monitor
3874: + -ts_monitor_lg_timestep_log - automatically sets line graph monitor
3875: . -ts_monitor_lg_solution - monitor the solution (or certain values of the solution by calling TSMonitorLGSetDisplayVariables() or TSMonitorLGCtxSetDisplayVariables())
3876: . -ts_monitor_lg_error - monitor the error
3877: . -ts_monitor_lg_ksp_iterations - monitor the number of KSP iterations needed for each timestep
3878: . -ts_monitor_lg_snes_iterations - monitor the number of SNES iterations needed for each timestep
3879: - -lg_use_markers <true,false> - mark the data points (at each time step) on the plot; default is true
3881: Notes:
3882: Use TSMonitorLGCtxDestroy() to destroy.
3884: One can provide a function that transforms the solution before plotting it with TSMonitorLGCtxSetTransform() or TSMonitorLGSetTransform()
3886: Many of the functions that control the monitoring have two forms: TSMonitorLGSet/GetXXXX() and TSMonitorLGCtxSet/GetXXXX() the first take a TS object as the
3887: first argument (if that TS object does not have a TSMonitorLGCtx associated with it the function call is ignored) and the second takes a TSMonitorLGCtx object
3888: as the first argument.
3890: One can control the names displayed for each solution or error variable with TSMonitorLGCtxSetVariableNames() or TSMonitorLGSetVariableNames()
3892: Level: intermediate
3894: .seealso: TSMonitorLGTimeStep(), TSMonitorSet(), TSMonitorLGSolution(), TSMonitorLGError(), TSMonitorDefault(), VecView(),
3895: TSMonitorLGCtxCreate(), TSMonitorLGCtxSetVariableNames(), TSMonitorLGCtxGetVariableNames(),
3896: TSMonitorLGSetVariableNames(), TSMonitorLGGetVariableNames(), TSMonitorLGSetDisplayVariables(), TSMonitorLGCtxSetDisplayVariables(),
3897: TSMonitorLGCtxSetTransform(), TSMonitorLGSetTransform(), TSMonitorLGError(), TSMonitorLGSNESIterations(), TSMonitorLGKSPIterations(),
3898: TSMonitorEnvelopeCtxCreate(), TSMonitorEnvelopeGetBounds(), TSMonitorEnvelopeCtxDestroy(), TSMonitorEnvelop()
3900: @*/
3901: PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscInt howoften,TSMonitorLGCtx *ctx)
3902: {
3903: PetscDraw draw;
3907: PetscNew(ctx);
3908: PetscDrawCreate(comm,host,label,x,y,m,n,&draw);
3909: PetscDrawSetFromOptions(draw);
3910: PetscDrawLGCreate(draw,1,&(*ctx)->lg);
3911: PetscDrawLGSetFromOptions((*ctx)->lg);
3912: PetscDrawDestroy(&draw);
3913: (*ctx)->howoften = howoften;
3914: return(0);
3915: }
3917: PetscErrorCode TSMonitorLGTimeStep(TS ts,PetscInt step,PetscReal ptime,Vec v,void *monctx)
3918: {
3919: TSMonitorLGCtx ctx = (TSMonitorLGCtx) monctx;
3920: PetscReal x = ptime,y;
3924: if (step < 0) return(0); /* -1 indicates an interpolated solution */
3925: if (!step) {
3926: PetscDrawAxis axis;
3927: const char *ylabel = ctx->semilogy ? "Log Time Step" : "Time Step";
3928: PetscDrawLGGetAxis(ctx->lg,&axis);
3929: PetscDrawAxisSetLabels(axis,"Timestep as function of time","Time",ylabel);
3930: PetscDrawLGReset(ctx->lg);
3931: }
3932: TSGetTimeStep(ts,&y);
3933: if (ctx->semilogy) y = PetscLog10Real(y);
3934: PetscDrawLGAddPoint(ctx->lg,&x,&y);
3935: if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
3936: PetscDrawLGDraw(ctx->lg);
3937: PetscDrawLGSave(ctx->lg);
3938: }
3939: return(0);
3940: }
3942: /*@C
3943: TSMonitorLGCtxDestroy - Destroys a line graph context that was created
3944: with TSMonitorLGCtxCreate().
3946: Collective on TSMonitorLGCtx
3948: Input Parameter:
3949: . ctx - the monitor context
3951: Level: intermediate
3953: .seealso: TSMonitorLGCtxCreate(), TSMonitorSet(), TSMonitorLGTimeStep();
3954: @*/
3955: PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx *ctx)
3956: {
3960: if ((*ctx)->transformdestroy) {
3961: ((*ctx)->transformdestroy)((*ctx)->transformctx);
3962: }
3963: PetscDrawLGDestroy(&(*ctx)->lg);
3964: PetscStrArrayDestroy(&(*ctx)->names);
3965: PetscStrArrayDestroy(&(*ctx)->displaynames);
3966: PetscFree((*ctx)->displayvariables);
3967: PetscFree((*ctx)->displayvalues);
3968: PetscFree(*ctx);
3969: return(0);
3970: }
3972: /*
3973:
3974: Creates a TS Monitor SPCtx for use with DM Swarm particle visualizations
3976: */
3977: PetscErrorCode TSMonitorSPCtxCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscInt howoften,TSMonitorSPCtx *ctx)
3978: {
3979: PetscDraw draw;
3983: PetscNew(ctx);
3984: PetscDrawCreate(comm,host,label,x,y,m,n,&draw);
3985: PetscDrawSetFromOptions(draw);
3986: PetscDrawSPCreate(draw,1,&(*ctx)->sp);
3987: PetscDrawDestroy(&draw);
3988: (*ctx)->howoften = howoften;
3989: return(0);
3991: }
3993: /*
3994: Destroys a TSMonitorSPCtx that was created with TSMonitorSPCtxCreate
3995: */
3996: PetscErrorCode TSMonitorSPCtxDestroy(TSMonitorSPCtx *ctx)
3997: {
4001:
4002: PetscDrawSPDestroy(&(*ctx)->sp);
4003: PetscFree(*ctx);
4004:
4005: return(0);
4007: }
4009: /*@
4010: TSGetTime - Gets the time of the most recently completed step.
4012: Not Collective
4014: Input Parameter:
4015: . ts - the TS context obtained from TSCreate()
4017: Output Parameter:
4018: . t - the current time. This time may not corresponds to the final time set with TSSetMaxTime(), use TSGetSolveTime().
4020: Level: beginner
4022: Note:
4023: When called during time step evaluation (e.g. during residual evaluation or via hooks set using TSSetPreStep(),
4024: TSSetPreStage(), TSSetPostStage(), or TSSetPostStep()), the time is the time at the start of the step being evaluated.
4026: .seealso: TSGetSolveTime(), TSSetTime(), TSGetTimeStep(), TSGetStepNumber()
4028: @*/
4029: PetscErrorCode TSGetTime(TS ts,PetscReal *t)
4030: {
4034: *t = ts->ptime;
4035: return(0);
4036: }
4038: /*@
4039: TSGetPrevTime - Gets the starting time of the previously completed step.
4041: Not Collective
4043: Input Parameter:
4044: . ts - the TS context obtained from TSCreate()
4046: Output Parameter:
4047: . t - the previous time
4049: Level: beginner
4051: .seealso: TSGetTime(), TSGetSolveTime(), TSGetTimeStep()
4053: @*/
4054: PetscErrorCode TSGetPrevTime(TS ts,PetscReal *t)
4055: {
4059: *t = ts->ptime_prev;
4060: return(0);
4061: }
4063: /*@
4064: TSSetTime - Allows one to reset the time.
4066: Logically Collective on TS
4068: Input Parameters:
4069: + ts - the TS context obtained from TSCreate()
4070: - time - the time
4072: Level: intermediate
4074: .seealso: TSGetTime(), TSSetMaxSteps()
4076: @*/
4077: PetscErrorCode TSSetTime(TS ts, PetscReal t)
4078: {
4082: ts->ptime = t;
4083: return(0);
4084: }
4086: /*@C
4087: TSSetOptionsPrefix - Sets the prefix used for searching for all
4088: TS options in the database.
4090: Logically Collective on TS
4092: Input Parameter:
4093: + ts - The TS context
4094: - prefix - The prefix to prepend to all option names
4096: Notes:
4097: A hyphen (-) must NOT be given at the beginning of the prefix name.
4098: The first character of all runtime options is AUTOMATICALLY the
4099: hyphen.
4101: Level: advanced
4103: .seealso: TSSetFromOptions()
4105: @*/
4106: PetscErrorCode TSSetOptionsPrefix(TS ts,const char prefix[])
4107: {
4109: SNES snes;
4113: PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);
4114: TSGetSNES(ts,&snes);
4115: SNESSetOptionsPrefix(snes,prefix);
4116: return(0);
4117: }
4119: /*@C
4120: TSAppendOptionsPrefix - Appends to the prefix used for searching for all
4121: TS options in the database.
4123: Logically Collective on TS
4125: Input Parameter:
4126: + ts - The TS context
4127: - prefix - The prefix to prepend to all option names
4129: Notes:
4130: A hyphen (-) must NOT be given at the beginning of the prefix name.
4131: The first character of all runtime options is AUTOMATICALLY the
4132: hyphen.
4134: Level: advanced
4136: .seealso: TSGetOptionsPrefix()
4138: @*/
4139: PetscErrorCode TSAppendOptionsPrefix(TS ts,const char prefix[])
4140: {
4142: SNES snes;
4146: PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);
4147: TSGetSNES(ts,&snes);
4148: SNESAppendOptionsPrefix(snes,prefix);
4149: return(0);
4150: }
4152: /*@C
4153: TSGetOptionsPrefix - Sets the prefix used for searching for all
4154: TS options in the database.
4156: Not Collective
4158: Input Parameter:
4159: . ts - The TS context
4161: Output Parameter:
4162: . prefix - A pointer to the prefix string used
4164: Notes:
4165: On the fortran side, the user should pass in a string 'prifix' of
4166: sufficient length to hold the prefix.
4168: Level: intermediate
4170: .seealso: TSAppendOptionsPrefix()
4171: @*/
4172: PetscErrorCode TSGetOptionsPrefix(TS ts,const char *prefix[])
4173: {
4179: PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);
4180: return(0);
4181: }
4183: /*@C
4184: TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
4186: Not Collective, but parallel objects are returned if TS is parallel
4188: Input Parameter:
4189: . ts - The TS context obtained from TSCreate()
4191: Output Parameters:
4192: + Amat - The (approximate) Jacobian J of G, where U_t = G(U,t) (or NULL)
4193: . Pmat - The matrix from which the preconditioner is constructed, usually the same as Amat (or NULL)
4194: . func - Function to compute the Jacobian of the RHS (or NULL)
4195: - ctx - User-defined context for Jacobian evaluation routine (or NULL)
4197: Notes:
4198: You can pass in NULL for any return argument you do not need.
4200: Level: intermediate
4202: .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetStepNumber()
4204: @*/
4205: PetscErrorCode TSGetRHSJacobian(TS ts,Mat *Amat,Mat *Pmat,TSRHSJacobian *func,void **ctx)
4206: {
4208: DM dm;
4211: if (Amat || Pmat) {
4212: SNES snes;
4213: TSGetSNES(ts,&snes);
4214: SNESSetUpMatrices(snes);
4215: SNESGetJacobian(snes,Amat,Pmat,NULL,NULL);
4216: }
4217: TSGetDM(ts,&dm);
4218: DMTSGetRHSJacobian(dm,func,ctx);
4219: return(0);
4220: }
4222: /*@C
4223: TSGetIJacobian - Returns the implicit Jacobian at the present timestep.
4225: Not Collective, but parallel objects are returned if TS is parallel
4227: Input Parameter:
4228: . ts - The TS context obtained from TSCreate()
4230: Output Parameters:
4231: + Amat - The (approximate) Jacobian of F(t,U,U_t)
4232: . Pmat - The matrix from which the preconditioner is constructed, often the same as Amat
4233: . f - The function to compute the matrices
4234: - ctx - User-defined context for Jacobian evaluation routine
4236: Notes:
4237: You can pass in NULL for any return argument you do not need.
4239: Level: advanced
4241: .seealso: TSGetTimeStep(), TSGetRHSJacobian(), TSGetMatrices(), TSGetTime(), TSGetStepNumber()
4243: @*/
4244: PetscErrorCode TSGetIJacobian(TS ts,Mat *Amat,Mat *Pmat,TSIJacobian *f,void **ctx)
4245: {
4247: DM dm;
4250: if (Amat || Pmat) {
4251: SNES snes;
4252: TSGetSNES(ts,&snes);
4253: SNESSetUpMatrices(snes);
4254: SNESGetJacobian(snes,Amat,Pmat,NULL,NULL);
4255: }
4256: TSGetDM(ts,&dm);
4257: DMTSGetIJacobian(dm,f,ctx);
4258: return(0);
4259: }
4261: /*@C
4262: TSMonitorDrawSolution - Monitors progress of the TS solvers by calling
4263: VecView() for the solution at each timestep
4265: Collective on TS
4267: Input Parameters:
4268: + ts - the TS context
4269: . step - current time-step
4270: . ptime - current time
4271: - dummy - either a viewer or NULL
4273: Options Database:
4274: . -ts_monitor_draw_solution_initial - show initial solution as well as current solution
4276: Notes:
4277: the initial solution and current solution are not display with a common axis scaling so generally the option -ts_monitor_draw_solution_initial
4278: will look bad
4280: Level: intermediate
4282: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
4283: @*/
4284: PetscErrorCode TSMonitorDrawSolution(TS ts,PetscInt step,PetscReal ptime,Vec u,void *dummy)
4285: {
4286: PetscErrorCode ierr;
4287: TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy;
4288: PetscDraw draw;
4291: if (!step && ictx->showinitial) {
4292: if (!ictx->initialsolution) {
4293: VecDuplicate(u,&ictx->initialsolution);
4294: }
4295: VecCopy(u,ictx->initialsolution);
4296: }
4297: if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) return(0);
4299: if (ictx->showinitial) {
4300: PetscReal pause;
4301: PetscViewerDrawGetPause(ictx->viewer,&pause);
4302: PetscViewerDrawSetPause(ictx->viewer,0.0);
4303: VecView(ictx->initialsolution,ictx->viewer);
4304: PetscViewerDrawSetPause(ictx->viewer,pause);
4305: PetscViewerDrawSetHold(ictx->viewer,PETSC_TRUE);
4306: }
4307: VecView(u,ictx->viewer);
4308: if (ictx->showtimestepandtime) {
4309: PetscReal xl,yl,xr,yr,h;
4310: char time[32];
4312: PetscViewerDrawGetDraw(ictx->viewer,0,&draw);
4313: PetscSNPrintf(time,32,"Timestep %d Time %g",(int)step,(double)ptime);
4314: PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);
4315: h = yl + .95*(yr - yl);
4316: PetscDrawStringCentered(draw,.5*(xl+xr),h,PETSC_DRAW_BLACK,time);
4317: PetscDrawFlush(draw);
4318: }
4320: if (ictx->showinitial) {
4321: PetscViewerDrawSetHold(ictx->viewer,PETSC_FALSE);
4322: }
4323: return(0);
4324: }
4326: /*@C
4327: TSMonitorDrawSolutionPhase - Monitors progress of the TS solvers by plotting the solution as a phase diagram
4329: Collective on TS
4331: Input Parameters:
4332: + ts - the TS context
4333: . step - current time-step
4334: . ptime - current time
4335: - dummy - either a viewer or NULL
4337: Level: intermediate
4339: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
4340: @*/
4341: PetscErrorCode TSMonitorDrawSolutionPhase(TS ts,PetscInt step,PetscReal ptime,Vec u,void *dummy)
4342: {
4343: PetscErrorCode ierr;
4344: TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy;
4345: PetscDraw draw;
4346: PetscDrawAxis axis;
4347: PetscInt n;
4348: PetscMPIInt size;
4349: PetscReal U0,U1,xl,yl,xr,yr,h;
4350: char time[32];
4351: const PetscScalar *U;
4354: MPI_Comm_size(PetscObjectComm((PetscObject)ts),&size);
4355: if (size != 1) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Only allowed for sequential runs");
4356: VecGetSize(u,&n);
4357: if (n != 2) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Only for ODEs with two unknowns");
4359: PetscViewerDrawGetDraw(ictx->viewer,0,&draw);
4360: PetscViewerDrawGetDrawAxis(ictx->viewer,0,&axis);
4361: PetscDrawAxisGetLimits(axis,&xl,&xr,&yl,&yr);
4362: if (!step) {
4363: PetscDrawClear(draw);
4364: PetscDrawAxisDraw(axis);
4365: }
4367: VecGetArrayRead(u,&U);
4368: U0 = PetscRealPart(U[0]);
4369: U1 = PetscRealPart(U[1]);
4370: VecRestoreArrayRead(u,&U);
4371: if ((U0 < xl) || (U1 < yl) || (U0 > xr) || (U1 > yr)) return(0);
4373: PetscDrawCollectiveBegin(draw);
4374: PetscDrawPoint(draw,U0,U1,PETSC_DRAW_BLACK);
4375: if (ictx->showtimestepandtime) {
4376: PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);
4377: PetscSNPrintf(time,32,"Timestep %d Time %g",(int)step,(double)ptime);
4378: h = yl + .95*(yr - yl);
4379: PetscDrawStringCentered(draw,.5*(xl+xr),h,PETSC_DRAW_BLACK,time);
4380: }
4381: PetscDrawCollectiveEnd(draw);
4382: PetscDrawFlush(draw);
4383: PetscDrawPause(draw);
4384: PetscDrawSave(draw);
4385: return(0);
4386: }
4388: /*@C
4389: TSMonitorDrawCtxDestroy - Destroys the monitor context for TSMonitorDrawSolution()
4391: Collective on TS
4393: Input Parameters:
4394: . ctx - the monitor context
4396: Level: intermediate
4398: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorDrawSolution(), TSMonitorDrawError()
4399: @*/
4400: PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx *ictx)
4401: {
4405: PetscViewerDestroy(&(*ictx)->viewer);
4406: VecDestroy(&(*ictx)->initialsolution);
4407: PetscFree(*ictx);
4408: return(0);
4409: }
4411: /*@C
4412: TSMonitorDrawCtxCreate - Creates the monitor context for TSMonitorDrawCtx
4414: Collective on TS
4416: Input Parameter:
4417: . ts - time-step context
4419: Output Patameter:
4420: . ctx - the monitor context
4422: Options Database:
4423: . -ts_monitor_draw_solution_initial - show initial solution as well as current solution
4425: Level: intermediate
4427: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorDrawCtx()
4428: @*/
4429: PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscInt howoften,TSMonitorDrawCtx *ctx)
4430: {
4431: PetscErrorCode ierr;
4434: PetscNew(ctx);
4435: PetscViewerDrawOpen(comm,host,label,x,y,m,n,&(*ctx)->viewer);
4436: PetscViewerSetFromOptions((*ctx)->viewer);
4438: (*ctx)->howoften = howoften;
4439: (*ctx)->showinitial = PETSC_FALSE;
4440: PetscOptionsGetBool(NULL,NULL,"-ts_monitor_draw_solution_initial",&(*ctx)->showinitial,NULL);
4442: (*ctx)->showtimestepandtime = PETSC_FALSE;
4443: PetscOptionsGetBool(NULL,NULL,"-ts_monitor_draw_solution_show_time",&(*ctx)->showtimestepandtime,NULL);
4444: return(0);
4445: }
4447: /*@C
4448: TSMonitorDrawSolutionFunction - Monitors progress of the TS solvers by calling
4449: VecView() for the solution provided by TSSetSolutionFunction() at each timestep
4451: Collective on TS
4453: Input Parameters:
4454: + ts - the TS context
4455: . step - current time-step
4456: . ptime - current time
4457: - dummy - either a viewer or NULL
4459: Options Database:
4460: . -ts_monitor_draw_solution_function - Monitor error graphically, requires user to have provided TSSetSolutionFunction()
4462: Level: intermediate
4464: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSSetSolutionFunction()
4465: @*/
4466: PetscErrorCode TSMonitorDrawSolutionFunction(TS ts,PetscInt step,PetscReal ptime,Vec u,void *dummy)
4467: {
4468: PetscErrorCode ierr;
4469: TSMonitorDrawCtx ctx = (TSMonitorDrawCtx)dummy;
4470: PetscViewer viewer = ctx->viewer;
4471: Vec work;
4474: if (!(((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason))) return(0);
4475: VecDuplicate(u,&work);
4476: TSComputeSolutionFunction(ts,ptime,work);
4477: VecView(work,viewer);
4478: VecDestroy(&work);
4479: return(0);
4480: }
4482: /*@C
4483: TSMonitorDrawError - Monitors progress of the TS solvers by calling
4484: VecView() for the error at each timestep
4486: Collective on TS
4488: Input Parameters:
4489: + ts - the TS context
4490: . step - current time-step
4491: . ptime - current time
4492: - dummy - either a viewer or NULL
4494: Options Database:
4495: . -ts_monitor_draw_error - Monitor error graphically, requires user to have provided TSSetSolutionFunction()
4497: Level: intermediate
4499: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSSetSolutionFunction()
4500: @*/
4501: PetscErrorCode TSMonitorDrawError(TS ts,PetscInt step,PetscReal ptime,Vec u,void *dummy)
4502: {
4503: PetscErrorCode ierr;
4504: TSMonitorDrawCtx ctx = (TSMonitorDrawCtx)dummy;
4505: PetscViewer viewer = ctx->viewer;
4506: Vec work;
4509: if (!(((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason))) return(0);
4510: VecDuplicate(u,&work);
4511: TSComputeSolutionFunction(ts,ptime,work);
4512: VecAXPY(work,-1.0,u);
4513: VecView(work,viewer);
4514: VecDestroy(&work);
4515: return(0);
4516: }
4518: #include <petsc/private/dmimpl.h>
4519: /*@
4520: TSSetDM - Sets the DM that may be used by some nonlinear solvers or preconditioners under the TS
4522: Logically Collective on ts
4524: Input Parameters:
4525: + ts - the ODE integrator object
4526: - dm - the dm, cannot be NULL
4528: Notes:
4529: A DM can only be used for solving one problem at a time because information about the problem is stored on the DM,
4530: even when not using interfaces like DMTSSetIFunction(). Use DMClone() to get a distinct DM when solving
4531: different problems using the same function space.
4533: Level: intermediate
4535: .seealso: TSGetDM(), SNESSetDM(), SNESGetDM()
4536: @*/
4537: PetscErrorCode TSSetDM(TS ts,DM dm)
4538: {
4540: SNES snes;
4541: DMTS tsdm;
4546: PetscObjectReference((PetscObject)dm);
4547: if (ts->dm) { /* Move the DMTS context over to the new DM unless the new DM already has one */
4548: if (ts->dm->dmts && !dm->dmts) {
4549: DMCopyDMTS(ts->dm,dm);
4550: DMGetDMTS(ts->dm,&tsdm);
4551: if (tsdm->originaldm == ts->dm) { /* Grant write privileges to the replacement DM */
4552: tsdm->originaldm = dm;
4553: }
4554: }
4555: DMDestroy(&ts->dm);
4556: }
4557: ts->dm = dm;
4559: TSGetSNES(ts,&snes);
4560: SNESSetDM(snes,dm);
4561: return(0);
4562: }
4564: /*@
4565: TSGetDM - Gets the DM that may be used by some preconditioners
4567: Not Collective
4569: Input Parameter:
4570: . ts - the preconditioner context
4572: Output Parameter:
4573: . dm - the dm
4575: Level: intermediate
4577: .seealso: TSSetDM(), SNESSetDM(), SNESGetDM()
4578: @*/
4579: PetscErrorCode TSGetDM(TS ts,DM *dm)
4580: {
4585: if (!ts->dm) {
4586: DMShellCreate(PetscObjectComm((PetscObject)ts),&ts->dm);
4587: if (ts->snes) {SNESSetDM(ts->snes,ts->dm);}
4588: }
4589: *dm = ts->dm;
4590: return(0);
4591: }
4593: /*@
4594: SNESTSFormFunction - Function to evaluate nonlinear residual
4596: Logically Collective on SNES
4598: Input Parameter:
4599: + snes - nonlinear solver
4600: . U - the current state at which to evaluate the residual
4601: - ctx - user context, must be a TS
4603: Output Parameter:
4604: . F - the nonlinear residual
4606: Notes:
4607: This function is not normally called by users and is automatically registered with the SNES used by TS.
4608: It is most frequently passed to MatFDColoringSetFunction().
4610: Level: advanced
4612: .seealso: SNESSetFunction(), MatFDColoringSetFunction()
4613: @*/
4614: PetscErrorCode SNESTSFormFunction(SNES snes,Vec U,Vec F,void *ctx)
4615: {
4616: TS ts = (TS)ctx;
4624: (ts->ops->snesfunction)(snes,U,F,ts);
4625: return(0);
4626: }
4628: /*@
4629: SNESTSFormJacobian - Function to evaluate the Jacobian
4631: Collective on SNES
4633: Input Parameter:
4634: + snes - nonlinear solver
4635: . U - the current state at which to evaluate the residual
4636: - ctx - user context, must be a TS
4638: Output Parameter:
4639: + A - the Jacobian
4640: . B - the preconditioning matrix (may be the same as A)
4641: - flag - indicates any structure change in the matrix
4643: Notes:
4644: This function is not normally called by users and is automatically registered with the SNES used by TS.
4646: Level: developer
4648: .seealso: SNESSetJacobian()
4649: @*/
4650: PetscErrorCode SNESTSFormJacobian(SNES snes,Vec U,Mat A,Mat B,void *ctx)
4651: {
4652: TS ts = (TS)ctx;
4663: (ts->ops->snesjacobian)(snes,U,A,B,ts);
4664: return(0);
4665: }
4667: /*@C
4668: TSComputeRHSFunctionLinear - Evaluate the right hand side via the user-provided Jacobian, for linear problems Udot = A U only
4670: Collective on TS
4672: Input Arguments:
4673: + ts - time stepping context
4674: . t - time at which to evaluate
4675: . U - state at which to evaluate
4676: - ctx - context
4678: Output Arguments:
4679: . F - right hand side
4681: Level: intermediate
4683: Notes:
4684: This function is intended to be passed to TSSetRHSFunction() to evaluate the right hand side for linear problems.
4685: The matrix (and optionally the evaluation context) should be passed to TSSetRHSJacobian().
4687: .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSJacobianConstant()
4688: @*/
4689: PetscErrorCode TSComputeRHSFunctionLinear(TS ts,PetscReal t,Vec U,Vec F,void *ctx)
4690: {
4692: Mat Arhs,Brhs;
4695: TSGetRHSMats_Private(ts,&Arhs,&Brhs);
4696: TSComputeRHSJacobian(ts,t,U,Arhs,Brhs);
4697: MatMult(Arhs,U,F);
4698: return(0);
4699: }
4701: /*@C
4702: TSComputeRHSJacobianConstant - Reuses a Jacobian that is time-independent.
4704: Collective on TS
4706: Input Arguments:
4707: + ts - time stepping context
4708: . t - time at which to evaluate
4709: . U - state at which to evaluate
4710: - ctx - context
4712: Output Arguments:
4713: + A - pointer to operator
4714: . B - pointer to preconditioning matrix
4715: - flg - matrix structure flag
4717: Level: intermediate
4719: Notes:
4720: This function is intended to be passed to TSSetRHSJacobian() to evaluate the Jacobian for linear time-independent problems.
4722: .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSFunctionLinear()
4723: @*/
4724: PetscErrorCode TSComputeRHSJacobianConstant(TS ts,PetscReal t,Vec U,Mat A,Mat B,void *ctx)
4725: {
4727: return(0);
4728: }
4730: /*@C
4731: TSComputeIFunctionLinear - Evaluate the left hand side via the user-provided Jacobian, for linear problems only
4733: Collective on TS
4735: Input Arguments:
4736: + ts - time stepping context
4737: . t - time at which to evaluate
4738: . U - state at which to evaluate
4739: . Udot - time derivative of state vector
4740: - ctx - context
4742: Output Arguments:
4743: . F - left hand side
4745: Level: intermediate
4747: Notes:
4748: 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
4749: user is required to write their own TSComputeIFunction.
4750: This function is intended to be passed to TSSetIFunction() to evaluate the left hand side for linear problems.
4751: The matrix (and optionally the evaluation context) should be passed to TSSetIJacobian().
4753: Note that using this function is NOT equivalent to using TSComputeRHSFunctionLinear() since that solves Udot = A U
4755: .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIJacobianConstant(), TSComputeRHSFunctionLinear()
4756: @*/
4757: PetscErrorCode TSComputeIFunctionLinear(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void *ctx)
4758: {
4760: Mat A,B;
4763: TSGetIJacobian(ts,&A,&B,NULL,NULL);
4764: TSComputeIJacobian(ts,t,U,Udot,1.0,A,B,PETSC_TRUE);
4765: MatMult(A,Udot,F);
4766: return(0);
4767: }
4769: /*@C
4770: TSComputeIJacobianConstant - Reuses a time-independent for a semi-implicit DAE or ODE
4772: Collective on TS
4774: Input Arguments:
4775: + ts - time stepping context
4776: . t - time at which to evaluate
4777: . U - state at which to evaluate
4778: . Udot - time derivative of state vector
4779: . shift - shift to apply
4780: - ctx - context
4782: Output Arguments:
4783: + A - pointer to operator
4784: . B - pointer to preconditioning matrix
4785: - flg - matrix structure flag
4787: Level: advanced
4789: Notes:
4790: This function is intended to be passed to TSSetIJacobian() to evaluate the Jacobian for linear time-independent problems.
4792: It is only appropriate for problems of the form
4794: $ M Udot = F(U,t)
4796: where M is constant and F is non-stiff. The user must pass M to TSSetIJacobian(). The current implementation only
4797: works with IMEX time integration methods such as TSROSW and TSARKIMEX, since there is no support for de-constructing
4798: an implicit operator of the form
4800: $ shift*M + J
4802: where J is the Jacobian of -F(U). Support may be added in a future version of PETSc, but for now, the user must store
4803: a copy of M or reassemble it when requested.
4805: .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIFunctionLinear()
4806: @*/
4807: PetscErrorCode TSComputeIJacobianConstant(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal shift,Mat A,Mat B,void *ctx)
4808: {
4812: MatScale(A, shift / ts->ijacobian.shift);
4813: ts->ijacobian.shift = shift;
4814: return(0);
4815: }
4817: /*@
4818: TSGetEquationType - Gets the type of the equation that TS is solving.
4820: Not Collective
4822: Input Parameter:
4823: . ts - the TS context
4825: Output Parameter:
4826: . equation_type - see TSEquationType
4828: Level: beginner
4830: .seealso: TSSetEquationType(), TSEquationType
4831: @*/
4832: PetscErrorCode TSGetEquationType(TS ts,TSEquationType *equation_type)
4833: {
4837: *equation_type = ts->equation_type;
4838: return(0);
4839: }
4841: /*@
4842: TSSetEquationType - Sets the type of the equation that TS is solving.
4844: Not Collective
4846: Input Parameter:
4847: + ts - the TS context
4848: - equation_type - see TSEquationType
4850: Level: advanced
4852: .seealso: TSGetEquationType(), TSEquationType
4853: @*/
4854: PetscErrorCode TSSetEquationType(TS ts,TSEquationType equation_type)
4855: {
4858: ts->equation_type = equation_type;
4859: return(0);
4860: }
4862: /*@
4863: TSGetConvergedReason - Gets the reason the TS iteration was stopped.
4865: Not Collective
4867: Input Parameter:
4868: . ts - the TS context
4870: Output Parameter:
4871: . reason - negative value indicates diverged, positive value converged, see TSConvergedReason or the
4872: manual pages for the individual convergence tests for complete lists
4874: Level: beginner
4876: Notes:
4877: Can only be called after the call to TSSolve() is complete.
4879: .seealso: TSSetConvergenceTest(), TSConvergedReason
4880: @*/
4881: PetscErrorCode TSGetConvergedReason(TS ts,TSConvergedReason *reason)
4882: {
4886: *reason = ts->reason;
4887: return(0);
4888: }
4890: /*@
4891: TSSetConvergedReason - Sets the reason for handling the convergence of TSSolve.
4893: Logically Collective; reason must contain common value
4895: Input Parameters:
4896: + ts - the TS context
4897: - reason - negative value indicates diverged, positive value converged, see TSConvergedReason or the
4898: manual pages for the individual convergence tests for complete lists
4900: Level: advanced
4902: Notes:
4903: Can only be called while TSSolve() is active.
4905: .seealso: TSConvergedReason
4906: @*/
4907: PetscErrorCode TSSetConvergedReason(TS ts,TSConvergedReason reason)
4908: {
4911: ts->reason = reason;
4912: return(0);
4913: }
4915: /*@
4916: TSGetSolveTime - Gets the time after a call to TSSolve()
4918: Not Collective
4920: Input Parameter:
4921: . ts - the TS context
4923: Output Parameter:
4924: . ftime - the final time. This time corresponds to the final time set with TSSetMaxTime()
4926: Level: beginner
4928: Notes:
4929: Can only be called after the call to TSSolve() is complete.
4931: .seealso: TSSetConvergenceTest(), TSConvergedReason
4932: @*/
4933: PetscErrorCode TSGetSolveTime(TS ts,PetscReal *ftime)
4934: {
4938: *ftime = ts->solvetime;
4939: return(0);
4940: }
4942: /*@
4943: TSGetSNESIterations - Gets the total number of nonlinear iterations
4944: used by the time integrator.
4946: Not Collective
4948: Input Parameter:
4949: . ts - TS context
4951: Output Parameter:
4952: . nits - number of nonlinear iterations
4954: Notes:
4955: This counter is reset to zero for each successive call to TSSolve().
4957: Level: intermediate
4959: .seealso: TSGetKSPIterations()
4960: @*/
4961: PetscErrorCode TSGetSNESIterations(TS ts,PetscInt *nits)
4962: {
4966: *nits = ts->snes_its;
4967: return(0);
4968: }
4970: /*@
4971: TSGetKSPIterations - Gets the total number of linear iterations
4972: used by the time integrator.
4974: Not Collective
4976: Input Parameter:
4977: . ts - TS context
4979: Output Parameter:
4980: . lits - number of linear iterations
4982: Notes:
4983: This counter is reset to zero for each successive call to TSSolve().
4985: Level: intermediate
4987: .seealso: TSGetSNESIterations(), SNESGetKSPIterations()
4988: @*/
4989: PetscErrorCode TSGetKSPIterations(TS ts,PetscInt *lits)
4990: {
4994: *lits = ts->ksp_its;
4995: return(0);
4996: }
4998: /*@
4999: TSGetStepRejections - Gets the total number of rejected steps.
5001: Not Collective
5003: Input Parameter:
5004: . ts - TS context
5006: Output Parameter:
5007: . rejects - number of steps rejected
5009: Notes:
5010: This counter is reset to zero for each successive call to TSSolve().
5012: Level: intermediate
5014: .seealso: TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetSNESFailures(), TSSetMaxSNESFailures(), TSSetErrorIfStepFails()
5015: @*/
5016: PetscErrorCode TSGetStepRejections(TS ts,PetscInt *rejects)
5017: {
5021: *rejects = ts->reject;
5022: return(0);
5023: }
5025: /*@
5026: TSGetSNESFailures - Gets the total number of failed SNES solves
5028: Not Collective
5030: Input Parameter:
5031: . ts - TS context
5033: Output Parameter:
5034: . fails - number of failed nonlinear solves
5036: Notes:
5037: This counter is reset to zero for each successive call to TSSolve().
5039: Level: intermediate
5041: .seealso: TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSSetMaxSNESFailures()
5042: @*/
5043: PetscErrorCode TSGetSNESFailures(TS ts,PetscInt *fails)
5044: {
5048: *fails = ts->num_snes_failures;
5049: return(0);
5050: }
5052: /*@
5053: TSSetMaxStepRejections - Sets the maximum number of step rejections before a step fails
5055: Not Collective
5057: Input Parameter:
5058: + ts - TS context
5059: - rejects - maximum number of rejected steps, pass -1 for unlimited
5061: Notes:
5062: The counter is reset to zero for each step
5064: Options Database Key:
5065: . -ts_max_reject - Maximum number of step rejections before a step fails
5067: Level: intermediate
5069: .seealso: TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxSNESFailures(), TSGetStepRejections(), TSGetSNESFailures(), TSSetErrorIfStepFails(), TSGetConvergedReason()
5070: @*/
5071: PetscErrorCode TSSetMaxStepRejections(TS ts,PetscInt rejects)
5072: {
5075: ts->max_reject = rejects;
5076: return(0);
5077: }
5079: /*@
5080: TSSetMaxSNESFailures - Sets the maximum number of failed SNES solves
5082: Not Collective
5084: Input Parameter:
5085: + ts - TS context
5086: - fails - maximum number of failed nonlinear solves, pass -1 for unlimited
5088: Notes:
5089: The counter is reset to zero for each successive call to TSSolve().
5091: Options Database Key:
5092: . -ts_max_snes_failures - Maximum number of nonlinear solve failures
5094: Level: intermediate
5096: .seealso: TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSGetSNESFailures(), SNESGetConvergedReason(), TSGetConvergedReason()
5097: @*/
5098: PetscErrorCode TSSetMaxSNESFailures(TS ts,PetscInt fails)
5099: {
5102: ts->max_snes_failures = fails;
5103: return(0);
5104: }
5106: /*@
5107: TSSetErrorIfStepFails - Error if no step succeeds
5109: Not Collective
5111: Input Parameter:
5112: + ts - TS context
5113: - err - PETSC_TRUE to error if no step succeeds, PETSC_FALSE to return without failure
5115: Options Database Key:
5116: . -ts_error_if_step_fails - Error if no step succeeds
5118: Level: intermediate
5120: .seealso: TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSGetSNESFailures(), TSSetErrorIfStepFails(), TSGetConvergedReason()
5121: @*/
5122: PetscErrorCode TSSetErrorIfStepFails(TS ts,PetscBool err)
5123: {
5126: ts->errorifstepfailed = err;
5127: return(0);
5128: }
5130: /*@C
5131: TSMonitorSolution - Monitors progress of the TS solvers by VecView() for the solution at each timestep. Normally the viewer is a binary file or a PetscDraw object
5133: Collective on TS
5135: Input Parameters:
5136: + ts - the TS context
5137: . step - current time-step
5138: . ptime - current time
5139: . u - current state
5140: - vf - viewer and its format
5142: Level: intermediate
5144: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
5145: @*/
5146: PetscErrorCode TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscViewerAndFormat *vf)
5147: {
5151: PetscViewerPushFormat(vf->viewer,vf->format);
5152: VecView(u,vf->viewer);
5153: PetscViewerPopFormat(vf->viewer);
5154: return(0);
5155: }
5157: /*@C
5158: TSMonitorSolutionVTK - Monitors progress of the TS solvers by VecView() for the solution at each timestep.
5160: Collective on TS
5162: Input Parameters:
5163: + ts - the TS context
5164: . step - current time-step
5165: . ptime - current time
5166: . u - current state
5167: - filenametemplate - string containing a format specifier for the integer time step (e.g. %03D)
5169: Level: intermediate
5171: Notes:
5172: 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.
5173: These are named according to the file name template.
5175: This function is normally passed as an argument to TSMonitorSet() along with TSMonitorSolutionVTKDestroy().
5177: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
5178: @*/
5179: PetscErrorCode TSMonitorSolutionVTK(TS ts,PetscInt step,PetscReal ptime,Vec u,void *filenametemplate)
5180: {
5182: char filename[PETSC_MAX_PATH_LEN];
5183: PetscViewer viewer;
5186: if (step < 0) return(0); /* -1 indicates interpolated solution */
5187: PetscSNPrintf(filename,sizeof(filename),(const char*)filenametemplate,step);
5188: PetscViewerVTKOpen(PetscObjectComm((PetscObject)ts),filename,FILE_MODE_WRITE,&viewer);
5189: VecView(u,viewer);
5190: PetscViewerDestroy(&viewer);
5191: return(0);
5192: }
5194: /*@C
5195: TSMonitorSolutionVTKDestroy - Destroy context for monitoring
5197: Collective on TS
5199: Input Parameters:
5200: . filenametemplate - string containing a format specifier for the integer time step (e.g. %03D)
5202: Level: intermediate
5204: Note:
5205: This function is normally passed to TSMonitorSet() along with TSMonitorSolutionVTK().
5207: .seealso: TSMonitorSet(), TSMonitorSolutionVTK()
5208: @*/
5209: PetscErrorCode TSMonitorSolutionVTKDestroy(void *filenametemplate)
5210: {
5214: PetscFree(*(char**)filenametemplate);
5215: return(0);
5216: }
5218: /*@
5219: TSGetAdapt - Get the adaptive controller context for the current method
5221: Collective on TS if controller has not been created yet
5223: Input Arguments:
5224: . ts - time stepping context
5226: Output Arguments:
5227: . adapt - adaptive controller
5229: Level: intermediate
5231: .seealso: TSAdapt, TSAdaptSetType(), TSAdaptChoose()
5232: @*/
5233: PetscErrorCode TSGetAdapt(TS ts,TSAdapt *adapt)
5234: {
5240: if (!ts->adapt) {
5241: TSAdaptCreate(PetscObjectComm((PetscObject)ts),&ts->adapt);
5242: PetscLogObjectParent((PetscObject)ts,(PetscObject)ts->adapt);
5243: PetscObjectIncrementTabLevel((PetscObject)ts->adapt,(PetscObject)ts,1);
5244: }
5245: *adapt = ts->adapt;
5246: return(0);
5247: }
5249: /*@
5250: TSSetTolerances - Set tolerances for local truncation error when using adaptive controller
5252: Logically Collective
5254: Input Arguments:
5255: + ts - time integration context
5256: . atol - scalar absolute tolerances, PETSC_DECIDE to leave current value
5257: . vatol - vector of absolute tolerances or NULL, used in preference to atol if present
5258: . rtol - scalar relative tolerances, PETSC_DECIDE to leave current value
5259: - vrtol - vector of relative tolerances or NULL, used in preference to atol if present
5261: Options Database keys:
5262: + -ts_rtol <rtol> - relative tolerance for local truncation error
5263: - -ts_atol <atol> Absolute tolerance for local truncation error
5265: Notes:
5266: With PETSc's implicit schemes for DAE problems, the calculation of the local truncation error
5267: (LTE) includes both the differential and the algebraic variables. If one wants the LTE to be
5268: computed only for the differential or the algebraic part then this can be done using the vector of
5269: tolerances vatol. For example, by setting the tolerance vector with the desired tolerance for the
5270: differential part and infinity for the algebraic part, the LTE calculation will include only the
5271: differential variables.
5273: Level: beginner
5275: .seealso: TS, TSAdapt, TSVecNormWRMS(), TSGetTolerances()
5276: @*/
5277: PetscErrorCode TSSetTolerances(TS ts,PetscReal atol,Vec vatol,PetscReal rtol,Vec vrtol)
5278: {
5282: if (atol != PETSC_DECIDE && atol != PETSC_DEFAULT) ts->atol = atol;
5283: if (vatol) {
5284: PetscObjectReference((PetscObject)vatol);
5285: VecDestroy(&ts->vatol);
5286: ts->vatol = vatol;
5287: }
5288: if (rtol != PETSC_DECIDE && rtol != PETSC_DEFAULT) ts->rtol = rtol;
5289: if (vrtol) {
5290: PetscObjectReference((PetscObject)vrtol);
5291: VecDestroy(&ts->vrtol);
5292: ts->vrtol = vrtol;
5293: }
5294: return(0);
5295: }
5297: /*@
5298: TSGetTolerances - Get tolerances for local truncation error when using adaptive controller
5300: Logically Collective
5302: Input Arguments:
5303: . ts - time integration context
5305: Output Arguments:
5306: + atol - scalar absolute tolerances, NULL to ignore
5307: . vatol - vector of absolute tolerances, NULL to ignore
5308: . rtol - scalar relative tolerances, NULL to ignore
5309: - vrtol - vector of relative tolerances, NULL to ignore
5311: Level: beginner
5313: .seealso: TS, TSAdapt, TSVecNormWRMS(), TSSetTolerances()
5314: @*/
5315: PetscErrorCode TSGetTolerances(TS ts,PetscReal *atol,Vec *vatol,PetscReal *rtol,Vec *vrtol)
5316: {
5318: if (atol) *atol = ts->atol;
5319: if (vatol) *vatol = ts->vatol;
5320: if (rtol) *rtol = ts->rtol;
5321: if (vrtol) *vrtol = ts->vrtol;
5322: return(0);
5323: }
5325: /*@
5326: TSErrorWeightedNorm2 - compute a weighted 2-norm of the difference between two state vectors
5328: Collective on TS
5330: Input Arguments:
5331: + ts - time stepping context
5332: . U - state vector, usually ts->vec_sol
5333: - Y - state vector to be compared to U
5335: Output Arguments:
5336: + norm - weighted norm, a value of 1.0 means that the error matches the tolerances
5337: . norma - weighted norm based on the absolute tolerance, a value of 1.0 means that the error matches the tolerances
5338: - normr - weighted norm based on the relative tolerance, a value of 1.0 means that the error matches the tolerances
5340: Level: developer
5342: .seealso: TSErrorWeightedNorm(), TSErrorWeightedNormInfinity()
5343: @*/
5344: PetscErrorCode TSErrorWeightedNorm2(TS ts,Vec U,Vec Y,PetscReal *norm,PetscReal *norma,PetscReal *normr)
5345: {
5346: PetscErrorCode ierr;
5347: PetscInt i,n,N,rstart;
5348: PetscInt n_loc,na_loc,nr_loc;
5349: PetscReal n_glb,na_glb,nr_glb;
5350: const PetscScalar *u,*y;
5351: PetscReal sum,suma,sumr,gsum,gsuma,gsumr,diff;
5352: PetscReal tol,tola,tolr;
5353: PetscReal err_loc[6],err_glb[6];
5365: if (U == Y) SETERRQ(PetscObjectComm((PetscObject)U),PETSC_ERR_ARG_IDN,"U and Y cannot be the same vector");
5367: VecGetSize(U,&N);
5368: VecGetLocalSize(U,&n);
5369: VecGetOwnershipRange(U,&rstart,NULL);
5370: VecGetArrayRead(U,&u);
5371: VecGetArrayRead(Y,&y);
5372: sum = 0.; n_loc = 0;
5373: suma = 0.; na_loc = 0;
5374: sumr = 0.; nr_loc = 0;
5375: if (ts->vatol && ts->vrtol) {
5376: const PetscScalar *atol,*rtol;
5377: VecGetArrayRead(ts->vatol,&atol);
5378: VecGetArrayRead(ts->vrtol,&rtol);
5379: for (i=0; i<n; i++) {
5380: SkipSmallValue(y[i],u[i],ts->adapt->ignore_max);
5381: diff = PetscAbsScalar(y[i] - u[i]);
5382: tola = PetscRealPart(atol[i]);
5383: if(tola>0.){
5384: suma += PetscSqr(diff/tola);
5385: na_loc++;
5386: }
5387: tolr = PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5388: if(tolr>0.){
5389: sumr += PetscSqr(diff/tolr);
5390: nr_loc++;
5391: }
5392: tol=tola+tolr;
5393: if(tol>0.){
5394: sum += PetscSqr(diff/tol);
5395: n_loc++;
5396: }
5397: }
5398: VecRestoreArrayRead(ts->vatol,&atol);
5399: VecRestoreArrayRead(ts->vrtol,&rtol);
5400: } else if (ts->vatol) { /* vector atol, scalar rtol */
5401: const PetscScalar *atol;
5402: VecGetArrayRead(ts->vatol,&atol);
5403: for (i=0; i<n; i++) {
5404: SkipSmallValue(y[i],u[i],ts->adapt->ignore_max);
5405: diff = PetscAbsScalar(y[i] - u[i]);
5406: tola = PetscRealPart(atol[i]);
5407: if(tola>0.){
5408: suma += PetscSqr(diff/tola);
5409: na_loc++;
5410: }
5411: tolr = ts->rtol * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5412: if(tolr>0.){
5413: sumr += PetscSqr(diff/tolr);
5414: nr_loc++;
5415: }
5416: tol=tola+tolr;
5417: if(tol>0.){
5418: sum += PetscSqr(diff/tol);
5419: n_loc++;
5420: }
5421: }
5422: VecRestoreArrayRead(ts->vatol,&atol);
5423: } else if (ts->vrtol) { /* scalar atol, vector rtol */
5424: const PetscScalar *rtol;
5425: VecGetArrayRead(ts->vrtol,&rtol);
5426: for (i=0; i<n; i++) {
5427: SkipSmallValue(y[i],u[i],ts->adapt->ignore_max);
5428: diff = PetscAbsScalar(y[i] - u[i]);
5429: tola = ts->atol;
5430: if(tola>0.){
5431: suma += PetscSqr(diff/tola);
5432: na_loc++;
5433: }
5434: tolr = PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5435: if(tolr>0.){
5436: sumr += PetscSqr(diff/tolr);
5437: nr_loc++;
5438: }
5439: tol=tola+tolr;
5440: if(tol>0.){
5441: sum += PetscSqr(diff/tol);
5442: n_loc++;
5443: }
5444: }
5445: VecRestoreArrayRead(ts->vrtol,&rtol);
5446: } else { /* scalar atol, scalar rtol */
5447: for (i=0; i<n; i++) {
5448: SkipSmallValue(y[i],u[i],ts->adapt->ignore_max);
5449: diff = PetscAbsScalar(y[i] - u[i]);
5450: tola = ts->atol;
5451: if(tola>0.){
5452: suma += PetscSqr(diff/tola);
5453: na_loc++;
5454: }
5455: tolr = ts->rtol * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5456: if(tolr>0.){
5457: sumr += PetscSqr(diff/tolr);
5458: nr_loc++;
5459: }
5460: tol=tola+tolr;
5461: if(tol>0.){
5462: sum += PetscSqr(diff/tol);
5463: n_loc++;
5464: }
5465: }
5466: }
5467: VecRestoreArrayRead(U,&u);
5468: VecRestoreArrayRead(Y,&y);
5470: err_loc[0] = sum;
5471: err_loc[1] = suma;
5472: err_loc[2] = sumr;
5473: err_loc[3] = (PetscReal)n_loc;
5474: err_loc[4] = (PetscReal)na_loc;
5475: err_loc[5] = (PetscReal)nr_loc;
5477: MPIU_Allreduce(err_loc,err_glb,6,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)ts));
5479: gsum = err_glb[0];
5480: gsuma = err_glb[1];
5481: gsumr = err_glb[2];
5482: n_glb = err_glb[3];
5483: na_glb = err_glb[4];
5484: nr_glb = err_glb[5];
5486: *norm = 0.;
5487: if(n_glb>0. ){*norm = PetscSqrtReal(gsum / n_glb );}
5488: *norma = 0.;
5489: if(na_glb>0.){*norma = PetscSqrtReal(gsuma / na_glb);}
5490: *normr = 0.;
5491: if(nr_glb>0.){*normr = PetscSqrtReal(gsumr / nr_glb);}
5493: if (PetscIsInfOrNanScalar(*norm)) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
5494: if (PetscIsInfOrNanScalar(*norma)) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_FP,"Infinite or not-a-number generated in norma");
5495: if (PetscIsInfOrNanScalar(*normr)) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_FP,"Infinite or not-a-number generated in normr");
5496: return(0);
5497: }
5499: /*@
5500: TSErrorWeightedNormInfinity - compute a weighted infinity-norm of the difference between two state vectors
5502: Collective on TS
5504: Input Arguments:
5505: + ts - time stepping context
5506: . U - state vector, usually ts->vec_sol
5507: - Y - state vector to be compared to U
5509: Output Arguments:
5510: + norm - weighted norm, a value of 1.0 means that the error matches the tolerances
5511: . norma - weighted norm based on the absolute tolerance, a value of 1.0 means that the error matches the tolerances
5512: - normr - weighted norm based on the relative tolerance, a value of 1.0 means that the error matches the tolerances
5514: Level: developer
5516: .seealso: TSErrorWeightedNorm(), TSErrorWeightedNorm2()
5517: @*/
5518: PetscErrorCode TSErrorWeightedNormInfinity(TS ts,Vec U,Vec Y,PetscReal *norm,PetscReal *norma,PetscReal *normr)
5519: {
5520: PetscErrorCode ierr;
5521: PetscInt i,n,N,rstart;
5522: const PetscScalar *u,*y;
5523: PetscReal max,gmax,maxa,gmaxa,maxr,gmaxr;
5524: PetscReal tol,tola,tolr,diff;
5525: PetscReal err_loc[3],err_glb[3];
5537: if (U == Y) SETERRQ(PetscObjectComm((PetscObject)U),PETSC_ERR_ARG_IDN,"U and Y cannot be the same vector");
5539: VecGetSize(U,&N);
5540: VecGetLocalSize(U,&n);
5541: VecGetOwnershipRange(U,&rstart,NULL);
5542: VecGetArrayRead(U,&u);
5543: VecGetArrayRead(Y,&y);
5545: max=0.;
5546: maxa=0.;
5547: maxr=0.;
5549: if (ts->vatol && ts->vrtol) { /* vector atol, vector rtol */
5550: const PetscScalar *atol,*rtol;
5551: VecGetArrayRead(ts->vatol,&atol);
5552: VecGetArrayRead(ts->vrtol,&rtol);
5554: for (i=0; i<n; i++) {
5555: SkipSmallValue(y[i],u[i],ts->adapt->ignore_max);
5556: diff = PetscAbsScalar(y[i] - u[i]);
5557: tola = PetscRealPart(atol[i]);
5558: tolr = PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5559: tol = tola+tolr;
5560: if(tola>0.){
5561: maxa = PetscMax(maxa,diff / tola);
5562: }
5563: if(tolr>0.){
5564: maxr = PetscMax(maxr,diff / tolr);
5565: }
5566: if(tol>0.){
5567: max = PetscMax(max,diff / tol);
5568: }
5569: }
5570: VecRestoreArrayRead(ts->vatol,&atol);
5571: VecRestoreArrayRead(ts->vrtol,&rtol);
5572: } else if (ts->vatol) { /* vector atol, scalar rtol */
5573: const PetscScalar *atol;
5574: VecGetArrayRead(ts->vatol,&atol);
5575: for (i=0; i<n; i++) {
5576: SkipSmallValue(y[i],u[i],ts->adapt->ignore_max);
5577: diff = PetscAbsScalar(y[i] - u[i]);
5578: tola = PetscRealPart(atol[i]);
5579: tolr = ts->rtol * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5580: tol = tola+tolr;
5581: if(tola>0.){
5582: maxa = PetscMax(maxa,diff / tola);
5583: }
5584: if(tolr>0.){
5585: maxr = PetscMax(maxr,diff / tolr);
5586: }
5587: if(tol>0.){
5588: max = PetscMax(max,diff / tol);
5589: }
5590: }
5591: VecRestoreArrayRead(ts->vatol,&atol);
5592: } else if (ts->vrtol) { /* scalar atol, vector rtol */
5593: const PetscScalar *rtol;
5594: VecGetArrayRead(ts->vrtol,&rtol);
5596: for (i=0; i<n; i++) {
5597: SkipSmallValue(y[i],u[i],ts->adapt->ignore_max);
5598: diff = PetscAbsScalar(y[i] - u[i]);
5599: tola = ts->atol;
5600: tolr = PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5601: tol = tola+tolr;
5602: if(tola>0.){
5603: maxa = PetscMax(maxa,diff / tola);
5604: }
5605: if(tolr>0.){
5606: maxr = PetscMax(maxr,diff / tolr);
5607: }
5608: if(tol>0.){
5609: max = PetscMax(max,diff / tol);
5610: }
5611: }
5612: VecRestoreArrayRead(ts->vrtol,&rtol);
5613: } else { /* scalar atol, scalar rtol */
5615: for (i=0; i<n; i++) {
5616: SkipSmallValue(y[i],u[i],ts->adapt->ignore_max);
5617: diff = PetscAbsScalar(y[i] - u[i]);
5618: tola = ts->atol;
5619: tolr = ts->rtol * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5620: tol = tola+tolr;
5621: if(tola>0.){
5622: maxa = PetscMax(maxa,diff / tola);
5623: }
5624: if(tolr>0.){
5625: maxr = PetscMax(maxr,diff / tolr);
5626: }
5627: if(tol>0.){
5628: max = PetscMax(max,diff / tol);
5629: }
5630: }
5631: }
5632: VecRestoreArrayRead(U,&u);
5633: VecRestoreArrayRead(Y,&y);
5634: err_loc[0] = max;
5635: err_loc[1] = maxa;
5636: err_loc[2] = maxr;
5637: MPIU_Allreduce(err_loc,err_glb,3,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)ts));
5638: gmax = err_glb[0];
5639: gmaxa = err_glb[1];
5640: gmaxr = err_glb[2];
5642: *norm = gmax;
5643: *norma = gmaxa;
5644: *normr = gmaxr;
5645: if (PetscIsInfOrNanScalar(*norm)) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
5646: if (PetscIsInfOrNanScalar(*norma)) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_FP,"Infinite or not-a-number generated in norma");
5647: if (PetscIsInfOrNanScalar(*normr)) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_FP,"Infinite or not-a-number generated in normr");
5648: return(0);
5649: }
5651: /*@
5652: TSErrorWeightedNorm - compute a weighted norm of the difference between two state vectors based on supplied absolute and relative tolerances
5654: Collective on TS
5656: Input Arguments:
5657: + ts - time stepping context
5658: . U - state vector, usually ts->vec_sol
5659: . Y - state vector to be compared to U
5660: - wnormtype - norm type, either NORM_2 or NORM_INFINITY
5662: Output Arguments:
5663: + norm - weighted norm, a value of 1.0 achieves a balance between absolute and relative tolerances
5664: . norma - weighted norm, a value of 1.0 means that the error meets the absolute tolerance set by the user
5665: - normr - weighted norm, a value of 1.0 means that the error meets the relative tolerance set by the user
5667: Options Database Keys:
5668: . -ts_adapt_wnormtype <wnormtype> - 2, INFINITY
5670: Level: developer
5672: .seealso: TSErrorWeightedNormInfinity(), TSErrorWeightedNorm2(), TSErrorWeightedENorm
5673: @*/
5674: PetscErrorCode TSErrorWeightedNorm(TS ts,Vec U,Vec Y,NormType wnormtype,PetscReal *norm,PetscReal *norma,PetscReal *normr)
5675: {
5679: if (wnormtype == NORM_2) {
5680: TSErrorWeightedNorm2(ts,U,Y,norm,norma,normr);
5681: } else if(wnormtype == NORM_INFINITY) {
5682: TSErrorWeightedNormInfinity(ts,U,Y,norm,norma,normr);
5683: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for norm type %s",NormTypes[wnormtype]);
5684: return(0);
5685: }
5688: /*@
5689: TSErrorWeightedENorm2 - compute a weighted 2 error norm based on supplied absolute and relative tolerances
5691: Collective on TS
5693: Input Arguments:
5694: + ts - time stepping context
5695: . E - error vector
5696: . U - state vector, usually ts->vec_sol
5697: - Y - state vector, previous time step
5699: Output Arguments:
5700: + norm - weighted norm, a value of 1.0 means that the error matches the tolerances
5701: . norma - weighted norm based on the absolute tolerance, a value of 1.0 means that the error matches the tolerances
5702: - normr - weighted norm based on the relative tolerance, a value of 1.0 means that the error matches the tolerances
5704: Level: developer
5706: .seealso: TSErrorWeightedENorm(), TSErrorWeightedENormInfinity()
5707: @*/
5708: PetscErrorCode TSErrorWeightedENorm2(TS ts,Vec E,Vec U,Vec Y,PetscReal *norm,PetscReal *norma,PetscReal *normr)
5709: {
5710: PetscErrorCode ierr;
5711: PetscInt i,n,N,rstart;
5712: PetscInt n_loc,na_loc,nr_loc;
5713: PetscReal n_glb,na_glb,nr_glb;
5714: const PetscScalar *e,*u,*y;
5715: PetscReal err,sum,suma,sumr,gsum,gsuma,gsumr;
5716: PetscReal tol,tola,tolr;
5717: PetscReal err_loc[6],err_glb[6];
5733: VecGetSize(E,&N);
5734: VecGetLocalSize(E,&n);
5735: VecGetOwnershipRange(E,&rstart,NULL);
5736: VecGetArrayRead(E,&e);
5737: VecGetArrayRead(U,&u);
5738: VecGetArrayRead(Y,&y);
5739: sum = 0.; n_loc = 0;
5740: suma = 0.; na_loc = 0;
5741: sumr = 0.; nr_loc = 0;
5742: if (ts->vatol && ts->vrtol) {
5743: const PetscScalar *atol,*rtol;
5744: VecGetArrayRead(ts->vatol,&atol);
5745: VecGetArrayRead(ts->vrtol,&rtol);
5746: for (i=0; i<n; i++) {
5747: SkipSmallValue(y[i],u[i],ts->adapt->ignore_max);
5748: err = PetscAbsScalar(e[i]);
5749: tola = PetscRealPart(atol[i]);
5750: if(tola>0.){
5751: suma += PetscSqr(err/tola);
5752: na_loc++;
5753: }
5754: tolr = PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5755: if(tolr>0.){
5756: sumr += PetscSqr(err/tolr);
5757: nr_loc++;
5758: }
5759: tol=tola+tolr;
5760: if(tol>0.){
5761: sum += PetscSqr(err/tol);
5762: n_loc++;
5763: }
5764: }
5765: VecRestoreArrayRead(ts->vatol,&atol);
5766: VecRestoreArrayRead(ts->vrtol,&rtol);
5767: } else if (ts->vatol) { /* vector atol, scalar rtol */
5768: const PetscScalar *atol;
5769: VecGetArrayRead(ts->vatol,&atol);
5770: for (i=0; i<n; i++) {
5771: SkipSmallValue(y[i],u[i],ts->adapt->ignore_max);
5772: err = PetscAbsScalar(e[i]);
5773: tola = PetscRealPart(atol[i]);
5774: if(tola>0.){
5775: suma += PetscSqr(err/tola);
5776: na_loc++;
5777: }
5778: tolr = ts->rtol * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5779: if(tolr>0.){
5780: sumr += PetscSqr(err/tolr);
5781: nr_loc++;
5782: }
5783: tol=tola+tolr;
5784: if(tol>0.){
5785: sum += PetscSqr(err/tol);
5786: n_loc++;
5787: }
5788: }
5789: VecRestoreArrayRead(ts->vatol,&atol);
5790: } else if (ts->vrtol) { /* scalar atol, vector rtol */
5791: const PetscScalar *rtol;
5792: VecGetArrayRead(ts->vrtol,&rtol);
5793: for (i=0; i<n; i++) {
5794: SkipSmallValue(y[i],u[i],ts->adapt->ignore_max);
5795: err = PetscAbsScalar(e[i]);
5796: tola = ts->atol;
5797: if(tola>0.){
5798: suma += PetscSqr(err/tola);
5799: na_loc++;
5800: }
5801: tolr = PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5802: if(tolr>0.){
5803: sumr += PetscSqr(err/tolr);
5804: nr_loc++;
5805: }
5806: tol=tola+tolr;
5807: if(tol>0.){
5808: sum += PetscSqr(err/tol);
5809: n_loc++;
5810: }
5811: }
5812: VecRestoreArrayRead(ts->vrtol,&rtol);
5813: } else { /* scalar atol, scalar rtol */
5814: for (i=0; i<n; i++) {
5815: SkipSmallValue(y[i],u[i],ts->adapt->ignore_max);
5816: err = PetscAbsScalar(e[i]);
5817: tola = ts->atol;
5818: if(tola>0.){
5819: suma += PetscSqr(err/tola);
5820: na_loc++;
5821: }
5822: tolr = ts->rtol * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5823: if(tolr>0.){
5824: sumr += PetscSqr(err/tolr);
5825: nr_loc++;
5826: }
5827: tol=tola+tolr;
5828: if(tol>0.){
5829: sum += PetscSqr(err/tol);
5830: n_loc++;
5831: }
5832: }
5833: }
5834: VecRestoreArrayRead(E,&e);
5835: VecRestoreArrayRead(U,&u);
5836: VecRestoreArrayRead(Y,&y);
5838: err_loc[0] = sum;
5839: err_loc[1] = suma;
5840: err_loc[2] = sumr;
5841: err_loc[3] = (PetscReal)n_loc;
5842: err_loc[4] = (PetscReal)na_loc;
5843: err_loc[5] = (PetscReal)nr_loc;
5845: MPIU_Allreduce(err_loc,err_glb,6,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)ts));
5847: gsum = err_glb[0];
5848: gsuma = err_glb[1];
5849: gsumr = err_glb[2];
5850: n_glb = err_glb[3];
5851: na_glb = err_glb[4];
5852: nr_glb = err_glb[5];
5854: *norm = 0.;
5855: if(n_glb>0. ){*norm = PetscSqrtReal(gsum / n_glb );}
5856: *norma = 0.;
5857: if(na_glb>0.){*norma = PetscSqrtReal(gsuma / na_glb);}
5858: *normr = 0.;
5859: if(nr_glb>0.){*normr = PetscSqrtReal(gsumr / nr_glb);}
5861: if (PetscIsInfOrNanScalar(*norm)) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
5862: if (PetscIsInfOrNanScalar(*norma)) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_FP,"Infinite or not-a-number generated in norma");
5863: if (PetscIsInfOrNanScalar(*normr)) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_FP,"Infinite or not-a-number generated in normr");
5864: return(0);
5865: }
5867: /*@
5868: TSErrorWeightedENormInfinity - compute a weighted infinity error norm based on supplied absolute and relative tolerances
5869: Collective on TS
5871: Input Arguments:
5872: + ts - time stepping context
5873: . E - error vector
5874: . U - state vector, usually ts->vec_sol
5875: - Y - state vector, previous time step
5877: Output Arguments:
5878: + norm - weighted norm, a value of 1.0 means that the error matches the tolerances
5879: . norma - weighted norm based on the absolute tolerance, a value of 1.0 means that the error matches the tolerances
5880: - normr - weighted norm based on the relative tolerance, a value of 1.0 means that the error matches the tolerances
5882: Level: developer
5884: .seealso: TSErrorWeightedENorm(), TSErrorWeightedENorm2()
5885: @*/
5886: PetscErrorCode TSErrorWeightedENormInfinity(TS ts,Vec E,Vec U,Vec Y,PetscReal *norm,PetscReal *norma,PetscReal *normr)
5887: {
5888: PetscErrorCode ierr;
5889: PetscInt i,n,N,rstart;
5890: const PetscScalar *e,*u,*y;
5891: PetscReal err,max,gmax,maxa,gmaxa,maxr,gmaxr;
5892: PetscReal tol,tola,tolr;
5893: PetscReal err_loc[3],err_glb[3];
5909: VecGetSize(E,&N);
5910: VecGetLocalSize(E,&n);
5911: VecGetOwnershipRange(E,&rstart,NULL);
5912: VecGetArrayRead(E,&e);
5913: VecGetArrayRead(U,&u);
5914: VecGetArrayRead(Y,&y);
5916: max=0.;
5917: maxa=0.;
5918: maxr=0.;
5920: if (ts->vatol && ts->vrtol) { /* vector atol, vector rtol */
5921: const PetscScalar *atol,*rtol;
5922: VecGetArrayRead(ts->vatol,&atol);
5923: VecGetArrayRead(ts->vrtol,&rtol);
5925: for (i=0; i<n; i++) {
5926: SkipSmallValue(y[i],u[i],ts->adapt->ignore_max);
5927: err = PetscAbsScalar(e[i]);
5928: tola = PetscRealPart(atol[i]);
5929: tolr = PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5930: tol = tola+tolr;
5931: if(tola>0.){
5932: maxa = PetscMax(maxa,err / tola);
5933: }
5934: if(tolr>0.){
5935: maxr = PetscMax(maxr,err / tolr);
5936: }
5937: if(tol>0.){
5938: max = PetscMax(max,err / tol);
5939: }
5940: }
5941: VecRestoreArrayRead(ts->vatol,&atol);
5942: VecRestoreArrayRead(ts->vrtol,&rtol);
5943: } else if (ts->vatol) { /* vector atol, scalar rtol */
5944: const PetscScalar *atol;
5945: VecGetArrayRead(ts->vatol,&atol);
5946: for (i=0; i<n; i++) {
5947: SkipSmallValue(y[i],u[i],ts->adapt->ignore_max);
5948: err = PetscAbsScalar(e[i]);
5949: tola = PetscRealPart(atol[i]);
5950: tolr = ts->rtol * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5951: tol = tola+tolr;
5952: if(tola>0.){
5953: maxa = PetscMax(maxa,err / tola);
5954: }
5955: if(tolr>0.){
5956: maxr = PetscMax(maxr,err / tolr);
5957: }
5958: if(tol>0.){
5959: max = PetscMax(max,err / tol);
5960: }
5961: }
5962: VecRestoreArrayRead(ts->vatol,&atol);
5963: } else if (ts->vrtol) { /* scalar atol, vector rtol */
5964: const PetscScalar *rtol;
5965: VecGetArrayRead(ts->vrtol,&rtol);
5967: for (i=0; i<n; i++) {
5968: SkipSmallValue(y[i],u[i],ts->adapt->ignore_max);
5969: err = PetscAbsScalar(e[i]);
5970: tola = ts->atol;
5971: tolr = PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5972: tol = tola+tolr;
5973: if(tola>0.){
5974: maxa = PetscMax(maxa,err / tola);
5975: }
5976: if(tolr>0.){
5977: maxr = PetscMax(maxr,err / tolr);
5978: }
5979: if(tol>0.){
5980: max = PetscMax(max,err / tol);
5981: }
5982: }
5983: VecRestoreArrayRead(ts->vrtol,&rtol);
5984: } else { /* scalar atol, scalar rtol */
5986: for (i=0; i<n; i++) {
5987: SkipSmallValue(y[i],u[i],ts->adapt->ignore_max);
5988: err = PetscAbsScalar(e[i]);
5989: tola = ts->atol;
5990: tolr = ts->rtol * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5991: tol = tola+tolr;
5992: if(tola>0.){
5993: maxa = PetscMax(maxa,err / tola);
5994: }
5995: if(tolr>0.){
5996: maxr = PetscMax(maxr,err / tolr);
5997: }
5998: if(tol>0.){
5999: max = PetscMax(max,err / tol);
6000: }
6001: }
6002: }
6003: VecRestoreArrayRead(E,&e);
6004: VecRestoreArrayRead(U,&u);
6005: VecRestoreArrayRead(Y,&y);
6006: err_loc[0] = max;
6007: err_loc[1] = maxa;
6008: err_loc[2] = maxr;
6009: MPIU_Allreduce(err_loc,err_glb,3,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)ts));
6010: gmax = err_glb[0];
6011: gmaxa = err_glb[1];
6012: gmaxr = err_glb[2];
6014: *norm = gmax;
6015: *norma = gmaxa;
6016: *normr = gmaxr;
6017: if (PetscIsInfOrNanScalar(*norm)) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
6018: if (PetscIsInfOrNanScalar(*norma)) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_FP,"Infinite or not-a-number generated in norma");
6019: if (PetscIsInfOrNanScalar(*normr)) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_FP,"Infinite or not-a-number generated in normr");
6020: return(0);
6021: }
6023: /*@
6024: TSErrorWeightedENorm - compute a weighted error norm based on supplied absolute and relative tolerances
6026: Collective on TS
6028: Input Arguments:
6029: + ts - time stepping context
6030: . E - error vector
6031: . U - state vector, usually ts->vec_sol
6032: . Y - state vector, previous time step
6033: - wnormtype - norm type, either NORM_2 or NORM_INFINITY
6035: Output Arguments:
6036: + norm - weighted norm, a value of 1.0 achieves a balance between absolute and relative tolerances
6037: . norma - weighted norm, a value of 1.0 means that the error meets the absolute tolerance set by the user
6038: - normr - weighted norm, a value of 1.0 means that the error meets the relative tolerance set by the user
6040: Options Database Keys:
6041: . -ts_adapt_wnormtype <wnormtype> - 2, INFINITY
6043: Level: developer
6045: .seealso: TSErrorWeightedENormInfinity(), TSErrorWeightedENorm2(), TSErrorWeightedNormInfinity(), TSErrorWeightedNorm2()
6046: @*/
6047: PetscErrorCode TSErrorWeightedENorm(TS ts,Vec E,Vec U,Vec Y,NormType wnormtype,PetscReal *norm,PetscReal *norma,PetscReal *normr)
6048: {
6052: if (wnormtype == NORM_2) {
6053: TSErrorWeightedENorm2(ts,E,U,Y,norm,norma,normr);
6054: } else if(wnormtype == NORM_INFINITY) {
6055: TSErrorWeightedENormInfinity(ts,E,U,Y,norm,norma,normr);
6056: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for norm type %s",NormTypes[wnormtype]);
6057: return(0);
6058: }
6061: /*@
6062: TSSetCFLTimeLocal - Set the local CFL constraint relative to forward Euler
6064: Logically Collective on TS
6066: Input Arguments:
6067: + ts - time stepping context
6068: - cfltime - maximum stable time step if using forward Euler (value can be different on each process)
6070: Note:
6071: After calling this function, the global CFL time can be obtained by calling TSGetCFLTime()
6073: Level: intermediate
6075: .seealso: TSGetCFLTime(), TSADAPTCFL
6076: @*/
6077: PetscErrorCode TSSetCFLTimeLocal(TS ts,PetscReal cfltime)
6078: {
6081: ts->cfltime_local = cfltime;
6082: ts->cfltime = -1.;
6083: return(0);
6084: }
6086: /*@
6087: TSGetCFLTime - Get the maximum stable time step according to CFL criteria applied to forward Euler
6089: Collective on TS
6091: Input Arguments:
6092: . ts - time stepping context
6094: Output Arguments:
6095: . cfltime - maximum stable time step for forward Euler
6097: Level: advanced
6099: .seealso: TSSetCFLTimeLocal()
6100: @*/
6101: PetscErrorCode TSGetCFLTime(TS ts,PetscReal *cfltime)
6102: {
6106: if (ts->cfltime < 0) {
6107: MPIU_Allreduce(&ts->cfltime_local,&ts->cfltime,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)ts));
6108: }
6109: *cfltime = ts->cfltime;
6110: return(0);
6111: }
6113: /*@
6114: TSVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu
6116: Input Parameters:
6117: + ts - the TS context.
6118: . xl - lower bound.
6119: - xu - upper bound.
6121: Notes:
6122: If this routine is not called then the lower and upper bounds are set to
6123: PETSC_NINFINITY and PETSC_INFINITY respectively during SNESSetUp().
6125: Level: advanced
6127: @*/
6128: PetscErrorCode TSVISetVariableBounds(TS ts, Vec xl, Vec xu)
6129: {
6131: SNES snes;
6134: TSGetSNES(ts,&snes);
6135: SNESVISetVariableBounds(snes,xl,xu);
6136: return(0);
6137: }
6139: /*@C
6140: TSMonitorLGSolution - Monitors progress of the TS solvers by plotting each component of the solution vector
6141: in a time based line graph
6143: Collective on TS
6145: Input Parameters:
6146: + ts - the TS context
6147: . step - current time-step
6148: . ptime - current time
6149: . u - current solution
6150: - dctx - the TSMonitorLGCtx object that contains all the options for the monitoring, this is created with TSMonitorLGCtxCreate()
6152: Options Database:
6153: . -ts_monitor_lg_solution_variables
6155: Level: intermediate
6157: Notes:
6158: Each process in a parallel run displays its component solutions in a separate window
6160: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorLGCtxCreate(), TSMonitorLGCtxSetVariableNames(), TSMonitorLGCtxGetVariableNames(),
6161: TSMonitorLGSetVariableNames(), TSMonitorLGGetVariableNames(), TSMonitorLGSetDisplayVariables(), TSMonitorLGCtxSetDisplayVariables(),
6162: TSMonitorLGCtxSetTransform(), TSMonitorLGSetTransform(), TSMonitorLGError(), TSMonitorLGSNESIterations(), TSMonitorLGKSPIterations(),
6163: TSMonitorEnvelopeCtxCreate(), TSMonitorEnvelopeGetBounds(), TSMonitorEnvelopeCtxDestroy(), TSMonitorEnvelop()
6164: @*/
6165: PetscErrorCode TSMonitorLGSolution(TS ts,PetscInt step,PetscReal ptime,Vec u,void *dctx)
6166: {
6167: PetscErrorCode ierr;
6168: TSMonitorLGCtx ctx = (TSMonitorLGCtx)dctx;
6169: const PetscScalar *yy;
6170: Vec v;
6173: if (step < 0) return(0); /* -1 indicates interpolated solution */
6174: if (!step) {
6175: PetscDrawAxis axis;
6176: PetscInt dim;
6177: PetscDrawLGGetAxis(ctx->lg,&axis);
6178: PetscDrawAxisSetLabels(axis,"Solution as function of time","Time","Solution");
6179: if (!ctx->names) {
6180: PetscBool flg;
6181: /* user provides names of variables to plot but no names has been set so assume names are integer values */
6182: PetscOptionsHasName(((PetscObject)ts)->options,((PetscObject)ts)->prefix,"-ts_monitor_lg_solution_variables",&flg);
6183: if (flg) {
6184: PetscInt i,n;
6185: char **names;
6186: VecGetSize(u,&n);
6187: PetscMalloc1(n+1,&names);
6188: for (i=0; i<n; i++) {
6189: PetscMalloc1(5,&names[i]);
6190: PetscSNPrintf(names[i],5,"%D",i);
6191: }
6192: names[n] = NULL;
6193: ctx->names = names;
6194: }
6195: }
6196: if (ctx->names && !ctx->displaynames) {
6197: char **displaynames;
6198: PetscBool flg;
6199: VecGetLocalSize(u,&dim);
6200: PetscCalloc1(dim+1,&displaynames);
6201: PetscOptionsGetStringArray(((PetscObject)ts)->options,((PetscObject)ts)->prefix,"-ts_monitor_lg_solution_variables",displaynames,&dim,&flg);
6202: if (flg) {
6203: TSMonitorLGCtxSetDisplayVariables(ctx,(const char *const *)displaynames);
6204: }
6205: PetscStrArrayDestroy(&displaynames);
6206: }
6207: if (ctx->displaynames) {
6208: PetscDrawLGSetDimension(ctx->lg,ctx->ndisplayvariables);
6209: PetscDrawLGSetLegend(ctx->lg,(const char *const *)ctx->displaynames);
6210: } else if (ctx->names) {
6211: VecGetLocalSize(u,&dim);
6212: PetscDrawLGSetDimension(ctx->lg,dim);
6213: PetscDrawLGSetLegend(ctx->lg,(const char *const *)ctx->names);
6214: } else {
6215: VecGetLocalSize(u,&dim);
6216: PetscDrawLGSetDimension(ctx->lg,dim);
6217: }
6218: PetscDrawLGReset(ctx->lg);
6219: }
6221: if (!ctx->transform) v = u;
6222: else {(*ctx->transform)(ctx->transformctx,u,&v);}
6223: VecGetArrayRead(v,&yy);
6224: if (ctx->displaynames) {
6225: PetscInt i;
6226: for (i=0; i<ctx->ndisplayvariables; i++)
6227: ctx->displayvalues[i] = PetscRealPart(yy[ctx->displayvariables[i]]);
6228: PetscDrawLGAddCommonPoint(ctx->lg,ptime,ctx->displayvalues);
6229: } else {
6230: #if defined(PETSC_USE_COMPLEX)
6231: PetscInt i,n;
6232: PetscReal *yreal;
6233: VecGetLocalSize(v,&n);
6234: PetscMalloc1(n,&yreal);
6235: for (i=0; i<n; i++) yreal[i] = PetscRealPart(yy[i]);
6236: PetscDrawLGAddCommonPoint(ctx->lg,ptime,yreal);
6237: PetscFree(yreal);
6238: #else
6239: PetscDrawLGAddCommonPoint(ctx->lg,ptime,yy);
6240: #endif
6241: }
6242: VecRestoreArrayRead(v,&yy);
6243: if (ctx->transform) {VecDestroy(&v);}
6245: if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
6246: PetscDrawLGDraw(ctx->lg);
6247: PetscDrawLGSave(ctx->lg);
6248: }
6249: return(0);
6250: }
6252: /*@C
6253: TSMonitorLGSetVariableNames - Sets the name of each component in the solution vector so that it may be displayed in the plot
6255: Collective on TS
6257: Input Parameters:
6258: + ts - the TS context
6259: - names - the names of the components, final string must be NULL
6261: Level: intermediate
6263: Notes:
6264: If the TS object does not have a TSMonitorLGCtx associated with it then this function is ignored
6266: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorLGSetDisplayVariables(), TSMonitorLGCtxSetVariableNames()
6267: @*/
6268: PetscErrorCode TSMonitorLGSetVariableNames(TS ts,const char * const *names)
6269: {
6270: PetscErrorCode ierr;
6271: PetscInt i;
6274: for (i=0; i<ts->numbermonitors; i++) {
6275: if (ts->monitor[i] == TSMonitorLGSolution) {
6276: TSMonitorLGCtxSetVariableNames((TSMonitorLGCtx)ts->monitorcontext[i],names);
6277: break;
6278: }
6279: }
6280: return(0);
6281: }
6283: /*@C
6284: TSMonitorLGCtxSetVariableNames - Sets the name of each component in the solution vector so that it may be displayed in the plot
6286: Collective on TS
6288: Input Parameters:
6289: + ts - the TS context
6290: - names - the names of the components, final string must be NULL
6292: Level: intermediate
6294: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorLGSetDisplayVariables(), TSMonitorLGSetVariableNames()
6295: @*/
6296: PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx ctx,const char * const *names)
6297: {
6298: PetscErrorCode ierr;
6301: PetscStrArrayDestroy(&ctx->names);
6302: PetscStrArrayallocpy(names,&ctx->names);
6303: return(0);
6304: }
6306: /*@C
6307: TSMonitorLGGetVariableNames - Gets the name of each component in the solution vector so that it may be displayed in the plot
6309: Collective on TS
6311: Input Parameter:
6312: . ts - the TS context
6314: Output Parameter:
6315: . names - the names of the components, final string must be NULL
6317: Level: intermediate
6319: Notes:
6320: If the TS object does not have a TSMonitorLGCtx associated with it then this function is ignored
6322: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorLGSetDisplayVariables()
6323: @*/
6324: PetscErrorCode TSMonitorLGGetVariableNames(TS ts,const char *const **names)
6325: {
6326: PetscInt i;
6329: *names = NULL;
6330: for (i=0; i<ts->numbermonitors; i++) {
6331: if (ts->monitor[i] == TSMonitorLGSolution) {
6332: TSMonitorLGCtx ctx = (TSMonitorLGCtx) ts->monitorcontext[i];
6333: *names = (const char *const *)ctx->names;
6334: break;
6335: }
6336: }
6337: return(0);
6338: }
6340: /*@C
6341: TSMonitorLGCtxSetDisplayVariables - Sets the variables that are to be display in the monitor
6343: Collective on TS
6345: Input Parameters:
6346: + ctx - the TSMonitorLG context
6347: - displaynames - the names of the components, final string must be NULL
6349: Level: intermediate
6351: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorLGSetVariableNames()
6352: @*/
6353: PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx ctx,const char * const *displaynames)
6354: {
6355: PetscInt j = 0,k;
6356: PetscErrorCode ierr;
6359: if (!ctx->names) return(0);
6360: PetscStrArrayDestroy(&ctx->displaynames);
6361: PetscStrArrayallocpy(displaynames,&ctx->displaynames);
6362: while (displaynames[j]) j++;
6363: ctx->ndisplayvariables = j;
6364: PetscMalloc1(ctx->ndisplayvariables,&ctx->displayvariables);
6365: PetscMalloc1(ctx->ndisplayvariables,&ctx->displayvalues);
6366: j = 0;
6367: while (displaynames[j]) {
6368: k = 0;
6369: while (ctx->names[k]) {
6370: PetscBool flg;
6371: PetscStrcmp(displaynames[j],ctx->names[k],&flg);
6372: if (flg) {
6373: ctx->displayvariables[j] = k;
6374: break;
6375: }
6376: k++;
6377: }
6378: j++;
6379: }
6380: return(0);
6381: }
6383: /*@C
6384: TSMonitorLGSetDisplayVariables - Sets the variables that are to be display in the monitor
6386: Collective on TS
6388: Input Parameters:
6389: + ts - the TS context
6390: - displaynames - the names of the components, final string must be NULL
6392: Notes:
6393: If the TS object does not have a TSMonitorLGCtx associated with it then this function is ignored
6395: Level: intermediate
6397: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorLGSetVariableNames()
6398: @*/
6399: PetscErrorCode TSMonitorLGSetDisplayVariables(TS ts,const char * const *displaynames)
6400: {
6401: PetscInt i;
6402: PetscErrorCode ierr;
6405: for (i=0; i<ts->numbermonitors; i++) {
6406: if (ts->monitor[i] == TSMonitorLGSolution) {
6407: TSMonitorLGCtxSetDisplayVariables((TSMonitorLGCtx)ts->monitorcontext[i],displaynames);
6408: break;
6409: }
6410: }
6411: return(0);
6412: }
6414: /*@C
6415: TSMonitorLGSetTransform - Solution vector will be transformed by provided function before being displayed
6417: Collective on TS
6419: Input Parameters:
6420: + ts - the TS context
6421: . transform - the transform function
6422: . destroy - function to destroy the optional context
6423: - ctx - optional context used by transform function
6425: Notes:
6426: If the TS object does not have a TSMonitorLGCtx associated with it then this function is ignored
6428: Level: intermediate
6430: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorLGSetVariableNames(), TSMonitorLGCtxSetTransform()
6431: @*/
6432: PetscErrorCode TSMonitorLGSetTransform(TS ts,PetscErrorCode (*transform)(void*,Vec,Vec*),PetscErrorCode (*destroy)(void*),void *tctx)
6433: {
6434: PetscInt i;
6435: PetscErrorCode ierr;
6438: for (i=0; i<ts->numbermonitors; i++) {
6439: if (ts->monitor[i] == TSMonitorLGSolution) {
6440: TSMonitorLGCtxSetTransform((TSMonitorLGCtx)ts->monitorcontext[i],transform,destroy,tctx);
6441: }
6442: }
6443: return(0);
6444: }
6446: /*@C
6447: TSMonitorLGCtxSetTransform - Solution vector will be transformed by provided function before being displayed
6449: Collective on TSLGCtx
6451: Input Parameters:
6452: + ts - the TS context
6453: . transform - the transform function
6454: . destroy - function to destroy the optional context
6455: - ctx - optional context used by transform function
6457: Level: intermediate
6459: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorLGSetVariableNames(), TSMonitorLGSetTransform()
6460: @*/
6461: PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx ctx,PetscErrorCode (*transform)(void*,Vec,Vec*),PetscErrorCode (*destroy)(void*),void *tctx)
6462: {
6464: ctx->transform = transform;
6465: ctx->transformdestroy = destroy;
6466: ctx->transformctx = tctx;
6467: return(0);
6468: }
6470: /*@C
6471: TSMonitorLGError - Monitors progress of the TS solvers by plotting each component of the error
6472: in a time based line graph
6474: Collective on TS
6476: Input Parameters:
6477: + ts - the TS context
6478: . step - current time-step
6479: . ptime - current time
6480: . u - current solution
6481: - dctx - TSMonitorLGCtx object created with TSMonitorLGCtxCreate()
6483: Level: intermediate
6485: Notes:
6486: Each process in a parallel run displays its component errors in a separate window
6488: The user must provide the solution using TSSetSolutionFunction() to use this monitor.
6490: Options Database Keys:
6491: . -ts_monitor_lg_error - create a graphical monitor of error history
6493: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSSetSolutionFunction()
6494: @*/
6495: PetscErrorCode TSMonitorLGError(TS ts,PetscInt step,PetscReal ptime,Vec u,void *dummy)
6496: {
6497: PetscErrorCode ierr;
6498: TSMonitorLGCtx ctx = (TSMonitorLGCtx)dummy;
6499: const PetscScalar *yy;
6500: Vec y;
6503: if (!step) {
6504: PetscDrawAxis axis;
6505: PetscInt dim;
6506: PetscDrawLGGetAxis(ctx->lg,&axis);
6507: PetscDrawAxisSetLabels(axis,"Error in solution as function of time","Time","Error");
6508: VecGetLocalSize(u,&dim);
6509: PetscDrawLGSetDimension(ctx->lg,dim);
6510: PetscDrawLGReset(ctx->lg);
6511: }
6512: VecDuplicate(u,&y);
6513: TSComputeSolutionFunction(ts,ptime,y);
6514: VecAXPY(y,-1.0,u);
6515: VecGetArrayRead(y,&yy);
6516: #if defined(PETSC_USE_COMPLEX)
6517: {
6518: PetscReal *yreal;
6519: PetscInt i,n;
6520: VecGetLocalSize(y,&n);
6521: PetscMalloc1(n,&yreal);
6522: for (i=0; i<n; i++) yreal[i] = PetscRealPart(yy[i]);
6523: PetscDrawLGAddCommonPoint(ctx->lg,ptime,yreal);
6524: PetscFree(yreal);
6525: }
6526: #else
6527: PetscDrawLGAddCommonPoint(ctx->lg,ptime,yy);
6528: #endif
6529: VecRestoreArrayRead(y,&yy);
6530: VecDestroy(&y);
6531: if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
6532: PetscDrawLGDraw(ctx->lg);
6533: PetscDrawLGSave(ctx->lg);
6534: }
6535: return(0);
6536: }
6538: /*@C
6539: TSMonitorSPSwarmSolution - Graphically displays phase plots of DMSwarm particles on a scatter plot
6541: Input Parameters:
6542: + ts - the TS context
6543: . step - current time-step
6544: . ptime - current time
6545: . u - current solution
6546: - dctx - the TSMonitorSPCtx object that contains all the options for the monitoring, this is created with TSMonitorSPCtxCreate()
6548: Options Database:
6549: . -ts_monitor_sp_swarm
6551: Level: intermediate
6553: @*/
6554: PetscErrorCode TSMonitorSPSwarmSolution(TS ts,PetscInt step,PetscReal ptime,Vec u,void *dctx)
6555: {
6556: PetscErrorCode ierr;
6557: TSMonitorSPCtx ctx = (TSMonitorSPCtx)dctx;
6558: const PetscScalar *yy;
6559: PetscReal *y,*x;
6560: PetscInt Np, p, dim=2;
6561: DM dm;
6564:
6565: if (step < 0) return(0); /* -1 indicates interpolated solution */
6566: if (!step) {
6567: PetscDrawAxis axis;
6568: PetscDrawSPGetAxis(ctx->sp,&axis);
6569: PetscDrawAxisSetLabels(axis,"Particles","X","Y");
6570: PetscDrawAxisSetLimits(axis, -5, 5, -5, 5);
6571: PetscDrawAxisSetHoldLimits(axis, PETSC_TRUE);
6572: TSGetDM(ts, &dm);
6573: DMGetDimension(dm, &dim);
6574: if(dim!=2) SETERRQ(PETSC_COMM_SELF, ierr, "Dimensions improper for monitor arguments! Current support: two dimensions.");
6575: VecGetLocalSize(u, &Np);
6576: Np /= 2*dim;
6577: PetscDrawSPSetDimension(ctx->sp, Np);
6578: PetscDrawSPReset(ctx->sp);
6579: }
6580:
6581: VecGetLocalSize(u, &Np);
6582: Np /= 2*dim;
6583: VecGetArrayRead(u,&yy);
6584: PetscMalloc2(Np, &x, Np, &y);
6585: /* get points from solution vector */
6586: for (p=0; p<Np; ++p){
6587: x[p] = PetscRealPart(yy[2*dim*p]);
6588: y[p] = PetscRealPart(yy[2*dim*p+1]);
6589: }
6590: VecRestoreArrayRead(u,&yy);
6591:
6592: if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
6593: PetscDrawSPAddPoint(ctx->sp,x,y);
6594: PetscDrawSPDraw(ctx->sp,PETSC_FALSE);
6595: PetscDrawSPSave(ctx->sp);
6596: }
6598: PetscFree2(x, y);
6600: return(0);
6601: }
6605: /*@C
6606: TSMonitorError - Monitors progress of the TS solvers by printing the 2 norm of the error at each timestep
6608: Collective on TS
6610: Input Parameters:
6611: + ts - the TS context
6612: . step - current time-step
6613: . ptime - current time
6614: . u - current solution
6615: - dctx - unused context
6617: Level: intermediate
6619: The user must provide the solution using TSSetSolutionFunction() to use this monitor.
6621: Options Database Keys:
6622: . -ts_monitor_error - create a graphical monitor of error history
6624: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSSetSolutionFunction()
6625: @*/
6626: PetscErrorCode TSMonitorError(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscViewerAndFormat *vf)
6627: {
6628: PetscErrorCode ierr;
6629: Vec y;
6630: PetscReal nrm;
6631: PetscBool flg;
6634: VecDuplicate(u,&y);
6635: TSComputeSolutionFunction(ts,ptime,y);
6636: VecAXPY(y,-1.0,u);
6637: PetscObjectTypeCompare((PetscObject)vf->viewer,PETSCVIEWERASCII,&flg);
6638: if (flg) {
6639: VecNorm(y,NORM_2,&nrm);
6640: PetscViewerASCIIPrintf(vf->viewer,"2-norm of error %g\n",(double)nrm);
6641: }
6642: PetscObjectTypeCompare((PetscObject)vf->viewer,PETSCVIEWERDRAW,&flg);
6643: if (flg) {
6644: VecView(y,vf->viewer);
6645: }
6646: VecDestroy(&y);
6647: return(0);
6648: }
6650: PetscErrorCode TSMonitorLGSNESIterations(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx)
6651: {
6652: TSMonitorLGCtx ctx = (TSMonitorLGCtx) monctx;
6653: PetscReal x = ptime,y;
6655: PetscInt its;
6658: if (n < 0) return(0); /* -1 indicates interpolated solution */
6659: if (!n) {
6660: PetscDrawAxis axis;
6661: PetscDrawLGGetAxis(ctx->lg,&axis);
6662: PetscDrawAxisSetLabels(axis,"Nonlinear iterations as function of time","Time","SNES Iterations");
6663: PetscDrawLGReset(ctx->lg);
6664: ctx->snes_its = 0;
6665: }
6666: TSGetSNESIterations(ts,&its);
6667: y = its - ctx->snes_its;
6668: PetscDrawLGAddPoint(ctx->lg,&x,&y);
6669: if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) {
6670: PetscDrawLGDraw(ctx->lg);
6671: PetscDrawLGSave(ctx->lg);
6672: }
6673: ctx->snes_its = its;
6674: return(0);
6675: }
6677: PetscErrorCode TSMonitorLGKSPIterations(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx)
6678: {
6679: TSMonitorLGCtx ctx = (TSMonitorLGCtx) monctx;
6680: PetscReal x = ptime,y;
6682: PetscInt its;
6685: if (n < 0) return(0); /* -1 indicates interpolated solution */
6686: if (!n) {
6687: PetscDrawAxis axis;
6688: PetscDrawLGGetAxis(ctx->lg,&axis);
6689: PetscDrawAxisSetLabels(axis,"Linear iterations as function of time","Time","KSP Iterations");
6690: PetscDrawLGReset(ctx->lg);
6691: ctx->ksp_its = 0;
6692: }
6693: TSGetKSPIterations(ts,&its);
6694: y = its - ctx->ksp_its;
6695: PetscDrawLGAddPoint(ctx->lg,&x,&y);
6696: if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) {
6697: PetscDrawLGDraw(ctx->lg);
6698: PetscDrawLGSave(ctx->lg);
6699: }
6700: ctx->ksp_its = its;
6701: return(0);
6702: }
6704: /*@
6705: TSComputeLinearStability - computes the linear stability function at a point
6707: Collective on TS
6709: Input Parameters:
6710: + ts - the TS context
6711: - xr,xi - real and imaginary part of input arguments
6713: Output Parameters:
6714: . yr,yi - real and imaginary part of function value
6716: Level: developer
6718: .seealso: TSSetRHSFunction(), TSComputeIFunction()
6719: @*/
6720: PetscErrorCode TSComputeLinearStability(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
6721: {
6726: if (!ts->ops->linearstability) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Linearized stability function not provided for this method");
6727: (*ts->ops->linearstability)(ts,xr,xi,yr,yi);
6728: return(0);
6729: }
6731: /* ------------------------------------------------------------------------*/
6732: /*@C
6733: TSMonitorEnvelopeCtxCreate - Creates a context for use with TSMonitorEnvelope()
6735: Collective on TS
6737: Input Parameters:
6738: . ts - the ODE solver object
6740: Output Parameter:
6741: . ctx - the context
6743: Level: intermediate
6745: .seealso: TSMonitorLGTimeStep(), TSMonitorSet(), TSMonitorLGSolution(), TSMonitorLGError()
6747: @*/
6748: PetscErrorCode TSMonitorEnvelopeCtxCreate(TS ts,TSMonitorEnvelopeCtx *ctx)
6749: {
6753: PetscNew(ctx);
6754: return(0);
6755: }
6757: /*@C
6758: TSMonitorEnvelope - Monitors the maximum and minimum value of each component of the solution
6760: Collective on TS
6762: Input Parameters:
6763: + ts - the TS context
6764: . step - current time-step
6765: . ptime - current time
6766: . u - current solution
6767: - dctx - the envelope context
6769: Options Database:
6770: . -ts_monitor_envelope
6772: Level: intermediate
6774: Notes:
6775: after a solve you can use TSMonitorEnvelopeGetBounds() to access the envelope
6777: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorEnvelopeGetBounds(), TSMonitorEnvelopeCtxCreate()
6778: @*/
6779: PetscErrorCode TSMonitorEnvelope(TS ts,PetscInt step,PetscReal ptime,Vec u,void *dctx)
6780: {
6781: PetscErrorCode ierr;
6782: TSMonitorEnvelopeCtx ctx = (TSMonitorEnvelopeCtx)dctx;
6785: if (!ctx->max) {
6786: VecDuplicate(u,&ctx->max);
6787: VecDuplicate(u,&ctx->min);
6788: VecCopy(u,ctx->max);
6789: VecCopy(u,ctx->min);
6790: } else {
6791: VecPointwiseMax(ctx->max,u,ctx->max);
6792: VecPointwiseMin(ctx->min,u,ctx->min);
6793: }
6794: return(0);
6795: }
6797: /*@C
6798: TSMonitorEnvelopeGetBounds - Gets the bounds for the components of the solution
6800: Collective on TS
6802: Input Parameter:
6803: . ts - the TS context
6805: Output Parameter:
6806: + max - the maximum values
6807: - min - the minimum values
6809: Notes:
6810: If the TS does not have a TSMonitorEnvelopeCtx associated with it then this function is ignored
6812: Level: intermediate
6814: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorLGSetDisplayVariables()
6815: @*/
6816: PetscErrorCode TSMonitorEnvelopeGetBounds(TS ts,Vec *max,Vec *min)
6817: {
6818: PetscInt i;
6821: if (max) *max = NULL;
6822: if (min) *min = NULL;
6823: for (i=0; i<ts->numbermonitors; i++) {
6824: if (ts->monitor[i] == TSMonitorEnvelope) {
6825: TSMonitorEnvelopeCtx ctx = (TSMonitorEnvelopeCtx) ts->monitorcontext[i];
6826: if (max) *max = ctx->max;
6827: if (min) *min = ctx->min;
6828: break;
6829: }
6830: }
6831: return(0);
6832: }
6834: /*@C
6835: TSMonitorEnvelopeCtxDestroy - Destroys a context that was created with TSMonitorEnvelopeCtxCreate().
6837: Collective on TSMonitorEnvelopeCtx
6839: Input Parameter:
6840: . ctx - the monitor context
6842: Level: intermediate
6844: .seealso: TSMonitorLGCtxCreate(), TSMonitorSet(), TSMonitorLGTimeStep()
6845: @*/
6846: PetscErrorCode TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx *ctx)
6847: {
6851: VecDestroy(&(*ctx)->min);
6852: VecDestroy(&(*ctx)->max);
6853: PetscFree(*ctx);
6854: return(0);
6855: }
6857: /*@
6858: TSRestartStep - Flags the solver to restart the next step
6860: Collective on TS
6862: Input Parameter:
6863: . ts - the TS context obtained from TSCreate()
6865: Level: advanced
6867: Notes:
6868: Multistep methods like BDF or Runge-Kutta methods with FSAL property require restarting the solver in the event of
6869: discontinuities. These discontinuities may be introduced as a consequence of explicitly modifications to the solution
6870: vector (which PETSc attempts to detect and handle) or problem coefficients (which PETSc is not able to detect). For
6871: the sake of correctness and maximum safety, users are expected to call TSRestart() whenever they introduce
6872: discontinuities in callback routines (e.g. prestep and poststep routines, or implicit/rhs function routines with
6873: discontinuous source terms).
6875: .seealso: TSSolve(), TSSetPreStep(), TSSetPostStep()
6876: @*/
6877: PetscErrorCode TSRestartStep(TS ts)
6878: {
6881: ts->steprestart = PETSC_TRUE;
6882: return(0);
6883: }
6885: /*@
6886: TSRollBack - Rolls back one time step
6888: Collective on TS
6890: Input Parameter:
6891: . ts - the TS context obtained from TSCreate()
6893: Level: advanced
6895: .seealso: TSCreate(), TSSetUp(), TSDestroy(), TSSolve(), TSSetPreStep(), TSSetPreStage(), TSInterpolate()
6896: @*/
6897: PetscErrorCode TSRollBack(TS ts)
6898: {
6903: if (ts->steprollback) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"TSRollBack already called");
6904: if (!ts->ops->rollback) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRollBack not implemented for type '%s'",((PetscObject)ts)->type_name);
6905: (*ts->ops->rollback)(ts);
6906: ts->time_step = ts->ptime - ts->ptime_prev;
6907: ts->ptime = ts->ptime_prev;
6908: ts->ptime_prev = ts->ptime_prev_rollback;
6909: ts->steps--;
6910: ts->steprollback = PETSC_TRUE;
6911: return(0);
6912: }
6914: /*@
6915: TSGetStages - Get the number of stages and stage values
6917: Input Parameter:
6918: . ts - the TS context obtained from TSCreate()
6920: Output Parameters:
6921: + ns - the number of stages
6922: - Y - the current stage vectors
6924: Level: advanced
6926: Notes: Both ns and Y can be NULL.
6928: .seealso: TSCreate()
6929: @*/
6930: PetscErrorCode TSGetStages(TS ts,PetscInt *ns,Vec **Y)
6931: {
6938: if (!ts->ops->getstages) {
6939: if (ns) *ns = 0;
6940: if (Y) *Y = NULL;
6941: } else {
6942: (*ts->ops->getstages)(ts,ns,Y);
6943: }
6944: return(0);
6945: }
6947: /*@C
6948: TSComputeIJacobianDefaultColor - Computes the Jacobian using finite differences and coloring to exploit matrix sparsity.
6950: Collective on SNES
6952: Input Parameters:
6953: + ts - the TS context
6954: . t - current timestep
6955: . U - state vector
6956: . Udot - time derivative of state vector
6957: . shift - shift to apply, see note below
6958: - ctx - an optional user context
6960: Output Parameters:
6961: + J - Jacobian matrix (not altered in this routine)
6962: - B - newly computed Jacobian matrix to use with preconditioner (generally the same as J)
6964: Level: intermediate
6966: Notes:
6967: If F(t,U,Udot)=0 is the DAE, the required Jacobian is
6969: dF/dU + shift*dF/dUdot
6971: Most users should not need to explicitly call this routine, as it
6972: is used internally within the nonlinear solvers.
6974: This will first try to get the coloring from the DM. If the DM type has no coloring
6975: routine, then it will try to get the coloring from the matrix. This requires that the
6976: matrix have nonzero entries precomputed.
6978: .seealso: TSSetIJacobian(), MatFDColoringCreate(), MatFDColoringSetFunction()
6979: @*/
6980: PetscErrorCode TSComputeIJacobianDefaultColor(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal shift,Mat J,Mat B,void *ctx)
6981: {
6982: SNES snes;
6983: MatFDColoring color;
6984: PetscBool hascolor, matcolor = PETSC_FALSE;
6988: PetscOptionsGetBool(((PetscObject)ts)->options,((PetscObject) ts)->prefix, "-ts_fd_color_use_mat", &matcolor, NULL);
6989: PetscObjectQuery((PetscObject) B, "TSMatFDColoring", (PetscObject *) &color);
6990: if (!color) {
6991: DM dm;
6992: ISColoring iscoloring;
6994: TSGetDM(ts, &dm);
6995: DMHasColoring(dm, &hascolor);
6996: if (hascolor && !matcolor) {
6997: DMCreateColoring(dm, IS_COLORING_GLOBAL, &iscoloring);
6998: MatFDColoringCreate(B, iscoloring, &color);
6999: MatFDColoringSetFunction(color, (PetscErrorCode (*)(void)) SNESTSFormFunction, (void *) ts);
7000: MatFDColoringSetFromOptions(color);
7001: MatFDColoringSetUp(B, iscoloring, color);
7002: ISColoringDestroy(&iscoloring);
7003: } else {
7004: MatColoring mc;
7006: MatColoringCreate(B, &mc);
7007: MatColoringSetDistance(mc, 2);
7008: MatColoringSetType(mc, MATCOLORINGSL);
7009: MatColoringSetFromOptions(mc);
7010: MatColoringApply(mc, &iscoloring);
7011: MatColoringDestroy(&mc);
7012: MatFDColoringCreate(B, iscoloring, &color);
7013: MatFDColoringSetFunction(color, (PetscErrorCode (*)(void)) SNESTSFormFunction, (void *) ts);
7014: MatFDColoringSetFromOptions(color);
7015: MatFDColoringSetUp(B, iscoloring, color);
7016: ISColoringDestroy(&iscoloring);
7017: }
7018: PetscObjectCompose((PetscObject) B, "TSMatFDColoring", (PetscObject) color);
7019: PetscObjectDereference((PetscObject) color);
7020: }
7021: TSGetSNES(ts, &snes);
7022: MatFDColoringApply(B, color, U, snes);
7023: if (J != B) {
7024: MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);
7025: MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);
7026: }
7027: return(0);
7028: }
7030: /*@
7031: TSSetFunctionDomainError - Set a function that tests if the current state vector is valid
7033: Input Parameters:
7034: + ts - the TS context
7035: - func - function called within TSFunctionDomainError
7037: Calling sequence of func:
7038: $ PetscErrorCode func(TS ts,PetscReal time,Vec state,PetscBool reject)
7040: + ts - the TS context
7041: . time - the current time (of the stage)
7042: . state - the state to check if it is valid
7043: - reject - (output parameter) PETSC_FALSE if the state is acceptable, PETSC_TRUE if not acceptable
7045: Level: intermediate
7047: Notes:
7048: If an implicit ODE solver is being used then, in addition to providing this routine, the
7049: user's code should call SNESSetFunctionDomainError() when domain errors occur during
7050: function evaluations where the functions are provided by TSSetIFunction() or TSSetRHSFunction().
7051: Use TSGetSNES() to obtain the SNES object
7053: Developer Notes:
7054: The naming of this function is inconsistent with the SNESSetFunctionDomainError()
7055: since one takes a function pointer and the other does not.
7057: .seealso: TSAdaptCheckStage(), TSFunctionDomainError(), SNESSetFunctionDomainError(), TSGetSNES()
7058: @*/
7060: PetscErrorCode TSSetFunctionDomainError(TS ts, PetscErrorCode (*func)(TS,PetscReal,Vec,PetscBool*))
7061: {
7064: ts->functiondomainerror = func;
7065: return(0);
7066: }
7068: /*@
7069: TSFunctionDomainError - Checks if the current state is valid
7071: Input Parameters:
7072: + ts - the TS context
7073: . stagetime - time of the simulation
7074: - Y - state vector to check.
7076: Output Parameter:
7077: . accept - Set to PETSC_FALSE if the current state vector is valid.
7079: Note:
7080: This function is called by the TS integration routines and calls the user provided function (set with TSSetFunctionDomainError())
7081: to check if the current state is valid.
7083: Level: developer
7085: .seealso: TSSetFunctionDomainError()
7086: @*/
7087: PetscErrorCode TSFunctionDomainError(TS ts,PetscReal stagetime,Vec Y,PetscBool* accept)
7088: {
7091: *accept = PETSC_TRUE;
7092: if (ts->functiondomainerror) {
7093: PetscStackCallStandard((*ts->functiondomainerror),(ts,stagetime,Y,accept));
7094: }
7095: return(0);
7096: }
7098: /*@C
7099: TSClone - This function clones a time step object.
7101: Collective
7103: Input Parameter:
7104: . tsin - The input TS
7106: Output Parameter:
7107: . tsout - The output TS (cloned)
7109: Notes:
7110: This function is used to create a clone of a TS object. It is used in ARKIMEX for initializing the slope for first stage explicit methods. It will likely be replaced in the future with a mechanism of switching methods on the fly.
7112: When using TSDestroy() on a clone the user has to first reset the correct TS reference in the embedded SNES object: e.g.: by running SNES snes_dup=NULL; TSGetSNES(ts,&snes_dup); TSSetSNES(ts,snes_dup);
7114: Level: developer
7116: .seealso: TSCreate(), TSSetType(), TSSetUp(), TSDestroy(), TSSetProblemType()
7117: @*/
7118: PetscErrorCode TSClone(TS tsin, TS *tsout)
7119: {
7120: TS t;
7122: SNES snes_start;
7123: DM dm;
7124: TSType type;
7128: *tsout = NULL;
7130: PetscHeaderCreate(t, TS_CLASSID, "TS", "Time stepping", "TS", PetscObjectComm((PetscObject)tsin), TSDestroy, TSView);
7132: /* General TS description */
7133: t->numbermonitors = 0;
7134: t->setupcalled = 0;
7135: t->ksp_its = 0;
7136: t->snes_its = 0;
7137: t->nwork = 0;
7138: t->rhsjacobian.time = PETSC_MIN_REAL;
7139: t->rhsjacobian.scale = 1.;
7140: t->ijacobian.shift = 1.;
7142: TSGetSNES(tsin,&snes_start);
7143: TSSetSNES(t,snes_start);
7145: TSGetDM(tsin,&dm);
7146: TSSetDM(t,dm);
7148: t->adapt = tsin->adapt;
7149: PetscObjectReference((PetscObject)t->adapt);
7151: t->trajectory = tsin->trajectory;
7152: PetscObjectReference((PetscObject)t->trajectory);
7154: t->event = tsin->event;
7155: if (t->event) t->event->refct++;
7157: t->problem_type = tsin->problem_type;
7158: t->ptime = tsin->ptime;
7159: t->ptime_prev = tsin->ptime_prev;
7160: t->time_step = tsin->time_step;
7161: t->max_time = tsin->max_time;
7162: t->steps = tsin->steps;
7163: t->max_steps = tsin->max_steps;
7164: t->equation_type = tsin->equation_type;
7165: t->atol = tsin->atol;
7166: t->rtol = tsin->rtol;
7167: t->max_snes_failures = tsin->max_snes_failures;
7168: t->max_reject = tsin->max_reject;
7169: t->errorifstepfailed = tsin->errorifstepfailed;
7171: TSGetType(tsin,&type);
7172: TSSetType(t,type);
7174: t->vec_sol = NULL;
7176: t->cfltime = tsin->cfltime;
7177: t->cfltime_local = tsin->cfltime_local;
7178: t->exact_final_time = tsin->exact_final_time;
7180: PetscMemcpy(t->ops,tsin->ops,sizeof(struct _TSOps));
7182: if (((PetscObject)tsin)->fortran_func_pointers) {
7183: PetscInt i;
7184: PetscMalloc((10)*sizeof(void(*)(void)),&((PetscObject)t)->fortran_func_pointers);
7185: for (i=0; i<10; i++) {
7186: ((PetscObject)t)->fortran_func_pointers[i] = ((PetscObject)tsin)->fortran_func_pointers[i];
7187: }
7188: }
7189: *tsout = t;
7190: return(0);
7191: }
7193: static PetscErrorCode RHSWrapperFunction_TSRHSJacobianTest(void* ctx,Vec x,Vec y)
7194: {
7196: TS ts = (TS) ctx;
7199: TSComputeRHSFunction(ts,0,x,y);
7200: return(0);
7201: }
7203: /*@
7204: TSRHSJacobianTest - Compares the multiply routine provided to the MATSHELL with differencing on the TS given RHS function.
7206: Logically Collective on TS
7208: Input Parameters:
7209: TS - the time stepping routine
7211: Output Parameter:
7212: . flg - PETSC_TRUE if the multiply is likely correct
7214: Options Database:
7215: . -ts_rhs_jacobian_test_mult -mat_shell_test_mult_view - run the test at each timestep of the integrator
7217: Level: advanced
7219: Notes:
7220: This only works for problems defined only the RHS function and Jacobian NOT IFunction and IJacobian
7222: .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMultTranspose(), TSRHSJacobianTestTranspose()
7223: @*/
7224: PetscErrorCode TSRHSJacobianTest(TS ts,PetscBool *flg)
7225: {
7226: Mat J,B;
7228: TSRHSJacobian func;
7229: void* ctx;
7232: TSGetRHSJacobian(ts,&J,&B,&func,&ctx);
7233: (*func)(ts,0.0,ts->vec_sol,J,B,ctx);
7234: MatShellTestMult(J,RHSWrapperFunction_TSRHSJacobianTest,ts->vec_sol,ts,flg);
7235: return(0);
7236: }
7238: /*@C
7239: TSRHSJacobianTestTranspose - Compares the multiply transpose routine provided to the MATSHELL with differencing on the TS given RHS function.
7241: Logically Collective on TS
7243: Input Parameters:
7244: TS - the time stepping routine
7246: Output Parameter:
7247: . flg - PETSC_TRUE if the multiply is likely correct
7249: Options Database:
7250: . -ts_rhs_jacobian_test_mult_transpose -mat_shell_test_mult_transpose_view - run the test at each timestep of the integrator
7252: Notes:
7253: This only works for problems defined only the RHS function and Jacobian NOT IFunction and IJacobian
7255: Level: advanced
7257: .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMultTranspose(), TSRHSJacobianTest()
7258: @*/
7259: PetscErrorCode TSRHSJacobianTestTranspose(TS ts,PetscBool *flg)
7260: {
7261: Mat J,B;
7263: void *ctx;
7264: TSRHSJacobian func;
7267: TSGetRHSJacobian(ts,&J,&B,&func,&ctx);
7268: (*func)(ts,0.0,ts->vec_sol,J,B,ctx);
7269: MatShellTestMultTranspose(J,RHSWrapperFunction_TSRHSJacobianTest,ts->vec_sol,ts,flg);
7270: return(0);
7271: }
7273: /*@
7274: TSSetUseSplitRHSFunction - Use the split RHSFunction when a multirate method is used.
7276: Logically collective
7278: Input Parameter:
7279: + ts - timestepping context
7280: - use_splitrhsfunction - PETSC_TRUE indicates that the split RHSFunction will be used
7282: Options Database:
7283: . -ts_use_splitrhsfunction - <true,false>
7285: Notes:
7286: This is only useful for multirate methods
7288: Level: intermediate
7290: .seealso: TSGetUseSplitRHSFunction()
7291: @*/
7292: PetscErrorCode TSSetUseSplitRHSFunction(TS ts, PetscBool use_splitrhsfunction)
7293: {
7296: ts->use_splitrhsfunction = use_splitrhsfunction;
7297: return(0);
7298: }
7300: /*@
7301: TSGetUseSplitRHSFunction - Gets whether to use the split RHSFunction when a multirate method is used.
7303: Not collective
7305: Input Parameter:
7306: . ts - timestepping context
7308: Output Parameter:
7309: . use_splitrhsfunction - PETSC_TRUE indicates that the split RHSFunction will be used
7311: Level: intermediate
7313: .seealso: TSSetUseSplitRHSFunction()
7314: @*/
7315: PetscErrorCode TSGetUseSplitRHSFunction(TS ts, PetscBool *use_splitrhsfunction)
7316: {
7319: *use_splitrhsfunction = ts->use_splitrhsfunction;
7320: return(0);
7321: }