cashocs.nonlinear_solvers.ts.TSPseudoSolver#

class cashocs.nonlinear_solvers.ts.TSPseudoSolver(nonlinear_form: ufl.Form, u: fenics.Function, bcs: fenics.DirichletBC | list[fenics.DirichletBC], derivative: ufl.Form | None = None, petsc_options: _typing.KspOption | None = None, shift: ufl.Form | None = None, rtol: float | None = None, atol: float | None = None, max_iter: int | None = None, A_tensor: fenics.PETScMatrix | None = None, b_tensor: fenics.PETScVector | None = None, preconditioner_form: ufl.Form | None = None, excluded_from_time_derivative: list[int] | None = None)[source]#

Bases: object

Interface for using PETSc’s TS solver for pseudo time stepping.

This class implements only pseudo time stepping for nonlinear equations, it should not be used for transient simulations.

Initializes the TS pseudo time stepping solver.

Parameters:
  • nonlinear_form (ufl.Form) – The variational form of the nonlinear problem to be solved by Newton’s method.

  • u (fenics.Function) – The sought solution / initial guess. It is not assumed that the initial guess satisfies the Dirichlet boundary conditions, they are applied automatically. The method overwrites / updates this Function.

  • bcs (fenics.DirichletBC | list[fenics.DirichletBC]) – A list of DirichletBCs for the nonlinear variational problem.

  • derivative (ufl.Form | None, optional) – The Jacobian of nonlinear_form, used for the Newton method. Default is None, and in this case the Jacobian is computed automatically with AD.

  • petsc_options (_typing.KspOption | None, optional) – The options for PETSc TS object. Defaults to None.

  • shift (ufl.Form | None, optional) – A shift term, if the right-hand side of the nonlinear problem is not zero, but shift. Defaults to None.

  • rtol (float | None, optional) – Relative tolerance of the solver. If this is set to a float, the float is used as relative tolerance. If this is set to None, then the relative tolerance of the SNES object is used, which can be defined with the petsc options snes_rtol rtol. Defaults to None.

  • atol (float | None, optional) – Absolute tolerance of the solver. If this is set to a float, the float is used as absolute tolerance. If this is set to None, then the absolute tolerance of the SNES object is used, which can be defined with the petsc options snes_atol atol. Defaults to None.

  • max_iter (int | None, optional) – Maximum number of iterations carried out by the method. Overrides the specification in the petsc_options. Defaults to None.

  • A_tensor (fenics.PETScMatrix | None, optional) – A fenics.PETScMatrix for storing the left-hand side of the linear sub-problem. Defaults to None.

  • b_tensor (fenics.PETScVector | None, optional) – A fenics.PETScVector for storing the right-hand side of the linear sub-problem. Defaults to None.

  • preconditioner_form (ufl.Form | None, optional) – A UFL form which defines the preconditioner matrix. Defaults to None.

  • excluded_from_time_derivative (list[int] | None, optional) – A list of indices for those components that are not included in the time derivative. Example: The pressure for incompressible Navier-Stokes. Default is None, so that all components are included for the time derivative.

Methods Summary

assemble_i_function(ts, t, u, u_dot, f)

Interface for PETSc TSSetIFunction - this describes the time derivative.

assemble_jacobian(ts, t, u, J, P)

Interface for PETSc TSSetRHSJacobian.

assemble_mass_matrix(ts, t, u, u_dot, sigma, ...)

Interface for PETSc TSSetIJacobian.

assemble_residual(ts, t, u, f)

Interface for PETSc TSSetRHSFunction.

compute_nonlinear_residual(u)

Computes the residual of the nonlinear equation.

monitor(ts, i, t, u)

The monitoring function for the pseudo time stepping.

solve()

Solves the (nonlinear) problem with pseudo time stepping.

Methods Documentation

assemble_i_function(ts: petsc4py.PETSc.TS, t: float, u: petsc4py.PETSc.Vec, u_dot: petsc4py.PETSc.Vec, f: petsc4py.PETSc.Vec) None[source]#

Interface for PETSc TSSetIFunction - this describes the time derivative.

Parameters:
  • ts (PETSc.TS) – The TS solver.

  • t (float) – The current time step.

  • u (PETSc.Vec) – The current iterate.

  • u_dot (PETSc.Vec) – The time derivative.

  • f (PETSc.Vec) – The vector for storing the IFunction.

Return type:

None

assemble_jacobian(ts: petsc4py.PETSc.TS, t: float, u: petsc4py.PETSc.Vec, J: petsc4py.PETSc.Mat, P: petsc4py.PETSc.Mat) None[source]#

Interface for PETSc TSSetRHSJacobian.

Parameters:
  • ts (PETSc.TS) – The TS solver.

  • t (float) – The current time step.

  • u (PETSc.Vec) – The current iterate.

  • J (PETSc.Mat) – The matrix for storing the Jacobian.

  • P (PETSc.Mat) – The matrix for storing the preconditioner.

Return type:

None

assemble_mass_matrix(ts: petsc4py.PETSc.TS, t: float, u: petsc4py.PETSc.Vec, u_dot: petsc4py.PETSc.Vec, sigma: float, A: petsc4py.PETSc.Mat, B: petsc4py.PETSc.Mat) None[source]#

Interface for PETSc TSSetIJacobian.

This describes the Jacobian of the time derivative. Here, the IJacobian and its preconditioner are given by the FEM mass matrices.

Parameters:
  • ts (PETSc.TS) – The TS solver.

  • t (float) – The current time step.

  • u (PETSc.Vec) – The current iterate.

  • u_dot (PETSc.Vec) – The time derivative of the current iterate.

  • sigma (float) – See PETSc documentation.

  • A (PETSc.Mat) – The matrix for storing the IJacobian.

  • B (PETSc.Mat) – The matrix for storing the preconditioner of the IJacobian.

Return type:

None

assemble_residual(ts: petsc4py.PETSc.TS, t: float, u: petsc4py.PETSc.Vec, f: petsc4py.PETSc.Vec) None[source]#

Interface for PETSc TSSetRHSFunction.

Parameters:
  • ts (PETSc.TS) – The TS solver.

  • t (float) – The current time step.

  • u (PETSc.Vec) – The current iterate.

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

Return type:

None

compute_nonlinear_residual(u: petsc4py.PETSc.Vec) float[source]#

Computes the residual of the nonlinear equation.

Parameters:

u (PETSc.Vec) – The current iterate.

Returns:

The norm of the nonlinear residual.

Return type:

float

monitor(ts: petsc4py.PETSc.TS, i: int, t: float, u: petsc4py.PETSc.Vec) None[source]#

The monitoring function for the pseudo time stepping.

Parameters:
  • ts (PETSc.TS) – The TS solver.

  • i (int) – The current iteration number.

  • t (float) – The current time step.

  • u (PETSc.Vec) – The current iterate.

Return type:

None

solve() fenics.Function[source]#

Solves the (nonlinear) problem with pseudo time stepping.

Returns:

The solution obtained by the solver.

Return type:

fenics.Function