Actual source code: ts.c

petsc-3.9.4 2018-09-11
Report Typos and Errors
  1:  #include <petsc/private/tsimpl.h>
  2:  #include <petscdmshell.h>
  3:  #include <petscdmda.h>
  4:  #include <petscviewer.h>
  5:  #include <petscdraw.h>

  7: /* Logging support */
  8: PetscClassId  TS_CLASSID, DMTS_CLASSID;
  9: PetscLogEvent TS_Step, TS_PseudoComputeTimeStep, TS_FunctionEval, TS_JacobianEval;

 11: const char *const TSExactFinalTimeOptions[] = {"UNSPECIFIED","STEPOVER","INTERPOLATE","MATCHSTEP","TSExactFinalTimeOption","TS_EXACTFINALTIME_",0};

 13: /*@C
 14:    TSMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user

 16:    Collective on TS

 18:    Input Parameters:
 19: +  ts - TS object you wish to monitor
 20: .  name - the monitor type one is seeking
 21: .  help - message indicating what monitoring is done
 22: .  manual - manual page for the monitor
 23: .  monitor - the monitor function
 24: -  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

 26:    Level: developer

 28: .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
 29:           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
 30:           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
 31:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
 32:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
 33:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
 34:           PetscOptionsFList(), PetscOptionsEList()
 35: @*/
 36: PetscErrorCode  TSMonitorSetFromOptions(TS ts,const char name[],const char help[], const char manual[],PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*),PetscErrorCode (*monitorsetup)(TS,PetscViewerAndFormat*))
 37: {
 38:   PetscErrorCode    ierr;
 39:   PetscViewer       viewer;
 40:   PetscViewerFormat format;
 41:   PetscBool         flg;

 44:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)ts),((PetscObject)ts)->prefix,name,&viewer,&format,&flg);
 45:   if (flg) {
 46:     PetscViewerAndFormat *vf;
 47:     PetscViewerAndFormatCreate(viewer,format,&vf);
 48:     PetscObjectDereference((PetscObject)viewer);
 49:     if (monitorsetup) {
 50:       (*monitorsetup)(ts,vf);
 51:     }
 52:     TSMonitorSet(ts,(PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);
 53:   }
 54:   return(0);
 55: }

 57: static PetscErrorCode TSAdaptSetDefaultType(TSAdapt adapt,TSAdaptType default_type)
 58: {

 64:   if (!((PetscObject)adapt)->type_name) {
 65:     TSAdaptSetType(adapt,default_type);
 66:   }
 67:   return(0);
 68: }

 70: /*@
 71:    TSSetFromOptions - Sets various TS parameters from user options.

 73:    Collective on TS

 75:    Input Parameter:
 76: .  ts - the TS context obtained from TSCreate()

 78:    Options Database Keys:
 79: +  -ts_type <type> - TSEULER, TSBEULER, TSSUNDIALS, TSPSEUDO, TSCN, TSRK, TSTHETA, TSALPHA, TSGLLE, TSSSP, TSGLEE
 80: .  -ts_save_trajectory - checkpoint the solution at each time-step
 81: .  -ts_max_time <time> - maximum time to compute to
 82: .  -ts_max_steps <steps> - maximum number of time-steps to take
 83: .  -ts_init_time <time> - initial time to start computation
 84: .  -ts_final_time <time> - final time to compute to
 85: .  -ts_dt <dt> - initial time step
 86: .  -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
 87: .  -ts_max_snes_failures <maxfailures> - Maximum number of nonlinear solve failures allowed
 88: .  -ts_max_reject <maxrejects> - Maximum number of step rejections before step fails
 89: .  -ts_error_if_step_fails <true,false> - Error if no step succeeds
 90: .  -ts_rtol <rtol> - relative tolerance for local truncation error
 91: .  -ts_atol <atol> Absolute tolerance for local truncation error
 92: .  -ts_rhs_jacobian_test_mult -mat_shell_test_mult_view - test the Jacobian at each iteration against finite difference with RHS function
 93: .  -ts_rhs_jacobian_test_mult_transpose -mat_shell_test_mult_transpose_view - test the Jacobian at each iteration against finite difference with RHS function
 94: .  -ts_adjoint_solve <yes,no> After solving the ODE/DAE solve the adjoint problem (requires -ts_save_trajectory)
 95: .  -ts_fd_color - Use finite differences with coloring to compute IJacobian
 96: .  -ts_monitor - print information at each timestep
 97: .  -ts_monitor_lg_solution - Monitor solution graphically
 98: .  -ts_monitor_lg_error - Monitor error graphically
 99: .  -ts_monitor_error - Monitors norm of error
100: .  -ts_monitor_lg_timestep - Monitor timestep size graphically
101: .  -ts_monitor_lg_timestep_log - Monitor log timestep size graphically
102: .  -ts_monitor_lg_snes_iterations - Monitor number nonlinear iterations for each timestep graphically
103: .  -ts_monitor_lg_ksp_iterations - Monitor number nonlinear iterations for each timestep graphically
104: .  -ts_monitor_sp_eig - Monitor eigenvalues of linearized operator graphically
105: .  -ts_monitor_draw_solution - Monitor solution graphically
106: .  -ts_monitor_draw_solution_phase  <xleft,yleft,xright,yright> - Monitor solution graphically with phase diagram, requires problem with exactly 2 degrees of freedom
107: .  -ts_monitor_draw_error - Monitor error graphically, requires use to have provided TSSetSolutionFunction()
108: .  -ts_monitor_solution [ascii binary draw][:filename][:viewerformat] - monitors the solution at each timestep
109: .  -ts_monitor_solution_vtk <filename.vts,filename.vtu> - Save each time step to a binary file, use filename-%%03D.vts (filename-%%03D.vtu)
110: .  -ts_monitor_envelope - determine maximum and minimum value of each component of the solution over the solution time

112:    Developer Note: We should unify all the -ts_monitor options in the way that -xxx_view has been unified

114:    Level: beginner

116: .keywords: TS, timestep, set, options, database

118: .seealso: TSGetType()
119: @*/
120: PetscErrorCode  TSSetFromOptions(TS ts)
121: {
122:   PetscBool              opt,flg,tflg;
123:   PetscErrorCode         ierr;
124:   char                   monfilename[PETSC_MAX_PATH_LEN];
125:   PetscReal              time_step;
126:   TSExactFinalTimeOption eftopt;
127:   char                   dir[16];
128:   TSIFunction            ifun;
129:   const char             *defaultType;
130:   char                   typeName[256];


135:   TSRegisterAll();
136:   TSGetIFunction(ts,NULL,&ifun,NULL);

138:   PetscObjectOptionsBegin((PetscObject)ts);
139:   if (((PetscObject)ts)->type_name)
140:     defaultType = ((PetscObject)ts)->type_name;
141:   else
142:     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:   PetscOptionsReal("-ts_max_time","Maximum time to run to","TSSetMaxTime",ts->max_time,&ts->max_time,NULL);
152:   PetscOptionsInt("-ts_max_steps","Maximum number of time steps","TSSetMaxSteps",ts->max_steps,&ts->max_steps,NULL);
153:   PetscOptionsReal("-ts_init_time","Initial time","TSSetTime",ts->ptime,&ts->ptime,NULL);
154:   PetscOptionsReal("-ts_final_time","Final time to run to","TSSetMaxTime",ts->max_time,&ts->max_time,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: #if defined(PETSC_HAVE_SAWS)
168:   {
169:   PetscBool set;
170:   flg  = PETSC_FALSE;
171:   PetscOptionsBool("-ts_saws_block","Block for SAWs memory snooper at end of TSSolve","PetscObjectSAWsBlock",((PetscObject)ts)->amspublishblock,&flg,&set);
172:   if (set) {
173:     PetscObjectSAWsSetBlock((PetscObject)ts,flg);
174:   }
175:   }
176: #endif

178:   /* Monitor options */
179:   TSMonitorSetFromOptions(ts,"-ts_monitor","Monitor time and timestep size","TSMonitorDefault",TSMonitorDefault,NULL);
180:   TSMonitorSetFromOptions(ts,"-ts_monitor_solution","View the solution at each timestep","TSMonitorSolution",TSMonitorSolution,NULL);

182:   PetscOptionsString("-ts_monitor_python","Use Python function","TSMonitorSet",0,monfilename,PETSC_MAX_PATH_LEN,&flg);
183:   if (flg) {PetscPythonMonitorSet((PetscObject)ts,monfilename);}

185:   PetscOptionsName("-ts_monitor_lg_solution","Monitor solution graphically","TSMonitorLGSolution",&opt);
186:   if (opt) {
187:     TSMonitorLGCtx ctx;
188:     PetscInt       howoften = 1;

190:     PetscOptionsInt("-ts_monitor_lg_solution","Monitor solution graphically","TSMonitorLGSolution",howoften,&howoften,NULL);
191:     TSMonitorLGCtxCreate(PETSC_COMM_SELF,0,0,PETSC_DECIDE,PETSC_DECIDE,400,300,howoften,&ctx);
192:     TSMonitorSet(ts,TSMonitorLGSolution,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);
193:   }

195:   PetscOptionsName("-ts_monitor_lg_error","Monitor error graphically","TSMonitorLGError",&opt);
196:   if (opt) {
197:     TSMonitorLGCtx ctx;
198:     PetscInt       howoften = 1;

200:     PetscOptionsInt("-ts_monitor_lg_error","Monitor error graphically","TSMonitorLGError",howoften,&howoften,NULL);
201:     TSMonitorLGCtxCreate(PETSC_COMM_SELF,0,0,PETSC_DECIDE,PETSC_DECIDE,400,300,howoften,&ctx);
202:     TSMonitorSet(ts,TSMonitorLGError,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);
203:   }
204:   TSMonitorSetFromOptions(ts,"-ts_monitor_error","View the error at each timestep","TSMonitorError",TSMonitorError,NULL);

206:   PetscOptionsName("-ts_monitor_lg_timestep","Monitor timestep size graphically","TSMonitorLGTimeStep",&opt);
207:   if (opt) {
208:     TSMonitorLGCtx ctx;
209:     PetscInt       howoften = 1;

211:     PetscOptionsInt("-ts_monitor_lg_timestep","Monitor timestep size graphically","TSMonitorLGTimeStep",howoften,&howoften,NULL);
212:     TSMonitorLGCtxCreate(PetscObjectComm((PetscObject)ts),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,howoften,&ctx);
213:     TSMonitorSet(ts,TSMonitorLGTimeStep,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);
214:   }
215:   PetscOptionsName("-ts_monitor_lg_timestep_log","Monitor log timestep size graphically","TSMonitorLGTimeStep",&opt);
216:   if (opt) {
217:     TSMonitorLGCtx ctx;
218:     PetscInt       howoften = 1;

220:     PetscOptionsInt("-ts_monitor_lg_timestep_log","Monitor log timestep size graphically","TSMonitorLGTimeStep",howoften,&howoften,NULL);
221:     TSMonitorLGCtxCreate(PetscObjectComm((PetscObject)ts),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,howoften,&ctx);
222:     TSMonitorSet(ts,TSMonitorLGTimeStep,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);
223:     ctx->semilogy = PETSC_TRUE;
224:   }

226:   PetscOptionsName("-ts_monitor_lg_snes_iterations","Monitor number nonlinear iterations for each timestep graphically","TSMonitorLGSNESIterations",&opt);
227:   if (opt) {
228:     TSMonitorLGCtx ctx;
229:     PetscInt       howoften = 1;

231:     PetscOptionsInt("-ts_monitor_lg_snes_iterations","Monitor number nonlinear iterations for each timestep graphically","TSMonitorLGSNESIterations",howoften,&howoften,NULL);
232:     TSMonitorLGCtxCreate(PetscObjectComm((PetscObject)ts),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,howoften,&ctx);
233:     TSMonitorSet(ts,TSMonitorLGSNESIterations,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);
234:   }
235:   PetscOptionsName("-ts_monitor_lg_ksp_iterations","Monitor number nonlinear iterations for each timestep graphically","TSMonitorLGKSPIterations",&opt);
236:   if (opt) {
237:     TSMonitorLGCtx ctx;
238:     PetscInt       howoften = 1;

240:     PetscOptionsInt("-ts_monitor_lg_ksp_iterations","Monitor number nonlinear iterations for each timestep graphically","TSMonitorLGKSPIterations",howoften,&howoften,NULL);
241:     TSMonitorLGCtxCreate(PetscObjectComm((PetscObject)ts),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,howoften,&ctx);
242:     TSMonitorSet(ts,TSMonitorLGKSPIterations,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);
243:   }
244:   PetscOptionsName("-ts_monitor_sp_eig","Monitor eigenvalues of linearized operator graphically","TSMonitorSPEig",&opt);
245:   if (opt) {
246:     TSMonitorSPEigCtx ctx;
247:     PetscInt          howoften = 1;

249:     PetscOptionsInt("-ts_monitor_sp_eig","Monitor eigenvalues of linearized operator graphically","TSMonitorSPEig",howoften,&howoften,NULL);
250:     TSMonitorSPEigCtxCreate(PETSC_COMM_SELF,0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx);
251:     TSMonitorSet(ts,TSMonitorSPEig,ctx,(PetscErrorCode (*)(void**))TSMonitorSPEigCtxDestroy);
252:   }
253:   opt  = PETSC_FALSE;
254:   PetscOptionsName("-ts_monitor_draw_solution","Monitor solution graphically","TSMonitorDrawSolution",&opt);
255:   if (opt) {
256:     TSMonitorDrawCtx ctx;
257:     PetscInt         howoften = 1;

259:     PetscOptionsInt("-ts_monitor_draw_solution","Monitor solution graphically","TSMonitorDrawSolution",howoften,&howoften,NULL);
260:     TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),0,"Computed Solution",PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx);
261:     TSMonitorSet(ts,TSMonitorDrawSolution,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);
262:   }
263:   opt  = PETSC_FALSE;
264:   PetscOptionsName("-ts_monitor_draw_solution_phase","Monitor solution graphically","TSMonitorDrawSolutionPhase",&opt);
265:   if (opt) {
266:     TSMonitorDrawCtx ctx;
267:     PetscReal        bounds[4];
268:     PetscInt         n = 4;
269:     PetscDraw        draw;
270:     PetscDrawAxis    axis;

272:     PetscOptionsRealArray("-ts_monitor_draw_solution_phase","Monitor solution graphically","TSMonitorDrawSolutionPhase",bounds,&n,NULL);
273:     if (n != 4) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONG,"Must provide bounding box of phase field");
274:     TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,1,&ctx);
275:     PetscViewerDrawGetDraw(ctx->viewer,0,&draw);
276:     PetscViewerDrawGetDrawAxis(ctx->viewer,0,&axis);
277:     PetscDrawAxisSetLimits(axis,bounds[0],bounds[2],bounds[1],bounds[3]);
278:     PetscDrawAxisSetLabels(axis,"Phase Diagram","Variable 1","Variable 2");
279:     TSMonitorSet(ts,TSMonitorDrawSolutionPhase,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);
280:   }
281:   opt  = PETSC_FALSE;
282:   PetscOptionsName("-ts_monitor_draw_error","Monitor error graphically","TSMonitorDrawError",&opt);
283:   if (opt) {
284:     TSMonitorDrawCtx ctx;
285:     PetscInt         howoften = 1;

287:     PetscOptionsInt("-ts_monitor_draw_error","Monitor error graphically","TSMonitorDrawError",howoften,&howoften,NULL);
288:     TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),0,"Error",PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx);
289:     TSMonitorSet(ts,TSMonitorDrawError,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);
290:   }
291:   opt  = PETSC_FALSE;
292:   PetscOptionsName("-ts_monitor_draw_solution_function","Monitor solution provided by TSMonitorSetSolutionFunction() graphically","TSMonitorDrawSolutionFunction",&opt);
293:   if (opt) {
294:     TSMonitorDrawCtx ctx;
295:     PetscInt         howoften = 1;

297:     PetscOptionsInt("-ts_monitor_draw_solution_function","Monitor solution provided by TSMonitorSetSolutionFunction() graphically","TSMonitorDrawSolutionFunction",howoften,&howoften,NULL);
298:     TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),0,"Solution provided by user function",PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx);
299:     TSMonitorSet(ts,TSMonitorDrawSolutionFunction,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);
300:   }

302:   opt  = PETSC_FALSE;
303:   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);
304:   if (flg) {
305:     const char *ptr,*ptr2;
306:     char       *filetemplate;
307:     if (!monfilename[0]) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"-ts_monitor_solution_vtk requires a file template, e.g. filename-%%03D.vts");
308:     /* Do some cursory validation of the input. */
309:     PetscStrstr(monfilename,"%",(char**)&ptr);
310:     if (!ptr) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"-ts_monitor_solution_vtk requires a file template, e.g. filename-%%03D.vts");
311:     for (ptr++; ptr && *ptr; ptr++) {
312:       PetscStrchr("DdiouxX",*ptr,(char**)&ptr2);
313:       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");
314:       if (ptr2) break;
315:     }
316:     PetscStrallocpy(monfilename,&filetemplate);
317:     TSMonitorSet(ts,TSMonitorSolutionVTK,filetemplate,(PetscErrorCode (*)(void**))TSMonitorSolutionVTKDestroy);
318:   }

320:   PetscOptionsString("-ts_monitor_dmda_ray","Display a ray of the solution","None","y=0",dir,16,&flg);
321:   if (flg) {
322:     TSMonitorDMDARayCtx *rayctx;
323:     int                  ray = 0;
324:     DMDADirection        ddir;
325:     DM                   da;
326:     PetscMPIInt          rank;

328:     if (dir[1] != '=') SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONG,"Unknown ray %s",dir);
329:     if (dir[0] == 'x') ddir = DMDA_X;
330:     else if (dir[0] == 'y') ddir = DMDA_Y;
331:     else SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONG,"Unknown ray %s",dir);
332:     sscanf(dir+2,"%d",&ray);

334:     PetscInfo2(((PetscObject)ts),"Displaying DMDA ray %c = %D\n",dir[0],ray);
335:     PetscNew(&rayctx);
336:     TSGetDM(ts,&da);
337:     DMDAGetRay(da,ddir,ray,&rayctx->ray,&rayctx->scatter);
338:     MPI_Comm_rank(PetscObjectComm((PetscObject)ts),&rank);
339:     if (!rank) {
340:       PetscViewerDrawOpen(PETSC_COMM_SELF,0,0,0,0,600,300,&rayctx->viewer);
341:     }
342:     rayctx->lgctx = NULL;
343:     TSMonitorSet(ts,TSMonitorDMDARay,rayctx,TSMonitorDMDARayDestroy);
344:   }
345:   PetscOptionsString("-ts_monitor_lg_dmda_ray","Display a ray of the solution","None","x=0",dir,16,&flg);
346:   if (flg) {
347:     TSMonitorDMDARayCtx *rayctx;
348:     int                 ray = 0;
349:     DMDADirection       ddir;
350:     DM                  da;
351:     PetscInt            howoften = 1;

353:     if (dir[1] != '=') SETERRQ1(PetscObjectComm((PetscObject) ts), PETSC_ERR_ARG_WRONG, "Malformed ray %s", dir);
354:     if      (dir[0] == 'x') ddir = DMDA_X;
355:     else if (dir[0] == 'y') ddir = DMDA_Y;
356:     else SETERRQ1(PetscObjectComm((PetscObject) ts), PETSC_ERR_ARG_WRONG, "Unknown ray direction %s", dir);
357:     sscanf(dir+2, "%d", &ray);

359:     PetscInfo2(((PetscObject) ts),"Displaying LG DMDA ray %c = %D\n", dir[0], ray);
360:     PetscNew(&rayctx);
361:     TSGetDM(ts, &da);
362:     DMDAGetRay(da, ddir, ray, &rayctx->ray, &rayctx->scatter);
363:     TSMonitorLGCtxCreate(PETSC_COMM_SELF,0,0,PETSC_DECIDE,PETSC_DECIDE,600,400,howoften,&rayctx->lgctx);
364:     TSMonitorSet(ts, TSMonitorLGDMDARay, rayctx, TSMonitorDMDARayDestroy);
365:   }

367:   PetscOptionsName("-ts_monitor_envelope","Monitor maximum and minimum value of each component of the solution","TSMonitorEnvelope",&opt);
368:   if (opt) {
369:     TSMonitorEnvelopeCtx ctx;

371:     TSMonitorEnvelopeCtxCreate(ts,&ctx);
372:     TSMonitorSet(ts,TSMonitorEnvelope,ctx,(PetscErrorCode (*)(void**))TSMonitorEnvelopeCtxDestroy);
373:   }

375:   flg  = PETSC_FALSE;
376:   PetscOptionsBool("-ts_fd_color", "Use finite differences with coloring to compute IJacobian", "TSComputeJacobianDefaultColor", flg, &flg, NULL);
377:   if (flg) {
378:     DM   dm;
379:     DMTS tdm;

381:     TSGetDM(ts, &dm);
382:     DMGetDMTS(dm, &tdm);
383:     tdm->ijacobianctx = NULL;
384:     TSSetIJacobian(ts, NULL, NULL, TSComputeIJacobianDefaultColor, 0);
385:     PetscInfo(ts, "Setting default finite difference coloring Jacobian matrix\n");
386:   }

388:   /* Handle specific TS options */
389:   if (ts->ops->setfromoptions) {
390:     (*ts->ops->setfromoptions)(PetscOptionsObject,ts);
391:   }

393:   /* Handle TSAdapt options */
394:   TSGetAdapt(ts,&ts->adapt);
395:   TSAdaptSetDefaultType(ts->adapt,ts->default_adapt_type);
396:   TSAdaptSetFromOptions(PetscOptionsObject,ts->adapt);

398:   /* TS trajectory must be set after TS, since it may use some TS options above */
399:   tflg = ts->trajectory ? PETSC_TRUE : PETSC_FALSE;
400:   PetscOptionsBool("-ts_save_trajectory","Save the solution at each timestep","TSSetSaveTrajectory",tflg,&tflg,NULL);
401:   if (tflg) {
402:     TSSetSaveTrajectory(ts);
403:   }

405:   TSAdjointSetFromOptions(PetscOptionsObject,ts);

407:   /* process any options handlers added with PetscObjectAddOptionsHandler() */
408:   PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)ts);
409:   PetscOptionsEnd();

411:   if (ts->trajectory) {
412:     TSTrajectorySetFromOptions(ts->trajectory,ts);
413:   }

415:   TSGetSNES(ts,&ts->snes);
416:   if (ts->problem_type == TS_LINEAR) {SNESSetType(ts->snes,SNESKSPONLY);}
417:   SNESSetFromOptions(ts->snes);
418:   return(0);
419: }

421: /*@
422:    TSGetTrajectory - Gets the trajectory from a TS if it exists

424:    Collective on TS

426:    Input Parameters:
427: .  ts - the TS context obtained from TSCreate()

429:    Output Parameters;
430: .  tr - the TSTrajectory object, if it exists

432:    Note: This routine should be called after all TS options have been set

434:    Level: advanced

436: .seealso: TSGetTrajectory(), TSAdjointSolve(), TSTrajectory, TSTrajectoryCreate()

438: .keywords: TS, set, checkpoint,
439: @*/
440: PetscErrorCode  TSGetTrajectory(TS ts,TSTrajectory *tr)
441: {
444:   *tr = ts->trajectory;
445:   return(0);
446: }

448: /*@
449:    TSSetSaveTrajectory - Causes the TS to save its solutions as it iterates forward in time in a TSTrajectory object

451:    Collective on TS

453:    Input Parameters:
454: .  ts - the TS context obtained from TSCreate()

456:    Options Database:
457: +  -ts_save_trajectory - saves the trajectory to a file
458: -  -ts_trajectory_type type

460: Note: This routine should be called after all TS options have been set

462:     The TSTRAJECTORYVISUALIZATION files can be loaded into Python with $PETSC_DIR/lib/petsc/bin/PetscBinaryIOTrajectory.py and
463:    MATLAB with $PETSC_DIR/share/petsc/matlab/PetscReadBinaryTrajectory.m

465:    Level: intermediate

467: .seealso: TSGetTrajectory(), TSAdjointSolve(), TSTrajectoryType, TSSetTrajectoryType()

469: .keywords: TS, set, checkpoint,
470: @*/
471: PetscErrorCode  TSSetSaveTrajectory(TS ts)
472: {

477:   if (!ts->trajectory) {
478:     TSTrajectoryCreate(PetscObjectComm((PetscObject)ts),&ts->trajectory);
479:   }
480:   return(0);
481: }

483: /*@
484:    TSComputeRHSJacobian - Computes the Jacobian matrix that has been
485:       set with TSSetRHSJacobian().

487:    Collective on TS and Vec

489:    Input Parameters:
490: +  ts - the TS context
491: .  t - current timestep
492: -  U - input vector

494:    Output Parameters:
495: +  A - Jacobian matrix
496: .  B - optional preconditioning matrix
497: -  flag - flag indicating matrix structure

499:    Notes:
500:    Most users should not need to explicitly call this routine, as it
501:    is used internally within the nonlinear solvers.

503:    See KSPSetOperators() for important information about setting the
504:    flag parameter.

506:    Level: developer

508: .keywords: SNES, compute, Jacobian, matrix

510: .seealso:  TSSetRHSJacobian(), KSPSetOperators()
511: @*/
512: PetscErrorCode  TSComputeRHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B)
513: {
514:   PetscErrorCode   ierr;
515:   PetscObjectState Ustate;
516:   PetscObjectId    Uid;
517:   DM               dm;
518:   DMTS             tsdm;
519:   TSRHSJacobian    rhsjacobianfunc;
520:   void             *ctx;
521:   TSIJacobian      ijacobianfunc;
522:   TSRHSFunction    rhsfunction;

528:   TSGetDM(ts,&dm);
529:   DMGetDMTS(dm,&tsdm);
530:   DMTSGetRHSJacobian(dm,&rhsjacobianfunc,&ctx);
531:   DMTSGetIJacobian(dm,&ijacobianfunc,NULL);
532:   DMTSGetRHSFunction(dm,&rhsfunction,&ctx);
533:   PetscObjectStateGet((PetscObject)U,&Ustate);
534:   PetscObjectGetId((PetscObject)U,&Uid);
535:   if (ts->rhsjacobian.time == t && (ts->problem_type == TS_LINEAR || (ts->rhsjacobian.Xid == Uid && ts->rhsjacobian.Xstate == Ustate)) && (rhsfunction != TSComputeRHSFunctionLinear)) {
536:     return(0);
537:   }

539:   if (!rhsjacobianfunc && !ijacobianfunc) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must call TSSetRHSJacobian() and / or TSSetIJacobian()");

541:   if (ts->rhsjacobian.reuse) {
542:     MatShift(A,-ts->rhsjacobian.shift);
543:     MatScale(A,1./ts->rhsjacobian.scale);
544:     if (B && A != B) {
545:       MatShift(B,-ts->rhsjacobian.shift);
546:       MatScale(B,1./ts->rhsjacobian.scale);
547:     }
548:     ts->rhsjacobian.shift = 0;
549:     ts->rhsjacobian.scale = 1.;
550:   }

552:   if (rhsjacobianfunc) {
553:     PetscBool missing;
554:     PetscLogEventBegin(TS_JacobianEval,ts,U,A,B);
555:     PetscStackPush("TS user Jacobian function");
556:     (*rhsjacobianfunc)(ts,t,U,A,B,ctx);
557:     PetscStackPop;
558:     PetscLogEventEnd(TS_JacobianEval,ts,U,A,B);
559:     if (A) {
560:       MatMissingDiagonal(A,&missing,NULL);
561:       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");
562:     }
563:     if (B && B != A) {
564:       MatMissingDiagonal(B,&missing,NULL);
565:       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");
566:     }
567:   } else {
568:     MatZeroEntries(A);
569:     if (A != B) {MatZeroEntries(B);}
570:   }
571:   ts->rhsjacobian.time       = t;
572:   PetscObjectGetId((PetscObject)U,&ts->rhsjacobian.Xid);
573:   PetscObjectStateGet((PetscObject)U,&ts->rhsjacobian.Xstate);
574:   return(0);
575: }

577: /*@
578:    TSComputeRHSFunction - Evaluates the right-hand-side function.

580:    Collective on TS and Vec

582:    Input Parameters:
583: +  ts - the TS context
584: .  t - current time
585: -  U - state vector

587:    Output Parameter:
588: .  y - right hand side

590:    Note:
591:    Most users should not need to explicitly call this routine, as it
592:    is used internally within the nonlinear solvers.

594:    Level: developer

596: .keywords: TS, compute

598: .seealso: TSSetRHSFunction(), TSComputeIFunction()
599: @*/
600: PetscErrorCode TSComputeRHSFunction(TS ts,PetscReal t,Vec U,Vec y)
601: {
603:   TSRHSFunction  rhsfunction;
604:   TSIFunction    ifunction;
605:   void           *ctx;
606:   DM             dm;

612:   TSGetDM(ts,&dm);
613:   DMTSGetRHSFunction(dm,&rhsfunction,&ctx);
614:   DMTSGetIFunction(dm,&ifunction,NULL);

616:   if (!rhsfunction && !ifunction) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must call TSSetRHSFunction() and / or TSSetIFunction()");

618:   PetscLogEventBegin(TS_FunctionEval,ts,U,y,0);
619:   if (rhsfunction) {
620:     PetscStackPush("TS user right-hand-side function");
621:     (*rhsfunction)(ts,t,U,y,ctx);
622:     PetscStackPop;
623:   } else {
624:     VecZeroEntries(y);
625:   }

627:   PetscLogEventEnd(TS_FunctionEval,ts,U,y,0);
628:   return(0);
629: }

631: /*@
632:    TSComputeSolutionFunction - Evaluates the solution function.

634:    Collective on TS and Vec

636:    Input Parameters:
637: +  ts - the TS context
638: -  t - current time

640:    Output Parameter:
641: .  U - the solution

643:    Note:
644:    Most users should not need to explicitly call this routine, as it
645:    is used internally within the nonlinear solvers.

647:    Level: developer

649: .keywords: TS, compute

651: .seealso: TSSetSolutionFunction(), TSSetRHSFunction(), TSComputeIFunction()
652: @*/
653: PetscErrorCode TSComputeSolutionFunction(TS ts,PetscReal t,Vec U)
654: {
655:   PetscErrorCode     ierr;
656:   TSSolutionFunction solutionfunction;
657:   void               *ctx;
658:   DM                 dm;

663:   TSGetDM(ts,&dm);
664:   DMTSGetSolutionFunction(dm,&solutionfunction,&ctx);

666:   if (solutionfunction) {
667:     PetscStackPush("TS user solution function");
668:     (*solutionfunction)(ts,t,U,ctx);
669:     PetscStackPop;
670:   }
671:   return(0);
672: }
673: /*@
674:    TSComputeForcingFunction - Evaluates the forcing function.

676:    Collective on TS and Vec

678:    Input Parameters:
679: +  ts - the TS context
680: -  t - current time

682:    Output Parameter:
683: .  U - the function value

685:    Note:
686:    Most users should not need to explicitly call this routine, as it
687:    is used internally within the nonlinear solvers.

689:    Level: developer

691: .keywords: TS, compute

693: .seealso: TSSetSolutionFunction(), TSSetRHSFunction(), TSComputeIFunction()
694: @*/
695: PetscErrorCode TSComputeForcingFunction(TS ts,PetscReal t,Vec U)
696: {
697:   PetscErrorCode     ierr, (*forcing)(TS,PetscReal,Vec,void*);
698:   void               *ctx;
699:   DM                 dm;

704:   TSGetDM(ts,&dm);
705:   DMTSGetForcingFunction(dm,&forcing,&ctx);

707:   if (forcing) {
708:     PetscStackPush("TS user forcing function");
709:     (*forcing)(ts,t,U,ctx);
710:     PetscStackPop;
711:   }
712:   return(0);
713: }

715: static PetscErrorCode TSGetRHSVec_Private(TS ts,Vec *Frhs)
716: {
717:   Vec            F;

721:   *Frhs = NULL;
722:   TSGetIFunction(ts,&F,NULL,NULL);
723:   if (!ts->Frhs) {
724:     VecDuplicate(F,&ts->Frhs);
725:   }
726:   *Frhs = ts->Frhs;
727:   return(0);
728: }

730: static PetscErrorCode TSGetRHSMats_Private(TS ts,Mat *Arhs,Mat *Brhs)
731: {
732:   Mat            A,B;
734:   TSIJacobian    ijacobian;

737:   if (Arhs) *Arhs = NULL;
738:   if (Brhs) *Brhs = NULL;
739:   TSGetIJacobian(ts,&A,&B,&ijacobian,NULL);
740:   if (Arhs) {
741:     if (!ts->Arhs) {
742:       if (ijacobian) {
743:         MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&ts->Arhs);
744:       } else {
745:         ts->Arhs = A;
746:         PetscObjectReference((PetscObject)A);
747:       }
748:     } else {
749:       PetscBool flg;
750:       SNESGetUseMatrixFree(ts->snes,NULL,&flg);
751:       /* Handle case where user provided only RHSJacobian and used -snes_mf_operator */
752:       if (flg && !ijacobian && ts->Arhs == ts->Brhs){
753:         PetscObjectDereference((PetscObject)ts->Arhs);
754:         ts->Arhs = A;
755:         PetscObjectReference((PetscObject)A);
756:       }
757:     }
758:     *Arhs = ts->Arhs;
759:   }
760:   if (Brhs) {
761:     if (!ts->Brhs) {
762:       if (A != B) {
763:         if (ijacobian) {
764:           MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&ts->Brhs);
765:         } else {
766:           ts->Brhs = B;
767:           PetscObjectReference((PetscObject)B);
768:         }
769:       } else {
770:         PetscObjectReference((PetscObject)ts->Arhs);
771:         ts->Brhs = ts->Arhs;
772:       }
773:     }
774:     *Brhs = ts->Brhs;
775:   }
776:   return(0);
777: }

779: /*@
780:    TSComputeIFunction - Evaluates the DAE residual written in implicit form F(t,U,Udot)=0

782:    Collective on TS and Vec

784:    Input Parameters:
785: +  ts - the TS context
786: .  t - current time
787: .  U - state vector
788: .  Udot - time derivative of state vector
789: -  imex - flag indicates if the method is IMEX so that the RHSFunction should be kept separate

791:    Output Parameter:
792: .  Y - right hand side

794:    Note:
795:    Most users should not need to explicitly call this routine, as it
796:    is used internally within the nonlinear solvers.

798:    If the user did did not write their equations in implicit form, this
799:    function recasts them in implicit form.

801:    Level: developer

803: .keywords: TS, compute

805: .seealso: TSSetIFunction(), TSComputeRHSFunction()
806: @*/
807: PetscErrorCode TSComputeIFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec Y,PetscBool imex)
808: {
810:   TSIFunction    ifunction;
811:   TSRHSFunction  rhsfunction;
812:   void           *ctx;
813:   DM             dm;


821:   TSGetDM(ts,&dm);
822:   DMTSGetIFunction(dm,&ifunction,&ctx);
823:   DMTSGetRHSFunction(dm,&rhsfunction,NULL);

825:   if (!rhsfunction && !ifunction) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must call TSSetRHSFunction() and / or TSSetIFunction()");

827:   PetscLogEventBegin(TS_FunctionEval,ts,U,Udot,Y);
828:   if (ifunction) {
829:     PetscStackPush("TS user implicit function");
830:     (*ifunction)(ts,t,U,Udot,Y,ctx);
831:     PetscStackPop;
832:   }
833:   if (imex) {
834:     if (!ifunction) {
835:       VecCopy(Udot,Y);
836:     }
837:   } else if (rhsfunction) {
838:     if (ifunction) {
839:       Vec Frhs;
840:       TSGetRHSVec_Private(ts,&Frhs);
841:       TSComputeRHSFunction(ts,t,U,Frhs);
842:       VecAXPY(Y,-1,Frhs);
843:     } else {
844:       TSComputeRHSFunction(ts,t,U,Y);
845:       VecAYPX(Y,-1,Udot);
846:     }
847:   }
848:   PetscLogEventEnd(TS_FunctionEval,ts,U,Udot,Y);
849:   return(0);
850: }

852: /*@
853:    TSComputeIJacobian - Evaluates the Jacobian of the DAE

855:    Collective on TS and Vec

857:    Input
858:       Input Parameters:
859: +  ts - the TS context
860: .  t - current timestep
861: .  U - state vector
862: .  Udot - time derivative of state vector
863: .  shift - shift to apply, see note below
864: -  imex - flag indicates if the method is IMEX so that the RHSJacobian should be kept separate

866:    Output Parameters:
867: +  A - Jacobian matrix
868: -  B - matrix from which the preconditioner is constructed; often the same as A

870:    Notes:
871:    If F(t,U,Udot)=0 is the DAE, the required Jacobian is

873:    dF/dU + shift*dF/dUdot

875:    Most users should not need to explicitly call this routine, as it
876:    is used internally within the nonlinear solvers.

878:    Level: developer

880: .keywords: TS, compute, Jacobian, matrix

882: .seealso:  TSSetIJacobian()
883: @*/
884: PetscErrorCode TSComputeIJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal shift,Mat A,Mat B,PetscBool imex)
885: {
887:   TSIJacobian    ijacobian;
888:   TSRHSJacobian  rhsjacobian;
889:   DM             dm;
890:   void           *ctx;


901:   TSGetDM(ts,&dm);
902:   DMTSGetIJacobian(dm,&ijacobian,&ctx);
903:   DMTSGetRHSJacobian(dm,&rhsjacobian,NULL);

905:   if (!rhsjacobian && !ijacobian) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must call TSSetRHSJacobian() and / or TSSetIJacobian()");

907:   PetscLogEventBegin(TS_JacobianEval,ts,U,A,B);
908:   if (ijacobian) {
909:     PetscBool missing;
910:     PetscStackPush("TS user implicit Jacobian");
911:     (*ijacobian)(ts,t,U,Udot,shift,A,B,ctx);
912:     PetscStackPop;
913:     MatMissingDiagonal(A,&missing,NULL);
914:     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");
915:     if (B != A) {
916:       MatMissingDiagonal(B,&missing,NULL);
917:       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");
918:     }
919:   }
920:   if (imex) {
921:     if (!ijacobian) {  /* system was written as Udot = G(t,U) */
922:       PetscBool assembled;
923:       MatZeroEntries(A);
924:       MatAssembled(A,&assembled);
925:       if (!assembled) {
926:         MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
927:         MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
928:       }
929:       MatShift(A,shift);
930:       if (A != B) {
931:         MatZeroEntries(B);
932:         MatAssembled(B,&assembled);
933:         if (!assembled) {
934:           MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
935:           MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
936:         }
937:         MatShift(B,shift);
938:       }
939:     }
940:   } else {
941:     Mat Arhs = NULL,Brhs = NULL;
942:     if (rhsjacobian) {
943:       TSGetRHSMats_Private(ts,&Arhs,&Brhs);
944:       TSComputeRHSJacobian(ts,t,U,Arhs,Brhs);
945:     }
946:     if (Arhs == A) {           /* No IJacobian, so we only have the RHS matrix */
947:       PetscBool flg;
948:       ts->rhsjacobian.scale = -1;
949:       ts->rhsjacobian.shift = shift;
950:       SNESGetUseMatrixFree(ts->snes,NULL,&flg);
951:       /* since -snes_mf_operator uses the full SNES function it does not need to be shifted or scaled here */
952:       if (!flg) {
953:         MatScale(A,-1);
954:         MatShift(A,shift);
955:       }
956:       if (A != B) {
957:         MatScale(B,-1);
958:         MatShift(B,shift);
959:       }
960:     } else if (Arhs) {          /* Both IJacobian and RHSJacobian */
961:       MatStructure axpy = DIFFERENT_NONZERO_PATTERN;
962:       if (!ijacobian) {         /* No IJacobian provided, but we have a separate RHS matrix */
963:         MatZeroEntries(A);
964:         MatShift(A,shift);
965:         if (A != B) {
966:           MatZeroEntries(B);
967:           MatShift(B,shift);
968:         }
969:       }
970:       MatAXPY(A,-1,Arhs,axpy);
971:       if (A != B) {
972:         MatAXPY(B,-1,Brhs,axpy);
973:       }
974:     }
975:   }
976:   PetscLogEventEnd(TS_JacobianEval,ts,U,A,B);
977:   return(0);
978: }

980: /*@C
981:     TSSetRHSFunction - Sets the routine for evaluating the function,
982:     where U_t = G(t,u).

984:     Logically Collective on TS

986:     Input Parameters:
987: +   ts - the TS context obtained from TSCreate()
988: .   r - vector to put the computed right hand side (or NULL to have it created)
989: .   f - routine for evaluating the right-hand-side function
990: -   ctx - [optional] user-defined context for private data for the
991:           function evaluation routine (may be NULL)

993:     Calling sequence of func:
994: $     func (TS ts,PetscReal t,Vec u,Vec F,void *ctx);

996: +   t - current timestep
997: .   u - input vector
998: .   F - function vector
999: -   ctx - [optional] user-defined function context

1001:     Level: beginner

1003:     Notes: You must call this function or TSSetIFunction() to define your ODE. You cannot use this function when solving a DAE.

1005: .keywords: TS, timestep, set, right-hand-side, function

1007: .seealso: TSSetRHSJacobian(), TSSetIJacobian(), TSSetIFunction()
1008: @*/
1009: PetscErrorCode  TSSetRHSFunction(TS ts,Vec r,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx)
1010: {
1012:   SNES           snes;
1013:   Vec            ralloc = NULL;
1014:   DM             dm;


1020:   TSGetDM(ts,&dm);
1021:   DMTSSetRHSFunction(dm,f,ctx);
1022:   TSGetSNES(ts,&snes);
1023:   if (!r && !ts->dm && ts->vec_sol) {
1024:     VecDuplicate(ts->vec_sol,&ralloc);
1025:     r = ralloc;
1026:   }
1027:   SNESSetFunction(snes,r,SNESTSFormFunction,ts);
1028:   VecDestroy(&ralloc);
1029:   return(0);
1030: }

1032: /*@C
1033:     TSSetSolutionFunction - Provide a function that computes the solution of the ODE or DAE

1035:     Logically Collective on TS

1037:     Input Parameters:
1038: +   ts - the TS context obtained from TSCreate()
1039: .   f - routine for evaluating the solution
1040: -   ctx - [optional] user-defined context for private data for the
1041:           function evaluation routine (may be NULL)

1043:     Calling sequence of func:
1044: $     func (TS ts,PetscReal t,Vec u,void *ctx);

1046: +   t - current timestep
1047: .   u - output vector
1048: -   ctx - [optional] user-defined function context

1050:     Options Database:
1051: +  -ts_monitor_lg_error - create a graphical monitor of error history, requires user to have provided TSSetSolutionFunction()
1052: -  -ts_monitor_draw_error - Monitor error graphically, requires user to have provided TSSetSolutionFunction()

1054:     Notes:
1055:     This routine is used for testing accuracy of time integration schemes when you already know the solution.
1056:     If analytic solutions are not known for your system, consider using the Method of Manufactured Solutions to
1057:     create closed-form solutions with non-physical forcing terms.

1059:     For low-dimensional problems solved in serial, such as small discrete systems, TSMonitorLGError() can be used to monitor the error history.

1061:     Level: beginner

1063: .keywords: TS, timestep, set, right-hand-side, function

1065: .seealso: TSSetRHSJacobian(), TSSetIJacobian(), TSComputeSolutionFunction(), TSSetForcingFunction(), TSSetSolution(), TSGetSolution(), TSMonitorLGError(), TSMonitorDrawError()
1066: @*/
1067: PetscErrorCode  TSSetSolutionFunction(TS ts,PetscErrorCode (*f)(TS,PetscReal,Vec,void*),void *ctx)
1068: {
1070:   DM             dm;

1074:   TSGetDM(ts,&dm);
1075:   DMTSSetSolutionFunction(dm,f,ctx);
1076:   return(0);
1077: }

1079: /*@C
1080:     TSSetForcingFunction - Provide a function that computes a forcing term for a ODE or PDE

1082:     Logically Collective on TS

1084:     Input Parameters:
1085: +   ts - the TS context obtained from TSCreate()
1086: .   func - routine for evaluating the forcing function
1087: -   ctx - [optional] user-defined context for private data for the
1088:           function evaluation routine (may be NULL)

1090:     Calling sequence of func:
1091: $     func (TS ts,PetscReal t,Vec f,void *ctx);

1093: +   t - current timestep
1094: .   f - output vector
1095: -   ctx - [optional] user-defined function context

1097:     Notes:
1098:     This routine is useful for testing accuracy of time integration schemes when using the Method of Manufactured Solutions to
1099:     create closed-form solutions with a non-physical forcing term. It allows you to use the Method of Manufactored Solution without directly editing the
1100:     definition of the problem you are solving and hence possibly introducing bugs.

1102:     This replaces the ODE F(u,u_t,t) = 0 the TS is solving with F(u,u_t,t) - func(t) = 0

1104:     This forcing function does not depend on the solution to the equations, it can only depend on spatial location, time, and possibly parameters, the
1105:     parameters can be passed in the ctx variable.

1107:     For low-dimensional problems solved in serial, such as small discrete systems, TSMonitorLGError() can be used to monitor the error history.

1109:     Level: beginner

1111: .keywords: TS, timestep, set, right-hand-side, function

1113: .seealso: TSSetRHSJacobian(), TSSetIJacobian(), TSComputeSolutionFunction(), TSSetSolutionFunction()
1114: @*/
1115: PetscErrorCode  TSSetForcingFunction(TS ts,TSForcingFunction func,void *ctx)
1116: {
1118:   DM             dm;

1122:   TSGetDM(ts,&dm);
1123:   DMTSSetForcingFunction(dm,func,ctx);
1124:   return(0);
1125: }

1127: /*@C
1128:    TSSetRHSJacobian - Sets the function to compute the Jacobian of G,
1129:    where U_t = G(U,t), as well as the location to store the matrix.

1131:    Logically Collective on TS

1133:    Input Parameters:
1134: +  ts  - the TS context obtained from TSCreate()
1135: .  Amat - (approximate) Jacobian matrix
1136: .  Pmat - matrix from which preconditioner is to be constructed (usually the same as Amat)
1137: .  f   - the Jacobian evaluation routine
1138: -  ctx - [optional] user-defined context for private data for the
1139:          Jacobian evaluation routine (may be NULL)

1141:    Calling sequence of f:
1142: $     func (TS ts,PetscReal t,Vec u,Mat A,Mat B,void *ctx);

1144: +  t - current timestep
1145: .  u - input vector
1146: .  Amat - (approximate) Jacobian matrix
1147: .  Pmat - matrix from which preconditioner is to be constructed (usually the same as Amat)
1148: -  ctx - [optional] user-defined context for matrix evaluation routine

1150:    Notes:
1151:    You must set all the diagonal entries of the matrices, if they are zero you must still set them with a zero value

1153:    The TS solver may modify the nonzero structure and the entries of the matrices Amat and Pmat between the calls to f()
1154:    You should not assume the values are the same in the next call to f() as you set them in the previous call.

1156:    Level: beginner

1158: .keywords: TS, timestep, set, right-hand-side, Jacobian

1160: .seealso: SNESComputeJacobianDefaultColor(), TSSetRHSFunction(), TSRHSJacobianSetReuse(), TSSetIJacobian()

1162: @*/
1163: PetscErrorCode  TSSetRHSJacobian(TS ts,Mat Amat,Mat Pmat,TSRHSJacobian f,void *ctx)
1164: {
1166:   SNES           snes;
1167:   DM             dm;
1168:   TSIJacobian    ijacobian;


1177:   TSGetDM(ts,&dm);
1178:   DMTSSetRHSJacobian(dm,f,ctx);
1179:   if (f == TSComputeRHSJacobianConstant) {
1180:     /* Handle this case automatically for the user; otherwise user should call themselves. */
1181:     TSRHSJacobianSetReuse(ts,PETSC_TRUE);
1182:   }
1183:   DMTSGetIJacobian(dm,&ijacobian,NULL);
1184:   TSGetSNES(ts,&snes);
1185:   if (!ijacobian) {
1186:     SNESSetJacobian(snes,Amat,Pmat,SNESTSFormJacobian,ts);
1187:   }
1188:   if (Amat) {
1189:     PetscObjectReference((PetscObject)Amat);
1190:     MatDestroy(&ts->Arhs);
1191:     ts->Arhs = Amat;
1192:   }
1193:   if (Pmat) {
1194:     PetscObjectReference((PetscObject)Pmat);
1195:     MatDestroy(&ts->Brhs);
1196:     ts->Brhs = Pmat;
1197:   }
1198:   return(0);
1199: }

1201: /*@C
1202:    TSSetIFunction - Set the function to compute F(t,U,U_t) where F() = 0 is the DAE to be solved.

1204:    Logically Collective on TS

1206:    Input Parameters:
1207: +  ts  - the TS context obtained from TSCreate()
1208: .  r   - vector to hold the residual (or NULL to have it created internally)
1209: .  f   - the function evaluation routine
1210: -  ctx - user-defined context for private data for the function evaluation routine (may be NULL)

1212:    Calling sequence of f:
1213: $  f(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx);

1215: +  t   - time at step/stage being solved
1216: .  u   - state vector
1217: .  u_t - time derivative of state vector
1218: .  F   - function vector
1219: -  ctx - [optional] user-defined context for matrix evaluation routine

1221:    Important:
1222:    The user MUST call either this routine or TSSetRHSFunction() to define the ODE.  When solving DAEs you must use this function.

1224:    Level: beginner

1226: .keywords: TS, timestep, set, DAE, Jacobian

1228: .seealso: TSSetRHSJacobian(), TSSetRHSFunction(), TSSetIJacobian()
1229: @*/
1230: PetscErrorCode  TSSetIFunction(TS ts,Vec r,TSIFunction f,void *ctx)
1231: {
1233:   SNES           snes;
1234:   Vec            ralloc = NULL;
1235:   DM             dm;


1241:   TSGetDM(ts,&dm);
1242:   DMTSSetIFunction(dm,f,ctx);

1244:   TSGetSNES(ts,&snes);
1245:   if (!r && !ts->dm && ts->vec_sol) {
1246:     VecDuplicate(ts->vec_sol,&ralloc);
1247:     r  = ralloc;
1248:   }
1249:   SNESSetFunction(snes,r,SNESTSFormFunction,ts);
1250:   VecDestroy(&ralloc);
1251:   return(0);
1252: }

1254: /*@C
1255:    TSGetIFunction - Returns the vector where the implicit residual is stored and the function/contex to compute it.

1257:    Not Collective

1259:    Input Parameter:
1260: .  ts - the TS context

1262:    Output Parameter:
1263: +  r - vector to hold residual (or NULL)
1264: .  func - the function to compute residual (or NULL)
1265: -  ctx - the function context (or NULL)

1267:    Level: advanced

1269: .keywords: TS, nonlinear, get, function

1271: .seealso: TSSetIFunction(), SNESGetFunction()
1272: @*/
1273: PetscErrorCode TSGetIFunction(TS ts,Vec *r,TSIFunction *func,void **ctx)
1274: {
1276:   SNES           snes;
1277:   DM             dm;

1281:   TSGetSNES(ts,&snes);
1282:   SNESGetFunction(snes,r,NULL,NULL);
1283:   TSGetDM(ts,&dm);
1284:   DMTSGetIFunction(dm,func,ctx);
1285:   return(0);
1286: }

1288: /*@C
1289:    TSGetRHSFunction - Returns the vector where the right hand side is stored and the function/context to compute it.

1291:    Not Collective

1293:    Input Parameter:
1294: .  ts - the TS context

1296:    Output Parameter:
1297: +  r - vector to hold computed right hand side (or NULL)
1298: .  func - the function to compute right hand side (or NULL)
1299: -  ctx - the function context (or NULL)

1301:    Level: advanced

1303: .keywords: TS, nonlinear, get, function

1305: .seealso: TSSetRHSFunction(), SNESGetFunction()
1306: @*/
1307: PetscErrorCode TSGetRHSFunction(TS ts,Vec *r,TSRHSFunction *func,void **ctx)
1308: {
1310:   SNES           snes;
1311:   DM             dm;

1315:   TSGetSNES(ts,&snes);
1316:   SNESGetFunction(snes,r,NULL,NULL);
1317:   TSGetDM(ts,&dm);
1318:   DMTSGetRHSFunction(dm,func,ctx);
1319:   return(0);
1320: }

1322: /*@C
1323:    TSSetIJacobian - Set the function to compute the matrix dF/dU + a*dF/dU_t where F(t,U,U_t) is the function
1324:         provided with TSSetIFunction().

1326:    Logically Collective on TS

1328:    Input Parameters:
1329: +  ts  - the TS context obtained from TSCreate()
1330: .  Amat - (approximate) Jacobian matrix
1331: .  Pmat - matrix used to compute preconditioner (usually the same as Amat)
1332: .  f   - the Jacobian evaluation routine
1333: -  ctx - user-defined context for private data for the Jacobian evaluation routine (may be NULL)

1335:    Calling sequence of f:
1336: $  f(TS ts,PetscReal t,Vec U,Vec U_t,PetscReal a,Mat Amat,Mat Pmat,void *ctx);

1338: +  t    - time at step/stage being solved
1339: .  U    - state vector
1340: .  U_t  - time derivative of state vector
1341: .  a    - shift
1342: .  Amat - (approximate) Jacobian of F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t
1343: .  Pmat - matrix used for constructing preconditioner, usually the same as Amat
1344: -  ctx  - [optional] user-defined context for matrix evaluation routine

1346:    Notes:
1347:    The matrices Amat and Pmat are exactly the matrices that are used by SNES for the nonlinear solve.

1349:    If you know the operator Amat has a null space you can use MatSetNullSpace() and MatSetTransposeNullSpace() to supply the null
1350:    space to Amat and the KSP solvers will automatically use that null space as needed during the solution process.

1352:    The matrix dF/dU + a*dF/dU_t you provide turns out to be
1353:    the Jacobian of F(t,U,W+a*U) where F(t,U,U_t) = 0 is the DAE to be solved.
1354:    The time integrator internally approximates U_t by W+a*U where the positive "shift"
1355:    a and vector W depend on the integration method, step size, and past states. For example with
1356:    the backward Euler method a = 1/dt and W = -a*U(previous timestep) so
1357:    W + a*U = a*(U - U(previous timestep)) = (U - U(previous timestep))/dt

1359:    You must set all the diagonal entries of the matrices, if they are zero you must still set them with a zero value

1361:    The TS solver may modify the nonzero structure and the entries of the matrices Amat and Pmat between the calls to f()
1362:    You should not assume the values are the same in the next call to f() as you set them in the previous call.

1364:    Level: beginner

1366: .keywords: TS, timestep, DAE, Jacobian

1368: .seealso: TSSetIFunction(), TSSetRHSJacobian(), SNESComputeJacobianDefaultColor(), SNESComputeJacobianDefault(), TSSetRHSFunction()

1370: @*/
1371: PetscErrorCode  TSSetIJacobian(TS ts,Mat Amat,Mat Pmat,TSIJacobian f,void *ctx)
1372: {
1374:   SNES           snes;
1375:   DM             dm;


1384:   TSGetDM(ts,&dm);
1385:   DMTSSetIJacobian(dm,f,ctx);

1387:   TSGetSNES(ts,&snes);
1388:   SNESSetJacobian(snes,Amat,Pmat,SNESTSFormJacobian,ts);
1389:   return(0);
1390: }

1392: /*@
1393:    TSRHSJacobianSetReuse - restore RHS Jacobian before re-evaluating.  Without this flag, TS will change the sign and
1394:    shift the RHS Jacobian for a finite-time-step implicit solve, in which case the user function will need to recompute
1395:    the entire Jacobian.  The reuse flag must be set if the evaluation function will assume that the matrix entries have
1396:    not been changed by the TS.

1398:    Logically Collective

1400:    Input Arguments:
1401: +  ts - TS context obtained from TSCreate()
1402: -  reuse - PETSC_TRUE if the RHS Jacobian

1404:    Level: intermediate

1406: .seealso: TSSetRHSJacobian(), TSComputeRHSJacobianConstant()
1407: @*/
1408: PetscErrorCode TSRHSJacobianSetReuse(TS ts,PetscBool reuse)
1409: {
1411:   ts->rhsjacobian.reuse = reuse;
1412:   return(0);
1413: }

1415: /*@C
1416:    TSSetI2Function - Set the function to compute F(t,U,U_t,U_tt) where F = 0 is the DAE to be solved.

1418:    Logically Collective on TS

1420:    Input Parameters:
1421: +  ts  - the TS context obtained from TSCreate()
1422: .  F   - vector to hold the residual (or NULL to have it created internally)
1423: .  fun - the function evaluation routine
1424: -  ctx - user-defined context for private data for the function evaluation routine (may be NULL)

1426:    Calling sequence of fun:
1427: $  fun(TS ts,PetscReal t,Vec U,Vec U_t,Vec U_tt,Vec F,ctx);

1429: +  t    - time at step/stage being solved
1430: .  U    - state vector
1431: .  U_t  - time derivative of state vector
1432: .  U_tt - second time derivative of state vector
1433: .  F    - function vector
1434: -  ctx  - [optional] user-defined context for matrix evaluation routine (may be NULL)

1436:    Level: beginner

1438: .keywords: TS, timestep, set, ODE, DAE, Function

1440: .seealso: TSSetI2Jacobian()
1441: @*/
1442: PetscErrorCode TSSetI2Function(TS ts,Vec F,TSI2Function fun,void *ctx)
1443: {
1444:   DM             dm;

1450:   TSSetIFunction(ts,F,NULL,NULL);
1451:   TSGetDM(ts,&dm);
1452:   DMTSSetI2Function(dm,fun,ctx);
1453:   return(0);
1454: }

1456: /*@C
1457:   TSGetI2Function - Returns the vector where the implicit residual is stored and the function/contex to compute it.

1459:   Not Collective

1461:   Input Parameter:
1462: . ts - the TS context

1464:   Output Parameter:
1465: + r - vector to hold residual (or NULL)
1466: . fun - the function to compute residual (or NULL)
1467: - ctx - the function context (or NULL)

1469:   Level: advanced

1471: .keywords: TS, nonlinear, get, function

1473: .seealso: TSSetI2Function(), SNESGetFunction()
1474: @*/
1475: PetscErrorCode TSGetI2Function(TS ts,Vec *r,TSI2Function *fun,void **ctx)
1476: {
1478:   SNES           snes;
1479:   DM             dm;

1483:   TSGetSNES(ts,&snes);
1484:   SNESGetFunction(snes,r,NULL,NULL);
1485:   TSGetDM(ts,&dm);
1486:   DMTSGetI2Function(dm,fun,ctx);
1487:   return(0);
1488: }

1490: /*@C
1491:    TSSetI2Jacobian - Set the function to compute the matrix dF/dU + v*dF/dU_t  + a*dF/dU_tt
1492:         where F(t,U,U_t,U_tt) is the function you provided with TSSetI2Function().

1494:    Logically Collective on TS

1496:    Input Parameters:
1497: +  ts  - the TS context obtained from TSCreate()
1498: .  J   - Jacobian matrix
1499: .  P   - preconditioning matrix for J (may be same as J)
1500: .  jac - the Jacobian evaluation routine
1501: -  ctx - user-defined context for private data for the Jacobian evaluation routine (may be NULL)

1503:    Calling sequence of jac:
1504: $  jac(TS ts,PetscReal t,Vec U,Vec U_t,Vec U_tt,PetscReal v,PetscReal a,Mat J,Mat P,void *ctx);

1506: +  t    - time at step/stage being solved
1507: .  U    - state vector
1508: .  U_t  - time derivative of state vector
1509: .  U_tt - second time derivative of state vector
1510: .  v    - shift for U_t
1511: .  a    - shift for U_tt
1512: .  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
1513: .  P    - preconditioning matrix for J, may be same as J
1514: -  ctx  - [optional] user-defined context for matrix evaluation routine

1516:    Notes:
1517:    The matrices J and P are exactly the matrices that are used by SNES for the nonlinear solve.

1519:    The matrix dF/dU + v*dF/dU_t + a*dF/dU_tt you provide turns out to be
1520:    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.
1521:    The time integrator internally approximates U_t by W+v*U and U_tt by W'+a*U  where the positive "shift"
1522:    parameters 'v' and 'a' and vectors W, W' depend on the integration method, step size, and past states.

1524:    Level: beginner

1526: .keywords: TS, timestep, set, ODE, DAE, Jacobian

1528: .seealso: TSSetI2Function()
1529: @*/
1530: PetscErrorCode TSSetI2Jacobian(TS ts,Mat J,Mat P,TSI2Jacobian jac,void *ctx)
1531: {
1532:   DM             dm;

1539:   TSSetIJacobian(ts,J,P,NULL,NULL);
1540:   TSGetDM(ts,&dm);
1541:   DMTSSetI2Jacobian(dm,jac,ctx);
1542:   return(0);
1543: }

1545: /*@C
1546:   TSGetI2Jacobian - Returns the implicit Jacobian at the present timestep.

1548:   Not Collective, but parallel objects are returned if TS is parallel

1550:   Input Parameter:
1551: . ts  - The TS context obtained from TSCreate()

1553:   Output Parameters:
1554: + J  - The (approximate) Jacobian of F(t,U,U_t,U_tt)
1555: . P - The matrix from which the preconditioner is constructed, often the same as J
1556: . jac - The function to compute the Jacobian matrices
1557: - ctx - User-defined context for Jacobian evaluation routine

1559:   Notes: You can pass in NULL for any return argument you do not need.

1561:   Level: advanced

1563: .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetStepNumber()

1565: .keywords: TS, timestep, get, matrix, Jacobian
1566: @*/
1567: PetscErrorCode  TSGetI2Jacobian(TS ts,Mat *J,Mat *P,TSI2Jacobian *jac,void **ctx)
1568: {
1570:   SNES           snes;
1571:   DM             dm;

1574:   TSGetSNES(ts,&snes);
1575:   SNESSetUpMatrices(snes);
1576:   SNESGetJacobian(snes,J,P,NULL,NULL);
1577:   TSGetDM(ts,&dm);
1578:   DMTSGetI2Jacobian(dm,jac,ctx);
1579:   return(0);
1580: }

1582: /*@
1583:   TSComputeI2Function - Evaluates the DAE residual written in implicit form F(t,U,U_t,U_tt) = 0

1585:   Collective on TS and Vec

1587:   Input Parameters:
1588: + ts - the TS context
1589: . t - current time
1590: . U - state vector
1591: . V - time derivative of state vector (U_t)
1592: - A - second time derivative of state vector (U_tt)

1594:   Output Parameter:
1595: . F - the residual vector

1597:   Note:
1598:   Most users should not need to explicitly call this routine, as it
1599:   is used internally within the nonlinear solvers.

1601:   Level: developer

1603: .keywords: TS, compute, function, vector

1605: .seealso: TSSetI2Function()
1606: @*/
1607: PetscErrorCode TSComputeI2Function(TS ts,PetscReal t,Vec U,Vec V,Vec A,Vec F)
1608: {
1609:   DM             dm;
1610:   TSI2Function   I2Function;
1611:   void           *ctx;
1612:   TSRHSFunction  rhsfunction;


1622:   TSGetDM(ts,&dm);
1623:   DMTSGetI2Function(dm,&I2Function,&ctx);
1624:   DMTSGetRHSFunction(dm,&rhsfunction,NULL);

1626:   if (!I2Function) {
1627:     TSComputeIFunction(ts,t,U,A,F,PETSC_FALSE);
1628:     return(0);
1629:   }

1631:   PetscLogEventBegin(TS_FunctionEval,ts,U,V,F);

1633:   PetscStackPush("TS user implicit function");
1634:   I2Function(ts,t,U,V,A,F,ctx);
1635:   PetscStackPop;

1637:   if (rhsfunction) {
1638:     Vec Frhs;
1639:     TSGetRHSVec_Private(ts,&Frhs);
1640:     TSComputeRHSFunction(ts,t,U,Frhs);
1641:     VecAXPY(F,-1,Frhs);
1642:   }

1644:   PetscLogEventEnd(TS_FunctionEval,ts,U,V,F);
1645:   return(0);
1646: }

1648: /*@
1649:   TSComputeI2Jacobian - Evaluates the Jacobian of the DAE

1651:   Collective on TS and Vec

1653:   Input Parameters:
1654: + ts - the TS context
1655: . t - current timestep
1656: . U - state vector
1657: . V - time derivative of state vector
1658: . A - second time derivative of state vector
1659: . shiftV - shift to apply, see note below
1660: - shiftA - shift to apply, see note below

1662:   Output Parameters:
1663: + J - Jacobian matrix
1664: - P - optional preconditioning matrix

1666:   Notes:
1667:   If F(t,U,V,A)=0 is the DAE, the required Jacobian is

1669:   dF/dU + shiftV*dF/dV + shiftA*dF/dA

1671:   Most users should not need to explicitly call this routine, as it
1672:   is used internally within the nonlinear solvers.

1674:   Level: developer

1676: .keywords: TS, compute, Jacobian, matrix

1678: .seealso:  TSSetI2Jacobian()
1679: @*/
1680: PetscErrorCode TSComputeI2Jacobian(TS ts,PetscReal t,Vec U,Vec V,Vec A,PetscReal shiftV,PetscReal shiftA,Mat J,Mat P)
1681: {
1682:   DM             dm;
1683:   TSI2Jacobian   I2Jacobian;
1684:   void           *ctx;
1685:   TSRHSJacobian  rhsjacobian;


1696:   TSGetDM(ts,&dm);
1697:   DMTSGetI2Jacobian(dm,&I2Jacobian,&ctx);
1698:   DMTSGetRHSJacobian(dm,&rhsjacobian,NULL);

1700:   if (!I2Jacobian) {
1701:     TSComputeIJacobian(ts,t,U,A,shiftA,J,P,PETSC_FALSE);
1702:     return(0);
1703:   }

1705:   PetscLogEventBegin(TS_JacobianEval,ts,U,J,P);

1707:   PetscStackPush("TS user implicit Jacobian");
1708:   I2Jacobian(ts,t,U,V,A,shiftV,shiftA,J,P,ctx);
1709:   PetscStackPop;

1711:   if (rhsjacobian) {
1712:     Mat Jrhs,Prhs; MatStructure axpy = DIFFERENT_NONZERO_PATTERN;
1713:     TSGetRHSMats_Private(ts,&Jrhs,&Prhs);
1714:     TSComputeRHSJacobian(ts,t,U,Jrhs,Prhs);
1715:     MatAXPY(J,-1,Jrhs,axpy);
1716:     if (P != J) {MatAXPY(P,-1,Prhs,axpy);}
1717:   }

1719:   PetscLogEventEnd(TS_JacobianEval,ts,U,J,P);
1720:   return(0);
1721: }

1723: /*@
1724:    TS2SetSolution - Sets the initial solution and time derivative vectors
1725:    for use by the TS routines handling second order equations.

1727:    Logically Collective on TS and Vec

1729:    Input Parameters:
1730: +  ts - the TS context obtained from TSCreate()
1731: .  u - the solution vector
1732: -  v - the time derivative vector

1734:    Level: beginner

1736: .keywords: TS, timestep, set, solution, initial conditions
1737: @*/
1738: PetscErrorCode  TS2SetSolution(TS ts,Vec u,Vec v)
1739: {

1746:   TSSetSolution(ts,u);
1747:   PetscObjectReference((PetscObject)v);
1748:   VecDestroy(&ts->vec_dot);
1749:   ts->vec_dot = v;
1750:   return(0);
1751: }

1753: /*@
1754:    TS2GetSolution - Returns the solution and time derivative at the present timestep
1755:    for second order equations. It is valid to call this routine inside the function
1756:    that you are evaluating in order to move to the new timestep. This vector not
1757:    changed until the solution at the next timestep has been calculated.

1759:    Not Collective, but Vec returned is parallel if TS is parallel

1761:    Input Parameter:
1762: .  ts - the TS context obtained from TSCreate()

1764:    Output Parameter:
1765: +  u - the vector containing the solution
1766: -  v - the vector containing the time derivative

1768:    Level: intermediate

1770: .seealso: TS2SetSolution(), TSGetTimeStep(), TSGetTime()

1772: .keywords: TS, timestep, get, solution
1773: @*/
1774: PetscErrorCode  TS2GetSolution(TS ts,Vec *u,Vec *v)
1775: {
1780:   if (u) *u = ts->vec_sol;
1781:   if (v) *v = ts->vec_dot;
1782:   return(0);
1783: }

1785: /*@C
1786:   TSLoad - Loads a KSP that has been stored in binary  with KSPView().

1788:   Collective on PetscViewer

1790:   Input Parameters:
1791: + newdm - the newly loaded TS, this needs to have been created with TSCreate() or
1792:            some related function before a call to TSLoad().
1793: - viewer - binary file viewer, obtained from PetscViewerBinaryOpen()

1795:    Level: intermediate

1797:   Notes:
1798:    The type is determined by the data in the file, any type set into the TS before this call is ignored.

1800:   Notes for advanced users:
1801:   Most users should not need to know the details of the binary storage
1802:   format, since TSLoad() and TSView() completely hide these details.
1803:   But for anyone who's interested, the standard binary matrix storage
1804:   format is
1805: .vb
1806:      has not yet been determined
1807: .ve

1809: .seealso: PetscViewerBinaryOpen(), TSView(), MatLoad(), VecLoad()
1810: @*/
1811: PetscErrorCode  TSLoad(TS ts, PetscViewer viewer)
1812: {
1814:   PetscBool      isbinary;
1815:   PetscInt       classid;
1816:   char           type[256];
1817:   DMTS           sdm;
1818:   DM             dm;

1823:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
1824:   if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");

1826:   PetscViewerBinaryRead(viewer,&classid,1,NULL,PETSC_INT);
1827:   if (classid != TS_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONG,"Not TS next in file");
1828:   PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);
1829:   TSSetType(ts, type);
1830:   if (ts->ops->load) {
1831:     (*ts->ops->load)(ts,viewer);
1832:   }
1833:   DMCreate(PetscObjectComm((PetscObject)ts),&dm);
1834:   DMLoad(dm,viewer);
1835:   TSSetDM(ts,dm);
1836:   DMCreateGlobalVector(ts->dm,&ts->vec_sol);
1837:   VecLoad(ts->vec_sol,viewer);
1838:   DMGetDMTS(ts->dm,&sdm);
1839:   DMTSLoad(sdm,viewer);
1840:   return(0);
1841: }

1843:  #include <petscdraw.h>
1844: #if defined(PETSC_HAVE_SAWS)
1845:  #include <petscviewersaws.h>
1846: #endif
1847: /*@C
1848:     TSView - Prints the TS data structure.

1850:     Collective on TS

1852:     Input Parameters:
1853: +   ts - the TS context obtained from TSCreate()
1854: -   viewer - visualization context

1856:     Options Database Key:
1857: .   -ts_view - calls TSView() at end of TSStep()

1859:     Notes:
1860:     The available visualization contexts include
1861: +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
1862: -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
1863:          output where only the first processor opens
1864:          the file.  All other processors send their
1865:          data to the first processor to print.

1867:     The user can open an alternative visualization context with
1868:     PetscViewerASCIIOpen() - output to a specified file.

1870:     Level: beginner

1872: .keywords: TS, timestep, view

1874: .seealso: PetscViewerASCIIOpen()
1875: @*/
1876: PetscErrorCode  TSView(TS ts,PetscViewer viewer)
1877: {
1879:   TSType         type;
1880:   PetscBool      iascii,isstring,isundials,isbinary,isdraw;
1881:   DMTS           sdm;
1882: #if defined(PETSC_HAVE_SAWS)
1883:   PetscBool      issaws;
1884: #endif

1888:   if (!viewer) {
1889:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ts),&viewer);
1890:   }

1894:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1895:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
1896:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
1897:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
1898: #if defined(PETSC_HAVE_SAWS)
1899:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSAWS,&issaws);
1900: #endif
1901:   if (iascii) {
1902:     PetscObjectPrintClassNamePrefixType((PetscObject)ts,viewer);
1903:     if (ts->ops->view) {
1904:       PetscViewerASCIIPushTab(viewer);
1905:       (*ts->ops->view)(ts,viewer);
1906:       PetscViewerASCIIPopTab(viewer);
1907:     }
1908:     if (ts->max_steps < PETSC_MAX_INT) {
1909:       PetscViewerASCIIPrintf(viewer,"  maximum steps=%D\n",ts->max_steps);
1910:     }
1911:     if (ts->max_time < PETSC_MAX_REAL) {
1912:       PetscViewerASCIIPrintf(viewer,"  maximum time=%g\n",(double)ts->max_time);
1913:     }
1914:     if (ts->usessnes) {
1915:       PetscBool lin;
1916:       if (ts->problem_type == TS_NONLINEAR) {
1917:         PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solver iterations=%D\n",ts->snes_its);
1918:       }
1919:       PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%D\n",ts->ksp_its);
1920:       PetscObjectTypeCompare((PetscObject)ts->snes,SNESKSPONLY,&lin);
1921:       PetscViewerASCIIPrintf(viewer,"  total number of %slinear solve failures=%D\n",lin ? "" : "non",ts->num_snes_failures);
1922:     }
1923:     PetscViewerASCIIPrintf(viewer,"  total number of rejected steps=%D\n",ts->reject);
1924:     if (ts->vrtol) {
1925:       PetscViewerASCIIPrintf(viewer,"  using vector of relative error tolerances, ");
1926:     } else {
1927:       PetscViewerASCIIPrintf(viewer,"  using relative error tolerance of %g, ",(double)ts->rtol);
1928:     }
1929:     if (ts->vatol) {
1930:       PetscViewerASCIIPrintf(viewer,"  using vector of absolute error tolerances\n");
1931:     } else {
1932:       PetscViewerASCIIPrintf(viewer,"  using absolute error tolerance of %g\n",(double)ts->atol);
1933:     }
1934:     PetscViewerASCIIPushTab(viewer);
1935:     TSAdaptView(ts->adapt,viewer);
1936:     PetscViewerASCIIPopTab(viewer);
1937:     if (ts->snes && ts->usessnes)  {
1938:       PetscViewerASCIIPushTab(viewer);
1939:       SNESView(ts->snes,viewer);
1940:       PetscViewerASCIIPopTab(viewer);
1941:     }
1942:     DMGetDMTS(ts->dm,&sdm);
1943:     DMTSView(sdm,viewer);
1944:   } else if (isstring) {
1945:     TSGetType(ts,&type);
1946:     PetscViewerStringSPrintf(viewer," %-7.7s",type);
1947:   } else if (isbinary) {
1948:     PetscInt    classid = TS_FILE_CLASSID;
1949:     MPI_Comm    comm;
1950:     PetscMPIInt rank;
1951:     char        type[256];

1953:     PetscObjectGetComm((PetscObject)ts,&comm);
1954:     MPI_Comm_rank(comm,&rank);
1955:     if (!rank) {
1956:       PetscViewerBinaryWrite(viewer,&classid,1,PETSC_INT,PETSC_FALSE);
1957:       PetscStrncpy(type,((PetscObject)ts)->type_name,256);
1958:       PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);
1959:     }
1960:     if (ts->ops->view) {
1961:       (*ts->ops->view)(ts,viewer);
1962:     }
1963:     if (ts->adapt) {TSAdaptView(ts->adapt,viewer);}
1964:     DMView(ts->dm,viewer);
1965:     VecView(ts->vec_sol,viewer);
1966:     DMGetDMTS(ts->dm,&sdm);
1967:     DMTSView(sdm,viewer);
1968:   } else if (isdraw) {
1969:     PetscDraw draw;
1970:     char      str[36];
1971:     PetscReal x,y,bottom,h;

1973:     PetscViewerDrawGetDraw(viewer,0,&draw);
1974:     PetscDrawGetCurrentPoint(draw,&x,&y);
1975:     PetscStrcpy(str,"TS: ");
1976:     PetscStrcat(str,((PetscObject)ts)->type_name);
1977:     PetscDrawStringBoxed(draw,x,y,PETSC_DRAW_BLACK,PETSC_DRAW_BLACK,str,NULL,&h);
1978:     bottom = y - h;
1979:     PetscDrawPushCurrentPoint(draw,x,bottom);
1980:     if (ts->ops->view) {
1981:       (*ts->ops->view)(ts,viewer);
1982:     }
1983:     if (ts->adapt) {TSAdaptView(ts->adapt,viewer);}
1984:     if (ts->snes)  {SNESView(ts->snes,viewer);}
1985:     PetscDrawPopCurrentPoint(draw);
1986: #if defined(PETSC_HAVE_SAWS)
1987:   } else if (issaws) {
1988:     PetscMPIInt rank;
1989:     const char  *name;

1991:     PetscObjectGetName((PetscObject)ts,&name);
1992:     MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
1993:     if (!((PetscObject)ts)->amsmem && !rank) {
1994:       char       dir[1024];

1996:       PetscObjectViewSAWs((PetscObject)ts,viewer);
1997:       PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/time_step",name);
1998:       PetscStackCallSAWs(SAWs_Register,(dir,&ts->steps,1,SAWs_READ,SAWs_INT));
1999:       PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/time",name);
2000:       PetscStackCallSAWs(SAWs_Register,(dir,&ts->ptime,1,SAWs_READ,SAWs_DOUBLE));
2001:     }
2002:     if (ts->ops->view) {
2003:       (*ts->ops->view)(ts,viewer);
2004:     }
2005: #endif
2006:   }

2008:   PetscViewerASCIIPushTab(viewer);
2009:   PetscObjectTypeCompare((PetscObject)ts,TSSUNDIALS,&isundials);
2010:   PetscViewerASCIIPopTab(viewer);
2011:   return(0);
2012: }

2014: /*@
2015:    TSSetApplicationContext - Sets an optional user-defined context for
2016:    the timesteppers.

2018:    Logically Collective on TS

2020:    Input Parameters:
2021: +  ts - the TS context obtained from TSCreate()
2022: -  usrP - optional user context

2024:    Fortran Notes: To use this from Fortran you must write a Fortran interface definition for this
2025:     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.

2027:    Level: intermediate

2029: .keywords: TS, timestep, set, application, context

2031: .seealso: TSGetApplicationContext()
2032: @*/
2033: PetscErrorCode  TSSetApplicationContext(TS ts,void *usrP)
2034: {
2037:   ts->user = usrP;
2038:   return(0);
2039: }

2041: /*@
2042:     TSGetApplicationContext - Gets the user-defined context for the
2043:     timestepper.

2045:     Not Collective

2047:     Input Parameter:
2048: .   ts - the TS context obtained from TSCreate()

2050:     Output Parameter:
2051: .   usrP - user context

2053:    Fortran Notes: To use this from Fortran you must write a Fortran interface definition for this
2054:     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.

2056:     Level: intermediate

2058: .keywords: TS, timestep, get, application, context

2060: .seealso: TSSetApplicationContext()
2061: @*/
2062: PetscErrorCode  TSGetApplicationContext(TS ts,void *usrP)
2063: {
2066:   *(void**)usrP = ts->user;
2067:   return(0);
2068: }

2070: /*@
2071:    TSGetStepNumber - Gets the number of steps completed.

2073:    Not Collective

2075:    Input Parameter:
2076: .  ts - the TS context obtained from TSCreate()

2078:    Output Parameter:
2079: .  steps - number of steps completed so far

2081:    Level: intermediate

2083: .keywords: TS, timestep, get, iteration, number
2084: .seealso: TSGetTime(), TSGetTimeStep(), TSSetPreStep(), TSSetPreStage(), TSSetPostStage(), TSSetPostStep()
2085: @*/
2086: PetscErrorCode TSGetStepNumber(TS ts,PetscInt *steps)
2087: {
2091:   *steps = ts->steps;
2092:   return(0);
2093: }

2095: /*@
2096:    TSSetStepNumber - Sets the number of steps completed.

2098:    Logically Collective on TS

2100:    Input Parameters:
2101: +  ts - the TS context
2102: -  steps - number of steps completed so far

2104:    Notes:
2105:    For most uses of the TS solvers the user need not explicitly call
2106:    TSSetStepNumber(), as the step counter is appropriately updated in
2107:    TSSolve()/TSStep()/TSRollBack(). Power users may call this routine to
2108:    reinitialize timestepping by setting the step counter to zero (and time
2109:    to the initial time) to solve a similar problem with different initial
2110:    conditions or parameters. Other possible use case is to continue
2111:    timestepping from a previously interrupted run in such a way that TS
2112:    monitors will be called with a initial nonzero step counter.

2114:    Level: advanced

2116: .keywords: TS, timestep, set, iteration, number
2117: .seealso: TSGetStepNumber(), TSSetTime(), TSSetTimeStep(), TSSetSolution()
2118: @*/
2119: PetscErrorCode TSSetStepNumber(TS ts,PetscInt steps)
2120: {
2124:   if (steps < 0) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Step number must be non-negative");
2125:   ts->steps = steps;
2126:   return(0);
2127: }

2129: /*@
2130:    TSSetTimeStep - Allows one to reset the timestep at any time,
2131:    useful for simple pseudo-timestepping codes.

2133:    Logically Collective on TS

2135:    Input Parameters:
2136: +  ts - the TS context obtained from TSCreate()
2137: -  time_step - the size of the timestep

2139:    Level: intermediate

2141: .seealso: TSGetTimeStep(), TSSetTime()

2143: .keywords: TS, set, timestep
2144: @*/
2145: PetscErrorCode  TSSetTimeStep(TS ts,PetscReal time_step)
2146: {
2150:   ts->time_step = time_step;
2151:   return(0);
2152: }

2154: /*@
2155:    TSSetExactFinalTime - Determines whether to adapt the final time step to
2156:      match the exact final time, interpolate solution to the exact final time,
2157:      or just return at the final time TS computed.

2159:   Logically Collective on TS

2161:    Input Parameter:
2162: +   ts - the time-step context
2163: -   eftopt - exact final time option

2165: $  TS_EXACTFINALTIME_STEPOVER    - Don't do anything if final time is exceeded
2166: $  TS_EXACTFINALTIME_INTERPOLATE - Interpolate back to final time
2167: $  TS_EXACTFINALTIME_MATCHSTEP - Adapt final time step to match the final time

2169:    Options Database:
2170: .   -ts_exact_final_time <stepover,interpolate,matchstep> - select the final step at runtime

2172:    Warning: If you use the option TS_EXACTFINALTIME_STEPOVER the solution may be at a very different time
2173:     then the final time you selected.

2175:    Level: beginner

2177: .seealso: TSExactFinalTimeOption, TSGetExactFinalTime()
2178: @*/
2179: PetscErrorCode TSSetExactFinalTime(TS ts,TSExactFinalTimeOption eftopt)
2180: {
2184:   ts->exact_final_time = eftopt;
2185:   return(0);
2186: }

2188: /*@
2189:    TSGetExactFinalTime - Gets the exact final time option.

2191:    Not Collective

2193:    Input Parameter:
2194: .  ts - the TS context

2196:    Output Parameter:
2197: .  eftopt - exact final time option

2199:    Level: beginner

2201: .seealso: TSExactFinalTimeOption, TSSetExactFinalTime()
2202: @*/
2203: PetscErrorCode TSGetExactFinalTime(TS ts,TSExactFinalTimeOption *eftopt)
2204: {
2208:   *eftopt = ts->exact_final_time;
2209:   return(0);
2210: }

2212: /*@
2213:    TSGetTimeStep - Gets the current timestep size.

2215:    Not Collective

2217:    Input Parameter:
2218: .  ts - the TS context obtained from TSCreate()

2220:    Output Parameter:
2221: .  dt - the current timestep size

2223:    Level: intermediate

2225: .seealso: TSSetTimeStep(), TSGetTime()

2227: .keywords: TS, get, timestep
2228: @*/
2229: PetscErrorCode  TSGetTimeStep(TS ts,PetscReal *dt)
2230: {
2234:   *dt = ts->time_step;
2235:   return(0);
2236: }

2238: /*@
2239:    TSGetSolution - Returns the solution at the present timestep. It
2240:    is valid to call this routine inside the function that you are evaluating
2241:    in order to move to the new timestep. This vector not changed until
2242:    the solution at the next timestep has been calculated.

2244:    Not Collective, but Vec returned is parallel if TS is parallel

2246:    Input Parameter:
2247: .  ts - the TS context obtained from TSCreate()

2249:    Output Parameter:
2250: .  v - the vector containing the solution

2252:    Note: If you used TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP); this does not return the solution at the requested
2253:    final time. It returns the solution at the next timestep.

2255:    Level: intermediate

2257: .seealso: TSGetTimeStep(), TSGetTime(), TSGetSolveTime(), TSGetSolutionComponents(), TSSetSolutionFunction()

2259: .keywords: TS, timestep, get, solution
2260: @*/
2261: PetscErrorCode  TSGetSolution(TS ts,Vec *v)
2262: {
2266:   *v = ts->vec_sol;
2267:   return(0);
2268: }

2270: /*@
2271:    TSGetSolutionComponents - Returns any solution components at the present
2272:    timestep, if available for the time integration method being used.
2273:    Solution components are quantities that share the same size and
2274:    structure as the solution vector.

2276:    Not Collective, but Vec returned is parallel if TS is parallel

2278:    Parameters :
2279: .  ts - the TS context obtained from TSCreate() (input parameter).
2280: .  n - If v is PETSC_NULL, then the number of solution components is
2281:        returned through n, else the n-th solution component is
2282:        returned in v.
2283: .  v - the vector containing the n-th solution component
2284:        (may be PETSC_NULL to use this function to find out
2285:         the number of solutions components).

2287:    Level: advanced

2289: .seealso: TSGetSolution()

2291: .keywords: TS, timestep, get, solution
2292: @*/
2293: PetscErrorCode  TSGetSolutionComponents(TS ts,PetscInt *n,Vec *v)
2294: {

2299:   if (!ts->ops->getsolutioncomponents) *n = 0;
2300:   else {
2301:     (*ts->ops->getsolutioncomponents)(ts,n,v);
2302:   }
2303:   return(0);
2304: }

2306: /*@
2307:    TSGetAuxSolution - Returns an auxiliary solution at the present
2308:    timestep, if available for the time integration method being used.

2310:    Not Collective, but Vec returned is parallel if TS is parallel

2312:    Parameters :
2313: .  ts - the TS context obtained from TSCreate() (input parameter).
2314: .  v - the vector containing the auxiliary solution

2316:    Level: intermediate

2318: .seealso: TSGetSolution()

2320: .keywords: TS, timestep, get, solution
2321: @*/
2322: PetscErrorCode  TSGetAuxSolution(TS ts,Vec *v)
2323: {

2328:   if (ts->ops->getauxsolution) {
2329:     (*ts->ops->getauxsolution)(ts,v);
2330:   } else {
2331:     VecZeroEntries(*v);
2332:   }
2333:   return(0);
2334: }

2336: /*@
2337:    TSGetTimeError - Returns the estimated error vector, if the chosen
2338:    TSType has an error estimation functionality.

2340:    Not Collective, but Vec returned is parallel if TS is parallel

2342:    Note: MUST call after TSSetUp()

2344:    Parameters :
2345: .  ts - the TS context obtained from TSCreate() (input parameter).
2346: .  n - current estimate (n=0) or previous one (n=-1)
2347: .  v - the vector containing the error (same size as the solution).

2349:    Level: intermediate

2351: .seealso: TSGetSolution(), TSSetTimeError()

2353: .keywords: TS, timestep, get, error
2354: @*/
2355: PetscErrorCode  TSGetTimeError(TS ts,PetscInt n,Vec *v)
2356: {

2361:   if (ts->ops->gettimeerror) {
2362:     (*ts->ops->gettimeerror)(ts,n,v);
2363:   } else {
2364:     VecZeroEntries(*v);
2365:   }
2366:   return(0);
2367: }

2369: /*@
2370:    TSSetTimeError - Sets the estimated error vector, if the chosen
2371:    TSType has an error estimation functionality. This can be used
2372:    to restart such a time integrator with a given error vector.

2374:    Not Collective, but Vec returned is parallel if TS is parallel

2376:    Parameters :
2377: .  ts - the TS context obtained from TSCreate() (input parameter).
2378: .  v - the vector containing the error (same size as the solution).

2380:    Level: intermediate

2382: .seealso: TSSetSolution(), TSGetTimeError)

2384: .keywords: TS, timestep, get, error
2385: @*/
2386: PetscErrorCode  TSSetTimeError(TS ts,Vec v)
2387: {

2392:   if (!ts->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetUp() first");
2393:   if (ts->ops->settimeerror) {
2394:     (*ts->ops->settimeerror)(ts,v);
2395:   }
2396:   return(0);
2397: }

2399: /* ----- Routines to initialize and destroy a timestepper ---- */
2400: /*@
2401:   TSSetProblemType - Sets the type of problem to be solved.

2403:   Not collective

2405:   Input Parameters:
2406: + ts   - The TS
2407: - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
2408: .vb
2409:          U_t - A U = 0      (linear)
2410:          U_t - A(t) U = 0   (linear)
2411:          F(t,U,U_t) = 0     (nonlinear)
2412: .ve

2414:    Level: beginner

2416: .keywords: TS, problem type
2417: .seealso: TSSetUp(), TSProblemType, TS
2418: @*/
2419: PetscErrorCode  TSSetProblemType(TS ts, TSProblemType type)
2420: {

2425:   ts->problem_type = type;
2426:   if (type == TS_LINEAR) {
2427:     SNES snes;
2428:     TSGetSNES(ts,&snes);
2429:     SNESSetType(snes,SNESKSPONLY);
2430:   }
2431:   return(0);
2432: }

2434: /*@C
2435:   TSGetProblemType - Gets the type of problem to be solved.

2437:   Not collective

2439:   Input Parameter:
2440: . ts   - The TS

2442:   Output Parameter:
2443: . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
2444: .vb
2445:          M U_t = A U
2446:          M(t) U_t = A(t) U
2447:          F(t,U,U_t)
2448: .ve

2450:    Level: beginner

2452: .keywords: TS, problem type
2453: .seealso: TSSetUp(), TSProblemType, TS
2454: @*/
2455: PetscErrorCode  TSGetProblemType(TS ts, TSProblemType *type)
2456: {
2460:   *type = ts->problem_type;
2461:   return(0);
2462: }

2464: /*@
2465:    TSSetUp - Sets up the internal data structures for the later use
2466:    of a timestepper.

2468:    Collective on TS

2470:    Input Parameter:
2471: .  ts - the TS context obtained from TSCreate()

2473:    Notes:
2474:    For basic use of the TS solvers the user need not explicitly call
2475:    TSSetUp(), since these actions will automatically occur during
2476:    the call to TSStep() or TSSolve().  However, if one wishes to control this
2477:    phase separately, TSSetUp() should be called after TSCreate()
2478:    and optional routines of the form TSSetXXX(), but before TSStep() and TSSolve().

2480:    Level: advanced

2482: .keywords: TS, timestep, setup

2484: .seealso: TSCreate(), TSStep(), TSDestroy(), TSSolve()
2485: @*/
2486: PetscErrorCode  TSSetUp(TS ts)
2487: {
2489:   DM             dm;
2490:   PetscErrorCode (*func)(SNES,Vec,Vec,void*);
2491:   PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*);
2492:   TSIFunction    ifun;
2493:   TSIJacobian    ijac;
2494:   TSI2Jacobian   i2jac;
2495:   TSRHSJacobian  rhsjac;
2496:   PetscBool      isnone;

2500:   if (ts->setupcalled) return(0);

2502:   if (!((PetscObject)ts)->type_name) {
2503:     TSGetIFunction(ts,NULL,&ifun,NULL);
2504:     TSSetType(ts,ifun ? TSBEULER : TSEULER);
2505:   }

2507:   if (!ts->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first");

2509:   if (ts->rhsjacobian.reuse) {
2510:     Mat Amat,Pmat;
2511:     SNES snes;
2512:     TSGetSNES(ts,&snes);
2513:     SNESGetJacobian(snes,&Amat,&Pmat,NULL,NULL);
2514:     /* Matching matrices implies that an IJacobian is NOT set, because if it had been set, the IJacobian's matrix would
2515:      * have displaced the RHS matrix */
2516:     if (Amat == ts->Arhs) {
2517:       /* 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 */
2518:       MatDuplicate(ts->Arhs,MAT_COPY_VALUES,&Amat);
2519:       SNESSetJacobian(snes,Amat,NULL,NULL,NULL);
2520:       MatDestroy(&Amat);
2521:     }
2522:     if (Pmat == ts->Brhs) {
2523:       MatDuplicate(ts->Brhs,MAT_COPY_VALUES,&Pmat);
2524:       SNESSetJacobian(snes,NULL,Pmat,NULL,NULL);
2525:       MatDestroy(&Pmat);
2526:     }
2527:   }

2529:   TSGetAdapt(ts,&ts->adapt);
2530:   TSAdaptSetDefaultType(ts->adapt,ts->default_adapt_type);

2532:   if (ts->ops->setup) {
2533:     (*ts->ops->setup)(ts);
2534:   }

2536:   /* Attempt to check/preset a default value for the exact final time option */
2537:   PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&isnone);
2538:   if (!isnone && ts->exact_final_time == TS_EXACTFINALTIME_UNSPECIFIED)
2539:     ts->exact_final_time = TS_EXACTFINALTIME_MATCHSTEP;

2541:   /* In the case where we've set a DMTSFunction or what have you, we need the default SNESFunction
2542:      to be set right but can't do it elsewhere due to the overreliance on ctx=ts.
2543:    */
2544:   TSGetDM(ts,&dm);
2545:   DMSNESGetFunction(dm,&func,NULL);
2546:   if (!func) {
2547:     DMSNESSetFunction(dm,SNESTSFormFunction,ts);
2548:   }
2549:   /* If the SNES doesn't have a jacobian set and the TS has an ijacobian or rhsjacobian set, set the SNES to use it.
2550:      Otherwise, the SNES will use coloring internally to form the Jacobian.
2551:    */
2552:   DMSNESGetJacobian(dm,&jac,NULL);
2553:   DMTSGetIJacobian(dm,&ijac,NULL);
2554:   DMTSGetI2Jacobian(dm,&i2jac,NULL);
2555:   DMTSGetRHSJacobian(dm,&rhsjac,NULL);
2556:   if (!jac && (ijac || i2jac || rhsjac)) {
2557:     DMSNESSetJacobian(dm,SNESTSFormJacobian,ts);
2558:   }

2560:   /* if time integration scheme has a starting method, call it */
2561:   if (ts->ops->startingmethod) {
2562:     (*ts->ops->startingmethod)(ts);
2563:   }

2565:   ts->setupcalled = PETSC_TRUE;
2566:   return(0);
2567: }

2569: /*@
2570:    TSReset - Resets a TS context and removes any allocated Vecs and Mats.

2572:    Collective on TS

2574:    Input Parameter:
2575: .  ts - the TS context obtained from TSCreate()

2577:    Level: beginner

2579: .keywords: TS, timestep, reset

2581: .seealso: TSCreate(), TSSetup(), TSDestroy()
2582: @*/
2583: PetscErrorCode  TSReset(TS ts)
2584: {


2590:   if (ts->ops->reset) {
2591:     (*ts->ops->reset)(ts);
2592:   }
2593:   if (ts->snes) {SNESReset(ts->snes);}
2594:   if (ts->adapt) {TSAdaptReset(ts->adapt);}

2596:   MatDestroy(&ts->Arhs);
2597:   MatDestroy(&ts->Brhs);
2598:   VecDestroy(&ts->Frhs);
2599:   VecDestroy(&ts->vec_sol);
2600:   VecDestroy(&ts->vec_dot);
2601:   VecDestroy(&ts->vatol);
2602:   VecDestroy(&ts->vrtol);
2603:   VecDestroyVecs(ts->nwork,&ts->work);

2605:   VecDestroyVecs(ts->numcost,&ts->vecs_drdy);
2606:   VecDestroyVecs(ts->numcost,&ts->vecs_drdp);

2608:   MatDestroy(&ts->Jacp);
2609:   VecDestroy(&ts->vec_costintegral);
2610:   VecDestroy(&ts->vec_costintegrand);
2611:   MatDestroy(&ts->mat_sensip);

2613:   ts->setupcalled = PETSC_FALSE;
2614:   return(0);
2615: }

2617: /*@
2618:    TSDestroy - Destroys the timestepper context that was created
2619:    with TSCreate().

2621:    Collective on TS

2623:    Input Parameter:
2624: .  ts - the TS context obtained from TSCreate()

2626:    Level: beginner

2628: .keywords: TS, timestepper, destroy

2630: .seealso: TSCreate(), TSSetUp(), TSSolve()
2631: @*/
2632: PetscErrorCode  TSDestroy(TS *ts)
2633: {

2637:   if (!*ts) return(0);
2639:   if (--((PetscObject)(*ts))->refct > 0) {*ts = 0; return(0);}

2641:   TSReset((*ts));

2643:   /* if memory was published with SAWs then destroy it */
2644:   PetscObjectSAWsViewOff((PetscObject)*ts);
2645:   if ((*ts)->ops->destroy) {(*(*ts)->ops->destroy)((*ts));}

2647:   TSTrajectoryDestroy(&(*ts)->trajectory);

2649:   TSAdaptDestroy(&(*ts)->adapt);
2650:   TSEventDestroy(&(*ts)->event);

2652:   SNESDestroy(&(*ts)->snes);
2653:   DMDestroy(&(*ts)->dm);
2654:   TSMonitorCancel((*ts));
2655:   TSAdjointMonitorCancel((*ts));

2657:   PetscHeaderDestroy(ts);
2658:   return(0);
2659: }

2661: /*@
2662:    TSGetSNES - Returns the SNES (nonlinear solver) associated with
2663:    a TS (timestepper) context. Valid only for nonlinear problems.

2665:    Not Collective, but SNES is parallel if TS is parallel

2667:    Input Parameter:
2668: .  ts - the TS context obtained from TSCreate()

2670:    Output Parameter:
2671: .  snes - the nonlinear solver context

2673:    Notes:
2674:    The user can then directly manipulate the SNES context to set various
2675:    options, etc.  Likewise, the user can then extract and manipulate the
2676:    KSP, KSP, and PC contexts as well.

2678:    TSGetSNES() does not work for integrators that do not use SNES; in
2679:    this case TSGetSNES() returns NULL in snes.

2681:    Level: beginner

2683: .keywords: timestep, get, SNES
2684: @*/
2685: PetscErrorCode  TSGetSNES(TS ts,SNES *snes)
2686: {

2692:   if (!ts->snes) {
2693:     SNESCreate(PetscObjectComm((PetscObject)ts),&ts->snes);
2694:     SNESSetFunction(ts->snes,NULL,SNESTSFormFunction,ts);
2695:     PetscLogObjectParent((PetscObject)ts,(PetscObject)ts->snes);
2696:     PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);
2697:     if (ts->dm) {SNESSetDM(ts->snes,ts->dm);}
2698:     if (ts->problem_type == TS_LINEAR) {
2699:       SNESSetType(ts->snes,SNESKSPONLY);
2700:     }
2701:   }
2702:   *snes = ts->snes;
2703:   return(0);
2704: }

2706: /*@
2707:    TSSetSNES - Set the SNES (nonlinear solver) to be used by the timestepping context

2709:    Collective

2711:    Input Parameter:
2712: +  ts - the TS context obtained from TSCreate()
2713: -  snes - the nonlinear solver context

2715:    Notes:
2716:    Most users should have the TS created by calling TSGetSNES()

2718:    Level: developer

2720: .keywords: timestep, set, SNES
2721: @*/
2722: PetscErrorCode TSSetSNES(TS ts,SNES snes)
2723: {
2725:   PetscErrorCode (*func)(SNES,Vec,Mat,Mat,void*);

2730:   PetscObjectReference((PetscObject)snes);
2731:   SNESDestroy(&ts->snes);

2733:   ts->snes = snes;

2735:   SNESSetFunction(ts->snes,NULL,SNESTSFormFunction,ts);
2736:   SNESGetJacobian(ts->snes,NULL,NULL,&func,NULL);
2737:   if (func == SNESTSFormJacobian) {
2738:     SNESSetJacobian(ts->snes,NULL,NULL,SNESTSFormJacobian,ts);
2739:   }
2740:   return(0);
2741: }

2743: /*@
2744:    TSGetKSP - Returns the KSP (linear solver) associated with
2745:    a TS (timestepper) context.

2747:    Not Collective, but KSP is parallel if TS is parallel

2749:    Input Parameter:
2750: .  ts - the TS context obtained from TSCreate()

2752:    Output Parameter:
2753: .  ksp - the nonlinear solver context

2755:    Notes:
2756:    The user can then directly manipulate the KSP context to set various
2757:    options, etc.  Likewise, the user can then extract and manipulate the
2758:    KSP and PC contexts as well.

2760:    TSGetKSP() does not work for integrators that do not use KSP;
2761:    in this case TSGetKSP() returns NULL in ksp.

2763:    Level: beginner

2765: .keywords: timestep, get, KSP
2766: @*/
2767: PetscErrorCode  TSGetKSP(TS ts,KSP *ksp)
2768: {
2770:   SNES           snes;

2775:   if (!((PetscObject)ts)->type_name) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"KSP is not created yet. Call TSSetType() first");
2776:   if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()");
2777:   TSGetSNES(ts,&snes);
2778:   SNESGetKSP(snes,ksp);
2779:   return(0);
2780: }

2782: /* ----------- Routines to set solver parameters ---------- */

2784: /*@
2785:    TSSetMaxSteps - Sets the maximum number of steps to use.

2787:    Logically Collective on TS

2789:    Input Parameters:
2790: +  ts - the TS context obtained from TSCreate()
2791: -  maxsteps - maximum number of steps to use

2793:    Options Database Keys:
2794: .  -ts_max_steps <maxsteps> - Sets maxsteps

2796:    Notes:
2797:    The default maximum number of steps is 5000

2799:    Level: intermediate

2801: .keywords: TS, timestep, set, maximum, steps

2803: .seealso: TSGetMaxSteps(), TSSetMaxTime(), TSSetExactFinalTime()
2804: @*/
2805: PetscErrorCode TSSetMaxSteps(TS ts,PetscInt maxsteps)
2806: {
2810:   if (maxsteps < 0 ) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of steps must be non-negative");
2811:   ts->max_steps = maxsteps;
2812:   return(0);
2813: }

2815: /*@
2816:    TSGetMaxSteps - Gets the maximum number of steps to use.

2818:    Not Collective

2820:    Input Parameters:
2821: .  ts - the TS context obtained from TSCreate()

2823:    Output Parameter:
2824: .  maxsteps - maximum number of steps to use

2826:    Level: advanced

2828: .keywords: TS, timestep, get, maximum, steps

2830: .seealso: TSSetMaxSteps(), TSGetMaxTime(), TSSetMaxTime()
2831: @*/
2832: PetscErrorCode TSGetMaxSteps(TS ts,PetscInt *maxsteps)
2833: {
2837:   *maxsteps = ts->max_steps;
2838:   return(0);
2839: }

2841: /*@
2842:    TSSetMaxTime - Sets the maximum (or final) time for timestepping.

2844:    Logically Collective on TS

2846:    Input Parameters:
2847: +  ts - the TS context obtained from TSCreate()
2848: -  maxtime - final time to step to

2850:    Options Database Keys:
2851: .  -ts_max_time <maxtime> - Sets maxtime

2853:    Notes:
2854:    The default maximum time is 5.0

2856:    Level: intermediate

2858: .keywords: TS, timestep, set, maximum, time

2860: .seealso: TSGetMaxTime(), TSSetMaxSteps(), TSSetExactFinalTime()
2861: @*/
2862: PetscErrorCode TSSetMaxTime(TS ts,PetscReal maxtime)
2863: {
2867:   ts->max_time = maxtime;
2868:   return(0);
2869: }

2871: /*@
2872:    TSGetMaxTime - Gets the maximum (or final) time for timestepping.

2874:    Not Collective

2876:    Input Parameters:
2877: .  ts - the TS context obtained from TSCreate()

2879:    Output Parameter:
2880: .  maxtime - final time to step to

2882:    Level: advanced

2884: .keywords: TS, timestep, get, maximum, time

2886: .seealso: TSSetMaxTime(), TSGetMaxSteps(), TSSetMaxSteps()
2887: @*/
2888: PetscErrorCode TSGetMaxTime(TS ts,PetscReal *maxtime)
2889: {
2893:   *maxtime = ts->max_time;
2894:   return(0);
2895: }

2897: /*@
2898:    TSSetInitialTimeStep - Deprecated, use TSSetTime() and TSSetTimeStep().

2900:    Level: deprecated

2902: @*/
2903: PetscErrorCode  TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step)
2904: {
2908:   TSSetTime(ts,initial_time);
2909:   TSSetTimeStep(ts,time_step);
2910:   return(0);
2911: }

2913: /*@
2914:    TSGetDuration - Deprecated, use TSGetMaxSteps() and TSGetMaxTime().

2916:    Level: deprecated

2918: @*/
2919: PetscErrorCode TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime)
2920: {
2923:   if (maxsteps) {
2925:     *maxsteps = ts->max_steps;
2926:   }
2927:   if (maxtime) {
2929:     *maxtime = ts->max_time;
2930:   }
2931:   return(0);
2932: }

2934: /*@
2935:    TSSetDuration - Deprecated, use TSSetMaxSteps() and TSSetMaxTime().

2937:    Level: deprecated

2939: @*/
2940: PetscErrorCode TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime)
2941: {
2946:   if (maxsteps >= 0) ts->max_steps = maxsteps;
2947:   if (maxtime != PETSC_DEFAULT) ts->max_time = maxtime;
2948:   return(0);
2949: }

2951: /*@
2952:    TSGetTimeStepNumber - Deprecated, use TSGetStepNumber().

2954:    Level: deprecated

2956: @*/
2957: PetscErrorCode TSGetTimeStepNumber(TS ts,PetscInt *steps) { return TSGetStepNumber(ts,steps); }

2959: /*@
2960:    TSGetTotalSteps - Deprecated, use TSGetStepNumber().

2962:    Level: deprecated

2964: @*/
2965: PetscErrorCode TSGetTotalSteps(TS ts,PetscInt *steps) { return TSGetStepNumber(ts,steps); }

2967: /*@
2968:    TSSetSolution - Sets the initial solution vector
2969:    for use by the TS routines.

2971:    Logically Collective on TS and Vec

2973:    Input Parameters:
2974: +  ts - the TS context obtained from TSCreate()
2975: -  u - the solution vector

2977:    Level: beginner

2979: .keywords: TS, timestep, set, solution, initial values

2981: .seealso: TSSetSolutionFunction(), TSGetSolution(), TSCreate()
2982: @*/
2983: PetscErrorCode  TSSetSolution(TS ts,Vec u)
2984: {
2986:   DM             dm;

2991:   PetscObjectReference((PetscObject)u);
2992:   VecDestroy(&ts->vec_sol);
2993:   ts->vec_sol = u;

2995:   TSGetDM(ts,&dm);
2996:   DMShellSetGlobalVector(dm,u);
2997:   return(0);
2998: }

3000: /*@C
3001:   TSSetPreStep - Sets the general-purpose function
3002:   called once at the beginning of each time step.

3004:   Logically Collective on TS

3006:   Input Parameters:
3007: + ts   - The TS context obtained from TSCreate()
3008: - func - The function

3010:   Calling sequence of func:
3011: . func (TS ts);

3013:   Level: intermediate

3015: .keywords: TS, timestep
3016: .seealso: TSSetPreStage(), TSSetPostStage(), TSSetPostStep(), TSStep(), TSRestartStep()
3017: @*/
3018: PetscErrorCode  TSSetPreStep(TS ts, PetscErrorCode (*func)(TS))
3019: {
3022:   ts->prestep = func;
3023:   return(0);
3024: }

3026: /*@
3027:   TSPreStep - Runs the user-defined pre-step function.

3029:   Collective on TS

3031:   Input Parameters:
3032: . ts   - The TS context obtained from TSCreate()

3034:   Notes:
3035:   TSPreStep() is typically used within time stepping implementations,
3036:   so most users would not generally call this routine themselves.

3038:   Level: developer

3040: .keywords: TS, timestep
3041: .seealso: TSSetPreStep(), TSPreStage(), TSPostStage(), TSPostStep()
3042: @*/
3043: PetscErrorCode  TSPreStep(TS ts)
3044: {

3049:   if (ts->prestep) {
3050:     Vec              U;
3051:     PetscObjectState sprev,spost;

3053:     TSGetSolution(ts,&U);
3054:     PetscObjectStateGet((PetscObject)U,&sprev);
3055:     PetscStackCallStandard((*ts->prestep),(ts));
3056:     PetscObjectStateGet((PetscObject)U,&spost);
3057:     if (sprev != spost) {TSRestartStep(ts);}
3058:   }
3059:   return(0);
3060: }

3062: /*@C
3063:   TSSetPreStage - Sets the general-purpose function
3064:   called once at the beginning of each stage.

3066:   Logically Collective on TS

3068:   Input Parameters:
3069: + ts   - The TS context obtained from TSCreate()
3070: - func - The function

3072:   Calling sequence of func:
3073: . PetscErrorCode func(TS ts, PetscReal stagetime);

3075:   Level: intermediate

3077:   Note:
3078:   There may be several stages per time step. If the solve for a given stage fails, the step may be rejected and retried.
3079:   The time step number being computed can be queried using TSGetStepNumber() and the total size of the step being
3080:   attempted can be obtained using TSGetTimeStep(). The time at the start of the step is available via TSGetTime().

3082: .keywords: TS, timestep
3083: .seealso: TSSetPostStage(), TSSetPreStep(), TSSetPostStep(), TSGetApplicationContext()
3084: @*/
3085: PetscErrorCode  TSSetPreStage(TS ts, PetscErrorCode (*func)(TS,PetscReal))
3086: {
3089:   ts->prestage = func;
3090:   return(0);
3091: }

3093: /*@C
3094:   TSSetPostStage - Sets the general-purpose function
3095:   called once at the end of each stage.

3097:   Logically Collective on TS

3099:   Input Parameters:
3100: + ts   - The TS context obtained from TSCreate()
3101: - func - The function

3103:   Calling sequence of func:
3104: . PetscErrorCode func(TS ts, PetscReal stagetime, PetscInt stageindex, Vec* Y);

3106:   Level: intermediate

3108:   Note:
3109:   There may be several stages per time step. If the solve for a given stage fails, the step may be rejected and retried.
3110:   The time step number being computed can be queried using TSGetStepNumber() and the total size of the step being
3111:   attempted can be obtained using TSGetTimeStep(). The time at the start of the step is available via TSGetTime().

3113: .keywords: TS, timestep
3114: .seealso: TSSetPreStage(), TSSetPreStep(), TSSetPostStep(), TSGetApplicationContext()
3115: @*/
3116: PetscErrorCode  TSSetPostStage(TS ts, PetscErrorCode (*func)(TS,PetscReal,PetscInt,Vec*))
3117: {
3120:   ts->poststage = func;
3121:   return(0);
3122: }

3124: /*@C
3125:   TSSetPostEvaluate - Sets the general-purpose function
3126:   called once at the end of each step evaluation.

3128:   Logically Collective on TS

3130:   Input Parameters:
3131: + ts   - The TS context obtained from TSCreate()
3132: - func - The function

3134:   Calling sequence of func:
3135: . PetscErrorCode func(TS ts);

3137:   Level: intermediate

3139:   Note:
3140:   Semantically, TSSetPostEvaluate() differs from TSSetPostStep() since the function it sets is called before event-handling
3141:   thus guaranteeing the same solution (computed by the time-stepper) will be passed to it. On the other hand, TSPostStep()
3142:   may be passed a different solution, possibly changed by the event handler. TSPostEvaluate() is called after the next step
3143:   solution is evaluated allowing to modify it, if need be. The solution can be obtained with TSGetSolution(), the time step
3144:   with TSGetTimeStep(), and the time at the start of the step is available via TSGetTime()

3146: .keywords: TS, timestep
3147: .seealso: TSSetPreStage(), TSSetPreStep(), TSSetPostStep(), TSGetApplicationContext()
3148: @*/
3149: PetscErrorCode  TSSetPostEvaluate(TS ts, PetscErrorCode (*func)(TS))
3150: {
3153:   ts->postevaluate = func;
3154:   return(0);
3155: }

3157: /*@
3158:   TSPreStage - Runs the user-defined pre-stage function set using TSSetPreStage()

3160:   Collective on TS

3162:   Input Parameters:
3163: . ts          - The TS context obtained from TSCreate()
3164:   stagetime   - The absolute time of the current stage

3166:   Notes:
3167:   TSPreStage() is typically used within time stepping implementations,
3168:   most users would not generally call this routine themselves.

3170:   Level: developer

3172: .keywords: TS, timestep
3173: .seealso: TSPostStage(), TSSetPreStep(), TSPreStep(), TSPostStep()
3174: @*/
3175: PetscErrorCode  TSPreStage(TS ts, PetscReal stagetime)
3176: {

3181:   if (ts->prestage) {
3182:     PetscStackCallStandard((*ts->prestage),(ts,stagetime));
3183:   }
3184:   return(0);
3185: }

3187: /*@
3188:   TSPostStage - Runs the user-defined post-stage function set using TSSetPostStage()

3190:   Collective on TS

3192:   Input Parameters:
3193: . ts          - The TS context obtained from TSCreate()
3194:   stagetime   - The absolute time of the current stage
3195:   stageindex  - Stage number
3196:   Y           - Array of vectors (of size = total number
3197:                 of stages) with the stage solutions

3199:   Notes:
3200:   TSPostStage() is typically used within time stepping implementations,
3201:   most users would not generally call this routine themselves.

3203:   Level: developer

3205: .keywords: TS, timestep
3206: .seealso: TSPreStage(), TSSetPreStep(), TSPreStep(), TSPostStep()
3207: @*/
3208: PetscErrorCode  TSPostStage(TS ts, PetscReal stagetime, PetscInt stageindex, Vec *Y)
3209: {

3214:   if (ts->poststage) {
3215:     PetscStackCallStandard((*ts->poststage),(ts,stagetime,stageindex,Y));
3216:   }
3217:   return(0);
3218: }

3220: /*@
3221:   TSPostEvaluate - Runs the user-defined post-evaluate function set using TSSetPostEvaluate()

3223:   Collective on TS

3225:   Input Parameters:
3226: . ts          - The TS context obtained from TSCreate()

3228:   Notes:
3229:   TSPostEvaluate() is typically used within time stepping implementations,
3230:   most users would not generally call this routine themselves.

3232:   Level: developer

3234: .keywords: TS, timestep
3235: .seealso: TSSetPostEvaluate(), TSSetPreStep(), TSPreStep(), TSPostStep()
3236: @*/
3237: PetscErrorCode  TSPostEvaluate(TS ts)
3238: {

3243:   if (ts->postevaluate) {
3244:     Vec              U;
3245:     PetscObjectState sprev,spost;

3247:     TSGetSolution(ts,&U);
3248:     PetscObjectStateGet((PetscObject)U,&sprev);
3249:     PetscStackCallStandard((*ts->postevaluate),(ts));
3250:     PetscObjectStateGet((PetscObject)U,&spost);
3251:     if (sprev != spost) {TSRestartStep(ts);}
3252:   }
3253:   return(0);
3254: }

3256: /*@C
3257:   TSSetPostStep - Sets the general-purpose function
3258:   called once at the end of each time step.

3260:   Logically Collective on TS

3262:   Input Parameters:
3263: + ts   - The TS context obtained from TSCreate()
3264: - func - The function

3266:   Calling sequence of func:
3267: $ func (TS ts);

3269:   Notes:
3270:   The function set by TSSetPostStep() is called after each successful step. The solution vector X
3271:   obtained by TSGetSolution() may be different than that computed at the step end if the event handler
3272:   locates an event and TSPostEvent() modifies it. Use TSSetPostEvaluate() if an unmodified solution is needed instead.

3274:   Level: intermediate

3276: .keywords: TS, timestep
3277: .seealso: TSSetPreStep(), TSSetPreStage(), TSSetPostEvaluate(), TSGetTimeStep(), TSGetStepNumber(), TSGetTime(), TSRestartStep()
3278: @*/
3279: PetscErrorCode  TSSetPostStep(TS ts, PetscErrorCode (*func)(TS))
3280: {
3283:   ts->poststep = func;
3284:   return(0);
3285: }

3287: /*@
3288:   TSPostStep - Runs the user-defined post-step function.

3290:   Collective on TS

3292:   Input Parameters:
3293: . ts   - The TS context obtained from TSCreate()

3295:   Notes:
3296:   TSPostStep() is typically used within time stepping implementations,
3297:   so most users would not generally call this routine themselves.

3299:   Level: developer

3301: .keywords: TS, timestep
3302: @*/
3303: PetscErrorCode  TSPostStep(TS ts)
3304: {

3309:   if (ts->poststep) {
3310:     Vec              U;
3311:     PetscObjectState sprev,spost;

3313:     TSGetSolution(ts,&U);
3314:     PetscObjectStateGet((PetscObject)U,&sprev);
3315:     PetscStackCallStandard((*ts->poststep),(ts));
3316:     PetscObjectStateGet((PetscObject)U,&spost);
3317:     if (sprev != spost) {TSRestartStep(ts);}
3318:   }
3319:   return(0);
3320: }

3322: /* ------------ Routines to set performance monitoring options ----------- */

3324: /*@C
3325:    TSMonitorSet - Sets an ADDITIONAL function that is to be used at every
3326:    timestep to display the iteration's  progress.

3328:    Logically Collective on TS

3330:    Input Parameters:
3331: +  ts - the TS context obtained from TSCreate()
3332: .  monitor - monitoring routine
3333: .  mctx - [optional] user-defined context for private data for the
3334:              monitor routine (use NULL if no context is desired)
3335: -  monitordestroy - [optional] routine that frees monitor context
3336:           (may be NULL)

3338:    Calling sequence of monitor:
3339: $    PetscErrorCode monitor(TS ts,PetscInt steps,PetscReal time,Vec u,void *mctx)

3341: +    ts - the TS context
3342: .    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)
3343: .    time - current time
3344: .    u - current iterate
3345: -    mctx - [optional] monitoring context

3347:    Notes:
3348:    This routine adds an additional monitor to the list of monitors that
3349:    already has been loaded.

3351:    Fortran notes: Only a single monitor function can be set for each TS object

3353:    Level: intermediate

3355: .keywords: TS, timestep, set, monitor

3357: .seealso: TSMonitorDefault(), TSMonitorCancel()
3358: @*/
3359: PetscErrorCode  TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void**))
3360: {
3362:   PetscInt       i;
3363:   PetscBool      identical;

3367:   for (i=0; i<ts->numbermonitors;i++) {
3368:     PetscMonitorCompare((PetscErrorCode (*)(void))monitor,mctx,mdestroy,(PetscErrorCode (*)(void))ts->monitor[i],ts->monitorcontext[i],ts->monitordestroy[i],&identical);
3369:     if (identical) return(0);
3370:   }
3371:   if (ts->numbermonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
3372:   ts->monitor[ts->numbermonitors]          = monitor;
3373:   ts->monitordestroy[ts->numbermonitors]   = mdestroy;
3374:   ts->monitorcontext[ts->numbermonitors++] = (void*)mctx;
3375:   return(0);
3376: }

3378: /*@C
3379:    TSMonitorCancel - Clears all the monitors that have been set on a time-step object.

3381:    Logically Collective on TS

3383:    Input Parameters:
3384: .  ts - the TS context obtained from TSCreate()

3386:    Notes:
3387:    There is no way to remove a single, specific monitor.

3389:    Level: intermediate

3391: .keywords: TS, timestep, set, monitor

3393: .seealso: TSMonitorDefault(), TSMonitorSet()
3394: @*/
3395: PetscErrorCode  TSMonitorCancel(TS ts)
3396: {
3398:   PetscInt       i;

3402:   for (i=0; i<ts->numbermonitors; i++) {
3403:     if (ts->monitordestroy[i]) {
3404:       (*ts->monitordestroy[i])(&ts->monitorcontext[i]);
3405:     }
3406:   }
3407:   ts->numbermonitors = 0;
3408:   return(0);
3409: }

3411: /*@C
3412:    TSMonitorDefault - The Default monitor, prints the timestep and time for each step

3414:    Level: intermediate

3416: .keywords: TS, set, monitor

3418: .seealso:  TSMonitorSet()
3419: @*/
3420: PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscViewerAndFormat *vf)
3421: {
3423:   PetscViewer    viewer =  vf->viewer;
3424:   PetscBool      iascii,ibinary;

3428:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
3429:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&ibinary);
3430:   PetscViewerPushFormat(viewer,vf->format);
3431:   if (iascii) {
3432:     PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);
3433:     if (step == -1){ /* this indicates it is an interpolated solution */
3434:       PetscViewerASCIIPrintf(viewer,"Interpolated solution at time %g between steps %D and %D\n",(double)ptime,ts->steps-1,ts->steps);
3435:     } else {
3436:       PetscViewerASCIIPrintf(viewer,"%D TS dt %g time %g%s",step,(double)ts->time_step,(double)ptime,ts->steprollback ? " (r)\n" : "\n");
3437:     }
3438:     PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);
3439:   } else if (ibinary) {
3440:     PetscMPIInt rank;
3441:     MPI_Comm_rank(PetscObjectComm((PetscObject)viewer),&rank);
3442:     if (!rank) {
3443:       PetscBool skipHeader;
3444:       PetscInt  classid = REAL_FILE_CLASSID;

3446:       PetscViewerBinaryGetSkipHeader(viewer,&skipHeader);
3447:       if (!skipHeader) {
3448:          PetscViewerBinaryWrite(viewer,&classid,1,PETSC_INT,PETSC_FALSE);
3449:        }
3450:       PetscRealView(1,&ptime,viewer);
3451:     } else {
3452:       PetscRealView(0,&ptime,viewer);
3453:     }
3454:   }
3455:   PetscViewerPopFormat(viewer);
3456:   return(0);
3457: }

3459: /*@
3460:    TSInterpolate - Interpolate the solution computed during the previous step to an arbitrary location in the interval

3462:    Collective on TS

3464:    Input Argument:
3465: +  ts - time stepping context
3466: -  t - time to interpolate to

3468:    Output Argument:
3469: .  U - state at given time

3471:    Level: intermediate

3473:    Developer Notes:
3474:    TSInterpolate() and the storing of previous steps/stages should be generalized to support delay differential equations and continuous adjoints.

3476: .keywords: TS, set

3478: .seealso: TSSetExactFinalTime(), TSSolve()
3479: @*/
3480: PetscErrorCode TSInterpolate(TS ts,PetscReal t,Vec U)
3481: {

3487:   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);
3488:   if (!ts->ops->interpolate) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"%s does not provide interpolation",((PetscObject)ts)->type_name);
3489:   (*ts->ops->interpolate)(ts,t,U);
3490:   return(0);
3491: }

3493: /*@
3494:    TSStep - Steps one time step

3496:    Collective on TS

3498:    Input Parameter:
3499: .  ts - the TS context obtained from TSCreate()

3501:    Level: developer

3503:    Notes:
3504:    The public interface for the ODE/DAE solvers is TSSolve(), you should almost for sure be using that routine and not this routine.

3506:    The hook set using TSSetPreStep() is called before each attempt to take the step. In general, the time step size may
3507:    be changed due to adaptive error controller or solve failures. Note that steps may contain multiple stages.

3509:    This may over-step the final time provided in TSSetMaxTime() depending on the time-step used. TSSolve() interpolates to exactly the
3510:    time provided in TSSetMaxTime(). One can use TSInterpolate() to determine an interpolated solution within the final timestep.

3512: .keywords: TS, timestep, solve

3514: .seealso: TSCreate(), TSSetUp(), TSDestroy(), TSSolve(), TSSetPreStep(), TSSetPreStage(), TSSetPostStage(), TSInterpolate()
3515: @*/
3516: PetscErrorCode  TSStep(TS ts)
3517: {
3518:   PetscErrorCode   ierr;
3519:   static PetscBool cite = PETSC_FALSE;
3520:   PetscReal        ptime;

3524:   PetscCitationsRegister("@techreport{tspaper,\n"
3525:                                 "  title       = {{PETSc/TS}: A Modern Scalable {DAE/ODE} Solver Library},\n"
3526:                                 "  author      = {Shrirang Abhyankar and Jed Brown and Emil Constantinescu and Debojyoti Ghosh and Barry F. Smith},\n"
3527:                                 "  type        = {Preprint},\n"
3528:                                 "  number      = {ANL/MCS-P5061-0114},\n"
3529:                                 "  institution = {Argonne National Laboratory},\n"
3530:                                 "  year        = {2014}\n}\n",&cite);

3532:   TSSetUp(ts);
3533:   TSTrajectorySetUp(ts->trajectory,ts);

3535:   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>");
3536:   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()");
3537:   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");

3539:   if (!ts->steps) ts->ptime_prev = ts->ptime;
3540:   ptime = ts->ptime; ts->ptime_prev_rollback = ts->ptime_prev;
3541:   ts->reason = TS_CONVERGED_ITERATING;
3542:   if (!ts->ops->step) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSStep not implemented for type '%s'",((PetscObject)ts)->type_name);
3543:   PetscLogEventBegin(TS_Step,ts,0,0,0);
3544:   (*ts->ops->step)(ts);
3545:   PetscLogEventEnd(TS_Step,ts,0,0,0);
3546:   ts->ptime_prev = ptime;
3547:   ts->steps++;
3548:   ts->steprollback = PETSC_FALSE;
3549:   ts->steprestart  = PETSC_FALSE;

3551:   if (ts->reason < 0) {
3552:     if (ts->errorifstepfailed) {
3553:       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]);
3554:       else SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s",TSConvergedReasons[ts->reason]);
3555:     }
3556:   } else if (!ts->reason) {
3557:     if (ts->steps >= ts->max_steps) ts->reason = TS_CONVERGED_ITS;
3558:     else if (ts->ptime >= ts->max_time) ts->reason = TS_CONVERGED_TIME;
3559:   }
3560:   return(0);
3561: }

3563: /*@
3564:    TSEvaluateWLTE - Evaluate the weighted local truncation error norm
3565:    at the end of a time step with a given order of accuracy.

3567:    Collective on TS

3569:    Input Arguments:
3570: +  ts - time stepping context
3571: .  wnormtype - norm type, either NORM_2 or NORM_INFINITY
3572: -  order - optional, desired order for the error evaluation or PETSC_DECIDE

3574:    Output Arguments:
3575: +  order - optional, the actual order of the error evaluation
3576: -  wlte - the weighted local truncation error norm

3578:    Level: advanced

3580:    Notes:
3581:    If the timestepper cannot evaluate the error in a particular step
3582:    (eg. in the first step or restart steps after event handling),
3583:    this routine returns wlte=-1.0 .

3585: .seealso: TSStep(), TSAdapt, TSErrorWeightedNorm()
3586: @*/
3587: PetscErrorCode TSEvaluateWLTE(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte)
3588: {

3598:   if (wnormtype != NORM_2 && wnormtype != NORM_INFINITY) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"No support for norm type %s",NormTypes[wnormtype]);
3599:   if (!ts->ops->evaluatewlte) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSEvaluateWLTE not implemented for type '%s'",((PetscObject)ts)->type_name);
3600:   (*ts->ops->evaluatewlte)(ts,wnormtype,order,wlte);
3601:   return(0);
3602: }

3604: /*@
3605:    TSEvaluateStep - Evaluate the solution at the end of a time step with a given order of accuracy.

3607:    Collective on TS

3609:    Input Arguments:
3610: +  ts - time stepping context
3611: .  order - desired order of accuracy
3612: -  done - whether the step was evaluated at this order (pass NULL to generate an error if not available)

3614:    Output Arguments:
3615: .  U - state at the end of the current step

3617:    Level: advanced

3619:    Notes:
3620:    This function cannot be called until all stages have been evaluated.
3621:    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.

3623: .seealso: TSStep(), TSAdapt
3624: @*/
3625: PetscErrorCode TSEvaluateStep(TS ts,PetscInt order,Vec U,PetscBool *done)
3626: {

3633:   if (!ts->ops->evaluatestep) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSEvaluateStep not implemented for type '%s'",((PetscObject)ts)->type_name);
3634:   (*ts->ops->evaluatestep)(ts,order,U,done);
3635:   return(0);
3636: }

3638: /*@
3639:    TSSolve - Steps the requested number of timesteps.

3641:    Collective on TS

3643:    Input Parameter:
3644: +  ts - the TS context obtained from TSCreate()
3645: -  u - the solution vector  (can be null if TSSetSolution() was used and TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP) was not used,
3646:                              otherwise must contain the initial conditions and will contain the solution at the final requested time

3648:    Level: beginner

3650:    Notes:
3651:    The final time returned by this function may be different from the time of the internally
3652:    held state accessible by TSGetSolution() and TSGetTime() because the method may have
3653:    stepped over the final time.

3655: .keywords: TS, timestep, solve

3657: .seealso: TSCreate(), TSSetSolution(), TSStep(), TSGetTime(), TSGetSolveTime()
3658: @*/
3659: PetscErrorCode TSSolve(TS ts,Vec u)
3660: {
3661:   Vec               solution;
3662:   PetscErrorCode    ierr;


3668:   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 */
3669:     if (!ts->vec_sol || u == ts->vec_sol) {
3670:       VecDuplicate(u,&solution);
3671:       TSSetSolution(ts,solution);
3672:       VecDestroy(&solution); /* grant ownership */
3673:     }
3674:     VecCopy(u,ts->vec_sol);
3675:     if (ts->forward_solve) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Sensitivity analysis does not support the mode TS_EXACTFINALTIME_INTERPOLATE");
3676:   } else if (u) {
3677:     TSSetSolution(ts,u);
3678:   }
3679:   TSSetUp(ts);
3680:   TSTrajectorySetUp(ts->trajectory,ts);

3682:   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>");
3683:   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()");
3684:   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");

3686:   if (ts->forward_solve) {
3687:     TSForwardSetUp(ts);
3688:   }

3690:   /* reset number of steps only when the step is not restarted. ARKIMEX
3691:      restarts the step after an event. Resetting these counters in such case causes
3692:      TSTrajectory to incorrectly save the output files
3693:   */
3694:   /* reset time step and iteration counters */
3695:   if (!ts->steps) {
3696:     ts->ksp_its           = 0;
3697:     ts->snes_its          = 0;
3698:     ts->num_snes_failures = 0;
3699:     ts->reject            = 0;
3700:     ts->steprestart       = PETSC_TRUE;
3701:     ts->steprollback      = PETSC_FALSE;
3702:   }
3703:   if (ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP && ts->ptime + ts->time_step > ts->max_time) ts->time_step = ts->max_time - ts->ptime;
3704:   ts->reason = TS_CONVERGED_ITERATING;

3706:   TSViewFromOptions(ts,NULL,"-ts_view_pre");

3708:   if (ts->ops->solve) { /* This private interface is transitional and should be removed when all implementations are updated. */
3709:     (*ts->ops->solve)(ts);
3710:     if (u) {VecCopy(ts->vec_sol,u);}
3711:     ts->solvetime = ts->ptime;
3712:     solution = ts->vec_sol;
3713:   } else { /* Step the requested number of timesteps. */
3714:     if (ts->steps >= ts->max_steps) ts->reason = TS_CONVERGED_ITS;
3715:     else if (ts->ptime >= ts->max_time) ts->reason = TS_CONVERGED_TIME;

3717:     if (!ts->steps) {
3718:       TSTrajectorySet(ts->trajectory,ts,ts->steps,ts->ptime,ts->vec_sol);
3719:       TSEventInitialize(ts->event,ts,ts->ptime,ts->vec_sol);
3720:     }

3722:     while (!ts->reason) {
3723:       TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);
3724:       if (!ts->steprollback) {
3725:         TSPreStep(ts);
3726:       }
3727:       TSStep(ts);
3728:       if (ts->testjacobian) {
3729:         TSRHSJacobianTest(ts,NULL);
3730:       }
3731:       if (ts->testjacobiantranspose) {
3732:         TSRHSJacobianTestTranspose(ts,NULL);
3733:       }
3734:       if (ts->vec_costintegral && ts->costintegralfwd) { /* Must evaluate the cost integral before event is handled. The cost integral value can also be rolled back. */
3735:         TSForwardCostIntegral(ts);
3736:       }
3737:       if (ts->forward_solve) { /* compute forward sensitivities before event handling because postevent() may change RHS and jump conditions may have to be applied */
3738:         TSForwardStep(ts);
3739:       }
3740:       TSPostEvaluate(ts);
3741:       TSEventHandler(ts); /* The right-hand side may be changed due to event. Be careful with Any computation using the RHS information after this point. */
3742:       if (ts->steprollback) {
3743:         TSPostEvaluate(ts);
3744:       }
3745:       if (!ts->steprollback) {
3746:         TSTrajectorySet(ts->trajectory,ts,ts->steps,ts->ptime,ts->vec_sol);
3747:         TSPostStep(ts);
3748:       }
3749:     }
3750:     TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);

3752:     if (ts->exact_final_time == TS_EXACTFINALTIME_INTERPOLATE && ts->ptime > ts->max_time) {
3753:       TSInterpolate(ts,ts->max_time,u);
3754:       ts->solvetime = ts->max_time;
3755:       solution = u;
3756:       TSMonitor(ts,-1,ts->solvetime,solution);
3757:     } else {
3758:       if (u) {VecCopy(ts->vec_sol,u);}
3759:       ts->solvetime = ts->ptime;
3760:       solution = ts->vec_sol;
3761:     }
3762:   }

3764:   TSViewFromOptions(ts,NULL,"-ts_view");
3765:   VecViewFromOptions(solution,NULL,"-ts_view_solution");
3766:   PetscObjectSAWsBlock((PetscObject)ts);
3767:   if (ts->adjoint_solve) {
3768:     TSAdjointSolve(ts);
3769:   }
3770:   return(0);
3771: }

3773: /*@C
3774:    TSMonitor - Runs all user-provided monitor routines set using TSMonitorSet()

3776:    Collective on TS

3778:    Input Parameters:
3779: +  ts - time stepping context obtained from TSCreate()
3780: .  step - step number that has just completed
3781: .  ptime - model time of the state
3782: -  u - state at the current model time

3784:    Notes:
3785:    TSMonitor() is typically used automatically within the time stepping implementations.
3786:    Users would almost never call this routine directly.

3788:    A step of -1 indicates that the monitor is being called on a solution obtained by interpolating from computed solutions

3790:    Level: developer

3792: .keywords: TS, timestep
3793: @*/
3794: PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec u)
3795: {
3796:   DM             dm;
3797:   PetscInt       i,n = ts->numbermonitors;


3804:   TSGetDM(ts,&dm);
3805:   DMSetOutputSequenceNumber(dm,step,ptime);

3807:   VecLockPush(u);
3808:   for (i=0; i<n; i++) {
3809:     (*ts->monitor[i])(ts,step,ptime,u,ts->monitorcontext[i]);
3810:   }
3811:   VecLockPop(u);
3812:   return(0);
3813: }

3815: /* ------------------------------------------------------------------------*/
3816: /*@C
3817:    TSMonitorLGCtxCreate - Creates a TSMonitorLGCtx context for use with
3818:    TS to monitor the solution process graphically in various ways

3820:    Collective on TS

3822:    Input Parameters:
3823: +  host - the X display to open, or null for the local machine
3824: .  label - the title to put in the title bar
3825: .  x, y - the screen coordinates of the upper left coordinate of the window
3826: .  m, n - the screen width and height in pixels
3827: -  howoften - if positive then determines the frequency of the plotting, if -1 then only at the final time

3829:    Output Parameter:
3830: .  ctx - the context

3832:    Options Database Key:
3833: +  -ts_monitor_lg_timestep - automatically sets line graph monitor
3834: +  -ts_monitor_lg_timestep_log - automatically sets line graph monitor
3835: .  -ts_monitor_lg_solution - monitor the solution (or certain values of the solution by calling TSMonitorLGSetDisplayVariables() or TSMonitorLGCtxSetDisplayVariables())
3836: .  -ts_monitor_lg_error -  monitor the error
3837: .  -ts_monitor_lg_ksp_iterations - monitor the number of KSP iterations needed for each timestep
3838: .  -ts_monitor_lg_snes_iterations - monitor the number of SNES iterations needed for each timestep
3839: -  -lg_use_markers <true,false> - mark the data points (at each time step) on the plot; default is true

3841:    Notes:
3842:    Use TSMonitorLGCtxDestroy() to destroy.

3844:    One can provide a function that transforms the solution before plotting it with TSMonitorLGCtxSetTransform() or TSMonitorLGSetTransform()

3846:    Many of the functions that control the monitoring have two forms: TSMonitorLGSet/GetXXXX() and TSMonitorLGCtxSet/GetXXXX() the first take a TS object as the
3847:    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
3848:    as the first argument.

3850:    One can control the names displayed for each solution or error variable with TSMonitorLGCtxSetVariableNames() or TSMonitorLGSetVariableNames()

3852:    Level: intermediate

3854: .keywords: TS, monitor, line graph, residual

3856: .seealso: TSMonitorLGTimeStep(), TSMonitorSet(), TSMonitorLGSolution(), TSMonitorLGError(), TSMonitorDefault(), VecView(),
3857:            TSMonitorLGCtxCreate(), TSMonitorLGCtxSetVariableNames(), TSMonitorLGCtxGetVariableNames(),
3858:            TSMonitorLGSetVariableNames(), TSMonitorLGGetVariableNames(), TSMonitorLGSetDisplayVariables(), TSMonitorLGCtxSetDisplayVariables(),
3859:            TSMonitorLGCtxSetTransform(), TSMonitorLGSetTransform(), TSMonitorLGError(), TSMonitorLGSNESIterations(), TSMonitorLGKSPIterations(),
3860:            TSMonitorEnvelopeCtxCreate(), TSMonitorEnvelopeGetBounds(), TSMonitorEnvelopeCtxDestroy(), TSMonitorEnvelop()

3862: @*/
3863: PetscErrorCode  TSMonitorLGCtxCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscInt howoften,TSMonitorLGCtx *ctx)
3864: {
3865:   PetscDraw      draw;

3869:   PetscNew(ctx);
3870:   PetscDrawCreate(comm,host,label,x,y,m,n,&draw);
3871:   PetscDrawSetFromOptions(draw);
3872:   PetscDrawLGCreate(draw,1,&(*ctx)->lg);
3873:   PetscDrawLGSetFromOptions((*ctx)->lg);
3874:   PetscDrawDestroy(&draw);
3875:   (*ctx)->howoften = howoften;
3876:   return(0);
3877: }

3879: PetscErrorCode TSMonitorLGTimeStep(TS ts,PetscInt step,PetscReal ptime,Vec v,void *monctx)
3880: {
3881:   TSMonitorLGCtx ctx = (TSMonitorLGCtx) monctx;
3882:   PetscReal      x   = ptime,y;

3886:   if (step < 0) return(0); /* -1 indicates an interpolated solution */
3887:   if (!step) {
3888:     PetscDrawAxis axis;
3889:     const char *ylabel = ctx->semilogy ? "Log Time Step" : "Time Step";
3890:     PetscDrawLGGetAxis(ctx->lg,&axis);
3891:     PetscDrawAxisSetLabels(axis,"Timestep as function of time","Time",ylabel);
3892:     PetscDrawLGReset(ctx->lg);
3893:   }
3894:   TSGetTimeStep(ts,&y);
3895:   if (ctx->semilogy) y = PetscLog10Real(y);
3896:   PetscDrawLGAddPoint(ctx->lg,&x,&y);
3897:   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
3898:     PetscDrawLGDraw(ctx->lg);
3899:     PetscDrawLGSave(ctx->lg);
3900:   }
3901:   return(0);
3902: }

3904: /*@C
3905:    TSMonitorLGCtxDestroy - Destroys a line graph context that was created
3906:    with TSMonitorLGCtxCreate().

3908:    Collective on TSMonitorLGCtx

3910:    Input Parameter:
3911: .  ctx - the monitor context

3913:    Level: intermediate

3915: .keywords: TS, monitor, line graph, destroy

3917: .seealso: TSMonitorLGCtxCreate(),  TSMonitorSet(), TSMonitorLGTimeStep();
3918: @*/
3919: PetscErrorCode  TSMonitorLGCtxDestroy(TSMonitorLGCtx *ctx)
3920: {

3924:   if ((*ctx)->transformdestroy) {
3925:     ((*ctx)->transformdestroy)((*ctx)->transformctx);
3926:   }
3927:   PetscDrawLGDestroy(&(*ctx)->lg);
3928:   PetscStrArrayDestroy(&(*ctx)->names);
3929:   PetscStrArrayDestroy(&(*ctx)->displaynames);
3930:   PetscFree((*ctx)->displayvariables);
3931:   PetscFree((*ctx)->displayvalues);
3932:   PetscFree(*ctx);
3933:   return(0);
3934: }

3936: /*@
3937:    TSGetTime - Gets the time of the most recently completed step.

3939:    Not Collective

3941:    Input Parameter:
3942: .  ts - the TS context obtained from TSCreate()

3944:    Output Parameter:
3945: .  t  - the current time. This time may not corresponds to the final time set with TSSetMaxTime(), use TSGetSolveTime().

3947:    Level: beginner

3949:    Note:
3950:    When called during time step evaluation (e.g. during residual evaluation or via hooks set using TSSetPreStep(),
3951:    TSSetPreStage(), TSSetPostStage(), or TSSetPostStep()), the time is the time at the start of the step being evaluated.

3953: .seealso:  TSGetSolveTime(), TSSetTime(), TSGetTimeStep()

3955: .keywords: TS, get, time
3956: @*/
3957: PetscErrorCode  TSGetTime(TS ts,PetscReal *t)
3958: {
3962:   *t = ts->ptime;
3963:   return(0);
3964: }

3966: /*@
3967:    TSGetPrevTime - Gets the starting time of the previously completed step.

3969:    Not Collective

3971:    Input Parameter:
3972: .  ts - the TS context obtained from TSCreate()

3974:    Output Parameter:
3975: .  t  - the previous time

3977:    Level: beginner

3979: .seealso: TSGetTime(), TSGetSolveTime(), TSGetTimeStep()

3981: .keywords: TS, get, time
3982: @*/
3983: PetscErrorCode  TSGetPrevTime(TS ts,PetscReal *t)
3984: {
3988:   *t = ts->ptime_prev;
3989:   return(0);
3990: }

3992: /*@
3993:    TSSetTime - Allows one to reset the time.

3995:    Logically Collective on TS

3997:    Input Parameters:
3998: +  ts - the TS context obtained from TSCreate()
3999: -  time - the time

4001:    Level: intermediate

4003: .seealso: TSGetTime(), TSSetMaxSteps()

4005: .keywords: TS, set, time
4006: @*/
4007: PetscErrorCode  TSSetTime(TS ts, PetscReal t)
4008: {
4012:   ts->ptime = t;
4013:   return(0);
4014: }

4016: /*@C
4017:    TSSetOptionsPrefix - Sets the prefix used for searching for all
4018:    TS options in the database.

4020:    Logically Collective on TS

4022:    Input Parameter:
4023: +  ts     - The TS context
4024: -  prefix - The prefix to prepend to all option names

4026:    Notes:
4027:    A hyphen (-) must NOT be given at the beginning of the prefix name.
4028:    The first character of all runtime options is AUTOMATICALLY the
4029:    hyphen.

4031:    Level: advanced

4033: .keywords: TS, set, options, prefix, database

4035: .seealso: TSSetFromOptions()

4037: @*/
4038: PetscErrorCode  TSSetOptionsPrefix(TS ts,const char prefix[])
4039: {
4041:   SNES           snes;

4045:   PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);
4046:   TSGetSNES(ts,&snes);
4047:   SNESSetOptionsPrefix(snes,prefix);
4048:   return(0);
4049: }

4051: /*@C
4052:    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
4053:    TS options in the database.

4055:    Logically Collective on TS

4057:    Input Parameter:
4058: +  ts     - The TS context
4059: -  prefix - The prefix to prepend to all option names

4061:    Notes:
4062:    A hyphen (-) must NOT be given at the beginning of the prefix name.
4063:    The first character of all runtime options is AUTOMATICALLY the
4064:    hyphen.

4066:    Level: advanced

4068: .keywords: TS, append, options, prefix, database

4070: .seealso: TSGetOptionsPrefix()

4072: @*/
4073: PetscErrorCode  TSAppendOptionsPrefix(TS ts,const char prefix[])
4074: {
4076:   SNES           snes;

4080:   PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);
4081:   TSGetSNES(ts,&snes);
4082:   SNESAppendOptionsPrefix(snes,prefix);
4083:   return(0);
4084: }

4086: /*@C
4087:    TSGetOptionsPrefix - Sets the prefix used for searching for all
4088:    TS options in the database.

4090:    Not Collective

4092:    Input Parameter:
4093: .  ts - The TS context

4095:    Output Parameter:
4096: .  prefix - A pointer to the prefix string used

4098:    Notes: On the fortran side, the user should pass in a string 'prifix' of
4099:    sufficient length to hold the prefix.

4101:    Level: intermediate

4103: .keywords: TS, get, options, prefix, database

4105: .seealso: TSAppendOptionsPrefix()
4106: @*/
4107: PetscErrorCode  TSGetOptionsPrefix(TS ts,const char *prefix[])
4108: {

4114:   PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);
4115:   return(0);
4116: }

4118: /*@C
4119:    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.

4121:    Not Collective, but parallel objects are returned if TS is parallel

4123:    Input Parameter:
4124: .  ts  - The TS context obtained from TSCreate()

4126:    Output Parameters:
4127: +  Amat - The (approximate) Jacobian J of G, where U_t = G(U,t)  (or NULL)
4128: .  Pmat - The matrix from which the preconditioner is constructed, usually the same as Amat  (or NULL)
4129: .  func - Function to compute the Jacobian of the RHS  (or NULL)
4130: -  ctx - User-defined context for Jacobian evaluation routine  (or NULL)

4132:    Notes: You can pass in NULL for any return argument you do not need.

4134:    Level: intermediate

4136: .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetStepNumber()

4138: .keywords: TS, timestep, get, matrix, Jacobian
4139: @*/
4140: PetscErrorCode  TSGetRHSJacobian(TS ts,Mat *Amat,Mat *Pmat,TSRHSJacobian *func,void **ctx)
4141: {
4143:   DM             dm;

4146:   if (Amat || Pmat) {
4147:     SNES snes;
4148:     TSGetSNES(ts,&snes);
4149:     SNESSetUpMatrices(snes);
4150:     SNESGetJacobian(snes,Amat,Pmat,NULL,NULL);
4151:   }
4152:   TSGetDM(ts,&dm);
4153:   DMTSGetRHSJacobian(dm,func,ctx);
4154:   return(0);
4155: }

4157: /*@C
4158:    TSGetIJacobian - Returns the implicit Jacobian at the present timestep.

4160:    Not Collective, but parallel objects are returned if TS is parallel

4162:    Input Parameter:
4163: .  ts  - The TS context obtained from TSCreate()

4165:    Output Parameters:
4166: +  Amat  - The (approximate) Jacobian of F(t,U,U_t)
4167: .  Pmat - The matrix from which the preconditioner is constructed, often the same as Amat
4168: .  f   - The function to compute the matrices
4169: - ctx - User-defined context for Jacobian evaluation routine

4171:    Notes: You can pass in NULL for any return argument you do not need.

4173:    Level: advanced

4175: .seealso: TSGetTimeStep(), TSGetRHSJacobian(), TSGetMatrices(), TSGetTime(), TSGetStepNumber()

4177: .keywords: TS, timestep, get, matrix, Jacobian
4178: @*/
4179: PetscErrorCode  TSGetIJacobian(TS ts,Mat *Amat,Mat *Pmat,TSIJacobian *f,void **ctx)
4180: {
4182:   DM             dm;

4185:   if (Amat || Pmat) {
4186:     SNES snes;
4187:     TSGetSNES(ts,&snes);
4188:     SNESSetUpMatrices(snes);
4189:     SNESGetJacobian(snes,Amat,Pmat,NULL,NULL);
4190:   }
4191:   TSGetDM(ts,&dm);
4192:   DMTSGetIJacobian(dm,f,ctx);
4193:   return(0);
4194: }

4196: /*@C
4197:    TSMonitorDrawSolution - Monitors progress of the TS solvers by calling
4198:    VecView() for the solution at each timestep

4200:    Collective on TS

4202:    Input Parameters:
4203: +  ts - the TS context
4204: .  step - current time-step
4205: .  ptime - current time
4206: -  dummy - either a viewer or NULL

4208:    Options Database:
4209: .   -ts_monitor_draw_solution_initial - show initial solution as well as current solution

4211:    Notes: the initial solution and current solution are not display with a common axis scaling so generally the option -ts_monitor_draw_solution_initial
4212:        will look bad

4214:    Level: intermediate

4216: .keywords: TS,  vector, monitor, view

4218: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
4219: @*/
4220: PetscErrorCode  TSMonitorDrawSolution(TS ts,PetscInt step,PetscReal ptime,Vec u,void *dummy)
4221: {
4222:   PetscErrorCode   ierr;
4223:   TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy;
4224:   PetscDraw        draw;

4227:   if (!step && ictx->showinitial) {
4228:     if (!ictx->initialsolution) {
4229:       VecDuplicate(u,&ictx->initialsolution);
4230:     }
4231:     VecCopy(u,ictx->initialsolution);
4232:   }
4233:   if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) return(0);

4235:   if (ictx->showinitial) {
4236:     PetscReal pause;
4237:     PetscViewerDrawGetPause(ictx->viewer,&pause);
4238:     PetscViewerDrawSetPause(ictx->viewer,0.0);
4239:     VecView(ictx->initialsolution,ictx->viewer);
4240:     PetscViewerDrawSetPause(ictx->viewer,pause);
4241:     PetscViewerDrawSetHold(ictx->viewer,PETSC_TRUE);
4242:   }
4243:   VecView(u,ictx->viewer);
4244:   if (ictx->showtimestepandtime) {
4245:     PetscReal xl,yl,xr,yr,h;
4246:     char      time[32];

4248:     PetscViewerDrawGetDraw(ictx->viewer,0,&draw);
4249:     PetscSNPrintf(time,32,"Timestep %d Time %g",(int)step,(double)ptime);
4250:     PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);
4251:     h    = yl + .95*(yr - yl);
4252:     PetscDrawStringCentered(draw,.5*(xl+xr),h,PETSC_DRAW_BLACK,time);
4253:     PetscDrawFlush(draw);
4254:   }

4256:   if (ictx->showinitial) {
4257:     PetscViewerDrawSetHold(ictx->viewer,PETSC_FALSE);
4258:   }
4259:   return(0);
4260: }

4262: /*@C
4263:    TSMonitorDrawSolutionPhase - Monitors progress of the TS solvers by plotting the solution as a phase diagram

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:    Level: intermediate

4275: .keywords: TS,  vector, monitor, view

4277: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
4278: @*/
4279: PetscErrorCode  TSMonitorDrawSolutionPhase(TS ts,PetscInt step,PetscReal ptime,Vec u,void *dummy)
4280: {
4281:   PetscErrorCode    ierr;
4282:   TSMonitorDrawCtx  ictx = (TSMonitorDrawCtx)dummy;
4283:   PetscDraw         draw;
4284:   PetscDrawAxis     axis;
4285:   PetscInt          n;
4286:   PetscMPIInt       size;
4287:   PetscReal         U0,U1,xl,yl,xr,yr,h;
4288:   char              time[32];
4289:   const PetscScalar *U;

4292:   MPI_Comm_size(PetscObjectComm((PetscObject)ts),&size);
4293:   if (size != 1) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Only allowed for sequential runs");
4294:   VecGetSize(u,&n);
4295:   if (n != 2) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Only for ODEs with two unknowns");

4297:   PetscViewerDrawGetDraw(ictx->viewer,0,&draw);
4298:   PetscViewerDrawGetDrawAxis(ictx->viewer,0,&axis);
4299:   PetscDrawAxisGetLimits(axis,&xl,&xr,&yl,&yr);
4300:   if (!step) {
4301:     PetscDrawClear(draw);
4302:     PetscDrawAxisDraw(axis);
4303:   }

4305:   VecGetArrayRead(u,&U);
4306:   U0 = PetscRealPart(U[0]);
4307:   U1 = PetscRealPart(U[1]);
4308:   VecRestoreArrayRead(u,&U);
4309:   if ((U0 < xl) || (U1 < yl) || (U0 > xr) || (U1 > yr)) return(0);

4311:   PetscDrawCollectiveBegin(draw);
4312:   PetscDrawPoint(draw,U0,U1,PETSC_DRAW_BLACK);
4313:   if (ictx->showtimestepandtime) {
4314:     PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);
4315:     PetscSNPrintf(time,32,"Timestep %d Time %g",(int)step,(double)ptime);
4316:     h    = yl + .95*(yr - yl);
4317:     PetscDrawStringCentered(draw,.5*(xl+xr),h,PETSC_DRAW_BLACK,time);
4318:   }
4319:   PetscDrawCollectiveEnd(draw);
4320:   PetscDrawFlush(draw);
4321:   PetscDrawPause(draw);
4322:   PetscDrawSave(draw);
4323:   return(0);
4324: }

4326: /*@C
4327:    TSMonitorDrawCtxDestroy - Destroys the monitor context for TSMonitorDrawSolution()

4329:    Collective on TS

4331:    Input Parameters:
4332: .    ctx - the monitor context

4334:    Level: intermediate

4336: .keywords: TS,  vector, monitor, view

4338: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorDrawSolution(), TSMonitorDrawError()
4339: @*/
4340: PetscErrorCode  TSMonitorDrawCtxDestroy(TSMonitorDrawCtx *ictx)
4341: {

4345:   PetscViewerDestroy(&(*ictx)->viewer);
4346:   VecDestroy(&(*ictx)->initialsolution);
4347:   PetscFree(*ictx);
4348:   return(0);
4349: }

4351: /*@C
4352:    TSMonitorDrawCtxCreate - Creates the monitor context for TSMonitorDrawCtx

4354:    Collective on TS

4356:    Input Parameter:
4357: .    ts - time-step context

4359:    Output Patameter:
4360: .    ctx - the monitor context

4362:    Options Database:
4363: .   -ts_monitor_draw_solution_initial - show initial solution as well as current solution

4365:    Level: intermediate

4367: .keywords: TS,  vector, monitor, view

4369: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorDrawCtx()
4370: @*/
4371: PetscErrorCode  TSMonitorDrawCtxCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscInt howoften,TSMonitorDrawCtx *ctx)
4372: {
4373:   PetscErrorCode   ierr;

4376:   PetscNew(ctx);
4377:   PetscViewerDrawOpen(comm,host,label,x,y,m,n,&(*ctx)->viewer);
4378:   PetscViewerSetFromOptions((*ctx)->viewer);

4380:   (*ctx)->howoften    = howoften;
4381:   (*ctx)->showinitial = PETSC_FALSE;
4382:   PetscOptionsGetBool(NULL,NULL,"-ts_monitor_draw_solution_initial",&(*ctx)->showinitial,NULL);

4384:   (*ctx)->showtimestepandtime = PETSC_FALSE;
4385:   PetscOptionsGetBool(NULL,NULL,"-ts_monitor_draw_solution_show_time",&(*ctx)->showtimestepandtime,NULL);
4386:   return(0);
4387: }

4389: /*@C
4390:    TSMonitorDrawSolutionFunction - Monitors progress of the TS solvers by calling
4391:    VecView() for the solution provided by TSSetSolutionFunction() at each timestep

4393:    Collective on TS

4395:    Input Parameters:
4396: +  ts - the TS context
4397: .  step - current time-step
4398: .  ptime - current time
4399: -  dummy - either a viewer or NULL

4401:    Options Database:
4402: .  -ts_monitor_draw_solution_function - Monitor error graphically, requires user to have provided TSSetSolutionFunction()

4404:    Level: intermediate

4406: .keywords: TS,  vector, monitor, view

4408: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSSetSolutionFunction()
4409: @*/
4410: PetscErrorCode  TSMonitorDrawSolutionFunction(TS ts,PetscInt step,PetscReal ptime,Vec u,void *dummy)
4411: {
4412:   PetscErrorCode   ierr;
4413:   TSMonitorDrawCtx ctx    = (TSMonitorDrawCtx)dummy;
4414:   PetscViewer      viewer = ctx->viewer;
4415:   Vec              work;

4418:   if (!(((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason))) return(0);
4419:   VecDuplicate(u,&work);
4420:   TSComputeSolutionFunction(ts,ptime,work);
4421:   VecView(work,viewer);
4422:   VecDestroy(&work);
4423:   return(0);
4424: }

4426: /*@C
4427:    TSMonitorDrawError - Monitors progress of the TS solvers by calling
4428:    VecView() for the error at each timestep

4430:    Collective on TS

4432:    Input Parameters:
4433: +  ts - the TS context
4434: .  step - current time-step
4435: .  ptime - current time
4436: -  dummy - either a viewer or NULL

4438:    Options Database:
4439: .  -ts_monitor_draw_error - Monitor error graphically, requires user to have provided TSSetSolutionFunction()

4441:    Level: intermediate

4443: .keywords: TS,  vector, monitor, view

4445: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSSetSolutionFunction()
4446: @*/
4447: PetscErrorCode  TSMonitorDrawError(TS ts,PetscInt step,PetscReal ptime,Vec u,void *dummy)
4448: {
4449:   PetscErrorCode   ierr;
4450:   TSMonitorDrawCtx ctx    = (TSMonitorDrawCtx)dummy;
4451:   PetscViewer      viewer = ctx->viewer;
4452:   Vec              work;

4455:   if (!(((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason))) return(0);
4456:   VecDuplicate(u,&work);
4457:   TSComputeSolutionFunction(ts,ptime,work);
4458:   VecAXPY(work,-1.0,u);
4459:   VecView(work,viewer);
4460:   VecDestroy(&work);
4461:   return(0);
4462: }

4464:  #include <petsc/private/dmimpl.h>
4465: /*@
4466:    TSSetDM - Sets the DM that may be used by some nonlinear solvers or preconditioners under the TS

4468:    Logically Collective on TS and DM

4470:    Input Parameters:
4471: +  ts - the ODE integrator object
4472: -  dm - the dm, cannot be NULL

4474:    Level: intermediate

4476: .seealso: TSGetDM(), SNESSetDM(), SNESGetDM()
4477: @*/
4478: PetscErrorCode  TSSetDM(TS ts,DM dm)
4479: {
4481:   SNES           snes;
4482:   DMTS           tsdm;

4487:   PetscObjectReference((PetscObject)dm);
4488:   if (ts->dm) {               /* Move the DMTS context over to the new DM unless the new DM already has one */
4489:     if (ts->dm->dmts && !dm->dmts) {
4490:       DMCopyDMTS(ts->dm,dm);
4491:       DMGetDMTS(ts->dm,&tsdm);
4492:       if (tsdm->originaldm == ts->dm) { /* Grant write privileges to the replacement DM */
4493:         tsdm->originaldm = dm;
4494:       }
4495:     }
4496:     DMDestroy(&ts->dm);
4497:   }
4498:   ts->dm = dm;

4500:   TSGetSNES(ts,&snes);
4501:   SNESSetDM(snes,dm);
4502:   return(0);
4503: }

4505: /*@
4506:    TSGetDM - Gets the DM that may be used by some preconditioners

4508:    Not Collective

4510:    Input Parameter:
4511: . ts - the preconditioner context

4513:    Output Parameter:
4514: .  dm - the dm

4516:    Level: intermediate

4518: .seealso: TSSetDM(), SNESSetDM(), SNESGetDM()
4519: @*/
4520: PetscErrorCode  TSGetDM(TS ts,DM *dm)
4521: {

4526:   if (!ts->dm) {
4527:     DMShellCreate(PetscObjectComm((PetscObject)ts),&ts->dm);
4528:     if (ts->snes) {SNESSetDM(ts->snes,ts->dm);}
4529:   }
4530:   *dm = ts->dm;
4531:   return(0);
4532: }

4534: /*@
4535:    SNESTSFormFunction - Function to evaluate nonlinear residual

4537:    Logically Collective on SNES

4539:    Input Parameter:
4540: + snes - nonlinear solver
4541: . U - the current state at which to evaluate the residual
4542: - ctx - user context, must be a TS

4544:    Output Parameter:
4545: . F - the nonlinear residual

4547:    Notes:
4548:    This function is not normally called by users and is automatically registered with the SNES used by TS.
4549:    It is most frequently passed to MatFDColoringSetFunction().

4551:    Level: advanced

4553: .seealso: SNESSetFunction(), MatFDColoringSetFunction()
4554: @*/
4555: PetscErrorCode  SNESTSFormFunction(SNES snes,Vec U,Vec F,void *ctx)
4556: {
4557:   TS             ts = (TS)ctx;

4565:   (ts->ops->snesfunction)(snes,U,F,ts);
4566:   return(0);
4567: }

4569: /*@
4570:    SNESTSFormJacobian - Function to evaluate the Jacobian

4572:    Collective on SNES

4574:    Input Parameter:
4575: + snes - nonlinear solver
4576: . U - the current state at which to evaluate the residual
4577: - ctx - user context, must be a TS

4579:    Output Parameter:
4580: + A - the Jacobian
4581: . B - the preconditioning matrix (may be the same as A)
4582: - flag - indicates any structure change in the matrix

4584:    Notes:
4585:    This function is not normally called by users and is automatically registered with the SNES used by TS.

4587:    Level: developer

4589: .seealso: SNESSetJacobian()
4590: @*/
4591: PetscErrorCode  SNESTSFormJacobian(SNES snes,Vec U,Mat A,Mat B,void *ctx)
4592: {
4593:   TS             ts = (TS)ctx;

4604:   (ts->ops->snesjacobian)(snes,U,A,B,ts);
4605:   return(0);
4606: }

4608: /*@C
4609:    TSComputeRHSFunctionLinear - Evaluate the right hand side via the user-provided Jacobian, for linear problems Udot = A U only

4611:    Collective on TS

4613:    Input Arguments:
4614: +  ts - time stepping context
4615: .  t - time at which to evaluate
4616: .  U - state at which to evaluate
4617: -  ctx - context

4619:    Output Arguments:
4620: .  F - right hand side

4622:    Level: intermediate

4624:    Notes:
4625:    This function is intended to be passed to TSSetRHSFunction() to evaluate the right hand side for linear problems.
4626:    The matrix (and optionally the evaluation context) should be passed to TSSetRHSJacobian().

4628: .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSJacobianConstant()
4629: @*/
4630: PetscErrorCode TSComputeRHSFunctionLinear(TS ts,PetscReal t,Vec U,Vec F,void *ctx)
4631: {
4633:   Mat            Arhs,Brhs;

4636:   TSGetRHSMats_Private(ts,&Arhs,&Brhs);
4637:   TSComputeRHSJacobian(ts,t,U,Arhs,Brhs);
4638:   MatMult(Arhs,U,F);
4639:   return(0);
4640: }

4642: /*@C
4643:    TSComputeRHSJacobianConstant - Reuses a Jacobian that is time-independent.

4645:    Collective on TS

4647:    Input Arguments:
4648: +  ts - time stepping context
4649: .  t - time at which to evaluate
4650: .  U - state at which to evaluate
4651: -  ctx - context

4653:    Output Arguments:
4654: +  A - pointer to operator
4655: .  B - pointer to preconditioning matrix
4656: -  flg - matrix structure flag

4658:    Level: intermediate

4660:    Notes:
4661:    This function is intended to be passed to TSSetRHSJacobian() to evaluate the Jacobian for linear time-independent problems.

4663: .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSFunctionLinear()
4664: @*/
4665: PetscErrorCode TSComputeRHSJacobianConstant(TS ts,PetscReal t,Vec U,Mat A,Mat B,void *ctx)
4666: {
4668:   return(0);
4669: }

4671: /*@C
4672:    TSComputeIFunctionLinear - Evaluate the left hand side via the user-provided Jacobian, for linear problems only

4674:    Collective on TS

4676:    Input Arguments:
4677: +  ts - time stepping context
4678: .  t - time at which to evaluate
4679: .  U - state at which to evaluate
4680: .  Udot - time derivative of state vector
4681: -  ctx - context

4683:    Output Arguments:
4684: .  F - left hand side

4686:    Level: intermediate

4688:    Notes:
4689:    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
4690:    user is required to write their own TSComputeIFunction.
4691:    This function is intended to be passed to TSSetIFunction() to evaluate the left hand side for linear problems.
4692:    The matrix (and optionally the evaluation context) should be passed to TSSetIJacobian().

4694:    Note that using this function is NOT equivalent to using TSComputeRHSFunctionLinear() since that solves Udot = A U

4696: .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIJacobianConstant(), TSComputeRHSFunctionLinear()
4697: @*/
4698: PetscErrorCode TSComputeIFunctionLinear(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void *ctx)
4699: {
4701:   Mat            A,B;

4704:   TSGetIJacobian(ts,&A,&B,NULL,NULL);
4705:   TSComputeIJacobian(ts,t,U,Udot,1.0,A,B,PETSC_TRUE);
4706:   MatMult(A,Udot,F);
4707:   return(0);
4708: }

4710: /*@C
4711:    TSComputeIJacobianConstant - Reuses a time-independent for a semi-implicit DAE or ODE

4713:    Collective on TS

4715:    Input Arguments:
4716: +  ts - time stepping context
4717: .  t - time at which to evaluate
4718: .  U - state at which to evaluate
4719: .  Udot - time derivative of state vector
4720: .  shift - shift to apply
4721: -  ctx - context

4723:    Output Arguments:
4724: +  A - pointer to operator
4725: .  B - pointer to preconditioning matrix
4726: -  flg - matrix structure flag

4728:    Level: advanced

4730:    Notes:
4731:    This function is intended to be passed to TSSetIJacobian() to evaluate the Jacobian for linear time-independent problems.

4733:    It is only appropriate for problems of the form

4735: $     M Udot = F(U,t)

4737:   where M is constant and F is non-stiff.  The user must pass M to TSSetIJacobian().  The current implementation only
4738:   works with IMEX time integration methods such as TSROSW and TSARKIMEX, since there is no support for de-constructing
4739:   an implicit operator of the form

4741: $    shift*M + J

4743:   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
4744:   a copy of M or reassemble it when requested.

4746: .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIFunctionLinear()
4747: @*/
4748: PetscErrorCode TSComputeIJacobianConstant(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal shift,Mat A,Mat B,void *ctx)
4749: {

4753:   MatScale(A, shift / ts->ijacobian.shift);
4754:   ts->ijacobian.shift = shift;
4755:   return(0);
4756: }

4758: /*@
4759:    TSGetEquationType - Gets the type of the equation that TS is solving.

4761:    Not Collective

4763:    Input Parameter:
4764: .  ts - the TS context

4766:    Output Parameter:
4767: .  equation_type - see TSEquationType

4769:    Level: beginner

4771: .keywords: TS, equation type

4773: .seealso: TSSetEquationType(), TSEquationType
4774: @*/
4775: PetscErrorCode  TSGetEquationType(TS ts,TSEquationType *equation_type)
4776: {
4780:   *equation_type = ts->equation_type;
4781:   return(0);
4782: }

4784: /*@
4785:    TSSetEquationType - Sets the type of the equation that TS is solving.

4787:    Not Collective

4789:    Input Parameter:
4790: +  ts - the TS context
4791: -  equation_type - see TSEquationType

4793:    Level: advanced

4795: .keywords: TS, equation type

4797: .seealso: TSGetEquationType(), TSEquationType
4798: @*/
4799: PetscErrorCode  TSSetEquationType(TS ts,TSEquationType equation_type)
4800: {
4803:   ts->equation_type = equation_type;
4804:   return(0);
4805: }

4807: /*@
4808:    TSGetConvergedReason - Gets the reason the TS iteration was stopped.

4810:    Not Collective

4812:    Input Parameter:
4813: .  ts - the TS context

4815:    Output Parameter:
4816: .  reason - negative value indicates diverged, positive value converged, see TSConvergedReason or the
4817:             manual pages for the individual convergence tests for complete lists

4819:    Level: beginner

4821:    Notes:
4822:    Can only be called after the call to TSSolve() is complete.

4824: .keywords: TS, nonlinear, set, convergence, test

4826: .seealso: TSSetConvergenceTest(), TSConvergedReason
4827: @*/
4828: PetscErrorCode  TSGetConvergedReason(TS ts,TSConvergedReason *reason)
4829: {
4833:   *reason = ts->reason;
4834:   return(0);
4835: }

4837: /*@
4838:    TSSetConvergedReason - Sets the reason for handling the convergence of TSSolve.

4840:    Not Collective

4842:    Input Parameter:
4843: +  ts - the TS context
4844: .  reason - negative value indicates diverged, positive value converged, see TSConvergedReason or the
4845:             manual pages for the individual convergence tests for complete lists

4847:    Level: advanced

4849:    Notes:
4850:    Can only be called during TSSolve() is active.

4852: .keywords: TS, nonlinear, set, convergence, test

4854: .seealso: TSConvergedReason
4855: @*/
4856: PetscErrorCode  TSSetConvergedReason(TS ts,TSConvergedReason reason)
4857: {
4860:   ts->reason = reason;
4861:   return(0);
4862: }

4864: /*@
4865:    TSGetSolveTime - Gets the time after a call to TSSolve()

4867:    Not Collective

4869:    Input Parameter:
4870: .  ts - the TS context

4872:    Output Parameter:
4873: .  ftime - the final time. This time corresponds to the final time set with TSSetMaxTime()

4875:    Level: beginner

4877:    Notes:
4878:    Can only be called after the call to TSSolve() is complete.

4880: .keywords: TS, nonlinear, set, convergence, test

4882: .seealso: TSSetConvergenceTest(), TSConvergedReason
4883: @*/
4884: PetscErrorCode  TSGetSolveTime(TS ts,PetscReal *ftime)
4885: {
4889:   *ftime = ts->solvetime;
4890:   return(0);
4891: }

4893: /*@
4894:    TSGetSNESIterations - Gets the total number of nonlinear iterations
4895:    used by the time integrator.

4897:    Not Collective

4899:    Input Parameter:
4900: .  ts - TS context

4902:    Output Parameter:
4903: .  nits - number of nonlinear iterations

4905:    Notes:
4906:    This counter is reset to zero for each successive call to TSSolve().

4908:    Level: intermediate

4910: .keywords: TS, get, number, nonlinear, iterations

4912: .seealso:  TSGetKSPIterations()
4913: @*/
4914: PetscErrorCode TSGetSNESIterations(TS ts,PetscInt *nits)
4915: {
4919:   *nits = ts->snes_its;
4920:   return(0);
4921: }

4923: /*@
4924:    TSGetKSPIterations - Gets the total number of linear iterations
4925:    used by the time integrator.

4927:    Not Collective

4929:    Input Parameter:
4930: .  ts - TS context

4932:    Output Parameter:
4933: .  lits - number of linear iterations

4935:    Notes:
4936:    This counter is reset to zero for each successive call to TSSolve().

4938:    Level: intermediate

4940: .keywords: TS, get, number, linear, iterations

4942: .seealso:  TSGetSNESIterations(), SNESGetKSPIterations()
4943: @*/
4944: PetscErrorCode TSGetKSPIterations(TS ts,PetscInt *lits)
4945: {
4949:   *lits = ts->ksp_its;
4950:   return(0);
4951: }

4953: /*@
4954:    TSGetStepRejections - Gets the total number of rejected steps.

4956:    Not Collective

4958:    Input Parameter:
4959: .  ts - TS context

4961:    Output Parameter:
4962: .  rejects - number of steps rejected

4964:    Notes:
4965:    This counter is reset to zero for each successive call to TSSolve().

4967:    Level: intermediate

4969: .keywords: TS, get, number

4971: .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetSNESFailures(), TSSetMaxSNESFailures(), TSSetErrorIfStepFails()
4972: @*/
4973: PetscErrorCode TSGetStepRejections(TS ts,PetscInt *rejects)
4974: {
4978:   *rejects = ts->reject;
4979:   return(0);
4980: }

4982: /*@
4983:    TSGetSNESFailures - Gets the total number of failed SNES solves

4985:    Not Collective

4987:    Input Parameter:
4988: .  ts - TS context

4990:    Output Parameter:
4991: .  fails - number of failed nonlinear solves

4993:    Notes:
4994:    This counter is reset to zero for each successive call to TSSolve().

4996:    Level: intermediate

4998: .keywords: TS, get, number

5000: .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSSetMaxSNESFailures()
5001: @*/
5002: PetscErrorCode TSGetSNESFailures(TS ts,PetscInt *fails)
5003: {
5007:   *fails = ts->num_snes_failures;
5008:   return(0);
5009: }

5011: /*@
5012:    TSSetMaxStepRejections - Sets the maximum number of step rejections before a step fails

5014:    Not Collective

5016:    Input Parameter:
5017: +  ts - TS context
5018: -  rejects - maximum number of rejected steps, pass -1 for unlimited

5020:    Notes:
5021:    The counter is reset to zero for each step

5023:    Options Database Key:
5024:  .  -ts_max_reject - Maximum number of step rejections before a step fails

5026:    Level: intermediate

5028: .keywords: TS, set, maximum, number

5030: .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxSNESFailures(), TSGetStepRejections(), TSGetSNESFailures(), TSSetErrorIfStepFails(), TSGetConvergedReason()
5031: @*/
5032: PetscErrorCode TSSetMaxStepRejections(TS ts,PetscInt rejects)
5033: {
5036:   ts->max_reject = rejects;
5037:   return(0);
5038: }

5040: /*@
5041:    TSSetMaxSNESFailures - Sets the maximum number of failed SNES solves

5043:    Not Collective

5045:    Input Parameter:
5046: +  ts - TS context
5047: -  fails - maximum number of failed nonlinear solves, pass -1 for unlimited

5049:    Notes:
5050:    The counter is reset to zero for each successive call to TSSolve().

5052:    Options Database Key:
5053:  .  -ts_max_snes_failures - Maximum number of nonlinear solve failures

5055:    Level: intermediate

5057: .keywords: TS, set, maximum, number

5059: .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSGetSNESFailures(), SNESGetConvergedReason(), TSGetConvergedReason()
5060: @*/
5061: PetscErrorCode TSSetMaxSNESFailures(TS ts,PetscInt fails)
5062: {
5065:   ts->max_snes_failures = fails;
5066:   return(0);
5067: }

5069: /*@
5070:    TSSetErrorIfStepFails - Error if no step succeeds

5072:    Not Collective

5074:    Input Parameter:
5075: +  ts - TS context
5076: -  err - PETSC_TRUE to error if no step succeeds, PETSC_FALSE to return without failure

5078:    Options Database Key:
5079:  .  -ts_error_if_step_fails - Error if no step succeeds

5081:    Level: intermediate

5083: .keywords: TS, set, error

5085: .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSGetSNESFailures(), TSSetErrorIfStepFails(), TSGetConvergedReason()
5086: @*/
5087: PetscErrorCode TSSetErrorIfStepFails(TS ts,PetscBool err)
5088: {
5091:   ts->errorifstepfailed = err;
5092:   return(0);
5093: }

5095: /*@C
5096:    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

5098:    Collective on TS

5100:    Input Parameters:
5101: +  ts - the TS context
5102: .  step - current time-step
5103: .  ptime - current time
5104: .  u - current state
5105: -  vf - viewer and its format

5107:    Level: intermediate

5109: .keywords: TS,  vector, monitor, view

5111: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
5112: @*/
5113: PetscErrorCode  TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscViewerAndFormat *vf)
5114: {

5118:   PetscViewerPushFormat(vf->viewer,vf->format);
5119:   VecView(u,vf->viewer);
5120:   PetscViewerPopFormat(vf->viewer);
5121:   return(0);
5122: }

5124: /*@C
5125:    TSMonitorSolutionVTK - Monitors progress of the TS solvers by VecView() for the solution at each timestep.

5127:    Collective on TS

5129:    Input Parameters:
5130: +  ts - the TS context
5131: .  step - current time-step
5132: .  ptime - current time
5133: .  u - current state
5134: -  filenametemplate - string containing a format specifier for the integer time step (e.g. %03D)

5136:    Level: intermediate

5138:    Notes:
5139:    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.
5140:    These are named according to the file name template.

5142:    This function is normally passed as an argument to TSMonitorSet() along with TSMonitorSolutionVTKDestroy().

5144: .keywords: TS,  vector, monitor, view

5146: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
5147: @*/
5148: PetscErrorCode TSMonitorSolutionVTK(TS ts,PetscInt step,PetscReal ptime,Vec u,void *filenametemplate)
5149: {
5151:   char           filename[PETSC_MAX_PATH_LEN];
5152:   PetscViewer    viewer;

5155:   if (step < 0) return(0); /* -1 indicates interpolated solution */
5156:   PetscSNPrintf(filename,sizeof(filename),(const char*)filenametemplate,step);
5157:   PetscViewerVTKOpen(PetscObjectComm((PetscObject)ts),filename,FILE_MODE_WRITE,&viewer);
5158:   VecView(u,viewer);
5159:   PetscViewerDestroy(&viewer);
5160:   return(0);
5161: }

5163: /*@C
5164:    TSMonitorSolutionVTKDestroy - Destroy context for monitoring

5166:    Collective on TS

5168:    Input Parameters:
5169: .  filenametemplate - string containing a format specifier for the integer time step (e.g. %03D)

5171:    Level: intermediate

5173:    Note:
5174:    This function is normally passed to TSMonitorSet() along with TSMonitorSolutionVTK().

5176: .keywords: TS,  vector, monitor, view

5178: .seealso: TSMonitorSet(), TSMonitorSolutionVTK()
5179: @*/
5180: PetscErrorCode TSMonitorSolutionVTKDestroy(void *filenametemplate)
5181: {

5185:   PetscFree(*(char**)filenametemplate);
5186:   return(0);
5187: }

5189: /*@
5190:    TSGetAdapt - Get the adaptive controller context for the current method

5192:    Collective on TS if controller has not been created yet

5194:    Input Arguments:
5195: .  ts - time stepping context

5197:    Output Arguments:
5198: .  adapt - adaptive controller

5200:    Level: intermediate

5202: .seealso: TSAdapt, TSAdaptSetType(), TSAdaptChoose()
5203: @*/
5204: PetscErrorCode TSGetAdapt(TS ts,TSAdapt *adapt)
5205: {

5211:   if (!ts->adapt) {
5212:     TSAdaptCreate(PetscObjectComm((PetscObject)ts),&ts->adapt);
5213:     PetscLogObjectParent((PetscObject)ts,(PetscObject)ts->adapt);
5214:     PetscObjectIncrementTabLevel((PetscObject)ts->adapt,(PetscObject)ts,1);
5215:   }
5216:   *adapt = ts->adapt;
5217:   return(0);
5218: }

5220: /*@
5221:    TSSetTolerances - Set tolerances for local truncation error when using adaptive controller

5223:    Logically Collective

5225:    Input Arguments:
5226: +  ts - time integration context
5227: .  atol - scalar absolute tolerances, PETSC_DECIDE to leave current value
5228: .  vatol - vector of absolute tolerances or NULL, used in preference to atol if present
5229: .  rtol - scalar relative tolerances, PETSC_DECIDE to leave current value
5230: -  vrtol - vector of relative tolerances or NULL, used in preference to atol if present

5232:    Options Database keys:
5233: +  -ts_rtol <rtol> - relative tolerance for local truncation error
5234: -  -ts_atol <atol> Absolute tolerance for local truncation error

5236:    Notes:
5237:    With PETSc's implicit schemes for DAE problems, the calculation of the local truncation error
5238:    (LTE) includes both the differential and the algebraic variables. If one wants the LTE to be
5239:    computed only for the differential or the algebraic part then this can be done using the vector of
5240:    tolerances vatol. For example, by setting the tolerance vector with the desired tolerance for the
5241:    differential part and infinity for the algebraic part, the LTE calculation will include only the
5242:    differential variables.

5244:    Level: beginner

5246: .seealso: TS, TSAdapt, TSVecNormWRMS(), TSGetTolerances()
5247: @*/
5248: PetscErrorCode TSSetTolerances(TS ts,PetscReal atol,Vec vatol,PetscReal rtol,Vec vrtol)
5249: {

5253:   if (atol != PETSC_DECIDE && atol != PETSC_DEFAULT) ts->atol = atol;
5254:   if (vatol) {
5255:     PetscObjectReference((PetscObject)vatol);
5256:     VecDestroy(&ts->vatol);
5257:     ts->vatol = vatol;
5258:   }
5259:   if (rtol != PETSC_DECIDE && rtol != PETSC_DEFAULT) ts->rtol = rtol;
5260:   if (vrtol) {
5261:     PetscObjectReference((PetscObject)vrtol);
5262:     VecDestroy(&ts->vrtol);
5263:     ts->vrtol = vrtol;
5264:   }
5265:   return(0);
5266: }

5268: /*@
5269:    TSGetTolerances - Get tolerances for local truncation error when using adaptive controller

5271:    Logically Collective

5273:    Input Arguments:
5274: .  ts - time integration context

5276:    Output Arguments:
5277: +  atol - scalar absolute tolerances, NULL to ignore
5278: .  vatol - vector of absolute tolerances, NULL to ignore
5279: .  rtol - scalar relative tolerances, NULL to ignore
5280: -  vrtol - vector of relative tolerances, NULL to ignore

5282:    Level: beginner

5284: .seealso: TS, TSAdapt, TSVecNormWRMS(), TSSetTolerances()
5285: @*/
5286: PetscErrorCode TSGetTolerances(TS ts,PetscReal *atol,Vec *vatol,PetscReal *rtol,Vec *vrtol)
5287: {
5289:   if (atol)  *atol  = ts->atol;
5290:   if (vatol) *vatol = ts->vatol;
5291:   if (rtol)  *rtol  = ts->rtol;
5292:   if (vrtol) *vrtol = ts->vrtol;
5293:   return(0);
5294: }

5296: /*@
5297:    TSErrorWeightedNorm2 - compute a weighted 2-norm of the difference between two state vectors

5299:    Collective on TS

5301:    Input Arguments:
5302: +  ts - time stepping context
5303: .  U - state vector, usually ts->vec_sol
5304: -  Y - state vector to be compared to U

5306:    Output Arguments:
5307: .  norm - weighted norm, a value of 1.0 means that the error matches the tolerances
5308: .  norma - weighted norm based on the absolute tolerance, a value of 1.0 means that the error matches the tolerances
5309: .  normr - weighted norm based on the relative tolerance, a value of 1.0 means that the error matches the tolerances

5311:    Level: developer

5313: .seealso: TSErrorWeightedNorm(), TSErrorWeightedNormInfinity()
5314: @*/
5315: PetscErrorCode TSErrorWeightedNorm2(TS ts,Vec U,Vec Y,PetscReal *norm,PetscReal *norma,PetscReal *normr)
5316: {
5317:   PetscErrorCode    ierr;
5318:   PetscInt          i,n,N,rstart;
5319:   PetscInt          n_loc,na_loc,nr_loc;
5320:   PetscReal         n_glb,na_glb,nr_glb;
5321:   const PetscScalar *u,*y;
5322:   PetscReal         sum,suma,sumr,gsum,gsuma,gsumr,diff;
5323:   PetscReal         tol,tola,tolr;
5324:   PetscReal         err_loc[6],err_glb[6];

5336:   if (U == Y) SETERRQ(PetscObjectComm((PetscObject)U),PETSC_ERR_ARG_IDN,"U and Y cannot be the same vector");

5338:   VecGetSize(U,&N);
5339:   VecGetLocalSize(U,&n);
5340:   VecGetOwnershipRange(U,&rstart,NULL);
5341:   VecGetArrayRead(U,&u);
5342:   VecGetArrayRead(Y,&y);
5343:   sum  = 0.; n_loc  = 0;
5344:   suma = 0.; na_loc = 0;
5345:   sumr = 0.; nr_loc = 0;
5346:   if (ts->vatol && ts->vrtol) {
5347:     const PetscScalar *atol,*rtol;
5348:     VecGetArrayRead(ts->vatol,&atol);
5349:     VecGetArrayRead(ts->vrtol,&rtol);
5350:     for (i=0; i<n; i++) {
5351:       diff = PetscAbsScalar(y[i] - u[i]);
5352:       tola = PetscRealPart(atol[i]);
5353:       if(tola>0.){
5354:         suma  += PetscSqr(diff/tola);
5355:         na_loc++;
5356:       }
5357:       tolr = PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5358:       if(tolr>0.){
5359:         sumr  += PetscSqr(diff/tolr);
5360:         nr_loc++;
5361:       }
5362:       tol=tola+tolr;
5363:       if(tol>0.){
5364:         sum  += PetscSqr(diff/tol);
5365:         n_loc++;
5366:       }
5367:     }
5368:     VecRestoreArrayRead(ts->vatol,&atol);
5369:     VecRestoreArrayRead(ts->vrtol,&rtol);
5370:   } else if (ts->vatol) {       /* vector atol, scalar rtol */
5371:     const PetscScalar *atol;
5372:     VecGetArrayRead(ts->vatol,&atol);
5373:     for (i=0; i<n; i++) {
5374:       diff = PetscAbsScalar(y[i] - u[i]);
5375:       tola = PetscRealPart(atol[i]);
5376:       if(tola>0.){
5377:         suma  += PetscSqr(diff/tola);
5378:         na_loc++;
5379:       }
5380:       tolr = ts->rtol * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5381:       if(tolr>0.){
5382:         sumr  += PetscSqr(diff/tolr);
5383:         nr_loc++;
5384:       }
5385:       tol=tola+tolr;
5386:       if(tol>0.){
5387:         sum  += PetscSqr(diff/tol);
5388:         n_loc++;
5389:       }
5390:     }
5391:     VecRestoreArrayRead(ts->vatol,&atol);
5392:   } else if (ts->vrtol) {       /* scalar atol, vector rtol */
5393:     const PetscScalar *rtol;
5394:     VecGetArrayRead(ts->vrtol,&rtol);
5395:     for (i=0; i<n; i++) {
5396:       diff = PetscAbsScalar(y[i] - u[i]);
5397:       tola = ts->atol;
5398:       if(tola>0.){
5399:         suma  += PetscSqr(diff/tola);
5400:         na_loc++;
5401:       }
5402:       tolr = PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5403:       if(tolr>0.){
5404:         sumr  += PetscSqr(diff/tolr);
5405:         nr_loc++;
5406:       }
5407:       tol=tola+tolr;
5408:       if(tol>0.){
5409:         sum  += PetscSqr(diff/tol);
5410:         n_loc++;
5411:       }
5412:     }
5413:     VecRestoreArrayRead(ts->vrtol,&rtol);
5414:   } else {                      /* scalar atol, scalar rtol */
5415:     for (i=0; i<n; i++) {
5416:       diff = PetscAbsScalar(y[i] - u[i]);
5417:      tola = ts->atol;
5418:       if(tola>0.){
5419:         suma  += PetscSqr(diff/tola);
5420:         na_loc++;
5421:       }
5422:       tolr = ts->rtol * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5423:       if(tolr>0.){
5424:         sumr  += PetscSqr(diff/tolr);
5425:         nr_loc++;
5426:       }
5427:       tol=tola+tolr;
5428:       if(tol>0.){
5429:         sum  += PetscSqr(diff/tol);
5430:         n_loc++;
5431:       }
5432:     }
5433:   }
5434:   VecRestoreArrayRead(U,&u);
5435:   VecRestoreArrayRead(Y,&y);

5437:   err_loc[0] = sum;
5438:   err_loc[1] = suma;
5439:   err_loc[2] = sumr;
5440:   err_loc[3] = (PetscReal)n_loc;
5441:   err_loc[4] = (PetscReal)na_loc;
5442:   err_loc[5] = (PetscReal)nr_loc;

5444:   MPIU_Allreduce(err_loc,err_glb,6,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)ts));

5446:   gsum   = err_glb[0];
5447:   gsuma  = err_glb[1];
5448:   gsumr  = err_glb[2];
5449:   n_glb  = err_glb[3];
5450:   na_glb = err_glb[4];
5451:   nr_glb = err_glb[5];

5453:   *norm  = 0.;
5454:   if(n_glb>0. ){*norm  = PetscSqrtReal(gsum  / n_glb );}
5455:   *norma = 0.;
5456:   if(na_glb>0.){*norma = PetscSqrtReal(gsuma / na_glb);}
5457:   *normr = 0.;
5458:   if(nr_glb>0.){*normr = PetscSqrtReal(gsumr / nr_glb);}

5460:   if (PetscIsInfOrNanScalar(*norm)) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
5461:   if (PetscIsInfOrNanScalar(*norma)) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_FP,"Infinite or not-a-number generated in norma");
5462:   if (PetscIsInfOrNanScalar(*normr)) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_FP,"Infinite or not-a-number generated in normr");
5463:   return(0);
5464: }

5466: /*@
5467:    TSErrorWeightedNormInfinity - compute a weighted infinity-norm of the difference between two state vectors

5469:    Collective on TS

5471:    Input Arguments:
5472: +  ts - time stepping context
5473: .  U - state vector, usually ts->vec_sol
5474: -  Y - state vector to be compared to U

5476:    Output Arguments:
5477: .  norm - weighted norm, a value of 1.0 means that the error matches the tolerances
5478: .  norma - weighted norm based on the absolute tolerance, a value of 1.0 means that the error matches the tolerances
5479: .  normr - weighted norm based on the relative tolerance, a value of 1.0 means that the error matches the tolerances

5481:    Level: developer

5483: .seealso: TSErrorWeightedNorm(), TSErrorWeightedNorm2()
5484: @*/
5485: PetscErrorCode TSErrorWeightedNormInfinity(TS ts,Vec U,Vec Y,PetscReal *norm,PetscReal *norma,PetscReal *normr)
5486: {
5487:   PetscErrorCode    ierr;
5488:   PetscInt          i,n,N,rstart;
5489:   const PetscScalar *u,*y;
5490:   PetscReal         max,gmax,maxa,gmaxa,maxr,gmaxr;
5491:   PetscReal         tol,tola,tolr,diff;
5492:   PetscReal         err_loc[3],err_glb[3];

5504:   if (U == Y) SETERRQ(PetscObjectComm((PetscObject)U),PETSC_ERR_ARG_IDN,"U and Y cannot be the same vector");

5506:   VecGetSize(U,&N);
5507:   VecGetLocalSize(U,&n);
5508:   VecGetOwnershipRange(U,&rstart,NULL);
5509:   VecGetArrayRead(U,&u);
5510:   VecGetArrayRead(Y,&y);

5512:   max=0.;
5513:   maxa=0.;
5514:   maxr=0.;

5516:   if (ts->vatol && ts->vrtol) {     /* vector atol, vector rtol */
5517:     const PetscScalar *atol,*rtol;
5518:     VecGetArrayRead(ts->vatol,&atol);
5519:     VecGetArrayRead(ts->vrtol,&rtol);

5521:     for (i=0; i<n; i++) {
5522:       diff = PetscAbsScalar(y[i] - u[i]);
5523:       tola = PetscRealPart(atol[i]);
5524:       tolr = PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5525:       tol  = tola+tolr;
5526:       if(tola>0.){
5527:         maxa = PetscMax(maxa,diff / tola);
5528:       }
5529:       if(tolr>0.){
5530:         maxr = PetscMax(maxr,diff / tolr);
5531:       }
5532:       if(tol>0.){
5533:         max = PetscMax(max,diff / tol);
5534:       }
5535:     }
5536:     VecRestoreArrayRead(ts->vatol,&atol);
5537:     VecRestoreArrayRead(ts->vrtol,&rtol);
5538:   } else if (ts->vatol) {       /* vector atol, scalar rtol */
5539:     const PetscScalar *atol;
5540:     VecGetArrayRead(ts->vatol,&atol);
5541:     for (i=0; i<n; i++) {
5542:       diff = PetscAbsScalar(y[i] - u[i]);
5543:       tola = PetscRealPart(atol[i]);
5544:       tolr = ts->rtol  * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5545:       tol  = tola+tolr;
5546:       if(tola>0.){
5547:         maxa = PetscMax(maxa,diff / tola);
5548:       }
5549:       if(tolr>0.){
5550:         maxr = PetscMax(maxr,diff / tolr);
5551:       }
5552:       if(tol>0.){
5553:         max = PetscMax(max,diff / tol);
5554:       }
5555:     }
5556:     VecRestoreArrayRead(ts->vatol,&atol);
5557:   } else if (ts->vrtol) {       /* scalar atol, vector rtol */
5558:     const PetscScalar *rtol;
5559:     VecGetArrayRead(ts->vrtol,&rtol);

5561:     for (i=0; i<n; i++) {
5562:       diff = PetscAbsScalar(y[i] - u[i]);
5563:       tola = ts->atol;
5564:       tolr = PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5565:       tol  = tola+tolr;
5566:       if(tola>0.){
5567:         maxa = PetscMax(maxa,diff / tola);
5568:       }
5569:       if(tolr>0.){
5570:         maxr = PetscMax(maxr,diff / tolr);
5571:       }
5572:       if(tol>0.){
5573:         max = PetscMax(max,diff / tol);
5574:       }
5575:     }
5576:     VecRestoreArrayRead(ts->vrtol,&rtol);
5577:   } else {                      /* scalar atol, scalar rtol */

5579:     for (i=0; i<n; i++) {
5580:       diff = PetscAbsScalar(y[i] - u[i]);
5581:       tola = ts->atol;
5582:       tolr = ts->rtol * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5583:       tol  = tola+tolr;
5584:       if(tola>0.){
5585:         maxa = PetscMax(maxa,diff / tola);
5586:       }
5587:       if(tolr>0.){
5588:         maxr = PetscMax(maxr,diff / tolr);
5589:       }
5590:       if(tol>0.){
5591:         max = PetscMax(max,diff / tol);
5592:       }
5593:     }
5594:   }
5595:   VecRestoreArrayRead(U,&u);
5596:   VecRestoreArrayRead(Y,&y);
5597:   err_loc[0] = max;
5598:   err_loc[1] = maxa;
5599:   err_loc[2] = maxr;
5600:   MPIU_Allreduce(err_loc,err_glb,3,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)ts));
5601:   gmax   = err_glb[0];
5602:   gmaxa  = err_glb[1];
5603:   gmaxr  = err_glb[2];

5605:   *norm = gmax;
5606:   *norma = gmaxa;
5607:   *normr = gmaxr;
5608:   if (PetscIsInfOrNanScalar(*norm)) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
5609:     if (PetscIsInfOrNanScalar(*norma)) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_FP,"Infinite or not-a-number generated in norma");
5610:     if (PetscIsInfOrNanScalar(*normr)) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_FP,"Infinite or not-a-number generated in normr");
5611:   return(0);
5612: }

5614: /*@
5615:    TSErrorWeightedNorm - compute a weighted norm of the difference between two state vectors based on supplied absolute and relative tolerances

5617:    Collective on TS

5619:    Input Arguments:
5620: +  ts - time stepping context
5621: .  U - state vector, usually ts->vec_sol
5622: .  Y - state vector to be compared to U
5623: -  wnormtype - norm type, either NORM_2 or NORM_INFINITY

5625:    Output Arguments:
5626: .  norm  - weighted norm, a value of 1.0 achieves a balance between absolute and relative tolerances
5627: .  norma - weighted norm, a value of 1.0 means that the error meets the absolute tolerance set by the user
5628: .  normr - weighted norm, a value of 1.0 means that the error meets the relative tolerance set by the user

5630:    Options Database Keys:
5631: .  -ts_adapt_wnormtype <wnormtype> - 2, INFINITY

5633:    Level: developer

5635: .seealso: TSErrorWeightedNormInfinity(), TSErrorWeightedNorm2(), TSErrorWeightedENorm
5636: @*/
5637: PetscErrorCode TSErrorWeightedNorm(TS ts,Vec U,Vec Y,NormType wnormtype,PetscReal *norm,PetscReal *norma,PetscReal *normr)
5638: {

5642:   if (wnormtype == NORM_2) {
5643:     TSErrorWeightedNorm2(ts,U,Y,norm,norma,normr);
5644:   } else if(wnormtype == NORM_INFINITY) {
5645:     TSErrorWeightedNormInfinity(ts,U,Y,norm,norma,normr);
5646:   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for norm type %s",NormTypes[wnormtype]);
5647:   return(0);
5648: }


5651: /*@
5652:    TSErrorWeightedENorm2 - compute a weighted 2 error norm based on supplied absolute and relative tolerances

5654:    Collective on TS

5656:    Input Arguments:
5657: +  ts - time stepping context
5658: .  E - error vector
5659: .  U - state vector, usually ts->vec_sol
5660: -  Y - state vector, previous time step

5662:    Output Arguments:
5663: .  norm - weighted norm, a value of 1.0 means that the error matches the tolerances
5664: .  norma - weighted norm based on the absolute tolerance, a value of 1.0 means that the error matches the tolerances
5665: .  normr - weighted norm based on the relative tolerance, a value of 1.0 means that the error matches the tolerances

5667:    Level: developer

5669: .seealso: TSErrorWeightedENorm(), TSErrorWeightedENormInfinity()
5670: @*/
5671: PetscErrorCode TSErrorWeightedENorm2(TS ts,Vec E,Vec U,Vec Y,PetscReal *norm,PetscReal *norma,PetscReal *normr)
5672: {
5673:   PetscErrorCode    ierr;
5674:   PetscInt          i,n,N,rstart;
5675:   PetscInt          n_loc,na_loc,nr_loc;
5676:   PetscReal         n_glb,na_glb,nr_glb;
5677:   const PetscScalar *e,*u,*y;
5678:   PetscReal         err,sum,suma,sumr,gsum,gsuma,gsumr;
5679:   PetscReal         tol,tola,tolr;
5680:   PetscReal         err_loc[6],err_glb[6];


5696:   VecGetSize(E,&N);
5697:   VecGetLocalSize(E,&n);
5698:   VecGetOwnershipRange(E,&rstart,NULL);
5699:   VecGetArrayRead(E,&e);
5700:   VecGetArrayRead(U,&u);
5701:   VecGetArrayRead(Y,&y);
5702:   sum  = 0.; n_loc  = 0;
5703:   suma = 0.; na_loc = 0;
5704:   sumr = 0.; nr_loc = 0;
5705:   if (ts->vatol && ts->vrtol) {
5706:     const PetscScalar *atol,*rtol;
5707:     VecGetArrayRead(ts->vatol,&atol);
5708:     VecGetArrayRead(ts->vrtol,&rtol);
5709:     for (i=0; i<n; i++) {
5710:       err = PetscAbsScalar(e[i]);
5711:       tola = PetscRealPart(atol[i]);
5712:       if(tola>0.){
5713:         suma  += PetscSqr(err/tola);
5714:         na_loc++;
5715:       }
5716:       tolr = PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5717:       if(tolr>0.){
5718:         sumr  += PetscSqr(err/tolr);
5719:         nr_loc++;
5720:       }
5721:       tol=tola+tolr;
5722:       if(tol>0.){
5723:         sum  += PetscSqr(err/tol);
5724:         n_loc++;
5725:       }
5726:     }
5727:     VecRestoreArrayRead(ts->vatol,&atol);
5728:     VecRestoreArrayRead(ts->vrtol,&rtol);
5729:   } else if (ts->vatol) {       /* vector atol, scalar rtol */
5730:     const PetscScalar *atol;
5731:     VecGetArrayRead(ts->vatol,&atol);
5732:     for (i=0; i<n; i++) {
5733:       err = PetscAbsScalar(e[i]);
5734:       tola = PetscRealPart(atol[i]);
5735:       if(tola>0.){
5736:         suma  += PetscSqr(err/tola);
5737:         na_loc++;
5738:       }
5739:       tolr = ts->rtol * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5740:       if(tolr>0.){
5741:         sumr  += PetscSqr(err/tolr);
5742:         nr_loc++;
5743:       }
5744:       tol=tola+tolr;
5745:       if(tol>0.){
5746:         sum  += PetscSqr(err/tol);
5747:         n_loc++;
5748:       }
5749:     }
5750:     VecRestoreArrayRead(ts->vatol,&atol);
5751:   } else if (ts->vrtol) {       /* scalar atol, vector rtol */
5752:     const PetscScalar *rtol;
5753:     VecGetArrayRead(ts->vrtol,&rtol);
5754:     for (i=0; i<n; i++) {
5755:       err = PetscAbsScalar(e[i]);
5756:       tola = ts->atol;
5757:       if(tola>0.){
5758:         suma  += PetscSqr(err/tola);
5759:         na_loc++;
5760:       }
5761:       tolr = PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5762:       if(tolr>0.){
5763:         sumr  += PetscSqr(err/tolr);
5764:         nr_loc++;
5765:       }
5766:       tol=tola+tolr;
5767:       if(tol>0.){
5768:         sum  += PetscSqr(err/tol);
5769:         n_loc++;
5770:       }
5771:     }
5772:     VecRestoreArrayRead(ts->vrtol,&rtol);
5773:   } else {                      /* scalar atol, scalar rtol */
5774:     for (i=0; i<n; i++) {
5775:       err = PetscAbsScalar(e[i]);
5776:      tola = ts->atol;
5777:       if(tola>0.){
5778:         suma  += PetscSqr(err/tola);
5779:         na_loc++;
5780:       }
5781:       tolr = ts->rtol * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5782:       if(tolr>0.){
5783:         sumr  += PetscSqr(err/tolr);
5784:         nr_loc++;
5785:       }
5786:       tol=tola+tolr;
5787:       if(tol>0.){
5788:         sum  += PetscSqr(err/tol);
5789:         n_loc++;
5790:       }
5791:     }
5792:   }
5793:   VecRestoreArrayRead(E,&e);
5794:   VecRestoreArrayRead(U,&u);
5795:   VecRestoreArrayRead(Y,&y);

5797:   err_loc[0] = sum;
5798:   err_loc[1] = suma;
5799:   err_loc[2] = sumr;
5800:   err_loc[3] = (PetscReal)n_loc;
5801:   err_loc[4] = (PetscReal)na_loc;
5802:   err_loc[5] = (PetscReal)nr_loc;

5804:   MPIU_Allreduce(err_loc,err_glb,6,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)ts));

5806:   gsum   = err_glb[0];
5807:   gsuma  = err_glb[1];
5808:   gsumr  = err_glb[2];
5809:   n_glb  = err_glb[3];
5810:   na_glb = err_glb[4];
5811:   nr_glb = err_glb[5];

5813:   *norm  = 0.;
5814:   if(n_glb>0. ){*norm  = PetscSqrtReal(gsum  / n_glb );}
5815:   *norma = 0.;
5816:   if(na_glb>0.){*norma = PetscSqrtReal(gsuma / na_glb);}
5817:   *normr = 0.;
5818:   if(nr_glb>0.){*normr = PetscSqrtReal(gsumr / nr_glb);}

5820:   if (PetscIsInfOrNanScalar(*norm)) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
5821:   if (PetscIsInfOrNanScalar(*norma)) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_FP,"Infinite or not-a-number generated in norma");
5822:   if (PetscIsInfOrNanScalar(*normr)) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_FP,"Infinite or not-a-number generated in normr");
5823:   return(0);
5824: }

5826: /*@
5827:    TSErrorWeightedENormInfinity - compute a weighted infinity error norm based on supplied absolute and relative tolerances
5828:    Collective on TS

5830:    Input Arguments:
5831: +  ts - time stepping context
5832: .  E - error vector
5833: .  U - state vector, usually ts->vec_sol
5834: -  Y - state vector, previous time step

5836:    Output Arguments:
5837: .  norm - weighted norm, a value of 1.0 means that the error matches the tolerances
5838: .  norma - weighted norm based on the absolute tolerance, a value of 1.0 means that the error matches the tolerances
5839: .  normr - weighted norm based on the relative tolerance, a value of 1.0 means that the error matches the tolerances

5841:    Level: developer

5843: .seealso: TSErrorWeightedENorm(), TSErrorWeightedENorm2()
5844: @*/
5845: PetscErrorCode TSErrorWeightedENormInfinity(TS ts,Vec E,Vec U,Vec Y,PetscReal *norm,PetscReal *norma,PetscReal *normr)
5846: {
5847:   PetscErrorCode    ierr;
5848:   PetscInt          i,n,N,rstart;
5849:   const PetscScalar *e,*u,*y;
5850:   PetscReal         err,max,gmax,maxa,gmaxa,maxr,gmaxr;
5851:   PetscReal         tol,tola,tolr;
5852:   PetscReal         err_loc[3],err_glb[3];


5868:   VecGetSize(E,&N);
5869:   VecGetLocalSize(E,&n);
5870:   VecGetOwnershipRange(E,&rstart,NULL);
5871:   VecGetArrayRead(E,&e);
5872:   VecGetArrayRead(U,&u);
5873:   VecGetArrayRead(Y,&y);

5875:   max=0.;
5876:   maxa=0.;
5877:   maxr=0.;

5879:   if (ts->vatol && ts->vrtol) {     /* vector atol, vector rtol */
5880:     const PetscScalar *atol,*rtol;
5881:     VecGetArrayRead(ts->vatol,&atol);
5882:     VecGetArrayRead(ts->vrtol,&rtol);

5884:     for (i=0; i<n; i++) {
5885:       err = PetscAbsScalar(e[i]);
5886:       tola = PetscRealPart(atol[i]);
5887:       tolr = PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5888:       tol  = tola+tolr;
5889:       if(tola>0.){
5890:         maxa = PetscMax(maxa,err / tola);
5891:       }
5892:       if(tolr>0.){
5893:         maxr = PetscMax(maxr,err / tolr);
5894:       }
5895:       if(tol>0.){
5896:         max = PetscMax(max,err / tol);
5897:       }
5898:     }
5899:     VecRestoreArrayRead(ts->vatol,&atol);
5900:     VecRestoreArrayRead(ts->vrtol,&rtol);
5901:   } else if (ts->vatol) {       /* vector atol, scalar rtol */
5902:     const PetscScalar *atol;
5903:     VecGetArrayRead(ts->vatol,&atol);
5904:     for (i=0; i<n; i++) {
5905:       err = PetscAbsScalar(e[i]);
5906:       tola = PetscRealPart(atol[i]);
5907:       tolr = ts->rtol  * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5908:       tol  = tola+tolr;
5909:       if(tola>0.){
5910:         maxa = PetscMax(maxa,err / tola);
5911:       }
5912:       if(tolr>0.){
5913:         maxr = PetscMax(maxr,err / tolr);
5914:       }
5915:       if(tol>0.){
5916:         max = PetscMax(max,err / tol);
5917:       }
5918:     }
5919:     VecRestoreArrayRead(ts->vatol,&atol);
5920:   } else if (ts->vrtol) {       /* scalar atol, vector rtol */
5921:     const PetscScalar *rtol;
5922:     VecGetArrayRead(ts->vrtol,&rtol);

5924:     for (i=0; i<n; i++) {
5925:       err = PetscAbsScalar(e[i]);
5926:       tola = ts->atol;
5927:       tolr = PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5928:       tol  = tola+tolr;
5929:       if(tola>0.){
5930:         maxa = PetscMax(maxa,err / tola);
5931:       }
5932:       if(tolr>0.){
5933:         maxr = PetscMax(maxr,err / tolr);
5934:       }
5935:       if(tol>0.){
5936:         max = PetscMax(max,err / tol);
5937:       }
5938:     }
5939:     VecRestoreArrayRead(ts->vrtol,&rtol);
5940:   } else {                      /* scalar atol, scalar rtol */

5942:     for (i=0; i<n; i++) {
5943:       err = PetscAbsScalar(e[i]);
5944:       tola = ts->atol;
5945:       tolr = ts->rtol * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5946:       tol  = tola+tolr;
5947:       if(tola>0.){
5948:         maxa = PetscMax(maxa,err / tola);
5949:       }
5950:       if(tolr>0.){
5951:         maxr = PetscMax(maxr,err / tolr);
5952:       }
5953:       if(tol>0.){
5954:         max = PetscMax(max,err / tol);
5955:       }
5956:     }
5957:   }
5958:   VecRestoreArrayRead(E,&e);
5959:   VecRestoreArrayRead(U,&u);
5960:   VecRestoreArrayRead(Y,&y);
5961:   err_loc[0] = max;
5962:   err_loc[1] = maxa;
5963:   err_loc[2] = maxr;
5964:   MPIU_Allreduce(err_loc,err_glb,3,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)ts));
5965:   gmax   = err_glb[0];
5966:   gmaxa  = err_glb[1];
5967:   gmaxr  = err_glb[2];

5969:   *norm = gmax;
5970:   *norma = gmaxa;
5971:   *normr = gmaxr;
5972:   if (PetscIsInfOrNanScalar(*norm)) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
5973:     if (PetscIsInfOrNanScalar(*norma)) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_FP,"Infinite or not-a-number generated in norma");
5974:     if (PetscIsInfOrNanScalar(*normr)) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_FP,"Infinite or not-a-number generated in normr");
5975:   return(0);
5976: }

5978: /*@
5979:    TSErrorWeightedENorm - compute a weighted error norm based on supplied absolute and relative tolerances

5981:    Collective on TS

5983:    Input Arguments:
5984: +  ts - time stepping context
5985: .  E - error vector
5986: .  U - state vector, usually ts->vec_sol
5987: .  Y - state vector, previous time step
5988: -  wnormtype - norm type, either NORM_2 or NORM_INFINITY

5990:    Output Arguments:
5991: .  norm  - weighted norm, a value of 1.0 achieves a balance between absolute and relative tolerances
5992: .  norma - weighted norm, a value of 1.0 means that the error meets the absolute tolerance set by the user
5993: .  normr - weighted norm, a value of 1.0 means that the error meets the relative tolerance set by the user

5995:    Options Database Keys:
5996: .  -ts_adapt_wnormtype <wnormtype> - 2, INFINITY

5998:    Level: developer

6000: .seealso: TSErrorWeightedENormInfinity(), TSErrorWeightedENorm2(), TSErrorWeightedNormInfinity(), TSErrorWeightedNorm2()
6001: @*/
6002: PetscErrorCode TSErrorWeightedENorm(TS ts,Vec E,Vec U,Vec Y,NormType wnormtype,PetscReal *norm,PetscReal *norma,PetscReal *normr)
6003: {

6007:   if (wnormtype == NORM_2) {
6008:     TSErrorWeightedENorm2(ts,E,U,Y,norm,norma,normr);
6009:   } else if(wnormtype == NORM_INFINITY) {
6010:     TSErrorWeightedENormInfinity(ts,E,U,Y,norm,norma,normr);
6011:   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for norm type %s",NormTypes[wnormtype]);
6012:   return(0);
6013: }


6016: /*@
6017:    TSSetCFLTimeLocal - Set the local CFL constraint relative to forward Euler

6019:    Logically Collective on TS

6021:    Input Arguments:
6022: +  ts - time stepping context
6023: -  cfltime - maximum stable time step if using forward Euler (value can be different on each process)

6025:    Note:
6026:    After calling this function, the global CFL time can be obtained by calling TSGetCFLTime()

6028:    Level: intermediate

6030: .seealso: TSGetCFLTime(), TSADAPTCFL
6031: @*/
6032: PetscErrorCode TSSetCFLTimeLocal(TS ts,PetscReal cfltime)
6033: {
6036:   ts->cfltime_local = cfltime;
6037:   ts->cfltime       = -1.;
6038:   return(0);
6039: }

6041: /*@
6042:    TSGetCFLTime - Get the maximum stable time step according to CFL criteria applied to forward Euler

6044:    Collective on TS

6046:    Input Arguments:
6047: .  ts - time stepping context

6049:    Output Arguments:
6050: .  cfltime - maximum stable time step for forward Euler

6052:    Level: advanced

6054: .seealso: TSSetCFLTimeLocal()
6055: @*/
6056: PetscErrorCode TSGetCFLTime(TS ts,PetscReal *cfltime)
6057: {

6061:   if (ts->cfltime < 0) {
6062:     MPIU_Allreduce(&ts->cfltime_local,&ts->cfltime,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)ts));
6063:   }
6064:   *cfltime = ts->cfltime;
6065:   return(0);
6066: }

6068: /*@
6069:    TSVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu

6071:    Input Parameters:
6072: .  ts   - the TS context.
6073: .  xl   - lower bound.
6074: .  xu   - upper bound.

6076:    Notes:
6077:    If this routine is not called then the lower and upper bounds are set to
6078:    PETSC_NINFINITY and PETSC_INFINITY respectively during SNESSetUp().

6080:    Level: advanced

6082: @*/
6083: PetscErrorCode TSVISetVariableBounds(TS ts, Vec xl, Vec xu)
6084: {
6086:   SNES           snes;

6089:   TSGetSNES(ts,&snes);
6090:   SNESVISetVariableBounds(snes,xl,xu);
6091:   return(0);
6092: }

6094: #if defined(PETSC_HAVE_MATLAB_ENGINE)
6095: #include <mex.h>

6097: typedef struct {char *funcname; mxArray *ctx;} TSMatlabContext;

6099: /*
6100:    TSComputeFunction_Matlab - Calls the function that has been set with
6101:                          TSSetFunctionMatlab().

6103:    Collective on TS

6105:    Input Parameters:
6106: +  snes - the TS context
6107: -  u - input vector

6109:    Output Parameter:
6110: .  y - function vector, as set by TSSetFunction()

6112:    Notes:
6113:    TSComputeFunction() is typically used within nonlinear solvers
6114:    implementations, so most users would not generally call this routine
6115:    themselves.

6117:    Level: developer

6119: .keywords: TS, nonlinear, compute, function

6121: .seealso: TSSetFunction(), TSGetFunction()
6122: */
6123: PetscErrorCode  TSComputeFunction_Matlab(TS snes,PetscReal time,Vec u,Vec udot,Vec y, void *ctx)
6124: {
6125:   PetscErrorCode  ierr;
6126:   TSMatlabContext *sctx = (TSMatlabContext*)ctx;
6127:   int             nlhs  = 1,nrhs = 7;
6128:   mxArray         *plhs[1],*prhs[7];
6129:   long long int   lx = 0,lxdot = 0,ly = 0,ls = 0;


6139:   PetscMemcpy(&ls,&snes,sizeof(snes));
6140:   PetscMemcpy(&lx,&u,sizeof(u));
6141:   PetscMemcpy(&lxdot,&udot,sizeof(udot));
6142:   PetscMemcpy(&ly,&y,sizeof(u));

6144:   prhs[0] =  mxCreateDoubleScalar((double)ls);
6145:   prhs[1] =  mxCreateDoubleScalar(time);
6146:   prhs[2] =  mxCreateDoubleScalar((double)lx);
6147:   prhs[3] =  mxCreateDoubleScalar((double)lxdot);
6148:   prhs[4] =  mxCreateDoubleScalar((double)ly);
6149:   prhs[5] =  mxCreateString(sctx->funcname);
6150:   prhs[6] =  sctx->ctx;
6151:    mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeFunctionInternal");
6152:    mxGetScalar(plhs[0]);
6153:   mxDestroyArray(prhs[0]);
6154:   mxDestroyArray(prhs[1]);
6155:   mxDestroyArray(prhs[2]);
6156:   mxDestroyArray(prhs[3]);
6157:   mxDestroyArray(prhs[4]);
6158:   mxDestroyArray(prhs[5]);
6159:   mxDestroyArray(plhs[0]);
6160:   return(0);
6161: }

6163: /*
6164:    TSSetFunctionMatlab - Sets the function evaluation routine and function
6165:    vector for use by the TS routines in solving ODEs
6166:    equations from MATLAB. Here the function is a string containing the name of a MATLAB function

6168:    Logically Collective on TS

6170:    Input Parameters:
6171: +  ts - the TS context
6172: -  func - function evaluation routine

6174:    Calling sequence of func:
6175: $    func (TS ts,PetscReal time,Vec u,Vec udot,Vec f,void *ctx);

6177:    Level: beginner

6179: .keywords: TS, nonlinear, set, function

6181: .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
6182: */
6183: PetscErrorCode  TSSetFunctionMatlab(TS ts,const char *func,mxArray *ctx)
6184: {
6185:   PetscErrorCode  ierr;
6186:   TSMatlabContext *sctx;

6189:   /* currently sctx is memory bleed */
6190:   PetscNew(&sctx);
6191:   PetscStrallocpy(func,&sctx->funcname);
6192:   /*
6193:      This should work, but it doesn't
6194:   sctx->ctx = ctx;
6195:   mexMakeArrayPersistent(sctx->ctx);
6196:   */
6197:   sctx->ctx = mxDuplicateArray(ctx);

6199:   TSSetIFunction(ts,NULL,TSComputeFunction_Matlab,sctx);
6200:   return(0);
6201: }

6203: /*
6204:    TSComputeJacobian_Matlab - Calls the function that has been set with
6205:                          TSSetJacobianMatlab().

6207:    Collective on TS

6209:    Input Parameters:
6210: +  ts - the TS context
6211: .  u - input vector
6212: .  A, B - the matrices
6213: -  ctx - user context

6215:    Level: developer

6217: .keywords: TS, nonlinear, compute, function

6219: .seealso: TSSetFunction(), TSGetFunction()
6220: @*/
6221: PetscErrorCode  TSComputeJacobian_Matlab(TS ts,PetscReal time,Vec u,Vec udot,PetscReal shift,Mat A,Mat B,void *ctx)
6222: {
6223:   PetscErrorCode  ierr;
6224:   TSMatlabContext *sctx = (TSMatlabContext*)ctx;
6225:   int             nlhs  = 2,nrhs = 9;
6226:   mxArray         *plhs[2],*prhs[9];
6227:   long long int   lx = 0,lxdot = 0,lA = 0,ls = 0, lB = 0;


6233:   /* call Matlab function in ctx with arguments u and y */

6235:   PetscMemcpy(&ls,&ts,sizeof(ts));
6236:   PetscMemcpy(&lx,&u,sizeof(u));
6237:   PetscMemcpy(&lxdot,&udot,sizeof(u));
6238:   PetscMemcpy(&lA,A,sizeof(u));
6239:   PetscMemcpy(&lB,B,sizeof(u));

6241:   prhs[0] =  mxCreateDoubleScalar((double)ls);
6242:   prhs[1] =  mxCreateDoubleScalar((double)time);
6243:   prhs[2] =  mxCreateDoubleScalar((double)lx);
6244:   prhs[3] =  mxCreateDoubleScalar((double)lxdot);
6245:   prhs[4] =  mxCreateDoubleScalar((double)shift);
6246:   prhs[5] =  mxCreateDoubleScalar((double)lA);
6247:   prhs[6] =  mxCreateDoubleScalar((double)lB);
6248:   prhs[7] =  mxCreateString(sctx->funcname);
6249:   prhs[8] =  sctx->ctx;
6250:    mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeJacobianInternal");
6251:    mxGetScalar(plhs[0]);
6252:   mxDestroyArray(prhs[0]);
6253:   mxDestroyArray(prhs[1]);
6254:   mxDestroyArray(prhs[2]);
6255:   mxDestroyArray(prhs[3]);
6256:   mxDestroyArray(prhs[4]);
6257:   mxDestroyArray(prhs[5]);
6258:   mxDestroyArray(prhs[6]);
6259:   mxDestroyArray(prhs[7]);
6260:   mxDestroyArray(plhs[0]);
6261:   mxDestroyArray(plhs[1]);
6262:   return(0);
6263: }

6265: /*
6266:    TSSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices
6267:    vector for use by the TS routines in solving ODEs from MATLAB. Here the function is a string containing the name of a MATLAB function

6269:    Logically Collective on TS

6271:    Input Parameters:
6272: +  ts - the TS context
6273: .  A,B - Jacobian matrices
6274: .  func - function evaluation routine
6275: -  ctx - user context

6277:    Calling sequence of func:
6278: $    flag = func (TS ts,PetscReal time,Vec u,Vec udot,Mat A,Mat B,void *ctx);

6280:    Level: developer

6282: .keywords: TS, nonlinear, set, function

6284: .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
6285: */
6286: PetscErrorCode  TSSetJacobianMatlab(TS ts,Mat A,Mat B,const char *func,mxArray *ctx)
6287: {
6288:   PetscErrorCode  ierr;
6289:   TSMatlabContext *sctx;

6292:   /* currently sctx is memory bleed */
6293:   PetscNew(&sctx);
6294:   PetscStrallocpy(func,&sctx->funcname);
6295:   /*
6296:      This should work, but it doesn't
6297:   sctx->ctx = ctx;
6298:   mexMakeArrayPersistent(sctx->ctx);
6299:   */
6300:   sctx->ctx = mxDuplicateArray(ctx);

6302:   TSSetIJacobian(ts,A,B,TSComputeJacobian_Matlab,sctx);
6303:   return(0);
6304: }

6306: /*
6307:    TSMonitor_Matlab - Calls the function that has been set with TSMonitorSetMatlab().

6309:    Collective on TS

6311: .seealso: TSSetFunction(), TSGetFunction()
6312: @*/
6313: PetscErrorCode  TSMonitor_Matlab(TS ts,PetscInt it, PetscReal time,Vec u, void *ctx)
6314: {
6315:   PetscErrorCode  ierr;
6316:   TSMatlabContext *sctx = (TSMatlabContext*)ctx;
6317:   int             nlhs  = 1,nrhs = 6;
6318:   mxArray         *plhs[1],*prhs[6];
6319:   long long int   lx = 0,ls = 0;


6325:   PetscMemcpy(&ls,&ts,sizeof(ts));
6326:   PetscMemcpy(&lx,&u,sizeof(u));

6328:   prhs[0] =  mxCreateDoubleScalar((double)ls);
6329:   prhs[1] =  mxCreateDoubleScalar((double)it);
6330:   prhs[2] =  mxCreateDoubleScalar((double)time);
6331:   prhs[3] =  mxCreateDoubleScalar((double)lx);
6332:   prhs[4] =  mxCreateString(sctx->funcname);
6333:   prhs[5] =  sctx->ctx;
6334:    mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSMonitorInternal");
6335:    mxGetScalar(plhs[0]);
6336:   mxDestroyArray(prhs[0]);
6337:   mxDestroyArray(prhs[1]);
6338:   mxDestroyArray(prhs[2]);
6339:   mxDestroyArray(prhs[3]);
6340:   mxDestroyArray(prhs[4]);
6341:   mxDestroyArray(plhs[0]);
6342:   return(0);
6343: }

6345: /*
6346:    TSMonitorSetMatlab - Sets the monitor function from Matlab

6348:    Level: developer

6350: .keywords: TS, nonlinear, set, function

6352: .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
6353: */
6354: PetscErrorCode  TSMonitorSetMatlab(TS ts,const char *func,mxArray *ctx)
6355: {
6356:   PetscErrorCode  ierr;
6357:   TSMatlabContext *sctx;

6360:   /* currently sctx is memory bleed */
6361:   PetscNew(&sctx);
6362:   PetscStrallocpy(func,&sctx->funcname);
6363:   /*
6364:      This should work, but it doesn't
6365:   sctx->ctx = ctx;
6366:   mexMakeArrayPersistent(sctx->ctx);
6367:   */
6368:   sctx->ctx = mxDuplicateArray(ctx);

6370:   TSMonitorSet(ts,TSMonitor_Matlab,sctx,NULL);
6371:   return(0);
6372: }
6373: #endif

6375: /*@C
6376:    TSMonitorLGSolution - Monitors progress of the TS solvers by plotting each component of the solution vector
6377:        in a time based line graph

6379:    Collective on TS

6381:    Input Parameters:
6382: +  ts - the TS context
6383: .  step - current time-step
6384: .  ptime - current time
6385: .  u - current solution
6386: -  dctx - the TSMonitorLGCtx object that contains all the options for the monitoring, this is created with TSMonitorLGCtxCreate()

6388:    Options Database:
6389: .   -ts_monitor_lg_solution_variables

6391:    Level: intermediate

6393:    Notes: Each process in a parallel run displays its component solutions in a separate window

6395: .keywords: TS,  vector, monitor, view

6397: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorLGCtxCreate(), TSMonitorLGCtxSetVariableNames(), TSMonitorLGCtxGetVariableNames(),
6398:            TSMonitorLGSetVariableNames(), TSMonitorLGGetVariableNames(), TSMonitorLGSetDisplayVariables(), TSMonitorLGCtxSetDisplayVariables(),
6399:            TSMonitorLGCtxSetTransform(), TSMonitorLGSetTransform(), TSMonitorLGError(), TSMonitorLGSNESIterations(), TSMonitorLGKSPIterations(),
6400:            TSMonitorEnvelopeCtxCreate(), TSMonitorEnvelopeGetBounds(), TSMonitorEnvelopeCtxDestroy(), TSMonitorEnvelop()
6401: @*/
6402: PetscErrorCode  TSMonitorLGSolution(TS ts,PetscInt step,PetscReal ptime,Vec u,void *dctx)
6403: {
6404:   PetscErrorCode    ierr;
6405:   TSMonitorLGCtx    ctx = (TSMonitorLGCtx)dctx;
6406:   const PetscScalar *yy;
6407:   Vec               v;

6410:   if (step < 0) return(0); /* -1 indicates interpolated solution */
6411:   if (!step) {
6412:     PetscDrawAxis axis;
6413:     PetscInt      dim;
6414:     PetscDrawLGGetAxis(ctx->lg,&axis);
6415:     PetscDrawAxisSetLabels(axis,"Solution as function of time","Time","Solution");
6416:     if (!ctx->names) {
6417:       PetscBool flg;
6418:       /* user provides names of variables to plot but no names has been set so assume names are integer values */
6419:       PetscOptionsHasName(((PetscObject)ts)->options,((PetscObject)ts)->prefix,"-ts_monitor_lg_solution_variables",&flg);
6420:       if (flg) {
6421:         PetscInt i,n;
6422:         char     **names;
6423:         VecGetSize(u,&n);
6424:         PetscMalloc1(n+1,&names);
6425:         for (i=0; i<n; i++) {
6426:           PetscMalloc1(5,&names[i]);
6427:           PetscSNPrintf(names[i],5,"%D",i);
6428:         }
6429:         names[n] = NULL;
6430:         ctx->names = names;
6431:       }
6432:     }
6433:     if (ctx->names && !ctx->displaynames) {
6434:       char      **displaynames;
6435:       PetscBool flg;
6436:       VecGetLocalSize(u,&dim);
6437:       PetscMalloc1(dim+1,&displaynames);
6438:       PetscMemzero(displaynames,(dim+1)*sizeof(char*));
6439:       PetscOptionsGetStringArray(((PetscObject)ts)->options,((PetscObject)ts)->prefix,"-ts_monitor_lg_solution_variables",displaynames,&dim,&flg);
6440:       if (flg) {
6441:         TSMonitorLGCtxSetDisplayVariables(ctx,(const char *const *)displaynames);
6442:       }
6443:       PetscStrArrayDestroy(&displaynames);
6444:     }
6445:     if (ctx->displaynames) {
6446:       PetscDrawLGSetDimension(ctx->lg,ctx->ndisplayvariables);
6447:       PetscDrawLGSetLegend(ctx->lg,(const char *const *)ctx->displaynames);
6448:     } else if (ctx->names) {
6449:       VecGetLocalSize(u,&dim);
6450:       PetscDrawLGSetDimension(ctx->lg,dim);
6451:       PetscDrawLGSetLegend(ctx->lg,(const char *const *)ctx->names);
6452:     } else {
6453:       VecGetLocalSize(u,&dim);
6454:       PetscDrawLGSetDimension(ctx->lg,dim);
6455:     }
6456:     PetscDrawLGReset(ctx->lg);
6457:   }

6459:   if (!ctx->transform) v = u;
6460:   else {(*ctx->transform)(ctx->transformctx,u,&v);}
6461:   VecGetArrayRead(v,&yy);
6462:   if (ctx->displaynames) {
6463:     PetscInt i;
6464:     for (i=0; i<ctx->ndisplayvariables; i++)
6465:       ctx->displayvalues[i] = PetscRealPart(yy[ctx->displayvariables[i]]);
6466:     PetscDrawLGAddCommonPoint(ctx->lg,ptime,ctx->displayvalues);
6467:   } else {
6468: #if defined(PETSC_USE_COMPLEX)
6469:     PetscInt  i,n;
6470:     PetscReal *yreal;
6471:     VecGetLocalSize(v,&n);
6472:     PetscMalloc1(n,&yreal);
6473:     for (i=0; i<n; i++) yreal[i] = PetscRealPart(yy[i]);
6474:     PetscDrawLGAddCommonPoint(ctx->lg,ptime,yreal);
6475:     PetscFree(yreal);
6476: #else
6477:     PetscDrawLGAddCommonPoint(ctx->lg,ptime,yy);
6478: #endif
6479:   }
6480:   VecRestoreArrayRead(v,&yy);
6481:   if (ctx->transform) {VecDestroy(&v);}

6483:   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
6484:     PetscDrawLGDraw(ctx->lg);
6485:     PetscDrawLGSave(ctx->lg);
6486:   }
6487:   return(0);
6488: }

6490: /*@C
6491:    TSMonitorLGSetVariableNames - Sets the name of each component in the solution vector so that it may be displayed in the plot

6493:    Collective on TS

6495:    Input Parameters:
6496: +  ts - the TS context
6497: -  names - the names of the components, final string must be NULL

6499:    Level: intermediate

6501:    Notes: If the TS object does not have a TSMonitorLGCtx associated with it then this function is ignored

6503: .keywords: TS,  vector, monitor, view

6505: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorLGSetDisplayVariables(), TSMonitorLGCtxSetVariableNames()
6506: @*/
6507: PetscErrorCode  TSMonitorLGSetVariableNames(TS ts,const char * const *names)
6508: {
6509:   PetscErrorCode    ierr;
6510:   PetscInt          i;

6513:   for (i=0; i<ts->numbermonitors; i++) {
6514:     if (ts->monitor[i] == TSMonitorLGSolution) {
6515:       TSMonitorLGCtxSetVariableNames((TSMonitorLGCtx)ts->monitorcontext[i],names);
6516:       break;
6517:     }
6518:   }
6519:   return(0);
6520: }

6522: /*@C
6523:    TSMonitorLGCtxSetVariableNames - Sets the name of each component in the solution vector so that it may be displayed in the plot

6525:    Collective on TS

6527:    Input Parameters:
6528: +  ts - the TS context
6529: -  names - the names of the components, final string must be NULL

6531:    Level: intermediate

6533: .keywords: TS,  vector, monitor, view

6535: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorLGSetDisplayVariables(), TSMonitorLGSetVariableNames()
6536: @*/
6537: PetscErrorCode  TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx ctx,const char * const *names)
6538: {
6539:   PetscErrorCode    ierr;

6542:   PetscStrArrayDestroy(&ctx->names);
6543:   PetscStrArrayallocpy(names,&ctx->names);
6544:   return(0);
6545: }

6547: /*@C
6548:    TSMonitorLGGetVariableNames - Gets the name of each component in the solution vector so that it may be displayed in the plot

6550:    Collective on TS

6552:    Input Parameter:
6553: .  ts - the TS context

6555:    Output Parameter:
6556: .  names - the names of the components, final string must be NULL

6558:    Level: intermediate

6560:    Notes: If the TS object does not have a TSMonitorLGCtx associated with it then this function is ignored

6562: .keywords: TS,  vector, monitor, view

6564: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorLGSetDisplayVariables()
6565: @*/
6566: PetscErrorCode  TSMonitorLGGetVariableNames(TS ts,const char *const **names)
6567: {
6568:   PetscInt       i;

6571:   *names = NULL;
6572:   for (i=0; i<ts->numbermonitors; i++) {
6573:     if (ts->monitor[i] == TSMonitorLGSolution) {
6574:       TSMonitorLGCtx  ctx = (TSMonitorLGCtx) ts->monitorcontext[i];
6575:       *names = (const char *const *)ctx->names;
6576:       break;
6577:     }
6578:   }
6579:   return(0);
6580: }

6582: /*@C
6583:    TSMonitorLGCtxSetDisplayVariables - Sets the variables that are to be display in the monitor

6585:    Collective on TS

6587:    Input Parameters:
6588: +  ctx - the TSMonitorLG context
6589: .  displaynames - the names of the components, final string must be NULL

6591:    Level: intermediate

6593: .keywords: TS,  vector, monitor, view

6595: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorLGSetVariableNames()
6596: @*/
6597: PetscErrorCode  TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx ctx,const char * const *displaynames)
6598: {
6599:   PetscInt          j = 0,k;
6600:   PetscErrorCode    ierr;

6603:   if (!ctx->names) return(0);
6604:   PetscStrArrayDestroy(&ctx->displaynames);
6605:   PetscStrArrayallocpy(displaynames,&ctx->displaynames);
6606:   while (displaynames[j]) j++;
6607:   ctx->ndisplayvariables = j;
6608:   PetscMalloc1(ctx->ndisplayvariables,&ctx->displayvariables);
6609:   PetscMalloc1(ctx->ndisplayvariables,&ctx->displayvalues);
6610:   j = 0;
6611:   while (displaynames[j]) {
6612:     k = 0;
6613:     while (ctx->names[k]) {
6614:       PetscBool flg;
6615:       PetscStrcmp(displaynames[j],ctx->names[k],&flg);
6616:       if (flg) {
6617:         ctx->displayvariables[j] = k;
6618:         break;
6619:       }
6620:       k++;
6621:     }
6622:     j++;
6623:   }
6624:   return(0);
6625: }

6627: /*@C
6628:    TSMonitorLGSetDisplayVariables - Sets the variables that are to be display in the monitor

6630:    Collective on TS

6632:    Input Parameters:
6633: +  ts - the TS context
6634: .  displaynames - the names of the components, final string must be NULL

6636:    Notes: If the TS object does not have a TSMonitorLGCtx associated with it then this function is ignored

6638:    Level: intermediate

6640: .keywords: TS,  vector, monitor, view

6642: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorLGSetVariableNames()
6643: @*/
6644: PetscErrorCode  TSMonitorLGSetDisplayVariables(TS ts,const char * const *displaynames)
6645: {
6646:   PetscInt          i;
6647:   PetscErrorCode    ierr;

6650:   for (i=0; i<ts->numbermonitors; i++) {
6651:     if (ts->monitor[i] == TSMonitorLGSolution) {
6652:       TSMonitorLGCtxSetDisplayVariables((TSMonitorLGCtx)ts->monitorcontext[i],displaynames);
6653:       break;
6654:     }
6655:   }
6656:   return(0);
6657: }

6659: /*@C
6660:    TSMonitorLGSetTransform - Solution vector will be transformed by provided function before being displayed

6662:    Collective on TS

6664:    Input Parameters:
6665: +  ts - the TS context
6666: .  transform - the transform function
6667: .  destroy - function to destroy the optional context
6668: -  ctx - optional context used by transform function

6670:    Notes: If the TS object does not have a TSMonitorLGCtx associated with it then this function is ignored

6672:    Level: intermediate

6674: .keywords: TS,  vector, monitor, view

6676: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorLGSetVariableNames(), TSMonitorLGCtxSetTransform()
6677: @*/
6678: PetscErrorCode  TSMonitorLGSetTransform(TS ts,PetscErrorCode (*transform)(void*,Vec,Vec*),PetscErrorCode (*destroy)(void*),void *tctx)
6679: {
6680:   PetscInt          i;
6681:   PetscErrorCode    ierr;

6684:   for (i=0; i<ts->numbermonitors; i++) {
6685:     if (ts->monitor[i] == TSMonitorLGSolution) {
6686:       TSMonitorLGCtxSetTransform((TSMonitorLGCtx)ts->monitorcontext[i],transform,destroy,tctx);
6687:     }
6688:   }
6689:   return(0);
6690: }

6692: /*@C
6693:    TSMonitorLGCtxSetTransform - Solution vector will be transformed by provided function before being displayed

6695:    Collective on TSLGCtx

6697:    Input Parameters:
6698: +  ts - the TS context
6699: .  transform - the transform function
6700: .  destroy - function to destroy the optional context
6701: -  ctx - optional context used by transform function

6703:    Level: intermediate

6705: .keywords: TS,  vector, monitor, view

6707: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorLGSetVariableNames(), TSMonitorLGSetTransform()
6708: @*/
6709: PetscErrorCode  TSMonitorLGCtxSetTransform(TSMonitorLGCtx ctx,PetscErrorCode (*transform)(void*,Vec,Vec*),PetscErrorCode (*destroy)(void*),void *tctx)
6710: {
6712:   ctx->transform    = transform;
6713:   ctx->transformdestroy = destroy;
6714:   ctx->transformctx = tctx;
6715:   return(0);
6716: }

6718: /*@C
6719:    TSMonitorLGError - Monitors progress of the TS solvers by plotting each component of the error
6720:        in a time based line graph

6722:    Collective on TS

6724:    Input Parameters:
6725: +  ts - the TS context
6726: .  step - current time-step
6727: .  ptime - current time
6728: .  u - current solution
6729: -  dctx - TSMonitorLGCtx object created with TSMonitorLGCtxCreate()

6731:    Level: intermediate

6733:    Notes: Each process in a parallel run displays its component errors in a separate window

6735:    The user must provide the solution using TSSetSolutionFunction() to use this monitor.

6737:    Options Database Keys:
6738: .  -ts_monitor_lg_error - create a graphical monitor of error history

6740: .keywords: TS,  vector, monitor, view

6742: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSSetSolutionFunction()
6743: @*/
6744: PetscErrorCode  TSMonitorLGError(TS ts,PetscInt step,PetscReal ptime,Vec u,void *dummy)
6745: {
6746:   PetscErrorCode    ierr;
6747:   TSMonitorLGCtx    ctx = (TSMonitorLGCtx)dummy;
6748:   const PetscScalar *yy;
6749:   Vec               y;

6752:   if (!step) {
6753:     PetscDrawAxis axis;
6754:     PetscInt      dim;
6755:     PetscDrawLGGetAxis(ctx->lg,&axis);
6756:     PetscDrawAxisSetLabels(axis,"Error in solution as function of time","Time","Error");
6757:     VecGetLocalSize(u,&dim);
6758:     PetscDrawLGSetDimension(ctx->lg,dim);
6759:     PetscDrawLGReset(ctx->lg);
6760:   }
6761:   VecDuplicate(u,&y);
6762:   TSComputeSolutionFunction(ts,ptime,y);
6763:   VecAXPY(y,-1.0,u);
6764:   VecGetArrayRead(y,&yy);
6765: #if defined(PETSC_USE_COMPLEX)
6766:   {
6767:     PetscReal *yreal;
6768:     PetscInt  i,n;
6769:     VecGetLocalSize(y,&n);
6770:     PetscMalloc1(n,&yreal);
6771:     for (i=0; i<n; i++) yreal[i] = PetscRealPart(yy[i]);
6772:     PetscDrawLGAddCommonPoint(ctx->lg,ptime,yreal);
6773:     PetscFree(yreal);
6774:   }
6775: #else
6776:   PetscDrawLGAddCommonPoint(ctx->lg,ptime,yy);
6777: #endif
6778:   VecRestoreArrayRead(y,&yy);
6779:   VecDestroy(&y);
6780:   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
6781:     PetscDrawLGDraw(ctx->lg);
6782:     PetscDrawLGSave(ctx->lg);
6783:   }
6784:   return(0);
6785: }

6787: /*@C
6788:    TSMonitorError - Monitors progress of the TS solvers by printing the 2 norm of the error at each timestep

6790:    Collective on TS

6792:    Input Parameters:
6793: +  ts - the TS context
6794: .  step - current time-step
6795: .  ptime - current time
6796: .  u - current solution
6797: -  dctx - unused context

6799:    Level: intermediate

6801:    The user must provide the solution using TSSetSolutionFunction() to use this monitor.

6803:    Options Database Keys:
6804: .  -ts_monitor_error - create a graphical monitor of error history

6806: .keywords: TS,  vector, monitor, view

6808: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSSetSolutionFunction()
6809: @*/
6810: PetscErrorCode  TSMonitorError(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscViewerAndFormat *vf)
6811: {
6812:   PetscErrorCode    ierr;
6813:   Vec               y;
6814:   PetscReal         nrm;
6815:   PetscBool         flg;

6818:   VecDuplicate(u,&y);
6819:   TSComputeSolutionFunction(ts,ptime,y);
6820:   VecAXPY(y,-1.0,u);
6821:   PetscObjectTypeCompare((PetscObject)vf->viewer,PETSCVIEWERASCII,&flg);
6822:   if (flg) {
6823:     VecNorm(y,NORM_2,&nrm);
6824:     PetscViewerASCIIPrintf(vf->viewer,"2-norm of error %g\n",(double)nrm);
6825:   }
6826:   PetscObjectTypeCompare((PetscObject)vf->viewer,PETSCVIEWERDRAW,&flg);
6827:   if (flg) {
6828:     VecView(y,vf->viewer);
6829:   }
6830:   VecDestroy(&y);
6831:   return(0);
6832: }

6834: PetscErrorCode TSMonitorLGSNESIterations(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx)
6835: {
6836:   TSMonitorLGCtx ctx = (TSMonitorLGCtx) monctx;
6837:   PetscReal      x   = ptime,y;
6839:   PetscInt       its;

6842:   if (n < 0) return(0); /* -1 indicates interpolated solution */
6843:   if (!n) {
6844:     PetscDrawAxis axis;
6845:     PetscDrawLGGetAxis(ctx->lg,&axis);
6846:     PetscDrawAxisSetLabels(axis,"Nonlinear iterations as function of time","Time","SNES Iterations");
6847:     PetscDrawLGReset(ctx->lg);
6848:     ctx->snes_its = 0;
6849:   }
6850:   TSGetSNESIterations(ts,&its);
6851:   y    = its - ctx->snes_its;
6852:   PetscDrawLGAddPoint(ctx->lg,&x,&y);
6853:   if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) {
6854:     PetscDrawLGDraw(ctx->lg);
6855:     PetscDrawLGSave(ctx->lg);
6856:   }
6857:   ctx->snes_its = its;
6858:   return(0);
6859: }

6861: PetscErrorCode TSMonitorLGKSPIterations(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx)
6862: {
6863:   TSMonitorLGCtx ctx = (TSMonitorLGCtx) monctx;
6864:   PetscReal      x   = ptime,y;
6866:   PetscInt       its;

6869:   if (n < 0) return(0); /* -1 indicates interpolated solution */
6870:   if (!n) {
6871:     PetscDrawAxis axis;
6872:     PetscDrawLGGetAxis(ctx->lg,&axis);
6873:     PetscDrawAxisSetLabels(axis,"Linear iterations as function of time","Time","KSP Iterations");
6874:     PetscDrawLGReset(ctx->lg);
6875:     ctx->ksp_its = 0;
6876:   }
6877:   TSGetKSPIterations(ts,&its);
6878:   y    = its - ctx->ksp_its;
6879:   PetscDrawLGAddPoint(ctx->lg,&x,&y);
6880:   if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) {
6881:     PetscDrawLGDraw(ctx->lg);
6882:     PetscDrawLGSave(ctx->lg);
6883:   }
6884:   ctx->ksp_its = its;
6885:   return(0);
6886: }

6888: /*@
6889:    TSComputeLinearStability - computes the linear stability function at a point

6891:    Collective on TS and Vec

6893:    Input Parameters:
6894: +  ts - the TS context
6895: -  xr,xi - real and imaginary part of input arguments

6897:    Output Parameters:
6898: .  yr,yi - real and imaginary part of function value

6900:    Level: developer

6902: .keywords: TS, compute

6904: .seealso: TSSetRHSFunction(), TSComputeIFunction()
6905: @*/
6906: PetscErrorCode TSComputeLinearStability(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
6907: {

6912:   if (!ts->ops->linearstability) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Linearized stability function not provided for this method");
6913:   (*ts->ops->linearstability)(ts,xr,xi,yr,yi);
6914:   return(0);
6915: }

6917: /* ------------------------------------------------------------------------*/
6918: /*@C
6919:    TSMonitorEnvelopeCtxCreate - Creates a context for use with TSMonitorEnvelope()

6921:    Collective on TS

6923:    Input Parameters:
6924: .  ts  - the ODE solver object

6926:    Output Parameter:
6927: .  ctx - the context

6929:    Level: intermediate

6931: .keywords: TS, monitor, line graph, residual, seealso

6933: .seealso: TSMonitorLGTimeStep(), TSMonitorSet(), TSMonitorLGSolution(), TSMonitorLGError()

6935: @*/
6936: PetscErrorCode  TSMonitorEnvelopeCtxCreate(TS ts,TSMonitorEnvelopeCtx *ctx)
6937: {

6941:   PetscNew(ctx);
6942:   return(0);
6943: }

6945: /*@C
6946:    TSMonitorEnvelope - Monitors the maximum and minimum value of each component of the solution

6948:    Collective on TS

6950:    Input Parameters:
6951: +  ts - the TS context
6952: .  step - current time-step
6953: .  ptime - current time
6954: .  u  - current solution
6955: -  dctx - the envelope context

6957:    Options Database:
6958: .  -ts_monitor_envelope

6960:    Level: intermediate

6962:    Notes: after a solve you can use TSMonitorEnvelopeGetBounds() to access the envelope

6964: .keywords: TS,  vector, monitor, view

6966: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorEnvelopeGetBounds(), TSMonitorEnvelopeCtxCreate()
6967: @*/
6968: PetscErrorCode  TSMonitorEnvelope(TS ts,PetscInt step,PetscReal ptime,Vec u,void *dctx)
6969: {
6970:   PetscErrorCode       ierr;
6971:   TSMonitorEnvelopeCtx ctx = (TSMonitorEnvelopeCtx)dctx;

6974:   if (!ctx->max) {
6975:     VecDuplicate(u,&ctx->max);
6976:     VecDuplicate(u,&ctx->min);
6977:     VecCopy(u,ctx->max);
6978:     VecCopy(u,ctx->min);
6979:   } else {
6980:     VecPointwiseMax(ctx->max,u,ctx->max);
6981:     VecPointwiseMin(ctx->min,u,ctx->min);
6982:   }
6983:   return(0);
6984: }

6986: /*@C
6987:    TSMonitorEnvelopeGetBounds - Gets the bounds for the components of the solution

6989:    Collective on TS

6991:    Input Parameter:
6992: .  ts - the TS context

6994:    Output Parameter:
6995: +  max - the maximum values
6996: -  min - the minimum values

6998:    Notes: If the TS does not have a TSMonitorEnvelopeCtx associated with it then this function is ignored

7000:    Level: intermediate

7002: .keywords: TS,  vector, monitor, view

7004: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorLGSetDisplayVariables()
7005: @*/
7006: PetscErrorCode  TSMonitorEnvelopeGetBounds(TS ts,Vec *max,Vec *min)
7007: {
7008:   PetscInt i;

7011:   if (max) *max = NULL;
7012:   if (min) *min = NULL;
7013:   for (i=0; i<ts->numbermonitors; i++) {
7014:     if (ts->monitor[i] == TSMonitorEnvelope) {
7015:       TSMonitorEnvelopeCtx  ctx = (TSMonitorEnvelopeCtx) ts->monitorcontext[i];
7016:       if (max) *max = ctx->max;
7017:       if (min) *min = ctx->min;
7018:       break;
7019:     }
7020:   }
7021:   return(0);
7022: }

7024: /*@C
7025:    TSMonitorEnvelopeCtxDestroy - Destroys a context that was created  with TSMonitorEnvelopeCtxCreate().

7027:    Collective on TSMonitorEnvelopeCtx

7029:    Input Parameter:
7030: .  ctx - the monitor context

7032:    Level: intermediate

7034: .keywords: TS, monitor, line graph, destroy

7036: .seealso: TSMonitorLGCtxCreate(),  TSMonitorSet(), TSMonitorLGTimeStep()
7037: @*/
7038: PetscErrorCode  TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx *ctx)
7039: {

7043:   VecDestroy(&(*ctx)->min);
7044:   VecDestroy(&(*ctx)->max);
7045:   PetscFree(*ctx);
7046:   return(0);
7047: }

7049: /*@
7050:    TSRestartStep - Flags the solver to restart the next step

7052:    Collective on TS

7054:    Input Parameter:
7055: .  ts - the TS context obtained from TSCreate()

7057:    Level: advanced

7059:    Notes:
7060:    Multistep methods like BDF or Runge-Kutta methods with FSAL property require restarting the solver in the event of
7061:    discontinuities. These discontinuities may be introduced as a consequence of explicitly modifications to the solution
7062:    vector (which PETSc attempts to detect and handle) or problem coefficients (which PETSc is not able to detect). For
7063:    the sake of correctness and maximum safety, users are expected to call TSRestart() whenever they introduce
7064:    discontinuities in callback routines (e.g. prestep and poststep routines, or implicit/rhs function routines with
7065:    discontinuous source terms).

7067: .keywords: TS, timestep, restart

7069: .seealso: TSSolve(), TSSetPreStep(), TSSetPostStep()
7070: @*/
7071: PetscErrorCode TSRestartStep(TS ts)
7072: {
7075:   ts->steprestart = PETSC_TRUE;
7076:   return(0);
7077: }

7079: /*@
7080:    TSRollBack - Rolls back one time step

7082:    Collective on TS

7084:    Input Parameter:
7085: .  ts - the TS context obtained from TSCreate()

7087:    Level: advanced

7089: .keywords: TS, timestep, rollback

7091: .seealso: TSCreate(), TSSetUp(), TSDestroy(), TSSolve(), TSSetPreStep(), TSSetPreStage(), TSInterpolate()
7092: @*/
7093: PetscErrorCode  TSRollBack(TS ts)
7094: {

7099:   if (ts->steprollback) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"TSRollBack already called");
7100:   if (!ts->ops->rollback) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRollBack not implemented for type '%s'",((PetscObject)ts)->type_name);
7101:   (*ts->ops->rollback)(ts);
7102:   ts->time_step = ts->ptime - ts->ptime_prev;
7103:   ts->ptime = ts->ptime_prev;
7104:   ts->ptime_prev = ts->ptime_prev_rollback;
7105:   ts->steps--;
7106:   ts->steprollback = PETSC_TRUE;
7107:   return(0);
7108: }

7110: /*@
7111:    TSGetStages - Get the number of stages and stage values

7113:    Input Parameter:
7114: .  ts - the TS context obtained from TSCreate()

7116:    Level: advanced

7118: .keywords: TS, getstages

7120: .seealso: TSCreate()
7121: @*/
7122: PetscErrorCode  TSGetStages(TS ts,PetscInt *ns,Vec **Y)
7123: {


7130:   if (!ts->ops->getstages) *ns=0;
7131:   else {
7132:     (*ts->ops->getstages)(ts,ns,Y);
7133:   }
7134:   return(0);
7135: }

7137: /*@C
7138:   TSComputeIJacobianDefaultColor - Computes the Jacobian using finite differences and coloring to exploit matrix sparsity.

7140:   Collective on SNES

7142:   Input Parameters:
7143: + ts - the TS context
7144: . t - current timestep
7145: . U - state vector
7146: . Udot - time derivative of state vector
7147: . shift - shift to apply, see note below
7148: - ctx - an optional user context

7150:   Output Parameters:
7151: + J - Jacobian matrix (not altered in this routine)
7152: - B - newly computed Jacobian matrix to use with preconditioner (generally the same as J)

7154:   Level: intermediate

7156:   Notes:
7157:   If F(t,U,Udot)=0 is the DAE, the required Jacobian is

7159:   dF/dU + shift*dF/dUdot

7161:   Most users should not need to explicitly call this routine, as it
7162:   is used internally within the nonlinear solvers.

7164:   This will first try to get the coloring from the DM.  If the DM type has no coloring
7165:   routine, then it will try to get the coloring from the matrix.  This requires that the
7166:   matrix have nonzero entries precomputed.

7168: .keywords: TS, finite differences, Jacobian, coloring, sparse
7169: .seealso: TSSetIJacobian(), MatFDColoringCreate(), MatFDColoringSetFunction()
7170: @*/
7171: PetscErrorCode TSComputeIJacobianDefaultColor(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal shift,Mat J,Mat B,void *ctx)
7172: {
7173:   SNES           snes;
7174:   MatFDColoring  color;
7175:   PetscBool      hascolor, matcolor = PETSC_FALSE;

7179:   PetscOptionsGetBool(((PetscObject)ts)->options,((PetscObject) ts)->prefix, "-ts_fd_color_use_mat", &matcolor, NULL);
7180:   PetscObjectQuery((PetscObject) B, "TSMatFDColoring", (PetscObject *) &color);
7181:   if (!color) {
7182:     DM         dm;
7183:     ISColoring iscoloring;

7185:     TSGetDM(ts, &dm);
7186:     DMHasColoring(dm, &hascolor);
7187:     if (hascolor && !matcolor) {
7188:       DMCreateColoring(dm, IS_COLORING_GLOBAL, &iscoloring);
7189:       MatFDColoringCreate(B, iscoloring, &color);
7190:       MatFDColoringSetFunction(color, (PetscErrorCode (*)(void)) SNESTSFormFunction, (void *) ts);
7191:       MatFDColoringSetFromOptions(color);
7192:       MatFDColoringSetUp(B, iscoloring, color);
7193:       ISColoringDestroy(&iscoloring);
7194:     } else {
7195:       MatColoring mc;

7197:       MatColoringCreate(B, &mc);
7198:       MatColoringSetDistance(mc, 2);
7199:       MatColoringSetType(mc, MATCOLORINGSL);
7200:       MatColoringSetFromOptions(mc);
7201:       MatColoringApply(mc, &iscoloring);
7202:       MatColoringDestroy(&mc);
7203:       MatFDColoringCreate(B, iscoloring, &color);
7204:       MatFDColoringSetFunction(color, (PetscErrorCode (*)(void)) SNESTSFormFunction, (void *) ts);
7205:       MatFDColoringSetFromOptions(color);
7206:       MatFDColoringSetUp(B, iscoloring, color);
7207:       ISColoringDestroy(&iscoloring);
7208:     }
7209:     PetscObjectCompose((PetscObject) B, "TSMatFDColoring", (PetscObject) color);
7210:     PetscObjectDereference((PetscObject) color);
7211:   }
7212:   TSGetSNES(ts, &snes);
7213:   MatFDColoringApply(B, color, U, snes);
7214:   if (J != B) {
7215:     MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);
7216:     MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);
7217:   }
7218:   return(0);
7219: }

7221: /*@
7222:     TSSetFunctionDomainError - Set the function testing if the current state vector is valid

7224:     Input Parameters:
7225:     ts - the TS context
7226:     func - function called within TSFunctionDomainError

7228:     Level: intermediate

7230: .keywords: TS, state, domain
7231: .seealso: TSAdaptCheckStage(), TSFunctionDomainError()
7232: @*/

7234: PetscErrorCode TSSetFunctionDomainError(TS ts, PetscErrorCode (*func)(TS,PetscReal,Vec,PetscBool*))
7235: {
7238:   ts->functiondomainerror = func;
7239:   return(0);
7240: }

7242: /*@
7243:     TSFunctionDomainError - Check if the current state is valid

7245:     Input Parameters:
7246:     ts - the TS context
7247:     stagetime - time of the simulation
7248:     Y - state vector to check.

7250:     Output Parameter:
7251:     accept - Set to PETSC_FALSE if the current state vector is valid.

7253:     Note:
7254:     This function should be used to ensure the state is in a valid part of the space.
7255:     For example, one can ensure here all values are positive.

7257:     Level: advanced
7258: @*/
7259: PetscErrorCode TSFunctionDomainError(TS ts,PetscReal stagetime,Vec Y,PetscBool* accept)
7260: {


7266:   *accept = PETSC_TRUE;
7267:   if (ts->functiondomainerror) {
7268:     PetscStackCallStandard((*ts->functiondomainerror),(ts,stagetime,Y,accept));
7269:   }
7270:   return(0);
7271: }

7273: /*@C
7274:   TSClone - This function clones a time step object.

7276:   Collective on MPI_Comm

7278:   Input Parameter:
7279: . tsin    - The input TS

7281:   Output Parameter:
7282: . tsout   - The output TS (cloned)

7284:   Notes:
7285:   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.

7287:   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);

7289:   Level: developer

7291: .keywords: TS, clone
7292: .seealso: TSCreate(), TSSetType(), TSSetUp(), TSDestroy(), TSSetProblemType()
7293: @*/
7294: PetscErrorCode  TSClone(TS tsin, TS *tsout)
7295: {
7296:   TS             t;
7298:   SNES           snes_start;
7299:   DM             dm;
7300:   TSType         type;

7304:   *tsout = NULL;

7306:   PetscHeaderCreate(t, TS_CLASSID, "TS", "Time stepping", "TS", PetscObjectComm((PetscObject)tsin), TSDestroy, TSView);

7308:   /* General TS description */
7309:   t->numbermonitors    = 0;
7310:   t->setupcalled       = 0;
7311:   t->ksp_its           = 0;
7312:   t->snes_its          = 0;
7313:   t->nwork             = 0;
7314:   t->rhsjacobian.time  = -1e20;
7315:   t->rhsjacobian.scale = 1.;
7316:   t->ijacobian.shift   = 1.;

7318:   TSGetSNES(tsin,&snes_start);
7319:   TSSetSNES(t,snes_start);

7321:   TSGetDM(tsin,&dm);
7322:   TSSetDM(t,dm);

7324:   t->adapt = tsin->adapt;
7325:   PetscObjectReference((PetscObject)t->adapt);

7327:   t->trajectory = tsin->trajectory;
7328:   PetscObjectReference((PetscObject)t->trajectory);

7330:   t->event = tsin->event;
7331:   if (t->event) t->event->refct++;

7333:   t->problem_type      = tsin->problem_type;
7334:   t->ptime             = tsin->ptime;
7335:   t->ptime_prev        = tsin->ptime_prev;
7336:   t->time_step         = tsin->time_step;
7337:   t->max_time          = tsin->max_time;
7338:   t->steps             = tsin->steps;
7339:   t->max_steps         = tsin->max_steps;
7340:   t->equation_type     = tsin->equation_type;
7341:   t->atol              = tsin->atol;
7342:   t->rtol              = tsin->rtol;
7343:   t->max_snes_failures = tsin->max_snes_failures;
7344:   t->max_reject        = tsin->max_reject;
7345:   t->errorifstepfailed = tsin->errorifstepfailed;

7347:   TSGetType(tsin,&type);
7348:   TSSetType(t,type);

7350:   t->vec_sol           = NULL;

7352:   t->cfltime          = tsin->cfltime;
7353:   t->cfltime_local    = tsin->cfltime_local;
7354:   t->exact_final_time = tsin->exact_final_time;

7356:   PetscMemcpy(t->ops,tsin->ops,sizeof(struct _TSOps));

7358:   if (((PetscObject)tsin)->fortran_func_pointers) {
7359:     PetscInt i;
7360:     PetscMalloc((10)*sizeof(void(*)(void)),&((PetscObject)t)->fortran_func_pointers);
7361:     for (i=0; i<10; i++) {
7362:       ((PetscObject)t)->fortran_func_pointers[i] = ((PetscObject)tsin)->fortran_func_pointers[i];
7363:     }
7364:   }
7365:   *tsout = t;
7366:   return(0);
7367: }

7369: static PetscErrorCode RHSWrapperFunction_TSRHSJacobianTest(void* ctx,Vec x,Vec y)
7370: {
7372:   TS             ts = (TS) ctx;

7375:   TSComputeRHSFunction(ts,0,x,y);
7376:   return(0);
7377: }

7379: /*@
7380:     TSRHSJacobianTest - Compares the multiply routine provided to the MATSHELL with differencing on the TS given RHS function.

7382:    Logically Collective on TS and Mat

7384:     Input Parameters:
7385:     TS - the time stepping routine

7387:    Output Parameter:
7388: .   flg - PETSC_TRUE if the multiply is likely correct

7390:    Options Database:
7391:  .   -ts_rhs_jacobian_test_mult -mat_shell_test_mult_view - run the test at each timestep of the integrator

7393:    Level: advanced

7395:    Notes: This only works for problems defined only the RHS function and Jacobian NOT IFunction and IJacobian

7397: .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMultTranspose(), TSRHSJacobianTestTranspose()
7398: @*/
7399: PetscErrorCode  TSRHSJacobianTest(TS ts,PetscBool *flg)
7400: {
7401:   Mat            J,B;
7403:   TSRHSJacobian  func;
7404:   void*          ctx;

7407:   TSGetRHSJacobian(ts,&J,&B,&func,&ctx);
7408:   (*func)(ts,0.0,ts->vec_sol,J,B,ctx);
7409:   MatShellTestMult(J,RHSWrapperFunction_TSRHSJacobianTest,ts->vec_sol,ts,flg);
7410:   return(0);
7411: }

7413: /*@C
7414:     TSRHSJacobianTestTranspose - Compares the multiply transpose routine provided to the MATSHELL with differencing on the TS given RHS function.

7416:    Logically Collective on TS and Mat

7418:     Input Parameters:
7419:     TS - the time stepping routine

7421:    Output Parameter:
7422: .   flg - PETSC_TRUE if the multiply is likely correct

7424:    Options Database:
7425: .   -ts_rhs_jacobian_test_mult_transpose -mat_shell_test_mult_transpose_view - run the test at each timestep of the integrator

7427:    Notes: This only works for problems defined only the RHS function and Jacobian NOT IFunction and IJacobian

7429:    Level: advanced

7431: .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMultTranspose(), TSRHSJacobianTest()
7432: @*/
7433: PetscErrorCode  TSRHSJacobianTestTranspose(TS ts,PetscBool *flg)
7434: {
7435:   Mat            J,B;
7437:   void           *ctx;
7438:   TSRHSJacobian  func;

7441:   TSGetRHSJacobian(ts,&J,&B,&func,&ctx);
7442:   (*func)(ts,0.0,ts->vec_sol,J,B,ctx);
7443:   MatShellTestMultTranspose(J,RHSWrapperFunction_TSRHSJacobianTest,ts->vec_sol,ts,flg);
7444:   return(0);
7445: }