.. DO NOT EDIT. .. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY. .. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE: .. "examples/Meshes/Mesh9.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_examples_Meshes_Mesh9.py: Mesh9 ===== Meshing of a perforated plate with a structured mesh. .. GENERATED FROM PYTHON SOURCE LINES 11-115 .. tab-set:: .. tab-item:: Static Scene .. image-sg:: /examples/Meshes/images/sphx_glr_Mesh9_001.png :alt: Mesh9 :srcset: /examples/Meshes/images/sphx_glr_Mesh9_001.png :class: sphx-glr-single-img .. tab-item:: Interactive Scene .. offlineviewer:: /home/docs/checkouts/readthedocs.org/user_builds/easyfea/checkouts/v1.7.0/docs/examples/Meshes/images/sphx_glr_Mesh9_001.vtksz .. tab-set:: .. tab-item:: Static Scene .. image-sg:: /examples/Meshes/images/sphx_glr_Mesh9_002.png :alt: Mesh9 :srcset: /examples/Meshes/images/sphx_glr_Mesh9_002.png :class: sphx-glr-single-img .. tab-item:: Interactive Scene .. offlineviewer:: /home/docs/checkouts/readthedocs.org/user_builds/easyfea/checkouts/v1.7.0/docs/examples/Meshes/images/sphx_glr_Mesh9_002.vtksz .. rst-class:: sphx-glr-script-out .. code-block:: none ============= Init GMSH interface ============= Meshing with gmsh (3.625 ms) Construct mesh object (31.034 ms) Element type: QUAD4 Ne = 200, Nn = 240 | .. code-block:: Python :lineno-start: 12 from EasyFEA import Display, np, Mesher, ElemType, PyVista from EasyFEA.Geoms import Point, Circle, Points, Line, CircleArc, Contour if __name__ == "__main__": dim = 2 if dim == 2: elemType = ElemType.QUAD4 else: elemType = ElemType.HEXA8 Display.Clear() # ---------------------------------------------- # Geom # ---------------------------------------------- H = 90 L = 45 D = 10 e = 20 N = 5 mS = (np.pi / 4 * D / 2) / N # PI for Points # pi for gmsh points PC = Point(L / 2, H / 2, 0) circle = Circle(PC, D, mS) P1 = Point() P2 = Point(L, 0) P3 = Point(L, H) P4 = Point(0, H) contour1 = Points( [ (P3 + P2) / 2, P3, (P3 + P4) / 2, P4, (P4 + P1) / 2, P1, (P1 + P2) / 2, P2, (P3 + P2) / 2, ], mS, ) # ---------------------------------------------- # Mesh # ---------------------------------------------- mesher = Mesher(False, True, True) factory = mesher._factory contours1: list[Contour] = [] for c in range(4): pc = circle.center contour = circle.Get_Contour() pc1 = contour.geoms[c].pt1 pc2 = contour.geoms[c].pt2 pc3 = contour.geoms[c].pt3 p1, p2, p3 = contour1.points[c * 2 : c * 2 + 3] cont1 = Contour( [Line(pc1, p1), Line(p1, p2), Line(p2, pc3), CircleArc(pc3, pc1, pc)] ) loop1, lines1, points1 = mesher._Loop_From_Geom(cont1) cont2 = Contour( [Line(pc3, p2), Line(p2, p3), Line(p3, pc2), CircleArc(pc2, pc3, pc)] ) loop2, lines2, points2 = mesher._Loop_From_Geom(cont2) surf1 = factory.addSurfaceFilling(loop1) surf2 = factory.addSurfaceFilling(loop2) mesher._Surfaces_Organize([surf1, surf2], elemType, True, [N] * 4) contours1.extend([cont1, cont2]) PyVista.Plot_Geoms(contours1).show() if dim == 3: for cont1 in contours1: cont2 = cont1.copy() cont2.Translate(dz=e) # cont2.rotate(np.pi/8, PC.coord) mesher._Link_Contours(cont1, cont2, elemType, 3, [N] * 4) mesher._Set_PhysicalGroups() mesher._Mesh_Generate(dim, elemType) mesh = mesher._Mesh_Get_Mesh() if len(mesh.orphanNodes) > 0: plotter = PyVista.Plot_Nodes(mesh, mesh.orphanNodes) plotter.add_title("Orphan nodes detected") plotter.show() PyVista.Plot_Mesh(mesh).show() .. rst-class:: sphx-glr-timing **Total running time of the script:** (0 minutes 0.906 seconds) .. _sphx_glr_download_examples_Meshes_Mesh9.py: .. only:: html .. container:: sphx-glr-footer sphx-glr-footer-example .. container:: sphx-glr-download sphx-glr-download-jupyter :download:`Download Jupyter notebook: Mesh9.ipynb ` .. container:: sphx-glr-download sphx-glr-download-python :download:`Download Python source code: Mesh9.py ` .. container:: sphx-glr-download sphx-glr-download-zip :download:`Download zipped: Mesh9.zip ` .. only:: html .. rst-class:: sphx-glr-signature `Gallery generated by Sphinx-Gallery `_