User problem

The key to using the FESolvers.jl package is to define your problem. This problem should support a set of functions in order for the solver to solve your problem. While some functions are always required, some are only required by certain solvers. Furthermore, a two-level API exist: Simple and advanced. The simple API does not expose which solver is used, while the advanced API requires you to dispatch on the type of solver.

Applicable to all solvers

FESolvers.update_to_next_step!Function
update_to_next_step!(problem, time)

Update prescribed values, external loads etc. for the given time.

This function is called in the beginning of each new time step. Note that for adaptive time stepping, it may be called with a lower time than the previous time if the solution did not converge.

source
FESolvers.update_problem!Function
update_problem!(problem, Δx, update_spec::UpdateSpec)

Update the unknowns of the problem by Δx according to update_spec. Note that

  • Some linear solvers may be inaccurate, and if a modified stiffness is used to enforce constraints on x, it is good the force Δx=0 on these components inside this function.
  • Δx=nothing in the first call after update_to_next_step! in which case, typically, no change of x should be made. Dirichlet boundary conditions are typically updated in update_to_next_step!.

The update_spec gives the information about what and how to update. See the documentation for UpdateSpec for further details. This feature is used by some nonlinear solvers to customize the iteration strategy to speed up or aid convergence. For basic cases when getting started, this can be ignored and a simple function definition would be

FESolvers.update_problem!(problem, Δx, _)
update_problem!(problem, Δx; update_residual::Bool, update_jacobian::Bool)

The old but now deprecated interface is still available without update_spec. The instructions are here:

  • Assemble the residual if update_residual=true
  • Assemble the jacobian if update_jacobian=true
source
FESolvers.handle_converged!Function
handle_converged!(problem)

Do necessary update operations once it is known that the problem has converged. E.g., update old values to the current. Only called after the problem has converged, after postprocess!

source
FESolvers.handle_notconverged!Function
handle_notconverged!(problem, solver)

Optional function to make changes to the problem in case it did not converge. If not implemented, this defaults to a no-op.

source
FESolvers.postprocess!Function
postprocess!(problem, solver)

Perform any postprocessing at the current time. Called at the beginning of the simulation, and directly after time step converged (right before handle_converged!).

source
FESolvers.close_problemFunction
close_problem(problem)

This function is called after solving the problem, even if the solution fails due to an error thrown, for example if the problem doesn't converge. Use to close any file streams etc. that are open and should be closed

source

Simple API

FESolvers.calculate_convergence_measureFunction
calculate_convergence_measure(problem, Δa, iter) -> AbstractFloat

Calculate a value to be compared with the tolerance of the nonlinear solver. A standard case when using Ferrite.jl is norm(getresidual(problem)[Ferrite.free_dofs(ch)]) where ch::Ferrite.ConstraintHandler. Δa is the update of the unknowns from the previous iteration. Note that iter=1 implies Δa=0

source
FESolvers.getjacobianFunction

getjacobian(problem)

Return the jacobian drdx, or approximations thereof.

Must be defined for NewtonSolver, but can also be defined by the advanced API alternative getsystemmatrix: getsystemmatrix(problem, ::NewtonSolver)

source
FESolvers.getdescentpreconditionerFunction

getdescentpreconditioner(problem)

Return a preconditioner K for calculating the descent direction p, considering solving r(x)=0 as a minimization problem of f(x) where r=∇f. The descent direction is then p = K⁻¹ ∇f

Used by the SteepestDescent solver, and defaults to I if not defined. The advanced API alternative is getsystemmatrix: getsystemmatrix(problem, ::SteepestDescent)

source
FESolvers.calculate_energyFunction
calculate_energy(problem,𝐮)

Return the energy of the system (a scalar) which is the integrated energy density over the domain Ω.

source

Advanced API

FESolvers.getsystemmatrixFunction
getsystemmatrix(problem,nlsolver)

Return the system matrix of the problem. For a Newton solver this method should return the Jacobian, while for a steepest descent method this can be a preconditioner as e.g., the L2 product of the gradients. By default the system matrix for the SteepestDescent method is the unity matrix and thus, renders a vanilla gradient descent solver.

source

Additional functions

These functions are usually not necessary to overload

FESolvers.setunknowns!Function
setunknowns!(problem, x)

Copy the given values x into the unknown values of problem. Defaults to copy!(getunknowns(problem), x), which works as long as getunknowns returns the Vector{<:Number} stored in the problem struct. If, e.g. the unknowns is a custom type or a nested vector, this function should be overloaded.

source