API

The main types below are exported. Remaining functions are not exported to avoid polluting the name space. Tip: To simplify calling the following functions it is possible to write import FerriteProblems as FP as is done in the examples.

Main types

FerriteProblems.FerriteProblemType
FerriteProblem(def::FEDefinition, post=nothing, io=nothing; kwargs...)
FerriteProblem(def::FEDefinition, post, savefolder::String; kwargs...)

Create a FerriteProblem from def. Postprocessing can be added as post, see FESolvers.postprocess!. File input/output using FerriteIO can be added with io. It is also possible to just give the savefolder, i.e. where to save the output when the default FerriteIO can be used.

Supported keyword arguments are

  • autodiffbuffer::Bool: Should FerriteAssembly.jl's AutoDiffCellBuffer be used? This will make the assembly faster if automatic differentiation is used, and can also be used without automatic differentiation (but with a slight extra computational overhead)
  • threading::Bool: Should threading be used?
source
FerriteProblems.FEDefinitionType
FEDefinition(
    domainspec::Union{DomainSpec,Dict{String,DomainSpec}};
    ch, lh=LoadHandler(dh),
    convergence_criterion=AbsoluteResidual(), 
    initial_conditions=NamedTuple(),
    )

Create the FEDefinition which can later be used to create a complete FerriteProblem. domainspec is the FerriteAssembly domain specification.

In addition, the following keyword arguments can be given

  • ch: The constraint handler, ConstraintHandler (Ferrite.jl)
  • lh: The external load handler, LoadHandler (FerriteAssembly.jl)
  • initial_conditions: NamedTuple with a function f(x) for each field that has a nonzero initial condition. Used by the Ferrite.jl's apply_analytical! function. Example: initial_conditions = (u = x -> Vec((x[1]/10, 0.0)), p = x -> -100*x[2]). For fields not given here, the initial condition is zeros everywhere.
  • convergence_criterion: Determines how to calculate the convergence measure including scaling. See ConvergenceCriterion
source
FerriteProblems.FerriteIOType
FerriteIO(
    folder::String, def::FEDefinition, post=nothing; 
    def_file="FEDefinition.jld2", 
    postfile="FEPost.jld2",
    T=Float64, 
    nsteps_per_file=typemax(Int), 
    switchsize=Inf
    )

Constructor for creating a FerriteIO when simulating.

source
FerriteIO(filepath::String)

Constructor for reading a FerriteIO that was saved during a simulation. The default filename is FerriteIO.jld2, located in the savefolder given to the FerriteProblem constructor when setting up the simulation.

Can be called either with the do-block (recommended),

FerriteIO(filepath) do io
    # Do whatever with io 
end

or by manually closing,

io = FerriteIO(filepath)
# Do whatever with io
close(io)
source

Convergence criteria

In normal usage, the following convergence criteria can be used

FerriteProblems.AbsoluteResidualType
AbsoluteResidual()

The default convergence criterion that calculates the convergence measure as √(sum([r[i]^2 for i ∈ free dofs]) without any scaling.

source
FerriteProblems.RelativeResidualElementScalingType
RelativeResidualElementScaling(;p = 2, minfactors::Union{AbstractFloat,NamedTuple}=eps())

Use Ferriteassembly.ElementResidualScaling with the exponent p to calculate the scaling for each field individually, based on the Lp-norm of each cell's residual. To avoid issues when all cells have zero residual (e.g. in the first time step), supply minfactors as the minimum scaling factor. The convergence measure is calculated with the following pseudo-code

val = 0.0
for field in Ferrite.getfieldnames(dh)
    factor = max(element_residual_scaling[field], minfactors[field])
    dofs = free_field_dofs[field]   # Get the non-constrained dofs for `field`
    val += sum(i->(r[i]/factor)^2, dofs)
end
return √val

where the same minfactor is used for all fields if only a scalar value is given.

source

To create custom convergence criteria, the following functions may require overloading.

FerriteProblems.TolScalingType
TolScaling(criterion::ConvergenceCriterion, def::FEDefinition)

The TolScaling type contains the criterion that determines how to scale the residuals to determine convergence. The constructor is specialized on typeof(criterion), creating the following fields:

  • assemscaling: scaling to be used by FerriteAssembly to give potential scaling contribution based on each element's residual
  • buffer: Used to pre-calculate values, such as dof-ranges for each field that is used when calculating the convergence measure.
source
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
FerriteProblems.make_assemscalingFunction
make_assemscaling(criterion, def)

Create the scaling for use in FerriteAssembly if required by the given criterion. The default, if not overloaded, returns nothing.

source

Access functions

FerriteAssembly.get_dofhandlerMethod
FerriteProblems.get_dofhandler(p::FerriteProblem)

Get dh::Ferrite.AbstractDofHandler from the FerriteProblem (Technically overloaded from FerriteAssembly, but accessible via FerriteProblems)

source
FerriteAssembly.get_materialMethod
FerriteProblems.get_material(p::FerriteProblem)
FerriteProblems.get_material(p::FerriteProblem, domain_name::String)

Get the material in p. For multiple domains, it is necessary to give the domain_name for where to get the material. Note that this is type-unstable and should be avoided in performance-critical code sections. This function belongs to FerriteAssembly.jl, but can be accessed via FerriteProblems.get_material.

source
FESolvers.getjacobianMethod
FerriteProblems.getjacobian(p::FerriteProblem)

Get the current jacobian matrix from p. Note that this function belongs to FESolvers.jl, but can be accessed via FerriteProblems.getjacobian

source
FESolvers.getunknownsMethod
FerriteProblems.getunknowns(p::FerriteProblem)

Get the current vector of unknowns from p. Note that this function belongs to FESolvers.jl, but can be accessed via FerriteProblems.getunknowns

source
FESolvers.getresidualMethod
FerriteProblems.getresidual(p::FerriteProblem)

Get the current residual vector from p. Note that this function belongs to FESolvers.jl, but can be accessed via FerriteProblems.getresidual

source
FerriteProblems.get_external_forceMethod
FerriteProblems.getneumannforce(p::FerriteProblem)

Get the current external force vector caused by Neumann boundary conditions. Note that this vector does not include external forces added during the cell assembly; only forces added with the NeumannHandler

source

Saving and loading data

FESolvers.handle_notconverged!Method
FESolvers.handle_notconverged!(post, p::FerriteProblem, solver)

Optional overload which is called when the problem doesn't converge for the current attempt. Allows for example to modify the problem or add special postprocessing to investigate convergence issues.

Note

If the problem is modified, and an adaptive time stepper is used, the adaptive time stepper will still consider the problem as not converged, and adapt the time-stepping accordingly. Currently, there is no interface to prevent this, and usage with a fixed time stepping might make more sense. However, by modifying the solver's state, it is possible to trick the adaptive time stepper to not modify the time step, but this requires using non-stable and non-public API.

source
FESolvers.postprocess!Method
FESolvers.postprocess!(p::FerriteProblem, solver)

When FESolvers call this function for p::FerriteProblem, the following function

FESolvers.postprocess!(post, p::FerriteProblem, solver)

is called where post=p.post (unless you define a different override). This allows you to easily define the dispatch on your postprocessing type as FESolvers.postprocess!(post::MyPostType, p, solver)

source
FerriteProblems.close_postprocessingFunction
close_postprocessing(post::MyPostType, p::FerriteProblem)

This function is called to close any open files manually created during the postprocessing with the custom postprocessing type MyPostType. Note that the file streams in p.io::FerriteIO are automatically closed and don't require any special handling.

source
FerriteProblems.addstep!Function
FerriteProblems.addstep!(io::FerriteIO, p::FerriteProblem)

Add a new step to be saved by io at the time gettime(p) Must be called before adding any new data

source
FerriteProblems.savedofdata!Function
FerriteProblems.savedofdata!(io::FerriteIO, vals, dt_order=0, field="dof")

Save data pertaining to each degree of freedom. Use a different field than "dof"` to save data located at each dof, but not the actual dof values (e.g. the residual vector)

source
FerriteProblems.savenodedata!Function
FerriteProblems.savenodedata!(io::FerriteIO, vals, field, dt_order=0)

Save data located at each node. By convention this should be indexed by the node numbers in the grid. (E.g. a Vector for all nodes or a Dict{Int} with keys the node numbers)

source
FerriteProblems.savecelldata!Function
FerriteProblems.savecelldata!(io::FerriteIO, vals, field, dt_order=0)

Save data for each cell. By convention this should be indexed by the cell numbers in the grid. (E.g. a Vector for all cells or a Dict{Int} with keys the cell indices)

source
FerriteProblems.saveipdata!Function
FerriteProblems.saveipdata!(io::FerriteIO, vals, field, dt_order=0)

Save data for each integration point in cells in the grid. By convention the data for each cell should be indexed by the cell numbers in the grid. (E.g. a Vector for all cells or a Dict{Int} with keys the cell indices) Note that it is on the user to know how the integration points are numbered, i.e. which QuadratureRule that was used.

source
FerriteProblems.saveglobaldata!Function
FerriteProblems.saveglobaldata!(io::FerriteIO, vals, field, dt_order=0)

Save data that is global to the entire simulation, i.e. global quantites such as reaction forces, total dissipation, etc.

source