petsc4py.PETSc.TS#
- class petsc4py.PETSc.TS#
Bases:
Object
ODE integrator.
TS is described in the
PETSc manual
.See also
Enumerations
The ARKIMEX subtype.
The reason the time step is converging.
The DIRK subtype.
Distinguishes among types of explicit and implicit equations.
The method for ending time stepping.
Distinguishes linear and nonlinear problems.
The RK subtype.
The time stepping method.
Methods Summary
Reset a
TS
, removing any allocated vectors and matrices.adjointSetSteps
(adjoint_steps)Set the number of steps the adjoint solver should take backward in time.
Set up the internal data structures for the later use of an adjoint solver.
Solve the discrete adjoint problem for an ODE/DAE.
Step one time step backward in the adjoint run.
appendOptionsPrefix
(prefix)Append to the prefix used for all the
TS
options.clone
()Return a shallow clone of the
TS
object.computeI2Function
(t, x, xdot, xdotdot, f)Evaluate the DAE residual in implicit form.
computeI2Jacobian
(t, x, xdot, xdotdot, v, a, J)Evaluate the Jacobian of the DAE.
computeIFunction
(t, x, xdot, f[, imex])Evaluate the DAE residual written in implicit form.
computeIJacobian
(t, x, xdot, a, J[, P, imex])Evaluate the Jacobian of the DAE.
computeIJacobianP
(t, x, xdot, a, J[, imex])Evaluate the Jacobian with respect to parameters.
computeRHSFunction
(t, x, f)Evaluate the right-hand side function.
computeRHSFunctionLinear
(t, x, f)Evaluate the right-hand side via the user-provided Jacobian.
computeRHSJacobian
(t, x, J[, P])Compute the Jacobian matrix that has been set with
setRHSJacobian
.computeRHSJacobianConstant
(t, x, J[, P])Reuse a Jacobian that is time-independent.
computeRHSJacobianP
(t, x, J)Run the user-defined JacobianP function.
create
([comm])Create an empty
TS
.createPython
([context, comm])Create an integrator of Python type.
createQuadratureTS
([forward])Create a sub
TS
that evaluates integrals over time.destroy
()Return the
Type.ARKIMEX
scheme.Return the algorithmic parameters for
Type.ALPHA
.Return the application context.
Return the reason the
TS
step was stopped.Return the cost gradients.
Return a vector of values of the integral term in the cost functions.
Return the
Type.DIRK
scheme.getDM
()Get the type of the equation that
TS
is solving.Return the vector and function which computes the residual.
Return the matrices and function which computes the Jacobian.
Return the vector and function which computes the implicit residual.
Return the matrices and function which computes the implicit Jacobian.
getKSP
()Return the total number of linear iterations used by the
TS
.Return the maximum number of steps to use.
Return the maximum (final) time.
Return the monitor.
Return the number of events.
Return the prefix used for all the
TS
options.Return the poststep function.
Return the prestep function.
Return the starting time of the previously completed step.
Return the type of problem to be solved.
Return the instance of the class implementing the required Python methods.
Return the fully qualified Python name of the class used by the solver.
Return the sub
TS
that evaluates integrals over time.Return the vector where the rhs is stored and the function used to compute it.
Return the Jacobian and the function used to compute them.
Return the
Type.RK
scheme.getSNES
()Return the total number of nonlinear iterations used by the
TS
.Return the solution at the present timestep.
Return the solution and time derivative at the present timestep.
Return the time after a call to
solve
.Return the minimum and maximum allowed time step sizes.
Return the number of time steps completed.
Return the total number of rejected steps.
getTheta
()Return the abscissa of the stage in
(0, 1]
forType.THETA
.Return whether the endpoint variable of
Type.THETA
is used.getTime
()Return the time of the most recently completed step.
Return the time span.
Return the solutions at the times in the time span.
Return the duration of the current timestep.
Return the tolerances for local truncation error.
getType
()Return the
TS
type.interpolate
(t, u)Interpolate the solution to a given time.
load
(viewer)monitor
(step, time[, u])Monitor the solve.
Clear all the monitors that have been set.
Remove the internal
TS
trajectory object.reset
()Reset the
TS
, removing any allocated vectors and matrices.Flag the solver to restart the next step.
rollBack
()Roll back one time step.
setARKIMEXFastSlowSplit
(flag)Use ARKIMEX for solving a fast-slow system.
setARKIMEXFullyImplicit
(flag)Solve both parts of the equation implicitly.
setARKIMEXType
(ts_type)Set the type of
Type.ARKIMEX
scheme.setAlphaParams
([alpha_m, alpha_f, gamma])Set the algorithmic parameters for
Type.ALPHA
.setAlphaRadius
(radius)Set the spectral radius for
Type.ALPHA
.setAppCtx
(appctx)Set the application context.
setConvergedReason
(reason)Set the reason for handling the convergence of
solve
.setCostGradients
(vl[, vm])Set the cost gradients.
setDIRKType
(ts_type)Set the type of
Type.DIRK
scheme.setDM
(dm)Set the DM that may be used by some nonlinear solvers or preconditioners.
setEquationType
(eqtype)Set the type of the equation that
TS
is solving.setErrorIfStepFails
([flag])Immediately error is no step succeeds.
setEventHandler
(direction, terminate, indicator)Set a function used for detecting events.
setEventTolerances
([tol, vtol])Set tolerances for event zero crossings when using event handler.
setExactFinalTime
(option)Set method of computing the final time step.
Set various
TS
parameters from user options.setI2Function
(function[, f, args, kargs])Set the function to compute the 2nd order DAE.
setI2Jacobian
(jacobian[, J, P, args, kargs])Set the function to compute the Jacobian of the 2nd order DAE.
setIFunction
(function[, f, args, kargs])Set the function representing the DAE to be solved.
setIJacobian
(jacobian[, J, P, args, kargs])Set the function to compute the Jacobian.
setIJacobianP
(jacobian[, J, args, kargs])Set the function that computes the Jacobian.
Set the maximum number of SNES solves failures allowed.
Set the maximum number of step rejections before a time step fails.
setMaxSteps
(max_steps)Set the maximum number of steps to use.
setMaxTime
(max_time)Set the maximum (final) time.
setMonitor
(monitor[, args, kargs])Set an additional monitor to the
TS
.setOptionsPrefix
(prefix)Set the prefix used for all the
TS
options.setPostStep
(poststep[, args, kargs])Set a function to be called at the end of each time step.
setPreStep
(prestep[, args, kargs])Set a function to be called at the beginning of each time step.
setProblemType
(ptype)Set the type of problem to be solved.
setPythonContext
(context)Set the instance of the class implementing the required Python methods.
setPythonType
(py_type)Set the fully qualified Python name of the class to be used.
setRHSFunction
(function[, f, args, kargs])Set the routine for evaluating the function
G
inU_t = G(t, u)
.setRHSJacobian
(jacobian[, J, P, args, kargs])Set the function to compute the Jacobian of
G
inU_t = G(U, t)
.setRHSJacobianP
(jacobianp[, A, args, kargs])Set the function that computes the Jacobian with respect to the parameters.
setRHSSplitIFunction
(splitname, function[, ...])Set the split implicit functions.
setRHSSplitIJacobian
(splitname, jacobian[, ...])Set the Jacobian for the split implicit function.
setRHSSplitIS
(splitname, iss)Set the index set for the specified split.
setRHSSplitRHSFunction
(splitname, function)Set the split right-hand-side functions.
setRKType
(ts_type)Set the type of the Runge-Kutta scheme.
Enable to save solutions as an internal
TS
trajectory.setSolution
(u)Set the initial solution vector.
setSolution2
(u, v)Set the initial solution and its time derivative.
setStepLimits
(hmin, hmax)Set the minimum and maximum allowed step sizes.
setStepNumber
(step_number)Set the number of steps completed.
setTheta
(theta)Set the abscissa of the stage in
(0, 1]
forType.THETA
.setThetaEndpoint
([flag])Set to use the endpoint variant of
Type.THETA
.setTime
(t)Set the time.
setTimeSpan
(tspan)Set the time span.
setTimeStep
(time_step)Set the duration of the timestep.
setTolerances
([rtol, atol])Set tolerances for local truncation error when using an adaptive controller.
setType
(ts_type)Set the method to be used as the
TS
solver.setUp
()Set up the internal data structures for the
TS
.solve
([u])Step the requested number of timesteps.
step
()Take one step.
view
([viewer])Print the
TS
object.Attributes Summary
Application context.
The absolute tolerance.
Indicates the
TS
has converged.Indicates the
TS
has stopped.The
DM
.The equation type.
Indicates the
TS
is still iterating.The
KSP
.The maximum number of steps.
The maximum time.
The problem type.
The converged reason.
The relative tolerance.
The
SNES
.The current step number.
The current time.
The current time step size.
The solution vector.
Methods Documentation
- adjointReset()#
Reset a
TS
, removing any allocated vectors and matrices.Collective.
See also
Source code at petsc4py/PETSc/TS.pyx:2849
- Return type:
- adjointSetSteps(adjoint_steps)#
Set the number of steps the adjoint solver should take backward in time.
Logically collective.
See also
- adjointSetUp()#
Set up the internal data structures for the later use of an adjoint solver.
Collective.
See also
Source code at petsc4py/PETSc/TS.pyx:2813
- Return type:
- adjointSolve()#
Solve the discrete adjoint problem for an ODE/DAE.
Collective.
See also
Source code at petsc4py/PETSc/TS.pyx:2825
- Return type:
- adjointStep()#
Step one time step backward in the adjoint run.
Collective.
See also
Source code at petsc4py/PETSc/TS.pyx:2837
- Return type:
- appendOptionsPrefix(prefix)#
Append to the prefix used for all the
TS
options.Logically collective.
Notes
A hyphen must not be given at the beginning of the prefix name.
- clone()#
Return a shallow clone of the
TS
object.Collective.
See also
Source code at petsc4py/PETSc/TS.pyx:243
- Return type:
- computeI2Function(t, x, xdot, xdotdot, f)#
Evaluate the DAE residual in implicit form.
Collective.
- Parameters:
- Return type:
See also
- computeI2Jacobian(t, x, xdot, xdotdot, v, a, J, P=None)#
Evaluate the Jacobian of the DAE.
Collective.
If
F(t, U, V, A)=0
is the DAE, the required Jacobian isdF/dU + v dF/dV + a dF/dA
.- Parameters:
t (float) – The current time.
x (Vec) – The state vector.
xdot (Vec) – The time derivative of the state vector.
xdotdot (Vec) – The second time derivative of the state vector.
v (float) – The shift to apply to the first derivative.
a (float) – The shift to apply to the second derivative.
J (Mat) – The matrix into which the Jacobian is computed.
P (Mat | None) – The optional matrix to use for building a preconditioner matrix.
- Return type:
See also
- computeIFunction(t, x, xdot, f, imex=False)#
Evaluate the DAE residual written in implicit form.
Collective.
- Parameters:
- Return type:
See also
- computeIJacobian(t, x, xdot, a, J, P=None, imex=False)#
Evaluate the Jacobian of the DAE.
Collective.
If
F(t, U, Udot)=0
is the DAE, the required Jacobian isdF/dU + shift*dF/dUdot
- Parameters:
t (float) – The current time.
x (Vec) – The state vector.
xdot (Vec) – The time derivative of the state vector.
a (float) – The shift to apply
J (Mat) – The matrix into which the Jacobian is computed.
P (Mat | None) – The optional matrix to use for building a preconditioner matrix.
imex (bool) – A flag which indicates if the RHS should be kept separate.
- Return type:
See also
- computeIJacobianP(t, x, xdot, a, J, imex=False)#
Evaluate the Jacobian with respect to parameters.
Collective.
- Parameters:
- Return type:
See also
- computeRHSFunction(t, x, f)#
Evaluate the right-hand side function.
Collective.
- Parameters:
- Return type:
See also
- computeRHSFunctionLinear(t, x, f)#
Evaluate the right-hand side via the user-provided Jacobian.
Collective.
- Parameters:
- Return type:
See also
- computeRHSJacobian(t, x, J, P=None)#
Compute the Jacobian matrix that has been set with
setRHSJacobian
.Collective.
- Parameters:
- Return type:
See also
- computeRHSJacobianConstant(t, x, J, P=None)#
Reuse a Jacobian that is time-independent.
Collective.
- Parameters:
- Return type:
See also
- computeRHSJacobianP(t, x, J)#
Run the user-defined JacobianP function.
Collective.
- Parameters:
- Return type:
See also
- create(comm=None)#
Create an empty
TS
.Collective.
The problem type can then be set with
setProblemType
and the type of solver can then be set withsetType
.- Parameters:
comm (Comm | None) – MPI communicator, defaults to
Sys.getDefaultComm
.- Return type:
See also
- createPython(context=None, comm=None)#
Create an integrator of Python type.
Collective.
- Parameters:
context (Any) – An instance of the Python class implementing the required methods.
comm (Comm | None) – MPI communicator, defaults to
Sys.getDefaultComm
.
- Return type:
- createQuadratureTS(forward=True)#
Create a sub
TS
that evaluates integrals over time.Collective.
See also
- destroy()#
Destroy the
TS
that was created withcreate
.Collective.
See also
Source code at petsc4py/PETSc/TS.pyx:206
- Return type:
- getARKIMEXType()#
Return the
Type.ARKIMEX
scheme.Not collective.
See also
Source code at petsc4py/PETSc/TS.pyx:389
- Return type:
- getAlphaParams()#
Return the algorithmic parameters for
Type.ALPHA
.Not collective.
See also
- getAppCtx()#
Return the application context.
Source code at petsc4py/PETSc/TS.pyx:589
- Return type:
- getConvergedReason()#
Return the reason the
TS
step was stopped.Not collective.
Can only be called once
solve
is complete.See also
Source code at petsc4py/PETSc/TS.pyx:2090
- Return type:
- getCostGradients()#
Return the cost gradients.
Not collective.
See also
- getCostIntegral()#
Return a vector of values of the integral term in the cost functions.
Not collective.
See also
Source code at petsc4py/PETSc/TS.pyx:2573
- Return type:
- getDIRKType()#
Return the
Type.DIRK
scheme.Not collective.
See also
Source code at petsc4py/PETSc/TS.pyx:426
- Return type:
- getDM()#
Return the
DM
associated with theTS
.Not collective.
Only valid if nonlinear solvers or preconditioners are used which use the
DM
.See also
Source code at petsc4py/PETSc/TS.pyx:1589
- Return type:
- getEquationType()#
Get the type of the equation that
TS
is solving.Not collective.
See also
Source code at petsc4py/PETSc/TS.pyx:488
- Return type:
- getI2Function()#
Return the vector and function which computes the residual.
Not collective.
See also
Source code at petsc4py/PETSc/TS.pyx:1218
- Return type:
- getI2Jacobian()#
Return the matrices and function which computes the Jacobian.
Not collective.
See also
Source code at petsc4py/PETSc/TS.pyx:1234
- Return type:
tuple[Mat, Mat, TSI2Jacobian]
- getIFunction()#
Return the vector and function which computes the implicit residual.
Not collective.
See also
Source code at petsc4py/PETSc/TS.pyx:1031
- Return type:
- getIJacobian()#
Return the matrices and function which computes the implicit Jacobian.
Not collective.
See also
Source code at petsc4py/PETSc/TS.pyx:1047
- Return type:
tuple[Mat, Mat, TSIJacobian]
- getKSP()#
Return the
KSP
associated with theTS
.Not collective.
See also
Source code at petsc4py/PETSc/TS.pyx:1572
- Return type:
- getKSPIterations()#
Return the total number of linear iterations used by the
TS
.Not collective.
This counter is reset to zero for each successive call to
solve
.See also
Source code at petsc4py/PETSc/TS.pyx:1857
- Return type:
- getMaxSteps()#
Return the maximum number of steps to use.
Not collective.
See also
Source code at petsc4py/PETSc/TS.pyx:1826
- Return type:
- getMaxTime()#
Return the maximum (final) time.
Not collective.
Defaults to
5
.See also
Source code at petsc4py/PETSc/TS.pyx:1790
- Return type:
- getMonitor()#
Return the monitor.
Not collective.
See also
- getNumEvents()#
Return the number of events.
Logically collective.
See also
Source code at petsc4py/PETSc/TS.pyx:2289
- Return type:
- getOptionsPrefix()#
Return the prefix used for all the
TS
options.Not collective.
See also
Source code at petsc4py/PETSc/TS.pyx:525
- Return type:
- getPostStep()#
Return the poststep function.
- getPreStep()#
Return the prestep function.
Not collective.
See also
- getPrevTime()#
Return the starting time of the previously completed step.
Not collective.
See also
Source code at petsc4py/PETSc/TS.pyx:1664
- Return type:
- getProblemType()#
Return the type of problem to be solved.
Not collective.
See also
Source code at petsc4py/PETSc/TS.pyx:457
- Return type:
- getPythonContext()#
Return the instance of the class implementing the required Python methods.
Not collective.
Source code at petsc4py/PETSc/TS.pyx:2900
- Return type:
- getPythonType()#
Return the fully qualified Python name of the class used by the solver.
Not collective.
Source code at petsc4py/PETSc/TS.pyx:2929
- Return type:
- getQuadratureTS()#
Return the sub
TS
that evaluates integrals over time.Not collective.
- Returns:
- Return type:
See also
- getRHSFunction()#
Return the vector where the rhs is stored and the function used to compute it.
Not collective.
See also
Source code at petsc4py/PETSc/TS.pyx:770
- Return type:
- getRHSJacobian()#
Return the Jacobian and the function used to compute them.
Not collective.
See also
Source code at petsc4py/PETSc/TS.pyx:786
- Return type:
tuple[Mat, Mat, TSRHSJacobian]
- getRKType()#
Return the
Type.RK
scheme.Not collective.
See also
Source code at petsc4py/PETSc/TS.pyx:375
- Return type:
- getSNES()#
Return the
SNES
associated with theTS
.Not collective.
See also
Source code at petsc4py/PETSc/TS.pyx:1557
- Return type:
- getSNESFailures()#
Return the total number of failed
SNES
solves in theTS
.Not collective.
This counter is reset to zero for each successive call to
solve
.See also
Source code at petsc4py/PETSc/TS.pyx:1931
- Return type:
- getSNESIterations()#
Return the total number of nonlinear iterations used by the
TS
.Not collective.
This counter is reset to zero for each successive call to
solve
.See also
Source code at petsc4py/PETSc/TS.pyx:1840
- Return type:
- getSolution()#
Return the solution at the present timestep.
Not collective.
It is valid to call this routine inside the function that you are evaluating in order to move to the new timestep. This vector is not changed until the solution at the next timestep has been calculated.
See also
Source code at petsc4py/PETSc/TS.pyx:1428
- Return type:
- getSolution2()#
Return the solution and time derivative at the present timestep.
Not collective.
It is valid to call this routine inside the function that you are evaluating in order to move to the new timestep. These vectors are not changed until the solution at the next timestep has been calculated.
See also
- getSolveTime()#
Return the time after a call to
solve
.Not collective.
This time corresponds to the final time set with
setMaxTime
.See also
Source code at petsc4py/PETSc/TS.pyx:1678
- Return type:
- getStepLimits()#
Return the minimum and maximum allowed time step sizes.
Not collective.
See also
- getStepNumber()#
Return the number of time steps completed.
Not collective.
See also
Source code at petsc4py/PETSc/TS.pyx:1754
- Return type:
- getStepRejections()#
Return the total number of rejected steps.
Not collective.
This counter is reset to zero for each successive call to
solve
.See also
Source code at petsc4py/PETSc/TS.pyx:1896
- Return type:
- getTheta()#
Return the abscissa of the stage in
(0, 1]
forType.THETA
.Not collective.
See also
Source code at petsc4py/PETSc/TS.pyx:2967
- Return type:
- getThetaEndpoint()#
Return whether the endpoint variable of
Type.THETA
is used.Not collective.
See also
Source code at petsc4py/PETSc/TS.pyx:2999
- Return type:
- getTime()#
Return the time of the most recently completed step.
Not collective.
When called during time step evaluation (e.g. during residual evaluation or via hooks set using
setPreStep
orsetPostStep
), the time returned is at the start of the step.See also
Source code at petsc4py/PETSc/TS.pyx:1646
- Return type:
- getTimeSpan()#
Return the time span.
Not collective.
See also
Source code at petsc4py/PETSc/TS.pyx:1521
- Return type:
- getTimeSpanSolutions()#
Return the solutions at the times in the time span.
Not collective.
See also
- getTimeStep()#
Return the duration of the current timestep.
Not collective.
See also
Source code at petsc4py/PETSc/TS.pyx:1713
- Return type:
- getTolerances()#
Return the tolerances for local truncation error.
Logically collective.
- Returns:
- Return type:
See also
- getType()#
Return the
TS
type.Not collective.
See also
Source code at petsc4py/PETSc/TS.pyx:361
- Return type:
- interpolate(t, u)#
Interpolate the solution to a given time.
Collective.
- Parameters:
- Return type:
See also
- monitor(step, time, u=None)#
Monitor the solve.
Collective.
- Parameters:
- Return type:
See also
- monitorCancel()#
Clear all the monitors that have been set.
Logically collective.
See also
Source code at petsc4py/PETSc/TS.pyx:2154
- Return type:
- removeTrajectory()#
Remove the internal
TS
trajectory object.Collective.
See also
Source code at petsc4py/PETSc/TS.pyx:2561
- Return type:
- reset()#
Reset the
TS
, removing any allocated vectors and matrices.Collective.
See also
Source code at petsc4py/PETSc/TS.pyx:2398
- Return type:
- restartStep()#
Flag the solver to restart the next step.
Collective.
Multistep methods like TSBDF or Runge-Kutta methods with FSAL property require restarting the solver in the event of discontinuities. These discontinuities may be introduced as a consequence of explicitly modifications to the solution vector (which PETSc attempts to detect and handle) or problem coefficients (which PETSc is not able to detect). For the sake of correctness and maximum safety, users are expected to call TSRestart() whenever they introduce discontinuities in callback routines (e.g. prestep and poststep routines, or implicit/rhs function routines with discontinuous source terms).
See also
Source code at petsc4py/PETSc/TS.pyx:2426
- Return type:
- rollBack()#
Roll back one time step.
Collective.
See also
Source code at petsc4py/PETSc/TS.pyx:2449
- Return type:
- setARKIMEXFastSlowSplit(flag)#
Use ARKIMEX for solving a fast-slow system.
Logically collective.
See also
- setARKIMEXFullyImplicit(flag)#
Solve both parts of the equation implicitly.
Logically collective.
See also
- setARKIMEXType(ts_type)#
Set the type of
Type.ARKIMEX
scheme.Logically collective.
- Parameters:
ts_type (ARKIMEXType | str) – The type of
Type.ARKIMEX
scheme.- Return type:
Notes
-ts_arkimex_type
sets scheme type from the commandline.See also
- setAlphaParams(alpha_m=None, alpha_f=None, gamma=None)#
Set the algorithmic parameters for
Type.ALPHA
.Logically collective.
Users should call
setAlphaRadius
.- Parameters:
- Return type:
See also
- setAlphaRadius(radius)#
Set the spectral radius for
Type.ALPHA
.Logically collective.
Notes
-ts_alpha_radius
can be used to set this from the commandline.See also
- setAppCtx(appctx)#
Set the application context.
Not collective.
- setConvergedReason(reason)#
Set the reason for handling the convergence of
solve
.Logically collective.
Can only be called when
solve
is active andreason
must contain common value.- Parameters:
reason (ConvergedReason) – The reason for convergence.
- Return type:
See also
- setCostGradients(vl, vm=None)#
Set the cost gradients.
Logically collective.
- Parameters:
vl (Vec | Sequence[Vec] | None) – gradients with respect to the initial condition variables, the dimension and parallel layout of these vectors is the same as the ODE solution vector
vm (Vec | Sequence[Vec] | None) – gradients with respect to the parameters, the number of entries in these vectors is the same as the number of parameters
- Return type:
See also
- setDIRKType(ts_type)#
Set the type of
Type.DIRK
scheme.Logically collective.
Notes
-ts_dirk_type
sets scheme type from the commandline.See also
- setDM(dm)#
Set the DM that may be used by some nonlinear solvers or preconditioners.
Logically collective.
See also
- setEquationType(eqtype)#
Set the type of the equation that
TS
is solving.Logically collective.
- Parameters:
eqtype (EquationType) – The type of equation.
- Return type:
See also
- setErrorIfStepFails(flag=True)#
Immediately error is no step succeeds.
Not collective.
Notes
-ts_error_if_step_fails
to enable from the commandline.See also
- setEventHandler(direction, terminate, indicator, postevent=None, args=None, kargs=None)#
Set a function used for detecting events.
Logically collective.
- Parameters:
direction (Sequence[int]) – Direction of zero crossing to be detected {-1, 0, +1}.
terminate (Sequence[bool]) – Flags for each event to indicate stepping should be terminated.
indicator (TSIndicatorFunction | None) – Function for defining the indicator-functions marking the events
postevent (TSPostEventFunction) – Function to execute after the event
args (tuple[Any, ...] | None) – Additional positional arguments for
indicator
.kargs (dict[str, Any] | None) – Additional keyword arguments for
indicator
.
- Return type:
See also
- setEventTolerances(tol=None, vtol=None)#
Set tolerances for event zero crossings when using event handler.
Logically collective.
setEventHandler
must have already been called.- Parameters:
- Return type:
Notes
-ts_event_tol
can be used to set values from the commandline.See also
- setExactFinalTime(option)#
Set method of computing the final time step.
Logically collective.
- Parameters:
option (ExactFinalTime) – The exact final time option
- Return type:
Notes
-ts_exact_final_time
may be used to specify from the commandline.See also
- setFromOptions()#
Set various
TS
parameters from user options.Collective.
See also
Source code at petsc4py/PETSc/TS.pyx:562
- Return type:
- setI2Function(function, f=None, args=None, kargs=None)#
Set the function to compute the 2nd order DAE.
Logically collective.
- Parameters:
- Return type:
See also
- setI2Jacobian(jacobian, J=None, P=None, args=None, kargs=None)#
Set the function to compute the Jacobian of the 2nd order DAE.
Logically collective.
- Parameters:
jacobian (TSI2Jacobian | None) – The function which computes the Jacobian.
J (Mat | None) – The matrix into which the Jacobian is computed.
P (Mat | None) – The optional matrix to use for building a preconditioner matrix.
args – Additional positional arguments for
jacobian
.kargs – Additional keyword arguments for
jacobian
.
- Return type:
See also
- setIFunction(function, f=None, args=None, kargs=None)#
Set the function representing the DAE to be solved.
Logically collective.
- Parameters:
- Return type:
See also
- setIJacobian(jacobian, J=None, P=None, args=None, kargs=None)#
Set the function to compute the Jacobian.
Logically collective.
Set the function to compute the matrix
dF/dU + a*dF/dU_t
whereF(t, U, U_t)
is the function provided withsetIFunction
.- Parameters:
jacobian (TSIJacobian | None) – The function which computes the Jacobian.
J (Mat | None) – The matrix into which the Jacobian is computed.
P (Mat | None) – The optional matrix to use for building a preconditioner matrix.
args (tuple[Any, ...] | None) – Additional positional arguments for
jacobian
.kargs (dict[str, Any] | None) – Additional keyword arguments for
jacobian
.
- Return type:
See also
- setIJacobianP(jacobian, J=None, args=None, kargs=None)#
Set the function that computes the Jacobian.
Logically collective.
Set the function that computes the Jacobian of
F
with respect to the parametersP
whereF(Udot, U, t) = G(U, P, t)
, as well as the location to store the matrix.- Parameters:
- Return type:
See also
- setMaxSNESFailures(n)#
Set the maximum number of SNES solves failures allowed.
Not collective.
- Parameters:
n (int) – The maximum number of failed nonlinear solver, use
-1
for unlimited.- Return type:
See also
- setMaxStepRejections(n)#
Set the maximum number of step rejections before a time step fails.
Not collective.
Notes
-ts_max_reject
can be used to set this from the commandlineSee also
- setMaxSteps(max_steps)#
Set the maximum number of steps to use.
Logically collective.
Defaults to
5000
.See also
- setMaxTime(max_time)#
Set the maximum (final) time.
Logically collective.
Notes
-ts_max_time
sets the max time from the commandlineSee also
- setMonitor(monitor, args=None, kargs=None)#
Set an additional monitor to the
TS
.Logically collective.
- Parameters:
- Return type:
See also
- setOptionsPrefix(prefix)#
Set the prefix used for all the
TS
options.Logically collective.
Notes
A hyphen must not be given at the beginning of the prefix name.
See also
- setPostStep(poststep, args=None, kargs=None)#
Set a function to be called at the end of each time step.
Logically collective.
- Parameters:
- Return type:
See also
- setPreStep(prestep, args=None, kargs=None)#
Set a function to be called at the beginning of each time step.
Logically collective.
- Parameters:
- Return type:
See also
- setProblemType(ptype)#
Set the type of problem to be solved.
Logically collective.
- Parameters:
ptype (ProblemType) – The type of problem of the forms.
- Return type:
See also
- setPythonContext(context)#
Set the instance of the class implementing the required Python methods.
Not collective.
- setPythonType(py_type)#
Set the fully qualified Python name of the class to be used.
Collective.
- setRHSFunction(function, f=None, args=None, kargs=None)#
Set the routine for evaluating the function
G
inU_t = G(t, u)
.Logically collective.
- Parameters:
- Return type:
See also
- setRHSJacobian(jacobian, J=None, P=None, args=None, kargs=None)#
Set the function to compute the Jacobian of
G
inU_t = G(U, t)
.Logically collective.
- Parameters:
jacobian (TSRHSJacobian | None) – The right-hand side function.
J (Mat | None) – The matrix into which the jacobian is computed.
P (Mat | None) – The matrix into which the preconditioner is computed.
args (tuple[Any, ...] | None) – Additional positional arguments for
jacobian
.kargs (dict[str, Any] | None) – Additional keyword arguments for
jacobian
.
- Return type:
See also
- setRHSJacobianP(jacobianp, A=None, args=None, kargs=None)#
Set the function that computes the Jacobian with respect to the parameters.
Logically collective.
- Parameters:
- Return type:
See also
- setRHSSplitIFunction(splitname, function, r=None, args=None, kargs=None)#
Set the split implicit functions.
Logically collective.
- Parameters:
- Return type:
See also
- setRHSSplitIJacobian(splitname, jacobian, J=None, P=None, args=None, kargs=None)#
Set the Jacobian for the split implicit function.
Logically collective.
- Parameters:
splitname (str) – Name of the split.
jacobian (TSRHSJacobian) – The Jacobian evaluation routine.
J (Mat | None) – Matrix to store Jacobian entries computed by
jacobian
.P (Mat | None) – Matrix used to compute preconditioner (usually the same as
J
).args (tuple[Any, ...] | None) – Additional positional arguments for
jacobian
.kargs (dict[str, Any] | None) – Additional keyword arguments for
jacobian
.
- Return type:
See also
- setRHSSplitIS(splitname, iss)#
Set the index set for the specified split.
Logically collective.
- Parameters:
- Return type:
See also
- setRHSSplitRHSFunction(splitname, function, r=None, args=None, kargs=None)#
Set the split right-hand-side functions.
Logically collective.
- Parameters:
- Return type:
See also
- setRKType(ts_type)#
Set the type of the Runge-Kutta scheme.
Logically collective.
Notes
-ts_rk_type
sets scheme type from the commandline.See also
- setSaveTrajectory()#
Enable to save solutions as an internal
TS
trajectory.Collective.
This routine should be called after all
TS
options have been set.Notes
-ts_save_trajectory
can be used to save a trajectory to a file.See also
Source code at petsc4py/PETSc/TS.pyx:2542
- Return type:
- setSolution(u)#
Set the initial solution vector.
Logically collective.
See also
- setSolution2(u, v)#
Set the initial solution and its time derivative.
Logically collective.
See also
- setStepLimits(hmin, hmax)#
Set the minimum and maximum allowed step sizes.
Logically collective.
- Parameters:
- Return type:
See also
- setStepNumber(step_number)#
Set the number of steps completed.
Logically collective.
For most uses of the
TS
solvers the user need not explicitly callsetStepNumber
, as the step counter is appropriately updated insolve
/step
/rollBack
. Power users may call this routine to reinitialize timestepping by setting the step counter to zero (and time to the initial time) to solve a similar problem with different initial conditions or parameters. It may also be used to continue timestepping from a previously interrupted run in such a way thatTS
monitors will be called with a initial nonzero step counter.See also
- setTheta(theta)#
Set the abscissa of the stage in
(0, 1]
forType.THETA
.Logically collective.
Notes
-ts_theta_theta
can be used to set a value from the commandline.See also
- setThetaEndpoint(flag=True)#
Set to use the endpoint variant of
Type.THETA
.Logically collective.
- Parameters:
flag – Enable to use the endpoint variant.
- Return type:
See also
- setTime(t)#
Set the time.
Logically collective.
See also
- setTimeSpan(tspan)#
Set the time span.
Collective.
The solution will be computed and stored for each time requested in the span. The times must be all increasing and correspond to the intermediate points for time integration.
ExactFinalTime.MATCHSTEP
must be used to make the last time step in each sub-interval match the intermediate points specified. The intermediate solutions are saved in a vector array that can be accessed withgetTimeSpanSolutions
.Notes
-ts_time_span <t0, ..., tf>
sets the time span from the commandlineSee also
- setTimeStep(time_step)#
Set the duration of the timestep.
Logically collective.
See also
- setTolerances(rtol=None, atol=None)#
Set tolerances for local truncation error when using an adaptive controller.
Logically collective.
- Parameters:
- Return type:
Notes
-ts_rtol
and-ts_atol
may be used to set values from the commandline.See also
- setType(ts_type)#
Set the method to be used as the
TS
solver.Collective.
Notes
-ts_type
sets the method from the commandlineSee also
- setUp()#
Set up the internal data structures for the
TS
.Collective.
See also
Source code at petsc4py/PETSc/TS.pyx:2386
- Return type:
- solve(u=None)#
Step the requested number of timesteps.
Collective.
See also
- step()#
Take one step.
Collective.
The preferred interface for the
TS
solvers issolve
. If you need to execute code at the beginning or ending of each step, usesetPreStep
andsetPostStep
respectively.See also
Source code at petsc4py/PETSc/TS.pyx:2410
- Return type:
- view(viewer=None)#
Print the
TS
object.Collective.
Notes
-ts_view
calls TSView at the end of TSStepSee also
Attributes Documentation
- appctx#
Application context.
- atol#
The absolute tolerance.
- equation_type#
The equation type.
- max_steps#
The maximum number of steps.
- max_time#
The maximum time.
- problem_type#
The problem type.
- reason#
The converged reason.
- rtol#
The relative tolerance.
- step_number#
The current step number.
- time#
The current time.
- time_step#
The current time step size.
- vec_sol#
The solution vector.