Internal API
Note that the internal API may change without being considered a breaking change!
FerriteAssembly.reinit_buffer!
— Functionreinit_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.
FerriteAssembly.create_states
— Functioncreate_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.
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
FerriteAssembly._create_cell_state
— Function_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.
FerriteAssembly.CellBuffer
— TypeCellBuffer(
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.
This constructor is normally not used, and is instead called from setup_domainbuffer
FerriteAssembly.AutoDiffCellBuffer
— TypeAutoDiffCellBuffer(cb::CellBuffer)
FerriteAssembly.FacetBuffer
— TypeFacetBuffer(
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.
This constructor is normally not used, and is instead called from setup_domainbuffer
FerriteAssembly.get_itembuffer
— Functionget_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).
FerriteAssembly.skip_this_domain
— Functionskip_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.
FerriteAssembly.can_thread
— Functioncan_thread(worker)::Bool
Does the worker support multithreaded work? Defaults to false
. If this returns true
, the worker must support the TaskLocals
interface.
FerriteAssembly.fast_getindex
— Functionfast_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.
FerriteAssembly.work_single_cell!
— Functionwork_single_cell!(worker, cellbuffer)
Each worker that supports a cellbuffer should overload this function.
FerriteAssembly.work_single_facet!
— Functionwork_single_facet!(worker, facetbuffer)
Each worker that supports for a facetbuffer should overload this function.
FerriteAssembly.autogenerate_cellvalues
— Functionautogenerate_cellvalues(cv::AbstractCellValues, 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)
FerriteAssembly.autogenerate_facetvalues
— Functionautogenerate_facetvalues(fv::AbstractFacetValues, 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)
Threading model
FerriteAssembly.TaskChunks
— TypeTaskChunks(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.