.. DO NOT EDIT. .. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY. .. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE: .. "gallery/pre-processing/example_build_mesh.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_pre-processing_example_build_mesh.py: Building a model from a pyvista grid ==================================== When you don't have a CDB on disk, the easiest way to build a :class:`femorph_solver.Model` from scratch is to author a :class:`pyvista.UnstructuredGrid` directly and hand it to :meth:`femorph_solver.Model.from_grid`. Here we build a 4-element HEX8 cantilever beam this way — no foreign deck required. .. GENERATED FROM PYTHON SOURCE LINES 11-21 .. code-block:: Python from __future__ import annotations import numpy as np import pyvista as pv from vtkmodules.util.vtkConstants import VTK_HEXAHEDRON import femorph_solver from femorph_solver import ELEMENTS .. GENERATED FROM PYTHON SOURCE LINES 22-24 Author the nodes ---------------- .. GENERATED FROM PYTHON SOURCE LINES 24-39 .. code-block:: Python LX, LY, LZ = 4.0, 1.0, 1.0 NX = 4 points: list[tuple[float, float, float]] = [] node_id: dict[tuple[int, int, int], int] = {} nn = 0 for i in range(NX + 1): for j in range(2): for k in range(2): points.append((i * LX / NX, j * LY, k * LZ)) node_id[(i, j, k)] = nn # 0-based pyvista index nn += 1 points_arr = np.asarray(points, dtype=float) .. GENERATED FROM PYTHON SOURCE LINES 40-42 Author the hex connectivity --------------------------- .. GENERATED FROM PYTHON SOURCE LINES 42-61 .. code-block:: Python cells: list[int] = [] for i in range(NX): cells.extend( [ 8, node_id[(i, 0, 0)], node_id[(i + 1, 0, 0)], node_id[(i + 1, 1, 0)], node_id[(i, 1, 0)], node_id[(i, 0, 1)], node_id[(i + 1, 0, 1)], node_id[(i + 1, 1, 1)], node_id[(i, 1, 1)], ] ) cell_types = np.full(NX, VTK_HEXAHEDRON, dtype=np.uint8) grid = pv.UnstructuredGrid(np.asarray(cells), cell_types, points_arr) .. GENERATED FROM PYTHON SOURCE LINES 62-64 Wrap as a Model and stamp material ---------------------------------- .. GENERATED FROM PYTHON SOURCE LINES 64-69 .. code-block:: Python m = femorph_solver.Model.from_grid(grid) m.assign(ELEMENTS.HEX8, material={"EX": 2.1e11, "PRXY": 0.30, "DENS": 7850.0}) print(f"model: {m.n_nodes} nodes, {m.n_elements} elements") .. rst-class:: sphx-glr-script-out .. code-block:: none model: 20 nodes, 4 elements .. GENERATED FROM PYTHON SOURCE LINES 70-72 Apply BCs and solve ------------------- .. GENERATED FROM PYTHON SOURCE LINES 72-93 .. code-block:: Python node_nums = np.asarray(grid.point_data["ansys_node_num"]) pts = np.asarray(grid.points) fixed = node_nums[pts[:, 0] < 1e-9] m.fix(nodes=fixed.tolist(), dof="ALL") tip_mask = pts[:, 0] > LX - 1e-9 tip_nodes = node_nums[tip_mask] for tip_nn in tip_nodes: m.apply_force(int(tip_nn), fz=-250.0) res = m.solve() dof_map = m.dof_map() tip_uz = np.array( [ res.displacement[int(np.where((dof_map[:, 0] == nn) & (dof_map[:, 1] == 2))[0][0])] for nn in tip_nodes ] ) print(f"tip UZ (mean over 4 corner nodes): {tip_uz.mean():.4e} m") .. rst-class:: sphx-glr-script-out .. code-block:: none tip UZ (mean over 4 corner nodes): -1.4761e-06 m .. GENERATED FROM PYTHON SOURCE LINES 94-96 Render ------ .. GENERATED FROM PYTHON SOURCE LINES 96-110 .. code-block:: Python deformed = femorph_solver.io.static_result_to_grid(m, res) plotter = pv.Plotter(off_screen=True) plotter.add_mesh(deformed, style="wireframe", color="gray", opacity=0.35, label="undeformed") plotter.add_mesh( deformed.warp_by_vector("displacement", factor=500.0), scalars="displacement_magnitude", show_edges=True, cmap="viridis", scalar_bar_args={"title": "|u| [m]"}, label="deformed (×500)", ) plotter.add_legend() plotter.add_axes() plotter.show() .. image-sg:: /gallery/pre-processing/images/sphx_glr_example_build_mesh_001.png :alt: example build mesh :srcset: /gallery/pre-processing/images/sphx_glr_example_build_mesh_001.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-timing **Total running time of the script:** (0 minutes 1.112 seconds) .. _sphx_glr_download_gallery_pre-processing_example_build_mesh.py: .. only:: html .. container:: sphx-glr-footer sphx-glr-footer-example .. container:: sphx-glr-download sphx-glr-download-jupyter :download:`Download Jupyter notebook: example_build_mesh.ipynb ` .. container:: sphx-glr-download sphx-glr-download-python :download:`Download Python source code: example_build_mesh.py ` .. container:: sphx-glr-download sphx-glr-download-zip :download:`Download zipped: example_build_mesh.zip ` .. only:: html .. rst-class:: sphx-glr-signature `Gallery generated by Sphinx-Gallery `_