xref: /petsc/src/ts/interface/ts.c (revision d763cef207eb7fa01b459cfa1b5fffdd1cf73419)
1 #ifdef PETSC_RCS_HEADER
2 
3 #endif
4 
5 #include "src/ts/tsimpl.h"        /*I "ts.h"  I*/
6 
7 #undef __FUNC__
8 #define __FUNC__ "TSComputeRHSFunction"
9 /*
10    TSComputeRHSFunction - Evaluates the right-hand-side function.
11 
12    Note: If the user did not provide a function but merely a matrix,
13    this routine applies the matrix.
14 */
15 int TSComputeRHSFunction(TS ts,double t,Vec x, Vec y)
16 {
17   int ierr;
18 
19   PetscFunctionBegin;
20   PetscValidHeaderSpecific(ts,TS_COOKIE);
21   PetscValidHeader(x);  PetscValidHeader(y);
22 
23   if (ts->rhsfunction) {
24     PetscStackPush("TS user right-hand-side function");
25     ierr = (*ts->rhsfunction)(ts,t,x,y,ts->funP);CHKERRQ(ierr);
26     PetscStackPop;
27     PetscFunctionReturn(0);
28   }
29 
30   if (ts->rhsmatrix) { /* assemble matrix for this timestep */
31     MatStructure flg;
32     PetscStackPush("TS user right-hand-side matrix function");
33     ierr = (*ts->rhsmatrix)(ts,t,&ts->A,&ts->B,&flg,ts->jacP); CHKERRQ(ierr);
34     PetscStackPop;
35   }
36   ierr = MatMult(ts->A,x,y); CHKERRQ(ierr);
37 
38   /* apply user-provided boundary conditions (only needed if these are time dependent) */
39   ierr = TSComputeRHSBoundaryConditions(ts,t,y); CHKERRQ(ierr);
40 
41   PetscFunctionReturn(0);
42 }
43 
44 #undef __FUNC__
45 #define __FUNC__ "TSSetRHSFunction"
46 /*@C
47     TSSetRHSFunction - Sets the routine for evaluating the function,
48     F(t,u), where U_t = F(t,u).
49 
50     Collective on TS
51 
52     Input Parameters:
53 +   ts - the TS context obtained from TSCreate()
54 .   f - routine for evaluating the right-hand-side function
55 -   ctx - [optional] user-defined context for private data for the
56           function evaluation routine (may be PETSC_NULL)
57 
58     Calling sequence of func:
59 $     func (TS ts,double t,Vec u,Vec F,void *ctx);
60 
61 +   t - current timestep
62 .   u - input vector
63 .   F - function vector
64 -   ctx - [optional] user-defined function context
65 
66     Important:
67     The user MUST call either this routine or TSSetRHSMatrix().
68 
69     Level: beginner
70 
71 .keywords: TS, timestep, set, right-hand-side, function
72 
73 .seealso: TSSetRHSMatrix()
74 @*/
75 int TSSetRHSFunction(TS ts,int (*f)(TS,double,Vec,Vec,void*),void *ctx)
76 {
77   PetscFunctionBegin;
78 
79   PetscValidHeaderSpecific(ts,TS_COOKIE);
80   if (ts->problem_type == TS_LINEAR) {
81     SETERRQ(PETSC_ERR_ARG_WRONG,0,"Cannot set function for linear problem");
82   }
83   ts->rhsfunction = f;
84   ts->funP        = ctx;
85   PetscFunctionReturn(0);
86 }
87 
88 #undef __FUNC__
89 #define __FUNC__ "TSSetRHSMatrix"
90 /*@C
91    TSSetRHSMatrix - Sets the function to compute the matrix A, where U_t = A(t) U.
92    Also sets the location to store A.
93 
94    Collective on TS
95 
96    Input Parameters:
97 +  ts  - the TS context obtained from TSCreate()
98 .  A   - matrix
99 .  B   - preconditioner matrix (usually same as A)
100 .  f   - the matrix evaluation routine; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran)
101          if A is not a function of t.
102 -  ctx - [optional] user-defined context for private data for the
103           matrix evaluation routine (may be PETSC_NULL)
104 
105    Calling sequence of func:
106 $     func (TS ts,double t,Mat *A,Mat *B,int *flag,void *ctx);
107 
108 +  t - current timestep
109 .  A - matrix A, where U_t = A(t) U
110 .  B - preconditioner matrix, usually the same as A
111 .  flag - flag indicating information about the preconditioner matrix
112           structure (same as flag in SLESSetOperators())
113 -  ctx - [optional] user-defined context for matrix evaluation routine
114 
115    Notes:
116    See SLESSetOperators() for important information about setting the flag
117    output parameter in the routine func().  Be sure to read this information!
118 
119    The routine func() takes Mat * as the matrix arguments rather than Mat.
120    This allows the matrix evaluation routine to replace A and/or B with a
121    completely new new matrix structure (not just different matrix elements)
122    when appropriate, for instance, if the nonzero structure is changing
123    throughout the global iterations.
124 
125    Important:
126    The user MUST call either this routine or TSSetRHSFunction().
127 
128    Level: beginner
129 
130 .keywords: TS, timestep, set, right-hand-side, matrix
131 
132 .seealso: TSSetRHSFunction()
133 @*/
134 int TSSetRHSMatrix(TS ts,Mat A, Mat B,int (*f)(TS,double,Mat*,Mat*,MatStructure*,void*),void *ctx)
135 {
136   PetscFunctionBegin;
137   PetscValidHeaderSpecific(ts,TS_COOKIE);
138   if (ts->problem_type == TS_NONLINEAR) {
139     SETERRQ(PETSC_ERR_ARG_WRONG,0,"Not for nonlinear problems; use TSSetRHSJacobian()");
140   }
141 
142   ts->rhsmatrix = f;
143   ts->jacP      = ctx;
144   ts->A         = A;
145   ts->B         = B;
146 
147   PetscFunctionReturn(0);
148 }
149 
150 #undef __FUNC__
151 #define __FUNC__ "TSSetRHSJacobian"
152 /*@C
153    TSSetRHSJacobian - Sets the function to compute the Jacobian of F,
154    where U_t = F(U,t), as well as the location to store the matrix.
155 
156    Collective on TS
157 
158    Input Parameters:
159 +  ts  - the TS context obtained from TSCreate()
160 .  A   - Jacobian matrix
161 .  B   - preconditioner matrix (usually same as A)
162 .  f   - the Jacobian evaluation routine
163 -  ctx - [optional] user-defined context for private data for the
164          Jacobian evaluation routine (may be PETSC_NULL)
165 
166    Calling sequence of func:
167 $     func (TS ts,double t,Vec u,Mat *A,Mat *B,int *flag,void *ctx);
168 
169 +  t - current timestep
170 .  u - input vector
171 .  A - matrix A, where U_t = A(t)u
172 .  B - preconditioner matrix, usually the same as A
173 .  flag - flag indicating information about the preconditioner matrix
174           structure (same as flag in SLESSetOperators())
175 -  ctx - [optional] user-defined context for matrix evaluation routine
176 
177    Notes:
178    See SLESSetOperators() for important information about setting the flag
179    output parameter in the routine func().  Be sure to read this information!
180 
181    The routine func() takes Mat * as the matrix arguments rather than Mat.
182    This allows the matrix evaluation routine to replace A and/or B with a
183    completely new new matrix structure (not just different matrix elements)
184    when appropriate, for instance, if the nonzero structure is changing
185    throughout the global iterations.
186 
187    Level: beginner
188 
189 .keywords: TS, timestep, set, right-hand-side, Jacobian
190 
191 .seealso: TSDefaultComputeJacobianColor(),
192           SNESDefaultComputeJacobianColor()
193 
194 @*/
195 int TSSetRHSJacobian(TS ts,Mat A, Mat B,int (*f)(TS,double,Vec,Mat*,Mat*,
196                      MatStructure*,void*),void *ctx)
197 {
198   PetscFunctionBegin;
199   PetscValidHeaderSpecific(ts,TS_COOKIE);
200   if (ts->problem_type != TS_NONLINEAR) {
201     SETERRQ(PETSC_ERR_ARG_WRONG,0,"Not for linear problems; use TSSetRHSMatrix()");
202   }
203 
204   ts->rhsjacobian = f;
205   ts->jacP        = ctx;
206   ts->A           = A;
207   ts->B           = B;
208   PetscFunctionReturn(0);
209 }
210 
211 #undef __FUNC__
212 #define __FUNC__ "TSComputeRHSBoundaryConditions"
213 /*
214    TSComputeRHSBoundaryConditions - Evaluates the boundary condition function.
215 
216    Note: If the user did not provide a function but merely a matrix,
217    this routine applies the matrix.
218 */
219 int TSComputeRHSBoundaryConditions(TS ts,double t,Vec x)
220 {
221   int ierr;
222 
223   PetscFunctionBegin;
224   PetscValidHeaderSpecific(ts,TS_COOKIE);
225   PetscValidHeader(x);
226 
227   if (ts->rhsbc) {
228     PetscStackPush("TS user boundary condition function");
229     ierr = (*ts->rhsbc)(ts,t,x,ts->bcP);CHKERRQ(ierr);
230     PetscStackPop;
231     PetscFunctionReturn(0);
232   }
233 
234   PetscFunctionReturn(0);
235 }
236 
237 #undef __FUNC__
238 #define __FUNC__ "TSSetRHSBoundaryConditions"
239 /*@C
240     TSSetRHSBoundaryConditions - Sets the routine for evaluating the function,
241     boundary conditions for the function F.
242 
243     Collective on TS
244 
245     Input Parameters:
246 +   ts - the TS context obtained from TSCreate()
247 .   f - routine for evaluating the boundary condition function
248 -   ctx - [optional] user-defined context for private data for the
249           function evaluation routine (may be PETSC_NULL)
250 
251     Calling sequence of func:
252 $     func (TS ts,double t,Vec F,void *ctx);
253 
254 +   t - current timestep
255 .   F - function vector
256 -   ctx - [optional] user-defined function context
257 
258     Level: intermediate
259 
260 .keywords: TS, timestep, set, boundary conditions, function
261 @*/
262 int TSSetRHSBoundaryConditions(TS ts,int (*f)(TS,double,Vec,void*),void *ctx)
263 {
264   PetscFunctionBegin;
265 
266   PetscValidHeaderSpecific(ts,TS_COOKIE);
267   if (ts->problem_type != TS_LINEAR) {
268     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For linear problems only");
269   }
270   ts->rhsbc = f;
271   ts->bcP   = ctx;
272   PetscFunctionReturn(0);
273 }
274 
275 #undef __FUNC__
276 #define __FUNC__ "TSView"
277 /*@
278     TSView - Prints the TS data structure.
279 
280     Collective on TS, unless Viewer is VIEWER_STDOUT_SELF
281 
282     Input Parameters:
283 +   ts - the TS context obtained from TSCreate()
284 -   viewer - visualization context
285 
286     Options Database Key:
287 .   -ts_view - calls TSView() at end of TSStep()
288 
289     Notes:
290     The available visualization contexts include
291 +     VIEWER_STDOUT_SELF - standard output (default)
292 -     VIEWER_STDOUT_WORLD - synchronized standard
293          output where only the first processor opens
294          the file.  All other processors send their
295          data to the first processor to print.
296 
297     The user can open an alternative visualization context with
298     ViewerASCIIOpen() - output to a specified file.
299 
300     Level: beginner
301 
302 .keywords: TS, timestep, view
303 
304 .seealso: ViewerASCIIOpen()
305 @*/
306 int TSView(TS ts,Viewer viewer)
307 {
308   int                 ierr;
309   char                *method;
310   ViewerType          vtype;
311 
312   PetscFunctionBegin;
313   PetscValidHeaderSpecific(ts,TS_COOKIE);
314   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
315   if (PetscTypeCompare(vtype,ASCII_VIEWER)) {
316     ViewerASCIIPrintf(viewer,"TS Object:\n");
317     ierr = TSGetType(ts,(TSType *)&method);CHKERRQ(ierr);
318     ViewerASCIIPrintf(viewer,"  method: %s\n",method);
319     if (ts->view) {
320       ierr = ViewerASCIIPushTab(viewer);CHKERRQ(ierr);
321       ierr = (*ts->view)(ts,viewer);CHKERRQ(ierr);
322       ierr = ViewerASCIIPopTab(viewer);CHKERRQ(ierr);
323     }
324     ViewerASCIIPrintf(viewer,"  maximum steps=%d\n",ts->max_steps);
325     ViewerASCIIPrintf(viewer,"  maximum time=%g\n",ts->max_time);
326     if (ts->problem_type == TS_NONLINEAR) {
327       ViewerASCIIPrintf(viewer,"  total number of nonlinear solver iterations=%d\n",ts->nonlinear_its);
328     }
329     ViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%d\n",ts->linear_its);
330   } else if (PetscTypeCompare(vtype,STRING_VIEWER)) {
331     ierr = TSGetType(ts,(TSType *)&method);CHKERRQ(ierr);
332     ierr = ViewerStringSPrintf(viewer," %-7.7s",method);CHKERRQ(ierr);
333   }
334   ierr = ViewerASCIIPushTab(viewer);CHKERRQ(ierr);
335   if (ts->sles) {ierr = SLESView(ts->sles,viewer);CHKERRQ(ierr);}
336   if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
337   ierr = ViewerASCIIPopTab(viewer);CHKERRQ(ierr);
338   PetscFunctionReturn(0);
339 }
340 
341 
342 #undef __FUNC__
343 #define __FUNC__ "TSSetApplicationContext"
344 /*@C
345    TSSetApplicationContext - Sets an optional user-defined context for
346    the timesteppers.
347 
348    Collective on TS
349 
350    Input Parameters:
351 +  ts - the TS context obtained from TSCreate()
352 -  usrP - optional user context
353 
354    Level: intermediate
355 
356 .keywords: TS, timestep, set, application, context
357 
358 .seealso: TSGetApplicationContext()
359 @*/
360 int TSSetApplicationContext(TS ts,void *usrP)
361 {
362   PetscFunctionBegin;
363   PetscValidHeaderSpecific(ts,TS_COOKIE);
364   ts->user = usrP;
365   PetscFunctionReturn(0);
366 }
367 
368 #undef __FUNC__
369 #define __FUNC__ "TSGetApplicationContext"
370 /*@C
371     TSGetApplicationContext - Gets the user-defined context for the
372     timestepper.
373 
374     Not Collective
375 
376     Input Parameter:
377 .   ts - the TS context obtained from TSCreate()
378 
379     Output Parameter:
380 .   usrP - user context
381 
382     Level: intermediate
383 
384 .keywords: TS, timestep, get, application, context
385 
386 .seealso: TSSetApplicationContext()
387 @*/
388 int TSGetApplicationContext( TS ts,  void **usrP )
389 {
390   PetscFunctionBegin;
391   PetscValidHeaderSpecific(ts,TS_COOKIE);
392   *usrP = ts->user;
393   PetscFunctionReturn(0);
394 }
395 
396 #undef __FUNC__
397 #define __FUNC__ "TSGetTimeStepNumber"
398 /*@
399    TSGetTimeStepNumber - Gets the current number of timesteps.
400 
401    Not Collective
402 
403    Input Parameter:
404 .  ts - the TS context obtained from TSCreate()
405 
406    Output Parameter:
407 .  iter - number steps so far
408 
409    Level: intermediate
410 
411 .keywords: TS, timestep, get, iteration, number
412 @*/
413 int TSGetTimeStepNumber(TS ts,int* iter)
414 {
415   PetscFunctionBegin;
416   PetscValidHeaderSpecific(ts,TS_COOKIE);
417   *iter = ts->steps;
418   PetscFunctionReturn(0);
419 }
420 
421 #undef __FUNC__
422 #define __FUNC__ "TSSetInitialTimeStep"
423 /*@
424    TSSetInitialTimeStep - Sets the initial timestep to be used,
425    as well as the initial time.
426 
427    Collective on TS
428 
429    Input Parameters:
430 +  ts - the TS context obtained from TSCreate()
431 .  initial_time - the initial time
432 -  time_step - the size of the timestep
433 
434    Level: intermediate
435 
436 .seealso: TSSetTimeStep(), TSGetTimeStep()
437 
438 .keywords: TS, set, initial, timestep
439 @*/
440 int TSSetInitialTimeStep(TS ts,double initial_time,double time_step)
441 {
442   PetscFunctionBegin;
443   PetscValidHeaderSpecific(ts,TS_COOKIE);
444   ts->time_step         = time_step;
445   ts->initial_time_step = time_step;
446   ts->ptime             = initial_time;
447   PetscFunctionReturn(0);
448 }
449 
450 #undef __FUNC__
451 #define __FUNC__ "TSSetTimeStep"
452 /*@
453    TSSetTimeStep - Allows one to reset the timestep at any time,
454    useful for simple pseudo-timestepping codes.
455 
456    Collective on TS
457 
458    Input Parameters:
459 +  ts - the TS context obtained from TSCreate()
460 -  time_step - the size of the timestep
461 
462    Level: intermediate
463 
464 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
465 
466 .keywords: TS, set, timestep
467 @*/
468 int TSSetTimeStep(TS ts,double time_step)
469 {
470   PetscFunctionBegin;
471   PetscValidHeaderSpecific(ts,TS_COOKIE);
472   ts->time_step = time_step;
473   PetscFunctionReturn(0);
474 }
475 
476 #undef __FUNC__
477 #define __FUNC__ "TSGetTimeStep"
478 /*@
479    TSGetTimeStep - Gets the current timestep size.
480 
481    Not Collective
482 
483    Input Parameter:
484 .  ts - the TS context obtained from TSCreate()
485 
486    Output Parameter:
487 .  dt - the current timestep size
488 
489    Level: intermediate
490 
491 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
492 
493 .keywords: TS, get, timestep
494 @*/
495 int TSGetTimeStep(TS ts,double* dt)
496 {
497   PetscFunctionBegin;
498   PetscValidHeaderSpecific(ts,TS_COOKIE);
499   *dt = ts->time_step;
500   PetscFunctionReturn(0);
501 }
502 
503 #undef __FUNC__
504 #define __FUNC__ "TSGetSolution"
505 /*@C
506    TSGetSolution - Returns the solution at the present timestep. It
507    is valid to call this routine inside the function that you are evaluating
508    in order to move to the new timestep. This vector not changed until
509    the solution at the next timestep has been calculated.
510 
511    Not Collective, but Vec returned is parallel if TS is parallel
512 
513    Input Parameter:
514 .  ts - the TS context obtained from TSCreate()
515 
516    Output Parameter:
517 .  v - the vector containing the solution
518 
519    Level: intermediate
520 
521 .seealso: TSGetTimeStep()
522 
523 .keywords: TS, timestep, get, solution
524 @*/
525 int TSGetSolution(TS ts,Vec *v)
526 {
527   PetscFunctionBegin;
528   PetscValidHeaderSpecific(ts,TS_COOKIE);
529   *v = ts->vec_sol_always;
530   PetscFunctionReturn(0);
531 }
532 
533 #undef __FUNC__
534 #define __FUNC__ "TSPublish_Petsc"
535 static int TSPublish_Petsc(PetscObject object)
536 {
537 #if defined(HAVE_AMS)
538   TS   v = (TS) object;
539   int  ierr;
540 
541   PetscFunctionBegin;
542 
543   /* if it is already published then return */
544   if (v->amem >=0 ) PetscFunctionReturn(0);
545 
546   ierr = PetscObjectPublishBaseBegin(object);CHKERRQ(ierr);
547   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Step",&v->steps,1,AMS_INT,AMS_READ,
548                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
549   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Time",&v->ptime,1,AMS_DOUBLE,AMS_READ,
550                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
551   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"CurrentTimeStep",&v->time_step,1,
552                                AMS_DOUBLE,AMS_READ,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
553   ierr = PetscObjectPublishBaseEnd(object);CHKERRQ(ierr);
554 #else
555   PetscFunctionBegin;
556 #endif
557   PetscFunctionReturn(0);
558 }
559 
560 /* -----------------------------------------------------------*/
561 
562 #undef __FUNC__
563 #define __FUNC__ "TSCreate"
564 /*@C
565    TSCreate - Creates a timestepper context.
566 
567    Collective on MPI_Comm
568 
569    Input Parameter:
570 +  comm - MPI communicator
571 -  type - One of  TS_LINEAR,TS_NONLINEAR
572    where these types refer to problems of the forms
573 .vb
574          U_t = A U
575          U_t = A(t) U
576          U_t = F(t,U)
577 .ve
578 
579    Output Parameter:
580 .  outts - the new TS context
581 
582    Level: beginner
583 
584 .keywords: TS, timestep, create, context
585 
586 .seealso: TSSetUp(), TSStep(), TSDestroy()
587 @*/
588 int TSCreate(MPI_Comm comm,TSProblemType problemtype,TS *outts)
589 {
590   TS   ts;
591 
592   PetscFunctionBegin;
593   *outts = 0;
594   PetscHeaderCreate(ts,_p_TS,int,TS_COOKIE,-1,"TS",comm,TSDestroy,TSView);
595   PLogObjectCreate(ts);
596   ts->bops->publish     = TSPublish_Petsc;
597   ts->max_steps         = 5000;
598   ts->max_time          = 5.0;
599   ts->time_step         = .1;
600   ts->initial_time_step = ts->time_step;
601   ts->steps             = 0;
602   ts->ptime             = 0.0;
603   ts->data              = 0;
604   ts->view              = 0;
605   ts->setupcalled       = 0;
606   ts->problem_type      = problemtype;
607   ts->numbermonitors    = 0;
608   ts->linear_its        = 0;
609   ts->nonlinear_its     = 0;
610 
611   *outts = ts;
612   PetscFunctionReturn(0);
613 }
614 
615 /* ----- Routines to initialize and destroy a timestepper ---- */
616 
617 #undef __FUNC__
618 #define __FUNC__ "TSSetUp"
619 /*@
620    TSSetUp - Sets up the internal data structures for the later use
621    of a timestepper.
622 
623    Collective on TS
624 
625    Input Parameter:
626 .  ts - the TS context obtained from TSCreate()
627 
628    Notes:
629    For basic use of the TS solvers the user need not explicitly call
630    TSSetUp(), since these actions will automatically occur during
631    the call to TSStep().  However, if one wishes to control this
632    phase separately, TSSetUp() should be called after TSCreate()
633    and optional routines of the form TSSetXXX(), but before TSStep().
634 
635    Level: advanced
636 
637 .keywords: TS, timestep, setup
638 
639 .seealso: TSCreate(), TSStep(), TSDestroy()
640 @*/
641 int TSSetUp(TS ts)
642 {
643   int ierr;
644 
645   PetscFunctionBegin;
646   PetscValidHeaderSpecific(ts,TS_COOKIE);
647   if (!ts->vec_sol) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call TSSetSolution() first");
648   if (!ts->type_name) {
649     ierr = TSSetType(ts,TS_EULER);CHKERRQ(ierr);
650   }
651   ierr = (*ts->setup)(ts); CHKERRQ(ierr);
652   ts->setupcalled = 1;
653   PetscFunctionReturn(0);
654 }
655 
656 #undef __FUNC__
657 #define __FUNC__ "TSDestroy"
658 /*@C
659    TSDestroy - Destroys the timestepper context that was created
660    with TSCreate().
661 
662    Collective on TS
663 
664    Input Parameter:
665 .  ts - the TS context obtained from TSCreate()
666 
667    Level: beginner
668 
669 .keywords: TS, timestepper, destroy
670 
671 .seealso: TSCreate(), TSSetUp(), TSSolve()
672 @*/
673 int TSDestroy(TS ts)
674 {
675   int ierr;
676 
677   PetscFunctionBegin;
678   PetscValidHeaderSpecific(ts,TS_COOKIE);
679   if (--ts->refct > 0) PetscFunctionReturn(0);
680 
681   if (ts->sles) {ierr = SLESDestroy(ts->sles); CHKERRQ(ierr);}
682   if (ts->snes) {ierr = SNESDestroy(ts->snes); CHKERRQ(ierr);}
683   ierr = (*(ts)->destroy)(ts); CHKERRQ(ierr);
684   PLogObjectDestroy((PetscObject)ts);
685   PetscHeaderDestroy((PetscObject)ts);
686   PetscFunctionReturn(0);
687 }
688 
689 #undef __FUNC__
690 #define __FUNC__ "TSGetSNES"
691 /*@C
692    TSGetSNES - Returns the SNES (nonlinear solver) associated with
693    a TS (timestepper) context. Valid only for nonlinear problems.
694 
695    Not Collective, but SNES is parallel if TS is parallel
696 
697    Input Parameter:
698 .  ts - the TS context obtained from TSCreate()
699 
700    Output Parameter:
701 .  snes - the nonlinear solver context
702 
703    Notes:
704    The user can then directly manipulate the SNES context to set various
705    options, etc.  Likewise, the user can then extract and manipulate the
706    SLES, KSP, and PC contexts as well.
707 
708    TSGetSNES() does not work for integrators that do not use SNES; in
709    this case TSGetSNES() returns PETSC_NULL in snes.
710 
711    Level: beginner
712 
713 .keywords: timestep, get, SNES
714 @*/
715 int TSGetSNES(TS ts,SNES *snes)
716 {
717   PetscFunctionBegin;
718   PetscValidHeaderSpecific(ts,TS_COOKIE);
719   if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Nonlinear only; use TSGetSLES()");
720   *snes = ts->snes;
721   PetscFunctionReturn(0);
722 }
723 
724 #undef __FUNC__
725 #define __FUNC__ "TSGetSLES"
726 /*@C
727    TSGetSLES - Returns the SLES (linear solver) associated with
728    a TS (timestepper) context.
729 
730    Not Collective, but SLES is parallel if TS is parallel
731 
732    Input Parameter:
733 .  ts - the TS context obtained from TSCreate()
734 
735    Output Parameter:
736 .  sles - the nonlinear solver context
737 
738    Notes:
739    The user can then directly manipulate the SLES context to set various
740    options, etc.  Likewise, the user can then extract and manipulate the
741    KSP and PC contexts as well.
742 
743    TSGetSLES() does not work for integrators that do not use SLES;
744    in this case TSGetSLES() returns PETSC_NULL in sles.
745 
746    Level: beginner
747 
748 .keywords: timestep, get, SLES
749 @*/
750 int TSGetSLES(TS ts,SLES *sles)
751 {
752   PetscFunctionBegin;
753   PetscValidHeaderSpecific(ts,TS_COOKIE);
754   if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Linear only; use TSGetSNES()");
755   *sles = ts->sles;
756   PetscFunctionReturn(0);
757 }
758 
759 /* ----------- Routines to set solver parameters ---------- */
760 
761 #undef __FUNC__
762 #define __FUNC__ "TSSetDuration"
763 /*@
764    TSSetDuration - Sets the maximum number of timesteps to use and
765    maximum time for iteration.
766 
767    Collective on TS
768 
769    Input Parameters:
770 +  ts - the TS context obtained from TSCreate()
771 .  maxsteps - maximum number of iterations to use
772 -  maxtime - final time to iterate to
773 
774    Options Database Keys:
775 .  -ts_max_steps <maxsteps> - Sets maxsteps
776 .  -ts_max_time <maxtime> - Sets maxtime
777 
778    Notes:
779    The default maximum number of iterations is 5000. Default time is 5.0
780 
781    Level: intermediate
782 
783 .keywords: TS, timestep, set, maximum, iterations
784 @*/
785 int TSSetDuration(TS ts,int maxsteps,double maxtime)
786 {
787   PetscFunctionBegin;
788   PetscValidHeaderSpecific(ts,TS_COOKIE);
789   ts->max_steps = maxsteps;
790   ts->max_time  = maxtime;
791   PetscFunctionReturn(0);
792 }
793 
794 #undef __FUNC__
795 #define __FUNC__ "TSSetSolution"
796 /*@
797    TSSetSolution - Sets the initial solution vector
798    for use by the TS routines.
799 
800    Collective on TS and Vec
801 
802    Input Parameters:
803 +  ts - the TS context obtained from TSCreate()
804 -  x - the solution vector
805 
806    Level: beginner
807 
808 .keywords: TS, timestep, set, solution, initial conditions
809 @*/
810 int TSSetSolution(TS ts,Vec x)
811 {
812   PetscFunctionBegin;
813   PetscValidHeaderSpecific(ts,TS_COOKIE);
814   ts->vec_sol        = ts->vec_sol_always = x;
815   PetscFunctionReturn(0);
816 }
817 
818 /* ------------ Routines to set performance monitoring options ----------- */
819 
820 #undef __FUNC__
821 #define __FUNC__ "TSSetMonitor"
822 /*@C
823    TSSetMonitor - Sets an ADDITIONAL function that is to be used at every
824    timestep to display the iteration's  progress.
825 
826    Collective on TS
827 
828    Input Parameters:
829 +  ts - the TS context obtained from TSCreate()
830 .  func - monitoring routine
831 -  mctx - [optional] user-defined context for private data for the
832           monitor routine (may be PETSC_NULL)
833 
834    Calling sequence of func:
835 $    int func(TS ts,int steps,double time,Vec x,void *mctx)
836 
837 +    ts - the TS context
838 .    steps - iteration number
839 .    time - current timestep
840 .    x - current iterate
841 -    mctx - [optional] monitoring context
842 
843    Notes:
844    This routine adds an additional monitor to the list of monitors that
845    already has been loaded.
846 
847    Level: intermediate
848 
849 .keywords: TS, timestep, set, monitor
850 
851 .seealso: TSDefaultMonitor(), TSClearMonitor()
852 @*/
853 int TSSetMonitor(TS ts, int (*monitor)(TS,int,double,Vec,void*), void *mctx )
854 {
855   PetscFunctionBegin;
856   PetscValidHeaderSpecific(ts,TS_COOKIE);
857   if (ts->numbermonitors >= MAXTSMONITORS) {
858     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many monitors set");
859   }
860   ts->monitor[ts->numbermonitors]           = monitor;
861   ts->monitorcontext[ts->numbermonitors++]  = (void*)mctx;
862   PetscFunctionReturn(0);
863 }
864 
865 #undef __FUNC__
866 #define __FUNC__ "TSClearMonitor"
867 /*@C
868    TSClearMonitor - Clears all the monitors that have been set on a time-step object.
869 
870    Collective on TS
871 
872    Input Parameters:
873 .  ts - the TS context obtained from TSCreate()
874 
875    Notes:
876    There is no way to remove a single, specific monitor.
877 
878    Level: intermediate
879 
880 .keywords: TS, timestep, set, monitor
881 
882 .seealso: TSDefaultMonitor(), TSSetMonitor()
883 @*/
884 int TSClearMonitor(TS ts)
885 {
886   PetscFunctionBegin;
887   PetscValidHeaderSpecific(ts,TS_COOKIE);
888   ts->numbermonitors = 0;
889   PetscFunctionReturn(0);
890 }
891 
892 #undef __FUNC__
893 #define __FUNC__ "TSDefaultMonitor"
894 int TSDefaultMonitor(TS ts, int step, double time,Vec v, void *ctx)
895 {
896   PetscFunctionBegin;
897   PetscPrintf(ts->comm,"timestep %d dt %g time %g\n",step,ts->time_step,time);
898   PetscFunctionReturn(0);
899 }
900 
901 #undef __FUNC__
902 #define __FUNC__ "TSStep"
903 /*@
904    TSStep - Steps the requested number of timesteps.
905 
906    Collective on TS
907 
908    Input Parameter:
909 .  ts - the TS context obtained from TSCreate()
910 
911    Output Parameters:
912 +  steps - number of iterations until termination
913 -  time - time until termination
914 
915    Level: beginner
916 
917 .keywords: TS, timestep, solve
918 
919 .seealso: TSCreate(), TSSetUp(), TSDestroy()
920 @*/
921 int TSStep(TS ts,int *steps,double *time)
922 {
923   int ierr,flg;
924 
925   PetscFunctionBegin;
926   PetscValidHeaderSpecific(ts,TS_COOKIE);
927   if (!ts->setupcalled) {ierr = TSSetUp(ts); CHKERRQ(ierr);}
928   PLogEventBegin(TS_Step,ts,0,0,0);
929   ierr = (*(ts)->step)(ts,steps,time); CHKERRQ(ierr);
930   PLogEventEnd(TS_Step,ts,0,0,0);
931   ierr = OptionsHasName(PETSC_NULL,"-ts_view",&flg); CHKERRQ(ierr);
932   if (flg) {
933     ierr = TSView(ts,VIEWER_STDOUT_WORLD); CHKERRQ(ierr);
934   }
935   PetscFunctionReturn(0);
936 }
937 
938 #undef __FUNC__
939 #define __FUNC__ "TSMonitor"
940 /*
941      Runs the user provided monitor routines, if they exists.
942 */
943 int TSMonitor(TS ts,int step,double time,Vec x)
944 {
945   int i,ierr,n = ts->numbermonitors;
946 
947   PetscFunctionBegin;
948   for ( i=0; i<n; i++ ) {
949     ierr = (*ts->monitor[i])(ts,step,time,x,ts->monitorcontext[i]);CHKERRQ(ierr);
950   }
951   PetscFunctionReturn(0);
952 }
953 
954 /* ------------------------------------------------------------------------*/
955 
956 /*@C
957    TSLGMonitorCreate - Creates a line graph context for use with
958    TS to monitor convergence of preconditioned residual norms.
959 
960    Collective on TS
961 
962    Input Parameters:
963 +  host - the X display to open, or null for the local machine
964 .  label - the title to put in the title bar
965 .  x, y - the screen coordinates of the upper left coordinate of
966           the window
967 -  m, n - the screen width and height in pixels
968 
969    Output Parameter:
970 .  draw - the drawing context
971 
972    Options Database Key:
973 .  -ts_xmonitor - automatically sets line graph monitor
974 
975    Notes:
976    Use TSLGMonitorDestroy() to destroy this line graph, not DrawLGDestroy().
977 
978    Level: intermediate
979 
980 .keywords: TS, monitor, line graph, residual, create
981 
982 .seealso: TSLGMonitorDestroy(), TSSetMonitor()
983 @*/
984 int TSLGMonitorCreate(char *host,char *label,int x,int y,int m,
985                        int n, DrawLG *draw)
986 {
987   Draw win;
988   int  ierr;
989 
990   PetscFunctionBegin;
991   ierr = DrawOpenX(PETSC_COMM_SELF,host,label,x,y,m,n,&win); CHKERRQ(ierr);
992   ierr = DrawLGCreate(win,1,draw); CHKERRQ(ierr);
993   ierr = DrawLGIndicateDataPoints(*draw); CHKERRQ(ierr);
994 
995   PLogObjectParent(*draw,win);
996   PetscFunctionReturn(0);
997 }
998 
999 #undef __FUNC__
1000 #define __FUNC__ "TSLGMonitor"
1001 int TSLGMonitor(TS ts,int n,double time,Vec v,void *monctx)
1002 {
1003   DrawLG lg = (DrawLG) monctx;
1004   double x,y = time;
1005   int    ierr;
1006 
1007   PetscFunctionBegin;
1008   if (!n) {ierr = DrawLGReset(lg);CHKERRQ(ierr);}
1009   x = (double) n;
1010   ierr = DrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1011   if (n < 20 || (n % 5)) {
1012     ierr = DrawLGDraw(lg);CHKERRQ(ierr);
1013   }
1014   PetscFunctionReturn(0);
1015 }
1016 
1017 #undef __FUNC__
1018 #define __FUNC__ "TSLGMonitorDestroy"
1019 /*@C
1020    TSLGMonitorDestroy - Destroys a line graph context that was created
1021    with TSLGMonitorCreate().
1022 
1023    Collective on DrawLG
1024 
1025    Input Parameter:
1026 .  draw - the drawing context
1027 
1028    Level: intermediate
1029 
1030 .keywords: TS, monitor, line graph, destroy
1031 
1032 .seealso: TSLGMonitorCreate(),  TSSetMonitor(), TSLGMonitor();
1033 @*/
1034 int TSLGMonitorDestroy(DrawLG drawlg)
1035 {
1036   Draw draw;
1037   int  ierr;
1038 
1039   PetscFunctionBegin;
1040   ierr = DrawLGGetDraw(drawlg,&draw);CHKERRQ(ierr);
1041   ierr = DrawDestroy(draw);CHKERRQ(ierr);
1042   ierr = DrawLGDestroy(drawlg);CHKERRQ(ierr);
1043   PetscFunctionReturn(0);
1044 }
1045 
1046 #undef __FUNC__
1047 #define __FUNC__ "TSGetTime"
1048 /*@
1049    TSGetTime - Gets the current time.
1050 
1051    Not Collective
1052 
1053    Input Parameter:
1054 .  ts - the TS context obtained from TSCreate()
1055 
1056    Output Parameter:
1057 .  t  - the current time
1058 
1059    Contributed by: Matthew Knepley
1060 
1061    Level: beginner
1062 
1063 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1064 
1065 .keywords: TS, get, time
1066 @*/
1067 int TSGetTime(TS ts, double* t)
1068 {
1069   PetscFunctionBegin;
1070   PetscValidHeaderSpecific(ts, TS_COOKIE);
1071   *t = ts->ptime;
1072   PetscFunctionReturn(0);
1073 }
1074 
1075 #undef __FUNC__
1076 #define __FUNC__ "TSGetProblemType"
1077 /*@C
1078    TSGetProblemType - Returns the problem type of a TS (timestepper) context.
1079 
1080    Not Collective
1081 
1082    Input Parameter:
1083 .  ts   - The TS context obtained from TSCreate()
1084 
1085    Output Parameter:
1086 .  type - The problem type, TS_LINEAR or TS_NONLINEAR
1087 
1088    Level: intermediate
1089 
1090    Contributed by: Matthew Knepley
1091 
1092 .keywords: ts, get, type
1093 
1094 @*/
1095 int TSGetProblemType(TS ts, TSProblemType *type)
1096 {
1097   PetscFunctionBegin;
1098   PetscValidHeaderSpecific(ts, TS_COOKIE);
1099   *type = ts->problem_type;
1100   PetscFunctionReturn(0);
1101 }
1102 
1103 #undef __FUNC__
1104 #define __FUNC__ "TSSetOptionsPrefix"
1105 /*@C
1106    TSSetOptionsPrefix - Sets the prefix used for searching for all
1107    TS options in the database.
1108 
1109    Collective on TS
1110 
1111    Input Parameter:
1112 +  ts     - The TS context
1113 -  prefix - The prefix to prepend to all option names
1114 
1115    Notes:
1116    A hyphen (-) must NOT be given at the beginning of the prefix name.
1117    The first character of all runtime options is AUTOMATICALLY the
1118    hyphen.
1119 
1120    Contributed by: Matthew Knepley
1121 
1122    Level: advanced
1123 
1124 .keywords: TS, set, options, prefix, database
1125 
1126 .seealso: TSSetFromOptions()
1127 
1128 @*/
1129 int TSSetOptionsPrefix(TS ts, char *prefix)
1130 {
1131   int ierr;
1132 
1133   PetscFunctionBegin;
1134   PetscValidHeaderSpecific(ts, TS_COOKIE);
1135   ierr = PetscObjectSetOptionsPrefix((PetscObject) ts, prefix); CHKERRQ(ierr);
1136   switch(ts->problem_type) {
1137     case TS_NONLINEAR:
1138       ierr = SNESSetOptionsPrefix(ts->snes, prefix);              CHKERRQ(ierr);
1139       break;
1140     case TS_LINEAR:
1141       ierr = SLESSetOptionsPrefix(ts->sles, prefix);              CHKERRQ(ierr);
1142       break;
1143   }
1144   PetscFunctionReturn(0);
1145 }
1146 
1147 
1148 #undef __FUNC__
1149 #define __FUNC__ "TSAppendOptionsPrefix"
1150 /*@C
1151    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
1152    TS options in the database.
1153 
1154    Collective on TS
1155 
1156    Input Parameter:
1157 +  ts     - The TS context
1158 -  prefix - The prefix to prepend to all option names
1159 
1160    Notes:
1161    A hyphen (-) must NOT be given at the beginning of the prefix name.
1162    The first character of all runtime options is AUTOMATICALLY the
1163    hyphen.
1164 
1165    Contributed by: Matthew Knepley
1166 
1167    Level: advanced
1168 
1169 .keywords: TS, append, options, prefix, database
1170 
1171 .seealso: TSGetOptionsPrefix()
1172 
1173 @*/
1174 int TSAppendOptionsPrefix(TS ts, char *prefix)
1175 {
1176   int ierr;
1177 
1178   PetscFunctionBegin;
1179   PetscValidHeaderSpecific(ts, TS_COOKIE);
1180   ierr = PetscObjectAppendOptionsPrefix((PetscObject) ts, prefix); CHKERRQ(ierr);
1181   switch(ts->problem_type) {
1182     case TS_NONLINEAR:
1183       ierr = SNESAppendOptionsPrefix(ts->snes, prefix);              CHKERRQ(ierr);
1184       break;
1185     case TS_LINEAR:
1186       ierr = SLESAppendOptionsPrefix(ts->sles, prefix);              CHKERRQ(ierr);
1187       break;
1188   }
1189   PetscFunctionReturn(0);
1190 }
1191 
1192 #undef __FUNC__
1193 #define __FUNC__ "TSGetOptionsPrefix"
1194 /*@C
1195    TSGetOptionsPrefix - Sets the prefix used for searching for all
1196    TS options in the database.
1197 
1198    Not Collective
1199 
1200    Input Parameter:
1201 .  ts - The TS context
1202 
1203    Output Parameter:
1204 .  prefix - A pointer to the prefix string used
1205 
1206    Contributed by: Matthew Knepley
1207 
1208    Notes: On the fortran side, the user should pass in a string 'prifix' of
1209    sufficient length to hold the prefix.
1210 
1211    Level: intermediate
1212 
1213 .keywords: TS, get, options, prefix, database
1214 
1215 .seealso: TSAppendOptionsPrefix()
1216 @*/
1217 int TSGetOptionsPrefix(TS ts, char **prefix)
1218 {
1219   int ierr;
1220 
1221   PetscFunctionBegin;
1222   PetscValidHeaderSpecific(ts, TS_COOKIE);
1223   ierr = PetscObjectGetOptionsPrefix((PetscObject) ts, prefix); CHKERRQ(ierr);
1224   PetscFunctionReturn(0);
1225 }
1226 
1227 #undef __FUNC__
1228 #define __FUNC__ "TSGetRHSMatrix"
1229 /*@C
1230    TSGetRHSMatrix - Returns the matrix A at the present timestep.
1231 
1232    Not Collective, but parallel objects are returned if TS is parallel
1233 
1234    Input Parameter:
1235 .  ts  - The TS context obtained from TSCreate()
1236 
1237    Output Parameters:
1238 +  A   - The matrix A, where U_t = A(t) U
1239 .  M   - The preconditioner matrix, usually the same as A
1240 -  ctx - User-defined context for matrix evaluation routine
1241 
1242    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1243 
1244    Contributed by: Matthew Knepley
1245 
1246    Level: intermediate
1247 
1248 .seealso: TSGetTimeStep(), TSGetTime(), TSGetTimeStepNumber(), TSGetRHSJacobian()
1249 
1250 .keywords: TS, timestep, get, matrix
1251 
1252 @*/
1253 int TSGetRHSMatrix(TS ts, Mat *A, Mat *M, void **ctx)
1254 {
1255   PetscFunctionBegin;
1256   PetscValidHeaderSpecific(ts, TS_COOKIE);
1257   if (A)   *A = ts->A;
1258   if (M)   *M = ts->B;
1259   if (ctx) *ctx = ts->jacP;
1260   PetscFunctionReturn(0);
1261 }
1262 
1263 #undef __FUNC__
1264 #define __FUNC__ "TSGetRHSJacobian"
1265 /*@C
1266    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
1267 
1268    Not Collective, but parallel objects are returned if TS is parallel
1269 
1270    Input Parameter:
1271 .  ts  - The TS context obtained from TSCreate()
1272 
1273    Output Parameters:
1274 +  J   - The Jacobian J of F, where U_t = F(U,t)
1275 .  M   - The preconditioner matrix, usually the same as J
1276 - ctx - User-defined context for Jacobian evaluation routine
1277 
1278    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1279 
1280    Contributed by: Matthew Knepley
1281 
1282    Level: intermediate
1283 
1284 .seealso: TSGetTimeStep(), TSGetRHSMatrix(), TSGetTime(), TSGetTimeStepNumber()
1285 
1286 .keywords: TS, timestep, get, matrix, Jacobian
1287 @*/
1288 int TSGetRHSJacobian(TS ts, Mat *J, Mat *M, void **ctx)
1289 {
1290   int ierr;
1291 
1292   PetscFunctionBegin;
1293   ierr = TSGetRHSMatrix(ts, J, M, ctx);CHKERRQ(ierr);
1294   PetscFunctionReturn(0);
1295 }
1296 
1297 /*MC
1298    TSRegister - Adds a method to the timestepping solver package.
1299 
1300    Synopsis:
1301 
1302    TSRegister(char *name_solver,char *path,char *name_create,int (*routine_create)(TS))
1303 
1304    Not collective
1305 
1306    Input Parameters:
1307 +  name_solver - name of a new user-defined solver
1308 .  path - path (either absolute or relative) the library containing this solver
1309 .  name_create - name of routine to create method context
1310 -  routine_create - routine to create method context
1311 
1312    Notes:
1313    TSRegister() may be called multiple times to add several user-defined solvers.
1314 
1315    If dynamic libraries are used, then the fourth input argument (routine_create)
1316    is ignored.
1317 
1318    Sample usage:
1319 .vb
1320    TSRegister("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
1321               "MySolverCreate",MySolverCreate);
1322 .ve
1323 
1324    Then, your solver can be chosen with the procedural interface via
1325 $     TSSetType(ts,"my_solver")
1326    or at runtime via the option
1327 $     -ts_type my_solver
1328 
1329    Level: advanced
1330 
1331    $PETSC_ARCH and $BOPT occuring in pathname will be replaced with appropriate values.
1332 
1333 .keywords: TS, register
1334 
1335 .seealso: TSRegisterAll(), TSRegisterDestroy()
1336 M*/
1337 
1338 #undef __FUNC__
1339 #define __FUNC__ "TSRegister_Private"
1340 int TSRegister_Private(char *sname,char *path,char *name,int (*function)(TS))
1341 {
1342   char fullname[256];
1343 
1344   PetscFunctionBegin;
1345   PetscStrcpy(fullname,path); PetscStrcat(fullname,":"); PetscStrcat(fullname,name);
1346   FListAdd_Private(&TSList,sname,fullname,        (int (*)(void*))function);
1347   PetscFunctionReturn(0);
1348 }
1349