Note
Go to the end to download the full example code
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.
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
First we get the STIFF and MASS matrices from the full file after running the input file from Verification Manual 153
Display size of the M and K matrices
print(m.shape)
print(k.shape)
(126, 126)
(126, 126)
Allocate an array to store the eigenshapes. where nev is the number of eigenvalues requested
Dense APDLMath Matrix (126, 10)
Perform the the modal analysis.
The algorithm is automatically chosen with respect to the matrices properties (e.g. scalar, storage, symmetry…)
Calling MAPDL to solve the eigenproblem...
Elapsed time to solve this problem: 0.04000568389892578
This is the vector of eigenfrequencies.
print(ev)
IHDCCA :
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
Verify the accuracy of eigenresults#
Check the residual error for the first eigenresult
First, we compute
Then we get the 1st Eigenshape
# shape
phi = a[0]
# APDL Command: *MULT,K,,Phi,,KPhi
kphi = k.dot(phi)
# APDL Command: *MULT,M,,Phi,,MPhi
mphi = m.dot(phi)
Next, compute the
# APDL Command: *MULT,K,,Phi,,KPhi
kphi = k.dot(phi)
# APDL Command: *NRM,KPhi,NRM2,KPhiNrm
kphinrm = kphi.norm()
Then we add these two vectors, using the
3.503338761293556e-11
This residual can be computed for all eigenmodes
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]}")
[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
Plot Accuracy of Eigenresults

Stop mapdl#
mapdl.exit()
Total running time of the script: ( 0 minutes 1.359 seconds)