petsc4py.PETSc.TS#

class petsc4py.PETSc.TS#

Bases: Object

ODE integrator.

TS is described in the PETSc manual.

See also

TS

Enumerations

ARKIMEXType

The ARKIMEX subtype.

ConvergedReason

The reason the time step is converging.

DIRKType

The DIRK subtype.

EquationType

Distinguishes among types of explicit and implicit equations.

ExactFinalTime

The method for ending time stepping.

ProblemType

Distinguishes linear and nonlinear problems.

RKType

The RK subtype.

Type

The time stepping method.

Methods Summary

adjointReset()

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.

adjointSetUp()

Set up the internal data structures for the later use of an adjoint solver.

adjointSolve()

Solve the discrete adjoint problem for an ODE/DAE.

adjointStep()

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()

Destroy the TS that was created with create.

getARKIMEXType()

Return the Type.ARKIMEX scheme.

getAlphaParams()

Return the algorithmic parameters for Type.ALPHA.

getAppCtx()

Return the application context.

getConvergedReason()

Return the reason the TS step was stopped.

getCostGradients()

Return the cost gradients.

getCostIntegral()

Return a vector of values of the integral term in the cost functions.

getDIRKType()

Return the Type.DIRK scheme.

getDM()

Return the DM associated with the TS.

getEquationType()

Get the type of the equation that TS is solving.

getI2Function()

Return the vector and function which computes the residual.

getI2Jacobian()

Return the matrices and function which computes the Jacobian.

getIFunction()

Return the vector and function which computes the implicit residual.

getIJacobian()

Return the matrices and function which computes the implicit Jacobian.

getKSP()

Return the KSP associated with the TS.

getKSPIterations()

Return the total number of linear iterations used by the TS.

getMaxSteps()

Return the maximum number of steps to use.

getMaxTime()

Return the maximum (final) time.

getMonitor()

Return the monitor.

getNumEvents()

Return the number of events.

getOptionsPrefix()

Return the prefix used for all the TS options.

getPostStep()

Return the poststep function.

getPreStep()

Return the prestep function.

getPrevTime()

Return the starting time of the previously completed step.

getProblemType()

Return the type of problem to be solved.

getPythonContext()

Return the instance of the class implementing the required Python methods.

getPythonType()

Return the fully qualified Python name of the class used by the solver.

getQuadratureTS()

Return the sub TS that evaluates integrals over time.

getRHSFunction()

Return the vector where the rhs is stored and the function used to compute it.

getRHSJacobian()

Return the Jacobian and the function used to compute them.

getRKType()

Return the Type.RK scheme.

getSNES()

Return the SNES associated with the TS.

getSNESFailures()

Return the total number of failed SNES solves in the TS.

getSNESIterations()

Return the total number of nonlinear iterations used by the TS.

getSolution()

Return the solution at the present timestep.

getSolution2()

Return the solution and time derivative at the present timestep.

getSolveTime()

Return the time after a call to solve.

getStepLimits()

Return the minimum and maximum allowed time step sizes.

getStepNumber()

Return the number of time steps completed.

getStepRejections()

Return the total number of rejected steps.

getTheta()

Return the abscissa of the stage in (0, 1] for Type.THETA.

getThetaEndpoint()

Return whether the endpoint variable of Type.THETA is used.

getTime()

Return the time of the most recently completed step.

getTimeSpan()

Return the time span.

getTimeSpanSolutions()

Return the solutions at the times in the time span.

getTimeStep()

Return the duration of the current timestep.

getTolerances()

Return the tolerances for local truncation error.

getType()

Return the TS type.

interpolate(t, u)

Interpolate the solution to a given time.

load(viewer)

Load a TS that has been stored in binary with view.

monitor(step, time[, u])

Monitor the solve.

monitorCancel()

Clear all the monitors that have been set.

removeTrajectory()

Remove the internal TS trajectory object.

reset()

Reset the TS, removing any allocated vectors and matrices.

restartStep()

Flag the solver to restart the next step.

rollBack()

Roll back one time step.

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.

setFromOptions()

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.

setMaxSNESFailures(n)

Set the maximum number of SNES solves failures allowed.

setMaxStepRejections(n)

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 in U_t = G(t, u).

setRHSJacobian(jacobian[, J, P, args, kargs])

Set the function to compute the Jacobian of G in U_t = G(U, t).

setRHSJacobianP(jacobianp[, A, args, kargs])

Set the function that computes the Jacobian with respect to the parameters.

setRKType(ts_type)

Set the type of the Runge-Kutta scheme.

setSaveTrajectory()

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] for Type.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

appctx

Application context.

atol

The absolute tolerance.

converged

Indicates the TS has converged.

diverged

Indicates the TS has stopped.

dm

The DM.

equation_type

The equation type.

iterating

Indicates the TS is still iterating.

ksp

The KSP.

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.

snes

The SNES.

step_number

The current step number.

time

The current time.

time_step

The current time step size.

vec_sol

The solution vector.

Methods Documentation

adjointReset()#

Reset a TS, removing any allocated vectors and matrices.

Collective.

See also

TSAdjointReset

Source code at petsc4py/PETSc/TS.pyx:2669

Return type:

None

adjointSetSteps(adjoint_steps)#

Set the number of steps the adjoint solver should take backward in time.

Logically collective.

Parameters:

adjoint_steps (int) – The number of steps to take.

Return type:

None

Source code at petsc4py/PETSc/TS.pyx:2615

adjointSetUp()#

Set up the internal data structures for the later use of an adjoint solver.

Collective.

See also

TSAdjointSetUp

Source code at petsc4py/PETSc/TS.pyx:2633

Return type:

None

adjointSolve()#

Solve the discrete adjoint problem for an ODE/DAE.

Collective.

See also

TSAdjointSolve

Source code at petsc4py/PETSc/TS.pyx:2645

Return type:

None

adjointStep()#

Step one time step backward in the adjoint run.

Collective.

See also

TSAdjointStep

Source code at petsc4py/PETSc/TS.pyx:2657

Return type:

None

appendOptionsPrefix(prefix)#

Append to the prefix used for all the TS options.

Logically collective.

Parameters:

prefix (str | None) – The prefix to append to the current prefix.

Return type:

None

Notes

A hyphen must not be given at the beginning of the prefix name.

Source code at petsc4py/PETSc/TS.pyx:522

clone()#

Return a shallow clone of the TS object.

Collective.

See also

TSClone

Source code at petsc4py/PETSc/TS.pyx:243

Return type:

TS

computeI2Function(t, x, xdot, xdotdot, f)#

Evaluate the DAE residual in implicit form.

Collective.

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.

  • f (Vec) – The vector into which the residual is stored.

Return type:

None

Source code at petsc4py/PETSc/TS.pyx:1125

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 is dF/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:

None

Source code at petsc4py/PETSc/TS.pyx:1152

computeIFunction(t, x, xdot, f, imex=False)#

Evaluate the DAE residual written in implicit form.

Collective.

Parameters:
  • t (float) – The current time.

  • x (Vec) – The state vector.

  • xdot (Vec) – The time derivative of the state vector.

  • f (Vec) – The vector into which the residual is stored.

  • imex (bool) – A flag which indicates if the RHS should be kept separate.

Return type:

None

Source code at petsc4py/PETSc/TS.pyx:910

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 is dF/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:

None

Source code at petsc4py/PETSc/TS.pyx:940

computeIJacobianP(t, x, xdot, a, J, imex=False)#

Evaluate the Jacobian with respect to parameters.

Collective.

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.

  • imex (bool) – A flag which indicates if the RHS should be kept separate.

Return type:

None

Source code at petsc4py/PETSc/TS.pyx:980

computeRHSFunction(t, x, f)#

Evaluate the right-hand side function.

Collective.

Parameters:
  • t (float) – The time at which to evaluate the RHS.

  • x (Vec) – The state vector.

  • f (Vec) – The Vec into which the RHS is computed.

Return type:

None

Source code at petsc4py/PETSc/TS.pyx:657

computeRHSFunctionLinear(t, x, f)#

Evaluate the right-hand side via the user-provided Jacobian.

Collective.

Parameters:
  • t (float) – The time at which to evaluate the RHS.

  • x (Vec) – The state vector.

  • f (Vec) – The Vec into which the RHS is computed.

Return type:

None

Source code at petsc4py/PETSc/TS.pyx:679

computeRHSJacobian(t, x, J, P=None)#

Compute the Jacobian matrix that has been set with setRHSJacobian.

Collective.

Parameters:
  • t (float) – The time at which to evaluate the Jacobian.

  • x (Vec) – The state vector.

  • 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:

None

Source code at petsc4py/PETSc/TS.pyx:701

computeRHSJacobianConstant(t, x, J, P=None)#

Reuse a Jacobian that is time-independent.

Collective.

Parameters:
  • t (float) – The time at which to evaluate the Jacobian.

  • x (Vec) – The state vector.

  • J (Mat) – A pointer to the stored Jacobian.

  • P (Mat | None) – An optional pointer to the preconditioner matrix.

Return type:

None

Source code at petsc4py/PETSc/TS.pyx:727

computeRHSJacobianP(t, x, J)#

Run the user-defined JacobianP function.

Collective.

Parameters:
  • t (float) – The time at which to compute the Jacobian.

  • x (Vec) – The solution at which to compute the Jacobian.

  • J (Mat) – The output Jacobian matrx.

Return type:

None

Source code at petsc4py/PETSc/TS.pyx:2593

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 with setType.

Parameters:

comm (Comm | None) – MPI communicator, defaults to Sys.getDefaultComm.

Return type:

Self

See also

TSCreate

Source code at petsc4py/PETSc/TS.pyx:219

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:

Self

Source code at petsc4py/PETSc/TS.pyx:2683

createQuadratureTS(forward=True)#

Create a sub TS that evaluates integrals over time.

Collective.

Parameters:

forward (bool) – Enable to evaluate forward in time.

Return type:

TS

Source code at petsc4py/PETSc/TS.pyx:2512

destroy()#

Destroy the TS that was created with create.

Collective.

See also

TSDestroy

Source code at petsc4py/PETSc/TS.pyx:206

Return type:

Self

getARKIMEXType()#

Return the Type.ARKIMEX scheme.

Not collective.

See also

TSARKIMEXGetType

Source code at petsc4py/PETSc/TS.pyx:372

Return type:

str

getAlphaParams()#

Return the algorithmic parameters for Type.ALPHA.

Not collective.

See also

TSAlphaGetParams

Source code at petsc4py/PETSc/TS.pyx:2890

Return type:

tuple[float, float, float]

getAppCtx()#

Return the application context.

Source code at petsc4py/PETSc/TS.pyx:572

Return type:

Any

getConvergedReason()#

Return the reason the TS step was stopped.

Not collective.

Can only be called once solve is complete.

Source code at petsc4py/PETSc/TS.pyx:1910

Return type:

ConvergedReason

getCostGradients()#

Return the cost gradients.

Not collective.

Source code at petsc4py/PETSc/TS.pyx:2454

Return type:

tuple[list[Vec], list[Vec]]

getCostIntegral()#

Return a vector of values of the integral term in the cost functions.

Not collective.

Source code at petsc4py/PETSc/TS.pyx:2393

Return type:

Vec

getDIRKType()#

Return the Type.DIRK scheme.

Not collective.

Source code at petsc4py/PETSc/TS.pyx:409

Return type:

str

getDM()#

Return the DM associated with the TS.

Not collective.

Only valid if nonlinear solvers or preconditioners are used which use the DM.

See also

TSGetDM

Source code at petsc4py/PETSc/TS.pyx:1413

Return type:

DM

getEquationType()#

Get the type of the equation that TS is solving.

Not collective.

Source code at petsc4py/PETSc/TS.pyx:471

Return type:

EquationType

getI2Function()#

Return the vector and function which computes the residual.

Not collective.

See also

TSGetI2Function

Source code at petsc4py/PETSc/TS.pyx:1201

Return type:

tuple[Vec, TSI2Function]

getI2Jacobian()#

Return the matrices and function which computes the Jacobian.

Not collective.

See also

TSGetI2Jacobian

Source code at petsc4py/PETSc/TS.pyx:1217

Return type:

tuple[Mat, Mat, TSI2Jacobian]

getIFunction()#

Return the vector and function which computes the implicit residual.

Not collective.

See also

TSGetIFunction

Source code at petsc4py/PETSc/TS.pyx:1014

Return type:

tuple[Vec, TSIFunction]

getIJacobian()#

Return the matrices and function which computes the implicit Jacobian.

Not collective.

See also

TSGetIJacobian

Source code at petsc4py/PETSc/TS.pyx:1030

Return type:

tuple[Mat, Mat, TSIJacobian]

getKSP()#

Return the KSP associated with the TS.

Not collective.

See also

TSGetKSP

Source code at petsc4py/PETSc/TS.pyx:1396

Return type:

KSP

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.

Source code at petsc4py/PETSc/TS.pyx:1681

Return type:

int

getMaxSteps()#

Return the maximum number of steps to use.

Not collective.

See also

TSGetMaxSteps

Source code at petsc4py/PETSc/TS.pyx:1650

Return type:

int

getMaxTime()#

Return the maximum (final) time.

Not collective.

Defaults to 5.

See also

TSGetMaxTime

Source code at petsc4py/PETSc/TS.pyx:1614

Return type:

float

getMonitor()#

Return the monitor.

Not collective.

See also

setMonitor

Source code at petsc4py/PETSc/TS.pyx:1962

Return type:

list[tuple[TSMonitorFunction, tuple[Any, …], dict[str, Any]]]

getNumEvents()#

Return the number of events.

Logically collective.

See also

TSGetNumEvents

Source code at petsc4py/PETSc/TS.pyx:2109

Return type:

int

getOptionsPrefix()#

Return the prefix used for all the TS options.

Not collective.

Source code at petsc4py/PETSc/TS.pyx:508

Return type:

str

getPostStep()#

Return the poststep function.

Source code at petsc4py/PETSc/TS.pyx:2202

Return type:

tuple[TSPostStepFunction, tuple[Any, …] | None, dict[str, Any] | None]

getPreStep()#

Return the prestep function.

Not collective.

See also

setPreStep

Source code at petsc4py/PETSc/TS.pyx:2158

Return type:

tuple[TSPreStepFunction, tuple[Any, …] | None, dict[str, Any] | None]

getPrevTime()#

Return the starting time of the previously completed step.

Not collective.

See also

TSGetPrevTime

Source code at petsc4py/PETSc/TS.pyx:1488

Return type:

float

getProblemType()#

Return the type of problem to be solved.

Not collective.

See also

TSGetProblemType

Source code at petsc4py/PETSc/TS.pyx:440

Return type:

ProblemType

getPythonContext()#

Return the instance of the class implementing the required Python methods.

Not collective.

Source code at petsc4py/PETSc/TS.pyx:2720

Return type:

Any

getPythonType()#

Return the fully qualified Python name of the class used by the solver.

Not collective.

Source code at petsc4py/PETSc/TS.pyx:2749

Return type:

str

getQuadratureTS()#

Return the sub TS that evaluates integrals over time.

Not collective.

Returns:

  • forward (bool) – True if evaluating the integral forward in time

  • qts (TS) – The sub TS

Return type:

tuple[bool, TS]

Source code at petsc4py/PETSc/TS.pyx:2533

getRHSFunction()#

Return the vector where the rhs is stored and the function used to compute it.

Not collective.

See also

TSGetRHSFunction

Source code at petsc4py/PETSc/TS.pyx:753

Return type:

tuple[Vec, TSRHSFunction]

getRHSJacobian()#

Return the Jacobian and the function used to compute them.

Not collective.

See also

TSGetRHSJacobian

Source code at petsc4py/PETSc/TS.pyx:769

Return type:

tuple[Mat, Mat, TSRHSJacobian]

getRKType()#

Return the Type.RK scheme.

Not collective.

See also

TSRKGetType

Source code at petsc4py/PETSc/TS.pyx:358

Return type:

str

getSNES()#

Return the SNES associated with the TS.

Not collective.

See also

TSGetSNES

Source code at petsc4py/PETSc/TS.pyx:1381

Return type:

SNES

getSNESFailures()#

Return the total number of failed SNES solves in the TS.

Not collective.

This counter is reset to zero for each successive call to solve.

Source code at petsc4py/PETSc/TS.pyx:1755

Return type:

int

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.

Source code at petsc4py/PETSc/TS.pyx:1664

Return type:

int

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

TSGetSolution

Source code at petsc4py/PETSc/TS.pyx:1252

Return type:

Vec

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

TS2GetSolution

Source code at petsc4py/PETSc/TS.pyx:1290

Return type:

tuple[Vec, Vec]

getSolveTime()#

Return the time after a call to solve.

Not collective.

This time corresponds to the final time set with setMaxTime.

See also

TSGetSolveTime

Source code at petsc4py/PETSc/TS.pyx:1502

Return type:

float

getStepLimits()#

Return the minimum and maximum allowed time step sizes.

Not collective.

Source code at petsc4py/PETSc/TS.pyx:2343

Return type:

tuple[float, float]

getStepNumber()#

Return the number of time steps completed.

Not collective.

See also

TSGetStepNumber

Source code at petsc4py/PETSc/TS.pyx:1578

Return type:

int

getStepRejections()#

Return the total number of rejected steps.

Not collective.

This counter is reset to zero for each successive call to solve.

Source code at petsc4py/PETSc/TS.pyx:1720

Return type:

int

getTheta()#

Return the abscissa of the stage in (0, 1] for Type.THETA.

Not collective.

See also

TSThetaGetTheta

Source code at petsc4py/PETSc/TS.pyx:2787

Return type:

float

getThetaEndpoint()#

Return whether the endpoint variable of Type.THETA is used.

Not collective.

Source code at petsc4py/PETSc/TS.pyx:2819

Return type:

bool

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 or setPostStep), the time returned is at the start of the step.

See also

TSGetTime

Source code at petsc4py/PETSc/TS.pyx:1470

Return type:

float

getTimeSpan()#

Return the time span.

Not collective.

See also

TSGetTimeSpan

Source code at petsc4py/PETSc/TS.pyx:1345

Return type:

ArrayReal

getTimeSpanSolutions()#

Return the solutions at the times in the time span.

Not collective.

Source code at petsc4py/PETSc/TS.pyx:1361

Return type:

list[Vec]

getTimeStep()#

Return the duration of the current timestep.

Not collective.

See also

TSGetTimeStep

Source code at petsc4py/PETSc/TS.pyx:1537

Return type:

float

getTolerances()#

Return the tolerances for local truncation error.

Logically collective.

Returns:

  • rtol (float) – the relative tolerance

  • atol (float) – the absolute tolerance

Return type:

tuple[float, float]

See also

TSGetTolerances

Source code at petsc4py/PETSc/TS.pyx:1833

getType()#

Return the TS type.

Not collective.

See also

TSGetType

Source code at petsc4py/PETSc/TS.pyx:344

Return type:

str

interpolate(t, u)#

Interpolate the solution to a given time.

Collective.

Parameters:
  • t (float) – The time to interpolate.

  • u (Vec) – The state vector to interpolate.

Return type:

None

See also

TSInterpolate

Source code at petsc4py/PETSc/TS.pyx:2300

load(viewer)#

Load a TS that has been stored in binary with view.

Collective.

Parameters:

viewer (Viewer) – The visualization context.

Return type:

None

See also

TSLoad

Source code at petsc4py/PETSc/TS.pyx:189

monitor(step, time, u=None)#

Monitor the solve.

Collective.

Parameters:
  • step (int) – The step number that has just completed.

  • time (float) – The model time of the state.

  • u (Vec | None) – The state at the current model time.

Return type:

None

See also

TSMonitor

Source code at petsc4py/PETSc/TS.pyx:1989

monitorCancel()#

Clear all the monitors that have been set.

Logically collective.

See also

TSMonitorCancel

Source code at petsc4py/PETSc/TS.pyx:1974

Return type:

None

removeTrajectory()#

Remove the internal TS trajectory object.

Collective.

Source code at petsc4py/PETSc/TS.pyx:2381

Return type:

None

reset()#

Reset the TS, removing any allocated vectors and matrices.

Collective.

See also

TSReset

Source code at petsc4py/PETSc/TS.pyx:2218

Return type:

None

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

TSRestartStep

Source code at petsc4py/PETSc/TS.pyx:2246

Return type:

None

rollBack()#

Roll back one time step.

Collective.

See also

TSRollBack

Source code at petsc4py/PETSc/TS.pyx:2269

Return type:

None

setARKIMEXFullyImplicit(flag)#

Solve both parts of the equation implicitly.

Logically collective.

Parameters:

flag (bool) – Set to True for fully implicit.

Return type:

None

Source code at petsc4py/PETSc/TS.pyx:326

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:

None

Notes

-ts_arkimex_type sets scheme type from the commandline.

See also

TSARKIMEXSetType

Source code at petsc4py/PETSc/TS.pyx:303

setAlphaParams(alpha_m=None, alpha_f=None, gamma=None)#

Set the algorithmic parameters for Type.ALPHA.

Logically collective.

Users should call setAlphaRadius.

Parameters:
  • alpha_m (float | None) – Parameter, leave None to keep current value.

  • alpha_f (float | None) – Parameter, leave None to keep current value.

  • gamma (float | None) – Parameter, leave None to keep current value.

Return type:

None

See also

TSAlphaSetParams

Source code at petsc4py/PETSc/TS.pyx:2857

setAlphaRadius(radius)#

Set the spectral radius for Type.ALPHA.

Logically collective.

Parameters:

radius (float) – the spectral radius

Return type:

None

Notes

-ts_alpha_radius can be used to set this from the commandline.

See also

TSAlphaSetRadius

Source code at petsc4py/PETSc/TS.pyx:2835

setAppCtx(appctx)#

Set the application context.

Not collective.

Parameters:

appctx (Any) – The application context.

Return type:

None

Source code at petsc4py/PETSc/TS.pyx:559

setConvergedReason(reason)#

Set the reason for handling the convergence of solve.

Logically collective.

Can only be called when solve is active and reason must contain common value.

Parameters:

reason (ConvergedReason) – The reason for convergence.

Return type:

None

Source code at petsc4py/PETSc/TS.pyx:1889

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:

None

Source code at petsc4py/PETSc/TS.pyx:2408

setDIRKType(ts_type)#

Set the type of Type.DIRK scheme.

Logically collective.

Parameters:

ts_type (DIRKType | str) – The type of Type.DIRK scheme.

Return type:

None

Notes

-ts_dirk_type sets scheme type from the commandline.

See also

TSDIRKSetType

Source code at petsc4py/PETSc/TS.pyx:386

setDM(dm)#

Set the DM that may be used by some nonlinear solvers or preconditioners.

Logically collective.

Parameters:

dm (DM) – The DM object.

Return type:

None

See also

TSSetDM

Source code at petsc4py/PETSc/TS.pyx:1433

setEquationType(eqtype)#

Set the type of the equation that TS is solving.

Logically collective.

Parameters:

eqtype (EquationType) – The type of equation.

Return type:

None

Source code at petsc4py/PETSc/TS.pyx:454

setErrorIfStepFails(flag=True)#

Immediately error is no step succeeds.

Not collective.

Parameters:

flag (bool) – Enable to error if no step succeeds.

Return type:

None

Notes

-ts_error_if_step_fails to enable from the commandline.

Source code at petsc4py/PETSc/TS.pyx:1772

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:

None

Source code at petsc4py/PETSc/TS.pyx:2018

setEventTolerances(tol=None, vtol=None)#

Set tolerances for event zero crossings when using event handler.

Logically collective.

setEventHandler must have already been called.

Parameters:
  • tol (float) – The scalar tolerance or None to leave at the current value

  • vtol (Sequence[float]) – A sequence of scalar tolerance for each event. Used in preference to tol if present. Set to None to leave at the current value.

Return type:

None

Notes

-ts_event_tol can be used to set values from the commandline.

Source code at petsc4py/PETSc/TS.pyx:2073

setExactFinalTime(option)#

Set method of computing the final time step.

Logically collective.

Parameters:

option (ExactFinalTime) – The exact final time option

Return type:

None

Notes

-ts_exact_final_time may be used to specify from the commandline.

Source code at petsc4py/PETSc/TS.pyx:1867

setFromOptions()#

Set various TS parameters from user options.

Collective.

Source code at petsc4py/PETSc/TS.pyx:545

Return type:

None

setI2Function(function, f=None, args=None, kargs=None)#

Set the function to compute the 2nd order DAE.

Logically collective.

Parameters:
  • function (TSI2Function | None) – The right-hand side function.

  • f (Vec | None) – The vector to store values or None to be created internally.

  • args (tuple[Any, ...] | None) – Additional positional arguments for function.

  • kargs (dict[str, Any] | None) – Additional keyword arguments for function.

Return type:

None

See also

TSSetI2Function

Source code at petsc4py/PETSc/TS.pyx:1046

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:

None

See also

TSSetI2Jacobian

Source code at petsc4py/PETSc/TS.pyx:1083

setIFunction(function, f=None, args=None, kargs=None)#

Set the function representing the DAE to be solved.

Logically collective.

Parameters:
  • function (TSIFunction | None) – The right-hand side function.

  • f (Vec | None) – The vector to store values or None to be created internally.

  • args (tuple[Any, ...] | None) – Additional positional arguments for function.

  • kargs (dict[str, Any] | None) – Additional keyword arguments for function.

Return type:

None

See also

TSSetIFunction

Source code at petsc4py/PETSc/TS.pyx:787

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 where F(t, U, U_t) is the function provided with setIFunction.

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:

None

See also

TSSetIJacobian

Source code at petsc4py/PETSc/TS.pyx:824

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 parameters P where F(Udot, U, t) = G(U, P, t), as well as the location to store the matrix.

Parameters:
  • jacobian – The function which computes the Jacobian.

  • J (Mat | None) – The matrix into which the Jacobian is computed.

  • args (tuple[Any, ...] | None) – Additional positional arguments for jacobian.

  • kargs (dict[str, Any] | None) – Additional keyword arguments for jacobian.

Return type:

None

See also

TSSetIJacobianP

Source code at petsc4py/PETSc/TS.pyx:869

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:

None

Source code at petsc4py/PETSc/TS.pyx:1737

setMaxStepRejections(n)#

Set the maximum number of step rejections before a time step fails.

Not collective.

Parameters:

n (int) – The maximum number of rejected steps, use -1 for unlimited.

Return type:

None

Notes

-ts_max_reject can be used to set this from the commandline

Source code at petsc4py/PETSc/TS.pyx:1698

setMaxSteps(max_steps)#

Set the maximum number of steps to use.

Logically collective.

Defaults to 5000.

Parameters:

max_steps (int) – The maximum number of steps to use.

Return type:

None

See also

TSSetMaxSteps

Source code at petsc4py/PETSc/TS.pyx:1630

setMaxTime(max_time)#

Set the maximum (final) time.

Logically collective.

Parameters:

max_time (float) – the final time

Return type:

None

Notes

-ts_max_time sets the max time from the commandline

See also

TSSetMaxTime

Source code at petsc4py/PETSc/TS.pyx:1592

setMonitor(monitor, args=None, kargs=None)#

Set an additional monitor to the TS.

Logically collective.

Parameters:
  • monitor (TSMonitorFunction | None) – The custom monitor function.

  • args (tuple[Any, ...] | None) – Additional positional arguments for monitor.

  • kargs (dict[str, Any] | None) – Additional keyword arguments for monitor.

Return type:

None

See also

TSMonitorSet

Source code at petsc4py/PETSc/TS.pyx:1928

setOptionsPrefix(prefix)#

Set the prefix used for all the TS options.

Logically collective.

Parameters:

prefix (str | None) – The prefix to prepend to all option names.

Return type:

None

Notes

A hyphen must not be given at the beginning of the prefix name.

Source code at petsc4py/PETSc/TS.pyx:485

setPostStep(poststep, args=None, kargs=None)#

Set a function to be called at the end of each time step.

Logically collective.

Parameters:
  • poststep (TSPostStepFunction | None) – The function to be called at the end of each step.

  • args (tuple[Any, ...] | None) – Additional positional arguments for poststep.

  • kargs (dict[str, Any] | None) – Additional keyword arguments for poststep.

Return type:

None

See also

TSSetPostStep

Source code at petsc4py/PETSc/TS.pyx:2170

setPreStep(prestep, args=None, kargs=None)#

Set a function to be called at the beginning of each time step.

Logically collective.

Parameters:
  • prestep (TSPreStepFunction | None) – The function to be called at the beginning of each step.

  • args (tuple[Any, ...] | None) – Additional positional arguments for prestep.

  • kargs (dict[str, Any] | None) – Additional keyword arguments for prestep.

Return type:

None

See also

TSSetPreStep

Source code at petsc4py/PETSc/TS.pyx:2125

setProblemType(ptype)#

Set the type of problem to be solved.

Logically collective.

Parameters:

ptype (ProblemType) – The type of problem of the forms.

Return type:

None

See also

TSSetProblemType

Source code at petsc4py/PETSc/TS.pyx:423

setPythonContext(context)#

Set the instance of the class implementing the required Python methods.

Not collective.

Source code at petsc4py/PETSc/TS.pyx:2708

Parameters:

context (Any)

Return type:

None

setPythonType(py_type)#

Set the fully qualified Python name of the class to be used.

Collective.

Source code at petsc4py/PETSc/TS.pyx:2735

Parameters:

py_type (str)

Return type:

None

setRHSFunction(function, f=None, args=None, kargs=None)#

Set the routine for evaluating the function G in U_t = G(t, u).

Logically collective.

Parameters:
  • function (TSRHSFunction | None) – The right-hand side function.

  • f (Vec | None) – The vector into which the right-hand side is computed.

  • args (tuple[Any, ...] | None) – Additional positional arguments for function.

  • kargs (dict[str, Any] | None) – Additional keyword arguments for function.

Return type:

None

See also

TSSetRHSFunction

Source code at petsc4py/PETSc/TS.pyx:578

setRHSJacobian(jacobian, J=None, P=None, args=None, kargs=None)#

Set the function to compute the Jacobian of G in U_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:

None

See also

TSSetRHSJacobian

Source code at petsc4py/PETSc/TS.pyx:615

setRHSJacobianP(jacobianp, A=None, args=None, kargs=None)#

Set the function that computes the Jacobian with respect to the parameters.

Logically collective.

Parameters:
  • jacobianp (TSRHSJacobianP | None) – The user-defined function.

  • A (Mat | None) – The matrix into which the Jacobian will be computed.

  • args (tuple[Any, ...] | None) – Additional positional arguments for jacobianp.

  • kargs (dict[str, Any] | None) – Additional keyword arguments for jacobianp.

Return type:

None

Source code at petsc4py/PETSc/TS.pyx:2475

setRKType(ts_type)#

Set the type of the Runge-Kutta scheme.

Logically collective.

Parameters:

ts_type (RKType | str) – The type of scheme.

Return type:

None

Notes

-ts_rk_type sets scheme type from the commandline.

See also

TSRKSetType

Source code at petsc4py/PETSc/TS.pyx:280

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.

Source code at petsc4py/PETSc/TS.pyx:2362

Return type:

None

setSolution(u)#

Set the initial solution vector.

Logically collective.

Parameters:

u (Vec) – The solution vector.

Return type:

None

See also

TSSetSolution

Source code at petsc4py/PETSc/TS.pyx:1235

setSolution2(u, v)#

Set the initial solution and its time derivative.

Logically collective.

Parameters:
  • u (Vec) – The solution vector.

  • v (Vec) – The time derivative vector.

Return type:

None

See also

TS2SetSolution

Source code at petsc4py/PETSc/TS.pyx:1271

setStepLimits(hmin, hmax)#

Set the minimum and maximum allowed step sizes.

Logically collective.

Parameters:
  • hmin (float) – the minimum step size

  • hmax (float) – the maximum step size

Return type:

None

Source code at petsc4py/PETSc/TS.pyx:2320

setStepNumber(step_number)#

Set the number of steps completed.

Logically collective.

For most uses of the TS solvers the user need not explicitly call setStepNumber, as the step counter is appropriately updated in solve/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 that TS monitors will be called with a initial nonzero step counter.

Parameters:

step_number (int) – the number of steps completed

Return type:

None

See also

TSSetStepNumber

Source code at petsc4py/PETSc/TS.pyx:1551

setTheta(theta)#

Set the abscissa of the stage in (0, 1] for Type.THETA.

Logically collective.

Parameters:

theta (float) – stage abscissa

Return type:

None

Notes

-ts_theta_theta can be used to set a value from the commandline.

See also

TSThetaSetTheta

Source code at petsc4py/PETSc/TS.pyx:2765

setThetaEndpoint(flag=True)#

Set to use the endpoint variant of Type.THETA.

Logically collective.

Parameters:

flag – Enable to use the endpoint variant.

Return type:

None

Source code at petsc4py/PETSc/TS.pyx:2801

setTime(t)#

Set the time.

Logically collective.

Parameters:

t (float) – The time.

Return type:

None

See also

TSSetTime

Source code at petsc4py/PETSc/TS.pyx:1452

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 with getTimeSpanSolutions.

Parameters:

tspan (Sequence[float]) – The sequence of time points.

Return type:

None

Notes

-ts_time_span <t0, ..., tf> sets the time span from the commandline

See also

TSSetTimeSpan

Source code at petsc4py/PETSc/TS.pyx:1313

setTimeStep(time_step)#

Set the duration of the timestep.

Logically collective.

Parameters:

time_step (float) – the duration of the timestep

Return type:

None

See also

TSSetTimeStep

Source code at petsc4py/PETSc/TS.pyx:1519

setTolerances(rtol=None, atol=None)#

Set tolerances for local truncation error when using an adaptive controller.

Logically collective.

Parameters:
  • rtol (float) – The relative tolerance or None to leave the current value.

  • atol (float) – The absolute tolerance or None to leave the current value.

Return type:

None

Notes

-ts_rtol and -ts_atol may be used to set values from the commandline.

See also

TSSetTolerances

Source code at petsc4py/PETSc/TS.pyx:1794

setType(ts_type)#

Set the method to be used as the TS solver.

Collective.

Parameters:

ts_type (Type | str) – The solver type.

Return type:

None

Notes

-ts_type sets the method from the commandline

See also

TSSetType

Source code at petsc4py/PETSc/TS.pyx:257

setUp()#

Set up the internal data structures for the TS.

Collective.

See also

TSSetUp

Source code at petsc4py/PETSc/TS.pyx:2206

Return type:

None

solve(u=None)#

Step the requested number of timesteps.

Collective.

Parameters:

u (Vec | None) – The solution vector. Can be None.

Return type:

None

See also

TSSolve

Source code at petsc4py/PETSc/TS.pyx:2281

step()#

Take one step.

Collective.

The preferred interface for the TS solvers is solve. If you need to execute code at the beginning or ending of each step, use setPreStep and setPostStep respectively.

See also

TSStep

Source code at petsc4py/PETSc/TS.pyx:2230

Return type:

None

view(viewer=None)#

Print the TS object.

Collective.

Parameters:

viewer (Viewer | None) – The visualization context.

Return type:

None

Notes

-ts_view calls TSView at the end of TSStep

See also

TSView

Source code at petsc4py/PETSc/TS.pyx:166

Attributes Documentation

appctx#

Application context.

Source code at petsc4py/PETSc/TS.pyx:2906

atol#

The absolute tolerance.

Source code at petsc4py/PETSc/TS.pyx:3009

converged#

Indicates the TS has converged.

Source code at petsc4py/PETSc/TS.pyx:3030

diverged#

Indicates the TS has stopped.

Source code at petsc4py/PETSc/TS.pyx:3035

dm#

The DM.

Source code at petsc4py/PETSc/TS.pyx:2916

equation_type#

The equation type.

Source code at petsc4py/PETSc/TS.pyx:2934

iterating#

Indicates the TS is still iterating.

Source code at petsc4py/PETSc/TS.pyx:3025

ksp#

The KSP.

Source code at petsc4py/PETSc/TS.pyx:2947

max_steps#

The maximum number of steps.

Source code at petsc4py/PETSc/TS.pyx:2991

max_time#

The maximum time.

Source code at petsc4py/PETSc/TS.pyx:2983

problem_type#

The problem type.

Source code at petsc4py/PETSc/TS.pyx:2926

reason#

The converged reason.

Source code at petsc4py/PETSc/TS.pyx:3017

rtol#

The relative tolerance.

Source code at petsc4py/PETSc/TS.pyx:3001

snes#

The SNES.

Source code at petsc4py/PETSc/TS.pyx:2942

step_number#

The current step number.

Source code at petsc4py/PETSc/TS.pyx:2975

time#

The current time.

Source code at petsc4py/PETSc/TS.pyx:2959

time_step#

The current time step size.

Source code at petsc4py/PETSc/TS.pyx:2967

vec_sol#

The solution vector.

Source code at petsc4py/PETSc/TS.pyx:2952