Actual source code: ts.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  2: #include <petsc-private/tsimpl.h>        /*I "petscts.h"  I*/
  3: #include <petscdmshell.h>
  4: #include <petscdmda.h>
  5: #include <petscviewer.h>
  6: #include <petscdraw.h>

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

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

 16: /*
 17:   TSSetTypeFromOptions - Sets the type of ts from user options.

 19:   Collective on TS

 21:   Input Parameter:
 22: . ts - The ts

 24:   Level: intermediate

 26: .keywords: TS, set, options, database, type
 27: .seealso: TSSetFromOptions(), TSSetType()
 28: */
 29: static PetscErrorCode TSSetTypeFromOptions(TS ts)
 30: {
 31:   PetscBool      opt;
 32:   const char     *defaultType;
 33:   char           typeName[256];

 37:   if (((PetscObject)ts)->type_name) defaultType = ((PetscObject)ts)->type_name;
 38:   else defaultType = TSEULER;

 40:   if (!TSRegisterAllCalled) {TSRegisterAll();}
 41:   PetscOptionsFList("-ts_type", "TS method"," TSSetType", TSList, defaultType, typeName, 256, &opt);
 42:   if (opt) {
 43:     TSSetType(ts, typeName);
 44:   } else {
 45:     TSSetType(ts, defaultType);
 46:   }
 47:   return(0);
 48: }

 50: struct _n_TSMonitorDrawCtx {
 51:   PetscViewer   viewer;
 52:   PetscDrawAxis axis;
 53:   Vec           initialsolution;
 54:   PetscBool     showinitial;
 55:   PetscInt      howoften;  /* when > 0 uses step % howoften, when negative only final solution plotted */
 56:   PetscBool     showtimestepandtime;
 57:   int           color;
 58: };

 62: /*@
 63:    TSSetFromOptions - Sets various TS parameters from user options.

 65:    Collective on TS

 67:    Input Parameter:
 68: .  ts - the TS context obtained from TSCreate()

 70:    Options Database Keys:
 71: +  -ts_type <type> - TSEULER, TSBEULER, TSSUNDIALS, TSPSEUDO, TSCN, TSRK, TSTHETA, TSGL, TSSSP
 72: .  -ts_max_steps maxsteps - maximum number of time-steps to take
 73: .  -ts_final_time time - maximum time to compute to
 74: .  -ts_dt dt - initial time step
 75: .  -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
 76: .  -ts_max_snes_failures <maxfailures> - Maximum number of nonlinear solve failures allowed
 77: .  -ts_max_reject <maxrejects> - Maximum number of step rejections before step fails
 78: .  -ts_error_if_step_fails <true,false> - Error if no step succeeds
 79: .  -ts_rtol <rtol> - relative tolerance for local truncation error
 80: .  -ts_atol <atol> Absolute tolerance for local truncation error
 81: .  -ts_monitor - print information at each timestep
 82: .  -ts_monitor_lg_timestep - Monitor timestep size graphically
 83: .  -ts_monitor_lg_solution - Monitor solution graphically
 84: .  -ts_monitor_lg_error - Monitor error graphically
 85: .  -ts_monitor_lg_snes_iterations - Monitor number nonlinear iterations for each timestep graphically
 86: .  -ts_monitor_lg_ksp_iterations - Monitor number nonlinear iterations for each timestep graphically
 87: .  -ts_monitor_sp_eig - Monitor eigenvalues of linearized operator graphically
 88: .  -ts_monitor_draw_solution - Monitor solution graphically
 89: .  -ts_monitor_draw_solution_phase - Monitor solution graphically with phase diagram
 90: .  -ts_monitor_draw_error - Monitor error graphically
 91: .  -ts_monitor_solution_binary <filename> - Save each solution to a binary file
 92: -  -ts_monitor_solution_vtk <filename.vts> - Save each time step to a binary file, use filename-%%03D.vts

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

 96:    Level: beginner

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

100: .seealso: TSGetType()
101: @*/
102: PetscErrorCode  TSSetFromOptions(TS ts)
103: {
104:   PetscBool              opt,flg;
105:   PetscErrorCode         ierr;
106:   PetscViewer            monviewer;
107:   char                   monfilename[PETSC_MAX_PATH_LEN];
108:   SNES                   snes;
109:   TSAdapt                adapt;
110:   PetscReal              time_step;
111:   TSExactFinalTimeOption eftopt;
112:   char                   dir[16];

116:   PetscObjectOptionsBegin((PetscObject)ts);
117:   /* Handle TS type options */
118:   TSSetTypeFromOptions(ts);

120:   /* Handle generic TS options */
121:   PetscOptionsInt("-ts_max_steps","Maximum number of time steps","TSSetDuration",ts->max_steps,&ts->max_steps,NULL);
122:   PetscOptionsReal("-ts_final_time","Time to run to","TSSetDuration",ts->max_time,&ts->max_time,NULL);
123:   PetscOptionsReal("-ts_init_time","Initial time","TSSetTime",ts->ptime,&ts->ptime,NULL);
124:   PetscOptionsReal("-ts_dt","Initial time step","TSSetTimeStep",ts->time_step,&time_step,&flg);
125:   if (flg) {
126:     TSSetTimeStep(ts,time_step);
127:   }
128:   PetscOptionsEnum("-ts_exact_final_time","Option for handling of final time step","TSSetExactFinalTime",TSExactFinalTimeOptions,(PetscEnum)ts->exact_final_time,(PetscEnum*)&eftopt,&flg);
129:   if (flg) {TSSetExactFinalTime(ts,eftopt);}
130:   PetscOptionsInt("-ts_max_snes_failures","Maximum number of nonlinear solve failures","TSSetMaxSNESFailures",ts->max_snes_failures,&ts->max_snes_failures,NULL);
131:   PetscOptionsInt("-ts_max_reject","Maximum number of step rejections before step fails","TSSetMaxStepRejections",ts->max_reject,&ts->max_reject,NULL);
132:   PetscOptionsBool("-ts_error_if_step_fails","Error if no step succeeds","TSSetErrorIfStepFails",ts->errorifstepfailed,&ts->errorifstepfailed,NULL);
133:   PetscOptionsReal("-ts_rtol","Relative tolerance for local truncation error","TSSetTolerances",ts->rtol,&ts->rtol,NULL);
134:   PetscOptionsReal("-ts_atol","Absolute tolerance for local truncation error","TSSetTolerances",ts->atol,&ts->atol,NULL);

136: #if defined(PETSC_HAVE_SAWS)
137:   {
138:   PetscBool set;
139:   flg  = PETSC_FALSE;
140:   PetscOptionsBool("-ts_saws_block","Block for SAWs memory snooper at end of TSSolve","PetscObjectSAWsBlock",((PetscObject)ts)->amspublishblock,&flg,&set);
141:   if (set) {
142:     PetscObjectSAWsSetBlock((PetscObject)ts,flg);
143:   }
144:   }
145: #endif

147:   /* Monitor options */
148:   PetscOptionsString("-ts_monitor","Monitor timestep size","TSMonitorDefault","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
149:   if (flg) {
150:     PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ts),monfilename,&monviewer);
151:     TSMonitorSet(ts,TSMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);
152:   }
153:   PetscOptionsString("-ts_monitor_python","Use Python function","TSMonitorSet",0,monfilename,PETSC_MAX_PATH_LEN,&flg);
154:   if (flg) {PetscPythonMonitorSet((PetscObject)ts,monfilename);}

156:   PetscOptionsName("-ts_monitor_lg_timestep","Monitor timestep size graphically","TSMonitorLGTimeStep",&opt);
157:   if (opt) {
158:     TSMonitorLGCtx ctx;
159:     PetscInt       howoften = 1;

161:     PetscOptionsInt("-ts_monitor_lg_timestep","Monitor timestep size graphically","TSMonitorLGTimeStep",howoften,&howoften,NULL);
162:     TSMonitorLGCtxCreate(PetscObjectComm((PetscObject)ts),0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx);
163:     TSMonitorSet(ts,TSMonitorLGTimeStep,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);
164:   }
165:   PetscOptionsName("-ts_monitor_lg_solution","Monitor solution graphically","TSMonitorLGSolution",&opt);
166:   if (opt) {
167:     TSMonitorLGCtx ctx;
168:     PetscInt       howoften = 1;

170:     PetscOptionsInt("-ts_monitor_lg_solution","Monitor solution graphically","TSMonitorLGSolution",howoften,&howoften,NULL);
171:     TSMonitorLGCtxCreate(PETSC_COMM_SELF,0,0,PETSC_DECIDE,PETSC_DECIDE,600,400,howoften,&ctx);
172:     TSMonitorSet(ts,TSMonitorLGSolution,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);
173:   }
174:   PetscOptionsName("-ts_monitor_lg_error","Monitor error graphically","TSMonitorLGError",&opt);
175:   if (opt) {
176:     TSMonitorLGCtx ctx;
177:     PetscInt       howoften = 1;

179:     PetscOptionsInt("-ts_monitor_lg_error","Monitor error graphically","TSMonitorLGError",howoften,&howoften,NULL);
180:     TSMonitorLGCtxCreate(PETSC_COMM_SELF,0,0,PETSC_DECIDE,PETSC_DECIDE,600,400,howoften,&ctx);
181:     TSMonitorSet(ts,TSMonitorLGError,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);
182:   }
183:   PetscOptionsName("-ts_monitor_lg_snes_iterations","Monitor number nonlinear iterations for each timestep graphically","TSMonitorLGSNESIterations",&opt);
184:   if (opt) {
185:     TSMonitorLGCtx ctx;
186:     PetscInt       howoften = 1;

188:     PetscOptionsInt("-ts_monitor_lg_snes_iterations","Monitor number nonlinear iterations for each timestep graphically","TSMonitorLGSNESIterations",howoften,&howoften,NULL);
189:     TSMonitorLGCtxCreate(PETSC_COMM_SELF,0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx);
190:     TSMonitorSet(ts,TSMonitorLGSNESIterations,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);
191:   }
192:   PetscOptionsName("-ts_monitor_lg_ksp_iterations","Monitor number nonlinear iterations for each timestep graphically","TSMonitorLGKSPIterations",&opt);
193:   if (opt) {
194:     TSMonitorLGCtx ctx;
195:     PetscInt       howoften = 1;

197:     PetscOptionsInt("-ts_monitor_lg_ksp_iterations","Monitor number nonlinear iterations for each timestep graphically","TSMonitorLGKSPIterations",howoften,&howoften,NULL);
198:     TSMonitorLGCtxCreate(PETSC_COMM_SELF,0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx);
199:     TSMonitorSet(ts,TSMonitorLGKSPIterations,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);
200:   }
201:   PetscOptionsName("-ts_monitor_sp_eig","Monitor eigenvalues of linearized operator graphically","TSMonitorSPEig",&opt);
202:   if (opt) {
203:     TSMonitorSPEigCtx ctx;
204:     PetscInt          howoften = 1;

206:     PetscOptionsInt("-ts_monitor_sp_eig","Monitor eigenvalues of linearized operator graphically","TSMonitorSPEig",howoften,&howoften,NULL);
207:     TSMonitorSPEigCtxCreate(PETSC_COMM_SELF,0,0,PETSC_DECIDE,PETSC_DECIDE,600,400,howoften,&ctx);
208:     TSMonitorSet(ts,TSMonitorSPEig,ctx,(PetscErrorCode (*)(void**))TSMonitorSPEigCtxDestroy);
209:   }
210:   opt  = PETSC_FALSE;
211:   PetscOptionsName("-ts_monitor_draw_solution","Monitor solution graphically","TSMonitorDrawSolution",&opt);
212:   if (opt) {
213:     TSMonitorDrawCtx ctx;
214:     PetscInt         howoften = 1;

216:     PetscOptionsInt("-ts_monitor_draw_solution","Monitor solution graphically","TSMonitorDrawSolution",howoften,&howoften,NULL);
217:     TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),0,0,PETSC_DECIDE,PETSC_DECIDE,600,400,howoften,&ctx);
218:     TSMonitorSet(ts,TSMonitorDrawSolution,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);
219:   }
220:   opt  = PETSC_FALSE;
221:   PetscOptionsName("-ts_monitor_draw_solution_phase","Monitor solution graphically","TSMonitorDrawSolutionPhase",&opt);
222:   if (opt) {
223:     TSMonitorDrawCtx ctx;
224:     PetscReal        bounds[4];
225:     PetscInt         n = 4;
226:     PetscDraw        draw;

228:     PetscOptionsRealArray("-ts_monitor_draw_solution_phase","Monitor solution graphically","TSMonitorDrawSolutionPhase",bounds,&n,NULL);
229:     if (n != 4) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONG,"Must provide bounding box of phase field");
230:     TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),0,0,PETSC_DECIDE,PETSC_DECIDE,600,400,1,&ctx);
231:     PetscViewerDrawGetDraw(ctx->viewer,0,&draw);
232:     PetscDrawClear(draw);
233:     PetscDrawAxisCreate(draw,&ctx->axis);
234:     PetscDrawAxisSetLimits(ctx->axis,bounds[0],bounds[2],bounds[1],bounds[3]);
235:     PetscDrawAxisSetLabels(ctx->axis,"Phase Diagram","Variable 1","Variable 2");
236:     PetscDrawAxisDraw(ctx->axis);
237:     /* PetscDrawSetCoordinates(draw,bounds[0],bounds[1],bounds[2],bounds[3]); */
238:     TSMonitorSet(ts,TSMonitorDrawSolutionPhase,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);
239:   }
240:   opt  = PETSC_FALSE;
241:   PetscOptionsName("-ts_monitor_draw_error","Monitor error graphically","TSMonitorDrawError",&opt);
242:   if (opt) {
243:     TSMonitorDrawCtx ctx;
244:     PetscInt         howoften = 1;

246:     PetscOptionsInt("-ts_monitor_draw_error","Monitor error graphically","TSMonitorDrawError",howoften,&howoften,NULL);
247:     TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),0,0,PETSC_DECIDE,PETSC_DECIDE,600,400,howoften,&ctx);
248:     TSMonitorSet(ts,TSMonitorDrawError,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);
249:   }
250:   opt  = PETSC_FALSE;
251:   PetscOptionsString("-ts_monitor_solution_binary","Save each solution to a binary file","TSMonitorSolutionBinary",0,monfilename,PETSC_MAX_PATH_LEN,&flg);
252:   if (flg) {
253:     PetscViewer ctx;
254:     if (monfilename[0]) {
255:       PetscViewerBinaryOpen(PetscObjectComm((PetscObject)ts),monfilename,FILE_MODE_WRITE,&ctx);
256:       TSMonitorSet(ts,TSMonitorSolutionBinary,ctx,(PetscErrorCode (*)(void**))PetscViewerDestroy);
257:     } else {
258:       ctx = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)ts));
259:       TSMonitorSet(ts,TSMonitorSolutionBinary,ctx,(PetscErrorCode (*)(void**))NULL);
260:     }
261:   }
262:   opt  = PETSC_FALSE;
263:   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);
264:   if (flg) {
265:     const char *ptr,*ptr2;
266:     char       *filetemplate;
267:     if (!monfilename[0]) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"-ts_monitor_solution_vtk requires a file template, e.g. filename-%%03D.vts");
268:     /* Do some cursory validation of the input. */
269:     PetscStrstr(monfilename,"%",(char**)&ptr);
270:     if (!ptr) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"-ts_monitor_solution_vtk requires a file template, e.g. filename-%%03D.vts");
271:     for (ptr++; ptr && *ptr; ptr++) {
272:       PetscStrchr("DdiouxX",*ptr,(char**)&ptr2);
273:       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");
274:       if (ptr2) break;
275:     }
276:     PetscStrallocpy(monfilename,&filetemplate);
277:     TSMonitorSet(ts,TSMonitorSolutionVTK,filetemplate,(PetscErrorCode (*)(void**))TSMonitorSolutionVTKDestroy);
278:   }

280:   PetscOptionsString("-ts_monitor_dmda_ray","Display a ray of the solution","None","y=0",dir,16,&flg);
281:   if (flg) {
282:     TSMonitorDMDARayCtx *rayctx;
283:     int                  ray = 0;
284:     DMDADirection        ddir;
285:     DM                   da;
286:     PetscMPIInt          rank;

288:     if (dir[1] != '=') SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONG,"Unknown ray %s",dir);
289:     if (dir[0] == 'x') ddir = DMDA_X;
290:     else if (dir[0] == 'y') ddir = DMDA_Y;
291:     else SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONG,"Unknown ray %s",dir);
292:     sscanf(dir+2,"%d",&ray);

294:     PetscInfo2(((PetscObject)ts),"Displaying DMDA ray %c = %D\n",dir[0],ray);
295:     PetscNew(&rayctx);
296:     TSGetDM(ts,&da);
297:     DMDAGetRay(da,ddir,ray,&rayctx->ray,&rayctx->scatter);
298:     MPI_Comm_rank(PetscObjectComm((PetscObject)ts),&rank);
299:     if (!rank) {
300:       PetscViewerDrawOpen(PETSC_COMM_SELF,0,0,0,0,600,300,&rayctx->viewer);
301:     }
302:     rayctx->lgctx = NULL;
303:     TSMonitorSet(ts,TSMonitorDMDARay,rayctx,TSMonitorDMDARayDestroy);
304:   }
305:   PetscOptionsString("-ts_monitor_lg_dmda_ray","Display a ray of the solution","None","x=0",dir,16,&flg);
306:   if (flg) {
307:     TSMonitorDMDARayCtx *rayctx;
308:     int                 ray = 0;
309:     DMDADirection       ddir;
310:     DM                  da;
311:     PetscInt            howoften = 1;

313:     if (dir[1] != '=') SETERRQ1(PetscObjectComm((PetscObject) ts), PETSC_ERR_ARG_WRONG, "Malformed ray %s", dir);
314:     if      (dir[0] == 'x') ddir = DMDA_X;
315:     else if (dir[0] == 'y') ddir = DMDA_Y;
316:     else SETERRQ1(PetscObjectComm((PetscObject) ts), PETSC_ERR_ARG_WRONG, "Unknown ray direction %s", dir);
317:     sscanf(dir+2, "%d", &ray);

319:     PetscInfo2(((PetscObject) ts),"Displaying LG DMDA ray %c = %D\n", dir[0], ray);
320:     PetscNew(&rayctx);
321:     TSGetDM(ts, &da);
322:     DMDAGetRay(da, ddir, ray, &rayctx->ray, &rayctx->scatter);
323:     TSMonitorLGCtxCreate(PETSC_COMM_SELF,0,0,PETSC_DECIDE,PETSC_DECIDE,600,400,howoften,&rayctx->lgctx);
324:     TSMonitorSet(ts, TSMonitorLGDMDARay, rayctx, TSMonitorDMDARayDestroy);
325:   }

327:   /*
328:      This code is all wrong. One is creating objects inside the TSSetFromOptions() so if run with the options gui
329:      will bleed memory. Also one is using a PetscOptionsBegin() inside a PetscOptionsBegin()
330:   */
331:   TSGetAdapt(ts,&adapt);
332:   TSAdaptSetFromOptions(adapt);

334:   TSGetSNES(ts,&snes);
335:   if (ts->problem_type == TS_LINEAR) {SNESSetType(snes,SNESKSPONLY);}

337:   /* Handle specific TS options */
338:   if (ts->ops->setfromoptions) {
339:     (*ts->ops->setfromoptions)(ts);
340:   }

342:   /* process any options handlers added with PetscObjectAddOptionsHandler() */
343:   PetscObjectProcessOptionsHandlers((PetscObject)ts);
344:   PetscOptionsEnd();
345:   return(0);
346: }

351: /*@
352:    TSComputeRHSJacobian - Computes the Jacobian matrix that has been
353:       set with TSSetRHSJacobian().

355:    Collective on TS and Vec

357:    Input Parameters:
358: +  ts - the TS context
359: .  t - current timestep
360: -  U - input vector

362:    Output Parameters:
363: +  A - Jacobian matrix
364: .  B - optional preconditioning matrix
365: -  flag - flag indicating matrix structure

367:    Notes:
368:    Most users should not need to explicitly call this routine, as it
369:    is used internally within the nonlinear solvers.

371:    See KSPSetOperators() for important information about setting the
372:    flag parameter.

374:    Level: developer

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

378: .seealso:  TSSetRHSJacobian(), KSPSetOperators()
379: @*/
380: PetscErrorCode  TSComputeRHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B)
381: {
383:   PetscObjectState Ustate;
384:   DM             dm;
385:   DMTS           tsdm;
386:   TSRHSJacobian  rhsjacobianfunc;
387:   void           *ctx;
388:   TSIJacobian    ijacobianfunc;

394:   TSGetDM(ts,&dm);
395:   DMGetDMTS(dm,&tsdm);
396:   DMTSGetRHSJacobian(dm,&rhsjacobianfunc,&ctx);
397:   DMTSGetIJacobian(dm,&ijacobianfunc,NULL);
398:   PetscObjectStateGet((PetscObject)U,&Ustate);
399:   if (ts->rhsjacobian.time == t && (ts->problem_type == TS_LINEAR || (ts->rhsjacobian.X == U && ts->rhsjacobian.Xstate == Ustate))) {
400:     return(0);
401:   }

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

405:   if (ts->rhsjacobian.reuse) {
406:     MatShift(A,-ts->rhsjacobian.shift);
407:     MatScale(A,1./ts->rhsjacobian.scale);
408:     if (A != B) {
409:       MatShift(B,-ts->rhsjacobian.shift);
410:       MatScale(B,1./ts->rhsjacobian.scale);
411:     }
412:     ts->rhsjacobian.shift = 0;
413:     ts->rhsjacobian.scale = 1.;
414:   }

416:   if (rhsjacobianfunc) {
417:     PetscLogEventBegin(TS_JacobianEval,ts,U,A,B);
418:     PetscStackPush("TS user Jacobian function");
419:     (*rhsjacobianfunc)(ts,t,U,A,B,ctx);
420:     PetscStackPop;
421:     PetscLogEventEnd(TS_JacobianEval,ts,U,A,B);
422:     /* make sure user returned a correct Jacobian and preconditioner */
425:   } else {
426:     MatZeroEntries(A);
427:     if (A != B) {MatZeroEntries(B);}
428:   }
429:   ts->rhsjacobian.time       = t;
430:   ts->rhsjacobian.X          = U;
431:   PetscObjectStateGet((PetscObject)U,&ts->rhsjacobian.Xstate);
432:   return(0);
433: }

437: /*@
438:    TSComputeRHSFunction - Evaluates the right-hand-side function.

440:    Collective on TS and Vec

442:    Input Parameters:
443: +  ts - the TS context
444: .  t - current time
445: -  U - state vector

447:    Output Parameter:
448: .  y - right hand side

450:    Note:
451:    Most users should not need to explicitly call this routine, as it
452:    is used internally within the nonlinear solvers.

454:    Level: developer

456: .keywords: TS, compute

458: .seealso: TSSetRHSFunction(), TSComputeIFunction()
459: @*/
460: PetscErrorCode TSComputeRHSFunction(TS ts,PetscReal t,Vec U,Vec y)
461: {
463:   TSRHSFunction  rhsfunction;
464:   TSIFunction    ifunction;
465:   void           *ctx;
466:   DM             dm;

472:   TSGetDM(ts,&dm);
473:   DMTSGetRHSFunction(dm,&rhsfunction,&ctx);
474:   DMTSGetIFunction(dm,&ifunction,NULL);

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

478:   PetscLogEventBegin(TS_FunctionEval,ts,U,y,0);
479:   if (rhsfunction) {
480:     PetscStackPush("TS user right-hand-side function");
481:     (*rhsfunction)(ts,t,U,y,ctx);
482:     PetscStackPop;
483:   } else {
484:     VecZeroEntries(y);
485:   }

487:   PetscLogEventEnd(TS_FunctionEval,ts,U,y,0);
488:   return(0);
489: }

493: /*@
494:    TSComputeSolutionFunction - Evaluates the solution function.

496:    Collective on TS and Vec

498:    Input Parameters:
499: +  ts - the TS context
500: -  t - current time

502:    Output Parameter:
503: .  U - the solution

505:    Note:
506:    Most users should not need to explicitly call this routine, as it
507:    is used internally within the nonlinear solvers.

509:    Level: developer

511: .keywords: TS, compute

513: .seealso: TSSetSolutionFunction(), TSSetRHSFunction(), TSComputeIFunction()
514: @*/
515: PetscErrorCode TSComputeSolutionFunction(TS ts,PetscReal t,Vec U)
516: {
517:   PetscErrorCode     ierr;
518:   TSSolutionFunction solutionfunction;
519:   void               *ctx;
520:   DM                 dm;

525:   TSGetDM(ts,&dm);
526:   DMTSGetSolutionFunction(dm,&solutionfunction,&ctx);

528:   if (solutionfunction) {
529:     PetscStackPush("TS user solution function");
530:     (*solutionfunction)(ts,t,U,ctx);
531:     PetscStackPop;
532:   }
533:   return(0);
534: }
537: /*@
538:    TSComputeForcingFunction - Evaluates the forcing function.

540:    Collective on TS and Vec

542:    Input Parameters:
543: +  ts - the TS context
544: -  t - current time

546:    Output Parameter:
547: .  U - the function value

549:    Note:
550:    Most users should not need to explicitly call this routine, as it
551:    is used internally within the nonlinear solvers.

553:    Level: developer

555: .keywords: TS, compute

557: .seealso: TSSetSolutionFunction(), TSSetRHSFunction(), TSComputeIFunction()
558: @*/
559: PetscErrorCode TSComputeForcingFunction(TS ts,PetscReal t,Vec U)
560: {
561:   PetscErrorCode     ierr, (*forcing)(TS,PetscReal,Vec,void*);
562:   void               *ctx;
563:   DM                 dm;

568:   TSGetDM(ts,&dm);
569:   DMTSGetForcingFunction(dm,&forcing,&ctx);

571:   if (forcing) {
572:     PetscStackPush("TS user forcing function");
573:     (*forcing)(ts,t,U,ctx);
574:     PetscStackPop;
575:   }
576:   return(0);
577: }

581: static PetscErrorCode TSGetRHSVec_Private(TS ts,Vec *Frhs)
582: {
583:   Vec            F;

587:   *Frhs = NULL;
588:   TSGetIFunction(ts,&F,NULL,NULL);
589:   if (!ts->Frhs) {
590:     VecDuplicate(F,&ts->Frhs);
591:   }
592:   *Frhs = ts->Frhs;
593:   return(0);
594: }

598: static PetscErrorCode TSGetRHSMats_Private(TS ts,Mat *Arhs,Mat *Brhs)
599: {
600:   Mat            A,B;

604:   if (Arhs) *Arhs = NULL;
605:   if (Brhs) *Brhs = NULL;
606:   TSGetIJacobian(ts,&A,&B,NULL,NULL);
607:   if (Arhs) {
608:     if (!ts->Arhs) {
609:       MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&ts->Arhs);
610:     }
611:     *Arhs = ts->Arhs;
612:   }
613:   if (Brhs) {
614:     if (!ts->Brhs) {
615:       if (A != B) {
616:         MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&ts->Brhs);
617:       } else {
618:         ts->Brhs = ts->Arhs;
619:         PetscObjectReference((PetscObject)ts->Arhs);
620:       }
621:     }
622:     *Brhs = ts->Brhs;
623:   }
624:   return(0);
625: }

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

632:    Collective on TS and Vec

634:    Input Parameters:
635: +  ts - the TS context
636: .  t - current time
637: .  U - state vector
638: .  Udot - time derivative of state vector
639: -  imex - flag indicates if the method is IMEX so that the RHSFunction should be kept separate

641:    Output Parameter:
642: .  Y - right hand side

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

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

651:    Level: developer

653: .keywords: TS, compute

655: .seealso: TSSetIFunction(), TSComputeRHSFunction()
656: @*/
657: PetscErrorCode TSComputeIFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec Y,PetscBool imex)
658: {
660:   TSIFunction    ifunction;
661:   TSRHSFunction  rhsfunction;
662:   void           *ctx;
663:   DM             dm;


671:   TSGetDM(ts,&dm);
672:   DMTSGetIFunction(dm,&ifunction,&ctx);
673:   DMTSGetRHSFunction(dm,&rhsfunction,NULL);

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

677:   PetscLogEventBegin(TS_FunctionEval,ts,U,Udot,Y);
678:   if (ifunction) {
679:     PetscStackPush("TS user implicit function");
680:     (*ifunction)(ts,t,U,Udot,Y,ctx);
681:     PetscStackPop;
682:   }
683:   if (imex) {
684:     if (!ifunction) {
685:       VecCopy(Udot,Y);
686:     }
687:   } else if (rhsfunction) {
688:     if (ifunction) {
689:       Vec Frhs;
690:       TSGetRHSVec_Private(ts,&Frhs);
691:       TSComputeRHSFunction(ts,t,U,Frhs);
692:       VecAXPY(Y,-1,Frhs);
693:     } else {
694:       TSComputeRHSFunction(ts,t,U,Y);
695:       VecAYPX(Y,-1,Udot);
696:     }
697:   }
698:   PetscLogEventEnd(TS_FunctionEval,ts,U,Udot,Y);
699:   return(0);
700: }

704: /*@
705:    TSComputeIJacobian - Evaluates the Jacobian of the DAE

707:    Collective on TS and Vec

709:    Input
710:       Input Parameters:
711: +  ts - the TS context
712: .  t - current timestep
713: .  U - state vector
714: .  Udot - time derivative of state vector
715: .  shift - shift to apply, see note below
716: -  imex - flag indicates if the method is IMEX so that the RHSJacobian should be kept separate

718:    Output Parameters:
719: +  A - Jacobian matrix
720: .  B - optional preconditioning matrix
721: -  flag - flag indicating matrix structure

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

726:    dF/dU + shift*dF/dUdot

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

731:    Level: developer

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

735: .seealso:  TSSetIJacobian()
736: @*/
737: PetscErrorCode TSComputeIJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal shift,Mat A,Mat B,PetscBool imex)
738: {
740:   TSIJacobian    ijacobian;
741:   TSRHSJacobian  rhsjacobian;
742:   DM             dm;
743:   void           *ctx;


754:   TSGetDM(ts,&dm);
755:   DMTSGetIJacobian(dm,&ijacobian,&ctx);
756:   DMTSGetRHSJacobian(dm,&rhsjacobian,NULL);

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

760:   PetscLogEventBegin(TS_JacobianEval,ts,U,A,B);
761:   if (ijacobian) {
762:     PetscStackPush("TS user implicit Jacobian");
763:     (*ijacobian)(ts,t,U,Udot,shift,A,B,ctx);
764:     PetscStackPop;
765:     /* make sure user returned a correct Jacobian and preconditioner */
768:   }
769:   if (imex) {
770:     if (!ijacobian) {  /* system was written as Udot = G(t,U) */
771:       MatZeroEntries(A);
772:       MatShift(A,shift);
773:       if (A != B) {
774:         MatZeroEntries(B);
775:         MatShift(B,shift);
776:       }
777:     }
778:   } else {
779:     Mat Arhs = NULL,Brhs = NULL;
780:     if (rhsjacobian) {
781:       if (ijacobian) {
782:         TSGetRHSMats_Private(ts,&Arhs,&Brhs);
783:       } else {
784:         TSGetIJacobian(ts,&Arhs,&Brhs,NULL,NULL);
785:       }
786:       TSComputeRHSJacobian(ts,t,U,Arhs,Brhs);
787:     }
788:     if (Arhs == A) {           /* No IJacobian, so we only have the RHS matrix */
789:       ts->rhsjacobian.scale = -1;
790:       ts->rhsjacobian.shift = shift;
791:       MatScale(A,-1);
792:       MatShift(A,shift);
793:       if (A != B) {
794:         MatScale(B,-1);
795:         MatShift(B,shift);
796:       }
797:     } else if (Arhs) {          /* Both IJacobian and RHSJacobian */
798:       MatStructure axpy = DIFFERENT_NONZERO_PATTERN;
799:       if (!ijacobian) {         /* No IJacobian provided, but we have a separate RHS matrix */
800:         MatZeroEntries(A);
801:         MatShift(A,shift);
802:         if (A != B) {
803:           MatZeroEntries(B);
804:           MatShift(B,shift);
805:         }
806:       }
807:       MatAXPY(A,-1,Arhs,axpy);
808:       if (A != B) {
809:         MatAXPY(B,-1,Brhs,axpy);
810:       }
811:     }
812:   }
813:   PetscLogEventEnd(TS_JacobianEval,ts,U,A,B);
814:   return(0);
815: }

819: /*@C
820:     TSSetRHSFunction - Sets the routine for evaluating the function,
821:     where U_t = G(t,u).

823:     Logically Collective on TS

825:     Input Parameters:
826: +   ts - the TS context obtained from TSCreate()
827: .   r - vector to put the computed right hand side (or NULL to have it created)
828: .   f - routine for evaluating the right-hand-side function
829: -   ctx - [optional] user-defined context for private data for the
830:           function evaluation routine (may be NULL)

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

835: +   t - current timestep
836: .   u - input vector
837: .   F - function vector
838: -   ctx - [optional] user-defined function context

840:     Level: beginner

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

844: .seealso: TSSetRHSJacobian(), TSSetIJacobian(), TSSetIFunction()
845: @*/
846: PetscErrorCode  TSSetRHSFunction(TS ts,Vec r,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx)
847: {
849:   SNES           snes;
850:   Vec            ralloc = NULL;
851:   DM             dm;


857:   TSGetDM(ts,&dm);
858:   DMTSSetRHSFunction(dm,f,ctx);
859:   TSGetSNES(ts,&snes);
860:   if (!r && !ts->dm && ts->vec_sol) {
861:     VecDuplicate(ts->vec_sol,&ralloc);
862:     r    = ralloc;
863:   }
864:   SNESSetFunction(snes,r,SNESTSFormFunction,ts);
865:   VecDestroy(&ralloc);
866:   return(0);
867: }

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

874:     Logically Collective on TS

876:     Input Parameters:
877: +   ts - the TS context obtained from TSCreate()
878: .   f - routine for evaluating the solution
879: -   ctx - [optional] user-defined context for private data for the
880:           function evaluation routine (may be NULL)

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

885: +   t - current timestep
886: .   u - output vector
887: -   ctx - [optional] user-defined function context

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

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

896:     Level: beginner

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

900: .seealso: TSSetRHSJacobian(), TSSetIJacobian(), TSComputeSolutionFunction(), TSSetForcingFunction()
901: @*/
902: PetscErrorCode  TSSetSolutionFunction(TS ts,PetscErrorCode (*f)(TS,PetscReal,Vec,void*),void *ctx)
903: {
905:   DM             dm;

909:   TSGetDM(ts,&dm);
910:   DMTSSetSolutionFunction(dm,f,ctx);
911:   return(0);
912: }

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

919:     Logically Collective on TS

921:     Input Parameters:
922: +   ts - the TS context obtained from TSCreate()
923: .   f - routine for evaluating the forcing function
924: -   ctx - [optional] user-defined context for private data for the
925:           function evaluation routine (may be NULL)

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

930: +   t - current timestep
931: .   u - output vector
932: -   ctx - [optional] user-defined function context

934:     Notes:
935:     This routine is useful for testing accuracy of time integration schemes when using the Method of Manufactured Solutions to
936:     create closed-form solutions with a non-physical forcing term.

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

940:     Level: beginner

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

944: .seealso: TSSetRHSJacobian(), TSSetIJacobian(), TSComputeSolutionFunction(), TSSetSolutionFunction()
945: @*/
946: PetscErrorCode  TSSetForcingFunction(TS ts,PetscErrorCode (*f)(TS,PetscReal,Vec,void*),void *ctx)
947: {
949:   DM             dm;

953:   TSGetDM(ts,&dm);
954:   DMTSSetForcingFunction(dm,f,ctx);
955:   return(0);
956: }

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

964:    Logically Collective on TS

966:    Input Parameters:
967: +  ts  - the TS context obtained from TSCreate()
968: .  Amat - (approximate) Jacobian matrix
969: .  Pmat - matrix from which preconditioner is to be constructed (usually the same as Amat)
970: .  f   - the Jacobian evaluation routine
971: -  ctx - [optional] user-defined context for private data for the
972:          Jacobian evaluation routine (may be NULL)

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

977: +  t - current timestep
978: .  u - input vector
979: .  Amat - (approximate) Jacobian matrix
980: .  Pmat - matrix from which preconditioner is to be constructed (usually the same as Amat)
981: -  ctx - [optional] user-defined context for matrix evaluation routine


984:    Level: beginner

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

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

990: @*/
991: PetscErrorCode  TSSetRHSJacobian(TS ts,Mat Amat,Mat Pmat,TSRHSJacobian f,void *ctx)
992: {
994:   SNES           snes;
995:   DM             dm;
996:   TSIJacobian    ijacobian;


1005:   TSGetDM(ts,&dm);
1006:   DMTSSetRHSJacobian(dm,f,ctx);
1007:   if (f == TSComputeRHSJacobianConstant) {
1008:     /* Handle this case automatically for the user; otherwise user should call themselves. */
1009:     TSRHSJacobianSetReuse(ts,PETSC_TRUE);
1010:   }
1011:   DMTSGetIJacobian(dm,&ijacobian,NULL);
1012:   TSGetSNES(ts,&snes);
1013:   if (!ijacobian) {
1014:     SNESSetJacobian(snes,Amat,Pmat,SNESTSFormJacobian,ts);
1015:   }
1016:   if (Amat) {
1017:     PetscObjectReference((PetscObject)Amat);
1018:     MatDestroy(&ts->Arhs);

1020:     ts->Arhs = Amat;
1021:   }
1022:   if (Pmat) {
1023:     PetscObjectReference((PetscObject)Pmat);
1024:     MatDestroy(&ts->Brhs);

1026:     ts->Brhs = Pmat;
1027:   }
1028:   return(0);
1029: }


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

1037:    Logically Collective on TS

1039:    Input Parameters:
1040: +  ts  - the TS context obtained from TSCreate()
1041: .  r   - vector to hold the residual (or NULL to have it created internally)
1042: .  f   - the function evaluation routine
1043: -  ctx - user-defined context for private data for the function evaluation routine (may be NULL)

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

1048: +  t   - time at step/stage being solved
1049: .  u   - state vector
1050: .  u_t - time derivative of state vector
1051: .  F   - function vector
1052: -  ctx - [optional] user-defined context for matrix evaluation routine

1054:    Important:
1055:    The user MUST call either this routine, TSSetRHSFunction().  This routine must be used when not solving an ODE, for example a DAE.

1057:    Level: beginner

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

1061: .seealso: TSSetRHSJacobian(), TSSetRHSFunction(), TSSetIJacobian()
1062: @*/
1063: PetscErrorCode  TSSetIFunction(TS ts,Vec res,TSIFunction f,void *ctx)
1064: {
1066:   SNES           snes;
1067:   Vec            resalloc = NULL;
1068:   DM             dm;


1074:   TSGetDM(ts,&dm);
1075:   DMTSSetIFunction(dm,f,ctx);

1077:   TSGetSNES(ts,&snes);
1078:   if (!res && !ts->dm && ts->vec_sol) {
1079:     VecDuplicate(ts->vec_sol,&resalloc);
1080:     res  = resalloc;
1081:   }
1082:   SNESSetFunction(snes,res,SNESTSFormFunction,ts);
1083:   VecDestroy(&resalloc);
1084:   return(0);
1085: }

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

1092:    Not Collective

1094:    Input Parameter:
1095: .  ts - the TS context

1097:    Output Parameter:
1098: +  r - vector to hold residual (or NULL)
1099: .  func - the function to compute residual (or NULL)
1100: -  ctx - the function context (or NULL)

1102:    Level: advanced

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

1106: .seealso: TSSetIFunction(), SNESGetFunction()
1107: @*/
1108: PetscErrorCode TSGetIFunction(TS ts,Vec *r,TSIFunction *func,void **ctx)
1109: {
1111:   SNES           snes;
1112:   DM             dm;

1116:   TSGetSNES(ts,&snes);
1117:   SNESGetFunction(snes,r,NULL,NULL);
1118:   TSGetDM(ts,&dm);
1119:   DMTSGetIFunction(dm,func,ctx);
1120:   return(0);
1121: }

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

1128:    Not Collective

1130:    Input Parameter:
1131: .  ts - the TS context

1133:    Output Parameter:
1134: +  r - vector to hold computed right hand side (or NULL)
1135: .  func - the function to compute right hand side (or NULL)
1136: -  ctx - the function context (or NULL)

1138:    Level: advanced

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

1142: .seealso: TSSetRhsfunction(), SNESGetFunction()
1143: @*/
1144: PetscErrorCode TSGetRHSFunction(TS ts,Vec *r,TSRHSFunction *func,void **ctx)
1145: {
1147:   SNES           snes;
1148:   DM             dm;

1152:   TSGetSNES(ts,&snes);
1153:   SNESGetFunction(snes,r,NULL,NULL);
1154:   TSGetDM(ts,&dm);
1155:   DMTSGetRHSFunction(dm,func,ctx);
1156:   return(0);
1157: }

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

1165:    Logically Collective on TS

1167:    Input Parameters:
1168: +  ts  - the TS context obtained from TSCreate()
1169: .  Amat - (approximate) Jacobian matrix
1170: .  Pmat - matrix used to compute preconditioner (usually the same as Amat)
1171: .  f   - the Jacobian evaluation routine
1172: -  ctx - user-defined context for private data for the Jacobian evaluation routine (may be NULL)

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

1177: +  t    - time at step/stage being solved
1178: .  U    - state vector
1179: .  U_t  - time derivative of state vector
1180: .  a    - shift
1181: .  Amat - (approximate) Jacobian of F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t
1182: .  Pmat - matrix used for constructing preconditioner, usually the same as Amat
1183: -  ctx  - [optional] user-defined context for matrix evaluation routine

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

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

1195:    Level: beginner

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

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

1201: @*/
1202: PetscErrorCode  TSSetIJacobian(TS ts,Mat Amat,Mat Pmat,TSIJacobian f,void *ctx)
1203: {
1205:   SNES           snes;
1206:   DM             dm;


1215:   TSGetDM(ts,&dm);
1216:   DMTSSetIJacobian(dm,f,ctx);

1218:   TSGetSNES(ts,&snes);
1219:   SNESSetJacobian(snes,Amat,Pmat,SNESTSFormJacobian,ts);
1220:   return(0);
1221: }

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

1231:    Logically Collective

1233:    Input Arguments:
1234: +  ts - TS context obtained from TSCreate()
1235: -  reuse - PETSC_TRUE if the RHS Jacobian

1237:    Level: intermediate

1239: .seealso: TSSetRHSJacobian(), TSComputeRHSJacobianConstant()
1240: @*/
1241: PetscErrorCode TSRHSJacobianSetReuse(TS ts,PetscBool reuse)
1242: {
1244:   ts->rhsjacobian.reuse = reuse;
1245:   return(0);
1246: }

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

1253:   Collective on PetscViewer

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

1260:    Level: intermediate

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

1265:   Notes for advanced users:
1266:   Most users should not need to know the details of the binary storage
1267:   format, since TSLoad() and TSView() completely hide these details.
1268:   But for anyone who's interested, the standard binary matrix storage
1269:   format is
1270: .vb
1271:      has not yet been determined
1272: .ve

1274: .seealso: PetscViewerBinaryOpen(), TSView(), MatLoad(), VecLoad()
1275: @*/
1276: PetscErrorCode  TSLoad(TS ts, PetscViewer viewer)
1277: {
1279:   PetscBool      isbinary;
1280:   PetscInt       classid;
1281:   char           type[256];
1282:   DMTS           sdm;
1283:   DM             dm;

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

1291:   PetscViewerBinaryRead(viewer,&classid,1,PETSC_INT);
1292:   if (classid != TS_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONG,"Not TS next in file");
1293:   PetscViewerBinaryRead(viewer,type,256,PETSC_CHAR);
1294:   TSSetType(ts, type);
1295:   if (ts->ops->load) {
1296:     (*ts->ops->load)(ts,viewer);
1297:   }
1298:   DMCreate(PetscObjectComm((PetscObject)ts),&dm);
1299:   DMLoad(dm,viewer);
1300:   TSSetDM(ts,dm);
1301:   DMCreateGlobalVector(ts->dm,&ts->vec_sol);
1302:   VecLoad(ts->vec_sol,viewer);
1303:   DMGetDMTS(ts->dm,&sdm);
1304:   DMTSLoad(sdm,viewer);
1305:   return(0);
1306: }

1308: #include <petscdraw.h>
1309: #if defined(PETSC_HAVE_SAWS)
1310: #include <petscviewersaws.h>
1311: #endif
1314: /*@C
1315:     TSView - Prints the TS data structure.

1317:     Collective on TS

1319:     Input Parameters:
1320: +   ts - the TS context obtained from TSCreate()
1321: -   viewer - visualization context

1323:     Options Database Key:
1324: .   -ts_view - calls TSView() at end of TSStep()

1326:     Notes:
1327:     The available visualization contexts include
1328: +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
1329: -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
1330:          output where only the first processor opens
1331:          the file.  All other processors send their
1332:          data to the first processor to print.

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

1337:     Level: beginner

1339: .keywords: TS, timestep, view

1341: .seealso: PetscViewerASCIIOpen()
1342: @*/
1343: PetscErrorCode  TSView(TS ts,PetscViewer viewer)
1344: {
1346:   TSType         type;
1347:   PetscBool      iascii,isstring,isundials,isbinary,isdraw;
1348:   DMTS           sdm;
1349: #if defined(PETSC_HAVE_SAWS)
1350:   PetscBool      isams;
1351: #endif

1355:   if (!viewer) {
1356:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ts),&viewer);
1357:   }

1361:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1362:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
1363:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
1364:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
1365: #if defined(PETSC_HAVE_SAWS)
1366:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSAWS,&isams);
1367: #endif
1368:   if (iascii) {
1369:     PetscObjectPrintClassNamePrefixType((PetscObject)ts,viewer);
1370:     PetscViewerASCIIPrintf(viewer,"  maximum steps=%D\n",ts->max_steps);
1371:     PetscViewerASCIIPrintf(viewer,"  maximum time=%g\n",(double)ts->max_time);
1372:     if (ts->problem_type == TS_NONLINEAR) {
1373:       PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solver iterations=%D\n",ts->snes_its);
1374:       PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solve failures=%D\n",ts->num_snes_failures);
1375:     }
1376:     PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%D\n",ts->ksp_its);
1377:     PetscViewerASCIIPrintf(viewer,"  total number of rejected steps=%D\n",ts->reject);
1378:     DMGetDMTS(ts->dm,&sdm);
1379:     DMTSView(sdm,viewer);
1380:     if (ts->ops->view) {
1381:       PetscViewerASCIIPushTab(viewer);
1382:       (*ts->ops->view)(ts,viewer);
1383:       PetscViewerASCIIPopTab(viewer);
1384:     }
1385:   } else if (isstring) {
1386:     TSGetType(ts,&type);
1387:     PetscViewerStringSPrintf(viewer," %-7.7s",type);
1388:   } else if (isbinary) {
1389:     PetscInt    classid = TS_FILE_CLASSID;
1390:     MPI_Comm    comm;
1391:     PetscMPIInt rank;
1392:     char        type[256];

1394:     PetscObjectGetComm((PetscObject)ts,&comm);
1395:     MPI_Comm_rank(comm,&rank);
1396:     if (!rank) {
1397:       PetscViewerBinaryWrite(viewer,&classid,1,PETSC_INT,PETSC_FALSE);
1398:       PetscStrncpy(type,((PetscObject)ts)->type_name,256);
1399:       PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);
1400:     }
1401:     if (ts->ops->view) {
1402:       (*ts->ops->view)(ts,viewer);
1403:     }
1404:     DMView(ts->dm,viewer);
1405:     VecView(ts->vec_sol,viewer);
1406:     DMGetDMTS(ts->dm,&sdm);
1407:     DMTSView(sdm,viewer);
1408:   } else if (isdraw) {
1409:     PetscDraw draw;
1410:     char      str[36];
1411:     PetscReal x,y,bottom,h;

1413:     PetscViewerDrawGetDraw(viewer,0,&draw);
1414:     PetscDrawGetCurrentPoint(draw,&x,&y);
1415:     PetscStrcpy(str,"TS: ");
1416:     PetscStrcat(str,((PetscObject)ts)->type_name);
1417:     PetscDrawBoxedString(draw,x,y,PETSC_DRAW_BLACK,PETSC_DRAW_BLACK,str,NULL,&h);
1418:     bottom = y - h;
1419:     PetscDrawPushCurrentPoint(draw,x,bottom);
1420:     if (ts->ops->view) {
1421:       (*ts->ops->view)(ts,viewer);
1422:     }
1423:     PetscDrawPopCurrentPoint(draw);
1424: #if defined(PETSC_HAVE_SAWS)
1425:   } else if (isams) {
1426:     PetscMPIInt rank;
1427:     const char  *name;

1429:     PetscObjectGetName((PetscObject)ts,&name);
1430:     MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
1431:     if (!((PetscObject)ts)->amsmem && !rank) {
1432:       char       dir[1024];

1434:       PetscObjectViewSAWs((PetscObject)ts,viewer);
1435:       PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/time_step",name);
1436:       PetscStackCallSAWs(SAWs_Register,(dir,&ts->steps,1,SAWs_READ,SAWs_INT));
1437:       PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/time",name);
1438:       PetscStackCallSAWs(SAWs_Register,(dir,&ts->ptime,1,SAWs_READ,SAWs_DOUBLE));
1439:     }
1440:     if (ts->ops->view) {
1441:       (*ts->ops->view)(ts,viewer);
1442:     }
1443: #endif
1444:   }

1446:   PetscViewerASCIIPushTab(viewer);
1447:   PetscObjectTypeCompare((PetscObject)ts,TSSUNDIALS,&isundials);
1448:   PetscViewerASCIIPopTab(viewer);
1449:   return(0);
1450: }


1455: /*@
1456:    TSSetApplicationContext - Sets an optional user-defined context for
1457:    the timesteppers.

1459:    Logically Collective on TS

1461:    Input Parameters:
1462: +  ts - the TS context obtained from TSCreate()
1463: -  usrP - optional user context

1465:    Level: intermediate

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

1469: .seealso: TSGetApplicationContext()
1470: @*/
1471: PetscErrorCode  TSSetApplicationContext(TS ts,void *usrP)
1472: {
1475:   ts->user = usrP;
1476:   return(0);
1477: }

1481: /*@
1482:     TSGetApplicationContext - Gets the user-defined context for the
1483:     timestepper.

1485:     Not Collective

1487:     Input Parameter:
1488: .   ts - the TS context obtained from TSCreate()

1490:     Output Parameter:
1491: .   usrP - user context

1493:     Level: intermediate

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

1497: .seealso: TSSetApplicationContext()
1498: @*/
1499: PetscErrorCode  TSGetApplicationContext(TS ts,void *usrP)
1500: {
1503:   *(void**)usrP = ts->user;
1504:   return(0);
1505: }

1509: /*@
1510:    TSGetTimeStepNumber - Gets the number of time steps completed.

1512:    Not Collective

1514:    Input Parameter:
1515: .  ts - the TS context obtained from TSCreate()

1517:    Output Parameter:
1518: .  iter - number of steps completed so far

1520:    Level: intermediate

1522: .keywords: TS, timestep, get, iteration, number
1523: .seealso: TSGetTime(), TSGetTimeStep(), TSSetPreStep(), TSSetPreStage(), TSSetPostStage(), TSSetPostStep()
1524: @*/
1525: PetscErrorCode  TSGetTimeStepNumber(TS ts,PetscInt *iter)
1526: {
1530:   *iter = ts->steps;
1531:   return(0);
1532: }

1536: /*@
1537:    TSSetInitialTimeStep - Sets the initial timestep to be used,
1538:    as well as the initial time.

1540:    Logically Collective on TS

1542:    Input Parameters:
1543: +  ts - the TS context obtained from TSCreate()
1544: .  initial_time - the initial time
1545: -  time_step - the size of the timestep

1547:    Level: intermediate

1549: .seealso: TSSetTimeStep(), TSGetTimeStep()

1551: .keywords: TS, set, initial, timestep
1552: @*/
1553: PetscErrorCode  TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step)
1554: {

1559:   TSSetTimeStep(ts,time_step);
1560:   TSSetTime(ts,initial_time);
1561:   return(0);
1562: }

1566: /*@
1567:    TSSetTimeStep - Allows one to reset the timestep at any time,
1568:    useful for simple pseudo-timestepping codes.

1570:    Logically Collective on TS

1572:    Input Parameters:
1573: +  ts - the TS context obtained from TSCreate()
1574: -  time_step - the size of the timestep

1576:    Level: intermediate

1578: .seealso: TSSetInitialTimeStep(), TSGetTimeStep()

1580: .keywords: TS, set, timestep
1581: @*/
1582: PetscErrorCode  TSSetTimeStep(TS ts,PetscReal time_step)
1583: {
1587:   ts->time_step      = time_step;
1588:   ts->time_step_orig = time_step;
1589:   return(0);
1590: }

1594: /*@
1595:    TSSetExactFinalTime - Determines whether to adapt the final time step to
1596:      match the exact final time, interpolate solution to the exact final time,
1597:      or just return at the final time TS computed.

1599:   Logically Collective on TS

1601:    Input Parameter:
1602: +   ts - the time-step context
1603: -   eftopt - exact final time option

1605:    Level: beginner

1607: .seealso: TSExactFinalTimeOption
1608: @*/
1609: PetscErrorCode  TSSetExactFinalTime(TS ts,TSExactFinalTimeOption eftopt)
1610: {
1614:   ts->exact_final_time = eftopt;
1615:   return(0);
1616: }

1620: /*@
1621:    TSGetTimeStep - Gets the current timestep size.

1623:    Not Collective

1625:    Input Parameter:
1626: .  ts - the TS context obtained from TSCreate()

1628:    Output Parameter:
1629: .  dt - the current timestep size

1631:    Level: intermediate

1633: .seealso: TSSetInitialTimeStep(), TSGetTimeStep()

1635: .keywords: TS, get, timestep
1636: @*/
1637: PetscErrorCode  TSGetTimeStep(TS ts,PetscReal *dt)
1638: {
1642:   *dt = ts->time_step;
1643:   return(0);
1644: }

1648: /*@
1649:    TSGetSolution - Returns the solution at the present timestep. It
1650:    is valid to call this routine inside the function that you are evaluating
1651:    in order to move to the new timestep. This vector not changed until
1652:    the solution at the next timestep has been calculated.

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

1656:    Input Parameter:
1657: .  ts - the TS context obtained from TSCreate()

1659:    Output Parameter:
1660: .  v - the vector containing the solution

1662:    Level: intermediate

1664: .seealso: TSGetTimeStep()

1666: .keywords: TS, timestep, get, solution
1667: @*/
1668: PetscErrorCode  TSGetSolution(TS ts,Vec *v)
1669: {
1673:   *v = ts->vec_sol;
1674:   return(0);
1675: }

1677: /* ----- Routines to initialize and destroy a timestepper ---- */
1680: /*@
1681:   TSSetProblemType - Sets the type of problem to be solved.

1683:   Not collective

1685:   Input Parameters:
1686: + ts   - The TS
1687: - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
1688: .vb
1689:          U_t - A U = 0      (linear)
1690:          U_t - A(t) U = 0   (linear)
1691:          F(t,U,U_t) = 0     (nonlinear)
1692: .ve

1694:    Level: beginner

1696: .keywords: TS, problem type
1697: .seealso: TSSetUp(), TSProblemType, TS
1698: @*/
1699: PetscErrorCode  TSSetProblemType(TS ts, TSProblemType type)
1700: {

1705:   ts->problem_type = type;
1706:   if (type == TS_LINEAR) {
1707:     SNES snes;
1708:     TSGetSNES(ts,&snes);
1709:     SNESSetType(snes,SNESKSPONLY);
1710:   }
1711:   return(0);
1712: }

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

1719:   Not collective

1721:   Input Parameter:
1722: . ts   - The TS

1724:   Output Parameter:
1725: . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
1726: .vb
1727:          M U_t = A U
1728:          M(t) U_t = A(t) U
1729:          F(t,U,U_t)
1730: .ve

1732:    Level: beginner

1734: .keywords: TS, problem type
1735: .seealso: TSSetUp(), TSProblemType, TS
1736: @*/
1737: PetscErrorCode  TSGetProblemType(TS ts, TSProblemType *type)
1738: {
1742:   *type = ts->problem_type;
1743:   return(0);
1744: }

1748: /*@
1749:    TSSetUp - Sets up the internal data structures for the later use
1750:    of a timestepper.

1752:    Collective on TS

1754:    Input Parameter:
1755: .  ts - the TS context obtained from TSCreate()

1757:    Notes:
1758:    For basic use of the TS solvers the user need not explicitly call
1759:    TSSetUp(), since these actions will automatically occur during
1760:    the call to TSStep().  However, if one wishes to control this
1761:    phase separately, TSSetUp() should be called after TSCreate()
1762:    and optional routines of the form TSSetXXX(), but before TSStep().

1764:    Level: advanced

1766: .keywords: TS, timestep, setup

1768: .seealso: TSCreate(), TSStep(), TSDestroy()
1769: @*/
1770: PetscErrorCode  TSSetUp(TS ts)
1771: {
1773:   DM             dm;
1774:   PetscErrorCode (*func)(SNES,Vec,Vec,void*);
1775:   PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*);
1776:   TSIJacobian    ijac;
1777:   TSRHSJacobian  rhsjac;

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

1783:   if (!((PetscObject)ts)->type_name) {
1784:     TSSetType(ts,TSEULER);
1785:   }

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

1789:   TSGetAdapt(ts,&ts->adapt);

1791:   if (ts->rhsjacobian.reuse) {
1792:     Mat Amat,Pmat;
1793:     SNES snes;
1794:     TSGetSNES(ts,&snes);
1795:     SNESGetJacobian(snes,&Amat,&Pmat,NULL,NULL);
1796:     /* Matching matrices implies that an IJacobian is NOT set, because if it had been set, the IJacobian's matrix would
1797:      * have displaced the RHS matrix */
1798:     if (Amat == ts->Arhs) {
1799:       MatDuplicate(ts->Arhs,MAT_DO_NOT_COPY_VALUES,&Amat);
1800:       SNESSetJacobian(snes,Amat,NULL,NULL,NULL);
1801:       MatDestroy(&Amat);
1802:     }
1803:     if (Pmat == ts->Brhs) {
1804:       MatDuplicate(ts->Brhs,MAT_DO_NOT_COPY_VALUES,&Pmat);
1805:       SNESSetJacobian(snes,NULL,Pmat,NULL,NULL);
1806:       MatDestroy(&Pmat);
1807:     }
1808:   }

1810:   if (ts->ops->setup) {
1811:     (*ts->ops->setup)(ts);
1812:   }

1814:   /* in the case where we've set a DMTSFunction or what have you, we need the default SNESFunction
1815:    to be set right but can't do it elsewhere due to the overreliance on ctx=ts.
1816:    */
1817:   TSGetDM(ts,&dm);
1818:   DMSNESGetFunction(dm,&func,NULL);
1819:   if (!func) {
1820:     ierr =DMSNESSetFunction(dm,SNESTSFormFunction,ts);
1821:   }
1822:   /* if the SNES doesn't have a jacobian set and the TS has an ijacobian or rhsjacobian set, set the SNES to use it.
1823:      Otherwise, the SNES will use coloring internally to form the Jacobian.
1824:    */
1825:   DMSNESGetJacobian(dm,&jac,NULL);
1826:   DMTSGetIJacobian(dm,&ijac,NULL);
1827:   DMTSGetRHSJacobian(dm,&rhsjac,NULL);
1828:   if (!jac && (ijac || rhsjac)) {
1829:     DMSNESSetJacobian(dm,SNESTSFormJacobian,ts);
1830:   }
1831:   ts->setupcalled = PETSC_TRUE;
1832:   return(0);
1833: }

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

1840:    Collective on TS

1842:    Input Parameter:
1843: .  ts - the TS context obtained from TSCreate()

1845:    Level: beginner

1847: .keywords: TS, timestep, reset

1849: .seealso: TSCreate(), TSSetup(), TSDestroy()
1850: @*/
1851: PetscErrorCode  TSReset(TS ts)
1852: {

1857:   if (ts->ops->reset) {
1858:     (*ts->ops->reset)(ts);
1859:   }
1860:   if (ts->snes) {SNESReset(ts->snes);}

1862:   MatDestroy(&ts->Arhs);
1863:   MatDestroy(&ts->Brhs);
1864:   VecDestroy(&ts->Frhs);
1865:   VecDestroy(&ts->vec_sol);
1866:   VecDestroy(&ts->vatol);
1867:   VecDestroy(&ts->vrtol);
1868:   VecDestroyVecs(ts->nwork,&ts->work);

1870:   ts->setupcalled = PETSC_FALSE;
1871:   return(0);
1872: }

1876: /*@
1877:    TSDestroy - Destroys the timestepper context that was created
1878:    with TSCreate().

1880:    Collective on TS

1882:    Input Parameter:
1883: .  ts - the TS context obtained from TSCreate()

1885:    Level: beginner

1887: .keywords: TS, timestepper, destroy

1889: .seealso: TSCreate(), TSSetUp(), TSSolve()
1890: @*/
1891: PetscErrorCode  TSDestroy(TS *ts)
1892: {

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

1900:   TSReset((*ts));

1902:   /* if memory was published with SAWs then destroy it */
1903:   PetscObjectSAWsViewOff((PetscObject)*ts);
1904:   if ((*ts)->ops->destroy) {(*(*ts)->ops->destroy)((*ts));}

1906:   TSAdaptDestroy(&(*ts)->adapt);
1907:   if ((*ts)->event) {
1908:     TSEventMonitorDestroy(&(*ts)->event);
1909:   }
1910:   SNESDestroy(&(*ts)->snes);
1911:   DMDestroy(&(*ts)->dm);
1912:   TSMonitorCancel((*ts));

1914:   PetscHeaderDestroy(ts);
1915:   return(0);
1916: }

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

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

1926:    Input Parameter:
1927: .  ts - the TS context obtained from TSCreate()

1929:    Output Parameter:
1930: .  snes - the nonlinear solver context

1932:    Notes:
1933:    The user can then directly manipulate the SNES context to set various
1934:    options, etc.  Likewise, the user can then extract and manipulate the
1935:    KSP, KSP, and PC contexts as well.

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

1940:    Level: beginner

1942: .keywords: timestep, get, SNES
1943: @*/
1944: PetscErrorCode  TSGetSNES(TS ts,SNES *snes)
1945: {

1951:   if (!ts->snes) {
1952:     SNESCreate(PetscObjectComm((PetscObject)ts),&ts->snes);
1953:     SNESSetFunction(ts->snes,NULL,SNESTSFormFunction,ts);
1954:     PetscLogObjectParent((PetscObject)ts,(PetscObject)ts->snes);
1955:     PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);
1956:     if (ts->dm) {SNESSetDM(ts->snes,ts->dm);}
1957:     if (ts->problem_type == TS_LINEAR) {
1958:       SNESSetType(ts->snes,SNESKSPONLY);
1959:     }
1960:   }
1961:   *snes = ts->snes;
1962:   return(0);
1963: }

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

1970:    Collective

1972:    Input Parameter:
1973: +  ts - the TS context obtained from TSCreate()
1974: -  snes - the nonlinear solver context

1976:    Notes:
1977:    Most users should have the TS created by calling TSGetSNES()

1979:    Level: developer

1981: .keywords: timestep, set, SNES
1982: @*/
1983: PetscErrorCode TSSetSNES(TS ts,SNES snes)
1984: {
1986:   PetscErrorCode (*func)(SNES,Vec,Mat,Mat,void*);

1991:   PetscObjectReference((PetscObject)snes);
1992:   SNESDestroy(&ts->snes);

1994:   ts->snes = snes;

1996:   SNESSetFunction(ts->snes,NULL,SNESTSFormFunction,ts);
1997:   SNESGetJacobian(ts->snes,NULL,NULL,&func,NULL);
1998:   if (func == SNESTSFormJacobian) {
1999:     SNESSetJacobian(ts->snes,NULL,NULL,SNESTSFormJacobian,ts);
2000:   }
2001:   return(0);
2002: }

2006: /*@
2007:    TSGetKSP - Returns the KSP (linear solver) associated with
2008:    a TS (timestepper) context.

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

2012:    Input Parameter:
2013: .  ts - the TS context obtained from TSCreate()

2015:    Output Parameter:
2016: .  ksp - the nonlinear solver context

2018:    Notes:
2019:    The user can then directly manipulate the KSP context to set various
2020:    options, etc.  Likewise, the user can then extract and manipulate the
2021:    KSP and PC contexts as well.

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

2026:    Level: beginner

2028: .keywords: timestep, get, KSP
2029: @*/
2030: PetscErrorCode  TSGetKSP(TS ts,KSP *ksp)
2031: {
2033:   SNES           snes;

2038:   if (!((PetscObject)ts)->type_name) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"KSP is not created yet. Call TSSetType() first");
2039:   if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()");
2040:   TSGetSNES(ts,&snes);
2041:   SNESGetKSP(snes,ksp);
2042:   return(0);
2043: }

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

2049: /*@
2050:    TSGetDuration - Gets the maximum number of timesteps to use and
2051:    maximum time for iteration.

2053:    Not Collective

2055:    Input Parameters:
2056: +  ts       - the TS context obtained from TSCreate()
2057: .  maxsteps - maximum number of iterations to use, or NULL
2058: -  maxtime  - final time to iterate to, or NULL

2060:    Level: intermediate

2062: .keywords: TS, timestep, get, maximum, iterations, time
2063: @*/
2064: PetscErrorCode  TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime)
2065: {
2068:   if (maxsteps) {
2070:     *maxsteps = ts->max_steps;
2071:   }
2072:   if (maxtime) {
2074:     *maxtime = ts->max_time;
2075:   }
2076:   return(0);
2077: }

2081: /*@
2082:    TSSetDuration - Sets the maximum number of timesteps to use and
2083:    maximum time for iteration.

2085:    Logically Collective on TS

2087:    Input Parameters:
2088: +  ts - the TS context obtained from TSCreate()
2089: .  maxsteps - maximum number of iterations to use
2090: -  maxtime - final time to iterate to

2092:    Options Database Keys:
2093: .  -ts_max_steps <maxsteps> - Sets maxsteps
2094: .  -ts_final_time <maxtime> - Sets maxtime

2096:    Notes:
2097:    The default maximum number of iterations is 5000. Default time is 5.0

2099:    Level: intermediate

2101: .keywords: TS, timestep, set, maximum, iterations

2103: .seealso: TSSetExactFinalTime()
2104: @*/
2105: PetscErrorCode  TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime)
2106: {
2111:   if (maxsteps >= 0) ts->max_steps = maxsteps;
2112:   if (maxtime != PETSC_DEFAULT) ts->max_time = maxtime;
2113:   return(0);
2114: }

2118: /*@
2119:    TSSetSolution - Sets the initial solution vector
2120:    for use by the TS routines.

2122:    Logically Collective on TS and Vec

2124:    Input Parameters:
2125: +  ts - the TS context obtained from TSCreate()
2126: -  u - the solution vector

2128:    Level: beginner

2130: .keywords: TS, timestep, set, solution, initial conditions
2131: @*/
2132: PetscErrorCode  TSSetSolution(TS ts,Vec u)
2133: {
2135:   DM             dm;

2140:   PetscObjectReference((PetscObject)u);
2141:   VecDestroy(&ts->vec_sol);

2143:   ts->vec_sol = u;

2145:   TSGetDM(ts,&dm);
2146:   DMShellSetGlobalVector(dm,u);
2147:   return(0);
2148: }

2152: /*@C
2153:   TSSetPreStep - Sets the general-purpose function
2154:   called once at the beginning of each time step.

2156:   Logically Collective on TS

2158:   Input Parameters:
2159: + ts   - The TS context obtained from TSCreate()
2160: - func - The function

2162:   Calling sequence of func:
2163: . func (TS ts);

2165:   Level: intermediate

2167:   Note:
2168:   If a step is rejected, TSStep() will call this routine again before each attempt.
2169:   The last completed time step number can be queried using TSGetTimeStepNumber(), the
2170:   size of the step being attempted can be obtained using TSGetTimeStep().

2172: .keywords: TS, timestep
2173: .seealso: TSSetPreStage(), TSSetPostStage(), TSSetPostStep(), TSStep()
2174: @*/
2175: PetscErrorCode  TSSetPreStep(TS ts, PetscErrorCode (*func)(TS))
2176: {
2179:   ts->prestep = func;
2180:   return(0);
2181: }

2185: /*@
2186:   TSPreStep - Runs the user-defined pre-step function.

2188:   Collective on TS

2190:   Input Parameters:
2191: . ts   - The TS context obtained from TSCreate()

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

2197:   Level: developer

2199: .keywords: TS, timestep
2200: .seealso: TSSetPreStep(), TSPreStage(), TSPostStage(), TSPostStep()
2201: @*/
2202: PetscErrorCode  TSPreStep(TS ts)
2203: {

2208:   if (ts->prestep) {
2209:     PetscStackCallStandard((*ts->prestep),(ts));
2210:   }
2211:   return(0);
2212: }

2216: /*@C
2217:   TSSetPreStage - Sets the general-purpose function
2218:   called once at the beginning of each stage.

2220:   Logically Collective on TS

2222:   Input Parameters:
2223: + ts   - The TS context obtained from TSCreate()
2224: - func - The function

2226:   Calling sequence of func:
2227: . PetscErrorCode func(TS ts, PetscReal stagetime);

2229:   Level: intermediate

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

2236: .keywords: TS, timestep
2237: .seealso: TSSetPostStage(), TSSetPreStep(), TSSetPostStep(), TSGetApplicationContext()
2238: @*/
2239: PetscErrorCode  TSSetPreStage(TS ts, PetscErrorCode (*func)(TS,PetscReal))
2240: {
2243:   ts->prestage = func;
2244:   return(0);
2245: }

2249: /*@C
2250:   TSSetPostStage - Sets the general-purpose function
2251:   called once at the end of each stage.

2253:   Logically Collective on TS

2255:   Input Parameters:
2256: + ts   - The TS context obtained from TSCreate()
2257: - func - The function

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

2262:   Level: intermediate

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

2269: .keywords: TS, timestep
2270: .seealso: TSSetPreStage(), TSSetPreStep(), TSSetPostStep(), TSGetApplicationContext()
2271: @*/
2272: PetscErrorCode  TSSetPostStage(TS ts, PetscErrorCode (*func)(TS,PetscReal,PetscInt,Vec*))
2273: {
2276:   ts->poststage = func;
2277:   return(0);
2278: }

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

2285:   Collective on TS

2287:   Input Parameters:
2288: . ts          - The TS context obtained from TSCreate()
2289:   stagetime   - The absolute time of the current stage

2291:   Notes:
2292:   TSPreStage() is typically used within time stepping implementations,
2293:   most users would not generally call this routine themselves.

2295:   Level: developer

2297: .keywords: TS, timestep
2298: .seealso: TSPostStage(), TSSetPreStep(), TSPreStep(), TSPostStep()
2299: @*/
2300: PetscErrorCode  TSPreStage(TS ts, PetscReal stagetime)
2301: {

2306:   if (ts->prestage) {
2307:     PetscStackCallStandard((*ts->prestage),(ts,stagetime));
2308:   }
2309:   return(0);
2310: }

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

2317:   Collective on TS

2319:   Input Parameters:
2320: . ts          - The TS context obtained from TSCreate()
2321:   stagetime   - The absolute time of the current stage
2322:   stageindex  - Stage number
2323:   Y           - Array of vectors (of size = total number
2324:                 of stages) with the stage solutions

2326:   Notes:
2327:   TSPostStage() is typically used within time stepping implementations,
2328:   most users would not generally call this routine themselves.

2330:   Level: developer

2332: .keywords: TS, timestep
2333: .seealso: TSPreStage(), TSSetPreStep(), TSPreStep(), TSPostStep()
2334: @*/
2335: PetscErrorCode  TSPostStage(TS ts, PetscReal stagetime, PetscInt stageindex, Vec *Y)
2336: {

2341:   if (ts->poststage) {
2342:     PetscStackCallStandard((*ts->poststage),(ts,stagetime,stageindex,Y));
2343:   }
2344:   return(0);
2345: }

2349: /*@C
2350:   TSSetPostStep - Sets the general-purpose function
2351:   called once at the end of each time step.

2353:   Logically Collective on TS

2355:   Input Parameters:
2356: + ts   - The TS context obtained from TSCreate()
2357: - func - The function

2359:   Calling sequence of func:
2360: $ func (TS ts);

2362:   Level: intermediate

2364: .keywords: TS, timestep
2365: .seealso: TSSetPreStep(), TSSetPreStage(), TSGetTimeStep(), TSGetTimeStepNumber(), TSGetTime()
2366: @*/
2367: PetscErrorCode  TSSetPostStep(TS ts, PetscErrorCode (*func)(TS))
2368: {
2371:   ts->poststep = func;
2372:   return(0);
2373: }

2377: /*@
2378:   TSPostStep - Runs the user-defined post-step function.

2380:   Collective on TS

2382:   Input Parameters:
2383: . ts   - The TS context obtained from TSCreate()

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

2389:   Level: developer

2391: .keywords: TS, timestep
2392: @*/
2393: PetscErrorCode  TSPostStep(TS ts)
2394: {

2399:   if (ts->poststep) {
2400:     PetscStackCallStandard((*ts->poststep),(ts));
2401:   }
2402:   return(0);
2403: }

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

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

2413:    Logically Collective on TS

2415:    Input Parameters:
2416: +  ts - the TS context obtained from TSCreate()
2417: .  monitor - monitoring routine
2418: .  mctx - [optional] user-defined context for private data for the
2419:              monitor routine (use NULL if no context is desired)
2420: -  monitordestroy - [optional] routine that frees monitor context
2421:           (may be NULL)

2423:    Calling sequence of monitor:
2424: $    int monitor(TS ts,PetscInt steps,PetscReal time,Vec u,void *mctx)

2426: +    ts - the TS context
2427: .    steps - iteration number (after the final time step the monitor routine is called with a step of -1, this is at the final time which may have
2428:                                been interpolated to)
2429: .    time - current time
2430: .    u - current iterate
2431: -    mctx - [optional] monitoring context

2433:    Notes:
2434:    This routine adds an additional monitor to the list of monitors that
2435:    already has been loaded.

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

2439:    Level: intermediate

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

2443: .seealso: TSMonitorDefault(), TSMonitorCancel()
2444: @*/
2445: PetscErrorCode  TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void**))
2446: {
2449:   if (ts->numbermonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
2450:   ts->monitor[ts->numbermonitors]          = monitor;
2451:   ts->monitordestroy[ts->numbermonitors]   = mdestroy;
2452:   ts->monitorcontext[ts->numbermonitors++] = (void*)mctx;
2453:   return(0);
2454: }

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

2461:    Logically Collective on TS

2463:    Input Parameters:
2464: .  ts - the TS context obtained from TSCreate()

2466:    Notes:
2467:    There is no way to remove a single, specific monitor.

2469:    Level: intermediate

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

2473: .seealso: TSMonitorDefault(), TSMonitorSet()
2474: @*/
2475: PetscErrorCode  TSMonitorCancel(TS ts)
2476: {
2478:   PetscInt       i;

2482:   for (i=0; i<ts->numbermonitors; i++) {
2483:     if (ts->monitordestroy[i]) {
2484:       (*ts->monitordestroy[i])(&ts->monitorcontext[i]);
2485:     }
2486:   }
2487:   ts->numbermonitors = 0;
2488:   return(0);
2489: }

2493: /*@
2494:    TSMonitorDefault - Sets the Default monitor

2496:    Level: intermediate

2498: .keywords: TS, set, monitor

2500: .seealso: TSMonitorDefault(), TSMonitorSet()
2501: @*/
2502: PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy)
2503: {
2505:   PetscViewer    viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)ts));

2508:   PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);
2509:   PetscViewerASCIIPrintf(viewer,"%D TS dt %g time %g\n",step,(double)ts->time_step,(double)ptime);
2510:   PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);
2511:   return(0);
2512: }

2516: /*@
2517:    TSSetRetainStages - Request that all stages in the upcoming step be stored so that interpolation will be available.

2519:    Logically Collective on TS

2521:    Input Argument:
2522: .  ts - time stepping context

2524:    Output Argument:
2525: .  flg - PETSC_TRUE or PETSC_FALSE

2527:    Level: intermediate

2529: .keywords: TS, set

2531: .seealso: TSInterpolate(), TSSetPostStep()
2532: @*/
2533: PetscErrorCode TSSetRetainStages(TS ts,PetscBool flg)
2534: {
2537:   ts->retain_stages = flg;
2538:   return(0);
2539: }

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

2546:    Collective on TS

2548:    Input Argument:
2549: +  ts - time stepping context
2550: -  t - time to interpolate to

2552:    Output Argument:
2553: .  U - state at given time

2555:    Notes:
2556:    The user should call TSSetRetainStages() before taking a step in which interpolation will be requested.

2558:    Level: intermediate

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

2563: .keywords: TS, set

2565: .seealso: TSSetRetainStages(), TSSetPostStep()
2566: @*/
2567: PetscErrorCode TSInterpolate(TS ts,PetscReal t,Vec U)
2568: {

2574:   if (t < ts->ptime - ts->time_step_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-ts->time_step_prev),(double)ts->ptime);
2575:   if (!ts->ops->interpolate) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"%s does not provide interpolation",((PetscObject)ts)->type_name);
2576:   (*ts->ops->interpolate)(ts,t,U);
2577:   return(0);
2578: }

2582: /*@
2583:    TSStep - Steps one time step

2585:    Collective on TS

2587:    Input Parameter:
2588: .  ts - the TS context obtained from TSCreate()

2590:    Level: intermediate

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

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

2599: .keywords: TS, timestep, solve

2601: .seealso: TSCreate(), TSSetUp(), TSDestroy(), TSSolve(), TSSetPreStep(), TSSetPreStage(), TSSetPostStage(), TSInterpolate()
2602: @*/
2603: PetscErrorCode  TSStep(TS ts)
2604: {
2605:   DM               dm;
2606:   PetscErrorCode   ierr;
2607:   static PetscBool cite = PETSC_FALSE;

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

2619:   TSGetDM(ts, &dm);
2620:   TSSetUp(ts);

2622:   ts->reason = TS_CONVERGED_ITERATING;
2623:   ts->ptime_prev = ts->ptime;
2624:   DMSetOutputSequenceNumber(dm, ts->steps, ts->ptime);
2625:   VecViewFromOptions(ts->vec_sol, ((PetscObject) ts)->prefix, "-ts_view_solution");

2627:   if (!ts->ops->step) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSStep not implemented for type '%s'",((PetscObject)ts)->type_name);
2628:   PetscLogEventBegin(TS_Step,ts,0,0,0);
2629:   (*ts->ops->step)(ts);
2630:   PetscLogEventEnd(TS_Step,ts,0,0,0);

2632:   ts->time_step_prev = ts->ptime - ts->ptime_prev;
2633:   DMSetOutputSequenceNumber(dm, ts->steps, ts->ptime);

2635:   if (ts->reason < 0) {
2636:     if (ts->errorifstepfailed) {
2637:       if (ts->reason == TS_DIVERGED_NONLINEAR_SOLVE) {
2638:         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]);
2639:       } else if (ts->reason == TS_DIVERGED_STEP_REJECTED) {
2640:         SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s, increase -ts_max_reject or make negative to attempt recovery",TSConvergedReasons[ts->reason]);
2641:       } else SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s",TSConvergedReasons[ts->reason]);
2642:     }
2643:   } else if (!ts->reason) {
2644:     if (ts->steps >= ts->max_steps)     ts->reason = TS_CONVERGED_ITS;
2645:     else if (ts->ptime >= ts->max_time) ts->reason = TS_CONVERGED_TIME;
2646:   }
2647:   return(0);
2648: }

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

2655:    Collective on TS

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

2662:    Output Arguments:
2663: .  U - state at the end of the current step

2665:    Level: advanced

2667:    Notes:
2668:    This function cannot be called until all stages have been evaluated.
2669:    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.

2671: .seealso: TSStep(), TSAdapt
2672: @*/
2673: PetscErrorCode TSEvaluateStep(TS ts,PetscInt order,Vec U,PetscBool *done)
2674: {

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

2688: /*@
2689:    TSSolve - Steps the requested number of timesteps.

2691:    Collective on TS

2693:    Input Parameter:
2694: +  ts - the TS context obtained from TSCreate()
2695: -  u - the solution vector  (can be null if TSSetSolution() was used, otherwise must contain the initial conditions)

2697:    Level: beginner

2699:    Notes:
2700:    The final time returned by this function may be different from the time of the internally
2701:    held state accessible by TSGetSolution() and TSGetTime() because the method may have
2702:    stepped over the final time.

2704: .keywords: TS, timestep, solve

2706: .seealso: TSCreate(), TSSetSolution(), TSStep()
2707: @*/
2708: PetscErrorCode TSSolve(TS ts,Vec u)
2709: {
2710:   Vec               solution;
2711:   PetscErrorCode    ierr;

2716:   if (ts->exact_final_time == TS_EXACTFINALTIME_INTERPOLATE) {   /* Need ts->vec_sol to be distinct so it is not overwritten when we interpolate at the end */
2718:     if (!ts->vec_sol || u == ts->vec_sol) {
2719:       VecDuplicate(u,&solution);
2720:       TSSetSolution(ts,solution);
2721:       VecDestroy(&solution); /* grant ownership */
2722:     }
2723:     VecCopy(u,ts->vec_sol);
2724:   } else if (u) {
2725:     TSSetSolution(ts,u);
2726:   }
2727:   TSSetUp(ts);
2728:   /* reset time step and iteration counters */
2729:   ts->steps             = 0;
2730:   ts->ksp_its           = 0;
2731:   ts->snes_its          = 0;
2732:   ts->num_snes_failures = 0;
2733:   ts->reject            = 0;
2734:   ts->reason            = TS_CONVERGED_ITERATING;

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

2738:   if (ts->ops->solve) {         /* This private interface is transitional and should be removed when all implementations are updated. */
2739:     (*ts->ops->solve)(ts);
2740:     VecCopy(ts->vec_sol,u);
2741:     ts->solvetime = ts->ptime;
2742:   } else {
2743:     /* steps the requested number of timesteps. */
2744:     if (ts->steps >= ts->max_steps)     ts->reason = TS_CONVERGED_ITS;
2745:     else if (ts->ptime >= ts->max_time) ts->reason = TS_CONVERGED_TIME;
2746:     while (!ts->reason) {
2747:       TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);
2748:       TSStep(ts);
2749:       if (ts->event) {
2750:         TSEventMonitor(ts);
2751:         if (ts->event->status != TSEVENT_PROCESSING) {
2752:           TSPostStep(ts);
2753:         }
2754:       } else {
2755:         TSPostStep(ts);
2756:       }
2757:     }
2758:     if (ts->exact_final_time == TS_EXACTFINALTIME_INTERPOLATE && ts->ptime > ts->max_time) {
2759:       TSInterpolate(ts,ts->max_time,u);
2760:       ts->solvetime = ts->max_time;
2761:       solution = u;
2762:     } else {
2763:       if (u) {VecCopy(ts->vec_sol,u);}
2764:       ts->solvetime = ts->ptime;
2765:       solution = ts->vec_sol;
2766:     }
2767:     TSMonitor(ts,ts->steps,ts->solvetime,solution);
2768:     VecViewFromOptions(u, ((PetscObject) ts)->prefix, "-ts_view_solution");
2769:   }
2770:   TSViewFromOptions(ts,NULL,"-ts_view");
2771:   PetscObjectSAWsBlock((PetscObject)ts);
2772:   return(0);
2773: }

2777: /*@
2778:    TSMonitor - Runs all user-provided monitor routines set using TSMonitorSet()

2780:    Collective on TS

2782:    Input Parameters:
2783: +  ts - time stepping context obtained from TSCreate()
2784: .  step - step number that has just completed
2785: .  ptime - model time of the state
2786: -  u - state at the current model time

2788:    Notes:
2789:    TSMonitor() is typically used within the time stepping implementations.
2790:    Users might call this function when using the TSStep() interface instead of TSSolve().

2792:    Level: advanced

2794: .keywords: TS, timestep
2795: @*/
2796: PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec u)
2797: {
2799:   PetscInt       i,n = ts->numbermonitors;

2804:   for (i=0; i<n; i++) {
2805:     (*ts->monitor[i])(ts,step,ptime,u,ts->monitorcontext[i]);
2806:   }
2807:   return(0);
2808: }

2810: /* ------------------------------------------------------------------------*/
2813: /*@C
2814:    TSMonitorLGCtxCreate - Creates a line graph context for use with
2815:    TS to monitor the solution process graphically in various ways

2817:    Collective on TS

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

2826:    Output Parameter:
2827: .  ctx - the context

2829:    Options Database Key:
2830: +  -ts_monitor_lg_timestep - automatically sets line graph monitor
2831: .  -ts_monitor_lg_solution -
2832: .  -ts_monitor_lg_error -
2833: .  -ts_monitor_lg_ksp_iterations -
2834: .  -ts_monitor_lg_snes_iterations -
2835: -  -lg_indicate_data_points <true,false> - indicate the data points (at each time step) on the plot; default is true

2837:    Notes:
2838:    Use TSMonitorLGCtxDestroy() to destroy.

2840:    Level: intermediate

2842: .keywords: TS, monitor, line graph, residual, seealso

2844: .seealso: TSMonitorLGTimeStep(), TSMonitorSet(), TSMonitorLGSolution(), TSMonitorLGError()

2846: @*/
2847: PetscErrorCode  TSMonitorLGCtxCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscInt howoften,TSMonitorLGCtx *ctx)
2848: {
2849:   PetscDraw      win;

2853:   PetscNew(ctx);
2854:   PetscDrawCreate(comm,host,label,x,y,m,n,&win);
2855:   PetscDrawSetFromOptions(win);
2856:   PetscDrawLGCreate(win,1,&(*ctx)->lg);
2857:   PetscLogObjectParent((PetscObject)(*ctx)->lg,(PetscObject)win);
2858:   PetscDrawLGIndicateDataPoints((*ctx)->lg,PETSC_TRUE);
2859:   PetscDrawLGSetFromOptions((*ctx)->lg);
2860:   (*ctx)->howoften = howoften;
2861:   return(0);
2862: }

2866: PetscErrorCode TSMonitorLGTimeStep(TS ts,PetscInt step,PetscReal ptime,Vec v,void *monctx)
2867: {
2868:   TSMonitorLGCtx ctx = (TSMonitorLGCtx) monctx;
2869:   PetscReal      x   = ptime,y;

2873:   if (!step) {
2874:     PetscDrawAxis axis;
2875:     PetscDrawLGGetAxis(ctx->lg,&axis);
2876:     PetscDrawAxisSetLabels(axis,"Timestep as function of time","Time","Time step");
2877:     PetscDrawLGReset(ctx->lg);
2878:     PetscDrawLGIndicateDataPoints(ctx->lg,PETSC_TRUE);
2879:   }
2880:   TSGetTimeStep(ts,&y);
2881:   PetscDrawLGAddPoint(ctx->lg,&x,&y);
2882:   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
2883:     PetscDrawLGDraw(ctx->lg);
2884:   }
2885:   return(0);
2886: }

2890: /*@C
2891:    TSMonitorLGCtxDestroy - Destroys a line graph context that was created
2892:    with TSMonitorLGCtxCreate().

2894:    Collective on TSMonitorLGCtx

2896:    Input Parameter:
2897: .  ctx - the monitor context

2899:    Level: intermediate

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

2903: .seealso: TSMonitorLGCtxCreate(),  TSMonitorSet(), TSMonitorLGTimeStep();
2904: @*/
2905: PetscErrorCode  TSMonitorLGCtxDestroy(TSMonitorLGCtx *ctx)
2906: {
2907:   PetscDraw      draw;

2911:   PetscDrawLGGetDraw((*ctx)->lg,&draw);
2912:   PetscDrawDestroy(&draw);
2913:   PetscDrawLGDestroy(&(*ctx)->lg);
2914:   PetscFree(*ctx);
2915:   return(0);
2916: }

2920: /*@
2921:    TSGetTime - Gets the time of the most recently completed step.

2923:    Not Collective

2925:    Input Parameter:
2926: .  ts - the TS context obtained from TSCreate()

2928:    Output Parameter:
2929: .  t  - the current time

2931:    Level: beginner

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

2937: .seealso: TSSetInitialTimeStep(), TSGetTimeStep()

2939: .keywords: TS, get, time
2940: @*/
2941: PetscErrorCode  TSGetTime(TS ts,PetscReal *t)
2942: {
2946:   *t = ts->ptime;
2947:   return(0);
2948: }

2952: /*@
2953:    TSSetTime - Allows one to reset the time.

2955:    Logically Collective on TS

2957:    Input Parameters:
2958: +  ts - the TS context obtained from TSCreate()
2959: -  time - the time

2961:    Level: intermediate

2963: .seealso: TSGetTime(), TSSetDuration()

2965: .keywords: TS, set, time
2966: @*/
2967: PetscErrorCode  TSSetTime(TS ts, PetscReal t)
2968: {
2972:   ts->ptime = t;
2973:   return(0);
2974: }

2978: /*@C
2979:    TSSetOptionsPrefix - Sets the prefix used for searching for all
2980:    TS options in the database.

2982:    Logically Collective on TS

2984:    Input Parameter:
2985: +  ts     - The TS context
2986: -  prefix - The prefix to prepend to all option names

2988:    Notes:
2989:    A hyphen (-) must NOT be given at the beginning of the prefix name.
2990:    The first character of all runtime options is AUTOMATICALLY the
2991:    hyphen.

2993:    Level: advanced

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

2997: .seealso: TSSetFromOptions()

2999: @*/
3000: PetscErrorCode  TSSetOptionsPrefix(TS ts,const char prefix[])
3001: {
3003:   SNES           snes;

3007:   PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);
3008:   TSGetSNES(ts,&snes);
3009:   SNESSetOptionsPrefix(snes,prefix);
3010:   return(0);
3011: }


3016: /*@C
3017:    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
3018:    TS options in the database.

3020:    Logically Collective on TS

3022:    Input Parameter:
3023: +  ts     - The TS context
3024: -  prefix - The prefix to prepend to all option names

3026:    Notes:
3027:    A hyphen (-) must NOT be given at the beginning of the prefix name.
3028:    The first character of all runtime options is AUTOMATICALLY the
3029:    hyphen.

3031:    Level: advanced

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

3035: .seealso: TSGetOptionsPrefix()

3037: @*/
3038: PetscErrorCode  TSAppendOptionsPrefix(TS ts,const char prefix[])
3039: {
3041:   SNES           snes;

3045:   PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);
3046:   TSGetSNES(ts,&snes);
3047:   SNESAppendOptionsPrefix(snes,prefix);
3048:   return(0);
3049: }

3053: /*@C
3054:    TSGetOptionsPrefix - Sets the prefix used for searching for all
3055:    TS options in the database.

3057:    Not Collective

3059:    Input Parameter:
3060: .  ts - The TS context

3062:    Output Parameter:
3063: .  prefix - A pointer to the prefix string used

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

3068:    Level: intermediate

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

3072: .seealso: TSAppendOptionsPrefix()
3073: @*/
3074: PetscErrorCode  TSGetOptionsPrefix(TS ts,const char *prefix[])
3075: {

3081:   PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);
3082:   return(0);
3083: }

3087: /*@C
3088:    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.

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

3092:    Input Parameter:
3093: .  ts  - The TS context obtained from TSCreate()

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

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

3103:    Level: intermediate

3105: .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()

3107: .keywords: TS, timestep, get, matrix, Jacobian
3108: @*/
3109: PetscErrorCode  TSGetRHSJacobian(TS ts,Mat *Amat,Mat *Pmat,TSRHSJacobian *func,void **ctx)
3110: {
3112:   SNES           snes;
3113:   DM             dm;

3116:   TSGetSNES(ts,&snes);
3117:   SNESGetJacobian(snes,Amat,Pmat,NULL,NULL);
3118:   TSGetDM(ts,&dm);
3119:   DMTSGetRHSJacobian(dm,func,ctx);
3120:   return(0);
3121: }

3125: /*@C
3126:    TSGetIJacobian - Returns the implicit Jacobian at the present timestep.

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

3130:    Input Parameter:
3131: .  ts  - The TS context obtained from TSCreate()

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

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

3141:    Level: advanced

3143: .seealso: TSGetTimeStep(), TSGetRHSJacobian(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()

3145: .keywords: TS, timestep, get, matrix, Jacobian
3146: @*/
3147: PetscErrorCode  TSGetIJacobian(TS ts,Mat *Amat,Mat *Pmat,TSIJacobian *f,void **ctx)
3148: {
3150:   SNES           snes;
3151:   DM             dm;

3154:   TSGetSNES(ts,&snes);
3155:   SNESSetUpMatrices(snes);
3156:   SNESGetJacobian(snes,Amat,Pmat,NULL,NULL);
3157:   TSGetDM(ts,&dm);
3158:   DMTSGetIJacobian(dm,f,ctx);
3159:   return(0);
3160: }


3165: /*@C
3166:    TSMonitorDrawSolution - Monitors progress of the TS solvers by calling
3167:    VecView() for the solution at each timestep

3169:    Collective on TS

3171:    Input Parameters:
3172: +  ts - the TS context
3173: .  step - current time-step
3174: .  ptime - current time
3175: -  dummy - either a viewer or NULL

3177:    Options Database:
3178: .   -ts_monitor_draw_solution_initial - show initial solution as well as current solution

3180:    Notes: the initial solution and current solution are not displayed with a common axis scaling so generally the option -ts_monitor_draw_solution_initial
3181:        will look bad

3183:    Level: intermediate

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

3187: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
3188: @*/
3189: PetscErrorCode  TSMonitorDrawSolution(TS ts,PetscInt step,PetscReal ptime,Vec u,void *dummy)
3190: {
3191:   PetscErrorCode   ierr;
3192:   TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy;
3193:   PetscDraw        draw;

3196:   if (!step && ictx->showinitial) {
3197:     if (!ictx->initialsolution) {
3198:       VecDuplicate(u,&ictx->initialsolution);
3199:     }
3200:     VecCopy(u,ictx->initialsolution);
3201:   }
3202:   if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) return(0);

3204:   if (ictx->showinitial) {
3205:     PetscReal pause;
3206:     PetscViewerDrawGetPause(ictx->viewer,&pause);
3207:     PetscViewerDrawSetPause(ictx->viewer,0.0);
3208:     VecView(ictx->initialsolution,ictx->viewer);
3209:     PetscViewerDrawSetPause(ictx->viewer,pause);
3210:     PetscViewerDrawSetHold(ictx->viewer,PETSC_TRUE);
3211:   }
3212:   VecView(u,ictx->viewer);
3213:   if (ictx->showtimestepandtime) {
3214:     PetscReal xl,yl,xr,yr,tw,w,h;
3215:     char      time[32];
3216:     size_t    len;

3218:     PetscViewerDrawGetDraw(ictx->viewer,0,&draw);
3219:     PetscSNPrintf(time,32,"Timestep %d Time %g",(int)step,(double)ptime);
3220:     PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);
3221:      PetscStrlen(time,&len);
3222:     PetscDrawStringGetSize(draw,&tw,NULL);
3223:     w    = xl + .5*(xr - xl) - .5*len*tw;
3224:     h    = yl + .95*(yr - yl);
3225:     PetscDrawString(draw,w,h,PETSC_DRAW_BLACK,time);
3226:     PetscDrawFlush(draw);
3227:   }

3229:   if (ictx->showinitial) {
3230:     PetscViewerDrawSetHold(ictx->viewer,PETSC_FALSE);
3231:   }
3232:   return(0);
3233: }

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

3240:    Collective on TS

3242:    Input Parameters:
3243: +  ts - the TS context
3244: .  step - current time-step
3245: .  ptime - current time
3246: -  dummy - either a viewer or NULL

3248:    Level: intermediate

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

3252: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
3253: @*/
3254: PetscErrorCode  TSMonitorDrawSolutionPhase(TS ts,PetscInt step,PetscReal ptime,Vec u,void *dummy)
3255: {
3256:   PetscErrorCode    ierr;
3257:   TSMonitorDrawCtx  ictx = (TSMonitorDrawCtx)dummy;
3258:   PetscDraw         draw;
3259:   MPI_Comm          comm;
3260:   PetscInt          n;
3261:   PetscMPIInt       size;
3262:   PetscReal         xl,yl,xr,yr,tw,w,h;
3263:   char              time[32];
3264:   size_t            len;
3265:   const PetscScalar *U;

3268:   PetscObjectGetComm((PetscObject)ts,&comm);
3269:   MPI_Comm_size(comm,&size);
3270:   if (size != 1) SETERRQ(comm,PETSC_ERR_SUP,"Only allowed for sequential runs");
3271:   VecGetSize(u,&n);
3272:   if (n != 2) SETERRQ(comm,PETSC_ERR_SUP,"Only for ODEs with two unknowns");

3274:   PetscViewerDrawGetDraw(ictx->viewer,0,&draw);

3276:   VecGetArrayRead(u,&U);
3277:   PetscDrawAxisGetLimits(ictx->axis,&xl,&xr,&yl,&yr);
3278:   if ((PetscRealPart(U[0]) < xl) || (PetscRealPart(U[1]) < yl) || (PetscRealPart(U[0]) > xr) || (PetscRealPart(U[1]) > yr)) {
3279:       VecRestoreArrayRead(u,&U);
3280:       return(0);
3281:   }
3282:   if (!step) ictx->color++;
3283:   PetscDrawPoint(draw,PetscRealPart(U[0]),PetscRealPart(U[1]),ictx->color);
3284:   VecRestoreArrayRead(u,&U);

3286:   if (ictx->showtimestepandtime) {
3287:     PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);
3288:     PetscSNPrintf(time,32,"Timestep %d Time %g",(int)step,(double)ptime);
3289:     PetscStrlen(time,&len);
3290:     PetscDrawStringGetSize(draw,&tw,NULL);
3291:     w    = xl + .5*(xr - xl) - .5*len*tw;
3292:     h    = yl + .95*(yr - yl);
3293:     PetscDrawString(draw,w,h,PETSC_DRAW_BLACK,time);
3294:   }
3295:   PetscDrawFlush(draw);
3296:   return(0);
3297: }


3302: /*@C
3303:    TSMonitorDrawCtxDestroy - Destroys the monitor context for TSMonitorDrawSolution()

3305:    Collective on TS

3307:    Input Parameters:
3308: .    ctx - the monitor context

3310:    Level: intermediate

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

3314: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorDrawSolution(), TSMonitorDrawError()
3315: @*/
3316: PetscErrorCode  TSMonitorDrawCtxDestroy(TSMonitorDrawCtx *ictx)
3317: {

3321:   PetscDrawAxisDestroy(&(*ictx)->axis);
3322:   PetscViewerDestroy(&(*ictx)->viewer);
3323:   VecDestroy(&(*ictx)->initialsolution);
3324:   PetscFree(*ictx);
3325:   return(0);
3326: }

3330: /*@C
3331:    TSMonitorDrawCtxCreate - Creates the monitor context for TSMonitorDrawCtx

3333:    Collective on TS

3335:    Input Parameter:
3336: .    ts - time-step context

3338:    Output Patameter:
3339: .    ctx - the monitor context

3341:    Options Database:
3342: .   -ts_monitor_draw_solution_initial - show initial solution as well as current solution

3344:    Level: intermediate

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

3348: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorDrawCtx()
3349: @*/
3350: PetscErrorCode  TSMonitorDrawCtxCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscInt howoften,TSMonitorDrawCtx *ctx)
3351: {
3352:   PetscErrorCode   ierr;

3355:   PetscNew(ctx);
3356:   PetscViewerDrawOpen(comm,host,label,x,y,m,n,&(*ctx)->viewer);
3357:   PetscViewerSetFromOptions((*ctx)->viewer);

3359:   (*ctx)->howoften    = howoften;
3360:   (*ctx)->showinitial = PETSC_FALSE;
3361:   PetscOptionsGetBool(NULL,"-ts_monitor_draw_solution_initial",&(*ctx)->showinitial,NULL);

3363:   (*ctx)->showtimestepandtime = PETSC_FALSE;
3364:   PetscOptionsGetBool(NULL,"-ts_monitor_draw_solution_show_time",&(*ctx)->showtimestepandtime,NULL);
3365:   (*ctx)->color = PETSC_DRAW_WHITE;
3366:   return(0);
3367: }

3371: /*@C
3372:    TSMonitorDrawError - Monitors progress of the TS solvers by calling
3373:    VecView() for the error at each timestep

3375:    Collective on TS

3377:    Input Parameters:
3378: +  ts - the TS context
3379: .  step - current time-step
3380: .  ptime - current time
3381: -  dummy - either a viewer or NULL

3383:    Level: intermediate

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

3387: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
3388: @*/
3389: PetscErrorCode  TSMonitorDrawError(TS ts,PetscInt step,PetscReal ptime,Vec u,void *dummy)
3390: {
3391:   PetscErrorCode   ierr;
3392:   TSMonitorDrawCtx ctx    = (TSMonitorDrawCtx)dummy;
3393:   PetscViewer      viewer = ctx->viewer;
3394:   Vec              work;

3397:   if (!(((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason))) return(0);
3398:   VecDuplicate(u,&work);
3399:   TSComputeSolutionFunction(ts,ptime,work);
3400:   VecAXPY(work,-1.0,u);
3401:   VecView(work,viewer);
3402:   VecDestroy(&work);
3403:   return(0);
3404: }

3406: #include <petsc-private/dmimpl.h>
3409: /*@
3410:    TSSetDM - Sets the DM that may be used by some preconditioners

3412:    Logically Collective on TS and DM

3414:    Input Parameters:
3415: +  ts - the preconditioner context
3416: -  dm - the dm

3418:    Level: intermediate


3421: .seealso: TSGetDM(), SNESSetDM(), SNESGetDM()
3422: @*/
3423: PetscErrorCode  TSSetDM(TS ts,DM dm)
3424: {
3426:   SNES           snes;
3427:   DMTS           tsdm;

3431:   PetscObjectReference((PetscObject)dm);
3432:   if (ts->dm) {               /* Move the DMTS context over to the new DM unless the new DM already has one */
3433:     if (ts->dm->dmts && !dm->dmts) {
3434:       DMCopyDMTS(ts->dm,dm);
3435:       DMGetDMTS(ts->dm,&tsdm);
3436:       if (tsdm->originaldm == ts->dm) { /* Grant write privileges to the replacement DM */
3437:         tsdm->originaldm = dm;
3438:       }
3439:     }
3440:     DMDestroy(&ts->dm);
3441:   }
3442:   ts->dm = dm;

3444:   TSGetSNES(ts,&snes);
3445:   SNESSetDM(snes,dm);
3446:   return(0);
3447: }

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

3454:    Not Collective

3456:    Input Parameter:
3457: . ts - the preconditioner context

3459:    Output Parameter:
3460: .  dm - the dm

3462:    Level: intermediate


3465: .seealso: TSSetDM(), SNESSetDM(), SNESGetDM()
3466: @*/
3467: PetscErrorCode  TSGetDM(TS ts,DM *dm)
3468: {

3473:   if (!ts->dm) {
3474:     DMShellCreate(PetscObjectComm((PetscObject)ts),&ts->dm);
3475:     if (ts->snes) {SNESSetDM(ts->snes,ts->dm);}
3476:   }
3477:   *dm = ts->dm;
3478:   return(0);
3479: }

3483: /*@
3484:    SNESTSFormFunction - Function to evaluate nonlinear residual

3486:    Logically Collective on SNES

3488:    Input Parameter:
3489: + snes - nonlinear solver
3490: . U - the current state at which to evaluate the residual
3491: - ctx - user context, must be a TS

3493:    Output Parameter:
3494: . F - the nonlinear residual

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

3500:    Level: advanced

3502: .seealso: SNESSetFunction(), MatFDColoringSetFunction()
3503: @*/
3504: PetscErrorCode  SNESTSFormFunction(SNES snes,Vec U,Vec F,void *ctx)
3505: {
3506:   TS             ts = (TS)ctx;

3514:   (ts->ops->snesfunction)(snes,U,F,ts);
3515:   return(0);
3516: }

3520: /*@
3521:    SNESTSFormJacobian - Function to evaluate the Jacobian

3523:    Collective on SNES

3525:    Input Parameter:
3526: + snes - nonlinear solver
3527: . U - the current state at which to evaluate the residual
3528: - ctx - user context, must be a TS

3530:    Output Parameter:
3531: + A - the Jacobian
3532: . B - the preconditioning matrix (may be the same as A)
3533: - flag - indicates any structure change in the matrix

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

3538:    Level: developer

3540: .seealso: SNESSetJacobian()
3541: @*/
3542: PetscErrorCode  SNESTSFormJacobian(SNES snes,Vec U,Mat A,Mat B,void *ctx)
3543: {
3544:   TS             ts = (TS)ctx;

3555:   (ts->ops->snesjacobian)(snes,U,A,B,ts);
3556:   return(0);
3557: }

3561: /*@C
3562:    TSComputeRHSFunctionLinear - Evaluate the right hand side via the user-provided Jacobian, for linear problems only

3564:    Collective on TS

3566:    Input Arguments:
3567: +  ts - time stepping context
3568: .  t - time at which to evaluate
3569: .  U - state at which to evaluate
3570: -  ctx - context

3572:    Output Arguments:
3573: .  F - right hand side

3575:    Level: intermediate

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

3581: .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSJacobianConstant()
3582: @*/
3583: PetscErrorCode TSComputeRHSFunctionLinear(TS ts,PetscReal t,Vec U,Vec F,void *ctx)
3584: {
3586:   Mat            Arhs,Brhs;

3589:   TSGetRHSMats_Private(ts,&Arhs,&Brhs);
3590:   TSComputeRHSJacobian(ts,t,U,Arhs,Brhs);
3591:   MatMult(Arhs,U,F);
3592:   return(0);
3593: }

3597: /*@C
3598:    TSComputeRHSJacobianConstant - Reuses a Jacobian that is time-independent.

3600:    Collective on TS

3602:    Input Arguments:
3603: +  ts - time stepping context
3604: .  t - time at which to evaluate
3605: .  U - state at which to evaluate
3606: -  ctx - context

3608:    Output Arguments:
3609: +  A - pointer to operator
3610: .  B - pointer to preconditioning matrix
3611: -  flg - matrix structure flag

3613:    Level: intermediate

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

3618: .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSFunctionLinear()
3619: @*/
3620: PetscErrorCode TSComputeRHSJacobianConstant(TS ts,PetscReal t,Vec U,Mat A,Mat B,void *ctx)
3621: {
3623:   return(0);
3624: }

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

3631:    Collective on TS

3633:    Input Arguments:
3634: +  ts - time stepping context
3635: .  t - time at which to evaluate
3636: .  U - state at which to evaluate
3637: .  Udot - time derivative of state vector
3638: -  ctx - context

3640:    Output Arguments:
3641: .  F - left hand side

3643:    Level: intermediate

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

3651: .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIJacobianConstant()
3652: @*/
3653: PetscErrorCode TSComputeIFunctionLinear(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void *ctx)
3654: {
3656:   Mat            A,B;

3659:   TSGetIJacobian(ts,&A,&B,NULL,NULL);
3660:   TSComputeIJacobian(ts,t,U,Udot,1.0,A,B,PETSC_TRUE);
3661:   MatMult(A,Udot,F);
3662:   return(0);
3663: }

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

3670:    Collective on TS

3672:    Input Arguments:
3673: +  ts - time stepping context
3674: .  t - time at which to evaluate
3675: .  U - state at which to evaluate
3676: .  Udot - time derivative of state vector
3677: .  shift - shift to apply
3678: -  ctx - context

3680:    Output Arguments:
3681: +  A - pointer to operator
3682: .  B - pointer to preconditioning matrix
3683: -  flg - matrix structure flag

3685:    Level: advanced

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

3690:    It is only appropriate for problems of the form

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

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

3698: $    shift*M + J

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

3703: .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIFunctionLinear()
3704: @*/
3705: PetscErrorCode TSComputeIJacobianConstant(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal shift,Mat A,Mat B,void *ctx)
3706: {

3710:   MatScale(A, shift / ts->ijacobian.shift);
3711:   ts->ijacobian.shift = shift;
3712:   return(0);
3713: }

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

3720:    Not Collective

3722:    Input Parameter:
3723: .  ts - the TS context

3725:    Output Parameter:
3726: .  equation_type - see TSEquationType

3728:    Level: beginner

3730: .keywords: TS, equation type

3732: .seealso: TSSetEquationType(), TSEquationType
3733: @*/
3734: PetscErrorCode  TSGetEquationType(TS ts,TSEquationType *equation_type)
3735: {
3739:   *equation_type = ts->equation_type;
3740:   return(0);
3741: }

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

3748:    Not Collective

3750:    Input Parameter:
3751: +  ts - the TS context
3752: .  equation_type - see TSEquationType

3754:    Level: advanced

3756: .keywords: TS, equation type

3758: .seealso: TSGetEquationType(), TSEquationType
3759: @*/
3760: PetscErrorCode  TSSetEquationType(TS ts,TSEquationType equation_type)
3761: {
3764:   ts->equation_type = equation_type;
3765:   return(0);
3766: }

3770: /*@
3771:    TSGetConvergedReason - Gets the reason the TS iteration was stopped.

3773:    Not Collective

3775:    Input Parameter:
3776: .  ts - the TS context

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

3782:    Level: beginner

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

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

3789: .seealso: TSSetConvergenceTest(), TSConvergedReason
3790: @*/
3791: PetscErrorCode  TSGetConvergedReason(TS ts,TSConvergedReason *reason)
3792: {
3796:   *reason = ts->reason;
3797:   return(0);
3798: }

3802: /*@
3803:    TSSetConvergedReason - Sets the reason for handling the convergence of TSSolve.

3805:    Not Collective

3807:    Input Parameter:
3808: +  ts - the TS context
3809: .  reason - negative value indicates diverged, positive value converged, see TSConvergedReason or the
3810:             manual pages for the individual convergence tests for complete lists

3812:    Level: advanced

3814:    Notes:
3815:    Can only be called during TSSolve() is active.

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

3819: .seealso: TSConvergedReason
3820: @*/
3821: PetscErrorCode  TSSetConvergedReason(TS ts,TSConvergedReason reason)
3822: {
3825:   ts->reason = reason;
3826:   return(0);
3827: }

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

3834:    Not Collective

3836:    Input Parameter:
3837: .  ts - the TS context

3839:    Output Parameter:
3840: .  ftime - the final time. This time should correspond to the final time set with TSSetDuration()

3842:    Level: beginner

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

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

3849: .seealso: TSSetConvergenceTest(), TSConvergedReason
3850: @*/
3851: PetscErrorCode  TSGetSolveTime(TS ts,PetscReal *ftime)
3852: {
3856:   *ftime = ts->solvetime;
3857:   return(0);
3858: }

3862: /*@
3863:    TSGetSNESIterations - Gets the total number of nonlinear iterations
3864:    used by the time integrator.

3866:    Not Collective

3868:    Input Parameter:
3869: .  ts - TS context

3871:    Output Parameter:
3872: .  nits - number of nonlinear iterations

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

3877:    Level: intermediate

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

3881: .seealso:  TSGetKSPIterations()
3882: @*/
3883: PetscErrorCode TSGetSNESIterations(TS ts,PetscInt *nits)
3884: {
3888:   *nits = ts->snes_its;
3889:   return(0);
3890: }

3894: /*@
3895:    TSGetKSPIterations - Gets the total number of linear iterations
3896:    used by the time integrator.

3898:    Not Collective

3900:    Input Parameter:
3901: .  ts - TS context

3903:    Output Parameter:
3904: .  lits - number of linear iterations

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

3909:    Level: intermediate

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

3913: .seealso:  TSGetSNESIterations(), SNESGetKSPIterations()
3914: @*/
3915: PetscErrorCode TSGetKSPIterations(TS ts,PetscInt *lits)
3916: {
3920:   *lits = ts->ksp_its;
3921:   return(0);
3922: }

3926: /*@
3927:    TSGetStepRejections - Gets the total number of rejected steps.

3929:    Not Collective

3931:    Input Parameter:
3932: .  ts - TS context

3934:    Output Parameter:
3935: .  rejects - number of steps rejected

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

3940:    Level: intermediate

3942: .keywords: TS, get, number

3944: .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetSNESFailures(), TSSetMaxSNESFailures(), TSSetErrorIfStepFails()
3945: @*/
3946: PetscErrorCode TSGetStepRejections(TS ts,PetscInt *rejects)
3947: {
3951:   *rejects = ts->reject;
3952:   return(0);
3953: }

3957: /*@
3958:    TSGetSNESFailures - Gets the total number of failed SNES solves

3960:    Not Collective

3962:    Input Parameter:
3963: .  ts - TS context

3965:    Output Parameter:
3966: .  fails - number of failed nonlinear solves

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

3971:    Level: intermediate

3973: .keywords: TS, get, number

3975: .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSSetMaxSNESFailures()
3976: @*/
3977: PetscErrorCode TSGetSNESFailures(TS ts,PetscInt *fails)
3978: {
3982:   *fails = ts->num_snes_failures;
3983:   return(0);
3984: }

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

3991:    Not Collective

3993:    Input Parameter:
3994: +  ts - TS context
3995: -  rejects - maximum number of rejected steps, pass -1 for unlimited

3997:    Notes:
3998:    The counter is reset to zero for each step

4000:    Options Database Key:
4001:  .  -ts_max_reject - Maximum number of step rejections before a step fails

4003:    Level: intermediate

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

4007: .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxSNESFailures(), TSGetStepRejections(), TSGetSNESFailures(), TSSetErrorIfStepFails(), TSGetConvergedReason()
4008: @*/
4009: PetscErrorCode TSSetMaxStepRejections(TS ts,PetscInt rejects)
4010: {
4013:   ts->max_reject = rejects;
4014:   return(0);
4015: }

4019: /*@
4020:    TSSetMaxSNESFailures - Sets the maximum number of failed SNES solves

4022:    Not Collective

4024:    Input Parameter:
4025: +  ts - TS context
4026: -  fails - maximum number of failed nonlinear solves, pass -1 for unlimited

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

4031:    Options Database Key:
4032:  .  -ts_max_snes_failures - Maximum number of nonlinear solve failures

4034:    Level: intermediate

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

4038: .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSGetSNESFailures(), SNESGetConvergedReason(), TSGetConvergedReason()
4039: @*/
4040: PetscErrorCode TSSetMaxSNESFailures(TS ts,PetscInt fails)
4041: {
4044:   ts->max_snes_failures = fails;
4045:   return(0);
4046: }

4050: /*@
4051:    TSSetErrorIfStepFails - Error if no step succeeds

4053:    Not Collective

4055:    Input Parameter:
4056: +  ts - TS context
4057: -  err - PETSC_TRUE to error if no step succeeds, PETSC_FALSE to return without failure

4059:    Options Database Key:
4060:  .  -ts_error_if_step_fails - Error if no step succeeds

4062:    Level: intermediate

4064: .keywords: TS, set, error

4066: .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSGetSNESFailures(), TSSetErrorIfStepFails(), TSGetConvergedReason()
4067: @*/
4068: PetscErrorCode TSSetErrorIfStepFails(TS ts,PetscBool err)
4069: {
4072:   ts->errorifstepfailed = err;
4073:   return(0);
4074: }

4078: /*@C
4079:    TSMonitorSolutionBinary - Monitors progress of the TS solvers by VecView() for the solution at each timestep. Normally the viewer is a binary file

4081:    Collective on TS

4083:    Input Parameters:
4084: +  ts - the TS context
4085: .  step - current time-step
4086: .  ptime - current time
4087: .  u - current state
4088: -  viewer - binary viewer

4090:    Level: intermediate

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

4094: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
4095: @*/
4096: PetscErrorCode  TSMonitorSolutionBinary(TS ts,PetscInt step,PetscReal ptime,Vec u,void *viewer)
4097: {
4099:   PetscViewer    v = (PetscViewer)viewer;

4102:   VecView(u,v);
4103:   return(0);
4104: }

4108: /*@C
4109:    TSMonitorSolutionVTK - Monitors progress of the TS solvers by VecView() for the solution at each timestep.

4111:    Collective on TS

4113:    Input Parameters:
4114: +  ts - the TS context
4115: .  step - current time-step
4116: .  ptime - current time
4117: .  u - current state
4118: -  filenametemplate - string containing a format specifier for the integer time step (e.g. %03D)

4120:    Level: intermediate

4122:    Notes:
4123:    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.
4124:    These are named according to the file name template.

4126:    This function is normally passed as an argument to TSMonitorSet() along with TSMonitorSolutionVTKDestroy().

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

4130: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
4131: @*/
4132: PetscErrorCode TSMonitorSolutionVTK(TS ts,PetscInt step,PetscReal ptime,Vec u,void *filenametemplate)
4133: {
4135:   char           filename[PETSC_MAX_PATH_LEN];
4136:   PetscViewer    viewer;

4139:   PetscSNPrintf(filename,sizeof(filename),(const char*)filenametemplate,step);
4140:   PetscViewerVTKOpen(PetscObjectComm((PetscObject)ts),filename,FILE_MODE_WRITE,&viewer);
4141:   VecView(u,viewer);
4142:   PetscViewerDestroy(&viewer);
4143:   return(0);
4144: }

4148: /*@C
4149:    TSMonitorSolutionVTKDestroy - Destroy context for monitoring

4151:    Collective on TS

4153:    Input Parameters:
4154: .  filenametemplate - string containing a format specifier for the integer time step (e.g. %03D)

4156:    Level: intermediate

4158:    Note:
4159:    This function is normally passed to TSMonitorSet() along with TSMonitorSolutionVTK().

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

4163: .seealso: TSMonitorSet(), TSMonitorSolutionVTK()
4164: @*/
4165: PetscErrorCode TSMonitorSolutionVTKDestroy(void *filenametemplate)
4166: {

4170:   PetscFree(*(char**)filenametemplate);
4171:   return(0);
4172: }

4176: /*@
4177:    TSGetAdapt - Get the adaptive controller context for the current method

4179:    Collective on TS if controller has not been created yet

4181:    Input Arguments:
4182: .  ts - time stepping context

4184:    Output Arguments:
4185: .  adapt - adaptive controller

4187:    Level: intermediate

4189: .seealso: TSAdapt, TSAdaptSetType(), TSAdaptChoose()
4190: @*/
4191: PetscErrorCode TSGetAdapt(TS ts,TSAdapt *adapt)
4192: {

4198:   if (!ts->adapt) {
4199:     TSAdaptCreate(PetscObjectComm((PetscObject)ts),&ts->adapt);
4200:     PetscLogObjectParent((PetscObject)ts,(PetscObject)ts->adapt);
4201:     PetscObjectIncrementTabLevel((PetscObject)ts->adapt,(PetscObject)ts,1);
4202:   }
4203:   *adapt = ts->adapt;
4204:   return(0);
4205: }

4209: /*@
4210:    TSSetTolerances - Set tolerances for local truncation error when using adaptive controller

4212:    Logically Collective

4214:    Input Arguments:
4215: +  ts - time integration context
4216: .  atol - scalar absolute tolerances, PETSC_DECIDE to leave current value
4217: .  vatol - vector of absolute tolerances or NULL, used in preference to atol if present
4218: .  rtol - scalar relative tolerances, PETSC_DECIDE to leave current value
4219: -  vrtol - vector of relative tolerances or NULL, used in preference to atol if present

4221:    Options Database keys:
4222: +  -ts_rtol <rtol> - relative tolerance for local truncation error
4223: -  -ts_atol <atol> Absolute tolerance for local truncation error

4225:    Level: beginner

4227: .seealso: TS, TSAdapt, TSVecNormWRMS(), TSGetTolerances()
4228: @*/
4229: PetscErrorCode TSSetTolerances(TS ts,PetscReal atol,Vec vatol,PetscReal rtol,Vec vrtol)
4230: {

4234:   if (atol != PETSC_DECIDE && atol != PETSC_DEFAULT) ts->atol = atol;
4235:   if (vatol) {
4236:     PetscObjectReference((PetscObject)vatol);
4237:     VecDestroy(&ts->vatol);

4239:     ts->vatol = vatol;
4240:   }
4241:   if (rtol != PETSC_DECIDE && rtol != PETSC_DEFAULT) ts->rtol = rtol;
4242:   if (vrtol) {
4243:     PetscObjectReference((PetscObject)vrtol);
4244:     VecDestroy(&ts->vrtol);

4246:     ts->vrtol = vrtol;
4247:   }
4248:   return(0);
4249: }

4253: /*@
4254:    TSGetTolerances - Get tolerances for local truncation error when using adaptive controller

4256:    Logically Collective

4258:    Input Arguments:
4259: .  ts - time integration context

4261:    Output Arguments:
4262: +  atol - scalar absolute tolerances, NULL to ignore
4263: .  vatol - vector of absolute tolerances, NULL to ignore
4264: .  rtol - scalar relative tolerances, NULL to ignore
4265: -  vrtol - vector of relative tolerances, NULL to ignore

4267:    Level: beginner

4269: .seealso: TS, TSAdapt, TSVecNormWRMS(), TSSetTolerances()
4270: @*/
4271: PetscErrorCode TSGetTolerances(TS ts,PetscReal *atol,Vec *vatol,PetscReal *rtol,Vec *vrtol)
4272: {
4274:   if (atol)  *atol  = ts->atol;
4275:   if (vatol) *vatol = ts->vatol;
4276:   if (rtol)  *rtol  = ts->rtol;
4277:   if (vrtol) *vrtol = ts->vrtol;
4278:   return(0);
4279: }

4283: /*@
4284:    TSErrorNormWRMS - compute a weighted norm of the difference between a vector and the current state

4286:    Collective on TS

4288:    Input Arguments:
4289: +  ts - time stepping context
4290: -  Y - state vector to be compared to ts->vec_sol

4292:    Output Arguments:
4293: .  norm - weighted norm, a value of 1.0 is considered small

4295:    Level: developer

4297: .seealso: TSSetTolerances()
4298: @*/
4299: PetscErrorCode TSErrorNormWRMS(TS ts,Vec Y,PetscReal *norm)
4300: {
4301:   PetscErrorCode    ierr;
4302:   PetscInt          i,n,N;
4303:   const PetscScalar *u,*y;
4304:   Vec               U;
4305:   PetscReal         sum,gsum;

4311:   U = ts->vec_sol;
4313:   if (U == Y) SETERRQ(PetscObjectComm((PetscObject)U),PETSC_ERR_ARG_IDN,"Y cannot be the TS solution vector");

4315:   VecGetSize(U,&N);
4316:   VecGetLocalSize(U,&n);
4317:   VecGetArrayRead(U,&u);
4318:   VecGetArrayRead(Y,&y);
4319:   sum  = 0.;
4320:   if (ts->vatol && ts->vrtol) {
4321:     const PetscScalar *atol,*rtol;
4322:     VecGetArrayRead(ts->vatol,&atol);
4323:     VecGetArrayRead(ts->vrtol,&rtol);
4324:     for (i=0; i<n; i++) {
4325:       PetscReal tol = PetscRealPart(atol[i]) + PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
4326:       sum += PetscSqr(PetscAbsScalar(y[i] - u[i]) / tol);
4327:     }
4328:     VecRestoreArrayRead(ts->vatol,&atol);
4329:     VecRestoreArrayRead(ts->vrtol,&rtol);
4330:   } else if (ts->vatol) {       /* vector atol, scalar rtol */
4331:     const PetscScalar *atol;
4332:     VecGetArrayRead(ts->vatol,&atol);
4333:     for (i=0; i<n; i++) {
4334:       PetscReal tol = PetscRealPart(atol[i]) + ts->rtol * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
4335:       sum += PetscSqr(PetscAbsScalar(y[i] - u[i]) / tol);
4336:     }
4337:     VecRestoreArrayRead(ts->vatol,&atol);
4338:   } else if (ts->vrtol) {       /* scalar atol, vector rtol */
4339:     const PetscScalar *rtol;
4340:     VecGetArrayRead(ts->vrtol,&rtol);
4341:     for (i=0; i<n; i++) {
4342:       PetscReal tol = ts->atol + PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
4343:       sum += PetscSqr(PetscAbsScalar(y[i] - u[i]) / tol);
4344:     }
4345:     VecRestoreArrayRead(ts->vrtol,&rtol);
4346:   } else {                      /* scalar atol, scalar rtol */
4347:     for (i=0; i<n; i++) {
4348:       PetscReal tol = ts->atol + ts->rtol * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
4349:       sum += PetscSqr(PetscAbsScalar(y[i] - u[i]) / tol);
4350:     }
4351:   }
4352:   VecRestoreArrayRead(U,&u);
4353:   VecRestoreArrayRead(Y,&y);

4355:   MPI_Allreduce(&sum,&gsum,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)ts));
4356:   *norm = PetscSqrtReal(gsum / N);
4357:   if (PetscIsInfOrNanScalar(*norm)) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
4358:   return(0);
4359: }

4363: /*@
4364:    TSSetCFLTimeLocal - Set the local CFL constraint relative to forward Euler

4366:    Logically Collective on TS

4368:    Input Arguments:
4369: +  ts - time stepping context
4370: -  cfltime - maximum stable time step if using forward Euler (value can be different on each process)

4372:    Note:
4373:    After calling this function, the global CFL time can be obtained by calling TSGetCFLTime()

4375:    Level: intermediate

4377: .seealso: TSGetCFLTime(), TSADAPTCFL
4378: @*/
4379: PetscErrorCode TSSetCFLTimeLocal(TS ts,PetscReal cfltime)
4380: {
4383:   ts->cfltime_local = cfltime;
4384:   ts->cfltime       = -1.;
4385:   return(0);
4386: }

4390: /*@
4391:    TSGetCFLTime - Get the maximum stable time step according to CFL criteria applied to forward Euler

4393:    Collective on TS

4395:    Input Arguments:
4396: .  ts - time stepping context

4398:    Output Arguments:
4399: .  cfltime - maximum stable time step for forward Euler

4401:    Level: advanced

4403: .seealso: TSSetCFLTimeLocal()
4404: @*/
4405: PetscErrorCode TSGetCFLTime(TS ts,PetscReal *cfltime)
4406: {

4410:   if (ts->cfltime < 0) {
4411:     MPI_Allreduce(&ts->cfltime_local,&ts->cfltime,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)ts));
4412:   }
4413:   *cfltime = ts->cfltime;
4414:   return(0);
4415: }

4419: /*@
4420:    TSVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu

4422:    Input Parameters:
4423: .  ts   - the TS context.
4424: .  xl   - lower bound.
4425: .  xu   - upper bound.

4427:    Notes:
4428:    If this routine is not called then the lower and upper bounds are set to
4429:    PETSC_NINFINITY and PETSC_INFINITY respectively during SNESSetUp().

4431:    Level: advanced

4433: @*/
4434: PetscErrorCode TSVISetVariableBounds(TS ts, Vec xl, Vec xu)
4435: {
4437:   SNES           snes;

4440:   TSGetSNES(ts,&snes);
4441:   SNESVISetVariableBounds(snes,xl,xu);
4442:   return(0);
4443: }

4445: #if defined(PETSC_HAVE_MATLAB_ENGINE)
4446: #include <mex.h>

4448: typedef struct {char *funcname; mxArray *ctx;} TSMatlabContext;

4452: /*
4453:    TSComputeFunction_Matlab - Calls the function that has been set with
4454:                          TSSetFunctionMatlab().

4456:    Collective on TS

4458:    Input Parameters:
4459: +  snes - the TS context
4460: -  u - input vector

4462:    Output Parameter:
4463: .  y - function vector, as set by TSSetFunction()

4465:    Notes:
4466:    TSComputeFunction() is typically used within nonlinear solvers
4467:    implementations, so most users would not generally call this routine
4468:    themselves.

4470:    Level: developer

4472: .keywords: TS, nonlinear, compute, function

4474: .seealso: TSSetFunction(), TSGetFunction()
4475: */
4476: PetscErrorCode  TSComputeFunction_Matlab(TS snes,PetscReal time,Vec u,Vec udot,Vec y, void *ctx)
4477: {
4478:   PetscErrorCode  ierr;
4479:   TSMatlabContext *sctx = (TSMatlabContext*)ctx;
4480:   int             nlhs  = 1,nrhs = 7;
4481:   mxArray         *plhs[1],*prhs[7];
4482:   long long int   lx = 0,lxdot = 0,ly = 0,ls = 0;


4492:   PetscMemcpy(&ls,&snes,sizeof(snes));
4493:   PetscMemcpy(&lx,&u,sizeof(u));
4494:   PetscMemcpy(&lxdot,&udot,sizeof(udot));
4495:   PetscMemcpy(&ly,&y,sizeof(u));

4497:   prhs[0] =  mxCreateDoubleScalar((double)ls);
4498:   prhs[1] =  mxCreateDoubleScalar(time);
4499:   prhs[2] =  mxCreateDoubleScalar((double)lx);
4500:   prhs[3] =  mxCreateDoubleScalar((double)lxdot);
4501:   prhs[4] =  mxCreateDoubleScalar((double)ly);
4502:   prhs[5] =  mxCreateString(sctx->funcname);
4503:   prhs[6] =  sctx->ctx;
4504:    mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeFunctionInternal");
4505:    mxGetScalar(plhs[0]);
4506:   mxDestroyArray(prhs[0]);
4507:   mxDestroyArray(prhs[1]);
4508:   mxDestroyArray(prhs[2]);
4509:   mxDestroyArray(prhs[3]);
4510:   mxDestroyArray(prhs[4]);
4511:   mxDestroyArray(prhs[5]);
4512:   mxDestroyArray(plhs[0]);
4513:   return(0);
4514: }


4519: /*
4520:    TSSetFunctionMatlab - Sets the function evaluation routine and function
4521:    vector for use by the TS routines in solving ODEs
4522:    equations from MATLAB. Here the function is a string containing the name of a MATLAB function

4524:    Logically Collective on TS

4526:    Input Parameters:
4527: +  ts - the TS context
4528: -  func - function evaluation routine

4530:    Calling sequence of func:
4531: $    func (TS ts,PetscReal time,Vec u,Vec udot,Vec f,void *ctx);

4533:    Level: beginner

4535: .keywords: TS, nonlinear, set, function

4537: .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
4538: */
4539: PetscErrorCode  TSSetFunctionMatlab(TS ts,const char *func,mxArray *ctx)
4540: {
4541:   PetscErrorCode  ierr;
4542:   TSMatlabContext *sctx;

4545:   /* currently sctx is memory bleed */
4546:   PetscMalloc(sizeof(TSMatlabContext),&sctx);
4547:   PetscStrallocpy(func,&sctx->funcname);
4548:   /*
4549:      This should work, but it doesn't
4550:   sctx->ctx = ctx;
4551:   mexMakeArrayPersistent(sctx->ctx);
4552:   */
4553:   sctx->ctx = mxDuplicateArray(ctx);

4555:   TSSetIFunction(ts,NULL,TSComputeFunction_Matlab,sctx);
4556:   return(0);
4557: }

4561: /*
4562:    TSComputeJacobian_Matlab - Calls the function that has been set with
4563:                          TSSetJacobianMatlab().

4565:    Collective on TS

4567:    Input Parameters:
4568: +  ts - the TS context
4569: .  u - input vector
4570: .  A, B - the matrices
4571: -  ctx - user context

4573:    Level: developer

4575: .keywords: TS, nonlinear, compute, function

4577: .seealso: TSSetFunction(), TSGetFunction()
4578: @*/
4579: PetscErrorCode  TSComputeJacobian_Matlab(TS ts,PetscReal time,Vec u,Vec udot,PetscReal shift,Mat A,Mat B,void *ctx)
4580: {
4581:   PetscErrorCode  ierr;
4582:   TSMatlabContext *sctx = (TSMatlabContext*)ctx;
4583:   int             nlhs  = 2,nrhs = 9;
4584:   mxArray         *plhs[2],*prhs[9];
4585:   long long int   lx = 0,lxdot = 0,lA = 0,ls = 0, lB = 0;


4591:   /* call Matlab function in ctx with arguments u and y */

4593:   PetscMemcpy(&ls,&ts,sizeof(ts));
4594:   PetscMemcpy(&lx,&u,sizeof(u));
4595:   PetscMemcpy(&lxdot,&udot,sizeof(u));
4596:   PetscMemcpy(&lA,A,sizeof(u));
4597:   PetscMemcpy(&lB,B,sizeof(u));

4599:   prhs[0] =  mxCreateDoubleScalar((double)ls);
4600:   prhs[1] =  mxCreateDoubleScalar((double)time);
4601:   prhs[2] =  mxCreateDoubleScalar((double)lx);
4602:   prhs[3] =  mxCreateDoubleScalar((double)lxdot);
4603:   prhs[4] =  mxCreateDoubleScalar((double)shift);
4604:   prhs[5] =  mxCreateDoubleScalar((double)lA);
4605:   prhs[6] =  mxCreateDoubleScalar((double)lB);
4606:   prhs[7] =  mxCreateString(sctx->funcname);
4607:   prhs[8] =  sctx->ctx;
4608:    mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeJacobianInternal");
4609:    mxGetScalar(plhs[0]);
4610:   mxDestroyArray(prhs[0]);
4611:   mxDestroyArray(prhs[1]);
4612:   mxDestroyArray(prhs[2]);
4613:   mxDestroyArray(prhs[3]);
4614:   mxDestroyArray(prhs[4]);
4615:   mxDestroyArray(prhs[5]);
4616:   mxDestroyArray(prhs[6]);
4617:   mxDestroyArray(prhs[7]);
4618:   mxDestroyArray(plhs[0]);
4619:   mxDestroyArray(plhs[1]);
4620:   return(0);
4621: }


4626: /*
4627:    TSSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices
4628:    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

4630:    Logically Collective on TS

4632:    Input Parameters:
4633: +  ts - the TS context
4634: .  A,B - Jacobian matrices
4635: .  func - function evaluation routine
4636: -  ctx - user context

4638:    Calling sequence of func:
4639: $    flag = func (TS ts,PetscReal time,Vec u,Vec udot,Mat A,Mat B,void *ctx);


4642:    Level: developer

4644: .keywords: TS, nonlinear, set, function

4646: .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
4647: */
4648: PetscErrorCode  TSSetJacobianMatlab(TS ts,Mat A,Mat B,const char *func,mxArray *ctx)
4649: {
4650:   PetscErrorCode  ierr;
4651:   TSMatlabContext *sctx;

4654:   /* currently sctx is memory bleed */
4655:   PetscMalloc(sizeof(TSMatlabContext),&sctx);
4656:   PetscStrallocpy(func,&sctx->funcname);
4657:   /*
4658:      This should work, but it doesn't
4659:   sctx->ctx = ctx;
4660:   mexMakeArrayPersistent(sctx->ctx);
4661:   */
4662:   sctx->ctx = mxDuplicateArray(ctx);

4664:   TSSetIJacobian(ts,A,B,TSComputeJacobian_Matlab,sctx);
4665:   return(0);
4666: }

4670: /*
4671:    TSMonitor_Matlab - Calls the function that has been set with TSMonitorSetMatlab().

4673:    Collective on TS

4675: .seealso: TSSetFunction(), TSGetFunction()
4676: @*/
4677: PetscErrorCode  TSMonitor_Matlab(TS ts,PetscInt it, PetscReal time,Vec u, void *ctx)
4678: {
4679:   PetscErrorCode  ierr;
4680:   TSMatlabContext *sctx = (TSMatlabContext*)ctx;
4681:   int             nlhs  = 1,nrhs = 6;
4682:   mxArray         *plhs[1],*prhs[6];
4683:   long long int   lx = 0,ls = 0;


4689:   PetscMemcpy(&ls,&ts,sizeof(ts));
4690:   PetscMemcpy(&lx,&u,sizeof(u));

4692:   prhs[0] =  mxCreateDoubleScalar((double)ls);
4693:   prhs[1] =  mxCreateDoubleScalar((double)it);
4694:   prhs[2] =  mxCreateDoubleScalar((double)time);
4695:   prhs[3] =  mxCreateDoubleScalar((double)lx);
4696:   prhs[4] =  mxCreateString(sctx->funcname);
4697:   prhs[5] =  sctx->ctx;
4698:    mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSMonitorInternal");
4699:    mxGetScalar(plhs[0]);
4700:   mxDestroyArray(prhs[0]);
4701:   mxDestroyArray(prhs[1]);
4702:   mxDestroyArray(prhs[2]);
4703:   mxDestroyArray(prhs[3]);
4704:   mxDestroyArray(prhs[4]);
4705:   mxDestroyArray(plhs[0]);
4706:   return(0);
4707: }


4712: /*
4713:    TSMonitorSetMatlab - Sets the monitor function from Matlab

4715:    Level: developer

4717: .keywords: TS, nonlinear, set, function

4719: .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
4720: */
4721: PetscErrorCode  TSMonitorSetMatlab(TS ts,const char *func,mxArray *ctx)
4722: {
4723:   PetscErrorCode  ierr;
4724:   TSMatlabContext *sctx;

4727:   /* currently sctx is memory bleed */
4728:   PetscMalloc(sizeof(TSMatlabContext),&sctx);
4729:   PetscStrallocpy(func,&sctx->funcname);
4730:   /*
4731:      This should work, but it doesn't
4732:   sctx->ctx = ctx;
4733:   mexMakeArrayPersistent(sctx->ctx);
4734:   */
4735:   sctx->ctx = mxDuplicateArray(ctx);

4737:   TSMonitorSet(ts,TSMonitor_Matlab,sctx,NULL);
4738:   return(0);
4739: }
4740: #endif



4746: /*@C
4747:    TSMonitorLGSolution - Monitors progress of the TS solvers by plotting each component of the solution vector
4748:        in a time based line graph

4750:    Collective on TS

4752:    Input Parameters:
4753: +  ts - the TS context
4754: .  step - current time-step
4755: .  ptime - current time
4756: -  lg - a line graph object

4758:    Level: intermediate

4760:     Notes: each process in a parallel run displays its component solutions in a separate window

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

4764: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
4765: @*/
4766: PetscErrorCode  TSMonitorLGSolution(TS ts,PetscInt step,PetscReal ptime,Vec u,void *dummy)
4767: {
4768:   PetscErrorCode    ierr;
4769:   TSMonitorLGCtx    ctx = (TSMonitorLGCtx)dummy;
4770:   const PetscScalar *yy;
4771:   PetscInt          dim;

4774:   if (!step) {
4775:     PetscDrawAxis axis;
4776:     PetscDrawLGGetAxis(ctx->lg,&axis);
4777:     PetscDrawAxisSetLabels(axis,"Solution as function of time","Time","Solution");
4778:     VecGetLocalSize(u,&dim);
4779:     PetscDrawLGSetDimension(ctx->lg,dim);
4780:     PetscDrawLGReset(ctx->lg);
4781:   }
4782:   VecGetArrayRead(u,&yy);
4783: #if defined(PETSC_USE_COMPLEX)
4784:   {
4785:     PetscReal *yreal;
4786:     PetscInt  i,n;
4787:     VecGetLocalSize(u,&n);
4788:     PetscMalloc1(n,&yreal);
4789:     for (i=0; i<n; i++) yreal[i] = PetscRealPart(yy[i]);
4790:     PetscDrawLGAddCommonPoint(ctx->lg,ptime,yreal);
4791:     PetscFree(yreal);
4792:   }
4793: #else
4794:   PetscDrawLGAddCommonPoint(ctx->lg,ptime,yy);
4795: #endif
4796:   VecRestoreArrayRead(u,&yy);
4797:   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
4798:     PetscDrawLGDraw(ctx->lg);
4799:   }
4800:   return(0);
4801: }

4805: /*@C
4806:    TSMonitorLGError - Monitors progress of the TS solvers by plotting each component of the solution vector
4807:        in a time based line graph

4809:    Collective on TS

4811:    Input Parameters:
4812: +  ts - the TS context
4813: .  step - current time-step
4814: .  ptime - current time
4815: -  lg - a line graph object

4817:    Level: intermediate

4819:    Notes:
4820:    Only for sequential solves.

4822:    The user must provide the solution using TSSetSolutionFunction() to use this monitor.

4824:    Options Database Keys:
4825: .  -ts_monitor_lg_error - create a graphical monitor of error history

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

4829: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSSetSolutionFunction()
4830: @*/
4831: PetscErrorCode  TSMonitorLGError(TS ts,PetscInt step,PetscReal ptime,Vec u,void *dummy)
4832: {
4833:   PetscErrorCode    ierr;
4834:   TSMonitorLGCtx    ctx = (TSMonitorLGCtx)dummy;
4835:   const PetscScalar *yy;
4836:   Vec               y;
4837:   PetscInt          dim;

4840:   if (!step) {
4841:     PetscDrawAxis axis;
4842:     PetscDrawLGGetAxis(ctx->lg,&axis);
4843:     PetscDrawAxisSetLabels(axis,"Error in solution as function of time","Time","Solution");
4844:     VecGetLocalSize(u,&dim);
4845:     PetscDrawLGSetDimension(ctx->lg,dim);
4846:     PetscDrawLGReset(ctx->lg);
4847:   }
4848:   VecDuplicate(u,&y);
4849:   TSComputeSolutionFunction(ts,ptime,y);
4850:   VecAXPY(y,-1.0,u);
4851:   VecGetArrayRead(y,&yy);
4852: #if defined(PETSC_USE_COMPLEX)
4853:   {
4854:     PetscReal *yreal;
4855:     PetscInt  i,n;
4856:     VecGetLocalSize(y,&n);
4857:     PetscMalloc1(n,&yreal);
4858:     for (i=0; i<n; i++) yreal[i] = PetscRealPart(yy[i]);
4859:     PetscDrawLGAddCommonPoint(ctx->lg,ptime,yreal);
4860:     PetscFree(yreal);
4861:   }
4862: #else
4863:   PetscDrawLGAddCommonPoint(ctx->lg,ptime,yy);
4864: #endif
4865:   VecRestoreArrayRead(y,&yy);
4866:   VecDestroy(&y);
4867:   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
4868:     PetscDrawLGDraw(ctx->lg);
4869:   }
4870:   return(0);
4871: }

4875: PetscErrorCode TSMonitorLGSNESIterations(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx)
4876: {
4877:   TSMonitorLGCtx ctx = (TSMonitorLGCtx) monctx;
4878:   PetscReal      x   = ptime,y;
4880:   PetscInt       its;

4883:   if (!n) {
4884:     PetscDrawAxis axis;

4886:     PetscDrawLGGetAxis(ctx->lg,&axis);
4887:     PetscDrawAxisSetLabels(axis,"Nonlinear iterations as function of time","Time","SNES Iterations");
4888:     PetscDrawLGReset(ctx->lg);

4890:     ctx->snes_its = 0;
4891:   }
4892:   TSGetSNESIterations(ts,&its);
4893:   y    = its - ctx->snes_its;
4894:   PetscDrawLGAddPoint(ctx->lg,&x,&y);
4895:   if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) {
4896:     PetscDrawLGDraw(ctx->lg);
4897:   }
4898:   ctx->snes_its = its;
4899:   return(0);
4900: }

4904: PetscErrorCode TSMonitorLGKSPIterations(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx)
4905: {
4906:   TSMonitorLGCtx ctx = (TSMonitorLGCtx) monctx;
4907:   PetscReal      x   = ptime,y;
4909:   PetscInt       its;

4912:   if (!n) {
4913:     PetscDrawAxis axis;

4915:     PetscDrawLGGetAxis(ctx->lg,&axis);
4916:     PetscDrawAxisSetLabels(axis,"Linear iterations as function of time","Time","KSP Iterations");
4917:     PetscDrawLGReset(ctx->lg);

4919:     ctx->ksp_its = 0;
4920:   }
4921:   TSGetKSPIterations(ts,&its);
4922:   y    = its - ctx->ksp_its;
4923:   PetscDrawLGAddPoint(ctx->lg,&x,&y);
4924:   if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) {
4925:     PetscDrawLGDraw(ctx->lg);
4926:   }
4927:   ctx->ksp_its = its;
4928:   return(0);
4929: }

4933: /*@
4934:    TSComputeLinearStability - computes the linear stability function at a point

4936:    Collective on TS and Vec

4938:    Input Parameters:
4939: +  ts - the TS context
4940: -  xr,xi - real and imaginary part of input arguments

4942:    Output Parameters:
4943: .  yr,yi - real and imaginary part of function value

4945:    Level: developer

4947: .keywords: TS, compute

4949: .seealso: TSSetRHSFunction(), TSComputeIFunction()
4950: @*/
4951: PetscErrorCode TSComputeLinearStability(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
4952: {

4957:   if (!ts->ops->linearstability) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Linearized stability function not provided for this method");
4958:   (*ts->ops->linearstability)(ts,xr,xi,yr,yi);
4959:   return(0);
4960: }

4964: /*@
4965:    TSRollBack - Rolls back one time step

4967:    Collective on TS

4969:    Input Parameter:
4970: .  ts - the TS context obtained from TSCreate()

4972:    Level: advanced

4974: .keywords: TS, timestep, rollback

4976: .seealso: TSCreate(), TSSetUp(), TSDestroy(), TSSolve(), TSSetPreStep(), TSSetPreStage(), TSInterpolate()
4977: @*/
4978: PetscErrorCode  TSRollBack(TS ts)
4979: {


4985:   if (!ts->ops->rollback) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRollBack not implemented for type '%s'",((PetscObject)ts)->type_name);
4986:   (*ts->ops->rollback)(ts);
4987:   ts->time_step = ts->ptime - ts->ptime_prev;
4988:   ts->ptime = ts->ptime_prev;
4989:   return(0);
4990: }