Assembly#

This module contains classes for the integration and assembly of (weak) integral forms into dense or sparse vectors and matrices. The integration algorithm switches automatically between general cartesion, plane strain or axisymmetric routines, dependent on the given fields.

Core

Take arrays for some pre-defined weak-forms and integrate them into dense or assembly them into sparse vectors or matrices.

IntegralForm(fun, v, dV[, u, grad_v, grad_u])

Form Expressions

Define weak-form expressions on-the-fly for flexible and general form expressions.

Form(v[, u, grad_v, grad_u, dx, args, ...])

A linear or bilinear form object as function decorator on a weak-form with methods for integration and assembly of vectors or sparse matrices.

Create an item out of bilinear and linear weak-form expressions for a Step.

FormItem(bilinearform[, linearform, sym, ...])

An item to be used in a felupe.Step with bilinear and optional linear form objects based on weak-forms with methods for integration and assembly of vectors / sparse matrices.

Detailed API Reference

class felupe.IntegralForm(fun, v, dV, u=None, grad_v=None, grad_u=None)[source]#
assemble(values=None, parallel=False, block=True)[source]#
integrate(parallel=False)[source]#
felupe.Form(v, u=None, grad_v=None, grad_u=None, dx=None, args=None, kwargs=None, parallel=False)#

A linear or bilinear form object as function decorator on a weak-form with methods for integration and assembly of vectors or sparse matrices.

Parameters:
  • v (FieldContainer) – A container for the v fields. May be updated during integration / assembly.

  • u (FieldContainer) – A container for the u fields. May be updated during integration / assembly.

  • grad_v (bool, optional) – Flag to use the gradient of v (default is None).

  • grad_u (bool, optional) – Flag to use the gradient of u (default is None).

  • dx (ndarray or None, optional) – Array with (numerical) differential volumes (default is None).

  • args (tuple or None, optional) – Tuple with initial optional weakform-arguments. May be updated during integration / assembly (default is None).

  • kwargs (dict or None, optional) – Dictionary with initial optional weakform-keyword-arguments. May be updated during integration / assembly (default is None).

Returns:

A form object based on LinearForm, LinearFormExpression, BilinearForm or BilinearFormExpression with methods for integration and assembly.

Return type:

FormExpression

Notes

Linear Form

\[L(v) = \int_\Omega f \cdot v \ dx\]

Bilinear Form

\[a(v, u) = \int_\Omega v \cdot f \cdot u \ dx\]

Examples

FElupe requires a pre-evaluated array for the definition of a bilinear felupe.IntegralForm object on interpolated field values or their gradients. While this has two benefits, namely a fast integration of the form is easy to code and the array may be computed in any programming language, sometimes numeric representations of analytic linear and bilinear form expressions may be easier in user-code and less error prone compared to the calculation of explicit second or fourth-order tensors. Therefore, FElupe provides a function decorator felupe.Form() as an easy-to-use high-level interface, similar to what scikit-fem offers. The felupe.Form() decorator handles a field container. The form class is similar, but not identical in its usage compared to felupe.IntegralForm. It requires a callable function (with optional arguments and keyword arguments) instead of a pre-computed array to be passed. The bilinear form of linear elasticity serves as a reference example for the demonstration on how to use this feature of FElupe. The stiffness matrix is assembled for a unit cube out of hexahedrons.

>>> import felupe as fem
>>> mesh = fem.Cube(n=11)
>>> region = fem.RegionHexahedron(mesh)
>>> displacement = fem.Field(region, dim=3)
>>> field = fem.FieldContainer([displacement])
>>> boundaries, loadcase = fem.dof.uniaxial(field, move=0.5, clamped=True)

The bilinear form of linear elasticity is defined as

\[a(v, u) = \int_\Omega 2 \mu \ \delta\boldsymbol{\varepsilon} : \boldsymbol{\varepsilon} + \lambda \ \text{tr}(\delta\boldsymbol{\varepsilon}) \ \text{tr}(\boldsymbol{\varepsilon}) \ dV\]

with

\[ \begin{align}\begin{aligned}\delta\boldsymbol{\varepsilon} &= \text{sym}(\text{grad}(\boldsymbol{v}))\\\boldsymbol{\varepsilon} &= \text{sym}(\text{grad}(\boldsymbol{u}))\end{aligned}\end{align} \]

and implemented in FElupe closely to the analytic expression. The first two arguments for the callable weak-form function of a bilinear form are always arrays of field (gradients) (v, u) followed by arguments and keyword arguments. Optionally, the integration/assembly may be performed in parallel (threaded). Please note that this is only faster for relatively large systems. The weak-form function is decorated by felupe.Form() where the appropriate fields are linked to v and u along with the gradient flags for both fields. Arguments as well as keyword arguments of the weak-form may be defined inside the decorator or as part of the assembly arguments.

>>> from felupe.math import ddot, trace, sym
>>> @fem.Form(
>>>     v=field, u=field, grad_v=[True], grad_u=[True], kwargs={"μ": 1.0, "λ": 2.0}
>>> )
>>> def bilinearform():
>>>     "A container for a bilinear form."
>>>
>>>     def linear_elasticity(gradv, gradu, μ, λ):
>>>         "Linear elasticity."
>>>
>>>         δε, ε = sym(gradv), sym(gradu)
>>>         return 2 * μ * ddot(δε, ε) + λ * trace(δε) * trace(ε)
>>>
>>>     return [linear_elasticity,]
>>> stiffness_matrix = bilinearform.assemble(v=field, u=field, parallel=False)
>>> system = fem.solve.partition(
>>>     field, stiffness_matrix, dof1=loadcase["dof1"], dof0=loadcase["dof0"]
>>> )
>>> field += fem.solve.solve(*system, ext0=loadcase["ext0"])
>>> field.plot("Principal Values of Logarithmic Strain").show()
class felupe.FormItem(bilinearform, linearform=None, sym=False, args=None, kwargs=None)[source]#

An item to be used in a felupe.Step with bilinear and optional linear form objects based on weak-forms with methods for integration and assembly of vectors / sparse matrices.

Parameters:
  • bilinearform (Form) – A bilinear form object.

  • linearform (Form or None, optional) – A linear form object (default is None). If None, the resulting vector will be filled with zeros.

  • sym (bool, optional) – Flag to active symmetric integration/assembly for bilinear forms (default is False).

  • args (tuple or None, optional) – Tuple with initial optional weakform-arguments (default is None).

  • kwargs (dict or None, optional) – Dictionary with initial optional weakform-keyword-arguments (default is None).

Examples

>>> import felupe as fem
>>> from felupe.math import ddot, sym, trace
>>> mesh = fem.Cube(n=11)
>>> region = fem.RegionHexahedron(mesh)
>>> field = fem.FieldContainer([fem.Field(region, dim=3)])
>>> boundaries, loadcase = fem.dof.uniaxial(field, clamped=True)
>>> @fem.Form(v=field, u=field, grad_v=[True], grad_u=[True])
>>> def bilinearform():
>>>     def a(gradv, gradu, μ=1.0, λ=2.0):
>>>         δε, ε = sym(gradv), sym(gradu)
>>>         return 2 * μ * ddot(δε, ε) + λ * trace(δε) * trace(ε)
>>>     return [a]
>>> item = fem.FormItem(bilinearform, linearform=None, sym=True)
>>> step = fem.Step(items=[item], boundaries=boundaries)
>>> fem.Job(steps=[step]).evaluate()

See also

felupe.Form

A function decorator for a linear- or bilinear-form object.

felupe.Step

A Step with multiple substeps.