Actual source code: ts.c
petsc-3.3-p0 2012-06-05
2: #include <petsc-private/tsimpl.h> /*I "petscts.h" I*/
3: #include <petscdmshell.h>
5: /* Logging support */
6: PetscClassId TS_CLASSID;
7: PetscLogEvent TS_Step, TS_PseudoComputeTimeStep, TS_FunctionEval, TS_JacobianEval;
11: /*
12: TSSetTypeFromOptions - Sets the type of ts from user options.
14: Collective on TS
16: Input Parameter:
17: . ts - The ts
19: Level: intermediate
21: .keywords: TS, set, options, database, type
22: .seealso: TSSetFromOptions(), TSSetType()
23: */
24: static PetscErrorCode TSSetTypeFromOptions(TS ts)
25: {
26: PetscBool opt;
27: const char *defaultType;
28: char typeName[256];
32: if (((PetscObject)ts)->type_name) {
33: defaultType = ((PetscObject)ts)->type_name;
34: } else {
35: defaultType = TSEULER;
36: }
38: if (!TSRegisterAllCalled) {TSRegisterAll(PETSC_NULL);}
39: PetscOptionsList("-ts_type", "TS method"," TSSetType", TSList, defaultType, typeName, 256, &opt);
40: if (opt) {
41: TSSetType(ts, typeName);
42: } else {
43: TSSetType(ts, defaultType);
44: }
45: return(0);
46: }
50: /*@
51: TSSetFromOptions - Sets various TS parameters from user options.
53: Collective on TS
55: Input Parameter:
56: . ts - the TS context obtained from TSCreate()
58: Options Database Keys:
59: + -ts_type <type> - TSEULER, TSBEULER, TSSUNDIALS, TSPSEUDO, TSCN, TSRK, TSTHETA, TSGL, TSSSP
60: . -ts_max_steps maxsteps - maximum number of time-steps to take
61: . -ts_final_time time - maximum time to compute to
62: . -ts_dt dt - initial time step
63: . -ts_monitor - print information at each timestep
64: - -ts_monitor_draw - plot information at each timestep
66: Level: beginner
68: .keywords: TS, timestep, set, options, database
70: .seealso: TSGetType()
71: @*/
72: PetscErrorCode TSSetFromOptions(TS ts)
73: {
74: PetscBool opt,flg;
76: PetscViewer monviewer;
77: char monfilename[PETSC_MAX_PATH_LEN];
78: SNES snes;
79: TSAdapt adapt;
83: PetscObjectOptionsBegin((PetscObject)ts);
84: /* Handle TS type options */
85: TSSetTypeFromOptions(ts);
87: /* Handle generic TS options */
88: PetscOptionsInt("-ts_max_steps","Maximum number of time steps","TSSetDuration",ts->max_steps,&ts->max_steps,PETSC_NULL);
89: PetscOptionsReal("-ts_final_time","Time to run to","TSSetDuration",ts->max_time,&ts->max_time,PETSC_NULL);
90: PetscOptionsReal("-ts_init_time","Initial time","TSSetTime",ts->ptime,&ts->ptime,PETSC_NULL);
91: PetscOptionsReal("-ts_dt","Initial time step","TSSetTimeStep",ts->time_step,&ts->time_step,PETSC_NULL);
92: opt = ts->exact_final_time == PETSC_DECIDE ? PETSC_FALSE : (PetscBool)ts->exact_final_time;
93: PetscOptionsBool("-ts_exact_final_time","Interpolate output to stop exactly at the final time","TSSetExactFinalTime",opt,&opt,&flg);
94: if (flg) {TSSetExactFinalTime(ts,opt);}
95: PetscOptionsInt("-ts_max_snes_failures","Maximum number of nonlinear solve failures","TSSetMaxSNESFailures",ts->max_snes_failures,&ts->max_snes_failures,PETSC_NULL);
96: PetscOptionsInt("-ts_max_reject","Maximum number of step rejections before step fails","TSSetMaxStepRejections",ts->max_reject,&ts->max_reject,PETSC_NULL);
97: PetscOptionsBool("-ts_error_if_step_fails","Error if no step succeeds","TSSetErrorIfStepFails",ts->errorifstepfailed,&ts->errorifstepfailed,PETSC_NULL);
98: PetscOptionsReal("-ts_rtol","Relative tolerance for local truncation error","TSSetTolerances",ts->rtol,&ts->rtol,PETSC_NULL);
99: PetscOptionsReal("-ts_atol","Absolute tolerance for local truncation error","TSSetTolerances",ts->atol,&ts->atol,PETSC_NULL);
101: /* Monitor options */
102: PetscOptionsString("-ts_monitor","Monitor timestep size","TSMonitorDefault","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
103: if (flg) {
104: PetscViewerASCIIOpen(((PetscObject)ts)->comm,monfilename,&monviewer);
105: TSMonitorSet(ts,TSMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);
106: }
107: PetscOptionsString("-ts_monitor_python","Use Python function","TSMonitorSet",0,monfilename,PETSC_MAX_PATH_LEN,&flg);
108: if (flg) {PetscPythonMonitorSet((PetscObject)ts,monfilename);}
110: opt = PETSC_FALSE;
111: PetscOptionsBool("-ts_monitor_draw","Monitor timestep size graphically","TSMonitorLG",opt,&opt,PETSC_NULL);
112: if (opt) {
113: TSMonitorSet(ts,TSMonitorLG,PETSC_NULL,PETSC_NULL);
114: }
115: opt = PETSC_FALSE;
116: PetscOptionsBool("-ts_monitor_solution","Monitor solution graphically","TSMonitorSolution",opt,&opt,PETSC_NULL);
117: if (opt) {
118: void *ctx;
119: TSMonitorSolutionCreate(ts,PETSC_NULL,&ctx);
120: TSMonitorSet(ts,TSMonitorSolution,ctx,TSMonitorSolutionDestroy);
121: }
122: opt = PETSC_FALSE;
123: PetscOptionsString("-ts_monitor_solution_binary","Save each solution to a binary file","TSMonitorSolutionBinary",0,monfilename,PETSC_MAX_PATH_LEN,&flg);
124: if (flg) {
125: PetscViewer ctx;
126: if (monfilename[0]) {
127: PetscViewerBinaryOpen(((PetscObject)ts)->comm,monfilename,FILE_MODE_WRITE,&ctx);
128: } else {
129: ctx = PETSC_VIEWER_BINARY_(((PetscObject)ts)->comm);
130: }
131: TSMonitorSet(ts,TSMonitorSolutionBinary,ctx,(PetscErrorCode (*)(void**))PetscViewerDestroy);
132: }
133: opt = PETSC_FALSE;
134: PetscOptionsString("-ts_monitor_solution_vtk","Save each time step to a binary file, use filename-%%03D.vts","TSMonitorSolutionVTK",0,monfilename,PETSC_MAX_PATH_LEN,&flg);
135: if (flg) {
136: const char *ptr,*ptr2;
137: char *filetemplate;
138: if (!monfilename[0]) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"-ts_monitor_solution_vtk requires a file template, e.g. filename-%%03D.vts");
139: /* Do some cursory validation of the input. */
140: PetscStrstr(monfilename,"%",(char**)&ptr);
141: if (!ptr) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"-ts_monitor_solution_vtk requires a file template, e.g. filename-%%03D.vts");
142: for (ptr++ ; ptr && *ptr; ptr++) {
143: PetscStrchr("DdiouxX",*ptr,(char**)&ptr2);
144: if (!ptr2 && (*ptr < '0' || '9' < *ptr)) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Invalid file template argument to -ts_monitor_solution_vtk, should look like filename-%%03D.vts");
145: if (ptr2) break;
146: }
147: PetscStrallocpy(monfilename,&filetemplate);
148: TSMonitorSet(ts,TSMonitorSolutionVTK,filetemplate,(PetscErrorCode (*)(void**))TSMonitorSolutionVTKDestroy);
149: }
151: TSGetAdapt(ts,&adapt);
152: TSAdaptSetFromOptions(adapt);
154: TSGetSNES(ts,&snes);
155: if (ts->problem_type == TS_LINEAR) {SNESSetType(snes,SNESKSPONLY);}
157: /* Handle specific TS options */
158: if (ts->ops->setfromoptions) {
159: (*ts->ops->setfromoptions)(ts);
160: }
162: /* process any options handlers added with PetscObjectAddOptionsHandler() */
163: PetscObjectProcessOptionsHandlers((PetscObject)ts);
164: PetscOptionsEnd();
165: return(0);
166: }
171: /*@
172: TSComputeRHSJacobian - Computes the Jacobian matrix that has been
173: set with TSSetRHSJacobian().
175: Collective on TS and Vec
177: Input Parameters:
178: + ts - the TS context
179: . t - current timestep
180: - x - input vector
182: Output Parameters:
183: + A - Jacobian matrix
184: . B - optional preconditioning matrix
185: - flag - flag indicating matrix structure
187: Notes:
188: Most users should not need to explicitly call this routine, as it
189: is used internally within the nonlinear solvers.
191: See KSPSetOperators() for important information about setting the
192: flag parameter.
194: Level: developer
196: .keywords: SNES, compute, Jacobian, matrix
198: .seealso: TSSetRHSJacobian(), KSPSetOperators()
199: @*/
200: PetscErrorCode TSComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg)
201: {
203: PetscInt Xstate;
209: PetscObjectStateQuery((PetscObject)X,&Xstate);
210: if (ts->rhsjacobian.time == t && (ts->problem_type == TS_LINEAR || (ts->rhsjacobian.X == X && ts->rhsjacobian.Xstate == Xstate))) {
211: *flg = ts->rhsjacobian.mstructure;
212: return(0);
213: }
215: if (!ts->userops->rhsjacobian && !ts->userops->ijacobian) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSJacobian() and / or TSSetIJacobian()");
217: if (ts->userops->rhsjacobian) {
218: PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);
219: *flg = DIFFERENT_NONZERO_PATTERN;
220: PetscStackPush("TS user Jacobian function");
221: (*ts->userops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);
222: PetscStackPop;
223: PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);
224: /* make sure user returned a correct Jacobian and preconditioner */
227: } else {
228: MatZeroEntries(*A);
229: if (*A != *B) {MatZeroEntries(*B);}
230: *flg = SAME_NONZERO_PATTERN;
231: }
232: ts->rhsjacobian.time = t;
233: ts->rhsjacobian.X = X;
234: PetscObjectStateQuery((PetscObject)X,&ts->rhsjacobian.Xstate);
235: ts->rhsjacobian.mstructure = *flg;
236: return(0);
237: }
241: /*@
242: TSComputeRHSFunction - Evaluates the right-hand-side function.
244: Collective on TS and Vec
246: Input Parameters:
247: + ts - the TS context
248: . t - current time
249: - x - state vector
251: Output Parameter:
252: . y - right hand side
254: Note:
255: Most users should not need to explicitly call this routine, as it
256: is used internally within the nonlinear solvers.
258: Level: developer
260: .keywords: TS, compute
262: .seealso: TSSetRHSFunction(), TSComputeIFunction()
263: @*/
264: PetscErrorCode TSComputeRHSFunction(TS ts,PetscReal t,Vec x,Vec y)
265: {
273: if (!ts->userops->rhsfunction && !ts->userops->ifunction) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSFunction() and / or TSSetIFunction()");
275: PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);
276: if (ts->userops->rhsfunction) {
277: PetscStackPush("TS user right-hand-side function");
278: (*ts->userops->rhsfunction)(ts,t,x,y,ts->funP);
279: PetscStackPop;
280: } else {
281: VecZeroEntries(y);
282: }
284: PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);
285: return(0);
286: }
290: static PetscErrorCode TSGetRHSVec_Private(TS ts,Vec *Frhs)
291: {
292: Vec F;
296: *Frhs = PETSC_NULL;
297: TSGetIFunction(ts,&F,PETSC_NULL,PETSC_NULL);
298: if (!ts->Frhs) {
299: VecDuplicate(F,&ts->Frhs);
300: }
301: *Frhs = ts->Frhs;
302: return(0);
303: }
307: static PetscErrorCode TSGetRHSMats_Private(TS ts,Mat *Arhs,Mat *Brhs)
308: {
309: Mat A,B;
313: TSGetIJacobian(ts,&A,&B,PETSC_NULL,PETSC_NULL);
314: if (Arhs) {
315: if (!ts->Arhs) {
316: MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&ts->Arhs);
317: }
318: *Arhs = ts->Arhs;
319: }
320: if (Brhs) {
321: if (!ts->Brhs) {
322: MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&ts->Brhs);
323: }
324: *Brhs = ts->Brhs;
325: }
326: return(0);
327: }
331: /*@
332: TSComputeIFunction - Evaluates the DAE residual written in implicit form F(t,X,Xdot)=0
334: Collective on TS and Vec
336: Input Parameters:
337: + ts - the TS context
338: . t - current time
339: . X - state vector
340: . Xdot - time derivative of state vector
341: - imex - flag indicates if the method is IMEX so that the RHSFunction should be kept separate
343: Output Parameter:
344: . Y - right hand side
346: Note:
347: Most users should not need to explicitly call this routine, as it
348: is used internally within the nonlinear solvers.
350: If the user did did not write their equations in implicit form, this
351: function recasts them in implicit form.
353: Level: developer
355: .keywords: TS, compute
357: .seealso: TSSetIFunction(), TSComputeRHSFunction()
358: @*/
359: PetscErrorCode TSComputeIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec Y,PetscBool imex)
360: {
369: if (!ts->userops->rhsfunction && !ts->userops->ifunction) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSFunction() and / or TSSetIFunction()");
371: PetscLogEventBegin(TS_FunctionEval,ts,X,Xdot,Y);
372: if (ts->userops->ifunction) {
373: PetscStackPush("TS user implicit function");
374: (*ts->userops->ifunction)(ts,t,X,Xdot,Y,ts->funP);
375: PetscStackPop;
376: }
377: if (imex) {
378: if (!ts->userops->ifunction) {
379: VecCopy(Xdot,Y);
380: }
381: } else if (ts->userops->rhsfunction) {
382: if (ts->userops->ifunction) {
383: Vec Frhs;
384: TSGetRHSVec_Private(ts,&Frhs);
385: TSComputeRHSFunction(ts,t,X,Frhs);
386: VecAXPY(Y,-1,Frhs);
387: } else {
388: TSComputeRHSFunction(ts,t,X,Y);
389: VecAYPX(Y,-1,Xdot);
390: }
391: }
392: PetscLogEventEnd(TS_FunctionEval,ts,X,Xdot,Y);
393: return(0);
394: }
398: /*@
399: TSComputeIJacobian - Evaluates the Jacobian of the DAE
401: Collective on TS and Vec
403: Input
404: Input Parameters:
405: + ts - the TS context
406: . t - current timestep
407: . X - state vector
408: . Xdot - time derivative of state vector
409: . shift - shift to apply, see note below
410: - imex - flag indicates if the method is IMEX so that the RHSJacobian should be kept separate
412: Output Parameters:
413: + A - Jacobian matrix
414: . B - optional preconditioning matrix
415: - flag - flag indicating matrix structure
417: Notes:
418: If F(t,X,Xdot)=0 is the DAE, the required Jacobian is
420: dF/dX + shift*dF/dXdot
422: Most users should not need to explicitly call this routine, as it
423: is used internally within the nonlinear solvers.
425: Level: developer
427: .keywords: TS, compute, Jacobian, matrix
429: .seealso: TSSetIJacobian()
430: @*/
431: PetscErrorCode TSComputeIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg,PetscBool imex)
432: {
433: PetscInt Xstate, Xdotstate;
445: PetscObjectStateQuery((PetscObject)X,&Xstate);
446: PetscObjectStateQuery((PetscObject)Xdot,&Xdotstate);
447: if (ts->ijacobian.time == t && (ts->problem_type == TS_LINEAR || (ts->ijacobian.X == X && ts->ijacobian.Xstate == Xstate && ts->ijacobian.Xdot == Xdot && ts->ijacobian.Xdotstate == Xdotstate && ts->ijacobian.imex == imex))) {
448: *flg = ts->ijacobian.mstructure;
449: MatScale(*A, shift / ts->ijacobian.shift);
450: return(0);
451: }
453: if (!ts->userops->rhsjacobian && !ts->userops->ijacobian) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSJacobian() and / or TSSetIJacobian()");
455: *flg = SAME_NONZERO_PATTERN; /* In case we're solving a linear problem in which case it wouldn't get initialized below. */
456: PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);
457: if (ts->userops->ijacobian) {
458: *flg = DIFFERENT_NONZERO_PATTERN;
459: PetscStackPush("TS user implicit Jacobian");
460: (*ts->userops->ijacobian)(ts,t,X,Xdot,shift,A,B,flg,ts->jacP);
461: PetscStackPop;
462: /* make sure user returned a correct Jacobian and preconditioner */
465: }
466: if (imex) {
467: if (!ts->userops->ijacobian) { /* system was written as Xdot = F(t,X) */
468: MatZeroEntries(*A);
469: MatShift(*A,shift);
470: if (*A != *B) {
471: MatZeroEntries(*B);
472: MatShift(*B,shift);
473: }
474: *flg = SAME_PRECONDITIONER;
475: }
476: } else {
477: if (!ts->userops->ijacobian) {
478: TSComputeRHSJacobian(ts,t,X,A,B,flg);
479: MatScale(*A,-1);
480: MatShift(*A,shift);
481: if (*A != *B) {
482: MatScale(*B,-1);
483: MatShift(*B,shift);
484: }
485: } else if (ts->userops->rhsjacobian) {
486: Mat Arhs,Brhs;
487: MatStructure axpy,flg2 = DIFFERENT_NONZERO_PATTERN;
488: TSGetRHSMats_Private(ts,&Arhs,&Brhs);
489: TSComputeRHSJacobian(ts,t,X,&Arhs,&Brhs,&flg2);
490: axpy = (*flg == flg2) ? SAME_NONZERO_PATTERN : DIFFERENT_NONZERO_PATTERN;
491: MatAXPY(*A,-1,Arhs,axpy);
492: if (*A != *B) {
493: MatAXPY(*B,-1,Brhs,axpy);
494: }
495: *flg = PetscMin(*flg,flg2);
496: }
497: }
499: ts->ijacobian.time = t;
500: ts->ijacobian.X = X;
501: ts->ijacobian.Xdot = Xdot;
502: PetscObjectStateQuery((PetscObject)X,&ts->ijacobian.Xstate);
503: PetscObjectStateQuery((PetscObject)Xdot,&ts->ijacobian.Xdotstate);
504: ts->ijacobian.shift = shift;
505: ts->ijacobian.imex = imex;
506: ts->ijacobian.mstructure = *flg;
507: PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);
508: return(0);
509: }
513: /*@C
514: TSSetRHSFunction - Sets the routine for evaluating the function,
515: F(t,u), where U_t = F(t,u).
517: Logically Collective on TS
519: Input Parameters:
520: + ts - the TS context obtained from TSCreate()
521: . r - vector to put the computed right hand side (or PETSC_NULL to have it created)
522: . f - routine for evaluating the right-hand-side function
523: - ctx - [optional] user-defined context for private data for the
524: function evaluation routine (may be PETSC_NULL)
526: Calling sequence of func:
527: $ func (TS ts,PetscReal t,Vec u,Vec F,void *ctx);
529: + t - current timestep
530: . u - input vector
531: . F - function vector
532: - ctx - [optional] user-defined function context
534: Level: beginner
536: .keywords: TS, timestep, set, right-hand-side, function
538: .seealso: TSSetRHSJacobian(), TSSetIJacobian()
539: @*/
540: PetscErrorCode TSSetRHSFunction(TS ts,Vec r,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx)
541: {
543: SNES snes;
544: Vec ralloc = PETSC_NULL;
549: if (f) ts->userops->rhsfunction = f;
550: if (ctx) ts->funP = ctx;
551: TSGetSNES(ts,&snes);
552: if (!r && !ts->dm && ts->vec_sol) {
553: VecDuplicate(ts->vec_sol,&ralloc);
554: r = ralloc;
555: }
556: SNESSetFunction(snes,r,SNESTSFormFunction,ts);
557: VecDestroy(&ralloc);
558: return(0);
559: }
563: /*@C
564: TSSetRHSJacobian - Sets the function to compute the Jacobian of F,
565: where U_t = F(U,t), as well as the location to store the matrix.
567: Logically Collective on TS
569: Input Parameters:
570: + ts - the TS context obtained from TSCreate()
571: . A - Jacobian matrix
572: . B - preconditioner matrix (usually same as A)
573: . f - the Jacobian evaluation routine
574: - ctx - [optional] user-defined context for private data for the
575: Jacobian evaluation routine (may be PETSC_NULL)
577: Calling sequence of func:
578: $ func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx);
580: + t - current timestep
581: . u - input vector
582: . A - matrix A, where U_t = A(t)u
583: . B - preconditioner matrix, usually the same as A
584: . flag - flag indicating information about the preconditioner matrix
585: structure (same as flag in KSPSetOperators())
586: - ctx - [optional] user-defined context for matrix evaluation routine
588: Notes:
589: See KSPSetOperators() for important information about setting the flag
590: output parameter in the routine func(). Be sure to read this information!
592: The routine func() takes Mat * as the matrix arguments rather than Mat.
593: This allows the matrix evaluation routine to replace A and/or B with a
594: completely new matrix structure (not just different matrix elements)
595: when appropriate, for instance, if the nonzero structure is changing
596: throughout the global iterations.
598: Level: beginner
599:
600: .keywords: TS, timestep, set, right-hand-side, Jacobian
602: .seealso: SNESDefaultComputeJacobianColor(), TSSetRHSFunction()
604: @*/
605: PetscErrorCode TSSetRHSJacobian(TS ts,Mat A,Mat B,TSRHSJacobian f,void *ctx)
606: {
608: SNES snes;
617: if (f) ts->userops->rhsjacobian = f;
618: if (ctx) ts->jacP = ctx;
619: TSGetSNES(ts,&snes);
620: if (!ts->userops->ijacobian) {
621: SNESSetJacobian(snes,A,B,SNESTSFormJacobian,ts);
622: }
623: if (A) {
624: PetscObjectReference((PetscObject)A);
625: MatDestroy(&ts->Arhs);
626: ts->Arhs = A;
627: }
628: if (B) {
629: PetscObjectReference((PetscObject)B);
630: MatDestroy(&ts->Brhs);
631: ts->Brhs = B;
632: }
633: return(0);
634: }
639: /*@C
640: TSSetIFunction - Set the function to compute F(t,U,U_t) where F = 0 is the DAE to be solved.
642: Logically Collective on TS
644: Input Parameters:
645: + ts - the TS context obtained from TSCreate()
646: . r - vector to hold the residual (or PETSC_NULL to have it created internally)
647: . f - the function evaluation routine
648: - ctx - user-defined context for private data for the function evaluation routine (may be PETSC_NULL)
650: Calling sequence of f:
651: $ f(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx);
653: + t - time at step/stage being solved
654: . u - state vector
655: . u_t - time derivative of state vector
656: . F - function vector
657: - ctx - [optional] user-defined context for matrix evaluation routine
659: Important:
660: The user MUST call either this routine, TSSetRHSFunction(). This routine must be used when not solving an ODE, for example a DAE.
662: Level: beginner
664: .keywords: TS, timestep, set, DAE, Jacobian
666: .seealso: TSSetRHSJacobian(), TSSetRHSFunction(), TSSetIJacobian()
667: @*/
668: PetscErrorCode TSSetIFunction(TS ts,Vec res,TSIFunction f,void *ctx)
669: {
671: SNES snes;
672: Vec resalloc = PETSC_NULL;
677: if (f) ts->userops->ifunction = f;
678: if (ctx) ts->funP = ctx;
679: TSGetSNES(ts,&snes);
680: if (!res && !ts->dm && ts->vec_sol) {
681: VecDuplicate(ts->vec_sol,&resalloc);
682: res = resalloc;
683: }
684: SNESSetFunction(snes,res,SNESTSFormFunction,ts);
685: VecDestroy(&resalloc);
686: return(0);
687: }
691: /*@C
692: TSGetIFunction - Returns the vector where the implicit residual is stored and the function/contex to compute it.
694: Not Collective
696: Input Parameter:
697: . ts - the TS context
699: Output Parameter:
700: + r - vector to hold residual (or PETSC_NULL)
701: . func - the function to compute residual (or PETSC_NULL)
702: - ctx - the function context (or PETSC_NULL)
704: Level: advanced
706: .keywords: TS, nonlinear, get, function
708: .seealso: TSSetIFunction(), SNESGetFunction()
709: @*/
710: PetscErrorCode TSGetIFunction(TS ts,Vec *r,TSIFunction *func,void **ctx)
711: {
713: SNES snes;
717: TSGetSNES(ts,&snes);
718: SNESGetFunction(snes,r,PETSC_NULL,PETSC_NULL);
719: if (func) *func = ts->userops->ifunction;
720: if (ctx) *ctx = ts->funP;
721: return(0);
722: }
726: /*@C
727: TSGetRHSFunction - Returns the vector where the right hand side is stored and the function/context to compute it.
729: Not Collective
731: Input Parameter:
732: . ts - the TS context
734: Output Parameter:
735: + r - vector to hold computed right hand side (or PETSC_NULL)
736: . func - the function to compute right hand side (or PETSC_NULL)
737: - ctx - the function context (or PETSC_NULL)
739: Level: advanced
741: .keywords: TS, nonlinear, get, function
743: .seealso: TSSetRhsfunction(), SNESGetFunction()
744: @*/
745: PetscErrorCode TSGetRHSFunction(TS ts,Vec *r,TSRHSFunction *func,void **ctx)
746: {
748: SNES snes;
752: TSGetSNES(ts,&snes);
753: SNESGetFunction(snes,r,PETSC_NULL,PETSC_NULL);
754: if (func) *func = ts->userops->rhsfunction;
755: if (ctx) *ctx = ts->funP;
756: return(0);
757: }
761: /*@C
762: TSSetIJacobian - Set the function to compute the matrix dF/dU + a*dF/dU_t where F(t,U,U_t) is the function
763: you provided with TSSetIFunction().
765: Logically Collective on TS
767: Input Parameters:
768: + ts - the TS context obtained from TSCreate()
769: . A - Jacobian matrix
770: . B - preconditioning matrix for A (may be same as A)
771: . f - the Jacobian evaluation routine
772: - ctx - user-defined context for private data for the Jacobian evaluation routine (may be PETSC_NULL)
774: Calling sequence of f:
775: $ f(TS ts,PetscReal t,Vec U,Vec U_t,PetscReal a,Mat *A,Mat *B,MatStructure *flag,void *ctx);
777: + t - time at step/stage being solved
778: . U - state vector
779: . U_t - time derivative of state vector
780: . a - shift
781: . A - Jacobian of G(U) = F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t
782: . B - preconditioning matrix for A, may be same as A
783: . flag - flag indicating information about the preconditioner matrix
784: structure (same as flag in KSPSetOperators())
785: - ctx - [optional] user-defined context for matrix evaluation routine
787: Notes:
788: The matrices A and B are exactly the matrices that are used by SNES for the nonlinear solve.
790: The matrix dF/dU + a*dF/dU_t you provide turns out to be
791: the Jacobian of G(U) = F(t,U,W+a*U) where F(t,U,U_t) = 0 is the DAE to be solved.
792: The time integrator internally approximates U_t by W+a*U where the positive "shift"
793: a and vector W depend on the integration method, step size, and past states. For example with
794: the backward Euler method a = 1/dt and W = -a*U(previous timestep) so
795: W + a*U = a*(U - U(previous timestep)) = (U - U(previous timestep))/dt
797: Level: beginner
799: .keywords: TS, timestep, DAE, Jacobian
801: .seealso: TSSetIFunction(), TSSetRHSJacobian(), SNESDefaultComputeJacobianColor(), SNESDefaultComputeJacobian()
803: @*/
804: PetscErrorCode TSSetIJacobian(TS ts,Mat A,Mat B,TSIJacobian f,void *ctx)
805: {
807: SNES snes;
815: if (f) ts->userops->ijacobian = f;
816: if (ctx) ts->jacP = ctx;
817: TSGetSNES(ts,&snes);
818: SNESSetJacobian(snes,A,B,SNESTSFormJacobian,ts);
819: return(0);
820: }
824: /*@C
825: TSView - Prints the TS data structure.
827: Collective on TS
829: Input Parameters:
830: + ts - the TS context obtained from TSCreate()
831: - viewer - visualization context
833: Options Database Key:
834: . -ts_view - calls TSView() at end of TSStep()
836: Notes:
837: The available visualization contexts include
838: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
839: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
840: output where only the first processor opens
841: the file. All other processors send their
842: data to the first processor to print.
844: The user can open an alternative visualization context with
845: PetscViewerASCIIOpen() - output to a specified file.
847: Level: beginner
849: .keywords: TS, timestep, view
851: .seealso: PetscViewerASCIIOpen()
852: @*/
853: PetscErrorCode TSView(TS ts,PetscViewer viewer)
854: {
856: const TSType type;
857: PetscBool iascii,isstring,isundials;
861: if (!viewer) {
862: PetscViewerASCIIGetStdout(((PetscObject)ts)->comm,&viewer);
863: }
867: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
868: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
869: if (iascii) {
870: PetscObjectPrintClassNamePrefixType((PetscObject)ts,viewer,"TS Object");
871: PetscViewerASCIIPrintf(viewer," maximum steps=%D\n",ts->max_steps);
872: PetscViewerASCIIPrintf(viewer," maximum time=%G\n",ts->max_time);
873: if (ts->problem_type == TS_NONLINEAR) {
874: PetscViewerASCIIPrintf(viewer," total number of nonlinear solver iterations=%D\n",ts->snes_its);
875: PetscViewerASCIIPrintf(viewer," total number of nonlinear solve failures=%D\n",ts->num_snes_failures);
876: }
877: PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%D\n",ts->ksp_its);
878: PetscViewerASCIIPrintf(viewer," total number of rejected steps=%D\n",ts->reject);
879: if (ts->ops->view) {
880: PetscViewerASCIIPushTab(viewer);
881: (*ts->ops->view)(ts,viewer);
882: PetscViewerASCIIPopTab(viewer);
883: }
884: } else if (isstring) {
885: TSGetType(ts,&type);
886: PetscViewerStringSPrintf(viewer," %-7.7s",type);
887: }
888: PetscViewerASCIIPushTab(viewer);
889: PetscObjectTypeCompare((PetscObject)ts,TSSUNDIALS,&isundials);
890: PetscViewerASCIIPopTab(viewer);
891: return(0);
892: }
897: /*@
898: TSSetApplicationContext - Sets an optional user-defined context for
899: the timesteppers.
901: Logically Collective on TS
903: Input Parameters:
904: + ts - the TS context obtained from TSCreate()
905: - usrP - optional user context
907: Level: intermediate
909: .keywords: TS, timestep, set, application, context
911: .seealso: TSGetApplicationContext()
912: @*/
913: PetscErrorCode TSSetApplicationContext(TS ts,void *usrP)
914: {
917: ts->user = usrP;
918: return(0);
919: }
923: /*@
924: TSGetApplicationContext - Gets the user-defined context for the
925: timestepper.
927: Not Collective
929: Input Parameter:
930: . ts - the TS context obtained from TSCreate()
932: Output Parameter:
933: . usrP - user context
935: Level: intermediate
937: .keywords: TS, timestep, get, application, context
939: .seealso: TSSetApplicationContext()
940: @*/
941: PetscErrorCode TSGetApplicationContext(TS ts,void *usrP)
942: {
945: *(void**)usrP = ts->user;
946: return(0);
947: }
951: /*@
952: TSGetTimeStepNumber - Gets the current number of timesteps.
954: Not Collective
956: Input Parameter:
957: . ts - the TS context obtained from TSCreate()
959: Output Parameter:
960: . iter - number steps so far
962: Level: intermediate
964: .keywords: TS, timestep, get, iteration, number
965: @*/
966: PetscErrorCode TSGetTimeStepNumber(TS ts,PetscInt* iter)
967: {
971: *iter = ts->steps;
972: return(0);
973: }
977: /*@
978: TSSetInitialTimeStep - Sets the initial timestep to be used,
979: as well as the initial time.
981: Logically Collective on TS
983: Input Parameters:
984: + ts - the TS context obtained from TSCreate()
985: . initial_time - the initial time
986: - time_step - the size of the timestep
988: Level: intermediate
990: .seealso: TSSetTimeStep(), TSGetTimeStep()
992: .keywords: TS, set, initial, timestep
993: @*/
994: PetscErrorCode TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step)
995: {
1000: TSSetTimeStep(ts,time_step);
1001: TSSetTime(ts,initial_time);
1002: return(0);
1003: }
1007: /*@
1008: TSSetTimeStep - Allows one to reset the timestep at any time,
1009: useful for simple pseudo-timestepping codes.
1011: Logically Collective on TS
1013: Input Parameters:
1014: + ts - the TS context obtained from TSCreate()
1015: - time_step - the size of the timestep
1017: Level: intermediate
1019: .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1021: .keywords: TS, set, timestep
1022: @*/
1023: PetscErrorCode TSSetTimeStep(TS ts,PetscReal time_step)
1024: {
1028: ts->time_step = time_step;
1029: return(0);
1030: }
1034: /*@
1035: TSSetExactFinalTime - Determines whether to interpolate solution to the
1036: exact final time requested by the user or just returns it at the final time
1037: it computed.
1039: Logically Collective on TS
1041: Input Parameter:
1042: + ts - the time-step context
1043: - ft - PETSC_TRUE if interpolates, else PETSC_FALSE
1045: Level: beginner
1047: .seealso: TSSetDuration()
1048: @*/
1049: PetscErrorCode TSSetExactFinalTime(TS ts,PetscBool flg)
1050: {
1055: ts->exact_final_time = flg;
1056: return(0);
1057: }
1061: /*@
1062: TSGetTimeStep - Gets the current timestep size.
1064: Not Collective
1066: Input Parameter:
1067: . ts - the TS context obtained from TSCreate()
1069: Output Parameter:
1070: . dt - the current timestep size
1072: Level: intermediate
1074: .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1076: .keywords: TS, get, timestep
1077: @*/
1078: PetscErrorCode TSGetTimeStep(TS ts,PetscReal* dt)
1079: {
1083: *dt = ts->time_step;
1084: return(0);
1085: }
1089: /*@
1090: TSGetSolution - Returns the solution at the present timestep. It
1091: is valid to call this routine inside the function that you are evaluating
1092: in order to move to the new timestep. This vector not changed until
1093: the solution at the next timestep has been calculated.
1095: Not Collective, but Vec returned is parallel if TS is parallel
1097: Input Parameter:
1098: . ts - the TS context obtained from TSCreate()
1100: Output Parameter:
1101: . v - the vector containing the solution
1103: Level: intermediate
1105: .seealso: TSGetTimeStep()
1107: .keywords: TS, timestep, get, solution
1108: @*/
1109: PetscErrorCode TSGetSolution(TS ts,Vec *v)
1110: {
1114: *v = ts->vec_sol;
1115: return(0);
1116: }
1118: /* ----- Routines to initialize and destroy a timestepper ---- */
1121: /*@
1122: TSSetProblemType - Sets the type of problem to be solved.
1124: Not collective
1126: Input Parameters:
1127: + ts - The TS
1128: - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
1129: .vb
1130: U_t = A U
1131: U_t = A(t) U
1132: U_t = F(t,U)
1133: .ve
1135: Level: beginner
1137: .keywords: TS, problem type
1138: .seealso: TSSetUp(), TSProblemType, TS
1139: @*/
1140: PetscErrorCode TSSetProblemType(TS ts, TSProblemType type)
1141: {
1146: ts->problem_type = type;
1147: if (type == TS_LINEAR) {
1148: SNES snes;
1149: TSGetSNES(ts,&snes);
1150: SNESSetType(snes,SNESKSPONLY);
1151: }
1152: return(0);
1153: }
1157: /*@C
1158: TSGetProblemType - Gets the type of problem to be solved.
1160: Not collective
1162: Input Parameter:
1163: . ts - The TS
1165: Output Parameter:
1166: . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
1167: .vb
1168: M U_t = A U
1169: M(t) U_t = A(t) U
1170: U_t = F(t,U)
1171: .ve
1173: Level: beginner
1175: .keywords: TS, problem type
1176: .seealso: TSSetUp(), TSProblemType, TS
1177: @*/
1178: PetscErrorCode TSGetProblemType(TS ts, TSProblemType *type)
1179: {
1183: *type = ts->problem_type;
1184: return(0);
1185: }
1189: /*@
1190: TSSetUp - Sets up the internal data structures for the later use
1191: of a timestepper.
1193: Collective on TS
1195: Input Parameter:
1196: . ts - the TS context obtained from TSCreate()
1198: Notes:
1199: For basic use of the TS solvers the user need not explicitly call
1200: TSSetUp(), since these actions will automatically occur during
1201: the call to TSStep(). However, if one wishes to control this
1202: phase separately, TSSetUp() should be called after TSCreate()
1203: and optional routines of the form TSSetXXX(), but before TSStep().
1205: Level: advanced
1207: .keywords: TS, timestep, setup
1209: .seealso: TSCreate(), TSStep(), TSDestroy()
1210: @*/
1211: PetscErrorCode TSSetUp(TS ts)
1212: {
1217: if (ts->setupcalled) return(0);
1219: if (!((PetscObject)ts)->type_name) {
1220: TSSetType(ts,TSEULER);
1221: }
1222: if (ts->exact_final_time == PETSC_DECIDE) ts->exact_final_time = PETSC_FALSE;
1224: if (!ts->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first");
1226: TSGetAdapt(ts,&ts->adapt);
1228: if (ts->ops->setup) {
1229: (*ts->ops->setup)(ts);
1230: }
1232: ts->setupcalled = PETSC_TRUE;
1233: return(0);
1234: }
1238: /*@
1239: TSReset - Resets a TS context and removes any allocated Vecs and Mats.
1241: Collective on TS
1243: Input Parameter:
1244: . ts - the TS context obtained from TSCreate()
1246: Level: beginner
1248: .keywords: TS, timestep, reset
1250: .seealso: TSCreate(), TSSetup(), TSDestroy()
1251: @*/
1252: PetscErrorCode TSReset(TS ts)
1253: {
1258: if (ts->ops->reset) {
1259: (*ts->ops->reset)(ts);
1260: }
1261: if (ts->snes) {SNESReset(ts->snes);}
1262: MatDestroy(&ts->Arhs);
1263: MatDestroy(&ts->Brhs);
1264: VecDestroy(&ts->Frhs);
1265: VecDestroy(&ts->vec_sol);
1266: VecDestroyVecs(ts->nwork,&ts->work);
1267: ts->setupcalled = PETSC_FALSE;
1268: return(0);
1269: }
1273: /*@
1274: TSDestroy - Destroys the timestepper context that was created
1275: with TSCreate().
1277: Collective on TS
1279: Input Parameter:
1280: . ts - the TS context obtained from TSCreate()
1282: Level: beginner
1284: .keywords: TS, timestepper, destroy
1286: .seealso: TSCreate(), TSSetUp(), TSSolve()
1287: @*/
1288: PetscErrorCode TSDestroy(TS *ts)
1289: {
1293: if (!*ts) return(0);
1295: if (--((PetscObject)(*ts))->refct > 0) {*ts = 0; return(0);}
1297: TSReset((*ts));
1299: /* if memory was published with AMS then destroy it */
1300: PetscObjectDepublish((*ts));
1301: if ((*ts)->ops->destroy) {(*(*ts)->ops->destroy)((*ts));}
1303: TSAdaptDestroy(&(*ts)->adapt);
1304: SNESDestroy(&(*ts)->snes);
1305: DMDestroy(&(*ts)->dm);
1306: TSMonitorCancel((*ts));
1308: PetscFree((*ts)->userops);
1310: PetscHeaderDestroy(ts);
1311: return(0);
1312: }
1316: /*@
1317: TSGetSNES - Returns the SNES (nonlinear solver) associated with
1318: a TS (timestepper) context. Valid only for nonlinear problems.
1320: Not Collective, but SNES is parallel if TS is parallel
1322: Input Parameter:
1323: . ts - the TS context obtained from TSCreate()
1325: Output Parameter:
1326: . snes - the nonlinear solver context
1328: Notes:
1329: The user can then directly manipulate the SNES context to set various
1330: options, etc. Likewise, the user can then extract and manipulate the
1331: KSP, KSP, and PC contexts as well.
1333: TSGetSNES() does not work for integrators that do not use SNES; in
1334: this case TSGetSNES() returns PETSC_NULL in snes.
1336: Level: beginner
1338: .keywords: timestep, get, SNES
1339: @*/
1340: PetscErrorCode TSGetSNES(TS ts,SNES *snes)
1341: {
1347: if (!ts->snes) {
1348: SNESCreate(((PetscObject)ts)->comm,&ts->snes);
1349: PetscLogObjectParent(ts,ts->snes);
1350: PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);
1351: if (ts->dm) {SNESSetDM(ts->snes,ts->dm);}
1352: if (ts->problem_type == TS_LINEAR) {
1353: SNESSetType(ts->snes,SNESKSPONLY);
1354: }
1355: }
1356: *snes = ts->snes;
1357: return(0);
1358: }
1362: /*@
1363: TSGetKSP - Returns the KSP (linear solver) associated with
1364: a TS (timestepper) context.
1366: Not Collective, but KSP is parallel if TS is parallel
1368: Input Parameter:
1369: . ts - the TS context obtained from TSCreate()
1371: Output Parameter:
1372: . ksp - the nonlinear solver context
1374: Notes:
1375: The user can then directly manipulate the KSP context to set various
1376: options, etc. Likewise, the user can then extract and manipulate the
1377: KSP and PC contexts as well.
1379: TSGetKSP() does not work for integrators that do not use KSP;
1380: in this case TSGetKSP() returns PETSC_NULL in ksp.
1382: Level: beginner
1384: .keywords: timestep, get, KSP
1385: @*/
1386: PetscErrorCode TSGetKSP(TS ts,KSP *ksp)
1387: {
1389: SNES snes;
1394: if (!((PetscObject)ts)->type_name) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"KSP is not created yet. Call TSSetType() first");
1395: if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()");
1396: TSGetSNES(ts,&snes);
1397: SNESGetKSP(snes,ksp);
1398: return(0);
1399: }
1401: /* ----------- Routines to set solver parameters ---------- */
1405: /*@
1406: TSGetDuration - Gets the maximum number of timesteps to use and
1407: maximum time for iteration.
1409: Not Collective
1411: Input Parameters:
1412: + ts - the TS context obtained from TSCreate()
1413: . maxsteps - maximum number of iterations to use, or PETSC_NULL
1414: - maxtime - final time to iterate to, or PETSC_NULL
1416: Level: intermediate
1418: .keywords: TS, timestep, get, maximum, iterations, time
1419: @*/
1420: PetscErrorCode TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime)
1421: {
1424: if (maxsteps) {
1426: *maxsteps = ts->max_steps;
1427: }
1428: if (maxtime) {
1430: *maxtime = ts->max_time;
1431: }
1432: return(0);
1433: }
1437: /*@
1438: TSSetDuration - Sets the maximum number of timesteps to use and
1439: maximum time for iteration.
1441: Logically Collective on TS
1443: Input Parameters:
1444: + ts - the TS context obtained from TSCreate()
1445: . maxsteps - maximum number of iterations to use
1446: - maxtime - final time to iterate to
1448: Options Database Keys:
1449: . -ts_max_steps <maxsteps> - Sets maxsteps
1450: . -ts_final_time <maxtime> - Sets maxtime
1452: Notes:
1453: The default maximum number of iterations is 5000. Default time is 5.0
1455: Level: intermediate
1457: .keywords: TS, timestep, set, maximum, iterations
1459: .seealso: TSSetExactFinalTime()
1460: @*/
1461: PetscErrorCode TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime)
1462: {
1467: if (maxsteps >= 0) ts->max_steps = maxsteps;
1468: if (maxtime != PETSC_DEFAULT) ts->max_time = maxtime;
1469: return(0);
1470: }
1474: /*@
1475: TSSetSolution - Sets the initial solution vector
1476: for use by the TS routines.
1478: Logically Collective on TS and Vec
1480: Input Parameters:
1481: + ts - the TS context obtained from TSCreate()
1482: - x - the solution vector
1484: Level: beginner
1486: .keywords: TS, timestep, set, solution, initial conditions
1487: @*/
1488: PetscErrorCode TSSetSolution(TS ts,Vec x)
1489: {
1491: DM dm;
1496: PetscObjectReference((PetscObject)x);
1497: VecDestroy(&ts->vec_sol);
1498: ts->vec_sol = x;
1499: TSGetDM(ts,&dm);
1500: DMShellSetGlobalVector(dm,x);
1501: return(0);
1502: }
1506: /*@C
1507: TSSetPreStep - Sets the general-purpose function
1508: called once at the beginning of each time step.
1510: Logically Collective on TS
1512: Input Parameters:
1513: + ts - The TS context obtained from TSCreate()
1514: - func - The function
1516: Calling sequence of func:
1517: . func (TS ts);
1519: Level: intermediate
1521: .keywords: TS, timestep
1522: @*/
1523: PetscErrorCode TSSetPreStep(TS ts, PetscErrorCode (*func)(TS))
1524: {
1527: ts->ops->prestep = func;
1528: return(0);
1529: }
1533: /*@
1534: TSPreStep - Runs the user-defined pre-step function.
1536: Collective on TS
1538: Input Parameters:
1539: . ts - The TS context obtained from TSCreate()
1541: Notes:
1542: TSPreStep() is typically used within time stepping implementations,
1543: so most users would not generally call this routine themselves.
1545: Level: developer
1547: .keywords: TS, timestep
1548: @*/
1549: PetscErrorCode TSPreStep(TS ts)
1550: {
1555: if (ts->ops->prestep) {
1556: PetscStackPush("TS PreStep function");
1557: (*ts->ops->prestep)(ts);
1558: PetscStackPop;
1559: }
1560: return(0);
1561: }
1565: /*@C
1566: TSSetPostStep - Sets the general-purpose function
1567: called once at the end of each time step.
1569: Logically Collective on TS
1571: Input Parameters:
1572: + ts - The TS context obtained from TSCreate()
1573: - func - The function
1575: Calling sequence of func:
1576: . func (TS ts);
1578: Level: intermediate
1580: .keywords: TS, timestep
1581: @*/
1582: PetscErrorCode TSSetPostStep(TS ts, PetscErrorCode (*func)(TS))
1583: {
1586: ts->ops->poststep = func;
1587: return(0);
1588: }
1592: /*@
1593: TSPostStep - Runs the user-defined post-step function.
1595: Collective on TS
1597: Input Parameters:
1598: . ts - The TS context obtained from TSCreate()
1600: Notes:
1601: TSPostStep() is typically used within time stepping implementations,
1602: so most users would not generally call this routine themselves.
1604: Level: developer
1606: .keywords: TS, timestep
1607: @*/
1608: PetscErrorCode TSPostStep(TS ts)
1609: {
1614: if (ts->ops->poststep) {
1615: PetscStackPush("TS PostStep function");
1616: (*ts->ops->poststep)(ts);
1617: PetscStackPop;
1618: }
1619: return(0);
1620: }
1622: /* ------------ Routines to set performance monitoring options ----------- */
1626: /*@C
1627: TSMonitorSet - Sets an ADDITIONAL function that is to be used at every
1628: timestep to display the iteration's progress.
1630: Logically Collective on TS
1632: Input Parameters:
1633: + ts - the TS context obtained from TSCreate()
1634: . monitor - monitoring routine
1635: . mctx - [optional] user-defined context for private data for the
1636: monitor routine (use PETSC_NULL if no context is desired)
1637: - monitordestroy - [optional] routine that frees monitor context
1638: (may be PETSC_NULL)
1640: Calling sequence of monitor:
1641: $ int monitor(TS ts,PetscInt steps,PetscReal time,Vec x,void *mctx)
1643: + ts - the TS context
1644: . steps - iteration number
1645: . time - current time
1646: . x - current iterate
1647: - mctx - [optional] monitoring context
1649: Notes:
1650: This routine adds an additional monitor to the list of monitors that
1651: already has been loaded.
1653: Fortran notes: Only a single monitor function can be set for each TS object
1655: Level: intermediate
1657: .keywords: TS, timestep, set, monitor
1659: .seealso: TSMonitorDefault(), TSMonitorCancel()
1660: @*/
1661: PetscErrorCode TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void**))
1662: {
1665: if (ts->numbermonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1666: ts->monitor[ts->numbermonitors] = monitor;
1667: ts->mdestroy[ts->numbermonitors] = mdestroy;
1668: ts->monitorcontext[ts->numbermonitors++] = (void*)mctx;
1669: return(0);
1670: }
1674: /*@C
1675: TSMonitorCancel - Clears all the monitors that have been set on a time-step object.
1677: Logically Collective on TS
1679: Input Parameters:
1680: . ts - the TS context obtained from TSCreate()
1682: Notes:
1683: There is no way to remove a single, specific monitor.
1685: Level: intermediate
1687: .keywords: TS, timestep, set, monitor
1689: .seealso: TSMonitorDefault(), TSMonitorSet()
1690: @*/
1691: PetscErrorCode TSMonitorCancel(TS ts)
1692: {
1694: PetscInt i;
1698: for (i=0; i<ts->numbermonitors; i++) {
1699: if (ts->mdestroy[i]) {
1700: (*ts->mdestroy[i])(&ts->monitorcontext[i]);
1701: }
1702: }
1703: ts->numbermonitors = 0;
1704: return(0);
1705: }
1709: /*@
1710: TSMonitorDefault - Sets the Default monitor
1712: Level: intermediate
1714: .keywords: TS, set, monitor
1716: .seealso: TSMonitorDefault(), TSMonitorSet()
1717: @*/
1718: PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy)
1719: {
1721: PetscViewer viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ts)->comm);
1724: PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);
1725: PetscViewerASCIIPrintf(viewer,"%D TS dt %g time %g\n",step,(double)ts->time_step,(double)ptime);
1726: PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);
1727: return(0);
1728: }
1732: /*@
1733: TSSetRetainStages - Request that all stages in the upcoming step be stored so that interpolation will be available.
1735: Logically Collective on TS
1737: Input Argument:
1738: . ts - time stepping context
1740: Output Argument:
1741: . flg - PETSC_TRUE or PETSC_FALSE
1743: Level: intermediate
1745: .keywords: TS, set
1747: .seealso: TSInterpolate(), TSSetPostStep()
1748: @*/
1749: PetscErrorCode TSSetRetainStages(TS ts,PetscBool flg)
1750: {
1754: ts->retain_stages = flg;
1755: return(0);
1756: }
1760: /*@
1761: TSInterpolate - Interpolate the solution computed during the previous step to an arbitrary location in the interval
1763: Collective on TS
1765: Input Argument:
1766: + ts - time stepping context
1767: - t - time to interpolate to
1769: Output Argument:
1770: . X - state at given time
1772: Notes:
1773: The user should call TSSetRetainStages() before taking a step in which interpolation will be requested.
1775: Level: intermediate
1777: Developer Notes:
1778: TSInterpolate() and the storing of previous steps/stages should be generalized to support delay differential equations and continuous adjoints.
1780: .keywords: TS, set
1782: .seealso: TSSetRetainStages(), TSSetPostStep()
1783: @*/
1784: PetscErrorCode TSInterpolate(TS ts,PetscReal t,Vec X)
1785: {
1790: if (t < ts->ptime - ts->time_step_prev || t > ts->ptime) SETERRQ3(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Requested time %G not in last time steps [%G,%G]",t,ts->ptime-ts->time_step_prev,ts->ptime);
1791: if (!ts->ops->interpolate) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"%s does not provide interpolation",((PetscObject)ts)->type_name);
1792: (*ts->ops->interpolate)(ts,t,X);
1793: return(0);
1794: }
1798: /*@
1799: TSStep - Steps one time step
1801: Collective on TS
1803: Input Parameter:
1804: . ts - the TS context obtained from TSCreate()
1806: Level: intermediate
1808: .keywords: TS, timestep, solve
1810: .seealso: TSCreate(), TSSetUp(), TSDestroy(), TSSolve()
1811: @*/
1812: PetscErrorCode TSStep(TS ts)
1813: {
1814: PetscReal ptime_prev;
1819: TSSetUp(ts);
1821: ts->reason = TS_CONVERGED_ITERATING;
1823: ptime_prev = ts->ptime;
1824: PetscLogEventBegin(TS_Step,ts,0,0,0);
1825: (*ts->ops->step)(ts);
1826: PetscLogEventEnd(TS_Step,ts,0,0,0);
1827: ts->time_step_prev = ts->ptime - ptime_prev;
1829: if (ts->reason < 0) {
1830: if (ts->errorifstepfailed) {
1831: if (ts->reason == TS_DIVERGED_NONLINEAR_SOLVE) {
1832: SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s, increase -ts_max_snes_failures or make negative to attempt recovery",TSConvergedReasons[ts->reason]);
1833: } else SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s",TSConvergedReasons[ts->reason]);
1834: }
1835: } else if (!ts->reason) {
1836: if (ts->steps >= ts->max_steps)
1837: ts->reason = TS_CONVERGED_ITS;
1838: else if (ts->ptime >= ts->max_time)
1839: ts->reason = TS_CONVERGED_TIME;
1840: }
1842: return(0);
1843: }
1847: /*@
1848: TSEvaluateStep - Evaluate the solution at the end of a time step with a given order of accuracy.
1850: Collective on TS
1852: Input Arguments:
1853: + ts - time stepping context
1854: . order - desired order of accuracy
1855: - done - whether the step was evaluated at this order (pass PETSC_NULL to generate an error if not available)
1857: Output Arguments:
1858: . X - state at the end of the current step
1860: Level: advanced
1862: Notes:
1863: This function cannot be called until all stages have been evaluated.
1864: It is normally called by adaptive controllers before a step has been accepted and may also be called by the user after TSStep() has returned.
1866: .seealso: TSStep(), TSAdapt
1867: @*/
1868: PetscErrorCode TSEvaluateStep(TS ts,PetscInt order,Vec X,PetscBool *done)
1869: {
1876: if (!ts->ops->evaluatestep) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSEvaluateStep not implemented for type '%s'",((PetscObject)ts)->type_name);
1877: (*ts->ops->evaluatestep)(ts,order,X,done);
1878: return(0);
1879: }
1883: /*@
1884: TSSolve - Steps the requested number of timesteps.
1886: Collective on TS
1888: Input Parameter:
1889: + ts - the TS context obtained from TSCreate()
1890: - x - the solution vector
1892: Output Parameter:
1893: . ftime - time of the state vector x upon completion
1895: Level: beginner
1897: Notes:
1898: The final time returned by this function may be different from the time of the internally
1899: held state accessible by TSGetSolution() and TSGetTime() because the method may have
1900: stepped over the final time.
1902: .keywords: TS, timestep, solve
1904: .seealso: TSCreate(), TSSetSolution(), TSStep()
1905: @*/
1906: PetscErrorCode TSSolve(TS ts,Vec x,PetscReal *ftime)
1907: {
1908: PetscBool flg;
1909: char filename[PETSC_MAX_PATH_LEN];
1910: PetscViewer viewer;
1916: if (ts->exact_final_time) { /* Need ts->vec_sol to be distinct so it is not overwritten when we interpolate at the end */
1917: if (!ts->vec_sol || x == ts->vec_sol) {
1918: Vec y;
1919: VecDuplicate(x,&y);
1920: TSSetSolution(ts,y);
1921: VecDestroy(&y); /* grant ownership */
1922: }
1923: VecCopy(x,ts->vec_sol);
1924: } else {
1925: TSSetSolution(ts,x);
1926: }
1927: TSSetUp(ts);
1928: /* reset time step and iteration counters */
1929: ts->steps = 0;
1930: ts->ksp_its = 0;
1931: ts->snes_its = 0;
1932: ts->num_snes_failures = 0;
1933: ts->reject = 0;
1934: ts->reason = TS_CONVERGED_ITERATING;
1936: if (ts->ops->solve) { /* This private interface is transitional and should be removed when all implementations are updated. */
1937: (*ts->ops->solve)(ts);
1938: VecCopy(ts->vec_sol,x);
1939: if (ftime) *ftime = ts->ptime;
1940: } else {
1941: /* steps the requested number of timesteps. */
1942: TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);
1943: if (ts->steps >= ts->max_steps)
1944: ts->reason = TS_CONVERGED_ITS;
1945: else if (ts->ptime >= ts->max_time)
1946: ts->reason = TS_CONVERGED_TIME;
1947: while (!ts->reason) {
1948: TSPreStep(ts);
1949: TSStep(ts);
1950: TSPostStep(ts);
1951: TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);
1952: }
1953: if (ts->exact_final_time && ts->ptime > ts->max_time) {
1954: TSInterpolate(ts,ts->max_time,x);
1955: if (ftime) *ftime = ts->max_time;
1956: } else {
1957: VecCopy(ts->vec_sol,x);
1958: if (ftime) *ftime = ts->ptime;
1959: }
1960: }
1961: PetscOptionsGetString(((PetscObject)ts)->prefix,"-ts_view",filename,PETSC_MAX_PATH_LEN,&flg);
1962: if (flg && !PetscPreLoadingOn) {
1963: PetscViewerASCIIOpen(((PetscObject)ts)->comm,filename,&viewer);
1964: TSView(ts,viewer);
1965: PetscViewerDestroy(&viewer);
1966: }
1967: return(0);
1968: }
1972: /*@
1973: TSMonitor - Runs all user-provided monitor routines set using TSMonitorSet()
1975: Collective on TS
1977: Input Parameters:
1978: + ts - time stepping context obtained from TSCreate()
1979: . step - step number that has just completed
1980: . ptime - model time of the state
1981: - x - state at the current model time
1983: Notes:
1984: TSMonitor() is typically used within the time stepping implementations.
1985: Users might call this function when using the TSStep() interface instead of TSSolve().
1987: Level: advanced
1989: .keywords: TS, timestep
1990: @*/
1991: PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x)
1992: {
1994: PetscInt i,n = ts->numbermonitors;
1997: for (i=0; i<n; i++) {
1998: (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);
1999: }
2000: return(0);
2001: }
2003: /* ------------------------------------------------------------------------*/
2007: /*@C
2008: TSMonitorLGCreate - Creates a line graph context for use with
2009: TS to monitor convergence of preconditioned residual norms.
2011: Collective on TS
2013: Input Parameters:
2014: + host - the X display to open, or null for the local machine
2015: . label - the title to put in the title bar
2016: . x, y - the screen coordinates of the upper left coordinate of the window
2017: - m, n - the screen width and height in pixels
2019: Output Parameter:
2020: . draw - the drawing context
2022: Options Database Key:
2023: . -ts_monitor_draw - automatically sets line graph monitor
2025: Notes:
2026: Use TSMonitorLGDestroy() to destroy this line graph, not PetscDrawLGDestroy().
2028: Level: intermediate
2030: .keywords: TS, monitor, line graph, residual, seealso
2032: .seealso: TSMonitorLGDestroy(), TSMonitorSet()
2034: @*/
2035: PetscErrorCode TSMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
2036: {
2037: PetscDraw win;
2041: PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);
2042: PetscDrawSetType(win,PETSC_DRAW_X);
2043: PetscDrawLGCreate(win,1,draw);
2044: PetscDrawLGIndicateDataPoints(*draw);
2046: PetscLogObjectParent(*draw,win);
2047: return(0);
2048: }
2052: PetscErrorCode TSMonitorLG(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx)
2053: {
2054: PetscDrawLG lg = (PetscDrawLG) monctx;
2055: PetscReal x,y = ptime;
2059: if (!monctx) {
2060: MPI_Comm comm;
2061: PetscViewer viewer;
2063: PetscObjectGetComm((PetscObject)ts,&comm);
2064: viewer = PETSC_VIEWER_DRAW_(comm);
2065: PetscViewerDrawGetDrawLG(viewer,0,&lg);
2066: }
2068: if (!n) {PetscDrawLGReset(lg);}
2069: x = (PetscReal)n;
2070: PetscDrawLGAddPoint(lg,&x,&y);
2071: if (n < 20 || (n % 5)) {
2072: PetscDrawLGDraw(lg);
2073: }
2074: return(0);
2075: }
2079: /*@C
2080: TSMonitorLGDestroy - Destroys a line graph context that was created
2081: with TSMonitorLGCreate().
2083: Collective on PetscDrawLG
2085: Input Parameter:
2086: . draw - the drawing context
2088: Level: intermediate
2090: .keywords: TS, monitor, line graph, destroy
2092: .seealso: TSMonitorLGCreate(), TSMonitorSet(), TSMonitorLG();
2093: @*/
2094: PetscErrorCode TSMonitorLGDestroy(PetscDrawLG *drawlg)
2095: {
2096: PetscDraw draw;
2100: PetscDrawLGGetDraw(*drawlg,&draw);
2101: PetscDrawDestroy(&draw);
2102: PetscDrawLGDestroy(drawlg);
2103: return(0);
2104: }
2108: /*@
2109: TSGetTime - Gets the current time.
2111: Not Collective
2113: Input Parameter:
2114: . ts - the TS context obtained from TSCreate()
2116: Output Parameter:
2117: . t - the current time
2119: Level: beginner
2121: .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
2123: .keywords: TS, get, time
2124: @*/
2125: PetscErrorCode TSGetTime(TS ts,PetscReal* t)
2126: {
2130: *t = ts->ptime;
2131: return(0);
2132: }
2136: /*@
2137: TSSetTime - Allows one to reset the time.
2139: Logically Collective on TS
2141: Input Parameters:
2142: + ts - the TS context obtained from TSCreate()
2143: - time - the time
2145: Level: intermediate
2147: .seealso: TSGetTime(), TSSetDuration()
2149: .keywords: TS, set, time
2150: @*/
2151: PetscErrorCode TSSetTime(TS ts, PetscReal t)
2152: {
2156: ts->ptime = t;
2157: return(0);
2158: }
2162: /*@C
2163: TSSetOptionsPrefix - Sets the prefix used for searching for all
2164: TS options in the database.
2166: Logically Collective on TS
2168: Input Parameter:
2169: + ts - The TS context
2170: - prefix - The prefix to prepend to all option names
2172: Notes:
2173: A hyphen (-) must NOT be given at the beginning of the prefix name.
2174: The first character of all runtime options is AUTOMATICALLY the
2175: hyphen.
2177: Level: advanced
2179: .keywords: TS, set, options, prefix, database
2181: .seealso: TSSetFromOptions()
2183: @*/
2184: PetscErrorCode TSSetOptionsPrefix(TS ts,const char prefix[])
2185: {
2187: SNES snes;
2191: PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);
2192: TSGetSNES(ts,&snes);
2193: SNESSetOptionsPrefix(snes,prefix);
2194: return(0);
2195: }
2200: /*@C
2201: TSAppendOptionsPrefix - Appends to the prefix used for searching for all
2202: TS options in the database.
2204: Logically Collective on TS
2206: Input Parameter:
2207: + ts - The TS context
2208: - prefix - The prefix to prepend to all option names
2210: Notes:
2211: A hyphen (-) must NOT be given at the beginning of the prefix name.
2212: The first character of all runtime options is AUTOMATICALLY the
2213: hyphen.
2215: Level: advanced
2217: .keywords: TS, append, options, prefix, database
2219: .seealso: TSGetOptionsPrefix()
2221: @*/
2222: PetscErrorCode TSAppendOptionsPrefix(TS ts,const char prefix[])
2223: {
2225: SNES snes;
2229: PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);
2230: TSGetSNES(ts,&snes);
2231: SNESAppendOptionsPrefix(snes,prefix);
2232: return(0);
2233: }
2237: /*@C
2238: TSGetOptionsPrefix - Sets the prefix used for searching for all
2239: TS options in the database.
2241: Not Collective
2243: Input Parameter:
2244: . ts - The TS context
2246: Output Parameter:
2247: . prefix - A pointer to the prefix string used
2249: Notes: On the fortran side, the user should pass in a string 'prifix' of
2250: sufficient length to hold the prefix.
2252: Level: intermediate
2254: .keywords: TS, get, options, prefix, database
2256: .seealso: TSAppendOptionsPrefix()
2257: @*/
2258: PetscErrorCode TSGetOptionsPrefix(TS ts,const char *prefix[])
2259: {
2265: PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);
2266: return(0);
2267: }
2271: /*@C
2272: TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
2274: Not Collective, but parallel objects are returned if TS is parallel
2276: Input Parameter:
2277: . ts - The TS context obtained from TSCreate()
2279: Output Parameters:
2280: + J - The Jacobian J of F, where U_t = F(U,t)
2281: . M - The preconditioner matrix, usually the same as J
2282: . func - Function to compute the Jacobian of the RHS
2283: - ctx - User-defined context for Jacobian evaluation routine
2285: Notes: You can pass in PETSC_NULL for any return argument you do not need.
2287: Level: intermediate
2289: .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()
2291: .keywords: TS, timestep, get, matrix, Jacobian
2292: @*/
2293: PetscErrorCode TSGetRHSJacobian(TS ts,Mat *J,Mat *M,TSRHSJacobian *func,void **ctx)
2294: {
2296: SNES snes;
2299: TSGetSNES(ts,&snes);
2300: SNESGetJacobian(snes,J,M,PETSC_NULL,PETSC_NULL);
2301: if (func) *func = ts->userops->rhsjacobian;
2302: if (ctx) *ctx = ts->jacP;
2303: return(0);
2304: }
2308: /*@C
2309: TSGetIJacobian - Returns the implicit Jacobian at the present timestep.
2311: Not Collective, but parallel objects are returned if TS is parallel
2313: Input Parameter:
2314: . ts - The TS context obtained from TSCreate()
2316: Output Parameters:
2317: + A - The Jacobian of F(t,U,U_t)
2318: . B - The preconditioner matrix, often the same as A
2319: . f - The function to compute the matrices
2320: - ctx - User-defined context for Jacobian evaluation routine
2322: Notes: You can pass in PETSC_NULL for any return argument you do not need.
2324: Level: advanced
2326: .seealso: TSGetTimeStep(), TSGetRHSJacobian(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()
2328: .keywords: TS, timestep, get, matrix, Jacobian
2329: @*/
2330: PetscErrorCode TSGetIJacobian(TS ts,Mat *A,Mat *B,TSIJacobian *f,void **ctx)
2331: {
2333: SNES snes;
2336: TSGetSNES(ts,&snes);
2337: SNESGetJacobian(snes,A,B,PETSC_NULL,PETSC_NULL);
2338: if (f) *f = ts->userops->ijacobian;
2339: if (ctx) *ctx = ts->jacP;
2340: return(0);
2341: }
2343: typedef struct {
2344: PetscViewer viewer;
2345: Vec initialsolution;
2346: PetscBool showinitial;
2347: } TSMonitorSolutionCtx;
2351: /*@C
2352: TSMonitorSolution - Monitors progress of the TS solvers by calling
2353: VecView() for the solution at each timestep
2355: Collective on TS
2357: Input Parameters:
2358: + ts - the TS context
2359: . step - current time-step
2360: . ptime - current time
2361: - dummy - either a viewer or PETSC_NULL
2363: Level: intermediate
2365: .keywords: TS, vector, monitor, view
2367: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
2368: @*/
2369: PetscErrorCode TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy)
2370: {
2371: PetscErrorCode ierr;
2372: TSMonitorSolutionCtx *ictx = (TSMonitorSolutionCtx*)dummy;
2375: if (!step && ictx->showinitial) {
2376: if (!ictx->initialsolution) {
2377: VecDuplicate(x,&ictx->initialsolution);
2378: }
2379: VecCopy(x,ictx->initialsolution);
2380: }
2381: if (ictx->showinitial) {
2382: PetscReal pause;
2383: PetscViewerDrawGetPause(ictx->viewer,&pause);
2384: PetscViewerDrawSetPause(ictx->viewer,0.0);
2385: VecView(ictx->initialsolution,ictx->viewer);
2386: PetscViewerDrawSetPause(ictx->viewer,pause);
2387: PetscViewerDrawSetHold(ictx->viewer,PETSC_TRUE);
2388: }
2389: VecView(x,ictx->viewer);
2390: if (ictx->showinitial) {
2391: PetscViewerDrawSetHold(ictx->viewer,PETSC_FALSE);
2392: }
2393: return(0);
2394: }
2399: /*@C
2400: TSMonitorSolutionDestroy - Destroys the monitor context for TSMonitorSolution
2402: Collective on TS
2404: Input Parameters:
2405: . ctx - the monitor context
2407: Level: intermediate
2409: .keywords: TS, vector, monitor, view
2411: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorSolution()
2412: @*/
2413: PetscErrorCode TSMonitorSolutionDestroy(void **ctx)
2414: {
2415: PetscErrorCode ierr;
2416: TSMonitorSolutionCtx *ictx = *(TSMonitorSolutionCtx**)ctx;
2417:
2419: PetscViewerDestroy(&ictx->viewer);
2420: VecDestroy(&ictx->initialsolution);
2421: PetscFree(ictx);
2422: return(0);
2423: }
2427: /*@C
2428: TSMonitorSolutionCreate - Creates the monitor context for TSMonitorSolution
2430: Collective on TS
2432: Input Parameter:
2433: . ts - time-step context
2435: Output Patameter:
2436: . ctx - the monitor context
2438: Level: intermediate
2440: .keywords: TS, vector, monitor, view
2442: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorSolution()
2443: @*/
2444: PetscErrorCode TSMonitorSolutionCreate(TS ts,PetscViewer viewer,void **ctx)
2445: {
2446: PetscErrorCode ierr;
2447: TSMonitorSolutionCtx *ictx;
2448:
2450: PetscNew(TSMonitorSolutionCtx,&ictx);
2451: *ctx = (void*)ictx;
2452: if (!viewer) {
2453: viewer = PETSC_VIEWER_DRAW_(((PetscObject)ts)->comm);
2454: }
2455: PetscObjectReference((PetscObject)viewer);
2456: ictx->viewer = viewer;
2457: ictx->showinitial = PETSC_FALSE;
2458: PetscOptionsGetBool(((PetscObject)ts)->prefix,"-ts_monitor_solution_initial",&ictx->showinitial,PETSC_NULL);
2459: return(0);
2460: }
2464: /*@
2465: TSSetDM - Sets the DM that may be used by some preconditioners
2467: Logically Collective on TS and DM
2469: Input Parameters:
2470: + ts - the preconditioner context
2471: - dm - the dm
2473: Level: intermediate
2476: .seealso: TSGetDM(), SNESSetDM(), SNESGetDM()
2477: @*/
2478: PetscErrorCode TSSetDM(TS ts,DM dm)
2479: {
2481: SNES snes;
2485: PetscObjectReference((PetscObject)dm);
2486: DMDestroy(&ts->dm);
2487: ts->dm = dm;
2488: TSGetSNES(ts,&snes);
2489: SNESSetDM(snes,dm);
2490: return(0);
2491: }
2495: /*@
2496: TSGetDM - Gets the DM that may be used by some preconditioners
2498: Not Collective
2500: Input Parameter:
2501: . ts - the preconditioner context
2503: Output Parameter:
2504: . dm - the dm
2506: Level: intermediate
2509: .seealso: TSSetDM(), SNESSetDM(), SNESGetDM()
2510: @*/
2511: PetscErrorCode TSGetDM(TS ts,DM *dm)
2512: {
2517: if (!ts->dm) {
2518: DMShellCreate(((PetscObject)ts)->comm,&ts->dm);
2519: if (ts->snes) {SNESSetDM(ts->snes,ts->dm);}
2520: }
2521: *dm = ts->dm;
2522: return(0);
2523: }
2527: /*@
2528: SNESTSFormFunction - Function to evaluate nonlinear residual
2530: Logically Collective on SNES
2532: Input Parameter:
2533: + snes - nonlinear solver
2534: . X - the current state at which to evaluate the residual
2535: - ctx - user context, must be a TS
2537: Output Parameter:
2538: . F - the nonlinear residual
2540: Notes:
2541: This function is not normally called by users and is automatically registered with the SNES used by TS.
2542: It is most frequently passed to MatFDColoringSetFunction().
2544: Level: advanced
2546: .seealso: SNESSetFunction(), MatFDColoringSetFunction()
2547: @*/
2548: PetscErrorCode SNESTSFormFunction(SNES snes,Vec X,Vec F,void *ctx)
2549: {
2550: TS ts = (TS)ctx;
2558: (ts->ops->snesfunction)(snes,X,F,ts);
2559: return(0);
2560: }
2564: /*@
2565: SNESTSFormJacobian - Function to evaluate the Jacobian
2567: Collective on SNES
2569: Input Parameter:
2570: + snes - nonlinear solver
2571: . X - the current state at which to evaluate the residual
2572: - ctx - user context, must be a TS
2574: Output Parameter:
2575: + A - the Jacobian
2576: . B - the preconditioning matrix (may be the same as A)
2577: - flag - indicates any structure change in the matrix
2579: Notes:
2580: This function is not normally called by users and is automatically registered with the SNES used by TS.
2582: Level: developer
2584: .seealso: SNESSetJacobian()
2585: @*/
2586: PetscErrorCode SNESTSFormJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flag,void *ctx)
2587: {
2588: TS ts = (TS)ctx;
2600: (ts->ops->snesjacobian)(snes,X,A,B,flag,ts);
2601: return(0);
2602: }
2606: /*@C
2607: TSComputeRHSFunctionLinear - Evaluate the right hand side via the user-provided Jacobian, for linear problems only
2609: Collective on TS
2611: Input Arguments:
2612: + ts - time stepping context
2613: . t - time at which to evaluate
2614: . X - state at which to evaluate
2615: - ctx - context
2617: Output Arguments:
2618: . F - right hand side
2620: Level: intermediate
2622: Notes:
2623: This function is intended to be passed to TSSetRHSFunction() to evaluate the right hand side for linear problems.
2624: The matrix (and optionally the evaluation context) should be passed to TSSetRHSJacobian().
2626: .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSJacobianConstant()
2627: @*/
2628: PetscErrorCode TSComputeRHSFunctionLinear(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
2629: {
2631: Mat Arhs,Brhs;
2632: MatStructure flg2;
2635: TSGetRHSMats_Private(ts,&Arhs,&Brhs);
2636: TSComputeRHSJacobian(ts,t,X,&Arhs,&Brhs,&flg2);
2637: MatMult(Arhs,X,F);
2638: return(0);
2639: }
2643: /*@C
2644: TSComputeRHSJacobianConstant - Reuses a Jacobian that is time-independent.
2646: Collective on TS
2648: Input Arguments:
2649: + ts - time stepping context
2650: . t - time at which to evaluate
2651: . X - state at which to evaluate
2652: - ctx - context
2654: Output Arguments:
2655: + A - pointer to operator
2656: . B - pointer to preconditioning matrix
2657: - flg - matrix structure flag
2659: Level: intermediate
2661: Notes:
2662: This function is intended to be passed to TSSetRHSJacobian() to evaluate the Jacobian for linear time-independent problems.
2664: .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSFunctionLinear()
2665: @*/
2666: PetscErrorCode TSComputeRHSJacobianConstant(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg,void *ctx)
2667: {
2670: *flg = SAME_PRECONDITIONER;
2671: return(0);
2672: }
2676: /*@C
2677: TSComputeIFunctionLinear - Evaluate the left hand side via the user-provided Jacobian, for linear problems only
2679: Collective on TS
2681: Input Arguments:
2682: + ts - time stepping context
2683: . t - time at which to evaluate
2684: . X - state at which to evaluate
2685: . Xdot - time derivative of state vector
2686: - ctx - context
2688: Output Arguments:
2689: . F - left hand side
2691: Level: intermediate
2693: Notes:
2694: The assumption here is that the left hand side is of the form A*Xdot (and not A*Xdot + B*X). For other cases, the
2695: user is required to write their own TSComputeIFunction.
2696: This function is intended to be passed to TSSetIFunction() to evaluate the left hand side for linear problems.
2697: The matrix (and optionally the evaluation context) should be passed to TSSetIJacobian().
2699: .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIJacobianConstant()
2700: @*/
2701: PetscErrorCode TSComputeIFunctionLinear(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
2702: {
2704: Mat A,B;
2705: MatStructure flg2;
2708: TSGetIJacobian(ts,&A,&B,PETSC_NULL,PETSC_NULL);
2709: TSComputeIJacobian(ts,t,X,Xdot,1.0,&A,&B,&flg2,PETSC_TRUE);
2710: MatMult(A,Xdot,F);
2711: return(0);
2712: }
2716: /*@C
2717: TSComputeIJacobianConstant - Reuses a Jacobian that is time-independent.
2719: Collective on TS
2721: Input Arguments:
2722: + ts - time stepping context
2723: . t - time at which to evaluate
2724: . X - state at which to evaluate
2725: . Xdot - time derivative of state vector
2726: . shift - shift to apply
2727: - ctx - context
2729: Output Arguments:
2730: + A - pointer to operator
2731: . B - pointer to preconditioning matrix
2732: - flg - matrix structure flag
2734: Level: intermediate
2736: Notes:
2737: This function is intended to be passed to TSSetIJacobian() to evaluate the Jacobian for linear time-independent problems.
2739: .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIFunctionLinear()
2740: @*/
2741: PetscErrorCode TSComputeIJacobianConstant(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg,void *ctx)
2742: {
2745: *flg = SAME_PRECONDITIONER;
2746: return(0);
2747: }
2752: /*@
2753: TSGetConvergedReason - Gets the reason the TS iteration was stopped.
2755: Not Collective
2757: Input Parameter:
2758: . ts - the TS context
2760: Output Parameter:
2761: . reason - negative value indicates diverged, positive value converged, see TSConvergedReason or the
2762: manual pages for the individual convergence tests for complete lists
2764: Level: intermediate
2766: Notes:
2767: Can only be called after the call to TSSolve() is complete.
2769: .keywords: TS, nonlinear, set, convergence, test
2771: .seealso: TSSetConvergenceTest(), TSConvergedReason
2772: @*/
2773: PetscErrorCode TSGetConvergedReason(TS ts,TSConvergedReason *reason)
2774: {
2778: *reason = ts->reason;
2779: return(0);
2780: }
2784: /*@
2785: TSGetSNESIterations - Gets the total number of nonlinear iterations
2786: used by the time integrator.
2788: Not Collective
2790: Input Parameter:
2791: . ts - TS context
2793: Output Parameter:
2794: . nits - number of nonlinear iterations
2796: Notes:
2797: This counter is reset to zero for each successive call to TSSolve().
2799: Level: intermediate
2801: .keywords: TS, get, number, nonlinear, iterations
2803: .seealso: TSGetKSPIterations()
2804: @*/
2805: PetscErrorCode TSGetSNESIterations(TS ts,PetscInt *nits)
2806: {
2810: *nits = ts->snes_its;
2811: return(0);
2812: }
2816: /*@
2817: TSGetKSPIterations - Gets the total number of linear iterations
2818: used by the time integrator.
2820: Not Collective
2822: Input Parameter:
2823: . ts - TS context
2825: Output Parameter:
2826: . lits - number of linear iterations
2828: Notes:
2829: This counter is reset to zero for each successive call to TSSolve().
2831: Level: intermediate
2833: .keywords: TS, get, number, linear, iterations
2835: .seealso: TSGetSNESIterations(), SNESGetKSPIterations()
2836: @*/
2837: PetscErrorCode TSGetKSPIterations(TS ts,PetscInt *lits)
2838: {
2842: *lits = ts->ksp_its;
2843: return(0);
2844: }
2848: /*@
2849: TSGetStepRejections - Gets the total number of rejected steps.
2851: Not Collective
2853: Input Parameter:
2854: . ts - TS context
2856: Output Parameter:
2857: . rejects - number of steps rejected
2859: Notes:
2860: This counter is reset to zero for each successive call to TSSolve().
2862: Level: intermediate
2864: .keywords: TS, get, number
2866: .seealso: TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetSNESFailures(), TSSetMaxSNESFailures(), TSSetErrorIfStepFails()
2867: @*/
2868: PetscErrorCode TSGetStepRejections(TS ts,PetscInt *rejects)
2869: {
2873: *rejects = ts->reject;
2874: return(0);
2875: }
2879: /*@
2880: TSGetSNESFailures - Gets the total number of failed SNES solves
2882: Not Collective
2884: Input Parameter:
2885: . ts - TS context
2887: Output Parameter:
2888: . fails - number of failed nonlinear solves
2890: Notes:
2891: This counter is reset to zero for each successive call to TSSolve().
2893: Level: intermediate
2895: .keywords: TS, get, number
2897: .seealso: TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSSetMaxSNESFailures()
2898: @*/
2899: PetscErrorCode TSGetSNESFailures(TS ts,PetscInt *fails)
2900: {
2904: *fails = ts->num_snes_failures;
2905: return(0);
2906: }
2910: /*@
2911: TSSetMaxStepRejections - Sets the maximum number of step rejections before a step fails
2913: Not Collective
2915: Input Parameter:
2916: + ts - TS context
2917: - rejects - maximum number of rejected steps, pass -1 for unlimited
2919: Notes:
2920: The counter is reset to zero for each step
2922: Options Database Key:
2923: . -ts_max_reject - Maximum number of step rejections before a step fails
2925: Level: intermediate
2927: .keywords: TS, set, maximum, number
2929: .seealso: TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxSNESFailures(), TSGetStepRejections(), TSGetSNESFailures(), TSSetErrorIfStepFails(), TSGetConvergedReason()
2930: @*/
2931: PetscErrorCode TSSetMaxStepRejections(TS ts,PetscInt rejects)
2932: {
2935: ts->max_reject = rejects;
2936: return(0);
2937: }
2941: /*@
2942: TSSetMaxSNESFailures - Sets the maximum number of failed SNES solves
2944: Not Collective
2946: Input Parameter:
2947: + ts - TS context
2948: - fails - maximum number of failed nonlinear solves, pass -1 for unlimited
2950: Notes:
2951: The counter is reset to zero for each successive call to TSSolve().
2953: Options Database Key:
2954: . -ts_max_snes_failures - Maximum number of nonlinear solve failures
2956: Level: intermediate
2958: .keywords: TS, set, maximum, number
2960: .seealso: TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSGetSNESFailures(), SNESGetConvergedReason(), TSGetConvergedReason()
2961: @*/
2962: PetscErrorCode TSSetMaxSNESFailures(TS ts,PetscInt fails)
2963: {
2966: ts->max_snes_failures = fails;
2967: return(0);
2968: }
2972: /*@
2973: TSSetErrorIfStepFails - Error if no step succeeds
2975: Not Collective
2977: Input Parameter:
2978: + ts - TS context
2979: - err - PETSC_TRUE to error if no step succeeds, PETSC_FALSE to return without failure
2981: Options Database Key:
2982: . -ts_error_if_step_fails - Error if no step succeeds
2984: Level: intermediate
2986: .keywords: TS, set, error
2988: .seealso: TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSGetSNESFailures(), TSSetErrorIfStepFails(), TSGetConvergedReason()
2989: @*/
2990: PetscErrorCode TSSetErrorIfStepFails(TS ts,PetscBool err)
2991: {
2994: ts->errorifstepfailed = err;
2995: return(0);
2996: }
3000: /*@C
3001: TSMonitorSolutionBinary - Monitors progress of the TS solvers by VecView() for the solution at each timestep. Normally the viewer is a binary file
3003: Collective on TS
3005: Input Parameters:
3006: + ts - the TS context
3007: . step - current time-step
3008: . ptime - current time
3009: . x - current state
3010: - viewer - binary viewer
3012: Level: intermediate
3014: .keywords: TS, vector, monitor, view
3016: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
3017: @*/
3018: PetscErrorCode TSMonitorSolutionBinary(TS ts,PetscInt step,PetscReal ptime,Vec x,void *viewer)
3019: {
3020: PetscErrorCode ierr;
3021: PetscViewer v = (PetscViewer)viewer;
3024: VecView(x,v);
3025: return(0);
3026: }
3030: /*@C
3031: TSMonitorSolutionVTK - Monitors progress of the TS solvers by VecView() for the solution at each timestep.
3033: Collective on TS
3035: Input Parameters:
3036: + ts - the TS context
3037: . step - current time-step
3038: . ptime - current time
3039: . x - current state
3040: - filenametemplate - string containing a format specifier for the integer time step (e.g. %03D)
3042: Level: intermediate
3044: Notes:
3045: The VTK format does not allow writing multiple time steps in the same file, therefore a different file will be written for each time step.
3046: These are named according to the file name template.
3048: This function is normally passed as an argument to TSMonitorSet() along with TSMonitorSolutionVTKDestroy().
3050: .keywords: TS, vector, monitor, view
3052: .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
3053: @*/
3054: PetscErrorCode TSMonitorSolutionVTK(TS ts,PetscInt step,PetscReal ptime,Vec x,void *filenametemplate)
3055: {
3057: char filename[PETSC_MAX_PATH_LEN];
3058: PetscViewer viewer;
3061: PetscSNPrintf(filename,sizeof filename,(const char*)filenametemplate,step);
3062: PetscViewerVTKOpen(((PetscObject)ts)->comm,filename,FILE_MODE_WRITE,&viewer);
3063: VecView(x,viewer);
3064: PetscViewerDestroy(&viewer);
3065: return(0);
3066: }
3070: /*@C
3071: TSMonitorSolutionVTKDestroy - Destroy context for monitoring
3073: Collective on TS
3075: Input Parameters:
3076: . filenametemplate - string containing a format specifier for the integer time step (e.g. %03D)
3078: Level: intermediate
3080: Note:
3081: This function is normally passed to TSMonitorSet() along with TSMonitorSolutionVTK().
3083: .keywords: TS, vector, monitor, view
3085: .seealso: TSMonitorSet(), TSMonitorSolutionVTK()
3086: @*/
3087: PetscErrorCode TSMonitorSolutionVTKDestroy(void *filenametemplate)
3088: {
3092: PetscFree(*(char**)filenametemplate);
3093: return(0);
3094: }
3098: /*@
3099: TSGetAdapt - Get the adaptive controller context for the current method
3101: Collective on TS if controller has not been created yet
3103: Input Arguments:
3104: . ts - time stepping context
3106: Output Arguments:
3107: . adapt - adaptive controller
3109: Level: intermediate
3111: .seealso: TSAdapt, TSAdaptSetType(), TSAdaptChoose()
3112: @*/
3113: PetscErrorCode TSGetAdapt(TS ts,TSAdapt *adapt)
3114: {
3120: if (!ts->adapt) {
3121: TSAdaptCreate(((PetscObject)ts)->comm,&ts->adapt);
3122: PetscLogObjectParent(ts,ts->adapt);
3123: PetscObjectIncrementTabLevel((PetscObject)ts->adapt,(PetscObject)ts,1);
3124: }
3125: *adapt = ts->adapt;
3126: return(0);
3127: }
3131: /*@
3132: TSSetTolerances - Set tolerances for local truncation error when using adaptive controller
3134: Logically Collective
3136: Input Arguments:
3137: + ts - time integration context
3138: . atol - scalar absolute tolerances, PETSC_DECIDE to leave current value
3139: . vatol - vector of absolute tolerances or PETSC_NULL, used in preference to atol if present
3140: . rtol - scalar relative tolerances, PETSC_DECIDE to leave current value
3141: - vrtol - vector of relative tolerances or PETSC_NULL, used in preference to atol if present
3143: Level: beginner
3145: .seealso: TS, TSAdapt, TSVecNormWRMS(), TSGetTolerances()
3146: @*/
3147: PetscErrorCode TSSetTolerances(TS ts,PetscReal atol,Vec vatol,PetscReal rtol,Vec vrtol)
3148: {
3152: if (atol != PETSC_DECIDE && atol != PETSC_DEFAULT) ts->atol = atol;
3153: if (vatol) {
3154: PetscObjectReference((PetscObject)vatol);
3155: VecDestroy(&ts->vatol);
3156: ts->vatol = vatol;
3157: }
3158: if (rtol != PETSC_DECIDE && rtol != PETSC_DEFAULT) ts->rtol = rtol;
3159: if (vrtol) {
3160: PetscObjectReference((PetscObject)vrtol);
3161: VecDestroy(&ts->vrtol);
3162: ts->vrtol = vrtol;
3163: }
3164: return(0);
3165: }
3169: /*@
3170: TSGetTolerances - Get tolerances for local truncation error when using adaptive controller
3172: Logically Collective
3174: Input Arguments:
3175: . ts - time integration context
3177: Output Arguments:
3178: + atol - scalar absolute tolerances, PETSC_NULL to ignore
3179: . vatol - vector of absolute tolerances, PETSC_NULL to ignore
3180: . rtol - scalar relative tolerances, PETSC_NULL to ignore
3181: - vrtol - vector of relative tolerances, PETSC_NULL to ignore
3183: Level: beginner
3185: .seealso: TS, TSAdapt, TSVecNormWRMS(), TSSetTolerances()
3186: @*/
3187: PetscErrorCode TSGetTolerances(TS ts,PetscReal *atol,Vec *vatol,PetscReal *rtol,Vec *vrtol)
3188: {
3191: if (atol) *atol = ts->atol;
3192: if (vatol) *vatol = ts->vatol;
3193: if (rtol) *rtol = ts->rtol;
3194: if (vrtol) *vrtol = ts->vrtol;
3195: return(0);
3196: }
3200: /*@
3201: TSErrorNormWRMS - compute a weighted norm of the difference between a vector and the current state
3203: Collective on TS
3205: Input Arguments:
3206: + ts - time stepping context
3207: - Y - state vector to be compared to ts->vec_sol
3209: Output Arguments:
3210: . norm - weighted norm, a value of 1.0 is considered small
3212: Level: developer
3214: .seealso: TSSetTolerances()
3215: @*/
3216: PetscErrorCode TSErrorNormWRMS(TS ts,Vec Y,PetscReal *norm)
3217: {
3219: PetscInt i,n,N;
3220: const PetscScalar *x,*y;
3221: Vec X;
3222: PetscReal sum,gsum;
3228: X = ts->vec_sol;
3230: if (X == Y) SETERRQ(((PetscObject)X)->comm,PETSC_ERR_ARG_IDN,"Y cannot be the TS solution vector");
3232: VecGetSize(X,&N);
3233: VecGetLocalSize(X,&n);
3234: VecGetArrayRead(X,&x);
3235: VecGetArrayRead(Y,&y);
3236: sum = 0.;
3237: if (ts->vatol && ts->vrtol) {
3238: const PetscScalar *atol,*rtol;
3239: VecGetArrayRead(ts->vatol,&atol);
3240: VecGetArrayRead(ts->vrtol,&rtol);
3241: for (i=0; i<n; i++) {
3242: PetscReal tol = PetscRealPart(atol[i]) + PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i]));
3243: sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol);
3244: }
3245: VecRestoreArrayRead(ts->vatol,&atol);
3246: VecRestoreArrayRead(ts->vrtol,&rtol);
3247: } else if (ts->vatol) { /* vector atol, scalar rtol */
3248: const PetscScalar *atol;
3249: VecGetArrayRead(ts->vatol,&atol);
3250: for (i=0; i<n; i++) {
3251: PetscReal tol = PetscRealPart(atol[i]) + ts->rtol * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i]));
3252: sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol);
3253: }
3254: VecRestoreArrayRead(ts->vatol,&atol);
3255: } else if (ts->vrtol) { /* scalar atol, vector rtol */
3256: const PetscScalar *rtol;
3257: VecGetArrayRead(ts->vrtol,&rtol);
3258: for (i=0; i<n; i++) {
3259: PetscReal tol = ts->atol + PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i]));
3260: sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol);
3261: }
3262: VecRestoreArrayRead(ts->vrtol,&rtol);
3263: } else { /* scalar atol, scalar rtol */
3264: for (i=0; i<n; i++) {
3265: PetscReal tol = ts->atol + ts->rtol * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i]));
3266: sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol);
3267: }
3268: }
3269: VecRestoreArrayRead(X,&x);
3270: VecRestoreArrayRead(Y,&y);
3272: MPI_Allreduce(&sum,&gsum,1,MPIU_REAL,MPIU_SUM,((PetscObject)ts)->comm);
3273: *norm = PetscSqrtReal(gsum / N);
3274: if (PetscIsInfOrNanScalar(*norm)) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
3275: return(0);
3276: }
3280: /*@
3281: TSSetCFLTimeLocal - Set the local CFL constraint relative to forward Euler
3283: Logically Collective on TS
3285: Input Arguments:
3286: + ts - time stepping context
3287: - cfltime - maximum stable time step if using forward Euler (value can be different on each process)
3289: Note:
3290: After calling this function, the global CFL time can be obtained by calling TSGetCFLTime()
3292: Level: intermediate
3294: .seealso: TSGetCFLTime(), TSADAPTCFL
3295: @*/
3296: PetscErrorCode TSSetCFLTimeLocal(TS ts,PetscReal cfltime)
3297: {
3301: ts->cfltime_local = cfltime;
3302: ts->cfltime = -1.;
3303: return(0);
3304: }
3308: /*@
3309: TSGetCFLTime - Get the maximum stable time step according to CFL criteria applied to forward Euler
3311: Collective on TS
3313: Input Arguments:
3314: . ts - time stepping context
3316: Output Arguments:
3317: . cfltime - maximum stable time step for forward Euler
3319: Level: advanced
3321: .seealso: TSSetCFLTimeLocal()
3322: @*/
3323: PetscErrorCode TSGetCFLTime(TS ts,PetscReal *cfltime)
3324: {
3328: if (ts->cfltime < 0) {
3329: MPI_Allreduce(&ts->cfltime_local,&ts->cfltime,1,MPIU_REAL,MPIU_MIN,((PetscObject)ts)->comm);
3330: }
3331: *cfltime = ts->cfltime;
3332: return(0);
3333: }
3337: /*@
3338: TSVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu
3340: Input Parameters:
3341: . ts - the TS context.
3342: . xl - lower bound.
3343: . xu - upper bound.
3345: Notes:
3346: If this routine is not called then the lower and upper bounds are set to
3347: SNES_VI_NINF and SNES_VI_INF respectively during SNESSetUp().
3349: Level: advanced
3351: @*/
3352: PetscErrorCode TSVISetVariableBounds(TS ts, Vec xl, Vec xu)
3353: {
3355: SNES snes;
3358: TSGetSNES(ts,&snes);
3359: SNESVISetVariableBounds(snes,xl,xu);
3360: return(0);
3361: }
3363: #if defined(PETSC_HAVE_MATLAB_ENGINE)
3364: #include <mex.h>
3366: typedef struct {char *funcname; mxArray *ctx;} TSMatlabContext;
3370: /*
3371: TSComputeFunction_Matlab - Calls the function that has been set with
3372: TSSetFunctionMatlab().
3374: Collective on TS
3376: Input Parameters:
3377: + snes - the TS context
3378: - x - input vector
3380: Output Parameter:
3381: . y - function vector, as set by TSSetFunction()
3383: Notes:
3384: TSComputeFunction() is typically used within nonlinear solvers
3385: implementations, so most users would not generally call this routine
3386: themselves.
3388: Level: developer
3390: .keywords: TS, nonlinear, compute, function
3392: .seealso: TSSetFunction(), TSGetFunction()
3393: */
3394: PetscErrorCode TSComputeFunction_Matlab(TS snes,PetscReal time,Vec x,Vec xdot,Vec y, void *ctx)
3395: {
3396: PetscErrorCode ierr;
3397: TSMatlabContext *sctx = (TSMatlabContext *)ctx;
3398: int nlhs = 1,nrhs = 7;
3399: mxArray *plhs[1],*prhs[7];
3400: long long int lx = 0,lxdot = 0,ly = 0,ls = 0;
3410: PetscMemcpy(&ls,&snes,sizeof(snes));
3411: PetscMemcpy(&lx,&x,sizeof(x));
3412: PetscMemcpy(&lxdot,&xdot,sizeof(xdot));
3413: PetscMemcpy(&ly,&y,sizeof(x));
3414: prhs[0] = mxCreateDoubleScalar((double)ls);
3415: prhs[1] = mxCreateDoubleScalar(time);
3416: prhs[2] = mxCreateDoubleScalar((double)lx);
3417: prhs[3] = mxCreateDoubleScalar((double)lxdot);
3418: prhs[4] = mxCreateDoubleScalar((double)ly);
3419: prhs[5] = mxCreateString(sctx->funcname);
3420: prhs[6] = sctx->ctx;
3421: mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeFunctionInternal");
3422: mxGetScalar(plhs[0]);
3423: mxDestroyArray(prhs[0]);
3424: mxDestroyArray(prhs[1]);
3425: mxDestroyArray(prhs[2]);
3426: mxDestroyArray(prhs[3]);
3427: mxDestroyArray(prhs[4]);
3428: mxDestroyArray(prhs[5]);
3429: mxDestroyArray(plhs[0]);
3430: return(0);
3431: }
3436: /*
3437: TSSetFunctionMatlab - Sets the function evaluation routine and function
3438: vector for use by the TS routines in solving ODEs
3439: equations from MATLAB. Here the function is a string containing the name of a MATLAB function
3441: Logically Collective on TS
3443: Input Parameters:
3444: + ts - the TS context
3445: - func - function evaluation routine
3447: Calling sequence of func:
3448: $ func (TS ts,PetscReal time,Vec x,Vec xdot,Vec f,void *ctx);
3450: Level: beginner
3452: .keywords: TS, nonlinear, set, function
3454: .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
3455: */
3456: PetscErrorCode TSSetFunctionMatlab(TS ts,const char *func,mxArray *ctx)
3457: {
3458: PetscErrorCode ierr;
3459: TSMatlabContext *sctx;
3462: /* currently sctx is memory bleed */
3463: PetscMalloc(sizeof(TSMatlabContext),&sctx);
3464: PetscStrallocpy(func,&sctx->funcname);
3465: /*
3466: This should work, but it doesn't
3467: sctx->ctx = ctx;
3468: mexMakeArrayPersistent(sctx->ctx);
3469: */
3470: sctx->ctx = mxDuplicateArray(ctx);
3471: TSSetIFunction(ts,PETSC_NULL,TSComputeFunction_Matlab,sctx);
3472: return(0);
3473: }
3477: /*
3478: TSComputeJacobian_Matlab - Calls the function that has been set with
3479: TSSetJacobianMatlab().
3481: Collective on TS
3483: Input Parameters:
3484: + ts - the TS context
3485: . x - input vector
3486: . A, B - the matrices
3487: - ctx - user context
3489: Output Parameter:
3490: . flag - structure of the matrix
3492: Level: developer
3494: .keywords: TS, nonlinear, compute, function
3496: .seealso: TSSetFunction(), TSGetFunction()
3497: @*/
3498: PetscErrorCode TSComputeJacobian_Matlab(TS ts,PetscReal time,Vec x,Vec xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flag, void *ctx)
3499: {
3500: PetscErrorCode ierr;
3501: TSMatlabContext *sctx = (TSMatlabContext *)ctx;
3502: int nlhs = 2,nrhs = 9;
3503: mxArray *plhs[2],*prhs[9];
3504: long long int lx = 0,lxdot = 0,lA = 0,ls = 0, lB = 0;
3510: /* call Matlab function in ctx with arguments x and y */
3512: PetscMemcpy(&ls,&ts,sizeof(ts));
3513: PetscMemcpy(&lx,&x,sizeof(x));
3514: PetscMemcpy(&lxdot,&xdot,sizeof(x));
3515: PetscMemcpy(&lA,A,sizeof(x));
3516: PetscMemcpy(&lB,B,sizeof(x));
3517: prhs[0] = mxCreateDoubleScalar((double)ls);
3518: prhs[1] = mxCreateDoubleScalar((double)time);
3519: prhs[2] = mxCreateDoubleScalar((double)lx);
3520: prhs[3] = mxCreateDoubleScalar((double)lxdot);
3521: prhs[4] = mxCreateDoubleScalar((double)shift);
3522: prhs[5] = mxCreateDoubleScalar((double)lA);
3523: prhs[6] = mxCreateDoubleScalar((double)lB);
3524: prhs[7] = mxCreateString(sctx->funcname);
3525: prhs[8] = sctx->ctx;
3526: mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeJacobianInternal");
3527: mxGetScalar(plhs[0]);
3528: *flag = (MatStructure) mxGetScalar(plhs[1]);
3529: mxDestroyArray(prhs[0]);
3530: mxDestroyArray(prhs[1]);
3531: mxDestroyArray(prhs[2]);
3532: mxDestroyArray(prhs[3]);
3533: mxDestroyArray(prhs[4]);
3534: mxDestroyArray(prhs[5]);
3535: mxDestroyArray(prhs[6]);
3536: mxDestroyArray(prhs[7]);
3537: mxDestroyArray(plhs[0]);
3538: mxDestroyArray(plhs[1]);
3539: return(0);
3540: }
3545: /*
3546: TSSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices
3547: vector for use by the TS routines in solving ODEs from MATLAB. Here the function is a string containing the name of a MATLAB function
3549: Logically Collective on TS
3551: Input Parameters:
3552: + ts - the TS context
3553: . A,B - Jacobian matrices
3554: . func - function evaluation routine
3555: - ctx - user context
3557: Calling sequence of func:
3558: $ flag = func (TS ts,PetscReal time,Vec x,Vec xdot,Mat A,Mat B,void *ctx);
3561: Level: developer
3563: .keywords: TS, nonlinear, set, function
3565: .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
3566: */
3567: PetscErrorCode TSSetJacobianMatlab(TS ts,Mat A,Mat B,const char *func,mxArray *ctx)
3568: {
3569: PetscErrorCode ierr;
3570: TSMatlabContext *sctx;
3573: /* currently sctx is memory bleed */
3574: PetscMalloc(sizeof(TSMatlabContext),&sctx);
3575: PetscStrallocpy(func,&sctx->funcname);
3576: /*
3577: This should work, but it doesn't
3578: sctx->ctx = ctx;
3579: mexMakeArrayPersistent(sctx->ctx);
3580: */
3581: sctx->ctx = mxDuplicateArray(ctx);
3582: TSSetIJacobian(ts,A,B,TSComputeJacobian_Matlab,sctx);
3583: return(0);
3584: }
3588: /*
3589: TSMonitor_Matlab - Calls the function that has been set with TSMonitorSetMatlab().
3591: Collective on TS
3593: .seealso: TSSetFunction(), TSGetFunction()
3594: @*/
3595: PetscErrorCode TSMonitor_Matlab(TS ts,PetscInt it, PetscReal time,Vec x, void *ctx)
3596: {
3597: PetscErrorCode ierr;
3598: TSMatlabContext *sctx = (TSMatlabContext *)ctx;
3599: int nlhs = 1,nrhs = 6;
3600: mxArray *plhs[1],*prhs[6];
3601: long long int lx = 0,ls = 0;
3607: PetscMemcpy(&ls,&ts,sizeof(ts));
3608: PetscMemcpy(&lx,&x,sizeof(x));
3609: prhs[0] = mxCreateDoubleScalar((double)ls);
3610: prhs[1] = mxCreateDoubleScalar((double)it);
3611: prhs[2] = mxCreateDoubleScalar((double)time);
3612: prhs[3] = mxCreateDoubleScalar((double)lx);
3613: prhs[4] = mxCreateString(sctx->funcname);
3614: prhs[5] = sctx->ctx;
3615: mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSMonitorInternal");
3616: mxGetScalar(plhs[0]);
3617: mxDestroyArray(prhs[0]);
3618: mxDestroyArray(prhs[1]);
3619: mxDestroyArray(prhs[2]);
3620: mxDestroyArray(prhs[3]);
3621: mxDestroyArray(prhs[4]);
3622: mxDestroyArray(plhs[0]);
3623: return(0);
3624: }
3629: /*
3630: TSMonitorSetMatlab - Sets the monitor function from Matlab
3632: Level: developer
3634: .keywords: TS, nonlinear, set, function
3636: .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
3637: */
3638: PetscErrorCode TSMonitorSetMatlab(TS ts,const char *func,mxArray *ctx)
3639: {
3640: PetscErrorCode ierr;
3641: TSMatlabContext *sctx;
3644: /* currently sctx is memory bleed */
3645: PetscMalloc(sizeof(TSMatlabContext),&sctx);
3646: PetscStrallocpy(func,&sctx->funcname);
3647: /*
3648: This should work, but it doesn't
3649: sctx->ctx = ctx;
3650: mexMakeArrayPersistent(sctx->ctx);
3651: */
3652: sctx->ctx = mxDuplicateArray(ctx);
3653: TSMonitorSet(ts,TSMonitor_Matlab,sctx,PETSC_NULL);
3654: return(0);
3655: }
3656: #endif