.. DO NOT EDIT.
.. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY.
.. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE:
.. "examples/gallery_examples/01-apdlmath-examples/eigen_solve.py"
.. LINE NUMBERS ARE GIVEN BELOW.

.. only:: html

    .. note::
        :class: sphx-glr-download-link-note

        :ref:`Go to the end <sphx_glr_download_examples_gallery_examples_01-apdlmath-examples_eigen_solve.py>`
        to download the full example code

.. rst-class:: sphx-glr-example-title

.. _sphx_glr_examples_gallery_examples_01-apdlmath-examples_eigen_solve.py:


.. _ref_mapdl_math_eigen_solve:

Using APDLMath to solve Eigenproblems
-------------------------------------

Use APDLMath to solve eigenproblems.

This example uses a verification manual input file, but you can use
your own sparse or dense matrices and solve those.

.. GENERATED FROM PYTHON SOURCE LINES 13-26

.. code-block:: default

    import time

    import matplotlib.pylab as plt
    import numpy as np

    from ansys.mapdl.core import launch_mapdl
    from ansys.mapdl.core.examples import vmfiles

    # Start MAPDL as a service and create an APDLMath object
    mapdl = launch_mapdl(loglevel="ERROR")
    mm = mapdl.math









.. GENERATED FROM PYTHON SOURCE LINES 27-30

First we get the `STIFF` and `MASS` matrices from the full file
after running the input file from Verification Manual 153


.. GENERATED FROM PYTHON SOURCE LINES 30-36

.. code-block:: default

    out = mapdl.input(vmfiles["vm153"])

    k = mm.stiff(fname="PRSMEMB.full")
    m = mm.mass(fname="PRSMEMB.full")









.. GENERATED FROM PYTHON SOURCE LINES 37-38

Display size of the M and K matrices

.. GENERATED FROM PYTHON SOURCE LINES 38-41

.. code-block:: default

    print(m.shape)
    print(k.shape)





.. rst-class:: sphx-glr-script-out

 .. code-block:: none

    (126, 126)
    (126, 126)




.. GENERATED FROM PYTHON SOURCE LINES 42-45

Allocate an array to store the eigenshapes.
where `nev` is the number of eigenvalues requested


.. GENERATED FROM PYTHON SOURCE LINES 45-49

.. code-block:: default

    nev = 10
    a = mm.mat(k.nrow, nev)
    a





.. rst-class:: sphx-glr-script-out

 .. code-block:: none


    Dense APDLMath Matrix (126, 10)



.. GENERATED FROM PYTHON SOURCE LINES 50-55

Perform the the modal analysis.

The algorithm is automatically chosen with respect to the matrices
properties (e.g. scalar, storage, symmetry...)


.. GENERATED FROM PYTHON SOURCE LINES 55-62

.. code-block:: default

    print("Calling MAPDL to solve the eigenproblem...")

    t1 = time.time()
    ev = mm.eigs(nev, k, m, phi=a)
    print(f"Elapsed time to solve this problem: {time.time() - t1}")






.. rst-class:: sphx-glr-script-out

 .. code-block:: none

    Calling MAPDL to solve the eigenproblem...
    Elapsed time to solve this problem: 0.028902530670166016




.. GENERATED FROM PYTHON SOURCE LINES 63-64

This is the vector of eigenfrequencies.

.. GENERATED FROM PYTHON SOURCE LINES 64-66

.. code-block:: default

    print(ev)





.. rst-class:: sphx-glr-script-out

 .. code-block:: none

    TGOEAE :
     Size : 10
      3.381e+02   3.381e+02   6.266e+02   6.266e+02   9.283e+02      <       5
      9.283e+02   1.250e+03   1.250e+03   1.424e+03   1.424e+03      <       10




.. GENERATED FROM PYTHON SOURCE LINES 67-73

Verify the accuracy of eigenresults
~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
Check the residual error for the first eigenresult
:math:`R_1=||(K-\lambda_1.M).\phi_1||_2`

First, we compute :math:`\lambda_1 = \omega_1^2 = (2.\pi.f_1)^2`

.. GENERATED FROM PYTHON SOURCE LINES 73-81

.. code-block:: default


    # Eigenfrequency (Hz)
    i = 0
    f = ev[0]
    omega = 2 * np.pi * f
    lam = omega * omega









.. GENERATED FROM PYTHON SOURCE LINES 82-84

Then we get the 1st Eigenshape :math:`\phi_1`, and compute
:math:`K.\phi_1` and :math:`M.\phi_1`

.. GENERATED FROM PYTHON SOURCE LINES 84-95

.. code-block:: default


    # shape
    phi = a[0]

    # APDL Command: *MULT,K,,Phi,,KPhi
    kphi = k.dot(phi)

    # APDL Command: *MULT,M,,Phi,,MPhi
    mphi = m.dot(phi)









.. GENERATED FROM PYTHON SOURCE LINES 96-98

Next, compute the :math:`||K.\phi_1||_2` quantity and normalize the
residual value.

.. GENERATED FROM PYTHON SOURCE LINES 98-107

.. code-block:: default


    # APDL Command: *MULT,K,,Phi,,KPhi
    kphi = k.dot(phi)


    # APDL Command: *NRM,KPhi,NRM2,KPhiNrm
    kphinrm = kphi.norm()









.. GENERATED FROM PYTHON SOURCE LINES 108-111

Then we add these two vectors, using the :math:`\lambda_1` scalar
factor and finally compute the normalized residual value
:math:`\frac{R_1}{||K.\phi_1||_2}`

.. GENERATED FROM PYTHON SOURCE LINES 111-120

.. code-block:: default


    # APDL Command: *AXPY,-lambda,,MPhi,1,,KPhi
    mphi *= lam
    kphi -= mphi

    # Compute the residual
    res = kphi.norm() / kphinrm
    print(res)





.. rst-class:: sphx-glr-script-out

 .. code-block:: none

    3.503338761293556e-11




.. GENERATED FROM PYTHON SOURCE LINES 121-123

This residual can be computed for all eigenmodes


.. GENERATED FROM PYTHON SOURCE LINES 123-163

.. code-block:: default



    def get_res(i):
        """Compute the residual for a given eigenmode"""
        # Eigenfrequency (Hz)
        f = ev[i]

        # omega = 2.pi.Frequency
        omega = 2 * np.pi * f

        # lambda = omega^2
        lam = omega * omega

        # i-th eigenshape
        phi = a[i]

        # K.Phi
        kphi = k.dot(phi)

        # M.Phi
        mphi = m.dot(phi)

        # Normalization scalar value
        kphinrm = kphi.norm()

        # (K-\lambda.M).Phi
        mphi *= lam
        kphi -= mphi

        # return the residual
        return kphi.norm() / kphinrm


    mapdl_acc = np.zeros(nev)

    for i in range(nev):
        f = ev[i]
        mapdl_acc[i] = get_res(i)
        print(f"[{i}] : Freq = {f}\t - Residual = {mapdl_acc[i]}")





.. rst-class:: sphx-glr-script-out

 .. code-block:: none

    [0] : Freq = 338.0666635506364   - Residual = 3.503338761293556e-11
    [1] : Freq = 338.06666355063646  - Residual = 5.095765571157227e-11
    [2] : Freq = 626.6450980927036   - Residual = 1.7852148894677666e-11
    [3] : Freq = 626.645098092704    - Residual = 2.17372918531887e-11
    [4] : Freq = 928.2598500574518   - Residual = 1.72025442075496e-11
    [5] : Freq = 928.2598500574528   - Residual = 6.747910207805831e-12
    [6] : Freq = 1249.8421074363514  - Residual = 2.815524037478371e-11
    [7] : Freq = 1249.8421074363525  - Residual = 7.546271373844357e-12
    [8] : Freq = 1423.993890941669   - Residual = 4.5372004789272195e-10
    [9] : Freq = 1423.9938909416699  - Residual = 1.2160010951553228e-09




.. GENERATED FROM PYTHON SOURCE LINES 164-165

Plot Accuracy of Eigenresults

.. GENERATED FROM PYTHON SOURCE LINES 165-177

.. code-block:: default


    fig = plt.figure(figsize=(12, 10))
    ax = plt.axes()
    x = np.linspace(1, nev, nev)
    plt.title("APDL Math Residual Error (%)")
    plt.yscale("log")
    plt.ylim([10e-13, 10e-7])
    plt.xlabel("Frequency #")
    plt.ylabel("Errors (%)")
    ax.bar(x, mapdl_acc, label="MAPDL Results")
    plt.show()




.. image-sg:: /examples/gallery_examples/01-apdlmath-examples/images/sphx_glr_eigen_solve_001.png
   :alt: APDL Math Residual Error (%)
   :srcset: /examples/gallery_examples/01-apdlmath-examples/images/sphx_glr_eigen_solve_001.png
   :class: sphx-glr-single-img





.. GENERATED FROM PYTHON SOURCE LINES 178-181

Stop mapdl
~~~~~~~~~~


.. GENERATED FROM PYTHON SOURCE LINES 181-182

.. code-block:: default

    mapdl.exit()








.. rst-class:: sphx-glr-timing

   **Total running time of the script:** ( 0 minutes  1.275 seconds)


.. _sphx_glr_download_examples_gallery_examples_01-apdlmath-examples_eigen_solve.py:

.. only:: html

  .. container:: sphx-glr-footer sphx-glr-footer-example




    .. container:: sphx-glr-download sphx-glr-download-python

      :download:`Download Python source code: eigen_solve.py <eigen_solve.py>`

    .. container:: sphx-glr-download sphx-glr-download-jupyter

      :download:`Download Jupyter notebook: eigen_solve.ipynb <eigen_solve.ipynb>`


.. only:: html

 .. rst-class:: sphx-glr-signature

    `Gallery generated by Sphinx-Gallery <https://sphinx-gallery.github.io>`_