Internal API

Note that the internal API may change without being considered a breaking change!

FerriteAssembly.reinit_buffer!Function
reinit_buffer!(cb::CellBuffer, db::AbstractDomainBuffer, cellnum::Int; a=nothing, aold=nothing)

Reinitialize the cb::CellBuffer for cell number cellnum. The global degree of freedom vectors a (current) and aold are used to update the cell degree of freedom vectors in c. If the global vectors are instead ::Nothing, the corresponding cell values are set to NaN The element stiffness, cb.Ke, and residual, cb.re, are also zeroed.

source
FerriteAssembly.create_statesFunction
create_states(sdh::SubDofHandler, material, cellvalues, a, cellset, dofrange)

Create a Dict of states for the cells in cellset, where the user should define the create_cell_state function for their material (and corresponding cellvalues) dofrange::NamedTuple is passed onto create_cell_state and contains the local dof ranges for each field.

source
FerriteAssembly._copydofs!Function
_copydofs!(edofs::Vector, gdofs::Vector, inds::Vector{Int})

Internal function for faster copying of global values into the element values. Equivalent to edofs .= gdofs[inds]

_copydofs!(edofs::Vector, gdofs::Nothing, inds::Vector{Int})

Fill edofs with NaN

source
FerriteAssembly._create_cell_stateFunction
_create_cell_state(cell, material, cellvalues, a, ae, dofrange, cellnr)

Internal function to reinit and extract the relevant quantities from the cell::CellCache, reinit cellvalues, update ae from a, and pass these into the create_cell_state function that the user should specify.

source
FerriteAssembly.CellBufferType
CellBuffer(
    numdofs::Int, numnodes::Int, ::Val{sdim}, 
    cellvalues, material, state, dofrange, user_data=nothing) -> CellBuffer

Create a cell cache for an element with numdofs degrees of freedom and numnodes nodes with dimension sdim. cellvalues are reinit!ed for each cell, and the state is updated to the old cell state. material will be passed as-is to the element. The given dofrange::NamedTuple, user_data::Any, and cache::Any are available to the element via the buffer input.

See setup_domainbuffer

This constructor is normally not used, and is instead called from setup_domainbuffer

source
FerriteAssembly.FacetBufferType
FacetBuffer(
    numdofs::Int, numnodes::Int, ::Val{sdim}, 
    facetvalues, material, state, dofrange, user_data=nothing)

Create a facetbuffer for an element with numdofs degrees of freedom and numnodes nodes with dimension sdim. facetvalues are reinit!ed for each cell, and the state is updated to the old cell state. material will be passed as-is to the element. The given dofrange::NamedTuple, user_data::Any, and cache::Any are available to the element via the buffer input.

See setup_domainbuffer

This constructor is normally not used, and is instead called from setup_domainbuffer

source
FerriteAssembly.get_itembufferFunction
get_itembuffer(dbs::Dict{String,AbstractDomainBuffer}, domain::String)
get_itembuffer(db::AbstractDomainBuffer)

Get the AbstractItemBuffer stored in db or dbs[domain]. This internal function might change, but currently the full TaskLocals is returned for a ThreadedDomainBuffer (used internally).

source
FerriteAssembly.skip_this_domainFunction
skip_this_domain(worker, name::String)

Should the domain with key name be skipped during work? Defaults to false. Can be used to e.g. only loop over parts of a domain.

source
FerriteAssembly.can_threadFunction
can_thread(worker)::Bool

Does the worker support multithreaded work? Defaults to false. If this returns true, the worker must support the TaskLocals interface.

source
FerriteAssembly.fast_getindexFunction
fast_getindex(collection)

If the output is known from the element type of the collection, this can be returned directly. Currently, this is implemented for AbstractDict with Nothing as type, which takes away overhead when state variables are not used.

source
FerriteAssembly.autogenerate_cellvaluesFunction
autogenerate_cellvalues(cv::CellValues, args...)

Just return the provided CellValues

autogenerate_cellvalues(order::Int, ip_fun::Interpolation{RefShape}, ip_geo::Interpolation{RefShape})

Using quadrature rule, qr = QuadratureRule{RefShape}(order), return CellValues(qr, ip_fun, ip_geo)

source
FerriteAssembly.autogenerate_facetvaluesFunction
autogenerate_facetvalues(fv::FacetValues, args...)

Just return the provided FacetValues

autogenerate_facetvalues(order::Int, ip_fun::Interpolation{RefShape}, ip_geo::Interpolation{RefShape})

Using quadrature rule, fqr = FacetQuadratureRule{RefShape}(order), create FacetValues(fqr, ip_fun, ip_geo)

source

Threading model

FerriteAssembly.TaskChunksType
TaskChunks(chunks)

This is strictly an internal implementation detail for now. It is used to provide a variant of Base.Channel, but where data is not removed, the index is just incremented.

source