.. DO NOT EDIT. .. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY. .. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE: .. "gallery/elements/solid187/example_solid187.py" .. LINE NUMBERS ARE GIVEN BELOW. .. only:: html .. note:: :class: sphx-glr-download-link-note :ref:`Go to the end ` to download the full example code. .. rst-class:: sphx-glr-example-title .. _sphx_glr_gallery_elements_solid187_example_solid187.py: .. _ref_solid187_example: TET10 — prescribed uniform strain on a 10-node tet ===================================================== Single TET10 (10-node quadratic tet) driven by a fully-prescribed displacement field that matches a uniform strain state. Under all-Dirichlet BCs the solve has a unique solution equal to the prescribed field, so the recovered strain tensor matches the analytical one to machine precision. Useful as a self-contained smoke test that ``TET10`` + :meth:`femorph_solver.Model.eel` are in lockstep. .. GENERATED FROM PYTHON SOURCE LINES 14-24 .. code-block:: Python from __future__ import annotations import numpy as np import pyvista as pv from vtkmodules.util.vtkConstants import VTK_QUADRATIC_TETRA import femorph_solver from femorph_solver import ELEMENTS .. GENERATED FROM PYTHON SOURCE LINES 25-29 Reference 10-node unit tet -------------------------- Corners at (0,0,0), (1,0,0), (0,1,0), (0,0,1); mid-edge nodes in VTK_QUADRATIC_TETRA order halfway along each edge. .. GENERATED FROM PYTHON SOURCE LINES 29-54 .. code-block:: Python E = 2.1e11 # Pa NU = 0.30 EPS_XX = 5.0e-4 EPS_YY = -NU * EPS_XX # Poisson contraction I = np.array([0.0, 0.0, 0.0]) J = np.array([1.0, 0.0, 0.0]) K = np.array([0.0, 1.0, 0.0]) L = np.array([0.0, 0.0, 1.0]) coords = np.array( [ I, J, K, L, 0.5 * (I + J), # 5 0.5 * (J + K), # 6 0.5 * (K + I), # 7 0.5 * (I + L), # 8 0.5 * (J + L), # 9 0.5 * (K + L), # 10 ], dtype=np.float64, ) .. GENERATED FROM PYTHON SOURCE LINES 55-58 Build the model + prescribed displacement ----------------------------------------- Apply ``u = (εxx·x, εyy·y, εyy·z)`` at every node. .. GENERATED FROM PYTHON SOURCE LINES 58-70 .. code-block:: Python cells = np.concatenate([[10], np.arange(10, dtype=np.int64)]) cell_types = np.array([VTK_QUADRATIC_TETRA], dtype=np.uint8) grid = pv.UnstructuredGrid(cells, cell_types, coords) m = femorph_solver.Model.from_grid(grid) m.assign(ELEMENTS.TET10, material={"EX": E, "PRXY": NU}) for i, (x, y, z) in enumerate(coords, start=1): m.fix(nodes=[i], dof="UX", value=EPS_XX * x) m.fix(nodes=[i], dof="UY", value=EPS_YY * y) m.fix(nodes=[i], dof="UZ", value=EPS_YY * z) .. GENERATED FROM PYTHON SOURCE LINES 71-73 Solve + verify strain recovery ------------------------------ .. GENERATED FROM PYTHON SOURCE LINES 73-82 .. code-block:: Python res = m.solve() eps = m.eel(res.displacement) print(f"εxx prescribed = {EPS_XX:.3e} | recovered (mean) = {eps[:, 0].mean():.3e}") print(f"εyy prescribed = {EPS_YY:.3e} | recovered (mean) = {eps[:, 1].mean():.3e}") np.testing.assert_allclose(eps[:, 0], EPS_XX, rtol=1e-10) np.testing.assert_allclose(eps[:, 1], EPS_YY, rtol=1e-10) np.testing.assert_allclose(eps[:, 2], EPS_YY, rtol=1e-10) .. rst-class:: sphx-glr-script-out .. code-block:: none εxx prescribed = 5.000e-04 | recovered (mean) = 5.000e-04 εyy prescribed = -1.500e-04 | recovered (mean) = -1.500e-04 .. GENERATED FROM PYTHON SOURCE LINES 83-85 Render the deformed tet coloured by εxx --------------------------------------- .. GENERATED FROM PYTHON SOURCE LINES 85-103 .. code-block:: Python grid = femorph_solver.io.static_result_to_grid(m, res) node_nums = m.node_numbers node_to_idx = {int(nn): i for i, nn in enumerate(node_nums)} grid.point_data["eps_xx"] = np.array( [eps[node_to_idx[int(nn)], 0] for nn in grid.point_data["ansys_node_num"]] ) plotter = pv.Plotter(off_screen=True) plotter.add_mesh(grid, style="wireframe", color="gray", line_width=1, opacity=0.4) plotter.add_mesh( grid.warp_by_vector("displacement", factor=200.0), scalars="eps_xx", show_edges=True, cmap="viridis", scalar_bar_args={"title": "εxx"}, ) plotter.add_axes() plotter.show() .. image-sg:: /gallery/elements/solid187/images/sphx_glr_example_solid187_001.png :alt: example solid187 :srcset: /gallery/elements/solid187/images/sphx_glr_example_solid187_001.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-timing **Total running time of the script:** (0 minutes 0.232 seconds) .. _sphx_glr_download_gallery_elements_solid187_example_solid187.py: .. only:: html .. container:: sphx-glr-footer sphx-glr-footer-example .. container:: sphx-glr-download sphx-glr-download-jupyter :download:`Download Jupyter notebook: example_solid187.ipynb ` .. container:: sphx-glr-download sphx-glr-download-python :download:`Download Python source code: example_solid187.py ` .. container:: sphx-glr-download sphx-glr-download-zip :download:`Download zipped: example_solid187.zip ` .. only:: html .. rst-class:: sphx-glr-signature `Gallery generated by Sphinx-Gallery `_