Adding Basic Components (v4.0)

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
from febio_python.feb import Feb40 # for type hinting

feb: Feb40 = Feb(version=4.0)
feb
[1]:
Feb40(4.0):
-> Module: Unknown
-> Control: 0
-> Material: 0
-> Globals: 0
-> Mesh: 0
-> MeshDomains: 0
-> MeshData: 0
-> Boundary: 0
-> Loads: 0
-> Discrete: 0
-> LoadData: 0
-> Output: 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, P=0, R=0, Fc=0) # default values
feb.setup_controls(analysis="STATIC") # here, you can change basic settings. See docs for more info.
feb.setup_output(variables=["displacement", "Lagrange strain", "stress"]) # default values
qn_method {'attr_type': 'BFGS', 'max_ups': 10, 'max_buffer_size': 0, 'cycle_buffer': 1, 'cmax': 100000}

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

[3]:
feb
[3]:
Feb40(4.0):
-> Module: solid
-> Control: 12
-> Material: 0
-> Globals: 1
-> Mesh: 0
-> MeshDomains: 0
-> MeshData: 0
-> Boundary: 0
-> Loads: 0
-> Discrete: 0
-> Output: 1
-> LoadData: 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_4_0_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="Part1",
                    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 ‘Mesh’ items:

[9]:
feb
[9]:
Feb40(4.0):
-> Module: solid
-> Control: 12
-> Material: 0
-> Globals: 1
-> Mesh: 2
--> Nodes: 28
--> Elements: 18
-> MeshDomains: 0
-> MeshData: 0
-> Boundary: 0
-> Loads: 0
-> Discrete: 0
-> Output: 1
-> LoadData: 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_4_0_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="Part1",
                    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]:
Feb40(4.0):
-> Module: solid
-> Control: 12
-> Material: 0
-> Globals: 1
-> Mesh: 0
-> MeshDomains: 0
-> MeshData: 0
-> Boundary: 0
-> Loads: 0
-> Discrete: 0
-> Output: 1
-> LoadData: 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_4_0_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="Part1",
                    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 Mesh Domains

[19]:
from febio_python.core import ShellDomain

# create a shell domain
shell = ShellDomain(
    id=1,
    name="plane_elements", # this must match one of the element names
    mat="plane material", # this must match one of the material names
    type="elastic-shell", # type of shell element
    shell_normal_nodal = 1, # normal to the shell
    shell_thickness=0.01, # thickness of the shell
    )
feb.add_mesh_domains([shell])

Adding Load

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

[20]:
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
[20]:
Feb40(4.0):
-> Module: solid
-> Control: 12
-> Material: 0
-> Globals: 1
-> Mesh: 3
--> Nodes: 231
--> Elements: 400
-> MeshDomains: 1
-> MeshData: 0
-> Boundary: 0
-> Loads: 1
-> Discrete: 0
-> Output: 1
-> LoadData: 1

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

[21]:
# 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_4_0_adding_basic_components_adding_basic_components_40_0.png
[21]:
Feb40(4.0):
-> Module: solid
-> Control: 12
-> Material: 0
-> Globals: 1
-> Mesh: 3
--> Nodes: 231
--> Elements: 400
-> MeshDomains: 1
-> MeshData: 1
-> Boundary: 0
-> Loads: 2
-> Discrete: 0
-> Output: 1
-> LoadData: 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.

[22]:
from febio_python.core import ZeroDisplacementCondition, ZeroShellDisplacementCondition

# 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 = ZeroDisplacementCondition(dof="x,y,z", node_set="left_edge")
left_shell_fix_condition = ZeroShellDisplacementCondition(dof="sx,sy,sz", node_set="left_edge")
# we will
all_fix_condition = ZeroDisplacementCondition(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, left_shell_fix_condition, all_fix_condition])

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",
    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:

[25]:
# from febio_python.feb import run
# # Save the FEB file
# feb.write("plane_mesh_v40.feb")
# run("plane_mesh_v40.feb")

Reading XPLT

[26]:
from febio_python import Xplt

xplt = Xplt("plane_mesh_v40.xplt")
xplt
[26]:
Xplt object [2887623865632]:
=== 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

[27]:
xplt.nodes # list of Nodes objects
[27]:
[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([  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]))]
[28]:
xplt.elements # list of Elements objects
[28]:
[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)]
[29]:
xplt.states # list of States objects
[29]:
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.00220876, -0.00138848,  0.        ],
        [-0.0044777 , -0.00291381,  0.        ],
        ...,
        [ 0.03501934, -0.07832458,  0.        ],
        [ 0.03535859, -0.08439115,  0.        ],
        [ 0.03562101, -0.09075736,  0.        ]],

       [[ 0.        ,  0.        ,  0.        ],
        [-0.0040143 , -0.00251072,  0.        ],
        [-0.00812449, -0.00527578,  0.        ],
        ...,
        [ 0.06155637, -0.14323959,  0.        ],
        [ 0.06199993, -0.15412821,  0.        ],
        [ 0.06229447, -0.16563064,  0.        ]],

       ...,

       [[ 0.        ,  0.        ,  0.        ],
        [-0.0102166 , -0.00648645,  0.        ],
        [-0.02036276, -0.01388767,  0.        ],
        ...,
        [ 0.14947748, -0.3811367 ,  0.        ],
        [ 0.14983469, -0.40797347,  0.        ],
        [ 0.14956348, -0.4375061 ,  0.        ]],

       [[ 0.        ,  0.        ,  0.        ],
        [-0.01080656, -0.00690494,  0.        ],
        [-0.02147301, -0.01484216,  0.        ],
        ...,
        [ 0.15917969, -0.40684003,  0.        ],
        [ 0.15957509, -0.43523568,  0.        ],
        [ 0.15924646, -0.46669516,  0.        ]],

       [[ 0.        ,  0.        ,  0.        ],
        [-0.01131956, -0.00728358,  0.        ],
        [-0.02242204, -0.01571946,  0.        ],
        ...,
        [ 0.16826567, -0.4302698 ,  0.        ],
        [ 0.16872025, -0.46006113,  0.        ],
        [ 0.16835007, -0.49328837,  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.17677094e-02,  0.00000000e+00,  3.74327949e-03,
         -6.93582371e-03, -4.37995950e-05, -2.12121537e-04],
        [-1.44304177e-02,  4.09667706e-03,  5.80157526e-03,
         -1.21883943e-03,  5.60197332e-05, -6.63803366e-05],
        [-2.22795270e-02,  4.09667706e-03,  7.83721078e-03,
         -3.98436794e-03,  3.43527390e-05,  1.65543141e-04],
        ...,
        [ 5.26422681e-03, -1.79072132e-03, -1.32321590e-03,
         -1.09698484e-03,  1.35720460e-04, -9.22642721e-05],
        [ 3.76291759e-03, -1.79072132e-03, -6.20350009e-04,
         -1.20366516e-03,  3.88481822e-05,  1.51849181e-05],
        [ 4.56744991e-03, -4.31843102e-03, -2.29693338e-04,
         -2.06525391e-03, -6.21005311e-04,  2.85787100e-04]],

       [[-3.90582196e-02,  0.00000000e+00,  6.69534458e-03,
         -1.25417085e-02, -7.91431958e-05, -3.80462821e-04],
        [-2.57335734e-02,  7.29539804e-03,  1.03692589e-02,
         -2.23168731e-03,  1.01221529e-04, -1.18771379e-04],
        [-3.98102477e-02,  7.29539804e-03,  1.40054580e-02,
         -7.24549964e-03,  6.26781766e-05,  2.96983548e-04],
        ...,
        [ 1.04291877e-02, -3.59612866e-03, -2.62246048e-03,
         -1.94583589e-03,  2.62564688e-04, -1.94546839e-04],
        [ 7.59910978e-03, -3.59612866e-03, -1.25865580e-03,
         -2.20169220e-03,  7.99803020e-05,  2.65620838e-05],
        [ 9.37132537e-03, -8.74294899e-03, -5.08905330e-04,
         -3.88899446e-03, -1.21698214e-03,  6.25088112e-04]],

       ...,

       [[-9.49234068e-02,  0.00000000e+00,  1.60137974e-02,
         -3.24016511e-02, -2.04027048e-04, -9.26311302e-04],
        [-5.90482466e-02,  1.65700149e-02,  2.44058762e-02,
         -6.17906218e-03,  2.60786124e-04, -2.74899823e-04],
        [-9.34232473e-02,  1.65700149e-02,  3.28529291e-02,
         -1.95991900e-02,  1.65617355e-04,  7.16031471e-04],
        ...,
        [ 3.98563445e-02, -1.43719194e-02, -9.90771782e-03,
         -4.43919096e-03,  9.21695610e-04, -8.83027446e-04],
        [ 3.06086857e-02, -1.43719194e-02, -5.02582826e-03,
         -6.17014943e-03,  3.49059643e-04,  5.59888649e-05],
        [ 3.98842469e-02, -3.57679091e-02, -2.52990308e-03,
         -1.22303627e-02, -4.48134402e-03,  3.11454432e-03]],

       [[-9.99256521e-02,  0.00000000e+00,  1.68171786e-02,
         -3.44921164e-02, -2.17095192e-04, -9.75765346e-04],
        [-6.14328049e-02,  1.72004160e-02,  2.55326405e-02,
         -6.64548576e-03,  2.77563930e-04, -2.86216731e-04],
        [-9.76665765e-02,  1.72004160e-02,  3.43420096e-02,
         -2.10359637e-02,  1.76507951e-04,  7.52359221e-04],
        ...,
        [ 4.45817485e-02, -1.61659475e-02, -1.10525191e-02,
         -4.63817036e-03,  1.02229952e-03, -1.00125966e-03],
        [ 3.43548208e-02, -1.61659475e-02, -5.61971823e-03,
         -6.68002805e-03,  3.96120304e-04,  5.75903759e-05],
        [ 4.50384207e-02, -4.02831174e-02, -2.87718838e-03,
         -1.33968135e-02, -4.99630766e-03,  3.55860265e-03]],

       [[-1.04221813e-01,  0.00000000e+00,  1.74986590e-02,
         -3.63835469e-02, -2.28882986e-04, -1.01843383e-03],
        [-6.32932782e-02,  1.76788010e-02,  2.64601391e-02,
         -7.07987975e-03,  2.92746816e-04, -2.95078993e-04],
        [-1.01138130e-01,  1.76788010e-02,  3.55594903e-02,
         -2.23720279e-02,  1.86323217e-04,  7.83163880e-04],
        ...,
        [ 4.92735021e-02, -1.79623552e-02, -1.21813323e-02,
         -4.80537629e-03,  1.12142053e-03, -1.11944962e-03],
        [ 3.80702429e-02, -1.79623552e-02, -6.20182930e-03,
         -7.17152003e-03,  4.43671306e-04,  5.86747301e-05],
        [ 5.01911044e-02, -4.47977670e-02, -3.22174374e-03,
         -1.45287961e-02, -5.50667616e-03,  4.00539767e-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.64484258e+04, -1.04529355e+04, -7.70807080e+03,
         -4.94026172e+03, -2.88648472e+01, -1.48734741e+02],
        [-1.33989375e+04,  5.57711304e+02,  1.88573279e+03,
         -7.94711731e+02,  4.41828728e+01, -4.59821739e+01],
        [-2.23616953e+04, -2.78293970e+03,  6.36981621e+01,
         -2.71356543e+03,  2.35810680e+01,  1.20599007e+02],
        ...,
        [ 5.21622900e+03, -1.42441969e+01,  2.21691864e+02,
         -1.16859033e+03,  1.08122070e+02, -6.51865768e+01],
        [ 3.57391211e+03, -4.66856659e+02,  3.02063782e+02,
         -1.17826062e+03,  2.89683361e+01,  1.33193684e+01],
        [ 3.34432422e+03, -3.05434204e+03, -1.65373367e+02,
         -1.99538770e+03, -4.88635071e+02,  1.92100922e+02]],

       [[-4.64621133e+04, -1.88682031e+04, -1.41482842e+04,
         -8.36959863e+03, -4.55878372e+01, -2.46629288e+02],
        [-2.34484590e+04,  1.02387891e+03,  3.44388110e+03,
         -1.25443335e+03,  8.13408356e+01, -7.52409744e+01],
        [-3.89257422e+04, -4.96432617e+03,  1.06509041e+02,
         -4.44692822e+03,  3.91465225e+01,  2.07021591e+02],
        ...,
        [ 1.01885801e+04,  1.16335869e+02,  4.08291443e+02,
         -2.64133667e+03,  2.15579498e+02, -1.29476761e+02],
        [ 7.06979736e+03, -7.04508606e+02,  6.12370789e+02,
         -2.58942432e+03,  5.83515320e+01,  2.63631763e+01],
        [ 6.61788672e+03, -5.68839160e+03, -3.21880890e+02,
         -4.44795850e+03, -9.77329163e+02,  3.83126587e+02]],

       ...,

       [[-1.04840055e+05, -4.70084375e+04, -3.73781641e+04,
         -1.69710664e+04, -6.83919601e+01, -4.50002289e+02],
        [-5.06672891e+04,  2.73838452e+03,  9.09365039e+03,
         -1.89116357e+03,  2.20604797e+02, -1.26405624e+02],
        [-8.42736016e+04, -1.19750898e+04, -8.18701248e+01,
         -8.20814648e+03,  7.18703384e+01,  4.27924744e+02],
        ...,
        [ 3.76055586e+04,  2.70066382e+03,  1.32220825e+03,
         -1.40359697e+04,  8.46136597e+02, -5.05699005e+02],
        [ 2.67870898e+04,  4.54238922e+02,  2.55717505e+03,
         -1.31917803e+04,  2.39443359e+02,  1.01728439e+02],
        [ 2.51188770e+04, -1.50789824e+04, -1.02254968e+03,
         -2.35351816e+04, -3.90821924e+03,  1.51564709e+03]],

       [[-1.09579477e+05, -4.95953594e+04, -3.96914570e+04,
         -1.76262246e+04, -6.85027390e+01, -4.60733307e+02],
        [-5.24544531e+04,  2.91616431e+03,  9.67994531e+03,
         -1.91455396e+03,  2.35612656e+02, -1.27949501e+02],
        [-8.75843594e+04, -1.25930576e+04, -1.67251129e+02,
         -8.49706055e+03,  7.33836594e+01,  4.42672363e+02],
        ...,
        [ 4.20180820e+04,  3.29068579e+03,  1.46246997e+03,
         -1.61456777e+04,  9.49143250e+02, -5.67585632e+02],
        [ 2.99093887e+04,  8.65932373e+02,  2.87275732e+03,
         -1.51393564e+04,  2.70139343e+02,  1.13898109e+02],
        [ 2.80943340e+04, -1.58814502e+04, -1.11977356e+03,
         -2.70906562e+04, -4.39748926e+03,  1.70320874e+03]],

       [[-1.13583438e+05, -5.18183477e+04, -4.17205117e+04,
         -1.81947441e+04, -6.83979645e+01, -4.69166656e+02],
        [-5.38241758e+04,  3.07517920e+03,  1.02040029e+04,
         -1.94041626e+03,  2.49166428e+02, -1.28910309e+02],
        [-9.02867969e+04, -1.31167275e+04, -2.62930817e+02,
         -8.76461328e+03,  7.44312668e+01,  4.54420135e+02],
        ...,
        [ 4.64248828e+04,  3.90854077e+03,  1.60259070e+03,
         -1.82951836e+04,  1.05164062e+03, -6.29355957e+02],
        [ 3.29966719e+04,  1.31243579e+03,  3.18418384e+03,
         -1.71187715e+04,  3.00957672e+02,  1.25958168e+02],
        [ 3.10599844e+04, -1.65141113e+04, -1.21385217e+03,
         -3.07127969e+04, -4.88739697e+03,  1.89066504e+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:

[30]:
xplt.states.nodes[0]
[30]:
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.00220876, -0.00138848,  0.        ],
        [-0.0044777 , -0.00291381,  0.        ],
        ...,
        [ 0.03501934, -0.07832458,  0.        ],
        [ 0.03535859, -0.08439115,  0.        ],
        [ 0.03562101, -0.09075736,  0.        ]],

       [[ 0.        ,  0.        ,  0.        ],
        [-0.0040143 , -0.00251072,  0.        ],
        [-0.00812449, -0.00527578,  0.        ],
        ...,
        [ 0.06155637, -0.14323959,  0.        ],
        [ 0.06199993, -0.15412821,  0.        ],
        [ 0.06229447, -0.16563064,  0.        ]],

       ...,

       [[ 0.        ,  0.        ,  0.        ],
        [-0.0102166 , -0.00648645,  0.        ],
        [-0.02036276, -0.01388767,  0.        ],
        ...,
        [ 0.14947748, -0.3811367 ,  0.        ],
        [ 0.14983469, -0.40797347,  0.        ],
        [ 0.14956348, -0.4375061 ,  0.        ]],

       [[ 0.        ,  0.        ,  0.        ],
        [-0.01080656, -0.00690494,  0.        ],
        [-0.02147301, -0.01484216,  0.        ],
        ...,
        [ 0.15917969, -0.40684003,  0.        ],
        [ 0.15957509, -0.43523568,  0.        ],
        [ 0.15924646, -0.46669516,  0.        ]],

       [[ 0.        ,  0.        ,  0.        ],
        [-0.01131956, -0.00728358,  0.        ],
        [-0.02242204, -0.01571946,  0.        ],
        ...,
        [ 0.16826567, -0.4302698 ,  0.        ],
        [ 0.16872025, -0.46006113,  0.        ],
        [ 0.16835007, -0.49328837,  0.        ]]], dtype=float32))
[31]:
xplt.states.elements[0]
[31]:
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.17677094e-02,  0.00000000e+00,  3.74327949e-03,
         -6.93582371e-03, -4.37995950e-05, -2.12121537e-04],
        [-1.44304177e-02,  4.09667706e-03,  5.80157526e-03,
         -1.21883943e-03,  5.60197332e-05, -6.63803366e-05],
        [-2.22795270e-02,  4.09667706e-03,  7.83721078e-03,
         -3.98436794e-03,  3.43527390e-05,  1.65543141e-04],
        ...,
        [ 5.26422681e-03, -1.79072132e-03, -1.32321590e-03,
         -1.09698484e-03,  1.35720460e-04, -9.22642721e-05],
        [ 3.76291759e-03, -1.79072132e-03, -6.20350009e-04,
         -1.20366516e-03,  3.88481822e-05,  1.51849181e-05],
        [ 4.56744991e-03, -4.31843102e-03, -2.29693338e-04,
         -2.06525391e-03, -6.21005311e-04,  2.85787100e-04]],

       [[-3.90582196e-02,  0.00000000e+00,  6.69534458e-03,
         -1.25417085e-02, -7.91431958e-05, -3.80462821e-04],
        [-2.57335734e-02,  7.29539804e-03,  1.03692589e-02,
         -2.23168731e-03,  1.01221529e-04, -1.18771379e-04],
        [-3.98102477e-02,  7.29539804e-03,  1.40054580e-02,
         -7.24549964e-03,  6.26781766e-05,  2.96983548e-04],
        ...,
        [ 1.04291877e-02, -3.59612866e-03, -2.62246048e-03,
         -1.94583589e-03,  2.62564688e-04, -1.94546839e-04],
        [ 7.59910978e-03, -3.59612866e-03, -1.25865580e-03,
         -2.20169220e-03,  7.99803020e-05,  2.65620838e-05],
        [ 9.37132537e-03, -8.74294899e-03, -5.08905330e-04,
         -3.88899446e-03, -1.21698214e-03,  6.25088112e-04]],

       ...,

       [[-9.49234068e-02,  0.00000000e+00,  1.60137974e-02,
         -3.24016511e-02, -2.04027048e-04, -9.26311302e-04],
        [-5.90482466e-02,  1.65700149e-02,  2.44058762e-02,
         -6.17906218e-03,  2.60786124e-04, -2.74899823e-04],
        [-9.34232473e-02,  1.65700149e-02,  3.28529291e-02,
         -1.95991900e-02,  1.65617355e-04,  7.16031471e-04],
        ...,
        [ 3.98563445e-02, -1.43719194e-02, -9.90771782e-03,
         -4.43919096e-03,  9.21695610e-04, -8.83027446e-04],
        [ 3.06086857e-02, -1.43719194e-02, -5.02582826e-03,
         -6.17014943e-03,  3.49059643e-04,  5.59888649e-05],
        [ 3.98842469e-02, -3.57679091e-02, -2.52990308e-03,
         -1.22303627e-02, -4.48134402e-03,  3.11454432e-03]],

       [[-9.99256521e-02,  0.00000000e+00,  1.68171786e-02,
         -3.44921164e-02, -2.17095192e-04, -9.75765346e-04],
        [-6.14328049e-02,  1.72004160e-02,  2.55326405e-02,
         -6.64548576e-03,  2.77563930e-04, -2.86216731e-04],
        [-9.76665765e-02,  1.72004160e-02,  3.43420096e-02,
         -2.10359637e-02,  1.76507951e-04,  7.52359221e-04],
        ...,
        [ 4.45817485e-02, -1.61659475e-02, -1.10525191e-02,
         -4.63817036e-03,  1.02229952e-03, -1.00125966e-03],
        [ 3.43548208e-02, -1.61659475e-02, -5.61971823e-03,
         -6.68002805e-03,  3.96120304e-04,  5.75903759e-05],
        [ 4.50384207e-02, -4.02831174e-02, -2.87718838e-03,
         -1.33968135e-02, -4.99630766e-03,  3.55860265e-03]],

       [[-1.04221813e-01,  0.00000000e+00,  1.74986590e-02,
         -3.63835469e-02, -2.28882986e-04, -1.01843383e-03],
        [-6.32932782e-02,  1.76788010e-02,  2.64601391e-02,
         -7.07987975e-03,  2.92746816e-04, -2.95078993e-04],
        [-1.01138130e-01,  1.76788010e-02,  3.55594903e-02,
         -2.23720279e-02,  1.86323217e-04,  7.83163880e-04],
        ...,
        [ 4.92735021e-02, -1.79623552e-02, -1.21813323e-02,
         -4.80537629e-03,  1.12142053e-03, -1.11944962e-03],
        [ 3.80702429e-02, -1.79623552e-02, -6.20182930e-03,
         -7.17152003e-03,  4.43671306e-04,  5.86747301e-05],
        [ 5.01911044e-02, -4.47977670e-02, -3.22174374e-03,
         -1.45287961e-02, -5.50667616e-03,  4.00539767e-03]]],
      dtype=float32))

FEBio Container

[32]:
from febio_python import FEBioContainer

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

Can also load from a file directly:

[33]:
container = FEBioContainer(feb="plane_mesh_v40.feb",
                           xplt="plane_mesh_v40.xplt")
feb:  plane_mesh_v40.feb
[34]:
container.feb
[34]:
Feb40(4.0):
-> Module: solid
-> Control: 12
-> Material: 1
-> Globals: 1
-> Mesh: 5
--> Nodes: 231
--> Elements: 400
-> MeshDomains: 1
-> MeshData: 2
-> Boundary: 3
-> Loads: 2
-> Discrete: 0
-> Output: 1
-> LoadData: 1
[35]:
container.xplt
[35]:
Xplt object [2887623868368]:
=== 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:
[36]:
container.nodes
[36]:
[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

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

grids_list = febio_to_pyvista(container)
len(grids_list)
C:\Users\igorp\AppData\Roaming\Python\Python310\site-packages\scipy\__init__.py:146: UserWarning: A NumPy version >=1.17.3 and <1.25.0 is required for this version of SciPy (detected version 1.26.3
  warnings.warn(f"A NumPy version >={np_minversion} and <{np_maxversion}"
[37]:
11
[38]:
grids_list[-1].plot(cpos='xy', show_edges=True, scalars="stress")
../../../../_images/modules_feb_spec_4_0_adding_basic_components_adding_basic_components_70_0.png

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

[39]:
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_4_0_adding_basic_components_adding_basic_components_72_0.png
[40]:
last_grid_as_nodal_data
[40]:
HeaderData Arrays
UnstructuredGridInformation
N Cells400
N Points231
X Bounds-1.000e+00, 1.168e+00
Y Bounds-9.488e-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+000.000e+00
fix_shellPointsint3230.000e+000.000e+00
displacementPointsfloat323-4.933e-011.687e-01
domainPointsfloat3210.000e+000.000e+00
Element thicknessPointsfloat6431.000e-021.000e-02
mat_parameters:1Pointsfloat6423.000e-011.000e+00
Lagrange strainPointsfloat326-1.042e-012.035e-01
stressPointsfloat326-1.136e+052.762e+05
stress-normedPointsfloat3216.978e+032.955e+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:

[41]:
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_4_0_adding_basic_components_adding_basic_components_75_0.png