View#

Visualization methods for meshes, mesh and field containers as well as solid bodies by PyVista.

view.Scene()

Base class for plotting a static scene.

ViewMesh(mesh[, point_data, cell_data, ...])

Provide Visualization methods for felupe.Mesh with optional given dicts of point- and cell-data items.

ViewField(field[, point_data, cell_data, ...])

Provide Visualization methods for felupe.FieldContainer.

ViewSolid(field[, solid, stress_type, ...])

Provide Visualization methods for felupe.FieldContainer and felupe.SolidBody.

ViewXdmf(filename[, time])

Provide Visualization methods for a XDMF file generated by Job.evaluate(filename="result.xdmf")().

Detailed API Reference

class felupe.view.Scene[ソース]#

Base class for plotting a static scene.

mesh#

A generalized Dataset with the mesh as well as point- and cell-data. This is not an instance of felupe.Mesh.

Type:

pyvista.UnstructuredGrid

サンプル

>>> import numpy as np
>>> import felupe as fem
>>>
>>> scene = fem.view.Scene()
>>> scene.mesh = fem.Cube(n=3).as_unstructured_grid()
>>> scene.mesh.point_data["Displacement"] = np.arange(81).reshape(27, 3) / 300
>>> scene.mesh.set_active_scalars(None)
>>>
>>> scene.plot("Displacement", component=None).show()
../_images/view-bcd3acedf710259f_00_00.png

参考

felupe.ViewMesh

Provide Visualization methods for a mesh with optional given dicts of point- and cell-data items.

felupe.ViewField

Provide Visualization methods for a field container.

felupe.ViewSolid

Provide Visualization methods for a field container or a solid body.

felupe.ViewMesh(mesh, point_data=None, cell_data=None, cell_type=None)[ソース]#

Provide Visualization methods for felupe.Mesh with optional given dicts of point- and cell-data items.

パラメータ:
  • mesh (felupe.Mesh) -- The mesh object.

  • 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).

felupe.mesh#

A generalized Dataset with the mesh as well as point- and cell-data. This is not an instance of felupe.Mesh.

Type:

pyvista.UnstructuredGrid

サンプル

>>> import numpy as np
>>> import felupe as fem
>>>
>>> mesh = fem.Cube(n=3)
>>> displacement = np.arange(81).reshape(27, 3) / 300
>>> view = fem.ViewMesh(mesh, point_data={"Displacement": displacement})
>>>
>>> view.plot("Displacement", component=None).show()
../_images/view-465046172fdae70e_00_00.png

参考

felupe.view.Scene

Base class for plotting a static scene.

felupe.ViewField

Provide Visualization methods for a field container.

felupe.ViewSolid

Provide Visualization methods for a field container or a solid body.

felupe.ViewField(field, point_data=None, cell_data=None, cell_type=None, project=None)[ソース]#

Provide Visualization methods for felupe.FieldContainer. The warped (deformed) mesh is created from the values of the first field (displacements). By default, the "Deformation Gradient" tensor, the "Logarithmic Strain" tensor and the "Principal Values of Logarithmic Strain" are evaluated as field-related items of the cell-data dict. Optional items of given point- and cell-data overwrite these default field-related cell-data items.

パラメータ:
  • field (felupe.FieldContainer) -- The field-container.

  • 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) -- Callable to project internal cell-data at quadrature-points to mesh-points (default is None). Valid callables are project or extrapolate.

felupe.mesh#

A generalized Dataset with the mesh as well as point- and cell-data. This is not an instance of felupe.Mesh.

Type:

pyvista.UnstructuredGrid

サンプル

参考

felupe.view.Scene

Base class for plotting a static scene.

felupe.ViewMesh

Provide Visualization methods for a mesh with optional given dicts of point- and cell-data items.

felupe.ViewSolid

Provide Visualization methods for a field container or a solid body.

felupe.project

Project given values at quadrature-points to mesh-points.

felupe.ViewSolid(field, solid=None, stress_type='Cauchy', point_data=None, cell_data=None, cell_type=None, project=None, **kwargs)[ソース]#

Provide Visualization methods for felupe.FieldContainer and felupe.SolidBody. The warped (deformed) mesh is created from the values of the first field (displacements). By default, the "Deformation Gradient" tensor, the "Logarithmic Strain" tensor and the "Principal Values of Logarithmic Strain" are evaluated as field-related items of the cell-data dict. Optional items of given point- and cell-data overwrite these default field-related cell-data items.

パラメータ:
  • field (felupe.FieldContainer) -- The field-container.

  • solid (felupe.SolidBody or felupe.SolidBodyIncompressible or None, optional) -- A solid body to evaluate the (Cauchy) stress (default is None).

  • stress_type (str or None, optional) -- The type of stress which is exported, either "Cauchy", "Kirchhoff" or None. If None, the first Piola-Kirchhoff stress (engineering stress in linear elasticity) is used. Default is "Cauchy".

  • 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) -- Callable to project stress at quadrature-points to mesh-points (default is None). Valid callables are project or extrapolate.

felupe.mesh#

A generalized Dataset with the mesh as well as point- and cell-data. This is not an instance of felupe.Mesh.

Type:

pyvista.UnstructuredGrid

サンプル

>>> import numpy as np
>>> import felupe as fem
>>>
>>> mesh = fem.Cube(n=3)
>>> region = fem.RegionHexahedron(mesh)
>>> u = np.sqrt(1 + np.arange(81)).reshape(27, 3) / 100
>>> field = fem.FieldContainer([fem.Field(region, values=u)])
>>> solid = fem.SolidBody(umat=fem.NeoHooke(mu=1, bulk=2), field=field)
>>>
>>> view = fem.ViewSolid(field, solid, project=fem.project)
>>> view.plot("Principal Values of Cauchy Stress").show()
../_images/view-c6744e2b779e845a_00_00.png
../_images/view-c6744e2b779e845a_00_01.png

参考

felupe.view.Scene

Base class for plotting a static scene.

felupe.ViewMesh

Provide Visualization methods for a mesh with optional given dicts of point- and cell-data items.

felupe.ViewField

Provide Visualization methods for a field container.

felupe.project

Project given values at quadrature-points to mesh-points.

felupe.ViewXdmf(filename, time=0)[ソース]#

Provide Visualization methods for a XDMF file generated by Job.evaluate(filename="result.xdmf")(). The warped (deformed) mesh is created from the values of the point-data "Displacement".

パラメータ:
  • filename (str) -- The filename of the XDMF file (including the extension).

  • time (float, optional) -- The time value at which the data is extracted (default is 0).

felupe.mesh#

A generalized Dataset with the mesh as well as point- and cell-data. This is not an instance of felupe.Mesh.

Type:

pyvista.UnstructuredGrid