xref: /petsc/src/ts/tutorials/ex3.c (revision d3c69ad0bb09681fda0ed612cc8275184ea14744)
1 
2 static char help[] ="Solves a simple time-dependent linear PDE (the heat equation).\n\
3 Input parameters include:\n\
4   -m <points>, where <points> = number of grid points\n\
5   -time_dependent_rhs : Treat the problem as having a time-dependent right-hand side\n\
6   -use_ifunc          : Use IFunction/IJacobian interface\n\
7   -debug              : Activate debugging printouts\n\
8   -nox                : Deactivate x-window graphics\n\n";
9 
10 /* ------------------------------------------------------------------------
11 
12    This program solves the one-dimensional heat equation (also called the
13    diffusion equation),
14        u_t = u_xx,
15    on the domain 0 <= x <= 1, with the boundary conditions
16        u(t,0) = 0, u(t,1) = 0,
17    and the initial condition
18        u(0,x) = sin(6*pi*x) + 3*sin(2*pi*x).
19    This is a linear, second-order, parabolic equation.
20 
21    We discretize the right-hand side using finite differences with
22    uniform grid spacing h:
23        u_xx = (u_{i+1} - 2u_{i} + u_{i-1})/(h^2)
24    We then demonstrate time evolution using the various TS methods by
25    running the program via
26        ex3 -ts_type <timestepping solver>
27 
28    We compare the approximate solution with the exact solution, given by
29        u_exact(x,t) = exp(-36*pi*pi*t) * sin(6*pi*x) +
30                       3*exp(-4*pi*pi*t) * sin(2*pi*x)
31 
32    Notes:
33    This code demonstrates the TS solver interface to two variants of
34    linear problems, u_t = f(u,t), namely
35      - time-dependent f:   f(u,t) is a function of t
36      - time-independent f: f(u,t) is simply f(u)
37 
38     The parallel version of this code is ts/tutorials/ex4.c
39 
40   ------------------------------------------------------------------------- */
41 
42 /*
43    Include "petscts.h" so that we can use TS solvers.  Note that this file
44    automatically includes:
45      petscsys.h       - base PETSc routines   petscvec.h  - vectors
46      petscmat.h  - matrices
47      petscis.h     - index sets            petscksp.h  - Krylov subspace methods
48      petscviewer.h - viewers               petscpc.h   - preconditioners
49      petscksp.h   - linear solvers        petscsnes.h - nonlinear solvers
50 */
51 
52 #include <petscts.h>
53 #include <petscdraw.h>
54 
55 /*
56    User-defined application context - contains data needed by the
57    application-provided call-back routines.
58 */
59 typedef struct {
60   Vec         solution;          /* global exact solution vector */
61   PetscInt    m;                 /* total number of grid points */
62   PetscReal   h;                 /* mesh width h = 1/(m-1) */
63   PetscBool   debug;             /* flag (1 indicates activation of debugging printouts) */
64   PetscViewer viewer1,viewer2;   /* viewers for the solution and error */
65   PetscReal   norm_2,norm_max;   /* error norms */
66   Mat         A;                 /* RHS mat, used with IFunction interface */
67   PetscReal   oshift;            /* old shift applied, prevent to recompute the IJacobian */
68 } AppCtx;
69 
70 /*
71    User-defined routines
72 */
73 extern PetscErrorCode InitialConditions(Vec,AppCtx*);
74 extern PetscErrorCode RHSMatrixHeat(TS,PetscReal,Vec,Mat,Mat,void*);
75 extern PetscErrorCode IFunctionHeat(TS,PetscReal,Vec,Vec,Vec,void*);
76 extern PetscErrorCode IJacobianHeat(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
77 extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*);
78 extern PetscErrorCode ExactSolution(PetscReal,Vec,AppCtx*);
79 
80 int main(int argc,char **argv)
81 {
82   AppCtx         appctx;                 /* user-defined application context */
83   TS             ts;                     /* timestepping context */
84   Mat            A;                      /* matrix data structure */
85   Vec            u;                      /* approximate solution vector */
86   PetscReal      time_total_max = 100.0; /* default max total time */
87   PetscInt       time_steps_max = 100;   /* default max timesteps */
88   PetscDraw      draw;                   /* drawing context */
89   PetscInt       steps,m;
90   PetscMPIInt    size;
91   PetscReal      dt;
92   PetscBool      flg,flg_string;
93 
94   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
95      Initialize program and set problem parameters
96      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
97 
98   PetscFunctionBeginUser;
99   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
100   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
101   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
102 
103   m    = 60;
104   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
105   PetscCall(PetscOptionsHasName(NULL,NULL,"-debug",&appctx.debug));
106   flg_string = PETSC_FALSE;
107   PetscCall(PetscOptionsGetBool(NULL,NULL,"-test_string_viewer",&flg_string,NULL));
108 
109   appctx.m        = m;
110   appctx.h        = 1.0/(m-1.0);
111   appctx.norm_2   = 0.0;
112   appctx.norm_max = 0.0;
113 
114   PetscCall(PetscPrintf(PETSC_COMM_SELF,"Solving a linear TS problem on 1 processor\n"));
115 
116   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
117      Create vector data structures
118      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
119 
120   /*
121      Create vector data structures for approximate and exact solutions
122   */
123   PetscCall(VecCreateSeq(PETSC_COMM_SELF,m,&u));
124   PetscCall(VecDuplicate(u,&appctx.solution));
125 
126   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127      Set up displays to show graphs of the solution and error
128      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
129 
130   PetscCall(PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,380,400,160,&appctx.viewer1));
131   PetscCall(PetscViewerDrawGetDraw(appctx.viewer1,0,&draw));
132   PetscCall(PetscDrawSetDoubleBuffer(draw));
133   PetscCall(PetscViewerDrawOpen(PETSC_COMM_SELF,0,"",80,0,400,160,&appctx.viewer2));
134   PetscCall(PetscViewerDrawGetDraw(appctx.viewer2,0,&draw));
135   PetscCall(PetscDrawSetDoubleBuffer(draw));
136 
137   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138      Create timestepping solver context
139      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
140 
141   PetscCall(TSCreate(PETSC_COMM_SELF,&ts));
142   PetscCall(TSSetProblemType(ts,TS_LINEAR));
143 
144   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
145      Set optional user-defined monitoring routine
146      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
147 
148   if (!flg_string) {
149     PetscCall(TSMonitorSet(ts,Monitor,&appctx,NULL));
150   }
151 
152   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
153 
154      Create matrix data structure; set matrix evaluation routine.
155      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
156 
157   PetscCall(MatCreate(PETSC_COMM_SELF,&A));
158   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,m));
159   PetscCall(MatSetFromOptions(A));
160   PetscCall(MatSetUp(A));
161 
162   flg  = PETSC_FALSE;
163   PetscCall(PetscOptionsGetBool(NULL,NULL,"-use_ifunc",&flg,NULL));
164   if (!flg) {
165     appctx.A = NULL;
166     PetscCall(PetscOptionsGetBool(NULL,NULL,"-time_dependent_rhs",&flg,NULL));
167     if (flg) {
168       /*
169          For linear problems with a time-dependent f(u,t) in the equation
170          u_t = f(u,t), the user provides the discretized right-hand-side
171          as a time-dependent matrix.
172       */
173       PetscCall(TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx));
174       PetscCall(TSSetRHSJacobian(ts,A,A,RHSMatrixHeat,&appctx));
175     } else {
176       /*
177          For linear problems with a time-independent f(u) in the equation
178          u_t = f(u), the user provides the discretized right-hand-side
179          as a matrix only once, and then sets the special Jacobian evaluation
180          routine TSComputeRHSJacobianConstant() which will NOT recompute the Jacobian.
181       */
182       PetscCall(RHSMatrixHeat(ts,0.0,u,A,A,&appctx));
183       PetscCall(TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx));
184       PetscCall(TSSetRHSJacobian(ts,A,A,TSComputeRHSJacobianConstant,&appctx));
185     }
186   } else {
187     Mat J;
188 
189     PetscCall(RHSMatrixHeat(ts,0.0,u,A,A,&appctx));
190     PetscCall(MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&J));
191     PetscCall(TSSetIFunction(ts,NULL,IFunctionHeat,&appctx));
192     PetscCall(TSSetIJacobian(ts,J,J,IJacobianHeat,&appctx));
193     PetscCall(MatDestroy(&J));
194 
195     PetscCall(PetscObjectReference((PetscObject)A));
196     appctx.A = A;
197     appctx.oshift = PETSC_MIN_REAL;
198   }
199   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
200      Set solution vector and initial timestep
201      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
202 
203   dt   = appctx.h*appctx.h/2.0;
204   PetscCall(TSSetTimeStep(ts,dt));
205 
206   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
207      Customize timestepping solver:
208        - Set the solution method to be the Backward Euler method.
209        - Set timestepping duration info
210      Then set runtime options, which can override these defaults.
211      For example,
212           -ts_max_steps <maxsteps> -ts_max_time <maxtime>
213      to override the defaults set by TSSetMaxSteps()/TSSetMaxTime().
214      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
215 
216   PetscCall(TSSetMaxSteps(ts,time_steps_max));
217   PetscCall(TSSetMaxTime(ts,time_total_max));
218   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
219   PetscCall(TSSetFromOptions(ts));
220 
221   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
222      Solve the problem
223      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
224 
225   /*
226      Evaluate initial conditions
227   */
228   PetscCall(InitialConditions(u,&appctx));
229 
230   /*
231      Run the timestepping solver
232   */
233   PetscCall(TSSolve(ts,u));
234   PetscCall(TSGetStepNumber(ts,&steps));
235 
236   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
237      View timestepping solver info
238      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
239 
240   PetscCall(PetscPrintf(PETSC_COMM_SELF,"avg. error (2 norm) = %g, avg. error (max norm) = %g\n",(double)(appctx.norm_2/steps),(double)(appctx.norm_max/steps)));
241   if (!flg_string) {
242     PetscCall(TSView(ts,PETSC_VIEWER_STDOUT_SELF));
243   } else {
244     PetscViewer stringviewer;
245     char        string[512];
246     const char  *outstring;
247 
248     PetscCall(PetscViewerStringOpen(PETSC_COMM_WORLD,string,sizeof(string),&stringviewer));
249     PetscCall(TSView(ts,stringviewer));
250     PetscCall(PetscViewerStringGetStringRead(stringviewer,&outstring,NULL));
251     PetscCheck((char*)outstring == (char*)string,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"String returned from viewer does not equal original string");
252     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Output from string viewer:%s\n",outstring));
253     PetscCall(PetscViewerDestroy(&stringviewer));
254   }
255 
256   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
257      Free work space.  All PETSc objects should be destroyed when they
258      are no longer needed.
259      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
260 
261   PetscCall(TSDestroy(&ts));
262   PetscCall(MatDestroy(&A));
263   PetscCall(VecDestroy(&u));
264   PetscCall(PetscViewerDestroy(&appctx.viewer1));
265   PetscCall(PetscViewerDestroy(&appctx.viewer2));
266   PetscCall(VecDestroy(&appctx.solution));
267   PetscCall(MatDestroy(&appctx.A));
268 
269   /*
270      Always call PetscFinalize() before exiting a program.  This routine
271        - finalizes the PETSc libraries as well as MPI
272        - provides summary and diagnostic information if certain runtime
273          options are chosen (e.g., -log_view).
274   */
275   PetscCall(PetscFinalize());
276   return 0;
277 }
278 /* --------------------------------------------------------------------- */
279 /*
280    InitialConditions - Computes the solution at the initial time.
281 
282    Input Parameter:
283    u - uninitialized solution vector (global)
284    appctx - user-defined application context
285 
286    Output Parameter:
287    u - vector with solution at initial time (global)
288 */
289 PetscErrorCode InitialConditions(Vec u,AppCtx *appctx)
290 {
291   PetscScalar    *u_localptr,h = appctx->h;
292   PetscInt       i;
293 
294   /*
295     Get a pointer to vector data.
296     - For default PETSc vectors, VecGetArray() returns a pointer to
297       the data array.  Otherwise, the routine is implementation dependent.
298     - You MUST call VecRestoreArray() when you no longer need access to
299       the array.
300     - Note that the Fortran interface to VecGetArray() differs from the
301       C version.  See the users manual for details.
302   */
303   PetscCall(VecGetArrayWrite(u,&u_localptr));
304 
305   /*
306      We initialize the solution array by simply writing the solution
307      directly into the array locations.  Alternatively, we could use
308      VecSetValues() or VecSetValuesLocal().
309   */
310   for (i=0; i<appctx->m; i++) u_localptr[i] = PetscSinScalar(PETSC_PI*i*6.*h) + 3.*PetscSinScalar(PETSC_PI*i*2.*h);
311 
312   /*
313      Restore vector
314   */
315   PetscCall(VecRestoreArrayWrite(u,&u_localptr));
316 
317   /*
318      Print debugging information if desired
319   */
320   if (appctx->debug) {
321     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Initial guess vector\n"));
322     PetscCall(VecView(u,PETSC_VIEWER_STDOUT_SELF));
323   }
324 
325   return 0;
326 }
327 /* --------------------------------------------------------------------- */
328 /*
329    ExactSolution - Computes the exact solution at a given time.
330 
331    Input Parameters:
332    t - current time
333    solution - vector in which exact solution will be computed
334    appctx - user-defined application context
335 
336    Output Parameter:
337    solution - vector with the newly computed exact solution
338 */
339 PetscErrorCode ExactSolution(PetscReal t,Vec solution,AppCtx *appctx)
340 {
341   PetscScalar    *s_localptr,h = appctx->h,ex1,ex2,sc1,sc2,tc = t;
342   PetscInt       i;
343 
344   /*
345      Get a pointer to vector data.
346   */
347   PetscCall(VecGetArrayWrite(solution,&s_localptr));
348 
349   /*
350      Simply write the solution directly into the array locations.
351      Alternatively, we culd use VecSetValues() or VecSetValuesLocal().
352   */
353   ex1 = PetscExpScalar(-36.*PETSC_PI*PETSC_PI*tc);
354   ex2 = PetscExpScalar(-4.*PETSC_PI*PETSC_PI*tc);
355   sc1 = PETSC_PI*6.*h;                 sc2 = PETSC_PI*2.*h;
356   for (i=0; i<appctx->m; i++) s_localptr[i] = PetscSinScalar(sc1*(PetscReal)i)*ex1 + 3.*PetscSinScalar(sc2*(PetscReal)i)*ex2;
357 
358   /*
359      Restore vector
360   */
361   PetscCall(VecRestoreArrayWrite(solution,&s_localptr));
362   return 0;
363 }
364 /* --------------------------------------------------------------------- */
365 /*
366    Monitor - User-provided routine to monitor the solution computed at
367    each timestep.  This example plots the solution and computes the
368    error in two different norms.
369 
370    This example also demonstrates changing the timestep via TSSetTimeStep().
371 
372    Input Parameters:
373    ts     - the timestep context
374    step   - the count of the current step (with 0 meaning the
375              initial condition)
376    time   - the current time
377    u      - the solution at this timestep
378    ctx    - the user-provided context for this monitoring routine.
379             In this case we use the application context which contains
380             information about the problem size, workspace and the exact
381             solution.
382 */
383 PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec u,void *ctx)
384 {
385   AppCtx         *appctx = (AppCtx*) ctx;   /* user-defined application context */
386   PetscReal      norm_2,norm_max,dt,dttol;
387 
388   /*
389      View a graph of the current iterate
390   */
391   PetscCall(VecView(u,appctx->viewer2));
392 
393   /*
394      Compute the exact solution
395   */
396   PetscCall(ExactSolution(time,appctx->solution,appctx));
397 
398   /*
399      Print debugging information if desired
400   */
401   if (appctx->debug) {
402     PetscCall(PetscPrintf(PETSC_COMM_SELF,"Computed solution vector\n"));
403     PetscCall(VecView(u,PETSC_VIEWER_STDOUT_SELF));
404     PetscCall(PetscPrintf(PETSC_COMM_SELF,"Exact solution vector\n"));
405     PetscCall(VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF));
406   }
407 
408   /*
409      Compute the 2-norm and max-norm of the error
410   */
411   PetscCall(VecAXPY(appctx->solution,-1.0,u));
412   PetscCall(VecNorm(appctx->solution,NORM_2,&norm_2));
413   norm_2 = PetscSqrtReal(appctx->h)*norm_2;
414   PetscCall(VecNorm(appctx->solution,NORM_MAX,&norm_max));
415 
416   PetscCall(TSGetTimeStep(ts,&dt));
417   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Timestep %3" PetscInt_FMT ": step size = %g, time = %g, 2-norm error = %g, max norm error = %g\n",step,(double)dt,(double)time,(double)norm_2,(double)norm_max));
418 
419   appctx->norm_2   += norm_2;
420   appctx->norm_max += norm_max;
421 
422   dttol = .0001;
423   PetscCall(PetscOptionsGetReal(NULL,NULL,"-dttol",&dttol,NULL));
424   if (dt < dttol) {
425     dt  *= .999;
426     PetscCall(TSSetTimeStep(ts,dt));
427   }
428 
429   /*
430      View a graph of the error
431   */
432   PetscCall(VecView(appctx->solution,appctx->viewer1));
433 
434   /*
435      Print debugging information if desired
436   */
437   if (appctx->debug) {
438     PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error vector\n"));
439     PetscCall(VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF));
440   }
441 
442   return 0;
443 }
444 /* --------------------------------------------------------------------- */
445 /*
446    RHSMatrixHeat - User-provided routine to compute the right-hand-side
447    matrix for the heat equation.
448 
449    Input Parameters:
450    ts - the TS context
451    t - current time
452    global_in - global input vector
453    dummy - optional user-defined context, as set by TSetRHSJacobian()
454 
455    Output Parameters:
456    AA - Jacobian matrix
457    BB - optionally different preconditioning matrix
458    str - flag indicating matrix structure
459 
460    Notes:
461    Recall that MatSetValues() uses 0-based row and column numbers
462    in Fortran as well as in C.
463 */
464 PetscErrorCode RHSMatrixHeat(TS ts,PetscReal t,Vec X,Mat AA,Mat BB,void *ctx)
465 {
466   Mat            A       = AA;                /* Jacobian matrix */
467   AppCtx         *appctx = (AppCtx*)ctx;     /* user-defined application context */
468   PetscInt       mstart  = 0;
469   PetscInt       mend    = appctx->m;
470   PetscInt       i,idx[3];
471   PetscScalar    v[3],stwo = -2./(appctx->h*appctx->h),sone = -.5*stwo;
472 
473   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
474      Compute entries for the locally owned part of the matrix
475      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
476   /*
477      Set matrix rows corresponding to boundary data
478   */
479 
480   mstart = 0;
481   v[0]   = 1.0;
482   PetscCall(MatSetValues(A,1,&mstart,1,&mstart,v,INSERT_VALUES));
483   mstart++;
484 
485   mend--;
486   v[0] = 1.0;
487   PetscCall(MatSetValues(A,1,&mend,1,&mend,v,INSERT_VALUES));
488 
489   /*
490      Set matrix rows corresponding to interior data.  We construct the
491      matrix one row at a time.
492   */
493   v[0] = sone; v[1] = stwo; v[2] = sone;
494   for (i=mstart; i<mend; i++) {
495     idx[0] = i-1; idx[1] = i; idx[2] = i+1;
496     PetscCall(MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES));
497   }
498 
499   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
500      Complete the matrix assembly process and set some options
501      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
502   /*
503      Assemble matrix, using the 2-step process:
504        MatAssemblyBegin(), MatAssemblyEnd()
505      Computations can be done while messages are in transition
506      by placing code between these two statements.
507   */
508   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
509   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
510 
511   /*
512      Set and option to indicate that we will never add a new nonzero location
513      to the matrix. If we do, it will generate an error.
514   */
515   PetscCall(MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE));
516 
517   return 0;
518 }
519 
520 PetscErrorCode IFunctionHeat(TS ts,PetscReal t,Vec X,Vec Xdot,Vec r,void *ctx)
521 {
522   AppCtx         *appctx = (AppCtx*)ctx;     /* user-defined application context */
523 
524   PetscCall(MatMult(appctx->A,X,r));
525   PetscCall(VecAYPX(r,-1.0,Xdot));
526   return 0;
527 }
528 
529 PetscErrorCode IJacobianHeat(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal s,Mat A,Mat B,void *ctx)
530 {
531   AppCtx         *appctx = (AppCtx*)ctx;     /* user-defined application context */
532 
533   if (appctx->oshift == s) return 0;
534   PetscCall(MatCopy(appctx->A,A,SAME_NONZERO_PATTERN));
535   PetscCall(MatScale(A,-1));
536   PetscCall(MatShift(A,s));
537   PetscCall(MatCopy(A,B,SAME_NONZERO_PATTERN));
538   appctx->oshift = s;
539   return 0;
540 }
541 
542 /*TEST
543 
544     test:
545       args: -nox -ts_type ssp -ts_dt 0.0005
546 
547     test:
548       suffix: 2
549       args: -nox -ts_type ssp -ts_dt 0.0005 -time_dependent_rhs 1
550 
551     test:
552       suffix: 3
553       args:  -nox -ts_type rosw -ts_max_steps 3 -ksp_converged_reason
554       filter: sed "s/ATOL/RTOL/g"
555       requires: !single
556 
557     test:
558       suffix: 4
559       args: -nox -ts_type beuler -ts_max_steps 3 -ksp_converged_reason
560       filter: sed "s/ATOL/RTOL/g"
561 
562     test:
563       suffix: 5
564       args: -nox -ts_type beuler -ts_max_steps 3 -ksp_converged_reason -time_dependent_rhs
565       filter: sed "s/ATOL/RTOL/g"
566 
567     test:
568       requires: !single
569       suffix: pod_guess
570       args: -nox -ts_type beuler -use_ifunc -ts_dt 0.0005 -ksp_guess_type pod -pc_type none -ksp_converged_reason
571 
572     test:
573       requires: !single
574       suffix: pod_guess_Ainner
575       args: -nox -ts_type beuler -use_ifunc -ts_dt 0.0005 -ksp_guess_type pod -ksp_guess_pod_Ainner -pc_type none -ksp_converged_reason
576 
577     test:
578       requires: !single
579       suffix: fischer_guess
580       args: -nox -ts_type beuler -use_ifunc -ts_dt 0.0005 -ksp_guess_type fischer -pc_type none -ksp_converged_reason
581 
582     test:
583       requires: !single
584       suffix: fischer_guess_2
585       args: -nox -ts_type beuler -use_ifunc -ts_dt 0.0005 -ksp_guess_type fischer -ksp_guess_fischer_model 2,10 -pc_type none -ksp_converged_reason
586 
587     test:
588       requires: !single
589       suffix: fischer_guess_3
590       args: -nox -ts_type beuler -use_ifunc -ts_dt 0.0005 -ksp_guess_type fischer -ksp_guess_fischer_model 3,10 -pc_type none -ksp_converged_reason
591 
592     test:
593       requires: !single
594       suffix: stringview
595       args: -nox -ts_type rosw -test_string_viewer
596 
597     test:
598       requires: !single
599       suffix: stringview_euler
600       args: -nox -ts_type euler -test_string_viewer
601 
602 TEST*/
603