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, **kwargs)[ソース]#
A container for fields which holds a list or a tuple of
Fieldinstances.- パラメータ:
fields (list or tuple of
Field, :class:``~felupe.FieldAxisymmetric`, :class:``~felupe.FieldPlaneStrain` orFieldContainer) -- List with fields. The region is linked to the first field.**kwargs (dict, optional) -- Extra class attributes for the field container.
- evaluate#
Methods to evaluate the deformation gradient and strain measures, see
EvaluateFieldContainerfor details on the provided methods.
サンプル
>>> 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
FieldContaineris also created by one of the logical-and combinations of aField,FieldAxisymmetric,FieldPlaneStrainorFieldContainer.>>> 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
参考
felupe.FieldField on points of a
Regionwith dimensiondimand initial pointvalues.felupe.FieldAxisymmetricAn axisymmetric
Fieldon points of a two dimensionalRegionwith dimensiondim(default is 2) and initial pointvalues(default is 0).felupe.FieldPlaneStrainA plane strain
Fieldon points of a two dimensionalRegionwith dimensiondim(default is 2) and initial pointvalues(default is 0).
- checkpoint()[ソース]#
Return a checkpoint of the field container.
- 戻り値:
A dict with the checkpoint array.
- 戻り値の型:
参考
felupe.FieldContainer.restoreRestore a checkpoint of a field container inplace.
- extract(grad=True, sym=False, add_identity=True, dtype=None, out=None, order='C')[ソース]#
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.
- パラメータ:
grad (bool or list of bool, optional) -- Flag(s) for gradient evaluation(s). A boolean value is applied 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).
orders (str or list of str, optional) -- Controls the memory layout of the outputs. 'C' means it should be C contiguous. 'F' means it should be Fortran contiguous, 'A' means it should be 'F' if the inputs are all 'F', 'C' otherwise. 'K' means it should be as close to the layout as the inputs as is possible, including arbitrarily permuted axes. Default is 'C'.
- 戻り値:
(Symmetric) gradient or interpolated field values evaluated at the integration points of each cell in the region.
- 戻り値の型:
tuple of ndarray
メモ
If the gradient is not requested, the interpolation method returns the field values evaluated at the numeric integration points
qfor each cellcin the region (so-called trailing axes).\[u_{i(qc)} = \hat{u}_{ai}\ h_{a(qc)}\]On the other hand, 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)}\]参考
felupe.Field.interpolateInterpolate field values located at mesh-points to the quadrature points in the region.
felupe.Field.gradGradient as partial derivative of field values w.r.t. undeformed coordinates.
- 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.
- merge(decimals=None)[ソース]#
Merge all fields and return a list of field containers as well as the top-level field container.
- パラメータ:
decimals (int or None, optional) -- Precision decimals for merging duplicated mesh points. Default is None.
- 戻り値:
メモ
注釈
This works only if all regions are template regions, like
RegionQuadorRegionHexahedron, which are supported byFieldDual.サンプル
>>> import felupe as fem >>> >>> mesh1 = fem.Rectangle(n=3) >>> field1 = fem.FieldAxisymmetric(fem.RegionQuad(mesh1), dim=2) >>> >>> mesh2 = fem.Rectangle(a=(1, 0), b=(2, 1), n=3) >>> field2 = fem.FieldAxisymmetric(fem.RegionQuad(mesh2), dim=2) >>> >>> fields, x0 = (field1 & field2).merge() >>> >>> umat = fem.NeoHookeCompressible(mu=1, lmbda=2) >>> solid1 = fem.SolidBody(umat, fields[0]) >>> solid2 = fem.SolidBody(umat, fields[1]) >>> >>> boundaries = fem.dof.uniaxial(x0, clamped=True, return_loadcase=False) >>> >>> step = fem.Step(items=[solid1, solid2], boundaries=boundaries) >>> job = fem.Job(steps=[step]).evaluate(x0=x0)
- plot(*args, project=None, **kwargs)[ソース]#
Plot the first field of the container.
参考
felupe.Scene.plotPlot method of a scene.
felupe.projectProject given values at quadrature-points to mesh-points.
felupe.topointsShift given values at quadrature-points to mesh-points.
- restore(checkpoint)[ソース]#
Restore a checkpoint inplace.
- パラメータ:
checkpoint (dict) -- A dict with checkpoint arrays / objects.
参考
felupe.FieldContainer.checkpointReturn a checkpoint of the field container.
- revolve(n=11, phi=180)[ソース]#
Return a revolved field container.
- パラメータ:
- 戻り値:
The revolved field container.
- 戻り値の型:
サンプル
First, create an axisymmetric field.
>>> import felupe as fem >>> >>> region = fem.RegionQuad(mesh=fem.Rectangle(n=6)) >>> field = fem.FieldContainer([fem.FieldAxisymmetric(region, dim=2)]) >>> field.plot().show()
The first field of the field container is now revolved around the x-axis.
>>> new_field = field.revolve(n=11, phi=180) >>> new_field.plot().show()
参考
SolidBody.revolveReturn a revolved solid body
SolidBodyNearlyIncompressible.revolveReturn a revolved solid body
- screenshot(*args, filename='field.png', transparent_background=None, scale=None, **kwargs)[ソース]#
Take a screenshot of the first field of the container.
参考
pyvista.Plotter.screenshotTake a screenshot of a PyVista plotter.
- 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.
- パラメータ:
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).
- 戻り値:
A object which provides visualization methods for
felupe.FieldContainer.- 戻り値の型:
felupe.ViewField
参考
felupe.ViewFieldVisualization methods for
felupe.FieldContainer.felupe.projectProject given values at quadrature-points to mesh-points.
felupe.topointsShift given values at quadrature-points to mesh-points.
- class felupe.field.EvaluateFieldContainer(field)[ソース]#
Methods to evaluate the deformation gradient and strain measures of a field container.
- パラメータ:
field (FieldContainer) -- A container for fields.
サンプル
>>> import felupe as fem >>> >>> mesh = fem.Rectangle(n=4) >>> region = fem.RegionQuad(mesh) >>> field = fem.FieldContainer([fem.FieldPlaneStrain(region, dim=2)]) >>> >>> evaluate = fem.field.EvaluateFieldContainer(field) >>> F = evaluate.deformation_gradient() >>> >>> F.shape # (3, 3, nquadraturepoints, ncells) (3, 3, 4, 9)
>>> F[..., 0, 0] # deformation gradient of first cell, first quadrature point array([[1., 0., 0.], [0., 1., 0.], [0., 0., 1.]])
参考
felupe.FieldContainerA container which holds one or multiple (mixed) fields.
- deformation_gradient()[ソース]#
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)[ソース]#
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\]- パラメータ:
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).
- 戻り値:
The strain tensor or its principal values.
- 戻り値の型:
ndarray of shape (N, N, ...) tensor, (N * (N + 1) / 2, ...) asvoigt or (N, ...) princ. values
参考
math.strainCompute a Lagrangian strain tensor.
math.strain_stretch_1dCompute the Seth-Hill strains.
- log_strain(tensor=True, asvoigt=False, n=0)[ソース]#
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\]- パラメータ:
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).
- 戻り値:
The strain tensor or its principal values.
- 戻り値の型:
ndarray of shape (N, N, ...) tensor, (N!, ...) asvoigt or (N, ...) princ. values
参考
math.strain_stretch_1dCompute the Seth-Hill strains.
math.strainCompute a Lagrangian strain tensor.
- right_cauchy_green_deformation()[ソース]#
Return the right Cauchy-Green deformation tensor.
\[ \begin{align}\begin{aligned}:label:right-cauchy-green-deformation-tensor\\\boldsymbol{F} &= \frac{\partial \boldsymbol{x}}{\partial \boldsymbol{X}}\\\boldsymbol{C} &= \boldsymbol{F}^T \boldsymbol{F}\\\boldsymbol{C} &= \sum_\alpha \lambda^2_\alpha \ \boldsymbol{N}_\alpha \otimes \boldsymbol{N}_\alpha\end{aligned}\end{align} \]
- strain(fun=<function strain_stretch_1d>, tensor=True, asvoigt=False, n=0, **kwargs)[ソース]#
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\]- パラメータ:
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.
- 戻り値:
The strain tensor or its principal values.
- 戻り値の型:
ndarray of shape (N, N, ...) tensor, (N!, ...) asvoigt or (N, ...) princ. values
参考
math.strainCompute a Lagrangian strain tensor.
math.strain_stretch_1dCompute the Seth-Hill strains.
- class felupe.Field(region, dim=1, values=0.0, dtype=None, **kwargs)[ソース]#
A Field on points of a
Regionwith dimensiondimand initial pointvalues.- パラメータ:
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.
メモ
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
qfor each cellcin 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)}\]サンプル
>>> 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}\]- as_container(**kwargs)[ソース]#
Create a
FieldContainerwith the field.
- extract(grad=True, sym=False, add_identity=True, dtype=None, out=None, order='C')[ソース]#
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.
- パラメータ:
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).
order ({'C', 'F', 'A', 'K'}, optional) -- Controls the memory layout of the output. 'C' means it should be C contiguous. 'F' means it should be Fortran contiguous, 'A' means it should be 'F' if the inputs are all 'F', 'C' otherwise. 'K' means it should be as close to the layout as the inputs as is possible, including arbitrarily permuted axes. Default is 'C'.
- 戻り値:
(Symmetric) gradient or interpolated field values evaluated at the integration points of each cell in the region.
- 戻り値の型:
ndarray
- classmethod from_mesh_container(mesh_container, dim=None, values=0.0)[ソース]#
Create a Field on a vertex mesh from a mesh container.
- grad(sym=False, dtype=None, out=None, order='C')[ソース]#
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)}\]- パラメータ:
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).
order ({'C', 'F', 'A', 'K'}, optional) -- Controls the memory layout of the output. 'C' means it should be C contiguous. 'F' means it should be Fortran contiguous, 'A' means it should be 'F' if the inputs are all 'F', 'C' otherwise. 'K' means it should be as close to the layout as the inputs as is possible, including arbitrarily permuted axes. Default is 'C'.
- 戻り値:
Gradient as partial derivative of field value components
iat points w.r.t. the undeformed coordinatesj, evaluated at the quadrature pointsqof cellscin the region.- 戻り値の型:
ndarray of shape (i, j, q, c)
- hess(dtype=None, out=None, order='C')[ソース]#
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)}\]- パラメータ:
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).
order ({'C', 'F', 'A', 'K'}, optional) -- Controls the memory layout of the output. 'C' means it should be C contiguous. 'F' means it should be Fortran contiguous, 'A' means it should be 'F' if the inputs are all 'F', 'C' otherwise. 'K' means it should be as close to the layout as the inputs as is possible, including arbitrarily permuted axes. Default is 'C'.
- 戻り値:
Hessian as partial derivative of field value components
iat points w.r.t. the undeformed coordinatesjandk, evaluated at the quadrature pointsqof cellscin the region.- 戻り値の型:
ndarray of shape (i, j, k, q, c)
- interpolate(dtype=None, out=None, order='C')[ソース]#
Interpolate field values located at mesh-points to the quadrature points
qof cellscin the region.\[u_{i(qc)} = \hat{u}_{ai}\ h_{a(qc)}\]- パラメータ:
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).
order ({'C', 'F', 'A', 'K'}, optional) -- Controls the memory layout of the output. 'C' means it should be C contiguous. 'F' means it should be Fortran contiguous, 'A' means it should be 'F' if the inputs are all 'F', 'C' otherwise. 'K' means it should be as close to the layout as the inputs as is possible, including arbitrarily permuted axes. Default is 'C'.
- 戻り値:
Interpolated field value components
i, evaluated at the quadrature pointsqof each cellcin the region.- 戻り値の型:
ndarray of shape (i, q, c)
- class felupe.FieldAxisymmetric(region, dim=1, values=0.0, dtype=None)[ソース]#
An axisymmetric
Fieldon points of a two-dimensionalRegionwith dimensiondim(default is 2) and initial pointvalues(default is 0).- パラメータ:
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.
メモ
component 1 = axial component
component 2 = radial component
x_2 (radial direction) ^ | _ | / \ --|-----------------> x_1 (axial rotation axis) \_^
This is a modified
Fieldin which the radial coordinates are evaluated at the numeric integration pointsqfor 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}\]参考
felupe.FieldField on points of a
Regionwith dimensiondimand initial pointvalues.
- as_container(**kwargs)#
Create a
FieldContainerwith the field.
- copy()#
Return a copy of the field.
- extract(grad=True, sym=False, add_identity=True, dtype=None, out=None, order='C')#
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.
- パラメータ:
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).
order ({'C', 'F', 'A', 'K'}, optional) -- Controls the memory layout of the output. 'C' means it should be C contiguous. 'F' means it should be Fortran contiguous, 'A' means it should be 'F' if the inputs are all 'F', 'C' otherwise. 'K' means it should be as close to the layout as the inputs as is possible, including arbitrarily permuted axes. Default is 'C'.
- 戻り値:
(Symmetric) gradient or interpolated field values evaluated at the integration points of each cell in the region.
- 戻り値の型:
ndarray
- fill(a)#
Fill all field values with a scalar value.
- classmethod from_mesh_container(mesh_container, dim=None, values=0.0)#
Create a Field on a vertex mesh from a mesh container.
- grad(sym=False, dtype=None, out=None, order='C')[ソース]#
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}\]- パラメータ:
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).
order ({'C', 'F', 'A', 'K'}, optional) -- Controls the memory layout of the output. 'C' means it should be C contiguous. 'F' means it should be Fortran contiguous, 'A' means it should be 'F' if the inputs are all 'F', 'C' otherwise. 'K' means it should be as close to the layout as the inputs as is possible, including arbitrarily permuted axes. Default is 'C'.
- 戻り値:
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.
- 戻り値の型:
ndarray
- hess(dtype=None, out=None, order='C')#
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)}\]- パラメータ:
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).
order ({'C', 'F', 'A', 'K'}, optional) -- Controls the memory layout of the output. 'C' means it should be C contiguous. 'F' means it should be Fortran contiguous, 'A' means it should be 'F' if the inputs are all 'F', 'C' otherwise. 'K' means it should be as close to the layout as the inputs as is possible, including arbitrarily permuted axes. Default is 'C'.
- 戻り値:
Hessian as partial derivative of field value components
iat points w.r.t. the undeformed coordinatesjandk, evaluated at the quadrature pointsqof cellscin the region.- 戻り値の型:
ndarray of shape (i, j, k, q, c)
- interpolate(dtype=None, out=None, order='C')[ソース]#
Interpolate field values located at mesh-points to the quadrature points
qof cellscin the region.\[u_{i(qc)} = \hat{u}_{ai}\ h_{a(qc)}\]- パラメータ:
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).
order ({'C', 'F', 'A', 'K'}, optional) -- Controls the memory layout of the output. 'C' means it should be C contiguous. 'F' means it should be Fortran contiguous, 'A' means it should be 'F' if the inputs are all 'F', 'C' otherwise. 'K' means it should be as close to the layout as the inputs as is possible, including arbitrarily permuted axes. Default is 'C'.
- 戻り値:
Interpolated field value components
i, evaluated at the quadrature pointsqof each cellcin the region.- 戻り値の型:
ndarray of shape (i, q, c)
- class felupe.FieldPlaneStrain(region, dim=2, values=0.0, dtype=None)[ソース]#
A plane strain
Fieldon points of a two-dimensionalRegionwith dimensiondim(default is 2) and initial pointvalues(default is 0).- パラメータ:
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.
メモ
This is a modified
Fieldfor 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}\]参考
felupe.FieldField on points of a
Regionwith dimensiondimand initial pointvalues.
- as_container(**kwargs)#
Create a
FieldContainerwith the field.
- copy()#
Return a copy of the field.
- extract(grad=True, sym=False, add_identity=True, dtype=None, out=None, order='C')#
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.
- パラメータ:
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).
order ({'C', 'F', 'A', 'K'}, optional) -- Controls the memory layout of the output. 'C' means it should be C contiguous. 'F' means it should be Fortran contiguous, 'A' means it should be 'F' if the inputs are all 'F', 'C' otherwise. 'K' means it should be as close to the layout as the inputs as is possible, including arbitrarily permuted axes. Default is 'C'.
- 戻り値:
(Symmetric) gradient or interpolated field values evaluated at the integration points of each cell in the region.
- 戻り値の型:
ndarray
- fill(a)#
Fill all field values with a scalar value.
- classmethod from_mesh_container(mesh_container, dim=None, values=0.0)#
Create a Field on a vertex mesh from a mesh container.
- grad(sym=False, dtype=None, out=None, order='C')[ソース]#
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 |
- パラメータ:
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).
order ({'C', 'F', 'A', 'K'}, optional) -- Controls the memory layout of the output. 'C' means it should be C contiguous. 'F' means it should be Fortran contiguous, 'A' means it should be 'F' if the inputs are all 'F', 'C' otherwise. 'K' means it should be as close to the layout as the inputs as is possible, including arbitrarily permuted axes. Default is 'C'.
- 戻り値:
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.
- 戻り値の型:
ndarray
- hess(dtype=None, out=None, order='C')[ソース]#
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.
- パラメータ:
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).
order ({'C', 'F', 'A', 'K'}, optional) -- Controls the memory layout of the output. 'C' means it should be C contiguous. 'F' means it should be Fortran contiguous, 'A' means it should be 'F' if the inputs are all 'F', 'C' otherwise. 'K' means it should be as close to the layout as the inputs as is possible, including arbitrarily permuted axes. Default is 'C'.
- 戻り値:
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.
- 戻り値の型:
ndarray
- interpolate(dtype=None, out=None, order='C')[ソース]#
Interpolate field values located at mesh-points to the quadrature points
qof cellscin the region.\[u_{i(qc)} = \hat{u}_{ai}\ h_{a(qc)}\]- パラメータ:
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).
order ({'C', 'F', 'A', 'K'}, optional) -- Controls the memory layout of the output. 'C' means it should be C contiguous. 'F' means it should be Fortran contiguous, 'A' means it should be 'F' if the inputs are all 'F', 'C' otherwise. 'K' means it should be as close to the layout as the inputs as is possible, including arbitrarily permuted axes. Default is 'C'.
- 戻り値:
Interpolated field value components
i, evaluated at the quadrature pointsqof each cellcin the region.- 戻り値の型:
ndarray of shape (i, q, c)
- class felupe.FieldDual(region, dim=1, values=0.0, offset=0, npoints=None, mesh=None, disconnect=None, **kwargs)[ソース]#
A dual field on points of a
Regionwith dimensiondimand initial pointvalues.- パラメータ:
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.
サンプル
>>> 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]) >>> field <felupe FieldContainer object> Number of fields: 2 Dimension of fields: Field: 3 FieldDual: 1
参考
felupe.FieldContainerA container which holds one or multiple (mixed) fields.
felupe.FieldField on points of a
Regionwith dimensiondimand initial pointvalues.felupe.FieldAxisymmetricAxisymmetric field on points of a
Regionwith dimensiondimand initial pointvalues.felupe.FieldPlaneStrainPlane strain field on points of a
Regionwith dimensiondimand initial pointvalues.felupe.mesh.dualCreate a dual
Mesh.
- as_container(**kwargs)#
Create a
FieldContainerwith the field.
- copy()#
Return a copy of the field.
- extract(grad=True, sym=False, add_identity=True, dtype=None, out=None, order='C')#
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.
- パラメータ:
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).
order ({'C', 'F', 'A', 'K'}, optional) -- Controls the memory layout of the output. 'C' means it should be C contiguous. 'F' means it should be Fortran contiguous, 'A' means it should be 'F' if the inputs are all 'F', 'C' otherwise. 'K' means it should be as close to the layout as the inputs as is possible, including arbitrarily permuted axes. Default is 'C'.
- 戻り値:
(Symmetric) gradient or interpolated field values evaluated at the integration points of each cell in the region.
- 戻り値の型:
ndarray
- fill(a)#
Fill all field values with a scalar value.
- classmethod from_mesh_container(mesh_container, dim=None, values=0.0)#
Create a Field on a vertex mesh from a mesh container.
- grad(sym=False, dtype=None, out=None, order='C')#
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)}\]- パラメータ:
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).
order ({'C', 'F', 'A', 'K'}, optional) -- Controls the memory layout of the output. 'C' means it should be C contiguous. 'F' means it should be Fortran contiguous, 'A' means it should be 'F' if the inputs are all 'F', 'C' otherwise. 'K' means it should be as close to the layout as the inputs as is possible, including arbitrarily permuted axes. Default is 'C'.
- 戻り値:
Gradient as partial derivative of field value components
iat points w.r.t. the undeformed coordinatesj, evaluated at the quadrature pointsqof cellscin the region.- 戻り値の型:
ndarray of shape (i, j, q, c)
- hess(dtype=None, out=None, order='C')#
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)}\]- パラメータ:
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).
order ({'C', 'F', 'A', 'K'}, optional) -- Controls the memory layout of the output. 'C' means it should be C contiguous. 'F' means it should be Fortran contiguous, 'A' means it should be 'F' if the inputs are all 'F', 'C' otherwise. 'K' means it should be as close to the layout as the inputs as is possible, including arbitrarily permuted axes. Default is 'C'.
- 戻り値:
Hessian as partial derivative of field value components
iat points w.r.t. the undeformed coordinatesjandk, evaluated at the quadrature pointsqof cellscin the region.- 戻り値の型:
ndarray of shape (i, j, k, q, c)
- interpolate(dtype=None, out=None, order='C')#
Interpolate field values located at mesh-points to the quadrature points
qof cellscin the region.\[u_{i(qc)} = \hat{u}_{ai}\ h_{a(qc)}\]- パラメータ:
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).
order ({'C', 'F', 'A', 'K'}, optional) -- Controls the memory layout of the output. 'C' means it should be C contiguous. 'F' means it should be Fortran contiguous, 'A' means it should be 'F' if the inputs are all 'F', 'C' otherwise. 'K' means it should be as close to the layout as the inputs as is possible, including arbitrarily permuted axes. Default is 'C'.
- 戻り値:
Interpolated field value components
i, evaluated at the quadrature pointsqof each cellcin the region.- 戻り値の型:
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)[ソース]#
A container with multiple (mixed)
Fieldsbased on aRegion.- パラメータ:
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
FieldAxisymmetricas the first field (default is False).planestrain (bool, optional) -- Flag to initiate a
FieldPlaneStrainas 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.
メモ
The dual region is chosen automatically, i.e. for a
RegionHexahedronthe dual region isRegionConstantHexahedron. A total number ofnfields 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.参考
felupe.FieldContainerA container which holds one or multiple (mixed) fields.
felupe.FieldField on points of a
Regionwith dimensiondimand initial pointvalues.felupe.FieldDualA dual field on points of a
Regionwith dimensiondimand initial pointvalues.felupe.FieldAxisymmetricAxisymmetric field on points of a
Regionwith dimensiondimand initial pointvalues.felupe.FieldPlaneStrainPlane strain field on points of a
Regionwith dimensiondimand initial pointvalues.felupe.mesh.dualCreate a dual
Mesh.
- checkpoint()#
Return a checkpoint of the field container.
- 戻り値:
A dict with the checkpoint array.
- 戻り値の型:
参考
felupe.FieldContainer.restoreRestore a checkpoint of a field container inplace.
- copy()#
Return a copy of the field.
- extract(grad=True, sym=False, add_identity=True, dtype=None, out=None, order='C')#
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.
- パラメータ:
grad (bool or list of bool, optional) -- Flag(s) for gradient evaluation(s). A boolean value is applied 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).
orders (str or list of str, optional) -- Controls the memory layout of the outputs. 'C' means it should be C contiguous. 'F' means it should be Fortran contiguous, 'A' means it should be 'F' if the inputs are all 'F', 'C' otherwise. 'K' means it should be as close to the layout as the inputs as is possible, including arbitrarily permuted axes. Default is 'C'.
- 戻り値:
(Symmetric) gradient or interpolated field values evaluated at the integration points of each cell in the region.
- 戻り値の型:
tuple of ndarray
メモ
If the gradient is not requested, the interpolation method returns the field values evaluated at the numeric integration points
qfor each cellcin the region (so-called trailing axes).\[u_{i(qc)} = \hat{u}_{ai}\ h_{a(qc)}\]On the other hand, 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)}\]参考
felupe.Field.interpolateInterpolate field values located at mesh-points to the quadrature points in the region.
felupe.Field.gradGradient as partial derivative of field values w.r.t. undeformed coordinates.
- 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=None)#
Link value array of other field.
- merge(decimals=None)#
Merge all fields and return a list of field containers as well as the top-level field container.
- パラメータ:
decimals (int or None, optional) -- Precision decimals for merging duplicated mesh points. Default is None.
- 戻り値:
メモ
注釈
This works only if all regions are template regions, like
RegionQuadorRegionHexahedron, which are supported byFieldDual.サンプル
>>> import felupe as fem >>> >>> mesh1 = fem.Rectangle(n=3) >>> field1 = fem.FieldAxisymmetric(fem.RegionQuad(mesh1), dim=2) >>> >>> mesh2 = fem.Rectangle(a=(1, 0), b=(2, 1), n=3) >>> field2 = fem.FieldAxisymmetric(fem.RegionQuad(mesh2), dim=2) >>> >>> fields, x0 = (field1 & field2).merge() >>> >>> umat = fem.NeoHookeCompressible(mu=1, lmbda=2) >>> solid1 = fem.SolidBody(umat, fields[0]) >>> solid2 = fem.SolidBody(umat, fields[1]) >>> >>> boundaries = fem.dof.uniaxial(x0, clamped=True, return_loadcase=False) >>> >>> step = fem.Step(items=[solid1, solid2], boundaries=boundaries) >>> job = fem.Job(steps=[step]).evaluate(x0=x0)
- plot(*args, project=None, **kwargs)#
Plot the first field of the container.
参考
felupe.Scene.plotPlot method of a scene.
felupe.projectProject given values at quadrature-points to mesh-points.
felupe.topointsShift given values at quadrature-points to mesh-points.
- restore(checkpoint)#
Restore a checkpoint inplace.
- パラメータ:
checkpoint (dict) -- A dict with checkpoint arrays / objects.
参考
felupe.FieldContainer.checkpointReturn a checkpoint of the field container.
- revolve(n=11, phi=180)#
Return a revolved field container.
- パラメータ:
- 戻り値:
The revolved field container.
- 戻り値の型:
サンプル
First, create an axisymmetric field.
>>> import felupe as fem >>> >>> region = fem.RegionQuad(mesh=fem.Rectangle(n=6)) >>> field = fem.FieldContainer([fem.FieldAxisymmetric(region, dim=2)]) >>> field.plot().show()
The first field of the field container is now revolved around the x-axis.
>>> new_field = field.revolve(n=11, phi=180) >>> new_field.plot().show()
参考
SolidBody.revolveReturn a revolved solid body
SolidBodyNearlyIncompressible.revolveReturn a revolved solid body
- screenshot(*args, filename='field.png', transparent_background=None, scale=None, **kwargs)#
Take a screenshot of the first field of the container.
参考
pyvista.Plotter.screenshotTake 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.
- パラメータ:
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).
- 戻り値:
A object which provides visualization methods for
felupe.FieldContainer.- 戻り値の型:
felupe.ViewField
参考
felupe.ViewFieldVisualization methods for
felupe.FieldContainer.felupe.projectProject given values at quadrature-points to mesh-points.
felupe.topointsShift given values at quadrature-points to mesh-points.