Actual source code: ts.c

  1: #include <petsc/private/tsimpl.h>
  2: #include <petscdmshell.h>
  3: #include <petscdmda.h>
  4: #include <petscviewer.h>
  5: #include <petscdraw.h>
  6: #include <petscconvest.h>

  8: #define SkipSmallValue(a,b,tol) if (PetscAbsScalar(a)< tol || PetscAbsScalar(b)< tol) continue;

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

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

 16: static PetscErrorCode TSAdaptSetDefaultType(TSAdapt adapt,TSAdaptType default_type)
 17: {
 20:   if (!((PetscObject)adapt)->type_name) {
 21:     TSAdaptSetType(adapt,default_type);
 22:   }
 23:   return 0;
 24: }

 26: /*@
 27:    TSSetFromOptions - Sets various TS parameters from user options.

 29:    Collective on TS

 31:    Input Parameter:
 32: .  ts - the TS context obtained from TSCreate()

 34:    Options Database Keys:
 35: +  -ts_type <type> - TSEULER, TSBEULER, TSSUNDIALS, TSPSEUDO, TSCN, TSRK, TSTHETA, TSALPHA, TSGLLE, TSSSP, TSGLEE, TSBSYMP, TSIRK
 36: .  -ts_save_trajectory - checkpoint the solution at each time-step
 37: .  -ts_max_time <time> - maximum time to compute to
 38: .  -ts_max_steps <steps> - maximum number of time-steps to take
 39: .  -ts_init_time <time> - initial time to start computation
 40: .  -ts_final_time <time> - final time to compute to (deprecated: use -ts_max_time)
 41: .  -ts_dt <dt> - initial time step
 42: .  -ts_exact_final_time <stepover,interpolate,matchstep> - whether to stop at the exact given final time and how to compute the solution at that time
 43: .  -ts_max_snes_failures <maxfailures> - Maximum number of nonlinear solve failures allowed
 44: .  -ts_max_reject <maxrejects> - Maximum number of step rejections before step fails
 45: .  -ts_error_if_step_fails <true,false> - Error if no step succeeds
 46: .  -ts_rtol <rtol> - relative tolerance for local truncation error
 47: .  -ts_atol <atol> - Absolute tolerance for local truncation error
 48: .  -ts_rhs_jacobian_test_mult -mat_shell_test_mult_view - test the Jacobian at each iteration against finite difference with RHS function
 49: .  -ts_rhs_jacobian_test_mult_transpose -mat_shell_test_mult_transpose_view - test the Jacobian at each iteration against finite difference with RHS function
 50: .  -ts_adjoint_solve <yes,no> - After solving the ODE/DAE solve the adjoint problem (requires -ts_save_trajectory)
 51: .  -ts_fd_color - Use finite differences with coloring to compute IJacobian
 52: .  -ts_monitor - print information at each timestep
 53: .  -ts_monitor_cancel - Cancel all monitors
 54: .  -ts_monitor_lg_solution - Monitor solution graphically
 55: .  -ts_monitor_lg_error - Monitor error graphically
 56: .  -ts_monitor_error - Monitors norm of error
 57: .  -ts_monitor_lg_timestep - Monitor timestep size graphically
 58: .  -ts_monitor_lg_timestep_log - Monitor log timestep size graphically
 59: .  -ts_monitor_lg_snes_iterations - Monitor number nonlinear iterations for each timestep graphically
 60: .  -ts_monitor_lg_ksp_iterations - Monitor number nonlinear iterations for each timestep graphically
 61: .  -ts_monitor_sp_eig - Monitor eigenvalues of linearized operator graphically
 62: .  -ts_monitor_draw_solution - Monitor solution graphically
 63: .  -ts_monitor_draw_solution_phase  <xleft,yleft,xright,yright> - Monitor solution graphically with phase diagram, requires problem with exactly 2 degrees of freedom
 64: .  -ts_monitor_draw_error - Monitor error graphically, requires use to have provided TSSetSolutionFunction()
 65: .  -ts_monitor_solution [ascii binary draw][:filename][:viewerformat] - monitors the solution at each timestep
 66: .  -ts_monitor_solution_vtk <filename.vts,filename.vtu> - Save each time step to a binary file, use filename-%%03D.vts (filename-%%03D.vtu)
 67: -  -ts_monitor_envelope - determine maximum and minimum value of each component of the solution over the solution time

 69:    Notes:
 70:      See SNESSetFromOptions() and KSPSetFromOptions() for how to control the nonlinear and linear solves used by the time-stepper.

 72:      Certain SNES options get reset for each new nonlinear solver, for example -snes_lag_jacobian <its> and -snes_lag_preconditioner <its>, in order
 73:      to retain them over the multiple nonlinear solves that TS uses you mush also provide -snes_lag_jacobian_persists true and
 74:      -snes_lag_preconditioner_persists true

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

 79:    Level: beginner

 81: .seealso: TSGetType()
 82: @*/
 83: PetscErrorCode  TSSetFromOptions(TS ts)
 84: {
 85:   PetscBool              opt,flg,tflg;
 86:   PetscErrorCode         ierr;
 87:   char                   monfilename[PETSC_MAX_PATH_LEN];
 88:   PetscReal              time_step;
 89:   TSExactFinalTimeOption eftopt;
 90:   char                   dir[16];
 91:   TSIFunction            ifun;
 92:   const char             *defaultType;
 93:   char                   typeName[256];


 97:   TSRegisterAll();
 98:   TSGetIFunction(ts,NULL,&ifun,NULL);

100:   PetscObjectOptionsBegin((PetscObject)ts);
101:   if (((PetscObject)ts)->type_name) defaultType = ((PetscObject)ts)->type_name;
102:   else defaultType = ifun ? TSBEULER : TSEULER;
103:   PetscOptionsFList("-ts_type","TS method","TSSetType",TSList,defaultType,typeName,256,&opt);
104:   if (opt) {
105:     TSSetType(ts,typeName);
106:   } else {
107:     TSSetType(ts,defaultType);
108:   }

110:   /* Handle generic TS options */
111:   PetscOptionsDeprecated("-ts_final_time","-ts_max_time","3.10",NULL);
112:   PetscOptionsReal("-ts_max_time","Maximum time to run to","TSSetMaxTime",ts->max_time,&ts->max_time,NULL);
113:   PetscOptionsInt("-ts_max_steps","Maximum number of time steps","TSSetMaxSteps",ts->max_steps,&ts->max_steps,NULL);
114:   PetscOptionsReal("-ts_init_time","Initial time","TSSetTime",ts->ptime,&ts->ptime,NULL);
115:   PetscOptionsReal("-ts_dt","Initial time step","TSSetTimeStep",ts->time_step,&time_step,&flg);
116:   if (flg) TSSetTimeStep(ts,time_step);
117:   PetscOptionsEnum("-ts_exact_final_time","Option for handling of final time step","TSSetExactFinalTime",TSExactFinalTimeOptions,(PetscEnum)ts->exact_final_time,(PetscEnum*)&eftopt,&flg);
118:   if (flg) TSSetExactFinalTime(ts,eftopt);
119:   PetscOptionsInt("-ts_max_snes_failures","Maximum number of nonlinear solve failures","TSSetMaxSNESFailures",ts->max_snes_failures,&ts->max_snes_failures,NULL);
120:   PetscOptionsInt("-ts_max_reject","Maximum number of step rejections before step fails","TSSetMaxStepRejections",ts->max_reject,&ts->max_reject,NULL);
121:   PetscOptionsBool("-ts_error_if_step_fails","Error if no step succeeds","TSSetErrorIfStepFails",ts->errorifstepfailed,&ts->errorifstepfailed,NULL);
122:   PetscOptionsReal("-ts_rtol","Relative tolerance for local truncation error","TSSetTolerances",ts->rtol,&ts->rtol,NULL);
123:   PetscOptionsReal("-ts_atol","Absolute tolerance for local truncation error","TSSetTolerances",ts->atol,&ts->atol,NULL);

125:   PetscOptionsBool("-ts_rhs_jacobian_test_mult","Test the RHS Jacobian for consistency with RHS at each solve ","None",ts->testjacobian,&ts->testjacobian,NULL);
126:   PetscOptionsBool("-ts_rhs_jacobian_test_mult_transpose","Test the RHS Jacobian transpose for consistency with RHS at each solve ","None",ts->testjacobiantranspose,&ts->testjacobiantranspose,NULL);
127:   PetscOptionsBool("-ts_use_splitrhsfunction","Use the split RHS function for multirate solvers ","TSSetUseSplitRHSFunction",ts->use_splitrhsfunction,&ts->use_splitrhsfunction,NULL);
128: #if defined(PETSC_HAVE_SAWS)
129:   {
130:     PetscBool set;
131:     flg  = PETSC_FALSE;
132:     PetscOptionsBool("-ts_saws_block","Block for SAWs memory snooper at end of TSSolve","PetscObjectSAWsBlock",((PetscObject)ts)->amspublishblock,&flg,&set);
133:     if (set) {
134:       PetscObjectSAWsSetBlock((PetscObject)ts,flg);
135:     }
136:   }
137: #endif

139:   /* Monitor options */
140:   PetscOptionsInt("-ts_monitor_frequency", "Number of time steps between monitor output", "TSMonitorSetFrequency", ts->monitorFrequency, &ts->monitorFrequency, NULL);
141:   TSMonitorSetFromOptions(ts,"-ts_monitor","Monitor time and timestep size","TSMonitorDefault",TSMonitorDefault,NULL);
142:   TSMonitorSetFromOptions(ts,"-ts_monitor_extreme","Monitor extreme values of the solution","TSMonitorExtreme",TSMonitorExtreme,NULL);
143:   TSMonitorSetFromOptions(ts,"-ts_monitor_solution","View the solution at each timestep","TSMonitorSolution",TSMonitorSolution,NULL);
144:   TSMonitorSetFromOptions(ts,"-ts_dmswarm_monitor_moments","Monitor moments of particle distribution","TSDMSwarmMonitorMoments",TSDMSwarmMonitorMoments,NULL);

146:   PetscOptionsString("-ts_monitor_python","Use Python function","TSMonitorSet",NULL,monfilename,sizeof(monfilename),&flg);
147:   if (flg) PetscPythonMonitorSet((PetscObject)ts,monfilename);

149:   PetscOptionsName("-ts_monitor_lg_solution","Monitor solution graphically","TSMonitorLGSolution",&opt);
150:   if (opt) {
151:     PetscInt       howoften = 1;
152:     DM             dm;
153:     PetscBool      net;

155:     PetscOptionsInt("-ts_monitor_lg_solution","Monitor solution graphically","TSMonitorLGSolution",howoften,&howoften,NULL);
156:     TSGetDM(ts,&dm);
157:     PetscObjectTypeCompare((PetscObject)dm,DMNETWORK,&net);
158:     if (net) {
159:       TSMonitorLGCtxNetwork ctx;
160:       TSMonitorLGCtxNetworkCreate(ts,NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,600,400,howoften,&ctx);
161:       TSMonitorSet(ts,TSMonitorLGCtxNetworkSolution,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxNetworkDestroy);
162:       PetscOptionsBool("-ts_monitor_lg_solution_semilogy","Plot the solution with a semi-log axis","",ctx->semilogy,&ctx->semilogy,NULL);
163:     } else {
164:       TSMonitorLGCtx ctx;
165:       TSMonitorLGCtxCreate(PETSC_COMM_SELF,NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,howoften,&ctx);
166:       TSMonitorSet(ts,TSMonitorLGSolution,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);
167:     }
168:   }

170:   PetscOptionsName("-ts_monitor_lg_error","Monitor error graphically","TSMonitorLGError",&opt);
171:   if (opt) {
172:     TSMonitorLGCtx ctx;
173:     PetscInt       howoften = 1;

175:     PetscOptionsInt("-ts_monitor_lg_error","Monitor error graphically","TSMonitorLGError",howoften,&howoften,NULL);
176:     TSMonitorLGCtxCreate(PETSC_COMM_SELF,NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,howoften,&ctx);
177:     TSMonitorSet(ts,TSMonitorLGError,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);
178:   }
179:   TSMonitorSetFromOptions(ts,"-ts_monitor_error","View the error at each timestep","TSMonitorError",TSMonitorError,NULL);

181:   PetscOptionsName("-ts_monitor_lg_timestep","Monitor timestep size graphically","TSMonitorLGTimeStep",&opt);
182:   if (opt) {
183:     TSMonitorLGCtx ctx;
184:     PetscInt       howoften = 1;

186:     PetscOptionsInt("-ts_monitor_lg_timestep","Monitor timestep size graphically","TSMonitorLGTimeStep",howoften,&howoften,NULL);
187:     TSMonitorLGCtxCreate(PetscObjectComm((PetscObject)ts),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,howoften,&ctx);
188:     TSMonitorSet(ts,TSMonitorLGTimeStep,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);
189:   }
190:   PetscOptionsName("-ts_monitor_lg_timestep_log","Monitor log timestep size graphically","TSMonitorLGTimeStep",&opt);
191:   if (opt) {
192:     TSMonitorLGCtx ctx;
193:     PetscInt       howoften = 1;

195:     PetscOptionsInt("-ts_monitor_lg_timestep_log","Monitor log timestep size graphically","TSMonitorLGTimeStep",howoften,&howoften,NULL);
196:     TSMonitorLGCtxCreate(PetscObjectComm((PetscObject)ts),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,howoften,&ctx);
197:     TSMonitorSet(ts,TSMonitorLGTimeStep,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);
198:     ctx->semilogy = PETSC_TRUE;
199:   }

201:   PetscOptionsName("-ts_monitor_lg_snes_iterations","Monitor number nonlinear iterations for each timestep graphically","TSMonitorLGSNESIterations",&opt);
202:   if (opt) {
203:     TSMonitorLGCtx ctx;
204:     PetscInt       howoften = 1;

206:     PetscOptionsInt("-ts_monitor_lg_snes_iterations","Monitor number nonlinear iterations for each timestep graphically","TSMonitorLGSNESIterations",howoften,&howoften,NULL);
207:     TSMonitorLGCtxCreate(PetscObjectComm((PetscObject)ts),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,howoften,&ctx);
208:     TSMonitorSet(ts,TSMonitorLGSNESIterations,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);
209:   }
210:   PetscOptionsName("-ts_monitor_lg_ksp_iterations","Monitor number nonlinear iterations for each timestep graphically","TSMonitorLGKSPIterations",&opt);
211:   if (opt) {
212:     TSMonitorLGCtx ctx;
213:     PetscInt       howoften = 1;

215:     PetscOptionsInt("-ts_monitor_lg_ksp_iterations","Monitor number nonlinear iterations for each timestep graphically","TSMonitorLGKSPIterations",howoften,&howoften,NULL);
216:     TSMonitorLGCtxCreate(PetscObjectComm((PetscObject)ts),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,howoften,&ctx);
217:     TSMonitorSet(ts,TSMonitorLGKSPIterations,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);
218:   }
219:   PetscOptionsName("-ts_monitor_sp_eig","Monitor eigenvalues of linearized operator graphically","TSMonitorSPEig",&opt);
220:   if (opt) {
221:     TSMonitorSPEigCtx ctx;
222:     PetscInt          howoften = 1;

224:     PetscOptionsInt("-ts_monitor_sp_eig","Monitor eigenvalues of linearized operator graphically","TSMonitorSPEig",howoften,&howoften,NULL);
225:     TSMonitorSPEigCtxCreate(PETSC_COMM_SELF,NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx);
226:     TSMonitorSet(ts,TSMonitorSPEig,ctx,(PetscErrorCode (*)(void**))TSMonitorSPEigCtxDestroy);
227:   }
228:   PetscOptionsName("-ts_monitor_sp_swarm","Display particle phase from the DMSwarm","TSMonitorSPSwarm",&opt);
229:   if (opt) {
230:     TSMonitorSPCtx  ctx;
231:     PetscInt        howoften = 1, retain = 0;
232:     PetscBool       phase = PETSC_TRUE;

234:     PetscOptionsInt("-ts_monitor_sp_swarm","Display particles phase from the DMSwarm", "TSMonitorSPSwarm", howoften, &howoften, NULL);
235:     PetscOptionsInt("-ts_monitor_sp_swarm_retain", "Retain n points plotted to show trajectory, -1 for all points", "TSMonitorSPSwarm", retain, &retain, NULL);
236:     PetscOptionsBool("-ts_monitor_sp_swarm_phase", "Plot in phase space rather than coordinate space", "TSMonitorSPSwarm", phase, &phase, NULL);
237:     TSMonitorSPCtxCreate(PetscObjectComm((PetscObject) ts), NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 300, 300, howoften, retain, phase, &ctx);
238:     TSMonitorSet(ts, TSMonitorSPSwarmSolution, ctx, (PetscErrorCode (*)(void**))TSMonitorSPCtxDestroy);
239:   }
240:   opt  = PETSC_FALSE;
241:   PetscOptionsName("-ts_monitor_draw_solution","Monitor solution graphically","TSMonitorDrawSolution",&opt);
242:   if (opt) {
243:     TSMonitorDrawCtx ctx;
244:     PetscInt         howoften = 1;

246:     PetscOptionsInt("-ts_monitor_draw_solution","Monitor solution graphically","TSMonitorDrawSolution",howoften,&howoften,NULL);
247:     TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),NULL,"Computed Solution",PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx);
248:     TSMonitorSet(ts,TSMonitorDrawSolution,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);
249:   }
250:   opt  = PETSC_FALSE;
251:   PetscOptionsName("-ts_monitor_draw_solution_phase","Monitor solution graphically","TSMonitorDrawSolutionPhase",&opt);
252:   if (opt) {
253:     TSMonitorDrawCtx ctx;
254:     PetscReal        bounds[4];
255:     PetscInt         n = 4;
256:     PetscDraw        draw;
257:     PetscDrawAxis    axis;

259:     PetscOptionsRealArray("-ts_monitor_draw_solution_phase","Monitor solution graphically","TSMonitorDrawSolutionPhase",bounds,&n,NULL);
261:     TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,300,300,1,&ctx);
262:     PetscViewerDrawGetDraw(ctx->viewer,0,&draw);
263:     PetscViewerDrawGetDrawAxis(ctx->viewer,0,&axis);
264:     PetscDrawAxisSetLimits(axis,bounds[0],bounds[2],bounds[1],bounds[3]);
265:     PetscDrawAxisSetLabels(axis,"Phase Diagram","Variable 1","Variable 2");
266:     TSMonitorSet(ts,TSMonitorDrawSolutionPhase,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);
267:   }
268:   opt  = PETSC_FALSE;
269:   PetscOptionsName("-ts_monitor_draw_error","Monitor error graphically","TSMonitorDrawError",&opt);
270:   if (opt) {
271:     TSMonitorDrawCtx ctx;
272:     PetscInt         howoften = 1;

274:     PetscOptionsInt("-ts_monitor_draw_error","Monitor error graphically","TSMonitorDrawError",howoften,&howoften,NULL);
275:     TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),NULL,"Error",PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx);
276:     TSMonitorSet(ts,TSMonitorDrawError,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);
277:   }
278:   opt  = PETSC_FALSE;
279:   PetscOptionsName("-ts_monitor_draw_solution_function","Monitor solution provided by TSMonitorSetSolutionFunction() graphically","TSMonitorDrawSolutionFunction",&opt);
280:   if (opt) {
281:     TSMonitorDrawCtx ctx;
282:     PetscInt         howoften = 1;

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

289:   opt  = PETSC_FALSE;
290:   PetscOptionsString("-ts_monitor_solution_vtk","Save each time step to a binary file, use filename-%%03D.vts","TSMonitorSolutionVTK",NULL,monfilename,sizeof(monfilename),&flg);
291:   if (flg) {
292:     const char *ptr,*ptr2;
293:     char       *filetemplate;
295:     /* Do some cursory validation of the input. */
296:     PetscStrstr(monfilename,"%",(char**)&ptr);
298:     for (ptr++; ptr && *ptr; ptr++) {
299:       PetscStrchr("DdiouxX",*ptr,(char**)&ptr2);
301:       if (ptr2) break;
302:     }
303:     PetscStrallocpy(monfilename,&filetemplate);
304:     TSMonitorSet(ts,TSMonitorSolutionVTK,filetemplate,(PetscErrorCode (*)(void**))TSMonitorSolutionVTKDestroy);
305:   }

307:   PetscOptionsString("-ts_monitor_dmda_ray","Display a ray of the solution","None","y=0",dir,sizeof(dir),&flg);
308:   if (flg) {
309:     TSMonitorDMDARayCtx *rayctx;
310:     int                  ray = 0;
311:     DMDirection          ddir;
312:     DM                   da;
313:     PetscMPIInt          rank;

316:     if (dir[0] == 'x') ddir = DM_X;
317:     else if (dir[0] == 'y') ddir = DM_Y;
318:     else SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONG,"Unknown ray %s",dir);
319:     sscanf(dir+2,"%d",&ray);

321:     PetscInfo(((PetscObject)ts),"Displaying DMDA ray %c = %d\n",dir[0],ray);
322:     PetscNew(&rayctx);
323:     TSGetDM(ts,&da);
324:     DMDAGetRay(da,ddir,ray,&rayctx->ray,&rayctx->scatter);
325:     MPI_Comm_rank(PetscObjectComm((PetscObject)ts),&rank);
326:     if (rank == 0) {
327:       PetscViewerDrawOpen(PETSC_COMM_SELF,NULL,NULL,0,0,600,300,&rayctx->viewer);
328:     }
329:     rayctx->lgctx = NULL;
330:     TSMonitorSet(ts,TSMonitorDMDARay,rayctx,TSMonitorDMDARayDestroy);
331:   }
332:   PetscOptionsString("-ts_monitor_lg_dmda_ray","Display a ray of the solution","None","x=0",dir,sizeof(dir),&flg);
333:   if (flg) {
334:     TSMonitorDMDARayCtx *rayctx;
335:     int                 ray = 0;
336:     DMDirection         ddir;
337:     DM                  da;
338:     PetscInt            howoften = 1;

341:     if      (dir[0] == 'x') ddir = DM_X;
342:     else if (dir[0] == 'y') ddir = DM_Y;
343:     else SETERRQ(PetscObjectComm((PetscObject) ts), PETSC_ERR_ARG_WRONG, "Unknown ray direction %s", dir);
344:     sscanf(dir+2, "%d", &ray);

346:     PetscInfo(((PetscObject) ts),"Displaying LG DMDA ray %c = %d\n", dir[0], ray);
347:     PetscNew(&rayctx);
348:     TSGetDM(ts, &da);
349:     DMDAGetRay(da, ddir, ray, &rayctx->ray, &rayctx->scatter);
350:     TSMonitorLGCtxCreate(PETSC_COMM_SELF,NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,600,400,howoften,&rayctx->lgctx);
351:     TSMonitorSet(ts, TSMonitorLGDMDARay, rayctx, TSMonitorDMDARayDestroy);
352:   }

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

358:     TSMonitorEnvelopeCtxCreate(ts,&ctx);
359:     TSMonitorSet(ts,TSMonitorEnvelope,ctx,(PetscErrorCode (*)(void**))TSMonitorEnvelopeCtxDestroy);
360:   }
361:   flg  = PETSC_FALSE;
362:   PetscOptionsBool("-ts_monitor_cancel","Remove all monitors","TSMonitorCancel",flg,&flg,&opt);
363:   if (opt && flg) TSMonitorCancel(ts);

365:   flg  = PETSC_FALSE;
366:   PetscOptionsBool("-ts_fd_color", "Use finite differences with coloring to compute IJacobian", "TSComputeJacobianDefaultColor", flg, &flg, NULL);
367:   if (flg) {
368:     DM   dm;
369:     DMTS tdm;

371:     TSGetDM(ts, &dm);
372:     DMGetDMTS(dm, &tdm);
373:     tdm->ijacobianctx = NULL;
374:     TSSetIJacobian(ts, NULL, NULL, TSComputeIJacobianDefaultColor, NULL);
375:     PetscInfo(ts, "Setting default finite difference coloring Jacobian matrix\n");
376:   }

378:   /* Handle specific TS options */
379:   if (ts->ops->setfromoptions) {
380:     (*ts->ops->setfromoptions)(PetscOptionsObject,ts);
381:   }

383:   /* Handle TSAdapt options */
384:   TSGetAdapt(ts,&ts->adapt);
385:   TSAdaptSetDefaultType(ts->adapt,ts->default_adapt_type);
386:   TSAdaptSetFromOptions(PetscOptionsObject,ts->adapt);

388:   /* TS trajectory must be set after TS, since it may use some TS options above */
389:   tflg = ts->trajectory ? PETSC_TRUE : PETSC_FALSE;
390:   PetscOptionsBool("-ts_save_trajectory","Save the solution at each timestep","TSSetSaveTrajectory",tflg,&tflg,NULL);
391:   if (tflg) {
392:     TSSetSaveTrajectory(ts);
393:   }

395:   TSAdjointSetFromOptions(PetscOptionsObject,ts);

397:   /* process any options handlers added with PetscObjectAddOptionsHandler() */
398:   PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)ts);
399:   PetscOptionsEnd();

401:   if (ts->trajectory) {
402:     TSTrajectorySetFromOptions(ts->trajectory,ts);
403:   }

405:   /* why do we have to do this here and not during TSSetUp? */
406:   TSGetSNES(ts,&ts->snes);
407:   if (ts->problem_type == TS_LINEAR) {
408:     PetscObjectTypeCompareAny((PetscObject)ts->snes,&flg,SNESKSPONLY,SNESKSPTRANSPOSEONLY,"");
409:     if (!flg) SNESSetType(ts->snes,SNESKSPONLY);
410:   }
411:   SNESSetFromOptions(ts->snes);
412:   return 0;
413: }

415: /*@
416:    TSGetTrajectory - Gets the trajectory from a TS if it exists

418:    Collective on TS

420:    Input Parameters:
421: .  ts - the TS context obtained from TSCreate()

423:    Output Parameters:
424: .  tr - the TSTrajectory object, if it exists

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

428:    Level: advanced

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

432: @*/
433: PetscErrorCode  TSGetTrajectory(TS ts,TSTrajectory *tr)
434: {
436:   *tr = ts->trajectory;
437:   return 0;
438: }

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

443:    Collective on TS

445:    Input Parameter:
446: .  ts - the TS context obtained from TSCreate()

448:    Options Database:
449: +  -ts_save_trajectory - saves the trajectory to a file
450: -  -ts_trajectory_type type - set trajectory type

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

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

457:    Level: intermediate

459: .seealso: TSGetTrajectory(), TSAdjointSolve()

461: @*/
462: PetscErrorCode  TSSetSaveTrajectory(TS ts)
463: {
465:   if (!ts->trajectory) {
466:     TSTrajectoryCreate(PetscObjectComm((PetscObject)ts),&ts->trajectory);
467:   }
468:   return 0;
469: }

471: /*@
472:    TSResetTrajectory - Destroys and recreates the internal TSTrajectory object

474:    Collective on TS

476:    Input Parameters:
477: .  ts - the TS context obtained from TSCreate()

479:    Level: intermediate

481: .seealso: TSGetTrajectory(), TSAdjointSolve(), TSRemoveTrajectory()

483: @*/
484: PetscErrorCode  TSResetTrajectory(TS ts)
485: {
487:   if (ts->trajectory) {
488:     TSTrajectoryDestroy(&ts->trajectory);
489:     TSTrajectoryCreate(PetscObjectComm((PetscObject)ts),&ts->trajectory);
490:   }
491:   return 0;
492: }

494: /*@
495:    TSRemoveTrajectory - Destroys and removes the internal TSTrajectory object from TS

497:    Collective on TS

499:    Input Parameters:
500: .  ts - the TS context obtained from TSCreate()

502:    Level: intermediate

504: .seealso: TSResetTrajectory(), TSAdjointSolve()

506: @*/
507: PetscErrorCode TSRemoveTrajectory(TS ts)
508: {
510:   if (ts->trajectory) {
511:     TSTrajectoryDestroy(&ts->trajectory);
512:   }
513:   return 0;
514: }

516: /*@
517:    TSComputeRHSJacobian - Computes the Jacobian matrix that has been
518:       set with TSSetRHSJacobian().

520:    Collective on TS

522:    Input Parameters:
523: +  ts - the TS context
524: .  t - current timestep
525: -  U - input vector

527:    Output Parameters:
528: +  A - Jacobian matrix
529: -  B - optional preconditioning matrix

531:    Notes:
532:    Most users should not need to explicitly call this routine, as it
533:    is used internally within the nonlinear solvers.

535:    Level: developer

537: .seealso:  TSSetRHSJacobian(), KSPSetOperators()
538: @*/
539: PetscErrorCode  TSComputeRHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B)
540: {
541:   PetscObjectState Ustate;
542:   PetscObjectId    Uid;
543:   DM               dm;
544:   DMTS             tsdm;
545:   TSRHSJacobian    rhsjacobianfunc;
546:   void             *ctx;
547:   TSRHSFunction    rhsfunction;

552:   TSGetDM(ts,&dm);
553:   DMGetDMTS(dm,&tsdm);
554:   DMTSGetRHSFunction(dm,&rhsfunction,NULL);
555:   DMTSGetRHSJacobian(dm,&rhsjacobianfunc,&ctx);
556:   PetscObjectStateGet((PetscObject)U,&Ustate);
557:   PetscObjectGetId((PetscObject)U,&Uid);

559:   if (ts->rhsjacobian.time == t && (ts->problem_type == TS_LINEAR || (ts->rhsjacobian.Xid == Uid && ts->rhsjacobian.Xstate == Ustate)) && (rhsfunction != TSComputeRHSFunctionLinear)) return 0;

562:   if (rhsjacobianfunc) {
563:     PetscLogEventBegin(TS_JacobianEval,ts,U,A,B);
564:     PetscStackPush("TS user Jacobian function");
565:     (*rhsjacobianfunc)(ts,t,U,A,B,ctx);
566:     PetscStackPop;
567:     ts->rhsjacs++;
568:     PetscLogEventEnd(TS_JacobianEval,ts,U,A,B);
569:   } else {
570:     MatZeroEntries(A);
571:     if (B && A != B) MatZeroEntries(B);
572:   }
573:   ts->rhsjacobian.time  = t;
574:   ts->rhsjacobian.shift = 0;
575:   ts->rhsjacobian.scale = 1.;
576:   PetscObjectGetId((PetscObject)U,&ts->rhsjacobian.Xid);
577:   PetscObjectStateGet((PetscObject)U,&ts->rhsjacobian.Xstate);
578:   return 0;
579: }

581: /*@
582:    TSComputeRHSFunction - Evaluates the right-hand-side function.

584:    Collective on TS

586:    Input Parameters:
587: +  ts - the TS context
588: .  t - current time
589: -  U - state vector

591:    Output Parameter:
592: .  y - right hand side

594:    Note:
595:    Most users should not need to explicitly call this routine, as it
596:    is used internally within the nonlinear solvers.

598:    Level: developer

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

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


618:   if (rhsfunction) {
619:     PetscLogEventBegin(TS_FunctionEval,ts,U,y,0);
620:     VecLockReadPush(U);
621:     PetscStackPush("TS user right-hand-side function");
622:     (*rhsfunction)(ts,t,U,y,ctx);
623:     PetscStackPop;
624:     VecLockReadPop(U);
625:     ts->rhsfuncs++;
626:     PetscLogEventEnd(TS_FunctionEval,ts,U,y,0);
627:   } else {
628:     VecZeroEntries(y);
629:   }
630:   return 0;
631: }

633: /*@
634:    TSComputeSolutionFunction - Evaluates the solution function.

636:    Collective on TS

638:    Input Parameters:
639: +  ts - the TS context
640: -  t - current time

642:    Output Parameter:
643: .  U - the solution

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

649:    Level: developer

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

661:   TSGetDM(ts,&dm);
662:   DMTSGetSolutionFunction(dm,&solutionfunction,&ctx);

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

674:    Collective on TS

676:    Input Parameters:
677: +  ts - the TS context
678: -  t - current time

680:    Output Parameter:
681: .  U - the function value

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

687:    Level: developer

689: .seealso: TSSetSolutionFunction(), TSSetRHSFunction(), TSComputeIFunction()
690: @*/
691: PetscErrorCode TSComputeForcingFunction(TS ts,PetscReal t,Vec U)
692: {
693:   void              *ctx;
694:   DM                 dm;
695:   TSForcingFunction  forcing;

699:   TSGetDM(ts,&dm);
700:   DMTSGetForcingFunction(dm,&forcing,&ctx);

702:   if (forcing) {
703:     PetscStackPush("TS user forcing function");
704:     (*forcing)(ts,t,U,ctx);
705:     PetscStackPop;
706:   }
707:   return 0;
708: }

710: static PetscErrorCode TSGetRHSVec_Private(TS ts,Vec *Frhs)
711: {
712:   Vec            F;

714:   *Frhs = NULL;
715:   TSGetIFunction(ts,&F,NULL,NULL);
716:   if (!ts->Frhs) {
717:     VecDuplicate(F,&ts->Frhs);
718:   }
719:   *Frhs = ts->Frhs;
720:   return 0;
721: }

723: PetscErrorCode TSGetRHSMats_Private(TS ts,Mat *Arhs,Mat *Brhs)
724: {
725:   Mat            A,B;
726:   TSIJacobian    ijacobian;

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

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

774:    Collective on TS

776:    Input Parameters:
777: +  ts - the TS context
778: .  t - current time
779: .  U - state vector
780: .  Udot - time derivative of state vector
781: -  imex - flag indicates if the method is IMEX so that the RHSFunction should be kept separate

783:    Output Parameter:
784: .  Y - right hand side

786:    Note:
787:    Most users should not need to explicitly call this routine, as it
788:    is used internally within the nonlinear solvers.

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

793:    Level: developer

795: .seealso: TSSetIFunction(), TSComputeRHSFunction()
796: @*/
797: PetscErrorCode TSComputeIFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec Y,PetscBool imex)
798: {
799:   TSIFunction    ifunction;
800:   TSRHSFunction  rhsfunction;
801:   void           *ctx;
802:   DM             dm;


809:   TSGetDM(ts,&dm);
810:   DMTSGetIFunction(dm,&ifunction,&ctx);
811:   DMTSGetRHSFunction(dm,&rhsfunction,NULL);


815:   PetscLogEventBegin(TS_FunctionEval,ts,U,Udot,Y);
816:   if (ifunction) {
817:     PetscStackPush("TS user implicit function");
818:     (*ifunction)(ts,t,U,Udot,Y,ctx);
819:     PetscStackPop;
820:     ts->ifuncs++;
821:   }
822:   if (imex) {
823:     if (!ifunction) {
824:       VecCopy(Udot,Y);
825:     }
826:   } else if (rhsfunction) {
827:     if (ifunction) {
828:       Vec Frhs;
829:       TSGetRHSVec_Private(ts,&Frhs);
830:       TSComputeRHSFunction(ts,t,U,Frhs);
831:       VecAXPY(Y,-1,Frhs);
832:     } else {
833:       TSComputeRHSFunction(ts,t,U,Y);
834:       VecAYPX(Y,-1,Udot);
835:     }
836:   }
837:   PetscLogEventEnd(TS_FunctionEval,ts,U,Udot,Y);
838:   return 0;
839: }

841: /*
842:    TSRecoverRHSJacobian - Recover the Jacobian matrix so that one can call TSComputeRHSJacobian() on it.

844:    Note:
845:    This routine is needed when one switches from TSComputeIJacobian() to TSComputeRHSJacobian() because the Jacobian matrix may be shifted or scaled in TSComputeIJacobian().

847: */
848: static PetscErrorCode TSRecoverRHSJacobian(TS ts,Mat A,Mat B)
849: {

854:   if (ts->rhsjacobian.shift) {
855:     MatShift(A,-ts->rhsjacobian.shift);
856:   }
857:   if (ts->rhsjacobian.scale == -1.) {
858:     MatScale(A,-1);
859:   }
860:   if (B && B == ts->Brhs && A != B) {
861:     if (ts->rhsjacobian.shift) {
862:       MatShift(B,-ts->rhsjacobian.shift);
863:     }
864:     if (ts->rhsjacobian.scale == -1.) {
865:       MatScale(B,-1);
866:     }
867:   }
868:   ts->rhsjacobian.shift = 0;
869:   ts->rhsjacobian.scale = 1.;
870:   return 0;
871: }

873: /*@
874:    TSComputeIJacobian - Evaluates the Jacobian of the DAE

876:    Collective on TS

878:    Input
879:       Input Parameters:
880: +  ts - the TS context
881: .  t - current timestep
882: .  U - state vector
883: .  Udot - time derivative of state vector
884: .  shift - shift to apply, see note below
885: -  imex - flag indicates if the method is IMEX so that the RHSJacobian should be kept separate

887:    Output Parameters:
888: +  A - Jacobian matrix
889: -  B - matrix from which the preconditioner is constructed; often the same as A

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

894:    dF/dU + shift*dF/dUdot

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

899:    Level: developer

901: .seealso:  TSSetIJacobian()
902: @*/
903: PetscErrorCode TSComputeIJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal shift,Mat A,Mat B,PetscBool imex)
904: {
905:   TSIJacobian    ijacobian;
906:   TSRHSJacobian  rhsjacobian;
907:   DM             dm;
908:   void           *ctx;


918:   TSGetDM(ts,&dm);
919:   DMTSGetIJacobian(dm,&ijacobian,&ctx);
920:   DMTSGetRHSJacobian(dm,&rhsjacobian,NULL);


924:   PetscLogEventBegin(TS_JacobianEval,ts,U,A,B);
925:   if (ijacobian) {
926:     PetscStackPush("TS user implicit Jacobian");
927:     (*ijacobian)(ts,t,U,Udot,shift,A,B,ctx);
928:     ts->ijacs++;
929:     PetscStackPop;
930:   }
931:   if (imex) {
932:     if (!ijacobian) {  /* system was written as Udot = G(t,U) */
933:       PetscBool assembled;
934:       if (rhsjacobian) {
935:         Mat Arhs = NULL;
936:         TSGetRHSMats_Private(ts,&Arhs,NULL);
937:         if (A == Arhs) {
939:           ts->rhsjacobian.time = PETSC_MIN_REAL;
940:         }
941:       }
942:       MatZeroEntries(A);
943:       MatAssembled(A,&assembled);
944:       if (!assembled) {
945:         MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
946:         MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
947:       }
948:       MatShift(A,shift);
949:       if (A != B) {
950:         MatZeroEntries(B);
951:         MatAssembled(B,&assembled);
952:         if (!assembled) {
953:           MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
954:           MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
955:         }
956:         MatShift(B,shift);
957:       }
958:     }
959:   } else {
960:     Mat Arhs = NULL,Brhs = NULL;
961:     if (rhsjacobian) { /* RHSJacobian needs to be converted to part of IJacobian if exists */
962:       TSGetRHSMats_Private(ts,&Arhs,&Brhs);
963:     }
964:     if (Arhs == A) { /* No IJacobian matrix, so we only have the RHS matrix */
965:       PetscObjectState Ustate;
966:       PetscObjectId    Uid;
967:       TSRHSFunction    rhsfunction;

969:       DMTSGetRHSFunction(dm,&rhsfunction,NULL);
970:       PetscObjectStateGet((PetscObject)U,&Ustate);
971:       PetscObjectGetId((PetscObject)U,&Uid);
972:       if ((rhsjacobian == TSComputeRHSJacobianConstant || (ts->rhsjacobian.time == t && (ts->problem_type == TS_LINEAR || (ts->rhsjacobian.Xid == Uid && ts->rhsjacobian.Xstate == Ustate)) && rhsfunction != TSComputeRHSFunctionLinear)) && ts->rhsjacobian.scale == -1.) { /* No need to recompute RHSJacobian */
973:         MatShift(A,shift-ts->rhsjacobian.shift); /* revert the old shift and add the new shift with a single call to MatShift */
974:         if (A != B) {
975:           MatShift(B,shift-ts->rhsjacobian.shift);
976:         }
977:       } else {
978:         PetscBool flg;

980:         if (ts->rhsjacobian.reuse) { /* Undo the damage */
981:           /* MatScale has a short path for this case.
982:              However, this code path is taken the first time TSComputeRHSJacobian is called
983:              and the matrices have not been assembled yet */
984:           TSRecoverRHSJacobian(ts,A,B);
985:         }
986:         TSComputeRHSJacobian(ts,t,U,A,B);
987:         SNESGetUseMatrixFree(ts->snes,NULL,&flg);
988:         /* since -snes_mf_operator uses the full SNES function it does not need to be shifted or scaled here */
989:         if (!flg) {
990:           MatScale(A,-1);
991:           MatShift(A,shift);
992:         }
993:         if (A != B) {
994:           MatScale(B,-1);
995:           MatShift(B,shift);
996:         }
997:       }
998:       ts->rhsjacobian.scale = -1;
999:       ts->rhsjacobian.shift = shift;
1000:     } else if (Arhs) {          /* Both IJacobian and RHSJacobian */
1001:       if (!ijacobian) {         /* No IJacobian provided, but we have a separate RHS matrix */
1002:         MatZeroEntries(A);
1003:         MatShift(A,shift);
1004:         if (A != B) {
1005:           MatZeroEntries(B);
1006:           MatShift(B,shift);
1007:         }
1008:       }
1009:       TSComputeRHSJacobian(ts,t,U,Arhs,Brhs);
1010:       MatAXPY(A,-1,Arhs,ts->axpy_pattern);
1011:       if (A != B) {
1012:         MatAXPY(B,-1,Brhs,ts->axpy_pattern);
1013:       }
1014:     }
1015:   }
1016:   PetscLogEventEnd(TS_JacobianEval,ts,U,A,B);
1017:   return 0;
1018: }

1020: /*@C
1021:     TSSetRHSFunction - Sets the routine for evaluating the function,
1022:     where U_t = G(t,u).

1024:     Logically Collective on TS

1026:     Input Parameters:
1027: +   ts - the TS context obtained from TSCreate()
1028: .   r - vector to put the computed right hand side (or NULL to have it created)
1029: .   f - routine for evaluating the right-hand-side function
1030: -   ctx - [optional] user-defined context for private data for the
1031:           function evaluation routine (may be NULL)

1033:     Calling sequence of f:
1034: $     PetscErrorCode f(TS ts,PetscReal t,Vec u,Vec F,void *ctx);

1036: +   ts - timestep context
1037: .   t - current timestep
1038: .   u - input vector
1039: .   F - function vector
1040: -   ctx - [optional] user-defined function context

1042:     Level: beginner

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

1047: .seealso: TSSetRHSJacobian(), TSSetIJacobian(), TSSetIFunction()
1048: @*/
1049: PetscErrorCode  TSSetRHSFunction(TS ts,Vec r,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx)
1050: {
1051:   SNES           snes;
1052:   Vec            ralloc = NULL;
1053:   DM             dm;


1058:   TSGetDM(ts,&dm);
1059:   DMTSSetRHSFunction(dm,f,ctx);
1060:   TSGetSNES(ts,&snes);
1061:   if (!r && !ts->dm && ts->vec_sol) {
1062:     VecDuplicate(ts->vec_sol,&ralloc);
1063:     r = ralloc;
1064:   }
1065:   SNESSetFunction(snes,r,SNESTSFormFunction,ts);
1066:   VecDestroy(&ralloc);
1067:   return 0;
1068: }

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

1073:     Logically Collective on TS

1075:     Input Parameters:
1076: +   ts - the TS context obtained from TSCreate()
1077: .   f - routine for evaluating the solution
1078: -   ctx - [optional] user-defined context for private data for the
1079:           function evaluation routine (may be NULL)

1081:     Calling sequence of f:
1082: $     PetscErrorCode f(TS ts,PetscReal t,Vec u,void *ctx);

1084: +   t - current timestep
1085: .   u - output vector
1086: -   ctx - [optional] user-defined function context

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

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

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

1099:     Level: beginner

1101: .seealso: TSSetRHSJacobian(), TSSetIJacobian(), TSComputeSolutionFunction(), TSSetForcingFunction(), TSSetSolution(), TSGetSolution(), TSMonitorLGError(), TSMonitorDrawError()
1102: @*/
1103: PetscErrorCode  TSSetSolutionFunction(TS ts,PetscErrorCode (*f)(TS,PetscReal,Vec,void*),void *ctx)
1104: {
1105:   DM             dm;

1108:   TSGetDM(ts,&dm);
1109:   DMTSSetSolutionFunction(dm,f,ctx);
1110:   return 0;
1111: }

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

1116:     Logically Collective on TS

1118:     Input Parameters:
1119: +   ts - the TS context obtained from TSCreate()
1120: .   func - routine for evaluating the forcing function
1121: -   ctx - [optional] user-defined context for private data for the
1122:           function evaluation routine (may be NULL)

1124:     Calling sequence of func:
1125: $     PetscErrorCode func (TS ts,PetscReal t,Vec f,void *ctx);

1127: +   t - current timestep
1128: .   f - output vector
1129: -   ctx - [optional] user-defined function context

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

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

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

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

1143:     Level: beginner

1145: .seealso: TSSetRHSJacobian(), TSSetIJacobian(), TSComputeSolutionFunction(), TSSetSolutionFunction()
1146: @*/
1147: PetscErrorCode  TSSetForcingFunction(TS ts,TSForcingFunction func,void *ctx)
1148: {
1149:   DM             dm;

1152:   TSGetDM(ts,&dm);
1153:   DMTSSetForcingFunction(dm,func,ctx);
1154:   return 0;
1155: }

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

1161:    Logically Collective on TS

1163:    Input Parameters:
1164: +  ts  - the TS context obtained from TSCreate()
1165: .  Amat - (approximate) Jacobian matrix
1166: .  Pmat - matrix from which preconditioner is to be constructed (usually the same as Amat)
1167: .  f   - the Jacobian evaluation routine
1168: -  ctx - [optional] user-defined context for private data for the
1169:          Jacobian evaluation routine (may be NULL)

1171:    Calling sequence of f:
1172: $     PetscErrorCode f(TS ts,PetscReal t,Vec u,Mat A,Mat B,void *ctx);

1174: +  t - current timestep
1175: .  u - input vector
1176: .  Amat - (approximate) Jacobian matrix
1177: .  Pmat - matrix from which preconditioner is to be constructed (usually the same as Amat)
1178: -  ctx - [optional] user-defined context for matrix evaluation routine

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

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

1186:    Level: beginner

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

1190: @*/
1191: PetscErrorCode  TSSetRHSJacobian(TS ts,Mat Amat,Mat Pmat,TSRHSJacobian f,void *ctx)
1192: {
1193:   SNES           snes;
1194:   DM             dm;
1195:   TSIJacobian    ijacobian;


1203:   TSGetDM(ts,&dm);
1204:   DMTSSetRHSJacobian(dm,f,ctx);
1205:   DMTSGetIJacobian(dm,&ijacobian,NULL);
1206:   TSGetSNES(ts,&snes);
1207:   if (!ijacobian) {
1208:     SNESSetJacobian(snes,Amat,Pmat,SNESTSFormJacobian,ts);
1209:   }
1210:   if (Amat) {
1211:     PetscObjectReference((PetscObject)Amat);
1212:     MatDestroy(&ts->Arhs);
1213:     ts->Arhs = Amat;
1214:   }
1215:   if (Pmat) {
1216:     PetscObjectReference((PetscObject)Pmat);
1217:     MatDestroy(&ts->Brhs);
1218:     ts->Brhs = Pmat;
1219:   }
1220:   return 0;
1221: }

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

1226:    Logically Collective on TS

1228:    Input Parameters:
1229: +  ts  - the TS context obtained from TSCreate()
1230: .  r   - vector to hold the residual (or NULL to have it created internally)
1231: .  f   - the function evaluation routine
1232: -  ctx - user-defined context for private data for the function evaluation routine (may be NULL)

1234:    Calling sequence of f:
1235: $     PetscErrorCode f(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx);

1237: +  t   - time at step/stage being solved
1238: .  u   - state vector
1239: .  u_t - time derivative of state vector
1240: .  F   - function vector
1241: -  ctx - [optional] user-defined context for matrix evaluation routine

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

1246:    Level: beginner

1248: .seealso: TSSetRHSJacobian(), TSSetRHSFunction(), TSSetIJacobian()
1249: @*/
1250: PetscErrorCode  TSSetIFunction(TS ts,Vec r,TSIFunction f,void *ctx)
1251: {
1252:   SNES           snes;
1253:   Vec            ralloc = NULL;
1254:   DM             dm;


1259:   TSGetDM(ts,&dm);
1260:   DMTSSetIFunction(dm,f,ctx);

1262:   TSGetSNES(ts,&snes);
1263:   if (!r && !ts->dm && ts->vec_sol) {
1264:     VecDuplicate(ts->vec_sol,&ralloc);
1265:     r  = ralloc;
1266:   }
1267:   SNESSetFunction(snes,r,SNESTSFormFunction,ts);
1268:   VecDestroy(&ralloc);
1269:   return 0;
1270: }

1272: /*@C
1273:    TSGetIFunction - Returns the vector where the implicit residual is stored and the function/context to compute it.

1275:    Not Collective

1277:    Input Parameter:
1278: .  ts - the TS context

1280:    Output Parameters:
1281: +  r - vector to hold residual (or NULL)
1282: .  func - the function to compute residual (or NULL)
1283: -  ctx - the function context (or NULL)

1285:    Level: advanced

1287: .seealso: TSSetIFunction(), SNESGetFunction()
1288: @*/
1289: PetscErrorCode TSGetIFunction(TS ts,Vec *r,TSIFunction *func,void **ctx)
1290: {
1291:   SNES           snes;
1292:   DM             dm;

1295:   TSGetSNES(ts,&snes);
1296:   SNESGetFunction(snes,r,NULL,NULL);
1297:   TSGetDM(ts,&dm);
1298:   DMTSGetIFunction(dm,func,ctx);
1299:   return 0;
1300: }

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

1305:    Not Collective

1307:    Input Parameter:
1308: .  ts - the TS context

1310:    Output Parameters:
1311: +  r - vector to hold computed right hand side (or NULL)
1312: .  func - the function to compute right hand side (or NULL)
1313: -  ctx - the function context (or NULL)

1315:    Level: advanced

1317: .seealso: TSSetRHSFunction(), SNESGetFunction()
1318: @*/
1319: PetscErrorCode TSGetRHSFunction(TS ts,Vec *r,TSRHSFunction *func,void **ctx)
1320: {
1321:   SNES           snes;
1322:   DM             dm;

1325:   TSGetSNES(ts,&snes);
1326:   SNESGetFunction(snes,r,NULL,NULL);
1327:   TSGetDM(ts,&dm);
1328:   DMTSGetRHSFunction(dm,func,ctx);
1329:   return 0;
1330: }

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

1336:    Logically Collective on TS

1338:    Input Parameters:
1339: +  ts  - the TS context obtained from TSCreate()
1340: .  Amat - (approximate) Jacobian matrix
1341: .  Pmat - matrix used to compute preconditioner (usually the same as Amat)
1342: .  f   - the Jacobian evaluation routine
1343: -  ctx - user-defined context for private data for the Jacobian evaluation routine (may be NULL)

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

1348: +  t    - time at step/stage being solved
1349: .  U    - state vector
1350: .  U_t  - time derivative of state vector
1351: .  a    - shift
1352: .  Amat - (approximate) Jacobian of F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t
1353: .  Pmat - matrix used for constructing preconditioner, usually the same as Amat
1354: -  ctx  - [optional] user-defined context for matrix evaluation routine

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

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

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

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

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

1374:    Level: beginner

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

1378: @*/
1379: PetscErrorCode  TSSetIJacobian(TS ts,Mat Amat,Mat Pmat,TSIJacobian f,void *ctx)
1380: {
1381:   SNES           snes;
1382:   DM             dm;


1390:   TSGetDM(ts,&dm);
1391:   DMTSSetIJacobian(dm,f,ctx);

1393:   TSGetSNES(ts,&snes);
1394:   SNESSetJacobian(snes,Amat,Pmat,SNESTSFormJacobian,ts);
1395:   return 0;
1396: }

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

1404:    Logically Collective

1406:    Input Parameters:
1407: +  ts - TS context obtained from TSCreate()
1408: -  reuse - PETSC_TRUE if the RHS Jacobian

1410:    Level: intermediate

1412: .seealso: TSSetRHSJacobian(), TSComputeRHSJacobianConstant()
1413: @*/
1414: PetscErrorCode TSRHSJacobianSetReuse(TS ts,PetscBool reuse)
1415: {
1416:   ts->rhsjacobian.reuse = reuse;
1417:   return 0;
1418: }

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

1423:    Logically Collective on TS

1425:    Input Parameters:
1426: +  ts  - the TS context obtained from TSCreate()
1427: .  F   - vector to hold the residual (or NULL to have it created internally)
1428: .  fun - the function evaluation routine
1429: -  ctx - user-defined context for private data for the function evaluation routine (may be NULL)

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

1434: +  t    - time at step/stage being solved
1435: .  U    - state vector
1436: .  U_t  - time derivative of state vector
1437: .  U_tt - second time derivative of state vector
1438: .  F    - function vector
1439: -  ctx  - [optional] user-defined context for matrix evaluation routine (may be NULL)

1441:    Level: beginner

1443: .seealso: TSSetI2Jacobian(), TSSetIFunction(), TSCreate(), TSSetRHSFunction()
1444: @*/
1445: PetscErrorCode TSSetI2Function(TS ts,Vec F,TSI2Function fun,void *ctx)
1446: {
1447:   DM             dm;

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

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

1460:   Not Collective

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

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

1470:   Level: advanced

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

1480:   TSGetSNES(ts,&snes);
1481:   SNESGetFunction(snes,r,NULL,NULL);
1482:   TSGetDM(ts,&dm);
1483:   DMTSGetI2Function(dm,fun,ctx);
1484:   return 0;
1485: }

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

1491:    Logically Collective on TS

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

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

1503: +  t    - time at step/stage being solved
1504: .  U    - state vector
1505: .  U_t  - time derivative of state vector
1506: .  U_tt - second time derivative of state vector
1507: .  v    - shift for U_t
1508: .  a    - shift for U_tt
1509: .  J    - Jacobian of G(U) = F(t,U,W+v*U,W'+a*U), equivalent to dF/dU + v*dF/dU_t  + a*dF/dU_tt
1510: .  P    - preconditioning matrix for J, may be same as J
1511: -  ctx  - [optional] user-defined context for matrix evaluation routine

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

1516:    The matrix dF/dU + v*dF/dU_t + a*dF/dU_tt you provide turns out to be
1517:    the Jacobian of G(U) = F(t,U,W+v*U,W'+a*U) where F(t,U,U_t,U_tt) = 0 is the DAE to be solved.
1518:    The time integrator internally approximates U_t by W+v*U and U_tt by W'+a*U  where the positive "shift"
1519:    parameters 'v' and 'a' and vectors W, W' depend on the integration method, step size, and past states.

1521:    Level: beginner

1523: .seealso: TSSetI2Function(), TSGetI2Jacobian()
1524: @*/
1525: PetscErrorCode TSSetI2Jacobian(TS ts,Mat J,Mat P,TSI2Jacobian jac,void *ctx)
1526: {
1527:   DM             dm;

1532:   TSSetIJacobian(ts,J,P,NULL,NULL);
1533:   TSGetDM(ts,&dm);
1534:   DMTSSetI2Jacobian(dm,jac,ctx);
1535:   return 0;
1536: }

1538: /*@C
1539:   TSGetI2Jacobian - Returns the implicit Jacobian at the present timestep.

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

1543:   Input Parameter:
1544: . ts  - The TS context obtained from TSCreate()

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

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

1555:   Level: advanced

1557: .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetStepNumber(), TSSetI2Jacobian(), TSGetI2Function(), TSCreate()

1559: @*/
1560: PetscErrorCode  TSGetI2Jacobian(TS ts,Mat *J,Mat *P,TSI2Jacobian *jac,void **ctx)
1561: {
1562:   SNES           snes;
1563:   DM             dm;

1565:   TSGetSNES(ts,&snes);
1566:   SNESSetUpMatrices(snes);
1567:   SNESGetJacobian(snes,J,P,NULL,NULL);
1568:   TSGetDM(ts,&dm);
1569:   DMTSGetI2Jacobian(dm,jac,ctx);
1570:   return 0;
1571: }

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

1576:   Collective on TS

1578:   Input Parameters:
1579: + ts - the TS context
1580: . t - current time
1581: . U - state vector
1582: . V - time derivative of state vector (U_t)
1583: - A - second time derivative of state vector (U_tt)

1585:   Output Parameter:
1586: . F - the residual vector

1588:   Note:
1589:   Most users should not need to explicitly call this routine, as it
1590:   is used internally within the nonlinear solvers.

1592:   Level: developer

1594: .seealso: TSSetI2Function(), TSGetI2Function()
1595: @*/
1596: PetscErrorCode TSComputeI2Function(TS ts,PetscReal t,Vec U,Vec V,Vec A,Vec F)
1597: {
1598:   DM             dm;
1599:   TSI2Function   I2Function;
1600:   void           *ctx;
1601:   TSRHSFunction  rhsfunction;


1609:   TSGetDM(ts,&dm);
1610:   DMTSGetI2Function(dm,&I2Function,&ctx);
1611:   DMTSGetRHSFunction(dm,&rhsfunction,NULL);

1613:   if (!I2Function) {
1614:     TSComputeIFunction(ts,t,U,A,F,PETSC_FALSE);
1615:     return 0;
1616:   }

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

1620:   PetscStackPush("TS user implicit function");
1621:   I2Function(ts,t,U,V,A,F,ctx);
1622:   PetscStackPop;

1624:   if (rhsfunction) {
1625:     Vec Frhs;
1626:     TSGetRHSVec_Private(ts,&Frhs);
1627:     TSComputeRHSFunction(ts,t,U,Frhs);
1628:     VecAXPY(F,-1,Frhs);
1629:   }

1631:   PetscLogEventEnd(TS_FunctionEval,ts,U,V,F);
1632:   return 0;
1633: }

1635: /*@
1636:   TSComputeI2Jacobian - Evaluates the Jacobian of the DAE

1638:   Collective on TS

1640:   Input Parameters:
1641: + ts - the TS context
1642: . t - current timestep
1643: . U - state vector
1644: . V - time derivative of state vector
1645: . A - second time derivative of state vector
1646: . shiftV - shift to apply, see note below
1647: - shiftA - shift to apply, see note below

1649:   Output Parameters:
1650: + J - Jacobian matrix
1651: - P - optional preconditioning matrix

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

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

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

1661:   Level: developer

1663: .seealso:  TSSetI2Jacobian()
1664: @*/
1665: PetscErrorCode TSComputeI2Jacobian(TS ts,PetscReal t,Vec U,Vec V,Vec A,PetscReal shiftV,PetscReal shiftA,Mat J,Mat P)
1666: {
1667:   DM             dm;
1668:   TSI2Jacobian   I2Jacobian;
1669:   void           *ctx;
1670:   TSRHSJacobian  rhsjacobian;


1679:   TSGetDM(ts,&dm);
1680:   DMTSGetI2Jacobian(dm,&I2Jacobian,&ctx);
1681:   DMTSGetRHSJacobian(dm,&rhsjacobian,NULL);

1683:   if (!I2Jacobian) {
1684:     TSComputeIJacobian(ts,t,U,A,shiftA,J,P,PETSC_FALSE);
1685:     return 0;
1686:   }

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

1690:   PetscStackPush("TS user implicit Jacobian");
1691:   I2Jacobian(ts,t,U,V,A,shiftV,shiftA,J,P,ctx);
1692:   PetscStackPop;

1694:   if (rhsjacobian) {
1695:     Mat Jrhs,Prhs;
1696:     TSGetRHSMats_Private(ts,&Jrhs,&Prhs);
1697:     TSComputeRHSJacobian(ts,t,U,Jrhs,Prhs);
1698:     MatAXPY(J,-1,Jrhs,ts->axpy_pattern);
1699:     if (P != J) MatAXPY(P,-1,Prhs,ts->axpy_pattern);
1700:   }

1702:   PetscLogEventEnd(TS_JacobianEval,ts,U,J,P);
1703:   return 0;
1704: }

1706: /*@C
1707:    TSSetTransientVariable - sets function to transform from state to transient variables

1709:    Logically Collective

1711:    Input Parameters:
1712: +  ts - time stepping context on which to change the transient variable
1713: .  tvar - a function that transforms to transient variables
1714: -  ctx - a context for tvar

1716:     Calling sequence of tvar:
1717: $     PetscErrorCode tvar(TS ts,Vec p,Vec c,void *ctx);

1719: +   ts - timestep context
1720: .   p - input vector (primitive form)
1721: .   c - output vector, transient variables (conservative form)
1722: -   ctx - [optional] user-defined function context

1724:    Level: advanced

1726:    Notes:
1727:    This is typically used to transform from primitive to conservative variables so that a time integrator (e.g., TSBDF)
1728:    can be conservative.  In this context, primitive variables P are used to model the state (e.g., because they lead to
1729:    well-conditioned formulations even in limiting cases such as low-Mach or zero porosity).  The transient variable is
1730:    C(P), specified by calling this function.  An IFunction thus receives arguments (P, Cdot) and the IJacobian must be
1731:    evaluated via the chain rule, as in

1733:      dF/dP + shift * dF/dCdot dC/dP.

1735: .seealso: DMTSSetTransientVariable(), DMTSGetTransientVariable(), TSSetIFunction(), TSSetIJacobian()
1736: @*/
1737: PetscErrorCode TSSetTransientVariable(TS ts,TSTransientVariable tvar,void *ctx)
1738: {
1739:   DM             dm;

1742:   TSGetDM(ts,&dm);
1743:   DMTSSetTransientVariable(dm,tvar,ctx);
1744:   return 0;
1745: }

1747: /*@
1748:    TSComputeTransientVariable - transforms state (primitive) variables to transient (conservative) variables

1750:    Logically Collective

1752:    Input Parameters:
1753: +  ts - TS on which to compute
1754: -  U - state vector to be transformed to transient variables

1756:    Output Parameters:
1757: .  C - transient (conservative) variable

1759:    Developer Notes:
1760:    If DMTSSetTransientVariable() has not been called, then C is not modified in this routine and C=NULL is allowed.
1761:    This makes it safe to call without a guard.  One can use TSHasTransientVariable() to check if transient variables are
1762:    being used.

1764:    Level: developer

1766: .seealso: DMTSSetTransientVariable(), TSComputeIFunction(), TSComputeIJacobian()
1767: @*/
1768: PetscErrorCode TSComputeTransientVariable(TS ts,Vec U,Vec C)
1769: {
1770:   DM             dm;
1771:   DMTS           dmts;

1775:   TSGetDM(ts,&dm);
1776:   DMGetDMTS(dm,&dmts);
1777:   if (dmts->ops->transientvar) {
1779:     (*dmts->ops->transientvar)(ts,U,C,dmts->transientvarctx);
1780:   }
1781:   return 0;
1782: }

1784: /*@
1785:    TSHasTransientVariable - determine whether transient variables have been set

1787:    Logically Collective

1789:    Input Parameters:
1790: .  ts - TS on which to compute

1792:    Output Parameters:
1793: .  has - PETSC_TRUE if transient variables have been set

1795:    Level: developer

1797: .seealso: DMTSSetTransientVariable(), TSComputeTransientVariable()
1798: @*/
1799: PetscErrorCode TSHasTransientVariable(TS ts,PetscBool *has)
1800: {
1801:   DM             dm;
1802:   DMTS           dmts;

1805:   TSGetDM(ts,&dm);
1806:   DMGetDMTS(dm,&dmts);
1807:   *has = dmts->ops->transientvar ? PETSC_TRUE : PETSC_FALSE;
1808:   return 0;
1809: }

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

1815:    Logically Collective on TS

1817:    Input Parameters:
1818: +  ts - the TS context obtained from TSCreate()
1819: .  u - the solution vector
1820: -  v - the time derivative vector

1822:    Level: beginner

1824: @*/
1825: PetscErrorCode  TS2SetSolution(TS ts,Vec u,Vec v)
1826: {
1830:   TSSetSolution(ts,u);
1831:   PetscObjectReference((PetscObject)v);
1832:   VecDestroy(&ts->vec_dot);
1833:   ts->vec_dot = v;
1834:   return 0;
1835: }

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

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

1845:    Input Parameter:
1846: .  ts - the TS context obtained from TSCreate()

1848:    Output Parameters:
1849: +  u - the vector containing the solution
1850: -  v - the vector containing the time derivative

1852:    Level: intermediate

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

1856: @*/
1857: PetscErrorCode  TS2GetSolution(TS ts,Vec *u,Vec *v)
1858: {
1862:   if (u) *u = ts->vec_sol;
1863:   if (v) *v = ts->vec_dot;
1864:   return 0;
1865: }

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

1870:   Collective on PetscViewer

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

1877:    Level: intermediate

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

1882:   Notes for advanced users:
1883:   Most users should not need to know the details of the binary storage
1884:   format, since TSLoad() and TSView() completely hide these details.
1885:   But for anyone who's interested, the standard binary matrix storage
1886:   format is
1887: .vb
1888:      has not yet been determined
1889: .ve

1891: .seealso: PetscViewerBinaryOpen(), TSView(), MatLoad(), VecLoad()
1892: @*/
1893: PetscErrorCode  TSLoad(TS ts, PetscViewer viewer)
1894: {
1895:   PetscBool      isbinary;
1896:   PetscInt       classid;
1897:   char           type[256];
1898:   DMTS           sdm;
1899:   DM             dm;

1903:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);

1906:   PetscViewerBinaryRead(viewer,&classid,1,NULL,PETSC_INT);
1908:   PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);
1909:   TSSetType(ts, type);
1910:   if (ts->ops->load) {
1911:     (*ts->ops->load)(ts,viewer);
1912:   }
1913:   DMCreate(PetscObjectComm((PetscObject)ts),&dm);
1914:   DMLoad(dm,viewer);
1915:   TSSetDM(ts,dm);
1916:   DMCreateGlobalVector(ts->dm,&ts->vec_sol);
1917:   VecLoad(ts->vec_sol,viewer);
1918:   DMGetDMTS(ts->dm,&sdm);
1919:   DMTSLoad(sdm,viewer);
1920:   return 0;
1921: }

1923: #include <petscdraw.h>
1924: #if defined(PETSC_HAVE_SAWS)
1925: #include <petscviewersaws.h>
1926: #endif

1928: /*@C
1929:    TSViewFromOptions - View from Options

1931:    Collective on TS

1933:    Input Parameters:
1934: +  A - the application ordering context
1935: .  obj - Optional object
1936: -  name - command line option

1938:    Level: intermediate
1939: .seealso:  TS, TSView, PetscObjectViewFromOptions(), TSCreate()
1940: @*/
1941: PetscErrorCode  TSViewFromOptions(TS A,PetscObject obj,const char name[])
1942: {
1944:   PetscObjectViewFromOptions((PetscObject)A,obj,name);
1945:   return 0;
1946: }

1948: /*@C
1949:     TSView - Prints the TS data structure.

1951:     Collective on TS

1953:     Input Parameters:
1954: +   ts - the TS context obtained from TSCreate()
1955: -   viewer - visualization context

1957:     Options Database Key:
1958: .   -ts_view - calls TSView() at end of TSStep()

1960:     Notes:
1961:     The available visualization contexts include
1962: +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
1963: -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
1964:          output where only the first processor opens
1965:          the file.  All other processors send their
1966:          data to the first processor to print.

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

1971:     In the debugger you can do "call TSView(ts,0)" to display the TS solver. (The same holds for any PETSc object viewer).

1973:     Level: beginner

1975: .seealso: PetscViewerASCIIOpen()
1976: @*/
1977: PetscErrorCode  TSView(TS ts,PetscViewer viewer)
1978: {
1979:   TSType         type;
1980:   PetscBool      iascii,isstring,isundials,isbinary,isdraw;
1981:   DMTS           sdm;
1982: #if defined(PETSC_HAVE_SAWS)
1983:   PetscBool      issaws;
1984: #endif

1987:   if (!viewer) {
1988:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ts),&viewer);
1989:   }

1993:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1994:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
1995:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
1996:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
1997: #if defined(PETSC_HAVE_SAWS)
1998:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSAWS,&issaws);
1999: #endif
2000:   if (iascii) {
2001:     PetscObjectPrintClassNamePrefixType((PetscObject)ts,viewer);
2002:     if (ts->ops->view) {
2003:       PetscViewerASCIIPushTab(viewer);
2004:       (*ts->ops->view)(ts,viewer);
2005:       PetscViewerASCIIPopTab(viewer);
2006:     }
2007:     if (ts->max_steps < PETSC_MAX_INT) {
2008:       PetscViewerASCIIPrintf(viewer,"  maximum steps=%D\n",ts->max_steps);
2009:     }
2010:     if (ts->max_time < PETSC_MAX_REAL) {
2011:       PetscViewerASCIIPrintf(viewer,"  maximum time=%g\n",(double)ts->max_time);
2012:     }
2013:     if (ts->ifuncs) {
2014:       PetscViewerASCIIPrintf(viewer,"  total number of I function evaluations=%D\n",ts->ifuncs);
2015:     }
2016:     if (ts->ijacs) {
2017:       PetscViewerASCIIPrintf(viewer,"  total number of I Jacobian evaluations=%D\n",ts->ijacs);
2018:     }
2019:     if (ts->rhsfuncs) {
2020:       PetscViewerASCIIPrintf(viewer,"  total number of RHS function evaluations=%D\n",ts->rhsfuncs);
2021:     }
2022:     if (ts->rhsjacs) {
2023:       PetscViewerASCIIPrintf(viewer,"  total number of RHS Jacobian evaluations=%D\n",ts->rhsjacs);
2024:     }
2025:     if (ts->usessnes) {
2026:       PetscBool lin;
2027:       if (ts->problem_type == TS_NONLINEAR) {
2028:         PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solver iterations=%D\n",ts->snes_its);
2029:       }
2030:       PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%D\n",ts->ksp_its);
2031:       PetscObjectTypeCompareAny((PetscObject)ts->snes,&lin,SNESKSPONLY,SNESKSPTRANSPOSEONLY,"");
2032:       PetscViewerASCIIPrintf(viewer,"  total number of %slinear solve failures=%D\n",lin ? "" : "non",ts->num_snes_failures);
2033:     }
2034:     PetscViewerASCIIPrintf(viewer,"  total number of rejected steps=%D\n",ts->reject);
2035:     if (ts->vrtol) {
2036:       PetscViewerASCIIPrintf(viewer,"  using vector of relative error tolerances, ");
2037:     } else {
2038:       PetscViewerASCIIPrintf(viewer,"  using relative error tolerance of %g, ",(double)ts->rtol);
2039:     }
2040:     if (ts->vatol) {
2041:       PetscViewerASCIIPrintf(viewer,"  using vector of absolute error tolerances\n");
2042:     } else {
2043:       PetscViewerASCIIPrintf(viewer,"  using absolute error tolerance of %g\n",(double)ts->atol);
2044:     }
2045:     PetscViewerASCIIPushTab(viewer);
2046:     TSAdaptView(ts->adapt,viewer);
2047:     PetscViewerASCIIPopTab(viewer);
2048:   } else if (isstring) {
2049:     TSGetType(ts,&type);
2050:     PetscViewerStringSPrintf(viewer," TSType: %-7.7s",type);
2051:     if (ts->ops->view) (*ts->ops->view)(ts,viewer);
2052:   } else if (isbinary) {
2053:     PetscInt    classid = TS_FILE_CLASSID;
2054:     MPI_Comm    comm;
2055:     PetscMPIInt rank;
2056:     char        type[256];

2058:     PetscObjectGetComm((PetscObject)ts,&comm);
2059:     MPI_Comm_rank(comm,&rank);
2060:     if (rank == 0) {
2061:       PetscViewerBinaryWrite(viewer,&classid,1,PETSC_INT);
2062:       PetscStrncpy(type,((PetscObject)ts)->type_name,256);
2063:       PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR);
2064:     }
2065:     if (ts->ops->view) {
2066:       (*ts->ops->view)(ts,viewer);
2067:     }
2068:     if (ts->adapt) TSAdaptView(ts->adapt,viewer);
2069:     DMView(ts->dm,viewer);
2070:     VecView(ts->vec_sol,viewer);
2071:     DMGetDMTS(ts->dm,&sdm);
2072:     DMTSView(sdm,viewer);
2073:   } else if (isdraw) {
2074:     PetscDraw draw;
2075:     char      str[36];
2076:     PetscReal x,y,bottom,h;

2078:     PetscViewerDrawGetDraw(viewer,0,&draw);
2079:     PetscDrawGetCurrentPoint(draw,&x,&y);
2080:     PetscStrcpy(str,"TS: ");
2081:     PetscStrcat(str,((PetscObject)ts)->type_name);
2082:     PetscDrawStringBoxed(draw,x,y,PETSC_DRAW_BLACK,PETSC_DRAW_BLACK,str,NULL,&h);
2083:     bottom = y - h;
2084:     PetscDrawPushCurrentPoint(draw,x,bottom);
2085:     if (ts->ops->view) {
2086:       (*ts->ops->view)(ts,viewer);
2087:     }
2088:     if (ts->adapt) TSAdaptView(ts->adapt,viewer);
2089:     if (ts->snes)  SNESView(ts->snes,viewer);
2090:     PetscDrawPopCurrentPoint(draw);
2091: #if defined(PETSC_HAVE_SAWS)
2092:   } else if (issaws) {
2093:     PetscMPIInt rank;
2094:     const char  *name;

2096:     PetscObjectGetName((PetscObject)ts,&name);
2097:     MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
2098:     if (!((PetscObject)ts)->amsmem && rank == 0) {
2099:       char       dir[1024];

2101:       PetscObjectViewSAWs((PetscObject)ts,viewer);
2102:       PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/time_step",name);
2103:       PetscStackCallSAWs(SAWs_Register,(dir,&ts->steps,1,SAWs_READ,SAWs_INT));
2104:       PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/time",name);
2105:       PetscStackCallSAWs(SAWs_Register,(dir,&ts->ptime,1,SAWs_READ,SAWs_DOUBLE));
2106:     }
2107:     if (ts->ops->view) {
2108:       (*ts->ops->view)(ts,viewer);
2109:     }
2110: #endif
2111:   }
2112:   if (ts->snes && ts->usessnes)  {
2113:     PetscViewerASCIIPushTab(viewer);
2114:     SNESView(ts->snes,viewer);
2115:     PetscViewerASCIIPopTab(viewer);
2116:   }
2117:   DMGetDMTS(ts->dm,&sdm);
2118:   DMTSView(sdm,viewer);

2120:   PetscViewerASCIIPushTab(viewer);
2121:   PetscObjectTypeCompare((PetscObject)ts,TSSUNDIALS,&isundials);
2122:   PetscViewerASCIIPopTab(viewer);
2123:   return 0;
2124: }

2126: /*@
2127:    TSSetApplicationContext - Sets an optional user-defined context for
2128:    the timesteppers.

2130:    Logically Collective on TS

2132:    Input Parameters:
2133: +  ts - the TS context obtained from TSCreate()
2134: -  usrP - optional user context

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

2140:    Level: intermediate

2142: .seealso: TSGetApplicationContext()
2143: @*/
2144: PetscErrorCode  TSSetApplicationContext(TS ts,void *usrP)
2145: {
2147:   ts->user = usrP;
2148:   return 0;
2149: }

2151: /*@
2152:     TSGetApplicationContext - Gets the user-defined context for the
2153:     timestepper.

2155:     Not Collective

2157:     Input Parameter:
2158: .   ts - the TS context obtained from TSCreate()

2160:     Output Parameter:
2161: .   usrP - user context

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

2167:     Level: intermediate

2169: .seealso: TSSetApplicationContext()
2170: @*/
2171: PetscErrorCode  TSGetApplicationContext(TS ts,void *usrP)
2172: {
2174:   *(void**)usrP = ts->user;
2175:   return 0;
2176: }

2178: /*@
2179:    TSGetStepNumber - Gets the number of steps completed.

2181:    Not Collective

2183:    Input Parameter:
2184: .  ts - the TS context obtained from TSCreate()

2186:    Output Parameter:
2187: .  steps - number of steps completed so far

2189:    Level: intermediate

2191: .seealso: TSGetTime(), TSGetTimeStep(), TSSetPreStep(), TSSetPreStage(), TSSetPostStage(), TSSetPostStep()
2192: @*/
2193: PetscErrorCode TSGetStepNumber(TS ts,PetscInt *steps)
2194: {
2197:   *steps = ts->steps;
2198:   return 0;
2199: }

2201: /*@
2202:    TSSetStepNumber - Sets the number of steps completed.

2204:    Logically Collective on TS

2206:    Input Parameters:
2207: +  ts - the TS context
2208: -  steps - number of steps completed so far

2210:    Notes:
2211:    For most uses of the TS solvers the user need not explicitly call
2212:    TSSetStepNumber(), as the step counter is appropriately updated in
2213:    TSSolve()/TSStep()/TSRollBack(). Power users may call this routine to
2214:    reinitialize timestepping by setting the step counter to zero (and time
2215:    to the initial time) to solve a similar problem with different initial
2216:    conditions or parameters. Other possible use case is to continue
2217:    timestepping from a previously interrupted run in such a way that TS
2218:    monitors will be called with a initial nonzero step counter.

2220:    Level: advanced

2222: .seealso: TSGetStepNumber(), TSSetTime(), TSSetTimeStep(), TSSetSolution()
2223: @*/
2224: PetscErrorCode TSSetStepNumber(TS ts,PetscInt steps)
2225: {
2229:   ts->steps = steps;
2230:   return 0;
2231: }

2233: /*@
2234:    TSSetTimeStep - Allows one to reset the timestep at any time,
2235:    useful for simple pseudo-timestepping codes.

2237:    Logically Collective on TS

2239:    Input Parameters:
2240: +  ts - the TS context obtained from TSCreate()
2241: -  time_step - the size of the timestep

2243:    Level: intermediate

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

2247: @*/
2248: PetscErrorCode  TSSetTimeStep(TS ts,PetscReal time_step)
2249: {
2252:   ts->time_step = time_step;
2253:   return 0;
2254: }

2256: /*@
2257:    TSSetExactFinalTime - Determines whether to adapt the final time step to
2258:      match the exact final time, interpolate solution to the exact final time,
2259:      or just return at the final time TS computed.

2261:   Logically Collective on TS

2263:    Input Parameters:
2264: +   ts - the time-step context
2265: -   eftopt - exact final time option

2267: $  TS_EXACTFINALTIME_STEPOVER    - Don't do anything if final time is exceeded
2268: $  TS_EXACTFINALTIME_INTERPOLATE - Interpolate back to final time
2269: $  TS_EXACTFINALTIME_MATCHSTEP - Adapt final time step to match the final time

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

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

2277:    Level: beginner

2279: .seealso: TSExactFinalTimeOption, TSGetExactFinalTime()
2280: @*/
2281: PetscErrorCode TSSetExactFinalTime(TS ts,TSExactFinalTimeOption eftopt)
2282: {
2285:   ts->exact_final_time = eftopt;
2286:   return 0;
2287: }

2289: /*@
2290:    TSGetExactFinalTime - Gets the exact final time option.

2292:    Not Collective

2294:    Input Parameter:
2295: .  ts - the TS context

2297:    Output Parameter:
2298: .  eftopt - exact final time option

2300:    Level: beginner

2302: .seealso: TSExactFinalTimeOption, TSSetExactFinalTime()
2303: @*/
2304: PetscErrorCode TSGetExactFinalTime(TS ts,TSExactFinalTimeOption *eftopt)
2305: {
2308:   *eftopt = ts->exact_final_time;
2309:   return 0;
2310: }

2312: /*@
2313:    TSGetTimeStep - Gets the current timestep size.

2315:    Not Collective

2317:    Input Parameter:
2318: .  ts - the TS context obtained from TSCreate()

2320:    Output Parameter:
2321: .  dt - the current timestep size

2323:    Level: intermediate

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

2327: @*/
2328: PetscErrorCode  TSGetTimeStep(TS ts,PetscReal *dt)
2329: {
2332:   *dt = ts->time_step;
2333:   return 0;
2334: }

2336: /*@
2337:    TSGetSolution - Returns the solution at the present timestep. It
2338:    is valid to call this routine inside the function that you are evaluating
2339:    in order to move to the new timestep. This vector not changed until
2340:    the solution at the next timestep has been calculated.

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

2344:    Input Parameter:
2345: .  ts - the TS context obtained from TSCreate()

2347:    Output Parameter:
2348: .  v - the vector containing the solution

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

2353:    Level: intermediate

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

2357: @*/
2358: PetscErrorCode  TSGetSolution(TS ts,Vec *v)
2359: {
2362:   *v = ts->vec_sol;
2363:   return 0;
2364: }

2366: /*@
2367:    TSGetSolutionComponents - Returns any solution components at the present
2368:    timestep, if available for the time integration method being used.
2369:    Solution components are quantities that share the same size and
2370:    structure as the solution vector.

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

2374:    Parameters :
2375: +  ts - the TS context obtained from TSCreate() (input parameter).
2376: .  n - If v is PETSC_NULL, then the number of solution components is
2377:        returned through n, else the n-th solution component is
2378:        returned in v.
2379: -  v - the vector containing the n-th solution component
2380:        (may be PETSC_NULL to use this function to find out
2381:         the number of solutions components).

2383:    Level: advanced

2385: .seealso: TSGetSolution()

2387: @*/
2388: PetscErrorCode  TSGetSolutionComponents(TS ts,PetscInt *n,Vec *v)
2389: {
2391:   if (!ts->ops->getsolutioncomponents) *n = 0;
2392:   else {
2393:     (*ts->ops->getsolutioncomponents)(ts,n,v);
2394:   }
2395:   return 0;
2396: }

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

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

2404:    Parameters :
2405: +  ts - the TS context obtained from TSCreate() (input parameter).
2406: -  v - the vector containing the auxiliary solution

2408:    Level: intermediate

2410: .seealso: TSGetSolution()

2412: @*/
2413: PetscErrorCode  TSGetAuxSolution(TS ts,Vec *v)
2414: {
2416:   if (ts->ops->getauxsolution) {
2417:     (*ts->ops->getauxsolution)(ts,v);
2418:   } else {
2419:     VecZeroEntries(*v);
2420:   }
2421:   return 0;
2422: }

2424: /*@
2425:    TSGetTimeError - Returns the estimated error vector, if the chosen
2426:    TSType has an error estimation functionality.

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

2430:    Note: MUST call after TSSetUp()

2432:    Parameters :
2433: +  ts - the TS context obtained from TSCreate() (input parameter).
2434: .  n - current estimate (n=0) or previous one (n=-1)
2435: -  v - the vector containing the error (same size as the solution).

2437:    Level: intermediate

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

2441: @*/
2442: PetscErrorCode  TSGetTimeError(TS ts,PetscInt n,Vec *v)
2443: {
2445:   if (ts->ops->gettimeerror) {
2446:     (*ts->ops->gettimeerror)(ts,n,v);
2447:   } else {
2448:     VecZeroEntries(*v);
2449:   }
2450:   return 0;
2451: }

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

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

2460:    Parameters :
2461: +  ts - the TS context obtained from TSCreate() (input parameter).
2462: -  v - the vector containing the error (same size as the solution).

2464:    Level: intermediate

2466: .seealso: TSSetSolution(), TSGetTimeError)

2468: @*/
2469: PetscErrorCode  TSSetTimeError(TS ts,Vec v)
2470: {
2473:   if (ts->ops->settimeerror) {
2474:     (*ts->ops->settimeerror)(ts,v);
2475:   }
2476:   return 0;
2477: }

2479: /* ----- Routines to initialize and destroy a timestepper ---- */
2480: /*@
2481:   TSSetProblemType - Sets the type of problem to be solved.

2483:   Not collective

2485:   Input Parameters:
2486: + ts   - The TS
2487: - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
2488: .vb
2489:          U_t - A U = 0      (linear)
2490:          U_t - A(t) U = 0   (linear)
2491:          F(t,U,U_t) = 0     (nonlinear)
2492: .ve

2494:    Level: beginner

2496: .seealso: TSSetUp(), TSProblemType, TS
2497: @*/
2498: PetscErrorCode  TSSetProblemType(TS ts, TSProblemType type)
2499: {
2501:   ts->problem_type = type;
2502:   if (type == TS_LINEAR) {
2503:     SNES snes;
2504:     TSGetSNES(ts,&snes);
2505:     SNESSetType(snes,SNESKSPONLY);
2506:   }
2507:   return 0;
2508: }

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

2513:   Not collective

2515:   Input Parameter:
2516: . ts   - The TS

2518:   Output Parameter:
2519: . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
2520: .vb
2521:          M U_t = A U
2522:          M(t) U_t = A(t) U
2523:          F(t,U,U_t)
2524: .ve

2526:    Level: beginner

2528: .seealso: TSSetUp(), TSProblemType, TS
2529: @*/
2530: PetscErrorCode  TSGetProblemType(TS ts, TSProblemType *type)
2531: {
2534:   *type = ts->problem_type;
2535:   return 0;
2536: }

2538: /*
2539:     Attempt to check/preset a default value for the exact final time option. This is needed at the beginning of TSSolve() and in TSSetUp()
2540: */
2541: static PetscErrorCode TSSetExactFinalTimeDefault(TS ts)
2542: {
2543:   PetscBool      isnone;

2545:   TSGetAdapt(ts,&ts->adapt);
2546:   TSAdaptSetDefaultType(ts->adapt,ts->default_adapt_type);

2548:   PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&isnone);
2549:   if (!isnone && ts->exact_final_time == TS_EXACTFINALTIME_UNSPECIFIED) {
2550:     ts->exact_final_time = TS_EXACTFINALTIME_MATCHSTEP;
2551:   } else if (ts->exact_final_time == TS_EXACTFINALTIME_UNSPECIFIED) {
2552:     ts->exact_final_time = TS_EXACTFINALTIME_INTERPOLATE;
2553:   }
2554:   return 0;
2555: }

2557: /*@
2558:    TSSetUp - Sets up the internal data structures for the later use of a timestepper.

2560:    Collective on TS

2562:    Input Parameter:
2563: .  ts - the TS context obtained from TSCreate()

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

2572:    Level: advanced

2574: .seealso: TSCreate(), TSStep(), TSDestroy(), TSSolve()
2575: @*/
2576: PetscErrorCode  TSSetUp(TS ts)
2577: {
2578:   DM             dm;
2579:   PetscErrorCode (*func)(SNES,Vec,Vec,void*);
2580:   PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*);
2581:   TSIFunction    ifun;
2582:   TSIJacobian    ijac;
2583:   TSI2Jacobian   i2jac;
2584:   TSRHSJacobian  rhsjac;

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

2589:   if (!((PetscObject)ts)->type_name) {
2590:     TSGetIFunction(ts,NULL,&ifun,NULL);
2591:     TSSetType(ts,ifun ? TSBEULER : TSEULER);
2592:   }

2594:   if (!ts->vec_sol) {
2595:     if (ts->dm) {
2596:       DMCreateGlobalVector(ts->dm,&ts->vec_sol);
2597:     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first");
2598:   }

2600:   if (!ts->Jacp && ts->Jacprhs) { /* IJacobianP shares the same matrix with RHSJacobianP if only RHSJacobianP is provided */
2601:     PetscObjectReference((PetscObject)ts->Jacprhs);
2602:     ts->Jacp = ts->Jacprhs;
2603:   }

2605:   if (ts->quadraturets) {
2606:     TSSetUp(ts->quadraturets);
2607:     VecDestroy(&ts->vec_costintegrand);
2608:     VecDuplicate(ts->quadraturets->vec_sol,&ts->vec_costintegrand);
2609:   }

2611:   TSGetRHSJacobian(ts,NULL,NULL,&rhsjac,NULL);
2612:   if (rhsjac == TSComputeRHSJacobianConstant) {
2613:     Mat Amat,Pmat;
2614:     SNES snes;
2615:     TSGetSNES(ts,&snes);
2616:     SNESGetJacobian(snes,&Amat,&Pmat,NULL,NULL);
2617:     /* Matching matrices implies that an IJacobian is NOT set, because if it had been set, the IJacobian's matrix would
2618:      * have displaced the RHS matrix */
2619:     if (Amat && Amat == ts->Arhs) {
2620:       /* we need to copy the values of the matrix because for the constant Jacobian case the user will never set the numerical values in this new location */
2621:       MatDuplicate(ts->Arhs,MAT_COPY_VALUES,&Amat);
2622:       SNESSetJacobian(snes,Amat,NULL,NULL,NULL);
2623:       MatDestroy(&Amat);
2624:     }
2625:     if (Pmat && Pmat == ts->Brhs) {
2626:       MatDuplicate(ts->Brhs,MAT_COPY_VALUES,&Pmat);
2627:       SNESSetJacobian(snes,NULL,Pmat,NULL,NULL);
2628:       MatDestroy(&Pmat);
2629:     }
2630:   }

2632:   TSGetAdapt(ts,&ts->adapt);
2633:   TSAdaptSetDefaultType(ts->adapt,ts->default_adapt_type);

2635:   if (ts->ops->setup) {
2636:     (*ts->ops->setup)(ts);
2637:   }

2639:   TSSetExactFinalTimeDefault(ts);

2641:   /* In the case where we've set a DMTSFunction or what have you, we need the default SNESFunction
2642:      to be set right but can't do it elsewhere due to the overreliance on ctx=ts.
2643:    */
2644:   TSGetDM(ts,&dm);
2645:   DMSNESGetFunction(dm,&func,NULL);
2646:   if (!func) {
2647:     DMSNESSetFunction(dm,SNESTSFormFunction,ts);
2648:   }
2649:   /* If the SNES doesn't have a jacobian set and the TS has an ijacobian or rhsjacobian set, set the SNES to use it.
2650:      Otherwise, the SNES will use coloring internally to form the Jacobian.
2651:    */
2652:   DMSNESGetJacobian(dm,&jac,NULL);
2653:   DMTSGetIJacobian(dm,&ijac,NULL);
2654:   DMTSGetI2Jacobian(dm,&i2jac,NULL);
2655:   DMTSGetRHSJacobian(dm,&rhsjac,NULL);
2656:   if (!jac && (ijac || i2jac || rhsjac)) {
2657:     DMSNESSetJacobian(dm,SNESTSFormJacobian,ts);
2658:   }

2660:   /* if time integration scheme has a starting method, call it */
2661:   if (ts->ops->startingmethod) {
2662:     (*ts->ops->startingmethod)(ts);
2663:   }

2665:   ts->setupcalled = PETSC_TRUE;
2666:   return 0;
2667: }

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

2672:    Collective on TS

2674:    Input Parameter:
2675: .  ts - the TS context obtained from TSCreate()

2677:    Level: beginner

2679: .seealso: TSCreate(), TSSetup(), TSDestroy()
2680: @*/
2681: PetscErrorCode  TSReset(TS ts)
2682: {
2683:   TS_RHSSplitLink ilink = ts->tsrhssplit,next;


2687:   if (ts->ops->reset) {
2688:     (*ts->ops->reset)(ts);
2689:   }
2690:   if (ts->snes) SNESReset(ts->snes);
2691:   if (ts->adapt) TSAdaptReset(ts->adapt);

2693:   MatDestroy(&ts->Arhs);
2694:   MatDestroy(&ts->Brhs);
2695:   VecDestroy(&ts->Frhs);
2696:   VecDestroy(&ts->vec_sol);
2697:   VecDestroy(&ts->vec_dot);
2698:   VecDestroy(&ts->vatol);
2699:   VecDestroy(&ts->vrtol);
2700:   VecDestroyVecs(ts->nwork,&ts->work);

2702:   MatDestroy(&ts->Jacprhs);
2703:   MatDestroy(&ts->Jacp);
2704:   if (ts->forward_solve) {
2705:     TSForwardReset(ts);
2706:   }
2707:   if (ts->quadraturets) {
2708:     TSReset(ts->quadraturets);
2709:     VecDestroy(&ts->vec_costintegrand);
2710:   }
2711:   while (ilink) {
2712:     next = ilink->next;
2713:     TSDestroy(&ilink->ts);
2714:     PetscFree(ilink->splitname);
2715:     ISDestroy(&ilink->is);
2716:     PetscFree(ilink);
2717:     ilink = next;
2718:   }
2719:   ts->tsrhssplit = NULL;
2720:   ts->num_rhs_splits = 0;
2721:   ts->setupcalled = PETSC_FALSE;
2722:   return 0;
2723: }

2725: /*@C
2726:    TSDestroy - Destroys the timestepper context that was created
2727:    with TSCreate().

2729:    Collective on TS

2731:    Input Parameter:
2732: .  ts - the TS context obtained from TSCreate()

2734:    Level: beginner

2736: .seealso: TSCreate(), TSSetUp(), TSSolve()
2737: @*/
2738: PetscErrorCode  TSDestroy(TS *ts)
2739: {
2740:   if (!*ts) return 0;
2742:   if (--((PetscObject)(*ts))->refct > 0) {*ts = NULL; return 0;}

2744:   TSReset(*ts);
2745:   TSAdjointReset(*ts);
2746:   if ((*ts)->forward_solve) {
2747:     TSForwardReset(*ts);
2748:   }
2749:   /* if memory was published with SAWs then destroy it */
2750:   PetscObjectSAWsViewOff((PetscObject)*ts);
2751:   if ((*ts)->ops->destroy) (*(*ts)->ops->destroy)((*ts));

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

2755:   TSAdaptDestroy(&(*ts)->adapt);
2756:   TSEventDestroy(&(*ts)->event);

2758:   SNESDestroy(&(*ts)->snes);
2759:   DMDestroy(&(*ts)->dm);
2760:   TSMonitorCancel((*ts));
2761:   TSAdjointMonitorCancel((*ts));

2763:   TSDestroy(&(*ts)->quadraturets);
2764:   PetscHeaderDestroy(ts);
2765:   return 0;
2766: }

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

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

2774:    Input Parameter:
2775: .  ts - the TS context obtained from TSCreate()

2777:    Output Parameter:
2778: .  snes - the nonlinear solver context

2780:    Notes:
2781:    The user can then directly manipulate the SNES context to set various
2782:    options, etc.  Likewise, the user can then extract and manipulate the
2783:    KSP, KSP, and PC contexts as well.

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

2788:    Level: beginner

2790: @*/
2791: PetscErrorCode  TSGetSNES(TS ts,SNES *snes)
2792: {
2795:   if (!ts->snes) {
2796:     SNESCreate(PetscObjectComm((PetscObject)ts),&ts->snes);
2797:     PetscObjectSetOptions((PetscObject)ts->snes,((PetscObject)ts)->options);
2798:     SNESSetFunction(ts->snes,NULL,SNESTSFormFunction,ts);
2799:     PetscLogObjectParent((PetscObject)ts,(PetscObject)ts->snes);
2800:     PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);
2801:     if (ts->dm) SNESSetDM(ts->snes,ts->dm);
2802:     if (ts->problem_type == TS_LINEAR) {
2803:       SNESSetType(ts->snes,SNESKSPONLY);
2804:     }
2805:   }
2806:   *snes = ts->snes;
2807:   return 0;
2808: }

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

2813:    Collective

2815:    Input Parameters:
2816: +  ts - the TS context obtained from TSCreate()
2817: -  snes - the nonlinear solver context

2819:    Notes:
2820:    Most users should have the TS created by calling TSGetSNES()

2822:    Level: developer

2824: @*/
2825: PetscErrorCode TSSetSNES(TS ts,SNES snes)
2826: {
2827:   PetscErrorCode (*func)(SNES,Vec,Mat,Mat,void*);

2831:   PetscObjectReference((PetscObject)snes);
2832:   SNESDestroy(&ts->snes);

2834:   ts->snes = snes;

2836:   SNESSetFunction(ts->snes,NULL,SNESTSFormFunction,ts);
2837:   SNESGetJacobian(ts->snes,NULL,NULL,&func,NULL);
2838:   if (func == SNESTSFormJacobian) {
2839:     SNESSetJacobian(ts->snes,NULL,NULL,SNESTSFormJacobian,ts);
2840:   }
2841:   return 0;
2842: }

2844: /*@
2845:    TSGetKSP - Returns the KSP (linear solver) associated with
2846:    a TS (timestepper) context.

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

2850:    Input Parameter:
2851: .  ts - the TS context obtained from TSCreate()

2853:    Output Parameter:
2854: .  ksp - the nonlinear solver context

2856:    Notes:
2857:    The user can then directly manipulate the KSP context to set various
2858:    options, etc.  Likewise, the user can then extract and manipulate the
2859:    KSP and PC contexts as well.

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

2864:    Level: beginner

2866: @*/
2867: PetscErrorCode  TSGetKSP(TS ts,KSP *ksp)
2868: {
2869:   SNES           snes;

2875:   TSGetSNES(ts,&snes);
2876:   SNESGetKSP(snes,ksp);
2877:   return 0;
2878: }

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

2882: /*@
2883:    TSSetMaxSteps - Sets the maximum number of steps to use.

2885:    Logically Collective on TS

2887:    Input Parameters:
2888: +  ts - the TS context obtained from TSCreate()
2889: -  maxsteps - maximum number of steps to use

2891:    Options Database Keys:
2892: .  -ts_max_steps <maxsteps> - Sets maxsteps

2894:    Notes:
2895:    The default maximum number of steps is 5000

2897:    Level: intermediate

2899: .seealso: TSGetMaxSteps(), TSSetMaxTime(), TSSetExactFinalTime()
2900: @*/
2901: PetscErrorCode TSSetMaxSteps(TS ts,PetscInt maxsteps)
2902: {
2906:   ts->max_steps = maxsteps;
2907:   return 0;
2908: }

2910: /*@
2911:    TSGetMaxSteps - Gets the maximum number of steps to use.

2913:    Not Collective

2915:    Input Parameters:
2916: .  ts - the TS context obtained from TSCreate()

2918:    Output Parameter:
2919: .  maxsteps - maximum number of steps to use

2921:    Level: advanced

2923: .seealso: TSSetMaxSteps(), TSGetMaxTime(), TSSetMaxTime()
2924: @*/
2925: PetscErrorCode TSGetMaxSteps(TS ts,PetscInt *maxsteps)
2926: {
2929:   *maxsteps = ts->max_steps;
2930:   return 0;
2931: }

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

2936:    Logically Collective on TS

2938:    Input Parameters:
2939: +  ts - the TS context obtained from TSCreate()
2940: -  maxtime - final time to step to

2942:    Options Database Keys:
2943: .  -ts_max_time <maxtime> - Sets maxtime

2945:    Notes:
2946:    The default maximum time is 5.0

2948:    Level: intermediate

2950: .seealso: TSGetMaxTime(), TSSetMaxSteps(), TSSetExactFinalTime()
2951: @*/
2952: PetscErrorCode TSSetMaxTime(TS ts,PetscReal maxtime)
2953: {
2956:   ts->max_time = maxtime;
2957:   return 0;
2958: }

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

2963:    Not Collective

2965:    Input Parameters:
2966: .  ts - the TS context obtained from TSCreate()

2968:    Output Parameter:
2969: .  maxtime - final time to step to

2971:    Level: advanced

2973: .seealso: TSSetMaxTime(), TSGetMaxSteps(), TSSetMaxSteps()
2974: @*/
2975: PetscErrorCode TSGetMaxTime(TS ts,PetscReal *maxtime)
2976: {
2979:   *maxtime = ts->max_time;
2980:   return 0;
2981: }

2983: /*@
2984:    TSSetInitialTimeStep - Deprecated, use TSSetTime() and TSSetTimeStep().

2986:    Level: deprecated

2988: @*/
2989: PetscErrorCode  TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step)
2990: {
2992:   TSSetTime(ts,initial_time);
2993:   TSSetTimeStep(ts,time_step);
2994:   return 0;
2995: }

2997: /*@
2998:    TSGetDuration - Deprecated, use TSGetMaxSteps() and TSGetMaxTime().

3000:    Level: deprecated

3002: @*/
3003: PetscErrorCode TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime)
3004: {
3006:   if (maxsteps) {
3008:     *maxsteps = ts->max_steps;
3009:   }
3010:   if (maxtime) {
3012:     *maxtime = ts->max_time;
3013:   }
3014:   return 0;
3015: }

3017: /*@
3018:    TSSetDuration - Deprecated, use TSSetMaxSteps() and TSSetMaxTime().

3020:    Level: deprecated

3022: @*/
3023: PetscErrorCode TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime)
3024: {
3028:   if (maxsteps >= 0) ts->max_steps = maxsteps;
3029:   if (maxtime != PETSC_DEFAULT) ts->max_time = maxtime;
3030:   return 0;
3031: }

3033: /*@
3034:    TSGetTimeStepNumber - Deprecated, use TSGetStepNumber().

3036:    Level: deprecated

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

3041: /*@
3042:    TSGetTotalSteps - Deprecated, use TSGetStepNumber().

3044:    Level: deprecated

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

3049: /*@
3050:    TSSetSolution - Sets the initial solution vector
3051:    for use by the TS routines.

3053:    Logically Collective on TS

3055:    Input Parameters:
3056: +  ts - the TS context obtained from TSCreate()
3057: -  u - the solution vector

3059:    Level: beginner

3061: .seealso: TSSetSolutionFunction(), TSGetSolution(), TSCreate()
3062: @*/
3063: PetscErrorCode  TSSetSolution(TS ts,Vec u)
3064: {
3065:   DM             dm;

3069:   PetscObjectReference((PetscObject)u);
3070:   VecDestroy(&ts->vec_sol);
3071:   ts->vec_sol = u;

3073:   TSGetDM(ts,&dm);
3074:   DMShellSetGlobalVector(dm,u);
3075:   return 0;
3076: }

3078: /*@C
3079:   TSSetPreStep - Sets the general-purpose function
3080:   called once at the beginning of each time step.

3082:   Logically Collective on TS

3084:   Input Parameters:
3085: + ts   - The TS context obtained from TSCreate()
3086: - func - The function

3088:   Calling sequence of func:
3089: .vb
3090:   PetscErrorCode func (TS ts);
3091: .ve

3093:   Level: intermediate

3095: .seealso: TSSetPreStage(), TSSetPostStage(), TSSetPostStep(), TSStep(), TSRestartStep()
3096: @*/
3097: PetscErrorCode  TSSetPreStep(TS ts, PetscErrorCode (*func)(TS))
3098: {
3100:   ts->prestep = func;
3101:   return 0;
3102: }

3104: /*@
3105:   TSPreStep - Runs the user-defined pre-step function.

3107:   Collective on TS

3109:   Input Parameters:
3110: . ts   - The TS context obtained from TSCreate()

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

3116:   Level: developer

3118: .seealso: TSSetPreStep(), TSPreStage(), TSPostStage(), TSPostStep()
3119: @*/
3120: PetscErrorCode  TSPreStep(TS ts)
3121: {
3123:   if (ts->prestep) {
3124:     Vec              U;
3125:     PetscObjectId    idprev;
3126:     PetscBool        sameObject;
3127:     PetscObjectState sprev,spost;

3129:     TSGetSolution(ts,&U);
3130:     PetscObjectGetId((PetscObject)U,&idprev);
3131:     PetscObjectStateGet((PetscObject)U,&sprev);
3132:     PetscStackCallStandard((*ts->prestep),ts);
3133:     TSGetSolution(ts,&U);
3134:     PetscObjectCompareId((PetscObject)U,idprev,&sameObject);
3135:     PetscObjectStateGet((PetscObject)U,&spost);
3136:     if (!sameObject || sprev != spost) TSRestartStep(ts);
3137:   }
3138:   return 0;
3139: }

3141: /*@C
3142:   TSSetPreStage - Sets the general-purpose function
3143:   called once at the beginning of each stage.

3145:   Logically Collective on TS

3147:   Input Parameters:
3148: + ts   - The TS context obtained from TSCreate()
3149: - func - The function

3151:   Calling sequence of func:
3152: .vb
3153:   PetscErrorCode func(TS ts, PetscReal stagetime);
3154: .ve

3156:   Level: intermediate

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

3163: .seealso: TSSetPostStage(), TSSetPreStep(), TSSetPostStep(), TSGetApplicationContext()
3164: @*/
3165: PetscErrorCode  TSSetPreStage(TS ts, PetscErrorCode (*func)(TS,PetscReal))
3166: {
3168:   ts->prestage = func;
3169:   return 0;
3170: }

3172: /*@C
3173:   TSSetPostStage - Sets the general-purpose function
3174:   called once at the end of each stage.

3176:   Logically Collective on TS

3178:   Input Parameters:
3179: + ts   - The TS context obtained from TSCreate()
3180: - func - The function

3182:   Calling sequence of func:
3183: .vb
3184:   PetscErrorCode func(TS ts, PetscReal stagetime, PetscInt stageindex, Vec* Y);
3185: .ve

3187:   Level: intermediate

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

3194: .seealso: TSSetPreStage(), TSSetPreStep(), TSSetPostStep(), TSGetApplicationContext()
3195: @*/
3196: PetscErrorCode  TSSetPostStage(TS ts, PetscErrorCode (*func)(TS,PetscReal,PetscInt,Vec*))
3197: {
3199:   ts->poststage = func;
3200:   return 0;
3201: }

3203: /*@C
3204:   TSSetPostEvaluate - Sets the general-purpose function
3205:   called once at the end of each step evaluation.

3207:   Logically Collective on TS

3209:   Input Parameters:
3210: + ts   - The TS context obtained from TSCreate()
3211: - func - The function

3213:   Calling sequence of func:
3214: .vb
3215:   PetscErrorCode func(TS ts);
3216: .ve

3218:   Level: intermediate

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

3227: .seealso: TSSetPreStage(), TSSetPreStep(), TSSetPostStep(), TSGetApplicationContext()
3228: @*/
3229: PetscErrorCode  TSSetPostEvaluate(TS ts, PetscErrorCode (*func)(TS))
3230: {
3232:   ts->postevaluate = func;
3233:   return 0;
3234: }

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

3239:   Collective on TS

3241:   Input Parameters:
3242: . ts          - The TS context obtained from TSCreate()
3243:   stagetime   - The absolute time of the current stage

3245:   Notes:
3246:   TSPreStage() is typically used within time stepping implementations,
3247:   most users would not generally call this routine themselves.

3249:   Level: developer

3251: .seealso: TSPostStage(), TSSetPreStep(), TSPreStep(), TSPostStep()
3252: @*/
3253: PetscErrorCode  TSPreStage(TS ts, PetscReal stagetime)
3254: {
3256:   if (ts->prestage) {
3257:     PetscStackCallStandard((*ts->prestage),ts,stagetime);
3258:   }
3259:   return 0;
3260: }

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

3265:   Collective on TS

3267:   Input Parameters:
3268: . ts          - The TS context obtained from TSCreate()
3269:   stagetime   - The absolute time of the current stage
3270:   stageindex  - Stage number
3271:   Y           - Array of vectors (of size = total number
3272:                 of stages) with the stage solutions

3274:   Notes:
3275:   TSPostStage() is typically used within time stepping implementations,
3276:   most users would not generally call this routine themselves.

3278:   Level: developer

3280: .seealso: TSPreStage(), TSSetPreStep(), TSPreStep(), TSPostStep()
3281: @*/
3282: PetscErrorCode  TSPostStage(TS ts, PetscReal stagetime, PetscInt stageindex, Vec *Y)
3283: {
3285:   if (ts->poststage) {
3286:     PetscStackCallStandard((*ts->poststage),ts,stagetime,stageindex,Y);
3287:   }
3288:   return 0;
3289: }

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

3294:   Collective on TS

3296:   Input Parameters:
3297: . ts          - The TS context obtained from TSCreate()

3299:   Notes:
3300:   TSPostEvaluate() is typically used within time stepping implementations,
3301:   most users would not generally call this routine themselves.

3303:   Level: developer

3305: .seealso: TSSetPostEvaluate(), TSSetPreStep(), TSPreStep(), TSPostStep()
3306: @*/
3307: PetscErrorCode  TSPostEvaluate(TS ts)
3308: {
3310:   if (ts->postevaluate) {
3311:     Vec              U;
3312:     PetscObjectState sprev,spost;

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

3323: /*@C
3324:   TSSetPostStep - Sets the general-purpose function
3325:   called once at the end of each time step.

3327:   Logically Collective on TS

3329:   Input Parameters:
3330: + ts   - The TS context obtained from TSCreate()
3331: - func - The function

3333:   Calling sequence of func:
3334: $ func (TS ts);

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

3341:   Level: intermediate

3343: .seealso: TSSetPreStep(), TSSetPreStage(), TSSetPostEvaluate(), TSGetTimeStep(), TSGetStepNumber(), TSGetTime(), TSRestartStep()
3344: @*/
3345: PetscErrorCode  TSSetPostStep(TS ts, PetscErrorCode (*func)(TS))
3346: {
3348:   ts->poststep = func;
3349:   return 0;
3350: }

3352: /*@
3353:   TSPostStep - Runs the user-defined post-step function.

3355:   Collective on TS

3357:   Input Parameters:
3358: . ts   - The TS context obtained from TSCreate()

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

3364:   Level: developer

3366: @*/
3367: PetscErrorCode  TSPostStep(TS ts)
3368: {
3370:   if (ts->poststep) {
3371:     Vec              U;
3372:     PetscObjectId    idprev;
3373:     PetscBool        sameObject;
3374:     PetscObjectState sprev,spost;

3376:     TSGetSolution(ts,&U);
3377:     PetscObjectGetId((PetscObject)U,&idprev);
3378:     PetscObjectStateGet((PetscObject)U,&sprev);
3379:     PetscStackCallStandard((*ts->poststep),ts);
3380:     TSGetSolution(ts,&U);
3381:     PetscObjectCompareId((PetscObject)U,idprev,&sameObject);
3382:     PetscObjectStateGet((PetscObject)U,&spost);
3383:     if (!sameObject || sprev != spost) TSRestartStep(ts);
3384:   }
3385:   return 0;
3386: }

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

3391:    Collective on TS

3393:    Input Parameters:
3394: +  ts - time stepping context
3395: -  t - time to interpolate to

3397:    Output Parameter:
3398: .  U - state at given time

3400:    Level: intermediate

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

3405: .seealso: TSSetExactFinalTime(), TSSolve()
3406: @*/
3407: PetscErrorCode TSInterpolate(TS ts,PetscReal t,Vec U)
3408: {
3413:   (*ts->ops->interpolate)(ts,t,U);
3414:   return 0;
3415: }

3417: /*@
3418:    TSStep - Steps one time step

3420:    Collective on TS

3422:    Input Parameter:
3423: .  ts - the TS context obtained from TSCreate()

3425:    Level: developer

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

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

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

3436: .seealso: TSCreate(), TSSetUp(), TSDestroy(), TSSolve(), TSSetPreStep(), TSSetPreStage(), TSSetPostStage(), TSInterpolate()
3437: @*/
3438: PetscErrorCode  TSStep(TS ts)
3439: {
3440:   PetscErrorCode   ierr;
3441:   static PetscBool cite = PETSC_FALSE;
3442:   PetscReal        ptime;

3445:   PetscCitationsRegister("@article{tspaper,\n"
3446:                                 "  title         = {{PETSc/TS}: A Modern Scalable {DAE/ODE} Solver Library},\n"
3447:                                 "  author        = {Abhyankar, Shrirang and Brown, Jed and Constantinescu, Emil and Ghosh, Debojyoti and Smith, Barry F. and Zhang, Hong},\n"
3448:                                 "  journal       = {arXiv e-preprints},\n"
3449:                                 "  eprint        = {1806.01437},\n"
3450:                                 "  archivePrefix = {arXiv},\n"
3451:                                 "  year          = {2018}\n}\n",&cite);

3453:   TSSetUp(ts);
3454:   TSTrajectorySetUp(ts->trajectory,ts);


3461:   if (!ts->steps) ts->ptime_prev = ts->ptime;
3462:   ptime = ts->ptime; ts->ptime_prev_rollback = ts->ptime_prev;
3463:   ts->reason = TS_CONVERGED_ITERATING;

3465:   PetscLogEventBegin(TS_Step,ts,0,0,0);
3466:   (*ts->ops->step)(ts);
3467:   PetscLogEventEnd(TS_Step,ts,0,0,0);

3469:   if (ts->reason >= 0) {
3470:     ts->ptime_prev = ptime;
3471:     ts->steps++;
3472:     ts->steprollback = PETSC_FALSE;
3473:     ts->steprestart  = PETSC_FALSE;
3474:   }

3476:   if (!ts->reason) {
3477:     if (ts->steps >= ts->max_steps) ts->reason = TS_CONVERGED_ITS;
3478:     else if (ts->ptime >= ts->max_time) ts->reason = TS_CONVERGED_TIME;
3479:   }

3483:   return 0;
3484: }

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

3490:    Collective on TS

3492:    Input Parameters:
3493: +  ts - time stepping context
3494: -  wnormtype - norm type, either NORM_2 or NORM_INFINITY

3496:    Input/Output Parameter:
3497: .  order - optional, desired order for the error evaluation or PETSC_DECIDE;
3498:            on output, the actual order of the error evaluation

3500:    Output Parameter:
3501: .  wlte - the weighted local truncation error norm

3503:    Level: advanced

3505:    Notes:
3506:    If the timestepper cannot evaluate the error in a particular step
3507:    (eg. in the first step or restart steps after event handling),
3508:    this routine returns wlte=-1.0 .

3510: .seealso: TSStep(), TSAdapt, TSErrorWeightedNorm()
3511: @*/
3512: PetscErrorCode TSEvaluateWLTE(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte)
3513: {
3522:   (*ts->ops->evaluatewlte)(ts,wnormtype,order,wlte);
3523:   return 0;
3524: }

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

3529:    Collective on TS

3531:    Input Parameters:
3532: +  ts - time stepping context
3533: .  order - desired order of accuracy
3534: -  done - whether the step was evaluated at this order (pass NULL to generate an error if not available)

3536:    Output Parameter:
3537: .  U - state at the end of the current step

3539:    Level: advanced

3541:    Notes:
3542:    This function cannot be called until all stages have been evaluated.
3543:    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.

3545: .seealso: TSStep(), TSAdapt
3546: @*/
3547: PetscErrorCode TSEvaluateStep(TS ts,PetscInt order,Vec U,PetscBool *done)
3548: {
3553:   (*ts->ops->evaluatestep)(ts,order,U,done);
3554:   return 0;
3555: }

3557: /*@C
3558:   TSGetComputeInitialCondition - Get the function used to automatically compute an initial condition for the timestepping.

3560:   Not collective

3562:   Input Parameter:
3563: . ts        - time stepping context

3565:   Output Parameter:
3566: . initConditions - The function which computes an initial condition

3568:    Level: advanced

3570:    Notes:
3571:    The calling sequence for the function is
3572: $ initCondition(TS ts, Vec u)
3573: $ ts - The timestepping context
3574: $ u  - The input vector in which the initial condition is stored

3576: .seealso: TSSetComputeInitialCondition(), TSComputeInitialCondition()
3577: @*/
3578: PetscErrorCode TSGetComputeInitialCondition(TS ts, PetscErrorCode (**initCondition)(TS, Vec))
3579: {
3582:   *initCondition = ts->ops->initcondition;
3583:   return 0;
3584: }

3586: /*@C
3587:   TSSetComputeInitialCondition - Set the function used to automatically compute an initial condition for the timestepping.

3589:   Logically collective on ts

3591:   Input Parameters:
3592: + ts        - time stepping context
3593: - initCondition - The function which computes an initial condition

3595:   Level: advanced

3597:   Calling sequence for initCondition:
3598: $ PetscErrorCode initCondition(TS ts, Vec u)

3600: + ts - The timestepping context
3601: - u  - The input vector in which the initial condition is to be stored

3603: .seealso: TSGetComputeInitialCondition(), TSComputeInitialCondition()
3604: @*/
3605: PetscErrorCode TSSetComputeInitialCondition(TS ts, PetscErrorCode (*initCondition)(TS, Vec))
3606: {
3609:   ts->ops->initcondition = initCondition;
3610:   return 0;
3611: }

3613: /*@
3614:   TSComputeInitialCondition - Compute an initial condition for the timestepping using the function previously set.

3616:   Collective on ts

3618:   Input Parameters:
3619: + ts - time stepping context
3620: - u  - The Vec to store the condition in which will be used in TSSolve()

3622:   Level: advanced

3624: .seealso: TSGetComputeInitialCondition(), TSSetComputeInitialCondition(), TSSolve()
3625: @*/
3626: PetscErrorCode TSComputeInitialCondition(TS ts, Vec u)
3627: {
3630:   if (ts->ops->initcondition) (*ts->ops->initcondition)(ts, u);
3631:   return 0;
3632: }

3634: /*@C
3635:   TSGetComputeExactError - Get the function used to automatically compute the exact error for the timestepping.

3637:   Not collective

3639:   Input Parameter:
3640: . ts         - time stepping context

3642:   Output Parameter:
3643: . exactError - The function which computes the solution error

3645:   Level: advanced

3647:   Calling sequence for exactError:
3648: $ PetscErrorCode exactError(TS ts, Vec u)

3650: + ts - The timestepping context
3651: . u  - The approximate solution vector
3652: - e  - The input vector in which the error is stored

3654: .seealso: TSGetComputeExactError(), TSComputeExactError()
3655: @*/
3656: PetscErrorCode TSGetComputeExactError(TS ts, PetscErrorCode (**exactError)(TS, Vec, Vec))
3657: {
3660:   *exactError = ts->ops->exacterror;
3661:   return 0;
3662: }

3664: /*@C
3665:   TSSetComputeExactError - Set the function used to automatically compute the exact error for the timestepping.

3667:   Logically collective on ts

3669:   Input Parameters:
3670: + ts         - time stepping context
3671: - exactError - The function which computes the solution error

3673:   Level: advanced

3675:   Calling sequence for exactError:
3676: $ PetscErrorCode exactError(TS ts, Vec u)

3678: + ts - The timestepping context
3679: . u  - The approximate solution vector
3680: - e  - The input vector in which the error is stored

3682: .seealso: TSGetComputeExactError(), TSComputeExactError()
3683: @*/
3684: PetscErrorCode TSSetComputeExactError(TS ts, PetscErrorCode (*exactError)(TS, Vec, Vec))
3685: {
3688:   ts->ops->exacterror = exactError;
3689:   return 0;
3690: }

3692: /*@
3693:   TSComputeExactError - Compute the solution error for the timestepping using the function previously set.

3695:   Collective on ts

3697:   Input Parameters:
3698: + ts - time stepping context
3699: . u  - The approximate solution
3700: - e  - The Vec used to store the error

3702:   Level: advanced

3704: .seealso: TSGetComputeInitialCondition(), TSSetComputeInitialCondition(), TSSolve()
3705: @*/
3706: PetscErrorCode TSComputeExactError(TS ts, Vec u, Vec e)
3707: {
3711:   if (ts->ops->exacterror) (*ts->ops->exacterror)(ts, u, e);
3712:   return 0;
3713: }

3715: /*@
3716:    TSSolve - Steps the requested number of timesteps.

3718:    Collective on TS

3720:    Input Parameters:
3721: +  ts - the TS context obtained from TSCreate()
3722: -  u - the solution vector  (can be null if TSSetSolution() was used and TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP) was not used,
3723:                              otherwise must contain the initial conditions and will contain the solution at the final requested time

3725:    Level: beginner

3727:    Notes:
3728:    The final time returned by this function may be different from the time of the internally
3729:    held state accessible by TSGetSolution() and TSGetTime() because the method may have
3730:    stepped over the final time.

3732: .seealso: TSCreate(), TSSetSolution(), TSStep(), TSGetTime(), TSGetSolveTime()
3733: @*/
3734: PetscErrorCode TSSolve(TS ts,Vec u)
3735: {
3736:   Vec               solution;


3741:   TSSetExactFinalTimeDefault(ts);
3742:   if (ts->exact_final_time == TS_EXACTFINALTIME_INTERPOLATE && u) {   /* Need ts->vec_sol to be distinct so it is not overwritten when we interpolate at the end */
3743:     if (!ts->vec_sol || u == ts->vec_sol) {
3744:       VecDuplicate(u,&solution);
3745:       TSSetSolution(ts,solution);
3746:       VecDestroy(&solution); /* grant ownership */
3747:     }
3748:     VecCopy(u,ts->vec_sol);
3750:   } else if (u) {
3751:     TSSetSolution(ts,u);
3752:   }
3753:   TSSetUp(ts);
3754:   TSTrajectorySetUp(ts->trajectory,ts);


3760:   if (ts->forward_solve) {
3761:     TSForwardSetUp(ts);
3762:   }

3764:   /* reset number of steps only when the step is not restarted. ARKIMEX
3765:      restarts the step after an event. Resetting these counters in such case causes
3766:      TSTrajectory to incorrectly save the output files
3767:   */
3768:   /* reset time step and iteration counters */
3769:   if (!ts->steps) {
3770:     ts->ksp_its           = 0;
3771:     ts->snes_its          = 0;
3772:     ts->num_snes_failures = 0;
3773:     ts->reject            = 0;
3774:     ts->steprestart       = PETSC_TRUE;
3775:     ts->steprollback      = PETSC_FALSE;
3776:     ts->rhsjacobian.time  = PETSC_MIN_REAL;
3777:   }

3779:   /* make sure initial time step does not overshoot final time */
3780:   if (ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) {
3781:     PetscReal maxdt = ts->max_time-ts->ptime;
3782:     PetscReal dt = ts->time_step;

3784:     ts->time_step = dt >= maxdt ? maxdt : (PetscIsCloseAtTol(dt,maxdt,10*PETSC_MACHINE_EPSILON,0) ? maxdt : dt);
3785:   }
3786:   ts->reason = TS_CONVERGED_ITERATING;

3788:   {
3789:     PetscViewer       viewer;
3790:     PetscViewerFormat format;
3791:     PetscBool         flg;
3792:     static PetscBool  incall = PETSC_FALSE;

3794:     if (!incall) {
3795:       /* Estimate the convergence rate of the time discretization */
3796:       PetscOptionsGetViewer(PetscObjectComm((PetscObject) ts),((PetscObject)ts)->options, ((PetscObject) ts)->prefix, "-ts_convergence_estimate", &viewer, &format, &flg);
3797:       if (flg) {
3798:         PetscConvEst conv;
3799:         DM           dm;
3800:         PetscReal   *alpha; /* Convergence rate of the solution error for each field in the L_2 norm */
3801:         PetscInt     Nf;
3802:         PetscBool    checkTemporal = PETSC_TRUE;

3804:         incall = PETSC_TRUE;
3805:         PetscOptionsGetBool(((PetscObject)ts)->options, ((PetscObject) ts)->prefix, "-ts_convergence_temporal", &checkTemporal, &flg);
3806:         TSGetDM(ts, &dm);
3807:         DMGetNumFields(dm, &Nf);
3808:         PetscCalloc1(PetscMax(Nf, 1), &alpha);
3809:         PetscConvEstCreate(PetscObjectComm((PetscObject) ts), &conv);
3810:         PetscConvEstUseTS(conv, checkTemporal);
3811:         PetscConvEstSetSolver(conv, (PetscObject) ts);
3812:         PetscConvEstSetFromOptions(conv);
3813:         PetscConvEstSetUp(conv);
3814:         PetscConvEstGetConvRate(conv, alpha);
3815:         PetscViewerPushFormat(viewer, format);
3816:         PetscConvEstRateView(conv, alpha, viewer);
3817:         PetscViewerPopFormat(viewer);
3818:         PetscViewerDestroy(&viewer);
3819:         PetscConvEstDestroy(&conv);
3820:         PetscFree(alpha);
3821:         incall = PETSC_FALSE;
3822:       }
3823:     }
3824:   }

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

3828:   if (ts->ops->solve) { /* This private interface is transitional and should be removed when all implementations are updated. */
3829:     (*ts->ops->solve)(ts);
3830:     if (u) VecCopy(ts->vec_sol,u);
3831:     ts->solvetime = ts->ptime;
3832:     solution = ts->vec_sol;
3833:   } else { /* Step the requested number of timesteps. */
3834:     if (ts->steps >= ts->max_steps) ts->reason = TS_CONVERGED_ITS;
3835:     else if (ts->ptime >= ts->max_time) ts->reason = TS_CONVERGED_TIME;

3837:     if (!ts->steps) {
3838:       TSTrajectorySet(ts->trajectory,ts,ts->steps,ts->ptime,ts->vec_sol);
3839:       TSEventInitialize(ts->event,ts,ts->ptime,ts->vec_sol);
3840:     }

3842:     while (!ts->reason) {
3843:       TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);
3844:       if (!ts->steprollback) {
3845:         TSPreStep(ts);
3846:       }
3847:       TSStep(ts);
3848:       if (ts->testjacobian) {
3849:         TSRHSJacobianTest(ts,NULL);
3850:       }
3851:       if (ts->testjacobiantranspose) {
3852:         TSRHSJacobianTestTranspose(ts,NULL);
3853:       }
3854:       if (ts->quadraturets && ts->costintegralfwd) { /* Must evaluate the cost integral before event is handled. The cost integral value can also be rolled back. */
3855:         if (ts->reason >= 0) ts->steps--; /* Revert the step number changed by TSStep() */
3856:         TSForwardCostIntegral(ts);
3857:         if (ts->reason >= 0) ts->steps++;
3858:       }
3859:       if (ts->forward_solve) { /* compute forward sensitivities before event handling because postevent() may change RHS and jump conditions may have to be applied */
3860:         if (ts->reason >= 0) ts->steps--; /* Revert the step number changed by TSStep() */
3861:         TSForwardStep(ts);
3862:         if (ts->reason >= 0) ts->steps++;
3863:       }
3864:       TSPostEvaluate(ts);
3865:       TSEventHandler(ts); /* The right-hand side may be changed due to event. Be careful with Any computation using the RHS information after this point. */
3866:       if (ts->steprollback) {
3867:         TSPostEvaluate(ts);
3868:       }
3869:       if (!ts->steprollback) {
3870:         TSTrajectorySet(ts->trajectory,ts,ts->steps,ts->ptime,ts->vec_sol);
3871:         TSPostStep(ts);
3872:       }
3873:     }
3874:     TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);

3876:     if (ts->exact_final_time == TS_EXACTFINALTIME_INTERPOLATE && ts->ptime > ts->max_time) {
3877:       TSInterpolate(ts,ts->max_time,u);
3878:       ts->solvetime = ts->max_time;
3879:       solution = u;
3880:       TSMonitor(ts,-1,ts->solvetime,solution);
3881:     } else {
3882:       if (u) VecCopy(ts->vec_sol,u);
3883:       ts->solvetime = ts->ptime;
3884:       solution = ts->vec_sol;
3885:     }
3886:   }

3888:   TSViewFromOptions(ts,NULL,"-ts_view");
3889:   VecViewFromOptions(solution,(PetscObject)ts,"-ts_view_solution");
3890:   PetscObjectSAWsBlock((PetscObject)ts);
3891:   if (ts->adjoint_solve) {
3892:     TSAdjointSolve(ts);
3893:   }
3894:   return 0;
3895: }

3897: /*@
3898:    TSGetTime - Gets the time of the most recently completed step.

3900:    Not Collective

3902:    Input Parameter:
3903: .  ts - the TS context obtained from TSCreate()

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

3908:    Level: beginner

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

3914: .seealso:  TSGetSolveTime(), TSSetTime(), TSGetTimeStep(), TSGetStepNumber()

3916: @*/
3917: PetscErrorCode  TSGetTime(TS ts,PetscReal *t)
3918: {
3921:   *t = ts->ptime;
3922:   return 0;
3923: }

3925: /*@
3926:    TSGetPrevTime - Gets the starting time of the previously completed step.

3928:    Not Collective

3930:    Input Parameter:
3931: .  ts - the TS context obtained from TSCreate()

3933:    Output Parameter:
3934: .  t  - the previous time

3936:    Level: beginner

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

3940: @*/
3941: PetscErrorCode  TSGetPrevTime(TS ts,PetscReal *t)
3942: {
3945:   *t = ts->ptime_prev;
3946:   return 0;
3947: }

3949: /*@
3950:    TSSetTime - Allows one to reset the time.

3952:    Logically Collective on TS

3954:    Input Parameters:
3955: +  ts - the TS context obtained from TSCreate()
3956: -  time - the time

3958:    Level: intermediate

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

3962: @*/
3963: PetscErrorCode  TSSetTime(TS ts, PetscReal t)
3964: {
3967:   ts->ptime = t;
3968:   return 0;
3969: }

3971: /*@C
3972:    TSSetOptionsPrefix - Sets the prefix used for searching for all
3973:    TS options in the database.

3975:    Logically Collective on TS

3977:    Input Parameters:
3978: +  ts     - The TS context
3979: -  prefix - The prefix to prepend to all option names

3981:    Notes:
3982:    A hyphen (-) must NOT be given at the beginning of the prefix name.
3983:    The first character of all runtime options is AUTOMATICALLY the
3984:    hyphen.

3986:    Level: advanced

3988: .seealso: TSSetFromOptions()

3990: @*/
3991: PetscErrorCode  TSSetOptionsPrefix(TS ts,const char prefix[])
3992: {
3993:   SNES           snes;

3996:   PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);
3997:   TSGetSNES(ts,&snes);
3998:   SNESSetOptionsPrefix(snes,prefix);
3999:   return 0;
4000: }

4002: /*@C
4003:    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
4004:    TS options in the database.

4006:    Logically Collective on TS

4008:    Input Parameters:
4009: +  ts     - The TS context
4010: -  prefix - The prefix to prepend to all option names

4012:    Notes:
4013:    A hyphen (-) must NOT be given at the beginning of the prefix name.
4014:    The first character of all runtime options is AUTOMATICALLY the
4015:    hyphen.

4017:    Level: advanced

4019: .seealso: TSGetOptionsPrefix()

4021: @*/
4022: PetscErrorCode  TSAppendOptionsPrefix(TS ts,const char prefix[])
4023: {
4024:   SNES           snes;

4027:   PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);
4028:   TSGetSNES(ts,&snes);
4029:   SNESAppendOptionsPrefix(snes,prefix);
4030:   return 0;
4031: }

4033: /*@C
4034:    TSGetOptionsPrefix - Sets the prefix used for searching for all
4035:    TS options in the database.

4037:    Not Collective

4039:    Input Parameter:
4040: .  ts - The TS context

4042:    Output Parameter:
4043: .  prefix - A pointer to the prefix string used

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

4049:    Level: intermediate

4051: .seealso: TSAppendOptionsPrefix()
4052: @*/
4053: PetscErrorCode  TSGetOptionsPrefix(TS ts,const char *prefix[])
4054: {
4057:   PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);
4058:   return 0;
4059: }

4061: /*@C
4062:    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.

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

4066:    Input Parameter:
4067: .  ts  - The TS context obtained from TSCreate()

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

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

4078:    Level: intermediate

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

4082: @*/
4083: PetscErrorCode  TSGetRHSJacobian(TS ts,Mat *Amat,Mat *Pmat,TSRHSJacobian *func,void **ctx)
4084: {
4085:   DM             dm;

4087:   if (Amat || Pmat) {
4088:     SNES snes;
4089:     TSGetSNES(ts,&snes);
4090:     SNESSetUpMatrices(snes);
4091:     SNESGetJacobian(snes,Amat,Pmat,NULL,NULL);
4092:   }
4093:   TSGetDM(ts,&dm);
4094:   DMTSGetRHSJacobian(dm,func,ctx);
4095:   return 0;
4096: }

4098: /*@C
4099:    TSGetIJacobian - Returns the implicit Jacobian at the present timestep.

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

4103:    Input Parameter:
4104: .  ts  - The TS context obtained from TSCreate()

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

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

4115:    Level: advanced

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

4119: @*/
4120: PetscErrorCode  TSGetIJacobian(TS ts,Mat *Amat,Mat *Pmat,TSIJacobian *f,void **ctx)
4121: {
4122:   DM             dm;

4124:   if (Amat || Pmat) {
4125:     SNES snes;
4126:     TSGetSNES(ts,&snes);
4127:     SNESSetUpMatrices(snes);
4128:     SNESGetJacobian(snes,Amat,Pmat,NULL,NULL);
4129:   }
4130:   TSGetDM(ts,&dm);
4131:   DMTSGetIJacobian(dm,f,ctx);
4132:   return 0;
4133: }

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

4139:    Logically Collective on ts

4141:    Input Parameters:
4142: +  ts - the ODE integrator object
4143: -  dm - the dm, cannot be NULL

4145:    Notes:
4146:    A DM can only be used for solving one problem at a time because information about the problem is stored on the DM,
4147:    even when not using interfaces like DMTSSetIFunction().  Use DMClone() to get a distinct DM when solving
4148:    different problems using the same function space.

4150:    Level: intermediate

4152: .seealso: TSGetDM(), SNESSetDM(), SNESGetDM()
4153: @*/
4154: PetscErrorCode  TSSetDM(TS ts,DM dm)
4155: {
4156:   SNES           snes;
4157:   DMTS           tsdm;

4161:   PetscObjectReference((PetscObject)dm);
4162:   if (ts->dm) {               /* Move the DMTS context over to the new DM unless the new DM already has one */
4163:     if (ts->dm->dmts && !dm->dmts) {
4164:       DMCopyDMTS(ts->dm,dm);
4165:       DMGetDMTS(ts->dm,&tsdm);
4166:       if (tsdm->originaldm == ts->dm) { /* Grant write privileges to the replacement DM */
4167:         tsdm->originaldm = dm;
4168:       }
4169:     }
4170:     DMDestroy(&ts->dm);
4171:   }
4172:   ts->dm = dm;

4174:   TSGetSNES(ts,&snes);
4175:   SNESSetDM(snes,dm);
4176:   return 0;
4177: }

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

4182:    Not Collective

4184:    Input Parameter:
4185: . ts - the preconditioner context

4187:    Output Parameter:
4188: .  dm - the dm

4190:    Level: intermediate

4192: .seealso: TSSetDM(), SNESSetDM(), SNESGetDM()
4193: @*/
4194: PetscErrorCode  TSGetDM(TS ts,DM *dm)
4195: {
4197:   if (!ts->dm) {
4198:     DMShellCreate(PetscObjectComm((PetscObject)ts),&ts->dm);
4199:     if (ts->snes) SNESSetDM(ts->snes,ts->dm);
4200:   }
4201:   *dm = ts->dm;
4202:   return 0;
4203: }

4205: /*@
4206:    SNESTSFormFunction - Function to evaluate nonlinear residual

4208:    Logically Collective on SNES

4210:    Input Parameters:
4211: + snes - nonlinear solver
4212: . U - the current state at which to evaluate the residual
4213: - ctx - user context, must be a TS

4215:    Output Parameter:
4216: . F - the nonlinear residual

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

4222:    Level: advanced

4224: .seealso: SNESSetFunction(), MatFDColoringSetFunction()
4225: @*/
4226: PetscErrorCode  SNESTSFormFunction(SNES snes,Vec U,Vec F,void *ctx)
4227: {
4228:   TS             ts = (TS)ctx;

4234:   (ts->ops->snesfunction)(snes,U,F,ts);
4235:   return 0;
4236: }

4238: /*@
4239:    SNESTSFormJacobian - Function to evaluate the Jacobian

4241:    Collective on SNES

4243:    Input Parameters:
4244: + snes - nonlinear solver
4245: . U - the current state at which to evaluate the residual
4246: - ctx - user context, must be a TS

4248:    Output Parameters:
4249: + A - the Jacobian
4250: - B - the preconditioning matrix (may be the same as A)

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

4255:    Level: developer

4257: .seealso: SNESSetJacobian()
4258: @*/
4259: PetscErrorCode  SNESTSFormJacobian(SNES snes,Vec U,Mat A,Mat B,void *ctx)
4260: {
4261:   TS             ts = (TS)ctx;

4270:   (ts->ops->snesjacobian)(snes,U,A,B,ts);
4271:   return 0;
4272: }

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

4277:    Collective on TS

4279:    Input Parameters:
4280: +  ts - time stepping context
4281: .  t - time at which to evaluate
4282: .  U - state at which to evaluate
4283: -  ctx - context

4285:    Output Parameter:
4286: .  F - right hand side

4288:    Level: intermediate

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

4294: .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSJacobianConstant()
4295: @*/
4296: PetscErrorCode TSComputeRHSFunctionLinear(TS ts,PetscReal t,Vec U,Vec F,void *ctx)
4297: {
4298:   Mat            Arhs,Brhs;

4300:   TSGetRHSMats_Private(ts,&Arhs,&Brhs);
4301:   /* undo the damage caused by shifting */
4302:   TSRecoverRHSJacobian(ts,Arhs,Brhs);
4303:   TSComputeRHSJacobian(ts,t,U,Arhs,Brhs);
4304:   MatMult(Arhs,U,F);
4305:   return 0;
4306: }

4308: /*@C
4309:    TSComputeRHSJacobianConstant - Reuses a Jacobian that is time-independent.

4311:    Collective on TS

4313:    Input Parameters:
4314: +  ts - time stepping context
4315: .  t - time at which to evaluate
4316: .  U - state at which to evaluate
4317: -  ctx - context

4319:    Output Parameters:
4320: +  A - pointer to operator
4321: -  B - pointer to preconditioning matrix

4323:    Level: intermediate

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

4328: .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSFunctionLinear()
4329: @*/
4330: PetscErrorCode TSComputeRHSJacobianConstant(TS ts,PetscReal t,Vec U,Mat A,Mat B,void *ctx)
4331: {
4332:   return 0;
4333: }

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

4338:    Collective on TS

4340:    Input Parameters:
4341: +  ts - time stepping context
4342: .  t - time at which to evaluate
4343: .  U - state at which to evaluate
4344: .  Udot - time derivative of state vector
4345: -  ctx - context

4347:    Output Parameter:
4348: .  F - left hand side

4350:    Level: intermediate

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

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

4360: .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIJacobianConstant(), TSComputeRHSFunctionLinear()
4361: @*/
4362: PetscErrorCode TSComputeIFunctionLinear(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void *ctx)
4363: {
4364:   Mat            A,B;

4366:   TSGetIJacobian(ts,&A,&B,NULL,NULL);
4367:   TSComputeIJacobian(ts,t,U,Udot,1.0,A,B,PETSC_TRUE);
4368:   MatMult(A,Udot,F);
4369:   return 0;
4370: }

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

4375:    Collective on TS

4377:    Input Parameters:
4378: +  ts - time stepping context
4379: .  t - time at which to evaluate
4380: .  U - state at which to evaluate
4381: .  Udot - time derivative of state vector
4382: .  shift - shift to apply
4383: -  ctx - context

4385:    Output Parameters:
4386: +  A - pointer to operator
4387: -  B - pointer to preconditioning matrix

4389:    Level: advanced

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

4394:    It is only appropriate for problems of the form

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

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

4402: $    shift*M + J

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

4407: .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIFunctionLinear()
4408: @*/
4409: PetscErrorCode TSComputeIJacobianConstant(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal shift,Mat A,Mat B,void *ctx)
4410: {
4411:   MatScale(A, shift / ts->ijacobian.shift);
4412:   ts->ijacobian.shift = shift;
4413:   return 0;
4414: }

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

4419:    Not Collective

4421:    Input Parameter:
4422: .  ts - the TS context

4424:    Output Parameter:
4425: .  equation_type - see TSEquationType

4427:    Level: beginner

4429: .seealso: TSSetEquationType(), TSEquationType
4430: @*/
4431: PetscErrorCode  TSGetEquationType(TS ts,TSEquationType *equation_type)
4432: {
4435:   *equation_type = ts->equation_type;
4436:   return 0;
4437: }

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

4442:    Not Collective

4444:    Input Parameters:
4445: +  ts - the TS context
4446: -  equation_type - see TSEquationType

4448:    Level: advanced

4450: .seealso: TSGetEquationType(), TSEquationType
4451: @*/
4452: PetscErrorCode  TSSetEquationType(TS ts,TSEquationType equation_type)
4453: {
4455:   ts->equation_type = equation_type;
4456:   return 0;
4457: }

4459: /*@
4460:    TSGetConvergedReason - Gets the reason the TS iteration was stopped.

4462:    Not Collective

4464:    Input Parameter:
4465: .  ts - the TS context

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

4471:    Level: beginner

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

4476: .seealso: TSSetConvergenceTest(), TSConvergedReason
4477: @*/
4478: PetscErrorCode  TSGetConvergedReason(TS ts,TSConvergedReason *reason)
4479: {
4482:   *reason = ts->reason;
4483:   return 0;
4484: }

4486: /*@
4487:    TSSetConvergedReason - Sets the reason for handling the convergence of TSSolve.

4489:    Logically Collective; reason must contain common value

4491:    Input Parameters:
4492: +  ts - the TS context
4493: -  reason - negative value indicates diverged, positive value converged, see TSConvergedReason or the
4494:             manual pages for the individual convergence tests for complete lists

4496:    Level: advanced

4498:    Notes:
4499:    Can only be called while TSSolve() is active.

4501: .seealso: TSConvergedReason
4502: @*/
4503: PetscErrorCode  TSSetConvergedReason(TS ts,TSConvergedReason reason)
4504: {
4506:   ts->reason = reason;
4507:   return 0;
4508: }

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

4513:    Not Collective

4515:    Input Parameter:
4516: .  ts - the TS context

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

4521:    Level: beginner

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

4526: .seealso: TSSetConvergenceTest(), TSConvergedReason
4527: @*/
4528: PetscErrorCode  TSGetSolveTime(TS ts,PetscReal *ftime)
4529: {
4532:   *ftime = ts->solvetime;
4533:   return 0;
4534: }

4536: /*@
4537:    TSGetSNESIterations - Gets the total number of nonlinear iterations
4538:    used by the time integrator.

4540:    Not Collective

4542:    Input Parameter:
4543: .  ts - TS context

4545:    Output Parameter:
4546: .  nits - number of nonlinear iterations

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

4551:    Level: intermediate

4553: .seealso:  TSGetKSPIterations()
4554: @*/
4555: PetscErrorCode TSGetSNESIterations(TS ts,PetscInt *nits)
4556: {
4559:   *nits = ts->snes_its;
4560:   return 0;
4561: }

4563: /*@
4564:    TSGetKSPIterations - Gets the total number of linear iterations
4565:    used by the time integrator.

4567:    Not Collective

4569:    Input Parameter:
4570: .  ts - TS context

4572:    Output Parameter:
4573: .  lits - number of linear iterations

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

4578:    Level: intermediate

4580: .seealso:  TSGetSNESIterations(), SNESGetKSPIterations()
4581: @*/
4582: PetscErrorCode TSGetKSPIterations(TS ts,PetscInt *lits)
4583: {
4586:   *lits = ts->ksp_its;
4587:   return 0;
4588: }

4590: /*@
4591:    TSGetStepRejections - Gets the total number of rejected steps.

4593:    Not Collective

4595:    Input Parameter:
4596: .  ts - TS context

4598:    Output Parameter:
4599: .  rejects - number of steps rejected

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

4604:    Level: intermediate

4606: .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetSNESFailures(), TSSetMaxSNESFailures(), TSSetErrorIfStepFails()
4607: @*/
4608: PetscErrorCode TSGetStepRejections(TS ts,PetscInt *rejects)
4609: {
4612:   *rejects = ts->reject;
4613:   return 0;
4614: }

4616: /*@
4617:    TSGetSNESFailures - Gets the total number of failed SNES solves

4619:    Not Collective

4621:    Input Parameter:
4622: .  ts - TS context

4624:    Output Parameter:
4625: .  fails - number of failed nonlinear solves

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

4630:    Level: intermediate

4632: .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSSetMaxSNESFailures()
4633: @*/
4634: PetscErrorCode TSGetSNESFailures(TS ts,PetscInt *fails)
4635: {
4638:   *fails = ts->num_snes_failures;
4639:   return 0;
4640: }

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

4645:    Not Collective

4647:    Input Parameters:
4648: +  ts - TS context
4649: -  rejects - maximum number of rejected steps, pass -1 for unlimited

4651:    Notes:
4652:    The counter is reset to zero for each step

4654:    Options Database Key:
4655: .  -ts_max_reject - Maximum number of step rejections before a step fails

4657:    Level: intermediate

4659: .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxSNESFailures(), TSGetStepRejections(), TSGetSNESFailures(), TSSetErrorIfStepFails(), TSGetConvergedReason()
4660: @*/
4661: PetscErrorCode TSSetMaxStepRejections(TS ts,PetscInt rejects)
4662: {
4664:   ts->max_reject = rejects;
4665:   return 0;
4666: }

4668: /*@
4669:    TSSetMaxSNESFailures - Sets the maximum number of failed SNES solves

4671:    Not Collective

4673:    Input Parameters:
4674: +  ts - TS context
4675: -  fails - maximum number of failed nonlinear solves, pass -1 for unlimited

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

4680:    Options Database Key:
4681: .  -ts_max_snes_failures - Maximum number of nonlinear solve failures

4683:    Level: intermediate

4685: .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSGetSNESFailures(), SNESGetConvergedReason(), TSGetConvergedReason()
4686: @*/
4687: PetscErrorCode TSSetMaxSNESFailures(TS ts,PetscInt fails)
4688: {
4690:   ts->max_snes_failures = fails;
4691:   return 0;
4692: }

4694: /*@
4695:    TSSetErrorIfStepFails - Error if no step succeeds

4697:    Not Collective

4699:    Input Parameters:
4700: +  ts - TS context
4701: -  err - PETSC_TRUE to error if no step succeeds, PETSC_FALSE to return without failure

4703:    Options Database Key:
4704: .  -ts_error_if_step_fails - Error if no step succeeds

4706:    Level: intermediate

4708: .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSGetSNESFailures(), TSSetErrorIfStepFails(), TSGetConvergedReason()
4709: @*/
4710: PetscErrorCode TSSetErrorIfStepFails(TS ts,PetscBool err)
4711: {
4713:   ts->errorifstepfailed = err;
4714:   return 0;
4715: }

4717: /*@
4718:    TSGetAdapt - Get the adaptive controller context for the current method

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

4722:    Input Parameter:
4723: .  ts - time stepping context

4725:    Output Parameter:
4726: .  adapt - adaptive controller

4728:    Level: intermediate

4730: .seealso: TSAdapt, TSAdaptSetType(), TSAdaptChoose()
4731: @*/
4732: PetscErrorCode TSGetAdapt(TS ts,TSAdapt *adapt)
4733: {
4736:   if (!ts->adapt) {
4737:     TSAdaptCreate(PetscObjectComm((PetscObject)ts),&ts->adapt);
4738:     PetscLogObjectParent((PetscObject)ts,(PetscObject)ts->adapt);
4739:     PetscObjectIncrementTabLevel((PetscObject)ts->adapt,(PetscObject)ts,1);
4740:   }
4741:   *adapt = ts->adapt;
4742:   return 0;
4743: }

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

4748:    Logically Collective

4750:    Input Parameters:
4751: +  ts - time integration context
4752: .  atol - scalar absolute tolerances, PETSC_DECIDE to leave current value
4753: .  vatol - vector of absolute tolerances or NULL, used in preference to atol if present
4754: .  rtol - scalar relative tolerances, PETSC_DECIDE to leave current value
4755: -  vrtol - vector of relative tolerances or NULL, used in preference to atol if present

4757:    Options Database keys:
4758: +  -ts_rtol <rtol> - relative tolerance for local truncation error
4759: -  -ts_atol <atol> - Absolute tolerance for local truncation error

4761:    Notes:
4762:    With PETSc's implicit schemes for DAE problems, the calculation of the local truncation error
4763:    (LTE) includes both the differential and the algebraic variables. If one wants the LTE to be
4764:    computed only for the differential or the algebraic part then this can be done using the vector of
4765:    tolerances vatol. For example, by setting the tolerance vector with the desired tolerance for the
4766:    differential part and infinity for the algebraic part, the LTE calculation will include only the
4767:    differential variables.

4769:    Level: beginner

4771: .seealso: TS, TSAdapt, TSErrorWeightedNorm(), TSGetTolerances()
4772: @*/
4773: PetscErrorCode TSSetTolerances(TS ts,PetscReal atol,Vec vatol,PetscReal rtol,Vec vrtol)
4774: {
4775:   if (atol != PETSC_DECIDE && atol != PETSC_DEFAULT) ts->atol = atol;
4776:   if (vatol) {
4777:     PetscObjectReference((PetscObject)vatol);
4778:     VecDestroy(&ts->vatol);
4779:     ts->vatol = vatol;
4780:   }
4781:   if (rtol != PETSC_DECIDE && rtol != PETSC_DEFAULT) ts->rtol = rtol;
4782:   if (vrtol) {
4783:     PetscObjectReference((PetscObject)vrtol);
4784:     VecDestroy(&ts->vrtol);
4785:     ts->vrtol = vrtol;
4786:   }
4787:   return 0;
4788: }

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

4793:    Logically Collective

4795:    Input Parameter:
4796: .  ts - time integration context

4798:    Output Parameters:
4799: +  atol - scalar absolute tolerances, NULL to ignore
4800: .  vatol - vector of absolute tolerances, NULL to ignore
4801: .  rtol - scalar relative tolerances, NULL to ignore
4802: -  vrtol - vector of relative tolerances, NULL to ignore

4804:    Level: beginner

4806: .seealso: TS, TSAdapt, TSErrorWeightedNorm(), TSSetTolerances()
4807: @*/
4808: PetscErrorCode TSGetTolerances(TS ts,PetscReal *atol,Vec *vatol,PetscReal *rtol,Vec *vrtol)
4809: {
4810:   if (atol)  *atol  = ts->atol;
4811:   if (vatol) *vatol = ts->vatol;
4812:   if (rtol)  *rtol  = ts->rtol;
4813:   if (vrtol) *vrtol = ts->vrtol;
4814:   return 0;
4815: }

4817: /*@
4818:    TSErrorWeightedNorm2 - compute a weighted 2-norm of the difference between two state vectors

4820:    Collective on TS

4822:    Input Parameters:
4823: +  ts - time stepping context
4824: .  U - state vector, usually ts->vec_sol
4825: -  Y - state vector to be compared to U

4827:    Output Parameters:
4828: +  norm - weighted norm, a value of 1.0 means that the error matches the tolerances
4829: .  norma - weighted norm based on the absolute tolerance, a value of 1.0 means that the error matches the tolerances
4830: -  normr - weighted norm based on the relative tolerance, a value of 1.0 means that the error matches the tolerances

4832:    Level: developer

4834: .seealso: TSErrorWeightedNorm(), TSErrorWeightedNormInfinity()
4835: @*/
4836: PetscErrorCode TSErrorWeightedNorm2(TS ts,Vec U,Vec Y,PetscReal *norm,PetscReal *norma,PetscReal *normr)
4837: {
4838:   PetscInt          i,n,N,rstart;
4839:   PetscInt          n_loc,na_loc,nr_loc;
4840:   PetscReal         n_glb,na_glb,nr_glb;
4841:   const PetscScalar *u,*y;
4842:   PetscReal         sum,suma,sumr,gsum,gsuma,gsumr,diff;
4843:   PetscReal         tol,tola,tolr;
4844:   PetscReal         err_loc[6],err_glb[6];


4857:   VecGetSize(U,&N);
4858:   VecGetLocalSize(U,&n);
4859:   VecGetOwnershipRange(U,&rstart,NULL);
4860:   VecGetArrayRead(U,&u);
4861:   VecGetArrayRead(Y,&y);
4862:   sum  = 0.; n_loc  = 0;
4863:   suma = 0.; na_loc = 0;
4864:   sumr = 0.; nr_loc = 0;
4865:   if (ts->vatol && ts->vrtol) {
4866:     const PetscScalar *atol,*rtol;
4867:     VecGetArrayRead(ts->vatol,&atol);
4868:     VecGetArrayRead(ts->vrtol,&rtol);
4869:     for (i=0; i<n; i++) {
4870:       SkipSmallValue(y[i],u[i],ts->adapt->ignore_max);
4871:       diff = PetscAbsScalar(y[i] - u[i]);
4872:       tola = PetscRealPart(atol[i]);
4873:       if (tola>0.) {
4874:         suma  += PetscSqr(diff/tola);
4875:         na_loc++;
4876:       }
4877:       tolr = PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
4878:       if (tolr>0.) {
4879:         sumr  += PetscSqr(diff/tolr);
4880:         nr_loc++;
4881:       }
4882:       tol=tola+tolr;
4883:       if (tol>0.) {
4884:         sum  += PetscSqr(diff/tol);
4885:         n_loc++;
4886:       }
4887:     }
4888:     VecRestoreArrayRead(ts->vatol,&atol);
4889:     VecRestoreArrayRead(ts->vrtol,&rtol);
4890:   } else if (ts->vatol) {       /* vector atol, scalar rtol */
4891:     const PetscScalar *atol;
4892:     VecGetArrayRead(ts->vatol,&atol);
4893:     for (i=0; i<n; i++) {
4894:       SkipSmallValue(y[i],u[i],ts->adapt->ignore_max);
4895:       diff = PetscAbsScalar(y[i] - u[i]);
4896:       tola = PetscRealPart(atol[i]);
4897:       if (tola>0.) {
4898:         suma  += PetscSqr(diff/tola);
4899:         na_loc++;
4900:       }
4901:       tolr = ts->rtol * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
4902:       if (tolr>0.) {
4903:         sumr  += PetscSqr(diff/tolr);
4904:         nr_loc++;
4905:       }
4906:       tol=tola+tolr;
4907:       if (tol>0.) {
4908:         sum  += PetscSqr(diff/tol);
4909:         n_loc++;
4910:       }
4911:     }
4912:     VecRestoreArrayRead(ts->vatol,&atol);
4913:   } else if (ts->vrtol) {       /* scalar atol, vector rtol */
4914:     const PetscScalar *rtol;
4915:     VecGetArrayRead(ts->vrtol,&rtol);
4916:     for (i=0; i<n; i++) {
4917:       SkipSmallValue(y[i],u[i],ts->adapt->ignore_max);
4918:       diff = PetscAbsScalar(y[i] - u[i]);
4919:       tola = ts->atol;
4920:       if (tola>0.) {
4921:         suma  += PetscSqr(diff/tola);
4922:         na_loc++;
4923:       }
4924:       tolr = PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
4925:       if (tolr>0.) {
4926:         sumr  += PetscSqr(diff/tolr);
4927:         nr_loc++;
4928:       }
4929:       tol=tola+tolr;
4930:       if (tol>0.) {
4931:         sum  += PetscSqr(diff/tol);
4932:         n_loc++;
4933:       }
4934:     }
4935:     VecRestoreArrayRead(ts->vrtol,&rtol);
4936:   } else {                      /* scalar atol, scalar rtol */
4937:     for (i=0; i<n; i++) {
4938:       SkipSmallValue(y[i],u[i],ts->adapt->ignore_max);
4939:       diff = PetscAbsScalar(y[i] - u[i]);
4940:       tola = ts->atol;
4941:       if (tola>0.) {
4942:         suma  += PetscSqr(diff/tola);
4943:         na_loc++;
4944:       }
4945:       tolr = ts->rtol * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
4946:       if (tolr>0.) {
4947:         sumr  += PetscSqr(diff/tolr);
4948:         nr_loc++;
4949:       }
4950:       tol=tola+tolr;
4951:       if (tol>0.) {
4952:         sum  += PetscSqr(diff/tol);
4953:         n_loc++;
4954:       }
4955:     }
4956:   }
4957:   VecRestoreArrayRead(U,&u);
4958:   VecRestoreArrayRead(Y,&y);

4960:   err_loc[0] = sum;
4961:   err_loc[1] = suma;
4962:   err_loc[2] = sumr;
4963:   err_loc[3] = (PetscReal)n_loc;
4964:   err_loc[4] = (PetscReal)na_loc;
4965:   err_loc[5] = (PetscReal)nr_loc;

4967:   MPIU_Allreduce(err_loc,err_glb,6,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)ts));

4969:   gsum   = err_glb[0];
4970:   gsuma  = err_glb[1];
4971:   gsumr  = err_glb[2];
4972:   n_glb  = err_glb[3];
4973:   na_glb = err_glb[4];
4974:   nr_glb = err_glb[5];

4976:   *norm  = 0.;
4977:   if (n_glb>0.) {*norm  = PetscSqrtReal(gsum  / n_glb);}
4978:   *norma = 0.;
4979:   if (na_glb>0.) {*norma = PetscSqrtReal(gsuma / na_glb);}
4980:   *normr = 0.;
4981:   if (nr_glb>0.) {*normr = PetscSqrtReal(gsumr / nr_glb);}

4986:   return 0;
4987: }

4989: /*@
4990:    TSErrorWeightedNormInfinity - compute a weighted infinity-norm of the difference between two state vectors

4992:    Collective on TS

4994:    Input Parameters:
4995: +  ts - time stepping context
4996: .  U - state vector, usually ts->vec_sol
4997: -  Y - state vector to be compared to U

4999:    Output Parameters:
5000: +  norm - weighted norm, a value of 1.0 means that the error matches the tolerances
5001: .  norma - weighted norm based on the absolute tolerance, a value of 1.0 means that the error matches the tolerances
5002: -  normr - weighted norm based on the relative tolerance, a value of 1.0 means that the error matches the tolerances

5004:    Level: developer

5006: .seealso: TSErrorWeightedNorm(), TSErrorWeightedNorm2()
5007: @*/
5008: PetscErrorCode TSErrorWeightedNormInfinity(TS ts,Vec U,Vec Y,PetscReal *norm,PetscReal *norma,PetscReal *normr)
5009: {
5010:   PetscInt          i,n,N,rstart;
5011:   const PetscScalar *u,*y;
5012:   PetscReal         max,gmax,maxa,gmaxa,maxr,gmaxr;
5013:   PetscReal         tol,tola,tolr,diff;
5014:   PetscReal         err_loc[3],err_glb[3];


5027:   VecGetSize(U,&N);
5028:   VecGetLocalSize(U,&n);
5029:   VecGetOwnershipRange(U,&rstart,NULL);
5030:   VecGetArrayRead(U,&u);
5031:   VecGetArrayRead(Y,&y);

5033:   max=0.;
5034:   maxa=0.;
5035:   maxr=0.;

5037:   if (ts->vatol && ts->vrtol) {     /* vector atol, vector rtol */
5038:     const PetscScalar *atol,*rtol;
5039:     VecGetArrayRead(ts->vatol,&atol);
5040:     VecGetArrayRead(ts->vrtol,&rtol);

5042:     for (i=0; i<n; i++) {
5043:       SkipSmallValue(y[i],u[i],ts->adapt->ignore_max);
5044:       diff = PetscAbsScalar(y[i] - u[i]);
5045:       tola = PetscRealPart(atol[i]);
5046:       tolr = PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5047:       tol  = tola+tolr;
5048:       if (tola>0.) {
5049:         maxa = PetscMax(maxa,diff / tola);
5050:       }
5051:       if (tolr>0.) {
5052:         maxr = PetscMax(maxr,diff / tolr);
5053:       }
5054:       if (tol>0.) {
5055:         max = PetscMax(max,diff / tol);
5056:       }
5057:     }
5058:     VecRestoreArrayRead(ts->vatol,&atol);
5059:     VecRestoreArrayRead(ts->vrtol,&rtol);
5060:   } else if (ts->vatol) {       /* vector atol, scalar rtol */
5061:     const PetscScalar *atol;
5062:     VecGetArrayRead(ts->vatol,&atol);
5063:     for (i=0; i<n; i++) {
5064:       SkipSmallValue(y[i],u[i],ts->adapt->ignore_max);
5065:       diff = PetscAbsScalar(y[i] - u[i]);
5066:       tola = PetscRealPart(atol[i]);
5067:       tolr = ts->rtol  * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5068:       tol  = tola+tolr;
5069:       if (tola>0.) {
5070:         maxa = PetscMax(maxa,diff / tola);
5071:       }
5072:       if (tolr>0.) {
5073:         maxr = PetscMax(maxr,diff / tolr);
5074:       }
5075:       if (tol>0.) {
5076:         max = PetscMax(max,diff / tol);
5077:       }
5078:     }
5079:     VecRestoreArrayRead(ts->vatol,&atol);
5080:   } else if (ts->vrtol) {       /* scalar atol, vector rtol */
5081:     const PetscScalar *rtol;
5082:     VecGetArrayRead(ts->vrtol,&rtol);

5084:     for (i=0; i<n; i++) {
5085:       SkipSmallValue(y[i],u[i],ts->adapt->ignore_max);
5086:       diff = PetscAbsScalar(y[i] - u[i]);
5087:       tola = ts->atol;
5088:       tolr = PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5089:       tol  = tola+tolr;
5090:       if (tola>0.) {
5091:         maxa = PetscMax(maxa,diff / tola);
5092:       }
5093:       if (tolr>0.) {
5094:         maxr = PetscMax(maxr,diff / tolr);
5095:       }
5096:       if (tol>0.) {
5097:         max = PetscMax(max,diff / tol);
5098:       }
5099:     }
5100:     VecRestoreArrayRead(ts->vrtol,&rtol);
5101:   } else {                      /* scalar atol, scalar rtol */

5103:     for (i=0; i<n; i++) {
5104:       SkipSmallValue(y[i],u[i],ts->adapt->ignore_max);
5105:       diff = PetscAbsScalar(y[i] - u[i]);
5106:       tola = ts->atol;
5107:       tolr = ts->rtol * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5108:       tol  = tola+tolr;
5109:       if (tola>0.) {
5110:         maxa = PetscMax(maxa,diff / tola);
5111:       }
5112:       if (tolr>0.) {
5113:         maxr = PetscMax(maxr,diff / tolr);
5114:       }
5115:       if (tol>0.) {
5116:         max = PetscMax(max,diff / tol);
5117:       }
5118:     }
5119:   }
5120:   VecRestoreArrayRead(U,&u);
5121:   VecRestoreArrayRead(Y,&y);
5122:   err_loc[0] = max;
5123:   err_loc[1] = maxa;
5124:   err_loc[2] = maxr;
5125:   MPIU_Allreduce(err_loc,err_glb,3,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)ts));
5126:   gmax   = err_glb[0];
5127:   gmaxa  = err_glb[1];
5128:   gmaxr  = err_glb[2];

5130:   *norm = gmax;
5131:   *norma = gmaxa;
5132:   *normr = gmaxr;
5136:   return 0;
5137: }

5139: /*@
5140:    TSErrorWeightedNorm - compute a weighted norm of the difference between two state vectors based on supplied absolute and relative tolerances

5142:    Collective on TS

5144:    Input Parameters:
5145: +  ts - time stepping context
5146: .  U - state vector, usually ts->vec_sol
5147: .  Y - state vector to be compared to U
5148: -  wnormtype - norm type, either NORM_2 or NORM_INFINITY

5150:    Output Parameters:
5151: +  norm  - weighted norm, a value of 1.0 achieves a balance between absolute and relative tolerances
5152: .  norma - weighted norm, a value of 1.0 means that the error meets the absolute tolerance set by the user
5153: -  normr - weighted norm, a value of 1.0 means that the error meets the relative tolerance set by the user

5155:    Options Database Keys:
5156: .  -ts_adapt_wnormtype <wnormtype> - 2, INFINITY

5158:    Level: developer

5160: .seealso: TSErrorWeightedNormInfinity(), TSErrorWeightedNorm2(), TSErrorWeightedENorm
5161: @*/
5162: PetscErrorCode TSErrorWeightedNorm(TS ts,Vec U,Vec Y,NormType wnormtype,PetscReal *norm,PetscReal *norma,PetscReal *normr)
5163: {
5164:   if (wnormtype == NORM_2) {
5165:     TSErrorWeightedNorm2(ts,U,Y,norm,norma,normr);
5166:   } else if (wnormtype == NORM_INFINITY) {
5167:     TSErrorWeightedNormInfinity(ts,U,Y,norm,norma,normr);
5168:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for norm type %s",NormTypes[wnormtype]);
5169:   return 0;
5170: }

5172: /*@
5173:    TSErrorWeightedENorm2 - compute a weighted 2 error norm based on supplied absolute and relative tolerances

5175:    Collective on TS

5177:    Input Parameters:
5178: +  ts - time stepping context
5179: .  E - error vector
5180: .  U - state vector, usually ts->vec_sol
5181: -  Y - state vector, previous time step

5183:    Output Parameters:
5184: +  norm - weighted norm, a value of 1.0 means that the error matches the tolerances
5185: .  norma - weighted norm based on the absolute tolerance, a value of 1.0 means that the error matches the tolerances
5186: -  normr - weighted norm based on the relative tolerance, a value of 1.0 means that the error matches the tolerances

5188:    Level: developer

5190: .seealso: TSErrorWeightedENorm(), TSErrorWeightedENormInfinity()
5191: @*/
5192: PetscErrorCode TSErrorWeightedENorm2(TS ts,Vec E,Vec U,Vec Y,PetscReal *norm,PetscReal *norma,PetscReal *normr)
5193: {
5194:   PetscInt          i,n,N,rstart;
5195:   PetscInt          n_loc,na_loc,nr_loc;
5196:   PetscReal         n_glb,na_glb,nr_glb;
5197:   const PetscScalar *e,*u,*y;
5198:   PetscReal         err,sum,suma,sumr,gsum,gsuma,gsumr;
5199:   PetscReal         tol,tola,tolr;
5200:   PetscReal         err_loc[6],err_glb[6];


5215:   VecGetSize(E,&N);
5216:   VecGetLocalSize(E,&n);
5217:   VecGetOwnershipRange(E,&rstart,NULL);
5218:   VecGetArrayRead(E,&e);
5219:   VecGetArrayRead(U,&u);
5220:   VecGetArrayRead(Y,&y);
5221:   sum  = 0.; n_loc  = 0;
5222:   suma = 0.; na_loc = 0;
5223:   sumr = 0.; nr_loc = 0;
5224:   if (ts->vatol && ts->vrtol) {
5225:     const PetscScalar *atol,*rtol;
5226:     VecGetArrayRead(ts->vatol,&atol);
5227:     VecGetArrayRead(ts->vrtol,&rtol);
5228:     for (i=0; i<n; i++) {
5229:       SkipSmallValue(y[i],u[i],ts->adapt->ignore_max);
5230:       err = PetscAbsScalar(e[i]);
5231:       tola = PetscRealPart(atol[i]);
5232:       if (tola>0.) {
5233:         suma  += PetscSqr(err/tola);
5234:         na_loc++;
5235:       }
5236:       tolr = PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5237:       if (tolr>0.) {
5238:         sumr  += PetscSqr(err/tolr);
5239:         nr_loc++;
5240:       }
5241:       tol=tola+tolr;
5242:       if (tol>0.) {
5243:         sum  += PetscSqr(err/tol);
5244:         n_loc++;
5245:       }
5246:     }
5247:     VecRestoreArrayRead(ts->vatol,&atol);
5248:     VecRestoreArrayRead(ts->vrtol,&rtol);
5249:   } else if (ts->vatol) {       /* vector atol, scalar rtol */
5250:     const PetscScalar *atol;
5251:     VecGetArrayRead(ts->vatol,&atol);
5252:     for (i=0; i<n; i++) {
5253:       SkipSmallValue(y[i],u[i],ts->adapt->ignore_max);
5254:       err = PetscAbsScalar(e[i]);
5255:       tola = PetscRealPart(atol[i]);
5256:       if (tola>0.) {
5257:         suma  += PetscSqr(err/tola);
5258:         na_loc++;
5259:       }
5260:       tolr = ts->rtol * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5261:       if (tolr>0.) {
5262:         sumr  += PetscSqr(err/tolr);
5263:         nr_loc++;
5264:       }
5265:       tol=tola+tolr;
5266:       if (tol>0.) {
5267:         sum  += PetscSqr(err/tol);
5268:         n_loc++;
5269:       }
5270:     }
5271:     VecRestoreArrayRead(ts->vatol,&atol);
5272:   } else if (ts->vrtol) {       /* scalar atol, vector rtol */
5273:     const PetscScalar *rtol;
5274:     VecGetArrayRead(ts->vrtol,&rtol);
5275:     for (i=0; i<n; i++) {
5276:       SkipSmallValue(y[i],u[i],ts->adapt->ignore_max);
5277:       err = PetscAbsScalar(e[i]);
5278:       tola = ts->atol;
5279:       if (tola>0.) {
5280:         suma  += PetscSqr(err/tola);
5281:         na_loc++;
5282:       }
5283:       tolr = PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5284:       if (tolr>0.) {
5285:         sumr  += PetscSqr(err/tolr);
5286:         nr_loc++;
5287:       }
5288:       tol=tola+tolr;
5289:       if (tol>0.) {
5290:         sum  += PetscSqr(err/tol);
5291:         n_loc++;
5292:       }
5293:     }
5294:     VecRestoreArrayRead(ts->vrtol,&rtol);
5295:   } else {                      /* scalar atol, scalar rtol */
5296:     for (i=0; i<n; i++) {
5297:       SkipSmallValue(y[i],u[i],ts->adapt->ignore_max);
5298:       err = PetscAbsScalar(e[i]);
5299:       tola = ts->atol;
5300:       if (tola>0.) {
5301:         suma  += PetscSqr(err/tola);
5302:         na_loc++;
5303:       }
5304:       tolr = ts->rtol * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5305:       if (tolr>0.) {
5306:         sumr  += PetscSqr(err/tolr);
5307:         nr_loc++;
5308:       }
5309:       tol=tola+tolr;
5310:       if (tol>0.) {
5311:         sum  += PetscSqr(err/tol);
5312:         n_loc++;
5313:       }
5314:     }
5315:   }
5316:   VecRestoreArrayRead(E,&e);
5317:   VecRestoreArrayRead(U,&u);
5318:   VecRestoreArrayRead(Y,&y);

5320:   err_loc[0] = sum;
5321:   err_loc[1] = suma;
5322:   err_loc[2] = sumr;
5323:   err_loc[3] = (PetscReal)n_loc;
5324:   err_loc[4] = (PetscReal)na_loc;
5325:   err_loc[5] = (PetscReal)nr_loc;

5327:   MPIU_Allreduce(err_loc,err_glb,6,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)ts));

5329:   gsum   = err_glb[0];
5330:   gsuma  = err_glb[1];
5331:   gsumr  = err_glb[2];
5332:   n_glb  = err_glb[3];
5333:   na_glb = err_glb[4];
5334:   nr_glb = err_glb[5];

5336:   *norm  = 0.;
5337:   if (n_glb>0.) {*norm  = PetscSqrtReal(gsum  / n_glb);}
5338:   *norma = 0.;
5339:   if (na_glb>0.) {*norma = PetscSqrtReal(gsuma / na_glb);}
5340:   *normr = 0.;
5341:   if (nr_glb>0.) {*normr = PetscSqrtReal(gsumr / nr_glb);}

5346:   return 0;
5347: }

5349: /*@
5350:    TSErrorWeightedENormInfinity - compute a weighted infinity error norm based on supplied absolute and relative tolerances
5351:    Collective on TS

5353:    Input Parameters:
5354: +  ts - time stepping context
5355: .  E - error vector
5356: .  U - state vector, usually ts->vec_sol
5357: -  Y - state vector, previous time step

5359:    Output Parameters:
5360: +  norm - weighted norm, a value of 1.0 means that the error matches the tolerances
5361: .  norma - weighted norm based on the absolute tolerance, a value of 1.0 means that the error matches the tolerances
5362: -  normr - weighted norm based on the relative tolerance, a value of 1.0 means that the error matches the tolerances

5364:    Level: developer

5366: .seealso: TSErrorWeightedENorm(), TSErrorWeightedENorm2()
5367: @*/
5368: PetscErrorCode TSErrorWeightedENormInfinity(TS ts,Vec E,Vec U,Vec Y,PetscReal *norm,PetscReal *norma,PetscReal *normr)
5369: {
5370:   PetscInt          i,n,N,rstart;
5371:   const PetscScalar *e,*u,*y;
5372:   PetscReal         err,max,gmax,maxa,gmaxa,maxr,gmaxr;
5373:   PetscReal         tol,tola,tolr;
5374:   PetscReal         err_loc[3],err_glb[3];


5389:   VecGetSize(E,&N);
5390:   VecGetLocalSize(E,&n);
5391:   VecGetOwnershipRange(E,&rstart,NULL);
5392:   VecGetArrayRead(E,&e);
5393:   VecGetArrayRead(U,&u);
5394:   VecGetArrayRead(Y,&y);

5396:   max=0.;
5397:   maxa=0.;
5398:   maxr=0.;

5400:   if (ts->vatol && ts->vrtol) {     /* vector atol, vector rtol */
5401:     const PetscScalar *atol,*rtol;
5402:     VecGetArrayRead(ts->vatol,&atol);
5403:     VecGetArrayRead(ts->vrtol,&rtol);

5405:     for (i=0; i<n; i++) {
5406:       SkipSmallValue(y[i],u[i],ts->adapt->ignore_max);
5407:       err = PetscAbsScalar(e[i]);
5408:       tola = PetscRealPart(atol[i]);
5409:       tolr = PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5410:       tol  = tola+tolr;
5411:       if (tola>0.) {
5412:         maxa = PetscMax(maxa,err / tola);
5413:       }
5414:       if (tolr>0.) {
5415:         maxr = PetscMax(maxr,err / tolr);
5416:       }
5417:       if (tol>0.) {
5418:         max = PetscMax(max,err / tol);
5419:       }
5420:     }
5421:     VecRestoreArrayRead(ts->vatol,&atol);
5422:     VecRestoreArrayRead(ts->vrtol,&rtol);
5423:   } else if (ts->vatol) {       /* vector atol, scalar rtol */
5424:     const PetscScalar *atol;
5425:     VecGetArrayRead(ts->vatol,&atol);
5426:     for (i=0; i<n; i++) {
5427:       SkipSmallValue(y[i],u[i],ts->adapt->ignore_max);
5428:       err = PetscAbsScalar(e[i]);
5429:       tola = PetscRealPart(atol[i]);
5430:       tolr = ts->rtol  * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5431:       tol  = tola+tolr;
5432:       if (tola>0.) {
5433:         maxa = PetscMax(maxa,err / tola);
5434:       }
5435:       if (tolr>0.) {
5436:         maxr = PetscMax(maxr,err / tolr);
5437:       }
5438:       if (tol>0.) {
5439:         max = PetscMax(max,err / tol);
5440:       }
5441:     }
5442:     VecRestoreArrayRead(ts->vatol,&atol);
5443:   } else if (ts->vrtol) {       /* scalar atol, vector rtol */
5444:     const PetscScalar *rtol;
5445:     VecGetArrayRead(ts->vrtol,&rtol);

5447:     for (i=0; i<n; i++) {
5448:       SkipSmallValue(y[i],u[i],ts->adapt->ignore_max);
5449:       err = PetscAbsScalar(e[i]);
5450:       tola = ts->atol;
5451:       tolr = PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5452:       tol  = tola+tolr;
5453:       if (tola>0.) {
5454:         maxa = PetscMax(maxa,err / tola);
5455:       }
5456:       if (tolr>0.) {
5457:         maxr = PetscMax(maxr,err / tolr);
5458:       }
5459:       if (tol>0.) {
5460:         max = PetscMax(max,err / tol);
5461:       }
5462:     }
5463:     VecRestoreArrayRead(ts->vrtol,&rtol);
5464:   } else {                      /* scalar atol, scalar rtol */

5466:     for (i=0; i<n; i++) {
5467:       SkipSmallValue(y[i],u[i],ts->adapt->ignore_max);
5468:       err = PetscAbsScalar(e[i]);
5469:       tola = ts->atol;
5470:       tolr = ts->rtol * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5471:       tol  = tola+tolr;
5472:       if (tola>0.) {
5473:         maxa = PetscMax(maxa,err / tola);
5474:       }
5475:       if (tolr>0.) {
5476:         maxr = PetscMax(maxr,err / tolr);
5477:       }
5478:       if (tol>0.) {
5479:         max = PetscMax(max,err / tol);
5480:       }
5481:     }
5482:   }
5483:   VecRestoreArrayRead(E,&e);
5484:   VecRestoreArrayRead(U,&u);
5485:   VecRestoreArrayRead(Y,&y);
5486:   err_loc[0] = max;
5487:   err_loc[1] = maxa;
5488:   err_loc[2] = maxr;
5489:   MPIU_Allreduce(err_loc,err_glb,3,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)ts));
5490:   gmax   = err_glb[0];
5491:   gmaxa  = err_glb[1];
5492:   gmaxr  = err_glb[2];

5494:   *norm = gmax;
5495:   *norma = gmaxa;
5496:   *normr = gmaxr;
5500:   return 0;
5501: }

5503: /*@
5504:    TSErrorWeightedENorm - compute a weighted error norm based on supplied absolute and relative tolerances

5506:    Collective on TS

5508:    Input Parameters:
5509: +  ts - time stepping context
5510: .  E - error vector
5511: .  U - state vector, usually ts->vec_sol
5512: .  Y - state vector, previous time step
5513: -  wnormtype - norm type, either NORM_2 or NORM_INFINITY

5515:    Output Parameters:
5516: +  norm  - weighted norm, a value of 1.0 achieves a balance between absolute and relative tolerances
5517: .  norma - weighted norm, a value of 1.0 means that the error meets the absolute tolerance set by the user
5518: -  normr - weighted norm, a value of 1.0 means that the error meets the relative tolerance set by the user

5520:    Options Database Keys:
5521: .  -ts_adapt_wnormtype <wnormtype> - 2, INFINITY

5523:    Level: developer

5525: .seealso: TSErrorWeightedENormInfinity(), TSErrorWeightedENorm2(), TSErrorWeightedNormInfinity(), TSErrorWeightedNorm2()
5526: @*/
5527: PetscErrorCode TSErrorWeightedENorm(TS ts,Vec E,Vec U,Vec Y,NormType wnormtype,PetscReal *norm,PetscReal *norma,PetscReal *normr)
5528: {
5529:   if (wnormtype == NORM_2) {
5530:     TSErrorWeightedENorm2(ts,E,U,Y,norm,norma,normr);
5531:   } else if (wnormtype == NORM_INFINITY) {
5532:     TSErrorWeightedENormInfinity(ts,E,U,Y,norm,norma,normr);
5533:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for norm type %s",NormTypes[wnormtype]);
5534:   return 0;
5535: }

5537: /*@
5538:    TSSetCFLTimeLocal - Set the local CFL constraint relative to forward Euler

5540:    Logically Collective on TS

5542:    Input Parameters:
5543: +  ts - time stepping context
5544: -  cfltime - maximum stable time step if using forward Euler (value can be different on each process)

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

5549:    Level: intermediate

5551: .seealso: TSGetCFLTime(), TSADAPTCFL
5552: @*/
5553: PetscErrorCode TSSetCFLTimeLocal(TS ts,PetscReal cfltime)
5554: {
5556:   ts->cfltime_local = cfltime;
5557:   ts->cfltime       = -1.;
5558:   return 0;
5559: }

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

5564:    Collective on TS

5566:    Input Parameter:
5567: .  ts - time stepping context

5569:    Output Parameter:
5570: .  cfltime - maximum stable time step for forward Euler

5572:    Level: advanced

5574: .seealso: TSSetCFLTimeLocal()
5575: @*/
5576: PetscErrorCode TSGetCFLTime(TS ts,PetscReal *cfltime)
5577: {
5578:   if (ts->cfltime < 0) {
5579:     MPIU_Allreduce(&ts->cfltime_local,&ts->cfltime,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)ts));
5580:   }
5581:   *cfltime = ts->cfltime;
5582:   return 0;
5583: }

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

5588:    Input Parameters:
5589: +  ts   - the TS context.
5590: .  xl   - lower bound.
5591: -  xu   - upper bound.

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

5597:    Level: advanced

5599: @*/
5600: PetscErrorCode TSVISetVariableBounds(TS ts, Vec xl, Vec xu)
5601: {
5602:   SNES           snes;

5604:   TSGetSNES(ts,&snes);
5605:   SNESVISetVariableBounds(snes,xl,xu);
5606:   return 0;
5607: }

5609: /*@
5610:    TSComputeLinearStability - computes the linear stability function at a point

5612:    Collective on TS

5614:    Input Parameters:
5615: +  ts - the TS context
5616: -  xr,xi - real and imaginary part of input arguments

5618:    Output Parameters:
5619: .  yr,yi - real and imaginary part of function value

5621:    Level: developer

5623: .seealso: TSSetRHSFunction(), TSComputeIFunction()
5624: @*/
5625: PetscErrorCode TSComputeLinearStability(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
5626: {
5629:   (*ts->ops->linearstability)(ts,xr,xi,yr,yi);
5630:   return 0;
5631: }

5633: /*@
5634:    TSRestartStep - Flags the solver to restart the next step

5636:    Collective on TS

5638:    Input Parameter:
5639: .  ts - the TS context obtained from TSCreate()

5641:    Level: advanced

5643:    Notes:
5644:    Multistep methods like BDF or Runge-Kutta methods with FSAL property require restarting the solver in the event of
5645:    discontinuities. These discontinuities may be introduced as a consequence of explicitly modifications to the solution
5646:    vector (which PETSc attempts to detect and handle) or problem coefficients (which PETSc is not able to detect). For
5647:    the sake of correctness and maximum safety, users are expected to call TSRestart() whenever they introduce
5648:    discontinuities in callback routines (e.g. prestep and poststep routines, or implicit/rhs function routines with
5649:    discontinuous source terms).

5651: .seealso: TSSolve(), TSSetPreStep(), TSSetPostStep()
5652: @*/
5653: PetscErrorCode TSRestartStep(TS ts)
5654: {
5656:   ts->steprestart = PETSC_TRUE;
5657:   return 0;
5658: }

5660: /*@
5661:    TSRollBack - Rolls back one time step

5663:    Collective on TS

5665:    Input Parameter:
5666: .  ts - the TS context obtained from TSCreate()

5668:    Level: advanced

5670: .seealso: TSCreate(), TSSetUp(), TSDestroy(), TSSolve(), TSSetPreStep(), TSSetPreStage(), TSInterpolate()
5671: @*/
5672: PetscErrorCode  TSRollBack(TS ts)
5673: {
5677:   (*ts->ops->rollback)(ts);
5678:   ts->time_step = ts->ptime - ts->ptime_prev;
5679:   ts->ptime = ts->ptime_prev;
5680:   ts->ptime_prev = ts->ptime_prev_rollback;
5681:   ts->steps--;
5682:   ts->steprollback = PETSC_TRUE;
5683:   return 0;
5684: }

5686: /*@
5687:    TSGetStages - Get the number of stages and stage values

5689:    Input Parameter:
5690: .  ts - the TS context obtained from TSCreate()

5692:    Output Parameters:
5693: +  ns - the number of stages
5694: -  Y - the current stage vectors

5696:    Level: advanced

5698:    Notes: Both ns and Y can be NULL.

5700: .seealso: TSCreate()
5701: @*/
5702: PetscErrorCode  TSGetStages(TS ts,PetscInt *ns,Vec **Y)
5703: {
5707:   if (!ts->ops->getstages) {
5708:     if (ns) *ns = 0;
5709:     if (Y) *Y = NULL;
5710:   } else {
5711:     (*ts->ops->getstages)(ts,ns,Y);
5712:   }
5713:   return 0;
5714: }

5716: /*@C
5717:   TSComputeIJacobianDefaultColor - Computes the Jacobian using finite differences and coloring to exploit matrix sparsity.

5719:   Collective on SNES

5721:   Input Parameters:
5722: + ts - the TS context
5723: . t - current timestep
5724: . U - state vector
5725: . Udot - time derivative of state vector
5726: . shift - shift to apply, see note below
5727: - ctx - an optional user context

5729:   Output Parameters:
5730: + J - Jacobian matrix (not altered in this routine)
5731: - B - newly computed Jacobian matrix to use with preconditioner (generally the same as J)

5733:   Level: intermediate

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

5738:   dF/dU + shift*dF/dUdot

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

5743:   This will first try to get the coloring from the DM.  If the DM type has no coloring
5744:   routine, then it will try to get the coloring from the matrix.  This requires that the
5745:   matrix have nonzero entries precomputed.

5747: .seealso: TSSetIJacobian(), MatFDColoringCreate(), MatFDColoringSetFunction()
5748: @*/
5749: PetscErrorCode TSComputeIJacobianDefaultColor(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal shift,Mat J,Mat B,void *ctx)
5750: {
5751:   SNES           snes;
5752:   MatFDColoring  color;
5753:   PetscBool      hascolor, matcolor = PETSC_FALSE;

5755:   PetscOptionsGetBool(((PetscObject)ts)->options,((PetscObject) ts)->prefix, "-ts_fd_color_use_mat", &matcolor, NULL);
5756:   PetscObjectQuery((PetscObject) B, "TSMatFDColoring", (PetscObject *) &color);
5757:   if (!color) {
5758:     DM         dm;
5759:     ISColoring iscoloring;

5761:     TSGetDM(ts, &dm);
5762:     DMHasColoring(dm, &hascolor);
5763:     if (hascolor && !matcolor) {
5764:       DMCreateColoring(dm, IS_COLORING_GLOBAL, &iscoloring);
5765:       MatFDColoringCreate(B, iscoloring, &color);
5766:       MatFDColoringSetFunction(color, (PetscErrorCode (*)(void)) SNESTSFormFunction, (void *) ts);
5767:       MatFDColoringSetFromOptions(color);
5768:       MatFDColoringSetUp(B, iscoloring, color);
5769:       ISColoringDestroy(&iscoloring);
5770:     } else {
5771:       MatColoring mc;

5773:       MatColoringCreate(B, &mc);
5774:       MatColoringSetDistance(mc, 2);
5775:       MatColoringSetType(mc, MATCOLORINGSL);
5776:       MatColoringSetFromOptions(mc);
5777:       MatColoringApply(mc, &iscoloring);
5778:       MatColoringDestroy(&mc);
5779:       MatFDColoringCreate(B, iscoloring, &color);
5780:       MatFDColoringSetFunction(color, (PetscErrorCode (*)(void)) SNESTSFormFunction, (void *) ts);
5781:       MatFDColoringSetFromOptions(color);
5782:       MatFDColoringSetUp(B, iscoloring, color);
5783:       ISColoringDestroy(&iscoloring);
5784:     }
5785:     PetscObjectCompose((PetscObject) B, "TSMatFDColoring", (PetscObject) color);
5786:     PetscObjectDereference((PetscObject) color);
5787:   }
5788:   TSGetSNES(ts, &snes);
5789:   MatFDColoringApply(B, color, U, snes);
5790:   if (J != B) {
5791:     MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);
5792:     MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);
5793:   }
5794:   return 0;
5795: }

5797: /*@
5798:     TSSetFunctionDomainError - Set a function that tests if the current state vector is valid

5800:     Input Parameters:
5801: +    ts - the TS context
5802: -    func - function called within TSFunctionDomainError

5804:     Calling sequence of func:
5805: $     PetscErrorCode func(TS ts,PetscReal time,Vec state,PetscBool reject)

5807: +   ts - the TS context
5808: .   time - the current time (of the stage)
5809: .   state - the state to check if it is valid
5810: -   reject - (output parameter) PETSC_FALSE if the state is acceptable, PETSC_TRUE if not acceptable

5812:     Level: intermediate

5814:     Notes:
5815:       If an implicit ODE solver is being used then, in addition to providing this routine, the
5816:       user's code should call SNESSetFunctionDomainError() when domain errors occur during
5817:       function evaluations where the functions are provided by TSSetIFunction() or TSSetRHSFunction().
5818:       Use TSGetSNES() to obtain the SNES object

5820:     Developer Notes:
5821:       The naming of this function is inconsistent with the SNESSetFunctionDomainError()
5822:       since one takes a function pointer and the other does not.

5824: .seealso: TSAdaptCheckStage(), TSFunctionDomainError(), SNESSetFunctionDomainError(), TSGetSNES()
5825: @*/

5827: PetscErrorCode TSSetFunctionDomainError(TS ts, PetscErrorCode (*func)(TS,PetscReal,Vec,PetscBool*))
5828: {
5830:   ts->functiondomainerror = func;
5831:   return 0;
5832: }

5834: /*@
5835:     TSFunctionDomainError - Checks if the current state is valid

5837:     Input Parameters:
5838: +    ts - the TS context
5839: .    stagetime - time of the simulation
5840: -    Y - state vector to check.

5842:     Output Parameter:
5843: .    accept - Set to PETSC_FALSE if the current state vector is valid.

5845:     Note:
5846:     This function is called by the TS integration routines and calls the user provided function (set with TSSetFunctionDomainError())
5847:     to check if the current state is valid.

5849:     Level: developer

5851: .seealso: TSSetFunctionDomainError()
5852: @*/
5853: PetscErrorCode TSFunctionDomainError(TS ts,PetscReal stagetime,Vec Y,PetscBool* accept)
5854: {
5856:   *accept = PETSC_TRUE;
5857:   if (ts->functiondomainerror) {
5858:     PetscStackCallStandard((*ts->functiondomainerror),ts,stagetime,Y,accept);
5859:   }
5860:   return 0;
5861: }

5863: /*@C
5864:   TSClone - This function clones a time step object.

5866:   Collective

5868:   Input Parameter:
5869: . tsin    - The input TS

5871:   Output Parameter:
5872: . tsout   - The output TS (cloned)

5874:   Notes:
5875:   This function is used to create a clone of a TS object. It is used in ARKIMEX for initializing the slope for first stage explicit methods. It will likely be replaced in the future with a mechanism of switching methods on the fly.

5877:   When using TSDestroy() on a clone the user has to first reset the correct TS reference in the embedded SNES object: e.g.: by running SNES snes_dup=NULL; TSGetSNES(ts,&snes_dup); TSSetSNES(ts,snes_dup);

5879:   Level: developer

5881: .seealso: TSCreate(), TSSetType(), TSSetUp(), TSDestroy(), TSSetProblemType()
5882: @*/
5883: PetscErrorCode  TSClone(TS tsin, TS *tsout)
5884: {
5885:   TS             t;
5886:   SNES           snes_start;
5887:   DM             dm;
5888:   TSType         type;

5891:   *tsout = NULL;

5893:   PetscHeaderCreate(t, TS_CLASSID, "TS", "Time stepping", "TS", PetscObjectComm((PetscObject)tsin), TSDestroy, TSView);

5895:   /* General TS description */
5896:   t->numbermonitors    = 0;
5897:   t->monitorFrequency  = 1;
5898:   t->setupcalled       = 0;
5899:   t->ksp_its           = 0;
5900:   t->snes_its          = 0;
5901:   t->nwork             = 0;
5902:   t->rhsjacobian.time  = PETSC_MIN_REAL;
5903:   t->rhsjacobian.scale = 1.;
5904:   t->ijacobian.shift   = 1.;

5906:   TSGetSNES(tsin,&snes_start);
5907:   TSSetSNES(t,snes_start);

5909:   TSGetDM(tsin,&dm);
5910:   TSSetDM(t,dm);

5912:   t->adapt = tsin->adapt;
5913:   PetscObjectReference((PetscObject)t->adapt);

5915:   t->trajectory = tsin->trajectory;
5916:   PetscObjectReference((PetscObject)t->trajectory);

5918:   t->event = tsin->event;
5919:   if (t->event) t->event->refct++;

5921:   t->problem_type      = tsin->problem_type;
5922:   t->ptime             = tsin->ptime;
5923:   t->ptime_prev        = tsin->ptime_prev;
5924:   t->time_step         = tsin->time_step;
5925:   t->max_time          = tsin->max_time;
5926:   t->steps             = tsin->steps;
5927:   t->max_steps         = tsin->max_steps;
5928:   t->equation_type     = tsin->equation_type;
5929:   t->atol              = tsin->atol;
5930:   t->rtol              = tsin->rtol;
5931:   t->max_snes_failures = tsin->max_snes_failures;
5932:   t->max_reject        = tsin->max_reject;
5933:   t->errorifstepfailed = tsin->errorifstepfailed;

5935:   TSGetType(tsin,&type);
5936:   TSSetType(t,type);

5938:   t->vec_sol           = NULL;

5940:   t->cfltime          = tsin->cfltime;
5941:   t->cfltime_local    = tsin->cfltime_local;
5942:   t->exact_final_time = tsin->exact_final_time;

5944:   PetscMemcpy(t->ops,tsin->ops,sizeof(struct _TSOps));

5946:   if (((PetscObject)tsin)->fortran_func_pointers) {
5947:     PetscInt i;
5948:     PetscMalloc((10)*sizeof(void(*)(void)),&((PetscObject)t)->fortran_func_pointers);
5949:     for (i=0; i<10; i++) {
5950:       ((PetscObject)t)->fortran_func_pointers[i] = ((PetscObject)tsin)->fortran_func_pointers[i];
5951:     }
5952:   }
5953:   *tsout = t;
5954:   return 0;
5955: }

5957: static PetscErrorCode RHSWrapperFunction_TSRHSJacobianTest(void* ctx,Vec x,Vec y)
5958: {
5959:   TS             ts = (TS) ctx;

5961:   TSComputeRHSFunction(ts,0,x,y);
5962:   return 0;
5963: }

5965: /*@
5966:     TSRHSJacobianTest - Compares the multiply routine provided to the MATSHELL with differencing on the TS given RHS function.

5968:    Logically Collective on TS

5970:     Input Parameters:
5971:     TS - the time stepping routine

5973:    Output Parameter:
5974: .   flg - PETSC_TRUE if the multiply is likely correct

5976:    Options Database:
5977:  .   -ts_rhs_jacobian_test_mult -mat_shell_test_mult_view - run the test at each timestep of the integrator

5979:    Level: advanced

5981:    Notes:
5982:     This only works for problems defined only the RHS function and Jacobian NOT IFunction and IJacobian

5984: .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMultTranspose(), TSRHSJacobianTestTranspose()
5985: @*/
5986: PetscErrorCode  TSRHSJacobianTest(TS ts,PetscBool *flg)
5987: {
5988:   Mat            J,B;
5989:   TSRHSJacobian  func;
5990:   void*          ctx;

5992:   TSGetRHSJacobian(ts,&J,&B,&func,&ctx);
5993:   (*func)(ts,0.0,ts->vec_sol,J,B,ctx);
5994:   MatShellTestMult(J,RHSWrapperFunction_TSRHSJacobianTest,ts->vec_sol,ts,flg);
5995:   return 0;
5996: }

5998: /*@C
5999:     TSRHSJacobianTestTranspose - Compares the multiply transpose routine provided to the MATSHELL with differencing on the TS given RHS function.

6001:    Logically Collective on TS

6003:     Input Parameters:
6004:     TS - the time stepping routine

6006:    Output Parameter:
6007: .   flg - PETSC_TRUE if the multiply is likely correct

6009:    Options Database:
6010: .   -ts_rhs_jacobian_test_mult_transpose -mat_shell_test_mult_transpose_view - run the test at each timestep of the integrator

6012:    Notes:
6013:     This only works for problems defined only the RHS function and Jacobian NOT IFunction and IJacobian

6015:    Level: advanced

6017: .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMultTranspose(), TSRHSJacobianTest()
6018: @*/
6019: PetscErrorCode  TSRHSJacobianTestTranspose(TS ts,PetscBool *flg)
6020: {
6021:   Mat            J,B;
6022:   void           *ctx;
6023:   TSRHSJacobian  func;

6025:   TSGetRHSJacobian(ts,&J,&B,&func,&ctx);
6026:   (*func)(ts,0.0,ts->vec_sol,J,B,ctx);
6027:   MatShellTestMultTranspose(J,RHSWrapperFunction_TSRHSJacobianTest,ts->vec_sol,ts,flg);
6028:   return 0;
6029: }

6031: /*@
6032:   TSSetUseSplitRHSFunction - Use the split RHSFunction when a multirate method is used.

6034:   Logically collective

6036:   Input Parameters:
6037: +  ts - timestepping context
6038: -  use_splitrhsfunction - PETSC_TRUE indicates that the split RHSFunction will be used

6040:   Options Database:
6041: .   -ts_use_splitrhsfunction - <true,false>

6043:   Notes:
6044:     This is only useful for multirate methods

6046:   Level: intermediate

6048: .seealso: TSGetUseSplitRHSFunction()
6049: @*/
6050: PetscErrorCode TSSetUseSplitRHSFunction(TS ts, PetscBool use_splitrhsfunction)
6051: {
6053:   ts->use_splitrhsfunction = use_splitrhsfunction;
6054:   return 0;
6055: }

6057: /*@
6058:   TSGetUseSplitRHSFunction - Gets whether to use the split RHSFunction when a multirate method is used.

6060:   Not collective

6062:   Input Parameter:
6063: .  ts - timestepping context

6065:   Output Parameter:
6066: .  use_splitrhsfunction - PETSC_TRUE indicates that the split RHSFunction will be used

6068:   Level: intermediate

6070: .seealso: TSSetUseSplitRHSFunction()
6071: @*/
6072: PetscErrorCode TSGetUseSplitRHSFunction(TS ts, PetscBool *use_splitrhsfunction)
6073: {
6075:   *use_splitrhsfunction = ts->use_splitrhsfunction;
6076:   return 0;
6077: }

6079: /*@
6080:     TSSetMatStructure - sets the relationship between the nonzero structure of the RHS Jacobian matrix to the IJacobian matrix.

6082:    Logically  Collective on ts

6084:    Input Parameters:
6085: +  ts - the time-stepper
6086: -  str - the structure (the default is UNKNOWN_NONZERO_PATTERN)

6088:    Level: intermediate

6090:    Notes:
6091:      When the relationship between the nonzero structures is known and supplied the solution process can be much faster

6093: .seealso: MatAXPY(), MatStructure
6094:  @*/
6095: PetscErrorCode TSSetMatStructure(TS ts,MatStructure str)
6096: {
6098:   ts->axpy_pattern = str;
6099:   return 0;
6100: }