Field#
A FieldContainer
with pre-defined fields is created with:
|
A field container is created with a list of one or more fields.
|
A container for fields which holds a list or a tuple of |
Available kinds of fields:
|
A Field on points of a |
|
An axisymmetric |
|
A plane strain |
|
A dual field on points of a |
Detailed API Reference
- class felupe.FieldContainer(fields)[source]#
A container for fields which holds a list or a tuple of
Field
instances.- Parameters:
fields (list or tuple of Field, FieldAxisymmetric or FieldPlaneStrain) – List with fields. The region is linked to the first field.
- evaluate#
Methods to evaluate the deformation gradient and strain measures, see
EvaluateFieldContainer
for details on the provided methods.
Examples
>>> import felupe as fem >>> >>> mesh = fem.Cube(n=3) >>> region = fem.RegionHexahedron(mesh) >>> region_dual = fem.RegionConstantHexahedron(mesh.dual(points_per_cell=1)) >>> displacement = fem.Field(region, dim=3) >>> pressure = fem.Field(region_dual) >>> field = fem.FieldContainer([displacement, pressure]) >>> field <felupe FieldContainer object> Number of fields: 2 Dimension of fields: Field: 3 Field: 1
A new
FieldContainer
is also created by one of the logical-and combinations of aField
,FieldAxisymmetric
,FieldPlaneStrain
orFieldContainer
.>>> displacement & pressure <felupe FieldContainer object> Number of fields: 2 Dimension of fields: Field: 3 Field: 1
>>> volume_ratio = fem.Field(region_dual) >>> field & volume_ratio # displacement & pressure & volume_ratio <felupe FieldContainer object> Number of fields: 3 Dimension of fields: Field: 3 Field: 1 Field: 1
See also
felupe.Field
Field on points of a
Region
with dimensiondim
and initial pointvalues
.felupe.FieldAxisymmetric
An axisymmetric
Field
on points of a two dimensionalRegion
with dimensiondim
(default is 2) and initial pointvalues
(default is 0).felupe.FieldPlaneStrain
A plane strain
Field
on points of a two dimensionalRegion
with dimensiondim
(default is 2) and initial pointvalues
(default is 0).
- extract(grad=True, sym=False, add_identity=True, dtype=None, out=None)[source]#
Generalized extraction method which evaluates either the gradient or the field values at the integration points of all cells in the region. Optionally, the symmetric part of the gradient is evaluated and/or the identity matrix is added to the gradient.
- Parameters:
grad (bool or list of bool, optional) – Flag(s) for gradient evaluation(s). A boolean value is appplied on the first field only and all other fields are extracted with
grad=False
. To enable or disable gradients per-field, use a list of boolean values instead (default is True).sym (bool, optional) – Flag for symmetric part if the gradient is evaluated (default is False).
add_identity (bool, optional) – Flag for the addition of the identity matrix if the gradient is evaluated (default is True).
dtype (data-type or None, optional) – If provided, forces the calculation to use the data type specified. Default is None.
out (None or ndarray, optional) – A location into which the result is stored. If provided, it must have a shape that the inputs broadcast to. If not provided or None, a freshly- allocated array is returned (default is None).
- Returns:
(Symmetric) gradient or interpolated field values evaluated at the integration points of each cell in the region.
- Return type:
tuple of ndarray
- imshow(*args, ax=None, dpi=None, **kwargs)[source]#
Take a screenshot of the first field of the container, show the image data in a figure and return the ax.
- plot(*args, project=None, **kwargs)[source]#
Plot the first field of the container.
See also
felupe.Scene.plot
Plot method of a scene.
felupe.project
Project given values at quadrature-points to mesh-points.
felupe.topoints
Shift given values at quadrature-points to mesh-points.
- screenshot(*args, filename='field.png', transparent_background=None, scale=None, **kwargs)[source]#
Take a screenshot of the first field of the container.
See also
pyvista.Plotter.screenshot
Take a screenshot of a PyVista plotter.
- view(point_data=None, cell_data=None, cell_type=None, project=None)[source]#
View the field with optional given dicts of point- and cell-data items.
- Parameters:
point_data (dict or None, optional) – Additional point-data dict (default is None).
cell_data (dict or None, optional) – Additional cell-data dict (default is None).
cell_type (pyvista.CellType or None, optional) – Cell-type of PyVista (default is None).
project (callable or None, optional) – Project internal cell-data at quadrature-points to mesh-points (default is None).
- Returns:
A object which provides visualization methods for
felupe.FieldContainer
.- Return type:
felupe.ViewField
See also
felupe.ViewField
Visualization methods for
felupe.FieldContainer
.felupe.project
Project given values at quadrature-points to mesh-points.
felupe.topoints
Shift given values at quadrature-points to mesh-points.
- class felupe.field.EvaluateFieldContainer(field)[source]#
Methods to evaluate the deformation gradient and strain measures of a field container.
- Parameters:
field (FieldContainer) – A container for fields.
- deformation_gradient()[source]#
Return the Deformation gradient tensor.
(1)#\[ \begin{align}\begin{aligned}\boldsymbol{F} &= \frac{\partial \boldsymbol{x}}{\partial \boldsymbol{X}}\\\boldsymbol{F} &= \sum_\alpha \lambda_\alpha \ \boldsymbol{n}_\alpha \otimes \boldsymbol{N}_\alpha\end{aligned}\end{align} \]
- green_lagrange_strain(tensor=True, asvoigt=False, n=0)[source]#
Return the Green-Lagrange Lagrangian strain tensor or its principal values.
(2)#\[\boldsymbol{E} = \sum_\alpha \frac{1}{2} \left( \lambda_\alpha^2 - 1 \right) \ \boldsymbol{N}_\alpha \otimes \boldsymbol{N}_\alpha\]- Parameters:
tensor (bool, optional) – Assemble and return the strain tensor if True or return its principal values only if False. Default is True.
asvoigt (bool, optional) – Return the symmetric strain tensor in reduced vector storage (default is False).
n (int, optional) – The index of the displacement field (default is 0).
- Returns:
The strain tensor or its principal values.
- Return type:
ndarray of shape (N, N, …) tensor, (N * (N + 1) / 2, …) asvoigt or (N, …) princ. values
See also
math.strain
Compute a Lagrangian strain tensor.
math.strain_stretch_1d
Compute the Seth-Hill strains.
- log_strain(tensor=True, asvoigt=False, n=0)[source]#
Return the logarithmic Lagrangian strain tensor or its principal values.
(3)#\[\boldsymbol{E} = \sum_\alpha \ln(\lambda_\alpha) \ \boldsymbol{N}_\alpha \otimes \boldsymbol{N}_\alpha\]- Parameters:
tensor (bool, optional) – Assemble and return the strain tensor if True or return its principal values only if False. Default is True.
asvoigt (bool, optional) – Return the symmetric strain tensor in reduced vector storage (default is False).
n (int, optional) – The index of the displacement field (default is 0).
- Returns:
The strain tensor or its principal values.
- Return type:
ndarray of shape (N, N, …) tensor, (N!, …) asvoigt or (N, …) princ. values
See also
math.strain_stretch_1d
Compute the Seth-Hill strains.
math.strain
Compute a Lagrangian strain tensor.
- strain(fun=<function strain_stretch_1d>, tensor=True, asvoigt=False, n=0, **kwargs)[source]#
Return the Lagrangian strain tensor or its principal values.
(4)#\[\boldsymbol{E} = \sum_\alpha f_\alpha \left( \lambda_\alpha \right) \ \boldsymbol{N}_\alpha \otimes \boldsymbol{N}_\alpha\]By default, the Seth-Hill strain-stretch relation with a strain exponent of zero is used, see Eq. (5).
(5)#\[\boldsymbol{E} = \sum_\alpha \frac{1}{k} \left( \lambda_\alpha^k - 1 \right) \ \boldsymbol{N}_\alpha \otimes \boldsymbol{N}_\alpha\]- Parameters:
fun (callable, optional) – A callable for the one-dimensional strain-stretch relation. Its Signature must be
lambda stretch, **kwargs: strain
(default is the log. strain,strain_stretch_1d()
withk=0
).tensor (bool, optional) – Assemble and return the strain tensor if True or return its principal values only if False. Default is True.
asvoigt (bool, optional) – Return the symmetric strain tensor in reduced vector storage (default is False).
n (int, optional) – The index of the displacement field (default is 0).
**kwargs (dict, optional) – Optional keyword-arguments are passed to the 1d strain-stretch relation.
- Returns:
The strain tensor or its principal values.
- Return type:
ndarray of shape (N, N, …) tensor, (N!, …) asvoigt or (N, …) princ. values
See also
math.strain
Compute a Lagrangian strain tensor.
math.strain_stretch_1d
Compute the Seth-Hill strains.
- class felupe.Field(region, dim=1, values=0.0, dtype=None, **kwargs)[source]#
A Field on points of a
Region
with dimensiondim
and initial pointvalues
.- Parameters:
region (Region) – The region on which the field will be created.
dim (int, optional) – The dimension of the field (default is 1).
values (float or array) – A single value for all components of the field or an array of shape (region.mesh.npoints, dim). Default is 0.0.
dtype (data-type or None, optional) – Data-type of the array containing the field values.
**kwargs (dict, optional) – Extra class attributes for the field.
Notes
A slice of this field directly accesses the field-values as 1d-array. The interpolation method returns the field values evaluated at the numeric integration points
q
for each cellc
in the region (so-called trailing axes).\[u_{i(qc)} = \hat{u}_{ai}\ h_{a(qc)}\]The gradient method returns the gradient of the field values w.r.t. the undeformed mesh point coordinates, evaluated at the integration points of all cells in the region.
\[\left( \frac{\partial u_i}{\partial X_J} \right)_{(qc)} = \hat{u}_{ai} \left( \frac{\partial h_a}{\partial X_J} \right)_{(qc)}\]Examples
>>> import felupe as fem
>>> mesh = fem.Cube(n=6) >>> region = fem.RegionHexahedron(mesh) >>> displacement = fem.Field(region, dim=3)
>>> u = displacement.interpolate() >>> dudX = displacement.grad()
To obtain deformation-related quantities like the right Cauchy-Green deformation tensor or the principal stretches, use the math-helpers from FElupe. These functions operate on arrays with trailing axes.
\[\boldsymbol{C} = \boldsymbol{F}^T \boldsymbol{F}\]>>> from felupe.math import dot, transpose, eigvalsh, sqrt
>>> F = displacement.extract(grad=True, add_identity=True) >>> C = dot(transpose(F), F) >>> λ = sqrt(eigvalsh(C))
- as_container()[source]#
Create a
FieldContainer
with the field.
- extract(grad=True, sym=False, add_identity=True, dtype=None, out=None)[source]#
Generalized extraction method which evaluates either the gradient or the field values at the integration points of all cells in the region. Optionally, the symmetric part of the gradient is evaluated and/or the identity matrix is added to the gradient.
- Parameters:
grad (bool, optional) – Flag for gradient evaluation (default is True).
sym (bool, optional) – Flag for symmetric part if the gradient is evaluated (default is False).
add_identity (bool, optional) – Flag for the addition of the identity matrix if the gradient is evaluated (default is True).
dtype (data-type or None, optional) – If provided, forces the calculation to use the data type specified. Default is None.
out (None or ndarray, optional) – A location into which the result is stored. If provided, it must have a shape that the inputs broadcast to. If not provided or None, a freshly- allocated array is returned (default is None).
- Returns:
(Symmetric) gradient or interpolated field values evaluated at the integration points of each cell in the region.
- Return type:
ndarray
See also
- grad(sym=False, dtype=None, out=None)[source]#
Gradient as partial derivative of field values w.r.t. undeformed coordinates, evaluated at the integration points of all cells in the region. Optionally, the symmetric part the gradient is evaluated.
\[\left( \frac{\partial u_i}{\partial X_J} \right)_{(qc)} = \hat{u}_{ai} \left( \frac{\partial h_a}{\partial X_J} \right)_{(qc)}\]- Parameters:
sym (bool, optional) – Calculate the symmetric part of the gradient (default is False).
dtype (data-type or None, optional) – If provided, forces the calculation to use the data type specified. Default is None.
out (None or ndarray, optional) – A location into which the result is stored. If provided, it must have a shape that the inputs broadcast to. If not provided or None, a freshly- allocated array is returned (default is None).
- Returns:
Gradient as partial derivative of field value components
i
at points w.r.t. the undeformed coordinatesj
, evaluated at the quadrature pointsq
of cellsc
in the region.- Return type:
ndarray of shape (i, j, q, c)
- hess(dtype=None, out=None)[source]#
Hessian as second partial derivative of field values w.r.t. undeformed coordinates, evaluated at the integration points of all cells in the region.
\[\left( \frac{\partial^2 u_i}{\partial X_J~\partial X_K} \right)_{(qc)} = \hat{u}_{ai} \left( \frac{\partial^2 h_a}{\partial X_J~\partial X_K} \right)_{(qc)}\]- Parameters:
dtype (data-type or None, optional) – If provided, forces the calculation to use the data type specified. Default is None.
out (None or ndarray, optional) – A location into which the result is stored. If provided, it must have a shape that the inputs broadcast to. If not provided or None, a freshly- allocated array is returned (default is None).
- Returns:
Hessian as partial derivative of field value components
i
at points w.r.t. the undeformed coordinatesj
andk
, evaluated at the quadrature pointsq
of cellsc
in the region.- Return type:
ndarray of shape (i, j, k, q, c)
- interpolate(dtype=None, out=None)[source]#
Interpolate field values located at mesh-points to the quadrature points
q
of cellsc
in the region.\[u_{i(qc)} = \hat{u}_{ai}\ h_{a(qc)}\]- Parameters:
dtype (data-type or None, optional) – If provided, forces the calculation to use the data type specified. Default is None.
out (None or ndarray, optional) – A location into which the result is stored. If provided, it must have a shape that the inputs broadcast to. If not provided or None, a freshly- allocated array is returned (default is None).
- Returns:
Interpolated field value components
i
, evaluated at the quadrature pointsq
of each cellc
in the region.- Return type:
ndarray of shape (i, q, c)
- class felupe.FieldAxisymmetric(region, dim=2, values=0.0, dtype=None)[source]#
An axisymmetric
Field
on points of a two-dimensionalRegion
with dimensiondim
(default is 2) and initial pointvalues
(default is 0).- Parameters:
region (Region) – The region on which the field will be created.
dim (int, optional) – The dimension of the field (default is 2).
values (float or array, optional) – A single value for all components of the field or an array of shape (region.mesh.npoints, dim)`. Default is 0.0.
dtype (data-type or None, optional) – Data-type of the array containing the field values.
Notes
component 1 = axial component
component 2 = radial component
x_2 (radial direction) ^ | _ | / \ --|-----------------> x_1 (axial rotation axis) \_^
This is a modified
Field
in which the radial coordinates are evaluated at the numeric integration pointsq
for each cellc
. Thegrad()
-method is modified in such a way that it does not only contain the in-plane 2d-gradient but also the circumferential stretch, see Eq. (6).(6)#\[\begin{split}\frac{\partial \boldsymbol{u}}{\partial \boldsymbol{X}} = \begin{bmatrix} \left( \frac{\partial \boldsymbol{u}}{\partial \boldsymbol{X}} \right)_{2d} & \boldsymbol{0} \\ \boldsymbol{0}^T & \frac{u_r}{R} \end{bmatrix}\end{split}\]See also
felupe.Field
Field on points of a
Region
with dimensiondim
and initial pointvalues
.
- as_container()#
Create a
FieldContainer
with the field.
- copy()#
Return a copy of the field.
- extract(grad=True, sym=False, add_identity=True, dtype=None, out=None)#
Generalized extraction method which evaluates either the gradient or the field values at the integration points of all cells in the region. Optionally, the symmetric part of the gradient is evaluated and/or the identity matrix is added to the gradient.
- Parameters:
grad (bool, optional) – Flag for gradient evaluation (default is True).
sym (bool, optional) – Flag for symmetric part if the gradient is evaluated (default is False).
add_identity (bool, optional) – Flag for the addition of the identity matrix if the gradient is evaluated (default is True).
dtype (data-type or None, optional) – If provided, forces the calculation to use the data type specified. Default is None.
out (None or ndarray, optional) – A location into which the result is stored. If provided, it must have a shape that the inputs broadcast to. If not provided or None, a freshly- allocated array is returned (default is None).
- Returns:
(Symmetric) gradient or interpolated field values evaluated at the integration points of each cell in the region.
- Return type:
ndarray
See also
- fill(a)#
Fill all field values with a scalar value.
- grad(sym=False, dtype=None, out=None)[source]#
3D-gradient as partial derivative of field values at points w.r.t. the undeformed coordinates, evaluated at the integration points of all cells in the region. Optionally, the symmetric part of the gradient is returned.
\[\begin{split}\frac{\partial \boldsymbol{u}}{\partial \boldsymbol{X}} = \begin{bmatrix} \left( \frac{\partial \boldsymbol{u}}{\partial \boldsymbol{X}} \right)_{2d} & \boldsymbol{0} \\ \boldsymbol{0}^T & \frac{u_r}{R} \end{bmatrix}\end{split}\]- Parameters:
sym (bool, optional) – Calculate the symmetric part of the gradient (default is False).
dtype (data-type or None, optional) – If provided, forces the calculation to use the data type specified. Default is None.
out (None or ndarray, optional) – A location into which the result is stored. If provided, it must have a shape that the inputs broadcast to. If not provided or None, a freshly- allocated array is returned (default is None).
- Returns:
Full 3D-gradient as partial derivative of field values at points w.r.t. undeformed coordinates, evaluated at the integration points of all cells in the region.
- Return type:
ndarray
- hess(dtype=None, out=None)#
Hessian as second partial derivative of field values w.r.t. undeformed coordinates, evaluated at the integration points of all cells in the region.
\[\left( \frac{\partial^2 u_i}{\partial X_J~\partial X_K} \right)_{(qc)} = \hat{u}_{ai} \left( \frac{\partial^2 h_a}{\partial X_J~\partial X_K} \right)_{(qc)}\]- Parameters:
dtype (data-type or None, optional) – If provided, forces the calculation to use the data type specified. Default is None.
out (None or ndarray, optional) – A location into which the result is stored. If provided, it must have a shape that the inputs broadcast to. If not provided or None, a freshly- allocated array is returned (default is None).
- Returns:
Hessian as partial derivative of field value components
i
at points w.r.t. the undeformed coordinatesj
andk
, evaluated at the quadrature pointsq
of cellsc
in the region.- Return type:
ndarray of shape (i, j, k, q, c)
- interpolate(dtype=None, out=None)[source]#
Interpolate field values located at mesh-points to the quadrature points
q
of cellsc
in the region.\[u_{i(qc)} = \hat{u}_{ai}\ h_{a(qc)}\]- Parameters:
dtype (data-type or None, optional) – If provided, forces the calculation to use the data type specified. Default is None.
out (None or ndarray, optional) – A location into which the result is stored. If provided, it must have a shape that the inputs broadcast to. If not provided or None, a freshly- allocated array is returned (default is None).
- Returns:
Interpolated field value components
i
, evaluated at the quadrature pointsq
of each cellc
in the region.- Return type:
ndarray of shape (i, q, c)
- class felupe.FieldPlaneStrain(region, dim=2, values=0.0, dtype=None)[source]#
A plane strain
Field
on points of a two-dimensionalRegion
with dimensiondim
(default is 2) and initial pointvalues
(default is 0).- Parameters:
region (Region) – The region on which the field will be created.
dim (int, optional) – The dimension of the field (default is 2).
values (float or array) – A single value for all components of the field or an array of shape (region.mesh.npoints, dim)`. Default is 0.0.
dtype (data-type or None, optional) – Data-type of the array containing the field values.
Notes
This is a modified
Field
for plane strain. Thegrad()
-method is modified in such a way that the in-plane 2d-gradient is embedded in 3d-space, see Eq. (7).(7)#\[\begin{split}\frac{\partial \boldsymbol{u}}{\partial \boldsymbol{X}} = \begin{bmatrix} \left( \frac{\partial \boldsymbol{u}}{\partial \boldsymbol{X}} \right)_{2d} & \boldsymbol{0} \\ \boldsymbol{0}^T & 0 \end{bmatrix}\end{split}\]See also
felupe.Field
Field on points of a
Region
with dimensiondim
and initial pointvalues
.
- as_container()#
Create a
FieldContainer
with the field.
- copy()#
Return a copy of the field.
- extract(grad=True, sym=False, add_identity=True, dtype=None, out=None)#
Generalized extraction method which evaluates either the gradient or the field values at the integration points of all cells in the region. Optionally, the symmetric part of the gradient is evaluated and/or the identity matrix is added to the gradient.
- Parameters:
grad (bool, optional) – Flag for gradient evaluation (default is True).
sym (bool, optional) – Flag for symmetric part if the gradient is evaluated (default is False).
add_identity (bool, optional) – Flag for the addition of the identity matrix if the gradient is evaluated (default is True).
dtype (data-type or None, optional) – If provided, forces the calculation to use the data type specified. Default is None.
out (None or ndarray, optional) – A location into which the result is stored. If provided, it must have a shape that the inputs broadcast to. If not provided or None, a freshly- allocated array is returned (default is None).
- Returns:
(Symmetric) gradient or interpolated field values evaluated at the integration points of each cell in the region.
- Return type:
ndarray
See also
- fill(a)#
Fill all field values with a scalar value.
- grad(sym=False, dtype=None, out=None)[source]#
3D-gradient as partial derivative of field values at points w.r.t. the undeformed coordinates, evaluated at the integration points of all cells in the region. Optionally, the symmetric part of the gradient is returned.
| dudX(2d) : 0 | dudX(planestrain) = | ..................| | 0 : 0 |
- Parameters:
sym (bool, optional) – Calculate the symmetric part of the gradient (default is False).
dtype (data-type or None, optional) – If provided, forces the calculation to use the data type specified. Default is None.
out (None or ndarray, optional) – A location into which the result is stored. If provided, it must have a shape that the inputs broadcast to. If not provided or None, a freshly- allocated array is returned (default is None).
- Returns:
Full 3D-gradient as partial derivative of field values at points w.r.t. undeformed coordinates, evaluated at the integration points of all cells in the region.
- Return type:
ndarray
- hess(dtype=None, out=None)[source]#
3D-Hessian as second partial derivative of field values at points w.r.t. the undeformed coordinates, evaluated at the integration points of all cells in the region. Optionally, the symmetric part of the gradient is returned.
- Parameters:
sym (bool, optional) – Calculate the symmetric part of the gradient (default is False).
dtype (data-type or None, optional) – If provided, forces the calculation to use the data type specified. Default is None.
out (None or ndarray, optional) – A location into which the result is stored. If provided, it must have a shape that the inputs broadcast to. If not provided or None, a freshly- allocated array is returned (default is None).
- Returns:
Full 3D-hessian as second partial derivative of field values at points w.r.t. undeformed coordinates, evaluated at the integration points of all cells in the region.
- Return type:
ndarray
- interpolate(dtype=None, out=None)[source]#
Interpolate field values located at mesh-points to the quadrature points
q
of cellsc
in the region.\[u_{i(qc)} = \hat{u}_{ai}\ h_{a(qc)}\]- Parameters:
dtype (data-type or None, optional) – If provided, forces the calculation to use the data type specified. Default is None.
out (None or ndarray, optional) – A location into which the result is stored. If provided, it must have a shape that the inputs broadcast to. If not provided or None, a freshly- allocated array is returned (default is None).
- Returns:
Interpolated field value components
i
, evaluated at the quadrature pointsq
of each cellc
in the region.- Return type:
ndarray of shape (i, q, c)
- class felupe.FieldDual(region, dim=1, values=0.0, offset=0, npoints=None, mesh=None, disconnect=None, **kwargs)[source]#
A dual field on points of a
Region
with dimensiondim
and initial pointvalues
.- Parameters:
region (Region) – The region on which the field will be created.
dim (int, optional) – The dimension of the field (default is 1).
values (float or array) – A single value for all components of the field or an array of shape (region.mesh.npoints, dim). Default is 0.0.
offset (int, optional) – Offset for cell connectivity of the dual mesh (default is 0).
npoints (int or None, optional) – Specified number of mesh points for the dual mesh (default is None).
mesh (Mesh or None, optional) – A mesh which is used for the dual region (default is None). If None, the mesh is taken from the region.
disconnect (bool or None, optional) – A flag to disconnect the dual mesh (default is None). If None, a disconnected mesh is used except for regions with quadratic-triangle or -tetra or MINI element formulations.
**kwargs (dict, optional) – Optional keyword arguments for the dual region.
Examples
>>> import felupe as fem >>> >>> mesh = fem.Cube(n=6) >>> region = fem.RegionHexahedron(mesh) >>> >>> displacement = fem.Field(region, dim=3) >>> pressure = fem.FieldDual(region) >>> >>> field = fem.FieldContainer([displacement, pressure])
See also
felupe.FieldContainer
A container which holds one or multiple (mixed) fields.
felupe.Field
Field on points of a
Region
with dimensiondim
and initial pointvalues
.felupe.FieldAxisymmetric
Axisymmetric field on points of a
Region
with dimensiondim
and initial pointvalues
.felupe.FieldPlaneStrain
Plane strain field on points of a
Region
with dimensiondim
and initial pointvalues
.felupe.mesh.dual
Create a dual
Mesh
.
- as_container()#
Create a
FieldContainer
with the field.
- copy()#
Return a copy of the field.
- extract(grad=True, sym=False, add_identity=True, dtype=None, out=None)#
Generalized extraction method which evaluates either the gradient or the field values at the integration points of all cells in the region. Optionally, the symmetric part of the gradient is evaluated and/or the identity matrix is added to the gradient.
- Parameters:
grad (bool, optional) – Flag for gradient evaluation (default is True).
sym (bool, optional) – Flag for symmetric part if the gradient is evaluated (default is False).
add_identity (bool, optional) – Flag for the addition of the identity matrix if the gradient is evaluated (default is True).
dtype (data-type or None, optional) – If provided, forces the calculation to use the data type specified. Default is None.
out (None or ndarray, optional) – A location into which the result is stored. If provided, it must have a shape that the inputs broadcast to. If not provided or None, a freshly- allocated array is returned (default is None).
- Returns:
(Symmetric) gradient or interpolated field values evaluated at the integration points of each cell in the region.
- Return type:
ndarray
See also
- fill(a)#
Fill all field values with a scalar value.
- grad(sym=False, dtype=None, out=None)#
Gradient as partial derivative of field values w.r.t. undeformed coordinates, evaluated at the integration points of all cells in the region. Optionally, the symmetric part the gradient is evaluated.
\[\left( \frac{\partial u_i}{\partial X_J} \right)_{(qc)} = \hat{u}_{ai} \left( \frac{\partial h_a}{\partial X_J} \right)_{(qc)}\]- Parameters:
sym (bool, optional) – Calculate the symmetric part of the gradient (default is False).
dtype (data-type or None, optional) – If provided, forces the calculation to use the data type specified. Default is None.
out (None or ndarray, optional) – A location into which the result is stored. If provided, it must have a shape that the inputs broadcast to. If not provided or None, a freshly- allocated array is returned (default is None).
- Returns:
Gradient as partial derivative of field value components
i
at points w.r.t. the undeformed coordinatesj
, evaluated at the quadrature pointsq
of cellsc
in the region.- Return type:
ndarray of shape (i, j, q, c)
- hess(dtype=None, out=None)#
Hessian as second partial derivative of field values w.r.t. undeformed coordinates, evaluated at the integration points of all cells in the region.
\[\left( \frac{\partial^2 u_i}{\partial X_J~\partial X_K} \right)_{(qc)} = \hat{u}_{ai} \left( \frac{\partial^2 h_a}{\partial X_J~\partial X_K} \right)_{(qc)}\]- Parameters:
dtype (data-type or None, optional) – If provided, forces the calculation to use the data type specified. Default is None.
out (None or ndarray, optional) – A location into which the result is stored. If provided, it must have a shape that the inputs broadcast to. If not provided or None, a freshly- allocated array is returned (default is None).
- Returns:
Hessian as partial derivative of field value components
i
at points w.r.t. the undeformed coordinatesj
andk
, evaluated at the quadrature pointsq
of cellsc
in the region.- Return type:
ndarray of shape (i, j, k, q, c)
- interpolate(dtype=None, out=None)#
Interpolate field values located at mesh-points to the quadrature points
q
of cellsc
in the region.\[u_{i(qc)} = \hat{u}_{ai}\ h_{a(qc)}\]- Parameters:
dtype (data-type or None, optional) – If provided, forces the calculation to use the data type specified. Default is None.
out (None or ndarray, optional) – A location into which the result is stored. If provided, it must have a shape that the inputs broadcast to. If not provided or None, a freshly- allocated array is returned (default is None).
- Returns:
Interpolated field value components
i
, evaluated at the quadrature pointsq
of each cellc
in the region.- Return type:
ndarray of shape (i, q, c)
- class felupe.FieldsMixed(region, n=3, values=(0.0, 0.0, 1.0, 0.0), axisymmetric=False, planestrain=False, offset=0, npoints=None, mesh=None, **kwargs)[source]#
A container with multiple (mixed)
Fields
based on aRegion
.- Parameters:
region (RegionHexahedron, RegionQuad, RegionQuadraticQuad, RegionBiQuadraticQuad, RegionQuadraticHexahedron, RegionTriQuadraticHexahedron, RegionQuadraticTetra, RegionQuadraticTriangle, RegionTetraMINI, RegionTriangleMINI or RegionLagrange) – A template region.
n (int, optional) – Number of fields where the first one is a vector field of mesh- dimension and the following fields are scalar-fields (default is 3).
values (tuple of float or tuple of ndarray, optional) – Initial field values (default is (0.0, 0.0, 1.0, 0.0)).
axisymmetric (bool, optional) – Flag to initiate a
FieldAxisymmetric
as the first field (default is False).planestrain (bool, optional) – Flag to initiate a
FieldPlaneStrain
as the first field (default is False).offset (int, optional) – Offset for cell connectivity of the dual mesh (default is 0).
npoints (int or None, optional) – Specified number of mesh points for the dual mesh (default is None).
mesh (Mesh or None, optional) – A mesh which is used for the dual region (default is None). If None, the mesh is taken from the region.
Notes
The dual region is chosen automatically, i.e. for a
RegionHexahedron
the dual region isRegionConstantHexahedron
. A total number ofn
fields are generated inside aFieldContainer
. For compatibility withThreeFieldVariation
, the third field is created with ones, all values of the other fields are initiated with zeros by default.See also
felupe.FieldContainer
A container which holds one or multiple (mixed) fields.
felupe.Field
Field on points of a
Region
with dimensiondim
and initial pointvalues
.felupe.FieldDual
A dual field on points of a
Region
with dimensiondim
and initial pointvalues
.felupe.FieldAxisymmetric
Axisymmetric field on points of a
Region
with dimensiondim
and initial pointvalues
.felupe.FieldPlaneStrain
Plane strain field on points of a
Region
with dimensiondim
and initial pointvalues
.felupe.mesh.dual
Create a dual
Mesh
.
- copy()#
Return a copy of the field.
- extract(grad=True, sym=False, add_identity=True, dtype=None, out=None)#
Generalized extraction method which evaluates either the gradient or the field values at the integration points of all cells in the region. Optionally, the symmetric part of the gradient is evaluated and/or the identity matrix is added to the gradient.
- Parameters:
grad (bool or list of bool, optional) – Flag(s) for gradient evaluation(s). A boolean value is appplied on the first field only and all other fields are extracted with
grad=False
. To enable or disable gradients per-field, use a list of boolean values instead (default is True).sym (bool, optional) – Flag for symmetric part if the gradient is evaluated (default is False).
add_identity (bool, optional) – Flag for the addition of the identity matrix if the gradient is evaluated (default is True).
dtype (data-type or None, optional) – If provided, forces the calculation to use the data type specified. Default is None.
out (None or ndarray, optional) – A location into which the result is stored. If provided, it must have a shape that the inputs broadcast to. If not provided or None, a freshly- allocated array is returned (default is None).
- Returns:
(Symmetric) gradient or interpolated field values evaluated at the integration points of each cell in the region.
- Return type:
tuple of ndarray
- imshow(*args, ax=None, dpi=None, **kwargs)#
Take a screenshot of the first field of the container, show the image data in a figure and return the ax.
- link(other_field)#
Link value array of other field.
- plot(*args, project=None, **kwargs)#
Plot the first field of the container.
See also
felupe.Scene.plot
Plot method of a scene.
felupe.project
Project given values at quadrature-points to mesh-points.
felupe.topoints
Shift given values at quadrature-points to mesh-points.
- screenshot(*args, filename='field.png', transparent_background=None, scale=None, **kwargs)#
Take a screenshot of the first field of the container.
See also
pyvista.Plotter.screenshot
Take a screenshot of a PyVista plotter.
- values()#
Return the field values.
- view(point_data=None, cell_data=None, cell_type=None, project=None)#
View the field with optional given dicts of point- and cell-data items.
- Parameters:
point_data (dict or None, optional) – Additional point-data dict (default is None).
cell_data (dict or None, optional) – Additional cell-data dict (default is None).
cell_type (pyvista.CellType or None, optional) – Cell-type of PyVista (default is None).
project (callable or None, optional) – Project internal cell-data at quadrature-points to mesh-points (default is None).
- Returns:
A object which provides visualization methods for
felupe.FieldContainer
.- Return type:
felupe.ViewField
See also
felupe.ViewField
Visualization methods for
felupe.FieldContainer
.felupe.project
Project given values at quadrature-points to mesh-points.
felupe.topoints
Shift given values at quadrature-points to mesh-points.