.. DO NOT EDIT. .. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY. .. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE: .. "gallery/analyses/cyclic/example_rotor_sector.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_analyses_cyclic_example_rotor_sector.py: Cyclic-symmetry modal on a rotor sector ======================================= A rotor with :math:`N` identical sectors spanning 360° has a spectrum that factors by *harmonic index* :math:`k \in \{0, 1, \dots, N/2\}`. Analysing one base sector plus the constraint .. math:: u_\text{high} = e^{i k \alpha}\, R(\alpha)\, u_\text{low}, \qquad \alpha = 2\pi / N, recovers the full-rotor spectrum at a fraction of the DOF count. This example builds a simple 4-sector rotor from a pyvista annular sector, wraps it in a :class:`~femorph_solver.CyclicModel`, and runs :meth:`~femorph_solver.CyclicModel.modal_solve` across all harmonic indices. ``CyclicModel`` auto-detects the low/high cyclic faces from the rotation axis and angle, so the caller never has to assemble a manual DOF pair list. .. GENERATED FROM PYTHON SOURCE LINES 23-32 .. code-block:: Python from __future__ import annotations import numpy as np import pyvista as pv import femorph_solver as fs from femorph_solver import ELEMENTS .. GENERATED FROM PYTHON SOURCE LINES 33-35 Build a 90° annular sector -------------------------- .. GENERATED FROM PYTHON SOURCE LINES 35-51 .. code-block:: Python N_SECTORS = 4 ALPHA = 2.0 * np.pi / N_SECTORS N_T, N_R, N_Z = 12, 6, 3 theta = np.linspace(0.0, ALPHA, N_T + 1) r = np.linspace(0.8, 1.0, N_R + 1) z = np.linspace(0.0, 0.1, N_Z + 1) # (r, θ, z) order gives a CCW-from-above hex winding in the (r, θ) # plane — swapping r and θ would flip the winding inside-out. R, T, Z = np.meshgrid(r, theta, z, indexing="ij") X = R * np.cos(T) Y = R * np.sin(T) grid = pv.StructuredGrid(X, Y, Z).cast_to_unstructured_grid() .. GENERATED FROM PYTHON SOURCE LINES 52-59 Wrap as a CyclicModel --------------------- ``CyclicModel.from_grid`` builds the underlying base-sector :class:`~femorph_solver.Model`, detects the cyclic face pairs from the rotation axis (default ``"z"``) + ``n_sectors``, and validates that every low-face point has a rotated counterpart on the high face within ``pair_tol``. .. GENERATED FROM PYTHON SOURCE LINES 59-65 .. code-block:: Python cm = fs.CyclicModel.from_grid(grid, n_sectors=N_SECTORS) # Per-sector physics — element kernel + isotropic material — is set # on the inner Model exactly as for a non-cyclic analysis. cm.assign(ELEMENTS.HEX8, material={"EX": 2.1e11, "PRXY": 0.30, "DENS": 7850.0}) .. GENERATED FROM PYTHON SOURCE LINES 66-71 Solve every harmonic index -------------------------- ``modal_solve`` defaults to all harmonics ``k = 0 … N/2``, returning one :class:`~femorph_solver.solvers.cyclic.CyclicModalResult` per harmonic. Each result carries ``frequency`` and ``mode_shapes``. .. GENERATED FROM PYTHON SOURCE LINES 71-78 .. code-block:: Python results = cm.modal_solve(n_modes=4) print(f"Cyclic modal on {N_SECTORS}-sector rotor:") for res in results: fs_str = ", ".join(f"{f:.2f}" for f in res.frequency[:4]) print(f" k = {res.harmonic_index}: f [Hz] = [{fs_str}]") .. rst-class:: sphx-glr-script-out .. code-block:: none Cyclic modal on 4-sector rotor: k = 0: f [Hz] = [0.00, 0.00, 503.32, 503.32] k = 1: f [Hz] = [0.00, 0.00, 262.09, 452.92] k = 2: f [Hz] = [91.75, 91.75, 164.68, 164.68] .. GENERATED FROM PYTHON SOURCE LINES 79-81 Visualise mode 0 of k = 1 ------------------------- .. GENERATED FROM PYTHON SOURCE LINES 81-94 .. code-block:: Python r1 = next(res for res in results if res.harmonic_index == 1) phi = r1.mode_shapes[:, 0].reshape(-1, 3) grid.point_data["mode_k1"] = np.abs(phi) plotter = pv.Plotter(off_screen=True) plotter.add_mesh( grid.warp_by_vector("mode_k1", factor=0.05 / (np.max(np.abs(phi)) + 1e-30)), scalars=np.linalg.norm(phi, axis=1), cmap="viridis", scalar_bar_args={"title": f"|mode| (k=1, f={r1.frequency[0]:.1f} Hz)"}, ) plotter.add_axes() plotter.show() .. image-sg:: /gallery/analyses/cyclic/images/sphx_glr_example_rotor_sector_001.png :alt: example rotor sector :srcset: /gallery/analyses/cyclic/images/sphx_glr_example_rotor_sector_001.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-timing **Total running time of the script:** (0 minutes 1.209 seconds) .. _sphx_glr_download_gallery_analyses_cyclic_example_rotor_sector.py: .. only:: html .. container:: sphx-glr-footer sphx-glr-footer-example .. container:: sphx-glr-download sphx-glr-download-jupyter :download:`Download Jupyter notebook: example_rotor_sector.ipynb ` .. container:: sphx-glr-download sphx-glr-download-python :download:`Download Python source code: example_rotor_sector.py ` .. container:: sphx-glr-download sphx-glr-download-zip :download:`Download zipped: example_rotor_sector.zip ` .. only:: html .. rst-class:: sphx-glr-signature `Gallery generated by Sphinx-Gallery `_