Actual source code: ts.c
1: /* $Id: ts.c,v 1.40 2001/04/05 21:09:19 balay Exp bsmith $ */
2: #include src/ts/tsimpl.h
4: /*@
5: TSComputeRHSJacobian - Computes the Jacobian matrix that has been
6: set with TSSetRHSJacobian().
8: Collective on TS and Vec
10: Input Parameters:
11: + ts - the SNES context
12: . t - current timestep
13: - x - input vector
15: Output Parameters:
16: + A - Jacobian matrix
17: . B - optional preconditioning matrix
18: - flag - flag indicating matrix structure
20: Notes:
21: Most users should not need to explicitly call this routine, as it
22: is used internally within the nonlinear solvers.
24: See SLESSetOperators() for important information about setting the
25: flag parameter.
27: TSComputeJacobian() is valid only for TS_NONLINEAR
29: Level: developer
31: .keywords: SNES, compute, Jacobian, matrix
33: .seealso: TSSetRHSJacobian(), SLESSetOperators()
34: @*/
35: int TSComputeRHSJacobian(TS ts,double t,Vec X,Mat *A,Mat *B,MatStructure *flg)
36: {
43: if (ts->problem_type != TS_NONLINEAR) {
44: SETERRQ(PETSC_ERR_ARG_WRONG,"For TS_NONLINEAR only");
45: }
46: if (!ts->rhsjacobian) return(0);
47: PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);
48: *flg = DIFFERENT_NONZERO_PATTERN;
49: PetscStackPush("TS user Jacobian function");
50: (*ts->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);
51: PetscStackPop;
52: PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);
53: /* make sure user returned a correct Jacobian and preconditioner */
56: return(0);
57: }
59: /*
60: TSComputeRHSFunction - Evaluates the right-hand-side function.
62: Note: If the user did not provide a function but merely a matrix,
63: this routine applies the matrix.
64: */
65: int TSComputeRHSFunction(TS ts,double t,Vec x,Vec y)
66: {
74: PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);
75: if (ts->rhsfunction) {
76: PetscStackPush("TS user right-hand-side function");
77: (*ts->rhsfunction)(ts,t,x,y,ts->funP);
78: PetscStackPop;
79: } else {
80: if (ts->rhsmatrix) { /* assemble matrix for this timestep */
81: MatStructure flg;
82: PetscStackPush("TS user right-hand-side matrix function");
83: (*ts->rhsmatrix)(ts,t,&ts->A,&ts->B,&flg,ts->jacP);
84: PetscStackPop;
85: }
86: MatMult(ts->A,x,y);
87: }
89: /* apply user-provided boundary conditions (only needed if these are time dependent) */
90: TSComputeRHSBoundaryConditions(ts,t,y);
91: PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);
93: return(0);
94: }
96: /*@C
97: TSSetRHSFunction - Sets the routine for evaluating the function,
98: F(t,u), where U_t = F(t,u).
100: Collective on TS
102: Input Parameters:
103: + ts - the TS context obtained from TSCreate()
104: . f - routine for evaluating the right-hand-side function
105: - ctx - [optional] user-defined context for private data for the
106: function evaluation routine (may be PETSC_NULL)
108: Calling sequence of func:
109: $ func (TS ts,double t,Vec u,Vec F,void *ctx);
111: + t - current timestep
112: . u - input vector
113: . F - function vector
114: - ctx - [optional] user-defined function context
116: Important:
117: The user MUST call either this routine or TSSetRHSMatrix().
119: Level: beginner
121: .keywords: TS, timestep, set, right-hand-side, function
123: .seealso: TSSetRHSMatrix()
124: @*/
125: int TSSetRHSFunction(TS ts,int (*f)(TS,double,Vec,Vec,void*),void *ctx)
126: {
130: if (ts->problem_type == TS_LINEAR) {
131: SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot set function for linear problem");
132: }
133: ts->rhsfunction = f;
134: ts->funP = ctx;
135: return(0);
136: }
138: /*@C
139: TSSetRHSMatrix - Sets the function to compute the matrix A, where U_t = A(t) U.
140: Also sets the location to store A.
142: Collective on TS
144: Input Parameters:
145: + ts - the TS context obtained from TSCreate()
146: . A - matrix
147: . B - preconditioner matrix (usually same as A)
148: . f - the matrix evaluation routine; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran)
149: if A is not a function of t.
150: - ctx - [optional] user-defined context for private data for the
151: matrix evaluation routine (may be PETSC_NULL)
153: Calling sequence of func:
154: $ func (TS ts,double t,Mat *A,Mat *B,int *flag,void *ctx);
156: + t - current timestep
157: . A - matrix A, where U_t = A(t) U
158: . B - preconditioner matrix, usually the same as A
159: . flag - flag indicating information about the preconditioner matrix
160: structure (same as flag in SLESSetOperators())
161: - ctx - [optional] user-defined context for matrix evaluation routine
163: Notes:
164: See SLESSetOperators() for important information about setting the flag
165: output parameter in the routine func(). Be sure to read this information!
167: The routine func() takes Mat * as the matrix arguments rather than Mat.
168: This allows the matrix evaluation routine to replace A and/or B with a
169: completely new new matrix structure (not just different matrix elements)
170: when appropriate, for instance, if the nonzero structure is changing
171: throughout the global iterations.
173: Important:
174: The user MUST call either this routine or TSSetRHSFunction().
176: Level: beginner
178: .keywords: TS, timestep, set, right-hand-side, matrix
180: .seealso: TSSetRHSFunction()
181: @*/
182: int TSSetRHSMatrix(TS ts,Mat A,Mat B,int (*f)(TS,double,Mat*,Mat*,MatStructure*,void*),void *ctx)
183: {
190: if (ts->problem_type == TS_NONLINEAR) {
191: SETERRQ(PETSC_ERR_ARG_WRONG,"Not for nonlinear problems; use TSSetRHSJacobian()");
192: }
194: ts->rhsmatrix = f;
195: ts->jacP = ctx;
196: ts->A = A;
197: ts->B = B;
199: return(0);
200: }
202: /*@C
203: TSSetRHSJacobian - Sets the function to compute the Jacobian of F,
204: where U_t = F(U,t), as well as the location to store the matrix.
206: Collective on TS
208: Input Parameters:
209: + ts - the TS context obtained from TSCreate()
210: . A - Jacobian matrix
211: . B - preconditioner matrix (usually same as A)
212: . f - the Jacobian evaluation routine
213: - ctx - [optional] user-defined context for private data for the
214: Jacobian evaluation routine (may be PETSC_NULL)
216: Calling sequence of func:
217: $ func (TS ts,double t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx);
219: + t - current timestep
220: . u - input vector
221: . A - matrix A, where U_t = A(t)u
222: . B - preconditioner matrix, usually the same as A
223: . flag - flag indicating information about the preconditioner matrix
224: structure (same as flag in SLESSetOperators())
225: - ctx - [optional] user-defined context for matrix evaluation routine
227: Notes:
228: See SLESSetOperators() for important information about setting the flag
229: output parameter in the routine func(). Be sure to read this information!
231: The routine func() takes Mat * as the matrix arguments rather than Mat.
232: This allows the matrix evaluation routine to replace A and/or B with a
233: completely new new matrix structure (not just different matrix elements)
234: when appropriate, for instance, if the nonzero structure is changing
235: throughout the global iterations.
237: Level: beginner
238:
239: .keywords: TS, timestep, set, right-hand-side, Jacobian
241: .seealso: TSDefaultComputeJacobianColor(),
242: SNESDefaultComputeJacobianColor()
244: @*/
245: int TSSetRHSJacobian(TS ts,Mat A,Mat B,int (*f)(TS,double,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
246: {
253: if (ts->problem_type != TS_NONLINEAR) {
254: SETERRQ(PETSC_ERR_ARG_WRONG,"Not for linear problems; use TSSetRHSMatrix()");
255: }
257: ts->rhsjacobian = f;
258: ts->jacP = ctx;
259: ts->A = A;
260: ts->B = B;
261: return(0);
262: }
264: /*
265: TSComputeRHSBoundaryConditions - Evaluates the boundary condition function.
267: Note: If the user did not provide a function but merely a matrix,
268: this routine applies the matrix.
269: */
270: int TSComputeRHSBoundaryConditions(TS ts,double t,Vec x)
271: {
279: if (ts->rhsbc) {
280: PetscStackPush("TS user boundary condition function");
281: (*ts->rhsbc)(ts,t,x,ts->bcP);
282: PetscStackPop;
283: return(0);
284: }
286: return(0);
287: }
289: /*@C
290: TSSetRHSBoundaryConditions - Sets the routine for evaluating the function,
291: boundary conditions for the function F.
293: Collective on TS
295: Input Parameters:
296: + ts - the TS context obtained from TSCreate()
297: . f - routine for evaluating the boundary condition function
298: - ctx - [optional] user-defined context for private data for the
299: function evaluation routine (may be PETSC_NULL)
301: Calling sequence of func:
302: $ func (TS ts,double t,Vec F,void *ctx);
304: + t - current timestep
305: . F - function vector
306: - ctx - [optional] user-defined function context
308: Level: intermediate
310: .keywords: TS, timestep, set, boundary conditions, function
311: @*/
312: int TSSetRHSBoundaryConditions(TS ts,int (*f)(TS,double,Vec,void*),void *ctx)
313: {
317: if (ts->problem_type != TS_LINEAR) {
318: SETERRQ(PETSC_ERR_ARG_WRONG,"For linear problems only");
319: }
320: ts->rhsbc = f;
321: ts->bcP = ctx;
322: return(0);
323: }
325: /*@C
326: TSView - Prints the TS data structure.
328: Collective on TS
330: Input Parameters:
331: + ts - the TS context obtained from TSCreate()
332: - viewer - visualization context
334: Options Database Key:
335: . -ts_view - calls TSView() at end of TSStep()
337: Notes:
338: The available visualization contexts include
339: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
340: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
341: output where only the first processor opens
342: the file. All other processors send their
343: data to the first processor to print.
345: The user can open an alternative visualization context with
346: PetscViewerASCIIOpen() - output to a specified file.
348: Level: beginner
350: .keywords: TS, timestep, view
352: .seealso: PetscViewerASCIIOpen()
353: @*/
354: int TSView(TS ts,PetscViewer viewer)
355: {
356: int ierr;
357: char *type;
358: PetscTruth isascii,isstring;
362: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(ts->comm);
366: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
367: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
368: if (isascii) {
369: PetscViewerASCIIPrintf(viewer,"TS Object:n");
370: TSGetType(ts,(TSType *)&type);
371: if (type) {
372: PetscViewerASCIIPrintf(viewer," type: %sn",type);
373: } else {
374: PetscViewerASCIIPrintf(viewer," type: not yet setn");
375: }
376: if (ts->view) {
377: PetscViewerASCIIPushTab(viewer);
378: (*ts->view)(ts,viewer);
379: PetscViewerASCIIPopTab(viewer);
380: }
381: PetscViewerASCIIPrintf(viewer," maximum steps=%dn",ts->max_steps);
382: PetscViewerASCIIPrintf(viewer," maximum time=%gn",ts->max_time);
383: if (ts->problem_type == TS_NONLINEAR) {
384: PetscViewerASCIIPrintf(viewer," total number of nonlinear solver iterations=%dn",ts->nonlinear_its);
385: }
386: PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%dn",ts->linear_its);
387: } else if (isstring) {
388: TSGetType(ts,(TSType *)&type);
389: PetscViewerStringSPrintf(viewer," %-7.7s",type);
390: }
391: PetscViewerASCIIPushTab(viewer);
392: if (ts->sles) {SLESView(ts->sles,viewer);}
393: if (ts->snes) {SNESView(ts->snes,viewer);}
394: PetscViewerASCIIPopTab(viewer);
395: return(0);
396: }
399: /*@C
400: TSSetApplicationContext - Sets an optional user-defined context for
401: the timesteppers.
403: Collective on TS
405: Input Parameters:
406: + ts - the TS context obtained from TSCreate()
407: - usrP - optional user context
409: Level: intermediate
411: .keywords: TS, timestep, set, application, context
413: .seealso: TSGetApplicationContext()
414: @*/
415: int TSSetApplicationContext(TS ts,void *usrP)
416: {
419: ts->user = usrP;
420: return(0);
421: }
423: /*@C
424: TSGetApplicationContext - Gets the user-defined context for the
425: timestepper.
427: Not Collective
429: Input Parameter:
430: . ts - the TS context obtained from TSCreate()
432: Output Parameter:
433: . usrP - user context
435: Level: intermediate
437: .keywords: TS, timestep, get, application, context
439: .seealso: TSSetApplicationContext()
440: @*/
441: int TSGetApplicationContext(TS ts,void **usrP)
442: {
445: *usrP = ts->user;
446: return(0);
447: }
449: /*@
450: TSGetTimeStepNumber - Gets the current number of timesteps.
452: Not Collective
454: Input Parameter:
455: . ts - the TS context obtained from TSCreate()
457: Output Parameter:
458: . iter - number steps so far
460: Level: intermediate
462: .keywords: TS, timestep, get, iteration, number
463: @*/
464: int TSGetTimeStepNumber(TS ts,int* iter)
465: {
468: *iter = ts->steps;
469: return(0);
470: }
472: /*@
473: TSSetInitialTimeStep - Sets the initial timestep to be used,
474: as well as the initial time.
476: Collective on TS
478: Input Parameters:
479: + ts - the TS context obtained from TSCreate()
480: . initial_time - the initial time
481: - time_step - the size of the timestep
483: Level: intermediate
485: .seealso: TSSetTimeStep(), TSGetTimeStep()
487: .keywords: TS, set, initial, timestep
488: @*/
489: int TSSetInitialTimeStep(TS ts,double initial_time,double time_step)
490: {
493: ts->time_step = time_step;
494: ts->initial_time_step = time_step;
495: ts->ptime = initial_time;
496: return(0);
497: }
499: /*@
500: TSSetTimeStep - Allows one to reset the timestep at any time,
501: useful for simple pseudo-timestepping codes.
503: Collective on TS
505: Input Parameters:
506: + ts - the TS context obtained from TSCreate()
507: - time_step - the size of the timestep
509: Level: intermediate
511: .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
513: .keywords: TS, set, timestep
514: @*/
515: int TSSetTimeStep(TS ts,double time_step)
516: {
519: ts->time_step = time_step;
520: return(0);
521: }
523: /*@
524: TSGetTimeStep - Gets the current timestep size.
526: Not Collective
528: Input Parameter:
529: . ts - the TS context obtained from TSCreate()
531: Output Parameter:
532: . dt - the current timestep size
534: Level: intermediate
536: .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
538: .keywords: TS, get, timestep
539: @*/
540: int TSGetTimeStep(TS ts,double* dt)
541: {
544: *dt = ts->time_step;
545: return(0);
546: }
548: /*@C
549: TSGetSolution - Returns the solution at the present timestep. It
550: is valid to call this routine inside the function that you are evaluating
551: in order to move to the new timestep. This vector not changed until
552: the solution at the next timestep has been calculated.
554: Not Collective, but Vec returned is parallel if TS is parallel
556: Input Parameter:
557: . ts - the TS context obtained from TSCreate()
559: Output Parameter:
560: . v - the vector containing the solution
562: Level: intermediate
564: .seealso: TSGetTimeStep()
566: .keywords: TS, timestep, get, solution
567: @*/
568: int TSGetSolution(TS ts,Vec *v)
569: {
572: *v = ts->vec_sol_always;
573: return(0);
574: }
576: static int TSPublish_Petsc(PetscObject obj)
577: {
578: #if defined(PETSC_HAVE_AMS)
579: TS v = (TS) obj;
580: int ierr;
581: #endif
585: #if defined(PETSC_HAVE_AMS)
586: /* if it is already published then return */
587: if (v->amem >=0) return(0);
589: PetscObjectPublishBaseBegin(obj);
590: AMS_Memory_add_field((AMS_Memory)v->amem,"Step",&v->steps,1,AMS_INT,AMS_READ,
591: AMS_COMMON,AMS_REDUCT_UNDEF);
592: AMS_Memory_add_field((AMS_Memory)v->amem,"Time",&v->ptime,1,AMS_DOUBLE,AMS_READ,
593: AMS_COMMON,AMS_REDUCT_UNDEF);
594: AMS_Memory_add_field((AMS_Memory)v->amem,"CurrentTimeStep",&v->time_step,1,
595: AMS_DOUBLE,AMS_READ,AMS_COMMON,AMS_REDUCT_UNDEF);
596: PetscObjectPublishBaseEnd(obj);
597: #endif
598: return(0);
599: }
601: /* -----------------------------------------------------------*/
603: /*@C
604: TSCreate - Creates a timestepper context.
606: Collective on MPI_Comm
608: Input Parameter:
609: + comm - MPI communicator
610: - type - One of TS_LINEAR,TS_NONLINEAR
611: where these types refer to problems of the forms
612: .vb
613: U_t = A U
614: U_t = A(t) U
615: U_t = F(t,U)
616: .ve
618: Output Parameter:
619: . outts - the new TS context
621: Level: beginner
623: .keywords: TS, timestep, create, context
625: .seealso: TSSetUp(), TSStep(), TSDestroy(), TSProblemType, TS
626: @*/
627: int TSCreate(MPI_Comm comm,TSProblemType problemtype,TS *outts)
628: {
629: TS ts;
632: *outts = 0;
633: PetscHeaderCreate(ts,_p_TS,int,TS_COOKIE,-1,"TS",comm,TSDestroy,TSView);
634: PetscLogObjectCreate(ts);
635: ts->bops->publish = TSPublish_Petsc;
636: ts->max_steps = 5000;
637: ts->max_time = 5.0;
638: ts->time_step = .1;
639: ts->initial_time_step = ts->time_step;
640: ts->steps = 0;
641: ts->ptime = 0.0;
642: ts->data = 0;
643: ts->view = 0;
644: ts->setupcalled = 0;
645: ts->problem_type = problemtype;
646: ts->numbermonitors = 0;
647: ts->linear_its = 0;
648: ts->nonlinear_its = 0;
650: *outts = ts;
651: return(0);
652: }
654: /* ----- Routines to initialize and destroy a timestepper ---- */
656: /*@
657: TSSetUp - Sets up the internal data structures for the later use
658: of a timestepper.
660: Collective on TS
662: Input Parameter:
663: . ts - the TS context obtained from TSCreate()
665: Notes:
666: For basic use of the TS solvers the user need not explicitly call
667: TSSetUp(), since these actions will automatically occur during
668: the call to TSStep(). However, if one wishes to control this
669: phase separately, TSSetUp() should be called after TSCreate()
670: and optional routines of the form TSSetXXX(), but before TSStep().
672: Level: advanced
674: .keywords: TS, timestep, setup
676: .seealso: TSCreate(), TSStep(), TSDestroy()
677: @*/
678: int TSSetUp(TS ts)
679: {
684: if (!ts->vec_sol) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first");
685: if (!ts->type_name) {
686: TSSetType(ts,TS_EULER);
687: }
688: (*ts->setup)(ts);
689: ts->setupcalled = 1;
690: return(0);
691: }
693: /*@C
694: TSDestroy - Destroys the timestepper context that was created
695: with TSCreate().
697: Collective on TS
699: Input Parameter:
700: . ts - the TS context obtained from TSCreate()
702: Level: beginner
704: .keywords: TS, timestepper, destroy
706: .seealso: TSCreate(), TSSetUp(), TSSolve()
707: @*/
708: int TSDestroy(TS ts)
709: {
710: int ierr,i;
714: if (--ts->refct > 0) return(0);
716: /* if memory was published with AMS then destroy it */
717: PetscObjectDepublish(ts);
719: if (ts->sles) {SLESDestroy(ts->sles);}
720: if (ts->snes) {SNESDestroy(ts->snes);}
721: (*(ts)->destroy)(ts);
722: for (i=0; i<ts->numbermonitors; i++) {
723: if (ts->mdestroy[i]) {
724: (*ts->mdestroy[i])(ts->monitorcontext[i]);
725: }
726: }
727: PetscLogObjectDestroy((PetscObject)ts);
728: PetscHeaderDestroy((PetscObject)ts);
729: return(0);
730: }
732: /*@C
733: TSGetSNES - Returns the SNES (nonlinear solver) associated with
734: a TS (timestepper) context. Valid only for nonlinear problems.
736: Not Collective, but SNES is parallel if TS is parallel
738: Input Parameter:
739: . ts - the TS context obtained from TSCreate()
741: Output Parameter:
742: . snes - the nonlinear solver context
744: Notes:
745: The user can then directly manipulate the SNES context to set various
746: options, etc. Likewise, the user can then extract and manipulate the
747: SLES, KSP, and PC contexts as well.
749: TSGetSNES() does not work for integrators that do not use SNES; in
750: this case TSGetSNES() returns PETSC_NULL in snes.
752: Level: beginner
754: .keywords: timestep, get, SNES
755: @*/
756: int TSGetSNES(TS ts,SNES *snes)
757: {
760: if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Nonlinear only; use TSGetSLES()");
761: *snes = ts->snes;
762: return(0);
763: }
765: /*@C
766: TSGetSLES - Returns the SLES (linear solver) associated with
767: a TS (timestepper) context.
769: Not Collective, but SLES is parallel if TS is parallel
771: Input Parameter:
772: . ts - the TS context obtained from TSCreate()
774: Output Parameter:
775: . sles - the nonlinear solver context
777: Notes:
778: The user can then directly manipulate the SLES context to set various
779: options, etc. Likewise, the user can then extract and manipulate the
780: KSP and PC contexts as well.
782: TSGetSLES() does not work for integrators that do not use SLES;
783: in this case TSGetSLES() returns PETSC_NULL in sles.
785: Level: beginner
787: .keywords: timestep, get, SLES
788: @*/
789: int TSGetSLES(TS ts,SLES *sles)
790: {
793: if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()");
794: *sles = ts->sles;
795: return(0);
796: }
798: /* ----------- Routines to set solver parameters ---------- */
800: /*@
801: TSSetDuration - Sets the maximum number of timesteps to use and
802: maximum time for iteration.
804: Collective on TS
806: Input Parameters:
807: + ts - the TS context obtained from TSCreate()
808: . maxsteps - maximum number of iterations to use
809: - maxtime - final time to iterate to
811: Options Database Keys:
812: . -ts_max_steps <maxsteps> - Sets maxsteps
813: . -ts_max_time <maxtime> - Sets maxtime
815: Notes:
816: The default maximum number of iterations is 5000. Default time is 5.0
818: Level: intermediate
820: .keywords: TS, timestep, set, maximum, iterations
821: @*/
822: int TSSetDuration(TS ts,int maxsteps,double maxtime)
823: {
826: ts->max_steps = maxsteps;
827: ts->max_time = maxtime;
828: return(0);
829: }
831: /*@
832: TSSetSolution - Sets the initial solution vector
833: for use by the TS routines.
835: Collective on TS and Vec
837: Input Parameters:
838: + ts - the TS context obtained from TSCreate()
839: - x - the solution vector
841: Level: beginner
843: .keywords: TS, timestep, set, solution, initial conditions
844: @*/
845: int TSSetSolution(TS ts,Vec x)
846: {
849: ts->vec_sol = ts->vec_sol_always = x;
850: return(0);
851: }
853: /* ------------ Routines to set performance monitoring options ----------- */
855: /*@C
856: TSSetMonitor - Sets an ADDITIONAL function that is to be used at every
857: timestep to display the iteration's progress.
859: Collective on TS
861: Input Parameters:
862: + ts - the TS context obtained from TSCreate()
863: . func - monitoring routine
864: . mctx - [optional] user-defined context for private data for the
865: monitor routine (use PETSC_NULL if no context is desired)
866: - monitordestroy - [optional] routine that frees monitor context
867: (may be PETSC_NULL)
869: Calling sequence of func:
870: $ int func(TS ts,int steps,double time,Vec x,void *mctx)
872: + ts - the TS context
873: . steps - iteration number
874: . time - current time
875: . x - current iterate
876: - mctx - [optional] monitoring context
878: Notes:
879: This routine adds an additional monitor to the list of monitors that
880: already has been loaded.
882: Level: intermediate
884: .keywords: TS, timestep, set, monitor
886: .seealso: TSDefaultMonitor(), TSClearMonitor()
887: @*/
888: int TSSetMonitor(TS ts,int (*monitor)(TS,int,double,Vec,void*),void *mctx,int (*mdestroy)(void*))
889: {
892: if (ts->numbermonitors >= MAXTSMONITORS) {
893: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
894: }
895: ts->monitor[ts->numbermonitors] = monitor;
896: ts->mdestroy[ts->numbermonitors] = mdestroy;
897: ts->monitorcontext[ts->numbermonitors++] = (void*)mctx;
898: return(0);
899: }
901: /*@C
902: TSClearMonitor - Clears all the monitors that have been set on a time-step object.
904: Collective on TS
906: Input Parameters:
907: . ts - the TS context obtained from TSCreate()
909: Notes:
910: There is no way to remove a single, specific monitor.
912: Level: intermediate
914: .keywords: TS, timestep, set, monitor
916: .seealso: TSDefaultMonitor(), TSSetMonitor()
917: @*/
918: int TSClearMonitor(TS ts)
919: {
922: ts->numbermonitors = 0;
923: return(0);
924: }
926: int TSDefaultMonitor(TS ts,int step,double ptime,Vec v,void *ctx)
927: {
931: PetscPrintf(ts->comm,"timestep %d dt %g time %gn",step,ts->time_step,ptime);
932: return(0);
933: }
935: /*@
936: TSStep - Steps the requested number of timesteps.
938: Collective on TS
940: Input Parameter:
941: . ts - the TS context obtained from TSCreate()
943: Output Parameters:
944: + steps - number of iterations until termination
945: - ptime - time until termination
947: Level: beginner
949: .keywords: TS, timestep, solve
951: .seealso: TSCreate(), TSSetUp(), TSDestroy()
952: @*/
953: int TSStep(TS ts,int *steps,double *ptime)
954: {
955: int ierr;
956: PetscTruth flg;
960: if (!ts->setupcalled) {TSSetUp(ts);}
961: PetscLogEventBegin(TS_Step,ts,0,0,0);
962: (*ts->step)(ts,steps,ptime);
963: PetscLogEventEnd(TS_Step,ts,0,0,0);
964: PetscOptionsHasName(ts->prefix,"-ts_view",&flg);
965: if (flg && !PetscPreLoadingOn) {TSView(ts,PETSC_VIEWER_STDOUT_WORLD);}
966: return(0);
967: }
969: /*
970: Runs the user provided monitor routines, if they exists.
971: */
972: int TSMonitor(TS ts,int step,double ptime,Vec x)
973: {
974: int i,ierr,n = ts->numbermonitors;
977: for (i=0; i<n; i++) {
978: (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);
979: }
980: return(0);
981: }
983: /* ------------------------------------------------------------------------*/
985: /*@C
986: TSLGMonitorCreate - Creates a line graph context for use with
987: TS to monitor convergence of preconditioned residual norms.
989: Collective on TS
991: Input Parameters:
992: + host - the X display to open, or null for the local machine
993: . label - the title to put in the title bar
994: . x, y - the screen coordinates of the upper left coordinate of the window
995: - m, n - the screen width and height in pixels
997: Output Parameter:
998: . draw - the drawing context
1000: Options Database Key:
1001: . -ts_xmonitor - automatically sets line graph monitor
1003: Notes:
1004: Use TSLGMonitorDestroy() to destroy this line graph, not PetscDrawLGDestroy().
1006: Level: intermediate
1008: .keywords: TS, monitor, line graph, residual, seealso
1010: .seealso: TSLGMonitorDestroy(), TSSetMonitor()
1012: @*/
1013: int TSLGMonitorCreate(char *host,char *label,int x,int y,int m,int n,PetscDrawLG *draw)
1014: {
1015: PetscDraw win;
1016: int ierr;
1019: PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);
1020: PetscDrawSetType(win,PETSC_DRAW_X);
1021: PetscDrawLGCreate(win,1,draw);
1022: PetscDrawLGIndicateDataPoints(*draw);
1024: PetscLogObjectParent(*draw,win);
1025: return(0);
1026: }
1028: int TSLGMonitor(TS ts,int n,double ptime,Vec v,void *monctx)
1029: {
1030: PetscDrawLG lg = (PetscDrawLG) monctx;
1031: double x,y = ptime;
1032: int ierr;
1035: if (!monctx) {
1036: MPI_Comm comm;
1037: PetscViewer viewer;
1039: ierr = PetscObjectGetComm((PetscObject)ts,&comm);
1040: viewer = PETSC_VIEWER_DRAW_(comm);
1041: ierr = PetscViewerDrawGetDrawLG(viewer,0,&lg);
1042: }
1044: if (!n) {PetscDrawLGReset(lg);}
1045: x = (double)n;
1046: PetscDrawLGAddPoint(lg,&x,&y);
1047: if (n < 20 || (n % 5)) {
1048: PetscDrawLGDraw(lg);
1049: }
1050: return(0);
1051: }
1053: /*@C
1054: TSLGMonitorDestroy - Destroys a line graph context that was created
1055: with TSLGMonitorCreate().
1057: Collective on PetscDrawLG
1059: Input Parameter:
1060: . draw - the drawing context
1062: Level: intermediate
1064: .keywords: TS, monitor, line graph, destroy
1066: .seealso: TSLGMonitorCreate(), TSSetMonitor(), TSLGMonitor();
1067: @*/
1068: int TSLGMonitorDestroy(PetscDrawLG drawlg)
1069: {
1070: PetscDraw draw;
1071: int ierr;
1074: PetscDrawLGGetDraw(drawlg,&draw);
1075: PetscDrawDestroy(draw);
1076: PetscDrawLGDestroy(drawlg);
1077: return(0);
1078: }
1080: /*@
1081: TSGetTime - Gets the current time.
1083: Not Collective
1085: Input Parameter:
1086: . ts - the TS context obtained from TSCreate()
1088: Output Parameter:
1089: . t - the current time
1091: Contributed by: Matthew Knepley
1093: Level: beginner
1095: .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1097: .keywords: TS, get, time
1098: @*/
1099: int TSGetTime(TS ts,double* t)
1100: {
1103: *t = ts->ptime;
1104: return(0);
1105: }
1107: /*@C
1108: TSGetProblemType - Returns the problem type of a TS (timestepper) context.
1110: Not Collective
1112: Input Parameter:
1113: . ts - The TS context obtained from TSCreate()
1115: Output Parameter:
1116: . type - The problem type, TS_LINEAR or TS_NONLINEAR
1118: Level: intermediate
1120: Contributed by: Matthew Knepley
1122: .keywords: ts, get, type
1124: @*/
1125: int TSGetProblemType(TS ts,TSProblemType *type)
1126: {
1129: *type = ts->problem_type;
1130: return(0);
1131: }
1133: /*@C
1134: TSSetOptionsPrefix - Sets the prefix used for searching for all
1135: TS options in the database.
1137: Collective on TS
1139: Input Parameter:
1140: + ts - The TS context
1141: - prefix - The prefix to prepend to all option names
1143: Notes:
1144: A hyphen (-) must NOT be given at the beginning of the prefix name.
1145: The first character of all runtime options is AUTOMATICALLY the
1146: hyphen.
1148: Contributed by: Matthew Knepley
1150: Level: advanced
1152: .keywords: TS, set, options, prefix, database
1154: .seealso: TSSetFromOptions()
1156: @*/
1157: int TSSetOptionsPrefix(TS ts,char *prefix)
1158: {
1163: PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);
1164: switch(ts->problem_type) {
1165: case TS_NONLINEAR:
1166: SNESSetOptionsPrefix(ts->snes,prefix);
1167: break;
1168: case TS_LINEAR:
1169: SLESSetOptionsPrefix(ts->sles,prefix);
1170: break;
1171: }
1172: return(0);
1173: }
1176: /*@C
1177: TSAppendOptionsPrefix - Appends to the prefix used for searching for all
1178: TS options in the database.
1180: Collective on TS
1182: Input Parameter:
1183: + ts - The TS context
1184: - prefix - The prefix to prepend to all option names
1186: Notes:
1187: A hyphen (-) must NOT be given at the beginning of the prefix name.
1188: The first character of all runtime options is AUTOMATICALLY the
1189: hyphen.
1191: Contributed by: Matthew Knepley
1193: Level: advanced
1195: .keywords: TS, append, options, prefix, database
1197: .seealso: TSGetOptionsPrefix()
1199: @*/
1200: int TSAppendOptionsPrefix(TS ts,char *prefix)
1201: {
1206: PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);
1207: switch(ts->problem_type) {
1208: case TS_NONLINEAR:
1209: SNESAppendOptionsPrefix(ts->snes,prefix);
1210: break;
1211: case TS_LINEAR:
1212: SLESAppendOptionsPrefix(ts->sles,prefix);
1213: break;
1214: }
1215: return(0);
1216: }
1218: /*@C
1219: TSGetOptionsPrefix - Sets the prefix used for searching for all
1220: TS options in the database.
1222: Not Collective
1224: Input Parameter:
1225: . ts - The TS context
1227: Output Parameter:
1228: . prefix - A pointer to the prefix string used
1230: Contributed by: Matthew Knepley
1232: Notes: On the fortran side, the user should pass in a string 'prifix' of
1233: sufficient length to hold the prefix.
1235: Level: intermediate
1237: .keywords: TS, get, options, prefix, database
1239: .seealso: TSAppendOptionsPrefix()
1240: @*/
1241: int TSGetOptionsPrefix(TS ts,char **prefix)
1242: {
1247: PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);
1248: return(0);
1249: }
1251: /*@C
1252: TSGetRHSMatrix - Returns the matrix A at the present timestep.
1254: Not Collective, but parallel objects are returned if TS is parallel
1256: Input Parameter:
1257: . ts - The TS context obtained from TSCreate()
1259: Output Parameters:
1260: + A - The matrix A, where U_t = A(t) U
1261: . M - The preconditioner matrix, usually the same as A
1262: - ctx - User-defined context for matrix evaluation routine
1264: Notes: You can pass in PETSC_NULL for any return argument you do not need.
1266: Contributed by: Matthew Knepley
1268: Level: intermediate
1270: .seealso: TSGetTimeStep(), TSGetTime(), TSGetTimeStepNumber(), TSGetRHSJacobian()
1272: .keywords: TS, timestep, get, matrix
1274: @*/
1275: int TSGetRHSMatrix(TS ts,Mat *A,Mat *M,void **ctx)
1276: {
1279: if (A) *A = ts->A;
1280: if (M) *M = ts->B;
1281: if (ctx) *ctx = ts->jacP;
1282: return(0);
1283: }
1285: /*@C
1286: TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
1288: Not Collective, but parallel objects are returned if TS is parallel
1290: Input Parameter:
1291: . ts - The TS context obtained from TSCreate()
1293: Output Parameters:
1294: + J - The Jacobian J of F, where U_t = F(U,t)
1295: . M - The preconditioner matrix, usually the same as J
1296: - ctx - User-defined context for Jacobian evaluation routine
1298: Notes: You can pass in PETSC_NULL for any return argument you do not need.
1300: Contributed by: Matthew Knepley
1302: Level: intermediate
1304: .seealso: TSGetTimeStep(), TSGetRHSMatrix(), TSGetTime(), TSGetTimeStepNumber()
1306: .keywords: TS, timestep, get, matrix, Jacobian
1307: @*/
1308: int TSGetRHSJacobian(TS ts,Mat *J,Mat *M,void **ctx)
1309: {
1313: TSGetRHSMatrix(ts,J,M,ctx);
1314: return(0);
1315: }
1317: /*MC
1318: TSRegisterDynamic - Adds a method to the timestepping solver package.
1320: Synopsis:
1322: int TSRegisterDynamic(char *name_solver,char *path,char *name_create,int (*routine_create)(TS))
1324: Not collective
1326: Input Parameters:
1327: + name_solver - name of a new user-defined solver
1328: . path - path (either absolute or relative) the library containing this solver
1329: . name_create - name of routine to create method context
1330: - routine_create - routine to create method context
1332: Notes:
1333: TSRegisterDynamic() may be called multiple times to add several user-defined solvers.
1335: If dynamic libraries are used, then the fourth input argument (routine_create)
1336: is ignored.
1338: Sample usage:
1339: .vb
1340: TSRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
1341: "MySolverCreate",MySolverCreate);
1342: .ve
1344: Then, your solver can be chosen with the procedural interface via
1345: $ TSSetType(ts,"my_solver")
1346: or at runtime via the option
1347: $ -ts_type my_solver
1349: Level: advanced
1351: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, ${BOPT},
1352: and others of the form ${any_environmental_variable} occuring in pathname will be
1353: replaced with appropriate values.
1355: .keywords: TS, register
1357: .seealso: TSRegisterAll(), TSRegisterDestroy()
1358: M*/
1360: int TSRegister(char *sname,char *path,char *name,int (*function)(TS))
1361: {
1362: char fullname[256];
1363: int ierr;
1366: PetscFListConcat(path,name,fullname);
1367: PetscFListAdd(&TSList,sname,fullname,(void (*)())function);
1368: return(0);
1369: }
1371: /*@C
1372: TSVecViewMonitor - Monitors progress of the TS solvers by calling
1373: VecView() for the solution at each timestep
1375: Collective on TS
1377: Input Parameters:
1378: + ts - the TS context
1379: . step - current time-step
1380: . ptime - current time
1381: - dummy - either a viewer or PETSC_NULL
1383: Level: intermediate
1385: .keywords: TS, vector, monitor, view
1387: .seealso: TSSetMonitor(), TSDefaultMonitor(), VecView()
1388: @*/
1389: int TSVecViewMonitor(TS ts,int step,PetscReal ptime,Vec x,void *dummy)
1390: {
1391: int ierr;
1392: PetscViewer viewer = (PetscViewer) dummy;
1395: if (!viewer) {
1396: MPI_Comm comm;
1397: ierr = PetscObjectGetComm((PetscObject)ts,&comm);
1398: viewer = PETSC_VIEWER_DRAW_(comm);
1399: }
1400: VecView(x,viewer);
1401: return(0);
1402: }