dolfin.fem.solving

Simple interface for solving variational problems

A small Python layer on top of the C++ VariationalProblem/Solver classes as well as the solve function.

Functions

solve(*args, **kwargs)

Solve variational problem a == L or F == 0.

dolfin.fem.solving.solve(*args, **kwargs)[source]

Solve variational problem a == L or F == 0.

The following list explains the various ways in which the solve() function can be used.

1. Solving linear variational problems

A linear variational problem a(u, v) = L(v) for all v may be solved by calling solve(a == L, u, …), where a is a bilinear form, L is a linear form, u is a Function (the solution). Optional arguments may be supplied to specify boundary conditions or solver parameters. Some examples are given below:

solve(a == L, u)
solve(a == L, u, bcs=bc)
solve(a == L, u, bcs=[bc1, bc2])

solve(a == L, u, bcs=bcs,
      solver_parameters={"linear_solver": "lu"},
      form_compiler_parameters={"optimize": True})

For available choices for the ‘solver_parameters’ kwarg, look at:

info(LinearVariationalSolver.default_parameters(), True)

2. Solving nonlinear variational problems

A nonlinear variational problem F(u; v) = 0 for all v may be solved by calling solve(F == 0, u, …), where the residual F is a linear form (linear in the test function v but possibly nonlinear in the unknown u) and u is a Function (the solution). Optional arguments may be supplied to specify boundary conditions, the Jacobian form or solver parameters. If the Jacobian is not supplied, it will be computed by automatic differentiation of the residual form. Some examples are given below:

solve(F == 0, u)
solve(F == 0, u, bcs=bc)
solve(F == 0, u, bcs=[bc1, bc2])

solve(F == 0, u, bcs, J=J,
      petsc_options={"ksp_type": "preonly", "pc_type": "lu"},
      form_compiler_parameters={"optimize": True})

For available choices for the ‘solver_parameters’ kwarg, look at:

info(NonlinearVariationalSolver.default_parameters(), True)

4. Solving linear/nonlinear variational problems adaptively

Linear and nonlinear variational problems maybe solved adaptively, with automated goal-oriented error control. The automated error control algorithm is based on adaptive mesh refinement in combination with automated generation of dual-weighted residual-based error estimates and error indicators.

An adaptive solve may be invoked by giving two additional arguments to the solve call, a numerical error tolerance and a goal functional (a Form).

M = u*dx()
tol = 1.e-6

# Linear variational problem
solve(a == L, u, bcs=bc, tol=tol, M=M)

# Nonlinear problem:
solve(F == 0, u, bcs=bc, tol=tol, M=M)