Note
Go to the end to download the full example code.
PyVista Mesh Integration#
Run a modal analysis on a mesh generated from pyvista within MAPDL.
import os
import tempfile
from ansys.mapdl.reader import save_as_archive
import pyvista as pv
from ansys.mapdl.core import launch_mapdl
# launch MAPDL as a service
mapdl = launch_mapdl(loglevel="ERROR")
# Create a simple plane mesh centered at (0, 0, 0) on the XY plane
mesh = pv.Plane(i_resolution=100, j_resolution=100)
mesh.plot(color="w", show_edges=True)
Write the mesh to a temporary file
archive_filename = os.path.join(tempfile.gettempdir(), "tmp.cdb")
save_as_archive(archive_filename, mesh)
# Read in the archive file
response = mapdl.cdread("db", archive_filename)
mapdl.prep7()
print(mapdl.shpp("SUMM"))
# specify shell thickness
mapdl.sectype(1, "shell")
mapdl.secdata(0.01)
mapdl.emodif("ALL", "SECNUM", 1)
# specify material properties
# using aprox values for AISI 5000 Series Steel
# http://www.matweb.com/search/datasheet.aspx?matguid=89d4b891eece40fbbe6b71f028b64e9e
mapdl.units("SI") # not necessary, but helpful for book keeping
mapdl.mp("EX", 1, 200e9) # Elastic moduli in Pa (kg/(m*s**2))
mapdl.mp("DENS", 1, 7800) # Density in kg/m3
mapdl.mp("NUXY", 1, 0.3) # Poissons Ratio
mapdl.emodif("ALL", "MAT", 1)
# Run an unconstrained modal analysis
# for the first 20 modes above 1 Hz
mapdl.modal_analysis(nmode=20, freqb=1)
# you could have also run:
# mapdl.run('/SOLU')
# mapdl.antype('MODAL') # default NEW
# mapdl.modopt('LANB', 20, 1)
# mapdl.solve()
SUMMARIZE SHAPE TESTING FOR ALL SELECTED ELEMENTS
------------------------------------------------------------------------------
<<<<<< SHAPE TESTING SUMMARY >>>>>>
<<<<<< FOR ALL SELECTED ELEMENTS >>>>>>
------------------------------------------------------------------------------
--------------------------------------
| Element count 10000 SHELL181 |
--------------------------------------
Test Number tested Warning count Error count Warn+Err %
---- ------------- ------------- ----------- ----------
Aspect Ratio 10000 0 0 0.00 %
Parallel Deviation 10000 0 0 0.00 %
Maximum Angle 10000 0 0 0.00 %
Jacobian Ratio 10000 0 0 0.00 %
Warping Factor 10000 0 0 0.00 %
Any 10000 0 0 0.00 %
------------------------------------------------------------------------------
***** MAPDL SOLVE COMMAND *****
*** NOTE *** CP = 0.000 TIME= 00:00:00
There is no title defined for this analysis.
*** SELECTION OF ELEMENT TECHNOLOGIES FOR APPLICABLE ELEMENTS ***
---GIVE SUGGESTIONS ONLY---
ELEMENT TYPE 4 IS SHELL181. IT IS ASSOCIATED WITH ELASTOPLASTIC
MATERIALS ONLY. KEYOPT(8)=2 IS SUGGESTED AND KEYOPT(3)=2 IS SUGGESTED FOR
HIGHER ACCURACY OF MEMBRANE STRESSES; OTHERWISE, KEYOPT(3)=0 IS SUGGESTED.
*****MAPDL VERIFICATION RUN ONLY*****
DO NOT USE RESULTS FOR PRODUCTION
S O L U T I O N O P T I O N S
PROBLEM DIMENSIONALITY. . . . . . . . . . . . .3-D
DEGREES OF FREEDOM. . . . . . UX UY UZ ROTX ROTY ROTZ
ANALYSIS TYPE . . . . . . . . . . . . . . . . .MODAL
EXTRACTION METHOD. . . . . . . . . . . . . .BLOCK LANCZOS
NUMBER OF MODES TO EXTRACT. . . . . . . . . . . 20
SHIFT POINT . . . . . . . . . . . . . . . . . . 1.0000
GLOBALLY ASSEMBLED MATRIX . . . . . . . . . . .SYMMETRIC
*** WARNING *** CP = 0.000 TIME= 00:00:00
No constraints have been defined using the D command.
*** NOTE *** CP = 0.000 TIME= 00:00:00
The step data was checked and warning messages were found.
Please review output or errors file ( ) for these warning messages.
*** NOTE *** CP = 0.000 TIME= 00:00:00
The conditions for direct assembly have been met. No .emat or .erot
files will be produced.
D I S T R I B U T E D D O M A I N D E C O M P O S E R
...Number of elements: 10000
...Number of nodes: 10201
...Decompose to 0 CPU domains
...Element load balance ratio = 0.000
L O A D S T E P O P T I O N S
LOAD STEP NUMBER. . . . . . . . . . . . . . . . 1
THERMAL STRAINS INCLUDED IN THE LOAD VECTOR . . YES
*********** PRECISE MASS SUMMARY ***********
TOTAL RIGID BODY MASS MATRIX ABOUT ORIGIN
Translational mass | Coupled translational/rotational mass
78.000 0.0000 0.0000 | 0.0000 0.0000 -0.93259E-14
0.0000 78.000 0.0000 | 0.0000 0.0000 -0.35527E-14
0.0000 0.0000 78.000 | 0.93259E-14 0.17764E-14 0.0000
------------------------------------------ | ------------------------------------------
| Rotational mass (inertia)
| 6.5006 -0.63664E-15 0.0000
| -0.26108E-15 6.5006 0.0000
| 0.0000 0.0000 12.999
TOTAL MASS = 78.000
The mass principal axes coincide with the global Cartesian axes
CENTER OF MASS (X,Y,Z)= -0.45548E-16 0.11956E-15 0.0000
TOTAL INERTIA ABOUT CENTER OF MASS
6.5006 -0.63664E-15 0.0000
-0.63664E-15 6.5006 0.0000
0.0000 0.0000 12.999
The inertia principal axes coincide with the global Cartesian axes
*** MASS SUMMARY BY ELEMENT TYPE ***
TYPE MASS
4 78.0000
Range of element maximum matrix coefficients in global coordinates
Maximum = 749504998 at element 0.
Minimum = 749503915 at element 0.
*** ELEMENT MATRIX FORMULATION TIMES
TYPE NUMBER ENAME TOTAL CP AVE CP
4 10000 SHELL181 0.000 0.000000
Time at end of element matrix formulation CP = 0.
BLOCK LANCZOS CALCULATION OF UP TO 20 EIGENVECTORS.
NUMBER OF EQUATIONS = 61206
MAXIMUM WAVEFRONT = 0
MAXIMUM MODES STORED = 20
MINIMUM EIGENVALUE = 0.10000E+01
MAXIMUM EIGENVALUE = 0.10000E+31
Memory available (MB) = 0.0 , Memory required (MB) = 0.0
*****MAPDL VERIFICATION RUN ONLY*****
DO NOT USE RESULTS FOR PRODUCTION
*** FREQUENCIES FROM BLOCK LANCZOS ITERATION ***
MODE FREQUENCY (HERTZ)
FREQUENCY RANGE REQUESTED= 20 MODES ABOVE 1.00000 HERTZ
1 32.74664192985
2 47.77804288733
3 59.17108298327
4 84.60790505105
5 84.60790505226
6 148.8927748644
7 148.8927748651
8 154.6519861742
9 168.3840411456
10 187.7790667334
11 255.9751634914
12 255.9751634921
13 285.2444050787
14 298.2943388368
15 319.8023401629
16 319.8023401636
17 370.6023621808
18 391.8457651384
19 409.1039880842
20 482.6516135444
Load the result file within pyansys
and plot the 8th mode.
result = mapdl.result
print(result)
result.plot_nodal_displacement(7, show_displacement=True, displacement_factor=0.4)
PyMAPDL Result
Units : User Defined
Version : 24.2
Cyclic : False
Result Sets : 20
Nodes : 10201
Elements : 10000
Available Results:
NSL : Nodal displacements
plot the 1st mode using contours
result.plot_nodal_displacement(
0, show_displacement=True, displacement_factor=0.4, n_colors=10
)
Animate a high frequency mode
Get a smoother plot by disabling movie_filename and increasing n_frames
.
Enable a continuous plot looping with `loop=True`
.
result.animate_nodal_displacement(
18,
loop=False,
add_text=False,
n_frames=30,
displacement_factor=0.4,
show_axes=False,
background="w",
movie_filename="plane_vib.gif",
)
[(1.5773502691896262, 1.5773502691896262, 1.5773502691896262),
(0.0, 0.0, 0.0),
(0.0, 0.0, 1.0)]
Stop mapdl#
mapdl.exit()
Total running time of the script: (0 minutes 8.234 seconds)