xref: /petsc/src/ts/interface/ts.c (revision fcca9d3d5c45103cafcd0f6bc9babe194d935ced)
1 /* $Id: ts.c,v 1.43 2001/09/07 20:12:01 bsmith Exp $ */
2 #include "src/ts/tsimpl.h"        /*I "petscts.h"  I*/
3 
4 /* Logging support */
5 int TS_COOKIE = 0;
6 int TS_Step = 0, TS_PseudoComputeTimeStep = 0, TS_FunctionEval = 0, TS_JacobianEval = 0;
7 
8 #undef __FUNCT__
9 #define __FUNCT__ "TSSetTypeFromOptions"
10 /*
11   TSSetTypeFromOptions - Sets the type of ts from user options.
12 
13   Collective on TS
14 
15   Input Parameter:
16 . ts - The ts
17 
18   Level: intermediate
19 
20 .keywords: TS, set, options, database, type
21 .seealso: TSSetFromOptions(), TSSetType()
22 */
23 static int TSSetTypeFromOptions(TS ts)
24 {
25   PetscTruth opt;
26   const char *defaultType;
27   char       typeName[256];
28   int        ierr;
29 
30   PetscFunctionBegin;
31   if (ts->type_name != PETSC_NULL) {
32     defaultType = ts->type_name;
33   } else {
34     defaultType = TS_EULER;
35   }
36 
37   if (!TSRegisterAllCalled) {
38     ierr = TSRegisterAll(PETSC_NULL);                                                                     CHKERRQ(ierr);
39   }
40   ierr = PetscOptionsList("-ts_type", "TS method"," TSSetType", TSList, defaultType, typeName, 256, &opt);CHKERRQ(ierr);
41   if (opt == PETSC_TRUE) {
42     ierr = TSSetType(ts, typeName);                                                                       CHKERRQ(ierr);
43   } else {
44     ierr = TSSetType(ts, defaultType);                                                                    CHKERRQ(ierr);
45   }
46   PetscFunctionReturn(0);
47 }
48 
49 #undef __FUNCT__
50 #define __FUNCT__ "TSSetFromOptions"
51 /*@
52    TSSetFromOptions - Sets various TS parameters from user options.
53 
54    Collective on TS
55 
56    Input Parameter:
57 .  ts - the TS context obtained from TSCreate()
58 
59    Options Database Keys:
60 +  -ts_type <type> - TS_EULER, TS_BEULER, TS_PVODE, TS_PSEUDO, TS_CRANK_NICHOLSON
61 .  -ts_max_steps maxsteps - maximum number of time-steps to take
62 .  -ts_max_time time - maximum time to compute to
63 .  -ts_dt dt - initial time step
64 .  -ts_monitor - print information at each timestep
65 -  -ts_xmonitor - plot information at each timestep
66 
67    Level: beginner
68 
69 .keywords: TS, timestep, set, options, database
70 
71 .seealso: TSGetType
72 @*/
73 int TSSetFromOptions(TS ts)
74 {
75   PetscReal  dt;
76   PetscTruth opt;
77   int        ierr;
78 
79   PetscFunctionBegin;
80   PetscValidHeaderSpecific(ts, TS_COOKIE);
81   ierr = PetscOptionsBegin(ts->comm, ts->prefix, "Time step options", "TS");                              CHKERRQ(ierr);
82 
83   /* Handle generic TS options */
84   ierr = PetscOptionsInt("-ts_max_steps","Maximum number of time steps","TSSetDuration",ts->max_steps,&ts->max_steps,PETSC_NULL);CHKERRQ(ierr);
85   ierr = PetscOptionsReal("-ts_max_time","Time to run to","TSSetDuration",ts->max_time,&ts->max_time,PETSC_NULL);CHKERRQ(ierr);
86   ierr = PetscOptionsReal("-ts_init_time","Initial time","TSSetInitialTime", ts->ptime, &ts->ptime, PETSC_NULL);CHKERRQ(ierr);
87   ierr = PetscOptionsReal("-ts_dt","Initial time step","TSSetInitialTimeStep",ts->initial_time_step,&dt,&opt);CHKERRQ(ierr);
88   if (opt == PETSC_TRUE) {
89     ts->initial_time_step = ts->time_step = dt;
90   }
91 
92   /* Monitor options */
93     ierr = PetscOptionsName("-ts_monitor","Monitor timestep size","TSDefaultMonitor",&opt);               CHKERRQ(ierr);
94     if (opt == PETSC_TRUE) {
95       ierr = TSSetMonitor(ts,TSDefaultMonitor,PETSC_NULL,PETSC_NULL);                                     CHKERRQ(ierr);
96     }
97     ierr = PetscOptionsName("-ts_xmonitor","Monitor timestep size graphically","TSLGMonitor",&opt);       CHKERRQ(ierr);
98     if (opt == PETSC_TRUE) {
99       ierr = TSSetMonitor(ts,TSLGMonitor,PETSC_NULL,PETSC_NULL);                                          CHKERRQ(ierr);
100     }
101     ierr = PetscOptionsName("-ts_vecmonitor","Monitor solution graphically","TSVecViewMonitor",&opt);     CHKERRQ(ierr);
102     if (opt == PETSC_TRUE) {
103       ierr = TSSetMonitor(ts,TSVecViewMonitor,PETSC_NULL,PETSC_NULL);                                     CHKERRQ(ierr);
104     }
105 
106   /* Handle TS type options */
107   ierr = TSSetTypeFromOptions(ts);                                                                        CHKERRQ(ierr);
108 
109   /* Handle specific TS options */
110   if (ts->ops->setfromoptions != PETSC_NULL) {
111     ierr = (*ts->ops->setfromoptions)(ts);                                                                CHKERRQ(ierr);
112   }
113   ierr = PetscOptionsEnd();                                                                               CHKERRQ(ierr);
114 
115   /* Handle subobject options */
116   switch(ts->problem_type) {
117     /* Should check for implicit/explicit */
118   case TS_LINEAR:
119     if (ts->ksp != PETSC_NULL) {
120       ierr = KSPSetFromOptions(ts->ksp);                                                                CHKERRQ(ierr);
121     }
122     break;
123   case TS_NONLINEAR:
124     if (ts->snes != PETSC_NULL) {
125       ierr = SNESSetFromOptions(ts->snes);                                                                CHKERRQ(ierr);
126     }
127     break;
128   default:
129     SETERRQ1(PETSC_ERR_ARG_WRONG, "Invalid problem type: %d", ts->problem_type);
130   }
131 
132   PetscFunctionReturn(0);
133 }
134 
135 #undef  __FUNCT__
136 #define __FUNCT__ "TSViewFromOptions"
137 /*@
138   TSViewFromOptions - This function visualizes the ts based upon user options.
139 
140   Collective on TS
141 
142   Input Parameter:
143 . ts - The ts
144 
145   Level: intermediate
146 
147 .keywords: TS, view, options, database
148 .seealso: TSSetFromOptions(), TSView()
149 @*/
150 int TSViewFromOptions(TS ts,const char title[])
151 {
152   PetscViewer viewer;
153   PetscDraw   draw;
154   PetscTruth  opt;
155   char        typeName[1024];
156   char        fileName[PETSC_MAX_PATH_LEN];
157   int         len;
158   int         ierr;
159 
160   PetscFunctionBegin;
161   ierr = PetscOptionsHasName(ts->prefix, "-ts_view", &opt);                                               CHKERRQ(ierr);
162   if (opt == PETSC_TRUE) {
163     ierr = PetscOptionsGetString(ts->prefix, "-ts_view", typeName, 1024, &opt);                           CHKERRQ(ierr);
164     ierr = PetscStrlen(typeName, &len);                                                                   CHKERRQ(ierr);
165     if (len > 0) {
166       ierr = PetscViewerCreate(ts->comm, &viewer);                                                        CHKERRQ(ierr);
167       ierr = PetscViewerSetType(viewer, typeName);                                                        CHKERRQ(ierr);
168       ierr = PetscOptionsGetString(ts->prefix, "-ts_view_file", fileName, 1024, &opt);                    CHKERRQ(ierr);
169       if (opt == PETSC_TRUE) {
170         ierr = PetscViewerSetFilename(viewer, fileName);                                                  CHKERRQ(ierr);
171       } else {
172         ierr = PetscViewerSetFilename(viewer, ts->name);                                                  CHKERRQ(ierr);
173       }
174       ierr = TSView(ts, viewer);                                                                          CHKERRQ(ierr);
175       ierr = PetscViewerFlush(viewer);                                                                    CHKERRQ(ierr);
176       ierr = PetscViewerDestroy(viewer);                                                                  CHKERRQ(ierr);
177     } else {
178       ierr = TSView(ts, PETSC_NULL);                                                                      CHKERRQ(ierr);
179     }
180   }
181   ierr = PetscOptionsHasName(ts->prefix, "-ts_view_draw", &opt);                                          CHKERRQ(ierr);
182   if (opt == PETSC_TRUE) {
183     ierr = PetscViewerDrawOpen(ts->comm, 0, 0, 0, 0, 300, 300, &viewer);                                  CHKERRQ(ierr);
184     ierr = PetscViewerDrawGetDraw(viewer, 0, &draw);                                                      CHKERRQ(ierr);
185     if (title != PETSC_NULL) {
186       ierr = PetscDrawSetTitle(draw, (char *)title);                                                      CHKERRQ(ierr);
187     } else {
188       ierr = PetscObjectName((PetscObject) ts);                                                           CHKERRQ(ierr) ;
189       ierr = PetscDrawSetTitle(draw, ts->name);                                                           CHKERRQ(ierr);
190     }
191     ierr = TSView(ts, viewer);                                                                            CHKERRQ(ierr);
192     ierr = PetscViewerFlush(viewer);                                                                      CHKERRQ(ierr);
193     ierr = PetscDrawPause(draw);                                                                          CHKERRQ(ierr);
194     ierr = PetscViewerDestroy(viewer);                                                                    CHKERRQ(ierr);
195   }
196   PetscFunctionReturn(0);
197 }
198 
199 #undef __FUNCT__
200 #define __FUNCT__ "TSComputeRHSJacobian"
201 /*@
202    TSComputeRHSJacobian - Computes the Jacobian matrix that has been
203       set with TSSetRHSJacobian().
204 
205    Collective on TS and Vec
206 
207    Input Parameters:
208 +  ts - the SNES context
209 .  t - current timestep
210 -  x - input vector
211 
212    Output Parameters:
213 +  A - Jacobian matrix
214 .  B - optional preconditioning matrix
215 -  flag - flag indicating matrix structure
216 
217    Notes:
218    Most users should not need to explicitly call this routine, as it
219    is used internally within the nonlinear solvers.
220 
221    See KSPSetOperators() for important information about setting the
222    flag parameter.
223 
224    TSComputeJacobian() is valid only for TS_NONLINEAR
225 
226    Level: developer
227 
228 .keywords: SNES, compute, Jacobian, matrix
229 
230 .seealso:  TSSetRHSJacobian(), KSPSetOperators()
231 @*/
232 int TSComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg)
233 {
234   int ierr;
235 
236   PetscFunctionBegin;
237   PetscValidHeaderSpecific(ts,TS_COOKIE);
238   PetscValidHeaderSpecific(X,VEC_COOKIE);
239   PetscCheckSameComm(ts,X);
240   if (ts->problem_type != TS_NONLINEAR) {
241     SETERRQ(PETSC_ERR_ARG_WRONG,"For TS_NONLINEAR only");
242   }
243   if (ts->ops->rhsjacobian) {
244     ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
245     *flg = DIFFERENT_NONZERO_PATTERN;
246     PetscStackPush("TS user Jacobian function");
247     ierr = (*ts->ops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr);
248     PetscStackPop;
249     ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
250     /* make sure user returned a correct Jacobian and preconditioner */
251     PetscValidHeaderSpecific(*A,MAT_COOKIE);
252     PetscValidHeaderSpecific(*B,MAT_COOKIE);
253   } else {
254     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
255     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
256     if (*A != *B) {
257       ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
258       ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
259     }
260   }
261   PetscFunctionReturn(0);
262 }
263 
264 #undef __FUNCT__
265 #define __FUNCT__ "TSComputeRHSFunction"
266 /*
267    TSComputeRHSFunction - Evaluates the right-hand-side function.
268 
269    Note: If the user did not provide a function but merely a matrix,
270    this routine applies the matrix.
271 */
272 int TSComputeRHSFunction(TS ts,PetscReal t,Vec x,Vec y)
273 {
274   int ierr;
275 
276   PetscFunctionBegin;
277   PetscValidHeaderSpecific(ts,TS_COOKIE);
278   PetscValidHeader(x);
279   PetscValidHeader(y);
280 
281   ierr = PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
282   if (ts->ops->rhsfunction) {
283     PetscStackPush("TS user right-hand-side function");
284     ierr = (*ts->ops->rhsfunction)(ts,t,x,y,ts->funP);CHKERRQ(ierr);
285     PetscStackPop;
286   } else {
287     if (ts->ops->rhsmatrix) { /* assemble matrix for this timestep */
288       MatStructure flg;
289       PetscStackPush("TS user right-hand-side matrix function");
290       ierr = (*ts->ops->rhsmatrix)(ts,t,&ts->A,&ts->B,&flg,ts->jacP);CHKERRQ(ierr);
291       PetscStackPop;
292     }
293     ierr = MatMult(ts->A,x,y);CHKERRQ(ierr);
294   }
295 
296   /* apply user-provided boundary conditions (only needed if these are time dependent) */
297   ierr = TSComputeRHSBoundaryConditions(ts,t,y);CHKERRQ(ierr);
298   ierr = PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
299 
300   PetscFunctionReturn(0);
301 }
302 
303 #undef __FUNCT__
304 #define __FUNCT__ "TSSetRHSFunction"
305 /*@C
306     TSSetRHSFunction - Sets the routine for evaluating the function,
307     F(t,u), where U_t = F(t,u).
308 
309     Collective on TS
310 
311     Input Parameters:
312 +   ts - the TS context obtained from TSCreate()
313 .   f - routine for evaluating the right-hand-side function
314 -   ctx - [optional] user-defined context for private data for the
315           function evaluation routine (may be PETSC_NULL)
316 
317     Calling sequence of func:
318 $     func (TS ts,PetscReal t,Vec u,Vec F,void *ctx);
319 
320 +   t - current timestep
321 .   u - input vector
322 .   F - function vector
323 -   ctx - [optional] user-defined function context
324 
325     Important:
326     The user MUST call either this routine or TSSetRHSMatrix().
327 
328     Level: beginner
329 
330 .keywords: TS, timestep, set, right-hand-side, function
331 
332 .seealso: TSSetRHSMatrix()
333 @*/
334 int TSSetRHSFunction(TS ts,int (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx)
335 {
336   PetscFunctionBegin;
337 
338   PetscValidHeaderSpecific(ts,TS_COOKIE);
339   if (ts->problem_type == TS_LINEAR) {
340     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot set function for linear problem");
341   }
342   ts->ops->rhsfunction = f;
343   ts->funP             = ctx;
344   PetscFunctionReturn(0);
345 }
346 
347 #undef __FUNCT__
348 #define __FUNCT__ "TSSetRHSMatrix"
349 /*@C
350    TSSetRHSMatrix - Sets the function to compute the matrix A, where U_t = A(t) U.
351    Also sets the location to store A.
352 
353    Collective on TS
354 
355    Input Parameters:
356 +  ts  - the TS context obtained from TSCreate()
357 .  A   - matrix
358 .  B   - preconditioner matrix (usually same as A)
359 .  f   - the matrix evaluation routine; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran)
360          if A is not a function of t.
361 -  ctx - [optional] user-defined context for private data for the
362           matrix evaluation routine (may be PETSC_NULL)
363 
364    Calling sequence of func:
365 $     func (TS ts,PetscReal t,Mat *A,Mat *B,int *flag,void *ctx);
366 
367 +  t - current timestep
368 .  A - matrix A, where U_t = A(t) U
369 .  B - preconditioner matrix, usually the same as A
370 .  flag - flag indicating information about the preconditioner matrix
371           structure (same as flag in KSPSetOperators())
372 -  ctx - [optional] user-defined context for matrix evaluation routine
373 
374    Notes:
375    See KSPSetOperators() for important information about setting the flag
376    output parameter in the routine func().  Be sure to read this information!
377 
378    The routine func() takes Mat * as the matrix arguments rather than Mat.
379    This allows the matrix evaluation routine to replace A and/or B with a
380    completely new new matrix structure (not just different matrix elements)
381    when appropriate, for instance, if the nonzero structure is changing
382    throughout the global iterations.
383 
384    Important:
385    The user MUST call either this routine or TSSetRHSFunction().
386 
387    Level: beginner
388 
389 .keywords: TS, timestep, set, right-hand-side, matrix
390 
391 .seealso: TSSetRHSFunction()
392 @*/
393 int TSSetRHSMatrix(TS ts,Mat A,Mat B,int (*f)(TS,PetscReal,Mat*,Mat*,MatStructure*,void*),void *ctx)
394 {
395   PetscFunctionBegin;
396   PetscValidHeaderSpecific(ts,TS_COOKIE);
397   PetscValidHeaderSpecific(A,MAT_COOKIE);
398   PetscValidHeaderSpecific(B,MAT_COOKIE);
399   PetscCheckSameComm(ts,A);
400   PetscCheckSameComm(ts,B);
401   if (ts->problem_type == TS_NONLINEAR) {
402     SETERRQ(PETSC_ERR_ARG_WRONG,"Not for nonlinear problems; use TSSetRHSJacobian()");
403   }
404 
405   ts->ops->rhsmatrix = f;
406   ts->jacP           = ctx;
407   ts->A              = A;
408   ts->B              = B;
409 
410   PetscFunctionReturn(0);
411 }
412 
413 #undef __FUNCT__
414 #define __FUNCT__ "TSSetRHSJacobian"
415 /*@C
416    TSSetRHSJacobian - Sets the function to compute the Jacobian of F,
417    where U_t = F(U,t), as well as the location to store the matrix.
418 
419    Collective on TS
420 
421    Input Parameters:
422 +  ts  - the TS context obtained from TSCreate()
423 .  A   - Jacobian matrix
424 .  B   - preconditioner matrix (usually same as A)
425 .  f   - the Jacobian evaluation routine
426 -  ctx - [optional] user-defined context for private data for the
427          Jacobian evaluation routine (may be PETSC_NULL)
428 
429    Calling sequence of func:
430 $     func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx);
431 
432 +  t - current timestep
433 .  u - input vector
434 .  A - matrix A, where U_t = A(t)u
435 .  B - preconditioner matrix, usually the same as A
436 .  flag - flag indicating information about the preconditioner matrix
437           structure (same as flag in KSPSetOperators())
438 -  ctx - [optional] user-defined context for matrix evaluation routine
439 
440    Notes:
441    See KSPSetOperators() for important information about setting the flag
442    output parameter in the routine func().  Be sure to read this information!
443 
444    The routine func() takes Mat * as the matrix arguments rather than Mat.
445    This allows the matrix evaluation routine to replace A and/or B with a
446    completely new new matrix structure (not just different matrix elements)
447    when appropriate, for instance, if the nonzero structure is changing
448    throughout the global iterations.
449 
450    Level: beginner
451 
452 .keywords: TS, timestep, set, right-hand-side, Jacobian
453 
454 .seealso: TSDefaultComputeJacobianColor(),
455           SNESDefaultComputeJacobianColor()
456 
457 @*/
458 int TSSetRHSJacobian(TS ts,Mat A,Mat B,int (*f)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
459 {
460   PetscFunctionBegin;
461   PetscValidHeaderSpecific(ts,TS_COOKIE);
462   PetscValidHeaderSpecific(A,MAT_COOKIE);
463   PetscValidHeaderSpecific(B,MAT_COOKIE);
464   PetscCheckSameComm(ts,A);
465   PetscCheckSameComm(ts,B);
466   if (ts->problem_type != TS_NONLINEAR) {
467     SETERRQ(PETSC_ERR_ARG_WRONG,"Not for linear problems; use TSSetRHSMatrix()");
468   }
469 
470   ts->ops->rhsjacobian = f;
471   ts->jacP             = ctx;
472   ts->A                = A;
473   ts->B                = B;
474   PetscFunctionReturn(0);
475 }
476 
477 #undef __FUNCT__
478 #define __FUNCT__ "TSComputeRHSBoundaryConditions"
479 /*
480    TSComputeRHSBoundaryConditions - Evaluates the boundary condition function.
481 
482    Note: If the user did not provide a function but merely a matrix,
483    this routine applies the matrix.
484 */
485 int TSComputeRHSBoundaryConditions(TS ts,PetscReal t,Vec x)
486 {
487   int ierr;
488 
489   PetscFunctionBegin;
490   PetscValidHeaderSpecific(ts,TS_COOKIE);
491   PetscValidHeader(x);
492   PetscCheckSameComm(ts,x);
493 
494   if (ts->ops->rhsbc) {
495     PetscStackPush("TS user boundary condition function");
496     ierr = (*ts->ops->rhsbc)(ts,t,x,ts->bcP);CHKERRQ(ierr);
497     PetscStackPop;
498     PetscFunctionReturn(0);
499   }
500 
501   PetscFunctionReturn(0);
502 }
503 
504 #undef __FUNCT__
505 #define __FUNCT__ "TSSetRHSBoundaryConditions"
506 /*@C
507     TSSetRHSBoundaryConditions - Sets the routine for evaluating the function,
508     boundary conditions for the function F.
509 
510     Collective on TS
511 
512     Input Parameters:
513 +   ts - the TS context obtained from TSCreate()
514 .   f - routine for evaluating the boundary condition function
515 -   ctx - [optional] user-defined context for private data for the
516           function evaluation routine (may be PETSC_NULL)
517 
518     Calling sequence of func:
519 $     func (TS ts,PetscReal t,Vec F,void *ctx);
520 
521 +   t - current timestep
522 .   F - function vector
523 -   ctx - [optional] user-defined function context
524 
525     Level: intermediate
526 
527 .keywords: TS, timestep, set, boundary conditions, function
528 @*/
529 int TSSetRHSBoundaryConditions(TS ts,int (*f)(TS,PetscReal,Vec,void*),void *ctx)
530 {
531   PetscFunctionBegin;
532 
533   PetscValidHeaderSpecific(ts,TS_COOKIE);
534   if (ts->problem_type != TS_LINEAR) {
535     SETERRQ(PETSC_ERR_ARG_WRONG,"For linear problems only");
536   }
537   ts->ops->rhsbc = f;
538   ts->bcP        = ctx;
539   PetscFunctionReturn(0);
540 }
541 
542 #undef __FUNCT__
543 #define __FUNCT__ "TSView"
544 /*@C
545     TSView - Prints the TS data structure.
546 
547     Collective on TS
548 
549     Input Parameters:
550 +   ts - the TS context obtained from TSCreate()
551 -   viewer - visualization context
552 
553     Options Database Key:
554 .   -ts_view - calls TSView() at end of TSStep()
555 
556     Notes:
557     The available visualization contexts include
558 +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
559 -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
560          output where only the first processor opens
561          the file.  All other processors send their
562          data to the first processor to print.
563 
564     The user can open an alternative visualization context with
565     PetscViewerASCIIOpen() - output to a specified file.
566 
567     Level: beginner
568 
569 .keywords: TS, timestep, view
570 
571 .seealso: PetscViewerASCIIOpen()
572 @*/
573 int TSView(TS ts,PetscViewer viewer)
574 {
575   int        ierr;
576   char       *type;
577   PetscTruth isascii,isstring;
578 
579   PetscFunctionBegin;
580   PetscValidHeaderSpecific(ts,TS_COOKIE);
581   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(ts->comm);
582   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE);
583   PetscCheckSameComm(ts,viewer);
584 
585   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
586   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
587   if (isascii) {
588     ierr = PetscViewerASCIIPrintf(viewer,"TS Object:\n");CHKERRQ(ierr);
589     ierr = TSGetType(ts,(TSType *)&type);CHKERRQ(ierr);
590     if (type) {
591       ierr = PetscViewerASCIIPrintf(viewer,"  type: %s\n",type);CHKERRQ(ierr);
592     } else {
593       ierr = PetscViewerASCIIPrintf(viewer,"  type: not yet set\n");CHKERRQ(ierr);
594     }
595     if (ts->ops->view) {
596       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
597       ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr);
598       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
599     }
600     ierr = PetscViewerASCIIPrintf(viewer,"  maximum steps=%d\n",ts->max_steps);CHKERRQ(ierr);
601     ierr = PetscViewerASCIIPrintf(viewer,"  maximum time=%g\n",ts->max_time);CHKERRQ(ierr);
602     if (ts->problem_type == TS_NONLINEAR) {
603       ierr = PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solver iterations=%d\n",ts->nonlinear_its);CHKERRQ(ierr);
604     }
605     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%d\n",ts->linear_its);CHKERRQ(ierr);
606   } else if (isstring) {
607     ierr = TSGetType(ts,(TSType *)&type);CHKERRQ(ierr);
608     ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr);
609   }
610   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
611   if (ts->ksp) {ierr = KSPView(ts->ksp,viewer);CHKERRQ(ierr);}
612   if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
613   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
614   PetscFunctionReturn(0);
615 }
616 
617 
618 #undef __FUNCT__
619 #define __FUNCT__ "TSSetApplicationContext"
620 /*@C
621    TSSetApplicationContext - Sets an optional user-defined context for
622    the timesteppers.
623 
624    Collective on TS
625 
626    Input Parameters:
627 +  ts - the TS context obtained from TSCreate()
628 -  usrP - optional user context
629 
630    Level: intermediate
631 
632 .keywords: TS, timestep, set, application, context
633 
634 .seealso: TSGetApplicationContext()
635 @*/
636 int TSSetApplicationContext(TS ts,void *usrP)
637 {
638   PetscFunctionBegin;
639   PetscValidHeaderSpecific(ts,TS_COOKIE);
640   ts->user = usrP;
641   PetscFunctionReturn(0);
642 }
643 
644 #undef __FUNCT__
645 #define __FUNCT__ "TSGetApplicationContext"
646 /*@C
647     TSGetApplicationContext - Gets the user-defined context for the
648     timestepper.
649 
650     Not Collective
651 
652     Input Parameter:
653 .   ts - the TS context obtained from TSCreate()
654 
655     Output Parameter:
656 .   usrP - user context
657 
658     Level: intermediate
659 
660 .keywords: TS, timestep, get, application, context
661 
662 .seealso: TSSetApplicationContext()
663 @*/
664 int TSGetApplicationContext(TS ts,void **usrP)
665 {
666   PetscFunctionBegin;
667   PetscValidHeaderSpecific(ts,TS_COOKIE);
668   *usrP = ts->user;
669   PetscFunctionReturn(0);
670 }
671 
672 #undef __FUNCT__
673 #define __FUNCT__ "TSGetTimeStepNumber"
674 /*@
675    TSGetTimeStepNumber - Gets the current number of timesteps.
676 
677    Not Collective
678 
679    Input Parameter:
680 .  ts - the TS context obtained from TSCreate()
681 
682    Output Parameter:
683 .  iter - number steps so far
684 
685    Level: intermediate
686 
687 .keywords: TS, timestep, get, iteration, number
688 @*/
689 int TSGetTimeStepNumber(TS ts,int* iter)
690 {
691   PetscFunctionBegin;
692   PetscValidHeaderSpecific(ts,TS_COOKIE);
693   *iter = ts->steps;
694   PetscFunctionReturn(0);
695 }
696 
697 #undef __FUNCT__
698 #define __FUNCT__ "TSSetInitialTimeStep"
699 /*@
700    TSSetInitialTimeStep - Sets the initial timestep to be used,
701    as well as the initial time.
702 
703    Collective on TS
704 
705    Input Parameters:
706 +  ts - the TS context obtained from TSCreate()
707 .  initial_time - the initial time
708 -  time_step - the size of the timestep
709 
710    Level: intermediate
711 
712 .seealso: TSSetTimeStep(), TSGetTimeStep()
713 
714 .keywords: TS, set, initial, timestep
715 @*/
716 int TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step)
717 {
718   PetscFunctionBegin;
719   PetscValidHeaderSpecific(ts,TS_COOKIE);
720   ts->time_step         = time_step;
721   ts->initial_time_step = time_step;
722   ts->ptime             = initial_time;
723   PetscFunctionReturn(0);
724 }
725 
726 #undef __FUNCT__
727 #define __FUNCT__ "TSSetTimeStep"
728 /*@
729    TSSetTimeStep - Allows one to reset the timestep at any time,
730    useful for simple pseudo-timestepping codes.
731 
732    Collective on TS
733 
734    Input Parameters:
735 +  ts - the TS context obtained from TSCreate()
736 -  time_step - the size of the timestep
737 
738    Level: intermediate
739 
740 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
741 
742 .keywords: TS, set, timestep
743 @*/
744 int TSSetTimeStep(TS ts,PetscReal time_step)
745 {
746   PetscFunctionBegin;
747   PetscValidHeaderSpecific(ts,TS_COOKIE);
748   ts->time_step = time_step;
749   PetscFunctionReturn(0);
750 }
751 
752 #undef __FUNCT__
753 #define __FUNCT__ "TSGetTimeStep"
754 /*@
755    TSGetTimeStep - Gets the current timestep size.
756 
757    Not Collective
758 
759    Input Parameter:
760 .  ts - the TS context obtained from TSCreate()
761 
762    Output Parameter:
763 .  dt - the current timestep size
764 
765    Level: intermediate
766 
767 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
768 
769 .keywords: TS, get, timestep
770 @*/
771 int TSGetTimeStep(TS ts,PetscReal* dt)
772 {
773   PetscFunctionBegin;
774   PetscValidHeaderSpecific(ts,TS_COOKIE);
775   *dt = ts->time_step;
776   PetscFunctionReturn(0);
777 }
778 
779 #undef __FUNCT__
780 #define __FUNCT__ "TSGetSolution"
781 /*@C
782    TSGetSolution - Returns the solution at the present timestep. It
783    is valid to call this routine inside the function that you are evaluating
784    in order to move to the new timestep. This vector not changed until
785    the solution at the next timestep has been calculated.
786 
787    Not Collective, but Vec returned is parallel if TS is parallel
788 
789    Input Parameter:
790 .  ts - the TS context obtained from TSCreate()
791 
792    Output Parameter:
793 .  v - the vector containing the solution
794 
795    Level: intermediate
796 
797 .seealso: TSGetTimeStep()
798 
799 .keywords: TS, timestep, get, solution
800 @*/
801 int TSGetSolution(TS ts,Vec *v)
802 {
803   PetscFunctionBegin;
804   PetscValidHeaderSpecific(ts,TS_COOKIE);
805   *v = ts->vec_sol_always;
806   PetscFunctionReturn(0);
807 }
808 
809 /* ----- Routines to initialize and destroy a timestepper ---- */
810 #undef __FUNCT__
811 #define __FUNCT__ "TSSetProblemType"
812 /*@C
813   TSSetProblemType - Sets the type of problem to be solved.
814 
815   Not collective
816 
817   Input Parameters:
818 + ts   - The TS
819 - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
820 .vb
821          U_t = A U
822          U_t = A(t) U
823          U_t = F(t,U)
824 .ve
825 
826    Level: beginner
827 
828 .keywords: TS, problem type
829 .seealso: TSSetUp(), TSProblemType, TS
830 @*/
831 int TSSetProblemType(TS ts, TSProblemType type) {
832   PetscFunctionBegin;
833   PetscValidHeaderSpecific(ts, TS_COOKIE);
834   ts->problem_type = type;
835   PetscFunctionReturn(0);
836 }
837 
838 #undef __FUNCT__
839 #define __FUNCT__ "TSGetProblemType"
840 /*@C
841   TSGetProblemType - Gets the type of problem to be solved.
842 
843   Not collective
844 
845   Input Parameter:
846 . ts   - The TS
847 
848   Output Parameter:
849 . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
850 .vb
851          U_t = A U
852          U_t = A(t) U
853          U_t = F(t,U)
854 .ve
855 
856    Level: beginner
857 
858 .keywords: TS, problem type
859 .seealso: TSSetUp(), TSProblemType, TS
860 @*/
861 int TSGetProblemType(TS ts, TSProblemType *type) {
862   PetscFunctionBegin;
863   PetscValidHeaderSpecific(ts, TS_COOKIE);
864   PetscValidPointer(type);
865   *type = ts->problem_type;
866   PetscFunctionReturn(0);
867 }
868 
869 #undef __FUNCT__
870 #define __FUNCT__ "TSSetUp"
871 /*@
872    TSSetUp - Sets up the internal data structures for the later use
873    of a timestepper.
874 
875    Collective on TS
876 
877    Input Parameter:
878 .  ts - the TS context obtained from TSCreate()
879 
880    Notes:
881    For basic use of the TS solvers the user need not explicitly call
882    TSSetUp(), since these actions will automatically occur during
883    the call to TSStep().  However, if one wishes to control this
884    phase separately, TSSetUp() should be called after TSCreate()
885    and optional routines of the form TSSetXXX(), but before TSStep().
886 
887    Level: advanced
888 
889 .keywords: TS, timestep, setup
890 
891 .seealso: TSCreate(), TSStep(), TSDestroy()
892 @*/
893 int TSSetUp(TS ts)
894 {
895   int ierr;
896 
897   PetscFunctionBegin;
898   PetscValidHeaderSpecific(ts,TS_COOKIE);
899   if (!ts->vec_sol) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first");
900   if (!ts->type_name) {
901     ierr = TSSetType(ts,TS_EULER);CHKERRQ(ierr);
902   }
903   ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr);
904   ts->setupcalled = 1;
905   PetscFunctionReturn(0);
906 }
907 
908 #undef __FUNCT__
909 #define __FUNCT__ "TSDestroy"
910 /*@C
911    TSDestroy - Destroys the timestepper context that was created
912    with TSCreate().
913 
914    Collective on TS
915 
916    Input Parameter:
917 .  ts - the TS context obtained from TSCreate()
918 
919    Level: beginner
920 
921 .keywords: TS, timestepper, destroy
922 
923 .seealso: TSCreate(), TSSetUp(), TSSolve()
924 @*/
925 int TSDestroy(TS ts)
926 {
927   int ierr,i;
928 
929   PetscFunctionBegin;
930   PetscValidHeaderSpecific(ts,TS_COOKIE);
931   if (--ts->refct > 0) PetscFunctionReturn(0);
932 
933   /* if memory was published with AMS then destroy it */
934   ierr = PetscObjectDepublish(ts);CHKERRQ(ierr);
935 
936   if (ts->ksp) {ierr = KSPDestroy(ts->ksp);CHKERRQ(ierr);}
937   if (ts->snes) {ierr = SNESDestroy(ts->snes);CHKERRQ(ierr);}
938   ierr = (*(ts)->ops->destroy)(ts);CHKERRQ(ierr);
939   for (i=0; i<ts->numbermonitors; i++) {
940     if (ts->mdestroy[i]) {
941       ierr = (*ts->mdestroy[i])(ts->monitorcontext[i]);CHKERRQ(ierr);
942     }
943   }
944   PetscLogObjectDestroy((PetscObject)ts);
945   PetscHeaderDestroy((PetscObject)ts);
946   PetscFunctionReturn(0);
947 }
948 
949 #undef __FUNCT__
950 #define __FUNCT__ "TSGetSNES"
951 /*@C
952    TSGetSNES - Returns the SNES (nonlinear solver) associated with
953    a TS (timestepper) context. Valid only for nonlinear problems.
954 
955    Not Collective, but SNES is parallel if TS is parallel
956 
957    Input Parameter:
958 .  ts - the TS context obtained from TSCreate()
959 
960    Output Parameter:
961 .  snes - the nonlinear solver context
962 
963    Notes:
964    The user can then directly manipulate the SNES context to set various
965    options, etc.  Likewise, the user can then extract and manipulate the
966    KSP, KSP, and PC contexts as well.
967 
968    TSGetSNES() does not work for integrators that do not use SNES; in
969    this case TSGetSNES() returns PETSC_NULL in snes.
970 
971    Level: beginner
972 
973 .keywords: timestep, get, SNES
974 @*/
975 int TSGetSNES(TS ts,SNES *snes)
976 {
977   PetscFunctionBegin;
978   PetscValidHeaderSpecific(ts,TS_COOKIE);
979   if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Nonlinear only; use TSGetKSP()");
980   *snes = ts->snes;
981   PetscFunctionReturn(0);
982 }
983 
984 #undef __FUNCT__
985 #define __FUNCT__ "TSGetKSP"
986 /*@C
987    TSGetKSP - Returns the KSP (linear solver) associated with
988    a TS (timestepper) context.
989 
990    Not Collective, but KSP is parallel if TS is parallel
991 
992    Input Parameter:
993 .  ts - the TS context obtained from TSCreate()
994 
995    Output Parameter:
996 .  ksp - the nonlinear solver context
997 
998    Notes:
999    The user can then directly manipulate the KSP context to set various
1000    options, etc.  Likewise, the user can then extract and manipulate the
1001    KSP and PC contexts as well.
1002 
1003    TSGetKSP() does not work for integrators that do not use KSP;
1004    in this case TSGetKSP() returns PETSC_NULL in ksp.
1005 
1006    Level: beginner
1007 
1008 .keywords: timestep, get, KSP
1009 @*/
1010 int TSGetKSP(TS ts,KSP *ksp)
1011 {
1012   PetscFunctionBegin;
1013   PetscValidHeaderSpecific(ts,TS_COOKIE);
1014   if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()");
1015   *ksp = ts->ksp;
1016   PetscFunctionReturn(0);
1017 }
1018 
1019 /* ----------- Routines to set solver parameters ---------- */
1020 
1021 #undef __FUNCT__
1022 #define __FUNCT__ "TSGetDuration"
1023 /*@
1024    TSGetDuration - Gets the maximum number of timesteps to use and
1025    maximum time for iteration.
1026 
1027    Collective on TS
1028 
1029    Input Parameters:
1030 +  ts       - the TS context obtained from TSCreate()
1031 .  maxsteps - maximum number of iterations to use, or PETSC_NULL
1032 -  maxtime  - final time to iterate to, or PETSC_NULL
1033 
1034    Level: intermediate
1035 
1036 .keywords: TS, timestep, get, maximum, iterations, time
1037 @*/
1038 int TSGetDuration(TS ts, int *maxsteps, PetscReal *maxtime)
1039 {
1040   PetscFunctionBegin;
1041   PetscValidHeaderSpecific(ts, TS_COOKIE);
1042   if (maxsteps != PETSC_NULL) {
1043     PetscValidIntPointer(maxsteps);
1044     *maxsteps = ts->max_steps;
1045   }
1046   if (maxtime  != PETSC_NULL) {
1047     PetscValidScalarPointer(maxtime);
1048     *maxtime  = ts->max_time;
1049   }
1050   PetscFunctionReturn(0);
1051 }
1052 
1053 #undef __FUNCT__
1054 #define __FUNCT__ "TSSetDuration"
1055 /*@
1056    TSSetDuration - Sets the maximum number of timesteps to use and
1057    maximum time for iteration.
1058 
1059    Collective on TS
1060 
1061    Input Parameters:
1062 +  ts - the TS context obtained from TSCreate()
1063 .  maxsteps - maximum number of iterations to use
1064 -  maxtime - final time to iterate to
1065 
1066    Options Database Keys:
1067 .  -ts_max_steps <maxsteps> - Sets maxsteps
1068 .  -ts_max_time <maxtime> - Sets maxtime
1069 
1070    Notes:
1071    The default maximum number of iterations is 5000. Default time is 5.0
1072 
1073    Level: intermediate
1074 
1075 .keywords: TS, timestep, set, maximum, iterations
1076 @*/
1077 int TSSetDuration(TS ts,int maxsteps,PetscReal maxtime)
1078 {
1079   PetscFunctionBegin;
1080   PetscValidHeaderSpecific(ts,TS_COOKIE);
1081   ts->max_steps = maxsteps;
1082   ts->max_time  = maxtime;
1083   PetscFunctionReturn(0);
1084 }
1085 
1086 #undef __FUNCT__
1087 #define __FUNCT__ "TSSetSolution"
1088 /*@
1089    TSSetSolution - Sets the initial solution vector
1090    for use by the TS routines.
1091 
1092    Collective on TS and Vec
1093 
1094    Input Parameters:
1095 +  ts - the TS context obtained from TSCreate()
1096 -  x - the solution vector
1097 
1098    Level: beginner
1099 
1100 .keywords: TS, timestep, set, solution, initial conditions
1101 @*/
1102 int TSSetSolution(TS ts,Vec x)
1103 {
1104   PetscFunctionBegin;
1105   PetscValidHeaderSpecific(ts,TS_COOKIE);
1106   ts->vec_sol        = ts->vec_sol_always = x;
1107   PetscFunctionReturn(0);
1108 }
1109 
1110 #undef __FUNCT__
1111 #define __FUNCT__ "TSSetRhsBC"
1112 /*@C
1113   TSSetRhsBC - Sets the function which applies boundary conditions
1114   to the Rhs of each system.
1115 
1116   Collective on TS
1117 
1118   Input Parameters:
1119 + ts   - The TS context obtained from TSCreate()
1120 - func - The function
1121 
1122   Calling sequence of func:
1123 . func (TS ts, Vec rhs, void *ctx);
1124 
1125 + rhs - The current rhs vector
1126 - ctx - The user-context
1127 
1128   Level: intermediate
1129 
1130 .keywords: TS, Rhs, boundary conditions
1131 @*/
1132 int TSSetRhsBC(TS ts, int (*func)(TS, Vec, void *))
1133 {
1134   PetscFunctionBegin;
1135   PetscValidHeaderSpecific(ts, TS_COOKIE);
1136   ts->ops->applyrhsbc = func;
1137   PetscFunctionReturn(0);
1138 }
1139 
1140 #undef __FUNCT__
1141 #define __FUNCT__ "TSDefaultRhsBC"
1142 /*@
1143   TSDefaultRhsBC - The default boundary condition function which does nothing.
1144 
1145   Collective on TS
1146 
1147   Input Parameters:
1148 + ts  - The TS context obtained from TSCreate()
1149 . rhs - The Rhs
1150 - ctx - The user-context
1151 
1152   Level: developer
1153 
1154 .keywords: TS, Rhs, boundary conditions
1155 @*/
1156 int TSDefaultRhsBC(TS ts,  Vec rhs, void *ctx)
1157 {
1158   PetscFunctionBegin;
1159   PetscFunctionReturn(0);
1160 }
1161 
1162 #undef __FUNCT__
1163 #define __FUNCT__ "TSSetSystemMatrixBC"
1164 /*@C
1165   TSSetSystemMatrixBC - Sets the function which applies boundary conditions
1166   to the system matrix and preconditioner of each system.
1167 
1168   Collective on TS
1169 
1170   Input Parameters:
1171 + ts   - The TS context obtained from TSCreate()
1172 - func - The function
1173 
1174   Calling sequence of func:
1175 . func (TS ts, Mat A, Mat B, void *ctx);
1176 
1177 + A   - The current system matrix
1178 . B   - The current preconditioner
1179 - ctx - The user-context
1180 
1181   Level: intermediate
1182 
1183 .keywords: TS, System matrix, boundary conditions
1184 @*/
1185 int TSSetSystemMatrixBC(TS ts, int (*func)(TS, Mat, Mat, void *))
1186 {
1187   PetscFunctionBegin;
1188   PetscValidHeaderSpecific(ts, TS_COOKIE);
1189   ts->ops->applymatrixbc = func;
1190   PetscFunctionReturn(0);
1191 }
1192 
1193 #undef __FUNCT__
1194 #define __FUNCT__ "TSDefaultSystemMatrixBC"
1195 /*@
1196   TSDefaultSystemMatrixBC - The default boundary condition function which
1197   does nothing.
1198 
1199   Collective on TS
1200 
1201   Input Parameters:
1202 + ts  - The TS context obtained from TSCreate()
1203 . A   - The system matrix
1204 . B   - The preconditioner
1205 - ctx - The user-context
1206 
1207   Level: developer
1208 
1209 .keywords: TS, System matrix, boundary conditions
1210 @*/
1211 int TSDefaultSystemMatrixBC(TS ts, Mat A, Mat B, void *ctx)
1212 {
1213   PetscFunctionBegin;
1214   PetscFunctionReturn(0);
1215 }
1216 
1217 #undef __FUNCT__
1218 #define __FUNCT__ "TSSetSolutionBC"
1219 /*@C
1220   TSSetSolutionBC - Sets the function which applies boundary conditions
1221   to the solution of each system. This is necessary in nonlinear systems
1222   which time dependent boundary conditions.
1223 
1224   Collective on TS
1225 
1226   Input Parameters:
1227 + ts   - The TS context obtained from TSCreate()
1228 - func - The function
1229 
1230   Calling sequence of func:
1231 . func (TS ts, Vec rsol, void *ctx);
1232 
1233 + sol - The current solution vector
1234 - ctx - The user-context
1235 
1236   Level: intermediate
1237 
1238 .keywords: TS, solution, boundary conditions
1239 @*/
1240 int TSSetSolutionBC(TS ts, int (*func)(TS, Vec, void *))
1241 {
1242   PetscFunctionBegin;
1243   PetscValidHeaderSpecific(ts, TS_COOKIE);
1244   ts->ops->applysolbc = func;
1245   PetscFunctionReturn(0);
1246 }
1247 
1248 #undef __FUNCT__
1249 #define __FUNCT__ "TSDefaultSolutionBC"
1250 /*@
1251   TSDefaultSolutionBC - The default boundary condition function which
1252   does nothing.
1253 
1254   Collective on TS
1255 
1256   Input Parameters:
1257 + ts  - The TS context obtained from TSCreate()
1258 . sol - The solution
1259 - ctx - The user-context
1260 
1261   Level: developer
1262 
1263 .keywords: TS, solution, boundary conditions
1264 @*/
1265 int TSDefaultSolutionBC(TS ts, Vec sol, void *ctx)
1266 {
1267   PetscFunctionBegin;
1268   PetscFunctionReturn(0);
1269 }
1270 
1271 #undef __FUNCT__
1272 #define __FUNCT__ "TSSetPreStep"
1273 /*@C
1274   TSSetPreStep - Sets the general-purpose function
1275   called once at the beginning of time stepping.
1276 
1277   Collective on TS
1278 
1279   Input Parameters:
1280 + ts   - The TS context obtained from TSCreate()
1281 - func - The function
1282 
1283   Calling sequence of func:
1284 . func (TS ts);
1285 
1286   Level: intermediate
1287 
1288 .keywords: TS, timestep
1289 @*/
1290 int TSSetPreStep(TS ts, int (*func)(TS))
1291 {
1292   PetscFunctionBegin;
1293   PetscValidHeaderSpecific(ts, TS_COOKIE);
1294   ts->ops->prestep = func;
1295   PetscFunctionReturn(0);
1296 }
1297 
1298 #undef __FUNCT__
1299 #define __FUNCT__ "TSDefaultPreStep"
1300 /*@
1301   TSDefaultPreStep - The default pre-stepping function which does nothing.
1302 
1303   Collective on TS
1304 
1305   Input Parameters:
1306 . ts  - The TS context obtained from TSCreate()
1307 
1308   Level: developer
1309 
1310 .keywords: TS, timestep
1311 @*/
1312 int TSDefaultPreStep(TS ts)
1313 {
1314   PetscFunctionBegin;
1315   PetscFunctionReturn(0);
1316 }
1317 
1318 #undef __FUNCT__
1319 #define __FUNCT__ "TSSetUpdate"
1320 /*@C
1321   TSSetUpdate - Sets the general-purpose update function called
1322   at the beginning of every time step. This function can change
1323   the time step.
1324 
1325   Collective on TS
1326 
1327   Input Parameters:
1328 + ts   - The TS context obtained from TSCreate()
1329 - func - The function
1330 
1331   Calling sequence of func:
1332 . func (TS ts, double t, double *dt);
1333 
1334 + t   - The current time
1335 - dt  - The current time step
1336 
1337   Level: intermediate
1338 
1339 .keywords: TS, update, timestep
1340 @*/
1341 int TSSetUpdate(TS ts, int (*func)(TS, PetscReal, PetscReal *))
1342 {
1343   PetscFunctionBegin;
1344   PetscValidHeaderSpecific(ts, TS_COOKIE);
1345   ts->ops->update = func;
1346   PetscFunctionReturn(0);
1347 }
1348 
1349 #undef __FUNCT__
1350 #define __FUNCT__ "TSDefaultUpdate"
1351 /*@
1352   TSDefaultUpdate - The default update function which does nothing.
1353 
1354   Collective on TS
1355 
1356   Input Parameters:
1357 + ts  - The TS context obtained from TSCreate()
1358 - t   - The current time
1359 
1360   Output Parameters:
1361 . dt  - The current time step
1362 
1363   Level: developer
1364 
1365 .keywords: TS, update, timestep
1366 @*/
1367 int TSDefaultUpdate(TS ts, PetscReal t, PetscReal *dt)
1368 {
1369   PetscFunctionBegin;
1370   PetscFunctionReturn(0);
1371 }
1372 
1373 #undef __FUNCT__
1374 #define __FUNCT__ "TSSetPostStep"
1375 /*@C
1376   TSSetPostStep - Sets the general-purpose function
1377   called once at the end of time stepping.
1378 
1379   Collective on TS
1380 
1381   Input Parameters:
1382 + ts   - The TS context obtained from TSCreate()
1383 - func - The function
1384 
1385   Calling sequence of func:
1386 . func (TS ts);
1387 
1388   Level: intermediate
1389 
1390 .keywords: TS, timestep
1391 @*/
1392 int TSSetPostStep(TS ts, int (*func)(TS))
1393 {
1394   PetscFunctionBegin;
1395   PetscValidHeaderSpecific(ts, TS_COOKIE);
1396   ts->ops->poststep = func;
1397   PetscFunctionReturn(0);
1398 }
1399 
1400 #undef __FUNCT__
1401 #define __FUNCT__ "TSDefaultPostStep"
1402 /*@
1403   TSDefaultPostStep - The default post-stepping function which does nothing.
1404 
1405   Collective on TS
1406 
1407   Input Parameters:
1408 . ts  - The TS context obtained from TSCreate()
1409 
1410   Level: developer
1411 
1412 .keywords: TS, timestep
1413 @*/
1414 int TSDefaultPostStep(TS ts)
1415 {
1416   PetscFunctionBegin;
1417   PetscFunctionReturn(0);
1418 }
1419 
1420 /* ------------ Routines to set performance monitoring options ----------- */
1421 
1422 #undef __FUNCT__
1423 #define __FUNCT__ "TSSetMonitor"
1424 /*@C
1425    TSSetMonitor - Sets an ADDITIONAL function that is to be used at every
1426    timestep to display the iteration's  progress.
1427 
1428    Collective on TS
1429 
1430    Input Parameters:
1431 +  ts - the TS context obtained from TSCreate()
1432 .  func - monitoring routine
1433 .  mctx - [optional] user-defined context for private data for the
1434              monitor routine (use PETSC_NULL if no context is desired)
1435 -  monitordestroy - [optional] routine that frees monitor context
1436           (may be PETSC_NULL)
1437 
1438    Calling sequence of func:
1439 $    int func(TS ts,int steps,PetscReal time,Vec x,void *mctx)
1440 
1441 +    ts - the TS context
1442 .    steps - iteration number
1443 .    time - current time
1444 .    x - current iterate
1445 -    mctx - [optional] monitoring context
1446 
1447    Notes:
1448    This routine adds an additional monitor to the list of monitors that
1449    already has been loaded.
1450 
1451    Level: intermediate
1452 
1453 .keywords: TS, timestep, set, monitor
1454 
1455 .seealso: TSDefaultMonitor(), TSClearMonitor()
1456 @*/
1457 int TSSetMonitor(TS ts,int (*monitor)(TS,int,PetscReal,Vec,void*),void *mctx,int (*mdestroy)(void*))
1458 {
1459   PetscFunctionBegin;
1460   PetscValidHeaderSpecific(ts,TS_COOKIE);
1461   if (ts->numbermonitors >= MAXTSMONITORS) {
1462     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1463   }
1464   ts->monitor[ts->numbermonitors]           = monitor;
1465   ts->mdestroy[ts->numbermonitors]          = mdestroy;
1466   ts->monitorcontext[ts->numbermonitors++]  = (void*)mctx;
1467   PetscFunctionReturn(0);
1468 }
1469 
1470 #undef __FUNCT__
1471 #define __FUNCT__ "TSClearMonitor"
1472 /*@C
1473    TSClearMonitor - Clears all the monitors that have been set on a time-step object.
1474 
1475    Collective on TS
1476 
1477    Input Parameters:
1478 .  ts - the TS context obtained from TSCreate()
1479 
1480    Notes:
1481    There is no way to remove a single, specific monitor.
1482 
1483    Level: intermediate
1484 
1485 .keywords: TS, timestep, set, monitor
1486 
1487 .seealso: TSDefaultMonitor(), TSSetMonitor()
1488 @*/
1489 int TSClearMonitor(TS ts)
1490 {
1491   PetscFunctionBegin;
1492   PetscValidHeaderSpecific(ts,TS_COOKIE);
1493   ts->numbermonitors = 0;
1494   PetscFunctionReturn(0);
1495 }
1496 
1497 #undef __FUNCT__
1498 #define __FUNCT__ "TSDefaultMonitor"
1499 int TSDefaultMonitor(TS ts,int step,PetscReal ptime,Vec v,void *ctx)
1500 {
1501   int ierr;
1502 
1503   PetscFunctionBegin;
1504   ierr = PetscPrintf(ts->comm,"timestep %d dt %g time %g\n",step,ts->time_step,ptime);CHKERRQ(ierr);
1505   PetscFunctionReturn(0);
1506 }
1507 
1508 #undef __FUNCT__
1509 #define __FUNCT__ "TSStep"
1510 /*@
1511    TSStep - Steps the requested number of timesteps.
1512 
1513    Collective on TS
1514 
1515    Input Parameter:
1516 .  ts - the TS context obtained from TSCreate()
1517 
1518    Output Parameters:
1519 +  steps - number of iterations until termination
1520 -  ptime - time until termination
1521 
1522    Level: beginner
1523 
1524 .keywords: TS, timestep, solve
1525 
1526 .seealso: TSCreate(), TSSetUp(), TSDestroy()
1527 @*/
1528 int TSStep(TS ts,int *steps,PetscReal *ptime)
1529 {
1530   int        ierr;
1531 
1532   PetscFunctionBegin;
1533   PetscValidHeaderSpecific(ts, TS_COOKIE);
1534   if (!ts->setupcalled) {
1535     ierr = TSSetUp(ts);                                                                                   CHKERRQ(ierr);
1536   }
1537 
1538   ierr = PetscLogEventBegin(TS_Step, ts, 0, 0, 0);                                                        CHKERRQ(ierr);
1539   ierr = (*ts->ops->prestep)(ts);                                                                         CHKERRQ(ierr);
1540   ierr = (*ts->ops->step)(ts, steps, ptime);                                                              CHKERRQ(ierr);
1541   ierr = (*ts->ops->poststep)(ts);                                                                        CHKERRQ(ierr);
1542   ierr = PetscLogEventEnd(TS_Step, ts, 0, 0, 0);                                                          CHKERRQ(ierr);
1543 
1544   if (!PetscPreLoadingOn) {
1545     ierr = TSViewFromOptions(ts,ts->name);CHKERRQ(ierr);
1546   }
1547   PetscFunctionReturn(0);
1548 }
1549 
1550 #undef __FUNCT__
1551 #define __FUNCT__ "TSMonitor"
1552 /*
1553      Runs the user provided monitor routines, if they exists.
1554 */
1555 int TSMonitor(TS ts,int step,PetscReal ptime,Vec x)
1556 {
1557   int i,ierr,n = ts->numbermonitors;
1558 
1559   PetscFunctionBegin;
1560   for (i=0; i<n; i++) {
1561     ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr);
1562   }
1563   PetscFunctionReturn(0);
1564 }
1565 
1566 /* ------------------------------------------------------------------------*/
1567 
1568 #undef __FUNCT__
1569 #define __FUNCT__ "TSLGMonitorCreate"
1570 /*@C
1571    TSLGMonitorCreate - Creates a line graph context for use with
1572    TS to monitor convergence of preconditioned residual norms.
1573 
1574    Collective on TS
1575 
1576    Input Parameters:
1577 +  host - the X display to open, or null for the local machine
1578 .  label - the title to put in the title bar
1579 .  x, y - the screen coordinates of the upper left coordinate of the window
1580 -  m, n - the screen width and height in pixels
1581 
1582    Output Parameter:
1583 .  draw - the drawing context
1584 
1585    Options Database Key:
1586 .  -ts_xmonitor - automatically sets line graph monitor
1587 
1588    Notes:
1589    Use TSLGMonitorDestroy() to destroy this line graph, not PetscDrawLGDestroy().
1590 
1591    Level: intermediate
1592 
1593 .keywords: TS, monitor, line graph, residual, seealso
1594 
1595 .seealso: TSLGMonitorDestroy(), TSSetMonitor()
1596 
1597 @*/
1598 int TSLGMonitorCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1599 {
1600   PetscDraw win;
1601   int       ierr;
1602 
1603   PetscFunctionBegin;
1604   ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr);
1605   ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr);
1606   ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr);
1607   ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr);
1608 
1609   PetscLogObjectParent(*draw,win);
1610   PetscFunctionReturn(0);
1611 }
1612 
1613 #undef __FUNCT__
1614 #define __FUNCT__ "TSLGMonitor"
1615 int TSLGMonitor(TS ts,int n,PetscReal ptime,Vec v,void *monctx)
1616 {
1617   PetscDrawLG lg = (PetscDrawLG) monctx;
1618   PetscReal      x,y = ptime;
1619   int         ierr;
1620 
1621   PetscFunctionBegin;
1622   if (!monctx) {
1623     MPI_Comm    comm;
1624     PetscViewer viewer;
1625 
1626     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
1627     viewer = PETSC_VIEWER_DRAW_(comm);
1628     ierr   = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr);
1629   }
1630 
1631   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
1632   x = (PetscReal)n;
1633   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1634   if (n < 20 || (n % 5)) {
1635     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1636   }
1637   PetscFunctionReturn(0);
1638 }
1639 
1640 #undef __FUNCT__
1641 #define __FUNCT__ "TSLGMonitorDestroy"
1642 /*@C
1643    TSLGMonitorDestroy - Destroys a line graph context that was created
1644    with TSLGMonitorCreate().
1645 
1646    Collective on PetscDrawLG
1647 
1648    Input Parameter:
1649 .  draw - the drawing context
1650 
1651    Level: intermediate
1652 
1653 .keywords: TS, monitor, line graph, destroy
1654 
1655 .seealso: TSLGMonitorCreate(),  TSSetMonitor(), TSLGMonitor();
1656 @*/
1657 int TSLGMonitorDestroy(PetscDrawLG drawlg)
1658 {
1659   PetscDraw draw;
1660   int       ierr;
1661 
1662   PetscFunctionBegin;
1663   ierr = PetscDrawLGGetDraw(drawlg,&draw);CHKERRQ(ierr);
1664   ierr = PetscDrawDestroy(draw);CHKERRQ(ierr);
1665   ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr);
1666   PetscFunctionReturn(0);
1667 }
1668 
1669 #undef __FUNCT__
1670 #define __FUNCT__ "TSGetTime"
1671 /*@
1672    TSGetTime - Gets the current time.
1673 
1674    Not Collective
1675 
1676    Input Parameter:
1677 .  ts - the TS context obtained from TSCreate()
1678 
1679    Output Parameter:
1680 .  t  - the current time
1681 
1682    Contributed by: Matthew Knepley
1683 
1684    Level: beginner
1685 
1686 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1687 
1688 .keywords: TS, get, time
1689 @*/
1690 int TSGetTime(TS ts,PetscReal* t)
1691 {
1692   PetscFunctionBegin;
1693   PetscValidHeaderSpecific(ts,TS_COOKIE);
1694   *t = ts->ptime;
1695   PetscFunctionReturn(0);
1696 }
1697 
1698 #undef __FUNCT__
1699 #define __FUNCT__ "TSSetOptionsPrefix"
1700 /*@C
1701    TSSetOptionsPrefix - Sets the prefix used for searching for all
1702    TS options in the database.
1703 
1704    Collective on TS
1705 
1706    Input Parameter:
1707 +  ts     - The TS context
1708 -  prefix - The prefix to prepend to all option names
1709 
1710    Notes:
1711    A hyphen (-) must NOT be given at the beginning of the prefix name.
1712    The first character of all runtime options is AUTOMATICALLY the
1713    hyphen.
1714 
1715    Contributed by: Matthew Knepley
1716 
1717    Level: advanced
1718 
1719 .keywords: TS, set, options, prefix, database
1720 
1721 .seealso: TSSetFromOptions()
1722 
1723 @*/
1724 int TSSetOptionsPrefix(TS ts,const char prefix[])
1725 {
1726   int ierr;
1727 
1728   PetscFunctionBegin;
1729   PetscValidHeaderSpecific(ts,TS_COOKIE);
1730   ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1731   switch(ts->problem_type) {
1732     case TS_NONLINEAR:
1733       ierr = SNESSetOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr);
1734       break;
1735     case TS_LINEAR:
1736       ierr = KSPSetOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr);
1737       break;
1738   }
1739   PetscFunctionReturn(0);
1740 }
1741 
1742 
1743 #undef __FUNCT__
1744 #define __FUNCT__ "TSAppendOptionsPrefix"
1745 /*@C
1746    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
1747    TS options in the database.
1748 
1749    Collective on TS
1750 
1751    Input Parameter:
1752 +  ts     - The TS context
1753 -  prefix - The prefix to prepend to all option names
1754 
1755    Notes:
1756    A hyphen (-) must NOT be given at the beginning of the prefix name.
1757    The first character of all runtime options is AUTOMATICALLY the
1758    hyphen.
1759 
1760    Contributed by: Matthew Knepley
1761 
1762    Level: advanced
1763 
1764 .keywords: TS, append, options, prefix, database
1765 
1766 .seealso: TSGetOptionsPrefix()
1767 
1768 @*/
1769 int TSAppendOptionsPrefix(TS ts,const char prefix[])
1770 {
1771   int ierr;
1772 
1773   PetscFunctionBegin;
1774   PetscValidHeaderSpecific(ts,TS_COOKIE);
1775   ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1776   switch(ts->problem_type) {
1777     case TS_NONLINEAR:
1778       ierr = SNESAppendOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr);
1779       break;
1780     case TS_LINEAR:
1781       ierr = KSPAppendOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr);
1782       break;
1783   }
1784   PetscFunctionReturn(0);
1785 }
1786 
1787 #undef __FUNCT__
1788 #define __FUNCT__ "TSGetOptionsPrefix"
1789 /*@C
1790    TSGetOptionsPrefix - Sets the prefix used for searching for all
1791    TS options in the database.
1792 
1793    Not Collective
1794 
1795    Input Parameter:
1796 .  ts - The TS context
1797 
1798    Output Parameter:
1799 .  prefix - A pointer to the prefix string used
1800 
1801    Contributed by: Matthew Knepley
1802 
1803    Notes: On the fortran side, the user should pass in a string 'prifix' of
1804    sufficient length to hold the prefix.
1805 
1806    Level: intermediate
1807 
1808 .keywords: TS, get, options, prefix, database
1809 
1810 .seealso: TSAppendOptionsPrefix()
1811 @*/
1812 int TSGetOptionsPrefix(TS ts,char *prefix[])
1813 {
1814   int ierr;
1815 
1816   PetscFunctionBegin;
1817   PetscValidHeaderSpecific(ts,TS_COOKIE);
1818   ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1819   PetscFunctionReturn(0);
1820 }
1821 
1822 #undef __FUNCT__
1823 #define __FUNCT__ "TSGetRHSMatrix"
1824 /*@C
1825    TSGetRHSMatrix - Returns the matrix A at the present timestep.
1826 
1827    Not Collective, but parallel objects are returned if TS is parallel
1828 
1829    Input Parameter:
1830 .  ts  - The TS context obtained from TSCreate()
1831 
1832    Output Parameters:
1833 +  A   - The matrix A, where U_t = A(t) U
1834 .  M   - The preconditioner matrix, usually the same as A
1835 -  ctx - User-defined context for matrix evaluation routine
1836 
1837    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1838 
1839    Contributed by: Matthew Knepley
1840 
1841    Level: intermediate
1842 
1843 .seealso: TSGetTimeStep(), TSGetTime(), TSGetTimeStepNumber(), TSGetRHSJacobian()
1844 
1845 .keywords: TS, timestep, get, matrix
1846 
1847 @*/
1848 int TSGetRHSMatrix(TS ts,Mat *A,Mat *M,void **ctx)
1849 {
1850   PetscFunctionBegin;
1851   PetscValidHeaderSpecific(ts,TS_COOKIE);
1852   if (A)   *A = ts->A;
1853   if (M)   *M = ts->B;
1854   if (ctx) *ctx = ts->jacP;
1855   PetscFunctionReturn(0);
1856 }
1857 
1858 #undef __FUNCT__
1859 #define __FUNCT__ "TSGetRHSJacobian"
1860 /*@C
1861    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
1862 
1863    Not Collective, but parallel objects are returned if TS is parallel
1864 
1865    Input Parameter:
1866 .  ts  - The TS context obtained from TSCreate()
1867 
1868    Output Parameters:
1869 +  J   - The Jacobian J of F, where U_t = F(U,t)
1870 .  M   - The preconditioner matrix, usually the same as J
1871 - ctx - User-defined context for Jacobian evaluation routine
1872 
1873    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1874 
1875    Contributed by: Matthew Knepley
1876 
1877    Level: intermediate
1878 
1879 .seealso: TSGetTimeStep(), TSGetRHSMatrix(), TSGetTime(), TSGetTimeStepNumber()
1880 
1881 .keywords: TS, timestep, get, matrix, Jacobian
1882 @*/
1883 int TSGetRHSJacobian(TS ts,Mat *J,Mat *M,void **ctx)
1884 {
1885   int ierr;
1886 
1887   PetscFunctionBegin;
1888   ierr = TSGetRHSMatrix(ts,J,M,ctx);CHKERRQ(ierr);
1889   PetscFunctionReturn(0);
1890 }
1891 
1892 #undef __FUNCT__
1893 #define __FUNCT__ "TSVecViewMonitor"
1894 /*@C
1895    TSVecViewMonitor - Monitors progress of the TS solvers by calling
1896    VecView() for the solution at each timestep
1897 
1898    Collective on TS
1899 
1900    Input Parameters:
1901 +  ts - the TS context
1902 .  step - current time-step
1903 .  ptime - current time
1904 -  dummy - either a viewer or PETSC_NULL
1905 
1906    Level: intermediate
1907 
1908 .keywords: TS,  vector, monitor, view
1909 
1910 .seealso: TSSetMonitor(), TSDefaultMonitor(), VecView()
1911 @*/
1912 int TSVecViewMonitor(TS ts,int step,PetscReal ptime,Vec x,void *dummy)
1913 {
1914   int         ierr;
1915   PetscViewer viewer = (PetscViewer) dummy;
1916 
1917   PetscFunctionBegin;
1918   if (!viewer) {
1919     MPI_Comm comm;
1920     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
1921     viewer = PETSC_VIEWER_DRAW_(comm);
1922   }
1923   ierr = VecView(x,viewer);CHKERRQ(ierr);
1924   PetscFunctionReturn(0);
1925 }
1926 
1927 
1928 
1929