.. DO NOT EDIT. .. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY. .. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE: .. "gallery/elements/solid186/example_solid186.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_solid186_example_solid186.py: .. _ref_solid186_example: HEX20 — uniaxial tension on a 20-node hex ============================================ Single HEX20 brick (a unit cube, 20 nodes including the mid-edge nodes of the serendipity family) loaded in uniaxial tension. The consistent-load vector for a serendipity 8-node face has corner weight :math:`-1/12` and mid-edge weight :math:`+1/3`; applying these produces a uniform σxx field, so ``εxx = σ/E`` and ``εyy = εzz = −ν·εxx`` exactly at every node. .. 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_HEXAHEDRON import femorph_solver from femorph_solver import ELEMENTS .. GENERATED FROM PYTHON SOURCE LINES 25-29 Reference 20-node unit cube --------------------------- Corners 1-8 in VTK_QUADRATIC_HEXAHEDRON order, then mid-edge nodes 9-20 on the bottom, top, and vertical edges. .. GENERATED FROM PYTHON SOURCE LINES 29-59 .. code-block:: Python E = 2.1e11 # Pa NU = 0.30 F_TOTAL = 4.2e4 # N coords = np.array( [ [0.0, 0.0, 0.0], [1.0, 0.0, 0.0], [1.0, 1.0, 0.0], [0.0, 1.0, 0.0], [0.0, 0.0, 1.0], [1.0, 0.0, 1.0], [1.0, 1.0, 1.0], [0.0, 1.0, 1.0], [0.5, 0.0, 0.0], [1.0, 0.5, 0.0], [0.5, 1.0, 0.0], [0.0, 0.5, 0.0], [0.5, 0.0, 1.0], [1.0, 0.5, 1.0], [0.5, 1.0, 1.0], [0.0, 0.5, 1.0], [0.0, 0.0, 0.5], [1.0, 0.0, 0.5], [1.0, 1.0, 0.5], [0.0, 1.0, 0.5], ], dtype=np.float64, ) .. GENERATED FROM PYTHON SOURCE LINES 60-62 Build the model --------------- .. GENERATED FROM PYTHON SOURCE LINES 62-87 .. code-block:: Python cells = np.concatenate([[20], np.arange(20, dtype=np.int64)]) cell_types = np.array([VTK_QUADRATIC_HEXAHEDRON], dtype=np.uint8) grid = pv.UnstructuredGrid(cells, cell_types, coords) m = femorph_solver.Model.from_grid(grid) m.assign(ELEMENTS.HEX20, material={"EX": E, "PRXY": NU}) # Symmetry BC: x=0 face UX, y=0 face UY, z=0 face UZ. x0 = [1, 4, 5, 8, 12, 16, 17, 20] y0 = [1, 2, 5, 6, 9, 13, 17, 18] z0 = [1, 2, 3, 4, 9, 10, 11, 12] for nn in x0: m.fix(nodes=[nn], dof="UX", value=0.0) for nn in y0: m.fix(nodes=[nn], dof="UY", value=0.0) for nn in z0: m.fix(nodes=[nn], dof="UZ", value=0.0) # Consistent-load on x=1 face: 4 corners × −F/12 + 4 mid-edges × +F/3 # (integrates to F_TOTAL exactly for a uniform traction). for nn in (2, 3, 6, 7): m.apply_force(nn, fx=-F_TOTAL / 12.0) for nn in (10, 14, 18, 19): m.apply_force(nn, fx=F_TOTAL / 3.0) .. GENERATED FROM PYTHON SOURCE LINES 88-90 Static solve and strain recovery -------------------------------- .. GENERATED FROM PYTHON SOURCE LINES 90-98 .. code-block:: Python res = m.solve() eps = m.eel(res.displacement) eps_xx_expected = F_TOTAL / E eps_yy_expected = -NU * eps_xx_expected print(f"εxx expected = {eps_xx_expected:.3e} | recovered (mean) = {eps[:, 0].mean():.3e}") print(f"εyy expected = {eps_yy_expected:.3e} | recovered (mean) = {eps[:, 1].mean():.3e}") .. rst-class:: sphx-glr-script-out .. code-block:: none εxx expected = 2.000e-07 | recovered (mean) = 2.000e-07 εyy expected = -6.000e-08 | recovered (mean) = -6.000e-08 .. GENERATED FROM PYTHON SOURCE LINES 99-101 Plot the deformed cube, coloured by εxx --------------------------------------- .. GENERATED FROM PYTHON SOURCE LINES 101-119 .. 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=2.0e5), scalars="eps_xx", show_edges=True, cmap="viridis", scalar_bar_args={"title": "εxx"}, ) plotter.add_axes() plotter.show() .. image-sg:: /gallery/elements/solid186/images/sphx_glr_example_solid186_001.png :alt: example solid186 :srcset: /gallery/elements/solid186/images/sphx_glr_example_solid186_001.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-timing **Total running time of the script:** (0 minutes 0.251 seconds) .. _sphx_glr_download_gallery_elements_solid186_example_solid186.py: .. only:: html .. container:: sphx-glr-footer sphx-glr-footer-example .. container:: sphx-glr-download sphx-glr-download-jupyter :download:`Download Jupyter notebook: example_solid186.ipynb ` .. container:: sphx-glr-download sphx-glr-download-python :download:`Download Python source code: example_solid186.py ` .. container:: sphx-glr-download sphx-glr-download-zip :download:`Download zipped: example_solid186.zip ` .. only:: html .. rst-class:: sphx-glr-signature `Gallery generated by Sphinx-Gallery `_