.. DO NOT EDIT. .. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY. .. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE: .. "gallery/elements/link180/example_truss2_reference_geometry.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_link180_example_truss2_reference_geometry.py: TRUSS2 reference geometry — linear axial bar ============================================= The 2-node 3D truss (axial-only spar) maps to the natural coordinate :math:`s \in [-1, +1]`. Linear shape functions interpolate the three translational DOFs at each node: .. math:: N_1(s) = \tfrac{1 - s}{2}, \qquad N_2(s) = \tfrac{1 + s}{2}. Axial strain :math:`\varepsilon_{xx} = (u_2 - u_1) / L` is constant across the element; the resulting 6×6 stiffness in the global frame is .. math:: K_e = \frac{E A}{L} \begin{bmatrix} +T & -T \\ -T & +T \end{bmatrix}, \qquad T_{ij} = d_i d_j with :math:`(d_1, d_2, d_3)` the unit direction cosines along the I→J axis. References ---------- * Cook, R. D., Malkus, D. S., Plesha, M. E., Witt, R. J. (2002) *Concepts and Applications of Finite Element Analysis*, 4th ed., Wiley, §3.2. * Bathe, K.-J. (2014) *Finite Element Procedures*, 2nd ed., Prentice Hall, §5.3.1. * Przemieniecki, J. S. (1968) *Theory of Matrix Structural Analysis*, McGraw-Hill, §5.1. Implementation: :class:`femorph_solver.elements.truss2.Truss2`. .. GENERATED FROM PYTHON SOURCE LINES 38-45 .. code-block:: Python from __future__ import annotations import matplotlib.pyplot as plt import numpy as np import pyvista as pv .. GENERATED FROM PYTHON SOURCE LINES 46-47 Linear shape functions on [-1, +1] .. GENERATED FROM PYTHON SOURCE LINES 47-66 .. code-block:: Python s = np.linspace(-1.0, 1.0, 200) N1 = 0.5 * (1.0 - s) N2 = 0.5 * (1.0 + s) fig, ax = plt.subplots(1, 1, figsize=(6.4, 3.6)) ax.plot(s, N1, label="$N_1(s) = (1 - s) / 2$", color="#1f77b4", lw=2) ax.plot(s, N2, label="$N_2(s) = (1 + s) / 2$", color="#ff7f0e", lw=2) ax.scatter([-1.0, 1.0], [1.0, 1.0], color="black", zorder=5) ax.set_xlabel(r"$s$") ax.set_ylabel("$N_i(s)$") ax.set_title("TRUSS2 linear shape functions") ax.legend(loc="lower center", ncol=2, fontsize=10) ax.grid(True, ls=":", alpha=0.5) ax.set_xlim(-1.0, 1.0) ax.set_ylim(-0.05, 1.1) fig.tight_layout() fig.show() .. image-sg:: /gallery/elements/link180/images/sphx_glr_example_truss2_reference_geometry_001.png :alt: TRUSS2 linear shape functions :srcset: /gallery/elements/link180/images/sphx_glr_example_truss2_reference_geometry_001.png :class: sphx-glr-single-img .. GENERATED FROM PYTHON SOURCE LINES 67-68 Reference element with the 2-point Gauss-Legendre rule .. GENERATED FROM PYTHON SOURCE LINES 68-95 .. code-block:: Python cells = np.array([2, 0, 1]) cell_types = np.array([pv.CellType.LINE], dtype=np.uint8) truss = pv.UnstructuredGrid(cells, cell_types, np.array([[-1.0, 0.0, 0.0], [+1.0, 0.0, 0.0]])) g = 1.0 / np.sqrt(3.0) gauss = np.array([[-g, 0.0, 0.0], [+g, 0.0, 0.0]]) plotter = pv.Plotter(off_screen=True, window_size=(640, 200)) plotter.add_mesh(truss, color="black", line_width=4) plotter.add_points( np.array([[-1.0, 0.0, 0.0], [+1.0, 0.0, 0.0]]), render_points_as_spheres=True, point_size=18, color="black", label="end nodes", ) plotter.add_points( gauss, render_points_as_spheres=True, point_size=14, color="#d62728", label="2-pt Gauss-Legendre", ) plotter.view_xy() plotter.camera.zoom(1.6) plotter.show() .. image-sg:: /gallery/elements/link180/images/sphx_glr_example_truss2_reference_geometry_002.png :alt: example truss2 reference geometry :srcset: /gallery/elements/link180/images/sphx_glr_example_truss2_reference_geometry_002.png :class: sphx-glr-single-img .. GENERATED FROM PYTHON SOURCE LINES 96-97 Sanity — partition of unity at every Gauss point .. GENERATED FROM PYTHON SOURCE LINES 97-101 .. code-block:: Python sums = N1[np.searchsorted(s, [-g, +g])] + N2[np.searchsorted(s, [-g, +g])] np.testing.assert_allclose(sums, 1.0, atol=1e-12) print("OK — N_1(s) + N_2(s) = 1 at every Gauss point.") .. rst-class:: sphx-glr-script-out .. code-block:: none OK — N_1(s) + N_2(s) = 1 at every Gauss point. .. rst-class:: sphx-glr-timing **Total running time of the script:** (0 minutes 0.258 seconds) .. _sphx_glr_download_gallery_elements_link180_example_truss2_reference_geometry.py: .. only:: html .. container:: sphx-glr-footer sphx-glr-footer-example .. container:: sphx-glr-download sphx-glr-download-jupyter :download:`Download Jupyter notebook: example_truss2_reference_geometry.ipynb ` .. container:: sphx-glr-download sphx-glr-download-python :download:`Download Python source code: example_truss2_reference_geometry.py ` .. container:: sphx-glr-download sphx-glr-download-zip :download:`Download zipped: example_truss2_reference_geometry.zip ` .. only:: html .. rst-class:: sphx-glr-signature `Gallery generated by Sphinx-Gallery `_