Finite Element Library¶
- class sigmaepsilon.solid.fem.cells.elem.FiniteElement(*args, model: Optional[ndarray] = None, activity: Optional[ndarray] = None, density: Optional[ndarray] = None, loads: Optional[ndarray] = None, fields: Optional[dict] = None, strain_loads: Optional[ndarray] = None, t: Optional[ndarray] = None, thickness: Optional[ndarray] = None, fixity: Optional[ndarray] = None, areas: Optional[ndarray] = None, **kwargs)[source]¶
Base mixin class for all element types. Functions of this class can be called on any class of SigmaEpsilon.
- body_load_vector(*, squeeze: bool = True, **kwargs)¶
Builds the equivalent discrete representation of body loads and returns it in either the global frame or cell-local frames.
- Parameters
values (numpy.ndarray, Optional) – Body load values for all cells. Default is None.
constant (bool, Optional) – Set this True if the input represents a constant load. Default is False.
assemble (bool, Optional) – If True, the values are returned with a matching shape to the total system. Default is False.
return_zeroes (bool, Optional) – Controls what happends if there are no strain loads provided. If True, a zero array is retured with correct shape, otherwise None. Default is False.
transform (bool, Optional) – If True, local matrices are transformed to the global frame. Default is True.
See also
- Returns
The nodal load vector for all load cases as a 2d numpy array of shape (nX, nRHS), where nX and nRHS are the total number of unknowns of the structure and the number of load cases.
- Return type
- condensate()[source]¶
Applies static condensation to account for cell fixity.
References
- 1
Duan Jin, Li-Yun-gui “About the Finite Element Analysis for Beam-Hinged Frame,” Advances in Engineering Research, vol. 143, pp. 231-235, 2017.
- consistent_mass_matrix(*, squeeze: bool = True, **kwargs)¶
Returns the stiffness-consistent mass matrix of the cells.
- Parameters
transform (bool, Optional) – If True, local matrices are transformed to the global frame. Default is True.
minval (float, Optional) – A minimal value for the entries in the main diagonal. Set it to a negative value to diable its effect. Default is 1e-12.
sparse (bool, Optional) – If True, the returned object is a sparse COO matrix. Default is False.
- Returns
A sparse SciPy matrix if ‘sparse’ is True, or a 3d numpy array, where the elements run along the first axis.
- Return type
- direction_cosine_matrix(*args, source: Optional[Union[ndarray, str]] = None, target: Optional[Union[ndarray, str]] = None, N: Optional[int] = None, **kwargs) ndarray[source]¶
Returns the DCM matrix from local to global for all elements in the block.
- Parameters
source (str or ReferenceFrame, Optional) – A source frame. Default is None.
target (str or ReferenceFrame, Optional) – A target frame. Default is None.
N (int, Optional) – Number of points. If not specified, the number of nodes is inferred from the class of the instance the function is called upon. Default is None.
- Returns
The dcm matrix for linear transformations from source to target.
- Return type
- dof_solution(*, squeeze: bool = True, **kwargs)¶
Returns nodal displacements for the cells, wrt. their local frames.
- Parameters
points (float or Iterable, Optional) – Points of evaluation. If provided, it is assumed that the given values are wrt. the range [0, 1], unless specified otherwise with the ‘rng’ parameter. If not provided, results are returned for the nodes of the selected elements. Default is None.
rng (Iterable, Optional) – Range where the points of evauation are understood. Default is [0, 1].
cells (int or Iterable[int], Optional) – Indices of cells. If not provided, results are returned for all cells. If cells are provided, the function returns a dictionary, with the cell indices being the keys. Default is None.
target (str or ReferenceFrame, Optional) – Reference frame for the output. A value of None or ‘local’ refers to the local system of the cells. Default is ‘local’.
- Returns
An array of shape (nE, nEVAB, nRHS) if ‘flatten’ is True else (nE, nNE, nDOF, nRHS).
- Return type
- elastic_stiffness_matrix(*, squeeze: bool = True, **kwargs)¶
Returns the elastic stiffness matrix of the cells.
- Parameters
transform (bool, Optional) – If True, local matrices are transformed to the global frame. Default is True.
minval (float, Optional) – A minimal value for the entries in the main diagonal. Set it to a negative value to diable its effect. Default is 1e-12.
sparse (bool, Optional) – If True, the returned object is a sparse COO matrix. Default is False.
- Returns
A sparse SciPy matrix if ‘sparse’ is True, or a 3d numpy array, where the elements run along the first axis.
- Return type
- external_forces(*, squeeze: bool = True, **kwargs)¶
Evaluates \(\mathbf{f}_e = \mathbf{K}_e @ \mathbf{u}_e\) for one or several cells and load cases.
- Parameters
cells (int or Iterable[int], Optional) – Indices of cells. If not provided, results are returned for all cells. Default is None.
target (Union[str, ReferenceFrame], Optional) – The target frame. Default is ‘local’, which means that the returned forces should be understood as coordinates of generalized vectors in the local frames of the cells.
flatten (bool, Optional) – Determines the shape of the resulting array. Default is True.
- Returns
An array of shape (nE, nP * nDOF, nRHS) if ‘flatten’ is True else (nE, nP, nDOF, nRHS).
- Return type
- global_dof_numbering(*args, **kwargs) Union[TopologyArray, ndarray][source]¶
Returns global numbering of the degrees of freedom of the cells.
- integrate_body_loads(values: ndarray) ndarray[source]¶
Returns nodal representation of body loads.
- Parameters
values (numpy.ndarray) – 2d or 3d numpy float array of material densities of shape (nE, nNE * nDOF, nRHS) or (nE, nNE * nDOF), where nE, nNE, nDOF and nRHS stand for the number of elements, nodes per element, number of degrees of freedom and number of load cases respectively.
- Returns
An array of shape (nE, nNE * 6, nRHS).
- Return type
Notes
1) The returned array is always 3 dimensional, even if there is only one load case. 2) Reimplemented for elements with Hermite basis functions.
See also
body_load_vector_bulk()
- internal_forces(*, squeeze: bool = True, **kwargs)¶
Returns internal forces for many cells and evaluation points.
- Parameters
points (float or Iterable[float], Optional) – Points of evaluation. If provided, it is assumed that the given values are wrt. the range [0, 1], unless specified otherwise with the ‘rng’ parameter. If not provided, results are returned for the nodes of the selected elements. Default is None.
rng (Iterable[float], Optional) – Range where the points of evauation are understood. Only for 1d cells. Default is [0, 1].
cells (int or Iterable[int], Optional) – Indices of cells. If not provided, results are returned for all cells. Default is None.
target (Union[str, ReferenceFrame], Optional) – The target frame. Default is ‘local’, which means that the returned forces should be understood as coordinates of generalized vectors in the local frames of the cells.
flatten (bool, Optional) – Determines the shape of the resulting array. Default is True.
- Returns
An array of shape (nE, nP * nDOF, nRHS) if ‘flatten’ is True else (nE, nP, nDOF, nRHS).
- Return type
- kinetic_strains(*, squeeze: bool = True, **kwargs)¶
Returns kinetic strains for one or more cells.
- Parameters
points (float or Iterable, Optional) – Points of evaluation. If provided, it is assumed that the given values are wrt. the range [0, 1], unless specified otherwise with the ‘rng’ parameter. If not provided, results are returned for the nodes of the selected elements. Default is None.
rng (Iterable, Optional) – Range where the points of evauation are understood. Default is [0, 1].
cells (int or Iterable[int], Optional) – Indices of cells. If not provided, results are returned for all cells. Default is None.
- Returns
An array of shape (nE, nSTRE, nRHS), where nE is the number of cells, nSTRE is the number of strain components and nRHS is the number of load cases.
- Return type
- load_vector(*, squeeze: bool = True, **kwargs)¶
Builds the equivalent nodal load vector from all sources and returns it in either the global frame or cell-local frames.
- Parameters
See also
- Returns
The nodal load vector for all load cases as a 2d numpy array of shape (nX, nRHS), where nX and nRHS are the total number of unknowns of the structure and the number of load cases.
- Return type
- strain_load_vector(*, squeeze: bool = True, **kwargs)¶
Generates a load vector from strain loads specified for all cells.
- Parameters
values (numpy.ndarray, Optional) – Strain loads as a 3d numpy array of shape (nE, nS, nRHS). The array must contain values for all cells (nE), strain components (nS) and load cases (nL).
return_zeroes (bool, Optional) – Controls what happends if there are no strain loads provided. If True, a zero array is retured with correct shape, otherwise None. Default is False.
transform (bool, Optional) – If True, local matrices are transformed to the global frame. Default is True.
assemble (bool, Optional) – If True, the values are returned with a matching shape to the total system. Default is False.
See also
- Returns
The equivalent load vector.
- Return type
- strains(*, squeeze: bool = True, **kwargs)¶
Returns strains for one or more cells.
- Parameters
points (float or Iterable[float], Optional) – Points of evaluation. If provided, it is assumed that the given values are wrt. the range [0, 1], unless specified otherwise with the ‘rng’ parameter. If not provided, results are returned for the nodes of the selected elements. Default is None.
rng (Iterable, Optional) – Range where the points of evauation are understood. Only for 1d cells. Default is [0, 1].
cells (int or Iterable[int], Optional) – Indices of cells. If not provided, results are returned for all cells. Default is None.
- Returns
An array of shape (nE, nP, nSTRE, nRHS).
- Return type
Bernoulli Beams¶
- class sigmaepsilon.solid.fem.cells.bernoulli.BernoulliBase(*args, **kwargs)[source]¶
Skeleton class for 1d finite elements, whose bending behaviour is governed by the Bernoulli theory. The implementations covered by this class are independent of the number of nodes of the cells. Specification of an actual finite element consists of specifying class level attributes only.
See also
sigmaepsilon.solid.fem.cells.bernoulli2.Bernoulli2,sigmaepsilon.solid.fem.cells.bernoulli3.Bernoulli3Note
Among normal circumstances, you do not have to interact with this class directly. Only use it if you know what you are doing and understand the meaning of the inputs precisely.
- integrate_body_loads(values: ndarray) ndarray[source]¶
Returns nodal representation of body loads.
- Parameters
values (numpy.ndarray) – 2d or 3d numpy float array of material densities of shape (nE, nNE * nDOF, nRHS) or (nE, nNE * nDOF), where nE, nNE, nDOF and nRHS stand for the number of elements, nodes per element, number of degrees of freedom and number of load cases respectively.
- Returns
An array of shape (nE, nNE * 6, nRHS).
- Return type
Notes
The returned array is always 3 dimensional, even if there is only one load case.
See also
body_load_vector_bulk()
- lumped_mass_matrix(*, squeeze: bool = True, **kwargs)¶
Returns the lumped mass matrix of the block.
- Parameters
alpha (float, Optional) – A nonnegative parameter, typically between 0 and 1/50 (see notes). Default is 1/20.
lumping (str, Optional) – Controls lumping. Currently only direct lumping is available. Default is ‘direct’.
frmt (str, Optional) – Possible values are ‘full’ or ‘diag’. If ‘diag’, only the diagonal entries are returned, if ‘full’ a full matrix is returned. Default is ‘full’.
- shape_function_derivatives(pcoords: Optional[Union[float, Iterable[float]]] = None, *args, rng: Optional[Iterable] = None, jac: Optional[ndarray] = None, dshp: Optional[ndarray] = None, lengths: Optional[Iterable[float]] = None, **kwargs) ndarray[source]¶
Evaluates the shape function derivatives (up to third) at the points specified by ‘pcoords’. The function either returns the derivatives on the master or the actual element, depending on the inputs.
Valid combination of inputs are:
‘pcoords’ and optionally ‘jac’ : this can be used to calculate the derivatives in both the local (‘jac’ is provided) and the parametric frame (‘jac’ is not provided).
‘dshp’ and ‘jac’ : this combination can only be used to return the derivatives wrt. to any frame.
- Parameters
pcoords (float or Iterable, Optional) – Locations of the evaluation points. Default is None.
rng (Iterable, Optional) – The range in which the locations ought to be understood, typically [0, 1] or [-1, 1]. Only if ‘pcoords’ is provided. Default is [0, 1].
lengths (Iterable, Optional) – The lengths of the beams in the block. Default is None.
jac (Iterable, Optional) – The jacobian matrix as a float array of shape (nE, nP, 1, 1), evaluated for each point in each cell. Default is None.
dshp (Iterable, Optional) – The shape function derivatives of the master element as a float array of shape (nP, nNE, nDOF=6, 3), evaluated at a ‘nP’ number of points. Default is None.
- Returns
The returned array has a shape of (nP, nNE=2, nDOF=6, 3), where nP, nNE and nDOF stand for the number of evaluation points, nodes per element and number of degrees of freedom, respectively. Number 3 refers to first, second and third derivatives.
- Return type
Notes
The returned array is always 4 dimensional, even if there is only one evaluation point.
- shape_function_matrix(pcoords: Optional[Union[float, Iterable[float]]] = None, *args, rng: Optional[Iterable] = None, lengths: Optional[Iterable[float]] = None, **kwargs) ndarray[source]¶
Evaluates the shape function matrix at the points specified by ‘pcoords’.
- Parameters
pcoords (float or Iterable) – Locations of the evaluation points.
rng (Iterable, Optional) – The range in which the locations ought to be understood, typically [0, 1] or [-1, 1]. Default is [0, 1].
omitted. (These parameters are for advanced users and can be) –
quantities. (They are here only to avoid repeated evaulation of common) –
lengths (Iterable, Optional) – The lengths of the beams in the block. Default is None.
- Returns
The returned array has a shape of (nE, nP, nDOF, nDOF * nNE), where nE, nP, nNE and nDOF stand for the number of elements, evaluation points, nodes per element and number of degrees of freedom respectively.
- Return type
Notes
The returned array is always 4 dimensional, even if there is only one evaluation point.
- shape_function_values(pcoords: Union[float, Iterable[float]], *args, rng: Optional[Iterable] = None, lengths: Optional[Iterable[float]] = None, **kwargs) ndarray[source]¶
Evaluates the shape functions at the points specified by ‘pcoords’.
- Parameters
- Returns
The returned array has a shape of (nP, nNE=2, nDOF=6), where nP, nNE and nDOF stand for the number of evaluation points, nodes per element and number of degrees of freedom, respectively.
- Return type
- class sigmaepsilon.solid.fem.cells.bernoulli2.Bernoulli2(*args, **kwargs)[source]¶
Finite element class to handle 2-noded Bernoulli beams.
- Model¶
alias of
BernoulliBase
- dshpfnc(L: ndarray)¶
Evaluates the derivatives of the shape functions at several points in the range [-1, 1].
- Parameters
x (numpy.ndarray) – 1d array of floats in the range [-1, -1].
L (numpy.ndarray) – Lengths of the beam elements. Default is None.
- Returns
5d float array of shape (nE, nP, nNE=2, nDOF=6, 3).
- Return type
- shpfnc(L: ndarray)¶
Evaluates the shape functions at several points in the range [-1, 1].
- Parameters
x (numpy.ndarray) – 1d array of floats in the range [-1, -1].
L (numpy.ndarray) – Lengths of the beam elements. Default is None.
- Returns
3d float array of shape (nP, nNE=2, nDOF=6).
- Return type
- class sigmaepsilon.solid.fem.cells.bernoulli3.Bernoulli3(*args, **kwargs)[source]¶
Finite element class to handle 3-noded Bernoulli beams.
- Model¶
alias of
BernoulliBase
- dshpfnc(L: ndarray)¶
Evaluates the derivatives of the shape functions at several points in the range [-1, 1].
- Parameters
x (numpy.ndarray) – 1d array of floats in the range [-1, -1].
L (numpy.ndarray) – Lengths of the beam elements. Default is None.
- Returns
5d float array of shape (nE, nP, nNE=3, nDOF=6, 3).
- Return type
- shpfnc(L: ndarray)¶
Evaluates the shape functions at several points in the range [-1, 1].
- Parameters
x (numpy.ndarray) – 1d array of floats in the range [-1, -1].
L (numpy.ndarray) – Lengths of the beam elements. Default is None.
- Returns
3d float array of shape (nP, nNE=3, nDOF=6).
- Return type