Degrees of Freedom
Degrees of freedom (dofs) are distributed by the DofHandler or the MixedDofHandler.
Ferrite.DofHandler — TypeDofHandler(grid::Grid)Construct a DofHandler based on grid.
Operates slightly faster than MixedDofHandler. Supports:
Grids with a single concrete cell type.- One or several fields on the whole domaine.
Ferrite.MixedDofHandler — TypeMixedDofHandler(grid::Grid)Construct a MixedDofHandler based on grid. Supports:
Grids with or without concrete element type (E.g. "mixed" grids with several different element types.)- One or several fields, which can live on the whole domain or on subsets of the
Grid.
Adding fields to the DofHandlers
Ferrite.add! — Methodadd!(dh::AbstractDofHandler, name::Symbol, dim::Int[, ip::Interpolation])Add a dim-dimensional Field called name which is approximated by ip to dh.
The field is added to all cells of the underlying grid. In case no interpolation ip is given, the default interpolation of the grid's celltype is used. If the grid uses several celltypes, add!(dh::MixedDofHandler, fh::FieldHandler) must be used instead.
Ferrite.add! — Methodadd!(dh::MixedDofHandler, fh::FieldHandler)Add all fields of the FieldHandler fh to dh.
Ferrite.Field — TypeField(name::Symbol, interpolation::Interpolation, dim::Int)Construct dim-dimensional Field called name which is approximated by interpolation.
The interpolation is used for distributing the degrees of freedom.
Ferrite.FieldHandler — TypeFieldHandler(fields::Vector{Field}, cellset::Set{Int})Construct a FieldHandler based on an array of Fields and assigns it a set of cells.
A FieldHandler must fulfill the following requirements:
- All
Cells incellsetare of the same type. - Each field only uses a single interpolation on the
cellset. - Each cell belongs only to a single
FieldHandler, i.e. all fields on a cell must be added within the sameFieldHandler.
Notice that a FieldHandler can hold several fields.
Ferrite.close! — Methodclose!(dh::AbstractDofHandler)Closes dh and creates degrees of freedom for each cell.
If there are several fields, the dofs are added in the following order: For a MixedDofHandler, go through each FieldHandler in the order they were added. For each field in the FieldHandler or in the DofHandler (again, in the order the fields were added), create dofs for the cell. This means that dofs on a particular cell, the dofs will be numbered according to the fields; first dofs for field 1, then field 2, etc.
Dof renumbering
Ferrite.renumber! — Functionrenumber!(dh::AbstractDofHandler, order)
renumber!(dh::AbstractDofHandler, ch::ConstraintHandler, order)Renumber the degrees of freedom in the DofHandler and/or ConstraintHandler according to the ordering order.
order can be given by one of the following options:
- A permutation vector
perm::AbstractVector{Int}such that dofiis renumbered toperm[i]. DofOrder.FieldWise()for renumbering dofs field wise.DofOrder.ComponentWise()for renumbering dofs component wise.DofOrder.Ext{T}for "external" renumber permutations, see documentation forDofOrder.Extfor details.
The dof numbering in the DofHandler and ConstraintHandler must always be consistent. It is therefore necessary to either renumber before creating the ConstraintHandler in the first place, or to renumber the DofHandler and the ConstraintHandler together.
Ferrite.DofOrder.FieldWise — TypeDofOrder.FieldWise()
DofOrder.FieldWise(target_blocks::Vector{Int})Dof order passed to renumber! to renumber global dofs field wise resulting in a globally blocked system.
The default behavior is to group dofs of each field into their own block, with the same order as in the DofHandler. This can be customized by passing a vector of the same length as the total number of fields in the DofHandler (see getfieldnames(dh)) that maps each field to a "target block": to renumber a DofHandler with three fields :u, :v, :w such that dofs for :u and :w end up in the first global block, and dofs for :v in the second global block use DofOrder.FieldWise([1, 2, 1]).
This renumbering is stable such that the original relative ordering of dofs within each target block is maintained.
Ferrite.DofOrder.ComponentWise — TypeDofOrder.ComponentWise()
DofOrder.ComponentWise(target_blocks::Vector{Int})Dof order passed to renumber! to renumber global dofs component wise resulting in a globally blocked system.
The default behavior is to group dofs of each component into their own block, with the same order as in the DofHandler. This can be customized by passing a vector of length ncomponents that maps each component to a "target block" (see DofOrder.FieldWise for details).
This renumbering is stable such that the original relative ordering of dofs within each target block is maintained.
Common methods
Ferrite.ndofs — Functionndofs(dh::AbstractDofHandler)Return the number of degrees of freedom in dh
Ferrite.ndofs_per_cell — Functionndofs_per_cell(dh::AbstractDofHandler[, cell::Int=1])Return the number of degrees of freedom for the cell with index cell.
See also ndofs.
Ferrite.dof_range — Functiondof_range(fh::FieldHandler, field_idx::Int)
dof_range(fh::FieldHandler, field_name::Symbol)
dof_range(dh:MixedDofHandler, field_name::Symbol)Return the local dof range for a given field. The field can be specified by its name or index, where field_idx represents the index of a field within a FieldHandler and field_idxs is a tuple of the FieldHandler-index within the MixedDofHandler and the field_idx.
The dof_range of a field can vary between different FieldHandlers. Therefore, it is advised to use the field_idxs or refer to a given FieldHandler directly in case several FieldHandlers exist. Using the field_name will always refer to the first occurrence of field within the MixedDofHandler.
Example:
julia> grid = generate_grid(Triangle, (3, 3))
Grid{2, Triangle, Float64} with 18 Triangle cells and 16 nodes
julia> dh = MixedDofHandler(grid); add!(dh, :u, 3); add!(dh, :p, 1); close!(dh);
julia> dof_range(dh, :u)
1:9
julia> dof_range(dh, :p)
10:12
julia> dof_range(dh, (1,1)) # field :u
1:9
julia> dof_range(dh.fieldhandlers[1], 2) # field :p
10:12dof_range(dh:DofHandler, field_name)Return the local dof range for field_name. Example:
julia> grid = generate_grid(Triangle, (3, 3))
Grid{2, Triangle, Float64} with 18 Triangle cells and 16 nodes
julia> dh = DofHandler(grid); add!(dh, :u, 3); add!(dh, :p, 1); close!(dh);
julia> dof_range(dh, :u)
1:9
julia> dof_range(dh, :p)
10:12Ferrite.celldofs — Functioncelldofs(dh::AbstractDofHandler, i::Int)Return a vector with the degrees of freedom that belong to cell i.
See also celldofs!.
Ferrite.celldofs! — Functioncelldofs!(global_dofs::Vector{Int}, dh::AbstractDofHandler, i::Int)Store the degrees of freedom that belong to cell i in global_dofs.
See also celldofs.
CellIterator
Ferrite.CellCache — TypeCellCache(grid::Grid)
CellCache(dh::AbstractDofHandler)Create a cache object with pre-allocated memory for the nodes, coordinates, and dofs of a cell. The cache is updated for a new cell by calling reinit!(cache, cellid) where cellid::Int is the cell id.
Struct fields of CellCache
cc.nodes :: Vector{Int}: global node idscc.coords :: Vector{<:Vec}: node coordinatescc.dofs :: Vector{Int}: global dof ids (empty when constructing the cache from a grid)
Methods with CellCache
reinit!(cc, i): reinitialize the cache for cellicellid(cc): get the cell id of the currently cached cellgetnodes(cc): get the global node ids of the cellgetcoordinates(cc): get the coordinates of the cellcelldofs(cc): get the global dof ids of the cellreinit!(fev, cc): reinitializeCellValuesorFaceValues
See also CellIterator.
Ferrite.CellIterator — TypeCellIterator(grid::Grid, cellset=1:getncells(grid))
CellIterator(dh::AbstractDofHandler, cellset=1:getncells(dh))Create a CellIterator to conveniently iterate over all, or a subset, of the cells in a grid. The elements of the iterator are CellCaches which are properly reinit!ialized. See CellCache for more details.
Looping over a CellIterator, i.e.:
for cc in CellIterator(grid, cellset)
# ...
endis thus simply convenience for the following equivalent snippet:
cc = CellCache(grid)
for idx in cellset
reinit!(cc, idx)
# ...
end