Adding Basic Components

Creating a new Feb object

The basics

We can create a new FEB file with basically no information. To do so, we can simply call Feb and provide it a version. Right now, we accept versions 2.5 and 3.0.

[1]:
from febio_python import Feb

feb = Feb(version=2.5)
feb
[1]:
Feb25(2.5):
-> Module: Unknown
-> Control: 0
-> Material: 0
-> Globals: 0
-> Geometry: 0
-> Boundary: 0
-> Loads: 0
-> Discrete: 0
-> LoadData: 0
-> Output: 0
-> MeshData: 0

You will notices that all contents of the feb file are zero. This means that there is nothing stored in the feb right now. Let’s start modifying it. We will begin with simple structures, like module, globals and controls:

[2]:
feb.setup_module(module_type="solid") # default values
feb.setup_globals(T=0, R=0, Fc=0) # default values
feb.setup_controls(analysis_type="static") # here, you can change basic settings. See docs for more info.
feb.setup_output(variables=["displacement", "Lagrange strain", "stress"]) # default values

We can check the FEB object. It now has some data:

[3]:
feb
[3]:
Feb25(2.5):
-> Module: solid
-> Control: 15
-> Material: 0
-> Globals: 1
-> Geometry: 0
-> Boundary: 0
-> Loads: 0
-> Discrete: 0
-> LoadData: 0
-> Output: 1
-> MeshData: 0

Adding Mesh

We will now add more interestind data to our feb. Let’s start by creating a simple mesh with pyvista:

[4]:
import pyvista as pv
pv.set_jupyter_backend('static')

# Create a simple mesh (plane mesh)
mesh = pv.Plane(direction=(0,0,1), i_size=2, j_size=1, i_resolution=6, j_resolution=3)
mesh = mesh.cast_to_unstructured_grid() # we will be using unstructured grid for FEBio
mesh.plot(show_edges=True, window_size=(800, 400), cpos='xy')
../../../../_images/modules_feb_spec_2_5_adding_basic_components_adding_basic_components_11_0.png

In Pyvista, an unstructured grid is defined as “points” (nodes) and “cells” (elements). We can access them like this:

[5]:
mesh.points
[5]:
pyvista_ndarray([[-1.        , -0.5       ,  0.        ],
                 [-0.6666667 , -0.5       ,  0.        ],
                 [-0.33333334, -0.5       ,  0.        ],
                 [ 0.        , -0.5       ,  0.        ],
                 [ 0.33333334, -0.5       ,  0.        ],
                 [ 0.6666667 , -0.5       ,  0.        ],
                 [ 1.        , -0.5       ,  0.        ],
                 [-1.        , -0.16666667,  0.        ],
                 [-0.6666667 , -0.16666667,  0.        ],
                 [-0.33333334, -0.16666667,  0.        ],
                 [ 0.        , -0.16666667,  0.        ],
                 [ 0.33333334, -0.16666667,  0.        ],
                 [ 0.6666667 , -0.16666667,  0.        ],
                 [ 1.        , -0.16666667,  0.        ],
                 [-1.        ,  0.16666667,  0.        ],
                 [-0.6666667 ,  0.16666667,  0.        ],
                 [-0.33333334,  0.16666667,  0.        ],
                 [ 0.        ,  0.16666667,  0.        ],
                 [ 0.33333334,  0.16666667,  0.        ],
                 [ 0.6666667 ,  0.16666667,  0.        ],
                 [ 1.        ,  0.16666667,  0.        ],
                 [-1.        ,  0.5       ,  0.        ],
                 [-0.6666667 ,  0.5       ,  0.        ],
                 [-0.33333334,  0.5       ,  0.        ],
                 [ 0.        ,  0.5       ,  0.        ],
                 [ 0.33333334,  0.5       ,  0.        ],
                 [ 0.6666667 ,  0.5       ,  0.        ],
                 [ 1.        ,  0.5       ,  0.        ]], dtype=float32)
[6]:
mesh.cells_dict
[6]:
{9: array([[ 0,  1,  8,  7],
        [ 1,  2,  9,  8],
        [ 2,  3, 10,  9],
        [ 3,  4, 11, 10],
        [ 4,  5, 12, 11],
        [ 5,  6, 13, 12],
        [ 7,  8, 15, 14],
        [ 8,  9, 16, 15],
        [ 9, 10, 17, 16],
        [10, 11, 18, 17],
        [11, 12, 19, 18],
        [12, 13, 20, 19],
        [14, 15, 22, 21],
        [15, 16, 23, 22],
        [16, 17, 24, 23],
        [17, 18, 25, 24],
        [18, 19, 26, 25],
        [19, 20, 27, 26]], dtype=int64)}
[7]:
print(f"Number of nodes: {mesh.points.shape[0]}")
print(f"Number of elements: {mesh.cells_dict[pv.CellType.QUAD].shape[0]}")
Number of nodes: 28
Number of elements: 18
[8]:
from febio_python.core import Nodes, Elements
# Create nodes
nodes = Nodes(name="plane", coordinates=mesh.points)
# Create elements
elements = Elements(name="plane_elements",
                    mat="1",
                    part=None,
                    type="QUAD",
                    connectivity=mesh.cells_dict[pv.CellType.QUAD.value],
                    )
# Add nodes and elements to the feb object (need to be lists)
feb.add_nodes([nodes])
feb.add_elements([elements])

We can now see that we have two ‘Geometry’ items:

[9]:
feb
[9]:
Feb25(2.5):
-> Module: solid
-> Control: 15
-> Material: 0
-> Globals: 1
-> Geometry: 2
-> Boundary: 0
-> Loads: 0
-> Discrete: 0
-> LoadData: 0
-> Output: 1
-> MeshData: 0

We can further inspect them:

[10]:
feb.inspect_nodes_and_elements()
--> Nodes 'plane': 28
--> Elements 'plane_elements': 18

Modfying nodes/elements

Oh no! We have defined a too coarse mesh! Let’s create a fine mesh and update our feb.

[11]:
# Create a simple mesh (plane mesh)
mesh = pv.Plane(direction=(0,0,1), i_size=2, j_size=1, i_resolution=10, j_resolution=5)
mesh = mesh.cast_to_unstructured_grid() # we will be using unstructured grid for FEBio
mesh.plot(show_edges=True, window_size=(800, 400), cpos='xy')
../../../../_images/modules_feb_spec_2_5_adding_basic_components_adding_basic_components_23_0.png
[12]:
print(f"Number of nodes: {mesh.points.shape[0]}")
print(f"Number of elements: {mesh.cells_dict[pv.CellType.QUAD].shape[0]}")
Number of nodes: 66
Number of elements: 50

Note that if we use the function “ADD” again, it will add the new nodes to the existing nodes in the same mesh. This is used when we are creating mesh in an iterative process. Let’s try:

[13]:
# Create nodes
nodes = Nodes(name="plane", coordinates=mesh.points)
# Create elements
elements = Elements(name="plane_elements",
                    mat="1",
                    part=None,
                    type="QUAD",
                    connectivity=mesh.cells_dict[pv.CellType.QUAD.value],
                    )
# Add nodes and elements to the feb object (need to be lists)
feb.add_nodes([nodes])
feb.add_elements([elements])
# Inspect nodes and elements
feb.inspect_nodes_and_elements()
--> Nodes 'plane': 94
--> Elements 'plane_elements': 68

You can see that it had inscrease the number of nodes/elements. But this is not what we are looking for. Let’s clean the feb data and add it again:

[14]:
feb.clear_nodes()
feb.clear_elements()
feb
[14]:
Feb25(2.5):
-> Module: solid
-> Control: 15
-> Material: 0
-> Globals: 1
-> Geometry: 0
-> Boundary: 0
-> Loads: 0
-> Discrete: 0
-> LoadData: 0
-> Output: 1
-> MeshData: 0
[15]:
feb.add_nodes([nodes])
feb.add_elements([elements])
feb.inspect_nodes_and_elements()
--> Nodes 'plane': 66
--> Elements 'plane_elements': 50

Is there another method? We need triangle elements!! Also, What if we have multiple nodes or elements and we do not want to delete all of them and just update the existing ones? Let’s create another mesh and update the ‘plane’ mesh in the feb file.

[16]:
# Create a simple mesh (plane mesh)
mesh = pv.Plane(direction=(0,0,1), i_size=2, j_size=1, i_resolution=20, j_resolution=10)
mesh = mesh.triangulate() # we will be using unstructured grid for FEBio
mesh = mesh.cast_to_unstructured_grid() # we will be using unstructured grid for FEBio
mesh.plot(show_edges=True, window_size=(800, 400), cpos='xy')
../../../../_images/modules_feb_spec_2_5_adding_basic_components_adding_basic_components_31_0.png
[17]:
print(f"Number of nodes: {mesh.points.shape[0]}")
print(f"Number of elements: {mesh.cells_dict[pv.CellType.TRIANGLE].shape[0]}")
Number of nodes: 231
Number of elements: 400
[18]:
# Create nodes
nodes = Nodes(name="plane", coordinates=mesh.points)
# Create elements
elements = Elements(name="plane_elements",
                    mat="1",
                    part=None,
                    type="TRIANGLE",
                    connectivity=mesh.cells_dict[pv.CellType.TRIANGLE.value],
                    )
# Add nodes and elements to the feb object (need to be lists)
feb.update_nodes([nodes])
feb.update_elements([elements])
# Inspect nodes and elements
feb.inspect_nodes_and_elements()
--> Nodes 'plane': 231
--> Elements 'plane_elements': 400

Adding Load

Let’s try adding a shear load to the mesh:

[19]:
import numpy as np
from febio_python.core import NodalLoad, LoadCurve, NodeSet

# First, let's select the nodeset to apply the load.
# we will add load to the nodes at the right edge of the mesh
x_values = mesh.points[:, 0]
selected_nodes = np.where(x_values == x_values.max())[0]

# Create a nodeset
nodeset = NodeSet(name="right_edge", ids=selected_nodes)
# Create a nodal load
shear_load = NodalLoad(node_set="right_edge",
                       dof="y",
                       scale=-25.0, # negative value means force is pointing in the negative direction
                       load_curve=1,
                       )
# Create a load curve
load_curve = LoadCurve(id=1, interpolate_type="smooth", data=[(0, 0), (1, 1)])

# Add nodeset, nodal load and load curve to the feb object
feb.add_node_sets([nodeset])
feb.add_nodal_loads([shear_load])
feb.add_load_curves([load_curve])
feb
[19]:
Feb25(2.5):
-> Module: solid
-> Control: 15
-> Material: 0
-> Globals: 1
-> Geometry: 3
-> Boundary: 0
-> Loads: 1
-> Discrete: 0
-> LoadData: 1
-> Output: 1
-> MeshData: 0

Let’s try to add a non-linear load now. Pointing in the x-direction

[20]:
# first, get the y position of the nodes
y_values = mesh.points[selected_nodes, 1]
# then, normalize the y values, so that we can use them as 'maps'
y_values = (y_values - y_values.min()) / (y_values.max() - y_values.min())
# now, create a normal distribution curve using the y values, centered at
# middle of the y values (0.5)
tensile_load_map = np.exp(-((y_values - 0.5) ** 2) / 0.1)
# plot the load curve
import matplotlib.pyplot as plt
plt.plot(y_values, tensile_load_map)
plt.xlabel("Normalized y position")
plt.ylabel("Load curve")
plt.title("Load curve for the shear load")
plt.show()

# Create a nodal load
tensile_load = NodalLoad(node_set="right_edge",
                       dof="x",
                       scale=100.0*tensile_load_map, #
                       load_curve=1,
                       )
feb.add_nodal_loads([tensile_load])
feb
../../../../_images/modules_feb_spec_2_5_adding_basic_components_adding_basic_components_38_0.png
[20]:
Feb25(2.5):
-> Module: solid
-> Control: 15
-> Material: 0
-> Globals: 1
-> Geometry: 3
-> Boundary: 0
-> Loads: 2
-> Discrete: 0
-> LoadData: 1
-> Output: 1
-> MeshData: 1

Add boundary condition

Now, let’s fix the left boundary of the mesh. We will be applying fixed condition to restrain the mesh in all coordinates. In addition, since this is a simple plane mesh, we will apply constraint in z (for all nodes) and rotation (shell constrain) for all the left nodes.

[21]:
from febio_python.core import FixCondition

# Fix the left edge of the mesh
x_values = mesh.points[:, 0]
selected_nodes = np.where(x_values == x_values.min())[0]
all_nodes = np.arange(mesh.points.shape[0]) # used to fix all nodes in the z direction
# Create nodesets
left_nodeset = NodeSet(name="left_edge", ids=selected_nodes)
all_nodeset = NodeSet(name="all_nodes", ids=all_nodes)
# Create a fix condition
left_fix_condition = FixCondition(dof="x,y,z,sx,sy", node_set="left_edge")
# we will
all_fix_condition = FixCondition(dof="z", node_set="all_nodes")

# Add nodeset and fix condition to the feb object
feb.add_node_sets([left_nodeset, all_nodeset])
feb.add_boundary_conditions([left_fix_condition, all_fix_condition])
[22]:
feb
[22]:
Feb25(2.5):
-> Module: solid
-> Control: 15
-> Material: 0
-> Globals: 1
-> Geometry: 5
-> Boundary: 2
-> Loads: 2
-> Discrete: 0
-> LoadData: 1
-> Output: 1
-> MeshData: 1

Adding Material

Now, let’s add a material. We will use simple Isotropic-Elastic.

[23]:
from febio_python.core import Material

mat = Material(
    id=1,
    type="isotropic elastic",
    name="plane material",
    parameters=dict(
        E=1e6,
        v=0.3,
        density=1,
    )
)

feb.add_materials([mat])

Add Element data

[24]:
from febio_python.core import ElementData

shell_thickness = ElementData(
    name="Element thickness",
    var="shell thickness",
    elem_set="plane_elements",
    data=np.full((mesh.n_cells, 3), 0.01),
    ids=np.arange(0, mesh.n_cells + 1),
)
feb.add_element_data([shell_thickness])

Running FEB

Uncomment and run the code below:

[27]:
# from febio_python.feb import run
# # Save the FEB file
# feb.write("plane_mesh_v25.feb")
# run("plane_mesh_v25.feb")
Starting process for: plane_mesh_v25.feb
Waiting for process to complete...
Completed process for: plane_mesh_v25.feb

Reading XPLT

[29]:
from febio_python import Xplt

xplt = Xplt("plane_mesh_v25.xplt")
xplt
[29]:
Xplt object [1691985673520]:
=== xplt_mesh: ===
-Nodes:
--->None: (231, 3)
-Elements:
--->None: (400, 3)
=== States: ===
-Nodes:
--->displacement: (11, 231, 3)
-Elements:
--->Lagrange strain: (11, 400, 6)
--->stress: (11, 400, 6)
-Surfaces:

Example properties

[30]:
xplt.nodes # list of Nodes objects
[30]:
[Nodes(name=None, coordinates=array([[-1.        , -0.5       ,  0.        ],
        [-0.89999998, -0.5       ,  0.        ],
        [-0.80000001, -0.5       ,  0.        ],
        [-0.69999999, -0.5       ,  0.        ],
        [-0.60000002, -0.5       ,  0.        ],
        [-0.5       , -0.5       ,  0.        ],
        [-0.40000001, -0.5       ,  0.        ],
        [-0.30000001, -0.5       ,  0.        ],
        [-0.2       , -0.5       ,  0.        ],
        [-0.1       , -0.5       ,  0.        ],
        [ 0.        , -0.5       ,  0.        ],
        [ 0.1       , -0.5       ,  0.        ],
        [ 0.2       , -0.5       ,  0.        ],
        [ 0.30000001, -0.5       ,  0.        ],
        [ 0.40000001, -0.5       ,  0.        ],
        [ 0.5       , -0.5       ,  0.        ],
        [ 0.60000002, -0.5       ,  0.        ],
        [ 0.69999999, -0.5       ,  0.        ],
        [ 0.80000001, -0.5       ,  0.        ],
        [ 0.89999998, -0.5       ,  0.        ],
        [ 1.        , -0.5       ,  0.        ],
        [-1.        , -0.40000001,  0.        ],
        [-0.89999998, -0.40000001,  0.        ],
        [-0.80000001, -0.40000001,  0.        ],
        [-0.69999999, -0.40000001,  0.        ],
        [-0.60000002, -0.40000001,  0.        ],
        [-0.5       , -0.40000001,  0.        ],
        [-0.40000001, -0.40000001,  0.        ],
        [-0.30000001, -0.40000001,  0.        ],
        [-0.2       , -0.40000001,  0.        ],
        [-0.1       , -0.40000001,  0.        ],
        [ 0.        , -0.40000001,  0.        ],
        [ 0.1       , -0.40000001,  0.        ],
        [ 0.2       , -0.40000001,  0.        ],
        [ 0.30000001, -0.40000001,  0.        ],
        [ 0.40000001, -0.40000001,  0.        ],
        [ 0.5       , -0.40000001,  0.        ],
        [ 0.60000002, -0.40000001,  0.        ],
        [ 0.69999999, -0.40000001,  0.        ],
        [ 0.80000001, -0.40000001,  0.        ],
        [ 0.89999998, -0.40000001,  0.        ],
        [ 1.        , -0.40000001,  0.        ],
        [-1.        , -0.30000001,  0.        ],
        [-0.89999998, -0.30000001,  0.        ],
        [-0.80000001, -0.30000001,  0.        ],
        [-0.69999999, -0.30000001,  0.        ],
        [-0.60000002, -0.30000001,  0.        ],
        [-0.5       , -0.30000001,  0.        ],
        [-0.40000001, -0.30000001,  0.        ],
        [-0.30000001, -0.30000001,  0.        ],
        [-0.2       , -0.30000001,  0.        ],
        [-0.1       , -0.30000001,  0.        ],
        [ 0.        , -0.30000001,  0.        ],
        [ 0.1       , -0.30000001,  0.        ],
        [ 0.2       , -0.30000001,  0.        ],
        [ 0.30000001, -0.30000001,  0.        ],
        [ 0.40000001, -0.30000001,  0.        ],
        [ 0.5       , -0.30000001,  0.        ],
        [ 0.60000002, -0.30000001,  0.        ],
        [ 0.69999999, -0.30000001,  0.        ],
        [ 0.80000001, -0.30000001,  0.        ],
        [ 0.89999998, -0.30000001,  0.        ],
        [ 1.        , -0.30000001,  0.        ],
        [-1.        , -0.2       ,  0.        ],
        [-0.89999998, -0.2       ,  0.        ],
        [-0.80000001, -0.2       ,  0.        ],
        [-0.69999999, -0.2       ,  0.        ],
        [-0.60000002, -0.2       ,  0.        ],
        [-0.5       , -0.2       ,  0.        ],
        [-0.40000001, -0.2       ,  0.        ],
        [-0.30000001, -0.2       ,  0.        ],
        [-0.2       , -0.2       ,  0.        ],
        [-0.1       , -0.2       ,  0.        ],
        [ 0.        , -0.2       ,  0.        ],
        [ 0.1       , -0.2       ,  0.        ],
        [ 0.2       , -0.2       ,  0.        ],
        [ 0.30000001, -0.2       ,  0.        ],
        [ 0.40000001, -0.2       ,  0.        ],
        [ 0.5       , -0.2       ,  0.        ],
        [ 0.60000002, -0.2       ,  0.        ],
        [ 0.69999999, -0.2       ,  0.        ],
        [ 0.80000001, -0.2       ,  0.        ],
        [ 0.89999998, -0.2       ,  0.        ],
        [ 1.        , -0.2       ,  0.        ],
        [-1.        , -0.1       ,  0.        ],
        [-0.89999998, -0.1       ,  0.        ],
        [-0.80000001, -0.1       ,  0.        ],
        [-0.69999999, -0.1       ,  0.        ],
        [-0.60000002, -0.1       ,  0.        ],
        [-0.5       , -0.1       ,  0.        ],
        [-0.40000001, -0.1       ,  0.        ],
        [-0.30000001, -0.1       ,  0.        ],
        [-0.2       , -0.1       ,  0.        ],
        [-0.1       , -0.1       ,  0.        ],
        [ 0.        , -0.1       ,  0.        ],
        [ 0.1       , -0.1       ,  0.        ],
        [ 0.2       , -0.1       ,  0.        ],
        [ 0.30000001, -0.1       ,  0.        ],
        [ 0.40000001, -0.1       ,  0.        ],
        [ 0.5       , -0.1       ,  0.        ],
        [ 0.60000002, -0.1       ,  0.        ],
        [ 0.69999999, -0.1       ,  0.        ],
        [ 0.80000001, -0.1       ,  0.        ],
        [ 0.89999998, -0.1       ,  0.        ],
        [ 1.        , -0.1       ,  0.        ],
        [-1.        ,  0.        ,  0.        ],
        [-0.89999998,  0.        ,  0.        ],
        [-0.80000001,  0.        ,  0.        ],
        [-0.69999999,  0.        ,  0.        ],
        [-0.60000002,  0.        ,  0.        ],
        [-0.5       ,  0.        ,  0.        ],
        [-0.40000001,  0.        ,  0.        ],
        [-0.30000001,  0.        ,  0.        ],
        [-0.2       ,  0.        ,  0.        ],
        [-0.1       ,  0.        ,  0.        ],
        [ 0.        ,  0.        ,  0.        ],
        [ 0.1       ,  0.        ,  0.        ],
        [ 0.2       ,  0.        ,  0.        ],
        [ 0.30000001,  0.        ,  0.        ],
        [ 0.40000001,  0.        ,  0.        ],
        [ 0.5       ,  0.        ,  0.        ],
        [ 0.60000002,  0.        ,  0.        ],
        [ 0.69999999,  0.        ,  0.        ],
        [ 0.80000001,  0.        ,  0.        ],
        [ 0.89999998,  0.        ,  0.        ],
        [ 1.        ,  0.        ,  0.        ],
        [-1.        ,  0.1       ,  0.        ],
        [-0.89999998,  0.1       ,  0.        ],
        [-0.80000001,  0.1       ,  0.        ],
        [-0.69999999,  0.1       ,  0.        ],
        [-0.60000002,  0.1       ,  0.        ],
        [-0.5       ,  0.1       ,  0.        ],
        [-0.40000001,  0.1       ,  0.        ],
        [-0.30000001,  0.1       ,  0.        ],
        [-0.2       ,  0.1       ,  0.        ],
        [-0.1       ,  0.1       ,  0.        ],
        [ 0.        ,  0.1       ,  0.        ],
        [ 0.1       ,  0.1       ,  0.        ],
        [ 0.2       ,  0.1       ,  0.        ],
        [ 0.30000001,  0.1       ,  0.        ],
        [ 0.40000001,  0.1       ,  0.        ],
        [ 0.5       ,  0.1       ,  0.        ],
        [ 0.60000002,  0.1       ,  0.        ],
        [ 0.69999999,  0.1       ,  0.        ],
        [ 0.80000001,  0.1       ,  0.        ],
        [ 0.89999998,  0.1       ,  0.        ],
        [ 1.        ,  0.1       ,  0.        ],
        [-1.        ,  0.2       ,  0.        ],
        [-0.89999998,  0.2       ,  0.        ],
        [-0.80000001,  0.2       ,  0.        ],
        [-0.69999999,  0.2       ,  0.        ],
        [-0.60000002,  0.2       ,  0.        ],
        [-0.5       ,  0.2       ,  0.        ],
        [-0.40000001,  0.2       ,  0.        ],
        [-0.30000001,  0.2       ,  0.        ],
        [-0.2       ,  0.2       ,  0.        ],
        [-0.1       ,  0.2       ,  0.        ],
        [ 0.        ,  0.2       ,  0.        ],
        [ 0.1       ,  0.2       ,  0.        ],
        [ 0.2       ,  0.2       ,  0.        ],
        [ 0.30000001,  0.2       ,  0.        ],
        [ 0.40000001,  0.2       ,  0.        ],
        [ 0.5       ,  0.2       ,  0.        ],
        [ 0.60000002,  0.2       ,  0.        ],
        [ 0.69999999,  0.2       ,  0.        ],
        [ 0.80000001,  0.2       ,  0.        ],
        [ 0.89999998,  0.2       ,  0.        ],
        [ 1.        ,  0.2       ,  0.        ],
        [-1.        ,  0.30000001,  0.        ],
        [-0.89999998,  0.30000001,  0.        ],
        [-0.80000001,  0.30000001,  0.        ],
        [-0.69999999,  0.30000001,  0.        ],
        [-0.60000002,  0.30000001,  0.        ],
        [-0.5       ,  0.30000001,  0.        ],
        [-0.40000001,  0.30000001,  0.        ],
        [-0.30000001,  0.30000001,  0.        ],
        [-0.2       ,  0.30000001,  0.        ],
        [-0.1       ,  0.30000001,  0.        ],
        [ 0.        ,  0.30000001,  0.        ],
        [ 0.1       ,  0.30000001,  0.        ],
        [ 0.2       ,  0.30000001,  0.        ],
        [ 0.30000001,  0.30000001,  0.        ],
        [ 0.40000001,  0.30000001,  0.        ],
        [ 0.5       ,  0.30000001,  0.        ],
        [ 0.60000002,  0.30000001,  0.        ],
        [ 0.69999999,  0.30000001,  0.        ],
        [ 0.80000001,  0.30000001,  0.        ],
        [ 0.89999998,  0.30000001,  0.        ],
        [ 1.        ,  0.30000001,  0.        ],
        [-1.        ,  0.40000001,  0.        ],
        [-0.89999998,  0.40000001,  0.        ],
        [-0.80000001,  0.40000001,  0.        ],
        [-0.69999999,  0.40000001,  0.        ],
        [-0.60000002,  0.40000001,  0.        ],
        [-0.5       ,  0.40000001,  0.        ],
        [-0.40000001,  0.40000001,  0.        ],
        [-0.30000001,  0.40000001,  0.        ],
        [-0.2       ,  0.40000001,  0.        ],
        [-0.1       ,  0.40000001,  0.        ],
        [ 0.        ,  0.40000001,  0.        ],
        [ 0.1       ,  0.40000001,  0.        ],
        [ 0.2       ,  0.40000001,  0.        ],
        [ 0.30000001,  0.40000001,  0.        ],
        [ 0.40000001,  0.40000001,  0.        ],
        [ 0.5       ,  0.40000001,  0.        ],
        [ 0.60000002,  0.40000001,  0.        ],
        [ 0.69999999,  0.40000001,  0.        ],
        [ 0.80000001,  0.40000001,  0.        ],
        [ 0.89999998,  0.40000001,  0.        ],
        [ 1.        ,  0.40000001,  0.        ],
        [-1.        ,  0.5       ,  0.        ],
        [-0.89999998,  0.5       ,  0.        ],
        [-0.80000001,  0.5       ,  0.        ],
        [-0.69999999,  0.5       ,  0.        ],
        [-0.60000002,  0.5       ,  0.        ],
        [-0.5       ,  0.5       ,  0.        ],
        [-0.40000001,  0.5       ,  0.        ],
        [-0.30000001,  0.5       ,  0.        ],
        [-0.2       ,  0.5       ,  0.        ],
        [-0.1       ,  0.5       ,  0.        ],
        [ 0.        ,  0.5       ,  0.        ],
        [ 0.1       ,  0.5       ,  0.        ],
        [ 0.2       ,  0.5       ,  0.        ],
        [ 0.30000001,  0.5       ,  0.        ],
        [ 0.40000001,  0.5       ,  0.        ],
        [ 0.5       ,  0.5       ,  0.        ],
        [ 0.60000002,  0.5       ,  0.        ],
        [ 0.69999999,  0.5       ,  0.        ],
        [ 0.80000001,  0.5       ,  0.        ],
        [ 0.89999998,  0.5       ,  0.        ],
        [ 1.        ,  0.5       ,  0.        ]]), ids=array([  0,   1,   2,   3,   4,   5,   6,   7,   8,   9,  10,  11,  12,
         13,  14,  15,  16,  17,  18,  19,  20,  21,  22,  23,  24,  25,
         26,  27,  28,  29,  30,  31,  32,  33,  34,  35,  36,  37,  38,
         39,  40,  41,  42,  43,  44,  45,  46,  47,  48,  49,  50,  51,
         52,  53,  54,  55,  56,  57,  58,  59,  60,  61,  62,  63,  64,
         65,  66,  67,  68,  69,  70,  71,  72,  73,  74,  75,  76,  77,
         78,  79,  80,  81,  82,  83,  84,  85,  86,  87,  88,  89,  90,
         91,  92,  93,  94,  95,  96,  97,  98,  99, 100, 101, 102, 103,
        104, 105, 106, 107, 108, 109, 110, 111, 112, 113, 114, 115, 116,
        117, 118, 119, 120, 121, 122, 123, 124, 125, 126, 127, 128, 129,
        130, 131, 132, 133, 134, 135, 136, 137, 138, 139, 140, 141, 142,
        143, 144, 145, 146, 147, 148, 149, 150, 151, 152, 153, 154, 155,
        156, 157, 158, 159, 160, 161, 162, 163, 164, 165, 166, 167, 168,
        169, 170, 171, 172, 173, 174, 175, 176, 177, 178, 179, 180, 181,
        182, 183, 184, 185, 186, 187, 188, 189, 190, 191, 192, 193, 194,
        195, 196, 197, 198, 199, 200, 201, 202, 203, 204, 205, 206, 207,
        208, 209, 210, 211, 212, 213, 214, 215, 216, 217, 218, 219, 220,
        221, 222, 223, 224, 225, 226, 227, 228, 229, 230]))]
[31]:
xplt.elements # list of Elements objects
[31]:
[Elements(name=None, type='TRIANGLE', connectivity=array([[  0,   1,  21],
        [ 22,  21,   1],
        [  1,   2,  22],
        ...,
        [229, 228, 208],
        [208, 209, 229],
        [230, 229, 209]]), ids=array([  1,   2,   3,   4,   5,   6,   7,   8,   9,  10,  11,  12,  13,
         14,  15,  16,  17,  18,  19,  20,  21,  22,  23,  24,  25,  26,
         27,  28,  29,  30,  31,  32,  33,  34,  35,  36,  37,  38,  39,
         40,  41,  42,  43,  44,  45,  46,  47,  48,  49,  50,  51,  52,
         53,  54,  55,  56,  57,  58,  59,  60,  61,  62,  63,  64,  65,
         66,  67,  68,  69,  70,  71,  72,  73,  74,  75,  76,  77,  78,
         79,  80,  81,  82,  83,  84,  85,  86,  87,  88,  89,  90,  91,
         92,  93,  94,  95,  96,  97,  98,  99, 100, 101, 102, 103, 104,
        105, 106, 107, 108, 109, 110, 111, 112, 113, 114, 115, 116, 117,
        118, 119, 120, 121, 122, 123, 124, 125, 126, 127, 128, 129, 130,
        131, 132, 133, 134, 135, 136, 137, 138, 139, 140, 141, 142, 143,
        144, 145, 146, 147, 148, 149, 150, 151, 152, 153, 154, 155, 156,
        157, 158, 159, 160, 161, 162, 163, 164, 165, 166, 167, 168, 169,
        170, 171, 172, 173, 174, 175, 176, 177, 178, 179, 180, 181, 182,
        183, 184, 185, 186, 187, 188, 189, 190, 191, 192, 193, 194, 195,
        196, 197, 198, 199, 200, 201, 202, 203, 204, 205, 206, 207, 208,
        209, 210, 211, 212, 213, 214, 215, 216, 217, 218, 219, 220, 221,
        222, 223, 224, 225, 226, 227, 228, 229, 230, 231, 232, 233, 234,
        235, 236, 237, 238, 239, 240, 241, 242, 243, 244, 245, 246, 247,
        248, 249, 250, 251, 252, 253, 254, 255, 256, 257, 258, 259, 260,
        261, 262, 263, 264, 265, 266, 267, 268, 269, 270, 271, 272, 273,
        274, 275, 276, 277, 278, 279, 280, 281, 282, 283, 284, 285, 286,
        287, 288, 289, 290, 291, 292, 293, 294, 295, 296, 297, 298, 299,
        300, 301, 302, 303, 304, 305, 306, 307, 308, 309, 310, 311, 312,
        313, 314, 315, 316, 317, 318, 319, 320, 321, 322, 323, 324, 325,
        326, 327, 328, 329, 330, 331, 332, 333, 334, 335, 336, 337, 338,
        339, 340, 341, 342, 343, 344, 345, 346, 347, 348, 349, 350, 351,
        352, 353, 354, 355, 356, 357, 358, 359, 360, 361, 362, 363, 364,
        365, 366, 367, 368, 369, 370, 371, 372, 373, 374, 375, 376, 377,
        378, 379, 380, 381, 382, 383, 384, 385, 386, 387, 388, 389, 390,
        391, 392, 393, 394, 395, 396, 397, 398, 399, 400]), mat=None, part=1)]
[32]:
xplt.states # list of States objects
[32]:
States(nodes=[StateData(name='displacement', dom=-1, data=array([[[ 0.        ,  0.        ,  0.        ],
        [ 0.        ,  0.        ,  0.        ],
        [ 0.        ,  0.        ,  0.        ],
        ...,
        [ 0.        ,  0.        ,  0.        ],
        [ 0.        ,  0.        ,  0.        ],
        [ 0.        ,  0.        ,  0.        ]],

       [[ 0.        ,  0.        ,  0.        ],
        [-0.00241185, -0.00141284,  0.        ],
        [-0.00460829, -0.00296836,  0.        ],
        ...,
        [ 0.03514386, -0.07879543,  0.        ],
        [ 0.03548158, -0.08488715,  0.        ],
        [ 0.03574249, -0.09127856,  0.        ]],

       [[ 0.        ,  0.        ,  0.        ],
        [-0.00438763, -0.00255033,  0.        ],
        [-0.00835898, -0.00537066,  0.        ],
        ...,
        [ 0.06173341, -0.1440275 ,  0.        ],
        [ 0.06217256, -0.154956  ,  0.        ],
        [ 0.06246274, -0.16649833,  0.        ]],

       ...,

       [[ 0.        ,  0.        ,  0.        ],
        [-0.01119141, -0.00654127,  0.        ],
        [-0.02091032, -0.0140982 ,  0.        ],
        ...,
        [ 0.14965743, -0.3825299 ,  0.        ],
        [ 0.14999932, -0.40942118,  0.        ],
        [ 0.14971283, -0.4390082 ,  0.        ]],

       [[ 0.        ,  0.        ,  0.        ],
        [-0.01183677, -0.00695707,  0.        ],
        [-0.02204275, -0.01506192,  0.        ],
        ...,
        [ 0.15935734, -0.40824923,  0.        ],
        [ 0.15973683, -0.43669814,  0.        ],
        [ 0.15939233, -0.46821073,  0.        ]],

       [[ 0.        ,  0.        ,  0.        ],
        [-0.01239656, -0.00733244,  0.        ],
        [-0.02300878, -0.01594689,  0.        ],
        ...,
        [ 0.16844292, -0.43168595,  0.        ],
        [ 0.16888118, -0.4615291 ,  0.        ],
        [ 0.1684947 , -0.49480793,  0.        ]]], dtype=float32))], elements=[StateData(name='Lagrange strain', dom=0, data=array([[[ 0.00000000e+00,  0.00000000e+00,  0.00000000e+00,
          0.00000000e+00,  0.00000000e+00,  0.00000000e+00],
        [ 0.00000000e+00,  0.00000000e+00,  0.00000000e+00,
          0.00000000e+00,  0.00000000e+00,  0.00000000e+00],
        [ 0.00000000e+00,  0.00000000e+00,  0.00000000e+00,
          0.00000000e+00,  0.00000000e+00,  0.00000000e+00],
        ...,
        [ 0.00000000e+00,  0.00000000e+00,  0.00000000e+00,
          0.00000000e+00,  0.00000000e+00,  0.00000000e+00],
        [ 0.00000000e+00,  0.00000000e+00,  0.00000000e+00,
          0.00000000e+00,  0.00000000e+00,  0.00000000e+00],
        [ 0.00000000e+00,  0.00000000e+00,  0.00000000e+00,
          0.00000000e+00,  0.00000000e+00,  0.00000000e+00]],

       [[-2.37153228e-02,  4.96507617e-08,  9.25132260e-03,
         -7.05668749e-03,  8.77791172e-05,  8.75434343e-05],
        [-1.55513035e-02,  4.19573579e-03,  6.64580334e-03,
         -8.76835082e-04, -4.10028315e-06,  3.52151233e-06],
        [-2.16151774e-02,  4.19573579e-03,  6.46888278e-03,
         -3.71682271e-03, -1.87471978e-05, -3.00115735e-05],
        ...,
        [ 5.26419561e-03, -1.79060502e-03, -1.32319576e-03,
         -1.09630043e-03,  1.35697308e-04, -9.22947293e-05],
        [ 3.76348314e-03, -1.79060502e-03, -6.20446925e-04,
         -1.20302581e-03,  3.88543667e-05,  1.51752083e-05],
        [ 4.56833979e-03, -4.31854278e-03, -2.29850310e-04,
         -2.06456636e-03, -6.20934181e-04,  2.85940798e-04]],

       [[-4.25659008e-02,  1.59901958e-07,  1.65927690e-02,
         -1.27380788e-02,  1.59318617e-04,  1.58670024e-04],
        [-2.77379584e-02,  7.46000186e-03,  1.18833957e-02,
         -1.59418711e-03, -7.71454506e-06,  5.81513541e-06],
        [-3.85503881e-02,  7.46000186e-03,  1.15433391e-02,
         -6.74896734e-03, -3.39396829e-05, -5.31097103e-05],
        ...,
        [ 1.04284454e-02, -3.59612191e-03, -2.62235897e-03,
         -1.94375543e-03,  2.62485846e-04, -1.94640495e-04],
        [ 7.60030141e-03, -3.59612191e-03, -1.25896640e-03,
         -2.19966378e-03,  7.99982663e-05,  2.65306644e-05],
        [ 9.37342457e-03, -8.74353014e-03, -5.09420410e-04,
         -3.88680655e-03, -1.21673197e-03,  6.25563203e-04]],

       ...,

       [[-1.03454977e-01,  1.02291756e-06,  4.00860123e-02,
         -3.26714478e-02,  4.17361414e-04,  4.09782020e-04],
        [-6.36420026e-02,  1.68143846e-02,  2.79301554e-02,
         -4.35827812e-03, -2.22580402e-05,  9.12787800e-06],
        [-8.96661356e-02,  1.68143846e-02,  2.68445406e-02,
         -1.82178933e-02, -8.80838707e-05, -1.21186145e-04],
        ...,
        [ 3.98490652e-02, -1.43709024e-02, -9.90649406e-03,
         -4.42765653e-03,  9.21171857e-04, -8.83480243e-04],
        [ 3.06131430e-02, -1.43709024e-02, -5.02734119e-03,
         -6.15815260e-03,  3.49141104e-04,  5.58151733e-05],
        [ 3.98919024e-02, -3.57674919e-02, -2.53269938e-03,
         -1.22170681e-02, -4.47947998e-03,  3.11689964e-03]],

       [[-1.08881332e-01,  1.15577188e-06,  4.21402082e-02,
         -3.47481929e-02,  4.44912352e-04,  4.35584836e-04],
        [-6.61997572e-02,  1.74349323e-02,  2.92007625e-02,
         -4.69052047e-03, -2.38736247e-05,  8.95932772e-06],
        [-9.36254486e-02,  1.74349323e-02,  2.80210096e-02,
         -1.95591226e-02, -9.38163721e-05, -1.26437721e-04],
        ...,
        [ 4.45734523e-02, -1.61647107e-02, -1.10511193e-02,
         -4.62550670e-03,  1.02171360e-03, -1.00175163e-03],
        [ 3.43595631e-02, -1.61647107e-02, -5.62138064e-03,
         -6.66673388e-03,  3.96208838e-04,  5.73999459e-05],
        [ 4.50464040e-02, -4.02821973e-02, -2.88027129e-03,
         -1.33820577e-02, -4.99419915e-03,  3.56116169e-03]],

       [[-1.13529623e-01,  1.28303191e-06,  4.38876860e-02,
         -3.66229638e-02,  4.69880557e-04,  4.58679424e-04],
        [-6.81892633e-02,  1.79007482e-02,  3.02395485e-02,
         -5.00305137e-03, -2.53279559e-05,  8.70800432e-06],
        [-9.68414694e-02,  1.79007482e-02,  2.89726872e-02,
         -2.08097380e-02, -9.90094777e-05, -1.30736560e-04],
        ...,
        [ 4.92642596e-02, -1.79609023e-02, -1.21797724e-02,
         -4.79171658e-03,  1.12077757e-03, -1.11997582e-03],
        [ 3.80752347e-02, -1.79609023e-02, -6.20363094e-03,
         -7.15704774e-03,  4.43766272e-04,  5.84696500e-05],
        [ 5.01992963e-02, -4.47962992e-02, -3.22508742e-03,
         -1.45127149e-02, -5.50433807e-03,  4.00813343e-03]]],
      dtype=float32)), StateData(name='stress', dom=0, data=array([[[ 0.00000000e+00,  0.00000000e+00,  0.00000000e+00,
          0.00000000e+00,  0.00000000e+00,  0.00000000e+00],
        [ 0.00000000e+00,  0.00000000e+00,  0.00000000e+00,
          0.00000000e+00,  0.00000000e+00,  0.00000000e+00],
        [ 0.00000000e+00,  0.00000000e+00,  0.00000000e+00,
          0.00000000e+00,  0.00000000e+00,  0.00000000e+00],
        ...,
        [ 0.00000000e+00,  0.00000000e+00,  0.00000000e+00,
          0.00000000e+00,  0.00000000e+00,  0.00000000e+00],
        [ 0.00000000e+00,  0.00000000e+00,  0.00000000e+00,
          0.00000000e+00,  0.00000000e+00,  0.00000000e+00],
        [ 0.00000000e+00,  0.00000000e+00,  0.00000000e+00,
          0.00000000e+00,  0.00000000e+00,  0.00000000e+00]],

       [[-2.57104258e+04, -8.32338770e+03, -1.26783472e+03,
         -5.00688428e+03,  6.47943573e+01,  5.92901344e+01],
        [-1.43045918e+04,  5.29054504e+02,  2.43950635e+03,
         -5.20961426e+02, -3.57376099e+00,  6.30108118e-01],
        [-2.22493047e+04, -3.06733716e+03, -1.37425366e+03,
         -2.51297241e+03, -1.46889915e+01, -2.43459587e+01],
        ...,
        [ 5.21573779e+03, -1.35725298e+01,  2.21767639e+02,
         -1.16937219e+03,  1.08122086e+02, -6.51837006e+01],
        [ 3.57415332e+03, -4.65898315e+02,  3.02325745e+02,
         -1.17880908e+03,  2.89701767e+01,  1.33194933e+01],
        [ 3.34443286e+03, -3.05312988e+03, -1.65135666e+02,
         -1.99652502e+03, -4.88636292e+02,  1.92099396e+02]],

       [[-4.48986719e+04, -1.49390605e+04, -2.35327466e+03,
         -8.44365625e+03,  1.14041161e+02,  9.63064270e+01],
        [-2.49900234e+04,  9.53104980e+02,  4.44329688e+03,
         -7.30678345e+02, -7.26328564e+00, -2.20842171e+00],
        [-3.87194727e+04, -5.47384131e+03, -2.50470239e+03,
         -4.08328662e+03, -2.69665909e+01, -4.50076408e+01],
        ...,
        [ 1.01858477e+04,  1.17725616e+02,  4.08006592e+02,
         -2.64370044e+03,  2.15578400e+02, -1.29467834e+02],
        [ 7.06949512e+03, -7.02262695e+02,  6.12641724e+02,
         -2.59109229e+03,  5.83572578e+01,  2.63640137e+01],
        [ 6.61690186e+03, -5.68505859e+03, -3.21698456e+02,
         -4.45147803e+03, -9.77332092e+02,  3.83120911e+02]],

       ...,

       [[-9.92679531e+04, -3.66051602e+04, -6.64452246e+03,
         -1.68790059e+04,  2.73736389e+02,  1.60902344e+02],
        [-5.36645586e+04,  2.33872754e+03,  1.14601699e+04,
         -3.11623627e+02, -2.48822613e+01, -3.28758545e+01],
        [-8.36302812e+04, -1.31645332e+04, -6.56354932e+03,
         -7.30365869e+03, -7.28974838e+01, -1.22733002e+02],
        ...,
        [ 3.75855078e+04,  2.70891431e+03,  1.32031433e+03,
         -1.40452988e+04,  8.46112244e+02, -5.05644989e+02],
        [ 2.67828867e+04,  4.67094666e+02,  2.55825366e+03,
         -1.31989033e+04,  2.39480713e+02,  1.01752251e+02],
        [ 2.51069902e+04, -1.50549033e+04, -1.02165302e+03,
         -2.35506895e+04, -3.90820166e+03,  1.51560632e+03]],

       [[-1.03544789e+05, -3.85584258e+04, -7.13044482e+03,
         -1.75093594e+04,  2.89711700e+02,  1.62880554e+02],
        [-5.55197031e+04,  2.46515479e+03,  1.21460254e+04,
         -2.14675766e+02, -2.70178871e+01, -3.74516640e+01],
        [-8.68775391e+04, -1.38355576e+04, -6.98250635e+03,
         -7.54512988e+03, -7.78906708e+01, -1.31013916e+02],
        ...,
        [ 4.19955117e+04,  3.29983154e+03,  1.46036755e+03,
         -1.61555713e+04,  9.49114014e+02, -5.67525269e+02],
        [ 2.99045957e+04,  8.80143616e+02,  2.87391748e+03,
         -1.51469600e+04,  2.70181335e+02,  1.13927521e+02],
        [ 2.80805938e+04, -1.58543672e+04, -1.11878113e+03,
         -2.71071465e+04, -4.39746387e+03,  1.70316260e+03]],

       [[-1.07134500e+05, -4.02281328e+04, -7.57302295e+03,
         -1.80556445e+04,  3.04073822e+02,  1.64106705e+02],
        [-5.69336562e+04,  2.57574097e+03,  1.27474004e+04,
         -1.31921478e+02, -2.89687786e+01, -4.16831512e+01],
        [-8.95182266e+04, -1.44029912e+04, -7.35732373e+03,
         -7.77286230e+03, -8.24211807e+01, -1.38445099e+02],
        ...,
        [ 4.63999648e+04,  3.91850513e+03,  1.60029834e+03,
         -1.83055371e+04,  1.05160632e+03, -6.29289734e+02],
        [ 3.29913398e+04,  1.32787097e+03,  3.18541235e+03,
         -1.71267695e+04,  3.01004028e+02,  1.25993431e+02],
        [ 3.10444766e+04, -1.64842617e+04, -1.21277368e+03,
         -3.07300488e+04, -4.88736279e+03,  1.89061377e+03]]],
      dtype=float32))], surfaces=[], timesteps=array([0. , 0.1, 0.2, 0.3, 0.4, 0.5, 0.6, 0.7, 0.8, 0.9, 1. ],
      dtype=float32))

States data:

[33]:
xplt.states.nodes[0]
[33]:
StateData(name='displacement', dom=-1, data=array([[[ 0.        ,  0.        ,  0.        ],
        [ 0.        ,  0.        ,  0.        ],
        [ 0.        ,  0.        ,  0.        ],
        ...,
        [ 0.        ,  0.        ,  0.        ],
        [ 0.        ,  0.        ,  0.        ],
        [ 0.        ,  0.        ,  0.        ]],

       [[ 0.        ,  0.        ,  0.        ],
        [-0.00241185, -0.00141284,  0.        ],
        [-0.00460829, -0.00296836,  0.        ],
        ...,
        [ 0.03514386, -0.07879543,  0.        ],
        [ 0.03548158, -0.08488715,  0.        ],
        [ 0.03574249, -0.09127856,  0.        ]],

       [[ 0.        ,  0.        ,  0.        ],
        [-0.00438763, -0.00255033,  0.        ],
        [-0.00835898, -0.00537066,  0.        ],
        ...,
        [ 0.06173341, -0.1440275 ,  0.        ],
        [ 0.06217256, -0.154956  ,  0.        ],
        [ 0.06246274, -0.16649833,  0.        ]],

       ...,

       [[ 0.        ,  0.        ,  0.        ],
        [-0.01119141, -0.00654127,  0.        ],
        [-0.02091032, -0.0140982 ,  0.        ],
        ...,
        [ 0.14965743, -0.3825299 ,  0.        ],
        [ 0.14999932, -0.40942118,  0.        ],
        [ 0.14971283, -0.4390082 ,  0.        ]],

       [[ 0.        ,  0.        ,  0.        ],
        [-0.01183677, -0.00695707,  0.        ],
        [-0.02204275, -0.01506192,  0.        ],
        ...,
        [ 0.15935734, -0.40824923,  0.        ],
        [ 0.15973683, -0.43669814,  0.        ],
        [ 0.15939233, -0.46821073,  0.        ]],

       [[ 0.        ,  0.        ,  0.        ],
        [-0.01239656, -0.00733244,  0.        ],
        [-0.02300878, -0.01594689,  0.        ],
        ...,
        [ 0.16844292, -0.43168595,  0.        ],
        [ 0.16888118, -0.4615291 ,  0.        ],
        [ 0.1684947 , -0.49480793,  0.        ]]], dtype=float32))
[34]:
xplt.states.elements[0]
[34]:
StateData(name='Lagrange strain', dom=0, data=array([[[ 0.00000000e+00,  0.00000000e+00,  0.00000000e+00,
          0.00000000e+00,  0.00000000e+00,  0.00000000e+00],
        [ 0.00000000e+00,  0.00000000e+00,  0.00000000e+00,
          0.00000000e+00,  0.00000000e+00,  0.00000000e+00],
        [ 0.00000000e+00,  0.00000000e+00,  0.00000000e+00,
          0.00000000e+00,  0.00000000e+00,  0.00000000e+00],
        ...,
        [ 0.00000000e+00,  0.00000000e+00,  0.00000000e+00,
          0.00000000e+00,  0.00000000e+00,  0.00000000e+00],
        [ 0.00000000e+00,  0.00000000e+00,  0.00000000e+00,
          0.00000000e+00,  0.00000000e+00,  0.00000000e+00],
        [ 0.00000000e+00,  0.00000000e+00,  0.00000000e+00,
          0.00000000e+00,  0.00000000e+00,  0.00000000e+00]],

       [[-2.37153228e-02,  4.96507617e-08,  9.25132260e-03,
         -7.05668749e-03,  8.77791172e-05,  8.75434343e-05],
        [-1.55513035e-02,  4.19573579e-03,  6.64580334e-03,
         -8.76835082e-04, -4.10028315e-06,  3.52151233e-06],
        [-2.16151774e-02,  4.19573579e-03,  6.46888278e-03,
         -3.71682271e-03, -1.87471978e-05, -3.00115735e-05],
        ...,
        [ 5.26419561e-03, -1.79060502e-03, -1.32319576e-03,
         -1.09630043e-03,  1.35697308e-04, -9.22947293e-05],
        [ 3.76348314e-03, -1.79060502e-03, -6.20446925e-04,
         -1.20302581e-03,  3.88543667e-05,  1.51752083e-05],
        [ 4.56833979e-03, -4.31854278e-03, -2.29850310e-04,
         -2.06456636e-03, -6.20934181e-04,  2.85940798e-04]],

       [[-4.25659008e-02,  1.59901958e-07,  1.65927690e-02,
         -1.27380788e-02,  1.59318617e-04,  1.58670024e-04],
        [-2.77379584e-02,  7.46000186e-03,  1.18833957e-02,
         -1.59418711e-03, -7.71454506e-06,  5.81513541e-06],
        [-3.85503881e-02,  7.46000186e-03,  1.15433391e-02,
         -6.74896734e-03, -3.39396829e-05, -5.31097103e-05],
        ...,
        [ 1.04284454e-02, -3.59612191e-03, -2.62235897e-03,
         -1.94375543e-03,  2.62485846e-04, -1.94640495e-04],
        [ 7.60030141e-03, -3.59612191e-03, -1.25896640e-03,
         -2.19966378e-03,  7.99982663e-05,  2.65306644e-05],
        [ 9.37342457e-03, -8.74353014e-03, -5.09420410e-04,
         -3.88680655e-03, -1.21673197e-03,  6.25563203e-04]],

       ...,

       [[-1.03454977e-01,  1.02291756e-06,  4.00860123e-02,
         -3.26714478e-02,  4.17361414e-04,  4.09782020e-04],
        [-6.36420026e-02,  1.68143846e-02,  2.79301554e-02,
         -4.35827812e-03, -2.22580402e-05,  9.12787800e-06],
        [-8.96661356e-02,  1.68143846e-02,  2.68445406e-02,
         -1.82178933e-02, -8.80838707e-05, -1.21186145e-04],
        ...,
        [ 3.98490652e-02, -1.43709024e-02, -9.90649406e-03,
         -4.42765653e-03,  9.21171857e-04, -8.83480243e-04],
        [ 3.06131430e-02, -1.43709024e-02, -5.02734119e-03,
         -6.15815260e-03,  3.49141104e-04,  5.58151733e-05],
        [ 3.98919024e-02, -3.57674919e-02, -2.53269938e-03,
         -1.22170681e-02, -4.47947998e-03,  3.11689964e-03]],

       [[-1.08881332e-01,  1.15577188e-06,  4.21402082e-02,
         -3.47481929e-02,  4.44912352e-04,  4.35584836e-04],
        [-6.61997572e-02,  1.74349323e-02,  2.92007625e-02,
         -4.69052047e-03, -2.38736247e-05,  8.95932772e-06],
        [-9.36254486e-02,  1.74349323e-02,  2.80210096e-02,
         -1.95591226e-02, -9.38163721e-05, -1.26437721e-04],
        ...,
        [ 4.45734523e-02, -1.61647107e-02, -1.10511193e-02,
         -4.62550670e-03,  1.02171360e-03, -1.00175163e-03],
        [ 3.43595631e-02, -1.61647107e-02, -5.62138064e-03,
         -6.66673388e-03,  3.96208838e-04,  5.73999459e-05],
        [ 4.50464040e-02, -4.02821973e-02, -2.88027129e-03,
         -1.33820577e-02, -4.99419915e-03,  3.56116169e-03]],

       [[-1.13529623e-01,  1.28303191e-06,  4.38876860e-02,
         -3.66229638e-02,  4.69880557e-04,  4.58679424e-04],
        [-6.81892633e-02,  1.79007482e-02,  3.02395485e-02,
         -5.00305137e-03, -2.53279559e-05,  8.70800432e-06],
        [-9.68414694e-02,  1.79007482e-02,  2.89726872e-02,
         -2.08097380e-02, -9.90094777e-05, -1.30736560e-04],
        ...,
        [ 4.92642596e-02, -1.79609023e-02, -1.21797724e-02,
         -4.79171658e-03,  1.12077757e-03, -1.11997582e-03],
        [ 3.80752347e-02, -1.79609023e-02, -6.20363094e-03,
         -7.15704774e-03,  4.43766272e-04,  5.84696500e-05],
        [ 5.01992963e-02, -4.47962992e-02, -3.22508742e-03,
         -1.45127149e-02, -5.50433807e-03,  4.00813343e-03]]],
      dtype=float32))

FEBio Container

[35]:
from febio_python import FEBioContainer

container = FEBioContainer(feb=feb,
                           xplt=xplt
                           )

Can also load from a file directly:

[36]:
container = FEBioContainer(feb="plane_mesh_v25.feb",
                           xplt="plane_mesh_v25.xplt"
                           )
[37]:
container.feb
[37]:
Feb25(2.5):
-> Module: solid
-> Control: 15
-> Material: 1
-> Globals: 1
-> Geometry: 5
-> Boundary: 2
-> Loads: 2
-> Discrete: 0
-> LoadData: 1
-> Output: 1
-> MeshData: 2
[38]:
container.xplt
[38]:
Xplt object [1692092081728]:
=== xplt_mesh: ===
-Nodes:
--->None: (231, 3)
-Elements:
--->None: (400, 3)
=== States: ===
-Nodes:
--->displacement: (11, 231, 3)
-Elements:
--->Lagrange strain: (11, 400, 6)
--->stress: (11, 400, 6)
-Surfaces:
[39]:
container.nodes
[39]:
[Nodes(name='plane', coordinates=array([[-1. , -0.5,  0. ],
        [-0.9, -0.5,  0. ],
        [-0.8, -0.5,  0. ],
        [-0.7, -0.5,  0. ],
        [-0.6, -0.5,  0. ],
        [-0.5, -0.5,  0. ],
        [-0.4, -0.5,  0. ],
        [-0.3, -0.5,  0. ],
        [-0.2, -0.5,  0. ],
        [-0.1, -0.5,  0. ],
        [ 0. , -0.5,  0. ],
        [ 0.1, -0.5,  0. ],
        [ 0.2, -0.5,  0. ],
        [ 0.3, -0.5,  0. ],
        [ 0.4, -0.5,  0. ],
        [ 0.5, -0.5,  0. ],
        [ 0.6, -0.5,  0. ],
        [ 0.7, -0.5,  0. ],
        [ 0.8, -0.5,  0. ],
        [ 0.9, -0.5,  0. ],
        [ 1. , -0.5,  0. ],
        [-1. , -0.4,  0. ],
        [-0.9, -0.4,  0. ],
        [-0.8, -0.4,  0. ],
        [-0.7, -0.4,  0. ],
        [-0.6, -0.4,  0. ],
        [-0.5, -0.4,  0. ],
        [-0.4, -0.4,  0. ],
        [-0.3, -0.4,  0. ],
        [-0.2, -0.4,  0. ],
        [-0.1, -0.4,  0. ],
        [ 0. , -0.4,  0. ],
        [ 0.1, -0.4,  0. ],
        [ 0.2, -0.4,  0. ],
        [ 0.3, -0.4,  0. ],
        [ 0.4, -0.4,  0. ],
        [ 0.5, -0.4,  0. ],
        [ 0.6, -0.4,  0. ],
        [ 0.7, -0.4,  0. ],
        [ 0.8, -0.4,  0. ],
        [ 0.9, -0.4,  0. ],
        [ 1. , -0.4,  0. ],
        [-1. , -0.3,  0. ],
        [-0.9, -0.3,  0. ],
        [-0.8, -0.3,  0. ],
        [-0.7, -0.3,  0. ],
        [-0.6, -0.3,  0. ],
        [-0.5, -0.3,  0. ],
        [-0.4, -0.3,  0. ],
        [-0.3, -0.3,  0. ],
        [-0.2, -0.3,  0. ],
        [-0.1, -0.3,  0. ],
        [ 0. , -0.3,  0. ],
        [ 0.1, -0.3,  0. ],
        [ 0.2, -0.3,  0. ],
        [ 0.3, -0.3,  0. ],
        [ 0.4, -0.3,  0. ],
        [ 0.5, -0.3,  0. ],
        [ 0.6, -0.3,  0. ],
        [ 0.7, -0.3,  0. ],
        [ 0.8, -0.3,  0. ],
        [ 0.9, -0.3,  0. ],
        [ 1. , -0.3,  0. ],
        [-1. , -0.2,  0. ],
        [-0.9, -0.2,  0. ],
        [-0.8, -0.2,  0. ],
        [-0.7, -0.2,  0. ],
        [-0.6, -0.2,  0. ],
        [-0.5, -0.2,  0. ],
        [-0.4, -0.2,  0. ],
        [-0.3, -0.2,  0. ],
        [-0.2, -0.2,  0. ],
        [-0.1, -0.2,  0. ],
        [ 0. , -0.2,  0. ],
        [ 0.1, -0.2,  0. ],
        [ 0.2, -0.2,  0. ],
        [ 0.3, -0.2,  0. ],
        [ 0.4, -0.2,  0. ],
        [ 0.5, -0.2,  0. ],
        [ 0.6, -0.2,  0. ],
        [ 0.7, -0.2,  0. ],
        [ 0.8, -0.2,  0. ],
        [ 0.9, -0.2,  0. ],
        [ 1. , -0.2,  0. ],
        [-1. , -0.1,  0. ],
        [-0.9, -0.1,  0. ],
        [-0.8, -0.1,  0. ],
        [-0.7, -0.1,  0. ],
        [-0.6, -0.1,  0. ],
        [-0.5, -0.1,  0. ],
        [-0.4, -0.1,  0. ],
        [-0.3, -0.1,  0. ],
        [-0.2, -0.1,  0. ],
        [-0.1, -0.1,  0. ],
        [ 0. , -0.1,  0. ],
        [ 0.1, -0.1,  0. ],
        [ 0.2, -0.1,  0. ],
        [ 0.3, -0.1,  0. ],
        [ 0.4, -0.1,  0. ],
        [ 0.5, -0.1,  0. ],
        [ 0.6, -0.1,  0. ],
        [ 0.7, -0.1,  0. ],
        [ 0.8, -0.1,  0. ],
        [ 0.9, -0.1,  0. ],
        [ 1. , -0.1,  0. ],
        [-1. ,  0. ,  0. ],
        [-0.9,  0. ,  0. ],
        [-0.8,  0. ,  0. ],
        [-0.7,  0. ,  0. ],
        [-0.6,  0. ,  0. ],
        [-0.5,  0. ,  0. ],
        [-0.4,  0. ,  0. ],
        [-0.3,  0. ,  0. ],
        [-0.2,  0. ,  0. ],
        [-0.1,  0. ,  0. ],
        [ 0. ,  0. ,  0. ],
        [ 0.1,  0. ,  0. ],
        [ 0.2,  0. ,  0. ],
        [ 0.3,  0. ,  0. ],
        [ 0.4,  0. ,  0. ],
        [ 0.5,  0. ,  0. ],
        [ 0.6,  0. ,  0. ],
        [ 0.7,  0. ,  0. ],
        [ 0.8,  0. ,  0. ],
        [ 0.9,  0. ,  0. ],
        [ 1. ,  0. ,  0. ],
        [-1. ,  0.1,  0. ],
        [-0.9,  0.1,  0. ],
        [-0.8,  0.1,  0. ],
        [-0.7,  0.1,  0. ],
        [-0.6,  0.1,  0. ],
        [-0.5,  0.1,  0. ],
        [-0.4,  0.1,  0. ],
        [-0.3,  0.1,  0. ],
        [-0.2,  0.1,  0. ],
        [-0.1,  0.1,  0. ],
        [ 0. ,  0.1,  0. ],
        [ 0.1,  0.1,  0. ],
        [ 0.2,  0.1,  0. ],
        [ 0.3,  0.1,  0. ],
        [ 0.4,  0.1,  0. ],
        [ 0.5,  0.1,  0. ],
        [ 0.6,  0.1,  0. ],
        [ 0.7,  0.1,  0. ],
        [ 0.8,  0.1,  0. ],
        [ 0.9,  0.1,  0. ],
        [ 1. ,  0.1,  0. ],
        [-1. ,  0.2,  0. ],
        [-0.9,  0.2,  0. ],
        [-0.8,  0.2,  0. ],
        [-0.7,  0.2,  0. ],
        [-0.6,  0.2,  0. ],
        [-0.5,  0.2,  0. ],
        [-0.4,  0.2,  0. ],
        [-0.3,  0.2,  0. ],
        [-0.2,  0.2,  0. ],
        [-0.1,  0.2,  0. ],
        [ 0. ,  0.2,  0. ],
        [ 0.1,  0.2,  0. ],
        [ 0.2,  0.2,  0. ],
        [ 0.3,  0.2,  0. ],
        [ 0.4,  0.2,  0. ],
        [ 0.5,  0.2,  0. ],
        [ 0.6,  0.2,  0. ],
        [ 0.7,  0.2,  0. ],
        [ 0.8,  0.2,  0. ],
        [ 0.9,  0.2,  0. ],
        [ 1. ,  0.2,  0. ],
        [-1. ,  0.3,  0. ],
        [-0.9,  0.3,  0. ],
        [-0.8,  0.3,  0. ],
        [-0.7,  0.3,  0. ],
        [-0.6,  0.3,  0. ],
        [-0.5,  0.3,  0. ],
        [-0.4,  0.3,  0. ],
        [-0.3,  0.3,  0. ],
        [-0.2,  0.3,  0. ],
        [-0.1,  0.3,  0. ],
        [ 0. ,  0.3,  0. ],
        [ 0.1,  0.3,  0. ],
        [ 0.2,  0.3,  0. ],
        [ 0.3,  0.3,  0. ],
        [ 0.4,  0.3,  0. ],
        [ 0.5,  0.3,  0. ],
        [ 0.6,  0.3,  0. ],
        [ 0.7,  0.3,  0. ],
        [ 0.8,  0.3,  0. ],
        [ 0.9,  0.3,  0. ],
        [ 1. ,  0.3,  0. ],
        [-1. ,  0.4,  0. ],
        [-0.9,  0.4,  0. ],
        [-0.8,  0.4,  0. ],
        [-0.7,  0.4,  0. ],
        [-0.6,  0.4,  0. ],
        [-0.5,  0.4,  0. ],
        [-0.4,  0.4,  0. ],
        [-0.3,  0.4,  0. ],
        [-0.2,  0.4,  0. ],
        [-0.1,  0.4,  0. ],
        [ 0. ,  0.4,  0. ],
        [ 0.1,  0.4,  0. ],
        [ 0.2,  0.4,  0. ],
        [ 0.3,  0.4,  0. ],
        [ 0.4,  0.4,  0. ],
        [ 0.5,  0.4,  0. ],
        [ 0.6,  0.4,  0. ],
        [ 0.7,  0.4,  0. ],
        [ 0.8,  0.4,  0. ],
        [ 0.9,  0.4,  0. ],
        [ 1. ,  0.4,  0. ],
        [-1. ,  0.5,  0. ],
        [-0.9,  0.5,  0. ],
        [-0.8,  0.5,  0. ],
        [-0.7,  0.5,  0. ],
        [-0.6,  0.5,  0. ],
        [-0.5,  0.5,  0. ],
        [-0.4,  0.5,  0. ],
        [-0.3,  0.5,  0. ],
        [-0.2,  0.5,  0. ],
        [-0.1,  0.5,  0. ],
        [ 0. ,  0.5,  0. ],
        [ 0.1,  0.5,  0. ],
        [ 0.2,  0.5,  0. ],
        [ 0.3,  0.5,  0. ],
        [ 0.4,  0.5,  0. ],
        [ 0.5,  0.5,  0. ],
        [ 0.6,  0.5,  0. ],
        [ 0.7,  0.5,  0. ],
        [ 0.8,  0.5,  0. ],
        [ 0.9,  0.5,  0. ],
        [ 1. ,  0.5,  0. ]], dtype=float32), ids=array([  0,   1,   2,   3,   4,   5,   6,   7,   8,   9,  10,  11,  12,
         13,  14,  15,  16,  17,  18,  19,  20,  21,  22,  23,  24,  25,
         26,  27,  28,  29,  30,  31,  32,  33,  34,  35,  36,  37,  38,
         39,  40,  41,  42,  43,  44,  45,  46,  47,  48,  49,  50,  51,
         52,  53,  54,  55,  56,  57,  58,  59,  60,  61,  62,  63,  64,
         65,  66,  67,  68,  69,  70,  71,  72,  73,  74,  75,  76,  77,
         78,  79,  80,  81,  82,  83,  84,  85,  86,  87,  88,  89,  90,
         91,  92,  93,  94,  95,  96,  97,  98,  99, 100, 101, 102, 103,
        104, 105, 106, 107, 108, 109, 110, 111, 112, 113, 114, 115, 116,
        117, 118, 119, 120, 121, 122, 123, 124, 125, 126, 127, 128, 129,
        130, 131, 132, 133, 134, 135, 136, 137, 138, 139, 140, 141, 142,
        143, 144, 145, 146, 147, 148, 149, 150, 151, 152, 153, 154, 155,
        156, 157, 158, 159, 160, 161, 162, 163, 164, 165, 166, 167, 168,
        169, 170, 171, 172, 173, 174, 175, 176, 177, 178, 179, 180, 181,
        182, 183, 184, 185, 186, 187, 188, 189, 190, 191, 192, 193, 194,
        195, 196, 197, 198, 199, 200, 201, 202, 203, 204, 205, 206, 207,
        208, 209, 210, 211, 212, 213, 214, 215, 216, 217, 218, 219, 220,
        221, 222, 223, 224, 225, 226, 227, 228, 229, 230], dtype=int64))]

Ploting Feb or Xplt

[40]:
from febio_python.utils.pyvista_utils import febio_to_pyvista

grids_list = febio_to_pyvista(container)
len(grids_list)
[40]:
11
[41]:
grids_list[-1].plot(cpos='xy', show_edges=True, scalars="stress")
../../../../_images/modules_feb_spec_2_5_adding_basic_components_adding_basic_components_69_0.png

You can also quickly convert all cell data to node data:

[42]:
last_grid_as_nodal_data = grids_list[-1].cell_data_to_point_data()
last_grid_as_nodal_data.plot(cpos='xy', show_edges=True, scalars="stress")
../../../../_images/modules_feb_spec_2_5_adding_basic_components_adding_basic_components_71_0.png
[43]:
last_grid_as_nodal_data
[43]:
HeaderData Arrays
UnstructuredGridInformation
N Cells400
N Points231
X Bounds-1.000e+00, 1.168e+00
Y Bounds-9.501e-01, 5.000e-01
Z Bounds0.000e+00, 0.000e+00
N Arrays19
NameFieldTypeN CompMinMax
right_edgePointsbool1nannan
left_edgePointsbool1nannan
all_nodesPointsbool1nannan
nodal_load_right_edge_scalePointsfloat6418.209e+001.000e+02
nodal_loadPointsfloat643-2.500e+011.000e+02
fixPointsint3230.000e+001.000e+00
fix_shellPointsint3230.000e+001.000e+00
displacementPointsfloat323-4.948e-011.689e-01
domainPointsfloat3210.000e+000.000e+00
Element thicknessPointsfloat6431.000e-021.000e-02
mat_parameters:1Pointsfloat6423.000e-011.000e+00
Lagrange strainPointsfloat326-1.135e-012.052e-01
stressPointsfloat326-1.071e+052.798e+05
stress-normedPointsfloat3217.009e+032.953e+05
mat_parameters:1Fields1nannan
mat_type:1Fields1nannan
mat_name:1Fields1nannan
timestepFieldsfloat3211.000e+001.000e+00
lc_1Fieldsfloat3211.000e+001.000e+00

When using the FEBio container, we have access to both Feb and Xplt data. This means that we can retrieve nodal load, boundary conditions, etc. Here is a cool plot that we can do:

[44]:
plotter = pv.Plotter()
strain_xx = last_grid_as_nodal_data["Lagrange strain"][:, 0]
fixed_nodes = last_grid_as_nodal_data["fix"].sum(1)

plotter.add_mesh(last_grid_as_nodal_data,
                 scalars=strain_xx,
                 cmap="coolwarm",
                 show_edges=True,
                 scalar_bar_args={"title": "Strain - XX"})
plotter.add_mesh(last_grid_as_nodal_data.points,
                 scalars=fixed_nodes,
                 cmap="viridis",
                 style="points",
                 point_size=10,
                 render_points_as_spheres=True, show_scalar_bar=False)
plotter.add_arrows(last_grid_as_nodal_data.points,
                   last_grid_as_nodal_data["nodal_load"],
                   mag=5e-3, # This controls the mag of the arrows. Not the actual load. There may be a better way to control this.
                   show_scalar_bar=False,
                   color="orange")
plotter.show(cpos="xy")

../../../../_images/modules_feb_spec_2_5_adding_basic_components_adding_basic_components_74_0.png