xref: /petsc/src/ts/interface/ts.c (revision fd16b1779ca341c090fa91c269e2c2dd5dd0b73d)
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   if (!viewer) viewer = VIEWER_STDOUT_SELF;
315 
316   ierr = ViewerGetType(viewer,&vtype);CHKERRQ(ierr);
317   if (PetscTypeCompare(vtype,ASCII_VIEWER)) {
318     ierr = ViewerASCIIPrintf(viewer,"TS Object:\n");CHKERRQ(ierr);
319     ierr = TSGetType(ts,(TSType *)&method);CHKERRQ(ierr);
320     ierr = ViewerASCIIPrintf(viewer,"  method: %s\n",method);CHKERRQ(ierr);
321     if (ts->view) {
322       ierr = ViewerASCIIPushTab(viewer);CHKERRQ(ierr);
323       ierr = (*ts->view)(ts,viewer);CHKERRQ(ierr);
324       ierr = ViewerASCIIPopTab(viewer);CHKERRQ(ierr);
325     }
326     ierr = ViewerASCIIPrintf(viewer,"  maximum steps=%d\n",ts->max_steps);CHKERRQ(ierr);
327     ierr = ViewerASCIIPrintf(viewer,"  maximum time=%g\n",ts->max_time);CHKERRQ(ierr);
328     if (ts->problem_type == TS_NONLINEAR) {
329       ierr = ViewerASCIIPrintf(viewer,"  total number of nonlinear solver iterations=%d\n",ts->nonlinear_its);CHKERRQ(ierr);
330     }
331     ierr = ViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%d\n",ts->linear_its);CHKERRQ(ierr);
332   } else if (PetscTypeCompare(vtype,STRING_VIEWER)) {
333     ierr = TSGetType(ts,(TSType *)&method);CHKERRQ(ierr);
334     ierr = ViewerStringSPrintf(viewer," %-7.7s",method);CHKERRQ(ierr);
335   }
336   ierr = ViewerASCIIPushTab(viewer);CHKERRQ(ierr);
337   if (ts->sles) {ierr = SLESView(ts->sles,viewer);CHKERRQ(ierr);}
338   if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
339   ierr = ViewerASCIIPopTab(viewer);CHKERRQ(ierr);
340   PetscFunctionReturn(0);
341 }
342 
343 
344 #undef __FUNC__
345 #define __FUNC__ "TSSetApplicationContext"
346 /*@C
347    TSSetApplicationContext - Sets an optional user-defined context for
348    the timesteppers.
349 
350    Collective on TS
351 
352    Input Parameters:
353 +  ts - the TS context obtained from TSCreate()
354 -  usrP - optional user context
355 
356    Level: intermediate
357 
358 .keywords: TS, timestep, set, application, context
359 
360 .seealso: TSGetApplicationContext()
361 @*/
362 int TSSetApplicationContext(TS ts,void *usrP)
363 {
364   PetscFunctionBegin;
365   PetscValidHeaderSpecific(ts,TS_COOKIE);
366   ts->user = usrP;
367   PetscFunctionReturn(0);
368 }
369 
370 #undef __FUNC__
371 #define __FUNC__ "TSGetApplicationContext"
372 /*@C
373     TSGetApplicationContext - Gets the user-defined context for the
374     timestepper.
375 
376     Not Collective
377 
378     Input Parameter:
379 .   ts - the TS context obtained from TSCreate()
380 
381     Output Parameter:
382 .   usrP - user context
383 
384     Level: intermediate
385 
386 .keywords: TS, timestep, get, application, context
387 
388 .seealso: TSSetApplicationContext()
389 @*/
390 int TSGetApplicationContext( TS ts,  void **usrP )
391 {
392   PetscFunctionBegin;
393   PetscValidHeaderSpecific(ts,TS_COOKIE);
394   *usrP = ts->user;
395   PetscFunctionReturn(0);
396 }
397 
398 #undef __FUNC__
399 #define __FUNC__ "TSGetTimeStepNumber"
400 /*@
401    TSGetTimeStepNumber - Gets the current number of timesteps.
402 
403    Not Collective
404 
405    Input Parameter:
406 .  ts - the TS context obtained from TSCreate()
407 
408    Output Parameter:
409 .  iter - number steps so far
410 
411    Level: intermediate
412 
413 .keywords: TS, timestep, get, iteration, number
414 @*/
415 int TSGetTimeStepNumber(TS ts,int* iter)
416 {
417   PetscFunctionBegin;
418   PetscValidHeaderSpecific(ts,TS_COOKIE);
419   *iter = ts->steps;
420   PetscFunctionReturn(0);
421 }
422 
423 #undef __FUNC__
424 #define __FUNC__ "TSSetInitialTimeStep"
425 /*@
426    TSSetInitialTimeStep - Sets the initial timestep to be used,
427    as well as the initial time.
428 
429    Collective on TS
430 
431    Input Parameters:
432 +  ts - the TS context obtained from TSCreate()
433 .  initial_time - the initial time
434 -  time_step - the size of the timestep
435 
436    Level: intermediate
437 
438 .seealso: TSSetTimeStep(), TSGetTimeStep()
439 
440 .keywords: TS, set, initial, timestep
441 @*/
442 int TSSetInitialTimeStep(TS ts,double initial_time,double time_step)
443 {
444   PetscFunctionBegin;
445   PetscValidHeaderSpecific(ts,TS_COOKIE);
446   ts->time_step         = time_step;
447   ts->initial_time_step = time_step;
448   ts->ptime             = initial_time;
449   PetscFunctionReturn(0);
450 }
451 
452 #undef __FUNC__
453 #define __FUNC__ "TSSetTimeStep"
454 /*@
455    TSSetTimeStep - Allows one to reset the timestep at any time,
456    useful for simple pseudo-timestepping codes.
457 
458    Collective on TS
459 
460    Input Parameters:
461 +  ts - the TS context obtained from TSCreate()
462 -  time_step - the size of the timestep
463 
464    Level: intermediate
465 
466 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
467 
468 .keywords: TS, set, timestep
469 @*/
470 int TSSetTimeStep(TS ts,double time_step)
471 {
472   PetscFunctionBegin;
473   PetscValidHeaderSpecific(ts,TS_COOKIE);
474   ts->time_step = time_step;
475   PetscFunctionReturn(0);
476 }
477 
478 #undef __FUNC__
479 #define __FUNC__ "TSGetTimeStep"
480 /*@
481    TSGetTimeStep - Gets the current timestep size.
482 
483    Not Collective
484 
485    Input Parameter:
486 .  ts - the TS context obtained from TSCreate()
487 
488    Output Parameter:
489 .  dt - the current timestep size
490 
491    Level: intermediate
492 
493 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
494 
495 .keywords: TS, get, timestep
496 @*/
497 int TSGetTimeStep(TS ts,double* dt)
498 {
499   PetscFunctionBegin;
500   PetscValidHeaderSpecific(ts,TS_COOKIE);
501   *dt = ts->time_step;
502   PetscFunctionReturn(0);
503 }
504 
505 #undef __FUNC__
506 #define __FUNC__ "TSGetSolution"
507 /*@C
508    TSGetSolution - Returns the solution at the present timestep. It
509    is valid to call this routine inside the function that you are evaluating
510    in order to move to the new timestep. This vector not changed until
511    the solution at the next timestep has been calculated.
512 
513    Not Collective, but Vec returned is parallel if TS is parallel
514 
515    Input Parameter:
516 .  ts - the TS context obtained from TSCreate()
517 
518    Output Parameter:
519 .  v - the vector containing the solution
520 
521    Level: intermediate
522 
523 .seealso: TSGetTimeStep()
524 
525 .keywords: TS, timestep, get, solution
526 @*/
527 int TSGetSolution(TS ts,Vec *v)
528 {
529   PetscFunctionBegin;
530   PetscValidHeaderSpecific(ts,TS_COOKIE);
531   *v = ts->vec_sol_always;
532   PetscFunctionReturn(0);
533 }
534 
535 #undef __FUNC__
536 #define __FUNC__ "TSPublish_Petsc"
537 static int TSPublish_Petsc(PetscObject object)
538 {
539 #if defined(PETSC_HAVE_AMS)
540   TS   v = (TS) object;
541   int  ierr;
542 #endif
543 
544   PetscFunctionBegin;
545 
546 #if defined(PETSC_HAVE_AMS)
547   /* if it is already published then return */
548   if (v->amem >=0 ) PetscFunctionReturn(0);
549 
550   ierr = PetscObjectPublishBaseBegin(object);CHKERRQ(ierr);
551   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Step",&v->steps,1,AMS_INT,AMS_READ,
552                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
553   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Time",&v->ptime,1,AMS_DOUBLE,AMS_READ,
554                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
555   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"CurrentTimeStep",&v->time_step,1,
556                                AMS_DOUBLE,AMS_READ,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
557   ierr = PetscObjectPublishBaseEnd(object);CHKERRQ(ierr);
558 #endif
559   PetscFunctionReturn(0);
560 }
561 
562 /* -----------------------------------------------------------*/
563 
564 #undef __FUNC__
565 #define __FUNC__ "TSCreate"
566 /*@C
567    TSCreate - Creates a timestepper context.
568 
569    Collective on MPI_Comm
570 
571    Input Parameter:
572 +  comm - MPI communicator
573 -  type - One of  TS_LINEAR,TS_NONLINEAR
574    where these types refer to problems of the forms
575 .vb
576          U_t = A U
577          U_t = A(t) U
578          U_t = F(t,U)
579 .ve
580 
581    Output Parameter:
582 .  outts - the new TS context
583 
584    Level: beginner
585 
586 .keywords: TS, timestep, create, context
587 
588 .seealso: TSSetUp(), TSStep(), TSDestroy()
589 @*/
590 int TSCreate(MPI_Comm comm,TSProblemType problemtype,TS *outts)
591 {
592   TS   ts;
593 
594   PetscFunctionBegin;
595   *outts = 0;
596   PetscHeaderCreate(ts,_p_TS,int,TS_COOKIE,-1,"TS",comm,TSDestroy,TSView);
597   PLogObjectCreate(ts);
598   ts->bops->publish     = TSPublish_Petsc;
599   ts->max_steps         = 5000;
600   ts->max_time          = 5.0;
601   ts->time_step         = .1;
602   ts->initial_time_step = ts->time_step;
603   ts->steps             = 0;
604   ts->ptime             = 0.0;
605   ts->data              = 0;
606   ts->view              = 0;
607   ts->setupcalled       = 0;
608   ts->problem_type      = problemtype;
609   ts->numbermonitors    = 0;
610   ts->linear_its        = 0;
611   ts->nonlinear_its     = 0;
612 
613   *outts = ts;
614   PetscFunctionReturn(0);
615 }
616 
617 /* ----- Routines to initialize and destroy a timestepper ---- */
618 
619 #undef __FUNC__
620 #define __FUNC__ "TSSetUp"
621 /*@
622    TSSetUp - Sets up the internal data structures for the later use
623    of a timestepper.
624 
625    Collective on TS
626 
627    Input Parameter:
628 .  ts - the TS context obtained from TSCreate()
629 
630    Notes:
631    For basic use of the TS solvers the user need not explicitly call
632    TSSetUp(), since these actions will automatically occur during
633    the call to TSStep().  However, if one wishes to control this
634    phase separately, TSSetUp() should be called after TSCreate()
635    and optional routines of the form TSSetXXX(), but before TSStep().
636 
637    Level: advanced
638 
639 .keywords: TS, timestep, setup
640 
641 .seealso: TSCreate(), TSStep(), TSDestroy()
642 @*/
643 int TSSetUp(TS ts)
644 {
645   int ierr;
646 
647   PetscFunctionBegin;
648   PetscValidHeaderSpecific(ts,TS_COOKIE);
649   if (!ts->vec_sol) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call TSSetSolution() first");
650   if (!ts->type_name) {
651     ierr = TSSetType(ts,TS_EULER);CHKERRQ(ierr);
652   }
653   ierr = (*ts->setup)(ts);CHKERRQ(ierr);
654   ts->setupcalled = 1;
655   PetscFunctionReturn(0);
656 }
657 
658 #undef __FUNC__
659 #define __FUNC__ "TSDestroy"
660 /*@C
661    TSDestroy - Destroys the timestepper context that was created
662    with TSCreate().
663 
664    Collective on TS
665 
666    Input Parameter:
667 .  ts - the TS context obtained from TSCreate()
668 
669    Level: beginner
670 
671 .keywords: TS, timestepper, destroy
672 
673 .seealso: TSCreate(), TSSetUp(), TSSolve()
674 @*/
675 int TSDestroy(TS ts)
676 {
677   int ierr;
678 
679   PetscFunctionBegin;
680   PetscValidHeaderSpecific(ts,TS_COOKIE);
681   if (--ts->refct > 0) PetscFunctionReturn(0);
682 
683   /* if memory was published with AMS then destroy it */
684   ierr = PetscAMSDestroy(ts);CHKERRQ(ierr);
685 
686   if (ts->sles) {ierr = SLESDestroy(ts->sles);CHKERRQ(ierr);}
687   if (ts->snes) {ierr = SNESDestroy(ts->snes);CHKERRQ(ierr);}
688   ierr = (*(ts)->destroy)(ts);CHKERRQ(ierr);
689   PLogObjectDestroy((PetscObject)ts);
690   PetscHeaderDestroy((PetscObject)ts);
691   PetscFunctionReturn(0);
692 }
693 
694 #undef __FUNC__
695 #define __FUNC__ "TSGetSNES"
696 /*@C
697    TSGetSNES - Returns the SNES (nonlinear solver) associated with
698    a TS (timestepper) context. Valid only for nonlinear problems.
699 
700    Not Collective, but SNES is parallel if TS is parallel
701 
702    Input Parameter:
703 .  ts - the TS context obtained from TSCreate()
704 
705    Output Parameter:
706 .  snes - the nonlinear solver context
707 
708    Notes:
709    The user can then directly manipulate the SNES context to set various
710    options, etc.  Likewise, the user can then extract and manipulate the
711    SLES, KSP, and PC contexts as well.
712 
713    TSGetSNES() does not work for integrators that do not use SNES; in
714    this case TSGetSNES() returns PETSC_NULL in snes.
715 
716    Level: beginner
717 
718 .keywords: timestep, get, SNES
719 @*/
720 int TSGetSNES(TS ts,SNES *snes)
721 {
722   PetscFunctionBegin;
723   PetscValidHeaderSpecific(ts,TS_COOKIE);
724   if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Nonlinear only; use TSGetSLES()");
725   *snes = ts->snes;
726   PetscFunctionReturn(0);
727 }
728 
729 #undef __FUNC__
730 #define __FUNC__ "TSGetSLES"
731 /*@C
732    TSGetSLES - Returns the SLES (linear solver) associated with
733    a TS (timestepper) context.
734 
735    Not Collective, but SLES is parallel if TS is parallel
736 
737    Input Parameter:
738 .  ts - the TS context obtained from TSCreate()
739 
740    Output Parameter:
741 .  sles - the nonlinear solver context
742 
743    Notes:
744    The user can then directly manipulate the SLES context to set various
745    options, etc.  Likewise, the user can then extract and manipulate the
746    KSP and PC contexts as well.
747 
748    TSGetSLES() does not work for integrators that do not use SLES;
749    in this case TSGetSLES() returns PETSC_NULL in sles.
750 
751    Level: beginner
752 
753 .keywords: timestep, get, SLES
754 @*/
755 int TSGetSLES(TS ts,SLES *sles)
756 {
757   PetscFunctionBegin;
758   PetscValidHeaderSpecific(ts,TS_COOKIE);
759   if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Linear only; use TSGetSNES()");
760   *sles = ts->sles;
761   PetscFunctionReturn(0);
762 }
763 
764 /* ----------- Routines to set solver parameters ---------- */
765 
766 #undef __FUNC__
767 #define __FUNC__ "TSSetDuration"
768 /*@
769    TSSetDuration - Sets the maximum number of timesteps to use and
770    maximum time for iteration.
771 
772    Collective on TS
773 
774    Input Parameters:
775 +  ts - the TS context obtained from TSCreate()
776 .  maxsteps - maximum number of iterations to use
777 -  maxtime - final time to iterate to
778 
779    Options Database Keys:
780 .  -ts_max_steps <maxsteps> - Sets maxsteps
781 .  -ts_max_time <maxtime> - Sets maxtime
782 
783    Notes:
784    The default maximum number of iterations is 5000. Default time is 5.0
785 
786    Level: intermediate
787 
788 .keywords: TS, timestep, set, maximum, iterations
789 @*/
790 int TSSetDuration(TS ts,int maxsteps,double maxtime)
791 {
792   PetscFunctionBegin;
793   PetscValidHeaderSpecific(ts,TS_COOKIE);
794   ts->max_steps = maxsteps;
795   ts->max_time  = maxtime;
796   PetscFunctionReturn(0);
797 }
798 
799 #undef __FUNC__
800 #define __FUNC__ "TSSetSolution"
801 /*@
802    TSSetSolution - Sets the initial solution vector
803    for use by the TS routines.
804 
805    Collective on TS and Vec
806 
807    Input Parameters:
808 +  ts - the TS context obtained from TSCreate()
809 -  x - the solution vector
810 
811    Level: beginner
812 
813 .keywords: TS, timestep, set, solution, initial conditions
814 @*/
815 int TSSetSolution(TS ts,Vec x)
816 {
817   PetscFunctionBegin;
818   PetscValidHeaderSpecific(ts,TS_COOKIE);
819   ts->vec_sol        = ts->vec_sol_always = x;
820   PetscFunctionReturn(0);
821 }
822 
823 /* ------------ Routines to set performance monitoring options ----------- */
824 
825 #undef __FUNC__
826 #define __FUNC__ "TSSetMonitor"
827 /*@C
828    TSSetMonitor - Sets an ADDITIONAL function that is to be used at every
829    timestep to display the iteration's  progress.
830 
831    Collective on TS
832 
833    Input Parameters:
834 +  ts - the TS context obtained from TSCreate()
835 .  func - monitoring routine
836 -  mctx - [optional] user-defined context for private data for the
837           monitor routine (may be PETSC_NULL)
838 
839    Calling sequence of func:
840 $    int func(TS ts,int steps,double time,Vec x,void *mctx)
841 
842 +    ts - the TS context
843 .    steps - iteration number
844 .    time - current timestep
845 .    x - current iterate
846 -    mctx - [optional] monitoring context
847 
848    Notes:
849    This routine adds an additional monitor to the list of monitors that
850    already has been loaded.
851 
852    Level: intermediate
853 
854 .keywords: TS, timestep, set, monitor
855 
856 .seealso: TSDefaultMonitor(), TSClearMonitor()
857 @*/
858 int TSSetMonitor(TS ts, int (*monitor)(TS,int,double,Vec,void*), void *mctx )
859 {
860   PetscFunctionBegin;
861   PetscValidHeaderSpecific(ts,TS_COOKIE);
862   if (ts->numbermonitors >= MAXTSMONITORS) {
863     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many monitors set");
864   }
865   ts->monitor[ts->numbermonitors]           = monitor;
866   ts->monitorcontext[ts->numbermonitors++]  = (void*)mctx;
867   PetscFunctionReturn(0);
868 }
869 
870 #undef __FUNC__
871 #define __FUNC__ "TSClearMonitor"
872 /*@C
873    TSClearMonitor - Clears all the monitors that have been set on a time-step object.
874 
875    Collective on TS
876 
877    Input Parameters:
878 .  ts - the TS context obtained from TSCreate()
879 
880    Notes:
881    There is no way to remove a single, specific monitor.
882 
883    Level: intermediate
884 
885 .keywords: TS, timestep, set, monitor
886 
887 .seealso: TSDefaultMonitor(), TSSetMonitor()
888 @*/
889 int TSClearMonitor(TS ts)
890 {
891   PetscFunctionBegin;
892   PetscValidHeaderSpecific(ts,TS_COOKIE);
893   ts->numbermonitors = 0;
894   PetscFunctionReturn(0);
895 }
896 
897 #undef __FUNC__
898 #define __FUNC__ "TSDefaultMonitor"
899 int TSDefaultMonitor(TS ts, int step, double time,Vec v, void *ctx)
900 {
901   int ierr;
902 
903   PetscFunctionBegin;
904   ierr = PetscPrintf(ts->comm,"timestep %d dt %g time %g\n",step,ts->time_step,time);CHKERRQ(ierr);
905   PetscFunctionReturn(0);
906 }
907 
908 #undef __FUNC__
909 #define __FUNC__ "TSStep"
910 /*@
911    TSStep - Steps the requested number of timesteps.
912 
913    Collective on TS
914 
915    Input Parameter:
916 .  ts - the TS context obtained from TSCreate()
917 
918    Output Parameters:
919 +  steps - number of iterations until termination
920 -  time - time until termination
921 
922    Level: beginner
923 
924 .keywords: TS, timestep, solve
925 
926 .seealso: TSCreate(), TSSetUp(), TSDestroy()
927 @*/
928 int TSStep(TS ts,int *steps,double *time)
929 {
930   int ierr,flg;
931 
932   PetscFunctionBegin;
933   PetscValidHeaderSpecific(ts,TS_COOKIE);
934   if (!ts->setupcalled) {ierr = TSSetUp(ts);CHKERRQ(ierr);}
935   PLogEventBegin(TS_Step,ts,0,0,0);
936   ierr = (*(ts)->step)(ts,steps,time);CHKERRQ(ierr);
937   PLogEventEnd(TS_Step,ts,0,0,0);
938   ierr = OptionsHasName(PETSC_NULL,"-ts_view",&flg);CHKERRQ(ierr);
939   if (flg) {
940     ierr = TSView(ts,VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
941   }
942   PetscFunctionReturn(0);
943 }
944 
945 #undef __FUNC__
946 #define __FUNC__ "TSMonitor"
947 /*
948      Runs the user provided monitor routines, if they exists.
949 */
950 int TSMonitor(TS ts,int step,double time,Vec x)
951 {
952   int i,ierr,n = ts->numbermonitors;
953 
954   PetscFunctionBegin;
955   for ( i=0; i<n; i++ ) {
956     ierr = (*ts->monitor[i])(ts,step,time,x,ts->monitorcontext[i]);CHKERRQ(ierr);
957   }
958   PetscFunctionReturn(0);
959 }
960 
961 /* ------------------------------------------------------------------------*/
962 
963 #undef __FUNC__
964 #define __FUNC__ "TSLGMonitorCreate"
965 /*@C
966    TSLGMonitorCreate - Creates a line graph context for use with
967    TS to monitor convergence of preconditioned residual norms.
968 
969    Collective on TS
970 
971    Input Parameters:
972 +  host - the X display to open, or null for the local machine
973 .  label - the title to put in the title bar
974 .  x, y - the screen coordinates of the upper left coordinate of
975           the window
976 -  m, n - the screen width and height in pixels
977 
978    Output Parameter:
979 .  draw - the drawing context
980 
981    Options Database Key:
982 .  -ts_xmonitor - automatically sets line graph monitor
983 
984    Notes:
985    Use TSLGMonitorDestroy() to destroy this line graph, not DrawLGDestroy().
986 
987    Level: intermediate
988 
989 .keywords: TS, monitor, line graph, residual, create
990 
991 .seealso: TSLGMonitorDestroy(), TSSetMonitor()
992 @*/
993 int TSLGMonitorCreate(char *host,char *label,int x,int y,int m,
994                        int n, DrawLG *draw)
995 {
996   Draw win;
997   int  ierr;
998 
999   PetscFunctionBegin;
1000   ierr = DrawOpenX(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr);
1001   ierr = DrawLGCreate(win,1,draw);CHKERRQ(ierr);
1002   ierr = DrawLGIndicateDataPoints(*draw);CHKERRQ(ierr);
1003 
1004   PLogObjectParent(*draw,win);
1005   PetscFunctionReturn(0);
1006 }
1007 
1008 #undef __FUNC__
1009 #define __FUNC__ "TSLGMonitor"
1010 int TSLGMonitor(TS ts,int n,double time,Vec v,void *monctx)
1011 {
1012   DrawLG lg = (DrawLG) monctx;
1013   double x,y = time;
1014   int    ierr;
1015 
1016   PetscFunctionBegin;
1017   if (!n) {ierr = DrawLGReset(lg);CHKERRQ(ierr);}
1018   x = (double) n;
1019   ierr = DrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1020   if (n < 20 || (n % 5)) {
1021     ierr = DrawLGDraw(lg);CHKERRQ(ierr);
1022   }
1023   PetscFunctionReturn(0);
1024 }
1025 
1026 #undef __FUNC__
1027 #define __FUNC__ "TSLGMonitorDestroy"
1028 /*@C
1029    TSLGMonitorDestroy - Destroys a line graph context that was created
1030    with TSLGMonitorCreate().
1031 
1032    Collective on DrawLG
1033 
1034    Input Parameter:
1035 .  draw - the drawing context
1036 
1037    Level: intermediate
1038 
1039 .keywords: TS, monitor, line graph, destroy
1040 
1041 .seealso: TSLGMonitorCreate(),  TSSetMonitor(), TSLGMonitor();
1042 @*/
1043 int TSLGMonitorDestroy(DrawLG drawlg)
1044 {
1045   Draw draw;
1046   int  ierr;
1047 
1048   PetscFunctionBegin;
1049   ierr = DrawLGGetDraw(drawlg,&draw);CHKERRQ(ierr);
1050   ierr = DrawDestroy(draw);CHKERRQ(ierr);
1051   ierr = DrawLGDestroy(drawlg);CHKERRQ(ierr);
1052   PetscFunctionReturn(0);
1053 }
1054 
1055 #undef __FUNC__
1056 #define __FUNC__ "TSGetTime"
1057 /*@
1058    TSGetTime - Gets the current time.
1059 
1060    Not Collective
1061 
1062    Input Parameter:
1063 .  ts - the TS context obtained from TSCreate()
1064 
1065    Output Parameter:
1066 .  t  - the current time
1067 
1068    Contributed by: Matthew Knepley
1069 
1070    Level: beginner
1071 
1072 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1073 
1074 .keywords: TS, get, time
1075 @*/
1076 int TSGetTime(TS ts, double* t)
1077 {
1078   PetscFunctionBegin;
1079   PetscValidHeaderSpecific(ts, TS_COOKIE);
1080   *t = ts->ptime;
1081   PetscFunctionReturn(0);
1082 }
1083 
1084 #undef __FUNC__
1085 #define __FUNC__ "TSGetProblemType"
1086 /*@C
1087    TSGetProblemType - Returns the problem type of a TS (timestepper) context.
1088 
1089    Not Collective
1090 
1091    Input Parameter:
1092 .  ts   - The TS context obtained from TSCreate()
1093 
1094    Output Parameter:
1095 .  type - The problem type, TS_LINEAR or TS_NONLINEAR
1096 
1097    Level: intermediate
1098 
1099    Contributed by: Matthew Knepley
1100 
1101 .keywords: ts, get, type
1102 
1103 @*/
1104 int TSGetProblemType(TS ts, TSProblemType *type)
1105 {
1106   PetscFunctionBegin;
1107   PetscValidHeaderSpecific(ts, TS_COOKIE);
1108   *type = ts->problem_type;
1109   PetscFunctionReturn(0);
1110 }
1111 
1112 #undef __FUNC__
1113 #define __FUNC__ "TSSetOptionsPrefix"
1114 /*@C
1115    TSSetOptionsPrefix - Sets the prefix used for searching for all
1116    TS options in the database.
1117 
1118    Collective on TS
1119 
1120    Input Parameter:
1121 +  ts     - The TS context
1122 -  prefix - The prefix to prepend to all option names
1123 
1124    Notes:
1125    A hyphen (-) must NOT be given at the beginning of the prefix name.
1126    The first character of all runtime options is AUTOMATICALLY the
1127    hyphen.
1128 
1129    Contributed by: Matthew Knepley
1130 
1131    Level: advanced
1132 
1133 .keywords: TS, set, options, prefix, database
1134 
1135 .seealso: TSSetFromOptions()
1136 
1137 @*/
1138 int TSSetOptionsPrefix(TS ts, char *prefix)
1139 {
1140   int ierr;
1141 
1142   PetscFunctionBegin;
1143   PetscValidHeaderSpecific(ts, TS_COOKIE);
1144   ierr = PetscObjectSetOptionsPrefix((PetscObject) ts, prefix);CHKERRQ(ierr);
1145   switch(ts->problem_type) {
1146     case TS_NONLINEAR:
1147       ierr = SNESSetOptionsPrefix(ts->snes, prefix);CHKERRQ(ierr);
1148       break;
1149     case TS_LINEAR:
1150       ierr = SLESSetOptionsPrefix(ts->sles, prefix);CHKERRQ(ierr);
1151       break;
1152   }
1153   PetscFunctionReturn(0);
1154 }
1155 
1156 
1157 #undef __FUNC__
1158 #define __FUNC__ "TSAppendOptionsPrefix"
1159 /*@C
1160    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
1161    TS options in the database.
1162 
1163    Collective on TS
1164 
1165    Input Parameter:
1166 +  ts     - The TS context
1167 -  prefix - The prefix to prepend to all option names
1168 
1169    Notes:
1170    A hyphen (-) must NOT be given at the beginning of the prefix name.
1171    The first character of all runtime options is AUTOMATICALLY the
1172    hyphen.
1173 
1174    Contributed by: Matthew Knepley
1175 
1176    Level: advanced
1177 
1178 .keywords: TS, append, options, prefix, database
1179 
1180 .seealso: TSGetOptionsPrefix()
1181 
1182 @*/
1183 int TSAppendOptionsPrefix(TS ts, char *prefix)
1184 {
1185   int ierr;
1186 
1187   PetscFunctionBegin;
1188   PetscValidHeaderSpecific(ts, TS_COOKIE);
1189   ierr = PetscObjectAppendOptionsPrefix((PetscObject) ts, prefix);CHKERRQ(ierr);
1190   switch(ts->problem_type) {
1191     case TS_NONLINEAR:
1192       ierr = SNESAppendOptionsPrefix(ts->snes, prefix);CHKERRQ(ierr);
1193       break;
1194     case TS_LINEAR:
1195       ierr = SLESAppendOptionsPrefix(ts->sles, prefix);CHKERRQ(ierr);
1196       break;
1197   }
1198   PetscFunctionReturn(0);
1199 }
1200 
1201 #undef __FUNC__
1202 #define __FUNC__ "TSGetOptionsPrefix"
1203 /*@C
1204    TSGetOptionsPrefix - Sets the prefix used for searching for all
1205    TS options in the database.
1206 
1207    Not Collective
1208 
1209    Input Parameter:
1210 .  ts - The TS context
1211 
1212    Output Parameter:
1213 .  prefix - A pointer to the prefix string used
1214 
1215    Contributed by: Matthew Knepley
1216 
1217    Notes: On the fortran side, the user should pass in a string 'prifix' of
1218    sufficient length to hold the prefix.
1219 
1220    Level: intermediate
1221 
1222 .keywords: TS, get, options, prefix, database
1223 
1224 .seealso: TSAppendOptionsPrefix()
1225 @*/
1226 int TSGetOptionsPrefix(TS ts, char **prefix)
1227 {
1228   int ierr;
1229 
1230   PetscFunctionBegin;
1231   PetscValidHeaderSpecific(ts, TS_COOKIE);
1232   ierr = PetscObjectGetOptionsPrefix((PetscObject) ts, prefix);CHKERRQ(ierr);
1233   PetscFunctionReturn(0);
1234 }
1235 
1236 #undef __FUNC__
1237 #define __FUNC__ "TSGetRHSMatrix"
1238 /*@C
1239    TSGetRHSMatrix - Returns the matrix A at the present timestep.
1240 
1241    Not Collective, but parallel objects are returned if TS is parallel
1242 
1243    Input Parameter:
1244 .  ts  - The TS context obtained from TSCreate()
1245 
1246    Output Parameters:
1247 +  A   - The matrix A, where U_t = A(t) U
1248 .  M   - The preconditioner matrix, usually the same as A
1249 -  ctx - User-defined context for matrix evaluation routine
1250 
1251    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1252 
1253    Contributed by: Matthew Knepley
1254 
1255    Level: intermediate
1256 
1257 .seealso: TSGetTimeStep(), TSGetTime(), TSGetTimeStepNumber(), TSGetRHSJacobian()
1258 
1259 .keywords: TS, timestep, get, matrix
1260 
1261 @*/
1262 int TSGetRHSMatrix(TS ts, Mat *A, Mat *M, void **ctx)
1263 {
1264   PetscFunctionBegin;
1265   PetscValidHeaderSpecific(ts, TS_COOKIE);
1266   if (A)   *A = ts->A;
1267   if (M)   *M = ts->B;
1268   if (ctx) *ctx = ts->jacP;
1269   PetscFunctionReturn(0);
1270 }
1271 
1272 #undef __FUNC__
1273 #define __FUNC__ "TSGetRHSJacobian"
1274 /*@C
1275    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
1276 
1277    Not Collective, but parallel objects are returned if TS is parallel
1278 
1279    Input Parameter:
1280 .  ts  - The TS context obtained from TSCreate()
1281 
1282    Output Parameters:
1283 +  J   - The Jacobian J of F, where U_t = F(U,t)
1284 .  M   - The preconditioner matrix, usually the same as J
1285 - ctx - User-defined context for Jacobian evaluation routine
1286 
1287    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1288 
1289    Contributed by: Matthew Knepley
1290 
1291    Level: intermediate
1292 
1293 .seealso: TSGetTimeStep(), TSGetRHSMatrix(), TSGetTime(), TSGetTimeStepNumber()
1294 
1295 .keywords: TS, timestep, get, matrix, Jacobian
1296 @*/
1297 int TSGetRHSJacobian(TS ts, Mat *J, Mat *M, void **ctx)
1298 {
1299   int ierr;
1300 
1301   PetscFunctionBegin;
1302   ierr = TSGetRHSMatrix(ts, J, M, ctx);CHKERRQ(ierr);
1303   PetscFunctionReturn(0);
1304 }
1305 
1306 /*MC
1307    TSRegister - Adds a method to the timestepping solver package.
1308 
1309    Synopsis:
1310 
1311    TSRegister(char *name_solver,char *path,char *name_create,int (*routine_create)(TS))
1312 
1313    Not collective
1314 
1315    Input Parameters:
1316 +  name_solver - name of a new user-defined solver
1317 .  path - path (either absolute or relative) the library containing this solver
1318 .  name_create - name of routine to create method context
1319 -  routine_create - routine to create method context
1320 
1321    Notes:
1322    TSRegister() may be called multiple times to add several user-defined solvers.
1323 
1324    If dynamic libraries are used, then the fourth input argument (routine_create)
1325    is ignored.
1326 
1327    Sample usage:
1328 .vb
1329    TSRegister("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
1330               "MySolverCreate",MySolverCreate);
1331 .ve
1332 
1333    Then, your solver can be chosen with the procedural interface via
1334 $     TSSetType(ts,"my_solver")
1335    or at runtime via the option
1336 $     -ts_type my_solver
1337 
1338    Level: advanced
1339 
1340    $PETSC_ARCH, $PETSC_DIR, $PETSC_LDIR, and $BOPT occuring in pathname will be replaced with appropriate values.
1341 
1342 .keywords: TS, register
1343 
1344 .seealso: TSRegisterAll(), TSRegisterDestroy()
1345 M*/
1346 
1347 #undef __FUNC__
1348 #define __FUNC__ "TSRegister_Private"
1349 int TSRegister_Private(char *sname,char *path,char *name,int (*function)(TS))
1350 {
1351   char fullname[256];
1352   int  ierr;
1353 
1354   PetscFunctionBegin;
1355   ierr = PetscStrcpy(fullname,path);CHKERRQ(ierr);
1356   PetscStrcat(fullname,":"); PetscStrcat(fullname,name);
1357   FListAdd_Private(&TSList,sname,fullname,        (int (*)(void*))function);
1358   PetscFunctionReturn(0);
1359 }
1360